unityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel > Class Template Reference

Eddy-diffusivity based energy gradient heat flux model for RAS or LES of turbulent flow. Specie fluxes are computed assuming a unity turbulent Lewis number. More...

Inheritance diagram for unityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >:
Collaboration diagram for unityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >:

Public Types

typedef TurbulenceThermophysicalTransportModel::alphaField alphaField
 
typedef TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel
 
typedef TurbulenceThermophysicalTransportModel::thermoModel thermoModel
 

Public Member Functions

 TypeName ("unityLewisEddyDiffusivity")
 Runtime type information. More...
 
 unityLewisEddyDiffusivity (const momentumTransportModel &momentumTransport, const thermoModel &thermo, const bool allowDefaultPrt=false)
 Construct from a momentum transport model and a thermo model. More...
 
 unityLewisEddyDiffusivity (const word &type, const momentumTransportModel &momentumTransport, const thermoModel &thermo, const bool allowDefaultPrt=false)
 Construct from a type name, a momentum transport model and a thermo. More...
 
virtual ~unityLewisEddyDiffusivity ()
 Destructor. More...
 
virtual bool read ()
 Read thermophysicalTransport dictionary. More...
 
virtual tmp< volScalarFieldalphat () const
 Turbulent thermal diffusivity for enthalpy [kg/m/s]. More...
 
virtual tmp< scalarFieldalphat (const label patchi) const
 Turbulent thermal diffusivity for enthalpy for a patch [kg/m/s]. More...
 
virtual tmp< volScalarFieldkappaEff () const
 Effective thermal turbulent conductivity. More...
 
virtual tmp< scalarFieldkappaEff (const label patchi) const
 Effective thermal turbulent conductivity. More...
 
virtual tmp< scalarFieldalphaEff (const label patchi) const
 Effective thermal turbulent diffusivity. More...
 
virtual tmp< volScalarFieldDEff (const volScalarField &Yi) const
 Effective mass diffusion coefficient. More...
 
virtual tmp< scalarFieldDEff (const volScalarField &Yi, const label patchi) const
 Effective mass diffusion coefficient. More...
 
virtual tmp< surfaceScalarFieldq () const
 Return the heat flux [W/m^2]. More...
 
virtual tmp< scalarFieldq (const label patchi) const
 Return the patch heat flux [W/m^2]. More...
 
virtual tmp< fvScalarMatrixdivq (volScalarField &he) const
 Return the source term for the energy equation. More...
 
virtual tmp< surfaceScalarFieldj (const volScalarField &Yi) const
 Return the specie flux for the given specie mass-fraction [kg/m^2/s]. More...
 
virtual tmp< scalarFieldj (const volScalarField &Yi, const label patchi) const
 Return the specie flux. More...
 
virtual tmp< fvScalarMatrixdivj (volScalarField &Yi) const
 Return the source term for the given specie mass-fraction equation. More...
 
virtual void predict ()
 Correct the unityLewisEddyDiffusivity viscosity. More...
 

Protected Member Functions

virtual void correctAlphat ()
 
tmp< volScalarFieldalphaEff () const
 

Protected Attributes

dimensionedScalar Prt_
 Turbulent Prandtl number []. More...
 
volScalarField alphat_
 Turbulent thermal diffusivity of enthalpy [kg/m/s]. More...
 

Detailed Description

template<class TurbulenceThermophysicalTransportModel>
class Foam::turbulenceThermophysicalTransportModels::unityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >

Eddy-diffusivity based energy gradient heat flux model for RAS or LES of turbulent flow. Specie fluxes are computed assuming a unity turbulent Lewis number.

Usage
LES
{
    model           unityLewisEddyDiffusivity;
    Prt             0.85;
}
Source files

Definition at line 60 of file unityLewisEddyDiffusivity.H.

Member Typedef Documentation

◆ alphaField

typedef TurbulenceThermophysicalTransportModel::alphaField alphaField

Definition at line 98 of file unityLewisEddyDiffusivity.H.

◆ momentumTransportModel

typedef TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel

Definition at line 102 of file unityLewisEddyDiffusivity.H.

◆ thermoModel

typedef TurbulenceThermophysicalTransportModel::thermoModel thermoModel

Definition at line 105 of file unityLewisEddyDiffusivity.H.

Constructor & Destructor Documentation

◆ unityLewisEddyDiffusivity() [1/2]

unityLewisEddyDiffusivity ( const momentumTransportModel momentumTransport,
const thermoModel thermo,
const bool  allowDefaultPrt = false 
)

Construct from a momentum transport model and a thermo model.

Definition at line 53 of file unityLewisEddyDiffusivity.C.

◆ unityLewisEddyDiffusivity() [2/2]

unityLewisEddyDiffusivity ( const word type,
const momentumTransportModel momentumTransport,
const thermoModel thermo,
const bool  allowDefaultPrt = false 
)

Construct from a type name, a momentum transport model and a thermo.

model, and whether to default the turbulent Prandtl number to one if it is not specified

Definition at line 72 of file unityLewisEddyDiffusivity.C.

◆ ~unityLewisEddyDiffusivity()

virtual ~unityLewisEddyDiffusivity ( )
inlinevirtual

Destructor.

Definition at line 135 of file unityLewisEddyDiffusivity.H.

