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