JSBSim Flight Dynamics Model 1.2.3 (07 Jun 2025)
An Open Source Flight Dynamics and Control Software Library in C++
Loading...
Searching...
No Matches
FGWinds.cpp
1/*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
2
3 Module: FGWinds.cpp
4 Author: Jon Berndt, Tony Peden, Andreas Gaeb
5 Date started: Extracted from FGAtmosphere, which originated in 1998
6 5/2011
7 Purpose: Models winds, gusts, turbulence, and other atmospheric
8 disturbances
9 Called by: FGFDMExec
10
11 ------------- Copyright (C) 2011 Jon S. Berndt (jon@jsbsim.org) -------------
12
13 This program is free software; you can redistribute it and/or modify it under
14 the terms of the GNU Lesser General Public License as published by the Free
15 Software Foundation; either version 2 of the License, or (at your option) any
16 later version.
17
18 This program is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20 FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more
21 details.
22
23 You should have received a copy of the GNU Lesser General Public License along
24 with this program; if not, write to the Free Software Foundation, Inc., 59
25 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
26
27 Further information about the GNU Lesser General Public License can also be
28 found on the world wide web at http://www.gnu.org.
29
30FUNCTIONAL DESCRIPTION
31--------------------------------------------------------------------------------
32
33HISTORY
34--------------------------------------------------------------------------------
35
36%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
37COMMENTS, REFERENCES, and NOTES
38%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
39[1] Anderson, John D. "Introduction to Flight, Third Edition", McGraw-Hill,
40 1989, ISBN 0-07-001641-0
41
42%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
43INCLUDES
44%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
45
46#include "FGWinds.h"
47#include "FGFDMExec.h"
48#include "math/FGTable.h"
49
50using namespace std;
51
52namespace JSBSim {
53
54/*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
55CLASS IMPLEMENTATION
56%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
57
58//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
59// square a value, but preserve the original sign
60
61/*
62static inline double square_signed (double value)
63{
64 if (value < 0)
65 return value * value * -1;
66 else
67 return value * value;
68}
69*/
70
72constexpr double sqr(double x) { return x*x; }
73
75 : FGModel(fdmex), generator(fdmex->GetRandomGenerator())
76{
77 Name = "FGWinds";
78
79 MagnitudedAccelDt = MagnitudeAccel = Magnitude = TurbDirection = 0.0;
80 SetTurbType( ttMilspec );
81 TurbGain = 1.0;
82 TurbRate = 10.0;
83 Rhythmicity = 0.1;
84 spike = target_time = strength = 0.0;
85 wind_from_clockwise = 0.0;
86 psiw = 0.0;
87
88 vGustNED.InitMatrix();
89 vTurbulenceNED.InitMatrix();
90 vCosineGust.InitMatrix();
91
92 // Milspec turbulence model
93 windspeed_at_20ft = 0.;
94 probability_of_exceedence_index = 0;
95 POE_Table = new FGTable(7,12);
96 // this is Figure 7 from p. 49 of MIL-F-8785C
97 // rows: probability of exceedance curve index, cols: altitude in ft
98 *POE_Table
99 << 500.0 << 1750.0 << 3750.0 << 7500.0 << 15000.0 << 25000.0 << 35000.0 << 45000.0 << 55000.0 << 65000.0 << 75000.0 << 80000.0
100 << 1 << 3.2 << 2.2 << 1.5 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0
101 << 2 << 4.2 << 3.6 << 3.3 << 1.6 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0
102 << 3 << 6.6 << 6.9 << 7.4 << 6.7 << 4.6 << 2.7 << 0.4 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0
103 << 4 << 8.6 << 9.6 << 10.6 << 10.1 << 8.0 << 6.6 << 5.0 << 4.2 << 2.7 << 0.0 << 0.0 << 0.0
104 << 5 << 11.8 << 13.0 << 16.0 << 15.1 << 11.6 << 9.7 << 8.1 << 8.2 << 7.9 << 4.9 << 3.2 << 2.1
105 << 6 << 15.6 << 17.6 << 23.0 << 23.6 << 22.1 << 20.0 << 16.0 << 15.1 << 12.1 << 7.9 << 6.2 << 5.1
106 << 7 << 18.7 << 21.5 << 28.4 << 30.2 << 30.7 << 31.0 << 25.2 << 23.1 << 17.5 << 10.7 << 8.4 << 7.2;
107
108 bind();
109 Debug(0);
110}
111
112//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
113
115{
116 delete(POE_Table);
117 Debug(1);
118}
119
120//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
121
122bool FGWinds::InitModel(void)
123{
124 if (!FGModel::InitModel()) return false;
125
126 psiw = 0.0;
127
128 vGustNED.InitMatrix();
129 vTurbulenceNED.InitMatrix();
130 vCosineGust.InitMatrix();
131
132 oneMinusCosineGust.gustProfile.Running = false;
133 oneMinusCosineGust.gustProfile.elapsedTime = 0.0;
134
135 xi_u_km1 = nu_u_km1 = 0;
136 xi_v_km1 = xi_v_km2 = nu_v_km1 = nu_v_km2 = 0;
137 xi_w_km1 = xi_w_km2 = nu_w_km1 = nu_w_km2 = 0;
138 xi_p_km1 = nu_p_km1 = 0;
139 xi_q_km1 = xi_r_km1 = 0;
140
141 return true;
142}
143
144//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
145
146bool FGWinds::Run(bool Holding)
147{
148 if (FGModel::Run(Holding)) return true;
149 if (Holding) return false;
150
151 if (turbType != ttNone)
152 Turbulence(in.AltitudeASL);
153 else
154 vTurbulenceNED.InitMatrix();
155 if (oneMinusCosineGust.gustProfile.Running) CosineGust();
156
157 vTotalWindNED = vWindNED + vGustNED + vCosineGust + vTurbulenceNED;
158
159 // psiw (Wind heading) is the direction the wind is blowing towards
160 if (vWindNED(eX) != 0.0) psiw = atan2( vWindNED(eY), vWindNED(eX) );
161 if (psiw < 0) psiw += 2*M_PI;
162
163 Debug(2);
164 return false;
165}
166
167//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
168//
169// psi is the angle that the wind is blowing *towards*
170
171void FGWinds::SetWindspeed(double speed)
172{
173 if (vWindNED.Magnitude() == 0.0) {
174 psiw = 0.0;
175 vWindNED(eNorth) = speed;
176 } else {
177 vWindNED(eNorth) = speed * cos(psiw);
178 vWindNED(eEast) = speed * sin(psiw);
179 vWindNED(eDown) = 0.0;
180 }
181}
182
183//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
184
185double FGWinds::GetWindspeed(void) const
186{
187 return vWindNED.Magnitude();
188}
189
190//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
191//
192// psi is the angle that the wind is blowing *towards*
193
194void FGWinds::SetWindPsi(double dir)
195{
196 double mag = GetWindspeed();
197 psiw = dir;
198 SetWindspeed(mag);
199}
200
201//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
202
203void FGWinds::Turbulence(double h)
204{
205 switch (turbType) {
206
207 case ttCulp: {
208
209 vTurbPQR(eP) = wind_from_clockwise;
210 if (TurbGain == 0.0) return;
211
212 // keep the inputs within allowable limts for this model
213 if (TurbGain < 0.0) TurbGain = 0.0;
214 if (TurbGain > 1.0) TurbGain = 1.0;
215 if (TurbRate < 0.0) TurbRate = 0.0;
216 if (TurbRate > 30.0) TurbRate = 30.0;
217 if (Rhythmicity < 0.0) Rhythmicity = 0.0;
218 if (Rhythmicity > 1.0) Rhythmicity = 1.0;
219
220 // generate a sine wave corresponding to turbulence rate in hertz
221 double time = FDMExec->GetSimTime();
222 double sinewave = sin( time * TurbRate * 6.283185307 );
223
224 double random = 0.0;
225 if (target_time == 0.0) {
226 strength = random = generator->GetUniformRandomNumber();
227 target_time = time + 0.71 + (random * 0.5);
228 }
229 if (time > target_time) {
230 spike = 1.0;
231 target_time = 0.0;
232 }
233
234 // max vertical wind speed in fps, corresponds to TurbGain = 1.0
235 double max_vs = 40;
236
237 vTurbulenceNED.InitMatrix();
238 double delta = strength * max_vs * TurbGain * (1-Rhythmicity) * spike;
239
240 // Vertical component of turbulence.
241 vTurbulenceNED(eDown) = sinewave * max_vs * TurbGain * Rhythmicity;
242 vTurbulenceNED(eDown)+= delta;
243 if (in.DistanceAGL/in.wingspan < 3.0)
244 vTurbulenceNED(eDown) *= in.DistanceAGL/in.wingspan * 0.3333;
245
246 // Yaw component of turbulence.
247 vTurbulenceNED(eNorth) = sin( delta * 3.0 );
248 vTurbulenceNED(eEast) = cos( delta * 3.0 );
249
250 // Roll component of turbulence. Clockwise vortex causes left roll.
251 vTurbPQR(eP) += delta * 0.04;
252
253 spike = spike * 0.9;
254 break;
255 }
256 case ttMilspec:
257 case ttTustin: {
258
259 // an index of zero means turbulence is disabled
260 // airspeed occurs as divisor in the code below
261 if (probability_of_exceedence_index == 0 || in.V == 0) {
262 vTurbulenceNED(eNorth) = vTurbulenceNED(eEast) = vTurbulenceNED(eDown) = 0.0;
263 vTurbPQR(eP) = vTurbPQR(eQ) = vTurbPQR(eR) = 0.0;
264 return;
265 }
266
267 // Turbulence model according to MIL-F-8785C (Flying Qualities of Piloted Aircraft)
268 double b_w = in.wingspan, L_u, L_w, sig_u, sig_w;
269
270 if (b_w == 0.) b_w = 30.;
271
272 // clip height functions at 10 ft
273 if (h <= 10.) h = 10;
274
275 // Scale lengths L and amplitudes sigma as function of height
276 if (h <= 1000) {
277 L_u = h/pow(0.177 + 0.000823*h, 1.2); // MIL-F-8785c, Fig. 10, p. 55
278 L_w = h;
279 sig_w = 0.1*windspeed_at_20ft;
280 sig_u = sig_w/pow(0.177 + 0.000823*h, 0.4); // MIL-F-8785c, Fig. 11, p. 56
281 } else if (h <= 2000) {
282 // linear interpolation between low altitude and high altitude models
283 L_u = L_w = 1000 + (h-1000.)/1000.*750.;
284 sig_u = sig_w = 0.1*windspeed_at_20ft
285 + (h-1000.)/1000.*(POE_Table->GetValue(probability_of_exceedence_index, h) - 0.1*windspeed_at_20ft);
286 } else {
287 L_u = L_w = 1750.; // MIL-F-8785c, Sec. 3.7.2.1, p. 48
288 sig_u = sig_w = POE_Table->GetValue(probability_of_exceedence_index, h);
289 }
290
291 double
292 T_V = in.totalDeltaT, // for compatibility of nomenclature
293 sig_p = 1.9/sqrt(L_w*b_w)*sig_w, // Yeager1998, eq. (8)
294 //sig_q = sqrt(M_PI/2/L_w/b_w), // eq. (14)
295 //sig_r = sqrt(2*M_PI/3/L_w/b_w), // eq. (17)
296 L_p = sqrt(L_w*b_w)/2.6, // eq. (10)
297 tau_u = L_u/in.V, // eq. (6)
298 tau_w = L_w/in.V, // eq. (3)
299 tau_p = L_p/in.V, // eq. (9)
300 tau_q = 4*b_w/M_PI/in.V, // eq. (13)
301 tau_r =3*b_w/M_PI/in.V, // eq. (17)
302 nu_u = generator->GetNormalRandomNumber(),
303 nu_v = generator->GetNormalRandomNumber(),
304 nu_w = generator->GetNormalRandomNumber(),
305 nu_p = generator->GetNormalRandomNumber(),
306 xi_u=0, xi_v=0, xi_w=0, xi_p=0, xi_q=0, xi_r=0;
307
308 // values of turbulence NED velocities
309
310 if (turbType == ttTustin) {
311 // the following is the Tustin formulation of Yeager's report
312 double
313 omega_w = in.V/L_w, // hidden in nomenclature p. 3
314 omega_v = in.V/L_u, // this is defined nowhere
315 C_BL = 1/tau_u/tan(T_V/2/tau_u), // eq. (19)
316 C_BLp = 1/tau_p/tan(T_V/2/tau_p), // eq. (22)
317 C_BLq = 1/tau_q/tan(T_V/2/tau_q), // eq. (24)
318 C_BLr = 1/tau_r/tan(T_V/2/tau_r); // eq. (26)
319
320 // all values calculated so far are strictly positive, except for
321 // the random numbers nu_*. This means that in the code below, all
322 // divisors are strictly positive, too, and no floating point
323 // exception should occur.
324 xi_u = -(1 - C_BL*tau_u)/(1 + C_BL*tau_u)*xi_u_km1
325 + sig_u*sqrt(2*tau_u/T_V)/(1 + C_BL*tau_u)*(nu_u + nu_u_km1); // eq. (18)
326 xi_v = -2*(sqr(omega_v) - sqr(C_BL))/sqr(omega_v + C_BL)*xi_v_km1
327 - sqr(omega_v - C_BL)/sqr(omega_v + C_BL) * xi_v_km2
328 + sig_u*sqrt(3*omega_v/T_V)/sqr(omega_v + C_BL)*(
329 (C_BL + omega_v/sqrt(3.))*nu_v
330 + 2/sqrt(3.)*omega_v*nu_v_km1
331 + (omega_v/sqrt(3.) - C_BL)*nu_v_km2); // eq. (20) for v
332 xi_w = -2*(sqr(omega_w) - sqr(C_BL))/sqr(omega_w + C_BL)*xi_w_km1
333 - sqr(omega_w - C_BL)/sqr(omega_w + C_BL) * xi_w_km2
334 + sig_w*sqrt(3*omega_w/T_V)/sqr(omega_w + C_BL)*(
335 (C_BL + omega_w/sqrt(3.))*nu_w
336 + 2/sqrt(3.)*omega_w*nu_w_km1
337 + (omega_w/sqrt(3.) - C_BL)*nu_w_km2); // eq. (20) for w
338 xi_p = -(1 - C_BLp*tau_p)/(1 + C_BLp*tau_p)*xi_p_km1
339 + sig_p*sqrt(2*tau_p/T_V)/(1 + C_BLp*tau_p) * (nu_p + nu_p_km1); // eq. (21)
340 xi_q = -(1 - 4*b_w*C_BLq/M_PI/in.V)/(1 + 4*b_w*C_BLq/M_PI/in.V) * xi_q_km1
341 + C_BLq/in.V/(1 + 4*b_w*C_BLq/M_PI/in.V) * (xi_w - xi_w_km1); // eq. (23)
342 xi_r = - (1 - 3*b_w*C_BLr/M_PI/in.V)/(1 + 3*b_w*C_BLr/M_PI/in.V) * xi_r_km1
343 + C_BLr/in.V/(1 + 3*b_w*C_BLr/M_PI/in.V) * (xi_v - xi_v_km1); // eq. (25)
344
345 } else if (turbType == ttMilspec) {
346 // the following is the MIL-STD-1797A formulation
347 // as cited in Yeager's report
348 xi_u = (1 - T_V/tau_u) *xi_u_km1 + sig_u*sqrt(2*T_V/tau_u)*nu_u; // eq. (30)
349 xi_v = (1 - 2*T_V/tau_u)*xi_v_km1 + sig_u*sqrt(4*T_V/tau_u)*nu_v; // eq. (31)
350 xi_w = (1 - 2*T_V/tau_w)*xi_w_km1 + sig_w*sqrt(4*T_V/tau_w)*nu_w; // eq. (32)
351 xi_p = (1 - T_V/tau_p) *xi_p_km1 + sig_p*sqrt(2*T_V/tau_p)*nu_p; // eq. (33)
352 xi_q = (1 - T_V/tau_q) *xi_q_km1 + M_PI/4/b_w*(xi_w - xi_w_km1); // eq. (34)
353 xi_r = (1 - T_V/tau_r) *xi_r_km1 + M_PI/3/b_w*(xi_v - xi_v_km1); // eq. (35)
354 }
355
356 // rotate by wind azimuth and assign the velocities
357 double cospsi = cos(psiw), sinpsi = sin(psiw);
358 vTurbulenceNED(eNorth) = cospsi*xi_u + sinpsi*xi_v;
359 vTurbulenceNED(eEast) = -sinpsi*xi_u + cospsi*xi_v;
360 vTurbulenceNED(eDown) = xi_w;
361
362 vTurbPQR(eP) = cospsi*xi_p + sinpsi*xi_q;
363 vTurbPQR(eQ) = -sinpsi*xi_p + cospsi*xi_q;
364 vTurbPQR(eR) = xi_r;
365
366 // vTurbPQR is in the body fixed frame, not NED
367 vTurbPQR = in.Tl2b*vTurbPQR;
368
369 // hand on the values for the next timestep
370 xi_u_km1 = xi_u; nu_u_km1 = nu_u;
371 xi_v_km2 = xi_v_km1; xi_v_km1 = xi_v; nu_v_km2 = nu_v_km1; nu_v_km1 = nu_v;
372 xi_w_km2 = xi_w_km1; xi_w_km1 = xi_w; nu_w_km2 = nu_w_km1; nu_w_km1 = nu_w;
373 xi_p_km1 = xi_p; nu_p_km1 = nu_p;
374 xi_q_km1 = xi_q;
375 xi_r_km1 = xi_r;
376
377 }
378 default:
379 break;
380 }
381
382 TurbDirection = atan2( vTurbulenceNED(eEast), vTurbulenceNED(eNorth))*radtodeg;
383
384}
385
386//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
387
388double FGWinds::CosineGustProfile(double startDuration, double steadyDuration, double endDuration, double elapsedTime)
389{
390 double factor = 0.0;
391 if (elapsedTime >= 0 && elapsedTime <= startDuration) {
392 factor = (1.0 - cos(M_PI*elapsedTime/startDuration))/2.0;
393 } else if (elapsedTime > startDuration && (elapsedTime <= (startDuration + steadyDuration))) {
394 factor = 1.0;
395 } else if (elapsedTime > (startDuration + steadyDuration) && elapsedTime <= (startDuration + steadyDuration + endDuration)) {
396 factor = (1-cos(M_PI*(1-(elapsedTime-(startDuration + steadyDuration))/endDuration)))/2.0;
397 } else {
398 factor = 0.0;
399 }
400
401 return factor;
402}
403
404//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
405
406void FGWinds::CosineGust()
407{
408 struct OneMinusCosineProfile& profile = oneMinusCosineGust.gustProfile;
409
410 double factor = CosineGustProfile( profile.startupDuration,
411 profile.steadyDuration,
412 profile.endDuration,
413 profile.elapsedTime);
414 // Normalize the gust wind vector
415 oneMinusCosineGust.vWind.Normalize();
416
417 if (oneMinusCosineGust.vWindTransformed.Magnitude() == 0.0) {
418 switch (oneMinusCosineGust.gustFrame) {
419 case gfBody:
420 oneMinusCosineGust.vWindTransformed = in.Tl2b.Inverse() * oneMinusCosineGust.vWind;
421 break;
422 case gfWind:
423 oneMinusCosineGust.vWindTransformed = in.Tl2b.Inverse() * in.Tw2b * oneMinusCosineGust.vWind;
424 break;
425 case gfLocal:
426 // this is the native frame - and the default.
427 oneMinusCosineGust.vWindTransformed = oneMinusCosineGust.vWind;
428 break;
429 default:
430 break;
431 }
432 }
433
434 vCosineGust = factor * oneMinusCosineGust.vWindTransformed * oneMinusCosineGust.magnitude;
435
436 profile.elapsedTime += in.totalDeltaT;
437
438 if (profile.elapsedTime > (profile.startupDuration + profile.steadyDuration + profile.endDuration)) {
439 profile.Running = false;
440 profile.elapsedTime = 0.0;
441 oneMinusCosineGust.vWindTransformed.InitMatrix();
442 vCosineGust.InitMatrix(0);
443 }
444}
445
446//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
447
448void FGWinds::NumberOfUpDownburstCells(int num)
449{
450 for (unsigned int i=0; i<UpDownBurstCells.size();i++) delete UpDownBurstCells[i];
451 UpDownBurstCells.clear();
452 if (num >= 0) {
453 for (int i=0; i<num; i++) UpDownBurstCells.push_back(new struct UpDownBurst);
454 }
455}
456
457//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
458// Calculates the distance between a specified point (where presumably the
459// Up/Downburst is centered) and the current vehicle location. The distance
460// here is calculated from the Haversine formula.
461
462double FGWinds::DistanceFromRingCenter(double lat, double lon)
463{
464 double deltaLat = in.latitude - lat;
465 double deltaLong = in.longitude - lon;
466 double dLat2 = deltaLat/2.0;
467 double dLong2 = deltaLong/2.0;
468 double a = sin(dLat2)*sin(dLat2)
469 + cos(lat)*cos(in.latitude)*sin(dLong2)*sin(dLong2);
470 double c = 2.0*atan2(sqrt(a), sqrt(1.0 - a));
471 double d = in.planetRadius*c;
472 return d;
473}
474
475//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
476
477void FGWinds::UpDownBurst()
478{
479
480 for (unsigned int i=0; i<UpDownBurstCells.size(); i++) {
481 /*double d =*/ DistanceFromRingCenter(UpDownBurstCells[i]->ringLatitude, UpDownBurstCells[i]->ringLongitude);
482
483 }
484}
485
486//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
487// User is supplying a wind specific random seed to override the random seed
488// used by FGFDMExec. So initialise a new random number generator with this
489// random seed rather than referencing the FGFDMExec random number generator.
490
491void FGWinds::SetRandomSeed(int sr)
492{
493 RandomSeed = sr;
494 generator = std::make_shared<RandomNumberGenerator>(*RandomSeed);
495}
496
497int FGWinds::GetRandomSeed(void) const {
498 if (RandomSeed)
499 return *RandomSeed;
500 else
501 return FDMExec->SRand();
502}
503
504//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
505
506void FGWinds::bind(void)
507{
508 typedef double (FGWinds::*PMF)(int) const;
509 typedef int (FGWinds::*PMFt)(void) const;
510 typedef void (FGWinds::*PMFd)(int,double);
511 typedef void (FGWinds::*PMFi)(int);
512 typedef double (FGWinds::*Ptr)(void) const;
513
514 // User-specified steady, constant, wind properties (local navigational/geographic frame: N-E-D)
515 PropertyManager->Tie("atmosphere/psiw-rad", this, &FGWinds::GetWindPsi, &FGWinds::SetWindPsi);
516 PropertyManager->Tie("atmosphere/wind-north-fps", this, eNorth, (PMF)&FGWinds::GetWindNED,
517 (PMFd)&FGWinds::SetWindNED);
518 PropertyManager->Tie("atmosphere/wind-east-fps", this, eEast, (PMF)&FGWinds::GetWindNED,
519 (PMFd)&FGWinds::SetWindNED);
520 PropertyManager->Tie("atmosphere/wind-down-fps", this, eDown, (PMF)&FGWinds::GetWindNED,
521 (PMFd)&FGWinds::SetWindNED);
522 PropertyManager->Tie("atmosphere/wind-mag-fps", this, &FGWinds::GetWindspeed,
523 &FGWinds::SetWindspeed);
524
525 // User-specifieded gust (local navigational/geographic frame: N-E-D)
526 PropertyManager->Tie("atmosphere/gust-north-fps", this, eNorth, (PMF)&FGWinds::GetGustNED,
527 (PMFd)&FGWinds::SetGustNED);
528 PropertyManager->Tie("atmosphere/gust-east-fps", this, eEast, (PMF)&FGWinds::GetGustNED,
529 (PMFd)&FGWinds::SetGustNED);
530 PropertyManager->Tie("atmosphere/gust-down-fps", this, eDown, (PMF)&FGWinds::GetGustNED,
531 (PMFd)&FGWinds::SetGustNED);
532
533 // User-specified 1 - cosine gust parameters (in specified frame)
534 PropertyManager->Tie("atmosphere/cosine-gust/startup-duration-sec", this, (Ptr)0L, &FGWinds::StartupGustDuration);
535 PropertyManager->Tie("atmosphere/cosine-gust/steady-duration-sec", this, (Ptr)0L, &FGWinds::SteadyGustDuration);
536 PropertyManager->Tie("atmosphere/cosine-gust/end-duration-sec", this, (Ptr)0L, &FGWinds::EndGustDuration);
537 PropertyManager->Tie("atmosphere/cosine-gust/magnitude-ft_sec", this, (Ptr)0L, &FGWinds::GustMagnitude);
538 PropertyManager->Tie("atmosphere/cosine-gust/frame", this, (PMFt)0L, (PMFi)&FGWinds::GustFrame);
539 PropertyManager->Tie("atmosphere/cosine-gust/X-velocity-ft_sec", this, (Ptr)0L, &FGWinds::GustXComponent);
540 PropertyManager->Tie("atmosphere/cosine-gust/Y-velocity-ft_sec", this, (Ptr)0L, &FGWinds::GustYComponent);
541 PropertyManager->Tie("atmosphere/cosine-gust/Z-velocity-ft_sec", this, (Ptr)0L, &FGWinds::GustZComponent);
542 PropertyManager->Tie("atmosphere/cosine-gust/start", this, static_cast<bool (FGWinds::*)(void) const>(nullptr), &FGWinds::StartGust);
543
544 // User-specified Up- Down-burst parameters
545 PropertyManager->Tie("atmosphere/updownburst/number-of-cells", this, (PMFt)0L, &FGWinds::NumberOfUpDownburstCells);
546// PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
547// PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
548// PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
549// PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
550// PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
551// PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
552// PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
553
554 // User-specified turbulence (local navigational/geographic frame: N-E-D)
555 PropertyManager->Tie("atmosphere/turb-north-fps", this, eNorth, (PMF)&FGWinds::GetTurbNED,
556 (PMFd)&FGWinds::SetTurbNED);
557 PropertyManager->Tie("atmosphere/turb-east-fps", this, eEast, (PMF)&FGWinds::GetTurbNED,
558 (PMFd)&FGWinds::SetTurbNED);
559 PropertyManager->Tie("atmosphere/turb-down-fps", this, eDown, (PMF)&FGWinds::GetTurbNED,
560 (PMFd)&FGWinds::SetTurbNED);
561 // Experimental turbulence parameters
562 PropertyManager->Tie("atmosphere/p-turb-rad_sec", this,1, (PMF)&FGWinds::GetTurbPQR);
563 PropertyManager->Tie("atmosphere/q-turb-rad_sec", this,2, (PMF)&FGWinds::GetTurbPQR);
564 PropertyManager->Tie("atmosphere/r-turb-rad_sec", this,3, (PMF)&FGWinds::GetTurbPQR);
565 PropertyManager->Tie("atmosphere/turb-type", this, (PMFt)&FGWinds::GetTurbType, (PMFi)&FGWinds::SetTurbType);
566 PropertyManager->Tie("atmosphere/turb-rate", this, &FGWinds::GetTurbRate, &FGWinds::SetTurbRate);
567 PropertyManager->Tie("atmosphere/turb-gain", this, &FGWinds::GetTurbGain, &FGWinds::SetTurbGain);
568 PropertyManager->Tie("atmosphere/turb-rhythmicity", this, &FGWinds::GetRhythmicity,
569 &FGWinds::SetRhythmicity);
570
571 // Parameters for milspec turbulence
572 PropertyManager->Tie("atmosphere/turbulence/milspec/windspeed_at_20ft_AGL-fps",
573 this, &FGWinds::GetWindspeed20ft,
574 &FGWinds::SetWindspeed20ft);
575 PropertyManager->Tie("atmosphere/turbulence/milspec/severity",
576 this, &FGWinds::GetProbabilityOfExceedence,
578
579 // Total, calculated winds (local navigational/geographic frame: N-E-D). Read only.
580 PropertyManager->Tie("atmosphere/total-wind-north-fps", this, eNorth, (PMF)&FGWinds::GetTotalWindNED);
581 PropertyManager->Tie("atmosphere/total-wind-east-fps", this, eEast, (PMF)&FGWinds::GetTotalWindNED);
582 PropertyManager->Tie("atmosphere/total-wind-down-fps", this, eDown, (PMF)&FGWinds::GetTotalWindNED);
583
584 // Allow user to specify a separate random seed independent of the FDMExec random seed
585 PropertyManager->Tie("atmosphere/randomseed", this, (PMFt)&FGWinds::GetRandomSeed, &FGWinds::SetRandomSeed);
586}
587
588//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
589// The bitmasked value choices are as follows:
590// unset: In this case (the default) JSBSim would only print
591// out the normally expected messages, essentially echoing
592// the config files as they are read. If the environment
593// variable is not set, debug_lvl is set to 1 internally
594// 0: This requests JSBSim not to output any messages
595// whatsoever.
596// 1: This value explicity requests the normal JSBSim
597// startup messages
598// 2: This value asks for a message to be printed out when
599// a class is instantiated
600// 4: When this value is set, a message is displayed when a
601// FGModel object executes its Run() method
602// 8: When this value is set, various runtime state variables
603// are printed out periodically
604// 16: When set various parameters are sanity checked and
605// a message is printed out when they go out of bounds
606
607void FGWinds::Debug(int from)
608{
609 if (debug_lvl <= 0) return;
610
611 if (debug_lvl & 1) { // Standard console startup message output
612 if (from == 0) { // Constructor
613 }
614 }
615 if (debug_lvl & 2 ) { // Instantiation/Destruction notification
616 if (from == 0) cout << "Instantiated: FGWinds" << endl;
617 if (from == 1) cout << "Destroyed: FGWinds" << endl;
618 }
619 if (debug_lvl & 4 ) { // Run() method entry print for FGModel-derived objects
620 }
621 if (debug_lvl & 8 ) { // Runtime state variables
622 }
623 if (debug_lvl & 16) { // Sanity checking
624 }
625 if (debug_lvl & 128) { //
626 }
627 if (debug_lvl & 64) {
628 if (from == 0) { // Constructor
629 }
630 }
631}
632
633} // namespace JSBSim
double Magnitude(void) const
Length of the vector.
FGColumnVector3 & Normalize(void)
Normalize.
Encapsulates the JSBSim simulation executive.
Definition FGFDMExec.h:184
double GetSimTime(void) const
Returns the cumulative simulation time in seconds.
Definition FGFDMExec.h:549
FGMatrix33 Inverse(void) const
Return the inverse of the matrix.
Base class for all scheduled JSBSim models.
Definition FGModel.h:70
virtual bool Run(bool Holding)
Runs the model; called by the Executive.
Definition FGModel.cpp:89
Lookup table class.
Definition FGTable.h:234
double GetValue(void) const
Get the current table value.
Definition FGTable.cpp:465
virtual const FGColumnVector3 & GetWindNED(void) const
Retrieves the wind components in NED frame.
Definition FGWinds.h:205
virtual void GustMagnitude(double mag)
Specifies the magnitude of the gust in feet/second.
Definition FGWinds.h:337
virtual const FGColumnVector3 & GetTotalWindNED(void) const
Retrieves the total wind components in NED frame.
Definition FGWinds.h:188
virtual void GustXComponent(double x)
Specifies the X component of velocity in the specified gust frame (ft/sec).
Definition FGWinds.h:347
virtual void SetTurbNED(int idx, double turb)
Sets a turbulence component in NED frame.
Definition FGWinds.h:234
virtual void SetWindPsi(double dir)
Sets the direction that the wind is coming from.
Definition FGWinds.cpp:194
virtual void StartupGustDuration(double dur)
Specifies the duration of the startup portion of the gust.
Definition FGWinds.h:331
bool Run(bool Holding) override
Runs the winds model; called by the Executive Can pass in a value indicating if the executive is dire...
Definition FGWinds.cpp:146
virtual void SetTurbType(tType tt)
Turbulence models available: ttNone, ttStandard, ttBerndt, ttCulp, ttMilspec, ttTustin.
Definition FGWinds.h:250
virtual double GetTurbNED(int idx) const
Retrieves a turbulence component in NED frame.
Definition FGWinds.h:243
virtual void SteadyGustDuration(double dur)
Specifies the length of time that the gust is at a steady, full strength.
Definition FGWinds.h:333
~FGWinds()
Destructor.
Definition FGWinds.cpp:114
virtual void EndGustDuration(double dur)
Specifies the length of time it takes for the gust to return to zero velocity.
Definition FGWinds.h:335
virtual void SetWindNED(double wN, double wE, double wD)
Sets the wind components in NED frame.
Definition FGWinds.h:196
virtual void GustZComponent(double z)
Specifies the Z component of velocity in the specified gust frame (ft/sec).
Definition FGWinds.h:351
virtual void SetProbabilityOfExceedence(int idx)
allowable range: 0-7, 3=light, 4=moderate, 6=severe turbulence
Definition FGWinds.h:271
virtual void GustFrame(eGustFrame gFrame)
Specifies the frame that the gust direction vector components are specified in.
Definition FGWinds.h:345
FGWinds(FGFDMExec *)
Constructor.
Definition FGWinds.cpp:74
virtual void GustYComponent(double y)
Specifies the Y component of velocity in the specified gust frame (ft/sec).
Definition FGWinds.h:349
virtual void StartGust(bool running)
Initiates the execution of the gust.
Definition FGWinds.h:329
virtual const FGColumnVector3 & GetGustNED(void) const
Retrieves the gust components in NED frame.
Definition FGWinds.h:246
virtual double GetWindPsi(void) const
Retrieves the direction that the wind is coming from.
Definition FGWinds.h:213
virtual void SetGustNED(int idx, double gust)
Sets a gust component in NED frame.
Definition FGWinds.h:231
struct OneMinusCosineProfile gustProfile
Definition FGWinds.h:300