3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fvlocalresidual.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 *****************************************************************************/
24#ifndef DUMUX_FV_LOCAL_RESIDUAL_HH
25#define DUMUX_FV_LOCAL_RESIDUAL_HH
26
27#include <dune/common/exceptions.hh>
28#include <dune/istl/bvector.hh>
29
35
36namespace Dumux {
37
44template<class TypeTag>
46{
51 using Element = typename GridView::template Codim<0>::Entity;
52 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
55 using SubControlVolume = typename GridGeometry::SubControlVolume;
56 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
57 using Extrusion = Extrusion_t<GridGeometry>;
60 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
62 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
65
66public:
69
71 FVLocalResidual(const Problem* problem,
72 const TimeLoop* timeLoop = nullptr)
73 : problem_(problem)
74 , timeLoop_(timeLoop)
75 {}
76
82 // \{
83
98 const Element &element,
99 const GridGeometry& gridGeometry,
100 const GridVariables& gridVariables,
101 const SolutionVector& sol) const
102 {
103 // make sure FVElementGeometry and volume variables are bound to the element
104 auto fvGeometry = localView(gridGeometry);
105 fvGeometry.bind(element);
106
107 auto elemVolVars = localView(gridVariables.curGridVolVars());
108 elemVolVars.bind(element, fvGeometry, sol);
109
110 ElementResidualVector storage(fvGeometry.numScv());
111
112 // calculate the amount of conservation each quantity inside
113 // all sub control volumes
114 for (auto&& scv : scvs(fvGeometry))
115 {
116 auto localScvIdx = scv.localDofIndex();
117 const auto& volVars = elemVolVars[scv];
118 storage[localScvIdx] = asImp().computeStorage(problem, scv, volVars);
119 storage[localScvIdx] *= Extrusion::volume(scv) * volVars.extrusionFactor();
120 }
121
122 return storage;
123 }
124
125 // \}
126
131 // \{
132
145 ElementResidualVector evalStorage(const Element& element,
146 const FVElementGeometry& fvGeometry,
147 const ElementVolumeVariables& prevElemVolVars,
148 const ElementVolumeVariables& curElemVolVars) const
149 {
150 assert(timeLoop_ && "no time loop set for storage term evaluation");
151
152 // initialize the residual vector for all scvs in this element
153 ElementResidualVector residual(fvGeometry.numScv());
154
155 // evaluate the volume terms (storage + source terms)
156 // forward to the local residual specialized for the discretization methods
157 for (auto&& scv : scvs(fvGeometry))
158 asImp().evalStorage(residual, this->problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
159
160 return residual;
161 }
162
164 const FVElementGeometry& fvGeometry,
165 const ElementVolumeVariables& elemVolVars,
166 const ElementFluxVariablesCache& elemFluxVarsCache,
167 const ElementBoundaryTypes &bcTypes) const
168 {
169 // initialize the residual vector for all scvs in this element
170 ElementResidualVector residual(fvGeometry.numScv());
171
172 // evaluate the volume terms (storage + source terms)
173 // forward to the local residual specialized for the discretization methods
174 for (auto&& scv : scvs(fvGeometry))
175 asImp().evalSource(residual, this->problem(), element, fvGeometry, elemVolVars, scv);
176
177 // forward to the local residual specialized for the discretization methods
178 for (auto&& scvf : scvfs(fvGeometry))
179 asImp().evalFlux(residual, this->problem(), element, fvGeometry, elemVolVars, bcTypes, elemFluxVarsCache, scvf);
180
181 return residual;
182 }
183
184 // \}
185
186
191 // \{
192
202 NumEqVector computeStorage(const Problem& problem,
203 const SubControlVolume& scv,
204 const VolumeVariables& volVars) const
205 {
206 DUNE_THROW(Dune::NotImplemented, "This model does not implement a storage method!");
207 }
208
222 NumEqVector computeSource(const Problem& problem,
223 const Element& element,
224 const FVElementGeometry& fvGeometry,
225 const ElementVolumeVariables& elemVolVars,
226 const SubControlVolume &scv) const
227 {
228 NumEqVector source(0.0);
229
230 // add contributions from volume flux sources
231 source += problem.source(element, fvGeometry, elemVolVars, scv);
232
233 // add contribution from possible point sources
234 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
235
236 return source;
237 }
238
253 NumEqVector computeFlux(const Problem& problem,
254 const Element& element,
255 const FVElementGeometry& fvGeometry,
256 const ElementVolumeVariables& elemVolVars,
257 const SubControlVolumeFace& scvf,
258 const ElementFluxVariablesCache& elemFluxVarsCache) const
259 {
260 DUNE_THROW(Dune::NotImplemented, "This model does not implement a flux method!");
261 }
262
263 // \}
264
269 // \{
270
287 const Problem& problem,
288 const Element& element,
289 const FVElementGeometry& fvGeometry,
290 const ElementVolumeVariables& prevElemVolVars,
291 const ElementVolumeVariables& curElemVolVars,
292 const SubControlVolume& scv) const
293 {
294 const auto& curVolVars = curElemVolVars[scv];
295 const auto& prevVolVars = prevElemVolVars[scv];
296
297 // mass balance within the element. this is the
298 // \f$\frac{m}{\partial t}\f$ term if using implicit or explicit
299 // euler as time discretization.
300 //
301 // TODO: We might need a more explicit way for
302 // doing the time discretization...
303
305 NumEqVector prevStorage = asImp().computeStorage(problem, scv, prevVolVars);
306 NumEqVector storage = asImp().computeStorage(problem, scv, curVolVars);
307
308 prevStorage *= prevVolVars.extrusionFactor();
309 storage *= curVolVars.extrusionFactor();
310
311 storage -= prevStorage;
312 storage *= Extrusion::volume(scv);
313 storage /= timeLoop_->timeStepSize();
314
315 residual[scv.localDofIndex()] += storage;
316 }
317
332 const Problem& problem,
333 const Element& element,
334 const FVElementGeometry& fvGeometry,
335 const ElementVolumeVariables& curElemVolVars,
336 const SubControlVolume& scv) const
337 {
339 const auto& curVolVars = curElemVolVars[scv];
340 NumEqVector source = asImp().computeSource(problem, element, fvGeometry, curElemVolVars, scv);
341 source *= Extrusion::volume(scv)*curVolVars.extrusionFactor();
342
344 residual[scv.localDofIndex()] -= source;
345 }
346
364 const Problem& problem,
365 const Element& element,
366 const FVElementGeometry& fvGeometry,
367 const ElementVolumeVariables& elemVolVars,
368 const ElementBoundaryTypes& elemBcTypes,
369 const ElementFluxVariablesCache& elemFluxVarsCache,
370 const SubControlVolumeFace& scvf) const {}
371
385 NumEqVector evalFlux(const Problem& problem,
386 const Element& element,
387 const FVElementGeometry& fvGeometry,
388 const ElementVolumeVariables& elemVolVars,
389 const ElementFluxVariablesCache& elemFluxVarsCache,
390 const SubControlVolumeFace& scvf) const
391 {
392 return asImp().evalFlux(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
393 }
394
395 //\}
396
400 // \{
401
403 template<class PartialDerivativeMatrix>
404 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
405 const Problem& problem,
406 const Element& element,
407 const FVElementGeometry& fvGeometry,
408 const VolumeVariables& curVolVars,
409 const SubControlVolume& scv) const
410 {
411 DUNE_THROW(Dune::NotImplemented, "analytic storage derivative");
412 }
413
415 template<class PartialDerivativeMatrix>
416 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
417 const Problem& problem,
418 const Element& element,
419 const FVElementGeometry& fvGeometry,
420 const VolumeVariables& curVolVars,
421 const SubControlVolume& scv) const
422 {
423 DUNE_THROW(Dune::NotImplemented, "analytic source derivative");
424 }
425
427 template<class PartialDerivativeMatrices, class T = TypeTag>
428 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod != DiscretizationMethod::box, void>
429 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
430 const Problem& problem,
431 const Element& element,
432 const FVElementGeometry& fvGeometry,
433 const ElementVolumeVariables& curElemVolVars,
434 const ElementFluxVariablesCache& elemFluxVarsCache,
435 const SubControlVolumeFace& scvf) const
436 {
437 DUNE_THROW(Dune::NotImplemented, "analytic flux derivative for cell-centered models");
438 }
439
441 template<class JacobianMatrix, class T = TypeTag>
442 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethod::box, void>
443 addFluxDerivatives(JacobianMatrix& A,
444 const Problem& problem,
445 const Element& element,
446 const FVElementGeometry& fvGeometry,
447 const ElementVolumeVariables& curElemVolVars,
448 const ElementFluxVariablesCache& elemFluxVarsCache,
449 const SubControlVolumeFace& scvf) const
450 {
451 DUNE_THROW(Dune::NotImplemented, "analytic flux derivative for box models");
452 }
453
455 template<class PartialDerivativeMatrices>
456 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
457 const Problem& problem,
458 const Element& element,
459 const FVElementGeometry& fvGeometry,
460 const ElementVolumeVariables& curElemVolVars,
461 const ElementFluxVariablesCache& elemFluxVarsCache,
462 const SubControlVolumeFace& scvf) const
463 {
464 DUNE_THROW(Dune::NotImplemented, "analytic Dirichlet flux derivative");
465 }
466
468 template<class PartialDerivativeMatrices>
469 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
470 const Problem& problem,
471 const Element& element,
472 const FVElementGeometry& fvGeometry,
473 const ElementVolumeVariables& curElemVolVars,
474 const ElementFluxVariablesCache& elemFluxVarsCache,
475 const SubControlVolumeFace& scvf) const
476 {
477 DUNE_THROW(Dune::NotImplemented, "analytic Robin flux derivative");
478 }
479
480 //\}
481
485 // \{
486
488 const Problem& problem() const
489 { return *problem_; }
490
493 const TimeLoop& timeLoop() const
494 { return *timeLoop_; }
495
497 bool isStationary() const
498 { return !timeLoop_; }
499
500 // \}
501protected:
502
503 Implementation &asImp()
504 { return *static_cast<Implementation*>(this); }
505
506 const Implementation &asImp() const
507 { return *static_cast<const Implementation*>(this); }
508
509private:
510 const Problem* problem_;
511 const TimeLoop* timeLoop_;
512};
513
514} // end namespace Dumux
515
516#endif
A arithmetic block vector type based on DUNE's reserved vector.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:46
Implementation & asImp()
Definition: fvlocalresidual.hh:503
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Calculate the source term of the equation.
Definition: fvlocalresidual.hh:222
void addRobinFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the derivative of Robin type boundary conditions ("solution dependent Neumann")
Definition: fvlocalresidual.hh:469
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Calculate the source term of the equation.
Definition: fvlocalresidual.hh:202
void evalFlux(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the source local residual, i.e. the deviation of the source term from zero and add to the res...
Definition: fvlocalresidual.hh:363
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod !=DiscretizationMethod::box, void > addFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the derivative of the flux residual.
Definition: fvlocalresidual.hh:429
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:488
void addSourceDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Compute the derivative of the source residual.
Definition: fvlocalresidual.hh:416
const TimeLoop & timeLoop() const
Definition: fvlocalresidual.hh:493
void addStorageDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Compute the derivative of the storage residual.
Definition: fvlocalresidual.hh:404
bool isStationary() const
returns true if the residual is stationary
Definition: fvlocalresidual.hh:497
NumEqVector evalFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the flux local residual, i.e. the deviation of the flux term from zero.
Definition: fvlocalresidual.hh:385
const Implementation & asImp() const
Definition: fvlocalresidual.hh:506
ElementResidualVector evalStorage(const Problem &problem, const Element &element, const GridGeometry &gridGeometry, const GridVariables &gridVariables, const SolutionVector &sol) const
Compute the storage term for the current solution.
Definition: fvlocalresidual.hh:97
FVLocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition: fvlocalresidual.hh:71
ElementResidualVector evalStorage(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars) const
Compute the storage local residual, i.e. the deviation of the storage term from zero for instationary...
Definition: fvlocalresidual.hh:145
void evalSource(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const SubControlVolume &scv) const
Compute the source local residual, i.e. the deviation of the source term from zero.
Definition: fvlocalresidual.hh:331
ElementResidualVector evalFluxAndSource(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &bcTypes) const
Definition: fvlocalresidual.hh:163
void evalStorage(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars, const SubControlVolume &scv) const
Compute the storage local residual, i.e. the deviation of the storage term from zero for instationary...
Definition: fvlocalresidual.hh:286
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethod::box, void > addFluxDerivatives(JacobianMatrix &A, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the derivative of the flux residual for the box method.
Definition: fvlocalresidual.hh:443
void addCCDirichletFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the derivative of the Dirichlet flux residual for cell-centered schemes.
Definition: fvlocalresidual.hh:456
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Calculate the source term of the equation.
Definition: fvlocalresidual.hh:253
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:38
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.