39namespace functionObjects
55 { ISENTROPIC,
"isentropic" },
56 { STATIC_COEFF,
"staticCoeff" },
57 { TOTAL_COEFF,
"totalCoeff" },
69 { SUBTRACT,
"subtract" },
75Foam::word Foam::functionObjects::pressure::resultName()
const
83 else if (mode_ &
TOTAL)
94 <<
"Unhandled calculation mode " <<
modeNames[mode_]
98 switch (hydrostaticMode_)
106 rName = rName +
"+rgh";
112 rName = rName +
"-rgh";
139 p.mesh().time().timeName(),
146 fvPatchField<scalar>::calculatedType()
150 if (!rhoInfInitialised_)
154 <<
"pressure identified as incompressible, but reference "
155 <<
"density is not set. Please set 'rho' to 'rhoInf', and "
156 <<
"set an appropriate value for 'rhoInf'"
167 const tmp<volScalarField>& tsf
172 return lookupObject<volScalarField>(rhoName_)*tsf;
179void Foam::functionObjects::pressure::addHydrostaticContribution
187 if (hydrostaticMode_ == NONE)
197 if (!hRefInitialised_)
204 (g_ & (
cmptMag(g_.value())/
mag(g_.value())))*hRef_
207 tmp<volScalarField> rgh = rhoScale(
p, (g_ & mesh_.C()) - ghRef);
209 switch (hydrostaticMode_)
230 const tmp<volScalarField>& tp
240 mesh_.time().timeName(),
250 addHydrostaticContribution(
p, result);
262 + rhoScale(
p, 0.5*
magSqr(lookupObject<volVectorField>(UName_)));
266 if (mode_ & ISENTROPIC)
268 const basicThermo* thermoPtr =
274 <<
"Isentropic pressure calculation requires a "
275 <<
"thermodynamics package"
282 mag(lookupObject<volVectorField>(UName_))
296 const tmp<volScalarField>& tp
301 tmp<volScalarField> tpCoeff(tp.ptr());
315 return std::move(tp);
321bool Foam::functionObjects::pressure::calc()
323 if (foundObject<volScalarField>(fieldName_))
332 p.mesh().time().timeName(),
337 coeff(calcPressure(
p, rhoScale(
p)))
340 return store(resultName_, tp);
358 hydrostaticMode_(NONE),
365 rhoInfInitialised_(false),
367 gInitialised_(false),
369 hRefInitialised_(false)
383 UName_ =
dict.getOrDefault<
word>(
"U",
"U");
384 rhoName_ =
dict.getOrDefault<
word>(
"rho",
"rho");
386 if (rhoName_ ==
"rhoInf")
388 dict.readEntry(
"rhoInf", rhoInf_);
389 rhoInfInitialised_ =
true;
392 if (!modeNames.readIfPresent(
"mode",
dict, mode_))
398 dict.getOrDefaultCompat<
bool>(
"mode", {{
"calcTotal", 1812}},
false);
400 dict.getOrDefaultCompat<
bool>(
"mode", {{
"calcCoeff", 1812}},
false);
413 mode_ =
static_cast<mode>(COEFF | mode_);
417 Info<<
" Operating mode: " << modeNames[mode_] <<
nl;
419 pRef_ =
dict.getOrDefault<scalar>(
"pRef", 0);
423 hydrostaticModeNames.readIfPresent
432 Info<<
" Hydrostatic mode: "
433 << hydrostaticModeNames[hydrostaticMode_]
435 gInitialised_ =
dict.readIfPresent(
"g", g_);
436 hRefInitialised_ =
dict.readIfPresent(
"hRef", hRef_);
440 Info<<
" Not including hydrostatic effects" <<
nl;
446 dict.readEntry(
"pInf", pInf_);
447 dict.readEntry(
"UInf", UInf_);
448 dict.readEntry(
"rhoInf", rhoInf_);
450 const scalar zeroCheck = 0.5*rhoInf_*
magSqr(UInf_) + pInf_;
452 if (
mag(zeroCheck) < ROOTVSMALL)
456 <<
"Coefficient calculation requested, but reference "
457 <<
"pressure level is zero. Please check the supplied "
458 <<
"values of pInf, UInf and rhoInf" <<
endl;
461 rhoInfInitialised_ =
true;
464 resultName_ =
dict.getOrDefault<
word>(
"result", resultName());
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
static const word dictName
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base-class for Time/database function objects.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
word fieldName_
Name of field to process.
Computes the magnitude of the square of an input field.
Computes the magnitude of an input field.
Provides several methods to convert an input pressure field into derived forms, including:
static const Enum< hydrostaticMode > hydrostaticModeNames
hydrostaticMode
Enumeration for hydrostatic contributions.
static const Enum< mode > modeNames
mode
Enumeration for pressure calculation mode.
@ COEFF
Coefficient manipulator.
@ ISENTROPIC
Isentropic pressure.
virtual bool read(const dictionary &)
Read the pressure data.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
@ NONE
No type, or default initialized type.
const dimensionSet dimPressure
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UniformDimensionedField< scalar > uniformDimensionedScalarField
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UniformDimensionedField< vector > uniformDimensionedVectorField
errorManip< error > abort(error &err)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimAcceleration
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)