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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2012-2014 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::targetCoeffTrim
29
30Description
31 Trim model where the operating characteristics of rotor
32 (e.g. blade pitch angle) can vary to reach a specified
33 target thrust and torque.
34
35 Solves:
36
37 \f[
38 c^{old} + J.d(\theta) = c^{target}
39 \f]
40
41 where
42 \vartable
43 n | Time level
44 c | Coefficient vector (thrust force, roll moment, pitch moment)
45 \theta | Pitch angle vector (collective, roll, pitch)
46 J | Jacobian [3x3] matrix
47 \endvartable
48
49 The trimmed pitch angles are found via solving the above with a
50 Newton-Raphson iterative method. The solver tolerance can be user-input,
51 using the 'tol' entry.
52
53 If coefficients are requested (useCoeffs = true),
54 the force and moments are normalised using:
55
56 \f[
57 c = \frac{F}{\alpha \pi \rho \omega^2 R^4}
58 \f]
59
60 and
61
62 \f[
63 c = \frac{M}{\alpha \pi \rho \omega^2 R^5}
64 \f]
65
66 where
67 \vartable
68 F | Force
69 M | Moment
70 \alpha | User-input conversion coefficient (default = 1)
71 \rho | Fluid desity
72 \omega | Rotor angular velocity
73 \pi | Mathematical pi
74 R | Rotor radius
75 \endvartable
76
77Usage
78 Minimal example by using \c constant/fvOptions:
79 rotorDiskSource1
80 {
81 // Mandatory/Optional (inherited) entries
82 ...
83
84 // Mandatory entries (runtime modifiable)
85 trimModel targetForce;
86
87 targetForceCoeffs
88 {
89 // Conditional mandatory entries (runtime modifiable)
90
91 // when trimModel=targetForce
92 target
93 {
94 // Mandatory entries (runtime modifiable)
95 thrust 0.5;
96 pitch 0.5;
97 roll 0.5;
98
99 // Optional entries (runtime modifiable)
100 useCoeffs true;
101 }
102
103 pitchAngles
104 {
105 theta0Ini 0.1;
106 theta1cIni 0.1;
107 theta1sIni 0.1;
108 }
109
110 calcFrequency 1;
111
112 // Optional entries (runtime modifiable)
113 nIter 50;
114 tol 1e-8;
115 relax 1;
116 dTheta 0.1;
117 alpha 1.0;
118 }
119 }
120
121 where the entries mean:
122 \table
123 Property | Description | Type | Reqd | Dflt
124 calcFrequency | Number of iterations between calls to 'correct' <!--
125 --> | label | yes | -
126 useCoeffs | Flag to indicate whether to solve coeffs <!--
127 --> (true) or forces (false) | bool | no | true
128 thrust | Target thrust (coefficient) | scalar | yes | -
129 pitch | Target pitch (coefficient) | scalar | yes | -
130 roll | Target roll moment (coefficient) | scalar | yes | -
131 theta0Ini | Initial pitch angle [deg] | scalar | yes | -
132 theta1cIni | Initial lateral pitch angle (cos coeff) [deg] <!--
133 --> | scalar | yes | -
134 theta1sIni | Initial longitudinal pitch angle (sin coeff) [deg] <!--
135 --> | scalar | yes | -
136 nIter | Maximum number of iterations in trim routine <!--
137 --> | label | no | 50
138 tol | Convergence tolerance in trim routine | scalar | no | 1e-8
139 relax | Relaxation factor in trim routine | scalar | no | 1
140 dTheta | Perturbation angle used to determine Jacobian [deg] <!--
141 --> | scalar | no | 0.1
142 alpha | Coefficient to allow for conversion between US <!--
143 --> and EU definitions | scalar | no | 1
144 \endtable
145
146See also
147 - Foam::fv::rotorDiskSource
148 - Foam::trimModel
149 - Foam::fixedTrim
150
151SourceFiles
152 targetCoeffTrim.C
153
154\*---------------------------------------------------------------------------*/
155
156#ifndef targetCoeffTrim_H
157#define targetCoeffTrim_H
158
159#include "trimModel.H"
160#include "tensor.H"
161#include "vector.H"
162
163// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164
165namespace Foam
166{
167
168/*---------------------------------------------------------------------------*\
169 Class targetCoeffTrim Declaration
170\*---------------------------------------------------------------------------*/
171
172class targetCoeffTrim
173:
174 public trimModel
175{
176protected:
177
178 // Protected Data
179
180 //- Number of iterations between calls to 'correct'
181 label calcFrequency_;
182
183 //- Flag to indicate whether to solve coeffs (true) or forces (false)
184 bool useCoeffs_;
185
186 //- Target coefficient vector (thrust force, roll moment, pitch moment)
188
189 //- Pitch angles (collective, roll, pitch) [rad]
191
192 //- Maximum number of iterations in trim routine
193 label nIter_;
194
195 //- Convergence tolerance
196 scalar tol_;
197
198 //- Under-relaxation coefficient
199 scalar relax_;
200
201 //- Perturbation angle used to determine jacobian
202 scalar dTheta_;
203
204 //- Coefficient to allow for conversion between US and EU definitions
205 scalar alpha_;
206
207
208 // Protected Member Functions
209
210 //- Calculate the rotor force and moment coefficients vector
211 template<class RhoFieldType>
213 (
214 const RhoFieldType& rho,
215 const vectorField& U,
216 const scalarField& alphag,
217 vectorField& force
218 ) const;
219
220 //- Correct the model
221 template<class RhoFieldType>
222 void correctTrim
223 (
224 const RhoFieldType& rho,
225 const vectorField& U,
226 vectorField& force
227 );
228
229
230public:
231
232 //- Run-time type information
233 TypeName("targetCoeffTrim");
234
235 // Constructors
236
237 //- Constructor from rotor and dictionary
239 (
240 const fv::rotorDiskSource& rotor,
241 const dictionary& dict
242 );
243
244 //- No copy construct
245 targetCoeffTrim(const targetCoeffTrim&) = delete;
246
247 //- No copy assignment
248 void operator=(const targetCoeffTrim&) = delete;
249
250
251 //- Destructor
252 virtual ~targetCoeffTrim() = default;
253
254
255 // Member Functions
256
257 //- Read
258 void read(const dictionary& dict);
259
260 //- Return the geometric angle of attack [rad]
261 virtual tmp<scalarField> thetag() const;
262
263 //- Correct the model
264 virtual void correct
265 (
266 const vectorField& U,
267 vectorField& force
268 );
269
270 //- Correct the model for compressible flow
271 virtual void correct
272 (
273 const volScalarField rho,
274 const vectorField& U,
275 vectorField& force
276 );
277};
278
279
280// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281
282} // End namespace Foam
283
284// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285
286#endif
287
288// ************************************************************************* //
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Applies cell-based momentum sources on velocity (i.e. U) within a specified cylindrical region to app...
Trim model where the operating characteristics of rotor (e.g. blade pitch angle) can vary to reach a ...
label nIter_
Maximum number of iterations in trim routine.
label calcFrequency_
Number of iterations between calls to 'correct'.
targetCoeffTrim(const fv::rotorDiskSource &rotor, const dictionary &dict)
Constructor from rotor and dictionary.
virtual tmp< scalarField > thetag() const
Return the geometric angle of attack [rad].
vector theta_
Pitch angles (collective, roll, pitch) [rad].
vector target_
Target coefficient vector (thrust force, roll moment, pitch moment)
scalar dTheta_
Perturbation angle used to determine jacobian.
void read(const dictionary &dict)
Read.
TypeName("targetCoeffTrim")
Run-time type information.
bool useCoeffs_
Flag to indicate whether to solve coeffs (true) or forces (false)
virtual ~targetCoeffTrim()=default
Destructor.
scalar relax_
Under-relaxation coefficient.
targetCoeffTrim(const targetCoeffTrim &)=delete
No copy construct.
scalar alpha_
Coefficient to allow for conversion between US and EU definitions.
void operator=(const targetCoeffTrim &)=delete
No copy assignment.
scalar tol_
Convergence tolerance.
void correctTrim(const RhoFieldType &rho, const vectorField &U, vectorField &force)
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.
Definition: tmp.H:65
Base class for trim models for handling blade characteristics and thrust-torque relations.
Definition: trimModel.H:112
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
U
Definition: pEqn.H:72
thermo correct()
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73