targetCoeffTrim.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2024 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 \*---------------------------------------------------------------------------*/
25 
26 #include "targetCoeffTrim.H"
27 #include "rotorDisk.H"
28 #include "geometricOneField.H"
30 
31 using namespace Foam::constant;
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
38 
40 }
41 
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
45 template<class RhoFieldType>
47 (
48  const RhoFieldType& rho,
49  const vectorField& U,
50  const scalarField& thetag,
51  vectorField& force
52 ) const
53 {
54  rotor_.calculate(rho, U, thetag, force, false, false);
55 
56  const labelUList cells = rotor_.set().cells();
57  const vectorField& C = rotor_.mesh().C();
58  const List<point>& x = rotor_.x();
59 
60  const vector& origin = rotor_.coordSys().origin();
61  const vector& rollAxis = rotor_.coordSys().R().e1();
62  const vector& pitchAxis = rotor_.coordSys().R().e2();
63  const vector& yawAxis = rotor_.coordSys().R().e3();
64 
65  scalar coeff1 = alpha_*sqr(rotor_.omega())*mathematical::pi;
66 
67  vector cf(Zero);
68  forAll(cells, i)
69  {
70  label celli = cells[i];
71 
72  vector fc = force[celli];
73  vector mc = fc^(C[celli] - origin);
74 
75  if (useCoeffs_)
76  {
77  scalar radius = x[i].x();
78  scalar coeff2 = rho[celli]*coeff1*pow4(radius);
79 
80  // add to coefficient vector
81  cf[0] += (fc & yawAxis)/(coeff2 + rootVSmall);
82  cf[1] += (mc & pitchAxis)/(coeff2*radius + rootVSmall);
83  cf[2] += (mc & rollAxis)/(coeff2*radius + rootVSmall);
84  }
85  else
86  {
87  cf[0] += fc & yawAxis;
88  cf[1] += mc & pitchAxis;
89  cf[2] += mc & rollAxis;
90  }
91  }
92 
93  reduce(cf, sumOp<vector>());
94 
95  return cf;
96 }
97 
98 
99 template<class RhoFieldType>
101 (
102  const RhoFieldType& rho,
103  const vectorField& U,
104  vectorField& force
105 ) const
106 {
107  if (rotor_.mesh().time().timeIndex() % calcFrequency_ == 0)
108  {
109  word calcType = "forces";
110  if (useCoeffs_)
111  {
112  calcType = "coefficients";
113  }
114 
115  Info<< type() << ":" << nl
116  << " solving for target trim " << calcType << nl;
117 
118  const scalar rhoRef = rotor_.rhoRef();
119 
120  // iterate to find new pitch angles to achieve target force
121  scalar err = great;
122  label iter = 0;
123  tensor J(Zero);
124 
125  vector old = Zero;
126  while ((err > tol_) && (iter < nIter_))
127  {
128  // cache initial theta vector
129  vector theta0(theta_);
130 
131  // set initial values
132  old = calcCoeffs(rho, U, thetag(), force);
133 
134  // construct Jacobian by perturbing the pitch angles
135  // by +/-(dTheta_/2)
136  for (label pitchI = 0; pitchI < 3; pitchI++)
137  {
138  theta_[pitchI] -= dTheta_/2.0;
139  vector cf0 = calcCoeffs(rho, U, thetag(), force);
140 
141  theta_[pitchI] += dTheta_;
142  vector cf1 = calcCoeffs(rho, U, thetag(), force);
143 
144  vector ddTheta = (cf1 - cf0)/dTheta_;
145 
146  J[pitchI + 0] = ddTheta[0];
147  J[pitchI + 3] = ddTheta[1];
148  J[pitchI + 6] = ddTheta[2];
149 
150  theta_ = theta0;
151  }
152 
153  // calculate the change in pitch angle vector
154  vector dt = inv(J) & (target_/rhoRef - old);
155 
156  // update pitch angles
157  vector thetaNew = theta_ + relax_*dt;
158 
159  // update error
160  err = mag(thetaNew - theta_);
161 
162  // update for next iteration
163  theta_ = thetaNew;
164  iter++;
165  }
166 
167  if (iter == nIter_)
168  {
169  Info<< " solution not converged in " << iter
170  << " iterations, final residual = " << err
171  << "(" << tol_ << ")" << endl;
172  }
173  else
174  {
175  Info<< " final residual = " << err << "(" << tol_
176  << "), iterations = " << iter << endl;
177  }
178 
179  Info<< " current and target " << calcType << nl
180  << " thrust = " << old[0]*rhoRef << ", " << target_[0] << nl
181  << " pitch = " << old[1]*rhoRef << ", " << target_[1] << nl
182  << " roll = " << old[2]*rhoRef << ", " << target_[2] << nl
183  << " new pitch angles [deg]:" << nl
184  << " theta0 = " << radToDeg(theta_[0]) << nl
185  << " theta1c = " << radToDeg(theta_[1]) << nl
186  << " theta1s = " << radToDeg(theta_[2]) << nl
187  << endl;
188  }
189 }
190 
191 
192 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
193 
195 (
196  const fv::rotorDisk& rotor,
197  const dictionary& dict
198 )
199 :
200  trimModel(rotor, dict, typeName),
201  calcFrequency_(-1),
202  useCoeffs_(true),
203  target_(Zero),
204  theta_(Zero),
205  nIter_(50),
206  tol_(1e-8),
207  relax_(1.0),
208  dTheta_(degToRad(0.1)),
209  alpha_(1.0)
210 {
211  read(dict);
212 }
213 
214 
215 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
216 
218 {}
219 
220 
221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222 
224 {
226 
227  const dictionary& targetDict(coeffs_.subDict("target"));
228  useCoeffs_ = targetDict.lookupOrDefault<bool>("useCoeffs", true);
229  word ext = "";
230  if (useCoeffs_)
231  {
232  ext = "Coeff";
233  }
234 
235  target_[0] = targetDict.lookup<scalar>("thrust" + ext);
236  target_[1] = targetDict.lookup<scalar>("pitch" + ext);
237  target_[2] = targetDict.lookup<scalar>("roll" + ext);
238 
239  const dictionary& pitchAngleDict(coeffs_.subDict("pitchAngles"));
240  theta_[0] = pitchAngleDict.lookup<scalar>("theta0Ini", unitDegrees);
241  theta_[1] = pitchAngleDict.lookup<scalar>("theta1cIni", unitDegrees);
242  theta_[2] = pitchAngleDict.lookup<scalar>("theta1sIni", unitDegrees);
243 
244  coeffs_.lookup("calcFrequency") >> calcFrequency_;
245 
246  coeffs_.readIfPresent("nIter", nIter_);
247  coeffs_.readIfPresent("tol", tol_);
248  coeffs_.readIfPresent("relax", relax_);
249 
250  coeffs_.readIfPresent("dTheta", unitDegrees, dTheta_);
251 
252  alpha_ = coeffs_.lookup<scalar>("alpha");
253 }
254 
255 
257 {
258  const List<vector>& x = rotor_.x();
259 
260  tmp<scalarField> ttheta(new scalarField(x.size()));
261  scalarField& t = ttheta.ref();
262 
263  forAll(t, i)
264  {
265  scalar psi = x[i].y();
266  t[i] = theta_[0] + theta_[1]*cos(psi) + theta_[2]*sin(psi);
267  }
268 
269  return ttheta;
270 }
271 
272 
274 (
275  const vectorField& U,
276  vectorField& force
277 ) const
278 {
279  correctTrim(geometricOneField(), U, force);
280 }
281 
282 
284 (
285  const volScalarField rho,
286  const vectorField& U,
287  vectorField& force
288 ) const
289 {
290  correctTrim(rho, U, force);
291 }
292 
293 
294 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Graphite solid properties.
Definition: C.H:51
C()
Construct null.
Definition: C.C:40
Generic GeometricField class.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
Cell based momentum source which approximates the mean effects of rotor forces on a cylindrical regio...
Definition: rotorDisk.H:125
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Target trim forces/coefficients.
virtual ~targetCoeffTrim()
Destructor.
virtual tmp< scalarField > thetag() const
Return the geometric angle of attack [rad].
targetCoeffTrim(const fv::rotorDisk &rotor, const dictionary &dict)
Constructor.
void read(const dictionary &dict)
Read.
virtual void correct(const vectorField &U, vectorField &force) const
Correct the model.
void correctTrim(const RhoFieldType &rho, 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.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Trim model base class.
Definition: trimModel.H:54
virtual void read(const dictionary &dict)
Read.
Definition: trimModel.C:61
A class for handling words, derived from string.
Definition: word.H:62
#define calcType(Type, GeoField)
const cellShapeList & cells
const volScalarField & psi
U
Definition: pEqn.H:72
Collection of constants.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:106
scalar radToDeg(const scalar rad)
Convert radians to degrees.
scalar degToRad(const scalar deg)
Convert degrees to radians.
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
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
static const char nl
Definition: Ostream.H:266
dimensionedScalar cos(const dimensionedScalar &ds)
const unitConversion unitDegrees
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict