25#ifndef DUMUX_CVFE_LOCAL_RESIDUAL_HH
26#define DUMUX_CVFE_LOCAL_RESIDUAL_HH
28#include <dune/common/std/type_traits.hh>
29#include <dune/geometry/type.hh>
30#include <dune/istl/matrix.hh>
39template<
class Imp,
class P,
class G,
class S,
class V>
41 std::declval<Imp>().computeStorage(
42 std::declval<P>(), std::declval<G>(), std::declval<S>(), std::declval<V>(),
true
46template<
class Imp,
class P,
class G,
class S,
class V>
48{
return Dune::Std::is_detected<TimeInfoInterfaceCVFEDetector, Imp, P, G, S, V>::value; }
52 std::declval<Imp>().isOverlapping()
57{
return Dune::Std::is_detected<SCVFIsOverlappingDetector, Imp>::value; }
70template<
class TypeTag>
78 using GridView =
typename GridGeometry::GridView;
79 using Element =
typename GridView::template Codim<0>::Entity;
81 using FVElementGeometry =
typename GridGeometry::LocalView;
83 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
84 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
85 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
86 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
94 using ParentType::ParentType;
99 const Element& element,
100 const FVElementGeometry& fvGeometry,
101 const ElementVolumeVariables& elemVolVars,
102 const ElementBoundaryTypes& elemBcTypes,
103 const ElementFluxVariablesCache& elemFluxVarsCache,
104 const SubControlVolumeFace& scvf)
const
106 const auto flux =
evalFlux(
problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
107 if (!scvf.boundary())
109 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
110 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
111 residual[insideScv.localDofIndex()] += flux;
114 if constexpr (Detail::hasScvfIsOverlapping<SubControlVolumeFace>())
116 if (!scvf.isOverlapping())
117 residual[outsideScv.localDofIndex()] -= flux;
120 residual[outsideScv.localDofIndex()] -= flux;
124 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
125 residual[insideScv.localDofIndex()] += flux;
131 const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const ElementVolumeVariables& elemVolVars,
134 const ElementBoundaryTypes& elemBcTypes,
135 const ElementFluxVariablesCache& elemFluxVarsCache,
136 const SubControlVolumeFace& scvf)
const
138 NumEqVector flux(0.0);
141 if (!scvf.boundary())
142 flux += this->
asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
147 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
148 const auto& bcTypes = elemBcTypes.get(fvGeometry, scv);
153 if (bcTypes.hasNeumann())
155 auto neumannFluxes =
problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
158 neumannFluxes *= Extrusion::area(fvGeometry, scvf)*elemVolVars[scv].extrusionFactor();
161 for (
int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
162 if (bcTypes.isNeumann(eqIdx))
163 flux[eqIdx] += neumannFluxes[eqIdx];
188 const Element& element,
189 const FVElementGeometry& fvGeometry,
190 const ElementVolumeVariables& prevElemVolVars,
191 const ElementVolumeVariables& curElemVolVars,
192 const SubControlVolume& scv)
const
194 const auto& curVolVars = curElemVolVars[scv];
195 const auto& prevVolVars = prevElemVolVars[scv];
200 NumEqVector prevStorage = computeStorageImpl_(
problem, fvGeometry, scv, prevVolVars,
true);
201 NumEqVector storage = computeStorageImpl_(
problem, fvGeometry, scv, curVolVars,
false);
203 prevStorage *= prevVolVars.extrusionFactor();
204 storage *= curVolVars.extrusionFactor();
206 storage -= prevStorage;
210 residual[scv.localDofIndex()] += storage;
215 const FVElementGeometry& fvGeometry,
216 const SubControlVolume& scv,
217 const VolumeVariables& volVars,
218 [[maybe_unused]]
bool isPreviousTimeStep)
const
220 if constexpr (Detail::hasTimeInfoInterfaceCVFE<Implementation, Problem, FVElementGeometry, SubControlVolume, VolumeVariables>())
221 return this->
asImp().computeStorage(problem, fvGeometry, scv, volVars, isPreviousTimeStep);
223 return this->
asImp().computeStorage(problem, scv, volVars);
The element-wise residual for finite volume schemes.
A helper to deduce a vector with the same size as numbers of equations.
Helper classes to compute the integration elements.
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
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
Distance implementation details.
Definition: cvfelocalresidual.hh:37
constexpr bool hasScvfIsOverlapping()
Definition: cvfelocalresidual.hh:56
constexpr bool hasTimeInfoInterfaceCVFE()
Definition: cvfelocalresidual.hh:47
decltype(std::declval< Imp >().computeStorage(std::declval< P >(), std::declval< G >(), std::declval< S >(), std::declval< V >(), true)) TimeInfoInterfaceCVFEDetector
Definition: cvfelocalresidual.hh:44
decltype(std::declval< Imp >().isOverlapping()) SCVFIsOverlappingDetector
Definition: cvfelocalresidual.hh:53
The element-wise residual for control-volume finite element schemes.
Definition: cvfelocalresidual.hh:72
NumEqVector evalFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
evaluate flux residuals for one sub control volume face
Definition: cvfelocalresidual.hh:130
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: cvfelocalresidual.hh:186
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
evaluate flux residuals for one sub control volume face and add to residual
Definition: cvfelocalresidual.hh:97
typename ParentType::ElementResidualVector ElementResidualVector
Definition: cvfelocalresidual.hh:93
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:47
Implementation & asImp()
Definition: fvlocalresidual.hh:499
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:484
ReservedBlockVector< NumEqVector, FVElementGeometry::maxNumElementScvs > ElementResidualVector
the container storing all element residuals
Definition: fvlocalresidual.hh:69
const TimeLoop & timeLoop() const
Definition: fvlocalresidual.hh:489
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
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
Declares all properties used in Dumux.