24#ifndef DUMUX_DISCRETIZATION_STAGGERED_ELEMENT_VOLUMEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_ELEMENT_VOLUMEVARIABLES_HH
32#include <dune/common/exceptions.hh>
42template<
class GVV,
bool cachingEnabled>
63 : gridVolVarsPtr_(&gridVolVars)
64 , numScv_(gridVolVars.problem().gridGeometry().numScv())
68 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
71 if (scv.dofIndex() < numScv_)
72 return gridVolVars().volVars(scv.dofIndex());
74 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();
142 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
143 void bind_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
144 const FVElementGeometry& fvGeometry,
145 const SolutionVector& sol)
148 bind_(element, fvGeometry, sol[FVElementGeometry::GridGeometry::cellCenterIdx()]);
153 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<!isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
154 void bind_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
155 const FVElementGeometry& fvGeometry,
156 const SolutionVector& sol)
158 if (!fvGeometry.hasBoundaryScvf())
162 boundaryVolVarIndices_.reserve(fvGeometry.numScvf());
163 boundaryVolumeVariables_.reserve(fvGeometry.numScvf());
166 for (
auto&& scvf : scvfs(fvGeometry))
169 if (!scvf.boundary())
172 const auto& problem = gridVolVars().problem();
173 auto boundaryPriVars = gridVolVars().getBoundaryPriVars(problem, sol, element, scvf);
174 const auto elemSol = elementSolution<FVElementGeometry>(std::move(boundaryPriVars));
175 auto&& scvI = fvGeometry.scv(scvf.insideScvIdx());
177 VolumeVariables volVars;
178 volVars.update(elemSol,
183 boundaryVolumeVariables_.emplace_back(std::move(volVars));
184 boundaryVolVarIndices_.push_back(scvf.outsideScvIdx());
189 template<
class FVElementGeometry,
class SolutionVector>
190 void bindElement_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
191 const FVElementGeometry& fvGeometry,
192 const SolutionVector& sol)
195 const GridVolumeVariables* gridVolVarsPtr_;
198 int getLocalIdx_(
const int volVarIdx)
const
200 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
201 assert(it != boundaryVolVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
205 std::vector<std::size_t> boundaryVolVarIndices_;
206 std::vector<VolumeVariables> boundaryVolumeVariables_;
207 const std::size_t numScv_;
219 using PrimaryVariables =
typename GVV::VolumeVariables::PrimaryVariables;
230 : gridVolVarsPtr_(&gridVolVars) {}
237 template<
class FVElementGeometry,
class SolutionVector>
239 const FVElementGeometry& fvGeometry,
240 const SolutionVector& sol) &&
242 this->bind_(element, fvGeometry, sol);
243 return std::move(*
this);
246 template<
class FVElementGeometry,
class SolutionVector>
247 void bind(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
248 const FVElementGeometry& fvGeometry,
249 const SolutionVector& sol) &
250 { this->bind_(element, fvGeometry, sol); }
257 template<
class FVElementGeometry,
class SolutionVector>
259 const FVElementGeometry& fvGeometry,
260 const SolutionVector& sol) &&
262 this->bindElement_(element, fvGeometry, sol);
263 return std::move(*
this);
266 template<
class FVElementGeometry,
class SolutionVector>
267 void bindElement(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
268 const FVElementGeometry& fvGeometry,
269 const SolutionVector& sol) &
270 { this->bindElement_(element, fvGeometry, sol); }
273 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
275 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
278 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
280 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
284 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
288 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
292 {
return *gridVolVarsPtr_; }
297 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
298 void bind_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
299 const FVElementGeometry& fvGeometry,
300 const SolutionVector& sol)
303 bind_(element, fvGeometry, sol[FVElementGeometry::GridGeometry::cellCenterIdx()]);
308 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<!isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
309 void bind_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
310 const FVElementGeometry& fvGeometry,
311 const SolutionVector& sol)
315 const auto& problem = gridVolVars().problem();
316 const auto& gridGeometry = fvGeometry.gridGeometry();
317 const auto globalI = gridGeometry.elementMapper().index(element);
318 const auto& map = gridGeometry.connectivityMap();
319 constexpr auto cellCenterIdx = FVElementGeometry::GridGeometry::cellCenterIdx();
320 const auto& connectivityMapI = map(cellCenterIdx, cellCenterIdx, globalI);
321 const auto numDofs = connectivityMapI.size();
323 auto&& scvI = fvGeometry.scv(globalI);
326 volumeVariables_.resize(numDofs+1);
327 volVarIndices_.resize(numDofs+1);
331 auto doVolVarUpdate = [&](
int globalJ)
333 const auto& elementJ = gridGeometry.element(globalJ);
334 auto&& scvJ = fvGeometry.scv(globalJ);
335 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PrimaryVariables>(sol[globalJ]);
336 volumeVariables_[localIdx].update(elemSol,
340 volVarIndices_[localIdx] = scvJ.dofIndex();
345 doVolVarUpdate(globalI);
348 for (
const auto& globalJ : connectivityMapI)
349 doVolVarUpdate(globalJ);
351 if (fvGeometry.hasBoundaryScvf())
354 for (
auto&& scvf : scvfs(fvGeometry))
357 if (!scvf.boundary())
360 volumeVariables_.resize(localIdx+1);
361 volVarIndices_.resize(localIdx+1);
363 auto boundaryPriVars = gridVolVars().getBoundaryPriVars(problem, sol, element, scvf);
364 auto elemSol = elementSolution<FVElementGeometry>(std::move(boundaryPriVars));
365 volumeVariables_[localIdx].update(elemSol,
369 volVarIndices_[localIdx] = scvf.outsideScvIdx();
377 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
378 void bindElement_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
379 const FVElementGeometry& fvGeometry,
380 const SolutionVector& sol)
383 bindElement_(element, fvGeometry, sol[FVElementGeometry::GridGeometry::cellCenterIdx()]);
388 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<!isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
389 void bindElement_(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
390 const FVElementGeometry& fvGeometry,
391 const SolutionVector& sol)
395 const auto globalI = fvGeometry.gridGeometry().elementMapper().index(element);
396 volumeVariables_.resize(1);
397 volVarIndices_.resize(1);
400 auto&& scv = fvGeometry.scv(globalI);
402 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PrimaryVariables>(sol[globalI]);
403 volumeVariables_[0].update(elemSol,
404 gridVolVars().problem(),
407 volVarIndices_[0] = scv.dofIndex();
413 volVarIndices_.clear();
414 volumeVariables_.clear();
417 const GridVolumeVariables* gridVolVarsPtr_;
419 int getLocalIdx_(
const int volVarIdx)
const
421 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
422 assert(it != volVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
426 std::vector<std::size_t> volVarIndices_;
427 std::vector<VolumeVariables> volumeVariables_;
Type traits to be used with vector types.
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
Base class for the element volume variables vector for the staggered model.
Definition: staggered/freeflow/elementvolumevariables.hh:44
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:59
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: staggered/freeflow/elementvolumevariables.hh:128
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: staggered/freeflow/elementvolumevariables.hh:122
GVV GridVolumeVariables
export type of the grid volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:56
StaggeredElementVolumeVariables 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: staggered/freeflow/elementvolumevariables.hh:113
StaggeredElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: staggered/freeflow/elementvolumevariables.hh:62
StaggeredElementVolumeVariables 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: staggered/freeflow/elementvolumevariables.hh:93
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: staggered/freeflow/elementvolumevariables.hh:102
GVV GridVolumeVariables
export type of the grid volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:223
StaggeredElementVolumeVariables 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: staggered/freeflow/elementvolumevariables.hh:238
StaggeredElementVolumeVariables 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: staggered/freeflow/elementvolumevariables.hh:258
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: staggered/freeflow/elementvolumevariables.hh:247
StaggeredElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: staggered/freeflow/elementvolumevariables.hh:229
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: staggered/freeflow/elementvolumevariables.hh:267
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:226
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: staggered/freeflow/elementvolumevariables.hh:291
The local element solution class for staggered methods.