3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
27#ifndef DUMUX_INTEGRATION_POINTSOURCE_HH
28#define DUMUX_INTEGRATION_POINTSOURCE_HH
29
30#include <type_traits>
31#include <dune/common/reservedvector.hh>
35
36namespace Dumux {
37
43template<class GlobalPosition, class SourceValues, typename IdType = std::size_t>
44class IntegrationPointSource : public IdPointSource<GlobalPosition, SourceValues, IdType>
45{
47 using Scalar = std::decay_t<decltype(std::declval<SourceValues>()[0])>;
48
49public:
52 Scalar qpweight, Scalar integrationElement,
53 std::size_t elementIndex)
54 : ParentType(pos, values, id)
55 , qpweight_(qpweight)
56 , integrationElement_(integrationElement)
57 , elementIndex_(elementIndex)
58 {}
59
63 Scalar qpweight, Scalar integrationElement,
64 std::size_t elementIndex)
65 : ParentType(pos, id)
66 , qpweight_(qpweight)
67 , integrationElement_(integrationElement)
68 , elementIndex_(elementIndex)
69 {}
70
71 Scalar quadratureWeight() const
72 { return qpweight_; }
73
74 Scalar integrationElement() const
75 { return integrationElement_; }
76
77 void setQuadratureWeight(const Scalar qpWeight)
78 { return qpweight_ = qpWeight; }
79
80 void setIntegrationElement(const Scalar ie)
81 { integrationElement_ = ie; }
82
84 std::size_t elementIndex() const
85 { return elementIndex_; }
86
89 {
91 return *this;
92 }
93
96 {
98 return *this;
99 }
100
101private:
102 Scalar qpweight_;
103 Scalar integrationElement_;
104 std::size_t elementIndex_;
105};
106
112{
113
114public:
116 template<class GridGeometry, class PointSource, class PointSourceMap>
117 static void computePointSourceMap(const GridGeometry& gridGeometry,
118 const std::vector<PointSource>& sources,
119 PointSourceMap& pointSourceMap,
120 const std::string& paramGroup = "")
121 {
122 for (const auto& source : sources)
123 {
124 // get the index of the element in which the point source falls
125 const auto eIdx = source.elementIndex();
126 if constexpr (GridGeometry::discMethod == DiscretizationMethods::box)
127 {
128 auto fvGeometry = localView(gridGeometry);
129 // check in which subcontrolvolume(s) we are
130 const auto element = gridGeometry.boundingBoxTree().entitySet().entity(eIdx);
131 fvGeometry.bindElement(element);
132 const auto globalPos = source.position();
133
134 static const bool boxPointSourceLumping = getParamFromGroup<bool>(paramGroup, "PointSource.EnableBoxLumping", true);
135 if (boxPointSourceLumping)
136 {
137 // loop over all sub control volumes and check if the point source is inside
138 constexpr int dim = GridGeometry::GridView::dimension;
139 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
140 for (const auto& scv : scvs(fvGeometry))
141 if (intersectsPointGeometry(globalPos, scv.geometry()))
142 scvIndices.push_back(scv.indexInElement());
143
144 // for all scvs that where tested positiv add the point sources
145 // to the element/scv to point source map
146 for (auto scvIdx : scvIndices)
147 {
148 const auto key = std::make_pair(eIdx, scvIdx);
149 if (pointSourceMap.count(key))
150 pointSourceMap.at(key).push_back(source);
151 else
152 pointSourceMap.insert({key, {source}});
153 // split equally on the number of matched scvs
154 auto& s = pointSourceMap.at(key).back();
155 s.setEmbeddings(scvIndices.size()*s.embeddings());
156 }
157 }
158
159 // distribute the sources using the basis function weights
160 else
161 {
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))
168 {
169 const auto key = std::make_pair(eIdx, scv.indexInElement());
170 if (pointSourceMap.count(key))
171 pointSourceMap.at(key).push_back(source);
172 else
173 pointSourceMap.insert({key, {source}});
174
175 // adjust the integration element
176 auto& s = pointSourceMap.at(key).back();
177 s.setIntegrationElement(shapeValues[scv.indexInElement()]*s.integrationElement());
178 }
179 }
180 }
181 else
182 {
183 // add the pointsource to the DOF map
184 const auto key = std::make_pair(eIdx, /*scvIdx=*/ 0);
185 if (pointSourceMap.count(key))
186 pointSourceMap.at(key).push_back(source);
187 else
188 pointSourceMap.insert({key, {source}});
189 }
190 }
191 }
192};
193
194} // end namespace Dumux
195
196#endif
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
Definition: adapt.hh:29
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 &paramGroup="")
calculate a DOF index to point source map from given vector of point sources
Definition: integrationpointsource.hh:117