3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
cellcentered/tpfa/elementvolumevariables.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 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_CCTPFA_ELEMENT_VOLUMEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_CCTPFA_ELEMENT_VOLUMEVARIABLES_HH
26
27#include <algorithm>
28#include <type_traits>
29#include <vector>
30
32
33namespace Dumux {
34
42template<class GVV, bool cachingEnabled>
44{};
45
51template<class GVV>
52class CCTpfaElementVolumeVariables<GVV, /*cachingEnabled*/true>
53{
54public:
57
59 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
60
63 : gridVolVarsPtr_(&gridVolVars)
64 , numScv_(gridVolVars.problem().gridGeometry().numScv())
65 {}
66
68 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
69 const VolumeVariables& operator [](const SubControlVolume& scv) const
70 {
71 if (scv.dofIndex() < numScv_)
72 return gridVolVars().volVars(scv.dofIndex());
73 else
74 return boundaryVolumeVariables_[getLocalIdx_(scv.dofIndex())];
75 }
76
78 const VolumeVariables& operator [](const std::size_t scvIdx) const
79 {
80 if (scvIdx < numScv_)
81 return gridVolVars().volVars(scvIdx);
82 else
83 return boundaryVolumeVariables_[getLocalIdx_(scvIdx)];
84 }
85
87 template<class FVElementGeometry, class SolutionVector>
88 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
89 const FVElementGeometry& fvGeometry,
90 const SolutionVector& sol)
91 {
92 if (!fvGeometry.hasBoundaryScvf())
93 return;
94
95 clear_();
96 boundaryVolVarIndices_.reserve(fvGeometry.numScvf());
97 boundaryVolumeVariables_.reserve(fvGeometry.numScvf());
98
99 for (const auto& scvf : scvfs(fvGeometry))
100 {
101 if (!scvf.boundary())
102 continue;
103
104 // check if boundary is a pure dirichlet boundary
105 const auto& problem = gridVolVars().problem();
106 const auto bcTypes = problem.boundaryTypes(element, scvf);
107 if (bcTypes.hasOnlyDirichlet())
108 {
109 const auto dirichletPriVars = elementSolution<FVElementGeometry>(problem.dirichlet(element, scvf));
110 auto&& scvI = fvGeometry.scv(scvf.insideScvIdx());
111
112 VolumeVariables volVars;
113 volVars.update(dirichletPriVars,
114 problem,
115 element,
116 scvI);
117
118 boundaryVolumeVariables_.emplace_back(std::move(volVars));
119 boundaryVolVarIndices_.push_back(scvf.outsideScvIdx());
120 }
121 }
122 }
123
125 template<class FVElementGeometry, class SolutionVector>
126 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
127 const FVElementGeometry& fvGeometry,
128 const SolutionVector& sol)
129 {}
130
133 { return *gridVolVarsPtr_; }
134
135private:
137 void clear_()
138 {
139 boundaryVolVarIndices_.clear();
140 boundaryVolumeVariables_.clear();
141 }
142
143 const GridVolumeVariables* gridVolVarsPtr_;
144
146 int getLocalIdx_(const int volVarIdx) const
147 {
148 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
149 assert(it != boundaryVolVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
150 return std::distance(boundaryVolVarIndices_.begin(), it);
151 }
152
153 std::vector<std::size_t> boundaryVolVarIndices_;
154 std::vector<VolumeVariables> boundaryVolumeVariables_;
155 const std::size_t numScv_;
156};
157
162template<class GVV>
163class CCTpfaElementVolumeVariables<GVV, /*cachingEnabled*/false>
164{
165public:
168
170 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
171
174 : gridVolVarsPtr_(&gridVolVars) {}
175
177 template<class FVElementGeometry, class SolutionVector>
178 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
179 const FVElementGeometry& fvGeometry,
180 const SolutionVector& sol)
181 {
182 clear_();
183
184 const auto& problem = gridVolVars().problem();
185 const auto& gridGeometry = fvGeometry.gridGeometry();
186 const auto globalI = gridGeometry.elementMapper().index(element);
187 const auto& connectivityMapI = gridGeometry.connectivityMap()[globalI];
188 const auto numDofs = connectivityMapI.size() + 1;
189
190 // resize local containers to the required size (for internal elements)
191 volumeVariables_.resize(numDofs);
192 volVarIndices_.resize(numDofs);
193 int localIdx = 0;
194
195 // update the volume variables of the element at hand
196 auto&& scvI = fvGeometry.scv(globalI);
197 volumeVariables_[localIdx].update(elementSolution(element, sol, gridGeometry),
198 problem,
199 element,
200 scvI);
201 volVarIndices_[localIdx] = scvI.dofIndex();
202 ++localIdx;
203
204 // Update the volume variables of the neighboring elements
205 for (const auto& dataJ : connectivityMapI)
206 {
207 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
208 auto&& scvJ = fvGeometry.scv(dataJ.globalJ);
209 volumeVariables_[localIdx].update(elementSolution(elementJ, sol, gridGeometry),
210 problem,
211 elementJ,
212 scvJ);
213 volVarIndices_[localIdx] = scvJ.dofIndex();
214 ++localIdx;
215 }
216
217 if (fvGeometry.hasBoundaryScvf())
218 {
219 // Update boundary volume variables
220 for (auto&& scvf : scvfs(fvGeometry))
221 {
222 // if we are not on a boundary, skip to the next scvf
223 if (!scvf.boundary())
224 continue;
225
226 // check if boundary is a pure dirichlet boundary
227 const auto bcTypes = problem.boundaryTypes(element, scvf);
228 if (bcTypes.hasOnlyDirichlet())
229 {
230 const auto dirichletPriVars = elementSolution<FVElementGeometry>(problem.dirichlet(element, scvf));
231
232 volumeVariables_.resize(localIdx+1);
233 volVarIndices_.resize(localIdx+1);
234 volumeVariables_[localIdx].update(dirichletPriVars,
235 problem,
236 element,
237 scvI);
238 volVarIndices_[localIdx] = scvf.outsideScvIdx();
239 ++localIdx;
240 }
241 }
242 }
243
246 // const auto& additionalDofDependencies = problem.getAdditionalDofDependencies(globalI);
247 // if (!additionalDofDependencies.empty())
248 // {
249 // volumeVariables_.resize(volumeVariables_.size() + additionalDofDependencies.size());
250 // volVarIndices_.resize(volVarIndices_.size() + additionalDofDependencies.size());
251 // for (auto globalJ : additionalDofDependencies)
252 // {
253 // const auto& elementJ = gridGeometry.element(globalJ);
254 // auto&& scvJ = fvGeometry.scv(globalJ);
255
256 // volumeVariables_[localIdx].update(elementSolution(elementJ, sol, gridGeometry),
257 // problem,
258 // elementJ,
259 // scvJ);
260 // volVarIndices_[localIdx] = scvJ.dofIndex();
261 // ++localIdx;
262 // }
263 // }
264 }
265
266 template<class FVElementGeometry, class SolutionVector>
267 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
268 const FVElementGeometry& fvGeometry,
269 const SolutionVector& sol)
270 {
271 clear_();
272
273 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
274 volumeVariables_.resize(1);
275 volVarIndices_.resize(1);
276
277 // update the volume variables of the element
278 auto&& scv = fvGeometry.scv(eIdx);
279 volumeVariables_[0].update(elementSolution(element, sol, fvGeometry.gridGeometry()),
280 gridVolVars().problem(),
281 element,
282 scv);
283 volVarIndices_[0] = scv.dofIndex();
284 }
285
287 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
288 const VolumeVariables& operator [](const SubControlVolume& scv) const
289 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
290
292 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
293 VolumeVariables& operator [](const SubControlVolume& scv)
294 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
295
297 const VolumeVariables& operator [](std::size_t scvIdx) const
298 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
299
301 VolumeVariables& operator [](std::size_t scvIdx)
302 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
303
306 { return *gridVolVarsPtr_; }
307
308private:
310 void clear_()
311 {
312 volVarIndices_.clear();
313 volumeVariables_.clear();
314 }
315
316 const GridVolumeVariables* gridVolVarsPtr_;
317
319 int getLocalIdx_(const int volVarIdx) const
320 {
321 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
322 assert(it != volVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
323 return std::distance(volVarIndices_.begin(), it);
324 }
325
326 std::vector<std::size_t> volVarIndices_;
327 std::vector<VolumeVariables> volumeVariables_;
328};
329
330} // end namespace Dumux
331
332#endif
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:115
Definition: adapt.hh:29
The local (stencil) volume variables class for cell centered tpfa models.
Definition: cellcentered/tpfa/elementvolumevariables.hh:44
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol)
precompute the volume variables of an element - do nothing: volVars are cached
Definition: cellcentered/tpfa/elementvolumevariables.hh:126
CCTpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/tpfa/elementvolumevariables.hh:62
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/tpfa/elementvolumevariables.hh:132
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:56
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol)
precompute all boundary volume variables in a stencil of an element, the remaining ones are cached
Definition: cellcentered/tpfa/elementvolumevariables.hh:88
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:59
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol)
Prepares the volume variables within the element stencil.
Definition: cellcentered/tpfa/elementvolumevariables.hh:178
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:170
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol)
Definition: cellcentered/tpfa/elementvolumevariables.hh:267
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:167
CCTpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/tpfa/elementvolumevariables.hh:173
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/tpfa/elementvolumevariables.hh:305
The local element solution class for cell-centered methods.