3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
integrationpointsource.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
27#ifndef DUMUX_INTEGRATION_POINTSOURCE_HH
28#define DUMUX_INTEGRATION_POINTSOURCE_HH
29
30#include <type_traits>
31#include <dune/common/reservedvector.hh>
35
36namespace Dumux {
37
43template<class GlobalPosition, class SourceValues, typename IdType = std::size_t>
44class IntegrationPointSource : public IdPointSource<GlobalPosition, SourceValues, IdType>
45{
47 using Scalar = std::decay_t<decltype(std::declval<SourceValues>()[0])>;
48
49public:
51 [[deprecated("Call constructor with a single element index instead! Will be removed after 3.3")]]
53 Scalar qpweight, Scalar integrationElement,
54 const std::vector<std::size_t>& elementIndices)
55 : ParentType(pos, values, id),
56 qpweight_(qpweight), integrationElement_(integrationElement),
57 elementIndex_(elementIndices[0]) {}
58
60 // value known at the time of initialization
61 [[deprecated("Call constructor with a single element index instead! Will be removed after 3.3")]]
63 Scalar qpweight, Scalar integrationElement,
64 const std::vector<std::size_t>& elementIndices)
65 : ParentType(pos, id),
66 qpweight_(qpweight), integrationElement_(integrationElement),
67 elementIndex_(elementIndices[0]) {}
68
71 Scalar qpweight, Scalar integrationElement,
72 std::size_t elementIndex)
73 : ParentType(pos, values, id)
74 , qpweight_(qpweight)
75 , integrationElement_(integrationElement)
76 , elementIndex_(elementIndex)
77 {}
78
82 Scalar qpweight, Scalar integrationElement,
83 std::size_t elementIndex)
84 : ParentType(pos, id)
85 , qpweight_(qpweight)
86 , integrationElement_(integrationElement)
87 , elementIndex_(elementIndex)
88 {}
89
90 Scalar quadratureWeight() const
91 { return qpweight_; }
92
93 Scalar integrationElement() const
94 { return integrationElement_; }
95
96 void setQuadratureWeight(const Scalar qpWeight)
97 { return qpweight_ = qpWeight; }
98
99 void setIntegrationElement(const Scalar ie)
100 { integrationElement_ = ie; }
101
102 [[deprecated("Use elementIndex() instead. Will be removed after 3.3")]]
103 std::vector<std::size_t> elementIndices() const
104 { return std::vector<std::size_t>({elementIndex_}); }
105
107 std::size_t elementIndex() const
108 { return elementIndex_; }
109
112 {
114 return *this;
115 }
116
119 {
121 return *this;
122 }
123
124private:
125 Scalar qpweight_;
126 Scalar integrationElement_;
127 std::size_t elementIndex_;
128};
129
135{
136
137public:
139 template<class GridGeometry, class PointSource, class PointSourceMap>
140 static void computePointSourceMap(const GridGeometry& gridGeometry,
141 const std::vector<PointSource>& sources,
142 PointSourceMap& pointSourceMap,
143 const std::string& paramGroup = "")
144 {
145 for (const auto& source : sources)
146 {
147 // get the index of the element in which the point source falls
148 const auto eIdx = source.elementIndex();
149 if constexpr (GridGeometry::discMethod == DiscretizationMethod::box)
150 {
151 // check in which subcontrolvolume(s) we are
152 const auto element = gridGeometry.boundingBoxTree().entitySet().entity(eIdx);
153 auto fvGeometry = localView(gridGeometry);
154 fvGeometry.bindElement(element);
155 const auto globalPos = source.position();
156
157 static const bool boxPointSourceLumping = getParamFromGroup<bool>(paramGroup, "PointSource.EnableBoxLumping", true);
158 if (boxPointSourceLumping)
159 {
160 // loop over all sub control volumes and check if the point source is inside
161 constexpr int dim = GridGeometry::GridView::dimension;
162 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
163 for (const auto& scv : scvs(fvGeometry))
164 if (intersectsPointGeometry(globalPos, scv.geometry()))
165 scvIndices.push_back(scv.indexInElement());
166
167 // for all scvs that where tested positiv add the point sources
168 // to the element/scv to point source map
169 for (auto scvIdx : scvIndices)
170 {
171 const auto key = std::make_pair(eIdx, scvIdx);
172 if (pointSourceMap.count(key))
173 pointSourceMap.at(key).push_back(source);
174 else
175 pointSourceMap.insert({key, {source}});
176 // split equally on the number of matched scvs
177 auto& s = pointSourceMap.at(key).back();
178 s.setEmbeddings(scvIndices.size()*s.embeddings());
179 }
180 }
181
182 // distribute the sources using the basis function weights
183 else
184 {
185 const auto& localBasis = fvGeometry.feLocalBasis();
186 const auto ipLocal = element.geometry().local(globalPos);
187 using Scalar = std::decay_t<decltype(source.values()[0])>;
188 std::vector<typename Dune::FieldVector<Scalar, 1>> shapeValues;
189 localBasis.evaluateFunction(ipLocal, shapeValues);
190 for (const auto& scv : scvs(fvGeometry))
191 {
192 const auto key = std::make_pair(eIdx, scv.indexInElement());
193 if (pointSourceMap.count(key))
194 pointSourceMap.at(key).push_back(source);
195 else
196 pointSourceMap.insert({key, {source}});
197
198 // adjust the integration element
199 auto& s = pointSourceMap.at(key).back();
200 s.setIntegrationElement(shapeValues[scv.indexInElement()]*s.integrationElement());
201 }
202 }
203 }
204 else
205 {
206 // add the pointsource to the DOF map
207 const auto key = std::make_pair(eIdx, /*scvIdx=*/ 0);
208 if (pointSourceMap.count(key))
209 pointSourceMap.at(key).push_back(source);
210 else
211 pointSourceMap.insert({key, {source}});
212 }
213 }
214 }
215};
216
217} // end namespace Dumux
218
219#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A point source class, i.e. sources located at a single point in space.
The available discretization methods in Dumux.
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: geometry/intersectspointgeometry.hh:38
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Definition: adapt.hh:29
Values values() const
return the source values
Definition: pointsource.hh:111
GlobalPosition GlobalPosition
Export the position type.
Definition: pointsource.hh:55
A point source class with an identifier to attach data.
Definition: pointsource.hh:167
std::size_t IdType
export the id type
Definition: pointsource.hh:173
IdType id() const
return the sources identifier
Definition: pointsource.hh:185
IdPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:189
An integration point source class with an identifier to attach data and a quadrature weight and integ...
Definition: integrationpointsource.hh:45
IntegrationPointSource(GlobalPosition pos, SourceValues values, IdType id, Scalar qpweight, Scalar integrationElement, const std::vector< std::size_t > &elementIndices)
Constructor for integration point sources.
Definition: integrationpointsource.hh:52
IntegrationPointSource(GlobalPosition pos, SourceValues values, IdType id, Scalar qpweight, Scalar integrationElement, std::size_t elementIndex)
Constructor for integration point sources.
Definition: integrationpointsource.hh:70
std::size_t elementIndex() const
The index of the element this intergration point source is associated with.
Definition: integrationpointsource.hh:107
Scalar quadratureWeight() const
Definition: integrationpointsource.hh:90
std::vector< std::size_t > elementIndices() const
Definition: integrationpointsource.hh:103
IntegrationPointSource(GlobalPosition pos, IdType id, Scalar qpweight, Scalar integrationElement, const std::vector< std::size_t > &elementIndices)
Constructor for integration point sources, when there is no.
Definition: integrationpointsource.hh:62
IntegrationPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: integrationpointsource.hh:111
void setQuadratureWeight(const Scalar qpWeight)
Definition: integrationpointsource.hh:96
void setIntegrationElement(const Scalar ie)
Definition: integrationpointsource.hh:99
Scalar integrationElement() const
Definition: integrationpointsource.hh:93
IntegrationPointSource(GlobalPosition pos, IdType id, Scalar qpweight, Scalar integrationElement, std::size_t elementIndex)
Definition: integrationpointsource.hh:81
A helper class calculating a DOF-index to point source map.
Definition: integrationpointsource.hh:135
static void computePointSourceMap(const GridGeometry &gridGeometry, const std::vector< PointSource > &sources, PointSourceMap &pointSourceMap, const std::string &paramGroup="")
calculate a DOF index to point source map from given vector of point sources
Definition: integrationpointsource.hh:140