List of all members
DAkEpsilon Class Reference
Inheritance diagram for DAkEpsilon:
Inheritance graph
[legend]
Collaboration diagram for DAkEpsilon:
Collaboration graph
[legend]

Protected Member Functions

SST functions
tmp< fvScalarMatrix > kSource () const
 
tmp< fvScalarMatrix > epsilonSource () const
 
tmp< volScalarField > DkEff () const
 
tmp< volScalarField > DepsilonEff () const
 
virtual tmp< volScalarField > k () const
 
virtual tmp< volScalarField > epsilon () const
 

Protected Attributes

SST parameters
dimensionedScalar Cmu_
 
dimensionedScalar C1_
 
dimensionedScalar C2_
 
dimensionedScalar C3_
 
dimensionedScalar sigmak_
 
dimensionedScalar sigmaEps_
 
- Protected Attributes inherited from DATurbulenceModel
const fvMesh & mesh_
 fvMesh More...
 
const DAOptiondaOption_
 DAOption object. More...
 
const dictionary & allOptions_
 all DAFoam options More...
 
DAGlobalVardaGlobalVar_
 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...
 

Augmented variables for adjoint residuals

volScalarField & epsilon_
 
volScalarField epsilonRes_
 
volScalarField & k_
 
volScalarField kRes_
 
autoPtr< volScalarField::Internal > GPtr_
 
volScalarField betaFIK_
 beta field for field inversion More...
 
volScalarField betaFIEpsilon_
 
scalarList epsilonNearWall_
 
label solveTurbState_ = 0
 whether to solve for turb states More...
 
 TypeName ("kEpsilon")
 
 DAkEpsilon (const word modelType, const fvMesh &mesh, const DAOption &daOption)
 
virtual ~DAkEpsilon ()
 
virtual void correctModelStates (wordList &modelStates) const
 update the turbulence state for DAStateInfo::regStates_ More...
 
virtual void correctNut ()
 update nut based on other turbulence variables and update the BCs More...
 
virtual void correctBoundaryConditions ()
 update turbulence variable boundary values More...
 
virtual void updateIntermediateVariables ()
 update any intermediate variables that are dependent on state variables and are used in calcResiduals More...
 
virtual void correctStateResidualModelCon (List< List< word >> &stateCon) const
 update the original variable connectivity for the adjoint state residuals in stateCon More...
 
virtual void addModelResidualCon (HashTable< List< List< word >>> &allCon) const
 add the model residual connectivity to stateCon More...
 
virtual void calcResiduals (const dictionary &options)
 compute the turbulence residuals More...
 
virtual void correct (label printToScreen)
 solve the residual equations and update the state More...
 
void saveEpsilonNearWall ()
 save near wall epsilon values to epsilonNearWall_ More...
 
void setEpsilonNearWall ()
 set epsilonNearWall_ to near wall epsilon values More...
 
void correctEpsilonBoundaryConditions ()
 specially treatment to correct epsilon BC 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 getTurbProdOverDestruct (volScalarField &PoD) const
 return the value of the 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...
 

Additional Inherited Members

- Public Member Functions inherited from DATurbulenceModel
 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 getTurbProdTerm (volScalarField &prodTerm) const
 return the value of the 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 solveAdjointFP (const word varName, const scalarList &rhs, scalarList &dPsi)
 solve the fvMatrixT field with given rhs and solution More...
 
- Static Public Member Functions inherited from DATurbulenceModel
static autoPtr< DATurbulenceModelNew (const word modelType, const fvMesh &mesh, const DAOption &daOption)
 

Detailed Description

Definition at line 48 of file DAkEpsilon.H.

Constructor & Destructor Documentation

◆ DAkEpsilon()

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

◆ ~DAkEpsilon()

virtual ~DAkEpsilon ( )
inlinevirtual

Definition at line 120 of file DAkEpsilon.H.

Member Function Documentation

◆ kSource()

tmp< fvScalarMatrix > kSource ( ) const
protected

Definition at line 160 of file DAkEpsilon.C.

References DAkEpsilon::k_, and DATurbulenceModel::rhoDimensions().

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::getFvMatrixFields().

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

◆ epsilonSource()

tmp< fvScalarMatrix > epsilonSource ( ) const
protected

