version 3.11-dev
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_POINTSOURCE_HH
15#define DUMUX_POINTSOURCE_HH
16
17#include <functional>
18
19#include <dune/common/exceptions.hh>
20#include <dune/common/reservedvector.hh>
27
29
30namespace Dumux {
31
38template<class PositionType, class ValueType>
40{
41public:
43 using Scalar = std::decay_t<decltype(std::declval<ValueType>()[0])>;
45 using GlobalPosition = PositionType;
47 using Values = ValueType;
48
51 : values_(values), pos_(pos), embeddings_(1) {}
52
54 // value known at the time of initialization
56 : values_(0.0), pos_(pos), embeddings_(1) {}
57
60 {
61 values_ += s;
62 return *this;
63 }
64
67 {
68 values_ -= s;
69 return *this;
70 }
71
74 {
75 values_ *= s;
76 return *this;
77 }
78
81 {
82 values_ /= s;
83 return *this;
84 }
85
88 {
90 return *this;
91 }
92
95 {
96 values_ = s;
97 return *this;
98 }
99
102 { return values_; }
103
106 { return pos_; }
107
109 // to the local residual in the problem in scvPointSources
110 // to be overloaded by derived classes
111 template<class Problem, class FVElementGeometry, class ElementVolumeVariables>
112 void update(const Problem &problem,
113 const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity &element,
114 const FVElementGeometry &fvGeometry,
115 const ElementVolumeVariables &elemVolVars,
116 const typename FVElementGeometry::SubControlVolume &scv)
117 {}
118
120 void setEmbeddings(std::size_t embeddings)
121 {
122 embeddings_ = embeddings;
123 }
124
136 std::size_t embeddings() const
137 {
138 return embeddings_;
139 }
140
141protected:
143private:
144 GlobalPosition pos_;
145 std::size_t embeddings_;
146};
147
155template<class GlobalPosition, class SourceValues, class I>
156class IdPointSource : public PointSource<GlobalPosition, SourceValues>
157{
159 using Scalar = typename ParentType::Scalar;
160
161public:
163 using IdType = I;
164
167 : ParentType(pos, values), id_(id) {}
168
170 // value known at the time of initialization
172 : ParentType(pos, SourceValues(0.0)), id_(id) {}
173
175 IdType id() const
176 { return id_; }
177
179 IdPointSource& operator= (const SourceValues& values)
180 {
182 return *this;
183 }
184
187 {
189 return *this;
190 }
191
192private:
193 IdType id_;
194};
195
200template<class TypeTag>
201class SolDependentPointSource : public PointSource<Dune::FieldVector<typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::ctype,
202 GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld>,
203 Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>>
204{
208 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
209 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
210 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
211 using Element = typename GridView::template Codim<0>::Entity;
212
213 static const int dimworld = GridView::dimensionworld;
214 using GlobalPosition = typename Dune::FieldVector<typename GridView::ctype, dimworld>;
215 // returns the PointSource values as PrimaryVariables
216 using ValueFunction = typename std::function<SourceValues(const Problem &problem,
217 const Element &element,
218 const FVElementGeometry &fvGeometry,
219 const ElementVolumeVariables &elemVolVars,
220 const SubControlVolume &scv)>;
221
223 using Scalar = typename ParentType::Scalar;
224
225public:
227 // value known at the time of initialization
228 SolDependentPointSource(GlobalPosition pos,
229 ValueFunction valueFunction)
230 : ParentType(pos, SourceValues(0.0)), valueFunction_(valueFunction) {}
231
233 // to the local residual in the problem in scvPointSources
234 // to be overloaded by derived classes
235 void update(const Problem &problem,
236 const Element &element,
237 const FVElementGeometry &fvGeometry,
238 const ElementVolumeVariables &elemVolVars,
239 const SubControlVolume &scv)
240 { this->values_ = valueFunction_(problem, element, fvGeometry, elemVolVars, scv); }
241
244 {
246 return *this;
247 }
248
251 {
253 return *this;
254 }
255
256private:
257 ValueFunction valueFunction_;
258};
259
267{
268public:
270 template<class GridGeometry, class PointSource, class PointSourceMap>
271 static void computePointSourceMap(const GridGeometry& gridGeometry,
272 const std::vector<PointSource>& sources,
273 PointSourceMap& pointSourceMap,
274 const std::string& paramGroup = "")
275 {
276 const auto& boundingBoxTree = gridGeometry.boundingBoxTree();
277
278 for (const auto& s : sources)
279 {
280 // compute in which elements the point source falls
281 const auto entities = intersectingEntities(s.position(), boundingBoxTree);
282
283 // continue with next point source if no intersection with the grid are found
284 if (entities.empty())
285 continue;
286
287 // make local copy of point source for the map
288 auto source = s;
289
290 // split the source values equally among all concerned entities
291 source.setEmbeddings(entities.size()*source.embeddings());
292
293 if constexpr (GridGeometry::discMethod == DiscretizationMethods::box
294 || GridGeometry::discMethod == DiscretizationMethods::fcdiamond)
295 {
296 // loop over all concerned elements
297 auto fvGeometry = localView(gridGeometry);
298 for (const auto eIdx : entities)
299 {
300 // check in which subcontrolvolume(s) we are
301 const auto element = boundingBoxTree.entitySet().entity(eIdx);
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 constexpr int dim = GridGeometry::GridView::dimension;
307 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
308 for (const auto& scv : scvs(fvGeometry))
309 if (intersectsPointGeometry(globalPos, fvGeometry.geometry(scv)))
310 scvIndices.push_back(scv.indexInElement());
311
312 // for all scvs that tested positive add the point sources
313 // to the element/scv to point source map
314 for (const 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 }
327 else if constexpr (GridGeometry::discMethod == DiscretizationMethods::cctpfa
328 || GridGeometry::discMethod == DiscretizationMethods::ccmpfa)
329 {
330 for (const auto eIdx : entities)
331 {
332 // add the pointsource to the DOF map
333 const auto key = std::make_pair(eIdx, /*scvIdx=*/ 0);
334 if (pointSourceMap.count(key))
335 pointSourceMap.at(key).push_back(source);
336 else
337 pointSourceMap.insert({key, {source}});
338 }
339 }
340 else
341 DUNE_THROW(Dune::NotImplemented, "BoundingBoxTreePointSourceHelper for discretization method " << GridGeometry::discMethod);
342 }
343 }
344};
345
346} // end namespace Dumux
347
348#endif
An axis-aligned bounding box volume hierarchy for dune grids.
A helper class calculating a sub control volume to point source map This class uses the bounding box ...
Definition: pointsource.hh:267
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:271
A point source class with an identifier to attach data.
Definition: pointsource.hh:157
I IdType
export the id type
Definition: pointsource.hh:163
IdPointSource(GlobalPosition pos, IdType id)
Constructor for sol dependent point sources, when there is no.
Definition: pointsource.hh:171
IdPointSource(GlobalPosition pos, SourceValues values, IdType id)
Constructor for constant point sources.
Definition: pointsource.hh:166
IdType id() const
return the sources identifier
Definition: pointsource.hh:175
IdPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:179
A point source base class.
Definition: pointsource.hh:40
void setEmbeddings(std::size_t embeddings)
set the number of embeddings for this point source
Definition: pointsource.hh:120
PointSource & operator*=(Scalar s)
Convenience *= operator overload modifying only the values.
Definition: pointsource.hh:73
ValueType Values
Export the value type.
Definition: pointsource.hh:47
PointSource & operator+=(Scalar s)
Convenience += operator overload modifying only the values.
Definition: pointsource.hh:59
PointSource(GlobalPosition pos)
Constructor for sol dependent point sources, when there is no.
Definition: pointsource.hh:55
std::decay_t< decltype(std::declval< ValueType >()[0])> Scalar
Export the scalar type.
Definition: pointsource.hh:43
Values values() const
return the source values
Definition: pointsource.hh:101
PositionType GlobalPosition
Export the position type.
Definition: pointsource.hh:45
Values values_
value of the point source for each equation
Definition: pointsource.hh:142
const GlobalPosition & position() const
return the source position
Definition: pointsource.hh:105
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:112
PointSource & operator-=(Scalar s)
Convenience -= operator overload modifying only the values.
Definition: pointsource.hh:66
PointSource(GlobalPosition pos, Values values)
Constructor for constant point sources.
Definition: pointsource.hh:50
PointSource & operator=(const Values &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:87
PointSource & operator/=(Scalar s)
Convenience /= operator overload modifying only the values.
Definition: pointsource.hh:80
std::size_t embeddings() const
get the number of embeddings for this point source
Definition: pointsource.hh:136
A point source class for time dependent point sources.
Definition: pointsource.hh:204
SolDependentPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:243
SolDependentPointSource(GlobalPosition pos, ValueFunction valueFunction)
Constructor for sol dependent point sources, when there is no.
Definition: pointsource.hh:228
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:235
Defines all properties used in Dumux.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:34
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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:28
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:102
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Algorithms that finds which geometric entities intersect.
Detect if a point intersects a geometry.
The available discretization methods in Dumux.
constexpr CCMpfa ccmpfa
Definition: method.hh:146
constexpr FCDiamond fcdiamond
Definition: method.hh:152
constexpr CCTpfa cctpfa
Definition: method.hh:145
constexpr Box box
Definition: method.hh:147
Definition: adapt.hh:17
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.