35 #ifndef momentumTransportModel_H
36 #define momentumTransportModel_H
86 template<
class MomentumTransportModel>
89 const typename MomentumTransportModel::alphaField&
alpha,
90 const typename MomentumTransportModel::rhoField&
rho,
127 virtual bool read() = 0;
183 return this->viscosity_.
nu();
Generic GeometricBoundaryField class.
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
bool & registerObject()
Register object created from this IOobject with registry if true.
static word group(const word &name)
Return group (extension part of name)
word group() const
Return group (extension part of name)
static word groupName(Name name, const word &group)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
const fileName & name() const
Return the dictionary name.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Mesh data needed to do the Finite Volume discretisation.
Abstract base class for turbulence models (RAS, LES and laminar).
const viscosity & viscosity_
const Time & time() const
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
word groupName(const word &name) const
const volScalarField::Boundary & yb() const
Return the near wall distance.
virtual void validate()
Validate the fields after construction.
const fvMesh & mesh() const
virtual const dictionary & coeffDict() const =0
Const access to the coefficients dictionary.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
word GName() const
Helper function to return the name of the turbulence G field.
const volVectorField & U() const
Access function to velocity field.
const surfaceScalarField & phi_
const class Foam::viscosity & properties() const
Access function to fluid properties.
virtual tmp< volScalarField > omega() const =0
Return the turbulence specific dissipation rate.
static typeIOobject< IOdictionary > readModelDict(const objectRegistry &obr, const word &group, bool registerObject=false)
virtual bool read()=0
Read model coefficients if they have changed.
const volVectorField & U_
momentumTransportModel(const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
void operator=(const momentumTransportModel &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual tmp< fvVectorMatrix > divDevTauCorr(const tmp< volTensorField > &devTauCorr, volVectorField &U) const
Return the explicit stress correction matrix.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
virtual ~momentumTransportModel()
Destructor.
virtual void correct()=0
Solve the momentum transport model equations.
virtual void predict()=0
Predict the momentum transport coefficients if possible.
const surfaceScalarField & alphaRhoPhi_
static autoPtr< MomentumTransportModel > New(const typename MomentumTransportModel::alphaField &alpha, const typename MomentumTransportModel::rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
const volScalarField & y() const
Return the wall distance.
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux field.
virtual tmp< volSymmTensorField > sigma() const =0
Return the stress tensor [m^2/s^2].
TypeName("momentumTransport")
Runtime type information.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
Registry of regIOobjects.
A class for managing temporary objects.
Templated form of IOobject providing type information for file reading and header type checking.
Abstract base class for all fluid physical properties.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
A class for handling words, derived from string.
Forward declarations of fvMatrix specialisations.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
word typedName(Name name)
Return the name of the object within the given type.
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.