24#ifndef DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_VOLUMEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_VOLUMEVARIABLES_HH
49 template<
class FVElementGeometry>
52 const auto& gridGeometry = fvGeometry.gridGeometry();
53 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
55 std::size_t numBoundaryVolVars = 0;
56 for (
const auto& scvf : scvfs(fvGeometry))
58 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
59 numBoundaryVolVars += gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
61 numBoundaryVolVars += gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
64 return numBoundaryVolVars;
81 template<
class VolumeVariables,
class IndexType,
class Problem,
class FVElemGeom,
class NodalIndexSet>
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)
89 if (nodalIndexSet.numBoundaryScvfs() == 0)
93 const auto boundElemIdx = fvGeometry.gridGeometry().elementMapper().index(element);
96 for (
auto scvfIdx : nodalIndexSet.gridScvfIndices())
98 const auto& ivScvf = fvGeometry.scvf(scvfIdx);
101 if (!ivScvf.boundary() || ivScvf.insideScvIdx() == boundElemIdx)
104 const auto insideScvIdx = ivScvf.insideScvIdx();
105 const auto insideElement = fvGeometry.gridGeometry().element(insideScvIdx);
106 const auto bcTypes = problem.boundaryTypes(insideElement, ivScvf);
110 if (bcTypes.hasOnlyDirichlet())
112 VolumeVariables dirichletVolVars;
113 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(insideElement, ivScvf)),
116 fvGeometry.scv(insideScvIdx));
118 volVars.emplace_back(std::move(dirichletVolVars));
119 volVarIndices.push_back(ivScvf.outsideScvIdx());
135 template<
class VolumeVariables,
class IndexType,
class Problem,
class FVElemGeom>
137 std::vector<IndexType>& volVarIndices,
138 const Problem& problem,
139 const typename FVElemGeom::GridGeometry::GridView::template Codim<0>::Entity& element,
140 const FVElemGeom& fvGeometry)
142 const auto& gridGeometry = fvGeometry.gridGeometry();
145 if (fvGeometry.hasBoundaryScvf())
147 const auto boundElemIdx = gridGeometry.elementMapper().index(element);
148 const auto& scvI = fvGeometry.scv(boundElemIdx);
150 for (
const auto& scvf : scvfs(fvGeometry))
152 if (!scvf.boundary())
157 if (problem.boundaryTypes(element, scvf).hasOnlyDirichlet())
159 VolumeVariables dirichletVolVars;
160 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(element, scvf)),
165 volVars.emplace_back(std::move(dirichletVolVars));
166 volVarIndices.push_back(scvf.outsideScvIdx());
172 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
173 for (
const auto& scvf : scvfs(fvGeometry))
175 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
177 gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet() );
180 gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet() );
192template<
class GVV,
bool cachingEnabled>
212 : gridVolVarsPtr_(&gridVolVars)
213 , numScv_(gridVolVars.problem().gridGeometry().numScv())
217 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
220 return scv.dofIndex() < numScv_ ? gridVolVars().volVars(scv.dofIndex())
221 : boundaryVolVars_[getLocalIdx_(scv.dofIndex())];
227 return scvIdx < numScv_ ? gridVolVars().volVars(scvIdx)
228 : boundaryVolVars_[getLocalIdx_(scvIdx)];
236 template<
class FVElementGeometry,
class SolutionVector>
238 const FVElementGeometry& fvGeometry,
239 const SolutionVector& sol) &&
241 this->bind_(element, fvGeometry, sol);
242 return std::move(*
this);
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); }
256 template<
class FVElementGeometry,
class SolutionVector>
258 const FVElementGeometry& fvGeometry,
259 const SolutionVector& sol) &&
261 this->bindElement_(element, fvGeometry, sol);
262 return std::move(*
this);
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); }
274 boundaryVolVarIndices_.clear();
275 boundaryVolVars_.clear();
280 {
return *gridVolVarsPtr_; }
284 int getLocalIdx_(
const int volVarIdx)
const
286 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
287 assert(it != boundaryVolVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
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)
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)
316 const GridVolumeVariables* gridVolVarsPtr_;
319 std::vector<std::size_t> boundaryVolVarIndices_;
320 std::vector<VolumeVariables> boundaryVolVars_;
340 : gridVolVarsPtr_(&gridVolVars) {}
343 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
345 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
348 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
350 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
354 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
358 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
365 template<
class FVElementGeometry,
class SolutionVector>
367 const FVElementGeometry& fvGeometry,
368 const SolutionVector& sol) &&
370 this->bind_(element, fvGeometry, sol);
371 return std::move(*
this);
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); }
385 template<
class FVElementGeometry,
class SolutionVector>
387 const FVElementGeometry& fvGeometry,
388 const SolutionVector& sol) &&
390 this->bindElement_(element, fvGeometry, sol);
391 return std::move(*
this);
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); }
402 {
return *gridVolVarsPtr_; }
407 volVarIndices_.clear();
408 volumeVariables_.clear();
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)
421 const auto& problem = gridVolVars().problem();
422 const auto& gridGeometry = fvGeometry.gridGeometry();
425 const auto globalI = gridGeometry.elementMapper().index(element);
426 const auto& assemblyMapI = gridGeometry.connectivityMap()[globalI];
427 const auto numVolVars = assemblyMapI.size() + 1;
434 VolumeVariables volVars;
435 const auto& scvI = fvGeometry.scv(globalI);
441 volVarIndices_.push_back(scvI.dofIndex());
442 volumeVariables_.emplace_back(std::move(volVars));
445 for (
auto&& dataJ : assemblyMapI)
447 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
448 const auto& scvJ = fvGeometry.scv(dataJ.globalJ);
449 VolumeVariables volVarsJ;
455 volVarIndices_.push_back(scvJ.dofIndex());
456 volumeVariables_.emplace_back(std::move(volVarsJ));
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)
495 const auto& gridGeometry = fvGeometry.gridGeometry();
496 auto eIdx = gridGeometry.elementMapper().index(element);
497 volumeVariables_.resize(1);
498 volVarIndices_.resize(1);
501 const auto& scv = fvGeometry.scv(eIdx);
502 volumeVariables_[0].update(
elementSolution(element, sol, gridGeometry),
503 gridVolVars().problem(),
506 volVarIndices_[0] = scv.dofIndex();
509 const GridVolumeVariables* gridVolVarsPtr_;
512 int getLocalIdx_(
const int volVarIdx)
const
514 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
515 assert(it != volVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
519 std::vector<std::size_t> volVarIndices_;
520 std::vector<VolumeVariables> volumeVariables_;
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.