|
|
|
RASModel Class ReferenceAbstract base class for turbulence models for compressible and combusting flows. More...
Inheritance diagram for RASModel:
![]()
Collaboration diagram for RASModel:
![]()
Detailed DescriptionAbstract base class for turbulence models for compressible and combusting flows.
Definition at line 67 of file RASModel.H. Constructor & Destructor Documentation
Construct from components.
Definition at line 49 of file RASModel.C. References RASModel::coeffDict_, Foam::endl(), Foam::Info, and Foam::type().
Here is the call graph for this function:
![]()
Member Function Documentation
Print model coefficients.
Definition at line 37 of file RASModel.C. Referenced by LaunderGibsonRSTM::devRhoReff(), LRR::devRhoReff(), RNGkEpsilon::R(), kOmegaSST::R(), realizableKE::R(), LaunderSharmaKE::R(), SpalartAllmaras::R(), and kEpsilon::R().
Here is the caller graph for this function:
![]()
Runtime type information.
Return a reference to the selected turbulence model.
Definition at line 92 of file RASModel.C. References turbulenceModel::mesh_.
Return the value of k0 which k is not allowed to be less than.
Definition at line 183 of file RASModel.H. References RASModel::k0(), and RASModel::k0_. Referenced by RASModel::k0().
Here is the call graph for this function:
![]()
Here is the caller graph for this function:
![]()
Return the value of epsilon0 which epsilon is not allowed to be. less than Definition at line 192 of file RASModel.H.
Return the value of epsilonSmall which is added to epsilon when. calculating nut Definition at line 201 of file RASModel.H. References RASModel::epsilonSmall(), and RASModel::epsilonSmall_. Referenced by RASModel::epsilonSmall().
Here is the call graph for this function:
![]()
Here is the caller graph for this function:
![]()
Return the value of omega0 which epsilon is not allowed to be. less than Definition at line 210 of file RASModel.H. References RASModel::omega0(), and RASModel::omega0_. Referenced by RASModel::omega0().
Here is the call graph for this function:
![]()
Here is the caller graph for this function:
![]()
Return the value of omegaSmall which is added to epsilon when. calculating nut Definition at line 219 of file RASModel.H. References RASModel::omegaSmall_.
Allow epsilon0 to be changed.
Definition at line 231 of file RASModel.H. References RASModel::epsilon0_.
Allow epsilonSmall to be changed.
Definition at line 237 of file RASModel.H. References RASModel::epsilonSmall_.
Allow omega0 to be changed.
Definition at line 243 of file RASModel.H. References RASModel::omega0_.
Allow omegaSmall to be changed.
Definition at line 249 of file RASModel.H. References RASModel::omegaSmall_.
Return the near wall distances.
Definition at line 255 of file RASModel.H. References RASModel::y_. Referenced by mutSpalartAllmarasWallFunctionFvPatchScalarField::calcUTau(), mutSpalartAllmarasStandardWallFunctionFvPatchScalarField::calcYPlus(), and mutSpalartAllmarasStandardRoughWallFunctionFvPatchScalarField::calcYPlus().
Here is the caller graph for this function:
![]()
Const access to the coefficients dictionary.
Definition at line 264 of file RASModel.H. References RASModel::coeffDict_.
Return the turbulence viscosity.
Implements turbulenceModel. Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon. Referenced by RASModel::muEff().
Here is the caller graph for this function:
![]()
Return the effective viscosity.
Implements turbulenceModel. Reimplemented in laminar. Definition at line 274 of file RASModel.H. References RASModel::mut().
Here is the call graph for this function:
![]()
Return the effective turbulent thermal diffusivity.
Implements turbulenceModel. Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.
Return the turbulence kinetic energy.
Implements turbulenceModel. Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.
Return the turbulence kinetic energy dissipation rate.
Implements turbulenceModel. Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.
Return the Reynolds stress tensor.
Implements turbulenceModel. Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.
Return the effective stress tensor including the laminar stress.
Implements turbulenceModel. Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.
Return the source term for the momentum equation.
Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon.
Return yPlus for the given patch.
Definition at line 161 of file RASModel.C. References Foam::log().
Here is the call graph for this function:
![]()
Solve the turbulence equations and correct the turbulence viscosity.
Implements turbulenceModel. Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon. Definition at line 193 of file RASModel.C.
Read RASProperties dictionary.
Implements turbulenceModel. Implemented in kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, realizableKE, RNGkEpsilon, SpalartAllmaras, and PDRkEpsilon. Definition at line 202 of file RASModel.C. References polyMesh::changing(), nearWallDist::correct(), turbulenceModel::mesh_, and RASModel::y_.
Here is the call graph for this function:
![]()
Member Data Documentation
Model coefficients dictionary.
Definition at line 84 of file RASModel.H. Referenced by RASModel::coeffDict(), and RASModel::RASModel().
Lower limit of epsilon.
Definition at line 93 of file RASModel.H. Referenced by RASModel::epsilon0().
Small epsilon value used to avoid divide by zero.
Definition at line 96 of file RASModel.H. Referenced by RASModel::epsilonSmall(), kEpsilon::kEpsilon(), and RNGkEpsilon::RNGkEpsilon().
Small omega value used to avoid divide by zero.
Definition at line 102 of file RASModel.H. Referenced by RASModel::omegaSmall().
Near wall distance boundary field.
Definition at line 105 of file RASModel.H. Referenced by RASModel::read(), and RASModel::y().
The documentation for this class was generated from the following files:
| ||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||