26#ifndef DUMUX_POINTSOURCE_HH
27#define DUMUX_POINTSOURCE_HH
31#include <dune/common/reservedvector.hh>
49template<
class PositionType,
class ValueType>
54 using Scalar = std::decay_t<decltype(std::declval<ValueType>()[0])>;
67 :
values_(0.0), pos_(pos), embeddings_(1) {}
122 template<
class Problem,
class FVElementGeometry,
class ElementVolumeVariables>
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)
156 std::size_t embeddings_;
166template<
class GlobalPosition,
class SourceValues,
class I>
211template<
class TypeTag>
213 GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld>,
214 Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>>
221 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
222 using Element =
typename GridView::template Codim<0>::Entity;
224 static const int dimworld = GridView::dimensionworld;
225 using GlobalPosition =
typename Dune::FieldVector<typename GridView::ctype, dimworld>;
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)>;
240 ValueFunction valueFunction)
241 :
ParentType(pos, SourceValues(0.0)), valueFunction_(valueFunction) {}
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); }
268 ValueFunction valueFunction_;
281 template<
class Gr
idGeometry,
class Po
intSource,
class Po
intSourceMap>
283 const std::vector<PointSource>& sources,
284 PointSourceMap& pointSourceMap,
285 const std::string& paramGroup =
"")
287 const auto& boundingBoxTree = gridGeometry.boundingBoxTree();
289 for (
const auto& s : sources)
295 if (entities.empty())
302 source.setEmbeddings(entities.size()*source.embeddings());
308 auto fvGeometry =
localView(gridGeometry);
309 for (
const auto eIdx : entities)
312 const auto element = boundingBoxTree.entitySet().entity(eIdx);
313 fvGeometry.bindElement(element);
315 const auto globalPos = source.position();
317 constexpr int dim = GridGeometry::GridView::dimension;
318 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
319 for (
const auto& scv : scvs(fvGeometry))
321 scvIndices.push_back(scv.indexInElement());
325 for (
const auto scvIdx : scvIndices)
327 const auto key = std::make_pair(eIdx, scvIdx);
328 if (pointSourceMap.count(key))
329 pointSourceMap.at(key).push_back(source);
331 pointSourceMap.insert({key, {source}});
333 auto& s = pointSourceMap.at(key).back();
334 s.setEmbeddings(scvIndices.size()*s.embeddings());
340 for (
const auto eIdx : entities)
343 const auto key = std::make_pair(eIdx, 0);
344 if (pointSourceMap.count(key))
345 pointSourceMap.at(key).push_back(source);
347 pointSourceMap.insert({key, {source}});
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.
Algorithms that finds which geometric entities intersect.
Detect if a point intersects a geometry.
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:40
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:114
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr FCDiamond fcdiamond
Definition: method.hh:141
constexpr Box box
Definition: method.hh:136
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 ¶mGroup="")
calculate a DOF index to point source map from given vector of point sources
Definition: pointsource.hh:282
Declares all properties used in Dumux.