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-------------------------------------------------------------------------------
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
27\*---------------------------------------------------------------------------*/
28
29#include "targetCoeffTrim.H"
30#include "geometricOneField.H"
32
33using namespace Foam::constant;
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
41}
42
43
44// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45
46template<class RhoFieldType>
48(
49 const RhoFieldType& rho,
50 const vectorField& U,
51 const scalarField& thetag,
52 vectorField& force
53) const
54{
55 rotor_.calculate(rho, U, thetag, force, false, false);
56
57 const labelList& cells = rotor_.cells();
58 const vectorField& C = rotor_.mesh().C();
59 const List<point>& x = rotor_.x();
60
61 const vector& origin = rotor_.coordSys().origin();
62 const vector& rollAxis = rotor_.coordSys().e1();
63 const vector& pitchAxis = rotor_.coordSys().e2();
64 const vector& yawAxis = rotor_.coordSys().e3();
65
66 scalar coeff1 = alpha_*sqr(rotor_.omega())*mathematical::pi;
67
68 vector cf(Zero);
69 forAll(cells, i)
70 {
71 label celli = cells[i];
72
73 vector fc = force[celli];
74 vector mc = fc^(C[celli] - origin);
75
76 if (useCoeffs_)
77 {
78 scalar radius = x[i].x();
79 scalar coeff2 = rho[celli]*coeff1*pow4(radius);
80
81 // add to coefficient vector
82 cf[0] += (fc & yawAxis)/(coeff2 + ROOTVSMALL);
83 cf[1] += (mc & pitchAxis)/(coeff2*radius + ROOTVSMALL);
84 cf[2] += (mc & rollAxis)/(coeff2*radius + ROOTVSMALL);
85 }
86 else
87 {
88 cf[0] += fc & yawAxis;
89 cf[1] += mc & pitchAxis;
90 cf[2] += mc & rollAxis;
91 }
92 }
93
94 reduce(cf, sumOp<vector>());
95
96 return cf;
97}
98
99
100template<class RhoFieldType>
102(
103 const RhoFieldType& rho,
104 const vectorField& U,
105 vectorField& force
106)
107{
108 if (rotor_.mesh().time().timeIndex() % calcFrequency_ == 0)
109 {
110 word calcType = "forces";
111 if (useCoeffs_)
112 {
113 calcType = "coefficients";
114 }
115
116 Info<< type() << ":" << nl
117 << " solving for target trim " << calcType << nl;
118
119 const scalar rhoRef = rotor_.rhoRef();
120
121 // iterate to find new pitch angles to achieve target force
122 scalar err = GREAT;
123 label iter = 0;
124 tensor J(Zero);
125
126 vector old = Zero;
127 while ((err > tol_) && (iter < nIter_))
128 {
129 // cache initial theta vector
130 vector theta0(theta_);
131
132 // set initial values
133 old = calcCoeffs(rho, U, thetag(), force);
134
135 // construct Jacobian by perturbing the pitch angles
136 // by +/-(dTheta_/2)
137 for (label pitchI = 0; pitchI < 3; pitchI++)
138 {
139 theta_[pitchI] -= dTheta_/2.0;
140 vector cf0 = calcCoeffs(rho, U, thetag(), force);
141
142 theta_[pitchI] += dTheta_;
143 vector cf1 = calcCoeffs(rho, U, thetag(), force);
144
145 vector ddTheta = (cf1 - cf0)/dTheta_;
146
147 J[pitchI + 0] = ddTheta[0];
148 J[pitchI + 3] = ddTheta[1];
149 J[pitchI + 6] = ddTheta[2];
150
151 theta_ = theta0;
152 }
153
154 // calculate the change in pitch angle vector
155 vector dt = inv(J) & (target_/rhoRef - old);
156
157 // update pitch angles
158 vector thetaNew = theta_ + relax_*dt;
159
160 // update error
161 err = mag(thetaNew - theta_);
162
163 // update for next iteration
164 theta_ = thetaNew;
165 iter++;
166 }
167
168 if (iter == nIter_)
169 {
170 Info<< " solution not converged in " << iter
171 << " iterations, final residual = " << err
172 << "(" << tol_ << ")" << endl;
173 }
174 else
175 {
176 Info<< " final residual = " << err << "(" << tol_
177 << "), iterations = " << iter << endl;
178 }
179
180 Info<< " current and target " << calcType << nl
181 << " thrust = " << old[0]*rhoRef << ", " << target_[0] << nl
182 << " pitch = " << old[1]*rhoRef << ", " << target_[1] << nl
183 << " roll = " << old[2]*rhoRef << ", " << target_[2] << nl
184 << " new pitch angles [deg]:" << nl
185 << " theta0 = " << radToDeg(theta_[0]) << nl
186 << " theta1c = " << radToDeg(theta_[1]) << nl
187 << " theta1s = " << radToDeg(theta_[2]) << nl
188 << endl;
189 }
190}
191
192
193// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
194
196(
197 const fv::rotorDiskSource& rotor,
198 const dictionary& dict
199)
200:
201 trimModel(rotor, dict, typeName),
202 calcFrequency_(-1),
203 useCoeffs_(true),
204 target_(Zero),
205 theta_(Zero),
206 nIter_(50),
207 tol_(1e-8),
208 relax_(1.0),
209 dTheta_(degToRad(0.1)),
210 alpha_(1.0)
211{
212 read(dict);
213}
214
215
216// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
217
219{
221
222 const dictionary& targetDict(coeffs_.subDict("target"));
223 useCoeffs_ = targetDict.getOrDefault("useCoeffs", true);
224 word ext = "";
225 if (useCoeffs_)
226 {
227 ext = "Coeff";
228 }
229
230 targetDict.readEntry("thrust" + ext, target_[0]);
231 targetDict.readEntry("pitch" + ext, target_[1]);
232 targetDict.readEntry("roll" + ext, target_[2]);
233
234 const dictionary& pitchAngleDict(coeffs_.subDict("pitchAngles"));
235 theta_[0] = degToRad(pitchAngleDict.get<scalar>("theta0Ini"));
236 theta_[1] = degToRad(pitchAngleDict.get<scalar>("theta1cIni"));
237 theta_[2] = degToRad(pitchAngleDict.get<scalar>("theta1sIni"));
238
239 coeffs_.readEntry("calcFrequency", calcFrequency_);
240
241 coeffs_.readIfPresent("nIter", nIter_);
242 coeffs_.readIfPresent("tol", tol_);
243 coeffs_.readIfPresent("relax", relax_);
244
245 if (coeffs_.readIfPresent("dTheta", dTheta_))
246 {
247 dTheta_ = degToRad(dTheta_);
248 }
249
250 coeffs_.readIfPresent("alpha", alpha_);
251}
252
253
255{
256 const List<vector>& x = rotor_.x();
257
258 auto ttheta = tmp<scalarField>::New(x.size());
259 auto& t = ttheta.ref();
260
261 forAll(t, i)
262 {
263 scalar psi = x[i].y();
264 t[i] = theta_[0] + theta_[1]*cos(psi) + theta_[2]*sin(psi);
265 }
266
267 return ttheta;
268}
269
270
272(
273 const vectorField& U,
274 vectorField& force
275)
276{
277 correctTrim(geometricOneField(), U, force);
278}
279
280
282(
283 const volScalarField rho,
284 const vectorField& U,
285 vectorField& force
286)
287{
288 correctTrim(rho, U, force);
289}
290
291
292// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Graphite solid properties.
Definition: C.H:53
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
virtual const vector e2() const
The local Cartesian y-axis in global coordinates.
virtual const point & origin() const
Return origin.
virtual const vector e1() const
The local Cartesian x-axis in global coordinates.
virtual const vector e3() const
The local Cartesian z-axis in global coordinates.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
const volVectorField & C() const
Return cell centres as volVectorField.
const labelList & cells() const noexcept
Return const access to the cell selection.
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition: fvOptionI.H:37
Applies cell-based momentum sources on velocity (i.e. U) within a specified cylindrical region to app...
void calculate(const RhoFieldType &rho, const vectorField &U, const scalarField &thetag, vectorField &force, const bool divideVolume=true, const bool output=true) const
Calculate forces.
const coordSystem::cylindrical & coordSys() const
Return the rotor coordinate system (r-theta-z)
const List< point > & x() const
Return the cell centre positions in local rotor frame.
scalar omega() const
Return the rotational speed [rad/s].
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Trim model where the operating characteristics of rotor (e.g. blade pitch angle) can vary to reach a ...
virtual tmp< scalarField > thetag() const
Return the geometric angle of attack [rad].
void read(const dictionary &dict)
Read.
bool useCoeffs_
Flag to indicate whether to solve coeffs (true) or forces (false)
scalar alpha_
Coefficient to allow for conversion between US and EU definitions.
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
const fv::rotorDiskSource & rotor_
Reference to the rotor source model.
Definition: trimModel.H:118
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
const volScalarField & psi
const cellShapeList & cells
constexpr scalar pi(M_PI)
Different types of constants.
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333