unityLewisFourier< BasicThermophysicalTransportModel > Class Template Reference

unityLewisFourier's energy gradient heat flux model for laminar flow. Specie fluxes are computed assuming a unity turbulent Lewis number. More...

Inheritance diagram for unityLewisFourier< BasicThermophysicalTransportModel >:
Collaboration diagram for unityLewisFourier< BasicThermophysicalTransportModel >:

Public Types

typedef BasicThermophysicalTransportModel::alphaField alphaField
 
typedef BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
 
typedef BasicThermophysicalTransportModel::thermoModel thermoModel
 
- Public Types inherited from laminarThermophysicalTransportModel< BasicThermophysicalTransportModel >
typedef BasicThermophysicalTransportModel::alphaField alphaField
 
typedef BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
 
typedef BasicThermophysicalTransportModel::thermoModel thermoModel
 

Public Member Functions

 TypeName ("unityLewisFourier")
 Runtime type information. More...
 
 unityLewisFourier (const momentumTransportModel &momentumTransport, const thermoModel &thermo)
 Construct from a momentum transport model and a thermo model. More...
 
 unityLewisFourier (const word &type, const momentumTransportModel &momentumTransport, const thermoModel &thermo)
 Construct from a type name, a momentum transport model and a thermo. More...
 
virtual ~unityLewisFourier ()
 Destructor. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
virtual bool read ()
 Read thermophysicalTransport dictionary. 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< 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< fvScalarMatrixdivj (volScalarField &Yi) const
 Return the source term for the given specie mass-fraction equation. More...
 
virtual void predict ()
 Correct the unityLewisFourier viscosity. More...
 
- Public Member Functions inherited from laminarThermophysicalTransportModel< BasicThermophysicalTransportModel >
 TypeName ("laminar")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, laminarThermophysicalTransportModel, dictionary,(const momentumTransportModel &momentumTransport, const thermoModel &thermo),(momentumTransport, thermo))
 
 laminarThermophysicalTransportModel (const word &type, const momentumTransportModel &momentumTransport, const thermoModel &thermo)
 Construct from components. More...
 
 laminarThermophysicalTransportModel (const laminarThermophysicalTransportModel &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~laminarThermophysicalTransportModel ()
 Destructor. More...
 
virtual tmp< volScalarFieldkappaEff () const
 Effective thermal turbulent conductivity. More...
 
virtual tmp< scalarFieldkappaEff (const label patchi) const
 Effective thermal turbulent conductivity. More...
 
virtual tmp< volScalarFieldalphaEff () const
 Effective thermal turbulent diffusivity of mixture [kg/m/s]. More...
 
virtual void correct ()
 Solve the thermophysical transport model equations. More...
 
void operator= (const laminarThermophysicalTransportModel &)=delete
 Disallow default bitwise assignment. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from laminarThermophysicalTransportModel< BasicThermophysicalTransportModel >
static autoPtr< laminarThermophysicalTransportModelNew (const momentumTransportModel &momentumTransport, const thermoModel &thermo)
 Return a reference to the selected laminar model. More...
 
- Protected Member Functions inherited from laminarThermophysicalTransportModel< BasicThermophysicalTransportModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 
- Protected Attributes inherited from laminarThermophysicalTransportModel< BasicThermophysicalTransportModel >
dictionary laminarDict_
 laminar coefficients dictionary More...
 
Switch printCoeffs_
 Flag to print the model coeffs at run-time. More...
 
dictionary coeffDict_
 Model coefficients dictionary. More...
 

Detailed Description

template<class BasicThermophysicalTransportModel>
class Foam::laminarThermophysicalTransportModels::unityLewisFourier< BasicThermophysicalTransportModel >

unityLewisFourier's energy gradient heat flux model for laminar flow. Specie fluxes are computed assuming a unity turbulent Lewis number.

Source files

Definition at line 52 of file unityLewisFourier.H.

Member Typedef Documentation

◆ alphaField

typedef BasicThermophysicalTransportModel::alphaField alphaField

Definition at line 63 of file unityLewisFourier.H.

◆ momentumTransportModel

typedef BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel

Definition at line 66 of file unityLewisFourier.H.

◆ thermoModel

typedef BasicThermophysicalTransportModel::thermoModel thermoModel

Definition at line 69 of file unityLewisFourier.H.

Constructor & Destructor Documentation

◆ unityLewisFourier() [1/2]

unityLewisFourier ( const momentumTransportModel momentumTransport,
const thermoModel thermo 
)

Construct from a momentum transport model and a thermo model.

Definition at line 41 of file unityLewisFourier.C.

◆ unityLewisFourier() [2/2]

unityLewisFourier ( const word type,
const momentumTransportModel momentumTransport,
const thermoModel thermo 
)

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

model

Definition at line 57 of file unityLewisFourier.C.

◆ ~unityLewisFourier()

virtual ~unityLewisFourier ( )
inlinevirtual

Destructor.

Definition at line 96 of file unityLewisFourier.H.

Member Function Documentation

◆ TypeName()

TypeName ( "unityLewisFourier< BasicThermophysicalTransportModel >"  )

Runtime type information.

◆ coeffDict()

const dictionary & coeffDict
virtual

Const access to the coefficients dictionary.

Reimplemented from laminarThermophysicalTransportModel< BasicThermophysicalTransportModel >.

Definition at line 77 of file unityLewisFourier.C.

References dictionary::null.

◆ read()

bool read
virtual

Read thermophysicalTransport dictionary.

Reimplemented from laminarThermophysicalTransportModel< BasicThermophysicalTransportModel >.

Definition at line 84 of file unityLewisFourier.C.

◆ 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]

Implements laminarThermophysicalTransportModel< BasicThermophysicalTransportModel >.

Definition at line 110 of file unityLewisFourier.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]

Implements laminarThermophysicalTransportModel< BasicThermophysicalTransportModel >.

Definition at line 121 of file unityLewisFourier.H.

◆ q()

◆ divq()

tmp< fvScalarMatrix > divq ( volScalarField he) const
virtual

Return the source term for the energy equation.

Definition at line 112 of file unityLewisFourier.C.

References alpha(), he(), Foam::constant::electromagnetic::kappa, Foam::fvm::laplacian(), GeometricField< Type, PatchField, GeoMesh >::New(), and thermo.

Here is the call graph for this function:

◆ j()

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 129 of file unityLewisFourier.C.

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

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 149 of file unityLewisFourier.C.

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

Here is the call graph for this function:

◆ predict()

void predict
virtual

Correct the unityLewisFourier viscosity.

Reimplemented from laminarThermophysicalTransportModel< BasicThermophysicalTransportModel >.

Definition at line 157 of file unityLewisFourier.C.


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