12#ifndef DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_VOLUMEVARIABLES_HH
13#define DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_VOLUMEVARIABLES_HH
37 template<
class FVElementGeometry>
40 const auto& gridGeometry = fvGeometry.gridGeometry();
41 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
43 std::size_t numBoundaryVolVars = 0;
44 for (
const auto& scvf : scvfs(fvGeometry))
46 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
47 numBoundaryVolVars += gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
49 numBoundaryVolVars += gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
52 return numBoundaryVolVars;
69 template<
class VolumeVariables,
class IndexType,
class Problem,
class FVElemGeom,
class NodalIndexSet>
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)
77 if (nodalIndexSet.numBoundaryScvfs() == 0)
81 const auto boundElemIdx = fvGeometry.gridGeometry().elementMapper().index(element);
84 for (
auto scvfIdx : nodalIndexSet.gridScvfIndices())
86 const auto& ivScvf = fvGeometry.scvf(scvfIdx);
89 if (!ivScvf.boundary() || ivScvf.insideScvIdx() == boundElemIdx)
92 const auto insideScvIdx = ivScvf.insideScvIdx();
93 const auto insideElement = fvGeometry.gridGeometry().element(insideScvIdx);
94 const auto bcTypes = problem.boundaryTypes(insideElement, ivScvf);
98 if (bcTypes.hasOnlyDirichlet())
100 VolumeVariables dirichletVolVars;
101 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(insideElement, ivScvf)),
104 fvGeometry.scv(insideScvIdx));
106 volVars.emplace_back(std::move(dirichletVolVars));
107 volVarIndices.push_back(ivScvf.outsideScvIdx());
123 template<
class VolumeVariables,
class IndexType,
class Problem,
class FVElemGeom>
125 std::vector<IndexType>& volVarIndices,
126 const Problem& problem,
127 const typename FVElemGeom::GridGeometry::GridView::template Codim<0>::Entity& element,
128 const FVElemGeom& fvGeometry)
130 const auto& gridGeometry = fvGeometry.gridGeometry();
133 if (fvGeometry.hasBoundaryScvf())
135 const auto boundElemIdx = gridGeometry.elementMapper().index(element);
136 const auto& scvI = fvGeometry.scv(boundElemIdx);
138 for (
const auto& scvf : scvfs(fvGeometry))
140 if (!scvf.boundary())
145 if (problem.boundaryTypes(element, scvf).hasOnlyDirichlet())
147 VolumeVariables dirichletVolVars;
148 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(element, scvf)),
153 volVars.emplace_back(std::move(dirichletVolVars));
154 volVarIndices.push_back(scvf.outsideScvIdx());
160 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
161 for (
const auto& scvf : scvfs(fvGeometry))
163 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
165 gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet() );
168 gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet() );
180template<
class GVV,
bool cachingEnabled>
200 : gridVolVarsPtr_(&gridVolVars)
201 , numScv_(gridVolVars.problem().gridGeometry().numScv())
205 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
208 return scv.dofIndex() < numScv_ ? gridVolVars().volVars(scv.dofIndex())
209 : boundaryVolVars_[getLocalIdx_(scv.dofIndex())];
215 return scvIdx < numScv_ ? gridVolVars().volVars(scvIdx)
216 : boundaryVolVars_[getLocalIdx_(scvIdx)];
224 template<
class FVElementGeometry,
class SolutionVector>
226 const FVElementGeometry& fvGeometry,
227 const SolutionVector& sol) &&
229 this->bind_(element, fvGeometry, sol);
230 return std::move(*
this);
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); }
244 template<
class FVElementGeometry,
class SolutionVector>
246 const FVElementGeometry& fvGeometry,
247 const SolutionVector& sol) &&
249 this->bindElement_(element, fvGeometry, sol);
250 return std::move(*
this);
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); }
262 boundaryVolVarIndices_.clear();
263 boundaryVolVars_.clear();
268 {
return *gridVolVarsPtr_; }
272 int getLocalIdx_(
const int volVarIdx)
const
274 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
275 assert(it != boundaryVolVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
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)
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)
304 const GridVolumeVariables* gridVolVarsPtr_;
307 std::vector<std::size_t> boundaryVolVarIndices_;
308 std::vector<VolumeVariables> boundaryVolVars_;
328 : gridVolVarsPtr_(&gridVolVars) {}
331 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
333 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
336 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
338 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
342 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
346 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
353 template<
class FVElementGeometry,
class SolutionVector>
355 const FVElementGeometry& fvGeometry,
356 const SolutionVector& sol) &&
358 this->bind_(element, fvGeometry, sol);
359 return std::move(*
this);
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); }
373 template<
class FVElementGeometry,
class SolutionVector>
375 const FVElementGeometry& fvGeometry,
376 const SolutionVector& sol) &&
378 this->bindElement_(element, fvGeometry, sol);
379 return std::move(*
this);
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); }
390 {
return *gridVolVarsPtr_; }
395 volVarIndices_.clear();
396 volumeVariables_.clear();
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)
409 const auto& problem = gridVolVars().problem();
410 const auto& gridGeometry = fvGeometry.gridGeometry();
413 const auto globalI = gridGeometry.elementMapper().index(element);
414 const auto& assemblyMapI = gridGeometry.connectivityMap()[globalI];
415 const auto numVolVars = assemblyMapI.size() + 1;
422 VolumeVariables volVars;
423 const auto& scvI = fvGeometry.scv(globalI);
429 volVarIndices_.push_back(scvI.dofIndex());
430 volumeVariables_.emplace_back(std::move(volVars));
433 for (
auto&& dataJ : assemblyMapI)
435 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
436 const auto& scvJ = fvGeometry.scv(dataJ.globalJ);
437 VolumeVariables volVarsJ;
443 volVarIndices_.push_back(scvJ.dofIndex());
444 volumeVariables_.emplace_back(std::move(volVarsJ));
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)
483 const auto& gridGeometry = fvGeometry.gridGeometry();
484 auto eIdx = gridGeometry.elementMapper().index(element);
485 volumeVariables_.resize(1);
486 volVarIndices_.resize(1);
489 const auto& scv = fvGeometry.scv(eIdx);
490 volumeVariables_[0].update(
elementSolution(element, sol, gridGeometry),
491 gridVolVars().problem(),
494 volVarIndices_[0] = scv.dofIndex();
497 const GridVolumeVariables* gridVolVarsPtr_;
500 int getLocalIdx_(
const int volVarIdx)
const
502 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
503 assert(it != volVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
507 std::vector<std::size_t> volVarIndices_;
508 std::vector<VolumeVariables> volumeVariables_;
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