JSBSim Flight Dynamics Model 1.0 (23 February 2013)
An Open Source Flight Dynamics and Control Software Library in C++

FGWinds.cpp

00001 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00002 
00003  Module:       FGWinds.cpp
00004  Author:       Jon Berndt, Tony Peden, Andreas Gaeb
00005  Date started: Extracted from FGAtmosphere, which originated in 1998
00006                5/2011
00007  Purpose:      Models winds, gusts, turbulence, and other atmospheric disturbances
00008  Called by:    FGFDMExec
00009 
00010  ------------- Copyright (C) 2011  Jon S. Berndt (jon@jsbsim.org) -------------
00011 
00012  This program is free software; you can redistribute it and/or modify it under
00013  the terms of the GNU Lesser General Public License as published by the Free Software
00014  Foundation; either version 2 of the License, or (at your option) any later
00015  version.
00016 
00017  This program is distributed in the hope that it will be useful, but WITHOUT
00018  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00019  FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public License for more
00020  details.
00021 
00022  You should have received a copy of the GNU Lesser General Public License along with
00023  this program; if not, write to the Free Software Foundation, Inc., 59 Temple
00024  Place - Suite 330, Boston, MA  02111-1307, USA.
00025 
00026  Further information about the GNU Lesser General Public License can also be found on
00027  the world wide web at http://www.gnu.org.
00028 
00029 FUNCTIONAL DESCRIPTION
00030 --------------------------------------------------------------------------------
00031 
00032 HISTORY
00033 --------------------------------------------------------------------------------
00034 
00035 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00036 COMMENTS, REFERENCES,  and NOTES
00037 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00038 [1]   Anderson, John D. "Introduction to Flight, Third Edition", McGraw-Hill,
00039       1989, ISBN 0-07-001641-0
00040 
00041 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00042 INCLUDES
00043 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
00044 
00045 #include <iostream>
00046 #include <cstdlib>
00047 #include "FGWinds.h"
00048 #include "FGFDMExec.h"
00049 
00050 using namespace std;
00051 
00052 namespace JSBSim {
00053 
00054 static const char *IdSrc = "$Id: FGWinds.cpp,v 1.8 2012/12/02 12:59:19 bcoconni Exp $";
00055 static const char *IdHdr = ID_WINDS;
00056 
00057 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00058 CLASS IMPLEMENTATION
00059 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
00060 
00061 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00062 // square a value, but preserve the original sign
00063 
00064 static inline double square_signed (double value)
00065 {
00066     if (value < 0)
00067         return value * value * -1;
00068     else
00069         return value * value;
00070 }
00071 
00073 static inline double sqr(double x) { return x*x; }
00074 
00075 FGWinds::FGWinds(FGFDMExec* fdmex) : FGModel(fdmex)
00076 {
00077   Name = "FGWinds";
00078 
00079   MagnitudedAccelDt = MagnitudeAccel = Magnitude = 0.0;
00080   SetTurbType( ttMilspec );
00081   TurbGain = 1.0;
00082   TurbRate = 10.0;
00083   Rhythmicity = 0.1;
00084   spike = target_time = strength = 0.0;
00085   wind_from_clockwise = 0.0;
00086   psiw = 0.0;
00087 
00088   vGustNED.InitMatrix();
00089   vTurbulenceNED.InitMatrix();
00090 
00091   // Milspec turbulence model
00092   windspeed_at_20ft = 0.;
00093   probability_of_exceedence_index = 0;
00094   POE_Table = new FGTable(7,12);
00095   // this is Figure 7 from p. 49 of MIL-F-8785C
00096   // rows: probability of exceedance curve index, cols: altitude in ft
00097   *POE_Table
00098            << 500.0 << 1750.0 << 3750.0 << 7500.0 << 15000.0 << 25000.0 << 35000.0 << 45000.0 << 55000.0 << 65000.0 << 75000.0 << 80000.0
00099     << 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
00100     << 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
00101     << 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
00102     << 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
00103     << 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
00104     << 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
00105     << 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;
00106 
00107   bind();
00108   Debug(0);
00109 }
00110 
00111 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00112 
00113 FGWinds::~FGWinds()
00114 {
00115   delete(POE_Table);
00116   Debug(1);
00117 }
00118 
00119 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00120 
00121 bool FGWinds::InitModel(void)
00122 {
00123   return true;
00124 }
00125 
00126 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00127 
00128 bool FGWinds::Run(bool Holding)
00129 {
00130   if (FGModel::Run(Holding)) return true;
00131   if (Holding) return false;
00132 
00133   if (turbType != ttNone) Turbulence(in.AltitudeASL);
00134   if (oneMinusCosineGust.gustProfile.Running) CosineGust();
00135 
00136   vTotalWindNED = vWindNED + vGustNED + vCosineGust + vTurbulenceNED;
00137 
00138    // psiw (Wind heading) is the direction the wind is blowing towards
00139   if (vWindNED(eX) != 0.0) psiw = atan2( vWindNED(eY), vWindNED(eX) );
00140   if (psiw < 0) psiw += 2*M_PI;
00141 
00142   Debug(2);
00143   return false;
00144 }
00145 
00146 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00147 //
00148 // psi is the angle that the wind is blowing *towards*
00149 
00150 void FGWinds::SetWindspeed(double speed)
00151 {
00152   if (vWindNED.Magnitude() == 0.0) {
00153     psiw = 0.0;
00154     vWindNED(eNorth) = speed;
00155   } else {
00156     vWindNED(eNorth) = speed * cos(psiw);
00157     vWindNED(eEast) = speed * sin(psiw);
00158     vWindNED(eDown) = 0.0;
00159   }
00160 }
00161 
00162 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00163 
00164 double FGWinds::GetWindspeed(void) const
00165 {
00166   return vWindNED.Magnitude();
00167 }
00168 
00169 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00170 //
00171 // psi is the angle that the wind is blowing *towards*
00172 
00173 void FGWinds::SetWindPsi(double dir)
00174 {
00175   double mag = GetWindspeed();
00176   psiw = dir;
00177   SetWindspeed(mag);
00178 }
00179 
00180 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00181 
00182 void FGWinds::Turbulence(double h)
00183 {
00184   switch (turbType) {
00185 
00186   case ttCulp: {
00187 
00188     vTurbPQR(eP) = wind_from_clockwise;
00189     if (TurbGain == 0.0) return;
00190 
00191     // keep the inputs within allowable limts for this model
00192     if (TurbGain < 0.0) TurbGain = 0.0;
00193     if (TurbGain > 1.0) TurbGain = 1.0;
00194     if (TurbRate < 0.0) TurbRate = 0.0;
00195     if (TurbRate > 30.0) TurbRate = 30.0;
00196     if (Rhythmicity < 0.0) Rhythmicity = 0.0;
00197     if (Rhythmicity > 1.0) Rhythmicity = 1.0;
00198 
00199     // generate a sine wave corresponding to turbulence rate in hertz
00200     double time = FDMExec->GetSimTime();
00201     double sinewave = sin( time * TurbRate * 6.283185307 );
00202 
00203     double random = 0.0;
00204     if (target_time == 0.0) {
00205       strength = random = 1 - 2.0*(double(rand())/double(RAND_MAX));
00206       target_time = time + 0.71 + (random * 0.5);
00207     }
00208     if (time > target_time) {
00209       spike = 1.0;
00210       target_time = 0.0;
00211     }
00212 
00213     // max vertical wind speed in fps, corresponds to TurbGain = 1.0
00214     double max_vs = 40;
00215 
00216     vTurbulenceNED(1) = vTurbulenceNED(2) = vTurbulenceNED(3) = 0.0;
00217     double delta = strength * max_vs * TurbGain * (1-Rhythmicity) * spike;
00218 
00219     // Vertical component of turbulence.
00220     vTurbulenceNED(3) = sinewave * max_vs * TurbGain * Rhythmicity;
00221     vTurbulenceNED(3)+= delta;
00222     if (in.DistanceAGL/in.wingspan < 3.0)
00223         vTurbulenceNED(3) *= in.DistanceAGL/in.wingspan * 0.3333;
00224 
00225     // Yaw component of turbulence.
00226     vTurbulenceNED(1) = sin( delta * 3.0 );
00227     vTurbulenceNED(2) = cos( delta * 3.0 );
00228 
00229     // Roll component of turbulence. Clockwise vortex causes left roll.
00230     vTurbPQR(eP) += delta * 0.04;
00231 
00232     spike = spike * 0.9;
00233     break;
00234   }
00235   case ttMilspec:
00236   case ttTustin: {
00237 
00238     // an index of zero means turbulence is disabled
00239     // airspeed occurs as divisor in the code below
00240     if (probability_of_exceedence_index == 0 || in.V == 0) {
00241       vTurbulenceNED(1) = vTurbulenceNED(2) = vTurbulenceNED(3) = 0.0;
00242       vTurbPQR(1) = vTurbPQR(2) = vTurbPQR(3) = 0.0;
00243       return;
00244     }
00245 
00246     // Turbulence model according to MIL-F-8785C (Flying Qualities of Piloted Aircraft)
00247     double b_w = in.wingspan, L_u, L_w, sig_u, sig_w;
00248 
00249       if (b_w == 0.) b_w = 30.;
00250 
00251     // clip height functions at 10 ft
00252     if (h <= 10.) h = 10;
00253 
00254     // Scale lengths L and amplitudes sigma as function of height
00255     if (h <= 1000) {
00256       L_u = h/pow(0.177 + 0.000823*h, 1.2); // MIL-F-8785c, Fig. 10, p. 55
00257       L_w = h;
00258       sig_w = 0.1*windspeed_at_20ft;
00259       sig_u = sig_w/pow(0.177 + 0.000823*h, 0.4); // MIL-F-8785c, Fig. 11, p. 56
00260     } else if (h <= 2000) {
00261       // linear interpolation between low altitude and high altitude models
00262       L_u = L_w = 1000 + (h-1000.)/1000.*750.;
00263       sig_u = sig_w = 0.1*windspeed_at_20ft
00264                     + (h-1000.)/1000.*(POE_Table->GetValue(probability_of_exceedence_index, h) - 0.1*windspeed_at_20ft);
00265     } else {
00266       L_u = L_w = 1750.; //  MIL-F-8785c, Sec. 3.7.2.1, p. 48
00267       sig_u = sig_w = POE_Table->GetValue(probability_of_exceedence_index, h);
00268     }
00269 
00270     // keep values from last timesteps
00271     // TODO maybe use deque?
00272     static double
00273       xi_u_km1 = 0, nu_u_km1 = 0,
00274       xi_v_km1 = 0, xi_v_km2 = 0, nu_v_km1 = 0, nu_v_km2 = 0,
00275       xi_w_km1 = 0, xi_w_km2 = 0, nu_w_km1 = 0, nu_w_km2 = 0,
00276       xi_p_km1 = 0, nu_p_km1 = 0,
00277       xi_q_km1 = 0, xi_r_km1 = 0;
00278 
00279 
00280     double
00281       T_V = in.totalDeltaT, // for compatibility of nomenclature
00282       sig_p = 1.9/sqrt(L_w*b_w)*sig_w, // Yeager1998, eq. (8)
00283       //sig_q = sqrt(M_PI/2/L_w/b_w), // eq. (14)
00284       //sig_r = sqrt(2*M_PI/3/L_w/b_w), // eq. (17)
00285       L_p = sqrt(L_w*b_w)/2.6, // eq. (10)
00286       tau_u = L_u/in.V, // eq. (6)
00287       tau_w = L_w/in.V, // eq. (3)
00288       tau_p = L_p/in.V, // eq. (9)
00289       tau_q = 4*b_w/M_PI/in.V, // eq. (13)
00290       tau_r =3*b_w/M_PI/in.V, // eq. (17)
00291       nu_u = GaussianRandomNumber(),
00292       nu_v = GaussianRandomNumber(),
00293       nu_w = GaussianRandomNumber(),
00294       nu_p = GaussianRandomNumber(),
00295       xi_u=0, xi_v=0, xi_w=0, xi_p=0, xi_q=0, xi_r=0;
00296 
00297     // values of turbulence NED velocities
00298 
00299     if (turbType == ttTustin) {
00300       // the following is the Tustin formulation of Yeager's report
00301       double
00302         omega_w = in.V/L_w, // hidden in nomenclature p. 3
00303         omega_v = in.V/L_u, // this is defined nowhere
00304         C_BL  = 1/tau_u/tan(T_V/2/tau_u), // eq. (19)
00305         C_BLp = 1/tau_p/tan(T_V/2/tau_p), // eq. (22)
00306         C_BLq = 1/tau_q/tan(T_V/2/tau_q), // eq. (24)
00307         C_BLr = 1/tau_r/tan(T_V/2/tau_r); // eq. (26)
00308 
00309       // all values calculated so far are strictly positive, except for
00310       // the random numbers nu_*. This means that in the code below, all
00311       // divisors are strictly positive, too, and no floating point
00312       // exception should occur.
00313       xi_u = -(1 - C_BL*tau_u)/(1 + C_BL*tau_u)*xi_u_km1
00314            + sig_u*sqrt(2*tau_u/T_V)/(1 + C_BL*tau_u)*(nu_u + nu_u_km1); // eq. (18)
00315       xi_v = -2*(sqr(omega_v) - sqr(C_BL))/sqr(omega_v + C_BL)*xi_v_km1
00316            - sqr(omega_v - C_BL)/sqr(omega_v + C_BL) * xi_v_km2
00317            + sig_u*sqrt(3*omega_v/T_V)/sqr(omega_v + C_BL)*(
00318                  (C_BL + omega_v/sqrt(3.))*nu_v
00319                + 2/sqrt(3.)*omega_v*nu_v_km1
00320                + (omega_v/sqrt(3.) - C_BL)*nu_v_km2); // eq. (20) for v
00321       xi_w = -2*(sqr(omega_w) - sqr(C_BL))/sqr(omega_w + C_BL)*xi_w_km1
00322            - sqr(omega_w - C_BL)/sqr(omega_w + C_BL) * xi_w_km2
00323            + sig_w*sqrt(3*omega_w/T_V)/sqr(omega_w + C_BL)*(
00324                  (C_BL + omega_w/sqrt(3.))*nu_w
00325                + 2/sqrt(3.)*omega_w*nu_w_km1
00326                + (omega_w/sqrt(3.) - C_BL)*nu_w_km2); // eq. (20) for w
00327       xi_p = -(1 - C_BLp*tau_p)/(1 + C_BLp*tau_p)*xi_p_km1
00328            + sig_p*sqrt(2*tau_p/T_V)/(1 + C_BLp*tau_p) * (nu_p + nu_p_km1); // eq. (21)
00329       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
00330            + C_BLq/in.V/(1 + 4*b_w*C_BLq/M_PI/in.V) * (xi_w - xi_w_km1); // eq. (23)
00331       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
00332            + C_BLr/in.V/(1 + 3*b_w*C_BLr/M_PI/in.V) * (xi_v - xi_v_km1); // eq. (25)
00333 
00334     } else if (turbType == ttMilspec) {
00335       // the following is the MIL-STD-1797A formulation
00336       // as cited in Yeager's report
00337       xi_u = (1 - T_V/tau_u)  *xi_u_km1 + sig_u*sqrt(2*T_V/tau_u)*nu_u;  // eq. (30)
00338       xi_v = (1 - 2*T_V/tau_u)*xi_v_km1 + sig_u*sqrt(4*T_V/tau_u)*nu_v;  // eq. (31)
00339       xi_w = (1 - 2*T_V/tau_w)*xi_w_km1 + sig_w*sqrt(4*T_V/tau_w)*nu_w;  // eq. (32)
00340       xi_p = (1 - T_V/tau_p)  *xi_p_km1 + sig_p*sqrt(2*T_V/tau_p)*nu_p;  // eq. (33)
00341       xi_q = (1 - T_V/tau_q)  *xi_q_km1 + M_PI/4/b_w*(xi_w - xi_w_km1);  // eq. (34)
00342       xi_r = (1 - T_V/tau_r)  *xi_r_km1 + M_PI/3/b_w*(xi_v - xi_v_km1);  // eq. (35)
00343     }
00344 
00345     // rotate by wind azimuth and assign the velocities
00346     double cospsi = cos(psiw), sinpsi = sin(psiw);
00347     vTurbulenceNED(1) =  cospsi*xi_u + sinpsi*xi_v;
00348     vTurbulenceNED(2) = -sinpsi*xi_u + cospsi*xi_v;
00349     vTurbulenceNED(3) = xi_w;
00350 
00351     vTurbPQR(1) =  cospsi*xi_p + sinpsi*xi_q;
00352     vTurbPQR(2) = -sinpsi*xi_p + cospsi*xi_q;
00353     vTurbPQR(3) = xi_r;
00354 
00355     // vTurbPQR is in the body fixed frame, not NED
00356     vTurbPQR = in.Tl2b*vTurbPQR;
00357 
00358     // hand on the values for the next timestep
00359     xi_u_km1 = xi_u; nu_u_km1 = nu_u;
00360     xi_v_km2 = xi_v_km1; xi_v_km1 = xi_v; nu_v_km2 = nu_v_km1; nu_v_km1 = nu_v;
00361     xi_w_km2 = xi_w_km1; xi_w_km1 = xi_w; nu_w_km2 = nu_w_km1; nu_w_km1 = nu_w;
00362     xi_p_km1 = xi_p; nu_p_km1 = nu_p;
00363     xi_q_km1 = xi_q;
00364     xi_r_km1 = xi_r;
00365 
00366   }
00367   default:
00368     break;
00369   }
00370 }
00371 
00372 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00373 
00374 double FGWinds::CosineGustProfile(double startDuration, double steadyDuration, double endDuration, double elapsedTime)
00375 {
00376   double factor = 0.0;
00377   if (elapsedTime >= 0 && elapsedTime <= startDuration) {
00378     factor = (1.0 - cos(M_PI*elapsedTime/startDuration))/2.0;
00379   } else if (elapsedTime > startDuration && (elapsedTime <= (startDuration + steadyDuration))) {
00380     factor = 1.0;
00381   } else if (elapsedTime > (startDuration + steadyDuration) && elapsedTime <= (startDuration + steadyDuration + endDuration)) {
00382     factor = (1-cos(M_PI*(1-(elapsedTime-(startDuration + steadyDuration))/endDuration)))/2.0;
00383   } else {
00384     factor = 0.0;
00385   }
00386 
00387   return factor;
00388 }
00389 
00390 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00391 
00392 void FGWinds::CosineGust()
00393 {
00394   struct OneMinusCosineProfile& profile = oneMinusCosineGust.gustProfile;
00395 
00396   double factor = CosineGustProfile( profile.startupDuration,
00397                                      profile.steadyDuration,
00398                                      profile.endDuration,
00399                                      profile.elapsedTime);
00400   // Normalize the gust wind vector
00401   oneMinusCosineGust.vWind.Normalize();
00402 
00403   if (oneMinusCosineGust.vWindTransformed.Magnitude() == 0.0) {
00404     switch (oneMinusCosineGust.gustFrame) {
00405     case gfBody:
00406       oneMinusCosineGust.vWindTransformed = in.Tl2b.Inverse() * oneMinusCosineGust.vWind;
00407       break;
00408     case gfWind:
00409       oneMinusCosineGust.vWindTransformed = in.Tl2b.Inverse() * in.Tw2b * oneMinusCosineGust.vWind;
00410       break;
00411     case gfLocal:
00412       // this is the native frame - and the default.
00413       oneMinusCosineGust.vWindTransformed = oneMinusCosineGust.vWind;
00414       break;
00415     default:
00416       break;
00417     }
00418   }
00419 
00420   vCosineGust = factor * oneMinusCosineGust.vWindTransformed * oneMinusCosineGust.magnitude;
00421 
00422   profile.elapsedTime += in.totalDeltaT;
00423 
00424   if (profile.elapsedTime > (profile.startupDuration + profile.steadyDuration + profile.endDuration)) {
00425     profile.Running = false;
00426     profile.elapsedTime = 0.0;
00427     oneMinusCosineGust.vWindTransformed.InitMatrix(0.0);
00428     vCosineGust.InitMatrix(0);
00429   }
00430 }
00431 
00432 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00433 
00434 void FGWinds::NumberOfUpDownburstCells(int num)
00435 {
00436   for (unsigned int i=0; i<UpDownBurstCells.size();i++) delete UpDownBurstCells[i];
00437   UpDownBurstCells.clear();
00438   if (num >= 0) {
00439     for (int i=0; i<num; i++) UpDownBurstCells.push_back(new struct UpDownBurst);
00440   }
00441 }
00442 
00443 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00444 // Calculates the distance between a specified point (where presumably the
00445 // Up/Downburst is centered) and the current vehicle location. The distance
00446 // here is calculated from the Haversine formula.
00447 
00448 double FGWinds::DistanceFromRingCenter(double lat, double lon)
00449 {
00450   double deltaLat = in.latitude - lat;
00451   double deltaLong = in.longitude - lon;
00452   double dLat2 = deltaLat/2.0;
00453   double dLong2 = deltaLong/2.0;
00454   double a = sin(dLat2)*sin(dLat2)
00455              + cos(lat)*cos(in.latitude)*sin(dLong2)*sin(dLong2);
00456   double c = 2.0*atan2(sqrt(a), sqrt(1.0 - a));
00457   double d = in.planetRadius*c;
00458   return d;
00459 }
00460 
00461 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00462 
00463 void FGWinds::UpDownBurst()
00464 {
00465 
00466   for (unsigned int i=0; i<UpDownBurstCells.size(); i++) {
00467     /*double d =*/ DistanceFromRingCenter(UpDownBurstCells[i]->ringLatitude, UpDownBurstCells[i]->ringLongitude);
00468 
00469   }
00470 }
00471 
00472 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00473 
00474 void FGWinds::bind(void)
00475 {
00476   typedef double (FGWinds::*PMF)(int) const;
00477   typedef int (FGWinds::*PMFt)(void) const;
00478   typedef void   (FGWinds::*PMFd)(int,double);
00479   typedef void   (FGWinds::*PMFi)(int);
00480   typedef double (FGWinds::*Ptr)(void) const;
00481 
00482   // User-specified steady, constant, wind properties (local navigational/geographic frame: N-E-D)
00483   PropertyManager->Tie("atmosphere/psiw-rad", this, &FGWinds::GetWindPsi, &FGWinds::SetWindPsi);
00484   PropertyManager->Tie("atmosphere/wind-north-fps", this, eNorth, (PMF)&FGWinds::GetWindNED,
00485                                                           (PMFd)&FGWinds::SetWindNED);
00486   PropertyManager->Tie("atmosphere/wind-east-fps",  this, eEast, (PMF)&FGWinds::GetWindNED,
00487                                                           (PMFd)&FGWinds::SetWindNED);
00488   PropertyManager->Tie("atmosphere/wind-down-fps",  this, eDown, (PMF)&FGWinds::GetWindNED,
00489                                                           (PMFd)&FGWinds::SetWindNED);
00490   PropertyManager->Tie("atmosphere/wind-mag-fps", this, &FGWinds::GetWindspeed,
00491                                                         &FGWinds::SetWindspeed);
00492 
00493   // User-specifieded gust (local navigational/geographic frame: N-E-D)
00494   PropertyManager->Tie("atmosphere/gust-north-fps", this, eNorth, (PMF)&FGWinds::GetGustNED,
00495                                                           (PMFd)&FGWinds::SetGustNED);
00496   PropertyManager->Tie("atmosphere/gust-east-fps",  this, eEast, (PMF)&FGWinds::GetGustNED,
00497                                                           (PMFd)&FGWinds::SetGustNED);
00498   PropertyManager->Tie("atmosphere/gust-down-fps",  this, eDown, (PMF)&FGWinds::GetGustNED,
00499                                                           (PMFd)&FGWinds::SetGustNED);
00500 
00501   // User-specified 1 - cosine gust parameters (in specified frame)
00502   PropertyManager->Tie("atmosphere/cosine-gust/startup-duration-sec", this, (Ptr)0L, &FGWinds::StartupGustDuration);
00503   PropertyManager->Tie("atmosphere/cosine-gust/steady-duration-sec", this, (Ptr)0L, &FGWinds::SteadyGustDuration);
00504   PropertyManager->Tie("atmosphere/cosine-gust/end-duration-sec", this, (Ptr)0L, &FGWinds::EndGustDuration);
00505   PropertyManager->Tie("atmosphere/cosine-gust/magnitude-ft_sec", this, (Ptr)0L, &FGWinds::GustMagnitude);
00506   PropertyManager->Tie("atmosphere/cosine-gust/frame", this, (PMFt)0L, (PMFi)&FGWinds::GustFrame);
00507   PropertyManager->Tie("atmosphere/cosine-gust/X-velocity-ft_sec", this, (Ptr)0L, &FGWinds::GustXComponent);
00508   PropertyManager->Tie("atmosphere/cosine-gust/Y-velocity-ft_sec", this, (Ptr)0L, &FGWinds::GustYComponent);
00509   PropertyManager->Tie("atmosphere/cosine-gust/Z-velocity-ft_sec", this, (Ptr)0L, &FGWinds::GustZComponent);
00510   PropertyManager->Tie("atmosphere/cosine-gust/start", this, (PMFt)0L, (PMFi)&FGWinds::StartGust);
00511 
00512   // User-specified Up- Down-burst parameters
00513   PropertyManager->Tie("atmosphere/updownburst/number-of-cells", this, (PMFt)0L, &FGWinds::NumberOfUpDownburstCells);
00514 //  PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
00515 //  PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
00516 //  PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
00517 //  PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
00518 //  PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
00519 //  PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
00520 //  PropertyManager->Tie("atmosphere/updownburst/", this, (Ptr)0L, &FGWinds::);
00521 
00522   // User-specified turbulence (local navigational/geographic frame: N-E-D)
00523   PropertyManager->Tie("atmosphere/turb-north-fps", this, eNorth, (PMF)&FGWinds::GetTurbNED,
00524                                                           (PMFd)&FGWinds::SetTurbNED);
00525   PropertyManager->Tie("atmosphere/turb-east-fps",  this, eEast, (PMF)&FGWinds::GetTurbNED,
00526                                                           (PMFd)&FGWinds::SetTurbNED);
00527   PropertyManager->Tie("atmosphere/turb-down-fps",  this, eDown, (PMF)&FGWinds::GetTurbNED,
00528                                                           (PMFd)&FGWinds::SetTurbNED);
00529   // Experimental turbulence parameters
00530   PropertyManager->Tie("atmosphere/p-turb-rad_sec", this,1, (PMF)&FGWinds::GetTurbPQR);
00531   PropertyManager->Tie("atmosphere/q-turb-rad_sec", this,2, (PMF)&FGWinds::GetTurbPQR);
00532   PropertyManager->Tie("atmosphere/r-turb-rad_sec", this,3, (PMF)&FGWinds::GetTurbPQR);
00533   PropertyManager->Tie("atmosphere/turb-type", this, (PMFt)&FGWinds::GetTurbType, (PMFi)&FGWinds::SetTurbType);
00534   PropertyManager->Tie("atmosphere/turb-rate", this, &FGWinds::GetTurbRate, &FGWinds::SetTurbRate);
00535   PropertyManager->Tie("atmosphere/turb-gain", this, &FGWinds::GetTurbGain, &FGWinds::SetTurbGain);
00536   PropertyManager->Tie("atmosphere/turb-rhythmicity", this, &FGWinds::GetRhythmicity,
00537                                                             &FGWinds::SetRhythmicity);
00538 
00539   // Parameters for milspec turbulence
00540   PropertyManager->Tie("atmosphere/turbulence/milspec/windspeed_at_20ft_AGL-fps",
00541                        this, &FGWinds::GetWindspeed20ft,
00542                              &FGWinds::SetWindspeed20ft);
00543   PropertyManager->Tie("atmosphere/turbulence/milspec/severity",
00544                        this, &FGWinds::GetProbabilityOfExceedence,
00545                              &FGWinds::SetProbabilityOfExceedence);
00546 
00547   // Total, calculated winds (local navigational/geographic frame: N-E-D). Read only.
00548   PropertyManager->Tie("atmosphere/total-wind-north-fps", this, eNorth, (PMF)&FGWinds::GetTotalWindNED);
00549   PropertyManager->Tie("atmosphere/total-wind-east-fps",  this, eEast,  (PMF)&FGWinds::GetTotalWindNED);
00550   PropertyManager->Tie("atmosphere/total-wind-down-fps",  this, eDown,  (PMF)&FGWinds::GetTotalWindNED);
00551 
00552 }
00553 
00554 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00555 //    The bitmasked value choices are as follows:
00556 //    unset: In this case (the default) JSBSim would only print
00557 //       out the normally expected messages, essentially echoing
00558 //       the config files as they are read. If the environment
00559 //       variable is not set, debug_lvl is set to 1 internally
00560 //    0: This requests JSBSim not to output any messages
00561 //       whatsoever.
00562 //    1: This value explicity requests the normal JSBSim
00563 //       startup messages
00564 //    2: This value asks for a message to be printed out when
00565 //       a class is instantiated
00566 //    4: When this value is set, a message is displayed when a
00567 //       FGModel object executes its Run() method
00568 //    8: When this value is set, various runtime state variables
00569 //       are printed out periodically
00570 //    16: When set various parameters are sanity checked and
00571 //       a message is printed out when they go out of bounds
00572 
00573 void FGWinds::Debug(int from)
00574 {
00575   if (debug_lvl <= 0) return;
00576 
00577   if (debug_lvl & 1) { // Standard console startup message output
00578     if (from == 0) { // Constructor
00579     }
00580   }
00581   if (debug_lvl & 2 ) { // Instantiation/Destruction notification
00582     if (from == 0) cout << "Instantiated: FGWinds" << endl;
00583     if (from == 1) cout << "Destroyed:    FGWinds" << endl;
00584   }
00585   if (debug_lvl & 4 ) { // Run() method entry print for FGModel-derived objects
00586   }
00587   if (debug_lvl & 8 ) { // Runtime state variables
00588   }
00589   if (debug_lvl & 16) { // Sanity checking
00590   }
00591   if (debug_lvl & 128) { //
00592   }
00593   if (debug_lvl & 64) {
00594     if (from == 0) { // Constructor
00595       cout << IdSrc << endl;
00596       cout << IdHdr << endl;
00597     }
00598   }
00599 }
00600 
00601 } // namespace JSBSim