Go to the documentation of this file.
43 { fanFlowDirection::ffdIn,
"in" },
44 { fanFlowDirection::ffdOut,
"out" },
59 nonDimensional_(
false),
74 fanCurve_(ptf.fanCurve_.clone()),
75 direction_(ptf.direction_),
76 nonDimensional_(ptf.nonDimensional_),
91 direction_(fanFlowDirectionNames_.get(
"direction",
dict)),
92 nonDimensional_(
dict.getOrDefault(
"nonDimensional",
false)),
97 if (
dict.found(
"file"))
114 dict.readEntry(
"rpm", rpm_);
115 dict.readEntry(
"dm", dm_);
126 fanCurve_(fppsf.fanCurve_.clone()),
127 direction_(fppsf.direction_),
128 nonDimensional_(fppsf.nonDimensional_),
141 fanCurve_(fppsf.fanCurve_.clone()),
142 direction_(fppsf.direction_),
143 nonDimensional_(fppsf.nonDimensional_),
163 const int dir = 2*direction_ - 1;
166 scalar volFlowRate = 0;
170 volFlowRate = dir*
gSum(phip);
176 volFlowRate = dir*
gSum(phip/rhop);
181 <<
"dimensions of phi are not correct\n"
182 <<
" on patch " <<
patch().name()
183 <<
" of field " << internalField().name()
184 <<
" in file " << internalField().objectPath() <<
nl
196 scalar pdFan = fanCurve_->value(
max(volFlowRate, scalar(0)));
207 patch().lookupPatchField<volVectorField, vector>(
UName())
215 fanCurve_->writeData(os);
216 os.
writeEntry(
"direction", fanFlowDirectionNames_[direction_]);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
This boundary condition can be applied to assign either a pressure inlet or outlet total pressure con...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
const dimensionSet dimVelocity
const dimensionSet dimDensity
virtual void write(Ostream &os) const
Write.
Type gSum(const FieldField< Field, Type > &f)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
fanPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
dimensionedScalar pow4(const dimensionedScalar &ds)
This boundary condition provides a total pressure condition. Four variants are possible:
const dimensionSet dimArea(sqr(dimLength))
dimensionedScalar pow3(const dimensionedScalar &ds)
const word & phiName() const
Return the name of the flux field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const word & UName() const
Return the name of the velocity field.
Templated table container function where data are read from file.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
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 word & rhoName() const
Return the name of the density field.
const scalarField & p0() const
Return the total pressure.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr scalar pi(M_PI)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const std::string patch
OpenFOAM patch number as a std::string.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
virtual void write(Ostream &) const
Write.
Foam::fvPatchFieldMapper.
static const Enum< fanFlowDirection > fanFlowDirectionNames_
Fan flow direction names.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
fanFlowDirection
Fan flow direction.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...