13#ifndef DUMUX_DISCRETIZATION_CVFE_VARIABLES_DEFLECTION_POLICY_HH
14#define DUMUX_DISCRETIZATION_CVFE_VARIABLES_DEFLECTION_POLICY_HH
16#include <dune/common/reservedvector.hh>
18#include <dumux/common/typetraits/localdofs_.hh>
23template<
class MutableVariablesView,
class FVElementGeometry>
26 using Variables =
typename MutableVariablesView::Variables;
27 static constexpr int maxNumLocalDofs = Detail::LocalDofs::maxNumLocalDofs<FVElementGeometry>();
31 const FVElementGeometry& fvGeometry,
32 bool deflectAllVariables)
33 : elementVariables_(elementVariables)
34 , fvGeometry_(fvGeometry)
35 , deflectAll_(deflectAllVariables)
38 for (
const auto& localDof :
localDofs(fvGeometry))
39 origVariables_.push_back(elementVariables_[localDof]);
42 template<
class LocalDof>
43 void store(
const LocalDof& localDof)
47 origVariables_.clear();
48 origVariables_.push_back(elementVariables_[localDof]);
52 template<
class ElementSolution,
class LocalDof,
class Problem>
53 void update(
const ElementSolution& elemSol,
54 const LocalDof& localDof,
55 const Problem& problem)
58 for (
const auto& localDofI :
localDofs(fvGeometry_))
59 elementVariables_[localDofI].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDofI));
62 elementVariables_[localDof].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDof));
66 template<
class LocalDof>
70 elementVariables_[localDof] = origVariables_[0];
74 for (
const auto& localDofI :
localDofs(fvGeometry_))
76 elementVariables_[localDofI] = origVariables_[idx];
83 MutableVariablesView elementVariables_;
84 const FVElementGeometry& fvGeometry_;
85 const bool deflectAll_;
86 Dune::ReservedVector<Variables, maxNumLocalDofs> origVariables_;
Definition: variablesdeflectionpolicy.hh:25
void restore(const LocalDof &localDof)
Definition: variablesdeflectionpolicy.hh:67
void update(const ElementSolution &elemSol, const LocalDof &localDof, const Problem &problem)
Definition: variablesdeflectionpolicy.hh:53
void store(const LocalDof &localDof)
Definition: variablesdeflectionpolicy.hh:43
VariablesDeflectionPolicy(MutableVariablesView elementVariables, const FVElementGeometry &fvGeometry, bool deflectAllVariables)
Definition: variablesdeflectionpolicy.hh:30
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