Line | Branch | Exec | Source |
---|---|---|---|
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-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup Experimental | ||
10 | * \ingroup Assembly | ||
11 | * \copydoc Dumux::FVLocalAssemblerBase | ||
12 | */ | ||
13 | #ifndef DUMUX_EXPERIMENTAL_FV_LOCAL_ASSEMBLER_BASE_HH | ||
14 | #define DUMUX_EXPERIMENTAL_FV_LOCAL_ASSEMBLER_BASE_HH | ||
15 | |||
16 | #include <dune/common/reservedvector.hh> | ||
17 | #include <dune/grid/common/gridenums.hh> // for GhostEntity | ||
18 | #include <dune/istl/matrixindexset.hh> | ||
19 | |||
20 | #include <dumux/common/reservedblockvector.hh> | ||
21 | #include <dumux/common/properties.hh> | ||
22 | #include <dumux/common/parameters.hh> | ||
23 | #include <dumux/discretization/extrusion.hh> | ||
24 | #include <dumux/assembly/diffmethod.hh> | ||
25 | |||
26 | namespace Dumux::Experimental { | ||
27 | |||
28 | /*! | ||
29 | * \ingroup Experimental | ||
30 | * \ingroup Assembly | ||
31 | * \brief A base class for all local assemblers | ||
32 | * \tparam TypeTag The TypeTag | ||
33 | * \tparam Assembler The assembler type | ||
34 | * \tparam Implementation The assembler implementation | ||
35 | */ | ||
36 | template<class TypeTag, class Assembler, class Implementation> | ||
37 |
2/8✓ Branch 0 taken 10368 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 3072 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10052992 | class FVLocalAssemblerBase |
38 | { | ||
39 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
40 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
41 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
42 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
43 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
44 | using SolutionVector = typename Assembler::SolutionVector; | ||
45 | using ElementBoundaryTypes = GetPropType<TypeTag, Properties::ElementBoundaryTypes>; | ||
46 | using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView; | ||
47 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
48 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
49 | using GridVolumeVariables = GetPropType<TypeTag, Properties::GridVolumeVariables>; | ||
50 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
51 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
52 | using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView; | ||
53 | using Element = typename GridView::template Codim<0>::Entity; | ||
54 | static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq(); | ||
55 | |||
56 | public: | ||
57 | using LocalResidual = std::decay_t<decltype(std::declval<Assembler>().localResidual())>; | ||
58 | using ElementResidualVector = typename LocalResidual::ElementResidualVector; | ||
59 | |||
60 | /*! | ||
61 | * \brief The constructor. Delegates to the general constructor. | ||
62 | */ | ||
63 | 10079104 | explicit FVLocalAssemblerBase(const Assembler& assembler, | |
64 | const Element& element, | ||
65 | const SolutionVector& curSol) | ||
66 | : FVLocalAssemblerBase(assembler, | ||
67 | element, | ||
68 | curSol, | ||
69 | 20158208 | localView(assembler.gridGeometry()), | |
70 |
1/2✓ Branch 1 taken 10065664 times.
✗ Branch 2 not taken.
|
10119424 | localView(assembler.gridVariables().curGridVolVars()), |
71 |
3/6✓ Branch 1 taken 10065664 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 5013056 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 5013056 times.
✗ Branch 8 not taken.
|
20145536 | localView(assembler.gridVariables().gridFluxVarsCache()), |
72 | assembler.localResidual(), | ||
73 | 20158208 | element.partitionType() == Dune::GhostEntity, | |
74 |
9/20✓ Branch 3 taken 10065664 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 10065664 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 10065664 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 10065664 times.
✗ Branch 13 not taken.
✓ Branch 15 taken 10065664 times.
✗ Branch 16 not taken.
✓ Branch 18 taken 5052608 times.
✗ Branch 19 not taken.
✓ Branch 21 taken 5052608 times.
✗ Branch 22 not taken.
✓ Branch 24 taken 5052608 times.
✗ Branch 25 not taken.
✗ Branch 27 not taken.
✓ Branch 28 taken 13440 times.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
|
10079104 | assembler.isImplicit()) |
75 | 10079104 | {} | |
76 | |||
77 | /*! | ||
78 | * \brief The constructor. General version explicitly expecting each argument. | ||
79 | */ | ||
80 | 10098048 | explicit FVLocalAssemblerBase(const Assembler& assembler, | |
81 | const Element& element, | ||
82 | const SolutionVector& curSol, | ||
83 | const FVElementGeometry& fvGeometry, | ||
84 | const ElementVolumeVariables& curElemVolVars, | ||
85 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
86 | const LocalResidual& localResidual, | ||
87 | const bool elementIsGhost, | ||
88 | const bool isImplicit) | ||
89 | : assembler_(assembler) | ||
90 | , element_(element) | ||
91 | , curSol_(curSol) | ||
92 | , fvGeometry_(fvGeometry) | ||
93 | , curElemVolVars_(curElemVolVars) | ||
94 | , elemFluxVarsCache_(elemFluxVarsCache) | ||
95 | , localOperator_(localResidual) | ||
96 | , elementIsGhost_(elementIsGhost) | ||
97 |
4/9✗ Branch 0 not taken.
✓ Branch 1 taken 10026112 times.
✓ Branch 2 taken 49024 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 10026112 times.
✓ Branch 5 taken 49024 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
|
15137984 | , isImplicit_(isImplicit) |
98 | 5071552 | {} | |
99 | |||
100 | /*! | ||
101 | * \brief Convenience function to evaluate the complete local residual for the current element. Automatically chooses the the appropriate | ||
102 | * element volume variables. | ||
103 | */ | ||
104 | 248832 | ElementResidualVector evalLocalResidual() const | |
105 | { | ||
106 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 328704 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 6912 times.
|
335616 | if (elementIsGhost()) |
107 | ✗ | return ElementResidualVector(0.0); | |
108 | |||
109 | 671232 | return evalLocalResidual(curElemVolVars()); | |
110 | } | ||
111 | |||
112 | /*! | ||
113 | * \brief Evaluates the complete local residual for the current element. | ||
114 | * \param elemVolVars The element volume variables | ||
115 | */ | ||
116 | 425472 | ElementResidualVector evalLocalResidual(const ElementVolumeVariables& elemVolVars) const | |
117 | { | ||
118 | 425472 | if (!assembler().isStationaryProblem()) | |
119 | { | ||
120 | 850944 | ElementResidualVector residual = evalLocalFluxAndSourceResidual(elemVolVars); | |
121 | 850944 | residual += evalStorage(); | |
122 | 425472 | return residual; | |
123 | } | ||
124 | else | ||
125 | return evalLocalFluxAndSourceResidual(elemVolVars); | ||
126 | } | ||
127 | |||
128 | /*! | ||
129 | * \brief Convenience function to evaluate the flux and source terms (i.e, the terms without a time derivative) | ||
130 | * of the local residual for the current element. Automatically chooses the the appropriate | ||
131 | * element volume variables. | ||
132 | */ | ||
133 | ElementResidualVector evalLocalFluxAndSourceResidual() const | ||
134 | { | ||
135 | return evalLocalFluxAndSourceResidual(curElemVolVars()); | ||
136 | } | ||
137 | |||
138 | /*! | ||
139 | * \brief Evaluates the flux and source terms (i.e, the terms without a time derivative) | ||
140 | * of the local residual for the current element. | ||
141 | * | ||
142 | * \param elemVolVars The element volume variables | ||
143 | */ | ||
144 | ElementResidualVector evalLocalFluxAndSourceResidual(const ElementVolumeVariables& elemVolVars) const | ||
145 | { | ||
146 |
0/2✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
5363520 | return localOperator_.evalFluxAndSource(element_, fvGeometry_, elemVolVars, elemFluxVarsCache_, elemBcTypes_); |
147 | } | ||
148 | |||
149 | /*! | ||
150 | * \brief Convenience function to evaluate storage term (i.e, the term with a time derivative) | ||
151 | * of the local residual for the current element. Automatically chooses the the appropriate | ||
152 | * element volume variables. | ||
153 | */ | ||
154 | ElementResidualVector evalStorage() const | ||
155 | { | ||
156 | 10335616 | return localOperator_.evalStorage(fvGeometry_, curElemVolVars_); | |
157 | } | ||
158 | |||
159 | /*! | ||
160 | * \brief Convenience function bind and prepare all relevant variables required for the | ||
161 | * evaluation of the local residual. | ||
162 | */ | ||
163 | 10079104 | void bindLocalViews() | |
164 | { | ||
165 | // get some references for convenience | ||
166 | 10079104 | const auto& element = this->element(); | |
167 | 10079104 | const auto& curSol = this->curSol(); | |
168 | 20158208 | auto&& fvGeometry = this->fvGeometry(); | |
169 | 20158208 | auto&& curElemVolVars = this->curElemVolVars(); | |
170 | 20158208 | auto&& elemFluxVarsCache = this->elemFluxVarsCache(); | |
171 | |||
172 | // bind the caches | ||
173 | 20158208 | fvGeometry.bind(element); | |
174 |
2/2✓ Branch 0 taken 10026112 times.
✓ Branch 1 taken 52992 times.
|
10079104 | if (std::abs(this->localResidual().spatialWeight()) < 1e-6) |
175 | ✗ | curElemVolVars.bindElement(element, fvGeometry, curSol); | |
176 | else | ||
177 | { | ||
178 | 10079104 | curElemVolVars.bind(element, fvGeometry, curSol); | |
179 | 5066048 | elemFluxVarsCache.bind(element, fvGeometry, curElemVolVars); | |
180 | } | ||
181 | 10079104 | } | |
182 | |||
183 | /*! | ||
184 | * \brief Enforces Dirichlet constraints if enabled in the problem | ||
185 | */ | ||
186 | template<typename ApplyFunction, class P = Problem, typename std::enable_if_t<P::enableInternalDirichletConstraints(), int> = 0> | ||
187 | void enforceInternalDirichletConstraints(const ApplyFunction& applyDirichlet) | ||
188 | { | ||
189 | // enforce Dirichlet constraints strongly by overwriting partial derivatives with 1 or 0 | ||
190 | // and set the residual to (privar - dirichletvalue) | ||
191 | for (const auto& scvI : scvs(this->fvGeometry())) | ||
192 | { | ||
193 | const auto internalDirichletConstraints = asImp_().problem().hasInternalDirichletConstraint(this->element(), scvI); | ||
194 | if (internalDirichletConstraints.any()) | ||
195 | { | ||
196 | const auto dirichletValues = asImp_().problem().internalDirichlet(this->element(), scvI); | ||
197 | // set the Dirichlet conditions in residual and jacobian | ||
198 | for (int eqIdx = 0; eqIdx < internalDirichletConstraints.size(); ++eqIdx) | ||
199 | if (internalDirichletConstraints[eqIdx]) | ||
200 | applyDirichlet(scvI, dirichletValues, eqIdx, eqIdx); | ||
201 | } | ||
202 | } | ||
203 | } | ||
204 | |||
205 | template<typename ApplyFunction, class P = Problem, typename std::enable_if_t<!P::enableInternalDirichletConstraints(), int> = 0> | ||
206 | ✗ | void enforceInternalDirichletConstraints(const ApplyFunction& applyDirichlet) | |
207 | ✗ | {} | |
208 | |||
209 | //! The problem | ||
210 | ✗ | const Problem& problem() const | |
211 |
14/33✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 4 times.
✗ Branch 18 not taken.
✓ Branch 20 taken 4 times.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✓ Branch 23 taken 4 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✓ Branch 26 taken 4 times.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✓ Branch 33 taken 2 times.
✗ Branch 34 not taken.
✓ Branch 36 taken 2 times.
✗ Branch 37 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✗ Branch 46 not taken.
✗ Branch 47 not taken.
|
23688164 | { return assembler_.problem(); } |
212 | |||
213 | //! The assembler | ||
214 | ✗ | const Assembler& assembler() const | |
215 | ✗ | { return assembler_; } | |
216 | |||
217 | //! The current element | ||
218 | ✗ | const Element& element() const | |
219 | ✗ | { return element_; } | |
220 | |||
221 | //! Returns if element is a ghost entity | ||
222 | ✗ | bool elementIsGhost() const | |
223 | ✗ | { return elementIsGhost_; } | |
224 | |||
225 | //! The current solution | ||
226 | ✗ | const SolutionVector& curSol() const | |
227 | ✗ | { return curSol_; } | |
228 | |||
229 | //! The global finite volume geometry | ||
230 | FVElementGeometry& fvGeometry() | ||
231 |
2/2✓ Branch 2 taken 5 times.
✓ Branch 3 taken 24187 times.
|
15213680 | { return fvGeometry_; } |
232 | |||
233 | //! The current element volume variables | ||
234 | ElementVolumeVariables& curElemVolVars() | ||
235 |
2/4✓ Branch 1 taken 5 times.
✓ Branch 2 taken 24187 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
10163280 | { return curElemVolVars_; } |
236 | |||
237 | //! The element flux variables cache | ||
238 | ElementFluxVariablesCache& elemFluxVarsCache() | ||
239 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 24187 times.
|
10163072 | { return elemFluxVarsCache_; } |
240 | |||
241 | //! The local residual for the current element | ||
242 | LocalResidual& localResidual() | ||
243 |
6/12✓ Branch 0 taken 10000000 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 3456 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 3456 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 26112 times.
✓ Branch 7 taken 26112 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 31616 times.
|
10111488 | { return localOperator_; } |
244 | |||
245 | //! The element's boundary types | ||
246 | ElementBoundaryTypes& elemBcTypes() | ||
247 |
4/4✓ Branch 0 taken 3888 times.
✓ Branch 1 taken 20304 times.
✓ Branch 3 taken 145152 times.
✓ Branch 4 taken 20736 times.
|
249344 | { return elemBcTypes_; } |
248 | |||
249 | //! The finite volume geometry | ||
250 | const FVElementGeometry& fvGeometry() const | ||
251 |
2/4✓ Branch 4 taken 3456 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 3456 times.
✗ Branch 8 not taken.
|
1986144 | { return fvGeometry_; } |
252 | |||
253 | //! The current element volume variables | ||
254 | const ElementVolumeVariables& curElemVolVars() const | ||
255 | 2314848 | { return curElemVolVars_; } | |
256 | |||
257 | //! The element flux variables cache | ||
258 | const ElementFluxVariablesCache& elemFluxVarsCache() const | ||
259 | 1979232 | { return elemFluxVarsCache_; } | |
260 | |||
261 | //! The element's boundary types | ||
262 | const ElementBoundaryTypes& elemBcTypes() const | ||
263 |
1/2✓ Branch 1 taken 3456 times.
✗ Branch 2 not taken.
|
369792 | { return elemBcTypes_; } |
264 | |||
265 | //! The local residual for the current element | ||
266 | const LocalResidual& localResidual() const | ||
267 |
1/2✓ Branch 2 taken 3456 times.
✗ Branch 3 not taken.
|
1982688 | { return localOperator_; } |
268 | |||
269 | //! If the time stepping scheme is implicit | ||
270 | ✗ | bool isImplicit() const | |
271 | ✗ | { return isImplicit_; } | |
272 | |||
273 | protected: | ||
274 | Implementation &asImp_() | ||
275 | { return *static_cast<Implementation*>(this); } | ||
276 | |||
277 | const Implementation &asImp_() const | ||
278 | { return *static_cast<const Implementation*>(this); } | ||
279 | |||
280 | template<class T = TypeTag, typename std::enable_if_t<!GetPropType<T, Properties::GridVariables>::GridVolumeVariables::cachingEnabled, int> = 0> | ||
281 | ✗ | VolumeVariables& getVolVarAccess(GridVolumeVariables& gridVolVars, ElementVolumeVariables& elemVolVars, const SubControlVolume& scv) | |
282 | 825216 | { return elemVolVars[scv]; } | |
283 | |||
284 | template<class T = TypeTag, typename std::enable_if_t<GetPropType<T, Properties::GridVariables>::GridVolumeVariables::cachingEnabled, int> = 0> | ||
285 | ✗ | VolumeVariables& getVolVarAccess(GridVolumeVariables& gridVolVars, ElementVolumeVariables& elemVolVars, const SubControlVolume& scv) | |
286 | 10454656 | { return gridVolVars.volVars(scv); } | |
287 | |||
288 | private: | ||
289 | |||
290 | const Assembler& assembler_; //!< access pointer to assembler instance | ||
291 | const Element& element_; //!< the element whose residual is assembled | ||
292 | const SolutionVector& curSol_; //!< the current solution | ||
293 | |||
294 | FVElementGeometry fvGeometry_; | ||
295 | ElementVolumeVariables curElemVolVars_; | ||
296 | ElementFluxVariablesCache elemFluxVarsCache_; | ||
297 | ElementBoundaryTypes elemBcTypes_; | ||
298 | |||
299 | LocalResidual localOperator_; //!< the local residual evaluating the equations per element | ||
300 | bool elementIsGhost_; //!< whether the element's partitionType is ghost | ||
301 | bool isImplicit_; //!< whether the time stepping scheme is implicit | ||
302 | }; | ||
303 | |||
304 | } // end namespace Dumux | ||
305 | |||
306 | #endif | ||
307 |