Go to the documentation of this file.
25 #ifndef DATurbulenceModel_H
26 #define DATurbulenceModel_H
28 #include "runTimeSelectionTables.H"
30 #include "nearWallDist.H"
35 #include "nutkWallFunctionFvPatchScalarField.H"
36 #include "wallFvPatch.H"
38 #ifdef IncompressibleFlow
39 #include "singlePhaseTransportModel.H"
40 #include "turbulentTransportModel.H"
45 #ifdef CompressibleFlow
46 #include "fluidThermo.H"
47 #include "turbulentFluidThermoModel.H"
97 #ifdef IncompressibleFlow
101 const singlePhaseTransportModel& laminarTransport_;
105 const incompressible::turbulenceModel& turbulence_;
110 #ifdef CompressibleFlow
114 const fluidThermo& thermo_;
118 const compressible::turbulenceModel& turbulence_;
120 volScalarField& rho_;
156 (
const word modelType,
165 const word modelType,
172 static autoPtr<DATurbulenceModel>
New(
173 const word modelType,
227 tmp<volScalarField>
nuEff()
const;
239 tmp<volScalarField>
nu()
const;
242 tmp<volScalarField>
getAlpha()
const;
266 FatalErrorIn(
"getPrt") << exit(FatalError);
271 #ifdef CompressibleFlow
272 const fluidThermo& getThermo()
const;
276 tmp<Foam::volScalarField>
getMu()
const;
286 const label printInterval)
const;
290 const volScalarField& mySource,
tmp< volScalarField > nuEff() const
return effective viscosity
virtual void constructPseudoNuTildaEqn()
surfaceScalarField & phi_
face flux
tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U)
divDev terms
tmp< fvVectorMatrix > divDevReff(volVectorField &U)
divDev terms
dimensionedScalar kMin_
Lower limit of k.
virtual void updateIntermediateVariables()=0
update any intermediate variables that are dependent on state variables and are used in calcResiduals
virtual void correctModelStates(wordList &modelStates) const =0
update the turbulence state for DAStateInfo::regStates_
TypeName("DATurbulenceModel")
DAOption daOption(mesh, pyOptions_)
volScalarField & pseudoNuTilda
dimensionedScalar nuTildaMin_
Lower limit for nuTilda.
scalar getPrt()
get the turbulent Prandtl number
virtual void correctNut()=0
update nut based on other turbulence variables and update the BCs
surfaceScalarField & phaseRhoPhi_
phase*phi*density field
const dictionary & allOptions_
all DAFoam options
static autoPtr< DATurbulenceModel > New(const word modelType, const fvMesh &mesh, const DAOption &daOption)
scalar Prt_
turbulent Prandtl number
tmp< volScalarField > getPhase()
get the phase field
virtual void correctBoundaryConditions()=0
update turbulence variable boundary values
virtual void calcLduResidualTurb(volScalarField &nuTildaRes)
calculate the turbulence residual using LDU matrix
declareRunTimeSelectionTable(autoPtr, DATurbulenceModel, dictionary,(const word modelType, const fvMesh &mesh, const DAOption &daOption),(modelType, mesh, daOption))
tmp< Foam::volScalarField > getMu() const
get mu
virtual void correct()=0
solve the residual equations and update the state
virtual void calcResiduals(const dictionary &options)=0
compute the turbulence residuals
tmp< volScalarField > getAlpha() const
get alpha field
tmp< volScalarField > getNut()
get the nut field
void correctWallDist()
update wall distance for d_. Note: y_ will be automatically updated in mesh_ object
volVectorField & U_
velocity
IOdictionary turbDict_
turbulence model property dict
volScalarField & nut_
turbulence viscosity
virtual void invTranProdNuTildaEqn(const volScalarField &mySource, volScalarField &pseudoNuTilda)
Inverse transpose product, M_nuTilda^(-T)
tmp< volSymmTensorField > devRhoReff() const
dev terms
virtual ~DATurbulenceModel()
tmp< volScalarField > getRho()
get the density field
virtual void rhsSolvePseudoNuTildaEqn(const volScalarField &nuTildaSource)
tmp< volScalarField > nu() const
get the nu field
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const =0
add the model residual connectivity to stateCon
void correctAlphat()
update alphat
const DAOption & daOption_
DAOption object.
virtual void getTurbProdTerm(scalarList &prodTerm) const
return the value of the production term from the Spalart Allmaras model
const fvMesh & mesh_
fvMesh
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const =0
update the original variable connectivity for the adjoint state residuals in stateCon
label isPrintTime(const Time &runTime, const label printInterval) const
tmp< volScalarField > alphaEff()
return effective thermal diffusivity
dimensionedScalar omegaMin_
Lower limit for omega.
void printYPlus() const
print the min max and mean yPlus to screen
dictionary coeffDict_
turbulence model parameters dict
scalar Pr_
Prandtl number.
volScalarField phase_
phase field
dimensionedScalar epsilonMin_
Lower limit of epsilon.
bool writeData(Ostream &os) const
this is a virtual function for regIOobject