46template<
class RhoFieldType>
49 const RhoFieldType&
rho,
71 label celli =
cells[i];
74 vector mc = fc^(
C[celli] - origin);
78 scalar radius =
x[i].x();
79 scalar coeff2 =
rho[celli]*coeff1*
pow4(radius);
82 cf[0] += (fc & yawAxis)/(coeff2 + ROOTVSMALL);
83 cf[1] += (mc & pitchAxis)/(coeff2*radius + ROOTVSMALL);
84 cf[2] += (mc & rollAxis)/(coeff2*radius + ROOTVSMALL);
88 cf[0] += fc & yawAxis;
89 cf[1] += mc & pitchAxis;
90 cf[2] += mc & rollAxis;
100template<
class RhoFieldType>
103 const RhoFieldType&
rho,
108 if (rotor_.mesh().time().timeIndex() % calcFrequency_ == 0)
110 word calcType =
"forces";
113 calcType =
"coefficients";
117 <<
" solving for target trim " << calcType <<
nl;
119 const scalar rhoRef = rotor_.rhoRef();
127 while ((err > tol_) && (iter < nIter_))
133 old = calcCoeffs(
rho,
U, thetag(), force);
137 for (label pitchI = 0; pitchI < 3; pitchI++)
139 theta_[pitchI] -= dTheta_/2.0;
140 vector cf0 = calcCoeffs(
rho,
U, thetag(), force);
142 theta_[pitchI] += dTheta_;
143 vector cf1 = calcCoeffs(
rho,
U, thetag(), force);
145 vector ddTheta = (cf1 - cf0)/dTheta_;
147 J[pitchI + 0] = ddTheta[0];
148 J[pitchI + 3] = ddTheta[1];
149 J[pitchI + 6] = ddTheta[2];
155 vector dt =
inv(J) & (target_/rhoRef - old);
158 vector thetaNew = theta_ + relax_*dt;
161 err =
mag(thetaNew - theta_);
170 Info<<
" solution not converged in " << iter
171 <<
" iterations, final residual = " << err
172 <<
"(" << tol_ <<
")" <<
endl;
176 Info<<
" final residual = " << err <<
"(" << tol_
177 <<
"), iterations = " << iter <<
endl;
180 Info<<
" current and target " << calcType <<
nl
181 <<
" thrust = " << old[0]*rhoRef <<
", " << target_[0] <<
nl
182 <<
" pitch = " << old[1]*rhoRef <<
", " << target_[1] <<
nl
183 <<
" roll = " << old[2]*rhoRef <<
", " << target_[2] <<
nl
184 <<
" new pitch angles [deg]:" <<
nl
222 const dictionary& targetDict(coeffs_.subDict(
"target"));
223 useCoeffs_ = targetDict.
getOrDefault(
"useCoeffs",
true);
230 targetDict.
readEntry(
"thrust" + ext, target_[0]);
231 targetDict.
readEntry(
"pitch" + ext, target_[1]);
232 targetDict.
readEntry(
"roll" + ext, target_[2]);
234 const dictionary& pitchAngleDict(coeffs_.subDict(
"pitchAngles"));
235 theta_[0] =
degToRad(pitchAngleDict.
get<scalar>(
"theta0Ini"));
236 theta_[1] =
degToRad(pitchAngleDict.
get<scalar>(
"theta1cIni"));
237 theta_[2] =
degToRad(pitchAngleDict.
get<scalar>(
"theta1sIni"));
239 coeffs_.readEntry(
"calcFrequency", calcFrequency_);
241 coeffs_.readIfPresent(
"nIter", nIter_);
242 coeffs_.readIfPresent(
"tol", tol_);
243 coeffs_.readIfPresent(
"relax", relax_);
245 if (coeffs_.readIfPresent(
"dTheta", dTheta_))
250 coeffs_.readIfPresent(
"alpha", alpha_);
259 auto& t = ttheta.ref();
263 scalar
psi =
x[i].y();
264 t[i] = theta_[0] + theta_[1]*
cos(
psi) + theta_[2]*
sin(
psi);
288 correctTrim(
rho,
U, force);
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Graphite solid properties.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
virtual const vector e2() const
The local Cartesian y-axis in global coordinates.
virtual const point & origin() const
Return origin.
virtual const vector e1() const
The local Cartesian x-axis in global coordinates.
virtual const vector e3() const
The local Cartesian z-axis in global coordinates.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
const volVectorField & C() const
Return cell centres as volVectorField.
const labelList & cells() const noexcept
Return const access to the cell selection.
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Applies cell-based momentum sources on velocity (i.e. U) within a specified cylindrical region to app...
void calculate(const RhoFieldType &rho, const vectorField &U, const scalarField &thetag, vectorField &force, const bool divideVolume=true, const bool output=true) const
Calculate forces.
const coordSystem::cylindrical & coordSys() const
Return the rotor coordinate system (r-theta-z)
const List< point > & x() const
Return the cell centre positions in local rotor frame.
scalar omega() const
Return the rotational speed [rad/s].
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Trim model where the operating characteristics of rotor (e.g. blade pitch angle) can vary to reach a ...
virtual tmp< scalarField > thetag() const
Return the geometric angle of attack [rad].
void read(const dictionary &dict)
Read.
bool useCoeffs_
Flag to indicate whether to solve coeffs (true) or forces (false)
scalar alpha_
Coefficient to allow for conversion between US and EU definitions.
void correctTrim(const RhoFieldType &rho, const vectorField &U, vectorField &force)
Correct the model.
vector calcCoeffs(const RhoFieldType &rho, const vectorField &U, const scalarField &alphag, vectorField &force) const
Calculate the rotor force and moment coefficients vector.
A class for managing temporary objects.
Base class for trim models for handling blade characteristics and thrust-torque relations.
const fv::rotorDiskSource & rotor_
Reference to the rotor source model.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const volScalarField & psi
constexpr scalar pi(M_PI)
Different types of constants.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Ostream & endl(Ostream &os)
Add newline and flush stream.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.