24#ifndef DUMUX_KEPSILON_PROBLEM_HH
25#define DUMUX_KEPSILON_PROBLEM_HH
47template<
class TypeTag>
55 using GridView =
typename GridGeometry::GridView;
56 using Element =
typename GridView::template Codim<0>::Entity;
59 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
60 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
61 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
62 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
65 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
74 using Indices =
typename ModelTraits::Indices;
76 static constexpr bool enableEnergyBalance = ModelTraits::enableEnergyBalance();
77 static constexpr bool isCompositional = ModelTraits::numFluidComponents() > 1;
80 static constexpr auto cellCenterOffset = ModelTraits::numEq() - CellCenterPrimaryVariables::dimension;
81 static_assert(cellCenterOffset == ModelTraits::dim(),
"cellCenterOffset must equal dim for staggered NavierStokes");
86 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
const std::string& paramGroup =
"")
87 : ParentType(gridGeometry, paramGroup)
95 if (!ParentType::isFlatWallBounded())
97 DUNE_THROW(Dune::NotImplemented,
"\n Due to grid/geometric concerns, k-epsilon models should only be used for "
98 <<
" wall bounded flows with flat channel geometries. "
99 <<
"\n If your geometry is a flat channel, please set the runtime parameter RANS.IsFlatWallBounded to true. \n");
102 ParentType::updateStaticWallProperties();
104 matchingPointIdx_.resize(this->gridGeometry().elementMapper().size(), 0);
105 storedDissipation_.resize(this->gridGeometry().elementMapper().size(), 0.0);
106 storedTurbulentKineticEnergy_.resize(this->gridGeometry().elementMapper().size(), 0.0);
107 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
108 zeroEqDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
118 ParentType::updateDynamicWallProperties(curSol);
121 for (
const auto& element : elements(this->gridGeometry().gridView()))
123 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
125 auto fvGeometry =
localView(this->gridGeometry());
126 fvGeometry.bindElement(element);
127 for (
auto&& scv : scvs(fvGeometry))
129 const int dofIdx = scv.dofIndex();
130 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
131 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
132 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
134 storedDissipation_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
135 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
137 VolumeVariables volVars;
138 volVars.update(elemSol, asImp_(), element, scv);
139 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity();
144 unsigned int numElementsInNearWallRegion = 0;
145 for (
const auto& element : elements(this->gridGeometry().gridView()))
147 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
148 unsigned int wallNormalAxis = asImp_().wallNormalAxis(elementIdx);
149 unsigned int neighborIndex0 = asImp_().neighborIndex(elementIdx, wallNormalAxis, 0);
150 unsigned int neighborIndex1 = asImp_().neighborIndex(elementIdx, wallNormalAxis, 1);
151 numElementsInNearWallRegion = inNearWallRegion(elementIdx)
152 ? numElementsInNearWallRegion + 1
153 : numElementsInNearWallRegion + 0;
154 if ((!inNearWallRegion(elementIdx) && (inNearWallRegion(neighborIndex0) || inNearWallRegion(neighborIndex1)))
155 || (!inNearWallRegion(elementIdx) && elementIdx == asImp_().wallElementIndex(elementIdx))
156 || (inNearWallRegion(elementIdx) && (asImp_().wallElementIndex(neighborIndex0) != asImp_().wallElementIndex(neighborIndex1))))
158 matchingPointIdx_[asImp_().wallElementIndex(elementIdx)] = elementIdx;
161 std::cout <<
"numElementsInNearWallRegion: " << numElementsInNearWallRegion << std::endl;
164 for (
const auto& element : elements(this->gridGeometry().gridView()))
166 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
167 zeroEqDynamicEddyViscosity_[elementIdx] = zeroEqEddyViscosityModel(elementIdx);
171 static const auto enableZeroEqScaling
172 = getParamFromGroup<bool>(this->paramGroup(),
"KEpsilon.EnableZeroEqScaling",
true);
173 if (enableZeroEqScaling)
175 for (
const auto& element : elements(this->gridGeometry().gridView()))
177 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
178 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
180 Scalar scalingFactor = storedDynamicEddyViscosity(matchingPointIndex)
181 / zeroEqDynamicEddyViscosity_[matchingPointIndex];
182 if (!isMatchingPoint(elementIdx)
183 && !std::isnan(scalingFactor) && !std::isinf(scalingFactor))
185 zeroEqDynamicEddyViscosity_[elementIdx] *= scalingFactor;
188 for (
const auto& element : elements(this->gridGeometry().gridView()))
190 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
191 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
192 if (isMatchingPoint(elementIdx))
194 zeroEqDynamicEddyViscosity_[matchingPointIndex] = storedDynamicEddyViscosity(matchingPointIndex);
205 unsigned int wallElementIdx = asImp_().wallElementIndex(elementIdx);
206 unsigned int matchingPointIndex = matchingPointIdx(wallElementIdx);
207 return (wallElementIdx == matchingPointIndex) ? yPlusNominal(elementIdx) < yPlusThreshold()
208 : yPlus(elementIdx) < yPlusThreshold();
215 {
return matchingPointIdx(asImp_().wallElementIndex(elementIdx)) == elementIdx; }
220 const Scalar
yPlus(
unsigned int elementIdx)
const
222 return asImp_().wallDistance(elementIdx) * uStar(elementIdx)
223 / asImp_().kinematicViscosity(elementIdx);
230 return asImp_().wallDistance(elementIdx) * uStarNominal(elementIdx)
231 / asImp_().kinematicViscosity(elementIdx);
244 Scalar yPlusValue = yPlus(elementIdx);
245 Scalar mixingLength = 0.0;
246 if (yPlusValue > 0.0)
248 mixingLength = asImp_().karmanConstant() * asImp_().wallDistance(elementIdx)
249 * (1.0 - exp(-yPlusValue / 26.0 ))
250 / sqrt(1.0 - exp(-0.26 * yPlusValue));
253 unsigned int wallNormalAxis = asImp_().wallNormalAxis(elementIdx);
254 unsigned int flowDirectionAxis = asImp_().flowDirectionAxis(elementIdx);
255 Scalar velocityGradient = asImp_().velocityGradient(elementIdx, flowDirectionAxis, wallNormalAxis);
256 return mixingLength * mixingLength * abs(velocityGradient) * asImp_().storedDensity(elementIdx);
260 const Scalar
uStar(
unsigned int elementIdx)
const
264 unsigned int wallElementIdx = asImp_().wallElementIndex(elementIdx);
265 unsigned int wallNormalAxis = asImp_().wallNormalAxis(elementIdx);
266 unsigned int flowDirectionAxis = asImp_().flowDirectionAxis(elementIdx);
267 return sqrt(asImp_().kinematicViscosity(wallElementIdx)
268 * abs(asImp_().velocityGradient(wallElementIdx, flowDirectionAxis, wallNormalAxis)));
276 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
277 return pow(cMu(), 0.25) * sqrt(storedTurbulentKineticEnergy(matchingPointIndex));
285 return uStarNominal(elementIdx) * uStarNominal(elementIdx) * uStarNominal(elementIdx)
286 / asImp_().karmanConstant() / asImp_().wallDistance(elementIdx);
294 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
295 return storedTurbulentKineticEnergy(matchingPointIndex);
302 Scalar velocityNominal = uStarNominal(elementIdx) * (1.0 / asImp_().karmanConstant() * log(yPlusNominal(elementIdx)) + 5.0);
303 return uStarNominal(elementIdx) * uStarNominal(elementIdx) * velocity / velocityNominal;
308 const SubControlVolumeFace& localSubFace,
309 const int& eqIdx)
const
311 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
312 auto bcTypes = asImp_().boundaryTypes(element, localSubFace);
313 return asImp_().isOnWall(localSubFace)
314 && bcTypes.isDirichlet(eqIdx)
315 && isMatchingPoint(elementIdx);
320 const FVElementGeometry& fvGeometry,
321 const ElementVolumeVariables& elemVolVars,
322 const ElementFaceVariables& elemFaceVars,
323 const SubControlVolumeFace& scvf,
324 const SubControlVolumeFace& localSubFace)
const
326 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
327 return FacePrimaryVariables(asImp_().tangentialMomentumWallFunction(elementIdx, elemFaceVars[scvf].velocitySelf())
328 * asImp_().storedDensity(elementIdx) );
332 template<
bool eB = enableEnergyBalance,
bool compositional = isCompositional,
333 typename std::enable_if_t<eB && compositional, int> = 0>
335 const FVElementGeometry& fvGeometry,
336 const ElementVolumeVariables& elemVolVars,
337 const ElementFaceVariables& elemFaceVars,
338 const SubControlVolumeFace& scvf)
const
340 return wallFunctionComponent(element, fvGeometry, elemVolVars, elemFaceVars, scvf)
341 + wallFunctionEnergy(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
345 template<
bool eB = enableEnergyBalance,
bool compositional = isCompositional,
346 typename std::enable_if_t<!eB && compositional, int> = 0>
348 const FVElementGeometry& fvGeometry,
349 const ElementVolumeVariables& elemVolVars,
350 const ElementFaceVariables& elemFaceVars,
351 const SubControlVolumeFace& scvf)
const
352 {
return wallFunctionComponent(element, fvGeometry, elemVolVars, elemFaceVars, scvf); }
355 template<
bool eB = enableEnergyBalance,
bool compositional = isCompositional,
356 typename std::enable_if_t<eB && !compositional, int> = 0>
358 const FVElementGeometry& fvGeometry,
359 const ElementVolumeVariables& elemVolVars,
360 const ElementFaceVariables& elemFaceVars,
361 const SubControlVolumeFace& scvf)
const
362 {
return wallFunctionEnergy(element, fvGeometry, elemVolVars, elemFaceVars, scvf); }
365 template<
bool eB = enableEnergyBalance,
bool compositional = isCompositional,
366 typename std::enable_if_t<!eB && !compositional, int> = 0>
368 const FVElementGeometry& fvGeometry,
369 const ElementVolumeVariables& elemVolVars,
370 const ElementFaceVariables& elemFaceVars,
371 const SubControlVolumeFace& scvf)
const
372 {
return CellCenterPrimaryVariables(0.0); }
376 const FVElementGeometry& fvGeometry,
377 const ElementVolumeVariables& elemVolVars,
378 const ElementFaceVariables& elemFaceVars,
379 const SubControlVolumeFace& scvf)
const
382 auto wallFunctionFlux = CellCenterPrimaryVariables(0.0);
383 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
386 for (
int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
388 if (ModelTraits::replaceCompEqIdx() == compIdx)
391 Scalar schmidtNumber = elemVolVars[scvf.insideScvIdx()].kinematicViscosity()
392 / elemVolVars[scvf.insideScvIdx()].diffusionCoefficient(0, 0, compIdx);
393 Scalar moleToMassConversionFactor = ModelTraits::useMoles()
394 ? 1.0 : FluidSystem::molarMass(compIdx);
395 wallFunctionFlux[compIdx] +=
396 -1.0 * (asImp_().dirichlet(element, scvf)[Indices::conti0EqIdx + compIdx]
397 - elemVolVars[scvf.insideScvIdx()].moleFraction(compIdx))
398 * elemVolVars[scvf.insideScvIdx()].molarDensity()
399 * moleToMassConversionFactor
400 * uStarNominal(elementIdx)
401 / asImp_().turbulentSchmidtNumber()
402 / (1. / asImp_().karmanConstant() * log(yPlusNominal(elementIdx) * 9.793)
403 + pFunction(schmidtNumber, asImp_().turbulentSchmidtNumber()));
406 if (ModelTraits::replaceCompEqIdx() < ModelTraits::numFluidComponents())
408 wallFunctionFlux[ModelTraits::replaceCompEqIdx()] =
409 -std::accumulate(wallFunctionFlux.begin(), wallFunctionFlux.end(), 0.0);
412 return wallFunctionFlux;
417 const FVElementGeometry& fvGeometry,
418 const ElementVolumeVariables& elemVolVars,
419 const ElementFaceVariables& elemFaceVars,
420 const SubControlVolumeFace& scvf)
const
423 auto wallFunctionFlux = CellCenterPrimaryVariables(0.0);
424 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
426 Scalar prandtlNumber = elemVolVars[scvf.insideScvIdx()].kinematicViscosity()
427 * elemVolVars[scvf.insideScvIdx()].density()
428 * elemVolVars[scvf.insideScvIdx()].heatCapacity()
429 / elemVolVars[scvf.insideScvIdx()].thermalConductivity();
430 wallFunctionFlux[Indices::energyEqIdx - cellCenterOffset] +=
431 -1.0 * (asImp_().dirichlet(element, scvf)[Indices::temperatureIdx]
432 - elemVolVars[scvf.insideScvIdx()].temperature())
433 * elemVolVars[scvf.insideScvIdx()].density()
434 * elemVolVars[scvf.insideScvIdx()].heatCapacity()
435 * uStarNominal(elementIdx)
436 / asImp_().turbulentPrandtlNumber()
437 / (1. / asImp_().karmanConstant() * log(yPlusNominal(elementIdx) * 9.793)
438 + pFunction(prandtlNumber, asImp_().turbulentPrandtlNumber()));
440 return wallFunctionFlux;
444 const Scalar
pFunction(Scalar molecularNumber, Scalar turbulentNumber)
const
449 * (pow(molecularNumber / turbulentNumber, 0.75) - 1.0)
450 * (1.0 + 0.28 * exp(-0.007 * molecularNumber / turbulentNumber));
459 static const Scalar yPlusThreshold = getParamFromGroup<Scalar>(this->paramGroup(),
"KEpsilon.YPlusThreshold", 30);
460 return yPlusThreshold;
465 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(),
"RANS.UseStoredEddyViscosity",
false);
466 return useStoredEddyViscosity;
470 {
return storedDissipation_[elementIdx]; }
473 {
return storedTurbulentKineticEnergy_[elementIdx]; }
476 {
return storedDynamicEddyViscosity_[elementIdx]; }
479 {
return zeroEqDynamicEddyViscosity_[elementIdx]; }
482 {
return matchingPointIdx_[elementIdx]; }
485 std::vector<unsigned int> matchingPointIdx_;
486 std::vector<Scalar> storedDissipation_;
487 std::vector<Scalar> storedTurbulentKineticEnergy_;
488 std::vector<Scalar> storedDynamicEddyViscosity_;
489 std::vector<Scalar> zeroEqDynamicEddyViscosity_;
492 Implementation &asImp_()
493 {
return *
static_cast<Implementation *
>(
this); }
496 const Implementation &asImp_()
const
497 {
return *
static_cast<const Implementation *
>(
this); }
Base class for all staggered fv problems.
Free function to get the local view of a grid cache object.
The available discretization methods in Dumux.
The available free flow turbulence models in Dumux.
TurbulenceModel
The available free flow turbulence models in Dumux.
Definition: turbulencemodel.hh:37
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
forward declare
Definition: freeflow/rans/problem.hh:43
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:59
const Scalar tangentialMomentumWallFunction(unsigned int elementIdx, Scalar velocity) const
Returns the nominal wall shear stress (accounts for poor approximation of viscous sublayer)
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:299
Scalar yPlusThreshold() const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:457
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:116
const Scalar yPlus(unsigned int elementIdx) const
Returns the value at an element center.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:220
bool useStoredEddyViscosity() const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:463
unsigned int matchingPointIdx(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:481
const Scalar zeroEqEddyViscosityModel(unsigned int elementIdx) const
Returns the kinematic eddy viscosity of a 0-Eq. model.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:237
const Scalar cMu() const
Returns the constant.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:454
bool inNearWallRegion(unsigned int elementIdx) const
Returns if an element is located in the near-wall region.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:203
CellCenterPrimaryVariables wallFunction(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const
Returns the flux for non-isothermal and compositional RANS models.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:334
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:472
Scalar storedDissipation(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:469
Scalar zeroEqDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:478
CellCenterPrimaryVariables wallFunctionEnergy(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const
Returns the energy wall-function flux.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:416
FacePrimaryVariables wallFunction(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, const SubControlVolumeFace &localSubFace) const
Returns an additional wall function momentum flux (only needed for RANS models)
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:319
const Scalar uStarNominal(unsigned int elementIdx) const
Returns the nominal wall shear stress velocity (accounts for poor approximation of viscous sublayer)
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:272
const Scalar uStar(unsigned int elementIdx) const
Returns the wall shear stress velocity.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:260
const Scalar dissipationWallFunction(unsigned int elementIdx) const
Returns the dissipation calculated from the wall function consideration.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:283
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
The constructor sets the gravity, if desired by the user.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:86
const Scalar yPlusNominal(unsigned int elementIdx) const
Returns the nominal value at an element center.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:228
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:475
const Scalar pFunction(Scalar molecularNumber, Scalar turbulentNumber) const
Returns the value of the P-function after Jayatilleke .
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:444
CellCenterPrimaryVariables wallFunctionComponent(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const
Returns the component wall-function flux.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:375
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:93
bool useWallFunction(const Element &element, const SubControlVolumeFace &localSubFace, const int &eqIdx) const
Checks whether a wall function should be used.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:307
bool isMatchingPoint(unsigned int elementIdx) const
Returns if an element is the matching point.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:214
const Scalar turbulentKineticEnergyWallFunction(unsigned int elementIdx) const
Returns the turbulentKineticEnergy calculated from the wall function consideration.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:292
Declares all properties used in Dumux.
Adaption of the fully implicit scheme to the tracer transport model.
The local element solution class for staggered methods.