13#ifndef DUMUX_FREEFLOW_NAVIERSTOKES_MOMENTUM_BOUNDARY_FLUXHELPER_HH
14#define DUMUX_FREEFLOW_NAVIERSTOKES_MOMENTUM_BOUNDARY_FLUXHELPER_HH
35template<
class DiscretizationMethod,
class SlipVelocityPolicy =
void>
48template<
class SlipVelocityPolicy>
65 template<
class Problem,
class FVElementGeometry,
class ElementVolumeVariables,
class ElementFluxVariablesCache,
class Scalar>
67 const FVElementGeometry& fvGeometry,
68 const typename FVElementGeometry::SubControlVolumeFace& scvf,
69 const ElementVolumeVariables& elemVolVars,
70 const ElementFluxVariablesCache& elemFluxVarsCache,
72 const bool zeroNormalVelocityGradient =
true)
75 using MomentumFlux =
typename Problem::MomentumFluxType;
76 constexpr std::size_t dim =
static_cast<std::size_t
>(FVElementGeometry::GridGeometry::GridView::dimension);
78 MomentumFlux::size() == dim,
79 "Expects momentum flux to have as many entries as dimensions."
82 const auto& element = fvGeometry.element();
85 std::decay_t<
decltype(problem.dirichlet(element, scvf))>::size() == dim,
86 "Expects problem.dirichlet to return an array with as many entries as dimensions."
91 MomentumFlux flux(0.0);
94 flux[scvf.normalAxis()] = (
pressure - problem.referencePressure(element, fvGeometry, scvf)) * scvf.directionSign();
97 if (problem.enableInertiaTerms())
99 const auto v = elemVolVars[scvf.insideScvIdx()].velocity();
100 flux[scvf.normalAxis()] += v*v * problem.density(element, fvGeometry, scvf) * scvf.directionSign();
105 else if (scvf.isLateral())
108 if (!zeroNormalVelocityGradient)
109 return MomentumFlux(0.0);
114 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
115 if (scv.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(scv)))
116 return slipVelocityMomentumFlux(problem, fvGeometry, scvf, elemVolVars, elemFluxVarsCache);
119 MomentumFlux flux(0.0);
120 const Scalar mu = problem.effectiveViscosity(element, fvGeometry, scvf);
125 * scvf.directionSign();
127 else if (scvf.boundary())
129 * scvf.directionSign();
132 if (problem.enableInertiaTerms())
134 const auto transportingVelocity = [&]()
136 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
137 const auto innerTransportingVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
139 static const bool useOldScheme = getParam<bool>(
"FreeFlow.UseOldTransportingVelocity",
true);
145 if (problem.boundaryTypes(element, scvf).isDirichlet(scvf.normalAxis()))
146 return problem.dirichlet(element, scvf)[scvf.normalAxis()];
149 return innerTransportingVelocity;
155 if (problem.boundaryTypes(element, scvf).isDirichlet(scvf.normalAxis()))
156 return 0.5*(problem.dirichlet(element, scvf)[scvf.normalAxis()] + innerTransportingVelocity);
159 if (orthogonalScvf.boundary())
161 if (problem.boundaryTypes(element, orthogonalScvf).isDirichlet(scvf.normalAxis()))
162 return 0.5*(problem.dirichlet(element, scvf)[scvf.normalAxis()] + innerTransportingVelocity);
164 return innerTransportingVelocity;
168 const auto insideVolume = fvGeometry.scv(orthogonalScvf.insideScvIdx()).
volume();
169 const auto outsideVolume = fvGeometry.scv(orthogonalScvf.outsideScvIdx()).
volume();
170 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
171 return (insideVolume*innerTransportingVelocity + outsideVolume*outerTransportingVelocity) / (insideVolume + outsideVolume);
177 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
178 const auto outerVelocity = elemVolVars[scvf.outsideScvIdx()].velocity();
180 const auto rho = problem.insideAndOutsideDensity(element, fvGeometry, scvf);
181 const bool selfIsUpstream = scvf.directionSign() ==
sign(transportingVelocity);
183 const auto insideMomentum = innerVelocity * rho.first;
184 const auto outsideMomentum = outerVelocity * rho.second;
186 static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
188 const auto transportedMomentum = selfIsUpstream ? (upwindWeight * insideMomentum + (1.0 - upwindWeight) * outsideMomentum)
189 : (upwindWeight * outsideMomentum + (1.0 - upwindWeight) * insideMomentum);
191 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
195 else if (scvf.boundary())
197 const auto insideDensity = problem.density(element, fvGeometry.scv(scvf.insideScvIdx()));
198 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
199 flux[scv.dofAxis()] += innerVelocity * transportingVelocity * insideDensity * scvf.directionSign();
206 DUNE_THROW(Dune::InvalidStateException,
"Face is neither lateral nor frontal.");
219 template<
class Problem,
class FVElementGeometry,
class ElementVolumeVariables,
class ElementFluxVariablesCache,
class TangentialVelocityGradient =
double>
221 const FVElementGeometry& fvGeometry,
222 const typename FVElementGeometry::SubControlVolumeFace& scvf,
223 const ElementVolumeVariables& elemVolVars,
224 const ElementFluxVariablesCache& elemFluxVarsCache,
225 const TangentialVelocityGradient& tangentialVelocityGradient = TangentialVelocityGradient(0.0))
227 using MomentumFlux =
typename Problem::MomentumFluxType;
228 constexpr std::size_t dim =
static_cast<std::size_t
>(FVElementGeometry::GridGeometry::GridView::dimension);
230 MomentumFlux::size() == dim,
231 "Expects momentum flux to have as many entries as dimensions."
234 MomentumFlux flux(0.0);
237 if (scvf.isFrontal())
240 using Scalar = std::decay_t<
decltype(flux[0])>;
241 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
242 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
243 const auto& orthogonalScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
245 if (scvf.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(orthogonalScv)))
263 if constexpr (!std::is_void_v<SlipVelocityPolicy>)
265 const Scalar velI = elemVolVars[scvf.insideScvIdx()].velocity();
268 const Scalar mu = problem.effectiveViscosity(fvGeometry.element(), fvGeometry, scvf);
269 const Scalar
distance = (scv.dofPosition()- scvf.ipGlobal()).two_norm();
270 const Scalar velGradJI = [&]
272 if (elemVolVars.hasVolVars(orthogonalScvf.outsideScvIdx()))
275 return tangentialVelocityGradient;
278 const Scalar slipVelocity = SlipVelocityPolicy::velocity(problem, fvGeometry, scvf, elemVolVars, velGradJI)[scv.dofAxis()];
279 const Scalar velGradIJ = (slipVelocity - velI) /
distance * scvf.directionSign();
281 flux[scv.dofAxis()] -= (mu * (velGradIJ + velGradJI))*scvf.directionSign();
284 if (problem.enableInertiaTerms())
287 const auto transportingVelocity = [&]()
289 const auto innerTransportingVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
291 if (!elemVolVars.hasVolVars(orthogonalScvf.outsideScvIdx()))
292 return innerTransportingVelocity;
294 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
297 if (orthogonalScvf.boundary())
298 return outerTransportingVelocity;
300 static const bool useOldScheme = getParam<bool>(
"FreeFlow.UseOldTransportingVelocity",
true);
302 return innerTransportingVelocity;
306 const auto insideVolume = fvGeometry.scv(orthogonalScvf.insideScvIdx()).volume();
307 const auto outsideVolume = fvGeometry.scv(orthogonalScvf.outsideScvIdx()).volume();
308 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
309 return (insideVolume*innerTransportingVelocity + outsideVolume*outerTransportingVelocity) / (insideVolume + outsideVolume);
315 const auto rho = problem.insideAndOutsideDensity(fvGeometry.element(), fvGeometry, scvf);
316 const auto transportedMomentum = slipVelocity * rho.second;
318 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
322 DUNE_THROW(Dune::InvalidStateException,
"SlipVelocityPolicy needs to be specified as template argument");
324 else if (scv.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(scv)))
347 if constexpr (!std::is_void_v<SlipVelocityPolicy>)
349 const Scalar velJ = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
352 const Scalar mu = problem.effectiveViscosity(fvGeometry.element(), fvGeometry, scvf);
353 const Scalar
distance = (fvGeometry.scv(orthogonalScvf.insideScvIdx()).dofPosition()- scvf.ipGlobal()).two_norm();
355 const Scalar velGradIJ = [&]
357 if (elemVolVars.hasVolVars(scvf.outsideScvIdx()))
360 return tangentialVelocityGradient;
363 const Scalar slipVelocity = SlipVelocityPolicy::velocity(problem, fvGeometry, orthogonalScvf, elemVolVars, velGradIJ)[scvf.normalAxis()];
364 const Scalar velGradJI = (slipVelocity - velJ) /
distance * orthogonalScvf.directionSign();
366 flux[scv.dofAxis()] -= (mu * (velGradIJ + velGradJI))*scvf.directionSign();
369 if (problem.enableInertiaTerms())
372 const auto transportingVelocity = slipVelocity;
375 if (scvf.boundary() && elemVolVars.hasVolVars(scvf.outsideScvIdx()))
377 flux[scv.dofAxis()] += problem.density(fvGeometry.element(), scv)
378 * elemVolVars[scvf.outsideScvIdx()].velocity() * transportingVelocity * scvf.directionSign();
382 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
383 const auto outerVelocity = [&]
385 if (!elemVolVars.hasVolVars(scvf.outsideScvIdx()))
386 return innerVelocity;
388 return elemVolVars[scvf.outsideScvIdx()].velocity();
391 const auto rho = problem.insideAndOutsideDensity(fvGeometry.element(), fvGeometry, scvf);
392 const bool selfIsUpstream = scvf.directionSign() ==
sign(transportingVelocity);
394 const auto insideMomentum = innerVelocity * rho.first;
395 const auto outsideMomentum = outerVelocity * rho.second;
397 static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
399 const auto transportedMomentum = selfIsUpstream ? (upwindWeight * insideMomentum + (1.0 - upwindWeight) * outsideMomentum)
400 : (upwindWeight * outsideMomentum + (1.0 - upwindWeight) * insideMomentum);
402 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
406 DUNE_THROW(Dune::InvalidStateException,
"SlipVelocityPolicy needs to be specified as template argument");
414 [[deprecated(
"Replace with implementation class `NavierStokesMomentumBoundaryFlux`with template arguments `DiscretizationMethod` and `SlipVelocityPolicy`. This will be removed after 3.9.")]]
Navier Stokes slip velocity policy.
Definition: slipcondition.hh:59
static auto velocityGradJI(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars)
Returns the velocity gradient in line with our current scvf.
Definition: momentum/velocitygradients.hh:222
static auto velocityGradIJ(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars)
Returns the velocity gradient perpendicular to the orientation of our current scvf.
Definition: momentum/velocitygradients.hh:181
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:658
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:282
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
Define some often used mathematical functions.
The available discretization methods in Dumux.
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:22
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Navier Stokes slip condition.
Definition: method.hh:130
static auto fixedPressureMomentumFlux(const Problem &problem, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const Scalar pressure, const bool zeroNormalVelocityGradient=true)
Returns the momentum flux a fixed-pressure boundary.
Definition: fluxhelper.hh:66
static auto slipVelocityMomentumFlux(const Problem &problem, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const TangentialVelocityGradient &tangentialVelocityGradient=TangentialVelocityGradient(0.0))
Returns the momentum flux for a boundary with a slip condition.
Definition: fluxhelper.hh:220
Class to compute the boundary flux for the momentum balance of the Navier-Stokes model This helper cl...
Definition: fluxhelper.hh:36