Go to the documentation of this file.
37 void Foam::swirlFanVelocityFvPatchField::calcFanJump()
41 const scalar rpm = rpm_->value(this->
db().time().timeOutputValue());
46 const fvPatchField<scalar>& pOwner =
53 const fvPatchField<scalar>& pSlave =
69 axisHat ^ (
patch().Cf() - origin_)
72 tanDir /= (
mag(tanDir) + SMALL);
82 const scalar rMag =
mag(pCf[i] - origin_);
84 if (rMag > rInner_ && rMag < rOuter_)
96 <<
"Effective radius rEff should be specified in the "<<
nl
97 <<
"dictionary." <<
nl
105 const vectorField tangentialVelocity(magTangU*tanDir);
107 this->
jump_ = tangentialVelocity;
130 useRealRadius_(
false)
157 this->cyclicPatch().owner()
161 rEff_(
dict.getOrDefault<scalar>(
"rEff", 0)),
162 fanEff_(
dict.getOrDefault<scalar>(
"fanEff", 1)),
163 rInner_(
dict.getOrDefault<scalar>(
"rInner", 0)),
164 rOuter_(
dict.getOrDefault<scalar>(
"rOuter", 0)),
165 useRealRadius_(
dict.getOrDefault(
"useRealRadius", false))
178 phiName_(ptf.phiName_),
180 rhoName_(ptf.rhoName_),
181 origin_(ptf.origin_),
182 rpm_(ptf.rpm_.clone()),
184 fanEff_(ptf.fanEff_),
185 rInner_(ptf.rInner_),
186 rOuter_(ptf.rOuter_),
187 useRealRadius_(ptf.useRealRadius_)
197 phiName_(ptf.phiName_),
199 rhoName_(ptf.rhoName_),
200 origin_(ptf.origin_),
201 rpm_(ptf.rpm_.clone()),
203 rInner_(ptf.rInner_),
204 rOuter_(ptf.rOuter_),
205 useRealRadius_(ptf.useRealRadius_)
216 phiName_(ptf.phiName_),
218 rhoName_(ptf.rhoName_),
219 origin_(ptf.origin_),
220 rpm_(ptf.rpm_.clone()),
222 rInner_(ptf.rInner_),
223 rOuter_(ptf.rOuter_),
224 useRealRadius_(ptf.useRealRadius_)
245 if (this->cyclicPatch().owner())
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
swirlFanVelocityFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A class for handling words, derived from Foam::string.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
static constexpr const zero Zero
Global zero (0)
Unit conversion functions.
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...
#define forAll(list, i)
Loop across all elements in list.
Field< vector > vectorField
Specialisation of Field<T> for vector.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
virtual void write(Ostream &) const
Write.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const cyclicFvPatch & cyclicPatch() const
Return local reference cast into the cyclic patch.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const Type & lookupObject(const word &name, const bool recursive=false) const
virtual void write(Ostream &os) const
Write.
const vectorField & Cf() const
Return face centres.
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
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 dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Vector< scalar > vector
A scalar version of the templated Vector.
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr scalar rpmToRads(const scalar rpm) noexcept
Conversion from revolutions/minute to radians/sec.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const objectRegistry & db() const
Return local objectRegistry.
const std::string patch
OpenFOAM patch number as a std::string.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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 label neighbPatchID() const
Return neighbour.
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
Return patch.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
This boundary condition provides a jump condition for U across a cyclic pressure jump condition and a...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Field< vector > jump_
"jump" field