12#ifndef DUMUX_BASE_LOCAL_RESIDUAL_HH
13#define DUMUX_BASE_LOCAL_RESIDUAL_HH
17#include <dune/common/exceptions.hh>
19#include <dumux/common/typetraits/localdofs_.hh>
37template<
class TypeTag>
44 using Element =
typename GridView::template Codim<0>::Entity;
51 using ElementVariables =
typename GridVariablesCache::LocalView;
82 const ElementDiscretization& elemDisc,
83 const ElementVariables& prevElemVars,
84 const ElementVariables& curElemVars)
const
86 assert(!this->
isStationary() &&
"no time loop set for storage term evaluation");
95 for (
const auto& scv :
scvs(elemDisc))
96 this->
asImp().evalStorage(residual, this->
problem(), element, elemDisc, prevElemVars, curElemVars, scv);
100 this->
asImp().addToElementStorageResidual(residual, this->
problem(), element, elemDisc, prevElemVars, curElemVars);
114 const ElementDiscretization& elemDisc,
115 const ElementVariables& elemVars)
const
124 for (
const auto& scv :
scvs(elemDisc))
125 this->
asImp().evalSource(residual, this->
problem(), element, elemDisc, elemVars, scv);
129 for (
auto&& scvf : scvfs(elemDisc))
131 this->
asImp().evalFlux(residual, this->
problem(), element, elemDisc, elemVars, scvf);
135 this->
asImp().addToElementFluxAndSourceResidual(residual, this->
problem(), element, elemDisc, elemVars);
138 this->
asImp().addBoundaryFluxIntegral(residual, this->
problem(), elemDisc, elemVars);
146 const Element& element,
147 const ElementDiscretization& elemDisc,
148 const ElementVariables& prevElemVars,
149 const ElementVariables& curElemVars)
const
155 const Element& element,
156 const ElementDiscretization& elemDisc,
157 const ElementVariables& curElemVars)
const
163 const ElementDiscretization& elemDisc,
164 const ElementVariables& elemVars)
const
166 if(!elemDisc.hasBoundaryFaces())
169 for (
const auto& boundaryFace : boundaryFaces(elemDisc))
171 const auto& bcTypes =
problem.boundaryTypes(elemDisc, boundaryFace);
172 if(!bcTypes.hasFluxBoundary())
175 problem.addFEBoundaryFluxIntegral(residual, elemDisc, elemVars, boundaryFace, bcTypes);
177 for(
const auto& scvf : scvfs(elemDisc, boundaryFace))
178 problem.addFVBoundaryFluxIntegral(residual, elemDisc, elemVars, scvf, bcTypes);
200 template<
class SubControlVolume>
202 const ElementVariables& elemVars,
203 const SubControlVolume& scv,
204 bool isPreviousTimeLevel)
const
206 DUNE_THROW(Dune::NotImplemented,
"This model does not implement a storageIntegral method!");
219 template<
class SubControlVolume>
221 const ElementVariables& elemVars,
222 const SubControlVolume& scv)
const
224 NumEqVector source(0.0);
228 source += qpData.weight() *
problem.source(elemDisc, elemVars, qpData.ipData());
230 source *= elemVars[scv].extrusionFactor();
233 const auto& pointSources =
problem.pointSources();
234 if (!pointSources.empty())
235 for (
const auto& context : pointSources.contexts(elemDisc, scv))
237 auto psValues = pointSources.eval(elemDisc, elemVars, context);
254 template<
class SubControlVolumeFace>
256 const ElementVariables& elemVars,
257 const SubControlVolumeFace& scvf)
const
259 DUNE_THROW(Dune::NotImplemented,
"This model does not implement a fluxIntegral method!");
283 template<
class SubControlVolume>
286 const Element& element,
287 const ElementDiscretization& elemDisc,
288 const ElementVariables& prevElemVars,
289 const ElementVariables& curElemVars,
290 const SubControlVolume& scv)
const
300 NumEqVector prevStorage = this->
asImp().storageIntegral(elemDisc, prevElemVars, scv,
true);
301 NumEqVector storage = this->
asImp().storageIntegral(elemDisc, curElemVars, scv,
false);
303 storage -= prevStorage;
304 storage /= timeLoop_->timeStepSize();
306 residual[scv.localDofIndex()] += storage;
321 template<
class SubControlVolume>
324 const Element& element,
325 const ElementDiscretization& elemDisc,
326 const ElementVariables& curElemVars,
327 const SubControlVolume& scv)
const
330 NumEqVector source = this->
asImp().sourceIntegral(elemDisc, curElemVars, scv);
332 residual[scv.localDofIndex()] -= source;
345 template<
class SubControlVolumeFace>
348 const Element& element,
349 const ElementDiscretization& elemDisc,
350 const ElementVariables& elemVars,
351 const SubControlVolumeFace& scvf)
const {}
363 template<
class SubControlVolumeFace>
365 const Element& element,
366 const ElementDiscretization& elemDisc,
367 const ElementVariables& elemVars,
368 const SubControlVolumeFace& scvf)
const
370 return asImp().evalFlux(
problem, element, elemDisc, elemVars, scvf);
382 {
return *problem_; }
387 {
return *timeLoop_; }
391 {
return !timeLoop_; }
397 {
return *
static_cast<Implementation*
>(
this); }
400 {
return *
static_cast<const Implementation*
>(
this); }
GVC GridVariablesCache
export type of the grid variables cache
Definition discretization/gridvariables.hh:33
ElementResidualVector evalStorage(const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &prevElemVars, const ElementVariables &curElemVars) const
Compute the storage local residual, i.e. the deviation of the storage term from zero for instationary...
Definition assembly/localresidual.hh:81
LocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition assembly/localresidual.hh:59
void evalSource(ElementResidualVector &residual, const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &curElemVars, const SubControlVolume &scv) const
Compute the source local residual, i.e. the deviation of the source term from zero.
Definition assembly/localresidual.hh:322
NumEqVector sourceIntegral(const ElementDiscretization &elemDisc, const ElementVariables &elemVars, const SubControlVolume &scv) const
Calculate the source term integral of the equation.
Definition assembly/localresidual.hh:220
NumEqVector evalFlux(const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &elemVars, const SubControlVolumeFace &scvf) const
Compute the fluxes of the local residual.
Definition assembly/localresidual.hh:364
ElementResidualVector evalFluxAndSource(const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &elemVars) const
Compute the flux and source.
Definition assembly/localresidual.hh:113
const Problem & problem() const
the problem
Definition assembly/localresidual.hh:381
void addToElementStorageResidual(ElementResidualVector &residual, const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &prevElemVars, const ElementVariables &curElemVars) const
add additional storage contributions (e.g. hybrid CVFE or FE schemes)
Definition assembly/localresidual.hh:144
const TimeLoop & timeLoop() const
Definition assembly/localresidual.hh:386
bool isStationary() const
returns true if the residual is stationary
Definition assembly/localresidual.hh:390
const Implementation & asImp() const
Definition assembly/localresidual.hh:399
void evalFlux(ElementResidualVector &residual, const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &elemVars, const SubControlVolumeFace &scvf) const
Compute the fluxes of the local residual.
Definition assembly/localresidual.hh:346
NumEqVector storageIntegral(const ElementDiscretization &elemDisc, const ElementVariables &elemVars, const SubControlVolume &scv, bool isPreviousTimeLevel) const
Calculate the source term integral of the equation.
Definition assembly/localresidual.hh:201
void addBoundaryFluxIntegral(ElementResidualVector &residual, const Problem &problem, const ElementDiscretization &elemDisc, const ElementVariables &elemVars) const
add boundary flux contributions
Definition assembly/localresidual.hh:161
void addToElementFluxAndSourceResidual(ElementResidualVector &residual, const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &curElemVars) const
add additional flux and source contributions (e.g. hybrid CVFE or FE schemes)
Definition assembly/localresidual.hh:153
Implementation & asImp()
Definition assembly/localresidual.hh:396
void evalStorage(ElementResidualVector &residual, const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &prevElemVars, const ElementVariables &curElemVars, const SubControlVolume &scv) const
Compute the storage local residual, i.e. the deviation of the storage term from zero for instationary...
Definition assembly/localresidual.hh:284
NumEqVector fluxIntegral(const ElementDiscretization &elemDisc, const ElementVariables &elemVars, const SubControlVolumeFace &scvf) const
Calculate the flux integral of the equation.
Definition assembly/localresidual.hh:255
ReservedBlockVector< NumEqVector, Dumux::Detail::LocalDofs::maxNumLocalDofs< ElementDiscretization >()> ElementResidualVector
the container storing all element residuals
Definition assembly/localresidual.hh:56
Base class for all standard finite volume or finite element problems.
Definition common/problem.hh:39
A arithmetic block vector type based on DUNE's reserved vector.
Definition reservedblockvector.hh:26
Manages the handling of time dependent problems.
Definition common/timeloop.hh:84
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.
Concept for finite-volume discretizations.
Definition concepts.hh:26
Concepts for discretization types.
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
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:296
auto quadratureRule(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolume &scv, QuadratureRules::MidpointQuadrature)
Midpoint quadrature for scv.
Definition quadraturerules.hh:159
Definition assembly/assembler.hh:44
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition localdof.hh:82
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:236
A helper to deduce a vector with the same size as numbers of equations.
Point source types and helpers for handling point sources.
Quadrature rules over sub-control volumes and sub-control volume faces.
A arithmetic block vector type based on DUNE's reserved vector.