Loading...
Searching...
No Matches
RASModelVariables Class Reference

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

#include <RASModelVariables.H>

Inheritance diagram for RASModelVariables:
Collaboration diagram for RASModelVariables:

Public Member Functions

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

Static Public Member Functions

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

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.

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() [2/2]

RASModelVariables ( const RASModelVariables & rmv)

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 203 of file RASModelVariables.C.

References cloneRefPtr(), distPtr_, mesh_, nutBaseName_, nutInitPtr_, nutMeanPtr_, nutPtr_, RASModelVariables(), solverControl_, TMVar1BaseName_, TMVar1InitPtr_, TMVar1MeanPtr_, TMVar1Ptr_, TMVar2BaseName_, TMVar2InitPtr_, TMVar2MeanPtr_, and TMVar2Ptr_.

Here is the call graph for this function:

◆ ~RASModelVariables()

virtual ~RASModelVariables ( )
virtualdefault

References laminarTransport().

Here is the call graph for this function:

Member Function Documentation

◆ allocateInitValues()

void allocateInitValues ( )
protectedvirtual

Definition at line 41 of file RASModelVariables.C.

References Foam::endl(), hasNut(), hasTMVar1(), hasTMVar2(), Foam::Info, name, nutInitPtr_, nutRefInst(), solverControl_, TMVar1InitPtr_, TMVar1Inst(), TMVar2InitPtr_, and TMVar2Inst().

Referenced by kEpsilon::kEpsilon(), kOmegaSST::kOmegaSST(), LaunderSharmaKE::LaunderSharmaKE(), and SpalartAllmaras::SpalartAllmaras().

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

◆ allocateMeanFields()

void allocateMeanFields ( )
protectedvirtual

◆ cloneRefPtr()

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

Definition at line 141 of file RASModelVariables.C.

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

Referenced by RASModelVariables().

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

◆ copyAndRename()

void copyAndRename ( volScalarField & f1,
volScalarField & f2 )
protected

Definition at line 156 of file RASModelVariables.C.

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

Referenced by 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.

References RASModelVariables().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "RASModelVariables" )

Runtime type information.

◆ declareRunTimeSelectionTable()

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

References mesh, and RASModelVariables().

Here is the call graph for this function:

◆ 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 238 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 TMVar1BaseName_.

◆ TMVar2BaseName()

const word & TMVar2BaseName ( ) const
inline

Definition at line 37 of file RASModelVariablesI.H.

References TMVar2BaseName_.

◆ nutBaseName()

const word & nutBaseName ( ) const
inline

Definition at line 43 of file RASModelVariablesI.H.

References nutBaseName_.

◆ hasTMVar1()

bool hasTMVar1 ( ) const
inlinevirtual

Bools to identify which turbulent fields are present.

Reimplemented in laminar, and LaunderSharmaKE.

Definition at line 49 of file RASModelVariablesI.H.

References TMVar1Ptr_.

Referenced by allocateInitValues(), allocateMeanFields(), computeMeanFields(), correctBoundaryConditions(), incompressibleVars::renameTurbulenceFields(), restoreInitValues(), TMVar1(), and transfer().

Here is the caller graph for this function:

◆ hasTMVar2()

bool hasTMVar2 ( ) const
inlinevirtual

◆ hasNut()

bool hasNut ( ) const
inlinevirtual

◆ hasDist()

bool hasDist ( ) const
inline

Definition at line 67 of file RASModelVariablesI.H.

References distPtr_.

Referenced by 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 solverControl_, TMVar1MeanPtr_, and TMVar1Ptr_.

Referenced by SpalartAllmaras::nutJacobianVar1(), and TMVar1().

Here is the caller graph for this function:

◆ TMVar1() [2/3]

volScalarField & TMVar1 ( )
inline

Definition at line 84 of file RASModelVariablesI.H.

References solverControl_, TMVar1MeanPtr_, and TMVar1Ptr_.

◆ TMVar2() [1/3]

const volScalarField & TMVar2 ( ) const
inline

Definition at line 95 of file RASModelVariablesI.H.

References solverControl_, TMVar2MeanPtr_, and TMVar2Ptr_.

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

Here is the caller graph for this function:

◆ TMVar2() [2/3]

volScalarField & TMVar2 ( )
inline

Definition at line 105 of file RASModelVariablesI.H.

References solverControl_, TMVar2MeanPtr_, and TMVar2Ptr_.

◆ nutRef() [1/2]

const volScalarField & nutRef ( ) const
inline

Definition at line 115 of file RASModelVariablesI.H.

References nutMeanPtr_, nutPtr_, and solverControl_.

Referenced by kOmegaSST::correctBoundaryConditions(), nut(), nut(), and nutPatchField().

Here is the caller graph for this function:

◆ nutRef() [2/2]

volScalarField & nutRef ( )
inline

Definition at line 126 of file RASModelVariablesI.H.

References nutMeanPtr_, nutPtr_, and solverControl_.

◆ nut() [1/2]

tmp< volScalarField > nut ( ) const
inline

Definition at line 137 of file RASModelVariablesI.H.

References Foam::dimViscosity, hasNut(), mesh_, GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, nutRef(), and Foam::Zero.

Referenced by devReff().

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

◆ d() [1/2]

const volScalarField & d ( ) const
inline

Definition at line 154 of file RASModelVariablesI.H.

References distPtr_.

Referenced by transfer().

Here is the caller graph for this function:

◆ d() [2/2]

volScalarField & d ( )
inline

Definition at line 160 of file RASModelVariablesI.H.

References distPtr_.

◆ TMVar1() [3/3]

tmp< scalarField > TMVar1 ( const label patchi) const
inline

