24#ifndef DUMUX_DISCRETIZATION_CCTPFA_ELEMENT_VOLUMEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_CCTPFA_ELEMENT_VOLUMEVARIABLES_HH
43template<
class GVV,
bool cachingEnabled>
64 : gridVolVarsPtr_(&gridVolVars)
65 , numScv_(gridVolVars.problem().gridGeometry().numScv())
69 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
72 if (scv.dofIndex() < numScv_)
73 return gridVolVars().volVars(scv.dofIndex());
75 return boundaryVolumeVariables_[getLocalIdx_(scv.dofIndex())];
82 return gridVolVars().volVars(scvIdx);
84 return boundaryVolumeVariables_[getLocalIdx_(scvIdx)];
92 template<
class FVElementGeometry,
class SolutionVector>
94 const FVElementGeometry& fvGeometry,
95 const SolutionVector& sol) &&
97 this->bind_(element, fvGeometry, sol);
98 return std::move(*
this);
101 template<
class FVElementGeometry,
class SolutionVector>
102 void bind(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
103 const FVElementGeometry& fvGeometry,
104 const SolutionVector& sol) &
105 { this->bind_(element, fvGeometry, sol); }
112 template<
class FVElementGeometry,
class SolutionVector>
114 const FVElementGeometry& fvGeometry,
115 const SolutionVector& sol) &&
117 this->bindElement_(element, fvGeometry, sol);
118 return std::move(*
this);
121 template<
class FVElementGeometry,
class SolutionVector>
122 void bindElement(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
123 const FVElementGeometry& fvGeometry,
124 const SolutionVector& sol) &
125 { this->bindElement_(element, fvGeometry, sol); }
129 {
return *gridVolVarsPtr_; }
136 boundaryVolVarIndices_.clear();
137 boundaryVolumeVariables_.clear();
141 template<
class FVElementGeometry,
class SolutionVector>
142 void bind_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
143 const FVElementGeometry& fvGeometry,
144 const SolutionVector& sol)
146 if (!fvGeometry.hasBoundaryScvf())
150 boundaryVolVarIndices_.reserve(fvGeometry.numScvf());
151 boundaryVolumeVariables_.reserve(fvGeometry.numScvf());
153 for (
const auto& scvf : scvfs(fvGeometry))
155 if (!scvf.boundary())
159 const auto& problem = gridVolVars().problem();
160 const auto bcTypes = problem.boundaryTypes(element, scvf);
161 if (bcTypes.hasOnlyDirichlet())
163 const auto dirichletPriVars = elementSolution<FVElementGeometry>(problem.dirichlet(element, scvf));
164 auto&& scvI = fvGeometry.scv(scvf.insideScvIdx());
166 VolumeVariables volVars;
167 volVars.update(dirichletPriVars,
172 boundaryVolumeVariables_.emplace_back(std::move(volVars));
173 boundaryVolVarIndices_.push_back(scvf.outsideScvIdx());
179 template<
class FVElementGeometry,
class SolutionVector>
180 void bindElement_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
181 const FVElementGeometry& fvGeometry,
182 const SolutionVector& sol)
185 const GridVolumeVariables* gridVolVarsPtr_;
188 int getLocalIdx_(
const int volVarIdx)
const
190 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
191 assert(it != boundaryVolVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
195 std::vector<std::size_t> boundaryVolVarIndices_;
196 std::vector<VolumeVariables> boundaryVolumeVariables_;
197 const std::size_t numScv_;
216 : gridVolVarsPtr_(&gridVolVars) {}
223 template<
class FVElementGeometry,
class SolutionVector>
225 const FVElementGeometry& fvGeometry,
226 const SolutionVector& sol) &&
228 this->bind_(element, fvGeometry, sol);
229 return std::move(*
this);
232 template<
class FVElementGeometry,
class SolutionVector>
233 void bind(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
234 const FVElementGeometry& fvGeometry,
235 const SolutionVector& sol) &
236 { this->bind_(element, fvGeometry, sol); }
243 template<
class FVElementGeometry,
class SolutionVector>
245 const FVElementGeometry& fvGeometry,
246 const SolutionVector& sol) &&
248 this->bindElement_(element, fvGeometry, sol);
249 return std::move(*
this);
252 template<
class FVElementGeometry,
class SolutionVector>
253 void bindElement(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
254 const FVElementGeometry& fvGeometry,
255 const SolutionVector& sol) &
256 { this->bindElement_(element, fvGeometry, sol); }
259 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
261 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
264 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
266 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
270 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
274 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
278 {
return *gridVolVarsPtr_; }
284 volVarIndices_.clear();
285 volumeVariables_.clear();
289 template<
class FVElementGeometry,
class SolutionVector>
290 void bind_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
291 const FVElementGeometry& fvGeometry,
292 const SolutionVector& sol)
296 const auto& problem = gridVolVars().problem();
297 const auto& gridGeometry = fvGeometry.gridGeometry();
298 const auto globalI = gridGeometry.elementMapper().index(element);
299 const auto& connectivityMapI = gridGeometry.connectivityMap()[globalI];
300 const auto numDofs = connectivityMapI.size() + 1;
303 volumeVariables_.resize(numDofs);
304 volVarIndices_.resize(numDofs);
308 auto&& scvI = fvGeometry.scv(globalI);
309 volumeVariables_[localIdx].update(
elementSolution(element, sol, gridGeometry),
313 volVarIndices_[localIdx] = scvI.dofIndex();
317 for (
const auto& dataJ : connectivityMapI)
319 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
320 auto&& scvJ = fvGeometry.scv(dataJ.globalJ);
321 volumeVariables_[localIdx].update(
elementSolution(elementJ, sol, gridGeometry),
325 volVarIndices_[localIdx] = scvJ.dofIndex();
329 if (fvGeometry.hasBoundaryScvf())
332 for (
auto&& scvf : scvfs(fvGeometry))
335 if (!scvf.boundary())
339 const auto bcTypes = problem.boundaryTypes(element, scvf);
340 if (bcTypes.hasOnlyDirichlet())
342 const auto dirichletPriVars = elementSolution<FVElementGeometry>(problem.dirichlet(element, scvf));
344 volumeVariables_.resize(localIdx+1);
345 volVarIndices_.resize(localIdx+1);
346 volumeVariables_[localIdx].update(dirichletPriVars,
350 volVarIndices_[localIdx] = scvf.outsideScvIdx();
378 template<
class FVElementGeometry,
class SolutionVector>
379 void bindElement_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
380 const FVElementGeometry& fvGeometry,
381 const SolutionVector& sol)
385 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
386 volumeVariables_.resize(1);
387 volVarIndices_.resize(1);
390 auto&& scv = fvGeometry.scv(eIdx);
391 volumeVariables_[0].update(
elementSolution(element, sol, fvGeometry.gridGeometry()),
392 gridVolVars().problem(),
395 volVarIndices_[0] = scv.dofIndex();
398 const GridVolumeVariables* gridVolVarsPtr_;
401 int getLocalIdx_(
const int volVarIdx)
const
403 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
404 assert(it != volVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
408 std::vector<std::size_t> volVarIndices_;
409 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:292
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
The local (stencil) volume variables class for cell centered tpfa models.
Definition: cellcentered/tpfa/elementvolumevariables.hh:45
CCTpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/tpfa/elementvolumevariables.hh:63
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:113
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:93
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:122
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:102
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/tpfa/elementvolumevariables.hh:128
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:57
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:60
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:233
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:253
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:212
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:244
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:209
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:224
CCTpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/tpfa/elementvolumevariables.hh:215
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/tpfa/elementvolumevariables.hh:277
The local element solution class for cell-centered methods.