27#ifndef DUMUX_INTEGRATION_POINTSOURCE_HH
28#define DUMUX_INTEGRATION_POINTSOURCE_HH
31#include <dune/common/reservedvector.hh>
42template<
class GlobalPosition,
class SourceValues,
typename IdType = std::
size_t>
46 using Scalar = std::decay_t<decltype(std::declval<SourceValues>()[0])>;
62 : ParentType(pos,
id),
74 return integrationElement_;
79 return elementIndices_;
98 Scalar integrationElement_;
99 std::vector<std::size_t> elementIndices_;
111 template<
class Gr
idGeometry,
class Po
intSource,
class Po
intSourceMap>
113 std::vector<PointSource>& sources,
114 PointSourceMap& pointSourceMap)
116 for (
auto&& source : sources)
119 const auto& entities = source.elementIndices();
121 source.setEmbeddings(source.embeddings()*entities.size());
123 for (
unsigned int eIdx : entities)
128 const auto element = gridGeometry.boundingBoxTree().entitySet().entity(eIdx);
129 auto fvGeometry =
localView(gridGeometry);
130 fvGeometry.bindElement(element);
131 const auto globalPos = source.position();
133 constexpr int dim = GridGeometry::GridView::dimension;
134 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
135 for (
auto&& scv : scvs(fvGeometry))
137 scvIndices.push_back(scv.indexInElement());
141 for (
auto scvIdx : scvIndices)
143 const auto key = std::make_pair(eIdx, scvIdx);
144 if (pointSourceMap.count(key))
145 pointSourceMap.at(key).push_back(source);
147 pointSourceMap.insert({key, {source}});
149 auto& s = pointSourceMap.at(key).back();
150 s.setEmbeddings(scvIndices.size()*s.embeddings());
156 const auto key = std::make_pair(eIdx, 0);
157 if (pointSourceMap.count(key))
158 pointSourceMap.at(key).push_back(source);
160 pointSourceMap.insert({key, {source}});
A point source class, i.e. sources located at a single point in space.
The available discretization methods in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition localview.hh:38
@ box
Definition method.hh:38
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
std::decay_t< decltype(std::declval< SourceValues >()[0])> Scalar
Definition pointsource.hh:52
Values values() const
Definition pointsource.hh:110
typename GridGeometry< i >::GlobalCoordinate GlobalPosition
Definition pointsource.hh:54
std::size_t IdType
Definition pointsource.hh:172
IdPointSource(GlobalPosition pos, SourceValues values, IdType id)
Definition pointsource.hh:175
std::size_t id() const
Definition pointsource.hh:184
IdPointSource & operator=(const SourceValues &values)
Definition pointsource.hh:188
Definition common/properties.hh:113
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:50
const std::vector< std::size_t > & elementIndices() const
Definition integrationpointsource.hh:77
Scalar quadratureWeight() const
Definition integrationpointsource.hh:67
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:59
IntegrationPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition integrationpointsource.hh:83
Scalar integrationElement() const
Definition integrationpointsource.hh:72
A helper class calculating a DOF-index to point source map.
Definition integrationpointsource.hh:107
static void computePointSourceMap(const GridGeometry &gridGeometry, std::vector< PointSource > &sources, PointSourceMap &pointSourceMap)
calculate a DOF index to point source map from given vector of point sources
Definition integrationpointsource.hh:112