Member Function Documentation

◆ correctAlphat()

void correctAlphat
protectedvirtual

Definition at line 40 of file unityLewisEddyDiffusivity.C.

◆ alphaEff() [1/2]

◆ TypeName()

TypeName ( "unityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Read thermophysicalTransport dictionary.

Reimplemented in nonUnityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >.

Definition at line 117 of file unityLewisEddyDiffusivity.C.

References Foam::blockMeshTools::read().

Here is the call graph for this function:

◆ alphat() [1/2]

virtual tmp<volScalarField> alphat ( ) const
inlinevirtual

◆ alphat() [2/2]

virtual tmp<scalarField> alphat ( const label  patchi) const
inlinevirtual

Turbulent thermal diffusivity for enthalpy for a patch [kg/m/s].

Definition at line 151 of file unityLewisEddyDiffusivity.H.

References unityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >::alphat(), and patchi.

Here is the call graph for this function:

◆ kappaEff() [1/2]

virtual tmp<volScalarField> kappaEff ( ) const
inlinevirtual

Effective thermal turbulent conductivity.

of mixture [W/m/K]

Definition at line 158 of file unityLewisEddyDiffusivity.H.

◆ kappaEff() [2/2]

virtual tmp<scalarField> kappaEff ( const label  patchi) const
inlinevirtual

Effective thermal turbulent conductivity.

of mixture for patch [W/m/K]

Definition at line 165 of file unityLewisEddyDiffusivity.H.

◆ alphaEff() [2/2]

virtual tmp<scalarField> alphaEff ( const label  patchi) const
inlinevirtual

Effective thermal turbulent diffusivity.

of mixture for a patch [kg/m/s]

Definition at line 174 of file unityLewisEddyDiffusivity.H.

◆ DEff() [1/2]

virtual tmp<volScalarField> DEff ( const volScalarField Yi) const
inlinevirtual

Effective mass diffusion coefficient.

for a given specie mass-fraction [kg/m/s]

Reimplemented in nonUnityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >.

Definition at line 184 of file unityLewisEddyDiffusivity.H.

◆ DEff() [2/2]

virtual tmp<scalarField> DEff ( const volScalarField Yi,
const label  patchi 
) const
inlinevirtual

Effective mass diffusion coefficient.

for a given specie mass-fraction for patch [kg/m/s]

Reimplemented in nonUnityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >.

Definition at line 195 of file unityLewisEddyDiffusivity.H.

◆ q() [1/2]

tmp< surfaceScalarField > q
virtual

Return the heat flux [W/m^2].

Reimplemented in nonUnityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >.

Definition at line 134 of file unityLewisEddyDiffusivity.C.

References Foam::constant::atomic::group, IOobject::groupName(), and GeometricField< Type, GeoMesh, PrimitiveField >::New().

Here is the call graph for this function:

◆ q() [2/2]

tmp< scalarField > q ( const label  patchi) const
virtual

Return the patch heat flux [W/m^2].

Reimplemented in nonUnityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >.

Definition at line 151 of file unityLewisEddyDiffusivity.C.

References alpha(), patchi, and thermo.

Here is the call graph for this function:

◆ divq()

tmp< fvScalarMatrix > divq ( volScalarField he) const
virtual

Return the source term for the energy equation.

Reimplemented in nonUnityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >.

Definition at line 167 of file unityLewisEddyDiffusivity.C.

References alpha(), he(), and Foam::fvm::laplacian().

Here is the call graph for this function:

◆ j() [1/2]

tmp< surfaceScalarField > j ( const volScalarField Yi) const
virtual

Return the specie flux for the given specie mass-fraction [kg/m^2/s].

Definition at line 178 of file unityLewisEddyDiffusivity.C.

References IOobject::groupName(), Foam::fvc::interpolate(), IOobject::name(), GeometricField< Type, GeoMesh, PrimitiveField >::New(), and Foam::fvc::snGrad().

Here is the call graph for this function:

◆ j() [2/2]

tmp< scalarField > j ( const volScalarField Yi,
const label  patchi 
) const
virtual

Return the specie flux.

for the given specie mass-fraction for patch [kg/m^2/s]

Definition at line 197 of file unityLewisEddyDiffusivity.C.

References alpha(), GeometricField< Type, GeoMesh, PrimitiveField >::boundaryField(), and patchi.

Here is the call graph for this function:

◆ divj()

tmp< fvScalarMatrix > divj ( volScalarField Yi) const
virtual

Return the source term for the given specie mass-fraction equation.

Definition at line 214 of file unityLewisEddyDiffusivity.C.

◆ predict()

void predict
virtual

Correct the unityLewisEddyDiffusivity viscosity.

Definition at line 224 of file unityLewisEddyDiffusivity.C.

Member Data Documentation

◆ Prt_

dimensionedScalar Prt_
protected

Turbulent Prandtl number [].

Definition at line 72 of file unityLewisEddyDiffusivity.H.

◆ alphat_

volScalarField alphat_
protected

Turbulent thermal diffusivity of enthalpy [kg/m/s].

Definition at line 77 of file unityLewisEddyDiffusivity.H.

Referenced by unityLewisEddyDiffusivity< TurbulenceThermophysicalTransportModel >::alphat().


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