Go to the documentation of this file.
76 phiName_(
dict.lookupOrDefault<
word>(
"phi",
"phi")),
77 rhoName_(
dict.lookupOrDefault<
word>(
"rho",
"rho")),
108 rpm_(ptf.
rpm_,
false)
121 const scalar t = this->db().time().timeOutputValue();
122 const scalar flowRate =
flowRate_->value(t);
123 const scalar rpm =
rpm_->value(t);
125 const scalar totArea =
gSum(patch().magSf());
126 const scalar avgU = -flowRate/totArea;
128 const vector avgCenter =
gSum(patch().Cf()*patch().magSf())/totArea;
129 const vector avgNormal =
gSum(patch().Sf())/totArea;
135 * (patch().Cf() - avgCenter) ^ avgNormal
159 <<
"dimensions of " <<
phiName_ <<
" are incorrect" <<
nl
160 <<
" on patch " << this->patch().name()
176 writeEntryIfDifferent<word>(os,
"phi",
"phi", phiName_);
177 writeEntryIfDifferent<word>(os,
"rho",
"rho", rhoName_);
178 flowRate_->writeData(os);
180 writeEntry(
"value", os);
autoPtr< DataEntry< scalar > > rpm_
Angular speed in revolutions per minute (RPM)
virtual void write(Ostream &) const
Write.
swirlFlowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
This boundary condition provides a volumetric- OR mass-flow normal vector boundary condition by its m...
A class for handling words, derived from string.
A class for managing temporary objects.
const dimensionSet dimVelocity
const dimensionSet dimDensity
rDeltaT dimensionedInternalField()
virtual void write(Ostream &) const
Write.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Type gSum(const FieldField< Field, Type > &f)
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const dimensionSet dimArea(sqr(dimLength))
const word phiName_
Name of the flux transporting the field.
autoPtr< DataEntry< scalar > > flowRate_
Inlet integral flow rate.
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.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const word rhoName_
Name of the density field used to normalize the mass flux.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...