Trim model where the operating characteristics of rotor (e.g. blade pitch angle) can vary to reach a specified target thrust and torque. More...
#include <targetCoeffTrim.H>


Public Member Functions | |
| TypeName ("targetCoeffTrim") | |
| Run-time type information. | |
| targetCoeffTrim (const fv::rotorDiskSource &rotor, const dictionary &dict) | |
| Constructor from rotor and dictionary. | |
| targetCoeffTrim (const targetCoeffTrim &)=delete | |
| No copy construct. | |
| void | operator= (const targetCoeffTrim &)=delete |
| No copy assignment. | |
| virtual | ~targetCoeffTrim ()=default |
| Destructor. | |
| void | read (const dictionary &dict) |
| Read. | |
| virtual tmp< scalarField > | thetag () const |
| Return the geometric angle of attack [rad]. | |
| virtual void | correct (const vectorField &U, vectorField &force) |
| Correct the model. | |
| virtual void | correct (const volScalarField rho, const vectorField &U, vectorField &force) |
| Correct the model for compressible flow. | |
| 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. | |
| declareRunTimeSelectionTable (autoPtr, trimModel, dictionary,(const fv::rotorDiskSource &rotor, const dictionary &dict),(rotor, dict)) | |
| trimModel (const fv::rotorDiskSource &rotor, const dictionary &dict, const word &name) | |
| Construct from components. | |
| virtual | ~trimModel ()=default |
| Destructor. | |
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. | |
| template<class RhoFieldType> | |
| void | correctTrim (const RhoFieldType &rho, const vectorField &U, vectorField &force) |
| Correct the model. | |
Protected Attributes | |
| label | calcFrequency_ |
| Number of iterations between calls to 'correct'. | |
| bool | useCoeffs_ |
| Flag to indicate whether to solve coeffs (true) or forces (false). | |
| vector | target_ |
| Target coefficient vector (thrust force, roll moment, pitch moment). | |
| vector | theta_ |
| Pitch angles (collective, roll, pitch) [rad]. | |
| label | nIter_ |
| Maximum number of iterations in trim routine. | |
| scalar | tol_ |
| Convergence tolerance. | |
| scalar | relax_ |
| Under-relaxation coefficient. | |
| scalar | dTheta_ |
| Perturbation angle used to determine jacobian. | |
| scalar | alpha_ |
| Coefficient to allow for conversion between US and EU definitions. | |
| Protected Attributes inherited from trimModel | |
| const fv::rotorDiskSource & | rotor_ |
| Reference to the rotor source model. | |
| const word | name_ |
| Name of model. | |
| dictionary | coeffs_ |
| Coefficients dictionary. | |
Additional Inherited Members | |
| Static Public Member Functions inherited from trimModel | |
| static autoPtr< trimModel > | New (const fv::rotorDiskSource &rotor, const dictionary &dict) |
| Return a reference to the selected trim model. | |
Trim model where the operating characteristics of rotor (e.g. blade pitch angle) can vary to reach a specified target thrust and torque.
Solves:
![\[ c^{old} + J.d(\theta) = c^{target}
\]](form_475.png)
where
![]() | = | Time level |
![]() | = | Coefficient vector (thrust force, roll moment, pitch moment) |
![]() | = | Pitch angle vector (collective, roll, pitch) |
![]() | = | 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:
![\[ c = \frac{F}{\alpha \pi \rho \omega^2 R^4}
\]](form_480.png)
and
![\[ c = \frac{M}{\alpha \pi \rho \omega^2 R^5}
\]](form_481.png)
where
![]() | = | Force |
![]() | = | Moment |
![]() | = | User-input conversion coefficient (default = 1) |
![]() | = | Fluid desity |
![]() | = | Rotor angular velocity |
![]() | = | Mathematical pi |
![]() | = | Rotor radius |
constant/fvOptions: rotorDiskSource1 { Mandatory/Optional (inherited) entries ...Mandatory entries (runtime modifiable) trimModel targetForce;
targetForceCoeffs { Conditional mandatory entries (runtime modifiable)
when trimModel=targetForce target { Mandatory entries (runtime modifiable) thrust 0.5; pitch 0.5; roll 0.5;
Optional entries (runtime modifiable) useCoeffs true; }
pitchAngles { theta0Ini 0.1; theta1cIni 0.1; theta1sIni 0.1; }
calcFrequency 1;
Optional entries (runtime modifiable) nIter 50; tol 1e-8; relax 1; dTheta 0.1; alpha 1.0; } }
where the entries mean:
| Property | Description | Type | Reqd | Dflt |
|---|---|---|---|---|
calcFrequency | Number of iterations between calls to 'correct' | label | yes | - |
useCoeffs | Flag to indicate whether to solve coeffs (true) or forces (false) | bool | no | true |
thrust | Target thrust (coefficient) | scalar | yes | - |
pitch | Target pitch (coefficient) | scalar | yes | - |
roll | Target roll moment (coefficient) | scalar | yes | - |
theta0Ini | Initial pitch angle [deg] | scalar | yes | - |
theta1cIni | Initial lateral pitch angle (cos coeff) [deg] | scalar | yes | - |
theta1sIni | Initial longitudinal pitch angle (sin coeff) [deg] | scalar | yes | - |
nIter | Maximum number of iterations in trim routine | label | no | 50 |
tol | Convergence tolerance in trim routine | scalar | no | 1e-8 |
relax | Relaxation factor in trim routine | scalar | no | 1 |
dTheta | Perturbation angle used to determine Jacobian [deg] | scalar | no | 0.1 |
alpha | Coefficient to allow for conversion between US and EU definitions | scalar | no | 1 |
Definition at line 309 of file targetCoeffTrim.H.
| targetCoeffTrim | ( | const fv::rotorDiskSource & | rotor, |
| const dictionary & | dict ) |
Constructor from rotor and dictionary.
Definition at line 188 of file targetCoeffTrim.C.
References alpha_, calcFrequency_, Foam::degToRad(), dict, dTheta_, e, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, nIter_, read(), relax_, target_, theta_, tol_, trimModel::trimModel(), useCoeffs_, and Foam::Zero.
Referenced by operator=(), and targetCoeffTrim().


