Abstract base class for objective functions. No point in making this runTime selectable since its children will have different constructors. More...
#include <RASModelVariables.H>


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< RASModelVariables > | clone () const |
| Clone. | |
| virtual | ~RASModelVariables ()=default |
| const word & | TMVar1BaseName () const |
| Turbulence field names. | |
| const word & | TMVar2BaseName () const |
| const word & | nutBaseName () 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 volScalarField & | TMVar1 () const |
| Return references to turbulence fields. | |
| volScalarField & | TMVar1 () |
| const volScalarField & | TMVar2 () const |
| volScalarField & | TMVar2 () |
| const volScalarField & | nutRef () const |
| volScalarField & | nutRef () |
| tmp< volScalarField > | nut () const |
| const volScalarField & | d () const |
| volScalarField & | d () |
| tmp< scalarField > | TMVar1 (const label patchi) const |
| tmp< scalarField > | TMVar2 (const label patchi) const |
| tmp< scalarField > | nut (const label patchi) const |
| tmp< fvPatchScalarField > | nutPatchField (const label patchi) const |
| const volScalarField & | TMVar1Inst () const |
| Return references to instantaneous turbulence fields. | |
| volScalarField & | TMVar1Inst () |
| const volScalarField & | TMVar2Inst () const |
| volScalarField & | TMVar2Inst () |
| const volScalarField & | nutRefInst () const |
| volScalarField & | nutRefInst () |
| virtual tmp< volScalarField > | nutJacobianVar1 (const singlePhaseTransportModel &laminarTransport) const |
| Return nut Jacobian wrt the TM vars. | |
| virtual tmp< volScalarField > | nutJacobianVar2 (const singlePhaseTransportModel &laminarTransport) const |
| virtual tmp< volScalarField::Internal > | G () |
| 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< volSymmTensorField > | devReff (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< RASModelVariables > | New (const fvMesh &mesh, const solverControl &SolverControl) |
| Return a reference to the selected turbulence model. | |
Protected Member Functions | |
| virtual void | allocateInitValues () |
| virtual void | allocateMeanFields () |
| refPtr< volScalarField > | cloneRefPtr (const refPtr< volScalarField > &obj) const |
| void | copyAndRename (volScalarField &f1, volScalarField &f2) |
| void | operator= (const RASModelVariables &)=delete |
| No copy assignment. | |
Protected Attributes | |
| const fvMesh & | mesh_ |
| const solverControl & | solverControl_ |
| word | TMVar1BaseName_ |
| word | TMVar2BaseName_ |
| word | nutBaseName_ |
| refPtr< volScalarField > | TMVar1Ptr_ |
| refPtr< volScalarField > | TMVar2Ptr_ |
| refPtr< volScalarField > | nutPtr_ |
| refPtr< volScalarField > | distPtr_ |
| refPtr< volScalarField > | TMVar1InitPtr_ |
| refPtr< volScalarField > | TMVar2InitPtr_ |
| refPtr< volScalarField > | nutInitPtr_ |
| refPtr< volScalarField > | TMVar1MeanPtr_ |
| refPtr< volScalarField > | TMVar2MeanPtr_ |
| refPtr< volScalarField > | nutMeanPtr_ |
Abstract base class for objective functions. No point in making this runTime selectable since its children will have different constructors.
Definition at line 55 of file RASModelVariables.H.
| RASModelVariables | ( | const fvMesh & | mesh, |
| const solverControl & | SolverControl ) |
Construct from components.
Definition at line 175 of file RASModelVariables.C.
References distPtr_, mesh, mesh_, nutBaseName_, nutInitPtr_, nutMeanPtr_, nutPtr_, solverControl_, TMVar1BaseName_, TMVar1InitPtr_, TMVar1MeanPtr_, TMVar1Ptr_, TMVar2BaseName_, TMVar2InitPtr_, TMVar2MeanPtr_, and TMVar2Ptr_.
Referenced by declareRunTimeSelectionTable(), kEpsilon::kEpsilon(), kOmegaSST::kOmegaSST(), laminar::laminar(), LaunderSharmaKE::LaunderSharmaKE(), operator=(), RASModelVariables(), SpalartAllmaras::SpalartAllmaras(), and transfer().

| 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_.

|
virtualdefault |
|
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().


|
protectedvirtual |
Reimplemented in kEpsilon, and kOmegaSST.
Definition at line 74 of file RASModelVariables.C.
References IOobjectOption::AUTO_WRITE, Foam::endl(), hasNut(), hasTMVar1(), hasTMVar2(), Foam::Info, mesh_, name, nutMeanPtr_, nutRefInst(), IOobjectOption::READ_IF_PRESENT, solverControl_, TMVar1Inst(), TMVar1MeanPtr_, TMVar2Inst(), and TMVar2MeanPtr_.
Referenced by kEpsilon::allocateMeanFields(), kOmegaSST::allocateMeanFields(), LaunderSharmaKE::LaunderSharmaKE(), and SpalartAllmaras::SpalartAllmaras().


|
protected |
Definition at line 141 of file RASModelVariables.C.
References mesh_, IOobject::name(), refPtr< T >::New(), and timeName.
Referenced by RASModelVariables().


|
protected |
Definition at line 156 of file RASModelVariables.C.
References IOobject::name(), and regIOobject::rename().
Referenced by transfer().


|
protecteddelete |
| TypeName | ( | "RASModelVariables" | ) |
Runtime type information.
| declareRunTimeSelectionTable | ( | autoPtr | , |
| RASModelVariables | , | ||
| dictionary | , | ||
| (const fvMesh &mesh, const solverControl &SolverControl) | , | ||
| (mesh, SolverControl) | ) |
| 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().

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


|
inline |
Turbulence field names.
Definition at line 31 of file RASModelVariablesI.H.
References TMVar1BaseName_.
|
inline |
Definition at line 37 of file RASModelVariablesI.H.
References TMVar2BaseName_.
|
inline |
Definition at line 43 of file RASModelVariablesI.H.
References nutBaseName_.
|
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().

|
inlinevirtual |
Reimplemented in laminar, and LaunderSharmaKE.
Definition at line 55 of file RASModelVariablesI.H.
References TMVar2Ptr_.
Referenced by allocateInitValues(), allocateMeanFields(), computeMeanFields(), correctBoundaryConditions(), incompressibleVars::renameTurbulenceFields(), restoreInitValues(), TMVar2(), and transfer().

|
inlinevirtual |
Reimplemented in laminar, and LaunderSharmaKE.
Definition at line 61 of file RASModelVariablesI.H.
References nutPtr_.
Referenced by allocateInitValues(), allocateMeanFields(), computeMeanFields(), correctBoundaryConditions(), nut(), nut(), nutPatchField(), incompressibleVars::renameTurbulenceFields(), restoreInitValues(), and transfer().

|
inline |
Definition at line 67 of file RASModelVariablesI.H.
References distPtr_.
Referenced by transfer().

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

|
inline |
Definition at line 84 of file RASModelVariablesI.H.
References solverControl_, TMVar1MeanPtr_, and TMVar1Ptr_.
|
inline |
Definition at line 95 of file RASModelVariablesI.H.
References solverControl_, TMVar2MeanPtr_, and TMVar2Ptr_.
Referenced by kEpsilon::computeG(), kOmegaSST::computeG(), and TMVar2().

|
inline |
Definition at line 105 of file RASModelVariablesI.H.
References solverControl_, TMVar2MeanPtr_, and TMVar2Ptr_.
|
inline |
Definition at line 115 of file RASModelVariablesI.H.
References nutMeanPtr_, nutPtr_, and solverControl_.
Referenced by kOmegaSST::correctBoundaryConditions(), nut(), nut(), and nutPatchField().

|
inline |
Definition at line 126 of file RASModelVariablesI.H.
References nutMeanPtr_, nutPtr_, and solverControl_.
|
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().


|
inline |
Definition at line 154 of file RASModelVariablesI.H.
References distPtr_.
Referenced by transfer().

|
inline |
Definition at line 160 of file RASModelVariablesI.H.
References distPtr_.
|
inline |
Definition at line 166 of file RASModelVariablesI.H.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), hasTMVar1(), mesh_, tmp< T >::New(), TMVar1(), and Foam::Zero.

