12#ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_ELEMENT_SOLUTION_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_ELEMENT_SOLUTION_HH
28template<
class FVElementGeometry,
class PV>
31 using GridGeometry =
typename FVElementGeometry::GridGeometry;
32 using GridView =
typename GridGeometry::GridView;
33 using Element =
typename GridView::template Codim<0>::Entity;
36 static constexpr auto numScvsPerElement = GridView::Grid::dimension * 2;
45 template<
class SolutionVector>
47 const SolutionVector& sol,
48 const GridGeometry& gridGeometry)
50 const auto fvGeometry =
localView(gridGeometry).bindElement(element);
52 for (
const auto& scv : scvs(fvGeometry))
53 priVars_[scv.indexInElement()] = sol[scv.dofIndex()];
57 template<
class ElementVolumeVariables>
59 const ElementVolumeVariables& elemVolVars,
60 const FVElementGeometry& fvGeometry)
62 for (
const auto& scv : scvs(fvGeometry))
63 priVars_ = elemVolVars[scv].priVars();
69 priVars_[0] = std::move(priVars);
70 for (
int i = 1; i < priVars_.size(); ++i)
71 priVars_[i] = priVars_[0];
77 priVars_[0] = priVars;
78 for (
int i = 1; i < priVars_.size(); ++i)
79 priVars_[i] = priVars_[0];
83 template<
class SolutionVector>
84 void update(
const Element& element,
const SolutionVector& sol,
85 const GridGeometry& gridGeometry)
87 const auto fvGeometry =
localView(gridGeometry).bindElement(element);
89 for (
const auto& scv : scvs(fvGeometry))
90 priVars_[scv.indexInElement()] = sol[scv.dofIndex()];
96 return priVars_[localScvIdx];
102 return priVars_[localScvIdx];
106 static constexpr std::size_t
size()
107 {
return numScvsPerElement; }
110 std::array<PrimaryVariables, numScvsPerElement> priVars_;
117template<
class Element,
class SolutionVector,
class Gr
idGeometry>
118auto elementSolution(
const Element& element,
const SolutionVector& sol,
const GridGeometry& gg)
121 FaceCenteredStaggeredElementSolution<
122 typename GridGeometry::LocalView,
123 std::decay_t<decltype(std::declval<SolutionVector>()[0])>
126{
return { element, sol, gg }; }
132template<
class Element,
class ElementVolumeVariables,
class FVElementGeometry>
133auto elementSolution(
const Element& element,
const ElementVolumeVariables& elemVolVars,
const FVElementGeometry& gg)
136 FaceCenteredStaggeredElementSolution<
138 typename ElementVolumeVariables::VolumeVariables::PrimaryVariables
141{
return {
element, elemVolVars, gg }; }
148template<
class FVElementGeometry,
class PrimaryVariables>
152 FaceCenteredStaggeredElementSolution<
154 std::decay_t<PrimaryVariables>
157{
return { std::forward<PrimaryVariables>(priVars) }; }
The global face variables class for staggered models.
Definition: facecentered/staggered/elementsolution.hh:30
FaceCenteredStaggeredElementSolution(PrimaryVariables &&priVars)
Constructor with a primary variable object.
Definition: facecentered/staggered/elementsolution.hh:67
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:58
FaceCenteredStaggeredElementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
Constructor with element, solution vector and grid geometry.
Definition: facecentered/staggered/elementsolution.hh:46
FaceCenteredStaggeredElementSolution()=default
PV PrimaryVariables
export the primary variables type
Definition: facecentered/staggered/elementsolution.hh:40
const PrimaryVariables & operator[](SmallLocalIndexType localScvIdx) const
bracket operator const access
Definition: facecentered/staggered/elementsolution.hh:94
static constexpr std::size_t size()
return the size of the element solution
Definition: facecentered/staggered/elementsolution.hh:106
FaceCenteredStaggeredElementSolution(const PrimaryVariables &priVars)
Constructor with a primary variable object.
Definition: facecentered/staggered/elementsolution.hh:75
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:84
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
The available discretization methods in Dumux.
constexpr FCStaggered fcstaggered
Definition: method.hh:151
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:29