37void Foam::swirlFanVelocityFvPatchField::calcFanJump()
41 const scalar rpm = rpm_->value(this->
db().time().timeOutputValue());
51 const auto& pSlave = nbrPatch.lookupPatchField<
volScalarField>(pName_);
64 axisHat ^ (
patch().Cf() - origin_)
67 tanDir /= (
mag(tanDir) + SMALL);
77 const scalar rMag =
mag(pCf[i] - origin_);
79 if (rMag > rInner_ && rMag < rOuter_)
98 <<
"Effective radius 'rEff' was ill-specified in the "
99 <<
"dictionary." <<
nl
114 const vectorField tangentialVelocity(magTangU*tanDir);
126 const DimensionedField<vector, volMesh>& iF
129 fixedJumpFvPatchField<
vector>(
p, iF),
139 useRealRadius_(false)
151 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi")),
152 pName_(
dict.getOrDefault<
word>(
"p",
"p")),
153 rhoName_(
dict.getOrDefault<
word>(
"rho",
"rho")),
166 this->cyclicPatch().owner()
170 fanEff_(
dict.getOrDefault<scalar>(
"fanEff", 1)),
171 rEff_(
dict.getOrDefault<scalar>(
"rEff", 0)),
172 rInner_(
dict.getOrDefault<scalar>(
"rInner", 0)),
173 rOuter_(
dict.getOrDefault<scalar>(
"rOuter", 0)),
174 useRealRadius_(
dict.getOrDefault(
"useRealRadius", false))
187 phiName_(
rhs.phiName_),
189 rhoName_(
rhs.rhoName_),
190 origin_(
rhs.origin_),
191 rpm_(
rhs.rpm_.clone()),
192 fanEff_(
rhs.fanEff_),
194 rInner_(
rhs.rInner_),
195 rOuter_(
rhs.rOuter_),
196 useRealRadius_(
rhs.useRealRadius_)
206 phiName_(
rhs.phiName_),
208 rhoName_(
rhs.rhoName_),
209 origin_(
rhs.origin_),
210 rpm_(
rhs.rpm_.clone()),
211 fanEff_(
rhs.fanEff_),
213 rInner_(
rhs.rInner_),
214 rOuter_(
rhs.rOuter_),
215 useRealRadius_(
rhs.useRealRadius_)
226 phiName_(
rhs.phiName_),
228 rhoName_(
rhs.rhoName_),
229 origin_(
rhs.origin_),
230 rpm_(
rhs.rpm_.clone()),
231 fanEff_(
rhs.fanEff_),
233 rInner_(
rhs.rInner_),
235 useRealRadius_(
rhs.useRealRadius_)
256 if (this->cyclicPatch().owner())
258 os.writeEntryIfDifferent<word>(
"phi",
"phi", phiName_);
259 os.writeEntryIfDifferent<word>(
"p",
"p", pName_);
260 os.writeEntryIfDifferent<word>(
"rho",
"rho", rhoName_);
261 os.writeEntry(
"origin", origin_);
291 swirlFanVelocityFvPatchField
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...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
label size() const noexcept
const cyclicFvPatch & cyclicPatch() const
virtual label neighbPatchID() const
Return neighbour.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
This boundary condition provides a jump condition, using the cyclic condition as a base.
virtual void write(Ostream &) const
Write.
virtual void setJump(const Field< vector > &jump)
fixedJumpFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const fvBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
const vectorField & Cf() const
Return face centres.
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup the named field from the local registry and return the patch field corresponding to this patch...
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
This boundary condition provides a jump condition for U across a cyclic pressure jump condition and a...
virtual tmp< fvPatchField< vector > > clone() const
Return a clone.
virtual void write(Ostream &os) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
swirlFanVelocityFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr scalar rpmToRads() noexcept
Multiplication factor for revolutions/minute to radians/sec.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
fvPatchField< vector > fvPatchVectorField
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.
Unit conversion functions.