26#ifndef DUMUX_POINTSOURCE_HH
27#define DUMUX_POINTSOURCE_HH
31#include <dune/common/reservedvector.hh>
48template<
class PositionType,
class ValueType>
53 using Scalar = std::decay_t<decltype(std::declval<ValueType>()[0])>;
66 :
values_(0.0), pos_(pos), embeddings_(1) {}
121 template<
class Problem,
class FVElementGeometry,
class ElementVolumeVariables>
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)
155 std::size_t embeddings_;
165template<
class GlobalPosition,
class SourceValues,
class I>
210template<
class TypeTag>
212 GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld>,
213 GetPropType<TypeTag, Properties::NumEqVector>>
220 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
221 using Element =
typename GridView::template Codim<0>::Entity;
223 static const int dimworld = GridView::dimensionworld;
224 using GlobalPosition =
typename Dune::FieldVector<typename GridView::ctype, dimworld>;
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)>;
239 ValueFunction valueFunction)
240 :
ParentType(pos, SourceValues(0.0)), valueFunction_(valueFunction) {}
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); }
267 ValueFunction valueFunction_;
280 template<
class Gr
idGeometry,
class Po
intSource,
class Po
intSourceMap>
282 const std::vector<PointSource>& sources,
283 PointSourceMap& pointSourceMap,
284 const std::string& paramGroup =
"")
286 const auto& boundingBoxTree = gridGeometry.boundingBoxTree();
288 for (
const auto& s : sources)
294 if (entities.empty())
301 source.setEmbeddings(entities.size()*source.embeddings());
303 for (
const auto eIdx : entities)
308 const auto element = boundingBoxTree.entitySet().entity(eIdx);
309 auto fvGeometry =
localView(gridGeometry);
310 fvGeometry.bindElement(element);
312 const auto globalPos = source.position();
314 constexpr int dim = GridGeometry::GridView::dimension;
315 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
316 for (
const auto& scv : scvs(fvGeometry))
318 scvIndices.push_back(scv.indexInElement());
322 for (
const auto scvIdx : scvIndices)
324 const auto key = std::make_pair(eIdx, scvIdx);
325 if (pointSourceMap.count(key))
326 pointSourceMap.at(key).push_back(source);
328 pointSourceMap.insert({key, {source}});
330 auto& s = pointSourceMap.at(key).back();
331 s.setEmbeddings(scvIndices.size()*s.embeddings());
337 const auto key = std::make_pair(eIdx, 0);
338 if (pointSourceMap.count(key))
339 pointSourceMap.at(key).push_back(source);
341 pointSourceMap.insert({key, {source}});
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
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 ¶mGroup="")
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.