Definition at line 169 of file DAkEpsilon.C.

References DAkEpsilon::epsilon_, and DATurbulenceModel::rhoDimensions().

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::getFvMatrixFields().

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

◆ DkEff()

tmp< volScalarField > DkEff ( ) const
protected

Definition at line 178 of file DAkEpsilon.C.

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

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::getFvMatrixFields().

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

◆ DepsilonEff()

tmp< volScalarField > DepsilonEff ( ) const
protected

Definition at line 186 of file DAkEpsilon.C.

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

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::getFvMatrixFields().

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

◆ k()

virtual tmp<volScalarField> k ( ) const
inlineprotectedvirtual

Definition at line 72 of file DAkEpsilon.H.

References DAkEpsilon::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlineprotectedvirtual

Definition at line 78 of file DAkEpsilon.H.

References DAkEpsilon::epsilon_.

◆ TypeName()

TypeName ( "kEpsilon"  )

◆ correctModelStates()

void correctModelStates ( wordList &  modelStates) const
virtual

update the turbulence state for DAStateInfo::regStates_

Implements DATurbulenceModel.

Definition at line 195 of file DAkEpsilon.C.

References forAll().

Here is the call graph for this function:

◆ correctNut()

void correctNut ( )
virtual

update nut based on other turbulence variables and update the BCs

Implements DATurbulenceModel.

Definition at line 231 of file DAkEpsilon.C.

References DAkEpsilon::Cmu_, DATurbulenceModel::correctAlphat(), DAkEpsilon::epsilon_, DAkEpsilon::k_, and DATurbulenceModel::nut_.

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::updateIntermediateVariables().

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

◆ correctBoundaryConditions()

void correctBoundaryConditions ( )
virtual

update turbulence variable boundary values

Implements DATurbulenceModel.

Definition at line 249 of file DAkEpsilon.C.

References DAkEpsilon::k_.

◆ updateIntermediateVariables()

void updateIntermediateVariables ( )
virtual

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

Implements DATurbulenceModel.

Definition at line 338 of file DAkEpsilon.C.

References DAkEpsilon::correctNut().

Here is the call graph for this function:

◆ correctStateResidualModelCon()

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

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

Implements DATurbulenceModel.

Definition at line 349 of file DAkEpsilon.C.

References forAll().

Here is the call graph for this function:

◆ addModelResidualCon()

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

add the model residual connectivity to stateCon

Implements DATurbulenceModel.

Definition at line 400 of file DAkEpsilon.C.

References DATurbulenceModel::mesh_, and DATurbulenceModel::turbModelType_.

◆ calcResiduals()

void calcResiduals ( const dictionary &  options)
virtual

◆ correct()

void correct ( label  printToScreen)
virtual

solve the residual equations and update the state

Implements DATurbulenceModel.

Definition at line 493 of file DAkEpsilon.C.

References DAkEpsilon::calcResiduals(), and DAkEpsilon::solveTurbState_.

Here is the call graph for this function:

◆ saveEpsilonNearWall()

void saveEpsilonNearWall ( )

save near wall epsilon values to epsilonNearWall_

Definition at line 288 of file DAkEpsilon.C.

References DAkEpsilon::epsilon_, DAkEpsilon::epsilonNearWall_, forAll(), and DATurbulenceModel::mesh_.

Referenced by DAkEpsilon::correctEpsilonBoundaryConditions().

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

◆ setEpsilonNearWall()

void setEpsilonNearWall ( )

set epsilonNearWall_ to near wall epsilon values

Definition at line 312 of file DAkEpsilon.C.

References DAkEpsilon::epsilon_, DAkEpsilon::epsilonNearWall_, forAll(), and DATurbulenceModel::mesh_.

Referenced by DAkEpsilon::calcResiduals(), DAkEpsilon::correctEpsilonBoundaryConditions(), and DAkEpsilon::getFvMatrixFields().

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

◆ correctEpsilonBoundaryConditions()

void correctEpsilonBoundaryConditions ( )

specially treatment to correct epsilon BC

Definition at line 261 of file DAkEpsilon.C.

References DAkEpsilon::epsilon_, DAkEpsilon::saveEpsilonNearWall(), and DAkEpsilon::setEpsilonNearWall().

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::getFvMatrixFields().

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

