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>
24#include <dumux/common/concepts/ipdata_.hh>
32template<
class Imp,
class P,
class G,
class S,
class V>
34 std::declval<Imp>().computeStorage(
35 std::declval<P>(), std::declval<G>(), std::declval<S>(), std::declval<V>(),
true
39template<
class Imp,
class P,
class G,
class S,
class V>
41{
return Dune::Std::is_detected<TimeInfoInterfaceCVFEDetector, Imp, P, G, S, V>::value; }
45 std::declval<Imp>().isOverlapping()
50{
return Dune::Std::is_detected<SCVFIsOverlappingDetector, Imp>::value; }
52template<
class Problem>
55template<
class Problem>
58 if constexpr (Dune::Std::is_detected<ProvidesIntegralInterfaceDetector, Problem>::value)
75template<
class TypeTag>
83 using GridView =
typename GridGeometry::GridView;
84 using Element =
typename GridView::template Codim<0>::Entity;
86 using FVElementGeometry =
typename GridGeometry::LocalView;
88 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
89 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
90 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
91 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
99 using ParentType::ParentType;
115 const FVElementGeometry& fvGeometry,
116 const ElementVolumeVariables& prevElemVolVars,
117 const ElementVolumeVariables& curElemVolVars)
const
119 assert(!this->
isStationary() &&
"no time loop set for storage term evaluation");
126 for (
const auto& scv :
scvs(fvGeometry))
127 this->
asImp().evalStorage(residual, this->
problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
130 this->
asImp().addToElementStorageResidual(residual, this->
problem(), element, fvGeometry, prevElemVolVars, curElemVolVars);
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& elemVolVars,
149 const ElementFluxVariablesCache& elemFluxVarsCache,
150 const ElementBoundaryTypes &bcTypes)
const
157 for (
const auto& scv :
scvs(fvGeometry))
158 this->
asImp().evalSource(residual, this->
problem(), element, fvGeometry, elemVolVars, scv);
161 for (
auto&& scvf : scvfs(fvGeometry))
162 this->
asImp().evalFlux(residual, this->
problem(), element, fvGeometry, elemVolVars, bcTypes, elemFluxVarsCache, scvf);
165 this->
asImp().addToElementFluxAndSourceResidual(residual, this->
problem(), element, fvGeometry, elemVolVars, elemFluxVarsCache, bcTypes);
173 const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& prevElemVolVars,
176 const ElementVolumeVariables& curElemVolVars)
const
182 const Element& element,
183 const FVElementGeometry& fvGeometry,
184 const ElementVolumeVariables& curElemVolVars,
185 const ElementFluxVariablesCache& elemFluxVarsCache,
186 const ElementBoundaryTypes &bcTypes)
const
192 const Element& element,
193 const FVElementGeometry& fvGeometry,
194 const ElementVolumeVariables& elemVolVars,
195 const ElementBoundaryTypes& elemBcTypes,
196 const ElementFluxVariablesCache& elemFluxVarsCache,
197 const SubControlVolumeFace& scvf)
const
199 const auto flux = this->
asImp().evalFlux(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
200 if (!scvf.boundary())
202 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
203 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
204 residual[insideScv.localDofIndex()] += flux;
207 if constexpr (Detail::hasScvfIsOverlapping<SubControlVolumeFace>())
209 if (!scvf.isOverlapping())
210 residual[outsideScv.localDofIndex()] -= flux;
213 residual[outsideScv.localDofIndex()] -= flux;
217 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
218 residual[insideScv.localDofIndex()] += flux;
224 const Element& element,
225 const FVElementGeometry& fvGeometry,
226 const ElementVolumeVariables& elemVolVars,
227 const ElementBoundaryTypes& elemBcTypes,
228 const ElementFluxVariablesCache& elemFluxVarsCache,
229 const SubControlVolumeFace& scvf)
const
231 NumEqVector flux(0.0);
233 if constexpr (Detail::providesIntegralInterface<Problem>())
236 flux += this->
asImp().fluxIntegral(fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
239 const auto& bcTypes = bcTypes_(
problem, fvGeometry, scvf, elemBcTypes);
244 if (bcTypes.hasNeumann())
246 NumEqVector boundaryFluxes =
problem.boundaryFluxIntegral(fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
249 for (
int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
250 if (bcTypes.isNeumann(eqIdx))
251 flux[eqIdx] += boundaryFluxes[eqIdx];
259 flux += this->
asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
263 const auto& bcTypes = bcTypes_(
problem, fvGeometry, scvf, elemBcTypes);
268 if (bcTypes.hasNeumann())
270 NumEqVector boundaryFluxes;
271 boundaryFluxes =
problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
273 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
274 boundaryFluxes *= Extrusion::area(fvGeometry, scvf) * elemVolVars[scv].extrusionFactor();
277 for (
int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
278 if (bcTypes.isNeumann(eqIdx))
279 flux[eqIdx] += boundaryFluxes[eqIdx];
304 const Element& element,
305 const FVElementGeometry& fvGeometry,
306 const ElementVolumeVariables& prevElemVolVars,
307 const ElementVolumeVariables& curElemVolVars,
308 const SubControlVolume& scv)
const
310 if constexpr (Detail::providesIntegralInterface<Problem>())
312 NumEqVector prevStorage = this->
asImp().storageIntegral(fvGeometry, prevElemVolVars, scv,
true);
313 NumEqVector storage = this->
asImp().storageIntegral(fvGeometry, curElemVolVars, scv,
false);
315 storage -= prevStorage;
318 residual[scv.localDofIndex()] += storage;
322 const auto& curVolVars = curElemVolVars[scv];
323 const auto& prevVolVars = prevElemVolVars[scv];
328 NumEqVector prevStorage = computeStorageImpl_(
problem, fvGeometry, scv, prevVolVars,
true);
329 NumEqVector storage = computeStorageImpl_(
problem, fvGeometry, scv, curVolVars,
false);
331 prevStorage *= prevVolVars.extrusionFactor();
332 storage *= curVolVars.extrusionFactor();
334 storage -= prevStorage;
338 residual[scv.localDofIndex()] += storage;
358 const Element& element,
359 const FVElementGeometry& fvGeometry,
360 const ElementVolumeVariables& curElemVolVars,
361 const SubControlVolume& scv)
const
363 NumEqVector source(0.0);
364 if constexpr (Detail::providesIntegralInterface<Problem>())
365 source = this->
asImp().sourceIntegral(fvGeometry, curElemVolVars, scv);
369 source = this->
asImp().computeSource(problem, element, fvGeometry, curElemVolVars, scv);
372 const auto& curVolVars = curElemVolVars[scv];
373 source *= curVolVars.extrusionFactor();
376 residual[scv.localDofIndex()] -= source;
381 const FVElementGeometry& fvGeometry,
382 const SubControlVolume& scv,
383 const VolumeVariables& volVars,
384 [[maybe_unused]]
bool isPreviousTimeStep)
const
386 if constexpr (Detail::hasTimeInfoInterfaceCVFE<Implementation, Problem, FVElementGeometry, SubControlVolume, VolumeVariables>())
387 return this->
asImp().computeStorage(problem, fvGeometry, scv, volVars, isPreviousTimeStep);
389 return this->
asImp().computeStorage(problem, scv, volVars);
392 auto bcTypes_(
const Problem&
problem,
393 const FVElementGeometry& fvGeometry,
394 const SubControlVolumeFace& scvf,
395 const ElementBoundaryTypes& elemBcTypes)
const
399 if constexpr (Detail::hasProblemBoundaryTypesForIntersectionFunction<Problem, FVElementGeometry, typename GridView::Intersection>())
400 return elemBcTypes.get(fvGeometry, scvf);
402 return elemBcTypes.get(fvGeometry, fvGeometry.scv(scvf.insideScvIdx()));
The element-wise residual for control-volume finite element schemes.
Definition: cvfelocalresidual.hh:77
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:114
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: cvfelocalresidual.hh:356
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:223
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:302
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:146
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:190
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:171
typename ParentType::ElementResidualVector ElementResidualVector
Definition: cvfelocalresidual.hh:98
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:180
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
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:319
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:30
decltype(Problem::providesIntegralInterface()) ProvidesIntegralInterfaceDetector
Definition: cvfelocalresidual.hh:53
constexpr bool hasScvfIsOverlapping()
Definition: cvfelocalresidual.hh:49
constexpr bool providesIntegralInterface()
Definition: cvfelocalresidual.hh:56
constexpr bool hasTimeInfoInterfaceCVFE()
Definition: cvfelocalresidual.hh:40
decltype(std::declval< Imp >().computeStorage(std::declval< P >(), std::declval< G >(), std::declval< S >(), std::declval< V >(), true)) TimeInfoInterfaceCVFEDetector
Definition: cvfelocalresidual.hh:37
decltype(std::declval< Imp >().isOverlapping()) SCVFIsOverlappingDetector
Definition: cvfelocalresidual.hh:46
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.
Quadrature rules over sub-control volumes and sub-control volume faces.