12#ifndef DUMUX_CVFE_ELEMENT_SOLUTION_HH
13#define DUMUX_CVFE_ELEMENT_SOLUTION_HH
16#include <dune/common/reservedvector.hh>
18#include <dumux/common/typetraits/localdofs_.hh>
26template<
class ElementVariables>
29 if constexpr (
requires {
typename ElementVariables::Variables::PrimaryVariables; })
30 return std::type_identity<typename ElementVariables::Variables::PrimaryVariables>{};
31 else if constexpr (
requires {
typename ElementVariables::VolumeVariables::PrimaryVariables; })
32 return std::type_identity<typename ElementVariables::VolumeVariables::PrimaryVariables>{};
34 return std::type_identity<void>{};
37template<
class ElementVariables>
46template<
class FVElementGeometry,
class PV>
49 using GridGeometry =
typename FVElementGeometry::GridGeometry;
50 using GridView =
typename GridGeometry::GridView;
51 using Element =
typename GridView::template Codim<0>::Entity;
53 static constexpr int dim = GridView::dimension;
56 static constexpr int maxNumLocalDofs = Detail::LocalDofs::maxNumLocalDofs<FVElementGeometry>();
66 template<
class SolutionVector>
68 const GridGeometry& gridGeometry)
70 update(element, sol, gridGeometry);
74 template<
class ElementVariables>
76 const FVElementGeometry& fvGeometry)
78 priVars_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry));
79 for (
const auto& localDof :
localDofs(fvGeometry))
80 priVars_[localDof.index()] = elemVars[localDof.index()].priVars();
84 template<
class SolutionVector>
85 void update(
const Element& element,
const SolutionVector& sol,
86 const GridGeometry& gridGeometry)
91 const auto& localCoeff = gridGeometry.feCache().get(element.type()).localCoefficients();
92 priVars_.resize(localCoeff.size());
93 for (
int localDofIdx = 0; localDofIdx < localCoeff.size(); ++localDofIdx)
95 const auto& localKey = localCoeff.localKey(localDofIdx);
96 priVars_[localDofIdx] = sol[
97 gridGeometry.dofMapper().subIndex(
98 element, localKey.subEntity(), localKey.codim()
105 template<
class SolutionVector>
106 void update(
const Element& element,
const SolutionVector& sol,
107 const FVElementGeometry& fvGeometry)
109 priVars_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry));
110 for (
const auto& localDof :
localDofs(fvGeometry))
111 priVars_[localDof.index()] = sol[localDof.dofIndex()];
116 {
return priVars_.size(); }
119 template<
typename IndexType>
121 {
return priVars_[i]; }
124 template<
typename IndexType>
126 {
return priVars_[i]; }
129 Dune::ReservedVector<PrimaryVariables, maxNumLocalDofs> priVars_;
133template<
class Element,
class SolutionVector,
class Gr
idGeometry>
134auto elementSolution(
const Element& element,
const SolutionVector& sol,
const GridGeometry& gg)
135-> std::enable_if_t<DiscretizationMethods::isCVFE<typename GridGeometry::DiscretizationMethod>,
136 CVFEElementSolution<
typename GridGeometry::LocalView,
137 std::decay_t<decltype(std::declval<SolutionVector>()[0])>>
140 using PrimaryVariables = std::decay_t<decltype(std::declval<SolutionVector>()[0])>;
141 return CVFEElementSolution<typename GridGeometry::LocalView, PrimaryVariables>(element, sol, gg);
145template<
class Element,
class ElementVariables,
class FVElementGeometry>
146auto elementSolution(
const Element& element,
const ElementVariables& elemVars,
const FVElementGeometry& gg)
147-> std::enable_if_t<DiscretizationMethods::isCVFE<typename FVElementGeometry::GridGeometry::DiscretizationMethod>
148 && !std::is_void_v<Detail::PrimaryVariables_t<ElementVariables>>,
The element solution vector.
Definition: cvfe/elementsolution.hh:48
CVFEElementSolution(const Element &element, const ElementVariables &elemVars, const FVElementGeometry &fvGeometry)
Constructor with element and elemVolVars and fvGeometry.
Definition: cvfe/elementsolution.hh:75
const PrimaryVariables & operator[](IndexType i) const
bracket operator const access
Definition: cvfe/elementsolution.hh:120
PV PrimaryVariables
export the primary variables type
Definition: cvfe/elementsolution.hh:60
CVFEElementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
Constructor with element and solution and grid geometry.
Definition: cvfe/elementsolution.hh:67
void update(const Element &element, const SolutionVector &sol, const FVElementGeometry &fvGeometry)
extract the element solution from the solution vector using a local fv geometry
Definition: cvfe/elementsolution.hh:106
std::size_t size() const
return the size of the element solution
Definition: cvfe/elementsolution.hh:115
void update(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
extract the element solution from the solution vector using a mapper
Definition: cvfe/elementsolution.hh:85
CVFEElementSolution()=default
Default constructor.
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
Class representing dofs on elements for control-volume finite element schemes.
The available discretization methods in Dumux.
constexpr auto primaryVariablesType()
Definition: cvfe/elementsolution.hh:27
typename decltype(primaryVariablesType< ElementVariables >())::type PrimaryVariables_t
Definition: cvfe/elementsolution.hh:38
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50