24#ifndef DUMUX_DISCRETIZATION_STAGGERED_ELEMENT_VOLUMEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_ELEMENT_VOLUMEVARIABLES_HH
31#include <dune/common/exceptions.hh>
41template<
class GVV,
bool cachingEnabled>
62 : gridVolVarsPtr_(&gridVolVars)
63 , numScv_(gridVolVars.problem().gridGeometry().numScv())
67 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
70 if (scv.dofIndex() < numScv_)
71 return gridVolVars().volVars(scv.dofIndex());
73 return boundaryVolumeVariables_[getLocalIdx_(scv.dofIndex())];
81 return gridVolVars().volVars(scvIdx);
83 return boundaryVolumeVariables_[getLocalIdx_(scvIdx)];
88 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
89 void bind(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
90 const FVElementGeometry& fvGeometry,
91 const SolutionVector& sol)
94 bind(element, fvGeometry, sol[FVElementGeometry::GridGeometry::cellCenterIdx()]);
99 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<!isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
100 void bind(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
101 const FVElementGeometry& fvGeometry,
102 const SolutionVector& sol)
104 if (!fvGeometry.hasBoundaryScvf())
108 boundaryVolVarIndices_.reserve(fvGeometry.numScvf());
109 boundaryVolumeVariables_.reserve(fvGeometry.numScvf());
112 for (
auto&& scvf : scvfs(fvGeometry))
115 if (!scvf.boundary())
118 const auto& problem = gridVolVars().problem();
119 auto boundaryPriVars = gridVolVars().getBoundaryPriVars(problem, sol, element, scvf);
120 const auto elemSol = elementSolution<FVElementGeometry>(std::move(boundaryPriVars));
121 auto&& scvI = fvGeometry.scv(scvf.insideScvIdx());
124 volVars.update(elemSol,
129 boundaryVolumeVariables_.emplace_back(std::move(volVars));
130 boundaryVolVarIndices_.push_back(scvf.outsideScvIdx());
135 template<
class FVElementGeometry,
class SolutionVector>
136 void bindElement(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
137 const FVElementGeometry& fvGeometry,
138 const SolutionVector& sol)
143 {
return *gridVolVarsPtr_; }
150 boundaryVolVarIndices_.clear();
151 boundaryVolumeVariables_.clear();
154 const GridVolumeVariables* gridVolVarsPtr_;
157 int getLocalIdx_(
const int volVarIdx)
const
159 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
160 assert(it != boundaryVolVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
164 std::vector<std::size_t> boundaryVolVarIndices_;
165 std::vector<VolumeVariables> boundaryVolumeVariables_;
166 const std::size_t numScv_;
178 using PrimaryVariables =
typename GVV::VolumeVariables::PrimaryVariables;
189 : gridVolVarsPtr_(&gridVolVars) {}
193 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
194 void bind(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
195 const FVElementGeometry& fvGeometry,
196 const SolutionVector& sol)
199 bind(element, fvGeometry, sol[FVElementGeometry::GridGeometry::cellCenterIdx()]);
204 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<!isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
205 void bind(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
206 const FVElementGeometry& fvGeometry,
207 const SolutionVector& sol)
211 const auto& problem = gridVolVars().problem();
212 const auto& gridGeometry = fvGeometry.gridGeometry();
213 const auto globalI = gridGeometry.elementMapper().index(element);
214 const auto& map = gridGeometry.connectivityMap();
215 constexpr auto cellCenterIdx = FVElementGeometry::GridGeometry::cellCenterIdx();
216 const auto& connectivityMapI = map(cellCenterIdx, cellCenterIdx, globalI);
217 const auto numDofs = connectivityMapI.size();
219 auto&& scvI = fvGeometry.scv(globalI);
222 volumeVariables_.resize(numDofs+1);
223 volVarIndices_.resize(numDofs+1);
227 auto doVolVarUpdate = [&](
int globalJ)
229 const auto& elementJ = gridGeometry.element(globalJ);
230 auto&& scvJ = fvGeometry.scv(globalJ);
231 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PrimaryVariables>(sol[globalJ]);
232 volumeVariables_[localIdx].update(elemSol,
236 volVarIndices_[localIdx] = scvJ.dofIndex();
241 doVolVarUpdate(globalI);
244 for (
const auto& globalJ : connectivityMapI)
245 doVolVarUpdate(globalJ);
247 if (fvGeometry.hasBoundaryScvf())
250 for (
auto&& scvf : scvfs(fvGeometry))
253 if (!scvf.boundary())
256 volumeVariables_.resize(localIdx+1);
257 volVarIndices_.resize(localIdx+1);
259 auto boundaryPriVars = gridVolVars().getBoundaryPriVars(problem, sol, element, scvf);
260 auto elemSol = elementSolution<FVElementGeometry>(std::move(boundaryPriVars));
261 volumeVariables_[localIdx].update(elemSol,
265 volVarIndices_[localIdx] = scvf.outsideScvIdx();
273 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
274 void bindElement(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
275 const FVElementGeometry& fvGeometry,
276 const SolutionVector& sol)
279 bindElement(element, fvGeometry, sol[FVElementGeometry::GridGeometry::cellCenterIdx()]);
284 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<!isMultiTypeBlockVector<SolutionVector>::value,
int> = 0>
285 void bindElement(
const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
286 const FVElementGeometry& fvGeometry,
287 const SolutionVector& sol)
291 const auto globalI = fvGeometry.gridGeometry().elementMapper().index(element);
292 volumeVariables_.resize(1);
293 volVarIndices_.resize(1);
296 auto&& scv = fvGeometry.scv(globalI);
298 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PrimaryVariables>(sol[globalI]);
299 volumeVariables_[0].update(elemSol,
300 gridVolVars().problem(),
303 volVarIndices_[0] = scv.dofIndex();
307 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
309 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
312 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value,
int> = 0>
314 {
return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
318 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
322 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
326 {
return *gridVolVarsPtr_; }
332 volVarIndices_.clear();
333 volumeVariables_.clear();
336 const GridVolumeVariables* gridVolVarsPtr_;
338 int getLocalIdx_(
const int volVarIdx)
const
340 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
341 assert(it != volVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
345 std::vector<std::size_t> volVarIndices_;
346 std::vector<VolumeVariables> volumeVariables_;
Type traits to be used with vector types.
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
Base class for the element volume variables vector for the staggered model.
Definition: staggered/freeflow/elementvolumevariables.hh:43
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:58
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: staggered/freeflow/elementvolumevariables.hh:142
GVV GridVolumeVariables
export type of the grid volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:55
StaggeredElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: staggered/freeflow/elementvolumevariables.hh:61
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol)
function to prepare the vol vars within the element
Definition: staggered/freeflow/elementvolumevariables.hh:136
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol)
Definition: staggered/freeflow/elementvolumevariables.hh:89
GVV GridVolumeVariables
export type of the grid volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:182
StaggeredElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: staggered/freeflow/elementvolumevariables.hh:188
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:185
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: staggered/freeflow/elementvolumevariables.hh:325
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol)
Definition: staggered/freeflow/elementvolumevariables.hh:194
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol)
Definition: staggered/freeflow/elementvolumevariables.hh:274
The local element solution class for staggered methods.