24#ifndef DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_VOLUMEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_VOLUMEVARIABLES_HH
48 template<
class FVElementGeometry>
51 const auto& gridGeometry = fvGeometry.gridGeometry();
52 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
54 std::size_t numBoundaryVolVars = 0;
55 for (
const auto& scvf : scvfs(fvGeometry))
57 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
58 numBoundaryVolVars += gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
60 numBoundaryVolVars += gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs();
63 return numBoundaryVolVars;
80 template<
class VolumeVariables,
class IndexType,
class Problem,
class FVElemGeom,
class NodalIndexSet>
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)
88 if (nodalIndexSet.numBoundaryScvfs() == 0)
92 const auto boundElemIdx = fvGeometry.gridGeometry().elementMapper().index(element);
95 for (
auto scvfIdx : nodalIndexSet.gridScvfIndices())
97 const auto& ivScvf = fvGeometry.scvf(scvfIdx);
100 if (!ivScvf.boundary() || ivScvf.insideScvIdx() == boundElemIdx)
103 const auto insideScvIdx = ivScvf.insideScvIdx();
104 const auto insideElement = fvGeometry.gridGeometry().element(insideScvIdx);
105 const auto bcTypes = problem.boundaryTypes(insideElement, ivScvf);
109 if (bcTypes.hasOnlyDirichlet())
111 VolumeVariables dirichletVolVars;
112 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(insideElement, ivScvf)),
115 fvGeometry.scv(insideScvIdx));
117 volVars.emplace_back(std::move(dirichletVolVars));
118 volVarIndices.push_back(ivScvf.outsideScvIdx());
134 template<
class VolumeVariables,
class IndexType,
class Problem,
class FVElemGeom>
136 std::vector<IndexType>& volVarIndices,
137 const Problem& problem,
138 const typename FVElemGeom::GridGeometry::GridView::template Codim<0>::Entity& element,
139 const FVElemGeom& fvGeometry)
141 const auto& gridGeometry = fvGeometry.gridGeometry();
144 if (fvGeometry.hasBoundaryScvf())
146 const auto boundElemIdx = gridGeometry.elementMapper().index(element);
147 const auto& scvI = fvGeometry.scv(boundElemIdx);
149 for (
const auto& scvf : scvfs(fvGeometry))
151 if (!scvf.boundary())
156 if (problem.boundaryTypes(element, scvf).hasOnlyDirichlet())
158 VolumeVariables dirichletVolVars;
159 dirichletVolVars.update(elementSolution<FVElemGeom>(problem.dirichlet(element, scvf)),
164 volVars.emplace_back(std::move(dirichletVolVars));
165 volVarIndices.push_back(scvf.outsideScvIdx());
171 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
172 for (
const auto& scvf : scvfs(fvGeometry))
174 if (!gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
176 gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet() );
179 gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet() );
191template<
class GVV,
bool cachingEnabled>
211 : gridVolVarsPtr_(&gridVolVars)
212 , numScv_(gridVolVars.problem().gridGeometry().numScv())
216 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
219 return scv.dofIndex() < numScv_ ? gridVolVars().volVars(scv.dofIndex())
220 : boundaryVolVars_[getLocalIdx_(scv.dofIndex())];
226 return scvIdx < numScv_ ? gridVolVars().volVars(scvIdx)
227 : boundaryVolVars_[getLocalIdx_(scvIdx)];
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)
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)
258 boundaryVolVarIndices_.clear();
259 boundaryVolVars_.clear();
264 {
return *gridVolVarsPtr_; }
268 int getLocalIdx_(
const int volVarIdx)
const
270 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
271 assert(it != boundaryVolVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
275 const GridVolumeVariables* gridVolVarsPtr_;
278 std::vector<std::size_t> boundaryVolVarIndices_;
279 std::vector<VolumeVariables> boundaryVolVars_;
299 : gridVolVarsPtr_(&gridVolVars) {}
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)
309 const auto& problem = gridVolVars().problem();
310 const auto& gridGeometry = fvGeometry.gridGeometry();
313 const auto globalI = gridGeometry.elementMapper().index(element);
314 const auto& assemblyMapI = gridGeometry.connectivityMap()[globalI];
315 const auto numVolVars = assemblyMapI.size() + 1;
323 const auto& scvI = fvGeometry.scv(globalI);
329 volVarIndices_.push_back(scvI.dofIndex());
330 volumeVariables_.emplace_back(std::move(volVars));
333 for (
auto&& dataJ : assemblyMapI)
335 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
336 const auto& scvJ = fvGeometry.scv(dataJ.globalJ);
343 volVarIndices_.push_back(scvJ.dofIndex());
344 volumeVariables_.emplace_back(std::move(volVarsJ));
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)
383 const auto& gridGeometry = fvGeometry.gridGeometry();
384 auto eIdx = gridGeometry.elementMapper().index(element);
385 volumeVariables_.resize(1);
386 volVarIndices_.resize(1);
389 const auto& scv = fvGeometry.scv(eIdx);
390 volumeVariables_[0].update(
elementSolution(element, sol, gridGeometry),
391 gridVolVars().problem(),
394 volVarIndices_[0] = scv.dofIndex();
398 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
400 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
403 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
405 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
409 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
413 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
417 {
return *gridVolVarsPtr_; }
422 volVarIndices_.clear();
423 volumeVariables_.clear();
427 const GridVolumeVariables* gridVolVarsPtr_;
430 int getLocalIdx_(
const int volVarIdx)
const
432 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
433 assert(it != volVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
437 std::vector<std::size_t> volVarIndices_;
438 std::vector<VolumeVariables> volumeVariables_;
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: geometry/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:115
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
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.