3.2-git
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
30
31namespace Dumux {
32
38template<class TypeTag>
40{
42
46 using Element = typename GridView::template Codim<0>::Entity;
48 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
49 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
50 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
51 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
52 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
53 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
54
57 using ElementFaceVariables = typename GetPropType<TypeTag, Properties::GridFaceVariables>::LocalView;
58
60
61public:
65
68 const TimeLoop* timeLoop = nullptr)
69 : problem_(problem)
70 , timeLoop_(timeLoop)
71 {}
72
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const ElementFaceVariables& elemFaceVars,
78 const ElementBoundaryTypes& bcTypes,
79 const ElementFluxVariablesCache& elemFluxVarsCache) const
80 {
81 CellCenterResidualValue residual(0.0);
82
83 // evaluate the source term
84 for (auto&& scv : scvs(fvGeometry))
85 asImp().evalSourceForCellCenter(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, scv);
86
87 // evaluate the flux term
88 for (auto&& scvf : scvfs(fvGeometry))
89 asImp().evalFluxForCellCenter(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, bcTypes, elemFluxVarsCache, scvf);
90
91 return residual;
92 }
93
96 const Problem& problem,
97 const Element& element,
98 const FVElementGeometry& fvGeometry,
99 const ElementVolumeVariables& elemVolVars,
100 const ElementFaceVariables& elemFaceVars,
101 const ElementBoundaryTypes& elemBcTypes,
102 const ElementFluxVariablesCache& elemFluxVarsCache,
103 const SubControlVolumeFace& scvf) const
104 {
105 if (!scvf.boundary())
106 residual += asImp_().computeFluxForCellCenter(problem, element, fvGeometry, elemVolVars, elemFaceVars, scvf, elemFluxVarsCache);
107 else
108 residual += asImp_().computeBoundaryFluxForCellCenter(problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, elemBcTypes, elemFluxVarsCache);
109 }
110
113 const Problem& problem,
114 const Element& element,
115 const FVElementGeometry& fvGeometry,
116 const ElementVolumeVariables& curElemVolVars,
117 const ElementFaceVariables& curElemFaceVars,
118 const SubControlVolume& scv) const
119 {
120 const auto curExtrusionFactor = curElemVolVars[scv].extrusionFactor();
121
122 // subtract the source term from the local rate
123 auto source = asImp_().computeSourceForCellCenter(problem, element, fvGeometry, curElemVolVars, curElemFaceVars, scv);
124 source *= scv.volume()*curExtrusionFactor;
125 residual -= source;
126 }
127
130 const FVElementGeometry& fvGeometry,
131 const ElementVolumeVariables& prevElemVolVars,
132 const ElementVolumeVariables& curElemVolVars) const
133 {
134 assert(timeLoop_ && "no time loop set for storage term evaluation");
135 CellCenterResidualValue storage(0.0);
136
137 for (auto&& scv : scvs(fvGeometry))
138 asImp().evalStorageForCellCenter(storage, problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
139
140 return storage;
141 }
142
145 const Problem& problem,
146 const Element &element,
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& prevElemVolVars,
149 const ElementVolumeVariables& curElemVolVars,
150 const SubControlVolume& scv) const
151 {
152 CellCenterResidualValue storage(0.0);
153 const auto& curVolVars = curElemVolVars[scv];
154 const auto& prevVolVars = prevElemVolVars[scv];
155
156 // mass balance within the element. this is the
157 // \f$\frac{m}{\partial t}\f$ term if using implicit
158 // euler as time discretization.
159
160 // We might need a more explicit way for
161 // doing the time discretization...
162 auto prevCCStorage = asImp_().computeStorageForCellCenter(problem, scv, prevVolVars);
163 auto curCCStorage = asImp_().computeStorageForCellCenter(problem, scv, curVolVars);
164
165 prevCCStorage *= prevVolVars.extrusionFactor();
166 curCCStorage *= curVolVars.extrusionFactor();
167
168 storage = std::move(curCCStorage);
169 storage -= std::move(prevCCStorage);
170 storage *= scv.volume();
171 storage /= timeLoop_->timeStepSize();
172
173 residual += storage;
174 }
175
177 template<class... Args>
179 {
180 return CellCenterResidualValue(0.0);
181 }
182
184 template<class... Args>
186 {
187 return CellCenterResidualValue(0.0);
188 }
189
195 // \{
196
199 const FVElementGeometry& fvGeometry,
200 const ElementVolumeVariables& elemVolVars,
201 const ElementFaceVariables& elemFaceVars,
202 const ElementBoundaryTypes& bcTypes,
203 const ElementFluxVariablesCache& elemFluxVarsCache,
204 const SubControlVolumeFace& scvf) const
205 {
206 FaceResidualValue residual(0.0);
207 asImp().evalSourceForFace(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, scvf);
208 asImp().evalFluxForFace(residual, this->problem(), element, fvGeometry, elemVolVars, elemFaceVars, bcTypes, elemFluxVarsCache, scvf);
209
210 return residual;
211 }
212
215 const Problem& problem,
216 const Element& element,
217 const FVElementGeometry& fvGeometry,
218 const ElementVolumeVariables& elemVolVars,
219 const ElementFaceVariables& elemFaceVars,
220 const ElementBoundaryTypes& elemBcTypes,
221 const ElementFluxVariablesCache& elemFluxVarsCache,
222 const SubControlVolumeFace& scvf) const
223 {
224 if (!scvf.boundary())
225 residual += asImp_().computeFluxForFace(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, elemFluxVarsCache);
226 else
227 residual += asImp_().computeBoundaryFluxForFace(problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, elemBcTypes, elemFluxVarsCache);
228 }
229
232 const Problem& problem,
233 const Element& element,
234 const FVElementGeometry& fvGeometry,
235 const ElementVolumeVariables& elemVolVars,
236 const ElementFaceVariables& elemFaceVars,
237 const SubControlVolumeFace& scvf) const
238 {
239 // the source term:
240 auto source = asImp_().computeSourceForFace(problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars);
241 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
242 const auto extrusionFactor = elemVolVars[scv].extrusionFactor();
243
244 // multiply by 0.5 because we only consider half of a staggered control volume here
245 source *= 0.5*scv.volume()*extrusionFactor;
246 residual -= source;
247 }
248
250 FaceResidualValue evalStorageForFace(const Element& element,
251 const FVElementGeometry& fvGeometry,
252 const ElementVolumeVariables& prevElemVolVars,
253 const ElementVolumeVariables& curElemVolVars,
254 const ElementFaceVariables& prevElemFaceVars,
255 const ElementFaceVariables& curElemFaceVars,
256 const SubControlVolumeFace& scvf) const
257 {
258 assert(timeLoop_ && "no time loop set for storage term evaluation");
259 FaceResidualValue storage(0.0);
260 asImp().evalStorageForFace(storage, problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, prevElemFaceVars, curElemFaceVars, scvf);
261 return storage;
262 }
263
266 const Problem& problem,
267 const Element& element,
268 const FVElementGeometry& fvGeometry,
269 const ElementVolumeVariables& prevElemVolVars,
270 const ElementVolumeVariables& curElemVolVars,
271 const ElementFaceVariables& prevElemFaceVars,
272 const ElementFaceVariables& curElemFaceVars,
273 const SubControlVolumeFace& scvf) const
274 {
275 FaceResidualValue storage(0.0);
276 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
277 auto prevFaceStorage = asImp_().computeStorageForFace(problem, scvf, prevElemVolVars[scv], prevElemFaceVars);
278 auto curFaceStorage = asImp_().computeStorageForFace(problem, scvf, curElemVolVars[scv], curElemFaceVars);
279
280 storage = std::move(curFaceStorage);
281 storage -= std::move(prevFaceStorage);
282
283 const auto extrusionFactor = curElemVolVars[scv].extrusionFactor();
284
285 // multiply by 0.5 because we only consider half of a staggered control volume here
286 storage *= 0.5*scv.volume()*extrusionFactor;
287 storage /= timeLoop_->timeStepSize();
288
289 residual += storage;
290 }
291
293 bool isStationary() const
294 { return !timeLoop_; }
295
297 const Problem& problem() const
298 { return *problem_; }
299
300protected:
301
302 Implementation &asImp_()
303 { return *static_cast<Implementation*>(this); }
304
305 const Implementation &asImp_() const
306 { return *static_cast<const Implementation*>(this); }
307
308
310 { return *timeLoop_; }
311
312 const TimeLoop& timeLoop() const
313 { return *timeLoop_; }
314
315 Implementation &asImp()
316 { return *static_cast<Implementation*>(this); }
317
318 const Implementation &asImp() const
319 { return *static_cast<const Implementation*>(this); }
320
321private:
322 const Problem* problem_;
323 const TimeLoop* timeLoop_;
324
325};
326
327} // end namespace Dumux
328
329#endif // DUMUX_CC_LOCAL_RESIDUAL_HH
Manages the handling of time dependent problems.
Definition: adapt.hh:29
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:40
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:129
GetPropType< TypeTag, Properties::CellCenterPrimaryVariables > CellCenterResidualValue
Definition: staggeredlocalresidual.hh:62
Implementation & asImp()
Definition: staggeredlocalresidual.hh:315
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:144
Implementation & asImp_()
Definition: staggeredlocalresidual.hh:302
const TimeLoop & timeLoop() const
Definition: staggeredlocalresidual.hh:312
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:198
GetPropType< TypeTag, Properties::FacePrimaryVariables > FaceResidualValue
Definition: staggeredlocalresidual.hh:63
const Problem & problem() const
the problem
Definition: staggeredlocalresidual.hh:297
const Implementation & asImp_() const
Definition: staggeredlocalresidual.hh:305
CellCenterResidualValue evalStorage(Args &&... args) const
for compatibility with FVLocalAssemblerBase
Definition: staggeredlocalresidual.hh:185
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:231
bool isStationary() const
If no solution has been set, we treat the problem as stationary.
Definition: staggeredlocalresidual.hh:293
TimeLoop & timeLoop()
Definition: staggeredlocalresidual.hh:309
StaggeredLocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition: staggeredlocalresidual.hh:67
CellCenterResidualValue evalFluxAndSource(Args &&... args) const
for compatibility with FVLocalAssemblerBase
Definition: staggeredlocalresidual.hh:178
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:265
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:74
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:214
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:112
CellCenterResidualValue ElementResidualVector
Definition: staggeredlocalresidual.hh:64
const Implementation & asImp() const
Definition: staggeredlocalresidual.hh:318
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:95
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:250
Manages the handling of time dependent problems.
Definition: timeloop.hh:72
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
The default time loop for instationary simulations.
Definition: timeloop.hh:117
Declares all properties used in Dumux.