Public Member Functions | |
TypeName ("DATurbulenceModel") | |
declareRunTimeSelectionTable (autoPtr, DATurbulenceModel, dictionary,(const word modelType, const fvMesh &mesh, const DAOption &daOption),(modelType, mesh, daOption)) | |
DATurbulenceModel (const word modelType, const fvMesh &mesh, const DAOption &daOption) | |
virtual | ~DATurbulenceModel () |
void | correctWallDist () |
update wall distance for d_. Note: y_ will be automatically updated in mesh_ object More... | |
void | correctAlphat () |
update alphat More... | |
virtual void | correctNut ()=0 |
update nut based on other turbulence variables and update the BCs More... | |
virtual void | correctModelStates (wordList &modelStates) const =0 |
update the turbulence state for DAStateInfo::regStates_ More... | |
virtual void | correctBoundaryConditions ()=0 |
update turbulence variable boundary values More... | |
virtual void | updateIntermediateVariables ()=0 |
update any intermediate variables that are dependent on state variables and are used in calcResiduals More... | |
virtual void | correctStateResidualModelCon (List< List< word >> &stateCon) const =0 |
update the original variable connectivity for the adjoint state residuals in stateCon More... | |
virtual void | addModelResidualCon (HashTable< List< List< word >>> &allCon) const =0 |
add the model residual connectivity to stateCon More... | |
virtual void | calcResiduals (const dictionary &options)=0 |
compute the turbulence residuals More... | |
virtual void | correct (label printToScreen)=0 |
solve the residual equations and update the state More... | |
virtual void | getTurbProdTerm (volScalarField &prodTerm) const |
return the value of the production term from the turbulence model More... | |
virtual void | getTurbProdOverDestruct (volScalarField &PoD) const |
return the ratio of the production over destruction term from the turbulence model More... | |
virtual void | getTurbConvOverProd (volScalarField &CoP) const |
return the value of the convective over production term from the turbulence model More... | |
tmp< volSymmTensorField > | devRhoReff () const |
dev terms More... | |
tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) |
divDev terms More... | |
tmp< fvVectorMatrix > | divDevReff (volVectorField &U) |
divDev terms More... | |
tmp< volScalarField > | nuEff () const |
return effective viscosity More... | |
tmp< volScalarField > | nut () |
get the nut field More... | |
tmp< volScalarField > | alphaEff () |
return effective thermal diffusivity More... | |
tmp< volScalarField > | nu () const |
get the nu field More... | |
tmp< volScalarField > | alpha () const |
get alpha field More... | |
tmp< volScalarField > | rho () const |
get the density field More... | |
dimensionSet | rhoDimensions () const |
return the dimension of rho More... | |
tmp< volScalarField > | phase () |
get the phase field More... | |
scalar | prt () |
get the turbulent Prandtl number More... | |
tmp< Foam::volScalarField > | mu () const |
get mu More... | |
bool | writeData (Ostream &os) const |
this is a virtual function for regIOobject More... | |
void | printYPlus (const label printToScreen) const |
print the min max and mean yPlus to screen More... | |
word | getTurbModelType () const |
label | isPrintTime (const Time &runTime, const label printInterval) const |
virtual void | invTranProdNuTildaEqn (const volScalarField &mySource, volScalarField &pseudoNuTilda) |
Inverse transpose product, M_nuTilda^(-T) More... | |
virtual void | constructPseudoNuTildaEqn () |
virtual void | rhsSolvePseudoNuTildaEqn (const volScalarField &nuTildaSource) |
virtual void | calcLduResidualTurb (volScalarField &nuTildaRes) |
calculate the turbulence residual using LDU matrix More... | |
virtual void | getFvMatrixFields (const word varName, scalarField &diag, scalarField &upper, scalarField &lower) |
return the diag(), upper(), and lower() scalarFields from the turbulence model's fvMatrix More... | |
virtual void | solveAdjointFP (const word varName, const scalarList &rhs, scalarList &dPsi) |
solve the fvMatrixT field with given rhs and solution More... | |
Static Public Member Functions | |
static autoPtr< DATurbulenceModel > | New (const word modelType, const fvMesh &mesh, const DAOption &daOption) |
Protected Attributes | |
const fvMesh & | mesh_ |
fvMesh More... | |
const DAOption & | daOption_ |
DAOption object. More... | |
const dictionary & | allOptions_ |
all DAFoam options More... | |
DAGlobalVar & | daGlobalVar_ |
global variable More... | |
volScalarField & | nut_ |
turbulence viscosity More... | |
volVectorField & | U_ |
velocity More... | |
surfaceScalarField & | phi_ |
face flux More... | |
volScalarField | phase_ |
phase field More... | |
surfaceScalarField & | phaseRhoPhi_ |
phase*phi*density field More... | |
volScalarField | rhoOne_ |
an uniform rho field with ones More... | |
word | turbModelType_ = "None" |
whether the turbulence model is incompressible or compressible More... | |
IOdictionary | turbDict_ |
turbulence model property dict More... | |
dictionary | coeffDict_ |
turbulence model parameters dict More... | |
dimensionedScalar | kMin_ |
Lower limit of k. More... | |
dimensionedScalar | epsilonMin_ |
Lower limit of epsilon. More... | |
dimensionedScalar | omegaMin_ |
Lower limit for omega. More... | |
dimensionedScalar | nuTildaMin_ |
Lower limit for nuTilda. More... | |
scalar | Pr_ |
Prandtl number. More... | |
scalar | Prt_ = -9999.0 |
turbulent Prandtl number More... | |
Definition at line 52 of file DATurbulenceModel.H.
DATurbulenceModel | ( | const word | modelType, |
const fvMesh & | mesh, | ||
const DAOption & | daOption | ||
) |
Definition at line 22 of file DATurbulenceModel.C.
References mesh, DATurbulenceModel::mesh_, DATurbulenceModel::Pr_, DATurbulenceModel::Prt_, and DATurbulenceModel::turbModelType_.
|
inlinevirtual |
Definition at line 152 of file DATurbulenceModel.H.
TypeName | ( | "DATurbulenceModel" | ) |
declareRunTimeSelectionTable | ( | autoPtr | , |
DATurbulenceModel | , | ||
dictionary | , | ||
(const word modelType, const fvMesh &mesh, const DAOption &daOption) | , | ||
(modelType, mesh, daOption) | |||
) |
|
static |
Definition at line 161 of file DATurbulenceModel.C.
References DAOption::getAllOptions(), and mesh.
Referenced by DATurboFoam::initSolver(), DARhoSimpleCFoam::initSolver(), DARhoPimpleFoam::initSolver(), DARhoSimpleFoam::initSolver(), DAPimpleFoam::initSolver(), DAPimpleDyMFoam::initSolver(), and DASimpleFoam::initSolver().
void correctWallDist | ( | ) |
update wall distance for d_. Note: y_ will be automatically updated in mesh_ object
Definition at line 422 of file DATurbulenceModel.C.
References DATurbulenceModel::correctBoundaryConditions().
void correctAlphat | ( | ) |
update alphat
Definition at line 203 of file DATurbulenceModel.C.
References alphat, DATurbulenceModel::mesh_, DATurbulenceModel::nut_, Prt, DATurbulenceModel::Prt_, and DATurbulenceModel::rho().
Referenced by DASpalartAllmaras::correctNut(), DAkOmega::correctNut(), DAkEpsilon::correctNut(), DASpalartAllmarasFv3::correctNut(), DAkOmegaSST::correctNut(), and DAkOmegaSSTLM::correctNut().
|
pure virtual |
update nut based on other turbulence variables and update the BCs
Implemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, DASpalartAllmaras, and DADummyTurbulenceModel.
|
pure virtual |
update the turbulence state for DAStateInfo::regStates_
Implemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, DASpalartAllmaras, and DADummyTurbulenceModel.
Referenced by DAModel::correctModelStates().
|
pure virtual |
update turbulence variable boundary values
Implemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, DASpalartAllmaras, and DADummyTurbulenceModel.
Referenced by DAModel::correctBoundaryConditions(), and DATurbulenceModel::correctWallDist().
|
pure virtual |
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Implemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, DASpalartAllmaras, and DADummyTurbulenceModel.
Referenced by DAModel::updateIntermediateVariables().
|
pure virtual |
update the original variable connectivity for the adjoint state residuals in stateCon
Implemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, DASpalartAllmaras, and DADummyTurbulenceModel.
Referenced by DAModel::correctStateResidualModelCon().
|
pure virtual |
add the model residual connectivity to stateCon
Implemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, DASpalartAllmaras, and DADummyTurbulenceModel.
Referenced by DAModel::addModelResidualCon().
|
pure virtual |
compute the turbulence residuals
Implemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, DASpalartAllmaras, and DADummyTurbulenceModel.
Referenced by DAModel::calcResiduals().
|
pure virtual |
solve the residual equations and update the state
Implemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, DASpalartAllmaras, and DADummyTurbulenceModel.
|
virtual |
return the value of the production term from the turbulence model
Definition at line 540 of file DATurbulenceModel.C.
Referenced by DAModel::getTurbProdTerm().
|
virtual |
return the ratio of the production over destruction term from the turbulence model
Reimplemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, and DASpalartAllmaras.
Definition at line 552 of file DATurbulenceModel.C.
Referenced by DAModel::getTurbProdOverDestruct().
|
virtual |
return the value of the convective over production term from the turbulence model
Reimplemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, and DASpalartAllmaras.
Definition at line 564 of file DATurbulenceModel.C.
Referenced by DAModel::getTurbConvOverProd().
tmp< volSymmTensorField > devRhoReff | ( | ) | const |
dev terms
Definition at line 371 of file DATurbulenceModel.C.
References DATurbulenceModel::mesh_, DATurbulenceModel::nuEff(), DATurbulenceModel::phase_, DATurbulenceModel::rho(), and DATurbulenceModel::U_.
Referenced by DAFunctionMoment::calcFunction(), DAFunctionForce::calcFunction(), DAFunctionVariance::calcFunction(), DAResidualTurboFoam::calcResiduals(), and DAOutputForceCoupling::run().
tmp< fvVectorMatrix > divDevRhoReff | ( | volVectorField & | U | ) |
divDev terms
Definition at line 389 of file DATurbulenceModel.C.
References divScheme, DATurbulenceModel::nuEff(), DATurbulenceModel::phase_, DATurbulenceModel::rho(), T, DATurbulenceModel::turbModelType_, and U.
Referenced by DAResidualRhoPimpleFoam::calcPCMatWithFvMatrix(), DAResidualRhoSimpleCFoam::calcResiduals(), DAResidualTurboFoam::calcResiduals(), DAResidualRhoPimpleFoam::calcResiduals(), DAResidualRhoSimpleFoam::calcResiduals(), and DATurbulenceModel::divDevReff().
tmp< fvVectorMatrix > divDevReff | ( | volVectorField & | U | ) |
divDev terms
Definition at line 411 of file DATurbulenceModel.C.
References DATurbulenceModel::divDevRhoReff(), and U.
Referenced by DAResidualPimpleDyMFoam::calcPCMatWithFvMatrix(), DAResidualPimpleFoam::calcPCMatWithFvMatrix(), DAResidualPimpleFoam::calcResiduals(), DAResidualPimpleDyMFoam::calcResiduals(), and DAResidualSimpleFoam::calcResiduals().
tmp< volScalarField > nuEff | ( | ) | const |
return effective viscosity
Definition at line 220 of file DATurbulenceModel.C.
References DATurbulenceModel::nu(), and DATurbulenceModel::nut_.
Referenced by DATurbulenceModel::devRhoReff(), DATurbulenceModel::divDevRhoReff(), DATurbulenceModel::printYPlus(), and DAPimpleFoam::solveAdjointFP().
|
inline |
get the nut field
Definition at line 210 of file DATurbulenceModel.H.
References DATurbulenceModel::nut_.
tmp< volScalarField > alphaEff | ( | ) |
return effective thermal diffusivity
Definition at line 233 of file DATurbulenceModel.C.
References DATurbulenceModel::alpha(), alphat, DATurbulenceModel::mesh_, thermo, and DATurbulenceModel::turbModelType_.
Referenced by DAFunctionWallHeatFlux::calcFunction(), DAFunctionVariance::calcFunction(), DAOutputThermalCoupling::run(), and DAInputThermalCoupling::run().
tmp< volScalarField > nu | ( | ) | const |
get the nu field
Definition at line 267 of file DATurbulenceModel.C.
References DATurbulenceModel::mesh_, DATurbulenceModel::mu(), DATurbulenceModel::rho(), and DATurbulenceModel::turbModelType_.
Referenced by DATurbulenceModel::alpha(), DARegression::calcInputFeatures(), DAResidualPimpleDyMFoam::calcPCMatWithFvMatrix(), DAResidualPimpleFoam::calcPCMatWithFvMatrix(), DAResidualPimpleFoam::calcResiduals(), DAResidualPimpleDyMFoam::calcResiduals(), DAResidualSimpleFoam::calcResiduals(), DAkOmegaSSTLM::calcResiduals(), DASpalartAllmarasFv3::chi(), DASpalartAllmaras::chi(), DAkEpsilon::DepsilonEff(), DAkOmegaSSTLM::DgammaIntEff(), DAkOmega::DkEff(), DAkEpsilon::DkEff(), DAkOmegaSST::DkEff(), DAkOmegaSSTLM::DkEff(), DASpalartAllmaras::DnuTildaEff(), DASpalartAllmarasFv3::DnuTildaEff(), DAkOmega::DomegaEff(), DAkOmegaSST::DomegaEff(), DAkOmegaSSTLM::DomegaEff(), DAkOmegaSSTLM::DReThetatEff(), DAkOmegaSST::F1(), DAkOmegaSSTLM::F1(), DAkOmegaSSTLM::F1SST(), DAkOmegaSST::F2(), DAkOmegaSSTLM::F2(), DAkOmegaSST::F3(), DAkOmegaSSTLM::F3(), DAkOmegaSSTLM::Flength(), DAkOmegaSSTLM::Fthetat(), DATurbulenceModel::mu(), DATurbulenceModel::nuEff(), DATurbulenceModel::printYPlus(), and DAkOmegaSSTLM::ReThetat0().
tmp< volScalarField > alpha | ( | ) | const |
get alpha field
Definition at line 340 of file DATurbulenceModel.C.
References DATurbulenceModel::nu(), and DATurbulenceModel::Pr_.
Referenced by DATurbulenceModel::alphaEff(), and DATurbulenceModel::rho().
tmp< volScalarField > rho | ( | ) | const |
get the density field
Definition at line 294 of file DATurbulenceModel.C.
References DATurbulenceModel::alpha(), DATurbulenceModel::mesh_, DATurbulenceModel::rhoOne_, and DATurbulenceModel::turbModelType_.
Referenced by DAFunctionMassFlowRate::calcFunction(), DASpalartAllmaras::calcResiduals(), DAkOmega::calcResiduals(), DAkEpsilon::calcResiduals(), DASpalartAllmarasFv3::calcResiduals(), DAkOmegaSST::calcResiduals(), DAkOmegaSSTLM::calcResiduals(), DATurbulenceModel::correctAlphat(), DATurbulenceModel::devRhoReff(), DATurbulenceModel::divDevRhoReff(), DASpalartAllmaras::getFvMatrixFields(), DAkOmega::getFvMatrixFields(), DAkEpsilon::getFvMatrixFields(), DASpalartAllmarasFv3::getFvMatrixFields(), DAkOmegaSST::getFvMatrixFields(), DASpalartAllmaras::getTurbConvOverProd(), DAkOmega::getTurbConvOverProd(), DAkEpsilon::getTurbConvOverProd(), DASpalartAllmarasFv3::getTurbConvOverProd(), DAkOmegaSST::getTurbConvOverProd(), DAkOmegaSSTLM::getTurbConvOverProd(), DASpalartAllmaras::getTurbProdOverDestruct(), DAkOmega::getTurbProdOverDestruct(), DAkEpsilon::getTurbProdOverDestruct(), DASpalartAllmarasFv3::getTurbProdOverDestruct(), DAkOmegaSST::getTurbProdOverDestruct(), DAkOmegaSSTLM::getTurbProdOverDestruct(), DATurbulenceModel::nu(), DATurbulenceModel::rhoDimensions(), and DASpalartAllmaras::solveAdjointFP().
dimensionSet rhoDimensions | ( | ) | const |
return the dimension of rho
Definition at line 321 of file DATurbulenceModel.C.
References DATurbulenceModel::mesh_, DATurbulenceModel::rho(), DATurbulenceModel::rhoOne_, and DATurbulenceModel::turbModelType_.
Referenced by DAkEpsilon::epsilonSource(), DAkEpsilon::kSource(), DAkOmegaSSTLM::kSource(), DAkOmegaSST::kSource(), DAkOmegaSSTLM::omegaSource(), DAkOmegaSST::omegaSource(), DAkOmegaSSTLM::Qsas(), and DAkOmegaSST::Qsas().
|
inline |
get the phase field
Definition at line 231 of file DATurbulenceModel.H.
References DATurbulenceModel::phase_.
|
inline |
get the turbulent Prandtl number
Definition at line 237 of file DATurbulenceModel.H.
References DATurbulenceModel::Prt_.
tmp< Foam::volScalarField > mu | ( | ) | const |
get mu
Definition at line 349 of file DATurbulenceModel.C.
References DATurbulenceModel::mesh_, DATurbulenceModel::nu(), and DATurbulenceModel::turbModelType_.
Referenced by DATurbulenceModel::nu().
bool writeData | ( | Ostream & | os | ) | const |
this is a virtual function for regIOobject
Definition at line 197 of file DATurbulenceModel.C.
void printYPlus | ( | const label | printToScreen | ) | const |
print the min max and mean yPlus to screen
Definition at line 438 of file DATurbulenceModel.C.
References forAll(), DATurbulenceModel::mesh_, DATurbulenceModel::nu(), DATurbulenceModel::nuEff(), DATurbulenceModel::nut_, and DATurbulenceModel::U_.
|
inline |
Definition at line 260 of file DATurbulenceModel.H.
References DATurbulenceModel::turbModelType_.
Referenced by DAFunctionWallHeatFlux::calcFunction(), DAFunctionVariance::calcFunction(), DAFunctionVariance::DAFunctionVariance(), DAFunctionWallHeatFlux::DAFunctionWallHeatFlux(), DAOutputThermalCoupling::run(), and DAInputThermalCoupling::run().
label isPrintTime | ( | const Time & | runTime, |
const label | printInterval | ||
) | const |
Definition at line 521 of file DATurbulenceModel.C.
References runTime.
|
virtual |
Inverse transpose product, M_nuTilda^(-T)
Reimplemented in DASpalartAllmarasFv3.
Definition at line 576 of file DATurbulenceModel.C.
|
virtual |
Reimplemented in DASpalartAllmarasFv3.
Definition at line 602 of file DATurbulenceModel.C.
|
virtual |
Reimplemented in DASpalartAllmarasFv3.
Definition at line 614 of file DATurbulenceModel.C.
|
virtual |
calculate the turbulence residual using LDU matrix
Reimplemented in DASpalartAllmarasFv3.
Definition at line 590 of file DATurbulenceModel.C.
|
virtual |
return the diag(), upper(), and lower() scalarFields from the turbulence model's fvMatrix
Reimplemented in DAkOmegaSSTLM, DAkOmegaSST, DASpalartAllmarasFv3, DAkEpsilon, DAkOmega, and DASpalartAllmaras.
Definition at line 627 of file DATurbulenceModel.C.
Referenced by DASolver::calcPCMatWithFvMatrix().
|
virtual |
solve the fvMatrixT field with given rhs and solution
Reimplemented in DASpalartAllmaras.
Definition at line 639 of file DATurbulenceModel.C.
|
protected |
fvMesh
Definition at line 65 of file DATurbulenceModel.H.
Referenced by DASpalartAllmaras::addModelResidualCon(), DAkOmega::addModelResidualCon(), DAkEpsilon::addModelResidualCon(), DASpalartAllmarasFv3::addModelResidualCon(), DAkOmegaSST::addModelResidualCon(), DAkOmegaSSTLM::addModelResidualCon(), DATurbulenceModel::alphaEff(), DATurbulenceModel::correctAlphat(), DAkOmegaSST::correctNut(), DAkOmegaSSTLM::correctNut(), DATurbulenceModel::DATurbulenceModel(), DATurbulenceModel::devRhoReff(), DAkOmegaSSTLM::Flength(), DATurbulenceModel::mu(), DATurbulenceModel::nu(), DATurbulenceModel::printYPlus(), DAkOmegaSSTLM::ReThetac(), DAkOmegaSSTLM::ReThetat0(), DATurbulenceModel::rho(), DATurbulenceModel::rhoDimensions(), DAkEpsilon::saveEpsilonNearWall(), DAkOmega::saveOmegaNearWall(), DAkOmegaSST::saveOmegaNearWall(), DAkOmegaSSTLM::saveOmegaNearWall(), DAkEpsilon::setEpsilonNearWall(), DAkOmega::setOmegaNearWall(), DAkOmegaSST::setOmegaNearWall(), DAkOmegaSSTLM::setOmegaNearWall(), and DASpalartAllmaras::solveAdjointFP().
|
protected |
DAOption object.
Definition at line 68 of file DATurbulenceModel.H.
|
protected |
all DAFoam options
Definition at line 71 of file DATurbulenceModel.H.
Referenced by DASpalartAllmaras::calcResiduals(), DAkOmega::calcResiduals(), DAkEpsilon::calcResiduals(), DASpalartAllmarasFv3::calcResiduals(), DAkOmegaSST::calcResiduals(), and DAkOmegaSSTLM::calcResiduals().
|
protected |
global variable
Definition at line 74 of file DATurbulenceModel.H.
Referenced by DASpalartAllmaras::calcResiduals(), DAkOmega::calcResiduals(), DAkEpsilon::calcResiduals(), DASpalartAllmarasFv3::calcResiduals(), DAkOmegaSST::calcResiduals(), and DAkOmegaSSTLM::calcResiduals().
|
protected |
turbulence viscosity
Definition at line 77 of file DATurbulenceModel.H.
Referenced by DAkOmega::calcResiduals(), DAkEpsilon::calcResiduals(), DAkOmegaSST::calcResiduals(), DAkOmegaSSTLM::calcResiduals(), DATurbulenceModel::correctAlphat(), DASpalartAllmaras::correctNut(), DAkOmega::correctNut(), DAkEpsilon::correctNut(), DASpalartAllmarasFv3::correctNut(), DAkOmegaSST::correctNut(), DAkOmegaSSTLM::correctNut(), DAkEpsilon::DAkEpsilon(), DAkOmega::DAkOmega(), DAkOmegaSST::DAkOmegaSST(), DAkOmegaSSTLM::DAkOmegaSSTLM(), DAkEpsilon::DepsilonEff(), DAkOmegaSSTLM::DgammaIntEff(), DAkOmega::DkEff(), DAkEpsilon::DkEff(), DAkOmegaSST::DkEff(), DAkOmegaSSTLM::DkEff(), DAkOmega::DomegaEff(), DAkOmegaSST::DomegaEff(), DAkOmegaSSTLM::DomegaEff(), DAkOmegaSSTLM::DReThetatEff(), DAkOmega::getFvMatrixFields(), DAkEpsilon::getFvMatrixFields(), DAkOmegaSST::getFvMatrixFields(), DAkOmega::getTurbConvOverProd(), DAkEpsilon::getTurbConvOverProd(), DAkOmegaSST::getTurbConvOverProd(), DAkOmegaSSTLM::getTurbConvOverProd(), DAkOmega::getTurbProdOverDestruct(), DAkEpsilon::getTurbProdOverDestruct(), DAkOmegaSST::getTurbProdOverDestruct(), DAkOmegaSSTLM::getTurbProdOverDestruct(), DATurbulenceModel::nuEff(), DATurbulenceModel::nut(), and DATurbulenceModel::printYPlus().
|
protected |
velocity
Definition at line 80 of file DATurbulenceModel.H.
Referenced by DASpalartAllmarasFv3::calcLduResidualTurb(), DAkOmega::calcResiduals(), DAkEpsilon::calcResiduals(), DASpalartAllmarasFv3::calcResiduals(), DAkOmegaSST::calcResiduals(), DAkOmegaSSTLM::calcResiduals(), DASpalartAllmarasFv3::constructPseudoNuTildaEqn(), DAkEpsilon::DAkEpsilon(), DAkOmega::DAkOmega(), DAkOmegaSST::DAkOmegaSST(), DAkOmegaSSTLM::DAkOmegaSSTLM(), DATurbulenceModel::devRhoReff(), DAkOmegaSSTLM::Flength(), DAkOmegaSSTLM::Fonset(), DAkOmegaSSTLM::Fthetat(), DAkOmega::getFvMatrixFields(), DAkEpsilon::getFvMatrixFields(), DASpalartAllmarasFv3::getFvMatrixFields(), DAkOmegaSST::getFvMatrixFields(), DAkOmega::getTurbConvOverProd(), DAkEpsilon::getTurbConvOverProd(), DASpalartAllmarasFv3::getTurbConvOverProd(), DAkOmegaSST::getTurbConvOverProd(), DAkOmegaSSTLM::getTurbConvOverProd(), DAkOmega::getTurbProdOverDestruct(), DAkEpsilon::getTurbProdOverDestruct(), DASpalartAllmarasFv3::getTurbProdOverDestruct(), DAkOmegaSST::getTurbProdOverDestruct(), DAkOmegaSSTLM::getTurbProdOverDestruct(), DATurbulenceModel::printYPlus(), DAkOmegaSSTLM::ReThetac(), DAkOmegaSSTLM::ReThetat0(), and DASpalartAllmaras::Stilda().
|
protected |
face flux
Definition at line 83 of file DATurbulenceModel.H.
Referenced by DASpalartAllmarasFv3::calcLduResidualTurb(), DAkOmega::calcResiduals(), DAkEpsilon::calcResiduals(), DAkOmegaSST::calcResiduals(), DAkOmegaSSTLM::calcResiduals(), DASpalartAllmarasFv3::constructPseudoNuTildaEqn(), DAkOmega::getFvMatrixFields(), DAkEpsilon::getFvMatrixFields(), and DAkOmegaSST::getFvMatrixFields().
|
protected |
phase field
Definition at line 86 of file DATurbulenceModel.H.
Referenced by DASpalartAllmaras::calcResiduals(), DAkOmega::calcResiduals(), DAkEpsilon::calcResiduals(), DASpalartAllmarasFv3::calcResiduals(), DAkOmegaSST::calcResiduals(), DAkOmegaSSTLM::calcResiduals(), DATurbulenceModel::devRhoReff(), DATurbulenceModel::divDevRhoReff(), DASpalartAllmaras::getFvMatrixFields(), DAkOmega::getFvMatrixFields(), DAkEpsilon::getFvMatrixFields(), DASpalartAllmarasFv3::getFvMatrixFields(), DAkOmegaSST::getFvMatrixFields(), DASpalartAllmaras::getTurbConvOverProd(), DAkOmega::getTurbConvOverProd(), DAkEpsilon::getTurbConvOverProd(), DASpalartAllmarasFv3::getTurbConvOverProd(), DAkOmegaSST::getTurbConvOverProd(), DAkOmegaSSTLM::getTurbConvOverProd(), DASpalartAllmaras::getTurbProdOverDestruct(), DAkOmega::getTurbProdOverDestruct(), DAkEpsilon::getTurbProdOverDestruct(), DASpalartAllmarasFv3::getTurbProdOverDestruct(), DAkOmegaSST::getTurbProdOverDestruct(), DAkOmegaSSTLM::getTurbProdOverDestruct(), DATurbulenceModel::phase(), and DASpalartAllmaras::solveAdjointFP().
|
protected |
phase*phi*density field
Definition at line 89 of file DATurbulenceModel.H.
Referenced by DASpalartAllmaras::calcResiduals(), DAkOmega::calcResiduals(), DAkEpsilon::calcResiduals(), DASpalartAllmarasFv3::calcResiduals(), DAkOmegaSST::calcResiduals(), DAkOmegaSSTLM::calcResiduals(), DASpalartAllmaras::getFvMatrixFields(), DAkOmega::getFvMatrixFields(), DAkEpsilon::getFvMatrixFields(), DASpalartAllmarasFv3::getFvMatrixFields(), DAkOmegaSST::getFvMatrixFields(), DASpalartAllmaras::getTurbConvOverProd(), DAkOmega::getTurbConvOverProd(), DAkEpsilon::getTurbConvOverProd(), DASpalartAllmarasFv3::getTurbConvOverProd(), DAkOmegaSST::getTurbConvOverProd(), DAkOmegaSSTLM::getTurbConvOverProd(), and DASpalartAllmaras::solveAdjointFP().
|
protected |
an uniform rho field with ones
Definition at line 92 of file DATurbulenceModel.H.
Referenced by DATurbulenceModel::rho(), and DATurbulenceModel::rhoDimensions().
|
protected |
whether the turbulence model is incompressible or compressible
Definition at line 95 of file DATurbulenceModel.H.
Referenced by DASpalartAllmaras::addModelResidualCon(), DAkOmega::addModelResidualCon(), DAkEpsilon::addModelResidualCon(), DASpalartAllmarasFv3::addModelResidualCon(), DAkOmegaSST::addModelResidualCon(), DAkOmegaSSTLM::addModelResidualCon(), DATurbulenceModel::alphaEff(), DAkEpsilon::DAkEpsilon(), DAkOmega::DAkOmega(), DAkOmegaSST::DAkOmegaSST(), DAkOmegaSSTLM::DAkOmegaSSTLM(), DASpalartAllmaras::DASpalartAllmaras(), DASpalartAllmarasFv3::DASpalartAllmarasFv3(), DATurbulenceModel::DATurbulenceModel(), DATurbulenceModel::divDevRhoReff(), DATurbulenceModel::getTurbModelType(), DATurbulenceModel::mu(), DATurbulenceModel::nu(), DATurbulenceModel::rho(), and DATurbulenceModel::rhoDimensions().
|
protected |
turbulence model property dict
Definition at line 98 of file DATurbulenceModel.H.
|
protected |
turbulence model parameters dict
Definition at line 101 of file DATurbulenceModel.H.
|
protected |
Lower limit of k.
Definition at line 104 of file DATurbulenceModel.H.
|
protected |
Lower limit of epsilon.
Definition at line 107 of file DATurbulenceModel.H.
|
protected |
Lower limit for omega.
Definition at line 110 of file DATurbulenceModel.H.
|
protected |
Lower limit for nuTilda.
Definition at line 113 of file DATurbulenceModel.H.
|
protected |
Prandtl number.
Definition at line 116 of file DATurbulenceModel.H.
Referenced by DATurbulenceModel::alpha(), and DATurbulenceModel::DATurbulenceModel().
|
protected |
turbulent Prandtl number
Definition at line 119 of file DATurbulenceModel.H.
Referenced by DATurbulenceModel::correctAlphat(), DATurbulenceModel::DATurbulenceModel(), and DATurbulenceModel::prt().