25#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_BOUNDARY_FLUXHELPER_HH
26#define DUMUX_NAVIERSTOKES_MOMENTUM_BOUNDARY_FLUXHELPER_HH
55 template<
class Problem,
class FVElementGeometry,
class ElementVolumeVariables,
class ElementFluxVariablesCache,
class Scalar>
57 const FVElementGeometry& fvGeometry,
58 const typename FVElementGeometry::SubControlVolumeFace& scvf,
59 const ElementVolumeVariables& elemVolVars,
60 const ElementFluxVariablesCache& elemFluxVarsCache,
62 const bool zeroNormalVelocityGradient =
true)
66 using MomentumFlux =
typename Problem::MomentumFluxType;
67 constexpr std::size_t dim =
static_cast<std::size_t
>(FVElementGeometry::GridGeometry::GridView::dimension);
69 MomentumFlux::size() == dim,
70 "Expects momentum flux to have as many entries as dimensions."
73 const auto& element = fvGeometry.element();
76 std::decay_t<
decltype(problem.dirichlet(element, scvf))>::size() == dim,
77 "Expects problem.dirichlet to return an array with as many entries as dimensions."
82 MomentumFlux flux(0.0);
85 flux[scvf.normalAxis()] = (
pressure - problem.referencePressure(element, fvGeometry, scvf)) * scvf.directionSign();
88 if (problem.enableInertiaTerms())
90 const auto v = elemVolVars[scvf.insideScvIdx()].velocity();
91 flux[scvf.normalAxis()] += v*v * problem.density(element, fvGeometry, scvf) * scvf.directionSign();
96 else if (scvf.isLateral())
99 if (!zeroNormalVelocityGradient)
100 return MomentumFlux(0.0);
105 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
106 if (scv.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(scv)))
110 MomentumFlux flux(0.0);
111 const Scalar mu = problem.effectiveViscosity(element, fvGeometry, scvf);
116 * scvf.directionSign();
118 else if (scvf.boundary())
120 * scvf.directionSign();
123 if (problem.enableInertiaTerms())
125 const auto transportingVelocity = [&]()
127 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
128 const auto innerTransportingVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
130 static const bool useOldScheme = getParam<bool>(
"FreeFlow.UseOldTransportingVelocity",
true);
136 if (problem.boundaryTypes(element, scvf).isDirichlet(scvf.normalAxis()))
137 return problem.dirichlet(element, scvf)[scvf.normalAxis()];
140 return innerTransportingVelocity;
146 if (problem.boundaryTypes(element, scvf).isDirichlet(scvf.normalAxis()))
147 return 0.5*(problem.dirichlet(element, scvf)[scvf.normalAxis()] + innerTransportingVelocity);
150 if (orthogonalScvf.boundary())
152 if (problem.boundaryTypes(element, orthogonalScvf).isDirichlet(scvf.normalAxis()))
153 return 0.5*(problem.dirichlet(element, scvf)[scvf.normalAxis()] + innerTransportingVelocity);
155 return innerTransportingVelocity;
159 const auto insideVolume = fvGeometry.scv(orthogonalScvf.insideScvIdx()).
volume();
160 const auto outsideVolume = fvGeometry.scv(orthogonalScvf.outsideScvIdx()).
volume();
161 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
162 return (insideVolume*innerTransportingVelocity + outsideVolume*outerTransportingVelocity) / (insideVolume + outsideVolume);
168 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
169 const auto outerVelocity = elemVolVars[scvf.outsideScvIdx()].velocity();
171 const auto rho = problem.insideAndOutsideDensity(element, fvGeometry, scvf);
172 const bool selfIsUpstream = scvf.directionSign() ==
sign(transportingVelocity);
174 const auto insideMomentum = innerVelocity * rho.first;
175 const auto outsideMomentum = outerVelocity * rho.second;
177 static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
179 const auto transportedMomentum = selfIsUpstream ? (upwindWeight * insideMomentum + (1.0 - upwindWeight) * outsideMomentum)
180 : (upwindWeight * outsideMomentum + (1.0 - upwindWeight) * insideMomentum);
182 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
186 else if (scvf.boundary())
188 const auto insideDensity = problem.density(element, fvGeometry.scv(scvf.insideScvIdx()));
189 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
190 flux[scv.dofAxis()] += innerVelocity * transportingVelocity * insideDensity * scvf.directionSign();
197 DUNE_THROW(Dune::InvalidStateException,
"Face is neither lateral nor frontal.");
210 template<
class Problem,
class FVElementGeometry,
class ElementVolumeVariables,
class ElementFluxVariablesCache,
class TangentialVelocityGradient =
double>
212 const FVElementGeometry& fvGeometry,
213 const typename FVElementGeometry::SubControlVolumeFace& scvf,
214 const ElementVolumeVariables& elemVolVars,
215 const ElementFluxVariablesCache& elemFluxVarsCache,
216 const TangentialVelocityGradient& tangentialVelocityGradient = TangentialVelocityGradient(0.0))
219 using MomentumFlux =
typename Problem::MomentumFluxType;
220 constexpr std::size_t dim =
static_cast<std::size_t
>(FVElementGeometry::GridGeometry::GridView::dimension);
222 MomentumFlux::size() == dim,
223 "Expects momentum flux to have as many entries as dimensions."
226 MomentumFlux flux(0.0);
229 if (scvf.isFrontal())
232 using Scalar = std::decay_t<
decltype(flux[0])>;
233 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
234 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
235 const auto& orthogonalScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
237 if (scvf.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(orthogonalScv)))
254 const Scalar velI = elemVolVars[scvf.insideScvIdx()].velocity();
257 const Scalar mu = problem.effectiveViscosity(fvGeometry.element(), fvGeometry, scvf);
258 const Scalar
distance = (scv.dofPosition()- scvf.ipGlobal()).two_norm();
259 const Scalar velGradJI = [&]
261 if (elemVolVars.hasVolVars(orthogonalScvf.outsideScvIdx()))
264 return tangentialVelocityGradient;
267 const Scalar slipVelocity = problem.beaversJosephVelocity(fvGeometry, scvf, elemVolVars, velGradJI)[scv.dofAxis()];
268 const Scalar velGradIJ = (slipVelocity - velI) /
distance * scvf.directionSign();
270 flux[scv.dofAxis()] -= (mu * (velGradIJ + velGradJI))*scvf.directionSign();
273 if (problem.enableInertiaTerms())
276 const auto transportingVelocity = [&]()
278 const auto innerTransportingVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
280 if (!elemVolVars.hasVolVars(orthogonalScvf.outsideScvIdx()))
281 return innerTransportingVelocity;
283 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
286 if (orthogonalScvf.boundary())
287 return outerTransportingVelocity;
289 static const bool useOldScheme = getParam<bool>(
"FreeFlow.UseOldTransportingVelocity",
true);
291 return innerTransportingVelocity;
295 const auto insideVolume = fvGeometry.scv(orthogonalScvf.insideScvIdx()).volume();
296 const auto outsideVolume = fvGeometry.scv(orthogonalScvf.outsideScvIdx()).volume();
297 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
298 return (insideVolume*innerTransportingVelocity + outsideVolume*outerTransportingVelocity) / (insideVolume + outsideVolume);
304 const auto rho = problem.insideAndOutsideDensity(fvGeometry.element(), fvGeometry, scvf);
305 const auto transportedMomentum = slipVelocity * rho.second;
307 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
310 else if (scv.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(scv)))
332 const Scalar velJ = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
335 const Scalar mu = problem.effectiveViscosity(fvGeometry.element(), fvGeometry, scvf);
336 const Scalar
distance = (fvGeometry.scv(orthogonalScvf.insideScvIdx()).dofPosition()- scvf.ipGlobal()).two_norm();
338 const Scalar velGradIJ = [&]
340 if (elemVolVars.hasVolVars(scvf.outsideScvIdx()))
343 return tangentialVelocityGradient;
346 const Scalar slipVelocity = problem.beaversJosephVelocity(fvGeometry, orthogonalScvf, elemVolVars, velGradIJ)[scvf.normalAxis()];
347 const Scalar velGradJI = (slipVelocity - velJ) /
distance * orthogonalScvf.directionSign();
349 flux[scv.dofAxis()] -= (mu * (velGradIJ + velGradJI))*scvf.directionSign();
352 if (problem.enableInertiaTerms())
355 const auto transportingVelocity = slipVelocity;
358 if (scvf.boundary() && elemVolVars.hasVolVars(scvf.outsideScvIdx()))
360 flux[scv.dofAxis()] += problem.density(fvGeometry.element(), scv)
361 * elemVolVars[scvf.outsideScvIdx()].velocity() * transportingVelocity * scvf.directionSign();
365 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
366 const auto outerVelocity = [&]
368 if (!elemVolVars.hasVolVars(scvf.outsideScvIdx()))
369 return innerVelocity;
371 return elemVolVars[scvf.outsideScvIdx()].velocity();
374 const auto rho = problem.insideAndOutsideDensity(fvGeometry.element(), fvGeometry, scvf);
375 const bool selfIsUpstream = scvf.directionSign() ==
sign(transportingVelocity);
377 const auto insideMomentum = innerVelocity * rho.first;
378 const auto outsideMomentum = outerVelocity * rho.second;
380 static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
382 const auto transportedMomentum = selfIsUpstream ? (upwindWeight * insideMomentum + (1.0 - upwindWeight) * outsideMomentum)
383 : (upwindWeight * outsideMomentum + (1.0 - upwindWeight) * insideMomentum);
385 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
The available discretization methods in Dumux.
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:294
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:641
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
constexpr FCStaggered fcstaggered
Definition: method.hh:140
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
Struct containing flux helper functions to be used in the momentum problem's Neumann function.
Definition: fluxhelper.hh:40
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:211
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:56
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:234
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:193