3.3.0
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
31#include <dune/common/reservedvector.hh>
37
39
40namespace Dumux {
41
48template<class PositionType, class ValueType>
50{
51public:
53 using Scalar = std::decay_t<decltype(std::declval<ValueType>()[0])>;
55 using GlobalPosition = PositionType;
57 using Values = ValueType;
58
61 : values_(values), pos_(pos), embeddings_(1) {}
62
64 // value known at the time of initialization
66 : values_(0.0), pos_(pos), embeddings_(1) {}
67
70 {
71 values_ += s;
72 return *this;
73 }
74
77 {
78 values_ -= s;
79 return *this;
80 }
81
84 {
85 values_ *= s;
86 return *this;
87 }
88
91 {
92 values_ /= s;
93 return *this;
94 }
95
98 {
100 return *this;
101 }
102
105 {
106 values_ = s;
107 return *this;
108 }
109
112 { return values_; }
113
116 { return pos_; }
117
119 // to the local residual in the problem in scvPointSources
120 // to be overloaded by derived classes
121 template<class Problem, class FVElementGeometry, class ElementVolumeVariables>
122 void update(const Problem &problem,
123 const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity &element,
124 const FVElementGeometry &fvGeometry,
125 const ElementVolumeVariables &elemVolVars,
126 const typename FVElementGeometry::SubControlVolume &scv)
127 {}
128
130 void setEmbeddings(std::size_t embeddings)
131 {
132 embeddings_ = embeddings;
133 }
134
146 std::size_t embeddings() const
147 {
148 return embeddings_;
149 }
150
151protected:
153private:
154 GlobalPosition pos_;
155 std::size_t embeddings_;
156};
157
165template<class GlobalPosition, class SourceValues, class I>
166class IdPointSource : public PointSource<GlobalPosition, SourceValues>
167{
169 using Scalar = typename ParentType::Scalar;
170
171public:
173 using IdType = I;
174
177 : ParentType(pos, values), id_(id) {}
178
180 // value known at the time of initialization
182 : ParentType(pos, SourceValues(0.0)), id_(id) {}
183
185 IdType id() const
186 { return id_; }
187
189 IdPointSource& operator= (const SourceValues& values)
190 {
192 return *this;
193 }
194
197 {
199 return *this;
200 }
201
202private:
203 IdType id_;
204};
205
210template<class TypeTag>
211class SolDependentPointSource : public PointSource<Dune::FieldVector<typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::ctype,
212 GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld>,
213 GetPropType<TypeTag, Properties::NumEqVector>>
214{
218 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
219 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
220 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
221 using Element = typename GridView::template Codim<0>::Entity;
222
223 static const int dimworld = GridView::dimensionworld;
224 using GlobalPosition = typename Dune::FieldVector<typename GridView::ctype, dimworld>;
225 // returns the PointSource values as PrimaryVariables
226 using ValueFunction = typename std::function<SourceValues(const Problem &problem,
227 const Element &element,
228 const FVElementGeometry &fvGeometry,
229 const ElementVolumeVariables &elemVolVars,
230 const SubControlVolume &scv)>;
231
233 using Scalar = typename ParentType::Scalar;
234
235public:
237 // value known at the time of initialization
238 SolDependentPointSource(GlobalPosition pos,
239 ValueFunction valueFunction)
240 : ParentType(pos, SourceValues(0.0)), valueFunction_(valueFunction) {}
241
243 // to the local residual in the problem in scvPointSources
244 // to be overloaded by derived classes
245 void update(const Problem &problem,
246 const Element &element,
247 const FVElementGeometry &fvGeometry,
248 const ElementVolumeVariables &elemVolVars,
249 const SubControlVolume &scv)
250 { this->values_ = valueFunction_(problem, element, fvGeometry, elemVolVars, scv); }
251
254 {
256 return *this;
257 }
258
261 {
263 return *this;
264 }
265
266private:
267 ValueFunction valueFunction_;
268};
269
277{
278public:
280 template<class GridGeometry, class PointSource, class PointSourceMap>
281 static void computePointSourceMap(const GridGeometry& gridGeometry,
282 const std::vector<PointSource>& sources,
283 PointSourceMap& pointSourceMap,
284 const std::string& paramGroup = "")
285 {
286 const auto& boundingBoxTree = gridGeometry.boundingBoxTree();
287
288 for (const auto& s : sources)
289 {
290 // compute in which elements the point source falls
291 const auto entities = intersectingEntities(s.position(), boundingBoxTree);
292
293 // continue with next point source if no intersection with the grid are found
294 if (entities.empty())
295 continue;
296
297 // make local copy of point source for the map
298 auto source = s;
299
300 // split the source values equally among all concerned entities
301 source.setEmbeddings(entities.size()*source.embeddings());
302 // loop over all concernes elements
303 for (const auto eIdx : entities)
304 {
305 if constexpr (GridGeometry::discMethod == DiscretizationMethod::box)
306 {
307 // check in which subcontrolvolume(s) we are
308 const auto element = boundingBoxTree.entitySet().entity(eIdx);
309 auto fvGeometry = localView(gridGeometry);
310 fvGeometry.bindElement(element);
311
312 const auto globalPos = source.position();
313 // loop over all sub control volumes and check if the point source is inside
314 constexpr int dim = GridGeometry::GridView::dimension;
315 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
316 for (const auto& scv : scvs(fvGeometry))
317 if (intersectsPointGeometry(globalPos, scv.geometry()))
318 scvIndices.push_back(scv.indexInElement());
319
320 // for all scvs that tested positive add the point sources
321 // to the element/scv to point source map
322 for (const auto scvIdx : scvIndices)
323 {
324 const auto key = std::make_pair(eIdx, scvIdx);
325 if (pointSourceMap.count(key))
326 pointSourceMap.at(key).push_back(source);
327 else
328 pointSourceMap.insert({key, {source}});
329 // split equally on the number of matched scvs
330 auto& s = pointSourceMap.at(key).back();
331 s.setEmbeddings(scvIndices.size()*s.embeddings());
332 }
333 }
334 else
335 {
336 // add the pointsource to the DOF map
337 const auto key = std::make_pair(eIdx, /*scvIdx=*/ 0);
338 if (pointSourceMap.count(key))
339 pointSourceMap.at(key).push_back(source);
340 else
341 pointSourceMap.insert({key, {source}});
342 }
343 }
344 }
345 }
346};
347
348} // end namespace Dumux
349
350#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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
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: geometry/intersectingentities.hh:100
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
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:50
void setEmbeddings(std::size_t embeddings)
set the number of embeddings for this point source
Definition: pointsource.hh:130
PointSource & operator*=(Scalar s)
Convenience *= operator overload modifying only the values.
Definition: pointsource.hh:83
ValueType Values
Export the value type.
Definition: pointsource.hh:57
PointSource & operator+=(Scalar s)
Convenience += operator overload modifying only the values.
Definition: pointsource.hh:69
PointSource(GlobalPosition pos)
Constructor for sol dependent point sources, when there is no.
Definition: pointsource.hh:65
std::decay_t< decltype(std::declval< ValueType >()[0])> Scalar
Export the scalar type.
Definition: pointsource.hh:53
Values values() const
return the source values
Definition: pointsource.hh:111
PositionType GlobalPosition
Export the position type.
Definition: pointsource.hh:55
Values values_
value of the point source for each equation
Definition: pointsource.hh:152
const GlobalPosition & position() const
return the source position
Definition: pointsource.hh:115
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:122
PointSource & operator-=(Scalar s)
Convenience -= operator overload modifying only the values.
Definition: pointsource.hh:76
PointSource(GlobalPosition pos, Values values)
Constructor for constant point sources.
Definition: pointsource.hh:60
PointSource & operator=(const Values &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:97
PointSource & operator/=(Scalar s)
Convenience /= operator overload modifying only the values.
Definition: pointsource.hh:90
std::size_t embeddings() const
get the number of embeddings for this point source
Definition: pointsource.hh:146
A point source class with an identifier to attach data.
Definition: pointsource.hh:167
I IdType
export the id type
Definition: pointsource.hh:173
IdPointSource(GlobalPosition pos, IdType id)
Constructor for sol dependent point sources, when there is no.
Definition: pointsource.hh:181
IdPointSource(GlobalPosition pos, SourceValues values, IdType id)
Constructor for constant point sources.
Definition: pointsource.hh:176
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
A point source class for time dependent point sources.
Definition: pointsource.hh:214
SolDependentPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:253
SolDependentPointSource(GlobalPosition pos, ValueFunction valueFunction)
Constructor for sol dependent point sources, when there is no.
Definition: pointsource.hh:238
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:245
A helper class calculating a sub control volume to point source map This class uses the bounding box ...
Definition: pointsource.hh:277
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: pointsource.hh:281
An axis-aligned bounding box volume hierarchy for dune grids.
Algorithms that finds which geometric entites intersect.
Detect if a point intersects a geometry.
Declares all properties used in Dumux.