15#ifndef DUMUX_INTEGRATION_POINTSOURCE_HH
16#define DUMUX_INTEGRATION_POINTSOURCE_HH
19#include <dune/common/reservedvector.hh>
31template<
class GlobalPosition,
class SourceValues,
typename IdType = std::
size_t>
35 using Scalar = std::decay_t<decltype(std::declval<SourceValues>()[0])>;
63 {
return integrationElement_; }
66 {
return qpweight_ = qpWeight; }
69 { integrationElement_ = ie; }
73 {
return elementIndex_; }
91 Scalar integrationElement_;
92 std::size_t elementIndex_;
104 template<
class Gr
idGeometry,
class Po
intSource,
class Po
intSourceMap>
106 const std::vector<PointSource>& sources,
107 PointSourceMap& pointSourceMap,
108 const std::string& paramGroup =
"")
110 for (
const auto& source : sources)
113 const auto eIdx = source.elementIndex();
117 auto fvGeometry =
localView(gridGeometry);
119 const auto element = gridGeometry.boundingBoxTree().entitySet().entity(eIdx);
120 fvGeometry.bindElement(element);
121 const auto globalPos = source.position();
123 static const bool boxPointSourceLumping = getParamFromGroup<bool>(paramGroup,
"PointSource.EnableBoxLumping",
true);
124 if (boxPointSourceLumping)
127 constexpr int dim = GridGeometry::GridView::dimension;
128 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
129 for (
const auto& scv : scvs(fvGeometry))
131 scvIndices.push_back(scv.indexInElement());
135 for (
auto scvIdx : scvIndices)
137 const auto key = std::make_pair(eIdx, scvIdx);
138 if (pointSourceMap.count(key))
139 pointSourceMap.at(key).push_back(source);
141 pointSourceMap.insert({key, {source}});
143 auto& s = pointSourceMap.at(key).back();
144 s.setEmbeddings(scvIndices.size()*s.embeddings());
151 const auto& localBasis = fvGeometry.feLocalBasis();
152 const auto ipLocal = element.geometry().local(globalPos);
153 using Scalar = std::decay_t<
decltype(source.values()[0])>;
154 std::vector<typename Dune::FieldVector<Scalar, 1>> shapeValues;
155 localBasis.evaluateFunction(ipLocal, shapeValues);
156 for (
const auto& scv : scvs(fvGeometry))
158 const auto key = std::make_pair(eIdx, scv.indexInElement());
159 if (pointSourceMap.count(key))
160 pointSourceMap.at(key).push_back(source);
162 pointSourceMap.insert({key, {source}});
165 auto& s = pointSourceMap.at(key).back();
166 s.setIntegrationElement(shapeValues[scv.indexInElement()]*s.integrationElement());
173 const auto key = std::make_pair(eIdx, 0);
174 if (pointSourceMap.count(key))
175 pointSourceMap.at(key).push_back(source);
177 pointSourceMap.insert({key, {source}});
A point source class with an identifier to attach data.
Definition: pointsource.hh:156
std::size_t IdType
export the id type
Definition: pointsource.hh:162
IdType id() const
return the sources identifier
Definition: pointsource.hh:174
IdPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:178
A helper class calculating a DOF-index to point source map.
Definition: integrationpointsource.hh:100
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: integrationpointsource.hh:105
An integration point source class with an identifier to attach data and a quadrature weight and integ...
Definition: integrationpointsource.hh:33
IntegrationPointSource(GlobalPosition pos, SourceValues values, IdType id, Scalar qpweight, Scalar integrationElement, std::size_t elementIndex)
Constructor for integration point sources.
Definition: integrationpointsource.hh:39
std::size_t elementIndex() const
The index of the element this integration point source is associated with.
Definition: integrationpointsource.hh:72
Scalar quadratureWeight() const
Definition: integrationpointsource.hh:59
IntegrationPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: integrationpointsource.hh:76
void setQuadratureWeight(const Scalar qpWeight)
Definition: integrationpointsource.hh:65
void setIntegrationElement(const Scalar ie)
Definition: integrationpointsource.hh:68
Scalar integrationElement() const
Definition: integrationpointsource.hh:62
IntegrationPointSource(GlobalPosition pos, IdType id, Scalar qpweight, Scalar integrationElement, std::size_t elementIndex)
Definition: integrationpointsource.hh:50
Values values() const
return the source values
Definition: pointsource.hh:100
GlobalPosition GlobalPosition
Export the position type.
Definition: pointsource.hh:44
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
The available discretization methods in Dumux.
constexpr FCDiamond fcdiamond
Definition: method.hh:152
constexpr Box box
Definition: method.hh:147
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A point source class, i.e. sources located at a single point in space.