14#ifndef DUMUX_POINTSOURCE_HH 
   15#define DUMUX_POINTSOURCE_HH 
   19#include <dune/common/exceptions.hh> 
   20#include <dune/common/reservedvector.hh> 
   38template<
class PositionType, 
class ValueType>
 
   43    using Scalar = std::decay_t<decltype(std::declval<ValueType>()[0])>;
 
   56      : 
values_(0.0), pos_(pos), embeddings_(1) {}
 
 
  111    template<
class Problem, 
class FVElementGeometry, 
class ElementVolumeVariables>
 
  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)
 
 
  145    std::size_t embeddings_; 
 
 
  155template<
class GlobalPosition, 
class SourceValues, 
class I>
 
  167      :  ParentType(pos, 
values), id_(
id) {}
 
 
  172      : ParentType(pos, SourceValues(0.0)), id_(
id) {}
 
 
 
  200template<
class TypeTag>
 
  202                                                   GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld>,
 
  203                                                   Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>>
 
  210    using SubControlVolume = 
typename FVElementGeometry::SubControlVolume;
 
  211    using Element = 
typename GridView::template Codim<0>::Entity;
 
  213    static const int dimworld = GridView::dimensionworld;
 
  214    using GlobalPosition = 
typename Dune::FieldVector<typename GridView::ctype, dimworld>;
 
  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)>;
 
  229                            ValueFunction valueFunction)
 
  230      : ParentType(pos, SourceValues(0.0)), valueFunction_(valueFunction) {}
 
 
  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); }
 
 
  257    ValueFunction valueFunction_;
 
 
  270    template<
class Gr
idGeometry, 
class Po
intSource, 
class Po
intSourceMap>
 
  272                                      const std::vector<PointSource>& sources,
 
  273                                      PointSourceMap& pointSourceMap,
 
  274                                      const std::string& paramGroup = 
"")
 
  276        const auto& boundingBoxTree = gridGeometry.boundingBoxTree();
 
  278        for (
const auto& s : sources)
 
  284            if (entities.empty())
 
  291            source.setEmbeddings(entities.size()*source.embeddings());
 
  297                auto fvGeometry = 
localView(gridGeometry);
 
  298                for (
const auto eIdx : entities)
 
  301                    const auto element = boundingBoxTree.entitySet().entity(eIdx);
 
  302                    fvGeometry.bindElement(element);
 
  304                    const auto globalPos = source.position();
 
  306                    constexpr int dim = GridGeometry::GridView::dimension;
 
  307                    Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
 
  308                    for (
const auto& scv : scvs(fvGeometry))
 
  310                            scvIndices.push_back(scv.indexInElement());
 
  314                    for (
const auto scvIdx : scvIndices)
 
  316                        const auto key = std::make_pair(eIdx, scvIdx);
 
  317                        if (pointSourceMap.count(key))
 
  318                            pointSourceMap.at(key).push_back(source);
 
  320                            pointSourceMap.insert({key, {source}});
 
  322                        auto& s = pointSourceMap.at(key).back();
 
  323                        s.setEmbeddings(scvIndices.size()*s.embeddings());
 
  330                for (
const auto eIdx : entities)
 
  333                    const auto key = std::make_pair(eIdx,  0);
 
  334                    if (pointSourceMap.count(key))
 
  335                        pointSourceMap.at(key).push_back(source);
 
  337                        pointSourceMap.insert({key, {source}});
 
  341                DUNE_THROW(Dune::NotImplemented, 
"BoundingBoxTreePointSourceHelper for discretization method " << GridGeometry::discMethod);
 
 
 
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 ¶mGroup="")
calculate a DOF index to point source map from given vector of point sources
Definition pointsource.hh:271
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
Definition pointsource.hh:175
IdPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition pointsource.hh:179
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
Definition pointsource.hh:101
PositionType GlobalPosition
Export the position type.
Definition pointsource.hh:45
Values values_
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
Definition pointsource.hh:136
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
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.