Target trim forces/coefficients. More...
Public Member Functions | |
TypeName ("targetCoeffTrim") | |
Run-time type information. More... | |
targetCoeffTrim (const fv::rotorDisk &rotor, const dictionary &dict) | |
Constructor. More... | |
virtual | ~targetCoeffTrim () |
Destructor. More... | |
void | read (const dictionary &dict) |
Read. More... | |
virtual tmp< scalarField > | thetag () const |
Return the geometric angle of attack [rad]. More... | |
virtual void | correct (const vectorField &U, vectorField &force) const |
Correct the model. More... | |
virtual void | correct (const volScalarField rho, const vectorField &U, vectorField &force) const |
Correct the model for compressible flow. More... | |
template<class RhoFieldType > | |
Foam::vector | calcCoeffs (const RhoFieldType &rho, const vectorField &U, const scalarField &thetag, vectorField &force) const |
Public Member Functions inherited from trimModel | |
TypeName ("trimModel") | |
Run-time type information. More... | |
declareRunTimeSelectionTable (autoPtr, trimModel, dictionary,(const fv::rotorDisk &rotor, const dictionary &dict),(rotor, dict)) | |
trimModel (const fv::rotorDisk &rotor, const dictionary &dict, const word &name) | |
Construct from components. More... | |
virtual | ~trimModel () |
Destructor. More... | |
Protected Member Functions | |
template<class RhoFieldType > | |
vector | calcCoeffs (const RhoFieldType &rho, const vectorField &U, const scalarField &alphag, vectorField &force) const |
Calculate the rotor force and moment coefficients vector. More... | |
template<class RhoFieldType > | |
void | correctTrim (const RhoFieldType &rho, const vectorField &U, vectorField &force) const |
Correct the model. More... | |
Protected Attributes | |
label | calcFrequency_ |
Number of iterations between calls to 'correct'. More... | |
bool | useCoeffs_ |
Flag to indicate whether to solve coeffs (true) or forces (false) More... | |
vector | target_ |
Target coefficient vector (thrust force, roll moment, pitch moment) More... | |
vector | theta_ |
Pitch angles (collective, roll, pitch) [rad]. More... | |
label | nIter_ |
Maximum number of iterations in trim routine. More... | |
scalar | tol_ |
Convergence tolerance. More... | |
scalar | relax_ |
Under-relaxation coefficient. More... | |
scalar | dTheta_ |
Perturbation angle used to determine jacobian. More... | |
scalar | alpha_ |
Coefficient to allow for conversion between US and EU definitions. More... | |
Protected Attributes inherited from trimModel | |
const fv::rotorDisk & | rotor_ |
Reference to the rotor source model. More... | |
const word | name_ |
Name of model. More... | |
dictionary | coeffs_ |
Coefficients dictionary. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from trimModel | |
static autoPtr< trimModel > | New (const fv::rotorDisk &rotor, const dictionary &dict) |
Return a reference to the selected trim model. More... | |
Target trim forces/coefficients.
Solves:
c^old + J.d(theta) = c^target
Where:
n = time level c = coefficient vector (thrust force, roll moment, pitch moment) theta = pitch angle vector (collective, roll, pitch) J = Jacobian [3x3] matrix
The trimmed pitch angles are found via solving the above with a Newton-Raphson iterative method. The solver tolerance can be user-input, using the 'tol' entry.
If coefficients are requested (useCoeffs = true), the force and moments are normalised using:
force c = ——————————— alpha*pi*rho*(omega^2)*(radius^4)
and
moment c = ——————————— alpha*pi*rho*(omega^2)*(radius^5)
Where:
alpha = user-input conversion coefficient (default = 1) rho = density omega = rotor angular velocity pi = mathematical pi
Definition at line 87 of file targetCoeffTrim.H.
targetCoeffTrim | ( | const fv::rotorDisk & | rotor, |
const dictionary & | dict | ||
) |
Constructor.
Definition at line 194 of file targetCoeffTrim.C.
References dict, and targetCoeffTrim::read().
|
virtual |
Destructor.
Definition at line 217 of file targetCoeffTrim.C.
|
protected |
Calculate the rotor force and moment coefficients vector.
|
protected |
Correct the model.
Definition at line 100 of file targetCoeffTrim.C.
References calcType, Foam::endl(), Foam::Info, Foam::inv(), Foam::mag(), Foam::nl, Foam::radToDeg(), rho, Foam::type(), U, and Foam::Zero.
TypeName | ( | "targetCoeffTrim" | ) |
Run-time type information.
|
virtual |
Read.
Reimplemented from trimModel.
Definition at line 223 of file targetCoeffTrim.C.
References dict, dictionary::lookup(), dictionary::lookupOrDefault(), trimModel::read(), and Foam::unitDegrees.
Referenced by targetCoeffTrim::targetCoeffTrim().
|
virtual |
Return the geometric angle of attack [rad].
Implements trimModel.
Definition at line 256 of file targetCoeffTrim.C.
References Foam::cos(), forAll, psi, tmp< T >::ref(), Foam::sin(), and x.
|
virtual |
Correct the model.
Implements trimModel.
Definition at line 273 of file targetCoeffTrim.C.
References U.
|
virtual |
Correct the model for compressible flow.
Implements trimModel.
Definition at line 283 of file targetCoeffTrim.C.
Foam::vector calcCoeffs | ( | const RhoFieldType & | rho, |
const vectorField & | U, | ||
const scalarField & | thetag, | ||
vectorField & | force | ||
) | const |
Definition at line 46 of file targetCoeffTrim.C.
References C::C(), cells, forAll, Foam::constant::mathematical::pi(), Foam::pow4(), Foam::reduce(), rho, Foam::sqr(), U, x, and Foam::Zero.
|
protected |
Number of iterations between calls to 'correct'.
Definition at line 97 of file targetCoeffTrim.H.
|
protected |
Flag to indicate whether to solve coeffs (true) or forces (false)
Definition at line 100 of file targetCoeffTrim.H.
|
protected |
Target coefficient vector (thrust force, roll moment, pitch moment)
Definition at line 103 of file targetCoeffTrim.H.
|
mutableprotected |
Pitch angles (collective, roll, pitch) [rad].
Definition at line 106 of file targetCoeffTrim.H.
|
protected |
Maximum number of iterations in trim routine.
Definition at line 109 of file targetCoeffTrim.H.
|
protected |
Convergence tolerance.
Definition at line 112 of file targetCoeffTrim.H.
|
protected |
Under-relaxation coefficient.
Definition at line 115 of file targetCoeffTrim.H.
|
protected |
Perturbation angle used to determine jacobian.
Definition at line 118 of file targetCoeffTrim.H.
|
protected |
Coefficient to allow for conversion between US and EU definitions.
Definition at line 121 of file targetCoeffTrim.H.