3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/rans/twoeq/kepsilon/problem.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_KEPSILON_PROBLEM_HH
25#define DUMUX_KEPSILON_PROBLEM_HH
26
27#include <numeric>
28
37
38#include "model.hh"
39
40namespace Dumux {
41
48template<class TypeTag>
49class RANSProblemImpl<TypeTag, TurbulenceModel::kepsilon> : public RANSProblemBase<TypeTag>
50{
52 using Implementation = GetPropType<TypeTag, Properties::Problem>;
54
56 using GridView = typename GridGeometry::GridView;
57 using Element = typename GridView::template Codim<0>::Entity;
58
60 using GridFaceVariables = typename GridVariables::GridFaceVariables;
61 using ElementFaceVariables = typename GridFaceVariables::LocalView;
62 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
63 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
64
65 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
66 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
67 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
68
72 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
76 using Indices = typename ModelTraits::Indices;
77
78 static constexpr bool enableEnergyBalance = ModelTraits::enableEnergyBalance();
79 static constexpr bool isCompositional = ModelTraits::numFluidComponents() > 1;
80
81 // account for the offset of the cell center privars within the PrimaryVariables container
82 static constexpr auto cellCenterOffset = ModelTraits::numEq() - CellCenterPrimaryVariables::dimension;
83 static_assert(cellCenterOffset == ModelTraits::dim(), "cellCenterOffset must equal dim for staggered NavierStokes");
84
85public:
86
88 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
89 : ParentType(gridGeometry, paramGroup)
90 { }
91
96 {
97 if (!ParentType::isFlatWallBounded())
98 {
99 DUNE_THROW(Dune::NotImplemented, "\n Due to grid/geometric concerns, k-epsilon models should only be used for "
100 << " wall bounded flows with flat channel geometries. "
101 << "\n If your geometry is a flat channel, please set the runtime parameter RANS.IsFlatWallBounded to true. \n");
102 }
103
104 ParentType::updateStaticWallProperties();
105 // update size and initial values of the global vectors
106 matchingPointIdx_.resize(this->gridGeometry().elementMapper().size(), 0);
107 storedDissipation_.resize(this->gridGeometry().elementMapper().size(), 0.0);
108 storedTurbulentKineticEnergy_.resize(this->gridGeometry().elementMapper().size(), 0.0);
109 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
110 zeroEqDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
111 }
112
118 void updateDynamicWallProperties(const SolutionVector& curSol)
119 {
120 ParentType::updateDynamicWallProperties(curSol);
121
122 // update the stored eddy viscosities
123 auto fvGeometry = localView(this->gridGeometry());
124 for (const auto& element : elements(this->gridGeometry().gridView()))
125 {
126 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
127 fvGeometry.bindElement(element);
128 for (auto&& scv : scvs(fvGeometry))
129 {
130 const int dofIdx = scv.dofIndex();
131 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
132 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
133 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
134 // NOTE: first update the turbulence quantities
135 storedDissipation_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
136 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
137 // NOTE: then update the volVars
138 VolumeVariables volVars;
139 volVars.update(elemSol, asImp_(), element, scv);
140 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity();
141 }
142 }
143
144 // get matching point for k-epsilon wall function
145 unsigned int numElementsInNearWallRegion = 0;
146 for (const auto& element : elements(this->gridGeometry().gridView()))
147 {
148 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
149 unsigned int wallNormalAxis = asImp_().wallNormalAxis(elementIdx);
150 unsigned int neighborIndex0 = asImp_().neighborIndex(elementIdx, wallNormalAxis, 0);
151 unsigned int neighborIndex1 = asImp_().neighborIndex(elementIdx, wallNormalAxis, 1);
152 numElementsInNearWallRegion = inNearWallRegion(elementIdx)
153 ? numElementsInNearWallRegion + 1
154 : numElementsInNearWallRegion + 0;
155 if ((!inNearWallRegion(elementIdx) && (inNearWallRegion(neighborIndex0) || inNearWallRegion(neighborIndex1)))
156 || (!inNearWallRegion(elementIdx) && elementIdx == asImp_().wallElementIndex(elementIdx))
157 || (inNearWallRegion(elementIdx) && (asImp_().wallElementIndex(neighborIndex0) != asImp_().wallElementIndex(neighborIndex1))))
158 {
159 matchingPointIdx_[asImp_().wallElementIndex(elementIdx)] = elementIdx;
160 }
161 }
162 std::cout << "numElementsInNearWallRegion: " << numElementsInNearWallRegion << std::endl;
163
164 // calculate the potential zeroeq eddy viscosities for two-layer model
165 for (const auto& element : elements(this->gridGeometry().gridView()))
166 {
167 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
168 zeroEqDynamicEddyViscosity_[elementIdx] = zeroEqEddyViscosityModel(elementIdx);
169 }
170
171 // then make them match at the matching point
172 static const auto enableZeroEqScaling
173 = getParamFromGroup<bool>(this->paramGroup(), "KEpsilon.EnableZeroEqScaling", true);
174 if (enableZeroEqScaling)
175 {
176 for (const auto& element : elements(this->gridGeometry().gridView()))
177 {
178 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
179 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
180
181 Scalar scalingFactor = storedDynamicEddyViscosity(matchingPointIndex)
182 / zeroEqDynamicEddyViscosity_[matchingPointIndex];
183 if (!isMatchingPoint(elementIdx)
184 && !std::isnan(scalingFactor) && !std::isinf(scalingFactor))
185 {
186 zeroEqDynamicEddyViscosity_[elementIdx] *= scalingFactor;
187 }
188 }
189 for (const auto& element : elements(this->gridGeometry().gridView()))
190 {
191 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
192 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
193 if (isMatchingPoint(elementIdx))
194 {
195 zeroEqDynamicEddyViscosity_[matchingPointIndex] = storedDynamicEddyViscosity(matchingPointIndex);
196 }
197 }
198 }
199 }
200
204 bool inNearWallRegion(unsigned int elementIdx) const
205 {
206 unsigned int wallElementIdx = asImp_().wallElementIndex(elementIdx);
207 unsigned int matchingPointIndex = matchingPointIdx(wallElementIdx);
208 return (wallElementIdx == matchingPointIndex) ? yPlusNominal(elementIdx) < yPlusThreshold()
209 : yPlus(elementIdx) < yPlusThreshold();
210 }
211
215 bool isMatchingPoint(unsigned int elementIdx) const
216 { return matchingPointIdx(asImp_().wallElementIndex(elementIdx)) == elementIdx; }
217
221 const Scalar yPlus(unsigned int elementIdx) const
222 {
223 return asImp_().wallDistance(elementIdx) * uStar(elementIdx)
224 / asImp_().kinematicViscosity(elementIdx);
225 }
229 const Scalar yPlusNominal(unsigned int elementIdx) const
230 {
231 return asImp_().wallDistance(elementIdx) * uStarNominal(elementIdx)
232 / asImp_().kinematicViscosity(elementIdx);
233 }
234
238 const Scalar zeroEqEddyViscosityModel(unsigned int elementIdx) const
239 {
240 using std::abs;
241 using std::exp;
242 using std::sqrt;
243
244 // use VanDriest's model
245 Scalar yPlusValue = yPlus(elementIdx);
246 Scalar mixingLength = 0.0;
247 if (yPlusValue > 0.0)
248 {
249 mixingLength = asImp_().karmanConstant() * asImp_().wallDistance(elementIdx)
250 * (1.0 - exp(-yPlusValue / 26.0 ))
251 / sqrt(1.0 - exp(-0.26 * yPlusValue));
252 }
253
254 unsigned int wallNormalAxis = asImp_().wallNormalAxis(elementIdx);
255 unsigned int flowDirectionAxis = asImp_().flowDirectionAxis(elementIdx);
256 Scalar velocityGradient = asImp_().velocityGradient(elementIdx, flowDirectionAxis, wallNormalAxis);
257 return mixingLength * mixingLength * abs(velocityGradient) * asImp_().storedDensity(elementIdx);
258 }
259
261 const Scalar uStar(unsigned int elementIdx) const
262 {
263 using std::abs;
264 using std::sqrt;
265 unsigned int wallElementIdx = asImp_().wallElementIndex(elementIdx);
266 unsigned int wallNormalAxis = asImp_().wallNormalAxis(elementIdx);
267 unsigned int flowDirectionAxis = asImp_().flowDirectionAxis(elementIdx);
268 return sqrt(asImp_().kinematicViscosity(wallElementIdx)
269 * abs(asImp_().velocityGradient(wallElementIdx, flowDirectionAxis, wallNormalAxis)));
270 }
271
273 const Scalar uStarNominal(unsigned int elementIdx) const
274 {
275 using std::pow;
276 using std::sqrt;
277 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
278 return pow(cMu(), 0.25) * sqrt(storedTurbulentKineticEnergy(matchingPointIndex));
279 }
280
284 const Scalar dissipationWallFunction(unsigned int elementIdx) const
285 {
286 return uStarNominal(elementIdx) * uStarNominal(elementIdx) * uStarNominal(elementIdx)
287 / asImp_().karmanConstant() / asImp_().wallDistance(elementIdx);
288 }
289
293 const Scalar turbulentKineticEnergyWallFunction(unsigned int elementIdx) const
294 {
295 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
296 return storedTurbulentKineticEnergy(matchingPointIndex);
297 }
298
300 const Scalar tangentialMomentumWallFunction(unsigned int elementIdx, Scalar velocity) const
301 {
302 using std::log;
303 Scalar velocityNominal = uStarNominal(elementIdx) * (1.0 / asImp_().karmanConstant() * log(yPlusNominal(elementIdx)) + 5.0);
304 return uStarNominal(elementIdx) * uStarNominal(elementIdx) * velocity / velocityNominal;
305 }
306
308 bool useWallFunction(const Element& element,
309 const SubControlVolumeFace& localSubFace,
310 const int& eqIdx) const
311 {
312 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
313 auto bcTypes = asImp_().boundaryTypes(element, localSubFace);
314
315 // Remove this check after release 3.5. IsOnWall Interface is deprecated
316 if constexpr (Deprecated::hasIsOnWall<Implementation, GlobalPosition>())
317 {
318 // Remove this part
319 return asImp_().isOnWall(localSubFace)
320 && bcTypes.isDirichlet(eqIdx)
321 && isMatchingPoint(elementIdx);
322 }
323 else
324 {
325 // Keep this part
326 return bcTypes.hasWall()
327 && bcTypes.isDirichlet(eqIdx)
328 && isMatchingPoint(elementIdx);
329 }
330 }
331
333 FacePrimaryVariables wallFunction(const Element& element,
334 const FVElementGeometry& fvGeometry,
335 const ElementVolumeVariables& elemVolVars,
336 const ElementFaceVariables& elemFaceVars,
337 const SubControlVolumeFace& scvf,
338 const SubControlVolumeFace& localSubFace) const
339 {
340 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
341 return FacePrimaryVariables(asImp_().tangentialMomentumWallFunction(elementIdx, elemFaceVars[scvf].velocitySelf())
342 * asImp_().storedDensity(elementIdx) );
343 }
344
346 template<bool eB = enableEnergyBalance, bool compositional = isCompositional,
347 typename std::enable_if_t<eB && compositional, int> = 0>
348 CellCenterPrimaryVariables wallFunction(const Element& element,
349 const FVElementGeometry& fvGeometry,
350 const ElementVolumeVariables& elemVolVars,
351 const ElementFaceVariables& elemFaceVars,
352 const SubControlVolumeFace& scvf) const
353 {
354 return wallFunctionComponent(element, fvGeometry, elemVolVars, elemFaceVars, scvf)
355 + wallFunctionEnergy(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
356 }
357
359 template<bool eB = enableEnergyBalance, bool compositional = isCompositional,
360 typename std::enable_if_t<!eB && compositional, int> = 0>
361 CellCenterPrimaryVariables wallFunction(const Element& element,
362 const FVElementGeometry& fvGeometry,
363 const ElementVolumeVariables& elemVolVars,
364 const ElementFaceVariables& elemFaceVars,
365 const SubControlVolumeFace& scvf) const
366 { return wallFunctionComponent(element, fvGeometry, elemVolVars, elemFaceVars, scvf); }
367
369 template<bool eB = enableEnergyBalance, bool compositional = isCompositional,
370 typename std::enable_if_t<eB && !compositional, int> = 0>
371 CellCenterPrimaryVariables wallFunction(const Element& element,
372 const FVElementGeometry& fvGeometry,
373 const ElementVolumeVariables& elemVolVars,
374 const ElementFaceVariables& elemFaceVars,
375 const SubControlVolumeFace& scvf) const
376 { return wallFunctionEnergy(element, fvGeometry, elemVolVars, elemFaceVars, scvf); }
377
379 template<bool eB = enableEnergyBalance, bool compositional = isCompositional,
380 typename std::enable_if_t<!eB && !compositional, int> = 0>
381 CellCenterPrimaryVariables wallFunction(const Element& element,
382 const FVElementGeometry& fvGeometry,
383 const ElementVolumeVariables& elemVolVars,
384 const ElementFaceVariables& elemFaceVars,
385 const SubControlVolumeFace& scvf) const
386 { return CellCenterPrimaryVariables(0.0); }
387
389 CellCenterPrimaryVariables wallFunctionComponent(const Element& element,
390 const FVElementGeometry& fvGeometry,
391 const ElementVolumeVariables& elemVolVars,
392 const ElementFaceVariables& elemFaceVars,
393 const SubControlVolumeFace& scvf) const
394 {
395 using std::log;
396 auto wallFunctionFlux = CellCenterPrimaryVariables(0.0);
397 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
398
399 // component mass fluxes
400 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
401 {
402 if (ModelTraits::replaceCompEqIdx() == compIdx)
403 continue;
404
405 Scalar schmidtNumber = elemVolVars[scvf.insideScvIdx()].kinematicViscosity()
406 / elemVolVars[scvf.insideScvIdx()].diffusionCoefficient(0, 0, compIdx);
407 Scalar moleToMassConversionFactor = ModelTraits::useMoles()
408 ? 1.0 : FluidSystem::molarMass(compIdx);
409 wallFunctionFlux[compIdx] +=
410 -1.0 * (asImp_().dirichlet(element, scvf)[Indices::conti0EqIdx + compIdx]
411 - elemVolVars[scvf.insideScvIdx()].moleFraction(compIdx))
412 * elemVolVars[scvf.insideScvIdx()].molarDensity()
413 * moleToMassConversionFactor
414 * uStarNominal(elementIdx)
415 / asImp_().turbulentSchmidtNumber()
416 / (1. / asImp_().karmanConstant() * log(yPlusNominal(elementIdx) * 9.793)
417 + pFunction(schmidtNumber, asImp_().turbulentSchmidtNumber()));
418 }
419
420 if (ModelTraits::replaceCompEqIdx() < ModelTraits::numFluidComponents())
421 {
422 wallFunctionFlux[ModelTraits::replaceCompEqIdx()] =
423 -std::accumulate(wallFunctionFlux.begin(), wallFunctionFlux.end(), 0.0);
424 }
425
426 return wallFunctionFlux;
427 }
428
430 CellCenterPrimaryVariables wallFunctionEnergy(const Element& element,
431 const FVElementGeometry& fvGeometry,
432 const ElementVolumeVariables& elemVolVars,
433 const ElementFaceVariables& elemFaceVars,
434 const SubControlVolumeFace& scvf) const
435 {
436 using std::log;
437 auto wallFunctionFlux = CellCenterPrimaryVariables(0.0);
438 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
439 // energy fluxes
440 Scalar prandtlNumber = elemVolVars[scvf.insideScvIdx()].kinematicViscosity()
441 * elemVolVars[scvf.insideScvIdx()].density()
442 * elemVolVars[scvf.insideScvIdx()].heatCapacity()
443 / elemVolVars[scvf.insideScvIdx()].thermalConductivity();
444 wallFunctionFlux[Indices::energyEqIdx - cellCenterOffset] +=
445 -1.0 * (asImp_().dirichlet(element, scvf)[Indices::temperatureIdx]
446 - elemVolVars[scvf.insideScvIdx()].temperature())
447 * elemVolVars[scvf.insideScvIdx()].density()
448 * elemVolVars[scvf.insideScvIdx()].heatCapacity()
449 * uStarNominal(elementIdx)
450 / asImp_().turbulentPrandtlNumber()
451 / (1. / asImp_().karmanConstant() * log(yPlusNominal(elementIdx) * 9.793)
452 + pFunction(prandtlNumber, asImp_().turbulentPrandtlNumber()));
453
454 return wallFunctionFlux;
455 }
456
458 const Scalar pFunction(Scalar molecularNumber, Scalar turbulentNumber) const
459 {
460 using std::pow;
461 using std::exp;
462 return 9.24
463 * (pow(molecularNumber / turbulentNumber, 0.75) - 1.0)
464 * (1.0 + 0.28 * exp(-0.007 * molecularNumber / turbulentNumber));
465 }
466
468 const Scalar cMu() const
469 { return 0.09; }
470
471 Scalar yPlusThreshold() const
472 {
473 static const Scalar yPlusThreshold = getParamFromGroup<Scalar>(this->paramGroup(), "KEpsilon.YPlusThreshold", 30.0);
474 return yPlusThreshold;
475 }
476
478 {
479 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(), "RANS.UseStoredEddyViscosity", false);
480 return useStoredEddyViscosity;
481 }
482
483 Scalar storedDissipation(const int elementIdx) const
484 { return storedDissipation_[elementIdx]; }
485
486 Scalar storedTurbulentKineticEnergy(const int elementIdx) const
487 { return storedTurbulentKineticEnergy_[elementIdx]; }
488
489 Scalar storedDynamicEddyViscosity(const int elementIdx) const
490 { return storedDynamicEddyViscosity_[elementIdx]; }
491
492 Scalar zeroEqDynamicEddyViscosity(const int elementIdx) const
493 { return zeroEqDynamicEddyViscosity_[elementIdx]; }
494
495 unsigned int matchingPointIdx(const int elementIdx) const
496 { return matchingPointIdx_[elementIdx]; }
497
498private:
499 std::vector<unsigned int> matchingPointIdx_;
500 std::vector<Scalar> storedDissipation_;
501 std::vector<Scalar> storedTurbulentKineticEnergy_;
502 std::vector<Scalar> storedDynamicEddyViscosity_;
503 std::vector<Scalar> zeroEqDynamicEddyViscosity_;
504
506 Implementation &asImp_()
507 { return *static_cast<Implementation *>(this); }
508
510 const Implementation &asImp_() const
511 { return *static_cast<const Implementation *>(this); }
512};
513
514} // end namespace Dumux
515
516#endif
Helpers for deprecation.
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:38
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Definition: freeflow/navierstokes/problem.hh:55
forward declare
Definition: freeflow/rans/problem.hh:44
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:60
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:300
Scalar yPlusThreshold() const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:471
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:118
const Scalar yPlus(unsigned int elementIdx) const
Returns the value at an element center.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:221
bool useStoredEddyViscosity() const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:477
unsigned int matchingPointIdx(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:495
const Scalar zeroEqEddyViscosityModel(unsigned int elementIdx) const
Returns the kinematic eddy viscosity of a 0-Eq. model.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:238
const Scalar cMu() const
Returns the constant.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:468
bool inNearWallRegion(unsigned int elementIdx) const
Returns if an element is located in the near-wall region.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:204
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:348
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:486
Scalar storedDissipation(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:483
Scalar zeroEqDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:492
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:430
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:333
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:273
const Scalar uStar(unsigned int elementIdx) const
Returns the wall shear stress velocity.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:261
const Scalar dissipationWallFunction(unsigned int elementIdx) const
Returns the dissipation calculated from the wall function consideration.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:284
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
The constructor sets the gravity, if desired by the user.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:88
const Scalar yPlusNominal(unsigned int elementIdx) const
Returns the nominal value at an element center.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:229
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:489
const Scalar pFunction(Scalar molecularNumber, Scalar turbulentNumber) const
Returns the value of the P-function after Jayatilleke .
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:458
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:389
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:95
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:308
bool isMatchingPoint(unsigned int elementIdx) const
Returns if an element is the matching point.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:215
const Scalar turbulentKineticEnergyWallFunction(unsigned int elementIdx) const
Returns the turbulentKineticEnergy calculated from the wall function consideration.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:293
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.