Go to the documentation of this file.
49 if (this->cyclicPatch().owner())
57 int dir = reversed_ ? -1 : 1;
60 scalar volFlowRate = 0;
64 volFlowRate = dir*
gSum(phip);
70 volFlowRate = dir*
gSum(phip/rhop);
74 FatalErrorIn(
"fanFvPatchField<Foam::scalar>::calcFanJump()")
75 <<
"dimensions of phi are not correct"
76 <<
"\n on patch " << patch().name()
83 this->jump_ = this->jumpTable_->value(
max(volFlowRate, scalar(0))) * dir;
99 phiName_(
dict.lookupOrDefault<
word>(
"phi",
"phi")),
100 rhoName_(
dict.lookupOrDefault<
word>(
"rho",
"none")),
103 if (this->cyclicPatch().owner())
109 is.
format(IOstream::ASCII);
115 if (
mag(
f[powI]) > VSMALL)
124 if (
mag(
f[powI]) > VSMALL)
130 this->jumpTable_.reset
141 reversed_.readIfPresent(
"reversed",
dict);
144 if (
dict.found(
"value"))
146 fvPatchScalarField::operator=
153 this->evaluate(Pstream::blocking);
fvPatchField< scalar > fvPatchScalarField
streamFormat format() const
Return current stream format.
PolynomialEntry container data entry for scalars. Items are stored in a list of Tuple2's....
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
#define forAll(list, i)
Loop across all elements in list.
const dimensionSet dimVelocity
const dimensionSet dimDensity
rDeltaT dimensionedInternalField()
void calcFanJump()
Calculate the fan pressure jump.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Type gSum(const FieldField< Field, Type > &f)
dimensioned< scalar > mag(const dimensioned< Type > &)
makeTemplatePatchTypeField(fvPatchScalarField, fanFvPatchScalarField)
const dimensionSet dimArea(sqr(dimLength))
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
A list of keyword definitions, which are a keyword followed by any number of values (e....
Macros for easy insertion into run-time selection tables.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
volScalarField scalarField(fieldObject, mesh)
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Generic GeometricField class.
fanFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...