44 template<
class RhoFieldType>
47 const RhoFieldType& rho,
53 rotor_.calculate(rho, U, thetag, force,
false,
false);
59 const vector& origin = rotor_.coordSys().origin();
60 const vector& rollAxis = rotor_.coordSys().R().e1();
61 const vector& pitchAxis = rotor_.coordSys().R().e2();
62 const vector& yawAxis = rotor_.coordSys().R().e3();
69 label celli = cells[i];
72 vector mc = fc^(C[celli] - origin);
76 scalar radius = x[i].x();
77 scalar coeff2 = rho[celli]*coeff1*
pow4(radius);
80 cf[0] += (fc & yawAxis)/(coeff2 + rootVSmall);
81 cf[1] += (mc & pitchAxis)/(coeff2*radius + rootVSmall);
82 cf[2] += (mc & rollAxis)/(coeff2*radius + rootVSmall);
86 cf[0] += fc & yawAxis;
87 cf[1] += mc & pitchAxis;
88 cf[2] += mc & rollAxis;
98 template<
class RhoFieldType>
101 const RhoFieldType& rho,
106 if (rotor_.mesh().time().timeIndex() % calcFrequency_ == 0)
111 calcType =
"coefficients";
115 <<
" solving for target trim " << calcType <<
nl;
117 const scalar rhoRef = rotor_.rhoRef();
125 while ((err > tol_) && (iter < nIter_))
131 old = calcCoeffs(rho, U, thetag(), force);
135 for (
label pitchI = 0; pitchI < 3; pitchI++)
137 theta_[pitchI] -= dTheta_/2.0;
138 vector cf0 = calcCoeffs(rho, U, thetag(), force);
140 theta_[pitchI] += dTheta_;
141 vector cf1 = calcCoeffs(rho, U, thetag(), force);
143 vector ddTheta = (cf1 - cf0)/dTheta_;
145 J[pitchI + 0] = ddTheta[0];
146 J[pitchI + 3] = ddTheta[1];
147 J[pitchI + 6] = ddTheta[2];
153 vector dt =
inv(J) & (target_/rhoRef - old);
156 vector thetaNew = theta_ + relax_*dt;
159 err =
mag(thetaNew - theta_);
168 Info<<
" solution not converged in " << iter
169 <<
" iterations, final residual = " << err
170 <<
"(" << tol_ <<
")" <<
endl;
174 Info<<
" final residual = " << err <<
"(" << tol_
175 <<
"), iterations = " << iter <<
endl;
178 Info<<
" current and target " << calcType << nl
179 <<
" thrust = " << old[0]*rhoRef <<
", " << target_[0] << nl
180 <<
" pitch = " << old[1]*rhoRef <<
", " << target_[1] << nl
181 <<
" roll = " << old[2]*rhoRef <<
", " << target_[2] << nl
182 <<
" new pitch angles [deg]:" << nl
183 <<
" theta0 = " <<
radToDeg(theta_[0]) << nl
184 <<
" theta1c = " <<
radToDeg(theta_[1]) << nl
185 <<
" theta1s = " <<
radToDeg(theta_[2]) << nl
226 const dictionary& targetDict(coeffs_.subDict(
"target"));
234 target_[0] = targetDict.lookup<scalar>(
"thrust" + ext);
235 target_[1] = targetDict.lookup<scalar>(
"pitch" + ext);
236 target_[2] = targetDict.lookup<scalar>(
"roll" + ext);
238 const dictionary& pitchAngleDict(coeffs_.subDict(
"pitchAngles"));
239 theta_[0] =
degToRad(pitchAngleDict.lookup<scalar>(
"theta0Ini"));
240 theta_[1] =
degToRad(pitchAngleDict.lookup<scalar>(
"theta1cIni"));
241 theta_[2] =
degToRad(pitchAngleDict.lookup<scalar>(
"theta1sIni"));
243 coeffs_.lookup(
"calcFrequency") >> calcFrequency_;
245 coeffs_.readIfPresent(
"nIter", nIter_);
246 coeffs_.readIfPresent(
"tol", tol_);
247 coeffs_.readIfPresent(
"relax", relax_);
249 if (coeffs_.readIfPresent(
"dTheta", dTheta_))
254 alpha_ = coeffs_.lookup<scalar>(
"alpha");
267 scalar
psi = x[i].y();
268 t[i] = theta_[0] + theta_[1]*
cos(psi) + theta_[2]*
sin(psi);
292 correctTrim(rho, U, force);
Graphite solid properties.
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
A list of keyword definitions, which are a keyword followed by any number of values (e...
void correctTrim(const RhoFieldType &rho, const vectorField &U, vectorField &force) const
Correct the model.
virtual void read(const dictionary &dict)
Read.
T & ref() const
Return non-const reference or generate a fatal error.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
void size(const label)
Override size to be inconsistent with allocated storage.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Macros for easy insertion into run-time selection tables.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Cell based momentum source which approximates the mean effects of rotor forces on a cylindrical regio...
bool read(const char *, int32_t &)
virtual ~targetCoeffTrim()
Destructor.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
dimensionedScalar cos(const dimensionedScalar &ds)
A class for handling words, derived from string.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void correct(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.
virtual tmp< scalarField > thetag() const
Return the geometric angle of attack [rad].
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross product operators.
dimensionedScalar sin(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
targetCoeffTrim(const fv::rotorDiskSource &rotor, const dictionary &dict)
Constructor.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
void read(const dictionary &dict)
Read.
const volScalarField & psi
const doubleScalar e
Elementary charge.
A class for managing temporary objects.