Definition at line 166 of file RASModelVariablesI.H.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), hasTMVar1(), mesh_, tmp< T >::New(), TMVar1(), and Foam::Zero.

Here is the call graph for this function:

◆ TMVar2() [3/3]

tmp< scalarField > TMVar2 ( const label patchi) const
inline

Definition at line 177 of file RASModelVariablesI.H.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), hasTMVar2(), mesh_, tmp< T >::New(), TMVar2(), and Foam::Zero.

Here is the call graph for this function:

◆ nut() [2/2]

tmp< scalarField > nut ( const label patchi) const
inline

Definition at line 188 of file RASModelVariablesI.H.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), hasNut(), mesh_, tmp< T >::New(), nutRef(), and Foam::Zero.

Here is the call graph for this function:

◆ nutPatchField()

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

Definition at line 200 of file RASModelVariablesI.H.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), hasNut(), mesh_, tmp< T >::New(), nutRef(), and Foam::Zero.

Here is the call graph for this function:

◆ TMVar1Inst() [1/2]

const volScalarField & TMVar1Inst ( ) const
inline

Return references to instantaneous turbulence fields.

Definition at line 213 of file RASModelVariablesI.H.

References TMVar1Ptr_.

Referenced by allocateInitValues(), allocateMeanFields(), computeMeanFields(), correctBoundaryConditions(), incompressibleVars::renameTurbulenceFields(), restoreInitValues(), and transfer().

Here is the caller graph for this function:

◆ TMVar1Inst() [2/2]

volScalarField & TMVar1Inst ( )
inline

Definition at line 219 of file RASModelVariablesI.H.

References TMVar1Ptr_.

◆ TMVar2Inst() [1/2]

◆ TMVar2Inst() [2/2]

volScalarField & TMVar2Inst ( )
inline

Definition at line 231 of file RASModelVariablesI.H.

References TMVar2Ptr_.

◆ nutRefInst() [1/2]

◆ nutRefInst() [2/2]

volScalarField & nutRefInst ( )
inline

Definition at line 243 of file RASModelVariablesI.H.

References nutPtr_.

◆ nutJacobianVar1()

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

Return nut Jacobian wrt the TM vars.

Reimplemented in SpalartAllmaras.

Definition at line 292 of file RASModelVariables.C.

References Foam::dimless, Foam::endl(), laminarTransport(), mesh_, GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, WarningInFunction, and Foam::Zero.

Here is the call graph for this function:

◆ nutJacobianVar2()

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

◆ 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()

void restoreInitValues ( )

Restore turbulent fields to their initial values.

Definition at line 330 of file RASModelVariables.C.

References hasNut(), hasTMVar1(), hasTMVar2(), nutInitPtr_, nutRefInst(), solverControl_, TMVar1InitPtr_, TMVar1Inst(), TMVar2InitPtr_, and TMVar2Inst().

Here is the call graph for this function:

◆ resetMeanFields()

void resetMeanFields ( )

Reset mean fields to zero.

Definition at line 350 of file RASModelVariables.C.

References Foam::endl(), Foam::Info, nutMeanPtr_, solverControl_, TMVar1MeanPtr_, TMVar2MeanPtr_, and Foam::Zero.

Here is the call graph for this function:

◆ computeMeanFields()

void computeMeanFields ( )
virtual

Compute mean fields on the fly.

Reimplemented in kEpsilon, and kOmegaSST.

Definition at line 373 of file RASModelVariables.C.

References hasNut(), hasTMVar1(), hasTMVar2(), nutMeanPtr_, nutRefInst(), solverControl_, TMVar1Inst(), TMVar1MeanPtr_, TMVar2Inst(), and TMVar2MeanPtr_.

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

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

◆ devReff()

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

Return stress tensor based on the mean flow variables.

Definition at line 401 of file RASModelVariables.C.

References Foam::devTwoSymm(), Foam::fvc::grad(), laminarTransport(), GeometricField< symmTensor, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, 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 416 of file RASModelVariables.C.

References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), hasNut(), hasTMVar1(), hasTMVar2(), nutMeanPtr_, nutRefInst(), solverControl_, TMVar1Inst(), TMVar1MeanPtr_, TMVar2Inst(), and TMVar2MeanPtr_.

Referenced by kOmegaSST::correctBoundaryConditions().

Here is the call graph for this function:
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 450 of file RASModelVariables.C.

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

Here is the call graph for this function:

Member Data Documentation

◆ mesh_

◆ solverControl_

◆ TMVar1BaseName_

◆ TMVar2BaseName_

◆ nutBaseName_

◆ TMVar1Ptr_

◆ TMVar2Ptr_

◆ nutPtr_

◆ distPtr_

◆ TMVar1InitPtr_

refPtr<volScalarField> TMVar1InitPtr_
protected

◆ TMVar2InitPtr_

refPtr<volScalarField> TMVar2InitPtr_
protected

◆ nutInitPtr_

refPtr<volScalarField> nutInitPtr_
protected

◆ TMVar1MeanPtr_

◆ TMVar2MeanPtr_

◆ nutMeanPtr_


The documentation for this class was generated from the following files:
  • src/optimisation/adjointOptimisation/adjoint/turbulenceModels/turbulenceModelVariables/RAS/RASModelVariables/RASModelVariables.H
  • src/optimisation/adjointOptimisation/adjoint/turbulenceModels/turbulenceModelVariables/RAS/RASModelVariables/RASModelVariables.C
  • src/optimisation/adjointOptimisation/adjoint/turbulenceModels/turbulenceModelVariables/RAS/RASModelVariables/RASModelVariablesI.H