12#ifndef DUMUX_RANS_PROBLEM_HH
13#define DUMUX_RANS_PROBLEM_HH
17#include <dune/common/fmatrix.hh>
30template<
class TypeTag, TurbulenceModel turbulenceModel>
34template<
class TypeTag>
45template<
class TypeTag>
54 using FVElementGeometry =
typename GridGeometry::LocalView;
55 using GridView =
typename GridGeometry::GridView;
56 using Element =
typename GridView::template Codim<0>::Entity;
57 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
58 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
60 using PrimaryVariables =
typename VolumeVariables::PrimaryVariables;
65 using GlobalPosition =
typename SubControlVolumeFace::GlobalPosition;
67 static constexpr auto dim = GridView::dimension;
68 static constexpr int numCorners = SubControlVolumeFace::numCornersPerFace;
69 using DimVector = GlobalPosition;
70 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
72 struct WallElementInformation
75 unsigned int wallElementIdx;
77 unsigned int wallFaceNormalAxis;
79 GlobalPosition wallFaceCenter;
80 std::array<GlobalPosition, numCorners> wallFaceCorners;
95 std::cout <<
"The parameter \"Rans.IsFlatWallBounded\" is not specified. \n"
96 <<
" -- Based on the grid and the boundary conditions specified by the user,"
97 <<
" this parameter is set to be "<< std::boolalpha <<
isFlatWallBounded() <<
"\n";
101 wallDistance_.resize(this->
gridGeometry().elementMapper().size(), std::numeric_limits<Scalar>::max());
102 neighborIdx_.resize(this->
gridGeometry().elementMapper().size());
103 velocity_.resize(this->
gridGeometry().elementMapper().size(), DimVector(0.0));
104 velocityGradients_.resize(this->
gridGeometry().elementMapper().size(), DimMatrix(0.0));
105 stressTensorScalarProduct_.resize(this->
gridGeometry().elementMapper().size(), 0.0);
106 vorticityTensorScalarProduct_.resize(this->
gridGeometry().elementMapper().size(), 0.0);
107 flowDirectionAxis_.resize(this->
gridGeometry().elementMapper().size(), fixedFlowDirectionAxis_);
108 storedViscosity_.resize(this->
gridGeometry().elementMapper().size(), 0.0);
109 storedDensity_.resize(this->
gridGeometry().elementMapper().size(), 0.0);
117 std::cout <<
"Update static wall properties. ";
121 findWallDistances_();
122 findNeighborIndices_();
130 template<
class SolutionVector>
133 std::cout <<
"Update dynamic wall properties." << std::endl;
135 DUNE_THROW(Dune::InvalidStateException,
136 "You have to call updateStaticWallProperties() once before you call updateDynamicWallProperties().");
138 calculateCCVelocities_(curSol);
139 calculateCCVelocityGradients_();
140 calculateMaxMinVelocities_();
141 calculateStressTensor_();
142 calculateVorticityTensor_();
143 storeViscosities_(curSol);
154 const SubControlVolumeFace& scvf,
155 const int& eqIdx)
const
161 template<
class ElementVolumeVariables,
class ElementFaceVariables>
163 const FVElementGeometry& fvGeometry,
164 const ElementVolumeVariables& elemVolVars,
165 const ElementFaceVariables& elemFaceVars,
166 const SubControlVolumeFace& scvf,
167 const SubControlVolumeFace& lateralBoundaryFace)
const
168 {
return FacePrimaryVariables(0.0); }
173 template<
class ElementVolumeVariables,
class ElementFaceVariables>
175 const FVElementGeometry& fvGeometry,
176 const ElementVolumeVariables& elemVolVars,
177 const ElementFaceVariables& elemFaceVars,
178 const SubControlVolumeFace& scvf)
const
179 {
return CellCenterPrimaryVariables(0.0); }
186 static const bool hasAlignedWalls = hasAlignedWalls_();
187 return hasAlignedWalls;
207 = getParamFromGroup<Scalar>(this->
paramGroup(),
"RANS.TurbulentPrandtlNumber", 1.0);
218 = getParamFromGroup<Scalar>(this->
paramGroup(),
"RANS.TurbulentSchmidtNumber", 1.0);
225 DUNE_THROW(Dune::NotImplemented,
"\n Due to grid/geometric concerns, models requiring a wallNormalAxis "
226 <<
"can only be used for flat wall bounded flows. "
227 <<
"\n If your geometry is a flat channel, "
228 <<
"please set the runtime parameter RANS.IsFlatWallBounded to true. \n");
229 return wallNormalAxis_[elementIdx];
235 DUNE_THROW(Dune::NotImplemented,
"\n Due to grid/geometric concerns, models requiring a flowDirectionAxis "
236 <<
"can only be used for flat wall bounded flows. "
237 <<
"\n If your geometry is a flat channel, "
238 <<
"please set the runtime parameter RANS.IsFlatWallBounded to true. \n");
239 return flowDirectionAxis_[elementIdx];
245 DUNE_THROW(Dune::NotImplemented,
"\n Due to grid/geometric concerns, models requiring a wallElementIndex "
246 <<
"can only be used for flat wall bounded flows. "
247 <<
"\n If your geometry is a flat channel, "
248 <<
"please set the runtime parameter RANS.IsFlatWallBounded to true. \n");
249 return wallElementIdx_[elementIdx];
254 {
return wallDistance_[elementIdx]; }
258 const auto& element = this->
gridGeometry().element(elementIdx);
259 return element.geometry().center();
262 unsigned int neighborIndex(
const int elementIdx,
const int axisIdx,
const int sideIdx)
const
263 {
return neighborIdx_[elementIdx][axisIdx][sideIdx];}
266 {
return velocity_[elementIdx]; }
268 Scalar
ccVelocity(
const int elementIdx,
const int axisIdx)
const
269 {
return velocity_[elementIdx][axisIdx]; }
272 {
return velocityMaximum_[elementIdx]; }
275 {
return velocityMinimum_[elementIdx]; }
278 {
return velocityGradients_[elementIdx]; }
281 {
return velocityGradients_[elementIdx][i][j]; }
284 {
return stressTensorScalarProduct_[elementIdx]; }
287 {
return vorticityTensorScalarProduct_[elementIdx]; }
290 {
return storedViscosity_[elementIdx]; }
293 {
return storedDensity_[elementIdx]; }
302 bool hasAlignedWalls_()
const
310 std::vector<int> wallFaceAxis;
311 wallFaceAxis.reserve(this->
gridGeometry().numBoundaryScvf());
315 for (
const auto& element : elements(gridView))
317 fvGeometry.bindElement(element);
318 for (
const auto& scvf : scvfs(fvGeometry))
319 if (!scvf.boundary() && asImp_().
boundaryTypes(element, scvf).hasWall())
320 wallFaceAxis.push_back(scvf.directionIndex());
324 return std::all_of(wallFaceAxis.begin(), wallFaceAxis.end(), [firstDir=wallFaceAxis[0]](
auto dir){ return (dir == firstDir);} ) ;
327 void checkForWalls_()
329 for (
const auto& element : elements(this->
gridGeometry().gridView()))
332 fvGeometry.bindElement(element);
333 for (
auto&& scvf : scvfs(fvGeometry))
338 DUNE_THROW(Dune::InvalidStateException,
"No walls are are specified with the setWall() function");
346 void findWallDistances_()
349 [
this] (
const FVElementGeometry& fvGeometry,
const SubControlVolumeFace& scvf)
350 {
return asImp_().boundaryTypes(fvGeometry.element(), scvf).hasWall(); });
351 wallDistance_ = wallInformation.wallDistance();
352 storeWallElementAndDirectionIndex_(wallInformation.wallData());
355 template <
class WallData>
356 void storeWallElementAndDirectionIndex_(
const WallData& wallData)
360 DUNE_THROW(Dune::NotImplemented,
"The wall direction Index can only be calculated for quadrilateral structured grids");
365 wallNormalAxis_.resize(wallData.size());
366 wallElementIdx_.resize(wallData.size());
368 for (
const auto& element : elements(this->
gridGeometry().gridView()))
370 unsigned int elementIdx = this->
gridGeometry().elementMapper().index(element);
371 wallElementIdx_[elementIdx] = wallData[elementIdx].eIdx;
372 if ( ! (
hasParam(
"RANS.WallNormalAxis")) )
374 GlobalPosition wallOuterNormal = wallData[elementIdx].scvfOuterNormal;
375 if constexpr (dim == 2)
376 wallNormalAxis_[elementIdx] = (wallOuterNormal[0] == 1) ? 0 : 1;
378 wallNormalAxis_[elementIdx] = (wallOuterNormal[0] == 1) ? 0 : ((wallOuterNormal[1] == 1) ? 1 : 2);
381 wallNormalAxis_[elementIdx] = fixedWallNormalAxis_;
389 void findNeighborIndices_()
392 for (
const auto& element : elements(this->
gridGeometry().gridView()))
394 unsigned int elementIdx = this->
gridGeometry().elementMapper().index(element);
395 for (
unsigned int axisIdx = 0; axisIdx < dim; ++axisIdx)
397 neighborIdx_[elementIdx][axisIdx][0] = elementIdx;
398 neighborIdx_[elementIdx][axisIdx][1] = elementIdx;
401 for (
const auto& intersection : intersections(this->
gridGeometry().gridView(), element))
403 if (intersection.boundary())
406 unsigned int neighborIdx = this->
gridGeometry().elementMapper().index(intersection.outside());
407 for (
unsigned int axisIdx = 0; axisIdx < dim; ++axisIdx)
412 neighborIdx_[elementIdx][axisIdx][0] = neighborIdx;
415 neighborIdx_[elementIdx][axisIdx][1] = neighborIdx;
422 template<
class SolutionVector>
423 void calculateCCVelocities_(
const SolutionVector& curSol)
427 for (
const auto& element : elements(this->
gridGeometry().gridView()))
429 fvGeometry.bindElement(element);
430 unsigned int elementIdx = this->
gridGeometry().elementMapper().index(element);
433 DimVector velocityTemp(0.0);
434 for (
auto&& scvf : scvfs(fvGeometry))
436 const int dofIdxFace = scvf.dofIndex();
437 const auto numericalSolutionFace = curSol[GridGeometry::faceIdx()][dofIdxFace][Indices::velocity(scvf.directionIndex())];
438 velocityTemp[scvf.directionIndex()] += numericalSolutionFace;
440 for (
unsigned int axisIdx = 0; axisIdx < dim; ++axisIdx)
441 velocity_[elementIdx][axisIdx] = velocityTemp[axisIdx] * 0.5;
446 void calculateCCVelocityGradients_()
452 for (
const auto& element : elements(this->
gridGeometry().gridView()))
454 const unsigned int elementIdx = this->
gridGeometry().elementMapper().index(element);
456 for (
unsigned int j = 0; j < dim; ++j)
458 for (
unsigned int i = 0; i < dim; ++i)
460 const unsigned int neighborIndex0 =
neighborIndex(elementIdx, j, 0);
461 const unsigned int neighborIndex1 =
neighborIndex(elementIdx, j, 1);
463 velocityGradients_[elementIdx][i][j]
468 velocityGradients_[elementIdx][i][j] = 0.0;
472 fvGeometry.bindElement(element);
473 for (
auto&& scvf : scvfs(fvGeometry))
476 unsigned int axisIdx = scvf.directionIndex();
479 for (
unsigned int velIdx = 0; velIdx < dim; ++velIdx)
481 if (!asImp_().boundaryTypes(element, scvf).isDirichlet(Indices::velocity(velIdx)))
484 Scalar dirichletVelocity = asImp_().dirichlet(element, scvf)[Indices::velocity(velIdx)];
486 unsigned int neighborIdx =
neighborIndex(elementIdx, axisIdx, 0);
487 if (scvf.center()[axisIdx] <
cellCenter(elementIdx)[axisIdx])
490 velocityGradients_[elementIdx][velIdx][axisIdx]
491 = (
ccVelocity(neighborIdx, velIdx) - dirichletVelocity)
492 / (
cellCenter(neighborIdx)[axisIdx] - scvf.center()[axisIdx]);
497 std::vector<int> bjsNumFaces(dim, 0);
498 std::vector<unsigned int> bjsNeighbor(dim, 0);
499 DimVector bjsVelocityAverage(0.0);
500 DimVector normalNormCoordinate(0.0);
501 unsigned int velCompIdx = Indices::velocity(scvf.directionIndex());
502 const int numSubFaces = scvf.pairData().size();
503 for(
int localSubFaceIdx = 0; localSubFaceIdx < numSubFaces; ++localSubFaceIdx)
505 const auto& lateralFace = fvGeometry.scvf(scvf.insideScvIdx(), scvf.pairData()[localSubFaceIdx].localLateralFaceIdx);
508 unsigned int lateralAxisIdx = lateralFace.directionIndex();
509 if (lateralFace.boundary() && (asImp_().
boundaryTypes(element, lateralFace).isBeaversJoseph(Indices::velocity(velCompIdx))))
511 unsigned int neighborIdx =
neighborIndex(elementIdx, lateralAxisIdx, 0);
512 if (lateralFace.center()[lateralAxisIdx] <
cellCenter(elementIdx)[lateralAxisIdx])
515 const SubControlVolume& scv = fvGeometry.scv(scvf.insideScvIdx());
517 if (bjsNumFaces[lateralAxisIdx] > 0 && neighborIdx != bjsNeighbor[lateralAxisIdx])
518 DUNE_THROW(Dune::InvalidStateException,
"Two different neighborIdx should not occur");
519 bjsNeighbor[lateralAxisIdx] = neighborIdx;
520 normalNormCoordinate[lateralAxisIdx] = lateralFace.center()[lateralAxisIdx];
521 bjsNumFaces[lateralAxisIdx]++;
524 for (
unsigned axisIdx = 0; axisIdx < dim; ++axisIdx)
526 if (bjsNumFaces[axisIdx] == 0)
529 unsigned int neighborIdx = bjsNeighbor[axisIdx];
530 bjsVelocityAverage[axisIdx] /= bjsNumFaces[axisIdx];
532 velocityGradients_[elementIdx][velCompIdx][axisIdx]
533 = (
ccVelocity(neighborIdx, velCompIdx) - bjsVelocityAverage[axisIdx])
534 / (
cellCenter(neighborIdx)[axisIdx] - normalNormCoordinate[axisIdx]);
540 void calculateMaxMinVelocities_()
549 velocityMaximum_.assign(this->
gridGeometry().elementMapper().size(), DimVector(1e-16));
550 velocityMinimum_.assign(this->
gridGeometry().elementMapper().size(), DimVector(std::numeric_limits<Scalar>::max()));
553 for (
const auto& element : elements(this->
gridGeometry().gridView()))
555 const unsigned int elementIdx = this->
gridGeometry().elementMapper().index(element);
556 Scalar maxVelocity = 0.0;
559 for (
unsigned int axisIdx = 0; axisIdx < dim; ++axisIdx)
561 if (abs(
ccVelocity(elementIdx, axisIdx)) > abs(velocityMaximum_[wallElementIdx][axisIdx]))
562 velocityMaximum_[wallElementIdx][axisIdx] =
ccVelocity(elementIdx, axisIdx);
564 if (abs(
ccVelocity(elementIdx, axisIdx)) < abs(velocityMinimum_[wallElementIdx][axisIdx]))
565 velocityMinimum_[wallElementIdx][axisIdx] =
ccVelocity(elementIdx, axisIdx);
568 if ((
hasParam(
"RANS.FlowDirectionAxis") != 1) && (maxVelocity) < abs(
ccVelocity(elementIdx, axisIdx)))
570 maxVelocity = abs(
ccVelocity(elementIdx, axisIdx));
571 flowDirectionAxis_[elementIdx] = axisIdx;
581 DimVector maxVelocity(0.0);
582 DimVector minVelocity(std::numeric_limits<Scalar>::max());
584 for (
const auto& element : elements(this->
gridGeometry().gridView()))
586 const unsigned int elementIdx = this->
gridGeometry().elementMapper().index(element);
588 for (
unsigned int axisIdx = 0; axisIdx < dim; ++axisIdx)
590 if (abs(
ccVelocity(elementIdx, axisIdx)) > abs(maxVelocity[axisIdx]))
591 maxVelocity[axisIdx] =
ccVelocity(elementIdx, axisIdx);
593 if (abs(
ccVelocity(elementIdx, axisIdx)) < abs(minVelocity[axisIdx]))
594 minVelocity[axisIdx] =
ccVelocity(elementIdx, axisIdx);
597 velocityMaximum_.assign(this->
gridGeometry().elementMapper().size(), maxVelocity);
598 velocityMinimum_.assign(this->
gridGeometry().elementMapper().size(), minVelocity);
602 void calculateStressTensor_()
604 for (
const auto& element : elements(this->
gridGeometry().gridView()))
606 unsigned int elementIdx = this->
gridGeometry().elementMapper().index(element);
607 Dune::FieldMatrix<Scalar, GridView::dimension, GridView::dimension> stressTensor(0.0);
608 for (
unsigned int j = 0; j < dim; ++j)
610 for (
unsigned int i = 0; i < dim; ++i)
616 stressTensorScalarProduct_[elementIdx] = 0.0;
617 for (
unsigned int j = 0; j < dim; ++j)
619 for (
unsigned int i = 0; i < dim; ++i)
621 stressTensorScalarProduct_[elementIdx] += stressTensor[j][i] * stressTensor[j][i];
627 void calculateVorticityTensor_()
629 for (
const auto& element : elements(this->
gridGeometry().gridView()))
631 unsigned int elementIdx = this->
gridGeometry().elementMapper().index(element);
632 Dune::FieldMatrix<Scalar, GridView::dimension, GridView::dimension> vorticityTensor(0.0);
633 for (
unsigned int j = 0; j < dim; ++j)
635 for (
unsigned int i = 0; i < dim; ++i)
641 vorticityTensorScalarProduct_[elementIdx] = 0.0;
642 for (
unsigned int j = 0; j < dim; ++j)
644 for (
unsigned int i = 0; i < dim; ++i)
646 vorticityTensorScalarProduct_[elementIdx] += vorticityTensor[j][i] * vorticityTensor[j][i];
652 template<
class SolutionVector>
653 void storeViscosities_(
const SolutionVector& curSol)
657 for (
const auto& element : elements(this->
gridGeometry().gridView()))
659 unsigned int elementIdx = this->
gridGeometry().elementMapper().index(element);
660 fvGeometry.bindElement(element);
661 for (
auto&& scv : scvs(fvGeometry))
663 const int dofIdx = scv.dofIndex();
665 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
666 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
667 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
669 VolumeVariables volVars;
670 volVars.update(elemSol, asImp_(), element, scv);
671 storedDensity_[elementIdx] = volVars.density();
672 storedViscosity_[elementIdx] = volVars.viscosity();
677 const int fixedFlowDirectionAxis_ = getParam<int>(
"RANS.FlowDirectionAxis", 0);
678 const int fixedWallNormalAxis_ = getParam<int>(
"RANS.WallNormalAxis", 1);
680 std::vector<unsigned int> wallNormalAxis_;
681 std::vector<unsigned int> flowDirectionAxis_;
682 std::vector<Scalar> wallDistance_;
683 std::vector<unsigned int> wallElementIdx_;
684 std::vector<std::array<std::array<unsigned int, 2>, dim>> neighborIdx_;
686 std::vector<DimVector> velocity_;
687 std::vector<DimVector> velocityMaximum_;
688 std::vector<DimVector> velocityMinimum_;
689 std::vector<DimMatrix> velocityGradients_;
691 std::vector<Scalar> stressTensorScalarProduct_;
692 std::vector<Scalar> vorticityTensorScalarProduct_;
694 std::vector<Scalar> storedDensity_;
695 std::vector<Scalar> storedViscosity_;
698 Implementation &asImp_()
699 {
return *
static_cast<Implementation *
>(
this); }
702 const Implementation &asImp_()
const
703 {
return *
static_cast<const Implementation *
>(
this); }
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:43
const std::string & paramGroup() const
The parameter group in which to retrieve runtime parameters.
Definition: common/fvproblem.hh:528
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvproblem.hh:524
auto boundaryTypes(const Element &element, const SubControlVolume &scv) const
Specifies which kind of boundary condition should be used for which equation on a given boundary segm...
Definition: common/fvproblem.hh:130
Navier-Stokes staggered problem base class.
Definition: freeflow/navierstokes/staggered/problem.hh:33
const Scalar beaversJosephVelocity(const Element &element, const SubControlVolume &scv, const SubControlVolumeFace &ownScvf, const SubControlVolumeFace &faceOnPorousBoundary, const Scalar velocitySelf, const Scalar tangentialVelocityGradient) const
Returns the slip velocity at a porous boundary based on the Beavers-Joseph(-Saffman) condition.
Definition: freeflow/navierstokes/staggered/problem.hh:189
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:47
unsigned int wallElementIndex(const int elementIdx) const
Definition: freeflow/rans/problem.hh:242
int flowDirectionAxis(const int elementIdx) const
Definition: freeflow/rans/problem.hh:232
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) turbulence parameters.
Definition: freeflow/rans/problem.hh:131
DimMatrix velocityGradientTensor(const int elementIdx) const
Definition: freeflow/rans/problem.hh:277
Scalar vorticityTensorScalarProduct(const int elementIdx) const
Definition: freeflow/rans/problem.hh:286
Scalar ccVelocity(const int elementIdx, const int axisIdx) const
Definition: freeflow/rans/problem.hh:268
void updateStaticWallProperties()
Update the static (solution independent) relations to the walls and neighbors.
Definition: freeflow/rans/problem.hh:115
GlobalPosition cellCenter(const int elementIdx) const
Definition: freeflow/rans/problem.hh:256
FacePrimaryVariables wallFunction(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, const SubControlVolumeFace &lateralBoundaryFace) const
Returns an additional wall function momentum flux.
Definition: freeflow/rans/problem.hh:162
Scalar velocityGradient(const int elementIdx, const int i, const int j) const
Definition: freeflow/rans/problem.hh:280
DimVector velocityMaximum(const int elementIdx) const
Definition: freeflow/rans/problem.hh:271
const Scalar betaOmega() const
Returns the constant.
Definition: freeflow/rans/problem.hh:197
DimVector ccVelocityVector(const int elementIdx) const
Definition: freeflow/rans/problem.hh:265
Scalar kinematicViscosity(const int elementIdx) const
Definition: freeflow/rans/problem.hh:295
Scalar turbulentSchmidtNumber() const
Return the turbulent Schmidt number which is used to convert the eddy viscosity to an eddy diffusivi...
Definition: freeflow/rans/problem.hh:215
CellCenterPrimaryVariables wallFunction(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const
Returns an additional wall function flux for cell-centered quantities.
Definition: freeflow/rans/problem.hh:174
DimVector velocityMinimum(const int elementIdx) const
Definition: freeflow/rans/problem.hh:274
Scalar wallDistance(const int elementIdx) const
Definition: freeflow/rans/problem.hh:253
Scalar storedDensity(const int elementIdx) const
Definition: freeflow/rans/problem.hh:292
Scalar turbulentPrandtlNumber() const
Return the turbulent Prandtl number which is used to convert the eddy viscosity to an eddy thermal c...
Definition: freeflow/rans/problem.hh:204
int wallNormalAxis(const int elementIdx) const
Definition: freeflow/rans/problem.hh:222
unsigned int neighborIndex(const int elementIdx, const int axisIdx, const int sideIdx) const
Definition: freeflow/rans/problem.hh:262
RANSProblemBase(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
The constructor.
Definition: freeflow/rans/problem.hh:90
bool isFlatWallBounded() const
Returns whether a given sub control volume face is on a wall.
Definition: freeflow/rans/problem.hh:184
Scalar stressTensorScalarProduct(const int elementIdx) const
Definition: freeflow/rans/problem.hh:283
bool calledUpdateStaticWallProperties
Definition: freeflow/rans/problem.hh:298
bool useWallFunction(const Element &element, const SubControlVolumeFace &scvf, const int &eqIdx) const
Returns whether a wall function should be used at a given face.
Definition: freeflow/rans/problem.hh:153
const Scalar karmanConstant() const
Returns the Karman constant.
Definition: freeflow/rans/problem.hh:193
Scalar storedViscosity(const int elementIdx) const
Definition: freeflow/rans/problem.hh:289
forward declare
Definition: freeflow/rans/problem.hh:31
static constexpr auto atElementCenters
Definition: walldistance.hh:98
Defines all properties used in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
bool hasParamInGroup(const std::string ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:165
bool hasParam(const std::string ¶m)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:157
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Free function to get the local view of a grid cache object.
The available discretization methods in Dumux.
constexpr Staggered staggered
Definition: method.hh:149
Adaption of the fully implicit scheme to the tracer transport model.
The local element solution class for staggered methods.
Base class for all staggered fv problems.