version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_VOLUMEVARIABLES_HH
13#define DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_VOLUMEVARIABLES_HH
14
15#include <algorithm>
16#include <type_traits>
17#include <utility>
18#include <vector>
19#include <utility>
20
22
23namespace Dumux {
24namespace CCMpfa {
25
37 template<class FVElementGeometry>
38 std::size_t maxNumBoundaryVolVars(const FVElementGeometry& fvGeometry)
39 {
40 const auto& gridGeometry = fvGeometry.gridGeometry();
41 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
42
43 std::size_t numBoundaryVolVars = 0;
44 for (const auto& scvf : scvfs(fvGeometry))
45 {
46 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
47 numBoundaryVolVars += gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
48 else
49 numBoundaryVolVars += gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
50 }
51
52 return numBoundaryVolVars;
53 }
54
69 template<class VolumeVariables, class IndexType, class Problem, class FVElemGeom, class NodalIndexSet>
70 void addBoundaryVolVarsAtNode(std::vector<VolumeVariables>& volVars,
71 std::vector<IndexType>& volVarIndices,
72 const Problem& problem,
73 const typename FVElemGeom::GridGeometry::GridView::template Codim<0>::Entity& element,
74 const FVElemGeom& fvGeometry,
75 const NodalIndexSet& nodalIndexSet)
76 {
77 if (nodalIndexSet.numBoundaryScvfs() == 0)
78 return;
79
80 // index of the element the fvGeometry was bound to
81 const auto boundElemIdx = fvGeometry.gridGeometry().elementMapper().index(element);
82
83 // check each scvf in the index set for boundary presence
84 for (auto scvfIdx : nodalIndexSet.gridScvfIndices())
85 {
86 const auto& ivScvf = fvGeometry.scvf(scvfIdx);
87
88 // only proceed for scvfs on the boundary and not in the bound element
89 if (!ivScvf.boundary() || ivScvf.insideScvIdx() == boundElemIdx)
90 continue;
91
92 const auto insideScvIdx = ivScvf.insideScvIdx();
93 const auto insideElement = fvGeometry.gridGeometry().element(insideScvIdx);
94 const auto bcTypes = problem.boundaryTypes(insideElement, ivScvf);
95
96 // Only proceed on dirichlet boundaries. On Neumann
97 // boundaries the "outside" vol vars cannot be properly defined.
98 if (bcTypes.hasOnlyDirichlet())
99 {
100 VolumeVariables dirichletVolVars;
101 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(insideElement, ivScvf)),
102 problem,
103 insideElement,
104 fvGeometry.scv(insideScvIdx));
105
106 volVars.emplace_back(std::move(dirichletVolVars));
107 volVarIndices.push_back(ivScvf.outsideScvIdx());
108 }
109 }
110 }
111
123 template<class VolumeVariables, class IndexType, class Problem, class FVElemGeom>
124 void addBoundaryVolVars(std::vector<VolumeVariables>& volVars,
125 std::vector<IndexType>& volVarIndices,
126 const Problem& problem,
127 const typename FVElemGeom::GridGeometry::GridView::template Codim<0>::Entity& element,
128 const FVElemGeom& fvGeometry)
129 {
130 const auto& gridGeometry = fvGeometry.gridGeometry();
131
132 // treat the BCs inside the element
133 if (fvGeometry.hasBoundaryScvf())
134 {
135 const auto boundElemIdx = gridGeometry.elementMapper().index(element);
136 const auto& scvI = fvGeometry.scv(boundElemIdx);
137
138 for (const auto& scvf : scvfs(fvGeometry))
139 {
140 if (!scvf.boundary())
141 continue;
142
143 // Only proceed on dirichlet boundaries. On Neumann
144 // boundaries the "outside" vol vars cannot be properly defined.
145 if (problem.boundaryTypes(element, scvf).hasOnlyDirichlet())
146 {
147 VolumeVariables dirichletVolVars;
148 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(element, scvf)),
149 problem,
150 element,
151 scvI);
152
153 volVars.emplace_back(std::move(dirichletVolVars));
154 volVarIndices.push_back(scvf.outsideScvIdx());
155 }
156 }
157 }
158
159 // Update boundary volume variables in the neighbors
160 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
161 for (const auto& scvf : scvfs(fvGeometry))
162 {
163 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
164 addBoundaryVolVarsAtNode( volVars, volVarIndices, problem, element, fvGeometry,
165 gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet() );
166 else
167 addBoundaryVolVarsAtNode( volVars, volVarIndices, problem, element, fvGeometry,
168 gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet() );
169 }
170 }
171} // end namespace CCMpfa
172
180template<class GVV, bool cachingEnabled>
182
188template<class GVV>
189class CCMpfaElementVolumeVariables<GVV, /*cachingEnabled*/true>
190{
191public:
194
196 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
197
200 : gridVolVarsPtr_(&gridVolVars)
201 , numScv_(gridVolVars.problem().gridGeometry().numScv())
202 {}
203
205 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
206 const VolumeVariables& operator [](const SubControlVolume& scv) const
207 {
208 return scv.dofIndex() < numScv_ ? gridVolVars().volVars(scv.dofIndex())
209 : boundaryVolVars_[getLocalIdx_(scv.dofIndex())];
210 }
211
213 const VolumeVariables& operator [](const std::size_t scvIdx) const
214 {
215 return scvIdx < numScv_ ? gridVolVars().volVars(scvIdx)
216 : boundaryVolVars_[getLocalIdx_(scvIdx)];
217 }
218
224 template<class FVElementGeometry, class SolutionVector>
225 CCMpfaElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
226 const FVElementGeometry& fvGeometry,
227 const SolutionVector& sol) &&
228 {
229 this->bind_(element, fvGeometry, sol);
230 return std::move(*this);
231 }
232
233 template<class FVElementGeometry, class SolutionVector>
234 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
235 const FVElementGeometry& fvGeometry,
236 const SolutionVector& sol) &
237 { this->bind_(element, fvGeometry, sol); }
238
244 template<class FVElementGeometry, class SolutionVector>
245 CCMpfaElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
246 const FVElementGeometry& fvGeometry,
247 const SolutionVector& sol) &&
248 {
249 this->bindElement_(element, fvGeometry, sol);
250 return std::move(*this);
251 }
252
253 template<class FVElementGeometry, class SolutionVector>
254 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
255 const FVElementGeometry& fvGeometry,
256 const SolutionVector& sol) &
257 { this->bindElement_(element, fvGeometry, sol); }
258
260 void clear()
261 {
262 boundaryVolVarIndices_.clear();
263 boundaryVolVars_.clear();
264 }
265
268 { return *gridVolVarsPtr_; }
269
270private:
272 int getLocalIdx_(const int volVarIdx) const
273 {
274 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
275 assert(it != boundaryVolVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
276 return std::distance(boundaryVolVarIndices_.begin(), it);
277 }
278
280 template<class FVElementGeometry, class SolutionVector>
281 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
282 const FVElementGeometry& fvGeometry,
283 const SolutionVector& sol)
284 {
285 clear();
286
287 // maybe prepare boundary volume variables
289 if (maxNumBoundaryVolVars > 0)
290 {
291 boundaryVolVars_.reserve(maxNumBoundaryVolVars);
292 boundaryVolVarIndices_.reserve(maxNumBoundaryVolVars);
293 CCMpfa::addBoundaryVolVars(boundaryVolVars_, boundaryVolVarIndices_, gridVolVars().problem(), element, fvGeometry);
294 }
295 }
296
298 template<class FVElementGeometry, class SolutionVector>
299 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
300 const FVElementGeometry& fvGeometry,
301 const SolutionVector& sol)
302 {}
303
304 const GridVolumeVariables* gridVolVarsPtr_;
305
306 std::size_t numScv_;
307 std::vector<std::size_t> boundaryVolVarIndices_;
308 std::vector<VolumeVariables> boundaryVolVars_;
309};
310
311
316template<class GVV>
317class CCMpfaElementVolumeVariables<GVV, /*cachingEnabled*/false>
318{
319public:
322
324 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
325
328 : gridVolVarsPtr_(&gridVolVars) {}
329
331 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
332 const VolumeVariables& operator [](const SubControlVolume& scv) const
333 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
334
336 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
337 VolumeVariables& operator [](const SubControlVolume& scv)
338 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
339
341 const VolumeVariables& operator [](std::size_t scvIdx) const
342 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
343
345 VolumeVariables& operator [](std::size_t scvIdx)
346 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
347
353 template<class FVElementGeometry, class SolutionVector>
354 CCMpfaElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
355 const FVElementGeometry& fvGeometry,
356 const SolutionVector& sol) &&
357 {
358 this->bind_(element, fvGeometry, sol);
359 return std::move(*this);
360 }
361
362 template<class FVElementGeometry, class SolutionVector>
363 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
364 const FVElementGeometry& fvGeometry,
365 const SolutionVector& sol) &
366 { this->bind_(element, fvGeometry, sol); }
367
373 template<class FVElementGeometry, class SolutionVector>
374 CCMpfaElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
375 const FVElementGeometry& fvGeometry,
376 const SolutionVector& sol) &&
377 {
378 this->bindElement_(element, fvGeometry, sol);
379 return std::move(*this);
380 }
381
382 template<class FVElementGeometry, class SolutionVector>
383 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
384 const FVElementGeometry& fvGeometry,
385 const SolutionVector& sol) &
386 { this->bindElement_(element, fvGeometry, sol); }
387
390 { return *gridVolVarsPtr_; }
391
393 void clear()
394 {
395 volVarIndices_.clear();
396 volumeVariables_.clear();
397 }
398
399private:
400
402 template<class FVElementGeometry, class SolutionVector>
403 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
404 const FVElementGeometry& fvGeometry,
405 const SolutionVector& sol)
406 {
407 clear();
408
409 const auto& problem = gridVolVars().problem();
410 const auto& gridGeometry = fvGeometry.gridGeometry();
411
412 // stencil information
413 const auto globalI = gridGeometry.elementMapper().index(element);
414 const auto& assemblyMapI = gridGeometry.connectivityMap()[globalI];
415 const auto numVolVars = assemblyMapI.size() + 1;
416
417 // resize local containers to the required size (for internal elements)
419 volumeVariables_.reserve(numVolVars+maxNumBoundaryVolVars);
420 volVarIndices_.reserve(numVolVars+maxNumBoundaryVolVars);
421
422 VolumeVariables volVars;
423 const auto& scvI = fvGeometry.scv(globalI);
424 volVars.update(elementSolution(element, sol, gridGeometry),
425 problem,
426 element,
427 scvI);
428
429 volVarIndices_.push_back(scvI.dofIndex());
430 volumeVariables_.emplace_back(std::move(volVars));
431
432 // Update the volume variables of the neighboring elements
433 for (auto&& dataJ : assemblyMapI)
434 {
435 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
436 const auto& scvJ = fvGeometry.scv(dataJ.globalJ);
437 VolumeVariables volVarsJ;
438 volVarsJ.update(elementSolution(elementJ, sol, gridGeometry),
439 problem,
440 elementJ,
441 scvJ);
442
443 volVarIndices_.push_back(scvJ.dofIndex());
444 volumeVariables_.emplace_back(std::move(volVarsJ));
445 }
446
447 // maybe prepare boundary volume variables
448 if (maxNumBoundaryVolVars > 0)
449 CCMpfa::addBoundaryVolVars(volumeVariables_, volVarIndices_, problem, element, fvGeometry);
450
451 // //! TODO Check if user added additional DOF dependencies, i.e. the residual of DOF globalI depends
452 // //! on additional DOFs not included in the discretization schemes' occupation pattern
453 // const auto& additionalDofDependencies = problem.getAdditionalDofDependencies(globalI);
454 // if (!additionalDofDependencies.empty())
455 // {
456 // volumeVariables_.reserve(volumeVariables_.size() + additionalDofDependencies.size());
457 // volVarIndices_.reserve(volVarIndices_.size() + additionalDofDependencies.size());
458 // for (auto globalJ : additionalDofDependencies)
459 // {
460 // const auto& elementJ = gridGeometry.element(globalJ);
461 // const auto& scvJ = fvGeometry.scv(globalJ);
462
463 // VolumeVariables additionalVolVars;
464 // additionalVolVars.update(elementSolution(elementJ, sol, gridGeometry),
465 // problem,
466 // elementJ,
467 // scvJ);
468
469 // volumeVariables_.emplace_back(std::move(additionalVolVars));
470 // volVarIndices_.push_back(globalJ);
471 // }
472 // }
473 }
474
476 template<class FVElementGeometry, class SolutionVector>
477 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
478 const FVElementGeometry& fvGeometry,
479 const SolutionVector& sol)
480 {
481 clear();
482
483 const auto& gridGeometry = fvGeometry.gridGeometry();
484 auto eIdx = gridGeometry.elementMapper().index(element);
485 volumeVariables_.resize(1);
486 volVarIndices_.resize(1);
487
488 // update the volume variables of the element
489 const auto& scv = fvGeometry.scv(eIdx);
490 volumeVariables_[0].update(elementSolution(element, sol, gridGeometry),
491 gridVolVars().problem(),
492 element,
493 scv);
494 volVarIndices_[0] = scv.dofIndex();
495 }
496
497 const GridVolumeVariables* gridVolVarsPtr_;
498
500 int getLocalIdx_(const int volVarIdx) const
501 {
502 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
503 assert(it != volVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
504 return std::distance(volVarIndices_.begin(), it);
505 }
506
507 std::vector<std::size_t> volVarIndices_;
508 std::vector<VolumeVariables> volumeVariables_;
509};
510
511} // end namespace Dumux
512
513#endif
The local element solution class for cell-centered methods.
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/mpfa/elementvolumevariables.hh:389
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:324
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/mpfa/elementvolumevariables.hh:383
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:374
CCMpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/mpfa/elementvolumevariables.hh:327
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/mpfa/elementvolumevariables.hh:363
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:354
void clear()
Clear all local storage.
Definition: cellcentered/mpfa/elementvolumevariables.hh:393
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:321
CCMpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/mpfa/elementvolumevariables.hh:199
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/mpfa/elementvolumevariables.hh:254
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:245
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:196
void clear()
Clear all local storage.
Definition: cellcentered/mpfa/elementvolumevariables.hh:260
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:225
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/mpfa/elementvolumevariables.hh:267
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/mpfa/elementvolumevariables.hh:193
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/mpfa/elementvolumevariables.hh:234
The local (stencil) volume variables class for cell centered mpfa models.
Definition: cellcentered/mpfa/elementvolumevariables.hh:181
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:124
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:38
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:70
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
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:282
Definition: adapt.hh:17