3.2-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>
34
35namespace Dumux {
36
42template<class GlobalPosition, class SourceValues, typename IdType = std::size_t>
43class IntegrationPointSource : public IdPointSource<GlobalPosition, SourceValues, IdType>
44{
46 using Scalar = std::decay_t<decltype(std::declval<SourceValues>()[0])>;
47
48public:
51 Scalar qpweight, Scalar integrationElement,
52 const std::vector<std::size_t>& elementIndices)
53 : ParentType(pos, values, id),
54 qpweight_(qpweight), integrationElement_(integrationElement),
55 elementIndices_(elementIndices) {}
56
58 // value known at the time of initialization
60 Scalar qpweight, Scalar integrationElement,
61 const std::vector<std::size_t>& elementIndices)
62 : ParentType(pos, id),
63 qpweight_(qpweight), integrationElement_(integrationElement),
64 elementIndices_(elementIndices) {}
65
66
67 Scalar quadratureWeight() const
68 {
69 return qpweight_;
70 }
71
72 Scalar integrationElement() const
73 {
74 return integrationElement_;
75 }
76
77 const std::vector<std::size_t>& elementIndices() const
78 {
79 return elementIndices_;
80 }
81
84 {
86 return *this;
87 }
88
91 {
93 return *this;
94 }
95
96private:
97 Scalar qpweight_;
98 Scalar integrationElement_;
99 std::vector<std::size_t> elementIndices_;
100};
101
107{
108
109public:
111 template<class GridGeometry, class PointSource, class PointSourceMap>
112 static void computePointSourceMap(const GridGeometry& gridGeometry,
113 std::vector<PointSource>& sources,
114 PointSourceMap& pointSourceMap)
115 {
116 for (auto&& source : sources)
117 {
118 // compute in which elements the point source falls
119 const auto& entities = source.elementIndices();
120 // split the source values equally among all concerned entities
121 source.setEmbeddings(source.embeddings()*entities.size());
122 // loop over all intersected elements
123 for (unsigned int eIdx : entities)
124 {
125 if (GridGeometry::discMethod == DiscretizationMethod::box)
126 {
127 // check in which subcontrolvolume(s) we are
128 const auto element = gridGeometry.boundingBoxTree().entitySet().entity(eIdx);
129 auto fvGeometry = localView(gridGeometry);
130 fvGeometry.bindElement(element);
131 const auto globalPos = source.position();
132 // loop over all sub control volumes and check if the point source is inside
133 constexpr int dim = GridGeometry::GridView::dimension;
134 Dune::ReservedVector<std::size_t, 1<<dim> scvIndices;
135 for (auto&& scv : scvs(fvGeometry))
136 if (intersectsPointGeometry(globalPos, scv.geometry()))
137 scvIndices.push_back(scv.indexInElement());
138
139 // for all scvs that where tested positiv add the point sources
140 // to the element/scv to point source map
141 for (auto scvIdx : scvIndices)
142 {
143 const auto key = std::make_pair(eIdx, scvIdx);
144 if (pointSourceMap.count(key))
145 pointSourceMap.at(key).push_back(source);
146 else
147 pointSourceMap.insert({key, {source}});
148 // split equally on the number of matched scvs
149 auto& s = pointSourceMap.at(key).back();
150 s.setEmbeddings(scvIndices.size()*s.embeddings());
151 }
152 }
153 else
154 {
155 // add the pointsource to the DOF map
156 const auto key = std::make_pair(eIdx, /*scvIdx=*/ 0);
157 if (pointSourceMap.count(key))
158 pointSourceMap.at(key).push_back(source);
159 else
160 pointSourceMap.insert({key, {source}});
161 }
162 }
163 }
164 }
165};
166
167} // end namespace Dumux
168
169#endif
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
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
Definition: adapt.hh:29
Values values() const
return the source values
Definition: pointsource.hh:110
GlobalPosition GlobalPosition
Export the position type.
Definition: pointsource.hh:54
A point source class with an identifier to attach data.
Definition: pointsource.hh:166
std::size_t IdType
export the id type
Definition: pointsource.hh:172
IdType id() const
return the sources identifier
Definition: pointsource.hh:184
IdPointSource & operator=(const SourceValues &values)
Convenience = operator overload modifying only the values.
Definition: pointsource.hh:188
An integration point source class with an identifier to attach data and a quadrature weight and integ...
Definition: integrationpointsource.hh:44
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