kOmegaSST Class Reference
Inheritance diagram for kOmegaSST:
[legend]
Collaboration diagram for kOmegaSST:
[legend]

Public Member Functions

 TypeName ("kOmegaSST")
 Runtime type information. More...
 
 kOmegaSST (const fvMesh &mesh, const solverControl &SolverControl)
 Construct from components. More...
 
virtual ~kOmegaSST ()=default
 Destructor. More...
 
virtual tmp< volScalarField::InternalG ()
 Return the turbulence production term. More...
 
virtual void computeMeanFields ()
 Compute mean fields on the fly. More...
 
virtual bool hasTMVar1 () const
 Bools to identify which turbulent fields are present. More...
 
virtual bool hasTMVar2 () const
 
virtual bool hasNut () const
 
virtual void correctBoundaryConditions (const incompressible::turbulenceModel &turbulence)
 Correct boundary conditions of turbulent fields. More...
 
- Public Member Functions inherited from RASModelVariables
 TypeName ("RASModelVariables")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, RASModelVariables, dictionary,(const fvMesh &mesh, const solverControl &SolverControl),(mesh, SolverControl))
 
 RASModelVariables (const fvMesh &mesh, const solverControl &SolverControl)
 Construct from components. More...
 
 RASModelVariables (const RASModelVariables &rmv)
 Copy constructor. More...
 
autoPtr< RASModelVariablesclone () const
 Clone. More...
 
virtual ~RASModelVariables ()=default
 
const wordTMVar1BaseName () const
 Turbulence field names. More...
 
const wordTMVar2BaseName () const
 
const wordnutBaseName () const
 
virtual bool hasTMVar1 () const
 Bools to identify which turbulent fields are present. More...
 
virtual bool hasTMVar2 () const
 
virtual bool hasNut () const
 
bool hasDist () const
 
const volScalarFieldTMVar1 () const
 Return references to turbulence fields. More...
 
volScalarFieldTMVar1 ()
 
const volScalarFieldTMVar2 () const
 
volScalarFieldTMVar2 ()
 
const volScalarFieldnutRef () const
 
volScalarFieldnutRef ()
 
const volScalarFieldd () const
 
volScalarFieldd ()
 
const volScalarFieldTMVar1Inst () const
 return references to instantaneous turbulence fields More...
 
volScalarFieldTMVar1Inst ()
 
const volScalarFieldTMVar2Inst () const
 
volScalarFieldTMVar2Inst ()
 
const volScalarFieldnutRefInst () const
 
volScalarFieldnutRefInst ()
 
virtual tmp< volScalarFieldnutJacobianVar1 (const singlePhaseTransportModel &laminarTransport) const
 Return nut Jacobian wrt the TM vars. More...
 
virtual tmp< volScalarFieldnutJacobianVar2 (const singlePhaseTransportModel &laminarTransport) const
 
virtual tmp< volScalarField::InternalG ()
 Return the turbulence production term. More...
 
void restoreInitValues ()
 Restore turbulent fields to their initial values. More...
 
void resetMeanFields ()
 Reset mean fields to zero. More...
 
virtual void computeMeanFields ()
 Compute mean fields on the fly. More...
 
tmp< volSymmTensorFielddevReff (const singlePhaseTransportModel &laminarTransport, const volVectorField &U) const
 Return stress tensor based on the mean flow variables. More...
 
virtual void correctBoundaryConditions (const incompressible::turbulenceModel &turbulence)
 correct bounday conditions of turbulent fields More...
 
virtual void transfer (RASModelVariables &rmv)
 Transfer turbulence fields from an another object. More...
 

Protected Member Functions

virtual void allocateMeanFields ()
 
tmp< volScalarField::InternalcomputeG ()
 
- Protected Member Functions inherited from RASModelVariables
virtual void allocateInitValues ()
 
virtual void allocateMeanFields ()
 
refPtr< volScalarFieldcloneRefPtr (const refPtr< volScalarField > &obj) const
 
void copyAndRename (volScalarField &f1, volScalarField &f2)
 
void operator= (const RASModelVariables &)=delete
 No copy assignment. More...
 

Protected Attributes

autoPtr< volScalarField::InternalGMean_
 Average of the production term. More...
 
- Protected Attributes inherited from RASModelVariables
const fvMeshmesh_
 
const solverControlsolverControl_
 
word TMVar1BaseName_
 
word TMVar2BaseName_
 
word nutBaseName_
 
refPtr< volScalarFieldTMVar1Ptr_
 
refPtr< volScalarFieldTMVar2Ptr_
 
refPtr< volScalarFieldnutPtr_
 
refPtr< volScalarFielddistPtr_
 
refPtr< volScalarFieldTMVar1InitPtr_
 
