3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
pointsource.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 *****************************************************************************/
26#ifndef DUMUX_POINTSOURCE_HH
27#define DUMUX_POINTSOURCE_HH
28
29#include <functional>
30
36
38
39namespace Dumux {
40
47template<class PositionType, class ValueType>
49{
50public:
52 using Scalar = std::decay_t<decltype(std::declval<ValueType>()[0])>;
54 using GlobalPosition = PositionType;
56 using Values = ValueType;
57
60 : values_(values), pos_(pos), embeddings_(1) {}
61
63 // value known at the time of initialization
65 : values_(0.0), pos_(pos), embeddings_(1) {}
66
69 {
70 values_ += s;
71 return *this;
72 }
73
76 {
77 values_ -= s;
78 return *this;
79 }
80
83 {
84 values_ *= s;
85 return *this;
86 }
87
90 {
91 values_ /= s;
92 return *this;
93 }
94
97 {
99 return *this;
100 }
101
104 {
105 values_ = s;
106 return *this;
107 }
108
111 { return values_; }
112
115 { return pos_; }
116
118 // to the local residual in the problem in scvPointSources
119 // to be overloaded by derived classes
120 template<class Problem, class FVElementGeometry, class ElementVolumeVariables>
121 void update(const Problem &problem,
122 const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity &element,
123 const FVElementGeometry &fvGeometry,
124 const ElementVolumeVariables &elemVolVars,
125 const typename FVElementGeometry::SubControlVolume &scv)
126 {}
127
129 void setEmbeddings(std::size_t embeddings)
130 {
131 embeddings_ = embeddings;
132 }
133
145 std::size_t embeddings() const
146 {
147 return embeddings_;
148 }
149
150protected:
152private:
153 GlobalPosition pos_;
154 std::size_t embeddings_;
155};
156
164template<class GlobalPosition, class SourceValues, class I>
165class IdPointSource : public PointSource<GlobalPosition, SourceValues>
166{
168 using Scalar = typename ParentType::Scalar;
169
170public:
172 using IdType = I;
173
176 : ParentType(pos, values), id_(id) {}
177
179 // value known at the time of initialization
181 : ParentType(pos, SourceValues(0.0)), id_(id) {}
182
184 IdType id() const
185 { return id_; }
186
188 IdPointSource& operator= (const SourceValues& values)
189 {
191 return *this;
192 }
193
196 {
198 return *this;
199 }
200
201private:
202 IdType id_;
203};
204
209template<class TypeTag>
210class SolDependentPointSource : public PointSource<Dune::FieldVector<typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::ctype,
211 GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld>,
212 GetPropType<TypeTag, Properties::NumEqVector>>
213{
217 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
218 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
219 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
220 using Element = typename GridView::template Codim<0>::Entity;
221
222 static const int dimworld = GridView::dimensionworld;
223 using GlobalPosition = typename Dune::FieldVector<typename GridView::ctype, dimworld>;
224 // returns the PointSource values as PrimaryVariables
225 using ValueFunction = typename std::function<SourceValues(const Problem &problem,
226 const Element &element,
227 const FVElementGeometry &fvGeometry,
228 const ElementVolumeVariables &elemVolVars,
229 const SubControlVolume &scv)>;
230
232 using Scalar = typename ParentType::Scalar;
233
234public:
236 // value known at the time of initialization
237 SolDependentPointSource(GlobalPosition pos,
238 ValueFunction valueFunction)
239 : ParentType(pos, SourceValues(0.0)), valueFunction_(valueFunction) {}
240
242 // to the local residual in the problem in scvPointSources
243 // to be overloaded by derived classes
244 void update(const Problem &problem,
245 const Element &element,
246 const FVElementGeometry &fvGeometry,
247 const ElementVolumeVariables &elemVolVars,
248 const SubControlVolume &scv)
249 { this->values_ = valueFunction_(problem, element, fvGeometry, elemVolVars, scv); }
250
253 {
255 return *this;
256 }
257
260 {
262 return *this;
263 }
264
265private:
266 ValueFunction valueFunction_;
267};
268
276{
277public:
279 template<class GridGeometry, class PointSource, class PointSourceMap>
280 static void computePointSourceMap(const GridGeometry& gridGeometry,
281 std::vector<PointSource>& sources,
282 PointSourceMap& pointSourceMap)
283 {
284 constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethod::box;
285
286 const auto& boundingBoxTree = gridGeometry.boundingBoxTree();
287
288 for (auto&& source : sources)
289 {
290 // compute in which elements the point source falls
291 const auto entities = intersectingEntities(source.position(), boundingBoxTree);
292 // split the source values equally among all concerned entities
293 source.setEmbeddings(entities.size()*source.embeddings());
294 // loop over all concernes elements
295 for (unsigned int eIdx : entities)
296 {
297 if(isBox)
298 {
299 // check in which subcontrolvolume(s) we are
300 const auto element = boundingBoxTree.entitySet().entity(eIdx);
301 auto fvGeometry = localView(gridGeometry);
302 fvGeometry.bindElement(element);
303
304 const auto globalPos = source.position();
305 // loop over all sub control volumes and check if the point source is inside
306 std::vector<unsigned int> scvIndices;
307 for (auto&& scv : scvs(fvGeometry))
308 {
309 if (intersectsPointGeometry(globalPos, scv.geometry()))
310 scvIndices.push_back(scv.indexInElement());
311 }
312 // for all scvs that where tested positiv add the point sources
313 // to the element/scv to point source map
314 for (auto scvIdx : scvIndices)
315 {
316 const auto key = std::make_pair(eIdx, scvIdx);
317 if (pointSourceMap.count(key))
318 pointSourceMap.at(key).push_back(source);
319 else
320 pointSourceMap.insert({key, {source}});
321 // split equally on the number of matched scvs
322 auto& s = pointSourceMap.at(key).back();
323 s.setEmbeddings(scvIndices.size()*s.embeddings());
324 }
325 }
326 else
327 {
328 // add the pointsource to the DOF map
329 const auto key = std::make_pair(eIdx, /*scvIdx=*/ 0);
330 if (pointSourceMap.count(key))
331 pointSourceMap.at(key).push_back(source);
332 else
333 pointSourceMap.insert({key, {source}});
334 }
335 }
336 }
337 }
338};
339
340} // end namespace Dumux
341
342#endif
An axis-aligned bounding box volume hierarchy for dune grids.
Detect if a point intersects a geometry.
Algorithms that finds which geometric entites intersect.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:38
std::vector< std::size_t > intersectingEntities(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, bool isCartesianGrid=false)
Compute all intersections between entities and a point.
Definition: intersectingentities.hh:99
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
A point source base class.
Definition: pointsource.hh:49
void setEmbeddings(std::size_t embeddings)
set the number of embeddings for this point source
Definition: pointsource.hh:129
PointSource & operator*=(Scalar s)
Convenience *= operator overload modifying only the values.
Definition: pointsource.hh:82
ValueType Values
Export the value type.
Definition: pointsource.hh:56
PointSource & operator+=(Scalar s)
Convenience += operator overload modifying only the values.
Definition: pointsource.hh:68
PointSource(GlobalPosition pos)
Constructor for sol dependent point sources, when there is no.
Definition: pointsource.hh:64
std::decay_t< decltype(std::declval< ValueType >()[0])> Scalar
Export the scalar type.
Definition: pointsource.hh:52
Values values() const
return the source values
Definition: pointsource.hh:110
PositionType GlobalPosition
Export the position type.
Definition: pointsource.hh:54
Values values_
value of the point source for each equation
Definition: pointsource.hh:151
const GlobalPosition & position() const
return the source position
Definition: pointsource.hh:114
void update(const Problem &problem, const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolume &scv)
an update function called before adding the value
Definition: pointsource.hh:121
PointSource & operator-=(Scalar s)
Convenience -= operator overload modifying only the values.
Definition: pointsource.hh:75
PointSource(GlobalPosition pos, Values values)
Constructor for constant point sources.
Definition: pointsource.hh:59
PointSource & operator=(const Values &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:96
PointSource & operator/=(Scalar s)
Convenience /= operator overload modifying only the values.
Definition: pointsource.hh:89
std::size_t embeddings() const
get the number of embeddings for this point source
Definition: pointsource.hh:145
A point source class with an identifier to attach data.
Definition: pointsource.hh:166
I IdType
export the id type
Definition: pointsource.hh:172
IdPointSource(GlobalPosition pos, IdType id)
Constructor for sol dependent point sources, when there is no.
Definition: pointsource.hh:180
IdPointSource(GlobalPosition pos, SourceValues values, IdType id)
Constructor for constant point sources.
Definition: pointsource.hh:175
IdType id() const
return the sources identifier
Definition: pointsource.hh:184
IdPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:188
A point source class for time dependent point sources.
Definition: pointsource.hh:213
SolDependentPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:252
SolDependentPointSource(GlobalPosition pos, ValueFunction valueFunction)
Constructor for sol dependent point sources, when there is no.
Definition: pointsource.hh:237
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv)
an update function called before adding the value
Definition: pointsource.hh:244
A helper class calculating a sub control volume to point source map This class uses the bounding box ...
Definition: pointsource.hh:276
static void computePointSourceMap(const GridGeometry &gridGeometry, std::vector< PointSource > &sources, PointSourceMap &pointSourceMap)
calculate a DOF index to point source map from given vector of point sources
Definition: pointsource.hh:280
Declares all properties used in Dumux.