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 <cstddef>
17#include <type_traits>
18#include <utility>
19
20#include <dune/common/reservedvector.hh>
21
22#include <dumux/common/typetraits/localdofs_.hh>
24
25namespace Dumux::Detail::CVFE {
26
27template<class MutableVariablesView, class FVElementGeometry>
29{
30 using Variables = typename MutableVariablesView::Variables;
31 static constexpr int maxNumLocalDofs = Detail::LocalDofs::maxNumLocalDofs<FVElementGeometry>();
32
33public:
34 VariablesDeflectionPolicy(MutableVariablesView elementVariables,
35 const FVElementGeometry& fvGeometry,
36 bool deflectAllVariables)
37 : elementVariables_(elementVariables)
38 , fvGeometry_(fvGeometry)
39 , deflectAll_(deflectAllVariables)
40 {
41 if (deflectAll_)
42 for (const auto& localDof : localDofs(fvGeometry))
43 origVariables_.push_back(elementVariables_[localDof]);
44 }
45
46 template<class LocalDof>
47 void store(const LocalDof& localDof)
48 {
49 if (!deflectAll_)
50 {
51 origVariables_.clear();
52 origVariables_.push_back(elementVariables_[localDof]);
53 }
54 }
55
56 template<class ElementSolution, class LocalDof, class Problem>
57 void update(const ElementSolution& elemSol,
58 const LocalDof& localDof,
59 const Problem& problem)
60 {
61 if (deflectAll_)
62 for (const auto& localDofI : localDofs(fvGeometry_))
63 elementVariables_[localDofI].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDofI));
64 else
65 {
66 elementVariables_[localDof].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDof));
67 }
68 }
69
70 template<class LocalDof>
71 void restore(const LocalDof& localDof)
72 {
73 if (!deflectAll_)
74 elementVariables_[localDof] = origVariables_[0];
75 else
76 {
77 int idx = 0;
78 for (const auto& localDofI : localDofs(fvGeometry_))
79 {
80 elementVariables_[localDofI] = origVariables_[idx];
81 idx++;
82 }
83 }
84 }
85
86private:
87 MutableVariablesView elementVariables_;
88 const FVElementGeometry& fvGeometry_;
89 const bool deflectAll_;
90 Dune::ReservedVector<Variables, maxNumLocalDofs> origVariables_;
91};
92
93template<class MutableVariablesView, class FVElementGeometry>
95{
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>();
99
100public:
101 VariablesDeflectionPolicyWithIpCacheUpdate(MutableVariablesView elementVariables,
102 const FVElementGeometry& fvGeometry,
103 bool deflectAllVariables)
104 : elementVariables_(elementVariables)
105 , fvGeometry_(fvGeometry)
106 , deflectAll_(deflectAllVariables)
107 {
108 if (deflectAll_)
109 {
110 for (const auto& localDof : localDofs(fvGeometry))
111 origVariables_.push_back(elementVariables_[localDof]);
112 }
113 origCache_ = elementVariables_.cache(elementIndex_());
114 }
115
116 template<class LocalDof>
117 void store(const LocalDof& localDof)
118 {
119 if (!deflectAll_)
120 {
121 origVariables_.clear();
122 origVariables_.push_back(elementVariables_[localDof]);
123 }
124 }
125
126 template<class ElementSolution, class LocalDof, class Problem>
127 void update(const ElementSolution& elemSol,
128 const LocalDof& localDof,
129 const Problem& problem)
130 {
131 if (deflectAll_)
132 for (const auto& localDofI : localDofs(fvGeometry_))
133 elementVariables_[localDofI].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDofI));
134 else
135 elementVariables_[localDof].update(elemSol, problem, fvGeometry_, ipData(fvGeometry_, localDof));
136
137 elementVariables_.cache(elementIndex_()).update(problem, fvGeometry_.element(), fvGeometry_, elementVariables_);
138 }
139
140 template<class LocalDof>
141 void restore(const LocalDof& localDof)
142 {
143 if (!deflectAll_)
144 elementVariables_[localDof] = origVariables_[0];
145 else
146 {
147 int idx = 0;
148 for (const auto& localDofI : localDofs(fvGeometry_))
149 {
150 elementVariables_[localDofI] = origVariables_[idx];
151 idx++;
152 }
153 }
154
155 elementVariables_.cache(elementIndex_()) = origCache_;
156 }
157
158private:
159 const auto elementIndex_() const
160 { return fvGeometry_.elementIndex(); }
161
162 MutableVariablesView elementVariables_;
163 const FVElementGeometry& fvGeometry_;
164 const bool deflectAll_;
165 Dune::ReservedVector<Variables, maxNumLocalDofs> origVariables_;
166 ElementCache origCache_;
167};
168
169} // end namespace Dumux::Detail::CVFE
170
171#endif
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