refPtr< volScalarFieldTMVar2InitPtr_
 
refPtr< volScalarFieldnutInitPtr_
 
refPtr< volScalarFieldTMVar1MeanPtr_
 
refPtr< volScalarFieldTMVar2MeanPtr_
 
refPtr< volScalarFieldnutMeanPtr_
 

Additional Inherited Members

- Static Public Member Functions inherited from RASModelVariables
static autoPtr< RASModelVariablesNew (const fvMesh &mesh, const solverControl &SolverControl)
 Return a reference to the selected turbulence model. More...
 

Detailed Description

Source files

Definition at line 56 of file kOmegaSST.H.

Constructor & Destructor Documentation

◆ kOmegaSST()

kOmegaSST ( const fvMesh mesh,
const solverControl SolverControl 
)

◆ ~kOmegaSST()

virtual ~kOmegaSST ( )
virtualdefault

Destructor.

Member Function Documentation

◆ allocateMeanFields()

void allocateMeanFields ( )
protectedvirtual

Reimplemented from RASModelVariables.

Definition at line 51 of file kOmegaSST.C.

References RASModelVariables::allocateMeanFields(), IOobject::AUTO_WRITE, solverControl::average(), Foam::dimArea, Foam::dimTime, kOmegaSST::GMean_, RASModelVariables::mesh_, Foam::pow3(), IOobject::READ_IF_PRESENT, RASModelVariables::solverControl_, fvMesh::time(), Time::timeName(), and Foam::Zero.

Referenced by kOmegaSST::kOmegaSST().

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

◆ computeG()

tmp< volScalarField::Internal > computeG ( )
protected

Definition at line 100 of file kOmegaSST.C.

References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dev(), turbulenceModel::GName(), Foam::fvc::grad(), IOobject::groupName(), objectRegistry::lookupObject(), RASModelVariables::mesh_, Time::New(), RASModelVariables::nutRefInst(), phaseSystem::propertiesName, RASModelVariables::TMVar2(), RASModelVariables::TMVar2Inst(), Foam::twoSymm(), Foam::type(), U, and turbulenceModel::U().

Referenced by kOmegaSST::computeMeanFields(), and kOmegaSST::G().

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

◆ TypeName()

TypeName ( "kOmegaSST"  )

Runtime type information.

◆ G()

tmp< volScalarField::Internal > G ( )
virtual

Return the turbulence production term.

Reimplemented from RASModelVariables.

Definition at line 135 of file kOmegaSST.C.

References kOmegaSST::computeG(), DebugInfo, Foam::endl(), kOmegaSST::GMean_, RASModelVariables::solverControl_, and solverControl::useAveragedFields().

Referenced by kOmegaSST::correctBoundaryConditions().

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

◆ computeMeanFields()

void computeMeanFields ( )
virtual

Compute mean fields on the fly.

Reimplemented from RASModelVariables.

Definition at line 149 of file kOmegaSST.C.

References solverControl::averageIter(), kOmegaSST::computeG(), RASModelVariables::computeMeanFields(), solverControl::doAverageIter(), kOmegaSST::GMean_, and RASModelVariables::solverControl_.

Here is the call graph for this function:

◆ hasTMVar1()

virtual bool hasTMVar1 ( ) const
inlinevirtual

Bools to identify which turbulent fields are present.

Apart from the distance pointer, all other pointers are allocated even if the the corresponding field does not exist. Hence, the pointer itself cannot be used to determine the existance of the field

Reimplemented from RASModelVariables.

Definition at line 105 of file kOmegaSST.H.

◆ hasTMVar2()

virtual bool hasTMVar2 ( ) const
inlinevirtual

Reimplemented from RASModelVariables.

Definition at line 109 of file kOmegaSST.H.

◆ hasNut()

virtual bool hasNut ( ) const
inlinevirtual

Reimplemented from RASModelVariables.

Definition at line 113 of file kOmegaSST.H.

◆ correctBoundaryConditions()

void correctBoundaryConditions ( const incompressible::turbulenceModel turbulence)
virtual

Correct boundary conditions of turbulent fields.

Reimplemented from RASModelVariables.

Definition at line 163 of file kOmegaSST.C.

References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), kOmegaSST::G(), Foam::fvc::grad(), Foam::magSqr(), RASModelVariables::nutRef(), Foam::symm(), turbulence, and U.

Here is the call graph for this function:

Member Data Documentation

◆ GMean_

autoPtr<volScalarField::Internal> GMean_
protected

Average of the production term.

Avaraged separetely due the bi-zonal treatment next to the wall

Definition at line 67 of file kOmegaSST.H.

Referenced by kOmegaSST::allocateMeanFields(), kOmegaSST::computeMeanFields(), and kOmegaSST::G().


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