27#ifndef DUMUX_INTEGRATION_POINTSOURCE_HH
28#define DUMUX_INTEGRATION_POINTSOURCE_HH
31#include <dune/common/reservedvector.hh>
43template<
class GlobalPosition,
class SourceValues,
typename IdType = std::
size_t>
47 using Scalar = std::decay_t<decltype(std::declval<SourceValues>()[0])>;
75 {
return integrationElement_; }
78 {
return qpweight_ = qpWeight; }
81 { integrationElement_ = ie; }
85 {
return elementIndex_; }
103 Scalar integrationElement_;
104 std::size_t elementIndex_;
116 template<
class Gr
idGeometry,
class Po
intSource,
class Po
intSourceMap>
118 const std::vector<PointSource>& sources,
119 PointSourceMap& pointSourceMap,
120 const std::string& paramGroup =
"")
122 for (
const auto& source : sources)
125 const auto eIdx = source.elementIndex();
128 auto fvGeometry =
localView(gridGeometry);
130 const auto element = gridGeometry.boundingBoxTree().entitySet().entity(eIdx);
131 fvGeometry.bindElement(element);
132 const auto globalPos = source.position();
134 static const bool boxPointSourceLumping = getParamFromGroup<bool>(paramGroup,
"PointSource.EnableBoxLumping",
true);
135 if (boxPointSourceLumping)
138 constexpr int dim = GridGeometry::GridView::dimension;
139 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
140 for (
const auto& scv : scvs(fvGeometry))
142 scvIndices.push_back(scv.indexInElement());
146 for (
auto scvIdx : scvIndices)
148 const auto key = std::make_pair(eIdx, scvIdx);
149 if (pointSourceMap.count(key))
150 pointSourceMap.at(key).push_back(source);
152 pointSourceMap.insert({key, {source}});
154 auto& s = pointSourceMap.at(key).back();
155 s.setEmbeddings(scvIndices.size()*s.embeddings());
162 const auto& localBasis = fvGeometry.feLocalBasis();
163 const auto ipLocal = element.geometry().local(globalPos);
164 using Scalar = std::decay_t<
decltype(source.values()[0])>;
165 std::vector<typename Dune::FieldVector<Scalar, 1>> shapeValues;
166 localBasis.evaluateFunction(ipLocal, shapeValues);
167 for (
const auto& scv : scvs(fvGeometry))
169 const auto key = std::make_pair(eIdx, scv.indexInElement());
170 if (pointSourceMap.count(key))
171 pointSourceMap.at(key).push_back(source);
173 pointSourceMap.insert({key, {source}});
176 auto& s = pointSourceMap.at(key).back();
177 s.setIntegrationElement(shapeValues[scv.indexInElement()]*s.integrationElement());
184 const auto key = std::make_pair(eIdx, 0);
185 if (pointSourceMap.count(key))
186 pointSourceMap.at(key).push_back(source);
188 pointSourceMap.insert({key, {source}});
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A point source class, i.e. sources located at a single point in space.
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: intersectspointgeometry.hh:38
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
constexpr Box box
Definition: method.hh:139
Values values() const
return the source values
Definition: pointsource.hh:112
GlobalPosition GlobalPosition
Export the position type.
Definition: pointsource.hh:56
A point source class with an identifier to attach data.
Definition: pointsource.hh:168
std::size_t IdType
export the id type
Definition: pointsource.hh:174
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
An integration point source class with an identifier to attach data and a quadrature weight and integ...
Definition: integrationpointsource.hh:45
IntegrationPointSource(GlobalPosition pos, SourceValues values, IdType id, Scalar qpweight, Scalar integrationElement, std::size_t elementIndex)
Constructor for integration point sources.
Definition: integrationpointsource.hh:51
std::size_t elementIndex() const
The index of the element this intergration point source is associated with.
Definition: integrationpointsource.hh:84
Scalar quadratureWeight() const
Definition: integrationpointsource.hh:71
IntegrationPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: integrationpointsource.hh:88
void setQuadratureWeight(const Scalar qpWeight)
Definition: integrationpointsource.hh:77
void setIntegrationElement(const Scalar ie)
Definition: integrationpointsource.hh:80
Scalar integrationElement() const
Definition: integrationpointsource.hh:74
IntegrationPointSource(GlobalPosition pos, IdType id, Scalar qpweight, Scalar integrationElement, std::size_t elementIndex)
Definition: integrationpointsource.hh:62
A helper class calculating a DOF-index to point source map.
Definition: integrationpointsource.hh:112
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:117