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])>;
51 [[deprecated(
"Call constructor with a single element index instead! Will be removed after 3.3")]]
61 [[deprecated(
"Call constructor with a single element index instead! Will be removed after 3.3")]]
94 {
return integrationElement_; }
97 {
return qpweight_ = qpWeight; }
100 { integrationElement_ = ie; }
102 [[deprecated(
"Use elementIndex() instead. Will be removed after 3.3")]]
104 {
return std::vector<std::size_t>({elementIndex_}); }
108 {
return elementIndex_; }
126 Scalar integrationElement_;
127 std::size_t elementIndex_;
139 template<
class Gr
idGeometry,
class Po
intSource,
class Po
intSourceMap>
141 const std::vector<PointSource>& sources,
142 PointSourceMap& pointSourceMap,
143 const std::string& paramGroup =
"")
145 for (
const auto& source : sources)
148 const auto eIdx = source.elementIndex();
152 const auto element = gridGeometry.boundingBoxTree().entitySet().entity(eIdx);
153 auto fvGeometry =
localView(gridGeometry);
154 fvGeometry.bindElement(element);
155 const auto globalPos = source.position();
157 static const bool boxPointSourceLumping = getParamFromGroup<bool>(paramGroup,
"PointSource.EnableBoxLumping",
true);
158 if (boxPointSourceLumping)
161 constexpr int dim = GridGeometry::GridView::dimension;
162 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
163 for (
const auto& scv : scvs(fvGeometry))
165 scvIndices.push_back(scv.indexInElement());
169 for (
auto scvIdx : scvIndices)
171 const auto key = std::make_pair(eIdx, scvIdx);
172 if (pointSourceMap.count(key))
173 pointSourceMap.at(key).push_back(source);
175 pointSourceMap.insert({key, {source}});
177 auto& s = pointSourceMap.at(key).back();
178 s.setEmbeddings(scvIndices.size()*s.embeddings());
185 const auto& localBasis = fvGeometry.feLocalBasis();
186 const auto ipLocal = element.geometry().local(globalPos);
187 using Scalar = std::decay_t<
decltype(source.values()[0])>;
188 std::vector<typename Dune::FieldVector<Scalar, 1>> shapeValues;
189 localBasis.evaluateFunction(ipLocal, shapeValues);
190 for (
const auto& scv : scvs(fvGeometry))
192 const auto key = std::make_pair(eIdx, scv.indexInElement());
193 if (pointSourceMap.count(key))
194 pointSourceMap.at(key).push_back(source);
196 pointSourceMap.insert({key, {source}});
199 auto& s = pointSourceMap.at(key).back();
200 s.setIntegrationElement(shapeValues[scv.indexInElement()]*s.integrationElement());
207 const auto key = std::make_pair(eIdx, 0);
208 if (pointSourceMap.count(key))
209 pointSourceMap.at(key).push_back(source);
211 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: geometry/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
Values values() const
return the source values
Definition: pointsource.hh:111
GlobalPosition GlobalPosition
Export the position type.
Definition: pointsource.hh:55
A point source class with an identifier to attach data.
Definition: pointsource.hh:167
std::size_t IdType
export the id type
Definition: pointsource.hh:173
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
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, const std::vector< std::size_t > &elementIndices)
Constructor for integration point sources.
Definition: integrationpointsource.hh:52
IntegrationPointSource(GlobalPosition pos, SourceValues values, IdType id, Scalar qpweight, Scalar integrationElement, std::size_t elementIndex)
Constructor for integration point sources.
Definition: integrationpointsource.hh:70
std::size_t elementIndex() const
The index of the element this intergration point source is associated with.
Definition: integrationpointsource.hh:107
Scalar quadratureWeight() const
Definition: integrationpointsource.hh:90
std::vector< std::size_t > elementIndices() const
Definition: integrationpointsource.hh:103
IntegrationPointSource(GlobalPosition pos, IdType id, Scalar qpweight, Scalar integrationElement, const std::vector< std::size_t > &elementIndices)
Constructor for integration point sources, when there is no.
Definition: integrationpointsource.hh:62
IntegrationPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: integrationpointsource.hh:111
void setQuadratureWeight(const Scalar qpWeight)
Definition: integrationpointsource.hh:96
void setIntegrationElement(const Scalar ie)
Definition: integrationpointsource.hh:99
Scalar integrationElement() const
Definition: integrationpointsource.hh:93
IntegrationPointSource(GlobalPosition pos, IdType id, Scalar qpweight, Scalar integrationElement, std::size_t elementIndex)
Definition: integrationpointsource.hh:81
A helper class calculating a DOF-index to point source map.
Definition: integrationpointsource.hh:135
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:140