12#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_CVFE_FE_LOCAL_RESIDUAL_HELPER_HH
13#define DUMUX_NAVIERSTOKES_MOMENTUM_CVFE_FE_LOCAL_RESIDUAL_HELPER_HH
15#include <dune/geometry/quadraturerules.hh>
16#include <dumux/common/typetraits/localdofs_.hh>
29template<
class Scalar,
class NumEqVector,
class LocalBasis,
class Extrusion>
32 using RangeType =
typename LocalBasis::Traits::RangeType;
45 template<
class Res
idualVector,
class Problem,
class FVElementGeometry,
class ElementVariables>
47 const Problem& problem,
48 const FVElementGeometry& fvGeometry,
49 const ElementVariables& prevElemVars,
50 const ElementVariables& curElemVars,
51 const Scalar timeStepSize)
53 if constexpr (Detail::LocalDofs::hasNonCVLocalDofsInterface<FVElementGeometry>())
56 if (nonCVLocalDofs(fvGeometry).empty())
59 const auto& localBasis = fvGeometry.feLocalBasis();
60 std::vector<RangeType> integralShapeFunctions(localBasis.size(), RangeType(0.0));
63 const auto& geometry = fvGeometry.elementGeometry();
64 const auto& element = fvGeometry.element();
65 using GlobalPosition =
typename FVElementGeometry::GridGeometry::GlobalCoordinate;
70 const auto& ipData = qpData.ipData();
72 FeIpData feIpData(geometry, ipData.local(), ipData.global(), localBasis);
75 for (
const auto& localDof : nonCVLocalDofs(fvGeometry))
76 integralShapeFunctions[localDof.index()] += qpData.weight() * feIpData.shapeValue(localDof.index());
79 for (
const auto& localDof : nonCVLocalDofs(fvGeometry))
81 const auto localDofIdx = localDof.index();
82 const auto& data = ipData(fvGeometry, localDof);
83 const auto curDensity = problem.density(element, fvGeometry, data,
false);
84 const auto prevDensity = problem.density(element, fvGeometry, data,
true);
85 const auto curVelocity = curElemVars[localDofIdx].velocity();
86 const auto prevVelocity = prevElemVars[localDofIdx].velocity();
87 auto timeDeriv = (curDensity*curVelocity - prevDensity*prevVelocity);
88 timeDeriv /= timeStepSize;
91 for (
int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
92 residual[localDofIdx][eqIdx] += integralShapeFunctions[localDofIdx]*timeDeriv[eqIdx];
107 template<
class ResidualVector,
class Problem,
class FVElementGeometry,
108 class ElementVariables,
class ElementFluxVariablesCache,
class ElementBoundaryTypes>
109 [[deprecated(
"This function is deprecated and will be removed after release 3.11. "
110 "Use addFluxAndSourceTerms(residual, problem, fvGeometry, elemVars, elemBcTypes) instead.")]]
112 const Problem& problem,
113 const FVElementGeometry& fvGeometry,
114 const ElementVariables& elemVars,
115 const ElementFluxVariablesCache& elemFluxVarsCache,
116 const ElementBoundaryTypes& elemBcTypes)
118 if constexpr (Detail::LocalDofs::hasNonCVLocalDofsInterface<FVElementGeometry>())
121 if (nonCVLocalDofs(fvGeometry).empty())
124 if (!problem.pointSourceMap().empty())
125 DUNE_THROW(Dune::NotImplemented,
"Point sources are not implemented for hybrid momentum schemes.");
127 static const bool enableUnsymmetrizedVelocityGradient
128 = getParamFromGroup<bool>(problem.paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
130 const auto& element = fvGeometry.element();
131 using FluxVariablesCache =
typename ElementFluxVariablesCache::FluxVariablesCache;
135 const auto& ipData = qpData.ipData();
137 const auto& cache = elemFluxVarsCache[ipData];
138 FluxFunctionContext context(problem, fvGeometry, elemVars, cache);
139 const auto& v = context.velocity();
140 const auto& gradV = context.gradVelocity();
143 const Scalar mu = problem.effectiveViscosity(element, fvGeometry, ipData);
145 const Scalar
density = problem.density(element, fvGeometry, ipData);
147 for (
const auto& localDof : nonCVLocalDofs(fvGeometry))
149 const auto localDofIdx = localDof.index();
152 if (problem.enableInertiaTerms())
153 fluxAndSourceTerm -=
density*(v*cache.gradN(localDofIdx))*v;
156 fluxAndSourceTerm += enableUnsymmetrizedVelocityGradient ?
157 mu*
mv(gradV, cache.gradN(localDofIdx))
161 fluxAndSourceTerm -= problem.pressure(element, fvGeometry, ipData) * cache.gradN(localDofIdx);
164 auto sourceAtIp = problem.source(fvGeometry, elemVars, ipData);
166 sourceAtIp +=
density * problem.gravity();
168 const auto& shapeValues = cache.shapeValues();
169 for (
int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
171 fluxAndSourceTerm[eqIdx] -= shapeValues[localDofIdx] * sourceAtIp[eqIdx];
172 residual[localDofIdx][eqIdx] += qpData.weight()*fluxAndSourceTerm[eqIdx];
177 if (elemBcTypes.hasNeumann())
178 addBoundaryFluxes(residual, problem, fvGeometry, elemVars, elemFluxVarsCache, elemBcTypes);
192 template<
class ResidualVector,
class Problem,
class FVElementGeometry,
193 class ElementVariables,
class ElementFluxVariablesCache,
class ElementBoundaryTypes>
194 [[deprecated(
"This function is deprecated and will be removed after release 3.11. "
195 "Use addBoundaryFluxes(residual, problem, fvGeometry, elemVars, elemBcTypes) instead.")]]
197 const Problem& problem,
198 const FVElementGeometry& fvGeometry,
199 const ElementVariables& elemVars,
200 const ElementFluxVariablesCache& elemFluxVarsCache,
201 const ElementBoundaryTypes& elemBcTypes)
203 ResidualVector flux(0.0);
205 const auto& element = fvGeometry.element();
206 for (
const auto& intersection : intersections(fvGeometry.gridGeometry().gridView(), element))
208 if (!intersection.boundary())
211 const auto& bcTypes = elemBcTypes.get(fvGeometry, intersection);
212 if (!bcTypes.hasNeumann())
215 problem.addBoundaryFluxIntegrals(flux, fvGeometry, elemVars, elemFluxVarsCache, intersection, bcTypes);
229 template<
class ResidualVector,
class Problem,
class FVElementGeometry,
230 class ElementVariables,
class ElementBoundaryTypes>
232 const Problem& problem,
233 const FVElementGeometry& fvGeometry,
234 const ElementVariables& elemVars,
235 const ElementBoundaryTypes& elemBcTypes)
237 if constexpr (Detail::LocalDofs::hasNonCVLocalDofsInterface<FVElementGeometry>())
240 if (nonCVLocalDofs(fvGeometry).empty())
243 if (!problem.pointSourceMap().empty())
244 DUNE_THROW(Dune::NotImplemented,
"Point sources are not implemented for hybrid momentum schemes.");
246 static const bool enableUnsymmetrizedVelocityGradient
247 = getParamFromGroup<bool>(problem.paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
249 const auto& element = fvGeometry.element();
250 using Cache =
typename ElementVariables::InterpolationPointData;
254 const auto& ipData = qpData.ipData();
256 const auto& ipCache = cache(elemVars, ipData);
257 FluxFunctionContext context(problem, fvGeometry, elemVars, ipCache);
258 const auto& v = context.velocity();
259 const auto& gradV = context.gradVelocity();
262 const Scalar mu = problem.effectiveViscosity(element, fvGeometry, ipData);
264 const Scalar
density = problem.density(element, fvGeometry, ipData);
266 for (
const auto& localDof : nonCVLocalDofs(fvGeometry))
268 const auto localDofIdx = localDof.index();
271 if (problem.enableInertiaTerms())
272 fluxAndSourceTerm -=
density*(v*ipCache.gradN(localDofIdx))*v;
275 fluxAndSourceTerm += enableUnsymmetrizedVelocityGradient ?
276 mu*
mv(gradV, ipCache.gradN(localDofIdx))
280 fluxAndSourceTerm -= problem.pressure(element, fvGeometry, ipData) * ipCache.gradN(localDofIdx);
283 auto sourceAtIp = problem.source(fvGeometry, elemVars, ipData);
285 sourceAtIp +=
density * problem.gravity();
287 const auto& shapeValues = ipCache.shapeValues();
288 for (
int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
290 fluxAndSourceTerm[eqIdx] -= shapeValues[localDofIdx] * sourceAtIp[eqIdx];
291 residual[localDofIdx][eqIdx] += qpData.weight()*fluxAndSourceTerm[eqIdx];
296 if (elemBcTypes.hasNeumann())
310 template<
class ResidualVector,
class Problem,
class FVElementGeometry,
311 class ElementVariables,
class ElementBoundaryTypes>
313 const Problem& problem,
314 const FVElementGeometry& fvGeometry,
315 const ElementVariables& elemVars,
316 const ElementBoundaryTypes& elemBcTypes)
318 ResidualVector flux(0.0);
320 const auto& element = fvGeometry.element();
321 for (
const auto& intersection : intersections(fvGeometry.gridGeometry().gridView(), element))
323 if (!intersection.boundary())
326 const auto& bcTypes = elemBcTypes.get(fvGeometry, intersection);
327 if (!bcTypes.hasNeumann())
330 problem.addBoundaryFluxIntegrals(flux, fvGeometry, elemVars, intersection, bcTypes);
Boundary flag to store e.g. in sub control volume faces.
Definition: fem/interpolationpointdata.hh:21
Helper class for evaluating FE-based local residuals.
Definition: felocalresidual.hh:31
static void addFluxAndSourceTerms(ResidualVector &residual, const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &elemBcTypes)
Add flux and source residual contribution for non-CV local dofs.
Definition: felocalresidual.hh:111
static void addBoundaryFluxes(ResidualVector &residual, const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &elemBcTypes)
Evaluate Neumann boundary contributions.
Definition: felocalresidual.hh:196
static void addFluxAndSourceTerms(ResidualVector &residual, const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementBoundaryTypes &elemBcTypes)
Add flux and source residual contribution for non-CV local dofs.
Definition: felocalresidual.hh:231
static void addBoundaryFluxes(ResidualVector &residual, const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementBoundaryTypes &elemBcTypes)
Evaluate Neumann boundary contributions.
Definition: felocalresidual.hh:312
static void addStorageTerms(ResidualVector &residual, const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &prevElemVars, const ElementVariables &curElemVars, const Scalar timeStepSize)
Add storage residual contribution for non-CV local dofs.
Definition: felocalresidual.hh:46
Context for interpolating data on interpolation points.
Definition: flux.hh:99
Helper classes to compute the integration elements.
Shape functions and gradients at an interpolation point.
Dune::DenseVector< V >::derived_type mv(const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V > &v)
Returns the result of the projection of a vector v with a Matrix M.
Definition: math.hh:829
Dune::FieldMatrix< Scalar, n, m > getTransposed(const Dune::FieldMatrix< Scalar, m, n > &M)
Transpose a FieldMatrix.
Definition: math.hh:712
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 quadratureRule(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolume &scv, QuadratureRules::MidpointQuadrature)
Midpoint quadrature for scv.
Definition: quadraturerules.hh:148
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53