RASModelVariables Class Reference

Abstract base class for objective functions. No point in making this runTime selectable since its children will have different constructors. More...

Inheritance diagram for RASModelVariables:
Collaboration diagram for RASModelVariables:

Public Member Functions

 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 ()
 
tmp< volScalarFieldnut () const
 
const volScalarFieldd () const
 
volScalarFieldd ()
 
tmp< scalarFieldTMVar1 (const label patchi) const
 
tmp< scalarFieldTMVar2 (const label patchi) const
 
tmp< scalarFieldnut (const label patchi) const
 
tmp< fvPatchScalarFieldnutPatchField (const label patchi) const
 
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...
 

Static Public Member Functions

static autoPtr< RASModelVariablesNew (const fvMesh &mesh, const solverControl &SolverControl)
 Return a reference to the selected turbulence model. More...
 

Protected Member Functions

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

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_
 

Detailed Description

Abstract base class for objective functions. No point in making this runTime selectable since its children will have different constructors.

Source files

Definition at line 55 of file RASModelVariables.H.

Constructor & Destructor Documentation

◆ RASModelVariables() [1/2]

RASModelVariables ( const fvMesh mesh,
const solverControl SolverControl 
)

Construct from components.

Definition at line 176 of file RASModelVariables.C.

◆ RASModelVariables() [2/2]

Copy constructor.

Will allocate new fields (instead of referencing the ones in the turbulence model), so cannot be used directly to access the fields of the turbulence model. Mainly used for checkpointing in unsteady adjoint

Definition at line 204 of file RASModelVariables.C.

◆ ~RASModelVariables()

virtual ~RASModelVariables ( )
virtualdefault

Member Function Documentation

◆ allocateInitValues()

◆ allocateMeanFields()

◆ cloneRefPtr()

Foam::refPtr< Foam::volScalarField > cloneRefPtr ( const refPtr< volScalarField > &  obj) const
protected

Definition at line 141 of file RASModelVariables.C.

References RASModelVariables::mesh_, IOobject::name(), refPtr< T >::New(), fvMesh::time(), timeName, and Time::timeName().

Here is the call graph for this function:

◆ copyAndRename()

void copyAndRename ( volScalarField f1,
volScalarField f2 
)
protected

Definition at line 157 of file RASModelVariables.C.

References IOobject::name(), and regIOobject::rename().

Referenced by RASModelVariables::transfer().

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

◆ operator=()

void operator= ( const RASModelVariables )
protecteddelete

No copy assignment.

◆ TypeName()

TypeName ( "RASModelVariables"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
RASModelVariables  ,
dictionary  ,
(const fvMesh &mesh, const solverControl &SolverControl)  ,
(mesh, SolverControl)   
)

◆ clone()

autoPtr< RASModelVariables > clone ( ) const

Clone.

Will allocate new fields (instead of referencing the ones in the turbulence model), so cannot be used directly to access the fields of the turbulence model. Mainly used for checkpointing in unsteady adjoint

Definition at line 230 of file RASModelVariables.C.

References autoPtr< T >::New().

Here is the call graph for this function:

◆ New()

autoPtr< RASModelVariables > New ( const fvMesh mesh,
const solverControl SolverControl 
)
static

Return a reference to the selected turbulence model.

Definition at line 239 of file RASModelVariables.C.

References Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, dictionary::findDict(), Foam::Info, mesh, IOobjectOption::MUST_READ, IOobjectOption::NO_REGISTER, IOobjectOption::NO_WRITE, dictionary::null, turbulenceModel::propertiesName, and dictionary::readCompat().

Referenced by incompressibleVars::setFields().

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

◆ TMVar1BaseName()

const word & TMVar1BaseName ( ) const
inline

Turbulence field names.

Definition at line 31 of file RASModelVariablesI.H.

References RASModelVariables::TMVar1BaseName_.

◆ TMVar2BaseName()

const word & TMVar2BaseName ( ) const
inline

Definition at line 37 of file RASModelVariablesI.H.

References RASModelVariables::TMVar2BaseName_.

◆ nutBaseName()

const word & nutBaseName ( ) const
inline

Definition at line 43 of file RASModelVariablesI.H.

References RASModelVariables::nutBaseName_.

◆ hasTMVar1()

bool hasTMVar1 ( ) const
inlinevirtual

◆ hasTMVar2()

◆ hasNut()

◆ hasDist()

bool hasDist ( ) const
inline

Definition at line 67 of file RASModelVariablesI.H.

References RASModelVariables::distPtr_.

Referenced by RASModelVariables::transfer().

Here is the caller graph for this function:

◆ TMVar1() [1/3]

const volScalarField & TMVar1 ( ) const
inline

Return references to turbulence fields.

will return the mean field if it exists, otherwise the instantaneous one

Definition at line 73 of file RASModelVariablesI.H.

References RASModelVariables::solverControl_, RASModelVariables::TMVar1MeanPtr_, RASModelVariables::TMVar1Ptr_, and solverControl::useAveragedFields().

Referenced by RASModelVariables::TMVar1().

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

◆ TMVar1() [2/3]

volScalarField & TMVar1 ( )
inline

◆ TMVar2() [1/3]

const volScalarField & TMVar2 ( ) const
inline

Definition at line 95 of file RASModelVariablesI.H.

References RASModelVariables::solverControl_, RASModelVariables::TMVar2MeanPtr_, RASModelVariables::TMVar2Ptr_, and solverControl::useAveragedFields().

Referenced by kEpsilon::computeG(), kOmegaSST::computeG(), and RASModelVariables::TMVar2().

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

◆ TMVar2() [2/3]

