13#ifndef DUMUX_CVFE_LOCAL_RESIDUAL_HH
14#define DUMUX_CVFE_LOCAL_RESIDUAL_HH
16#include <dune/common/std/type_traits.hh>
17#include <dune/geometry/type.hh>
18#include <dune/istl/matrix.hh>
20#include <dumux/common/typetraits/localdofs_.hh>
21#include <dumux/common/typetraits/boundary_.hh>
30template<
class Imp,
class P,
class G,
class S,
class V>
32 std::declval<Imp>().computeStorage(
33 std::declval<P>(), std::declval<G>(), std::declval<S>(), std::declval<V>(),
true
37template<
class Imp,
class P,
class G,
class S,
class V>
39{
return Dune::Std::is_detected<TimeInfoInterfaceCVFEDetector, Imp, P, G, S, V>::value; }
43 std::declval<Imp>().isOverlapping()
48{
return Dune::Std::is_detected<SCVFIsOverlappingDetector, Imp>::value; }
61template<
class TypeTag>
69 using GridView =
typename GridGeometry::GridView;
70 using Element =
typename GridView::template Codim<0>::Entity;
72 using FVElementGeometry =
typename GridGeometry::LocalView;
74 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
75 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
76 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
77 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
84 typename Element::Geometry::GlobalCoordinate,
85 typename SubControlVolumeFace::Traits::LocalIndexType>;
89 using ParentType::ParentType;
105 const FVElementGeometry& fvGeometry,
106 const ElementVolumeVariables& prevElemVolVars,
107 const ElementVolumeVariables& curElemVolVars)
const
109 assert(!this->
isStationary() &&
"no time loop set for storage term evaluation");
116 for (
const auto& scv :
scvs(fvGeometry))
117 this->
asImp().evalStorage(residual, this->
problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
120 this->
asImp().addToElementStorageResidual(residual, this->
problem(), element, fvGeometry, prevElemVolVars, curElemVolVars);
137 const FVElementGeometry& fvGeometry,
138 const ElementVolumeVariables& elemVolVars,
139 const ElementFluxVariablesCache& elemFluxVarsCache,
140 const ElementBoundaryTypes &bcTypes)
const
147 for (
const auto& scv :
scvs(fvGeometry))
148 this->
asImp().evalSource(residual, this->
problem(), element, fvGeometry, elemVolVars, scv);
151 for (
auto&& scvf : scvfs(fvGeometry))
152 this->
asImp().evalFlux(residual, this->
problem(), element, fvGeometry, elemVolVars, bcTypes, elemFluxVarsCache, scvf);
155 this->
asImp().addToElementFluxAndSourceResidual(residual, this->
problem(), element, fvGeometry, elemVolVars, elemFluxVarsCache, bcTypes);
163 const Element& element,
164 const FVElementGeometry& fvGeometry,
165 const ElementVolumeVariables& prevElemVolVars,
166 const ElementVolumeVariables& curElemVolVars)
const
172 const Element& element,
173 const FVElementGeometry& fvGeometry,
174 const ElementVolumeVariables& curElemVolVars,
175 const ElementFluxVariablesCache& elemFluxVarsCache,
176 const ElementBoundaryTypes &bcTypes)
const
182 const Element& element,
183 const FVElementGeometry& fvGeometry,
184 const ElementVolumeVariables& elemVolVars,
185 const ElementBoundaryTypes& elemBcTypes,
186 const ElementFluxVariablesCache& elemFluxVarsCache,
187 const SubControlVolumeFace& scvf)
const
189 const auto flux = this->
asImp().evalFlux(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
190 if (!scvf.boundary())
192 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
193 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
194 residual[insideScv.localDofIndex()] += flux;
197 if constexpr (Detail::hasScvfIsOverlapping<SubControlVolumeFace>())
199 if (!scvf.isOverlapping())
200 residual[outsideScv.localDofIndex()] -= flux;
203 residual[outsideScv.localDofIndex()] -= flux;
207 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
208 residual[insideScv.localDofIndex()] += flux;
214 const Element& element,
215 const FVElementGeometry& fvGeometry,
216 const ElementVolumeVariables& elemVolVars,
217 const ElementBoundaryTypes& elemBcTypes,
218 const ElementFluxVariablesCache& elemFluxVarsCache,
219 const SubControlVolumeFace& scvf)
const
221 NumEqVector flux(0.0);
224 if (!scvf.boundary())
225 flux += this->
asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
230 const auto& bcTypes = bcTypes_(
problem, fvGeometry, scvf, elemBcTypes);
235 if (bcTypes.hasNeumann())
237 NumEqVector boundaryFluxes;
238 if constexpr (Dumux::Detail::hasProblemBoundaryFluxFunction<Problem, FVElementGeometry, ElementVolumeVariables, ElementFluxVariablesCache, FaceIpData>())
240 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
241 boundaryFluxes =
problem.boundaryFlux(fvGeometry, elemVolVars, elemFluxVarsCache,
242 FaceIpData(fluxVarsCache.ipLocal(), fluxVarsCache.ipGlobal(), scvf.unitOuterNormal(), scvf.index()));
245 boundaryFluxes =
problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
248 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
249 boundaryFluxes *= Extrusion::area(fvGeometry, scvf)*elemVolVars[scv].extrusionFactor();
252 for (
int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
253 if (bcTypes.isNeumann(eqIdx))
254 flux[eqIdx] += boundaryFluxes[eqIdx];
278 const Element& element,
279 const FVElementGeometry& fvGeometry,
280 const ElementVolumeVariables& prevElemVolVars,
281 const ElementVolumeVariables& curElemVolVars,
282 const SubControlVolume& scv)
const
284 const auto& curVolVars = curElemVolVars[scv];
285 const auto& prevVolVars = prevElemVolVars[scv];
290 NumEqVector prevStorage = computeStorageImpl_(
problem, fvGeometry, scv, prevVolVars,
true);
291 NumEqVector storage = computeStorageImpl_(
problem, fvGeometry, scv, curVolVars,
false);
293 prevStorage *= prevVolVars.extrusionFactor();
294 storage *= curVolVars.extrusionFactor();
296 storage -= prevStorage;
300 residual[scv.localDofIndex()] += storage;
305 const FVElementGeometry& fvGeometry,
306 const SubControlVolume& scv,
307 const VolumeVariables& volVars,
308 [[maybe_unused]]
bool isPreviousTimeStep)
const
310 if constexpr (Detail::hasTimeInfoInterfaceCVFE<Implementation, Problem, FVElementGeometry, SubControlVolume, VolumeVariables>())
311 return this->
asImp().computeStorage(problem, fvGeometry, scv, volVars, isPreviousTimeStep);
313 return this->
asImp().computeStorage(problem, scv, volVars);
316 auto bcTypes_(
const Problem&
problem,
317 const FVElementGeometry& fvGeometry,
318 const SubControlVolumeFace& scvf,
319 const ElementBoundaryTypes& elemBcTypes)
const
323 if constexpr (Detail::hasProblemBoundaryTypesForIntersectionFunction<Problem, FVElementGeometry, typename GridView::Intersection>())
324 return elemBcTypes.get(fvGeometry, scvf);
326 return elemBcTypes.get(fvGeometry, fvGeometry.scv(scvf.insideScvIdx()));
An interpolation point related to a face of an element.
Definition: cvfe/interpolationpointdata.hh:99
The element-wise residual for control-volume finite element schemes.
Definition: cvfelocalresidual.hh:63
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: cvfelocalresidual.hh:104
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:213
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:276
ElementResidualVector evalFluxAndSource(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &bcTypes) const
Compute the flux and source.
Definition: cvfelocalresidual.hh:136
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:180
void addToElementStorageResidual(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars) const
add additional storage contributions (e.g. hybrid CVFE schemes)
Definition: cvfelocalresidual.hh:161
typename ParentType::ElementResidualVector ElementResidualVector
Definition: cvfelocalresidual.hh:88
void addToElementFluxAndSourceResidual(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &bcTypes) const
add additional flux and source contributions (e.g. hybrid CVFE schemes)
Definition: cvfelocalresidual.hh:170
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:36
Implementation & asImp()
Definition: fvlocalresidual.hh:489
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:474
const TimeLoop & timeLoop() const
Definition: fvlocalresidual.hh:479
bool isStationary() const
returns true if the residual is stationary
Definition: fvlocalresidual.hh:483
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:87
ReservedBlockVector< NumEqVector, Detail::LocalDofs::maxNumLocalDofs< FVElementGeometry >()> ElementResidualVector
the container storing all element residuals
Definition: fvlocalresidual.hh:58
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
Defines all properties used in Dumux.
Classes representing interpolation point data for control-volume finite element schemes.
Helper classes to compute the integration elements.
The element-wise residual for finite volume schemes.
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
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
Definition: cvfelocalresidual.hh:28
constexpr bool hasScvfIsOverlapping()
Definition: cvfelocalresidual.hh:47
constexpr bool hasTimeInfoInterfaceCVFE()
Definition: cvfelocalresidual.hh:38
decltype(std::declval< Imp >().computeStorage(std::declval< P >(), std::declval< G >(), std::declval< S >(), std::declval< V >(), true)) TimeInfoInterfaceCVFEDetector
Definition: cvfelocalresidual.hh:35
decltype(std::declval< Imp >().isOverlapping()) SCVFIsOverlappingDetector
Definition: cvfelocalresidual.hh:44
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
A helper to deduce a vector with the same size as numbers of equations.