31 template<
class BasicMomentumTransportModel>
39 Info<< coeffDict_.dictName() << coeffDict_ <<
endl;
46 template<
class BasicMomentumTransportModel>
58 BasicMomentumTransportModel
69 laminarDict_(this->subOrEmptyDict(
"laminar")),
70 printCoeffs_(laminarDict_.lookupOrDefault<
Switch>(
"printCoeffs", false)),
71 coeffDict_(laminarDict_.optionalSubDict(
type +
"Coeffs"))
75 this->mesh_.deltaCoeffs();
81 template<
class BasicMomentumTransportModel>
102 if (modelDict.
found(
"laminar"))
104 const word modelType =
107 {
"model",
"laminarModel"}
110 Info<<
"Selecting laminar stress model " << modelType <<
endl;
112 typename dictionaryConstructorTable::iterator cstrIter =
113 dictionaryConstructorTablePtr_->find(modelType);
115 if (cstrIter == dictionaryConstructorTablePtr_->end())
118 <<
"Unknown laminarModel type "
119 << modelType <<
nl <<
nl
120 <<
"Valid laminarModel types:" <<
endl
121 << dictionaryConstructorTablePtr_->sortedToc()
140 Info<<
"Selecting laminar stress model "
162 template<
class BasicMomentumTransportModel>
167 laminarDict_ <<= this->subDict(
"laminar");
169 coeffDict_ <<= laminarDict_.optionalSubDict(
type() +
"Coeffs");
180 template<
class BasicMomentumTransportModel>
186 this->groupName(
"nut"),
193 template<
class BasicMomentumTransportModel>
207 template<
class BasicMomentumTransportModel>
213 this->groupName(
"k"),
220 template<
class BasicMomentumTransportModel>
226 this->groupName(
"epsilon"),
233 template<
class BasicMomentumTransportModel>
239 this->groupName(
"omega"),
246 template<
class BasicMomentumTransportModel>
252 this->groupName(
"sigma"),
259 template<
class BasicMomentumTransportModel>
262 BasicMomentumTransportModel::predict();
266 template<
class BasicMomentumTransportModel>
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
static word group(const word &name)
Return group (extension part of name)
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
ITstream & lookupBackwardsCompatible(const wordList &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream, trying a list of keywords.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
BasicMomentumTransportModel::alphaField alphaField
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate,.
virtual void printCoeffs(const word &type)
Print model coefficients.
static autoPtr< laminarModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Return a reference to the selected laminar model.
virtual void correct()
Predict the laminar viscosity.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for laminar flow.
virtual void predict()
Predict the laminar viscosity.
laminarModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2], i.e. 0 for laminar flow.
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for laminar flow.
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate,.
virtual bool read()
Read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
Momentum transport model for Stokes flow.
static typeIOobject< IOdictionary > readModelDict(const objectRegistry &obr, const word &group, bool registerObject=false)
A class for managing temporary objects.
Abstract base class for all fluid physical properties.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
errorManipArg< error, int > exit(error &err, const int errNo=1)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
const dimensionSet dimKinematicViscosity
const dimensionSet dimless
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTime
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
const dimensionSet dimVelocity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.