13#ifndef DUMUX_DISCRETIZATION_CVFE_VARIABLES_DEFLECTION_POLICY_HH
14#define DUMUX_DISCRETIZATION_CVFE_VARIABLES_DEFLECTION_POLICY_HH
20#include <dune/common/reservedvector.hh>
22#include <dumux/common/typetraits/localdofs_.hh>
27template<
class MutableVariablesView,
class FVElementGeometry>
30 using Variables =
typename MutableVariablesView::Variables;
31 static constexpr int maxNumLocalDofs = Detail::LocalDofs::maxNumLocalDofs<FVElementGeometry>();
35 const FVElementGeometry& fvGeometry,
36 bool deflectAllVariables)
37 : elementVariables_(elementVariables)
38 , fvGeometry_(fvGeometry)
39 , deflectAll_(deflectAllVariables)
42 for (
const auto& localDof :
localDofs(fvGeometry))
43 origVariables_.push_back(elementVariables_[localDof]);
46 template<
class LocalDof>
47 void store(
const LocalDof& localDof)
51 origVariables_.clear();
52 origVariables_.push_back(elementVariables_[localDof]);
56 template<
class ElementSolution,
class LocalDof,
class Problem>
57 void update(
const ElementSolution& elemSol,
58 const LocalDof& localDof,
59 const Problem& problem)
62 for (
const auto& localDofI :
localDofs(fvGeometry_))
63 elementVariables_[localDofI].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDofI));
66 elementVariables_[localDof].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDof));
70 template<
class LocalDof>
74 elementVariables_[localDof] = origVariables_[0];
78 for (
const auto& localDofI :
localDofs(fvGeometry_))
80 elementVariables_[localDofI] = origVariables_[idx];
87 MutableVariablesView elementVariables_;
88 const FVElementGeometry& fvGeometry_;
89 const bool deflectAll_;
90 Dune::ReservedVector<Variables, maxNumLocalDofs> origVariables_;
93template<
class MutableVariablesView,
class FVElementGeometry>
96 using Variables =
typename MutableVariablesView::Variables;
97 using ElementCache = std::remove_cvref_t<decltype(std::declval<MutableVariablesView&>().cache(std::size_t{}))>;
98 static constexpr int maxNumLocalDofs = Detail::LocalDofs::maxNumLocalDofs<FVElementGeometry>();
102 const FVElementGeometry& fvGeometry,
103 bool deflectAllVariables)
104 : elementVariables_(elementVariables)
105 , fvGeometry_(fvGeometry)
106 , deflectAll_(deflectAllVariables)
110 for (
const auto& localDof :
localDofs(fvGeometry))
111 origVariables_.push_back(elementVariables_[localDof]);
113 origCache_ = elementVariables_.cache(elementIndex_());
116 template<
class LocalDof>
117 void store(
const LocalDof& localDof)
121 origVariables_.clear();
122 origVariables_.push_back(elementVariables_[localDof]);
126 template<
class ElementSolution,
class LocalDof,
class Problem>
127 void update(
const ElementSolution& elemSol,
128 const LocalDof& localDof,
129 const Problem& problem)
132 for (
const auto& localDofI :
localDofs(fvGeometry_))
133 elementVariables_[localDofI].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDofI));
135 elementVariables_[localDof].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDof));
137 elementVariables_.cache(elementIndex_()).update(problem, fvGeometry_.element(), fvGeometry_, elementVariables_);
140 template<
class LocalDof>
144 elementVariables_[localDof] = origVariables_[0];
148 for (
const auto& localDofI :
localDofs(fvGeometry_))
150 elementVariables_[localDofI] = origVariables_[idx];
155 elementVariables_.cache(elementIndex_()) = origCache_;
159 const auto elementIndex_()
const
160 {
return fvGeometry_.elementIndex(); }
162 MutableVariablesView elementVariables_;
163 const FVElementGeometry& fvGeometry_;
164 const bool deflectAll_;
165 Dune::ReservedVector<Variables, maxNumLocalDofs> origVariables_;
166 ElementCache origCache_;
Definition: variablesdeflectionpolicy.hh:29
void restore(const LocalDof &localDof)
Definition: variablesdeflectionpolicy.hh:71
void update(const ElementSolution &elemSol, const LocalDof &localDof, const Problem &problem)
Definition: variablesdeflectionpolicy.hh:57
void store(const LocalDof &localDof)
Definition: variablesdeflectionpolicy.hh:47
VariablesDeflectionPolicy(MutableVariablesView elementVariables, const FVElementGeometry &fvGeometry, bool deflectAllVariables)
Definition: variablesdeflectionpolicy.hh:34
Definition: variablesdeflectionpolicy.hh:95
void restore(const LocalDof &localDof)
Definition: variablesdeflectionpolicy.hh:141
void store(const LocalDof &localDof)
Definition: variablesdeflectionpolicy.hh:117
VariablesDeflectionPolicyWithIpCacheUpdate(MutableVariablesView elementVariables, const FVElementGeometry &fvGeometry, bool deflectAllVariables)
Definition: variablesdeflectionpolicy.hh:101
void update(const ElementSolution &elemSol, const LocalDof &localDof, const Problem &problem)
Definition: variablesdeflectionpolicy.hh:127
Class representing dofs on elements for control-volume finite element schemes.
Definition: elementvariables.hh:24
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50