version 3.11-dev
variablesdeflectionpolicy.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_DISCRETIZATION_CVFE_VARIABLES_DEFLECTION_POLICY_HH
14#define DUMUX_DISCRETIZATION_CVFE_VARIABLES_DEFLECTION_POLICY_HH
15
16#include <dune/common/reservedvector.hh>
17
18#include <dumux/common/typetraits/localdofs_.hh>
20
21namespace Dumux::Detail::CVFE {
22
23template<class MutableVariablesView, class FVElementGeometry>
25{
26 using Variables = typename MutableVariablesView::Variables;
27 static constexpr int maxNumLocalDofs = Detail::LocalDofs::maxNumLocalDofs<FVElementGeometry>();
28
29public:
30 VariablesDeflectionPolicy(MutableVariablesView elementVariables,
31 const FVElementGeometry& fvGeometry,
32 bool deflectAllVariables)
33 : elementVariables_(elementVariables)
34 , fvGeometry_(fvGeometry)
35 , deflectAll_(deflectAllVariables)
36 {
37 if (deflectAll_)
38 for (const auto& localDof : localDofs(fvGeometry))
39 origVariables_.push_back(elementVariables_[localDof]);
40 }
41
42 template<class LocalDof>
43 void store(const LocalDof& localDof)
44 {
45 if (!deflectAll_)
46 {
47 origVariables_.clear();
48 origVariables_.push_back(elementVariables_[localDof]);
49 }
50 }
51
52 template<class ElementSolution, class LocalDof, class Problem>
53 void update(const ElementSolution& elemSol,
54 const LocalDof& localDof,
55 const Problem& problem)
56 {
57 if (deflectAll_)
58 for (const auto& localDofI : localDofs(fvGeometry_))
59 elementVariables_[localDofI].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDofI));
60 else
61 {
62 elementVariables_[localDof].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDof));
63 }
64 }
65
66 template<class LocalDof>
67 void restore(const LocalDof& localDof)
68 {
69 if (!deflectAll_)
70 elementVariables_[localDof] = origVariables_[0];
71 else
72 {
73 int idx = 0;
74 for (const auto& localDofI : localDofs(fvGeometry_))
75 {
76 elementVariables_[localDofI] = origVariables_[idx];
77 idx++;
78 }
79 }
80 }
81
82private:
83 MutableVariablesView elementVariables_;
84 const FVElementGeometry& fvGeometry_;
85 const bool deflectAll_;
86 Dune::ReservedVector<Variables, maxNumLocalDofs> origVariables_;
87};
88
89} // end namespace Dumux::Detail::CVFE
90
91#endif
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