Go to the documentation of this file.
37 template<
class RhoFieldType>
40 const RhoFieldType&
rho,
44 const bool divideVolume,
53 scalar AOAmin = GREAT;
54 scalar AOAmax = -GREAT;
57 const bool hasCache =
bool(Rcyl_);
61 if (area_[i] > ROOTVSMALL)
63 const label celli = cells_[i];
65 const scalar radius = x_[i].x();
71 : coordSys_.R(mesh_.C()[celli])
84 Uc.
y() = radius*omega_ - Uc.
y();
93 blade_.interpolate(radius, twist, chord, i1, i2, invDr);
96 scalar alphaGeom = thetag[i] + twist;
117 const label profile1 = blade_.profileID()[i1];
118 const label profile2 = blade_.profileID()[i2];
122 profiles_[profile1].Cdl(
alphaEff, Cd1, Cl1);
126 profiles_[profile2].Cdl(
alphaEff, Cd2, Cl2);
128 scalar Cd = invDr*(Cd2 - Cd1) + Cd1;
129 scalar Cl = invDr*(Cl2 - Cl1) + Cl1;
132 scalar tipFactor =
neg(radius/rMax_ - tipEffect_);
141 dragEff += rhoRef_*localForce.
y();
142 liftEff += rhoRef_*localForce.
z();
148 force[celli] =
transform(Rcyl, localForce);
152 force[celli] /= V[celli];
165 <<
" min/max(AOA) = " <<
radToDeg(AOAmin) <<
", "
167 <<
" Effective drag = " << dragEff <<
nl
168 <<
" Effective lift = " << liftEff <<
endl;
183 if (mesh_.time().writeTime() || writeNow)
190 mesh_.time().timeName(),
199 auto&
field = tfield.ref().primitiveFieldRef();
201 if (cells_.size() !=
values.size())
209 const label celli = cells_[i];
void calculate(const RhoFieldType &rho, const vectorField &U, const scalarField &thetag, vectorField &force, const bool divideVolume=true, const bool output=true) const
Calculate forces.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const Cmpt & x() const
Access to the vector x component.
static Ostream & output(Ostream &os, const IntRange< T > &range)
A class for handling words, derived from Foam::string.
volScalarField pDyn(IOobject("pDyn", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimPressure, Zero))
static constexpr const zero Zero
Global zero (0)
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Different types of constants.
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
Unit conversion functions.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
const Cmpt & z() const
Access to the vector z component.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
#define forAll(list, i)
Loop across all elements in list.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
constexpr scalar twoPi(2 *M_PI)
Generic dimensioned Type class.
void writeField(const word &name, const List< Type > &values, const bool writeNow=false) const
Helper function to write rotor values.
errorManip< error > abort(error &err)
Vector< scalar > vector
A scalar version of the templated Vector.
const Cmpt & y() const
Access to the vector y component.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
constexpr scalar pi(M_PI)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
dimensionSet invTransform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.