3.4
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
36
37namespace Dumux {
38
45template<class TypeTag>
47{
52 using Element = typename GridView::template Codim<0>::Entity;
53 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
56 using SubControlVolume = typename GridGeometry::SubControlVolume;
57 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
58 using Extrusion = Extrusion_t<GridGeometry>;
61 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
63 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
66
67public:
70
72 FVLocalResidual(const Problem* problem,
73 const TimeLoop* timeLoop = nullptr)
74 : problem_(problem)
75 , timeLoop_(timeLoop)
76 {}
77
83 // \{
84
99 const Element &element,
100 const GridGeometry& gridGeometry,
101 const GridVariables& gridVariables,
102 const SolutionVector& sol) const
103 {
104 // make sure FVElementGeometry and volume variables are bound to the element
105 auto fvGeometry = localView(gridGeometry);
106 fvGeometry.bind(element);
107
108 auto elemVolVars = localView(gridVariables.curGridVolVars());
109 elemVolVars.bind(element, fvGeometry, sol);
110
111 ElementResidualVector storage(fvGeometry.numScv());
112
113 // calculate the amount of conservation each quantity inside
114 // all sub control volumes
115 for (auto&& scv : scvs(fvGeometry))
116 {
117 auto localScvIdx = scv.localDofIndex();
118 const auto& volVars = elemVolVars[scv];
119 storage[localScvIdx] = asImp().computeStorage(problem, scv, volVars);
120 storage[localScvIdx] *= Extrusion::volume(scv) * volVars.extrusionFactor();
121 }
122
123 return storage;
124 }
125
126 // \}
127
132 // \{
133
146 ElementResidualVector evalStorage(const Element& element,
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& prevElemVolVars,
149 const ElementVolumeVariables& curElemVolVars) const
150 {
151 assert(timeLoop_ && "no time loop set for storage term evaluation");
152
153 // initialize the residual vector for all scvs in this element
154 ElementResidualVector residual(fvGeometry.numScv());
155
156 // evaluate the volume terms (storage + source terms)
157 // forward to the local residual specialized for the discretization methods
158 for (auto&& scv : scvs(fvGeometry))
159 asImp().evalStorage(residual, this->problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
160
161 return residual;
162 }
163
165 const FVElementGeometry& fvGeometry,
166 const ElementVolumeVariables& elemVolVars,
167 const ElementFluxVariablesCache& elemFluxVarsCache,
168 const ElementBoundaryTypes &bcTypes) const
169 {
170 // initialize the residual vector for all scvs in this element
171 ElementResidualVector residual(fvGeometry.numScv());
172
173 // evaluate the volume terms (storage + source terms)
174 // forward to the local residual specialized for the discretization methods
175 for (auto&& scv : scvs(fvGeometry))
176 asImp().evalSource(residual, this->problem(), element, fvGeometry, elemVolVars, scv);
177
178 // forward to the local residual specialized for the discretization methods
179 for (auto&& scvf : scvfs(fvGeometry))
180 asImp().evalFlux(residual, this->problem(), element, fvGeometry, elemVolVars, bcTypes, elemFluxVarsCache, scvf);
181
182 return residual;
183 }
184
185 // \}
186
187
192 // \{
193
203 NumEqVector computeStorage(const Problem& problem,
204 const SubControlVolume& scv,
205 const VolumeVariables& volVars) const
206 {
207 DUNE_THROW(Dune::NotImplemented, "This model does not implement a storage method!");
208 }
209
223 NumEqVector computeSource(const Problem& problem,
224 const Element& element,
225 const FVElementGeometry& fvGeometry,
226 const ElementVolumeVariables& elemVolVars,
227 const SubControlVolume &scv) const
228 {
229 NumEqVector source(0.0);
230
231 // add contributions from volume flux sources
232 source += problem.source(element, fvGeometry, elemVolVars, scv);
233
234 // add contribution from possible point sources
235 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
236
237 return source;
238 }
239
254 NumEqVector computeFlux(const Problem& problem,
255 const Element& element,
256 const FVElementGeometry& fvGeometry,
257 const ElementVolumeVariables& elemVolVars,
258 const SubControlVolumeFace& scvf,
259 const ElementFluxVariablesCache& elemFluxVarsCache) const
260 {
261 DUNE_THROW(Dune::NotImplemented, "This model does not implement a flux method!");
262 }
263
264 // \}
265
270 // \{
271
288 const Problem& problem,
289 const Element& element,
290 const FVElementGeometry& fvGeometry,
291 const ElementVolumeVariables& prevElemVolVars,
292 const ElementVolumeVariables& curElemVolVars,
293 const SubControlVolume& scv) const
294 {
295 const auto& curVolVars = curElemVolVars[scv];
296 const auto& prevVolVars = prevElemVolVars[scv];
297
298 // mass balance within the element. this is the
299 // \f$\frac{m}{\partial t}\f$ term if using implicit or explicit
300 // euler as time discretization.
301 //
302 // TODO: We might need a more explicit way for
303 // doing the time discretization...
304
306 NumEqVector prevStorage = asImp().computeStorage(problem, scv, prevVolVars);
307 NumEqVector storage = asImp().computeStorage(problem, scv, curVolVars);
308
309 prevStorage *= prevVolVars.extrusionFactor();
310 storage *= curVolVars.extrusionFactor();
311
312 storage -= prevStorage;
313 storage *= Extrusion::volume(scv);
314 storage /= timeLoop_->timeStepSize();
315
316 residual[scv.localDofIndex()] += storage;
317 }
318
333 const Problem& problem,
334 const Element& element,
335 const FVElementGeometry& fvGeometry,
336 const ElementVolumeVariables& curElemVolVars,
337 const SubControlVolume& scv) const
338 {
340 const auto& curVolVars = curElemVolVars[scv];
341 NumEqVector source = asImp().computeSource(problem, element, fvGeometry, curElemVolVars, scv);
342 source *= Extrusion::volume(scv)*curVolVars.extrusionFactor();
343
345 residual[scv.localDofIndex()] -= source;
346 }
347
365 const Problem& problem,
366 const Element& element,
367 const FVElementGeometry& fvGeometry,
368 const ElementVolumeVariables& elemVolVars,
369 const ElementBoundaryTypes& elemBcTypes,
370 const ElementFluxVariablesCache& elemFluxVarsCache,
371 const SubControlVolumeFace& scvf) const {}
372
386 NumEqVector evalFlux(const Problem& problem,
387 const Element& element,
388 const FVElementGeometry& fvGeometry,
389 const ElementVolumeVariables& elemVolVars,
390 const ElementFluxVariablesCache& elemFluxVarsCache,
391 const SubControlVolumeFace& scvf) const
392 {
393 return asImp().evalFlux(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
394 }
395
396 //\}
397
401 // \{
402
404 template<class PartialDerivativeMatrix>
405 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
406 const Problem& problem,
407 const Element& element,
408 const FVElementGeometry& fvGeometry,
409 const VolumeVariables& curVolVars,
410 const SubControlVolume& scv) const
411 {
412 DUNE_THROW(Dune::NotImplemented, "analytic storage derivative");
413 }
414
416 template<class PartialDerivativeMatrix>
417 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
418 const Problem& problem,
419 const Element& element,
420 const FVElementGeometry& fvGeometry,
421 const VolumeVariables& curVolVars,
422 const SubControlVolume& scv) const
423 {
424 DUNE_THROW(Dune::NotImplemented, "analytic source derivative");
425 }
426
428 template<class PartialDerivativeMatrices, class T = TypeTag>
429 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod != DiscretizationMethod::box, void>
430 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
431 const Problem& problem,
432 const Element& element,
433 const FVElementGeometry& fvGeometry,
434 const ElementVolumeVariables& curElemVolVars,
435 const ElementFluxVariablesCache& elemFluxVarsCache,
436 const SubControlVolumeFace& scvf) const
437 {
438 DUNE_THROW(Dune::NotImplemented, "analytic flux derivative for cell-centered models");
439 }
440
442 template<class JacobianMatrix, class T = TypeTag>
443 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethod::box, void>
444 addFluxDerivatives(JacobianMatrix& A,
445 const Problem& problem,
446 const Element& element,
447 const FVElementGeometry& fvGeometry,
448 const ElementVolumeVariables& curElemVolVars,
449 const ElementFluxVariablesCache& elemFluxVarsCache,
450 const SubControlVolumeFace& scvf) const
451 {
452 DUNE_THROW(Dune::NotImplemented, "analytic flux derivative for box models");
453 }
454
456 template<class PartialDerivativeMatrices>
457 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
458 const Problem& problem,
459 const Element& element,
460 const FVElementGeometry& fvGeometry,
461 const ElementVolumeVariables& curElemVolVars,
462 const ElementFluxVariablesCache& elemFluxVarsCache,
463 const SubControlVolumeFace& scvf) const
464 {
465 DUNE_THROW(Dune::NotImplemented, "analytic Dirichlet flux derivative");
466 }
467
469 template<class PartialDerivativeMatrices>
470 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
471 const Problem& problem,
472 const Element& element,
473 const FVElementGeometry& fvGeometry,
474 const ElementVolumeVariables& curElemVolVars,
475 const ElementFluxVariablesCache& elemFluxVarsCache,
476 const SubControlVolumeFace& scvf) const
477 {
478 DUNE_THROW(Dune::NotImplemented, "analytic Robin flux derivative");
479 }
480
481 //\}
482
486 // \{
487
489 const Problem& problem() const
490 { return *problem_; }
491
494 const TimeLoop& timeLoop() const
495 { return *timeLoop_; }
496
498 bool isStationary() const
499 { return !timeLoop_; }
500
501 // \}
502protected:
503
504 Implementation &asImp()
505 { return *static_cast<Implementation*>(this); }
506
507 const Implementation &asImp() const
508 { return *static_cast<const Implementation*>(this); }
509
510private:
511 const Problem* problem_;
512 const TimeLoop* timeLoop_;
513};
514
515} // end namespace Dumux
516
517#endif
A arithmetic block vector type based on DUNE's reserved vector.
A helper to deduce a vector with the same size as numbers of equations.
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
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
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
Definition: propertysystem.hh:150
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:47
Implementation & asImp()
Definition: fvlocalresidual.hh:504
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:223
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:470
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Calculate the source term of the equation.
Definition: fvlocalresidual.hh:203
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:364
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:430
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:489
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:417
const TimeLoop & timeLoop() const
Definition: fvlocalresidual.hh:494
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:405
bool isStationary() const
returns true if the residual is stationary
Definition: fvlocalresidual.hh:498
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:386
const Implementation & asImp() const
Definition: fvlocalresidual.hh:507
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:98
FVLocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition: fvlocalresidual.hh:72
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:146
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:332
ElementResidualVector evalFluxAndSource(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &bcTypes) const
Definition: fvlocalresidual.hh:164
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:287
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:444
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:457
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:254
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.