3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_STAGGERED_LOCAL_RESIDUAL_HH
26#define DUMUX_STAGGERED_LOCAL_RESIDUAL_HH
27
31
32namespace Dumux {
33
39template<class TypeTag>
41{
43
47 using Element = typename GridView::template Codim<0>::Entity;
49 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
50 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
52 using FVElementGeometry = typename GridGeometry::LocalView;
53 using SubControlVolume = typename GridGeometry::SubControlVolume;
54 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
55 using FaceSubControlVolume = typename GridGeometry::Traits::FaceSubControlVolume;
56 using Extrusion = Extrusion_t<GridGeometry>;
57 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
58
61 using ElementFaceVariables = typename GetPropType<TypeTag, Properties::GridFaceVariables>::LocalView;
62
64
65public:
69
72 const TimeLoop* timeLoop = nullptr)
73 : problem_(problem)
74 , timeLoop_(timeLoop)
75 {}
76
79 const FVElementGeometry& fvGeometry,
80 const ElementVolumeVariables& elemVolVars,
81 const ElementFaceVariables& elemFaceVars,
82 const ElementBoundaryTypes& bcTypes,
83 const ElementFluxVariablesCache& elemFluxVarsCache) const
84 {
85 CellCenterResidualValue residual(0.0);
86
87 // evaluate the source term
88 for (auto&& scv : scvs(fvGeometry))
89 asImp().evalSourceForCellCenter(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, scv);
90
91 // evaluate the flux term
92 for (auto&& scvf : scvfs(fvGeometry))
93 asImp().evalFluxForCellCenter(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, bcTypes, elemFluxVarsCache, scvf);
94
95 return residual;
96 }
97
100 const Problem& problem,
101 const Element& element,
102 const FVElementGeometry& fvGeometry,
103 const ElementVolumeVariables& elemVolVars,
104 const ElementFaceVariables& elemFaceVars,
105 const ElementBoundaryTypes& elemBcTypes,
106 const ElementFluxVariablesCache& elemFluxVarsCache,
107 const SubControlVolumeFace& scvf) const
108 {
109 if (!scvf.boundary())
110 residual += asImp_().computeFluxForCellCenter(problem, element, fvGeometry, elemVolVars, elemFaceVars, scvf, elemFluxVarsCache);
111 else
112 residual += asImp_().computeBoundaryFluxForCellCenter(problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, elemBcTypes, elemFluxVarsCache);
113 }
114
117 const Problem& problem,
118 const Element& element,
119 const FVElementGeometry& fvGeometry,
120 const ElementVolumeVariables& curElemVolVars,
121 const ElementFaceVariables& curElemFaceVars,
122 const SubControlVolume& scv) const
123 {
124 const auto curExtrusionFactor = curElemVolVars[scv].extrusionFactor();
125
126 // subtract the source term from the local rate
127 auto source = asImp_().computeSourceForCellCenter(problem, element, fvGeometry, curElemVolVars, curElemFaceVars, scv);
128 source *= Extrusion::volume(scv)*curExtrusionFactor;
129 residual -= source;
130 }
131
134 const FVElementGeometry& fvGeometry,
135 const ElementVolumeVariables& prevElemVolVars,
136 const ElementVolumeVariables& curElemVolVars) const
137 {
138 assert(timeLoop_ && "no time loop set for storage term evaluation");
139 CellCenterResidualValue storage(0.0);
140
141 for (auto&& scv : scvs(fvGeometry))
142 asImp().evalStorageForCellCenter(storage, problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
143
144 return storage;
145 }
146
149 const Problem& problem,
150 const Element &element,
151 const FVElementGeometry& fvGeometry,
152 const ElementVolumeVariables& prevElemVolVars,
153 const ElementVolumeVariables& curElemVolVars,
154 const SubControlVolume& scv) const
155 {
156 CellCenterResidualValue storage(0.0);
157 const auto& curVolVars = curElemVolVars[scv];
158 const auto& prevVolVars = prevElemVolVars[scv];
159
160 // mass balance within the element. this is the
161 // \f$\frac{m}{\partial t}\f$ term if using implicit
162 // euler as time discretization.
163
164 // We might need a more explicit way for
165 // doing the time discretization...
166 auto prevCCStorage = asImp_().computeStorageForCellCenter(problem, scv, prevVolVars);
167 auto curCCStorage = asImp_().computeStorageForCellCenter(problem, scv, curVolVars);
168
169 prevCCStorage *= prevVolVars.extrusionFactor();
170 curCCStorage *= curVolVars.extrusionFactor();
171
172 storage = std::move(curCCStorage);
173 storage -= std::move(prevCCStorage);
174 storage *= Extrusion::volume(scv);
175 storage /= timeLoop_->timeStepSize();
176
177 residual += storage;
178 }
179
181 template<class... Args>
183 {
184 return CellCenterResidualValue(0.0);
185 }
186
188 template<class... Args>
190 {
191 return CellCenterResidualValue(0.0);
192 }
193
199 // \{
200
203 const FVElementGeometry& fvGeometry,
204 const ElementVolumeVariables& elemVolVars,
205 const ElementFaceVariables& elemFaceVars,
206 const ElementBoundaryTypes& bcTypes,
207 const ElementFluxVariablesCache& elemFluxVarsCache,
208 const SubControlVolumeFace& scvf) const
209 {
210 FaceResidualValue residual(0.0);
211 asImp().evalSourceForFace(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, scvf);
212 asImp().evalFluxForFace(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, bcTypes, elemFluxVarsCache, scvf);
213
214 return residual;
215 }
216
219 const Problem& problem,
220 const Element& element,
221 const FVElementGeometry& fvGeometry,
222 const ElementVolumeVariables& elemVolVars,
223 const ElementFaceVariables& elemFaceVars,
224 const ElementBoundaryTypes& elemBcTypes,
225 const ElementFluxVariablesCache& elemFluxVarsCache,
226 const SubControlVolumeFace& scvf) const
227 {
228 if (!scvf.boundary())
229 residual += asImp_().computeFluxForFace(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, elemFluxVarsCache);
230 else
231 residual += asImp_().computeBoundaryFluxForFace(problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, elemBcTypes, elemFluxVarsCache);
232 }
233
236 const Problem& problem,
237 const Element& element,
238 const FVElementGeometry& fvGeometry,
239 const ElementVolumeVariables& elemVolVars,
240 const ElementFaceVariables& elemFaceVars,
241 const SubControlVolumeFace& scvf) const
242 {
243 // the source term:
244 auto source = asImp_().computeSourceForFace(problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars);
245 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
246 const auto extrusionFactor = elemVolVars[scv].extrusionFactor();
247
248 // construct staggered scv (half of the element)
249 auto faceScvCenter = scvf.center() + scv.center();
250 faceScvCenter *= 0.5;
251 FaceSubControlVolume faceScv(faceScvCenter, 0.5*scv.volume());
252
253 source *= Extrusion::volume(faceScv)*extrusionFactor;
254 residual -= source;
255 }
256
258 FaceResidualValue evalStorageForFace(const Element& element,
259 const FVElementGeometry& fvGeometry,
260 const ElementVolumeVariables& prevElemVolVars,
261 const ElementVolumeVariables& curElemVolVars,
262 const ElementFaceVariables& prevElemFaceVars,
263 const ElementFaceVariables& curElemFaceVars,
264 const SubControlVolumeFace& scvf) const
265 {
266 assert(timeLoop_ && "no time loop set for storage term evaluation");
267 FaceResidualValue storage(0.0);
268 asImp().evalStorageForFace(storage, problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, prevElemFaceVars, curElemFaceVars, scvf);
269 return storage;
270 }
271
274 const Problem& problem,
275 const Element& element,
276 const FVElementGeometry& fvGeometry,
277 const ElementVolumeVariables& prevElemVolVars,
278 const ElementVolumeVariables& curElemVolVars,
279 const ElementFaceVariables& prevElemFaceVars,
280 const ElementFaceVariables& curElemFaceVars,
281 const SubControlVolumeFace& scvf) const
282 {
283 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
284
285 auto storage = asImp_().computeStorageForFace(problem, scvf, curElemVolVars[scv], curElemFaceVars);
286 storage -= asImp_().computeStorageForFace(problem, scvf, prevElemVolVars[scv], prevElemFaceVars);
287
288 const auto extrusionFactor = curElemVolVars[scv].extrusionFactor();
289
290 // construct staggered scv (half of the element)
291 auto faceScvCenter = scvf.center() + scv.center();
292 faceScvCenter *= 0.5;
293 FaceSubControlVolume faceScv(faceScvCenter, 0.5*scv.volume());
294
295 storage *= Extrusion::volume(faceScv)*extrusionFactor;
296 storage /= timeLoop_->timeStepSize();
297
298 residual += storage;
299 }
300
302 bool isStationary() const
303 { return !timeLoop_; }
304
306 const Problem& problem() const
307 { return *problem_; }
308
309protected:
310
311 Implementation &asImp_()
312 { return *static_cast<Implementation*>(this); }
313
314 const Implementation &asImp_() const
315 { return *static_cast<const Implementation*>(this); }
316
317
319 { return *timeLoop_; }
320
321 const TimeLoop& timeLoop() const
322 { return *timeLoop_; }
323
324 Implementation &asImp()
325 { return *static_cast<Implementation*>(this); }
326
327 const Implementation &asImp() const
328 { return *static_cast<const Implementation*>(this); }
329
330private:
331 const Problem* problem_;
332 const TimeLoop* timeLoop_;
333
334};
335
336} // end namespace Dumux
337
338#endif
Helper classes to compute the integration elements.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Calculates the element-wise residual for the staggered FV scheme.
Definition: staggeredlocalresidual.hh:41
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:133
GetPropType< TypeTag, Properties::CellCenterPrimaryVariables > CellCenterResidualValue
Definition: staggeredlocalresidual.hh:66
Implementation & asImp()
Definition: staggeredlocalresidual.hh:324
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:148
Implementation & asImp_()
Definition: staggeredlocalresidual.hh:311
const TimeLoop & timeLoop() const
Definition: staggeredlocalresidual.hh:321
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:202
GetPropType< TypeTag, Properties::FacePrimaryVariables > FaceResidualValue
Definition: staggeredlocalresidual.hh:67
const Problem & problem() const
the problem
Definition: staggeredlocalresidual.hh:306
const Implementation & asImp_() const
Definition: staggeredlocalresidual.hh:314
CellCenterResidualValue evalStorage(Args &&... args) const
for compatibility with FVLocalAssemblerBase
Definition: staggeredlocalresidual.hh:189
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:235
bool isStationary() const
If no solution has been set, we treat the problem as stationary.
Definition: staggeredlocalresidual.hh:302
TimeLoop & timeLoop()
Definition: staggeredlocalresidual.hh:318
StaggeredLocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition: staggeredlocalresidual.hh:71
CellCenterResidualValue evalFluxAndSource(Args &&... args) const
for compatibility with FVLocalAssemblerBase
Definition: staggeredlocalresidual.hh:182
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:273
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:78
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:218
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:116
CellCenterResidualValue ElementResidualVector
Definition: staggeredlocalresidual.hh:68
const Implementation & asImp() const
Definition: staggeredlocalresidual.hh:327
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:99
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:258
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
The default time loop for instationary simulations.
Definition: common/timeloop.hh:114
Declares all properties used in Dumux.
Manages the handling of time dependent problems.