24#ifndef DUMUX_FV_LOCAL_RESIDUAL_HH
25#define DUMUX_FV_LOCAL_RESIDUAL_HH
27#include <dune/common/exceptions.hh>
28#include <dune/istl/bvector.hh>
45template<
class TypeTag>
52 using Element =
typename GridView::template Codim<0>::Entity;
56 using SubControlVolume =
typename GridGeometry::SubControlVolume;
57 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
99 const Element &element,
100 const GridGeometry& gridGeometry,
101 const GridVariables& gridVariables,
102 const SolutionVector& sol)
const
105 auto fvGeometry =
localView(gridGeometry);
106 fvGeometry.bind(element);
108 auto elemVolVars =
localView(gridVariables.curGridVolVars());
109 elemVolVars.bind(element, fvGeometry, sol);
115 for (
auto&& scv : scvs(fvGeometry))
117 auto localScvIdx = scv.localDofIndex();
118 const auto& volVars = elemVolVars[scv];
119 storage[localScvIdx] =
asImp().computeStorage(
problem, scv, volVars);
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& prevElemVolVars,
149 const ElementVolumeVariables& curElemVolVars)
const
151 assert(timeLoop_ &&
"no time loop set for storage term evaluation");
158 for (
auto&& scv : scvs(fvGeometry))
159 asImp().evalStorage(residual, this->
problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
165 const FVElementGeometry& fvGeometry,
166 const ElementVolumeVariables& elemVolVars,
167 const ElementFluxVariablesCache& elemFluxVarsCache,
168 const ElementBoundaryTypes &bcTypes)
const
175 for (
auto&& scv : scvs(fvGeometry))
176 asImp().evalSource(residual, this->
problem(), element, fvGeometry, elemVolVars, scv);
179 for (
auto&& scvf : scvfs(fvGeometry))
180 asImp().evalFlux(residual, this->
problem(), element, fvGeometry, elemVolVars, bcTypes, elemFluxVarsCache, scvf);
204 const SubControlVolume& scv,
205 const VolumeVariables& volVars)
const
207 DUNE_THROW(Dune::NotImplemented,
"This model does not implement a storage method!");
224 const Element& element,
225 const FVElementGeometry& fvGeometry,
226 const ElementVolumeVariables& elemVolVars,
227 const SubControlVolume &scv)
const
229 NumEqVector source(0.0);
232 source +=
problem.source(element, fvGeometry, elemVolVars, scv);
235 source +=
problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
255 const Element& element,
256 const FVElementGeometry& fvGeometry,
257 const ElementVolumeVariables& elemVolVars,
258 const SubControlVolumeFace& scvf,
259 const ElementFluxVariablesCache& elemFluxVarsCache)
const
261 DUNE_THROW(Dune::NotImplemented,
"This model does not implement a flux method!");
289 const Element& element,
290 const FVElementGeometry& fvGeometry,
291 const ElementVolumeVariables& prevElemVolVars,
292 const ElementVolumeVariables& curElemVolVars,
293 const SubControlVolume& scv)
const
295 const auto& curVolVars = curElemVolVars[scv];
296 const auto& prevVolVars = prevElemVolVars[scv];
306 NumEqVector prevStorage =
asImp().computeStorage(
problem, scv, prevVolVars);
307 NumEqVector storage =
asImp().computeStorage(
problem, scv, curVolVars);
309 prevStorage *= prevVolVars.extrusionFactor();
310 storage *= curVolVars.extrusionFactor();
312 storage -= prevStorage;
316 residual[scv.localDofIndex()] += storage;
334 const Element& element,
335 const FVElementGeometry& fvGeometry,
336 const ElementVolumeVariables& curElemVolVars,
337 const SubControlVolume& scv)
const
340 const auto& curVolVars = curElemVolVars[scv];
341 NumEqVector source =
asImp().computeSource(
problem, element, fvGeometry, curElemVolVars, scv);
345 residual[scv.localDofIndex()] -= source;
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 {}
387 const Element& element,
388 const FVElementGeometry& fvGeometry,
389 const ElementVolumeVariables& elemVolVars,
390 const ElementFluxVariablesCache& elemFluxVarsCache,
391 const SubControlVolumeFace& scvf)
const
393 return asImp().evalFlux(
problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
404 template<
class PartialDerivativeMatrix>
407 const Element& element,
408 const FVElementGeometry& fvGeometry,
409 const VolumeVariables& curVolVars,
410 const SubControlVolume& scv)
const
412 DUNE_THROW(Dune::NotImplemented,
"analytic storage derivative");
416 template<
class PartialDerivativeMatrix>
419 const Element& element,
420 const FVElementGeometry& fvGeometry,
421 const VolumeVariables& curVolVars,
422 const SubControlVolume& scv)
const
424 DUNE_THROW(Dune::NotImplemented,
"analytic source derivative");
428 template<
class PartialDerivativeMatrices,
class T = TypeTag>
432 const Element& element,
433 const FVElementGeometry& fvGeometry,
434 const ElementVolumeVariables& curElemVolVars,
435 const ElementFluxVariablesCache& elemFluxVarsCache,
436 const SubControlVolumeFace& scvf)
const
438 DUNE_THROW(Dune::NotImplemented,
"analytic flux derivative for cell-centered models");
442 template<
class JacobianMatrix,
class T = TypeTag>
446 const Element& element,
447 const FVElementGeometry& fvGeometry,
448 const ElementVolumeVariables& curElemVolVars,
449 const ElementFluxVariablesCache& elemFluxVarsCache,
450 const SubControlVolumeFace& scvf)
const
452 DUNE_THROW(Dune::NotImplemented,
"analytic flux derivative for box models");
456 template<
class PartialDerivativeMatrices>
459 const Element& element,
460 const FVElementGeometry& fvGeometry,
461 const ElementVolumeVariables& curElemVolVars,
462 const ElementFluxVariablesCache& elemFluxVarsCache,
463 const SubControlVolumeFace& scvf)
const
465 DUNE_THROW(Dune::NotImplemented,
"analytic Dirichlet flux derivative");
469 template<
class PartialDerivativeMatrices>
472 const Element& element,
473 const FVElementGeometry& fvGeometry,
474 const ElementVolumeVariables& curElemVolVars,
475 const ElementFluxVariablesCache& elemFluxVarsCache,
476 const SubControlVolumeFace& scvf)
const
478 DUNE_THROW(Dune::NotImplemented,
"analytic Robin flux derivative");
490 {
return *problem_; }
495 {
return *timeLoop_; }
499 {
return !timeLoop_; }
505 {
return *
static_cast<Implementation*
>(
this); }
508 {
return *
static_cast<const Implementation*
>(
this); }
511 const Problem* problem_;
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
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.