40#ifndef RASModelVariables_H
41#define RASModelVariables_H
53namespace incompressible
122 (
mesh, SolverControl)
185 inline virtual bool hasNut()
const;
Templated abstract base class for single-phase incompressible turbulence models.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Mesh data needed to do the Finite Volume discretisation.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
virtual void transfer(RASModelVariables &rmv)
Transfer turbulence fields from an another object.
refPtr< volScalarField > TMVar1Ptr_
static autoPtr< RASModelVariables > New(const fvMesh &mesh, const solverControl &SolverControl)
Return a reference to the selected turbulence model.
const volScalarField & TMVar2() const
const volScalarField & TMVar1Inst() const
return references to instantaneous turbulence fields
declareRunTimeSelectionTable(autoPtr, RASModelVariables, dictionary,(const fvMesh &mesh, const solverControl &SolverControl),(mesh, SolverControl))
refPtr< volScalarField > distPtr_
virtual void allocateMeanFields()
refPtr< volScalarField > nutInitPtr_
virtual void allocateInitValues()
autoPtr< RASModelVariables > clone() const
Clone.
void restoreInitValues()
Restore turbulent fields to their initial values.
virtual ~RASModelVariables()=default
const volScalarField & d() const
virtual bool hasNut() const
virtual bool hasTMVar2() const
refPtr< volScalarField > cloneRefPtr(const refPtr< volScalarField > &obj) const
refPtr< volScalarField > nutMeanPtr_
virtual tmp< volScalarField > nutJacobianVar1(const singlePhaseTransportModel &laminarTransport) const
Return nut Jacobian wrt the TM vars.
refPtr< volScalarField > TMVar1InitPtr_
virtual tmp< volScalarField > nutJacobianVar2(const singlePhaseTransportModel &laminarTransport) const
const word & nutBaseName() const
refPtr< volScalarField > TMVar2Ptr_
void copyAndRename(volScalarField &f1, volScalarField &f2)
const volScalarField & nutRefInst() const
virtual void computeMeanFields()
Compute mean fields on the fly.
const word & TMVar1BaseName() const
Turbulence field names.
tmp< volSymmTensorField > devReff(const singlePhaseTransportModel &laminarTransport, const volVectorField &U) const
Return stress tensor based on the mean flow variables.
void resetMeanFields()
Reset mean fields to zero.
const solverControl & solverControl_
const volScalarField & TMVar1() const
Return references to turbulence fields.
refPtr< volScalarField > TMVar1MeanPtr_
TypeName("RASModelVariables")
Runtime type information.
virtual bool hasTMVar1() const
Bools to identify which turbulent fields are present.
const word & TMVar2BaseName() const
refPtr< volScalarField > TMVar2MeanPtr_
virtual tmp< volScalarField::Internal > G()
Return the turbulence production term.
const volScalarField & nutRef() const
const volScalarField & TMVar2Inst() const
refPtr< volScalarField > TMVar2InitPtr_
void operator=(const RASModelVariables &)=delete
No copy assignment.
refPtr< volScalarField > nutPtr_
A class for managing references or pointers (no reference counting)
A simple single-phase transport model based on viscosityModel.
Base class for solver control classes.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
cellMask correctBoundaryConditions()
singlePhaseTransportModel laminarTransport(U, phi)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.