|
inline |
Definition at line 177 of file RASModelVariablesI.H.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), hasTMVar2(), mesh_, tmp< T >::New(), TMVar2(), and Foam::Zero.

|
inline |
Definition at line 188 of file RASModelVariablesI.H.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), hasNut(), mesh_, tmp< T >::New(), nutRef(), and Foam::Zero.

|
inline |
Definition at line 200 of file RASModelVariablesI.H.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), hasNut(), mesh_, tmp< T >::New(), nutRef(), and Foam::Zero.

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

|
inline |
Definition at line 219 of file RASModelVariablesI.H.
References TMVar1Ptr_.
|
inline |
Definition at line 225 of file RASModelVariablesI.H.
References TMVar2Ptr_.
Referenced by allocateInitValues(), allocateMeanFields(), kEpsilon::computeG(), kOmegaSST::computeG(), computeMeanFields(), correctBoundaryConditions(), incompressibleVars::renameTurbulenceFields(), restoreInitValues(), and transfer().

|
inline |
Definition at line 231 of file RASModelVariablesI.H.
References TMVar2Ptr_.
|
inline |
Definition at line 237 of file RASModelVariablesI.H.
References nutPtr_.
Referenced by allocateInitValues(), allocateMeanFields(), kEpsilon::computeG(), kOmegaSST::computeG(), computeMeanFields(), correctBoundaryConditions(), incompressibleVars::renameTurbulenceFields(), restoreInitValues(), and transfer().

|
inline |
Definition at line 243 of file RASModelVariablesI.H.
References nutPtr_.
|
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.

|
virtual |
Definition at line 311 of file RASModelVariables.C.
References Foam::dimless, Foam::endl(), laminarTransport(), mesh_, GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, WarningInFunction, and Foam::Zero.

|
inlinevirtual |
Return the turbulence production term.
Reimplemented in kEpsilon, and kOmegaSST.
Definition at line 245 of file RASModelVariables.H.
References NotImplemented.
| 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().

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

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


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

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


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

