48 #include "FGFDMExec.h" 54 IDENT(IdSrc,
"$Id: FGWinds.cpp,v 1.15 2015/02/27 20:49:36 bcoconni Exp $");
55 IDENT(IdHdr,ID_WINDS);
75 static inline double sqr(
double x) {
return x*x; }
81 MagnitudedAccelDt = MagnitudeAccel = Magnitude = TurbDirection = 0.0;
86 spike = target_time = strength = 0.0;
87 wind_from_clockwise = 0.0;
90 vGustNED.InitMatrix();
91 vTurbulenceNED.InitMatrix();
92 vCosineGust.InitMatrix();
95 windspeed_at_20ft = 0.;
96 probability_of_exceedence_index = 0;
101 << 500.0 << 1750.0 << 3750.0 << 7500.0 << 15000.0 << 25000.0 << 35000.0 << 45000.0 << 55000.0 << 65000.0 << 75000.0 << 80000.0
102 << 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
103 << 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
104 << 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
105 << 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
106 << 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
107 << 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
108 << 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;
124 bool FGWinds::InitModel(
void)
126 if (!FGModel::InitModel())
return false;
130 vGustNED.InitMatrix();
131 vTurbulenceNED.InitMatrix();
132 vCosineGust.InitMatrix();
145 if (Holding)
return false;
147 if (turbType != ttNone) Turbulence(in.AltitudeASL);
150 vTotalWindNED = vWindNED + vGustNED + vCosineGust + vTurbulenceNED;
153 if (vWindNED(eX) != 0.0) psiw = atan2( vWindNED(eY), vWindNED(eX) );
154 if (psiw < 0) psiw += 2*M_PI;
164 void FGWinds::SetWindspeed(
double speed)
168 vWindNED(eNorth) = speed;
170 vWindNED(eNorth) = speed * cos(psiw);
171 vWindNED(eEast) = speed * sin(psiw);
172 vWindNED(eDown) = 0.0;
178 double FGWinds::GetWindspeed(
void)
const 189 double mag = GetWindspeed();
196 void FGWinds::Turbulence(
double h)
202 vTurbPQR(eP) = wind_from_clockwise;
203 if (TurbGain == 0.0)
return;
206 if (TurbGain < 0.0) TurbGain = 0.0;
207 if (TurbGain > 1.0) TurbGain = 1.0;
208 if (TurbRate < 0.0) TurbRate = 0.0;
209 if (TurbRate > 30.0) TurbRate = 30.0;
210 if (Rhythmicity < 0.0) Rhythmicity = 0.0;
211 if (Rhythmicity > 1.0) Rhythmicity = 1.0;
215 double sinewave = sin( time * TurbRate * 6.283185307 );
218 if (target_time == 0.0) {
219 strength = random = 1 - 2.0*(double(rand())/double(RAND_MAX));
220 target_time = time + 0.71 + (random * 0.5);
222 if (time > target_time) {
230 vTurbulenceNED.InitMatrix();
231 double delta = strength * max_vs * TurbGain * (1-Rhythmicity) * spike;
234 vTurbulenceNED(eDown) = sinewave * max_vs * TurbGain * Rhythmicity;
235 vTurbulenceNED(eDown)+= delta;
236 if (in.DistanceAGL/in.wingspan < 3.0)
237 vTurbulenceNED(eDown) *= in.DistanceAGL/in.wingspan * 0.3333;
240 vTurbulenceNED(eNorth) = sin( delta * 3.0 );
241 vTurbulenceNED(eEast) = cos( delta * 3.0 );
244 vTurbPQR(eP) += delta * 0.04;
254 if (probability_of_exceedence_index == 0 || in.V == 0) {
255 vTurbulenceNED(eNorth) = vTurbulenceNED(eEast) = vTurbulenceNED(eDown) = 0.0;
256 vTurbPQR(eP) = vTurbPQR(eQ) = vTurbPQR(eR) = 0.0;
261 double b_w = in.wingspan, L_u, L_w, sig_u, sig_w;
263 if (b_w == 0.) b_w = 30.;
266 if (h <= 10.) h = 10;
270 L_u = h/pow(0.177 + 0.000823*h, 1.2);
272 sig_w = 0.1*windspeed_at_20ft;
273 sig_u = sig_w/pow(0.177 + 0.000823*h, 0.4);
274 }
else if (h <= 2000) {
276 L_u = L_w = 1000 + (h-1000.)/1000.*750.;
277 sig_u = sig_w = 0.1*windspeed_at_20ft
278 + (h-1000.)/1000.*(POE_Table->GetValue(probability_of_exceedence_index, h) - 0.1*windspeed_at_20ft);
281 sig_u = sig_w = POE_Table->GetValue(probability_of_exceedence_index, h);
287 xi_u_km1 = 0, nu_u_km1 = 0,
288 xi_v_km1 = 0, xi_v_km2 = 0, nu_v_km1 = 0, nu_v_km2 = 0,
289 xi_w_km1 = 0, xi_w_km2 = 0, nu_w_km1 = 0, nu_w_km2 = 0,
290 xi_p_km1 = 0, nu_p_km1 = 0,
291 xi_q_km1 = 0, xi_r_km1 = 0;
295 T_V = in.totalDeltaT,
296 sig_p = 1.9/sqrt(L_w*b_w)*sig_w,
299 L_p = sqrt(L_w*b_w)/2.6,
303 tau_q = 4*b_w/M_PI/in.V,
304 tau_r =3*b_w/M_PI/in.V,
305 nu_u = GaussianRandomNumber(),
306 nu_v = GaussianRandomNumber(),
307 nu_w = GaussianRandomNumber(),
308 nu_p = GaussianRandomNumber(),
309 xi_u=0, xi_v=0, xi_w=0, xi_p=0, xi_q=0, xi_r=0;
313 if (turbType == ttTustin) {
318 C_BL = 1/tau_u/tan(T_V/2/tau_u),
319 C_BLp = 1/tau_p/tan(T_V/2/tau_p),
320 C_BLq = 1/tau_q/tan(T_V/2/tau_q),
321 C_BLr = 1/tau_r/tan(T_V/2/tau_r);
327 xi_u = -(1 - C_BL*tau_u)/(1 + C_BL*tau_u)*xi_u_km1
328 + sig_u*sqrt(2*tau_u/T_V)/(1 + C_BL*tau_u)*(nu_u + nu_u_km1);
329 xi_v = -2*(sqr(omega_v) - sqr(C_BL))/sqr(omega_v + C_BL)*xi_v_km1
330 - sqr(omega_v - C_BL)/sqr(omega_v + C_BL) * xi_v_km2
331 + sig_u*sqrt(3*omega_v/T_V)/sqr(omega_v + C_BL)*(
332 (C_BL + omega_v/sqrt(3.))*nu_v
333 + 2/sqrt(3.)*omega_v*nu_v_km1
334 + (omega_v/sqrt(3.) - C_BL)*nu_v_km2);
335 xi_w = -2*(sqr(omega_w) - sqr(C_BL))/sqr(omega_w + C_BL)*xi_w_km1
336 - sqr(omega_w - C_BL)/sqr(omega_w + C_BL) * xi_w_km2
337 + sig_w*sqrt(3*omega_w/T_V)/sqr(omega_w + C_BL)*(
338 (C_BL + omega_w/sqrt(3.))*nu_w
339 + 2/sqrt(3.)*omega_w*nu_w_km1
340 + (omega_w/sqrt(3.) - C_BL)*nu_w_km2);
341 xi_p = -(1 - C_BLp*tau_p)/(1 + C_BLp*tau_p)*xi_p_km1
342 + sig_p*sqrt(2*tau_p/T_V)/(1 + C_BLp*tau_p) * (nu_p + nu_p_km1);
343 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
344 + C_BLq/in.V/(1 + 4*b_w*C_BLq/M_PI/in.V) * (xi_w - xi_w_km1);
345 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
346 + C_BLr/in.V/(1 + 3*b_w*C_BLr/M_PI/in.V) * (xi_v - xi_v_km1);
348 }
else if (turbType == ttMilspec) {
351 xi_u = (1 - T_V/tau_u) *xi_u_km1 + sig_u*sqrt(2*T_V/tau_u)*nu_u;
352 xi_v = (1 - 2*T_V/tau_u)*xi_v_km1 + sig_u*sqrt(4*T_V/tau_u)*nu_v;
353 xi_w = (1 - 2*T_V/tau_w)*xi_w_km1 + sig_w*sqrt(4*T_V/tau_w)*nu_w;
354 xi_p = (1 - T_V/tau_p) *xi_p_km1 + sig_p*sqrt(2*T_V/tau_p)*nu_p;
355 xi_q = (1 - T_V/tau_q) *xi_q_km1 + M_PI/4/b_w*(xi_w - xi_w_km1);
356 xi_r = (1 - T_V/tau_r) *xi_r_km1 + M_PI/3/b_w*(xi_v - xi_v_km1);
360 double cospsi = cos(psiw), sinpsi = sin(psiw);
361 vTurbulenceNED(eNorth) = cospsi*xi_u + sinpsi*xi_v;
362 vTurbulenceNED(eEast) = -sinpsi*xi_u + cospsi*xi_v;
363 vTurbulenceNED(eDown) = xi_w;
365 vTurbPQR(eP) = cospsi*xi_p + sinpsi*xi_q;
366 vTurbPQR(eQ) = -sinpsi*xi_p + cospsi*xi_q;
370 vTurbPQR = in.Tl2b*vTurbPQR;
373 xi_u_km1 = xi_u; nu_u_km1 = nu_u;
374 xi_v_km2 = xi_v_km1; xi_v_km1 = xi_v; nu_v_km2 = nu_v_km1; nu_v_km1 = nu_v;
375 xi_w_km2 = xi_w_km1; xi_w_km1 = xi_w; nu_w_km2 = nu_w_km1; nu_w_km1 = nu_w;
376 xi_p_km1 = xi_p; nu_p_km1 = nu_p;
385 TurbDirection = atan2( vTurbulenceNED(eEast), vTurbulenceNED(eNorth))*radtodeg;
391 double FGWinds::CosineGustProfile(
double startDuration,
double steadyDuration,
double endDuration,
double elapsedTime)
394 if (elapsedTime >= 0 && elapsedTime <= startDuration) {
395 factor = (1.0 - cos(M_PI*elapsedTime/startDuration))/2.0;
396 }
else if (elapsedTime > startDuration && (elapsedTime <= (startDuration + steadyDuration))) {
398 }
else if (elapsedTime > (startDuration + steadyDuration) && elapsedTime <= (startDuration + steadyDuration + endDuration)) {
399 factor = (1-cos(M_PI*(1-(elapsedTime-(startDuration + steadyDuration))/endDuration)))/2.0;
409 void FGWinds::CosineGust()
445 vCosineGust.InitMatrix(0);
451 void FGWinds::NumberOfUpDownburstCells(
int num)
453 for (
unsigned int i=0; i<UpDownBurstCells.size();i++)
delete UpDownBurstCells[i];
454 UpDownBurstCells.clear();
456 for (
int i=0; i<num; i++) UpDownBurstCells.push_back(
new struct UpDownBurst);
465 double FGWinds::DistanceFromRingCenter(
double lat,
double lon)
467 double deltaLat = in.latitude - lat;
468 double deltaLong = in.longitude - lon;
469 double dLat2 = deltaLat/2.0;
470 double dLong2 = deltaLong/2.0;
471 double a = sin(dLat2)*sin(dLat2)
472 + cos(lat)*cos(in.latitude)*sin(dLong2)*sin(dLong2);
473 double c = 2.0*atan2(sqrt(a), sqrt(1.0 - a));
474 double d = in.planetRadius*c;
483 for (
unsigned int i=0; i<UpDownBurstCells.size(); i++) {
484 DistanceFromRingCenter(UpDownBurstCells[i]->ringLatitude, UpDownBurstCells[i]->ringLongitude);
491 void FGWinds::bind(
void)
493 typedef double (
FGWinds::*PMF)(int)
const;
494 typedef int (
FGWinds::*PMFt)(void)
const;
495 typedef void (
FGWinds::*PMFd)(int,double);
496 typedef void (
FGWinds::*PMFi)(int);
497 typedef double (
FGWinds::*Ptr)(void)
const;
503 PropertyManager->
Tie(
"atmosphere/wind-east-fps",
this, eEast, (PMF)&FGWinds::GetWindNED,
504 (PMFd)&FGWinds::SetWindNED);
505 PropertyManager->
Tie(
"atmosphere/wind-down-fps",
this, eDown, (PMF)&FGWinds::GetWindNED,
506 (PMFd)&FGWinds::SetWindNED);
507 PropertyManager->
Tie(
"atmosphere/wind-mag-fps",
this, &FGWinds::GetWindspeed,
508 &FGWinds::SetWindspeed);
513 PropertyManager->
Tie(
"atmosphere/gust-east-fps",
this, eEast, (PMF)&FGWinds::GetGustNED,
514 (PMFd)&FGWinds::SetGustNED);
515 PropertyManager->
Tie(
"atmosphere/gust-down-fps",
this, eDown, (PMF)&FGWinds::GetGustNED,
516 (PMFd)&FGWinds::SetGustNED);
530 PropertyManager->
Tie(
"atmosphere/updownburst/number-of-cells",
this, (PMFt)0L, &FGWinds::NumberOfUpDownburstCells);
542 PropertyManager->
Tie(
"atmosphere/turb-east-fps",
this, eEast, (PMF)&FGWinds::GetTurbNED,
543 (PMFd)&FGWinds::SetTurbNED);
544 PropertyManager->
Tie(
"atmosphere/turb-down-fps",
this, eDown, (PMF)&FGWinds::GetTurbNED,
545 (PMFd)&FGWinds::SetTurbNED);
547 PropertyManager->
Tie(
"atmosphere/p-turb-rad_sec",
this,1, (PMF)&FGWinds::GetTurbPQR);
548 PropertyManager->
Tie(
"atmosphere/q-turb-rad_sec",
this,2, (PMF)&FGWinds::GetTurbPQR);
549 PropertyManager->
Tie(
"atmosphere/r-turb-rad_sec",
this,3, (PMF)&FGWinds::GetTurbPQR);
551 PropertyManager->
Tie(
"atmosphere/turb-rate",
this, &FGWinds::GetTurbRate, &FGWinds::SetTurbRate);
552 PropertyManager->
Tie(
"atmosphere/turb-gain",
this, &FGWinds::GetTurbGain, &FGWinds::SetTurbGain);
553 PropertyManager->
Tie(
"atmosphere/turb-rhythmicity",
this, &FGWinds::GetRhythmicity,
554 &FGWinds::SetRhythmicity);
557 PropertyManager->
Tie(
"atmosphere/turbulence/milspec/windspeed_at_20ft_AGL-fps",
558 this, &FGWinds::GetWindspeed20ft,
559 &FGWinds::SetWindspeed20ft);
560 PropertyManager->
Tie(
"atmosphere/turbulence/milspec/severity",
561 this, &FGWinds::GetProbabilityOfExceedence,
566 PropertyManager->
Tie(
"atmosphere/total-wind-east-fps",
this, eEast, (PMF)&FGWinds::GetTotalWindNED);
567 PropertyManager->
Tie(
"atmosphere/total-wind-down-fps",
this, eDown, (PMF)&FGWinds::GetTotalWindNED);
590 void FGWinds::Debug(
int from)
592 if (debug_lvl <= 0)
return;
598 if (debug_lvl & 2 ) {
599 if (from == 0) cout <<
"Instantiated: FGWinds" << endl;
600 if (from == 1) cout <<
"Destroyed: FGWinds" << endl;
602 if (debug_lvl & 4 ) {
604 if (debug_lvl & 8 ) {
606 if (debug_lvl & 16) {
608 if (debug_lvl & 128) {
610 if (debug_lvl & 64) {
612 cout << IdSrc << endl;
613 cout << IdHdr << endl;
virtual void SteadyGustDuration(double dur)
Specifies the length of time that the gust is at a steady, full strength.
virtual void SetWindPsi(double dir)
Sets the direction that the wind is coming from.
virtual void StartGust(bool running)
Initiates the execution of the gust.
virtual void SetGustNED(int idx, double gust)
Sets a gust component in NED frame.
virtual void StartupGustDuration(double dur)
Specifies the duration of the startup portion of the gust.
virtual const FGColumnVector3 & GetTotalWindNED(void) const
Retrieves the total wind components in NED frame.
FGMatrix33 Inverse(void) const
Return the inverse of the matrix.
virtual void SetWindNED(double wN, double wE, double wD)
Sets the wind components in NED frame.
virtual void SetTurbNED(int idx, double turb)
Sets a turbulence component in NED frame.
virtual void EndGustDuration(double dur)
Specifies the length of time it takes for the gust to return to zero velocity.
FGColumnVector3 & Normalize(void)
Normalize.
virtual void GustYComponent(double y)
Specifies the Y component of velocity in the specified gust frame (ft/sec).
FGColumnVector3 vWindTransformed
Models atmospheric disturbances: winds, gusts, turbulence, downbursts, etc.
virtual bool Run(bool Holding)
Runs the model; called by the Executive.
void Tie(const std::string &name, bool *pointer, bool useDefault=true)
Tie a property to an external bool variable.
virtual const FGColumnVector3 & GetGustNED(void) const
Retrieves the gust components in NED frame.
Stores information about a specified Up- or Down-burst.
struct OneMinusCosineProfile gustProfile
virtual double GetTurbNED(int idx) const
Retrieves a turbulence component in NED frame.
virtual double GetWindPsi(void) const
Retrieves the direction that the wind is coming from.
Base class for all scheduled JSBSim models.
virtual void GustMagnitude(double mag)
Specifies the magnitude of the gust in feet/second.
virtual void SetTurbType(tType tt)
Turbulence models available: ttNone, ttStandard, ttBerndt, ttCulp, ttMilspec, ttTustin.
virtual const FGColumnVector3 & GetWindNED(void) const
Retrieves the wind components in NED frame.
double GetSimTime(void) const
Returns the cumulative simulation time in seconds.
virtual void GustFrame(eGustFrame gFrame)
Specifies the frame that the gust direction vector components are specified in.
virtual void GustXComponent(double x)
Specifies the X component of velocity in the specified gust frame (ft/sec).
virtual void GustZComponent(double z)
Specifies the Z component of velocity in the specified gust frame (ft/sec).
virtual void SetProbabilityOfExceedence(int idx)
allowable range: 0-7, 3=light, 4=moderate, 6=severe turbulence
double Magnitude(void) const
Length of the vector.
Encapsulates the JSBSim simulation executive.
bool Run(bool Holding)
Runs the winds model; called by the Executive Can pass in a value indicating if the executive is dire...