3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
cellcentered/mpfa/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_CCMPFA_ELEMENT_VOLUMEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_VOLUMEVARIABLES_HH
26
27#include <algorithm>
28#include <type_traits>
29#include <utility>
30#include <vector>
31
33
34namespace Dumux {
35namespace CCMpfa {
36
48 template<class FVElementGeometry>
49 std::size_t maxNumBoundaryVolVars(const FVElementGeometry& fvGeometry)
50 {
51 const auto& gridGeometry = fvGeometry.gridGeometry();
52 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
53
54 std::size_t numBoundaryVolVars = 0;
55 for (const auto& scvf : scvfs(fvGeometry))
56 {
57 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
58 numBoundaryVolVars += gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
59 else
60 numBoundaryVolVars += gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
61 }
62
63 return numBoundaryVolVars;
64 }
65
80 template<class VolumeVariables, class IndexType, class Problem, class FVElemGeom, class NodalIndexSet>
81 void addBoundaryVolVarsAtNode(std::vector<VolumeVariables>& volVars,
82 std::vector<IndexType>& volVarIndices,
83 const Problem& problem,
84 const typename FVElemGeom::GridGeometry::GridView::template Codim<0>::Entity& element,
85 const FVElemGeom& fvGeometry,
86 const NodalIndexSet& nodalIndexSet)
87 {
88 if (nodalIndexSet.numBoundaryScvfs() == 0)
89 return;
90
91 // index of the element the fvGeometry was bound to
92 const auto boundElemIdx = fvGeometry.gridGeometry().elementMapper().index(element);
93
94 // check each scvf in the index set for boundary presence
95 for (auto scvfIdx : nodalIndexSet.gridScvfIndices())
96 {
97 const auto& ivScvf = fvGeometry.scvf(scvfIdx);
98
99 // only proceed for scvfs on the boundary and not in the bound element
100 if (!ivScvf.boundary() || ivScvf.insideScvIdx() == boundElemIdx)
101 continue;
102
103 const auto insideScvIdx = ivScvf.insideScvIdx();
104 const auto insideElement = fvGeometry.gridGeometry().element(insideScvIdx);
105 const auto bcTypes = problem.boundaryTypes(insideElement, ivScvf);
106
107 // Only proceed on dirichlet boundaries. On Neumann
108 // boundaries the "outside" vol vars cannot be properly defined.
109 if (bcTypes.hasOnlyDirichlet())
110 {
111 VolumeVariables dirichletVolVars;
112 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(insideElement, ivScvf)),
113 problem,
114 insideElement,
115 fvGeometry.scv(insideScvIdx));
116
117 volVars.emplace_back(std::move(dirichletVolVars));
118 volVarIndices.push_back(ivScvf.outsideScvIdx());
119 }
120 }
121 }
122
134 template<class VolumeVariables, class IndexType, class Problem, class FVElemGeom>
135 void addBoundaryVolVars(std::vector<VolumeVariables>& volVars,
136 std::vector<IndexType>& volVarIndices,
137 const Problem& problem,
138 const typename FVElemGeom::GridGeometry::GridView::template Codim<0>::Entity& element,
139 const FVElemGeom& fvGeometry)
140 {
141 const auto& gridGeometry = fvGeometry.gridGeometry();
142
143 // treat the BCs inside the element
144 if (fvGeometry.hasBoundaryScvf())
145 {
146 const auto boundElemIdx = gridGeometry.elementMapper().index(element);
147 const auto& scvI = fvGeometry.scv(boundElemIdx);
148
149 for (const auto& scvf : scvfs(fvGeometry))
150 {
151 if (!scvf.boundary())
152 continue;
153
154 // Only proceed on dirichlet boundaries. On Neumann
155 // boundaries the "outside" vol vars cannot be properly defined.
156 if (problem.boundaryTypes(element, scvf).hasOnlyDirichlet())
157 {
158 VolumeVariables dirichletVolVars;
159 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(element, scvf)),
160 problem,
161 element,
162 scvI);
163
164 volVars.emplace_back(std::move(dirichletVolVars));
165 volVarIndices.push_back(scvf.outsideScvIdx());
166 }
167 }
168 }
169
170 // Update boundary volume variables in the neighbors
171 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
172 for (const auto& scvf : scvfs(fvGeometry))
173 {
174 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
175 addBoundaryVolVarsAtNode( volVars, volVarIndices, problem, element, fvGeometry,
176 gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet() );
177 else
178 addBoundaryVolVarsAtNode( volVars, volVarIndices, problem, element, fvGeometry,
179 gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet() );
180 }
181 }
182} // end namespace CCMpfa
183
191template<class GVV, bool cachingEnabled>
193
199template<class GVV>
200class CCMpfaElementVolumeVariables<GVV, /*cachingEnabled*/true>
201{
202public:
205
207 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
208
211 : gridVolVarsPtr_(&gridVolVars)
212 , numScv_(gridVolVars.problem().gridGeometry().numScv())
213 {}
214
216 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
217 const VolumeVariables& operator [](const SubControlVolume& scv) const
218 {
219 return scv.dofIndex() < numScv_ ? gridVolVars().volVars(scv.dofIndex())
220 : boundaryVolVars_[getLocalIdx_(scv.dofIndex())];
221 }
222
224 const VolumeVariables& operator [](const std::size_t scvIdx) const
225 {
226 return scvIdx < numScv_ ? gridVolVars().volVars(scvIdx)
227 : boundaryVolVars_[getLocalIdx_(scvIdx)];
228 }
229
231 template<class FVElementGeometry, class SolutionVector>
232 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
233 const FVElementGeometry& fvGeometry,
234 const SolutionVector& sol)
235 {
236 clear();
237
238 // maybe prepare boundary volume variables
240 if (maxNumBoundaryVolVars > 0)
241 {
242 boundaryVolVars_.reserve(maxNumBoundaryVolVars);
243 boundaryVolVarIndices_.reserve(maxNumBoundaryVolVars);
244 CCMpfa::addBoundaryVolVars(boundaryVolVars_, boundaryVolVarIndices_, gridVolVars().problem(), element, fvGeometry);
245 }
246 }
247
249 template<class FVElementGeometry, class SolutionVector>
250 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
251 const FVElementGeometry& fvGeometry,
252 const SolutionVector& sol)
253 {}
254
256 void clear()
257 {
258 boundaryVolVarIndices_.clear();
259 boundaryVolVars_.clear();
260 }
261
264 { return *gridVolVarsPtr_; }
265
266private:
268 int getLocalIdx_(const int volVarIdx) const
269 {
270 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
271 assert(it != boundaryVolVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
272 return std::distance(boundaryVolVarIndices_.begin(), it);
273 }
274
275 const GridVolumeVariables* gridVolVarsPtr_;
276
277 std::size_t numScv_;
278 std::vector<std::size_t> boundaryVolVarIndices_;
279 std::vector<VolumeVariables> boundaryVolVars_;
280};
281
282
287template<class GVV>
288class CCMpfaElementVolumeVariables<GVV, /*cachingEnabled*/false>
289{
290public:
293
295 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
296
299 : gridVolVarsPtr_(&gridVolVars) {}
300
302 template<class FVElementGeometry, class SolutionVector>
303 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
304 const FVElementGeometry& fvGeometry,
305 const SolutionVector& sol)
306 {
307 clear();
308
309 const auto& problem = gridVolVars().problem();
310 const auto& gridGeometry = fvGeometry.gridGeometry();
311
312 // stencil information
313 const auto globalI = gridGeometry.elementMapper().index(element);
314 const auto& assemblyMapI = gridGeometry.connectivityMap()[globalI];
315 const auto numVolVars = assemblyMapI.size() + 1;
316
317 // resize local containers to the required size (for internal elements)
319 volumeVariables_.reserve(numVolVars+maxNumBoundaryVolVars);
320 volVarIndices_.reserve(numVolVars+maxNumBoundaryVolVars);
321
322 VolumeVariables volVars;
323 const auto& scvI = fvGeometry.scv(globalI);
324 volVars.update(elementSolution(element, sol, gridGeometry),
325 problem,
326 element,
327 scvI);
328
329 volVarIndices_.push_back(scvI.dofIndex());
330 volumeVariables_.emplace_back(std::move(volVars));
331
332 // Update the volume variables of the neighboring elements
333 for (auto&& dataJ : assemblyMapI)
334 {
335 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
336 const auto& scvJ = fvGeometry.scv(dataJ.globalJ);
337 VolumeVariables volVarsJ;
338 volVarsJ.update(elementSolution(elementJ, sol, gridGeometry),
339 problem,
340 elementJ,
341 scvJ);
342
343 volVarIndices_.push_back(scvJ.dofIndex());
344 volumeVariables_.emplace_back(std::move(volVarsJ));
345 }
346
347 // maybe prepare boundary volume variables
348 if (maxNumBoundaryVolVars > 0)
349 CCMpfa::addBoundaryVolVars(volumeVariables_, volVarIndices_, problem, element, fvGeometry);
350
351 // //! TODO Check if user added additional DOF dependencies, i.e. the residual of DOF globalI depends
352 // //! on additional DOFs not included in the discretization schemes' occupation pattern
353 // const auto& additionalDofDependencies = problem.getAdditionalDofDependencies(globalI);
354 // if (!additionalDofDependencies.empty())
355 // {
356 // volumeVariables_.reserve(volumeVariables_.size() + additionalDofDependencies.size());
357 // volVarIndices_.reserve(volVarIndices_.size() + additionalDofDependencies.size());
358 // for (auto globalJ : additionalDofDependencies)
359 // {
360 // const auto& elementJ = gridGeometry.element(globalJ);
361 // const auto& scvJ = fvGeometry.scv(globalJ);
362
363 // VolumeVariables additionalVolVars;
364 // additionalVolVars.update(elementSolution(elementJ, sol, gridGeometry),
365 // problem,
366 // elementJ,
367 // scvJ);
368
369 // volumeVariables_.emplace_back(std::move(additionalVolVars));
370 // volVarIndices_.push_back(globalJ);
371 // }
372 // }
373 }
374
376 template<class FVElementGeometry, class SolutionVector>
377 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
378 const FVElementGeometry& fvGeometry,
379 const SolutionVector& sol)
380 {
381 clear();
382
383 const auto& gridGeometry = fvGeometry.gridGeometry();
384 auto eIdx = gridGeometry.elementMapper().index(element);
385 volumeVariables_.resize(1);
386 volVarIndices_.resize(1);
387
388 // update the volume variables of the element
389 const auto& scv = fvGeometry.scv(eIdx);
390 volumeVariables_[0].update(elementSolution(element, sol, gridGeometry),
391 gridVolVars().problem(),
392 element,
393 scv);
394 volVarIndices_[0] = scv.dofIndex();
395 }
396
398 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
399 const VolumeVariables& operator [](const SubControlVolume& scv) const
400 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
401
403 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
404 VolumeVariables& operator [](const SubControlVolume& scv)
405 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
406
408 const VolumeVariables& operator [](std::size_t scvIdx) const
409 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
410
412 VolumeVariables& operator [](std::size_t scvIdx)
413 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
414
417 { return *gridVolVarsPtr_; }
418
420 void clear()
421 {
422 volVarIndices_.clear();
423 volumeVariables_.clear();
424 }
425
426private:
427 const GridVolumeVariables* gridVolVarsPtr_;
428
430 int getLocalIdx_(const int volVarIdx) const
431 {
432 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
433 assert(it != volVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
434 return std::distance(volVarIndices_.begin(), it);
435 }
436
437 std::vector<std::size_t> volVarIndices_;
438 std::vector<VolumeVariables> volumeVariables_;
439};
440
441} // end namespace Dumux
442
443#endif
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:138
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:118
void addBoundaryVolVars(std::vector< VolumeVariables > &volVars, std::vector< IndexType > &volVarIndices, const Problem &problem, const typename FVElemGeom::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElemGeom &fvGeometry)
Adds the boundary volume variables found within the stencil to the provided containers and stores the...
Definition: cellcentered/mpfa/elementvolumevariables.hh:135
std::size_t maxNumBoundaryVolVars(const FVElementGeometry &fvGeometry)
Computes how many boundary vol vars come into play for flux calculations on an element (for a given e...
Definition: cellcentered/mpfa/elementvolumevariables.hh:49
void addBoundaryVolVarsAtNode(std::vector< VolumeVariables > &volVars, std::vector< IndexType > &volVarIndices, const Problem &problem, const typename FVElemGeom::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElemGeom &fvGeometry, const NodalIndexSet &nodalIndexSet)
Adds the boundary volume variables found within a given nodal index set into the provided containers ...
Definition: cellcentered/mpfa/elementvolumevariables.hh:81
Definition: adapt.hh:29
The local (stencil) volume variables class for cell centered mpfa models.
Definition: cellcentered/mpfa/elementvolumevariables.hh:192
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol)
precompute all volume variables in a stencil of an element - bind Dirichlet vol vars in the stencil
Definition: cellcentered/mpfa/elementvolumevariables.hh:232
CCMpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/mpfa/elementvolumevariables.hh:210
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/mpfa/elementvolumevariables.hh:250
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:207
void clear()
Clear all local storage.
Definition: cellcentered/mpfa/elementvolumevariables.hh:256
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/mpfa/elementvolumevariables.hh:263
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:204
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/mpfa/elementvolumevariables.hh:416
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:295
CCMpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/mpfa/elementvolumevariables.hh:298
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/mpfa/elementvolumevariables.hh:303
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol)
Prepares the volume variables of an element.
Definition: cellcentered/mpfa/elementvolumevariables.hh:377
void clear()
Clear all local storage.
Definition: cellcentered/mpfa/elementvolumevariables.hh:420
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:292
The local element solution class for cell-centered methods.