3.6-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
36
37#include "model.hh"
38
39namespace Dumux {
40
47template<class TypeTag>
48class RANSProblemImpl<TypeTag, TurbulenceModel::kepsilon> : public RANSProblemBase<TypeTag>
49{
51 using Implementation = GetPropType<TypeTag, Properties::Problem>;
53
55 using GridView = typename GridGeometry::GridView;
56 using Element = typename GridView::template Codim<0>::Entity;
57
59 using GridFaceVariables = typename GridVariables::GridFaceVariables;
60 using ElementFaceVariables = typename GridFaceVariables::LocalView;
61 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
62 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
63
64 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
65 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
66 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
67
70 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
74 using Indices = typename ModelTraits::Indices;
75
76 static constexpr bool enableEnergyBalance = ModelTraits::enableEnergyBalance();
77 static constexpr bool isCompositional = ModelTraits::numFluidComponents() > 1;
78
79 // account for the offset of the cell center privars within the PrimaryVariables container
80 static constexpr auto cellCenterOffset = ModelTraits::numEq() - CellCenterPrimaryVariables::dimension;
81 static_assert(cellCenterOffset == ModelTraits::dim(), "cellCenterOffset must equal dim for staggered NavierStokes");
82
83public:
84
86 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
87 : ParentType(gridGeometry, paramGroup)
88 { }
89
94 {
95 if (!ParentType::isFlatWallBounded())
96 {
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");
100 }
101
102 ParentType::updateStaticWallProperties();
103 // update size and initial values of the global vectors
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);
109 }
110
116 template<class SolutionVector>
117 void updateDynamicWallProperties(const SolutionVector& curSol)
118 {
119 ParentType::updateDynamicWallProperties(curSol);
120
121 // update the stored eddy viscosities
122 auto fvGeometry = localView(this->gridGeometry());
123 for (const auto& element : elements(this->gridGeometry().gridView()))
124 {
125 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
126 fvGeometry.bindElement(element);
127 for (auto&& scv : scvs(fvGeometry))
128 {
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));
133 // NOTE: first update the turbulence quantities
134 storedDissipation_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
135 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
136 // NOTE: then update the volVars
137 VolumeVariables volVars;
138 volVars.update(elemSol, asImp_(), element, scv);
139 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity();
140 }
141 }
142
143 // get matching point for k-epsilon wall function
144 unsigned int numElementsInNearWallRegion = 0;
145 for (const auto& element : elements(this->gridGeometry().gridView()))
146 {
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))))
157 {
158 matchingPointIdx_[asImp_().wallElementIndex(elementIdx)] = elementIdx;
159 }
160 }
161 std::cout << "numElementsInNearWallRegion: " << numElementsInNearWallRegion << std::endl;
162
163 // calculate the potential zeroeq eddy viscosities for two-layer model
164 for (const auto& element : elements(this->gridGeometry().gridView()))
165 {
166 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
167 zeroEqDynamicEddyViscosity_[elementIdx] = zeroEqEddyViscosityModel(elementIdx);
168 }
169
170 // then make them match at the matching point
171 static const auto enableZeroEqScaling
172 = getParamFromGroup<bool>(this->paramGroup(), "KEpsilon.EnableZeroEqScaling", true);
173 if (enableZeroEqScaling)
174 {
175 for (const auto& element : elements(this->gridGeometry().gridView()))
176 {
177 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
178 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
179
180 Scalar scalingFactor = storedDynamicEddyViscosity(matchingPointIndex)
181 / zeroEqDynamicEddyViscosity_[matchingPointIndex];
182 if (!isMatchingPoint(elementIdx)
183 && !std::isnan(scalingFactor) && !std::isinf(scalingFactor))
184 {
185 zeroEqDynamicEddyViscosity_[elementIdx] *= scalingFactor;
186 }
187 }
188 for (const auto& element : elements(this->gridGeometry().gridView()))
189 {
190 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
191 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
192 if (isMatchingPoint(elementIdx))
193 {
194 zeroEqDynamicEddyViscosity_[matchingPointIndex] = storedDynamicEddyViscosity(matchingPointIndex);
195 }
196 }
197 }
198 }
199
203 bool inNearWallRegion(unsigned int elementIdx) const
204 {
205 unsigned int wallElementIdx = asImp_().wallElementIndex(elementIdx);
206 unsigned int matchingPointIndex = matchingPointIdx(wallElementIdx);
207 return (wallElementIdx == matchingPointIndex) ? yPlusNominal(elementIdx) < yPlusThreshold()
208 : yPlus(elementIdx) < yPlusThreshold();
209 }
210
214 bool isMatchingPoint(unsigned int elementIdx) const
215 { return matchingPointIdx(asImp_().wallElementIndex(elementIdx)) == elementIdx; }
216
220 const Scalar yPlus(unsigned int elementIdx) const
221 {
222 return asImp_().wallDistance(elementIdx) * uStar(elementIdx)
223 / asImp_().kinematicViscosity(elementIdx);
224 }
228 const Scalar yPlusNominal(unsigned int elementIdx) const
229 {
230 return asImp_().wallDistance(elementIdx) * uStarNominal(elementIdx)
231 / asImp_().kinematicViscosity(elementIdx);
232 }
233
237 const Scalar zeroEqEddyViscosityModel(unsigned int elementIdx) const
238 {
239 using std::abs;
240 using std::exp;
241 using std::sqrt;
242
243 // use VanDriest's model
244 Scalar yPlusValue = yPlus(elementIdx);
245 Scalar mixingLength = 0.0;
246 if (yPlusValue > 0.0)
247 {
248 mixingLength = asImp_().karmanConstant() * asImp_().wallDistance(elementIdx)
249 * (1.0 - exp(-yPlusValue / 26.0 ))
250 / sqrt(1.0 - exp(-0.26 * yPlusValue));
251 }
252
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);
257 }
258
260 const Scalar uStar(unsigned int elementIdx) const
261 {
262 using std::abs;
263 using std::sqrt;
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)));
269 }
270
272 const Scalar uStarNominal(unsigned int elementIdx) const
273 {
274 using std::pow;
275 using std::sqrt;
276 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
277 return pow(cMu(), 0.25) * sqrt(storedTurbulentKineticEnergy(matchingPointIndex));
278 }
279
283 const Scalar dissipationWallFunction(unsigned int elementIdx) const
284 {
285 return uStarNominal(elementIdx) * uStarNominal(elementIdx) * uStarNominal(elementIdx)
286 / asImp_().karmanConstant() / asImp_().wallDistance(elementIdx);
287 }
288
292 const Scalar turbulentKineticEnergyWallFunction(unsigned int elementIdx) const
293 {
294 unsigned int matchingPointIndex = matchingPointIdx(asImp_().wallElementIndex(elementIdx));
295 return storedTurbulentKineticEnergy(matchingPointIndex);
296 }
297
299 const Scalar tangentialMomentumWallFunction(unsigned int elementIdx, Scalar velocity) const
300 {
301 using std::log;
302 Scalar velocityNominal = uStarNominal(elementIdx) * (1.0 / asImp_().karmanConstant() * log(yPlusNominal(elementIdx)) + 5.0);
303 return uStarNominal(elementIdx) * uStarNominal(elementIdx) * velocity / velocityNominal;
304 }
305
307 bool useWallFunction(const Element& element,
308 const SubControlVolumeFace& localSubFace,
309 const int& eqIdx) const
310 {
311 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
312 auto bcTypes = asImp_().boundaryTypes(element, localSubFace);
313
314 return bcTypes.hasWall()
315 && bcTypes.isDirichlet(eqIdx)
316 && isMatchingPoint(elementIdx);
317 }
318
320 FacePrimaryVariables wallFunction(const Element& element,
321 const FVElementGeometry& fvGeometry,
322 const ElementVolumeVariables& elemVolVars,
323 const ElementFaceVariables& elemFaceVars,
324 const SubControlVolumeFace& scvf,
325 const SubControlVolumeFace& localSubFace) const
326 {
327 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
328 return FacePrimaryVariables(asImp_().tangentialMomentumWallFunction(elementIdx, elemFaceVars[scvf].velocitySelf())
329 * asImp_().storedDensity(elementIdx) );
330 }
331
333 template<bool eB = enableEnergyBalance, bool compositional = isCompositional,
334 typename std::enable_if_t<eB && compositional, int> = 0>
335 CellCenterPrimaryVariables wallFunction(const Element& element,
336 const FVElementGeometry& fvGeometry,
337 const ElementVolumeVariables& elemVolVars,
338 const ElementFaceVariables& elemFaceVars,
339 const SubControlVolumeFace& scvf) const
340 {
341 return wallFunctionComponent(element, fvGeometry, elemVolVars, elemFaceVars, scvf)
342 + wallFunctionEnergy(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
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 { return wallFunctionComponent(element, fvGeometry, elemVolVars, elemFaceVars, scvf); }
354
356 template<bool eB = enableEnergyBalance, bool compositional = isCompositional,
357 typename std::enable_if_t<eB && !compositional, int> = 0>
358 CellCenterPrimaryVariables wallFunction(const Element& element,
359 const FVElementGeometry& fvGeometry,
360 const ElementVolumeVariables& elemVolVars,
361 const ElementFaceVariables& elemFaceVars,
362 const SubControlVolumeFace& scvf) const
363 { return wallFunctionEnergy(element, fvGeometry, elemVolVars, elemFaceVars, scvf); }
364
366 template<bool eB = enableEnergyBalance, bool compositional = isCompositional,
367 typename std::enable_if_t<!eB && !compositional, int> = 0>
368 CellCenterPrimaryVariables wallFunction(const Element& element,
369 const FVElementGeometry& fvGeometry,
370 const ElementVolumeVariables& elemVolVars,
371 const ElementFaceVariables& elemFaceVars,
372 const SubControlVolumeFace& scvf) const
373 { return CellCenterPrimaryVariables(0.0); }
374
376 CellCenterPrimaryVariables wallFunctionComponent(const Element& element,
377 const FVElementGeometry& fvGeometry,
378 const ElementVolumeVariables& elemVolVars,
379 const ElementFaceVariables& elemFaceVars,
380 const SubControlVolumeFace& scvf) const
381 {
382 using std::log;
383 auto wallFunctionFlux = CellCenterPrimaryVariables(0.0);
384 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
385
386 // component mass fluxes
387 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
388 {
389 if (ModelTraits::replaceCompEqIdx() == compIdx)
390 continue;
391
392 Scalar schmidtNumber = elemVolVars[scvf.insideScvIdx()].kinematicViscosity()
393 / elemVolVars[scvf.insideScvIdx()].diffusionCoefficient(0, 0, compIdx);
394 Scalar moleToMassConversionFactor = ModelTraits::useMoles()
395 ? 1.0 : FluidSystem::molarMass(compIdx);
396 wallFunctionFlux[compIdx] +=
397 -1.0 * (asImp_().dirichlet(element, scvf)[Indices::conti0EqIdx + compIdx]
398 - elemVolVars[scvf.insideScvIdx()].moleFraction(compIdx))
399 * elemVolVars[scvf.insideScvIdx()].molarDensity()
400 * moleToMassConversionFactor
401 * uStarNominal(elementIdx)
402 / asImp_().turbulentSchmidtNumber()
403 / (1. / asImp_().karmanConstant() * log(yPlusNominal(elementIdx) * 9.793)
404 + pFunction(schmidtNumber, asImp_().turbulentSchmidtNumber()));
405 }
406
407 if (ModelTraits::replaceCompEqIdx() < ModelTraits::numFluidComponents())
408 {
409 wallFunctionFlux[ModelTraits::replaceCompEqIdx()] =
410 -std::accumulate(wallFunctionFlux.begin(), wallFunctionFlux.end(), 0.0);
411 }
412
413 return wallFunctionFlux;
414 }
415
417 CellCenterPrimaryVariables wallFunctionEnergy(const Element& element,
418 const FVElementGeometry& fvGeometry,
419 const ElementVolumeVariables& elemVolVars,
420 const ElementFaceVariables& elemFaceVars,
421 const SubControlVolumeFace& scvf) const
422 {
423 using std::log;
424 auto wallFunctionFlux = CellCenterPrimaryVariables(0.0);
425 unsigned int elementIdx = asImp_().gridGeometry().elementMapper().index(element);
426 // energy fluxes
427 Scalar prandtlNumber = elemVolVars[scvf.insideScvIdx()].kinematicViscosity()
428 * elemVolVars[scvf.insideScvIdx()].density()
429 * elemVolVars[scvf.insideScvIdx()].heatCapacity()
430 / elemVolVars[scvf.insideScvIdx()].thermalConductivity();
431 wallFunctionFlux[Indices::energyEqIdx - cellCenterOffset] +=
432 -1.0 * (asImp_().dirichlet(element, scvf)[Indices::temperatureIdx]
433 - elemVolVars[scvf.insideScvIdx()].temperature())
434 * elemVolVars[scvf.insideScvIdx()].density()
435 * elemVolVars[scvf.insideScvIdx()].heatCapacity()
436 * uStarNominal(elementIdx)
437 / asImp_().turbulentPrandtlNumber()
438 / (1. / asImp_().karmanConstant() * log(yPlusNominal(elementIdx) * 9.793)
439 + pFunction(prandtlNumber, asImp_().turbulentPrandtlNumber()));
440
441 return wallFunctionFlux;
442 }
443
445 const Scalar pFunction(Scalar molecularNumber, Scalar turbulentNumber) const
446 {
447 using std::pow;
448 using std::exp;
449 return 9.24
450 * (pow(molecularNumber / turbulentNumber, 0.75) - 1.0)
451 * (1.0 + 0.28 * exp(-0.007 * molecularNumber / turbulentNumber));
452 }
453
455 const Scalar cMu() const
456 { return 0.09; }
457
458 Scalar yPlusThreshold() const
459 {
460 static const Scalar yPlusThreshold = getParamFromGroup<Scalar>(this->paramGroup(), "KEpsilon.YPlusThreshold", 30.0);
461 return yPlusThreshold;
462 }
463
465 {
466 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(), "RANS.UseStoredEddyViscosity", false);
467 return useStoredEddyViscosity;
468 }
469
470 Scalar storedDissipation(const int elementIdx) const
471 { return storedDissipation_[elementIdx]; }
472
473 Scalar storedTurbulentKineticEnergy(const int elementIdx) const
474 { return storedTurbulentKineticEnergy_[elementIdx]; }
475
476 Scalar storedDynamicEddyViscosity(const int elementIdx) const
477 { return storedDynamicEddyViscosity_[elementIdx]; }
478
479 Scalar zeroEqDynamicEddyViscosity(const int elementIdx) const
480 { return zeroEqDynamicEddyViscosity_[elementIdx]; }
481
482 unsigned int matchingPointIdx(const int elementIdx) const
483 { return matchingPointIdx_[elementIdx]; }
484
485private:
486 std::vector<unsigned int> matchingPointIdx_;
487 std::vector<Scalar> storedDissipation_;
488 std::vector<Scalar> storedTurbulentKineticEnergy_;
489 std::vector<Scalar> storedDynamicEddyViscosity_;
490 std::vector<Scalar> zeroEqDynamicEddyViscosity_;
491
493 Implementation &asImp_()
494 { return *static_cast<Implementation *>(this); }
495
497 const Implementation &asImp_() const
498 { return *static_cast<const Implementation *>(this); }
499};
500
501} // end namespace Dumux
502
503#endif
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:55
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:458
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:117
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:464
unsigned int matchingPointIdx(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:482
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:455
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:335
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:473
Scalar storedDissipation(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:470
Scalar zeroEqDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:479
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:417
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:320
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 &paramGroup="")
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:476
const Scalar pFunction(Scalar molecularNumber, Scalar turbulentNumber) const
Returns the value of the P-function after Jayatilleke .
Definition: freeflow/rans/twoeq/kepsilon/problem.hh:445
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:376
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.