|
delete |
|
virtualdefault |
|
protected |
Calculate the rotor force and moment coefficients vector.
Referenced by correctTrim().

|
protected |
Correct the model.
Definition at line 94 of file targetCoeffTrim.C.
References calcCoeffs(), calcFrequency_, dTheta_, Foam::endl(), Foam::Info, Foam::inv(), Foam::mag(), nIter_, Foam::nl, Foam::radToDeg(), relax_, rho, trimModel::rotor_, target_, theta_, thetag(), tol_, Foam::type(), U, useCoeffs_, and Foam::Zero.
Referenced by correct(), and correct().


| TypeName | ( | "targetCoeffTrim" | ) |
Run-time type information.
References dict.
|
delete |
|
virtual |
Read.
Reimplemented from trimModel.
Definition at line 211 of file targetCoeffTrim.C.
References alpha_, calcFrequency_, trimModel::coeffs_, Foam::degToRad(), dict, dTheta_, dictionary::get(), dictionary::getOrDefault(), nIter_, trimModel::read(), dictionary::readEntry(), relax_, target_, theta_, tol_, and useCoeffs_.
Referenced by targetCoeffTrim().


|
virtual |
Return the geometric angle of attack [rad].
Implements trimModel.
Definition at line 247 of file targetCoeffTrim.C.
References Foam::cos(), forAll, tmp< T >::New(), psi, trimModel::rotor_, Foam::sin(), theta_, and x.
Referenced by calcCoeffs(), and correctTrim().


|
virtual |
Correct the model.
Implements trimModel.
Definition at line 264 of file targetCoeffTrim.C.
References correctTrim(), and U.

|
virtual |
Correct the model for compressible flow.
Implements trimModel.
Definition at line 274 of file targetCoeffTrim.C.
References correctTrim(), rho, and U.

| Foam::vector calcCoeffs | ( | const RhoFieldType & | rho, |
| const vectorField & | U, | ||
| const scalarField & | thetag, | ||
| vectorField & | force ) const |
Definition at line 40 of file targetCoeffTrim.C.
References alpha_, cells, forAll, Foam::constant::mathematical::pi(), Foam::pow4(), Foam::reduce(), rho, trimModel::rotor_, Foam::sqr(), thetag(), U, useCoeffs_, x, and Foam::Zero.

|
protected |
Number of iterations between calls to 'correct'.
Definition at line 320 of file targetCoeffTrim.H.
Referenced by correctTrim(), read(), and targetCoeffTrim().
|
protected |
Flag to indicate whether to solve coeffs (true) or forces (false).
Definition at line 325 of file targetCoeffTrim.H.
Referenced by calcCoeffs(), correctTrim(), read(), and targetCoeffTrim().
|
protected |
Target coefficient vector (thrust force, roll moment, pitch moment).
Definition at line 330 of file targetCoeffTrim.H.
Referenced by correctTrim(), read(), and targetCoeffTrim().
|
protected |
Pitch angles (collective, roll, pitch) [rad].
Definition at line 335 of file targetCoeffTrim.H.
Referenced by correctTrim(), read(), targetCoeffTrim(), and thetag().
|
protected |
Maximum number of iterations in trim routine.
Definition at line 340 of file targetCoeffTrim.H.
Referenced by correctTrim(), read(), and targetCoeffTrim().
|
protected |
Convergence tolerance.
Definition at line 345 of file targetCoeffTrim.H.
Referenced by correctTrim(), read(), and targetCoeffTrim().
|
protected |
Under-relaxation coefficient.
Definition at line 350 of file targetCoeffTrim.H.
Referenced by correctTrim(), read(), and targetCoeffTrim().
|
protected |
Perturbation angle used to determine jacobian.
Definition at line 355 of file targetCoeffTrim.H.
Referenced by correctTrim(), read(), and targetCoeffTrim().
|
protected |
Coefficient to allow for conversion between US and EU definitions.
Definition at line 360 of file targetCoeffTrim.H.
Referenced by calcCoeffs(), read(), and targetCoeffTrim().