12#ifndef DUMUX_FV_LOCAL_RESIDUAL_HH
13#define DUMUX_FV_LOCAL_RESIDUAL_HH
15#include <dune/common/exceptions.hh>
16#include <dune/istl/bvector.hh>
33template<
class TypeTag>
40 using Element =
typename GridView::template Codim<0>::Entity;
44 using SubControlVolume =
typename GridGeometry::SubControlVolume;
45 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
87 const Element &element,
88 const GridGeometry& gridGeometry,
89 const GridVariables& gridVariables,
90 const SolutionVector& sol)
const
93 const auto fvGeometry =
localView(gridGeometry).bind(element);
94 const auto elemVolVars =
localView(gridVariables.curGridVolVars()).bind(element, fvGeometry, sol);
100 for (
auto&& scv : scvs(fvGeometry))
102 auto localScvIdx = scv.localDofIndex();
103 const auto& volVars = elemVolVars[scv];
104 storage[localScvIdx] =
asImp().computeStorage(
problem, scv, volVars);
105 storage[localScvIdx] *=
Extrusion::volume(fvGeometry, scv) * volVars.extrusionFactor();
132 const FVElementGeometry& fvGeometry,
133 const ElementVolumeVariables& prevElemVolVars,
134 const ElementVolumeVariables& curElemVolVars)
const
136 assert(timeLoop_ &&
"no time loop set for storage term evaluation");
143 for (
auto&& scv : scvs(fvGeometry))
144 asImp().evalStorage(residual, this->
problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
150 const FVElementGeometry& fvGeometry,
151 const ElementVolumeVariables& elemVolVars,
152 const ElementFluxVariablesCache& elemFluxVarsCache,
153 const ElementBoundaryTypes &bcTypes)
const
160 for (
auto&& scv : scvs(fvGeometry))
161 asImp().evalSource(residual, this->
problem(), element, fvGeometry, elemVolVars, scv);
164 for (
auto&& scvf : scvfs(fvGeometry))
165 asImp().evalFlux(residual, this->
problem(), element, fvGeometry, elemVolVars, bcTypes, elemFluxVarsCache, scvf);
189 const SubControlVolume& scv,
190 const VolumeVariables& volVars)
const
192 DUNE_THROW(Dune::NotImplemented,
"This model does not implement a storage method!");
209 const Element& element,
210 const FVElementGeometry& fvGeometry,
211 const ElementVolumeVariables& elemVolVars,
212 const SubControlVolume &scv)
const
214 NumEqVector source(0.0);
217 source +=
problem.source(element, fvGeometry, elemVolVars, scv);
220 source +=
problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
240 const Element& element,
241 const FVElementGeometry& fvGeometry,
242 const ElementVolumeVariables& elemVolVars,
243 const SubControlVolumeFace& scvf,
244 const ElementFluxVariablesCache& elemFluxVarsCache)
const
246 DUNE_THROW(Dune::NotImplemented,
"This model does not implement a flux method!");
274 const Element& element,
275 const FVElementGeometry& fvGeometry,
276 const ElementVolumeVariables& prevElemVolVars,
277 const ElementVolumeVariables& curElemVolVars,
278 const SubControlVolume& scv)
const
280 const auto& curVolVars = curElemVolVars[scv];
281 const auto& prevVolVars = prevElemVolVars[scv];
291 NumEqVector prevStorage =
asImp().computeStorage(
problem, scv, prevVolVars);
292 NumEqVector storage =
asImp().computeStorage(
problem, scv, curVolVars);
294 prevStorage *= prevVolVars.extrusionFactor();
295 storage *= curVolVars.extrusionFactor();
297 storage -= prevStorage;
301 residual[scv.localDofIndex()] += storage;
319 const Element& element,
320 const FVElementGeometry& fvGeometry,
321 const ElementVolumeVariables& curElemVolVars,
322 const SubControlVolume& scv)
const
325 const auto& curVolVars = curElemVolVars[scv];
326 NumEqVector source =
asImp().computeSource(
problem, element, fvGeometry, curElemVolVars, scv);
330 residual[scv.localDofIndex()] -= source;
349 const Element& element,
350 const FVElementGeometry& fvGeometry,
351 const ElementVolumeVariables& elemVolVars,
352 const ElementBoundaryTypes& elemBcTypes,
353 const ElementFluxVariablesCache& elemFluxVarsCache,
354 const SubControlVolumeFace& scvf)
const {}
370 const Element& element,
371 const FVElementGeometry& fvGeometry,
372 const ElementVolumeVariables& elemVolVars,
373 const ElementFluxVariablesCache& elemFluxVarsCache,
374 const SubControlVolumeFace& scvf)
const
376 return asImp().evalFlux(
problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
387 template<
class PartialDerivativeMatrix>
390 const Element& element,
391 const FVElementGeometry& fvGeometry,
392 const VolumeVariables& curVolVars,
393 const SubControlVolume& scv)
const
395 DUNE_THROW(Dune::NotImplemented,
"analytic storage derivative");
399 template<
class PartialDerivativeMatrix>
402 const Element& element,
403 const FVElementGeometry& fvGeometry,
404 const VolumeVariables& curVolVars,
405 const SubControlVolume& scv)
const
407 DUNE_THROW(Dune::NotImplemented,
"analytic source derivative");
411 template<
class PartialDerivativeMatrices,
class T = TypeTag>
415 const Element& element,
416 const FVElementGeometry& fvGeometry,
417 const ElementVolumeVariables& curElemVolVars,
418 const ElementFluxVariablesCache& elemFluxVarsCache,
419 const SubControlVolumeFace& scvf)
const
421 DUNE_THROW(Dune::NotImplemented,
"analytic flux derivative for cell-centered models");
425 template<
class JacobianMatrix,
class T = TypeTag>
429 const Element& element,
430 const FVElementGeometry& fvGeometry,
431 const ElementVolumeVariables& curElemVolVars,
432 const ElementFluxVariablesCache& elemFluxVarsCache,
433 const SubControlVolumeFace& scvf)
const
435 DUNE_THROW(Dune::NotImplemented,
"analytic flux derivative for box models");
439 template<
class PartialDerivativeMatrices>
442 const Element& element,
443 const FVElementGeometry& fvGeometry,
444 const ElementVolumeVariables& curElemVolVars,
445 const ElementFluxVariablesCache& elemFluxVarsCache,
446 const SubControlVolumeFace& scvf)
const
448 DUNE_THROW(Dune::NotImplemented,
"analytic Dirichlet flux derivative");
452 template<
class PartialDerivativeMatrices>
455 const Element& element,
456 const FVElementGeometry& fvGeometry,
457 const ElementVolumeVariables& curElemVolVars,
458 const ElementFluxVariablesCache& elemFluxVarsCache,
459 const SubControlVolumeFace& scvf)
const
461 DUNE_THROW(Dune::NotImplemented,
"analytic Robin flux derivative");
473 {
return *problem_; }
478 {
return *timeLoop_; }
482 {
return !timeLoop_; }
488 {
return *
static_cast<Implementation*
>(
this); }
491 {
return *
static_cast<const Implementation*
>(
this); }
494 const Problem* problem_;
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:35
Implementation & asImp()
Definition: fvlocalresidual.hh:487
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:208
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:453
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Calculate the source term of the equation.
Definition: fvlocalresidual.hh:188
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:347
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:472
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:400
const TimeLoop & timeLoop() const
Definition: fvlocalresidual.hh:477
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:388
bool isStationary() const
returns true if the residual is stationary
Definition: fvlocalresidual.hh:481
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:369
const Implementation & asImp() const
Definition: fvlocalresidual.hh:490
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:86
FVLocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition: fvlocalresidual.hh:60
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:131
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:413
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:317
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:427
ElementResidualVector evalFluxAndSource(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &bcTypes) const
Definition: fvlocalresidual.hh:149
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:272
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:440
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:239
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:26
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
The default time loop for instationary simulations.
Definition: common/timeloop.hh:139
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
Helper classes to compute the integration elements.
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:34
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
constexpr Box box
Definition: method.hh:147
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
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.