version 3.8
integrationpointsource.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
15#ifndef DUMUX_INTEGRATION_POINTSOURCE_HH
16#define DUMUX_INTEGRATION_POINTSOURCE_HH
17
18#include <type_traits>
19#include <dune/common/reservedvector.hh>
23
24namespace Dumux {
25
31template<class GlobalPosition, class SourceValues, typename IdType = std::size_t>
32class IntegrationPointSource : public IdPointSource<GlobalPosition, SourceValues, IdType>
33{
35 using Scalar = std::decay_t<decltype(std::declval<SourceValues>()[0])>;
36
37public:
40 Scalar qpweight, Scalar integrationElement,
41 std::size_t elementIndex)
42 : ParentType(pos, values, id)
43 , qpweight_(qpweight)
44 , integrationElement_(integrationElement)
45 , elementIndex_(elementIndex)
46 {}
47
51 Scalar qpweight, Scalar integrationElement,
52 std::size_t elementIndex)
53 : ParentType(pos, id)
54 , qpweight_(qpweight)
55 , integrationElement_(integrationElement)
56 , elementIndex_(elementIndex)
57 {}
58
59 Scalar quadratureWeight() const
60 { return qpweight_; }
61
62 Scalar integrationElement() const
63 { return integrationElement_; }
64
65 void setQuadratureWeight(const Scalar qpWeight)
66 { return qpweight_ = qpWeight; }
67
68 void setIntegrationElement(const Scalar ie)
69 { integrationElement_ = ie; }
70
72 std::size_t elementIndex() const
73 { return elementIndex_; }
74
77 {
79 return *this;
80 }
81
84 {
86 return *this;
87 }
88
89private:
90 Scalar qpweight_;
91 Scalar integrationElement_;
92 std::size_t elementIndex_;
93};
94
100{
101
102public:
104 template<class GridGeometry, class PointSource, class PointSourceMap>
105 static void computePointSourceMap(const GridGeometry& gridGeometry,
106 const std::vector<PointSource>& sources,
107 PointSourceMap& pointSourceMap,
108 const std::string& paramGroup = "")
109 {
110 for (const auto& source : sources)
111 {
112 // get the index of the element in which the point source falls
113 const auto eIdx = source.elementIndex();
114 if constexpr (GridGeometry::discMethod == DiscretizationMethods::box
115 || GridGeometry::discMethod == DiscretizationMethods::fcdiamond)
116 {
117 auto fvGeometry = localView(gridGeometry);
118 // check in which subcontrolvolume(s) we are
119 const auto element = gridGeometry.boundingBoxTree().entitySet().entity(eIdx);
120 fvGeometry.bindElement(element);
121 const auto globalPos = source.position();
122
123 static const bool boxPointSourceLumping = getParamFromGroup<bool>(paramGroup, "PointSource.EnableBoxLumping", true);
124 if (boxPointSourceLumping)
125 {
126 // loop over all sub control volumes and check if the point source is inside
127 constexpr int dim = GridGeometry::GridView::dimension;
128 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
129 for (const auto& scv : scvs(fvGeometry))
130 if (intersectsPointGeometry(globalPos, fvGeometry.geometry(scv)))
131 scvIndices.push_back(scv.indexInElement());
132
133 // for all scvs that where tested positive add the point sources
134 // to the element/scv to point source map
135 for (auto scvIdx : scvIndices)
136 {
137 const auto key = std::make_pair(eIdx, scvIdx);
138 if (pointSourceMap.count(key))
139 pointSourceMap.at(key).push_back(source);
140 else
141 pointSourceMap.insert({key, {source}});
142 // split equally on the number of matched scvs
143 auto& s = pointSourceMap.at(key).back();
144 s.setEmbeddings(scvIndices.size()*s.embeddings());
145 }
146 }
147
148 // distribute the sources using the basis function weights
149 else
150 {
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))
157 {
158 const auto key = std::make_pair(eIdx, scv.indexInElement());
159 if (pointSourceMap.count(key))
160 pointSourceMap.at(key).push_back(source);
161 else
162 pointSourceMap.insert({key, {source}});
163
164 // adjust the integration element
165 auto& s = pointSourceMap.at(key).back();
166 s.setIntegrationElement(shapeValues[scv.indexInElement()]*s.integrationElement());
167 }
168 }
169 }
170 else
171 {
172 // add the pointsource to the DOF map
173 const auto key = std::make_pair(eIdx, /*scvIdx=*/ 0);
174 if (pointSourceMap.count(key))
175 pointSourceMap.at(key).push_back(source);
176 else
177 pointSourceMap.insert({key, {source}});
178 }
179 }
180 }
181};
182
183} // end namespace Dumux
184
185#endif
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 &paramGroup="")
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
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A point source class, i.e. sources located at a single point in space.