24#ifndef DUMUX_DISCRETIZATION_FACECENTERED_ELEMENTVOLUMEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_ELEMENTVOLUMEVARIABLES_HH
40template<
class Gr
idVolumeVariables,
bool cachingEnabled>
51 using GridGeometry = std::decay_t<decltype(std::declval<GVV>().problem().gridGeometry())>;
52 using FVElementGeometry =
typename GridGeometry::LocalView;
53 using SubControlVolume =
typename GridGeometry::SubControlVolume;
63 : gridVolumeVariablesPtr_(&gridVolumeVariables)
64 , numScv_(gridVolumeVariables.problem().gridGeometry().numScv())
70 if (scv.index() < numScv_)
71 return gridVolVars().volVars(scv.index());
73 return boundaryVolumeVariables_[getLocalIdx_(scv.index())];
81 return gridVolVars().volVars(scvIdx);
83 return boundaryVolumeVariables_[getLocalIdx_(scvIdx)];
91 template<
class SolutionVector>
93 const FVElementGeometry& fvGeometry,
94 const SolutionVector& sol) &&
96 this->bind(element, fvGeometry, sol);
97 return std::move(*
this);
102 template<
class SolutionVector>
103 void bind(
const typename FVElementGeometry::Element& element,
104 const FVElementGeometry& fvGeometry,
105 const SolutionVector& sol) &
107 if (!fvGeometry.hasBoundaryScvf())
111 boundaryVolVarIndices_.reserve(fvGeometry.gridGeometry().numBoundaryScv());
112 boundaryVolumeVariables_.reserve(fvGeometry.gridGeometry().numBoundaryScv());
114 for (
const auto& scvf : scvfs(fvGeometry))
116 if (!scvf.boundary() || scvf.isFrontal() || scvf.processorBoundary())
120 const auto& problem = gridVolVars().problem();
121 const auto bcTypes = problem.boundaryTypes(element, scvf);
125 const auto& scvI = fvGeometry.scv(scvFace.insideScvIdx());
126 typename VolumeVariables::PrimaryVariables pv(
127 problem.dirichlet(element, scvFace)[scvI.dofAxis()]
129 const auto dirichletPriVars = elementSolution<FVElementGeometry>(pv);
132 volVars.update(dirichletPriVars, problem, element, scvI);
134 boundaryVolumeVariables_.emplace_back(std::move(volVars));
135 boundaryVolVarIndices_.push_back(scvFace.outsideScvIdx());
138 if (bcTypes.hasDirichlet())
145 if (
const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); orthogonalScvf.boundary())
146 if (problem.boundaryTypes(element, orthogonalScvf).hasDirichlet())
151 assert(boundaryVolumeVariables_.size() == boundaryVolVarIndices_.size());
159 template<
class SolutionVector>
161 const FVElementGeometry& fvGeometry,
162 const SolutionVector& sol) &&
164 this->bindElement(element, fvGeometry, sol);
165 return std::move(*
this);
169 template<
class SolutionVector>
170 void bindElement(
const typename FVElementGeometry::Element& element,
171 const FVElementGeometry& fvGeometry,
172 const SolutionVector& sol) &
178 {
return *gridVolumeVariablesPtr_; }
183 if (scvIdx < numScv_)
187 const auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), scvIdx);
188 return it != boundaryVolVarIndices_.end();
196 boundaryVolVarIndices_.clear();
197 boundaryVolumeVariables_.clear();
201 int getLocalIdx_(
const std::size_t volVarIdx)
const
203 const auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
204 assert(it != boundaryVolVarIndices_.end() &&
"Could not find the current volume variables for volVarIdx!");
208 const GridVolumeVariables* gridVolumeVariablesPtr_;
209 const std::size_t numScv_;
210 std::vector<std::size_t> boundaryVolVarIndices_;
211 std::vector<VolumeVariables> boundaryVolumeVariables_;
222 using GridGeometry = std::decay_t<decltype(std::declval<GVV>().problem().gridGeometry())>;
223 using FVElementGeometry =
typename GridGeometry::LocalView;
224 using SubControlVolume =
typename GridGeometry::SubControlVolume;
226 static constexpr auto dim = GridGeometry::GridView::dimension;
227 static constexpr auto numInsideVolVars = dim * 2;
228 static constexpr auto numOutsideVolVars = numInsideVolVars * 2 * (dim - 1);
238 : gridVolumeVariablesPtr_(&globalFacesVars) {}
242 {
return volumeVariables_[getLocalIdx_(scv.index())]; }
246 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
250 {
return volumeVariables_[getLocalIdx_(scv.index())]; }
254 {
return volumeVariables_[getLocalIdx_(scvIdx)]; }
261 template<
class SolutionVector>
263 const FVElementGeometry& fvGeometry,
264 const SolutionVector& sol) &&
266 this->bind_(element, fvGeometry, sol);
267 return std::move(*
this);
270 template<
class SolutionVector>
271 void bind(
const typename FVElementGeometry::Element& element,
272 const FVElementGeometry& fvGeometry,
273 const SolutionVector& sol) &
275 this->bind_(element, fvGeometry, sol);
283 template<
class SolutionVector>
285 const FVElementGeometry& fvGeometry,
286 const SolutionVector& sol) &&
288 this->bindElement_(element, fvGeometry, sol);
289 return std::move(*
this);
292 template<
class SolutionVector>
293 void bindElement(
const typename FVElementGeometry::Element& element,
294 const FVElementGeometry& fvGeometry,
295 const SolutionVector& sol) &
296 { this->bindElement_(element, fvGeometry, sol); }
300 {
return *gridVolumeVariablesPtr_; }
304 {
return volVarsInserted_(scvIdx); }
309 template<
class SolutionVector>
310 void bind_(
const typename FVElementGeometry::Element& element,
311 const FVElementGeometry& fvGeometry,
312 const SolutionVector& sol)
316 const auto& problem = gridVolVars().problem();
317 const auto& gridGeometry = fvGeometry.gridGeometry();
319 volVarIndices_.reserve(numInsideVolVars + numInsideVolVars);
320 volumeVariables_.reserve(numInsideVolVars + numInsideVolVars);
322 for (
const auto& scv : scvs(fvGeometry))
324 for (
const auto otherScvIdx : gridGeometry.connectivityMap()[scv.index()])
326 if (!volVarsInserted_(otherScvIdx))
328 const auto& otherScv = fvGeometry.scv(otherScvIdx);
329 volVarIndices_.push_back(otherScvIdx);
330 volumeVariables_.emplace_back();
331 const auto& otherElement = gridGeometry.element(otherScv.elementIndex());
332 volumeVariables_.back().update(
334 problem, otherElement, otherScv
340 if (fvGeometry.hasBoundaryScvf())
342 for (
const auto& scvf : scvfs(fvGeometry))
344 if (!scvf.boundary() || scvf.isFrontal())
348 const auto& problem = gridVolVars().problem();
349 const auto bcTypes = problem.boundaryTypes(element, scvf);
353 const auto& scvI = fvGeometry.scv(scvFace.insideScvIdx());
354 typename VolumeVariables::PrimaryVariables pv(
355 problem.dirichlet(element, scvFace)[scvI.dofAxis()]
357 const auto dirichletPriVars = elementSolution<FVElementGeometry>(pv);
359 VolumeVariables volVars;
360 volVars.update(dirichletPriVars,
365 volumeVariables_.emplace_back(std::move(volVars));
366 volVarIndices_.push_back(scvFace.outsideScvIdx());
369 if (bcTypes.hasDirichlet())
376 if (
const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); orthogonalScvf.boundary())
377 if (problem.boundaryTypes(element, orthogonalScvf).hasDirichlet())
386 template<
class SolutionVector>
387 void bindElement_(
const typename FVElementGeometry::Element& element,
388 const FVElementGeometry& fvGeometry,
389 const SolutionVector& sol)
392 const auto& problem = gridVolVars().problem();
393 const auto& gridGeometry = fvGeometry.gridGeometry();
394 volVarIndices_.reserve(numInsideVolVars);
396 for (
const auto& scv : scvs(fvGeometry))
398 volVarIndices_.push_back(scv.index());
399 volumeVariables_.emplace_back();
400 volumeVariables_.back().update(
402 problem, element, scv
410 volVarIndices_.clear();
411 volumeVariables_.clear();
414 bool volVarsInserted_(
const std::size_t scvIdx)
const
416 return std::find(volVarIndices_.begin(), volVarIndices_.end(), scvIdx) != volVarIndices_.end();
419 int getLocalIdx_(
const int scvfIdx)
const
421 const auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), scvfIdx);
422 assert(it != volVarIndices_.end() &&
"Could not find the current face variables for scvfIdx!");
426 const GridVolumeVariables* gridVolumeVariablesPtr_;
427 std::vector<std::size_t> volVarIndices_;
428 std::vector<VolumeVariables> volumeVariables_;
Element solution classes and factory functions.
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
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:136
Base class for the face variables vector.
Definition: facecentered/staggered/elementvolumevariables.hh:41
Class for the face variables vector. Specialization for the case of storing the face variables global...
Definition: facecentered/staggered/elementvolumevariables.hh:49
void bindElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:170
GVV GridVolumeVariables
export type of the grid volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:57
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: facecentered/staggered/elementvolumevariables.hh:177
ThisType bindElement(const typename FVElementGeometry::Element &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: facecentered/staggered/elementvolumevariables.hh:160
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:60
bool hasVolVars(const std::size_t scvIdx) const
Returns true if volVars exist for the given scv index.
Definition: facecentered/staggered/elementvolumevariables.hh:181
FaceCenteredStaggeredElementVolumeVariables(const GridVolumeVariables &gridVolumeVariables)
Definition: facecentered/staggered/elementvolumevariables.hh:62
void bind(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:103
ThisType bind(const typename FVElementGeometry::Element &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: facecentered/staggered/elementvolumevariables.hh:92
Class for the face variables vector. Specialization for the case of not storing the face variables gl...
Definition: facecentered/staggered/elementvolumevariables.hh:220
GVV GridVolumeVariables
export type of the grid volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:232
ThisType bindElement(const typename FVElementGeometry::Element &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: facecentered/staggered/elementvolumevariables.hh:284
void bind(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:271
ThisType bind(const typename FVElementGeometry::Element &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: facecentered/staggered/elementvolumevariables.hh:262
bool hasVolVars(const std::size_t scvIdx) const
Returns true if volVars exist for the given scv index.
Definition: facecentered/staggered/elementvolumevariables.hh:303
void bindElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:293
FaceCenteredStaggeredElementVolumeVariables(const GridVolumeVariables &globalFacesVars)
Definition: facecentered/staggered/elementvolumevariables.hh:237
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:235
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: facecentered/staggered/elementvolumevariables.hh:299