Go to the documentation of this file.
31 #include "surfaceInterpolate.H"
38 namespace functionObjects
62 bool Foam::functionObjects::PecletNo::calc()
64 if (foundObject<surfaceScalarField>(fieldName_))
66 tmp<volScalarField> nuEff;
70 lookupObject<turbulenceModel>
75 nuEff = model.nuEff();
77 else if (mesh_.foundObject<dictionary>(
"transportProperties"))
79 const dictionary& model =
80 mesh_.lookupObject<dictionary>(
"transportProperties");
87 mesh_.time().timeName(),
99 <<
"Unable to determine the viscosity"
113 *mesh_.surfaceInterpolation::deltaCoeffs()
135 setResultName(
"Pe",
"phi");
144 rhoName_ =
dict.getOrDefault<
word>(
"rho",
"rho");
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
bool read(const char *buf, int32_t &val)
Same as readInt32.
static const word propertiesName
Default name of the turbulence properties dictionary.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
PecletNo(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
word name(const complex &c)
Return string representation of complex.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimViscosity
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
errorManipArg< error, int > exit(error &err, const int errNo=1)
virtual bool read(const dictionary &)
Read the PecletNo data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
defineTypeNameAndDebug(ObukhovLength, 0)