Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
DATurbulenceModel Class Referenceabstract
Inheritance diagram for DATurbulenceModel:
Inheritance graph
[legend]
Collaboration diagram for DATurbulenceModel:
Collaboration graph
[legend]

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 ()=0
 solve the residual equations and update the state More...
 
virtual void getTurbProdTerm (scalarList &prodTerm) const
 return the value of the production term from the Spalart Allmaras 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 > getNut ()
 get the nut field More...
 
tmp< volScalarField > alphaEff ()
 return effective thermal diffusivity More...
 
tmp< volScalarField > nu () const
 get the nu field More...
 
tmp< volScalarField > getAlpha () const
 get alpha field More...
 
tmp< volScalarField > getRho ()
 get the density field More...
 
tmp< volScalarField > getPhase ()
 get the phase field More...
 
scalar getPrt ()
 get the turbulent Prandtl number More...
 
tmp< Foam::volScalarField > getMu () const
 get mu More...
 
bool writeData (Ostream &os) const
 this is a virtual function for regIOobject More...
 
void printYPlus () const
 print the min max and mean yPlus to screen More...
 
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...
 

Static Public Member Functions

static autoPtr< DATurbulenceModelNew (const word modelType, const fvMesh &mesh, const DAOption &daOption)
 

Protected Attributes

const fvMesh & mesh_
 fvMesh More...
 
const DAOptiondaOption_
 DAOption object. More...
 
const dictionary & allOptions_
 all DAFoam options 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...
 
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...
 

Detailed Description

Definition at line 61 of file DATurbulenceModel.H.

Constructor & Destructor Documentation

◆ DATurbulenceModel()

DATurbulenceModel ( const word  modelType,
const fvMesh &  mesh,
const DAOption daOption 
)

◆ ~DATurbulenceModel()

virtual ~DATurbulenceModel ( )
inlinevirtual

Definition at line 178 of file DATurbulenceModel.H.

Member Function Documentation

◆ TypeName()

