]> git.mxchange.org Git - flightgear.git/blob - src/FDM/JSBSim/models/FGAtmosphere.cpp
Sync. with JSBSim CVS
[flightgear.git] / src / FDM / JSBSim / models / FGAtmosphere.cpp
1 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
2
3  Module:       FGAtmosphere.cpp
4  Author:       Jon Berndt
5                Implementation of 1959 Standard Atmosphere added by Tony Peden
6  Date started: 11/24/98
7  Purpose:      Models the atmosphere
8  Called by:    FGSimExec
9
10  ------------- Copyright (C) 1999  Jon S. Berndt (jon@jsbsim.org) -------------
11
12  This program is free software; you can redistribute it and/or modify it under
13  the terms of the GNU Lesser General Public License as published by the Free Software
14  Foundation; either version 2 of the License, or (at your option) any later
15  version.
16
17  This program is distributed in the hope that it will be useful, but WITHOUT
18  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19  FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public License for more
20  details.
21
22  You should have received a copy of the GNU Lesser General Public License along with
23  this program; if not, write to the Free Software Foundation, Inc., 59 Temple
24  Place - Suite 330, Boston, MA  02111-1307, USA.
25
26  Further information about the GNU Lesser General Public License can also be found on
27  the world wide web at http://www.gnu.org.
28
29 FUNCTIONAL DESCRIPTION
30 --------------------------------------------------------------------------------
31 Models the atmosphere. The equation used below was determined by a third order
32 curve fit using Excel. The data is from the ICAO atmosphere model.
33
34 HISTORY
35 --------------------------------------------------------------------------------
36 11/24/98   JSB   Created
37 07/23/99   TP    Added implementation of 1959 Standard Atmosphere
38                  Moved calculation of Mach number to FGPropagate
39                  Later updated to '76 model
40 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
41 COMMENTS, REFERENCES,  and NOTES
42 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
43 [1]   Anderson, John D. "Introduction to Flight, Third Edition", McGraw-Hill,
44       1989, ISBN 0-07-001641-0
45
46 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
47 INCLUDES
48 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
49
50 #include "FGAtmosphere.h"
51 #include "FGAircraft.h"
52 #include "FGPropagate.h"
53 #include "FGInertial.h"
54 #include "FGAuxiliary.h"
55 #include "FGFDMExec.h"
56 #include "input_output/FGPropertyManager.h"
57 #include <iostream>
58 #include <cstdlib>
59
60 using namespace std;
61
62 namespace JSBSim {
63
64 static const char *IdSrc = "$Id: FGAtmosphere.cpp,v 1.41 2010/11/30 12:19:57 jberndt Exp $";
65 static const char *IdHdr = ID_ATMOSPHERE;
66
67 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
68 CLASS IMPLEMENTATION
69 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
70
71 FGAtmosphere::FGAtmosphere(FGFDMExec* fdmex) : FGModel(fdmex)
72 {
73   Name = "FGAtmosphere";
74   lastIndex = 0;
75   h = 0.0;
76   psiw = 0.0;
77   htab[0]=0;
78   htab[1]= 36089.0;
79   htab[2]= 65617.0;
80   htab[3]=104987.0;
81   htab[4]=154199.0;
82   htab[5]=167322.0;
83   htab[6]=232940.0;
84   htab[7]=278385.0; //ft.
85
86   MagnitudedAccelDt = MagnitudeAccel = Magnitude = 0.0;
87 //  SetTurbType( ttCulp );
88   SetTurbType( ttNone );
89   TurbGain = 1.0;
90   TurbRate = 10.0;
91   Rhythmicity = 0.1;
92   spike = target_time = strength = 0.0;
93   wind_from_clockwise = 0.0;
94   SutherlandConstant = 198.72; // deg Rankine
95   Beta = 2.269690E-08; // slug/(sec ft R^0.5)
96
97   T_dev_sl = T_dev = delta_T = 0.0;
98   StandardTempOnly = false;
99   first_pass = true;
100   vGustNED.InitMatrix();
101   vTurbulenceNED.InitMatrix();
102
103   // Milspec turbulence model
104   windspeed_at_20ft = 0.;
105   probability_of_exceedence_index = 0;
106   POE_Table = new FGTable(7,12);
107   // this is Figure 7 from p. 49 of MIL-F-8785C
108   // rows: probability of exceedance curve index, cols: altitude in ft
109   *POE_Table
110            << 500.0 << 1750.0 << 3750.0 << 7500.0 << 15000.0 << 25000.0 << 35000.0 << 45000.0 << 55000.0 << 65000.0 << 75000.0 << 80000.0
111     << 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
112     << 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
113     << 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
114     << 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
115     << 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
116     << 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
117     << 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;
118
119   bind();
120   Debug(0);
121 }
122
123 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
124
125 FGAtmosphere::~FGAtmosphere()
126 {
127   Debug(1);
128 }
129
130 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
131
132 bool FGAtmosphere::InitModel(void)
133 {
134   if (!FGModel::InitModel()) return false;
135
136   UseInternal();  // this is the default
137
138   Calculate(h);
139   StdSLtemperature = SLtemperature = 518.67;
140   StdSLpressure    = SLpressure = 2116.22;
141   StdSLdensity     = SLdensity = 0.00237767;
142   StdSLsoundspeed  = SLsoundspeed = sqrt(SHRatio*Reng*StdSLtemperature);
143   rSLtemperature = 1.0/StdSLtemperature;
144   rSLpressure    = 1.0/StdSLpressure;
145   rSLdensity     = 1.0/StdSLdensity;
146   rSLsoundspeed  = 1.0/StdSLsoundspeed;
147
148   return true;
149 }
150
151 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
152
153 bool FGAtmosphere::Run(void)
154 {
155   if (FGModel::Run()) return true;
156   if (FDMExec->Holding()) return false;
157
158   RunPreFunctions();
159
160   T_dev = 0.0;
161   h = FDMExec->GetPropagate()->GetAltitudeASL();
162
163   if (!useExternal) {
164     Calculate(h);
165     CalculateDerived();
166   } else {
167     CalculateDerived();
168   }
169
170   RunPostFunctions();
171
172   Debug(2);
173   return false;
174 }
175
176 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
177 //
178 // See reference 1
179
180 void FGAtmosphere::Calculate(double altitude)
181 {
182   double slope, reftemp, refpress;
183   int i = lastIndex;
184
185   if (altitude < htab[lastIndex]) {
186     if (altitude <= 0) {
187       i = 0;
188       altitude=0;
189     } else {
190        i = lastIndex-1;
191        while (htab[i] > altitude) i--;
192     }
193   } else if (altitude > htab[lastIndex+1]) {
194     if (altitude >= htab[7]) {
195       i = 7;
196       altitude = htab[7];
197     } else {
198       i = lastIndex+1;
199       while (htab[i+1] < altitude) i++;
200     }
201   }
202
203   switch(i) {
204   case 0: // Sea level
205     slope     = -0.00356616; // R/ft.
206     reftemp   = 518.67;   // in degrees Rankine, 288.15 Kelvin
207     refpress  = 2116.22;    // psf
208     //refdens   = 0.00237767;  // slugs/cubic ft.
209     break;
210   case 1:     // 36089 ft. or 11 km
211     slope     = 0;
212     reftemp   = 389.97; // in degrees Rankine, 216.65 Kelvin
213     refpress  = 472.763;
214     //refdens   = 0.000706032;
215     break;
216   case 2:     // 65616 ft. or 20 km
217     slope     = 0.00054864;
218     reftemp   = 389.97; // in degrees Rankine, 216.65 Kelvin
219     refpress  = 114.636;
220     //refdens   = 0.000171306;
221     break;
222   case 3:     // 104986 ft. or 32 km
223     slope     = 0.001536192;
224     reftemp   = 411.57; // in degrees Rankine, 228.65 Kelvin
225     refpress  = 18.128;
226     //refdens   = 1.18422e-05;
227     break;
228   case 4:     // 154199 ft. 47 km
229     slope     = 0;
230     reftemp   = 487.17; // in degrees Rankine, 270.65 Kelvin
231     refpress  = 2.316;
232     //refdens   = 4.00585e-7;
233     break;
234   case 5:     // 167322 ft. or 51 km
235     slope     = -0.001536192;
236     reftemp   = 487.17; // in degrees Rankine, 270.65 Kelvin
237     refpress  = 1.398;
238     //refdens   = 8.17102e-7;
239     break;
240   case 6:     // 232940 ft. or 71 km
241     slope     = -0.00109728;
242     reftemp   = 386.368; // in degrees Rankine, 214.649 Kelvin
243     refpress  = 0.0826;
244     //refdens   = 8.77702e-9;
245     break;
246   case 7:     // 278385 ft. or 84.8520 km
247     slope     = 0;
248     reftemp   = 336.5; // in degrees Rankine, 186.94 Kelvin
249     refpress  = 0.00831;
250     //refdens   = 2.19541e-10;
251     break;
252   default:     // sea level
253     slope     = -0.00356616; // R/ft.
254     reftemp   = 518.67;   // in degrees Rankine, 288.15 Kelvin
255     refpress  = 2116.22;    // psf
256     //refdens   = 0.00237767;  // slugs/cubic ft.
257     break;
258
259   }
260
261   // If delta_T is set, then that is our temperature deviation at any altitude.
262   // If not, then we'll estimate a deviation based on the sea level deviation (if set).
263
264   if(!StandardTempOnly) {
265     T_dev = 0.0;
266     if (delta_T != 0.0) {
267       T_dev = delta_T;
268     } else {
269       if ((altitude < 36089.239) && (T_dev_sl != 0.0)) {
270         T_dev = T_dev_sl * ( 1.0 - (altitude/36089.239));
271       }
272     }
273     reftemp+=T_dev;
274   }
275
276   if (slope == 0) {
277     intTemperature = reftemp;
278     intPressure = refpress*exp(-FDMExec->GetInertial()->SLgravity()/(reftemp*Reng)*(altitude-htab[i]));
279     intDensity = intPressure/(Reng*intTemperature);
280   } else {
281     intTemperature = reftemp+slope*(altitude-htab[i]);
282     intPressure = refpress*pow(intTemperature/reftemp,-FDMExec->GetInertial()->SLgravity()/(slope*Reng));
283     intDensity = intPressure/(Reng*intTemperature);
284   }
285   
286   lastIndex=i;
287 }
288
289 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
290 // Calculate parameters derived from T, P and rho
291 // Sum gust and turbulence values in NED frame into the wind vector.
292
293 void FGAtmosphere::CalculateDerived(void)
294 {
295   T_dev = (*temperature) - GetTemperature(h);
296
297   if (T_dev == 0.0) density_altitude = h;
298   else              density_altitude = 518.67/0.00356616 * (1.0 - pow(GetDensityRatio(),0.235));
299
300   if (turbType != ttNone) Turbulence();
301
302   vTotalWindNED = vWindNED + vGustNED + vTurbulenceNED;
303
304    // psiw (Wind heading) is the direction the wind is blowing towards
305   if (vWindNED(eX) != 0.0) psiw = atan2( vWindNED(eY), vWindNED(eX) );
306   if (psiw < 0) psiw += 2*M_PI;
307
308   soundspeed = sqrt(SHRatio*Reng*(*temperature));
309
310   intViscosity = Beta * pow(intTemperature, 1.5) / (SutherlandConstant + intTemperature);
311   intKinematicViscosity = intViscosity / intDensity;
312 }
313
314
315 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
316 // Get the standard atmospheric properties at a specified altitude
317
318 void FGAtmosphere::GetStdAtmosphere(double altitude) {
319   StandardTempOnly = true;
320   Calculate(altitude);
321   StandardTempOnly = false;
322   atmosphere.Temperature = intTemperature;
323   atmosphere.Pressure = intPressure;
324   atmosphere.Density = intDensity;
325
326   // Reset the internal atmospheric state
327   Calculate(h);
328 }
329
330 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
331 // Get the standard pressure at a specified altitude
332
333 double FGAtmosphere::GetPressure(double altitude) {
334   GetStdAtmosphere(altitude);
335   return atmosphere.Pressure;
336 }
337
338 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
339 // Get the standard temperature at a specified altitude
340
341 double FGAtmosphere::GetTemperature(double altitude) {
342   GetStdAtmosphere(altitude);
343   return atmosphere.Temperature;
344 }
345
346 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
347 // Get the standard density at a specified altitude
348
349 double FGAtmosphere::GetDensity(double altitude) {
350   GetStdAtmosphere(altitude);
351   return atmosphere.Density;
352 }
353
354
355 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
356 // square a value, but preserve the original sign
357
358 static inline double square_signed (double value)
359 {
360     if (value < 0)
361         return value * value * -1;
362     else
363         return value * value;
364 }
365
366 /// simply square a value
367 static inline double sqr(double x) { return x*x; }
368
369 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
370 //
371 // psi is the angle that the wind is blowing *towards*
372
373 void FGAtmosphere::SetWindspeed(double speed)
374 {
375   if (vWindNED.Magnitude() == 0.0) {
376     psiw = 0.0;
377     vWindNED(eNorth) = speed;
378   } else {
379     vWindNED(eNorth) = speed * cos(psiw);
380     vWindNED(eEast) = speed * sin(psiw);
381     vWindNED(eDown) = 0.0;
382   }
383 }
384
385 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
386
387 double FGAtmosphere::GetWindspeed(void) const
388 {
389   return vWindNED.Magnitude();
390 }
391
392 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
393 //
394 // psi is the angle that the wind is blowing *towards*
395
396 void FGAtmosphere::SetWindPsi(double dir)
397 {
398   double mag = GetWindspeed();
399   psiw = dir;
400   SetWindspeed(mag);  
401 }
402
403 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
404
405 void FGAtmosphere::Turbulence(void)
406 {
407   const double DeltaT = rate*FDMExec->GetDeltaT();
408   const double wingspan = FDMExec->GetAircraft()->GetWingSpan();
409   const double HOverBMAC = FDMExec->GetAuxiliary()->GetHOverBMAC();
410   const FGMatrix33& Tl2b = FDMExec->GetPropagate()->GetTl2b();
411   const double HTailArm = FDMExec->GetAircraft()->GetHTailArm();
412   const double VTailArm = FDMExec->GetAircraft()->GetVTailArm();
413
414   switch (turbType) {
415   case ttStandard: {
416     // TurbGain = TurbGain * TurbGain * 100.0; // what is this!?
417
418     vDirectiondAccelDt(eX) = 1 - 2.0*(double(rand())/double(RAND_MAX));
419     vDirectiondAccelDt(eY) = 1 - 2.0*(double(rand())/double(RAND_MAX));
420     vDirectiondAccelDt(eZ) = 1 - 2.0*(double(rand())/double(RAND_MAX));
421
422     MagnitudedAccelDt = 1 - 2.0*(double(rand())/double(RAND_MAX)) - Magnitude;
423                                 // Scale the magnitude so that it moves
424                                 // away from the peaks
425     MagnitudedAccelDt = ((MagnitudedAccelDt - Magnitude) /
426                          (1 + fabs(Magnitude)));
427     MagnitudeAccel    += MagnitudedAccelDt*TurbRate*DeltaT;
428     Magnitude         += MagnitudeAccel*DeltaT;
429     Magnitude          = fabs(Magnitude);
430
431     vDirectiondAccelDt.Normalize();
432
433                                 // deemphasise non-vertical forces
434     vDirectiondAccelDt(eX) = square_signed(vDirectiondAccelDt(eX));
435     vDirectiondAccelDt(eY) = square_signed(vDirectiondAccelDt(eY));
436
437     vDirectionAccel += vDirectiondAccelDt*TurbRate*DeltaT;
438     vDirectionAccel.Normalize();
439     vDirection      += vDirectionAccel*DeltaT;
440
441     vDirection.Normalize();
442
443                                 // Diminish turbulence within three wingspans
444                                 // of the ground
445     vTurbulenceNED = TurbGain * Magnitude * vDirection;
446     if (HOverBMAC < 3.0)
447         vTurbulenceNED *= (HOverBMAC / 3.0) * (HOverBMAC / 3.0);
448
449     // I don't believe these next two statements calculate the proper gradient over
450     // the aircraft body. One reason is because this has no relationship with the
451     // orientation or velocity of the aircraft, which it must have. What is vTurbulenceGrad
452     // supposed to represent? And the direction and magnitude of the turbulence can change,
453     // so both accelerations need to be accounted for, no?
454
455     // Need to determine the turbulence change in body axes between two time points.
456
457     vTurbulenceGrad = TurbGain*MagnitudeAccel * vDirection;
458     vBodyTurbGrad = Tl2b*vTurbulenceGrad;
459
460     if (wingspan > 0) {
461       vTurbPQR(eP) = vBodyTurbGrad(eY)/wingspan;
462     } else {
463       vTurbPQR(eP) = vBodyTurbGrad(eY)/30.0;
464     }
465 //     if (HTailArm != 0.0)
466 //       vTurbPQR(eQ) = vBodyTurbGrad(eZ)/HTailArm;
467 //     else
468 //       vTurbPQR(eQ) = vBodyTurbGrad(eZ)/10.0;
469
470     if (VTailArm > 0)
471       vTurbPQR(eR) = vBodyTurbGrad(eX)/VTailArm;
472     else
473       vTurbPQR(eR) = vBodyTurbGrad(eX)/10.0;
474
475                                 // Clear the horizontal forces
476                                 // actually felt by the plane, now
477                                 // that we've used them to calculate
478                                 // moments.
479                                 // Why? (JSB)
480 //    vTurbulenceNED(eX) = 0.0;
481 //    vTurbulenceNED(eY) = 0.0;
482
483     break;
484   }
485   case ttBerndt: { // This is very experimental and incomplete at the moment.
486
487     vDirectiondAccelDt(eX) = GaussianRandomNumber();
488     vDirectiondAccelDt(eY) = GaussianRandomNumber();
489     vDirectiondAccelDt(eZ) = GaussianRandomNumber();
490 /*
491     MagnitudedAccelDt = GaussianRandomNumber();
492     MagnitudeAccel    += MagnitudedAccelDt * DeltaT;
493     Magnitude         += MagnitudeAccel * DeltaT;
494 */
495     Magnitude         += GaussianRandomNumber() * DeltaT;
496
497     vDirectiondAccelDt.Normalize();
498     vDirectionAccel += TurbRate * vDirectiondAccelDt * DeltaT;
499     vDirectionAccel.Normalize();
500     vDirection      += vDirectionAccel*DeltaT;
501
502     // Diminish z-vector within two wingspans of the ground
503     if (HOverBMAC < 2.0) vDirection(eZ) *= HOverBMAC / 2.0;
504
505     vDirection.Normalize();
506
507     vTurbulenceNED = TurbGain*Magnitude * vDirection;
508     vTurbulenceGrad = TurbGain*MagnitudeAccel * vDirection;
509
510     vBodyTurbGrad = Tl2b * vTurbulenceGrad;
511     vTurbPQR(eP) = vBodyTurbGrad(eY) / wingspan;
512     if (HTailArm > 0)
513       vTurbPQR(eQ) = vBodyTurbGrad(eZ) / HTailArm;
514     else
515       vTurbPQR(eQ) = vBodyTurbGrad(eZ) / 10.0;
516
517     if (VTailArm > 0)
518       vTurbPQR(eR) = vBodyTurbGrad(eX) / VTailArm;
519     else
520       vTurbPQR(eR) = vBodyTurbGrad(eX)/10.0;
521
522     break;
523   }
524   case ttCulp: { 
525
526     vTurbPQR(eP) = wind_from_clockwise;
527     if (TurbGain == 0.0) return;
528   
529     // keep the inputs within allowable limts for this model
530     if (TurbGain < 0.0) TurbGain = 0.0;
531     if (TurbGain > 1.0) TurbGain = 1.0;
532     if (TurbRate < 0.0) TurbRate = 0.0;
533     if (TurbRate > 30.0) TurbRate = 30.0;
534     if (Rhythmicity < 0.0) Rhythmicity = 0.0;
535     if (Rhythmicity > 1.0) Rhythmicity = 1.0;
536
537     // generate a sine wave corresponding to turbulence rate in hertz
538     double time = FDMExec->GetSimTime();
539     double sinewave = sin( time * TurbRate * 6.283185307 );
540
541     double random = 0.0;
542     if (target_time == 0.0) {
543       strength = random = 1 - 2.0*(double(rand())/double(RAND_MAX));
544       target_time = time + 0.71 + (random * 0.5);
545     }
546     if (time > target_time) {
547       spike = 1.0;
548       target_time = 0.0;
549     }    
550
551     // max vertical wind speed in fps, corresponds to TurbGain = 1.0
552     double max_vs = 40;
553
554     vTurbulenceNED(1) = vTurbulenceNED(2) = vTurbulenceNED(3) = 0.0;
555     double delta = strength * max_vs * TurbGain * (1-Rhythmicity) * spike;
556
557     // Vertical component of turbulence.
558     vTurbulenceNED(3) = sinewave * max_vs * TurbGain * Rhythmicity;
559     vTurbulenceNED(3)+= delta;
560     if (HOverBMAC < 3.0)
561         vTurbulenceNED(3) *= HOverBMAC * 0.3333;
562  
563     // Yaw component of turbulence.
564     vTurbulenceNED(1) = sin( delta * 3.0 );
565     vTurbulenceNED(2) = cos( delta * 3.0 );
566
567     // Roll component of turbulence. Clockwise vortex causes left roll.
568     vTurbPQR(eP) += delta * 0.04;
569
570     spike = spike * 0.9;
571     break;
572   }
573   case ttMilspec:
574   case ttTustin: {
575     double V = FDMExec->GetAuxiliary()->GetVt(); // true airspeed in ft/s
576
577     // an index of zero means turbulence is disabled
578     // airspeed occurs as divisor in the code below
579     if (probability_of_exceedence_index == 0 || V == 0) {
580       vTurbulenceNED(1) = vTurbulenceNED(2) = vTurbulenceNED(3) = 0.0;
581       vTurbPQR(1) = vTurbPQR(2) = vTurbPQR(3) = 0.0;
582       return;
583     }
584
585     // Turbulence model according to MIL-F-8785C (Flying Qualities of Piloted Aircraft)
586     double
587       h = FDMExec->GetPropagate()->GetDistanceAGL(),
588       b_w = wingspan,
589       L_u, L_w, sig_u, sig_w;
590
591       if (b_w == 0.) b_w = 30.;
592
593     // clip height functions at 10 ft
594     if (h <= 10.) h = 10;
595
596     // Scale lengths L and amplitudes sigma as function of height
597     if (h <= 1000) {
598       L_u = h/pow(0.177 + 0.000823*h, 1.2); // MIL-F-8785c, Fig. 10, p. 55
599       L_w = h;
600       sig_w = 0.1*windspeed_at_20ft;
601       sig_u = sig_w/pow(0.177 + 0.000823*h, 0.4); // MIL-F-8785c, Fig. 11, p. 56
602     } else if (h <= 2000) {
603       // linear interpolation between low altitude and high altitude models
604       L_u = L_w = 1000 + (h-1000.)/1000.*750.;
605       sig_u = sig_w = 0.1*windspeed_at_20ft
606                     + (h-1000.)/1000.*(POE_Table->GetValue(probability_of_exceedence_index, h) - 0.1*windspeed_at_20ft);
607     } else {
608       L_u = L_w = 1750.; //  MIL-F-8785c, Sec. 3.7.2.1, p. 48
609       sig_u = sig_w = POE_Table->GetValue(probability_of_exceedence_index, h);
610     }
611
612     // keep values from last timesteps
613     // TODO maybe use deque?
614     static double
615       xi_u_km1 = 0, nu_u_km1 = 0,
616       xi_v_km1 = 0, xi_v_km2 = 0, nu_v_km1 = 0, nu_v_km2 = 0,
617       xi_w_km1 = 0, xi_w_km2 = 0, nu_w_km1 = 0, nu_w_km2 = 0,
618       xi_p_km1 = 0, nu_p_km1 = 0,
619       xi_q_km1 = 0, xi_r_km1 = 0;
620
621
622     double
623       T_V = DeltaT, // for compatibility of nomenclature
624       sig_p = 1.9/sqrt(L_w*b_w)*sig_w, // Yeager1998, eq. (8)
625       sig_q = sqrt(M_PI/2/L_w/b_w), // eq. (14)
626       sig_r = sqrt(2*M_PI/3/L_w/b_w), // eq. (17)
627       L_p = sqrt(L_w*b_w)/2.6, // eq. (10)
628       tau_u = L_u/V, // eq. (6)
629       tau_w = L_w/V, // eq. (3)
630       tau_p = L_p/V, // eq. (9)
631       tau_q = 4*b_w/M_PI/V, // eq. (13)
632       tau_r =3*b_w/M_PI/V, // eq. (17)
633       nu_u = GaussianRandomNumber(),
634       nu_v = GaussianRandomNumber(),
635       nu_w = GaussianRandomNumber(),
636       nu_p = GaussianRandomNumber(),
637       xi_u=0, xi_v=0, xi_w=0, xi_p=0, xi_q=0, xi_r=0;
638
639     // values of turbulence NED velocities
640
641     if (turbType == ttTustin) {
642       // the following is the Tustin formulation of Yeager's report
643       double
644         omega_w = V/L_w, // hidden in nomenclature p. 3
645         omega_v = V/L_u, // this is defined nowhere
646         C_BL  = 1/tau_u/tan(T_V/2/tau_u), // eq. (19)
647         C_BLp = 1/tau_p/tan(T_V/2/tau_p), // eq. (22)
648         C_BLq = 1/tau_q/tan(T_V/2/tau_q), // eq. (24)
649         C_BLr = 1/tau_r/tan(T_V/2/tau_r); // eq. (26)
650
651       // all values calculated so far are strictly positive, except for
652       // the random numbers nu_*. This means that in the code below, all
653       // divisors are strictly positive, too, and no floating point
654       // exception should occur.
655       xi_u = -(1 - C_BL*tau_u)/(1 + C_BL*tau_u)*xi_u_km1
656            + sig_u*sqrt(2*tau_u/T_V)/(1 + C_BL*tau_u)*(nu_u + nu_u_km1); // eq. (18)
657       xi_v = -2*(sqr(omega_v) - sqr(C_BL))/sqr(omega_v + C_BL)*xi_v_km1
658            - sqr(omega_v - C_BL)/sqr(omega_v + C_BL) * xi_v_km2
659            + sig_u*sqrt(3*omega_v/T_V)/sqr(omega_v + C_BL)*(
660                  (C_BL + omega_v/sqrt(3.))*nu_v
661                + 2/sqrt(3.)*omega_v*nu_v_km1
662                + (omega_v/sqrt(3.) - C_BL)*nu_v_km2); // eq. (20) for v
663       xi_w = -2*(sqr(omega_w) - sqr(C_BL))/sqr(omega_w + C_BL)*xi_w_km1
664            - sqr(omega_w - C_BL)/sqr(omega_w + C_BL) * xi_w_km2
665            + sig_w*sqrt(3*omega_w/T_V)/sqr(omega_w + C_BL)*(
666                  (C_BL + omega_w/sqrt(3.))*nu_w
667                + 2/sqrt(3.)*omega_w*nu_w_km1
668                + (omega_w/sqrt(3.) - C_BL)*nu_w_km2); // eq. (20) for w
669       xi_p = -(1 - C_BLp*tau_p)/(1 + C_BLp*tau_p)*xi_p_km1
670            + sig_p*sqrt(2*tau_p/T_V)/(1 + C_BLp*tau_p) * (nu_p + nu_p_km1); // eq. (21)
671       xi_q = -(1 - 4*b_w*C_BLq/M_PI/V)/(1 + 4*b_w*C_BLq/M_PI/V) * xi_q_km1
672            + C_BLq/V/(1 + 4*b_w*C_BLq/M_PI/V) * (xi_w - xi_w_km1); // eq. (23)
673       xi_r = - (1 - 3*b_w*C_BLr/M_PI/V)/(1 + 3*b_w*C_BLr/M_PI/V) * xi_r_km1
674            + C_BLr/V/(1 + 3*b_w*C_BLr/M_PI/V) * (xi_v - xi_v_km1); // eq. (25)
675
676     } else if (turbType == ttMilspec) {
677       // the following is the MIL-STD-1797A formulation
678       // as cited in Yeager's report
679       xi_u = (1 - T_V/tau_u)  *xi_u_km1 + sig_u*sqrt(2*T_V/tau_u)*nu_u;  // eq. (30)
680       xi_v = (1 - 2*T_V/tau_u)*xi_v_km1 + sig_u*sqrt(4*T_V/tau_u)*nu_v;  // eq. (31)
681       xi_w = (1 - 2*T_V/tau_w)*xi_w_km1 + sig_w*sqrt(4*T_V/tau_w)*nu_w;  // eq. (32)
682       xi_p = (1 - T_V/tau_p)  *xi_p_km1 + sig_p*sqrt(2*T_V/tau_p)*nu_p;  // eq. (33)
683       xi_q = (1 - T_V/tau_q)  *xi_q_km1 + M_PI/4/b_w*(xi_w - xi_w_km1);  // eq. (34)
684       xi_r = (1 - T_V/tau_r)  *xi_r_km1 + M_PI/3/b_w*(xi_v - xi_v_km1);  // eq. (35)
685     }
686
687     // rotate by wind azimuth and assign the velocities
688     double cospsi = cos(psiw), sinpsi = sin(psiw);
689     vTurbulenceNED(1) =  cospsi*xi_u + sinpsi*xi_v;
690     vTurbulenceNED(2) = -sinpsi*xi_u + cospsi*xi_v;
691     vTurbulenceNED(3) = xi_w;
692
693     vTurbPQR(1) =  cospsi*xi_p + sinpsi*xi_q;
694     vTurbPQR(2) = -sinpsi*xi_p + cospsi*xi_q;
695     vTurbPQR(3) = xi_r;
696
697     // vTurbPQR is in the body fixed frame, not NED
698     vTurbPQR = Tl2b*vTurbPQR;
699
700     // hand on the values for the next timestep
701     xi_u_km1 = xi_u; nu_u_km1 = nu_u;
702     xi_v_km2 = xi_v_km1; xi_v_km1 = xi_v; nu_v_km2 = nu_v_km1; nu_v_km1 = nu_v;
703     xi_w_km2 = xi_w_km1; xi_w_km1 = xi_w; nu_w_km2 = nu_w_km1; nu_w_km1 = nu_w;
704     xi_p_km1 = xi_p; nu_p_km1 = nu_p;
705     xi_q_km1 = xi_q;
706     xi_r_km1 = xi_r;
707
708   }
709   default:
710     break;
711   }
712 }
713
714 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
715
716 void FGAtmosphere::UseExternal(void)
717 {
718   temperature=&exTemperature;
719   pressure=&exPressure;
720   density=&exDensity;
721   useExternal=true;
722 }
723
724 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
725
726 void FGAtmosphere::UseInternal(void)
727 {
728   temperature=&intTemperature;
729   pressure=&intPressure;
730   density=&intDensity;
731   useExternal=false;
732 }
733
734 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
735
736 void FGAtmosphere::bind(void)
737 {
738   typedef double (FGAtmosphere::*PMF)(int) const;
739   typedef double (FGAtmosphere::*PMFv)(void) const;
740   typedef int (FGAtmosphere::*PMFt)(void) const;
741   typedef void   (FGAtmosphere::*PMFd)(int,double);
742   typedef void   (FGAtmosphere::*PMFi)(int);
743   PropertyManager->Tie("atmosphere/T-R", this, (PMFv)&FGAtmosphere::GetTemperature);
744   PropertyManager->Tie("atmosphere/rho-slugs_ft3", this, (PMFv)&FGAtmosphere::GetDensity);
745   PropertyManager->Tie("atmosphere/P-psf", this, (PMFv)&FGAtmosphere::GetPressure);
746   PropertyManager->Tie("atmosphere/a-fps", this, &FGAtmosphere::GetSoundSpeed);
747   PropertyManager->Tie("atmosphere/T-sl-R", this, &FGAtmosphere::GetTemperatureSL);
748   PropertyManager->Tie("atmosphere/rho-sl-slugs_ft3", this, &FGAtmosphere::GetDensitySL);
749   PropertyManager->Tie("atmosphere/P-sl-psf", this, &FGAtmosphere::GetPressureSL);
750   PropertyManager->Tie("atmosphere/a-sl-fps", this, &FGAtmosphere::GetSoundSpeedSL);
751   PropertyManager->Tie("atmosphere/theta", this, &FGAtmosphere::GetTemperatureRatio);
752   PropertyManager->Tie("atmosphere/sigma", this, &FGAtmosphere::GetDensityRatio);
753   PropertyManager->Tie("atmosphere/delta", this, &FGAtmosphere::GetPressureRatio);
754   PropertyManager->Tie("atmosphere/a-ratio", this, &FGAtmosphere::GetSoundSpeedRatio);
755   PropertyManager->Tie("atmosphere/psiw-rad", this, &FGAtmosphere::GetWindPsi, &FGAtmosphere::SetWindPsi);
756   PropertyManager->Tie("atmosphere/delta-T", this, &FGAtmosphere::GetDeltaT, &FGAtmosphere::SetDeltaT);
757   PropertyManager->Tie("atmosphere/T-sl-dev-F", this, &FGAtmosphere::GetSLTempDev, &FGAtmosphere::SetSLTempDev);
758   PropertyManager->Tie("atmosphere/density-altitude", this, &FGAtmosphere::GetDensityAltitude);
759
760   PropertyManager->Tie("atmosphere/wind-north-fps", this, eNorth, (PMF)&FGAtmosphere::GetWindNED,
761                                                           (PMFd)&FGAtmosphere::SetWindNED);
762   PropertyManager->Tie("atmosphere/wind-east-fps",  this, eEast, (PMF)&FGAtmosphere::GetWindNED,
763                                                           (PMFd)&FGAtmosphere::SetWindNED);
764   PropertyManager->Tie("atmosphere/wind-down-fps",  this, eDown, (PMF)&FGAtmosphere::GetWindNED,
765                                                           (PMFd)&FGAtmosphere::SetWindNED);
766   PropertyManager->Tie("atmosphere/wind-mag-fps", this, &FGAtmosphere::GetWindspeed,
767                                                         &FGAtmosphere::SetWindspeed);
768   PropertyManager->Tie("atmosphere/total-wind-north-fps", this, eNorth, (PMF)&FGAtmosphere::GetTotalWindNED);
769   PropertyManager->Tie("atmosphere/total-wind-east-fps",  this, eEast,  (PMF)&FGAtmosphere::GetTotalWindNED);
770   PropertyManager->Tie("atmosphere/total-wind-down-fps",  this, eDown,  (PMF)&FGAtmosphere::GetTotalWindNED);
771
772   PropertyManager->Tie("atmosphere/gust-north-fps", this, eNorth, (PMF)&FGAtmosphere::GetGustNED,
773                                                           (PMFd)&FGAtmosphere::SetGustNED);
774   PropertyManager->Tie("atmosphere/gust-east-fps",  this, eEast, (PMF)&FGAtmosphere::GetGustNED,
775                                                           (PMFd)&FGAtmosphere::SetGustNED);
776   PropertyManager->Tie("atmosphere/gust-down-fps",  this, eDown, (PMF)&FGAtmosphere::GetGustNED,
777                                                           (PMFd)&FGAtmosphere::SetGustNED);
778
779   PropertyManager->Tie("atmosphere/turb-north-fps", this, eNorth, (PMF)&FGAtmosphere::GetTurbNED,
780                                                           (PMFd)&FGAtmosphere::SetTurbNED);
781   PropertyManager->Tie("atmosphere/turb-east-fps",  this, eEast, (PMF)&FGAtmosphere::GetTurbNED,
782                                                           (PMFd)&FGAtmosphere::SetTurbNED);
783   PropertyManager->Tie("atmosphere/turb-down-fps",  this, eDown, (PMF)&FGAtmosphere::GetTurbNED,
784                                                           (PMFd)&FGAtmosphere::SetTurbNED);
785
786   PropertyManager->Tie("atmosphere/p-turb-rad_sec", this,1, (PMF)&FGAtmosphere::GetTurbPQR);
787   PropertyManager->Tie("atmosphere/q-turb-rad_sec", this,2, (PMF)&FGAtmosphere::GetTurbPQR);
788   PropertyManager->Tie("atmosphere/r-turb-rad_sec", this,3, (PMF)&FGAtmosphere::GetTurbPQR);
789   PropertyManager->Tie("atmosphere/turb-type", this, (PMFt)&FGAtmosphere::GetTurbType, (PMFi)&FGAtmosphere::SetTurbType);
790   PropertyManager->Tie("atmosphere/turb-rate", this, &FGAtmosphere::GetTurbRate, &FGAtmosphere::SetTurbRate);
791   PropertyManager->Tie("atmosphere/turb-gain", this, &FGAtmosphere::GetTurbGain, &FGAtmosphere::SetTurbGain);
792   PropertyManager->Tie("atmosphere/turb-rhythmicity", this, &FGAtmosphere::GetRhythmicity,
793                                                             &FGAtmosphere::SetRhythmicity);
794
795   PropertyManager->Tie("atmosphere/turbulence/milspec/windspeed_at_20ft_AGL-fps",
796                        this, &FGAtmosphere::GetWindspeed20ft,
797                              &FGAtmosphere::SetWindspeed20ft);
798   PropertyManager->Tie("atmosphere/turbulence/milspec/severity",
799                        this, &FGAtmosphere::GetProbabilityOfExceedence,
800                              &FGAtmosphere::SetProbabilityOfExceedence);
801
802 }
803
804 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
805 //    The bitmasked value choices are as follows:
806 //    unset: In this case (the default) JSBSim would only print
807 //       out the normally expected messages, essentially echoing
808 //       the config files as they are read. If the environment
809 //       variable is not set, debug_lvl is set to 1 internally
810 //    0: This requests JSBSim not to output any messages
811 //       whatsoever.
812 //    1: This value explicity requests the normal JSBSim
813 //       startup messages
814 //    2: This value asks for a message to be printed out when
815 //       a class is instantiated
816 //    4: When this value is set, a message is displayed when a
817 //       FGModel object executes its Run() method
818 //    8: When this value is set, various runtime state variables
819 //       are printed out periodically
820 //    16: When set various parameters are sanity checked and
821 //       a message is printed out when they go out of bounds
822
823 void FGAtmosphere::Debug(int from)
824 {
825   if (debug_lvl <= 0) return;
826
827   if (debug_lvl & 1) { // Standard console startup message output
828     if (from == 0) { // Constructor
829     }
830   }
831   if (debug_lvl & 2 ) { // Instantiation/Destruction notification
832     if (from == 0) cout << "Instantiated: FGAtmosphere" << endl;
833     if (from == 1) cout << "Destroyed:    FGAtmosphere" << endl;
834   }
835   if (debug_lvl & 4 ) { // Run() method entry print for FGModel-derived objects
836   }
837   if (debug_lvl & 8 ) { // Runtime state variables
838   }
839   if (debug_lvl & 16) { // Sanity checking
840   }
841   if (debug_lvl & 128) { // Turbulence
842     if (first_pass && from == 2) {
843       first_pass = false;
844       cout << "vTurbulenceNED(X), vTurbulenceNED(Y), vTurbulenceNED(Z), "
845            << "vTurbulenceGrad(X), vTurbulenceGrad(Y), vTurbulenceGrad(Z), "
846            << "vDirection(X), vDirection(Y), vDirection(Z), "
847            << "Magnitude, "
848            << "vTurbPQR(P), vTurbPQR(Q), vTurbPQR(R), " << endl;
849     } 
850     if (from == 2) {
851       cout << vTurbulenceNED << ", " << vTurbulenceGrad << ", " << vDirection << ", " << Magnitude << ", " << vTurbPQR << endl;
852     }
853   }
854   if (debug_lvl & 64) {
855     if (from == 0) { // Constructor
856       cout << IdSrc << endl;
857       cout << IdHdr << endl;
858     }
859   }
860 }
861
862 } // namespace JSBSim