24#ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_ELEMENT_SOLUTION_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_ELEMENT_SOLUTION_HH
40template<
class FVElementGeometry,
class PV>
43 using GridGeometry =
typename FVElementGeometry::GridGeometry;
44 using GridView =
typename GridGeometry::GridView;
45 using Element =
typename GridView::template Codim<0>::Entity;
48 static constexpr auto numScvsPerElement = GridView::Grid::dimension * 2;
57 template<
class SolutionVector>
59 const SolutionVector& sol,
60 const GridGeometry& gridGeometry)
62 const auto fvGeometry =
localView(gridGeometry).bindElement(element);
64 for (
const auto& scv : scvs(fvGeometry))
65 priVars_[scv.indexInElement()] = sol[scv.dofIndex()];
69 template<
class ElementVolumeVariables>
71 const ElementVolumeVariables& elemVolVars,
72 const FVElementGeometry& fvGeometry)
74 for (
const auto& scv : scvs(fvGeometry))
75 priVars_ = elemVolVars[scv].priVars();
81 priVars_[0] = std::move(priVars);
82 for (
int i = 1; i < priVars_.size(); ++i)
83 priVars_[i] = priVars_[0];
89 priVars_[0] = priVars;
90 for (
int i = 1; i < priVars_.size(); ++i)
91 priVars_[i] = priVars_[0];
95 template<
class SolutionVector>
96 void update(
const Element& element,
const SolutionVector& sol,
97 const GridGeometry& gridGeometry)
99 const auto fvGeometry =
localView(gridGeometry).bindElement(element);
101 for (
const auto& scv : scvs(fvGeometry))
102 priVars_[scv.indexInElement()] = sol[scv.dofIndex()];
108 return priVars_[localScvIdx];
114 return priVars_[localScvIdx];
118 static constexpr std::size_t
size()
119 {
return numScvsPerElement; }
122 std::array<PrimaryVariables, numScvsPerElement> priVars_;
129template<
class Element,
class SolutionVector,
class Gr
idGeometry>
130auto elementSolution(
const Element& element,
const SolutionVector& sol,
const GridGeometry& gg)
133 FaceCenteredStaggeredElementSolution<
134 typename GridGeometry::LocalView,
135 std::decay_t<decltype(std::declval<SolutionVector>()[0])>
138{
return { element, sol, gg }; }
144template<
class Element,
class ElementVolumeVariables,
class FVElementGeometry>
145auto elementSolution(
const Element& element,
const ElementVolumeVariables& elemVolVars,
const FVElementGeometry& gg)
148 FaceCenteredStaggeredElementSolution<
150 typename ElementVolumeVariables::VolumeVariables::PrimaryVariables
153{
return {
element, elemVolVars, gg }; }
160template<
class FVElementGeometry,
class PrimaryVariables>
164 FaceCenteredStaggeredElementSolution<
166 std::decay_t<PrimaryVariables>
169{
return { std::forward<PrimaryVariables>(priVars) }; }
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
constexpr FCStaggered fcstaggered
Definition: method.hh:142
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
The global face variables class for staggered models.
Definition: facecentered/staggered/elementsolution.hh:42
FaceCenteredStaggeredElementSolution(PrimaryVariables &&priVars)
Constructor with a primary variable object.
Definition: facecentered/staggered/elementsolution.hh:79
FaceCenteredStaggeredElementSolution(const Element &element, const ElementVolumeVariables &elemVolVars, const FVElementGeometry &fvGeometry)
Constructor with element, element volume variables and fv element geometry.
Definition: facecentered/staggered/elementsolution.hh:70
FaceCenteredStaggeredElementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
Constructor with element, solution vector and grid geometry.
Definition: facecentered/staggered/elementsolution.hh:58
FaceCenteredStaggeredElementSolution()=default
PV PrimaryVariables
export the primary variables type
Definition: facecentered/staggered/elementsolution.hh:52
const PrimaryVariables & operator[](SmallLocalIndexType localScvIdx) const
bracket operator const access
Definition: facecentered/staggered/elementsolution.hh:106
static constexpr std::size_t size()
return the size of the element solution
Definition: facecentered/staggered/elementsolution.hh:118
FaceCenteredStaggeredElementSolution(const PrimaryVariables &priVars)
Constructor with a primary variable object.
Definition: facecentered/staggered/elementsolution.hh:87
void update(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
extract the element solution from the solution vector using a mapper
Definition: facecentered/staggered/elementsolution.hh:96