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-FileCopyrightText: 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 Assembly | ||
10 | * \ingroup CVFEDiscretization | ||
11 | * \brief Calculates the element-wise residual for control-volume finite element schemes | ||
12 | */ | ||
13 | #ifndef DUMUX_CVFE_LOCAL_RESIDUAL_HH | ||
14 | #define DUMUX_CVFE_LOCAL_RESIDUAL_HH | ||
15 | |||
16 | #include <dune/common/std/type_traits.hh> | ||
17 | #include <dune/geometry/type.hh> | ||
18 | #include <dune/istl/matrix.hh> | ||
19 | |||
20 | #include <dumux/common/typetraits/localdofs_.hh> | ||
21 | #include <dumux/common/typetraits/boundary_.hh> | ||
22 | #include <dumux/common/properties.hh> | ||
23 | #include <dumux/common/numeqvector.hh> | ||
24 | #include <dumux/assembly/fvlocalresidual.hh> | ||
25 | #include <dumux/discretization/extrusion.hh> | ||
26 | #include <dumux/discretization/cvfe/integrationpointdata.hh> | ||
27 | |||
28 | namespace Dumux::Detail { | ||
29 | |||
30 | template<class Imp, class P, class G, class S, class V> | ||
31 | using TimeInfoInterfaceCVFEDetector = decltype( | ||
32 | std::declval<Imp>().computeStorage( | ||
33 | std::declval<P>(), std::declval<G>(), std::declval<S>(), std::declval<V>(), true | ||
34 | ) | ||
35 | ); | ||
36 | |||
37 | template<class Imp, class P, class G, class S, class V> | ||
38 | constexpr inline bool hasTimeInfoInterfaceCVFE() | ||
39 | { return Dune::Std::is_detected<TimeInfoInterfaceCVFEDetector, Imp, P, G, S, V>::value; } | ||
40 | |||
41 | template<class Imp> | ||
42 | using SCVFIsOverlappingDetector = decltype( | ||
43 | std::declval<Imp>().isOverlapping() | ||
44 | ); | ||
45 | |||
46 | template<class Imp> | ||
47 | constexpr inline bool hasScvfIsOverlapping() | ||
48 | { return Dune::Std::is_detected<SCVFIsOverlappingDetector, Imp>::value; } | ||
49 | |||
50 | } // end namespace Dumux::Detail | ||
51 | |||
52 | |||
53 | namespace Dumux { | ||
54 | |||
55 | /*! | ||
56 | * \ingroup Assembly | ||
57 | * \ingroup CVFEDiscretization | ||
58 | * \brief The element-wise residual for control-volume finite element schemes | ||
59 | * \tparam TypeTag the TypeTag | ||
60 | */ | ||
61 | template<class TypeTag> | ||
62 | class CVFELocalResidual : public FVLocalResidual<TypeTag> | ||
63 | { | ||
64 | using ParentType = FVLocalResidual<TypeTag>; | ||
65 | using Implementation = GetPropType<TypeTag, Properties::LocalResidual>; | ||
66 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
67 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
68 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
69 | using GridView = typename GridGeometry::GridView; | ||
70 | using Element = typename GridView::template Codim<0>::Entity; | ||
71 | using ElementBoundaryTypes = GetPropType<TypeTag, Properties::ElementBoundaryTypes>; | ||
72 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
73 | using GridVolumeVariables = GetPropType<TypeTag, Properties::GridVolumeVariables>; | ||
74 | using ElementVolumeVariables = typename GridVolumeVariables::LocalView; | ||
75 | using VolumeVariables = typename GridVolumeVariables::VolumeVariables; | ||
76 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
77 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
78 | using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView; | ||
79 | using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>; | ||
80 | using NumEqVector = Dumux::NumEqVector<PrimaryVariables>; | ||
81 | using Extrusion = Extrusion_t<GridGeometry>; | ||
82 | |||
83 | using FaceIpData = Dumux::CVFE::FaceIntegrationPointData<typename Element::Geometry::GlobalCoordinate, | ||
84 | typename SubControlVolumeFace::Traits::LocalIndexType>; | ||
85 | |||
86 | public: | ||
87 | using ElementResidualVector = typename ParentType::ElementResidualVector; | ||
88 |
2/4✓ Branch 1 taken 8431842 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 50002 times.
✗ Branch 5 not taken.
|
14793341 | using ParentType::ParentType; |
89 | |||
90 | using ParentType::evalStorage; | ||
91 | /*! | ||
92 | * \brief Compute the storage local residual, i.e. the deviation of the | ||
93 | * storage term from zero for instationary problems. | ||
94 | * | ||
95 | * \param element The DUNE Codim<0> entity for which the residual | ||
96 | * ought to be calculated | ||
97 | * \param fvGeometry The finite-volume geometry of the element | ||
98 | * \param prevElemVolVars The volume averaged variables for all | ||
99 | * sub-control volumes of the element at the previous time level | ||
100 | * \param curElemVolVars The volume averaged variables for all | ||
101 | * sub-control volumes of the element at the current time level | ||
102 | */ | ||
103 | 86487556 | ElementResidualVector evalStorage(const Element& element, | |
104 | const FVElementGeometry& fvGeometry, | ||
105 | const ElementVolumeVariables& prevElemVolVars, | ||
106 | const ElementVolumeVariables& curElemVolVars) const | ||
107 | { | ||
108 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 86287556 times.
|
86487556 | assert(!this->isStationary() && "no time loop set for storage term evaluation"); |
109 | |||
110 | // initialize the residual vector for all scvs in this element | ||
111 |
0/2✗ Branch 0 not taken.
✗ Branch 1 not taken.
|
93274225 | ElementResidualVector residual(Detail::LocalDofs::numLocalDofs(fvGeometry)); |
112 | |||
113 | // evaluate the volume terms (storage + source terms) | ||
114 | // forward to the local residual specialized for the discretization methods | ||
115 |
2/2✓ Branch 0 taken 313829364 times.
✓ Branch 1 taken 86287556 times.
|
400716920 | for (const auto& scv : scvs(fvGeometry)) |
116 | 314229364 | this->asImp().evalStorage(residual, this->problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv); | |
117 | |||
118 | // allow for additional contributions (e.g. hybrid CVFE schemes) | ||
119 | 86487556 | this->asImp().addToElementStorageResidual(residual, this->problem(), element, fvGeometry, prevElemVolVars, curElemVolVars); | |
120 | |||
121 | 86487556 | return residual; | |
122 | } | ||
123 | |||
124 | /*! | ||
125 | * \brief Compute the flux and source | ||
126 | * | ||
127 | * \param element The DUNE Codim<0> entity for which the residual | ||
128 | * ought to be calculated | ||
129 | * \param fvGeometry The finite-volume geometry of the element | ||
130 | * \param elemVolVars The volume averaged variables for all | ||
131 | * sub-control volumes of the element at the current time level | ||
132 | * \param elemFluxVarsCache The element flux variables cache | ||
133 | * \param bcTypes The element boundary types | ||
134 | */ | ||
135 | 118618617 | ElementResidualVector evalFluxAndSource(const Element& element, | |
136 | const FVElementGeometry& fvGeometry, | ||
137 | const ElementVolumeVariables& elemVolVars, | ||
138 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
139 | const ElementBoundaryTypes &bcTypes) const | ||
140 | { | ||
141 | // initialize the residual vector for all scvs in this element | ||
142 | 142672307 | ElementResidualVector residual(Detail::LocalDofs::numLocalDofs(fvGeometry)); | |
143 | |||
144 | // evaluate the volume terms (storage + source terms) | ||
145 | // forward to the local residual specialized for the discretization methods | ||
146 |
2/2✓ Branch 0 taken 384507162 times.
✓ Branch 1 taken 104894521 times.
|
549305705 | for (const auto& scv : scvs(fvGeometry)) |
147 | 430687088 | this->asImp().evalSource(residual, this->problem(), element, fvGeometry, elemVolVars, scv); | |
148 | |||
149 | // forward to the local residual specialized for the discretization methods | ||
150 |
2/2✓ Branch 0 taken 442351232 times.
✓ Branch 1 taken 104894521 times.
|
616489362 | for (auto&& scvf : scvfs(fvGeometry)) |
151 | 497870745 | this->asImp().evalFlux(residual, this->problem(), element, fvGeometry, elemVolVars, bcTypes, elemFluxVarsCache, scvf); | |
152 | |||
153 | // allow for additional contributions (e.g. hybrid CVFE schemes) | ||
154 | 118618617 | this->asImp().addToElementFluxAndSourceResidual(residual, this->problem(), element, fvGeometry, elemVolVars, elemFluxVarsCache, bcTypes); | |
155 | |||
156 | 118618617 | return residual; | |
157 | } | ||
158 | |||
159 | //! add additional storage contributions (e.g. hybrid CVFE schemes) | ||
160 | void addToElementStorageResidual(ElementResidualVector& residual, | ||
161 | const Problem& problem, | ||
162 | const Element& element, | ||
163 | const FVElementGeometry& fvGeometry, | ||
164 | const ElementVolumeVariables& prevElemVolVars, | ||
165 | const ElementVolumeVariables& curElemVolVars) const | ||
166 | {} | ||
167 | |||
168 | //! add additional flux and source contributions (e.g. hybrid CVFE schemes) | ||
169 | void addToElementFluxAndSourceResidual(ElementResidualVector& residual, | ||
170 | const Problem& problem, | ||
171 | const Element& element, | ||
172 | const FVElementGeometry& fvGeometry, | ||
173 | const ElementVolumeVariables& curElemVolVars, | ||
174 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
175 | const ElementBoundaryTypes &bcTypes) const | ||
176 | {} | ||
177 | |||
178 | //! evaluate flux residuals for one sub control volume face and add to residual | ||
179 | 497870745 | void evalFlux(ElementResidualVector& residual, | |
180 | const Problem& problem, | ||
181 | const Element& element, | ||
182 | const FVElementGeometry& fvGeometry, | ||
183 | const ElementVolumeVariables& elemVolVars, | ||
184 | const ElementBoundaryTypes& elemBcTypes, | ||
185 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
186 | const SubControlVolumeFace& scvf) const | ||
187 | { | ||
188 |
2/2✓ Branch 1 taken 31107188 times.
✓ Branch 2 taken 3895064 times.
|
505221644 | const auto flux = this->asImp().evalFlux(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf); |
189 |
2/2✓ Branch 0 taken 397891916 times.
✓ Branch 1 taken 39589921 times.
|
489360957 | if (!scvf.boundary()) |
190 | { | ||
191 |
2/2✓ Branch 0 taken 7400 times.
✓ Branch 1 taken 7400 times.
|
457244808 | const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx()); |
192 | 457244808 | const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx()); | |
193 |
2/2✓ Branch 0 taken 18452304 times.
✓ Branch 1 taken 18231704 times.
|
526260024 | residual[insideScv.localDofIndex()] += flux; |
194 | |||
195 | // for control-volume finite element schemes with overlapping control volumes | ||
196 | if constexpr (Detail::hasScvfIsOverlapping<SubControlVolumeFace>()) | ||
197 | { | ||
198 |
2/2✓ Branch 0 taken 18452304 times.
✓ Branch 1 taken 18231704 times.
|
69030016 | if (!scvf.isOverlapping()) |
199 | 34735608 | residual[outsideScv.localDofIndex()] -= flux; | |
200 | } | ||
201 | else | ||
202 | 388214792 | residual[outsideScv.localDofIndex()] -= flux; | |
203 | } | ||
204 | else | ||
205 | { | ||
206 | 40625937 | const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx()); | |
207 | 40625937 | residual[insideScv.localDofIndex()] += flux; | |
208 | } | ||
209 | 497870745 | } | |
210 | |||
211 | //! evaluate flux residuals for one sub control volume face | ||
212 | 494230352 | NumEqVector evalFlux(const Problem& problem, | |
213 | const Element& element, | ||
214 | const FVElementGeometry& fvGeometry, | ||
215 | const ElementVolumeVariables& elemVolVars, | ||
216 | const ElementBoundaryTypes& elemBcTypes, | ||
217 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
218 | const SubControlVolumeFace& scvf) const | ||
219 | { | ||
220 |
2/2✓ Branch 0 taken 397891916 times.
✓ Branch 1 taken 39589921 times.
|
494230352 | NumEqVector flux(0.0); |
221 | |||
222 | // inner faces | ||
223 |
2/2✓ Branch 0 taken 397891916 times.
✓ Branch 1 taken 39589921 times.
|
489360957 | if (!scvf.boundary()) |
224 | 843272219 | flux += this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache); | |
225 | |||
226 | // boundary faces | ||
227 | else | ||
228 | { | ||
229 |
2/2✓ Branch 1 taken 8781642 times.
✓ Branch 2 taken 2819466 times.
|
40625937 | const auto& bcTypes = bcTypes_(problem, fvGeometry, scvf, elemBcTypes); |
230 | |||
231 | // Treat Neumann and Robin ("solution dependent Neumann") boundary conditions. | ||
232 | // For Dirichlet there is no addition to the residual here but they | ||
233 | // are enforced strongly by replacing the residual entry afterwards. | ||
234 |
2/2✓ Branch 0 taken 32178306 times.
✓ Branch 1 taken 7411615 times.
|
40625937 | if (bcTypes.hasNeumann()) |
235 | { | ||
236 |
2/2✓ Branch 0 taken 360657 times.
✓ Branch 1 taken 3192327 times.
|
14372826 | NumEqVector boundaryFluxes; |
237 | if constexpr (Dumux::Detail::hasProblemBoundaryFluxFunction<Problem, FVElementGeometry, ElementVolumeVariables, ElementFluxVariablesCache, FaceIpData>()) | ||
238 | boundaryFluxes = problem.boundaryFlux(fvGeometry, elemVolVars, elemFluxVarsCache, FaceIpData(scvf.ipGlobal(), scvf.unitOuterNormal(), scvf.index())); | ||
239 | else | ||
240 |
2/2✓ Branch 0 taken 362628 times.
✓ Branch 1 taken 9322047 times.
|
42666747 | boundaryFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf); |
241 | |||
242 | // multiply neumann fluxes with the area and the extrusion factor | ||
243 | 32689270 | const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); | |
244 | 32689270 | boundaryFluxes *= Extrusion::area(fvGeometry, scvf)*elemVolVars[scv].extrusionFactor(); | |
245 | |||
246 | // only add fluxes to equations for which Neumann is set | ||
247 |
2/2✓ Branch 0 taken 96135363 times.
✓ Branch 1 taken 32178306 times.
|
129371873 | for (int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx) |
248 |
3/4✓ Branch 0 taken 44048 times.
✓ Branch 1 taken 96091315 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 44048 times.
|
96682603 | if (bcTypes.isNeumann(eqIdx)) |
249 | 96638555 | flux[eqIdx] += boundaryFluxes[eqIdx]; | |
250 | } | ||
251 | } | ||
252 | |||
253 | 490496685 | return flux; | |
254 | } | ||
255 | |||
256 | /*! | ||
257 | * \brief Compute the storage local residual, i.e. the deviation of the | ||
258 | * storage term from zero for instationary problems. | ||
259 | * | ||
260 | * \param residual The residual vector to fill | ||
261 | * \param problem The problem to solve | ||
262 | * \param element The DUNE Codim<0> entity for which the residual | ||
263 | * ought to be calculated | ||
264 | * \param fvGeometry The finite-volume geometry of the element | ||
265 | * \param prevElemVolVars The volume averaged variables for all | ||
266 | * sub-control volumes of the element at the previous time level | ||
267 | * \param curElemVolVars The volume averaged variables for all | ||
268 | * sub-control volumes of the element at the current time level | ||
269 | * \param scv The sub control volume the storage term is integrated over | ||
270 | */ | ||
271 |
1/2✓ Branch 0 taken 400000 times.
✗ Branch 1 not taken.
|
302760592 | void evalStorage(ElementResidualVector& residual, |
272 | const Problem& problem, | ||
273 | const Element& element, | ||
274 | const FVElementGeometry& fvGeometry, | ||
275 | const ElementVolumeVariables& prevElemVolVars, | ||
276 | const ElementVolumeVariables& curElemVolVars, | ||
277 | const SubControlVolume& scv) const | ||
278 | { | ||
279 | 302760592 | const auto& curVolVars = curElemVolVars[scv]; | |
280 | 302760592 | const auto& prevVolVars = prevElemVolVars[scv]; | |
281 | |||
282 | // Compute storage with the model specific storage residual | ||
283 | // This addresses issues #792/#940 in ad-hoc way by additionally providing crude time level information (previous or current) | ||
284 | // to the low-level interfaces if this is supported by the LocalResidual implementation | ||
285 | 302760592 | NumEqVector prevStorage = computeStorageImpl_(problem, fvGeometry, scv, prevVolVars, /*previous time level?*/true); | |
286 | 302760592 | NumEqVector storage = computeStorageImpl_(problem, fvGeometry, scv, curVolVars, /*previous time level?*/false); | |
287 | |||
288 | 461581620 | prevStorage *= prevVolVars.extrusionFactor(); | |
289 |
2/2✓ Branch 0 taken 789356022 times.
✓ Branch 1 taken 283090896 times.
|
1387499926 | storage *= curVolVars.extrusionFactor(); |
290 | |||
291 | 302760592 | storage -= prevStorage; | |
292 | 302760592 | storage *= Extrusion::volume(fvGeometry, scv); | |
293 | 302760592 | storage /= this->timeLoop().timeStepSize(); | |
294 | |||
295 | 302760592 | residual[scv.localDofIndex()] += storage; | |
296 | 302760592 | } | |
297 | |||
298 | private: | ||
299 | 501255654 | NumEqVector computeStorageImpl_(const Problem& problem, | |
300 | const FVElementGeometry& fvGeometry, | ||
301 | const SubControlVolume& scv, | ||
302 | const VolumeVariables& volVars, | ||
303 | [[maybe_unused]] bool isPreviousTimeStep) const | ||
304 | { | ||
305 | if constexpr (Detail::hasTimeInfoInterfaceCVFE<Implementation, Problem, FVElementGeometry, SubControlVolume, VolumeVariables>()) | ||
306 | ✗ | return this->asImp().computeStorage(problem, fvGeometry, scv, volVars, isPreviousTimeStep); | |
307 | else | ||
308 |
2/4✓ Branch 0 taken 400000 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 400000 times.
✗ Branch 3 not taken.
|
451610980 | return this->asImp().computeStorage(problem, scv, volVars); |
309 | } | ||
310 | |||
311 | 39589921 | auto bcTypes_(const Problem& problem, | |
312 | const FVElementGeometry& fvGeometry, | ||
313 | const SubControlVolumeFace& scvf, | ||
314 | const ElementBoundaryTypes& elemBcTypes) const | ||
315 | { | ||
316 | // Check if problem supports the new boundaryTypes function for element intersections | ||
317 | // then we can always get bcTypes for intersections and the associated scvfs | ||
318 | if constexpr (Detail::hasProblemBoundaryTypesForIntersectionFunction<Problem, FVElementGeometry, typename GridView::Intersection>()) | ||
319 | return elemBcTypes.get(fvGeometry, scvf); | ||
320 | else | ||
321 |
3/4✓ Branch 1 taken 23282364 times.
✓ Branch 2 taken 3835101 times.
✓ Branch 4 taken 4000 times.
✗ Branch 5 not taken.
|
39589921 | return elemBcTypes.get(fvGeometry, fvGeometry.scv(scvf.insideScvIdx())); |
322 | } | ||
323 | |||
324 | }; | ||
325 | |||
326 | } // end namespace Dumux | ||
327 | |||
328 | #endif | ||
329 |