volScalarField & TMVar2 ( )
inline

◆ nutRef() [1/2]

const volScalarField & nutRef ( ) const
inline

Definition at line 115 of file RASModelVariablesI.H.

References RASModelVariables::nutMeanPtr_, RASModelVariables::nutPtr_, RASModelVariables::solverControl_, and solverControl::useAveragedFields().

Referenced by RASModelVariables::nut(), and RASModelVariables::nutPatchField().

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

◆ nutRef() [2/2]

volScalarField & nutRef ( )
inline

◆ nut() [1/2]

◆ d() [1/2]

const volScalarField & d ( ) const
inline

Definition at line 161 of file RASModelVariablesI.H.

References RASModelVariables::distPtr_.

Referenced by RASModelVariables::transfer().

Here is the caller graph for this function:

◆ d() [2/2]

volScalarField & d ( )
inline

Definition at line 167 of file RASModelVariablesI.H.

References RASModelVariables::distPtr_.

◆ TMVar1() [3/3]

◆ TMVar2() [3/3]

◆ nut() [2/2]

◆ nutPatchField()

tmp< fvPatchScalarField > nutPatchField ( const label  patchi) const
inline

◆ TMVar1Inst() [1/2]

◆ TMVar1Inst() [2/2]

volScalarField & TMVar1Inst ( )
inline

Definition at line 226 of file RASModelVariablesI.H.

References RASModelVariables::TMVar1Ptr_.

◆ TMVar2Inst() [1/2]

◆ TMVar2Inst() [2/2]

volScalarField & TMVar2Inst ( )
inline

Definition at line 238 of file RASModelVariablesI.H.

References RASModelVariables::TMVar2Ptr_.

◆ nutRefInst() [1/2]

◆ nutRefInst() [2/2]

volScalarField & nutRefInst ( )
inline

Definition at line 250 of file RASModelVariablesI.H.

References RASModelVariables::nutPtr_.

◆ nutJacobianVar1()

tmp< volScalarField > nutJacobianVar1 ( const singlePhaseTransportModel laminarTransport) const
virtual

Return nut Jacobian wrt the TM vars.

Reimplemented in SpalartAllmaras.

Definition at line 293 of file RASModelVariables.C.

References Foam::dimless, Foam::endl(), tmp< T >::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, WarningInFunction, and Foam::Zero.

Here is the call graph for this function:

◆ nutJacobianVar2()

tmp< volScalarField > nutJacobianVar2 ( const singlePhaseTransportModel laminarTransport) const
virtual

Definition at line 318 of file RASModelVariables.C.

References Foam::dimless, Foam::endl(), tmp< T >::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, WarningInFunction, and Foam::Zero.

Here is the call graph for this function:

◆ G()

virtual tmp<volScalarField::Internal> G ( )
inlinevirtual

Return the turbulence production term.

Reimplemented in kEpsilon, and kOmegaSST.

Definition at line 245 of file RASModelVariables.H.

References NotImplemented.

◆ restoreInitValues()

◆ resetMeanFields()

◆ computeMeanFields()

◆ devReff()

tmp< volSymmTensorField > devReff ( const singlePhaseTransportModel laminarTransport,
const volVectorField U 
) const

Return stress tensor based on the mean flow variables.

Definition at line 417 of file RASModelVariables.C.

References Foam::devTwoSymm(), Foam::fvc::grad(), laminarTransport(), tmp< T >::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, nut, and U.

Here is the call graph for this function:

◆ correctBoundaryConditions()

void correctBoundaryConditions ( const incompressible::turbulenceModel turbulence)
virtual

correct bounday conditions of turbulent fields

Reimplemented in kOmegaSST.

Definition at line 438 of file RASModelVariables.C.

Referenced by kOmegaSST::correctBoundaryConditions().

Here is the caller graph for this function:

◆ transfer()

void transfer ( RASModelVariables rmv)
virtual

Transfer turbulence fields from an another object.

Copies values since the ownership of the original fields is held by the turbulence model

Definition at line 471 of file RASModelVariables.C.

References RASModelVariables::copyAndRename(), RASModelVariables::d(), RASModelVariables::hasDist(), RASModelVariables::hasNut(), RASModelVariables::hasTMVar1(), RASModelVariables::hasTMVar2(), RASModelVariables::nutRefInst(), RASModelVariables::TMVar1Inst(), and RASModelVariables::TMVar2Inst().

Here is the call graph for this function:

Member Data Documentation

◆ mesh_

◆ solverControl_

◆ TMVar1BaseName_

word TMVar1BaseName_
protected

Definition at line 65 of file RASModelVariables.H.

Referenced by RASModelVariables::TMVar1BaseName().

◆ TMVar2BaseName_

word TMVar2BaseName_
protected

Definition at line 66 of file RASModelVariables.H.

Referenced by RASModelVariables::TMVar2BaseName().

◆ nutBaseName_

word nutBaseName_
protected

Definition at line 67 of file RASModelVariables.H.

Referenced by RASModelVariables::nutBaseName().

◆ TMVar1Ptr_

◆ TMVar2Ptr_

◆ nutPtr_

◆ distPtr_

refPtr<volScalarField> distPtr_
protected

Definition at line 72 of file RASModelVariables.H.

Referenced by RASModelVariables::d(), and RASModelVariables::hasDist().

◆ TMVar1InitPtr_

refPtr<volScalarField> TMVar1InitPtr_
protected

◆ TMVar2InitPtr_

refPtr<volScalarField> TMVar2InitPtr_
protected

◆ nutInitPtr_

◆ TMVar1MeanPtr_

◆ TMVar2MeanPtr_

◆ nutMeanPtr_


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