targetCoeffTrim.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2012-2014 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::targetCoeffTrim
26 
27 Description
28  Target trim forces/coefficients
29 
30  Solves:
31 
32  c^old + J.d(theta) = c^target
33 
34  Where:
35 
36  n = time level
37  c = coefficient vector (thrust force, roll moment, pitch moment)
38  theta = pitch angle vector (collective, roll, pitch)
39  J = Jacobian [3x3] matrix
40 
41 
42  The trimmed pitch angles are found via solving the above with a
43  Newton-Raphson iterative method. The solver tolerance can be user-input,
44  using the 'tol' entry.
45 
46  If coefficients are requested (useCoeffs = true), the force and moments
47  are normalised using:
48 
49  force
50  c = ---------------------------------
51  alpha*pi*rho*(omega^2)*(radius^4)
52 
53  and
54 
55  moment
56  c = ---------------------------------
57  alpha*pi*rho*(omega^2)*(radius^5)
58 
59  Where:
60 
61  alpha = user-input conversion coefficient (default = 1)
62  rho = desity
63  omega = rotor angulr velocity
64  pi = mathematical pi
65 
66 
67 SourceFiles
68  targetCoeffTrim.C
69 
70 \*---------------------------------------------------------------------------*/
71 
72 #ifndef targetCoeffTrim_H
73 #define targetCoeffTrim_H
74 
75 #include "trimModel.H"
76 #include "tensor.H"
77 #include "vector.H"
78 
79 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80 
81 namespace Foam
82 {
83 
84 /*---------------------------------------------------------------------------*\
85  Class targetCoeffTrim Declaration
86 \*---------------------------------------------------------------------------*/
87 
88 class targetCoeffTrim
89 :
90  public trimModel
91 {
92 
93 protected:
94 
95  // Protected data
96 
97  //- Number of iterations between calls to 'correct'
99 
100  //- Flag to indicate whether to solve coeffs (true) or forces (false)
101  bool useCoeffs_;
102 
103  //- Target coefficient vector (thrust force, roll moment, pitch moment)
104  vector target_;
105 
106  //- Pitch angles (collective, roll, pitch) [rad]
107  vector theta_;
108 
109  //- Maximum number of iterations in trim routine
110  label nIter_;
111 
112  //- Convergence tolerance
113  scalar tol_;
114 
115  //- Under-relaxation coefficient
116  scalar relax_;
117 
118  //- Perturbation angle used to determine jacobian
119  scalar dTheta_;
120 
121  //- Coefficient to allow for conversion between US and EU definitions
122  scalar alpha_;
123 
124 
125  // Protected member functions
126 
127  //- Calculate the rotor force and moment coefficients vector
128  template<class RhoFieldType>
130  (
131  const RhoFieldType& rho,
132  const vectorField& U,
133  const scalarField& alphag,
134  vectorField& force
135  ) const;
136 
137  //- Correct the model
138  template<class RhoFieldType>
139  void correctTrim
140  (
141  const RhoFieldType& rho,
142  const vectorField& U,
143  vectorField& force
144  );
145 
146 
147 public:
148 
149  //- Run-time type information
150  TypeName("targetCoeffTrim");
151 
152  //- Constructor
153  targetCoeffTrim(const fv::rotorDiskSource& rotor, const dictionary& dict);
154 
155  //- Destructor
156  virtual ~targetCoeffTrim();
157 
158 
159  // Member functions
160 
161  //- Read
162  void read(const dictionary& dict);
163 
164  //- Return the geometric angle of attack [rad]
165  virtual tmp<scalarField> thetag() const;
166 
167  //- Correct the model
168  virtual void correct
169  (
170  const vectorField& U,
171  vectorField& force
172  );
173 
174  //- Correct the model for compressible flow
175  virtual void correct
176  (
177  const volScalarField rho,
178  const vectorField& U,
179  vectorField& force
180  );
181 };
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 } // End namespace Foam
187 
188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 
190 #endif
191 
192 // ************************************************************************* //
dictionary dict
U
Definition: pEqn.H:83
virtual void correct(const vectorField &U, vectorField &force)
Correct the model.
scalar alpha_
Coefficient to allow for conversion between US and EU definitions.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void correctTrim(const RhoFieldType &rho, const vectorField &U, vectorField &force)
Correct the model.
Trim model base class.
Definition: trimModel.H:50
scalar tol_
Convergence tolerance.
label nIter_
Maximum number of iterations in trim routine.
vector theta_
Pitch angles (collective, roll, pitch) [rad].
vector calcCoeffs(const RhoFieldType &rho, const vectorField &U, const scalarField &alphag, vectorField &force) const
Calculate the rotor force and moment coefficients vector.
Rotor disk source.
virtual tmp< scalarField > thetag() const
Return the geometric angle of attack [rad].
virtual ~targetCoeffTrim()
Destructor.
Target trim forces/coefficients.
vector target_
Target coefficient vector (thrust force, roll moment, pitch moment)
label calcFrequency_
Number of iterations between calls to &#39;correct&#39;.
TypeName("targetCoeffTrim")
Run-time type information.
scalar dTheta_
Perturbation angle used to determine jacobian.
targetCoeffTrim(const fv::rotorDiskSource &rotor, const dictionary &dict)
Constructor.
void read(const dictionary &dict)
Read.
scalar relax_
Under-relaxation coefficient.
A class for managing temporary objects.
Definition: PtrList.H:54
bool useCoeffs_
Flag to indicate whether to solve coeffs (true) or forces (false)
Namespace for OpenFOAM.