TypeName ( "DATurbulenceModel"  )

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
DATurbulenceModel  ,
dictionary  ,
(const word modelType, const fvMesh &mesh, const DAOption &daOption ,
(modelType, mesh, daOption  
)

◆ New()

autoPtr< DATurbulenceModel > New ( const word  modelType,
const fvMesh &  mesh,
const DAOption daOption 
)
static

Definition at line 164 of file DATurbulenceModel.C.

References daOption(), and mesh.

Here is the call graph for this function:

◆ correctWallDist()

void correctWallDist ( )

update wall distance for d_. Note: y_ will be automatically updated in mesh_ object

Definition at line 358 of file DATurbulenceModel.C.

References DATurbulenceModel::correctBoundaryConditions().

Here is the call graph for this function:

◆ correctAlphat()

void correctAlphat ( )

◆ correctNut()

virtual void correctNut ( )
pure virtual

◆ correctModelStates()

virtual void correctModelStates ( wordList &  modelStates) const
pure virtual

update the turbulence state for DAStateInfo::regStates_

Implemented in DAkOmegaSSTLM, DAkOmegaSSTFIML, DAkOmegaSSTFieldInversion, DAkOmegaSST, DAkOmegaFieldInversionOmega, DASpalartAllmarasFv3FieldInversion, DAkOmega, DASpalartAllmarasFv3, DAkEpsilon, DASpalartAllmaras, and DADummyTurbulenceModel.

Referenced by DAModel::correctModelStates().

Here is the caller graph for this function:

◆ correctBoundaryConditions()

virtual void correctBoundaryConditions ( )
pure virtual

◆ updateIntermediateVariables()

virtual void updateIntermediateVariables ( )
pure virtual

update any intermediate variables that are dependent on state variables and are used in calcResiduals

Implemented in DAkOmegaSSTLM, DAkOmegaSSTFIML, DAkOmegaSSTFieldInversion, DAkOmegaSST, DAkOmegaFieldInversionOmega, DASpalartAllmarasFv3FieldInversion, DAkOmega, DASpalartAllmarasFv3, DAkEpsilon, DASpalartAllmaras, and DADummyTurbulenceModel.

Referenced by DAModel::updateIntermediateVariables().

Here is the caller graph for this function:

◆ correctStateResidualModelCon()

virtual void correctStateResidualModelCon ( List< List< word >> &  stateCon) const
pure virtual

update the original variable connectivity for the adjoint state residuals in stateCon

Implemented in DAkOmegaSSTLM, DAkOmegaSSTFIML, DAkOmegaSSTFieldInversion, DAkOmegaSST, DAkOmegaFieldInversionOmega, DASpalartAllmarasFv3FieldInversion, DAkOmega, DASpalartAllmarasFv3, DAkEpsilon, DASpalartAllmaras, and DADummyTurbulenceModel.

Referenced by DAModel::correctStateResidualModelCon().

Here is the caller graph for this function:

◆ addModelResidualCon()

virtual void addModelResidualCon ( HashTable< List< List< word >>> &  allCon) const
pure virtual

◆ calcResiduals()

virtual void calcResiduals ( const dictionary &  options)
pure virtual

◆ correct()

virtual void correct ( )
pure virtual

◆ getTurbProdTerm()

void getTurbProdTerm ( scalarList &  prodTerm) const
virtual

return the value of the production term from the Spalart Allmaras model

Reimplemented in DASpalartAllmarasFv3FieldInversion.

Definition at line 483 of file DATurbulenceModel.C.

Referenced by DAModel::getTurbProdTerm().

Here is the caller graph for this function:

◆ devRhoReff()

tmp< volSymmTensorField > devRhoReff ( ) const

◆ divDevRhoReff()

tmp< fvVectorMatrix > divDevRhoReff ( volVectorField &  U)

divDev terms

Definition at line 324 of file DATurbulenceModel.C.

References divScheme, DATurbulenceModel::nuEff(), DATurbulenceModel::phase_, rho, T, and U.

Referenced by DAResidualRhoSimpleCFoam::calcResiduals(), DAResidualTurboFoam::calcResiduals(), DAResidualRhoSimpleFoam::calcResiduals(), and DATurbulenceModel::divDevReff().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ divDevReff()

tmp< fvVectorMatrix > divDevReff ( volVectorField &  U)

divDev terms

Definition at line 347 of file DATurbulenceModel.C.

References DATurbulenceModel::divDevRhoReff(), and U.

Referenced by DAResidualPimpleFoam::calcResiduals(), DAResidualPisoFoam::calcResiduals(), DAResidualSimpleFoam::calcResiduals(), and DAResidualSimpleTFoam::calcResiduals().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ nuEff()

tmp< volScalarField > nuEff ( ) const

return effective viscosity

Definition at line 223 of file DATurbulenceModel.C.

References DATurbulenceModel::nu(), and DATurbulenceModel::nut_.

Referenced by DATurbulenceModel::devRhoReff(), DATurbulenceModel::divDevRhoReff(), and DATurbulenceModel::printYPlus().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getNut()

tmp<volScalarField> getNut ( )
inline

get the nut field

Definition at line 230 of file DATurbulenceModel.H.

References DATurbulenceModel::nut_.

Referenced by DAResidualSimpleTFoam::updateIntermediateVariables().

Here is the caller graph for this function:

◆ alphaEff()

tmp< volScalarField > alphaEff ( )

return effective thermal diffusivity

Definition at line 236 of file DATurbulenceModel.C.

References alphat, DATurbulenceModel::getAlpha(), and DATurbulenceModel::mesh_.

Referenced by DAObjFuncWallHeatFlux::calcObjFunc(), DASolver::getThermal(), and DASolver::setThermal().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ nu()

tmp< volScalarField > nu ( ) const

◆ getAlpha()

tmp< volScalarField > getAlpha ( ) const

get alpha field

Definition at line 280 of file DATurbulenceModel.C.

References DATurbulenceModel::nu(), and DATurbulenceModel::Pr_.

Referenced by DATurbulenceModel::alphaEff().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getRho()

tmp<volScalarField> getRho ( )
inline

get the density field

Definition at line 245 of file DATurbulenceModel.H.

◆ getPhase()

tmp<volScalarField> getPhase ( )
inline

get the phase field

Definition at line 251 of file DATurbulenceModel.H.

References DATurbulenceModel::phase_.

◆ getPrt()

scalar getPrt ( )
inline

get the turbulent Prandtl number

Definition at line 257 of file DATurbulenceModel.H.

References DATurbulenceModel::Prt_.

◆ getMu()

tmp< Foam::volScalarField > getMu ( ) const

get mu

Definition at line 289 of file DATurbulenceModel.C.

References DATurbulenceModel::nut_.

◆ writeData()

bool writeData ( Ostream &  os) const

this is a virtual function for regIOobject

Definition at line 200 of file DATurbulenceModel.C.

◆ printYPlus()

void printYPlus ( ) const

print the min max and mean yPlus to screen

Definition at line 385 of file DATurbulenceModel.C.

References forAll(), DATurbulenceModel::mesh_, DATurbulenceModel::nu(), DATurbulenceModel::nuEff(), DATurbulenceModel::nut_, and DATurbulenceModel::U_.

Here is the call graph for this function:

◆ isPrintTime()

label isPrintTime ( const Time &  runTime,
const label  printInterval 
) const

◆ invTranProdNuTildaEqn()

void invTranProdNuTildaEqn ( const volScalarField &  mySource,
volScalarField &  pseudoNuTilda 
)
virtual

Inverse transpose product, M_nuTilda^(-T)

Reimplemented in DASpalartAllmarasFv3.

Definition at line 495 of file DATurbulenceModel.C.

◆ constructPseudoNuTildaEqn()

void constructPseudoNuTildaEqn ( )
virtual

Reimplemented in DASpalartAllmarasFv3.

Definition at line 521 of file DATurbulenceModel.C.

◆ rhsSolvePseudoNuTildaEqn()

void rhsSolvePseudoNuTildaEqn ( const volScalarField &  nuTildaSource)
virtual

Reimplemented in DASpalartAllmarasFv3.

Definition at line 533 of file DATurbulenceModel.C.

◆ calcLduResidualTurb()

void calcLduResidualTurb ( volScalarField &  nuTildaRes)
virtual

calculate the turbulence residual using LDU matrix

Reimplemented in DASpalartAllmarasFv3.

Definition at line 509 of file DATurbulenceModel.C.

Member Data Documentation

◆ mesh_

const fvMesh& mesh_
protected

fvMesh

Definition at line 74 of file DATurbulenceModel.H.

Referenced by DASpalartAllmaras::addModelResidualCon(), DAkEpsilon::addModelResidualCon(), DASpalartAllmarasFv3::addModelResidualCon(), DAkOmega::addModelResidualCon(), DASpalartAllmarasFv3FieldInversion::addModelResidualCon(), DAkOmegaFieldInversionOmega::addModelResidualCon(), DAkOmegaSST::addModelResidualCon(), DAkOmegaSSTFieldInversion::addModelResidualCon(), DAkOmegaSSTFIML::addModelResidualCon(), DAkOmegaSSTLM::addModelResidualCon(), DATurbulenceModel::alphaEff(), DAkOmegaSSTFIML::calcBetaField(), DASpalartAllmaras::calcResiduals(), DAkEpsilon::calcResiduals(), DASpalartAllmarasFv3::calcResiduals(), DAkOmega::calcResiduals(), DASpalartAllmarasFv3FieldInversion::calcResiduals(), DAkOmegaFieldInversionOmega::calcResiduals(), DAkOmegaSST::calcResiduals(), DAkOmegaSSTFieldInversion::calcResiduals(), DAkOmegaSSTFIML::calcResiduals(), DAkOmegaSSTLM::calcResiduals(), DATurbulenceModel::correctAlphat(), DAkOmegaSST::correctNut(), DAkOmegaSSTFieldInversion::correctNut(), DAkOmegaSSTFIML::correctNut(), DAkOmegaSSTLM::correctNut(), DATurbulenceModel::DATurbulenceModel(), DATurbulenceModel::devRhoReff(), DAkOmegaSSTLM::Flength(), DATurbulenceModel::printYPlus(), DAkOmegaSSTLM::ReThetac(), DAkOmegaSSTLM::ReThetat0(), DAkEpsilon::saveEpsilonNearWall(), DAkOmega::saveOmegaNearWall(), DAkOmegaFieldInversionOmega::saveOmegaNearWall(), DAkOmegaSST::saveOmegaNearWall(), DAkOmegaSSTFieldInversion::saveOmegaNearWall(), DAkOmegaSSTFIML::saveOmegaNearWall(), DAkOmegaSSTLM::saveOmegaNearWall(), DAkEpsilon::setEpsilonNearWall(), DAkOmega::setOmegaNearWall(), DAkOmegaFieldInversionOmega::setOmegaNearWall(), DAkOmegaSST::setOmegaNearWall(), DAkOmegaSSTFieldInversion::setOmegaNearWall(), DAkOmegaSSTFIML::setOmegaNearWall(), and DAkOmegaSSTLM::setOmegaNearWall().

◆ daOption_

const DAOption& daOption_
protected

DAOption object.

Definition at line 77 of file DATurbulenceModel.H.

◆ allOptions_

const dictionary& allOptions_
protected

◆ nut_

volScalarField& nut_
protected

◆ U_

volVectorField& U_
protected

◆ phi_

surfaceScalarField& phi_
protected

◆ phase_

volScalarField phase_
protected

◆ phaseRhoPhi_

surfaceScalarField& phaseRhoPhi_
protected

◆ turbDict_

IOdictionary turbDict_
protected

turbulence model property dict

Definition at line 124 of file DATurbulenceModel.H.

◆ coeffDict_

dictionary coeffDict_
protected

turbulence model parameters dict

Definition at line 127 of file DATurbulenceModel.H.

◆ kMin_

dimensionedScalar kMin_
protected

Lower limit of k.

Definition at line 130 of file DATurbulenceModel.H.

◆ epsilonMin_

dimensionedScalar epsilonMin_
protected

Lower limit of epsilon.

Definition at line 133 of file DATurbulenceModel.H.

◆ omegaMin_

dimensionedScalar omegaMin_
protected

Lower limit for omega.

Definition at line 136 of file DATurbulenceModel.H.

◆ nuTildaMin_

dimensionedScalar nuTildaMin_
protected

Lower limit for nuTilda.

Definition at line 139 of file DATurbulenceModel.H.

◆ Pr_

scalar Pr_
protected

Prandtl number.

Definition at line 142 of file DATurbulenceModel.H.

Referenced by DATurbulenceModel::DATurbulenceModel(), and DATurbulenceModel::getAlpha().

◆ Prt_

scalar Prt_ = -9999.0
protected

The documentation for this class was generated from the following files: