3.6-git
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 const auto fvGeometry = localView(gridGeometry).bind(element);
106 const auto elemVolVars = localView(gridVariables.curGridVolVars()).bind(element, fvGeometry, sol);
107
108 ElementResidualVector storage(fvGeometry.numScv());
109
110 // calculate the amount of conservation each quantity inside
111 // all sub control volumes
112 for (auto&& scv : scvs(fvGeometry))
113 {
114 auto localScvIdx = scv.localDofIndex();
115 const auto& volVars = elemVolVars[scv];
116 storage[localScvIdx] = asImp().computeStorage(problem, scv, volVars);
117 storage[localScvIdx] *= Extrusion::volume(fvGeometry, scv) * volVars.extrusionFactor();
118 }
119
120 return storage;
121 }
122
123 // \}
124
129 // \{
130
143 ElementResidualVector evalStorage(const Element& element,
144 const FVElementGeometry& fvGeometry,
145 const ElementVolumeVariables& prevElemVolVars,
146 const ElementVolumeVariables& curElemVolVars) const
147 {
148 assert(timeLoop_ && "no time loop set for storage term evaluation");
149
150 // initialize the residual vector for all scvs in this element
151 ElementResidualVector residual(fvGeometry.numScv());
152
153 // evaluate the volume terms (storage + source terms)
154 // forward to the local residual specialized for the discretization methods
155 for (auto&& scv : scvs(fvGeometry))
156 asImp().evalStorage(residual, this->problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
157
158 return residual;
159 }
160
162 const FVElementGeometry& fvGeometry,
163 const ElementVolumeVariables& elemVolVars,
164 const ElementFluxVariablesCache& elemFluxVarsCache,
165 const ElementBoundaryTypes &bcTypes) const
166 {
167 // initialize the residual vector for all scvs in this element
168 ElementResidualVector residual(fvGeometry.numScv());
169
170 // evaluate the volume terms (storage + source terms)
171 // forward to the local residual specialized for the discretization methods
172 for (auto&& scv : scvs(fvGeometry))
173 asImp().evalSource(residual, this->problem(), element, fvGeometry, elemVolVars, scv);
174
175 // forward to the local residual specialized for the discretization methods
176 for (auto&& scvf : scvfs(fvGeometry))
177 asImp().evalFlux(residual, this->problem(), element, fvGeometry, elemVolVars, bcTypes, elemFluxVarsCache, scvf);
178
179 return residual;
180 }
181
182 // \}
183
184
189 // \{
190
200 NumEqVector computeStorage(const Problem& problem,
201 const SubControlVolume& scv,
202 const VolumeVariables& volVars) const
203 {
204 DUNE_THROW(Dune::NotImplemented, "This model does not implement a storage method!");
205 }
206
220 NumEqVector computeSource(const Problem& problem,
221 const Element& element,
222 const FVElementGeometry& fvGeometry,
223 const ElementVolumeVariables& elemVolVars,
224 const SubControlVolume &scv) const
225 {
226 NumEqVector source(0.0);
227
228 // add contributions from volume flux sources
229 source += problem.source(element, fvGeometry, elemVolVars, scv);
230
231 // add contribution from possible point sources
232 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
233
234 return source;
235 }
236
251 NumEqVector computeFlux(const Problem& problem,
252 const Element& element,
253 const FVElementGeometry& fvGeometry,
254 const ElementVolumeVariables& elemVolVars,
255 const SubControlVolumeFace& scvf,
256 const ElementFluxVariablesCache& elemFluxVarsCache) const
257 {
258 DUNE_THROW(Dune::NotImplemented, "This model does not implement a flux method!");
259 }
260
261 // \}
262
267 // \{
268
285 const Problem& problem,
286 const Element& element,
287 const FVElementGeometry& fvGeometry,
288 const ElementVolumeVariables& prevElemVolVars,
289 const ElementVolumeVariables& curElemVolVars,
290 const SubControlVolume& scv) const
291 {
292 const auto& curVolVars = curElemVolVars[scv];
293 const auto& prevVolVars = prevElemVolVars[scv];
294
295 // mass balance within the element. this is the
296 // \f$\frac{m}{\partial t}\f$ term if using implicit or explicit
297 // euler as time discretization.
298 //
299 // TODO: We might need a more explicit way for
300 // doing the time discretization...
301
303 NumEqVector prevStorage = asImp().computeStorage(problem, scv, prevVolVars);
304 NumEqVector storage = asImp().computeStorage(problem, scv, curVolVars);
305
306 prevStorage *= prevVolVars.extrusionFactor();
307 storage *= curVolVars.extrusionFactor();
308
309 storage -= prevStorage;
310 storage *= Extrusion::volume(fvGeometry, scv);
311 storage /= timeLoop_->timeStepSize();
312
313 residual[scv.localDofIndex()] += storage;
314 }
315
330 const Problem& problem,
331 const Element& element,
332 const FVElementGeometry& fvGeometry,
333 const ElementVolumeVariables& curElemVolVars,
334 const SubControlVolume& scv) const
335 {
337 const auto& curVolVars = curElemVolVars[scv];
338 NumEqVector source = asImp().computeSource(problem, element, fvGeometry, curElemVolVars, scv);
339 source *= Extrusion::volume(fvGeometry, scv)*curVolVars.extrusionFactor();
340
342 residual[scv.localDofIndex()] -= source;
343 }
344
360 const Problem& problem,
361 const Element& element,
362 const FVElementGeometry& fvGeometry,
363 const ElementVolumeVariables& elemVolVars,
364 const ElementBoundaryTypes& elemBcTypes,
365 const ElementFluxVariablesCache& elemFluxVarsCache,
366 const SubControlVolumeFace& scvf) const {}
367
381 NumEqVector evalFlux(const Problem& problem,
382 const Element& element,
383 const FVElementGeometry& fvGeometry,
384 const ElementVolumeVariables& elemVolVars,
385 const ElementFluxVariablesCache& elemFluxVarsCache,
386 const SubControlVolumeFace& scvf) const
387 {
388 return asImp().evalFlux(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
389 }
390
391 //\}
392
396 // \{
397
399 template<class PartialDerivativeMatrix>
400 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
401 const Problem& problem,
402 const Element& element,
403 const FVElementGeometry& fvGeometry,
404 const VolumeVariables& curVolVars,
405 const SubControlVolume& scv) const
406 {
407 DUNE_THROW(Dune::NotImplemented, "analytic storage derivative");
408 }
409
411 template<class PartialDerivativeMatrix>
412 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
413 const Problem& problem,
414 const Element& element,
415 const FVElementGeometry& fvGeometry,
416 const VolumeVariables& curVolVars,
417 const SubControlVolume& scv) const
418 {
419 DUNE_THROW(Dune::NotImplemented, "analytic source derivative");
420 }
421
423 template<class PartialDerivativeMatrices, class T = TypeTag>
424 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod != DiscretizationMethods::box, void>
425 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
426 const Problem& problem,
427 const Element& element,
428 const FVElementGeometry& fvGeometry,
429 const ElementVolumeVariables& curElemVolVars,
430 const ElementFluxVariablesCache& elemFluxVarsCache,
431 const SubControlVolumeFace& scvf) const
432 {
433 DUNE_THROW(Dune::NotImplemented, "analytic flux derivative for cell-centered models");
434 }
435
437 template<class JacobianMatrix, class T = TypeTag>
438 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
439 addFluxDerivatives(JacobianMatrix& A,
440 const Problem& problem,
441 const Element& element,
442 const FVElementGeometry& fvGeometry,
443 const ElementVolumeVariables& curElemVolVars,
444 const ElementFluxVariablesCache& elemFluxVarsCache,
445 const SubControlVolumeFace& scvf) const
446 {
447 DUNE_THROW(Dune::NotImplemented, "analytic flux derivative for box models");
448 }
449
451 template<class PartialDerivativeMatrices>
452 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
453 const Problem& problem,
454 const Element& element,
455 const FVElementGeometry& fvGeometry,
456 const ElementVolumeVariables& curElemVolVars,
457 const ElementFluxVariablesCache& elemFluxVarsCache,
458 const SubControlVolumeFace& scvf) const
459 {
460 DUNE_THROW(Dune::NotImplemented, "analytic Dirichlet flux derivative");
461 }
462
464 template<class PartialDerivativeMatrices>
465 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
466 const Problem& problem,
467 const Element& element,
468 const FVElementGeometry& fvGeometry,
469 const ElementVolumeVariables& curElemVolVars,
470 const ElementFluxVariablesCache& elemFluxVarsCache,
471 const SubControlVolumeFace& scvf) const
472 {
473 DUNE_THROW(Dune::NotImplemented, "analytic Robin flux derivative");
474 }
475
476 //\}
477
481 // \{
482
484 const Problem& problem() const
485 { return *problem_; }
486
489 const TimeLoop& timeLoop() const
490 { return *timeLoop_; }
491
493 bool isStationary() const
494 { return !timeLoop_; }
495
496 // \}
497protected:
498
499 Implementation &asImp()
500 { return *static_cast<Implementation*>(this); }
501
502 const Implementation &asImp() const
503 { return *static_cast<const Implementation*>(this); }
504
505private:
506 const Problem* problem_;
507 const TimeLoop* timeLoop_;
508};
509
510} // end namespace Dumux
511
512#endif
A helper to deduce a vector with the same size as numbers of equations.
A arithmetic block vector type based on DUNE's reserved vector.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr Box box
Definition: method.hh:136
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:47
Implementation & asImp()
Definition: fvlocalresidual.hh:499
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:220
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:465
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Calculate the source term of the equation.
Definition: fvlocalresidual.hh:200
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 flux local residual, i.e. the deviation of the flux term from zero.
Definition: fvlocalresidual.hh:359
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:484
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:412
const TimeLoop & timeLoop() const
Definition: fvlocalresidual.hh:489
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:400
bool isStationary() const
returns true if the residual is stationary
Definition: fvlocalresidual.hh:493
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:381
const Implementation & asImp() const
Definition: fvlocalresidual.hh:502
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:143
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod !=DiscretizationMethods::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:425
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:329
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::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:439
ElementResidualVector evalFluxAndSource(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &bcTypes) const
Definition: fvlocalresidual.hh:161
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:284
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:452
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Calculate the flux term of the equation.
Definition: fvlocalresidual.hh:251
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:68
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
The default time loop for instationary simulations.
Definition: common/timeloop.hh:113
Declares all properties used in Dumux.
Manages the handling of time dependent problems.