33 void Foam::fanPressureJumpFvPatchScalarField::calcFanJump()
35 const fvsPatchField<scalar>& phip =
38 const scalar
sign = reverse_ ? -1 : 1;
40 if (fanCurve_.valid())
44 scalar volFlowRate = 0;
46 if (phip.internalField().dimensions() ==
dimFlux)
48 volFlowRate =
gSum(phip);
58 volFlowRate =
gSum(phip/rhop);
63 <<
"dimensions of phi are not correct"
64 <<
"\n on patch " << patch().name()
65 <<
" of field " << internalField().name()
66 <<
" in file " << internalField().objectPath() <<
nl
78 if (phip.internalField().dimensions() ==
dimFlux)
87 const fvPatchField<scalar>& rhop =
95 <<
"dimensions of phi are not correct"
96 <<
"\n on patch " << patch().name()
97 <<
" of field " << internalField().name()
98 <<
" in file " << internalField().objectPath() <<
nl
102 jump_ =
sign*
max(jumpTable_->value(Un), scalar(0));
116 fixedJumpFvPatchScalarField(
p, iF),
119 reverse_(
dict.lookupOrDefault<
Switch>(
"reverse", false)),
120 phiName_(
dict.lookupOrDefault<
word>(
"phi",
"phi")),
121 rhoName_(
dict.lookupOrDefault<
word>(
"rho",
"rho"))
123 if (cyclicPatch().owner())
125 if (
dict.found(
"fanCurve"))
130 else if (
dict.found(
"jumpTable"))
142 if (
dict.found(
"value"))
164 fixedJumpFvPatchScalarField(ptf,
p, iF, mapper),
165 fanCurve_(ptf.fanCurve_, false),
166 jumpTable_(ptf.jumpTable_, false),
167 reverse_(ptf.reverse_),
168 phiName_(ptf.phiName_),
169 rhoName_(ptf.rhoName_)
179 fixedJumpFvPatchScalarField(ptf, iF),
180 fanCurve_(ptf.fanCurve_, false),
181 jumpTable_(ptf.jumpTable_, false),
182 reverse_(ptf.reverse_),
183 phiName_(ptf.phiName_),
184 rhoName_(ptf.rhoName_)
197 if (cyclicPatch().owner())
202 fixedJumpFvPatchScalarField::updateCoeffs();
210 if (cyclicPatch().owner())
212 if (fanCurve_.valid())
222 writeEntryIfDifferent<Switch>(os,
"reverse",
false, reverse_);
223 writeEntryIfDifferent<word>(os,
"phi",
"phi", phiName_);
224 writeEntryIfDifferent<word>(os,
"rho",
"rho", rhoName_);
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
A list of keyword definitions, which are a keyword followed by any number of values (e....
This boundary condition provides a pressure jump condition, using the cyclic condition as a base.
virtual void write(Ostream &) const
Write.
fanPressureJumpFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar sign(const dimensionedScalar &ds)
const dimensionSet dimMassFlux
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
VolField< scalar > volScalarField
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const dimensionSet dimFlux