3.6-git
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#include <utility>
32
34
35namespace Dumux {
36namespace CCMpfa {
37
49 template<class FVElementGeometry>
50 std::size_t maxNumBoundaryVolVars(const FVElementGeometry& fvGeometry)
51 {
52 const auto& gridGeometry = fvGeometry.gridGeometry();
53 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
54
55 std::size_t numBoundaryVolVars = 0;
56 for (const auto& scvf : scvfs(fvGeometry))
57 {
58 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
59 numBoundaryVolVars += gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
60 else
61 numBoundaryVolVars += gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
62 }
63
64 return numBoundaryVolVars;
65 }
66
81 template<class VolumeVariables, class IndexType, class Problem, class FVElemGeom, class NodalIndexSet>
82 void addBoundaryVolVarsAtNode(std::vector<VolumeVariables>& volVars,
83 std::vector<IndexType>& volVarIndices,
84 const Problem& problem,
85 const typename FVElemGeom::GridGeometry::GridView::template Codim<0>::Entity& element,
86 const FVElemGeom& fvGeometry,
87 const NodalIndexSet& nodalIndexSet)
88 {
89 if (nodalIndexSet.numBoundaryScvfs() == 0)
90 return;
91
92 // index of the element the fvGeometry was bound to
93 const auto boundElemIdx = fvGeometry.gridGeometry().elementMapper().index(element);
94
95 // check each scvf in the index set for boundary presence
96 for (auto scvfIdx : nodalIndexSet.gridScvfIndices())
97 {
98 const auto& ivScvf = fvGeometry.scvf(scvfIdx);
99
100 // only proceed for scvfs on the boundary and not in the bound element
101 if (!ivScvf.boundary() || ivScvf.insideScvIdx() == boundElemIdx)
102 continue;
103
104 const auto insideScvIdx = ivScvf.insideScvIdx();
105 const auto insideElement = fvGeometry.gridGeometry().element(insideScvIdx);
106 const auto bcTypes = problem.boundaryTypes(insideElement, ivScvf);
107
108 // Only proceed on dirichlet boundaries. On Neumann
109 // boundaries the "outside" vol vars cannot be properly defined.
110 if (bcTypes.hasOnlyDirichlet())
111 {
112 VolumeVariables dirichletVolVars;
113 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(insideElement, ivScvf)),
114 problem,
115 insideElement,
116 fvGeometry.scv(insideScvIdx));
117
118 volVars.emplace_back(std::move(dirichletVolVars));
119 volVarIndices.push_back(ivScvf.outsideScvIdx());
120 }
121 }
122 }
123
135 template<class VolumeVariables, class IndexType, class Problem, class FVElemGeom>
136 void addBoundaryVolVars(std::vector<VolumeVariables>& volVars,
137 std::vector<IndexType>& volVarIndices,
138 const Problem& problem,
139 const typename FVElemGeom::GridGeometry::GridView::template Codim<0>::Entity& element,
140 const FVElemGeom& fvGeometry)
141 {
142 const auto& gridGeometry = fvGeometry.gridGeometry();
143
144 // treat the BCs inside the element
145 if (fvGeometry.hasBoundaryScvf())
146 {
147 const auto boundElemIdx = gridGeometry.elementMapper().index(element);
148 const auto& scvI = fvGeometry.scv(boundElemIdx);
149
150 for (const auto& scvf : scvfs(fvGeometry))
151 {
152 if (!scvf.boundary())
153 continue;
154
155 // Only proceed on dirichlet boundaries. On Neumann
156 // boundaries the "outside" vol vars cannot be properly defined.
157 if (problem.boundaryTypes(element, scvf).hasOnlyDirichlet())
158 {
159 VolumeVariables dirichletVolVars;
160 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(element, scvf)),
161 problem,
162 element,
163 scvI);
164
165 volVars.emplace_back(std::move(dirichletVolVars));
166 volVarIndices.push_back(scvf.outsideScvIdx());
167 }
168 }
169 }
170
171 // Update boundary volume variables in the neighbors
172 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
173 for (const auto& scvf : scvfs(fvGeometry))
174 {
175 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
176 addBoundaryVolVarsAtNode( volVars, volVarIndices, problem, element, fvGeometry,
177 gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet() );
178 else
179 addBoundaryVolVarsAtNode( volVars, volVarIndices, problem, element, fvGeometry,
180 gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet() );
181 }
182 }
183} // end namespace CCMpfa
184
192template<class GVV, bool cachingEnabled>
194
200template<class GVV>
201class CCMpfaElementVolumeVariables<GVV, /*cachingEnabled*/true>
202{
203public:
206
208 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
209
212 : gridVolVarsPtr_(&gridVolVars)
213 , numScv_(gridVolVars.problem().gridGeometry().numScv())
214 {}
215
217 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
218 const VolumeVariables& operator [](const SubControlVolume& scv) const
219 {
220 return scv.dofIndex() < numScv_ ? gridVolVars().volVars(scv.dofIndex())
221 : boundaryVolVars_[getLocalIdx_(scv.dofIndex())];
222 }
223
225 const VolumeVariables& operator [](const std::size_t scvIdx) const
226 {
227 return scvIdx < numScv_ ? gridVolVars().volVars(scvIdx)
228 : boundaryVolVars_[getLocalIdx_(scvIdx)];
229 }
230
236 template<class FVElementGeometry, class SolutionVector>
237 CCMpfaElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
238 const FVElementGeometry& fvGeometry,
239 const SolutionVector& sol) &&
240 {
241 this->bind_(element, fvGeometry, sol);
242 return std::move(*this);
243 }
244
245 template<class FVElementGeometry, class SolutionVector>
246 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
247 const FVElementGeometry& fvGeometry,
248 const SolutionVector& sol) &
249 { this->bind_(element, fvGeometry, sol); }
250
256 template<class FVElementGeometry, class SolutionVector>
257 CCMpfaElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
258 const FVElementGeometry& fvGeometry,
259 const SolutionVector& sol) &&
260 {
261 this->bindElement_(element, fvGeometry, sol);
262 return std::move(*this);
263 }
264
265 template<class FVElementGeometry, class SolutionVector>
266 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
267 const FVElementGeometry& fvGeometry,
268 const SolutionVector& sol) &
269 { this->bindElement_(element, fvGeometry, sol); }
270
272 void clear()
273 {
274 boundaryVolVarIndices_.clear();
275 boundaryVolVars_.clear();
276 }
277
280 { return *gridVolVarsPtr_; }
281
282private:
284 int getLocalIdx_(const int volVarIdx) const
285 {
286 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
287 assert(it != boundaryVolVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
288 return std::distance(boundaryVolVarIndices_.begin(), it);
289 }
290
292 template<class FVElementGeometry, class SolutionVector>
293 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
294 const FVElementGeometry& fvGeometry,
295 const SolutionVector& sol)
296 {
297 clear();
298
299 // maybe prepare boundary volume variables
301 if (maxNumBoundaryVolVars > 0)
302 {
303 boundaryVolVars_.reserve(maxNumBoundaryVolVars);
304 boundaryVolVarIndices_.reserve(maxNumBoundaryVolVars);
305 CCMpfa::addBoundaryVolVars(boundaryVolVars_, boundaryVolVarIndices_, gridVolVars().problem(), element, fvGeometry);
306 }
307 }
308
310 template<class FVElementGeometry, class SolutionVector>
311 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
312 const FVElementGeometry& fvGeometry,
313 const SolutionVector& sol)
314 {}
315
316 const GridVolumeVariables* gridVolVarsPtr_;
317
318 std::size_t numScv_;
319 std::vector<std::size_t> boundaryVolVarIndices_;
320 std::vector<VolumeVariables> boundaryVolVars_;
321};
322
323
328template<class GVV>
329class CCMpfaElementVolumeVariables<GVV, /*cachingEnabled*/false>
330{
331public:
334
336 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
337
340 : gridVolVarsPtr_(&gridVolVars) {}
341
343 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
344 const VolumeVariables& operator [](const SubControlVolume& scv) const
345 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
346
348 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
349 VolumeVariables& operator [](const SubControlVolume& scv)
350 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
351
353 const VolumeVariables& operator [](std::size_t scvIdx) const
354 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
355
357 VolumeVariables& operator [](std::size_t scvIdx)
358 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
359
365 template<class FVElementGeometry, class SolutionVector>
366 CCMpfaElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
367 const FVElementGeometry& fvGeometry,
368 const SolutionVector& sol) &&
369 {
370 this->bind_(element, fvGeometry, sol);
371 return std::move(*this);
372 }
373
374 template<class FVElementGeometry, class SolutionVector>
375 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
376 const FVElementGeometry& fvGeometry,
377 const SolutionVector& sol) &
378 { this->bind_(element, fvGeometry, sol); }
379
385 template<class FVElementGeometry, class SolutionVector>
386 CCMpfaElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
387 const FVElementGeometry& fvGeometry,
388 const SolutionVector& sol) &&
389 {
390 this->bindElement_(element, fvGeometry, sol);
391 return std::move(*this);
392 }
393
394 template<class FVElementGeometry, class SolutionVector>
395 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
396 const FVElementGeometry& fvGeometry,
397 const SolutionVector& sol) &
398 { this->bindElement_(element, fvGeometry, sol); }
399
402 { return *gridVolVarsPtr_; }
403
405 void clear()
406 {
407 volVarIndices_.clear();
408 volumeVariables_.clear();
409 }
410
411private:
412
414 template<class FVElementGeometry, class SolutionVector>
415 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
416 const FVElementGeometry& fvGeometry,
417 const SolutionVector& sol)
418 {
419 clear();
420
421 const auto& problem = gridVolVars().problem();
422 const auto& gridGeometry = fvGeometry.gridGeometry();
423
424 // stencil information
425 const auto globalI = gridGeometry.elementMapper().index(element);
426 const auto& assemblyMapI = gridGeometry.connectivityMap()[globalI];
427 const auto numVolVars = assemblyMapI.size() + 1;
428
429 // resize local containers to the required size (for internal elements)
431 volumeVariables_.reserve(numVolVars+maxNumBoundaryVolVars);
432 volVarIndices_.reserve(numVolVars+maxNumBoundaryVolVars);
433
434 VolumeVariables volVars;
435 const auto& scvI = fvGeometry.scv(globalI);
436 volVars.update(elementSolution(element, sol, gridGeometry),
437 problem,
438 element,
439 scvI);
440
441 volVarIndices_.push_back(scvI.dofIndex());
442 volumeVariables_.emplace_back(std::move(volVars));
443
444 // Update the volume variables of the neighboring elements
445 for (auto&& dataJ : assemblyMapI)
446 {
447 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
448 const auto& scvJ = fvGeometry.scv(dataJ.globalJ);
449 VolumeVariables volVarsJ;
450 volVarsJ.update(elementSolution(elementJ, sol, gridGeometry),
451 problem,
452 elementJ,
453 scvJ);
454
455 volVarIndices_.push_back(scvJ.dofIndex());
456 volumeVariables_.emplace_back(std::move(volVarsJ));
457 }
458
459 // maybe prepare boundary volume variables
460 if (maxNumBoundaryVolVars > 0)
461 CCMpfa::addBoundaryVolVars(volumeVariables_, volVarIndices_, problem, element, fvGeometry);
462
463 // //! TODO Check if user added additional DOF dependencies, i.e. the residual of DOF globalI depends
464 // //! on additional DOFs not included in the discretization schemes' occupation pattern
465 // const auto& additionalDofDependencies = problem.getAdditionalDofDependencies(globalI);
466 // if (!additionalDofDependencies.empty())
467 // {
468 // volumeVariables_.reserve(volumeVariables_.size() + additionalDofDependencies.size());
469 // volVarIndices_.reserve(volVarIndices_.size() + additionalDofDependencies.size());
470 // for (auto globalJ : additionalDofDependencies)
471 // {
472 // const auto& elementJ = gridGeometry.element(globalJ);
473 // const auto& scvJ = fvGeometry.scv(globalJ);
474
475 // VolumeVariables additionalVolVars;
476 // additionalVolVars.update(elementSolution(elementJ, sol, gridGeometry),
477 // problem,
478 // elementJ,
479 // scvJ);
480
481 // volumeVariables_.emplace_back(std::move(additionalVolVars));
482 // volVarIndices_.push_back(globalJ);
483 // }
484 // }
485 }
486
488 template<class FVElementGeometry, class SolutionVector>
489 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
490 const FVElementGeometry& fvGeometry,
491 const SolutionVector& sol)
492 {
493 clear();
494
495 const auto& gridGeometry = fvGeometry.gridGeometry();
496 auto eIdx = gridGeometry.elementMapper().index(element);
497 volumeVariables_.resize(1);
498 volVarIndices_.resize(1);
499
500 // update the volume variables of the element
501 const auto& scv = fvGeometry.scv(eIdx);
502 volumeVariables_[0].update(elementSolution(element, sol, gridGeometry),
503 gridVolVars().problem(),
504 element,
505 scv);
506 volVarIndices_[0] = scv.dofIndex();
507 }
508
509 const GridVolumeVariables* gridVolVarsPtr_;
510
512 int getLocalIdx_(const int volVarIdx) const
513 {
514 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
515 assert(it != volVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
516 return std::distance(volVarIndices_.begin(), it);
517 }
518
519 std::vector<std::size_t> volVarIndices_;
520 std::vector<VolumeVariables> volumeVariables_;
521};
522
523} // end namespace Dumux
524
525#endif
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:294
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::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:136
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:50
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:82
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
The local (stencil) volume variables class for cell centered mpfa models.
Definition: cellcentered/mpfa/elementvolumevariables.hh:193
CCMpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/mpfa/elementvolumevariables.hh:211
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/mpfa/elementvolumevariables.hh:266
CCMpfaElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: cellcentered/mpfa/elementvolumevariables.hh:257
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:208
void clear()
Clear all local storage.
Definition: cellcentered/mpfa/elementvolumevariables.hh:272
CCMpfaElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: cellcentered/mpfa/elementvolumevariables.hh:237
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/mpfa/elementvolumevariables.hh:279
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:205
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/mpfa/elementvolumevariables.hh:246
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/mpfa/elementvolumevariables.hh:401
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:336
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/mpfa/elementvolumevariables.hh:395
CCMpfaElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: cellcentered/mpfa/elementvolumevariables.hh:386
CCMpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/mpfa/elementvolumevariables.hh:339
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/mpfa/elementvolumevariables.hh:375
CCMpfaElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: cellcentered/mpfa/elementvolumevariables.hh:366
void clear()
Clear all local storage.
Definition: cellcentered/mpfa/elementvolumevariables.hh:405
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:333
The local element solution class for cell-centered methods.