|
protected |
Definition at line 61 of file RASModelVariables.H.
Referenced by allocateMeanFields(), kEpsilon::allocateMeanFields(), kOmegaSST::allocateMeanFields(), cloneRefPtr(), kEpsilon::computeG(), kOmegaSST::computeG(), kEpsilon::kEpsilon(), kOmegaSST::kOmegaSST(), LaunderSharmaKE::LaunderSharmaKE(), nut(), nut(), nutJacobianVar1(), SpalartAllmaras::nutJacobianVar1(), nutJacobianVar2(), nutPatchField(), RASModelVariables(), RASModelVariables(), SpalartAllmaras::SpalartAllmaras(), TMVar1(), and TMVar2().
|
protected |
Definition at line 62 of file RASModelVariables.H.
Referenced by allocateInitValues(), allocateMeanFields(), kEpsilon::allocateMeanFields(), kOmegaSST::allocateMeanFields(), computeMeanFields(), kEpsilon::computeMeanFields(), kOmegaSST::computeMeanFields(), correctBoundaryConditions(), kEpsilon::G(), kOmegaSST::G(), nutRef(), nutRef(), RASModelVariables(), RASModelVariables(), resetMeanFields(), restoreInitValues(), TMVar1(), TMVar1(), TMVar2(), and TMVar2().
|
protected |
Definition at line 65 of file RASModelVariables.H.
Referenced by kEpsilon::kEpsilon(), kOmegaSST::kOmegaSST(), LaunderSharmaKE::LaunderSharmaKE(), RASModelVariables(), RASModelVariables(), SpalartAllmaras::SpalartAllmaras(), and TMVar1BaseName().
|
protected |
Definition at line 66 of file RASModelVariables.H.
Referenced by kEpsilon::kEpsilon(), kOmegaSST::kOmegaSST(), LaunderSharmaKE::LaunderSharmaKE(), RASModelVariables(), RASModelVariables(), and TMVar2BaseName().
|
protected |
Definition at line 67 of file RASModelVariables.H.
Referenced by kEpsilon::kEpsilon(), kOmegaSST::kOmegaSST(), LaunderSharmaKE::LaunderSharmaKE(), nutBaseName(), RASModelVariables(), RASModelVariables(), and SpalartAllmaras::SpalartAllmaras().
|
protected |
Definition at line 69 of file RASModelVariables.H.
Referenced by hasTMVar1(), kEpsilon::kEpsilon(), kOmegaSST::kOmegaSST(), LaunderSharmaKE::LaunderSharmaKE(), RASModelVariables(), RASModelVariables(), SpalartAllmaras::SpalartAllmaras(), TMVar1(), TMVar1(), TMVar1Inst(), and TMVar1Inst().
|
protected |
Definition at line 70 of file RASModelVariables.H.
Referenced by hasTMVar2(), kEpsilon::kEpsilon(), kOmegaSST::kOmegaSST(), LaunderSharmaKE::LaunderSharmaKE(), RASModelVariables(), RASModelVariables(), TMVar2(), TMVar2(), TMVar2Inst(), and TMVar2Inst().
|
protected |
Definition at line 71 of file RASModelVariables.H.
Referenced by hasNut(), kEpsilon::kEpsilon(), kOmegaSST::kOmegaSST(), LaunderSharmaKE::LaunderSharmaKE(), nutRef(), nutRef(), nutRefInst(), nutRefInst(), RASModelVariables(), RASModelVariables(), and SpalartAllmaras::SpalartAllmaras().
|
protected |
Definition at line 72 of file RASModelVariables.H.
Referenced by d(), d(), hasDist(), kOmegaSST::kOmegaSST(), RASModelVariables(), RASModelVariables(), and SpalartAllmaras::SpalartAllmaras().
|
protected |
Definition at line 76 of file RASModelVariables.H.
Referenced by allocateInitValues(), RASModelVariables(), RASModelVariables(), and restoreInitValues().
|
protected |
Definition at line 77 of file RASModelVariables.H.
Referenced by allocateInitValues(), RASModelVariables(), RASModelVariables(), and restoreInitValues().
|
protected |
Definition at line 78 of file RASModelVariables.H.
Referenced by allocateInitValues(), RASModelVariables(), RASModelVariables(), and restoreInitValues().
|
protected |
Definition at line 81 of file RASModelVariables.H.
Referenced by allocateMeanFields(), computeMeanFields(), correctBoundaryConditions(), RASModelVariables(), RASModelVariables(), resetMeanFields(), TMVar1(), and TMVar1().
|
protected |
Definition at line 82 of file RASModelVariables.H.
Referenced by allocateMeanFields(), computeMeanFields(), correctBoundaryConditions(), RASModelVariables(), RASModelVariables(), resetMeanFields(), TMVar2(), and TMVar2().
|
protected |
Definition at line 83 of file RASModelVariables.H.
Referenced by allocateMeanFields(), computeMeanFields(), correctBoundaryConditions(), nutRef(), nutRef(), RASModelVariables(), RASModelVariables(), and resetMeanFields().