45 template<
class RhoFieldType>
48 const RhoFieldType&
rho,
54 rotor_.calculate(
rho,
U, thetag, force,
false,
false);
60 const vector& origin = rotor_.coordSys().origin();
61 const vector& rollAxis = rotor_.coordSys().R().e1();
62 const vector& pitchAxis = rotor_.coordSys().R().e2();
63 const vector& yawAxis = rotor_.coordSys().R().e3();
73 vector mc = fc^(
C[celli] - origin);
77 scalar radius =
x[i].x();
78 scalar coeff2 =
rho[celli]*coeff1*
pow4(radius);
81 cf[0] += (fc & yawAxis)/(coeff2 + rootVSmall);
82 cf[1] += (mc & pitchAxis)/(coeff2*radius + rootVSmall);
83 cf[2] += (mc & rollAxis)/(coeff2*radius + rootVSmall);
87 cf[0] += fc & yawAxis;
88 cf[1] += mc & pitchAxis;
89 cf[2] += mc & rollAxis;
99 template<
class RhoFieldType>
102 const RhoFieldType&
rho,
107 if (rotor_.mesh().time().timeIndex() % calcFrequency_ == 0)
116 <<
" solving for target trim " <<
calcType <<
nl;
118 const scalar rhoRef = rotor_.rhoRef();
126 while ((err > tol_) && (iter < nIter_))
132 old = calcCoeffs(
rho,
U, thetag(), force);
136 for (
label pitchI = 0; pitchI < 3; pitchI++)
138 theta_[pitchI] -= dTheta_/2.0;
139 vector cf0 = calcCoeffs(
rho,
U, thetag(), force);
141 theta_[pitchI] += dTheta_;
142 vector cf1 = calcCoeffs(
rho,
U, thetag(), force);
144 vector ddTheta = (cf1 - cf0)/dTheta_;
146 J[pitchI + 0] = ddTheta[0];
147 J[pitchI + 3] = ddTheta[1];
148 J[pitchI + 6] = ddTheta[2];
154 vector dt =
inv(J) & (target_/rhoRef - old);
157 vector thetaNew = theta_ + relax_*dt;
160 err =
mag(thetaNew - theta_);
169 Info<<
" solution not converged in " << iter
170 <<
" iterations, final residual = " << err
171 <<
"(" << tol_ <<
")" <<
endl;
175 Info<<
" final residual = " << err <<
"(" << tol_
176 <<
"), iterations = " << iter <<
endl;
180 <<
" thrust = " << old[0]*rhoRef <<
", " << target_[0] <<
nl
181 <<
" pitch = " << old[1]*rhoRef <<
", " << target_[1] <<
nl
182 <<
" roll = " << old[2]*rhoRef <<
", " << target_[2] <<
nl
183 <<
" new pitch angles [deg]:" <<
nl
227 const dictionary& targetDict(coeffs_.subDict(
"target"));
235 target_[0] = targetDict.
lookup<scalar>(
"thrust" + ext);
236 target_[1] = targetDict.
lookup<scalar>(
"pitch" + ext);
237 target_[2] = targetDict.
lookup<scalar>(
"roll" + ext);
239 const dictionary& pitchAngleDict(coeffs_.subDict(
"pitchAngles"));
240 theta_[0] = pitchAngleDict.lookup<scalar>(
"theta0Ini",
unitDegrees);
241 theta_[1] = pitchAngleDict.lookup<scalar>(
"theta1cIni",
unitDegrees);
242 theta_[2] = pitchAngleDict.lookup<scalar>(
"theta1sIni",
unitDegrees);
244 coeffs_.lookup(
"calcFrequency") >> calcFrequency_;
246 coeffs_.readIfPresent(
"nIter", nIter_);
247 coeffs_.readIfPresent(
"tol", tol_);
248 coeffs_.readIfPresent(
"relax", relax_);
250 coeffs_.readIfPresent(
"dTheta",
unitDegrees, dTheta_);
252 alpha_ = coeffs_.lookup<scalar>(
"alpha");
265 scalar
psi =
x[i].y();
266 t[i] = theta_[0] + theta_[1]*
cos(
psi) + theta_[2]*
sin(
psi);
290 correctTrim(
rho,
U, force);
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
Graphite solid properties.
Generic GeometricField class.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A list of keyword definitions, which are a keyword followed by any number of values (e....
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Cell based momentum source which approximates the mean effects of rotor forces on a cylindrical regio...
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Target trim forces/coefficients.
virtual ~targetCoeffTrim()
Destructor.
virtual tmp< scalarField > thetag() const
Return the geometric angle of attack [rad].
targetCoeffTrim(const fv::rotorDisk &rotor, const dictionary &dict)
Constructor.
void read(const dictionary &dict)
Read.
virtual void correct(const vectorField &U, vectorField &force) const
Correct the model.
void correctTrim(const RhoFieldType &rho, const vectorField &U, vectorField &force) const
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.
T & ref() const
Return non-const reference or generate a fatal error.
virtual void read(const dictionary &dict)
Read.
A class for handling words, derived from string.
const volScalarField & psi
scalar radToDeg(const scalar rad)
Convert radians to degrees.
scalar degToRad(const scalar deg)
Convert degrees to radians.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sin(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
const unitConversion unitDegrees
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.