version 3.8
staggeredlocalresidual.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-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_STAGGERED_LOCAL_RESIDUAL_HH
14#define DUMUX_STAGGERED_LOCAL_RESIDUAL_HH
15
19
20namespace Dumux {
21
27template<class TypeTag>
29{
31
35 using Element = typename GridView::template Codim<0>::Entity;
37 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
38 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
40 using FVElementGeometry = typename GridGeometry::LocalView;
41 using SubControlVolume = typename GridGeometry::SubControlVolume;
42 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
43 using FaceSubControlVolume = typename GridGeometry::Traits::FaceSubControlVolume;
44 using Extrusion = Extrusion_t<GridGeometry>;
45 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
46
49 using ElementFaceVariables = typename GetPropType<TypeTag, Properties::GridFaceVariables>::LocalView;
50
52
53public:
57
60 const TimeLoop* timeLoop = nullptr)
61 : problem_(problem)
62 , timeLoop_(timeLoop)
63 {}
64
67 const FVElementGeometry& fvGeometry,
68 const ElementVolumeVariables& elemVolVars,
69 const ElementFaceVariables& elemFaceVars,
70 const ElementBoundaryTypes& bcTypes,
71 const ElementFluxVariablesCache& elemFluxVarsCache) const
72 {
73 CellCenterResidualValue residual(0.0);
74
75 // evaluate the source term
76 for (auto&& scv : scvs(fvGeometry))
77 asImp().evalSourceForCellCenter(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, scv);
78
79 // evaluate the flux term
80 for (auto&& scvf : scvfs(fvGeometry))
81 asImp().evalFluxForCellCenter(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, bcTypes, elemFluxVarsCache, scvf);
82
83 return residual;
84 }
85
88 const Problem& problem,
89 const Element& element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const ElementFaceVariables& elemFaceVars,
93 const ElementBoundaryTypes& elemBcTypes,
94 const ElementFluxVariablesCache& elemFluxVarsCache,
95 const SubControlVolumeFace& scvf) const
96 {
97 if (!scvf.boundary())
98 residual += asImp_().computeFluxForCellCenter(problem, element, fvGeometry, elemVolVars, elemFaceVars, scvf, elemFluxVarsCache);
99 else
100 residual += asImp_().computeBoundaryFluxForCellCenter(problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, elemBcTypes, elemFluxVarsCache);
101 }
102
105 const Problem& problem,
106 const Element& element,
107 const FVElementGeometry& fvGeometry,
108 const ElementVolumeVariables& curElemVolVars,
109 const ElementFaceVariables& curElemFaceVars,
110 const SubControlVolume& scv) const
111 {
112 const auto curExtrusionFactor = curElemVolVars[scv].extrusionFactor();
113
114 // subtract the source term from the local rate
115 auto source = asImp_().computeSourceForCellCenter(problem, element, fvGeometry, curElemVolVars, curElemFaceVars, scv);
116 source *= Extrusion::volume(fvGeometry, scv)*curExtrusionFactor;
117 residual -= source;
118 }
119
122 const FVElementGeometry& fvGeometry,
123 const ElementVolumeVariables& prevElemVolVars,
124 const ElementVolumeVariables& curElemVolVars) const
125 {
126 assert(timeLoop_ && "no time loop set for storage term evaluation");
127 CellCenterResidualValue storage(0.0);
128
129 for (auto&& scv : scvs(fvGeometry))
130 asImp().evalStorageForCellCenter(storage, problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
131
132 return storage;
133 }
134
137 const Problem& problem,
138 const Element &element,
139 const FVElementGeometry& fvGeometry,
140 const ElementVolumeVariables& prevElemVolVars,
141 const ElementVolumeVariables& curElemVolVars,
142 const SubControlVolume& scv) const
143 {
144 CellCenterResidualValue storage(0.0);
145 const auto& curVolVars = curElemVolVars[scv];
146 const auto& prevVolVars = prevElemVolVars[scv];
147
148 // mass balance within the element. this is the
149 // \f$\frac{m}{\partial t}\f$ term if using implicit
150 // euler as time discretization.
151
152 // We might need a more explicit way for
153 // doing the time discretization...
154 auto prevCCStorage = asImp_().computeStorageForCellCenter(problem, scv, prevVolVars);
155 auto curCCStorage = asImp_().computeStorageForCellCenter(problem, scv, curVolVars);
156
157 prevCCStorage *= prevVolVars.extrusionFactor();
158 curCCStorage *= curVolVars.extrusionFactor();
159
160 storage = std::move(curCCStorage);
161 storage -= std::move(prevCCStorage);
162 storage *= Extrusion::volume(fvGeometry, scv);
163 storage /= timeLoop_->timeStepSize();
164
165 residual += storage;
166 }
167
169 template<class... Args>
171 {
172 return CellCenterResidualValue(0.0);
173 }
174
176 template<class... Args>
178 {
179 return CellCenterResidualValue(0.0);
180 }
181
187 // \{
188
191 const FVElementGeometry& fvGeometry,
192 const ElementVolumeVariables& elemVolVars,
193 const ElementFaceVariables& elemFaceVars,
194 const ElementBoundaryTypes& bcTypes,
195 const ElementFluxVariablesCache& elemFluxVarsCache,
196 const SubControlVolumeFace& scvf) const
197 {
198 FaceResidualValue residual(0.0);
199 asImp().evalSourceForFace(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, scvf);
200 asImp().evalFluxForFace(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, bcTypes, elemFluxVarsCache, scvf);
201
202 return residual;
203 }
204
207 const Problem& problem,
208 const Element& element,
209 const FVElementGeometry& fvGeometry,
210 const ElementVolumeVariables& elemVolVars,
211 const ElementFaceVariables& elemFaceVars,
212 const ElementBoundaryTypes& elemBcTypes,
213 const ElementFluxVariablesCache& elemFluxVarsCache,
214 const SubControlVolumeFace& scvf) const
215 {
216 if (!scvf.boundary())
217 residual += asImp_().computeFluxForFace(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, elemFluxVarsCache);
218 else
219 residual += asImp_().computeBoundaryFluxForFace(problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, elemBcTypes, elemFluxVarsCache);
220 }
221
224 const Problem& problem,
225 const Element& element,
226 const FVElementGeometry& fvGeometry,
227 const ElementVolumeVariables& elemVolVars,
228 const ElementFaceVariables& elemFaceVars,
229 const SubControlVolumeFace& scvf) const
230 {
231 // the source term:
232 auto source = asImp_().computeSourceForFace(problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars);
233 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
234 const auto extrusionFactor = elemVolVars[scv].extrusionFactor();
235
236 // construct staggered scv (half of the element)
237 auto faceScvCenter = scvf.center() + scv.center();
238 faceScvCenter *= 0.5;
239 FaceSubControlVolume faceScv(faceScvCenter, 0.5*scv.volume());
240
241 source *= Extrusion::volume(fvGeometry, faceScv)*extrusionFactor;
242 residual -= source;
243 }
244
246 FaceResidualValue evalStorageForFace(const Element& element,
247 const FVElementGeometry& fvGeometry,
248 const ElementVolumeVariables& prevElemVolVars,
249 const ElementVolumeVariables& curElemVolVars,
250 const ElementFaceVariables& prevElemFaceVars,
251 const ElementFaceVariables& curElemFaceVars,
252 const SubControlVolumeFace& scvf) const
253 {
254 assert(timeLoop_ && "no time loop set for storage term evaluation");
255 FaceResidualValue storage(0.0);
256 asImp().evalStorageForFace(storage, problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, prevElemFaceVars, curElemFaceVars, scvf);
257 return storage;
258 }
259
262 const Problem& problem,
263 const Element& element,
264 const FVElementGeometry& fvGeometry,
265 const ElementVolumeVariables& prevElemVolVars,
266 const ElementVolumeVariables& curElemVolVars,
267 const ElementFaceVariables& prevElemFaceVars,
268 const ElementFaceVariables& curElemFaceVars,
269 const SubControlVolumeFace& scvf) const
270 {
271 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
272
273 auto storage = asImp_().computeStorageForFace(problem, scvf, curElemVolVars[scv], curElemFaceVars);
274 storage -= asImp_().computeStorageForFace(problem, scvf, prevElemVolVars[scv], prevElemFaceVars);
275
276 const auto extrusionFactor = curElemVolVars[scv].extrusionFactor();
277
278 // construct staggered scv (half of the element)
279 auto faceScvCenter = scvf.center() + scv.center();
280 faceScvCenter *= 0.5;
281 FaceSubControlVolume faceScv(faceScvCenter, 0.5*scv.volume());
282
283 storage *= Extrusion::volume(fvGeometry, faceScv)*extrusionFactor;
284 storage /= timeLoop_->timeStepSize();
285
286 residual += storage;
287 }
288
290 bool isStationary() const
291 { return !timeLoop_; }
292
294 const Problem& problem() const
295 { return *problem_; }
296
297protected:
298
299 Implementation &asImp_()
300 { return *static_cast<Implementation*>(this); }
301
302 const Implementation &asImp_() const
303 { return *static_cast<const Implementation*>(this); }
304
305
307 { return *timeLoop_; }
308
309 const TimeLoop& timeLoop() const
310 { return *timeLoop_; }
311
312 Implementation &asImp()
313 { return *static_cast<Implementation*>(this); }
314
315 const Implementation &asImp() const
316 { return *static_cast<const Implementation*>(this); }
317
318private:
319 const Problem* problem_;
320 const TimeLoop* timeLoop_;
321
322};
323
324} // end namespace Dumux
325
326#endif
Calculates the element-wise residual for the staggered FV scheme.
Definition: staggeredlocalresidual.hh:29
CellCenterResidualValue evalStorageForCellCenter(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars) const
Evaluate the storage terms for a cell center residual.
Definition: staggeredlocalresidual.hh:121
GetPropType< TypeTag, Properties::CellCenterPrimaryVariables > CellCenterResidualValue
Definition: staggeredlocalresidual.hh:54
Implementation & asImp()
Definition: staggeredlocalresidual.hh:312
void evalStorageForCellCenter(CellCenterResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars, const SubControlVolume &scv) const
Evaluate the storage terms for a cell center residual.
Definition: staggeredlocalresidual.hh:136
Implementation & asImp_()
Definition: staggeredlocalresidual.hh:299
const TimeLoop & timeLoop() const
Definition: staggeredlocalresidual.hh:309
FaceResidualValue evalFluxAndSourceForFace(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &bcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Convenience function to evaluate the flux and source terms for the face residual.
Definition: staggeredlocalresidual.hh:190
GetPropType< TypeTag, Properties::FacePrimaryVariables > FaceResidualValue
Definition: staggeredlocalresidual.hh:55
const Problem & problem() const
the problem
Definition: staggeredlocalresidual.hh:294
const Implementation & asImp_() const
Definition: staggeredlocalresidual.hh:302
CellCenterResidualValue evalStorage(Args &&... args) const
for compatibility with FVLocalAssemblerBase
Definition: staggeredlocalresidual.hh:177
void evalSourceForFace(FaceResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const
Evaluate the source terms for a face residual.
Definition: staggeredlocalresidual.hh:223
bool isStationary() const
If no solution has been set, we treat the problem as stationary.
Definition: staggeredlocalresidual.hh:290
TimeLoop & timeLoop()
Definition: staggeredlocalresidual.hh:306
StaggeredLocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition: staggeredlocalresidual.hh:59
CellCenterResidualValue evalFluxAndSource(Args &&... args) const
for compatibility with FVLocalAssemblerBase
Definition: staggeredlocalresidual.hh:170
void evalStorageForFace(FaceResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars, const ElementFaceVariables &prevElemFaceVars, const ElementFaceVariables &curElemFaceVars, const SubControlVolumeFace &scvf) const
Evaluate the storage terms for a face residual.
Definition: staggeredlocalresidual.hh:261
CellCenterResidualValue evalFluxAndSourceForCellCenter(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &bcTypes, const ElementFluxVariablesCache &elemFluxVarsCache) const
Convenience function to evaluate the flux and source terms for the cell center residual.
Definition: staggeredlocalresidual.hh:66
void evalFluxForFace(FaceResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Evaluate the flux terms for a face residual.
Definition: staggeredlocalresidual.hh:206
void evalSourceForCellCenter(CellCenterResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFaceVariables &curElemFaceVars, const SubControlVolume &scv) const
Evaluate the source terms for a cell center residual.
Definition: staggeredlocalresidual.hh:104
CellCenterResidualValue ElementResidualVector
Definition: staggeredlocalresidual.hh:56
const Implementation & asImp() const
Definition: staggeredlocalresidual.hh:315
void evalFluxForCellCenter(CellCenterResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Evaluate the flux terms for a cell center residual.
Definition: staggeredlocalresidual.hh:87
FaceResidualValue evalStorageForFace(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars, const ElementFaceVariables &prevElemFaceVars, const ElementFaceVariables &curElemFaceVars, const SubControlVolumeFace &scvf) const
Evaluate the storage terms for a face residual.
Definition: staggeredlocalresidual.hh:246
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
The default time loop for instationary simulations.
Definition: common/timeloop.hh:139
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
Helper classes to compute the integration elements.
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166