targetCoeffTrim.C
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-2016 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "targetCoeffTrim.H"
30 #include "geometricOneField.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(targetCoeffTrim, 0);
40 
41  addToRunTimeSelectionTable(trimModel, targetCoeffTrim, dictionary);
42 }
43 
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
47 template<class RhoFieldType>
49 (
50  const RhoFieldType& rho,
51  const vectorField& U,
52  const scalarField& thetag,
53  vectorField& force
54 ) const
55 {
56  rotor_.calculate(rho, U, thetag, force, false, false);
57 
58  const labelList& cells = rotor_.cells();
59  const vectorField& C = rotor_.mesh().C();
60  const List<point>& x = rotor_.x();
61 
62  const vector& origin = rotor_.coordSys().origin();
63  const vector& rollAxis = rotor_.coordSys().e1();
64  const vector& pitchAxis = rotor_.coordSys().e2();
65  const vector& yawAxis = rotor_.coordSys().e3();
66 
67  scalar coeff1 = alpha_*sqr(rotor_.omega())*mathematical::pi;
68 
69  vector cf(Zero);
70  forAll(cells, i)
71  {
72  label celli = cells[i];
73 
74  vector fc = force[celli];
75  vector mc = fc^(C[celli] - origin);
76 
77  if (useCoeffs_)
78  {
79  scalar radius = x[i].x();
80  scalar coeff2 = rho[celli]*coeff1*pow4(radius);
81 
82  // add to coefficient vector
83  cf[0] += (fc & yawAxis)/(coeff2 + ROOTVSMALL);
84  cf[1] += (mc & pitchAxis)/(coeff2*radius + ROOTVSMALL);
85  cf[2] += (mc & rollAxis)/(coeff2*radius + ROOTVSMALL);
86  }
87  else
88  {
89  cf[0] += fc & yawAxis;
90  cf[1] += mc & pitchAxis;
91  cf[2] += mc & rollAxis;
92  }
93  }
94 
95  reduce(cf, sumOp<vector>());
96 
97  return cf;
98 }
99 
100 
101 template<class RhoFieldType>
103 (
104  const RhoFieldType& rho,
105  const vectorField& U,
106  vectorField& force
107 )
108 {
109  if (rotor_.mesh().time().timeIndex() % calcFrequency_ == 0)
110  {
111  word calcType = "forces";
112  if (useCoeffs_)
113  {
114  calcType = "coefficients";
115  }
116 
117  Info<< type() << ":" << nl
118  << " solving for target trim " << calcType << nl;
119 
120  const scalar rhoRef = rotor_.rhoRef();
121 
122  // iterate to find new pitch angles to achieve target force
123  scalar err = GREAT;
124  label iter = 0;
125  tensor J(Zero);
126 
127  vector old = Zero;
128  while ((err > tol_) && (iter < nIter_))
129  {
130  // cache initial theta vector
131  vector theta0(theta_);
132 
133  // set initial values
134  old = calcCoeffs(rho, U, thetag(), force);
135 
136  // construct Jacobian by perturbing the pitch angles
137  // by +/-(dTheta_/2)
138  for (label pitchI = 0; pitchI < 3; pitchI++)
139  {
140  theta_[pitchI] -= dTheta_/2.0;
141  vector cf0 = calcCoeffs(rho, U, thetag(), force);
142 
143  theta_[pitchI] += dTheta_;
144  vector cf1 = calcCoeffs(rho, U, thetag(), force);
145 
146  vector ddTheta = (cf1 - cf0)/dTheta_;
147 
148  J[pitchI + 0] = ddTheta[0];
149  J[pitchI + 3] = ddTheta[1];
150  J[pitchI + 6] = ddTheta[2];
151 
152  theta_ = theta0;
153  }
154 
155  // calculate the change in pitch angle vector
156  vector dt = inv(J) & (target_/rhoRef - old);
157 
158  // update pitch angles
159  vector thetaNew = theta_ + relax_*dt;
160 
161  // update error
162  err = mag(thetaNew - theta_);
163 
164  // update for next iteration
165  theta_ = thetaNew;
166  iter++;
167  }
168 
169  if (iter == nIter_)
170  {
171  Info<< " solution not converged in " << iter
172  << " iterations, final residual = " << err
173  << "(" << tol_ << ")" << endl;
174  }
175  else
176  {
177  Info<< " final residual = " << err << "(" << tol_
178  << "), iterations = " << iter << endl;
179  }
180 
181  Info<< " current and target " << calcType << nl
182  << " thrust = " << old[0]*rhoRef << ", " << target_[0] << nl
183  << " pitch = " << old[1]*rhoRef << ", " << target_[1] << nl
184  << " roll = " << old[2]*rhoRef << ", " << target_[2] << nl
185  << " new pitch angles [deg]:" << nl
186  << " theta0 = " << radToDeg(theta_[0]) << nl
187  << " theta1c = " << radToDeg(theta_[1]) << nl
188  << " theta1s = " << radToDeg(theta_[2]) << nl
189  << endl;
190  }
191 }
192 
193 
194 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
195 
197 (
198  const fv::rotorDiskSource& rotor,
199  const dictionary& dict
200 )
201 :
202  trimModel(rotor, dict, typeName),
203  calcFrequency_(-1),
204  useCoeffs_(true),
205  target_(Zero),
206  theta_(Zero),
207  nIter_(50),
208  tol_(1e-8),
209  relax_(1.0),
210  dTheta_(degToRad(0.1)),
211  alpha_(1.0)
212 {
213  read(dict);
214 }
215 
216 
217 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218 
220 {
222 
223  const dictionary& targetDict(coeffs_.subDict("target"));
224  useCoeffs_ = targetDict.getOrDefault("useCoeffs", true);
225  word ext = "";
226  if (useCoeffs_)
227  {
228  ext = "Coeff";
229  }
230 
231  targetDict.readEntry("thrust" + ext, target_[0]);
232  targetDict.readEntry("pitch" + ext, target_[1]);
233  targetDict.readEntry("roll" + ext, target_[2]);
234 
235  const dictionary& pitchAngleDict(coeffs_.subDict("pitchAngles"));
236  theta_[0] = degToRad(pitchAngleDict.get<scalar>("theta0Ini"));
237  theta_[1] = degToRad(pitchAngleDict.get<scalar>("theta1cIni"));
238  theta_[2] = degToRad(pitchAngleDict.get<scalar>("theta1sIni"));
239 
240  coeffs_.readEntry("calcFrequency", calcFrequency_);
241 
242  coeffs_.readIfPresent("nIter", nIter_);
243  coeffs_.readIfPresent("tol", tol_);
244  coeffs_.readIfPresent("relax", relax_);
245 
246  if (coeffs_.readIfPresent("dTheta", dTheta_))
247  {
248  dTheta_ = degToRad(dTheta_);
249  }
250 
251  coeffs_.readEntry("alpha", alpha_);
252 }
253 
254 
256 {
257  const List<vector>& x = rotor_.x();
258 
259  auto ttheta = tmp<scalarField>::New(x.size());
260  auto& t = ttheta.ref();
261 
262  forAll(t, i)
263  {
264  scalar psi = x[i].y();
265  t[i] = theta_[0] + theta_[1]*cos(psi) + theta_[2]*sin(psi);
266  }
267 
268  return ttheta;
269 }
270 
271 
273 (
274  const vectorField& U,
275  vectorField& force
276 )
277 {
278  correctTrim(geometricOneField(), U, force);
279 }
280 
281 
283 (
284  const volScalarField rho,
285  const vectorField& U,
286  vectorField& force
287 )
288 {
289  correctTrim(rho, U, force);
290 }
291 
292 
293 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::Tensor< scalar >
Foam::targetCoeffTrim::read
void read(const dictionary &dict)
Read.
Definition: targetCoeffTrim.C:219
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:107
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:54
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
rho
rho
Definition: readInitialConditions.H:88
Foam::targetCoeffTrim::correctTrim
void correctTrim(const RhoFieldType &rho, const vectorField &U, vectorField &force)
Correct the model.
Definition: targetCoeffTrim.C:103
Foam::targetCoeffTrim::calcCoeffs
vector calcCoeffs(const RhoFieldType &rho, const vectorField &U, const scalarField &alphag, vectorField &force) const
Calculate the rotor force and moment coefficients vector.
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::trimModel::read
virtual void read(const dictionary &dict)
Read.
Definition: trimModel.C:58
Foam::Field< vector >
Foam::targetCoeffTrim::thetag
virtual tmp< scalarField > thetag() const
Return the geometric angle of attack [rad].
Definition: targetCoeffTrim.C:255
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fv::rotorDiskSource
Rotor disk source.
Definition: rotorDiskSource.H:127
Foam::radToDeg
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
Definition: unitConversion.H:54
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::C::C
C()
Construct null.
Definition: C.C:43
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::trimModel
Trim model base class.
Definition: trimModel.H:52
U
U
Definition: pEqn.H:72
targetCoeffTrim.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
geometricOneField.H
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::targetCoeffTrim::targetCoeffTrim
targetCoeffTrim(const fv::rotorDiskSource &rotor, const dictionary &dict)
Constructor.
Definition: targetCoeffTrim.C:197
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::GeometricField< scalar, fvPatchField, volMesh >
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::targetCoeffTrim::correct
virtual void correct(const vectorField &U, vectorField &force)
Correct the model.
Definition: targetCoeffTrim.C:273
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265