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 const Problem& problem,
108 const FVElementGeometry& fvGeometry,
109 const ElementVariables& elemVars)
111 if constexpr (Detail::LocalDofs::hasNonCVLocalDofsInterface<FVElementGeometry>())
114 if (nonCVLocalDofs(fvGeometry).empty())
117 if constexpr (
requires { problem.pointSources(); })
119 if (!problem.pointSources().empty())
120 DUNE_THROW(Dune::NotImplemented,
"Point sources are not implemented for hybrid momentum schemes.");
123 static const bool enableUnsymmetrizedVelocityGradient
126 const auto& element = fvGeometry.element();
127 using Cache =
typename ElementVariables::InterpolationPointData;
131 const auto& ipData = qpData.ipData();
133 const auto& ipCache = cache(elemVars, ipData);
134 FluxFunctionContext context(problem, fvGeometry, elemVars, ipCache);
135 const auto& v = context.velocity();
136 const auto& gradV = context.gradVelocity();
139 const Scalar mu = problem.effectiveViscosity(element, fvGeometry, ipData);
141 const Scalar density = problem.density(element, fvGeometry, ipData);
143 for (
const auto& localDof : nonCVLocalDofs(fvGeometry))
145 const auto localDofIdx = localDof.index();
148 if (problem.enableInertiaTerms())
149 fluxAndSourceTerm -= density*(v*ipCache.gradN(localDofIdx))*v;
152 fluxAndSourceTerm += enableUnsymmetrizedVelocityGradient ?
153 mu*
mv(gradV, ipCache.gradN(localDofIdx))
157 fluxAndSourceTerm -= problem.pressure(element, fvGeometry, ipData) * ipCache.gradN(localDofIdx);
160 auto sourceAtIp = problem.source(fvGeometry, elemVars, ipData);
162 sourceAtIp += density * problem.gravity();
164 const auto& shapeValues = ipCache.shapeValues();
165 for (
int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
167 fluxAndSourceTerm[eqIdx] -= shapeValues[localDofIdx] * sourceAtIp[eqIdx];
168 residual[localDofIdx][eqIdx] += qpData.weight()*fluxAndSourceTerm[eqIdx];