24#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_LOCAL_RESIDUAL_HH
25#define DUMUX_NAVIERSTOKES_MOMENTUM_LOCAL_RESIDUAL_HH
32#include <dune/common/hybridutilities.hh>
38static constexpr bool isRotationalExtrusion =
false;
40template<
int radialAxis>
41static constexpr bool isRotationalExtrusion<RotationalExtrusion<radialAxis>> =
true;
48template<
class TypeTag>
57 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
58 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
59 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
61 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
62 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
68 using FVElementGeometry =
typename GridGeometry::LocalView;
69 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
70 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
71 using GridView =
typename GridGeometry::GridView;
72 using Element =
typename GridView::template Codim<0>::Entity;
82 static constexpr auto dim = GridView::dimension;
86 using ParentType::ParentType;
99 const SubControlVolume& scv,
100 const VolumeVariables& volVars,
101 const bool isPreviousStorage =
false)
const
103 const auto& element =
problem.gridGeometry().element(scv.elementIndex());
104 return problem.density(element, scv, isPreviousStorage) * volVars.velocity();
121 const Element& element,
122 const FVElementGeometry& fvGeometry,
123 const ElementVolumeVariables& elemVolVars,
124 const SubControlVolume& scv)
const
127 source +=
problem.gravity()[scv.dofAxis()] *
problem.density(element, scv);
132 if constexpr (dim == 2 && Impl::isRotationalExtrusion<Extrusion>)
134 if (scv.dofAxis() == Extrusion::radialAxis)
136 const auto r = scv.center()[scv.dofAxis()] - fvGeometry.gridGeometry().bBoxMin()[scv.dofAxis()];
137 const auto& scvf = (*scvfs(fvGeometry, scv).begin());
140 source -= -2.0*
problem.effectiveViscosity(element, fvGeometry, scvf) * elemVolVars[scv].velocity() / (r*r);
145 source +=
problem.pressure(element, fvGeometry, scvf)/r;
164 const Element& element,
165 const FVElementGeometry& fvGeometry,
166 const ElementVolumeVariables& elemVolVars,
167 const SubControlVolumeFace& scvf,
168 const ElementFluxVariablesCache& elemFluxVarsCache,
169 const ElementBoundaryTypes& elemBcTypes)
const
171 FluxVariables fluxVars(
problem, element, fvGeometry, scvf, elemVolVars, elemFluxVarsCache, elemBcTypes);
173 NumEqVector flux(0.0);
174 flux += fluxVars.advectiveMomentumFlux();
175 flux += fluxVars.diffusiveMomentumFlux();
176 flux += fluxVars.pressureContribution();
183 const Element& element,
184 const FVElementGeometry& fvGeometry,
185 const ElementVolumeVariables& elemVolVars,
186 const ElementBoundaryTypes& elemBcTypes,
187 const ElementFluxVariablesCache& elemFluxVarsCache,
188 const SubControlVolumeFace& scvf)
const
190 if (
const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); scv.boundary())
192 const auto& frontalScvfOnBoundary = fvGeometry.frontalScvfOnBoundary(scv);
193 const auto& bcTypes = elemBcTypes[frontalScvfOnBoundary.localIndex()];
194 if (bcTypes.isDirichlet(scv.dofAxis()))
200 if (elemBcTypes.hasNeumann())
201 return this->
asImp().maybeHandleNeumannBoundary(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
203 return this->
asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, elemBcTypes);
209 const Element& element,
210 const FVElementGeometry& fvGeometry,
211 const ElementVolumeVariables& elemVolVars,
212 const ElementBoundaryTypes& elemBcTypes,
213 const ElementFluxVariablesCache& elemFluxVarsCache,
214 const SubControlVolumeFace& scvf)
const
218 std::decay_t<
decltype(
219 problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf)
220 )>::size() == ModelTraits::dim(),
221 "The momentum model expects problem.neumann to return a vector of size dim. "
222 "When in doubt you should be able to use 'using NumEqVector = typename ParentType::NumEqVector;'."
224 assert(elemBcTypes.hasNeumann());
226 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
245 if (scvf.boundary() && scvf.isFrontal())
247 const auto& bcTypes = elemBcTypes[scvf.localIndex()];
248 if (bcTypes.hasNeumann() && bcTypes.isNeumann(scv.dofAxis()))
250 const auto neumannFluxes =
problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
251 return neumannFluxes[scv.dofAxis()] * Extrusion::area(scvf) * elemVolVars[scv].extrusionFactor();
254 else if (scvf.isLateral())
278 const auto& frontalScvfOnBoundary = fvGeometry.frontalScvfOnBoundary(scv);
279 assert(frontalScvfOnBoundary.isFrontal() && frontalScvfOnBoundary.boundary());
280 const auto& bcTypes = elemBcTypes[frontalScvfOnBoundary.localIndex()];
281 if (bcTypes.hasNeumann() && bcTypes.isNeumann(scvf.normalAxis()))
283 const auto neumannFluxes =
problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
284 return neumannFluxes[scv.dofAxis()] * Extrusion::area(scvf) * elemVolVars[scv].extrusionFactor();
288 else if (scvf.boundary())
309 assert(scvf.isLateral());
310 const auto& bcTypes = elemBcTypes[scvf.localIndex()];
311 if (bcTypes.hasNeumann() && bcTypes.isNeumann(scv.dofAxis()))
313 const auto neumannFluxes =
problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
314 return neumannFluxes[scv.dofAxis()] * Extrusion::area(scvf) * elemVolVars[scv].extrusionFactor();
319 return this->
asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, elemBcTypes);
324 {
return *
static_cast<Implementation *
>(
this); }
328 {
return *
static_cast<const Implementation *
>(
this); }
Calculates the element-wise residual for the box scheme.
A helper to deduce a vector with the same size as numbers of equations.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
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
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
constexpr FCStaggered fcstaggered
Definition: method.hh:142
The element-wise residual for the box scheme.
Definition: fclocalresidual.hh:46
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: fclocalresidual.hh:112
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:47
Implementation & asImp()
Definition: fvlocalresidual.hh:501
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:486
Element-wise calculation of the Navier-Stokes residual for models using the staggered discretization.
Definition: freeflow/navierstokes/momentum/localresidual.hh:51
NumEqVector maybeHandleNeumannBoundary(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Definition: freeflow/navierstokes/momentum/localresidual.hh:208
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: freeflow/navierstokes/momentum/localresidual.hh:120
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &elemBcTypes) const
Evaluates the mass flux over a face of a sub control volume.
Definition: freeflow/navierstokes/momentum/localresidual.hh:163
Implementation & asImp_()
Returns the implementation of the problem (i.e. static polymorphism)
Definition: freeflow/navierstokes/momentum/localresidual.hh:323
const Implementation & asImp_() const
Returns the implementation of the problem (i.e. static polymorphism)
Definition: freeflow/navierstokes/momentum/localresidual.hh:327
NumEqVector maybeHandleDirichletBoundary(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Definition: freeflow/navierstokes/momentum/localresidual.hh:182
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars, const bool isPreviousStorage=false) const
Calculate the source term of the equation.
Definition: freeflow/navierstokes/momentum/localresidual.hh:98
Declares all properties used in Dumux.