◆ getFvMatrixFields()

void getFvMatrixFields ( const word  varName,
scalarField &  diag,
scalarField &  upper,
scalarField &  lower 
)
virtual

◆ getTurbProdOverDestruct()

void getTurbProdOverDestruct ( volScalarField &  PoD) const
virtual

return the value of the destruction term from the turbulence model

Reimplemented from DATurbulenceModel.

Definition at line 722 of file DAkEpsilon.C.

References D, DAkEpsilon::epsilon_, forAll(), DAkEpsilon::GPtr_, DATurbulenceModel::nut_, DATurbulenceModel::phase_, DATurbulenceModel::rho(), and DATurbulenceModel::U_.

Here is the call graph for this function:

◆ getTurbConvOverProd()

void getTurbConvOverProd ( volScalarField &  CoP) const
virtual

return the value of the convective over production term from the turbulence model

Reimplemented from DATurbulenceModel.

Definition at line 743 of file DAkEpsilon.C.

References forAll(), DAkEpsilon::GPtr_, DAkEpsilon::k_, DATurbulenceModel::nut_, DATurbulenceModel::phase_, DATurbulenceModel::phaseRhoPhi_, DATurbulenceModel::rho(), and DATurbulenceModel::U_.

Here is the call graph for this function:

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 55 of file DAkEpsilon.H.

Referenced by DAkEpsilon::correctNut().

◆ C1_

dimensionedScalar C1_
protected

Definition at line 56 of file DAkEpsilon.H.

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::getFvMatrixFields().

◆ C2_

dimensionedScalar C2_
protected

Definition at line 57 of file DAkEpsilon.H.

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::getFvMatrixFields().

◆ C3_

dimensionedScalar C3_
protected

Definition at line 58 of file DAkEpsilon.H.

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::getFvMatrixFields().

◆ sigmak_

dimensionedScalar sigmak_
protected

Definition at line 59 of file DAkEpsilon.H.

Referenced by DAkEpsilon::DkEff().

◆ sigmaEps_

dimensionedScalar sigmaEps_
protected

Definition at line 60 of file DAkEpsilon.H.

Referenced by DAkEpsilon::DepsilonEff().

◆ epsilon_

volScalarField& epsilon_
protected

◆ epsilonRes_

volScalarField epsilonRes_
protected

Definition at line 87 of file DAkEpsilon.H.

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::DAkEpsilon().

◆ k_

volScalarField& k_
protected

◆ kRes_

volScalarField kRes_
protected

Definition at line 89 of file DAkEpsilon.H.

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::DAkEpsilon().

◆ GPtr_

autoPtr<volScalarField::Internal> GPtr_
protected

we need to make the G field a class variable and register it to the mesh.db such that the omegaWallFunction BC can find it

Definition at line 94 of file DAkEpsilon.H.

Referenced by DAkEpsilon::calcResiduals(), DAkEpsilon::DAkEpsilon(), DAkEpsilon::getFvMatrixFields(), DAkEpsilon::getTurbConvOverProd(), and DAkEpsilon::getTurbProdOverDestruct().

◆ betaFIK_

volScalarField betaFIK_
protected

beta field for field inversion

Definition at line 97 of file DAkEpsilon.H.

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::getFvMatrixFields().

◆ betaFIEpsilon_

volScalarField betaFIEpsilon_
protected

Definition at line 98 of file DAkEpsilon.H.

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::getFvMatrixFields().

◆ epsilonNearWall_

scalarList epsilonNearWall_
protected

cell-center epsilon values near the wall, this is to fix the issue that the epsilonWallFunction will try to modify epsilon values for the cells near walls this will cause issue for FD-based partial derivatives, so here we basically implement a zeroGradient BC for near wall epsilon.

Definition at line 104 of file DAkEpsilon.H.

Referenced by DAkEpsilon::DAkEpsilon(), DAkEpsilon::saveEpsilonNearWall(), and DAkEpsilon::setEpsilonNearWall().

◆ solveTurbState_

label solveTurbState_ = 0
protected

whether to solve for turb states

Definition at line 107 of file DAkEpsilon.H.

Referenced by DAkEpsilon::calcResiduals(), and DAkEpsilon::correct().


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