12#ifndef DUMUX_DISCRETIZATION_CCTPFA_ELEMENT_VOLUMEVARIABLES_HH
13#define DUMUX_DISCRETIZATION_CCTPFA_ELEMENT_VOLUMEVARIABLES_HH
31template<
class GVV,
bool cachingEnabled>
52 : gridVolVarsPtr_(&gridVolVars)
53 , numScv_(gridVolVars.problem().gridGeometry().numScv())
57 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
60 if (scv.dofIndex() < numScv_)
61 return gridVolVars().volVars(scv.dofIndex());
63 return boundaryVolumeVariables_[getLocalIdx_(scv.dofIndex())];
70 return gridVolVars().volVars(scvIdx);
72 return boundaryVolumeVariables_[getLocalIdx_(scvIdx)];
80 template<
class FVElementGeometry,
class SolutionVector>
82 const FVElementGeometry& fvGeometry,
83 const SolutionVector& sol) &&
85 this->bind_(element, fvGeometry, sol);
86 return std::move(*
this);
89 template<
class FVElementGeometry,
class SolutionVector>
90 void bind(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
91 const FVElementGeometry& fvGeometry,
92 const SolutionVector& sol) &
93 { this->bind_(element, fvGeometry, sol); }
100 template<
class FVElementGeometry,
class SolutionVector>
102 const FVElementGeometry& fvGeometry,
103 const SolutionVector& sol) &&
105 this->bindElement_(element, fvGeometry, sol);
106 return std::move(*
this);
109 template<
class FVElementGeometry,
class SolutionVector>
110 void bindElement(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
111 const FVElementGeometry& fvGeometry,
112 const SolutionVector& sol) &
113 { this->bindElement_(element, fvGeometry, sol); }
117 {
return *gridVolVarsPtr_; }
124 boundaryVolVarIndices_.clear();
125 boundaryVolumeVariables_.clear();
129 template<
class FVElementGeometry,
class SolutionVector>
130 void bind_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
131 const FVElementGeometry& fvGeometry,
132 const SolutionVector& sol)
134 if (!fvGeometry.hasBoundaryScvf())
138 boundaryVolVarIndices_.reserve(fvGeometry.numScvf());
139 boundaryVolumeVariables_.reserve(fvGeometry.numScvf());
141 for (
const auto& scvf : scvfs(fvGeometry))
143 if (!scvf.boundary())
147 const auto& problem = gridVolVars().problem();
148 const auto bcTypes = problem.boundaryTypes(element, scvf);
149 if (bcTypes.hasOnlyDirichlet())
151 const auto dirichletPriVars = elementSolution<FVElementGeometry>(problem.dirichlet(element, scvf));
152 auto&& scvI = fvGeometry.scv(scvf.insideScvIdx());
154 VolumeVariables volVars;
155 volVars.update(dirichletPriVars,
160 boundaryVolumeVariables_.emplace_back(std::move(volVars));
161 boundaryVolVarIndices_.push_back(scvf.outsideScvIdx());
167 template<
class FVElementGeometry,
class SolutionVector>
168 void bindElement_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
169 const FVElementGeometry& fvGeometry,
170 const SolutionVector& sol)
173 const GridVolumeVariables* gridVolVarsPtr_;
176 int getLocalIdx_(
const int volVarIdx)
const
178 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
179 assert(it != boundaryVolVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
183 std::vector<std::size_t> boundaryVolVarIndices_;
184 std::vector<VolumeVariables> boundaryVolumeVariables_;
185 const std::size_t numScv_;
204 : gridVolVarsPtr_(&gridVolVars) {}
211 template<
class FVElementGeometry,
class SolutionVector>
213 const FVElementGeometry& fvGeometry,
214 const SolutionVector& sol) &&
216 this->bind_(element, fvGeometry, sol);
217 return std::move(*
this);
220 template<
class FVElementGeometry,
class SolutionVector>
221 void bind(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
222 const FVElementGeometry& fvGeometry,
223 const SolutionVector& sol) &
224 { this->bind_(element, fvGeometry, sol); }
231 template<
class FVElementGeometry,
class SolutionVector>
233 const FVElementGeometry& fvGeometry,
234 const SolutionVector& sol) &&
236 this->bindElement_(element, fvGeometry, sol);
237 return std::move(*
this);
240 template<
class FVElementGeometry,
class SolutionVector>
241 void bindElement(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
242 const FVElementGeometry& fvGeometry,
243 const SolutionVector& sol) &
244 { this->bindElement_(element, fvGeometry, sol); }
247 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
249 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
252 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
254 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
258 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
262 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
266 {
return *gridVolVarsPtr_; }
272 volVarIndices_.clear();
273 volumeVariables_.clear();
277 template<
class FVElementGeometry,
class SolutionVector>
278 void bind_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
279 const FVElementGeometry& fvGeometry,
280 const SolutionVector& sol)
284 const auto& problem = gridVolVars().problem();
285 const auto& gridGeometry = fvGeometry.gridGeometry();
286 const auto globalI = gridGeometry.elementMapper().index(element);
287 const auto& connectivityMapI = gridGeometry.connectivityMap()[globalI];
288 const auto numDofs = connectivityMapI.size() + 1;
291 volumeVariables_.resize(numDofs);
292 volVarIndices_.resize(numDofs);
296 auto&& scvI = fvGeometry.scv(globalI);
297 volumeVariables_[localIdx].update(
elementSolution(element, sol, gridGeometry),
301 volVarIndices_[localIdx] = scvI.dofIndex();
305 for (
const auto& dataJ : connectivityMapI)
307 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
308 auto&& scvJ = fvGeometry.scv(dataJ.globalJ);
309 volumeVariables_[localIdx].update(
elementSolution(elementJ, sol, gridGeometry),
313 volVarIndices_[localIdx] = scvJ.dofIndex();
317 if (fvGeometry.hasBoundaryScvf())
320 for (
auto&& scvf : scvfs(fvGeometry))
323 if (!scvf.boundary())
327 const auto bcTypes = problem.boundaryTypes(element, scvf);
328 if (bcTypes.hasOnlyDirichlet())
330 const auto dirichletPriVars = elementSolution<FVElementGeometry>(problem.dirichlet(element, scvf));
332 volumeVariables_.resize(localIdx+1);
333 volVarIndices_.resize(localIdx+1);
334 volumeVariables_[localIdx].update(dirichletPriVars,
338 volVarIndices_[localIdx] = scvf.outsideScvIdx();
366 template<
class FVElementGeometry,
class SolutionVector>
367 void bindElement_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
368 const FVElementGeometry& fvGeometry,
369 const SolutionVector& sol)
373 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
374 volumeVariables_.resize(1);
375 volVarIndices_.resize(1);
378 auto&& scv = fvGeometry.scv(eIdx);
379 volumeVariables_[0].update(
elementSolution(element, sol, fvGeometry.gridGeometry()),
380 gridVolVars().problem(),
383 volVarIndices_[0] = scv.dofIndex();
386 const GridVolumeVariables* gridVolVarsPtr_;
389 int getLocalIdx_(
const int volVarIdx)
const
391 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
392 assert(it != volVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
396 std::vector<std::size_t> volVarIndices_;
397 std::vector<VolumeVariables> volumeVariables_;
The local element solution class for cell-centered methods.
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:221
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:241
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:200
CCTpfaElementVolumeVariables 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/tpfa/elementvolumevariables.hh:232
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:197
CCTpfaElementVolumeVariables 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/tpfa/elementvolumevariables.hh:212
CCTpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/tpfa/elementvolumevariables.hh:203
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/tpfa/elementvolumevariables.hh:265
CCTpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/tpfa/elementvolumevariables.hh:51
CCTpfaElementVolumeVariables 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/tpfa/elementvolumevariables.hh:101
CCTpfaElementVolumeVariables 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/tpfa/elementvolumevariables.hh:81
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:110
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:90
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/tpfa/elementvolumevariables.hh:116
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:45
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:48
The local (stencil) volume variables class for cell centered tpfa models.
Definition: cellcentered/tpfa/elementvolumevariables.hh:33
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