48#include "math/FGTable.h"
49#include "input_output/FGLog.h"
73constexpr double sqr(
double x) {
return x*x; }
76 :
FGModel(fdmex), generator(fdmex->GetRandomGenerator())
80 MagnitudedAccelDt = MagnitudeAccel = Magnitude = TurbDirection = 0.0;
85 spike = target_time = strength = 0.0;
86 wind_from_clockwise = 0.0;
89 vGustNED.InitMatrix();
90 vTurbulenceNED.InitMatrix();
91 vCosineGust.InitMatrix();
94 windspeed_at_20ft = 0.;
95 probability_of_exceedence_index = 0;
100 << 500.0 << 1750.0 << 3750.0 << 7500.0 << 15000.0 << 25000.0 << 35000.0 << 45000.0 << 55000.0 << 65000.0 << 75000.0 << 80000.0
101 << 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
102 << 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
103 << 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
104 << 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
105 << 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
106 << 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
107 << 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;
123bool FGWinds::InitModel(
void)
125 if (!FGModel::InitModel())
return false;
129 vGustNED.InitMatrix();
130 vTurbulenceNED.InitMatrix();
131 vCosineGust.InitMatrix();
136 xi_u_km1 = nu_u_km1 = 0;
137 xi_v_km1 = xi_v_km2 = nu_v_km1 = nu_v_km2 = 0;
138 xi_w_km1 = xi_w_km2 = nu_w_km1 = nu_w_km2 = 0;
139 xi_p_km1 = nu_p_km1 = 0;
140 xi_q_km1 = xi_r_km1 = 0;
150 if (Holding)
return false;
152 if (turbType != ttNone)
153 Turbulence(in.AltitudeASL);
155 vTurbulenceNED.InitMatrix();
158 vTotalWindNED = vWindNED + vGustNED + vCosineGust + vTurbulenceNED;
161 if (vWindNED(eX) != 0.0) psiw = atan2( vWindNED(eY), vWindNED(eX) );
162 if (psiw < 0) psiw += 2*M_PI;
172void FGWinds::SetWindspeed(
double speed)
176 vWindNED(eNorth) = speed;
178 vWindNED(eNorth) = speed * cos(psiw);
179 vWindNED(eEast) = speed * sin(psiw);
180 vWindNED(eDown) = 0.0;
186double FGWinds::GetWindspeed(
void)
const
197 double mag = GetWindspeed();
204void FGWinds::Turbulence(
double h)
210 vTurbPQR(eP) = wind_from_clockwise;
211 if (TurbGain == 0.0)
return;
214 if (TurbGain < 0.0) TurbGain = 0.0;
215 if (TurbGain > 1.0) TurbGain = 1.0;
216 if (TurbRate < 0.0) TurbRate = 0.0;
217 if (TurbRate > 30.0) TurbRate = 30.0;
218 if (Rhythmicity < 0.0) Rhythmicity = 0.0;
219 if (Rhythmicity > 1.0) Rhythmicity = 1.0;
223 double sinewave = sin( time * TurbRate * 6.283185307 );
226 if (target_time == 0.0) {
227 strength = random = generator->GetUniformRandomNumber();
228 target_time = time + 0.71 + (random * 0.5);
230 if (time > target_time) {
238 vTurbulenceNED.InitMatrix();
239 double delta = strength * max_vs * TurbGain * (1-Rhythmicity) * spike;
242 vTurbulenceNED(eDown) = sinewave * max_vs * TurbGain * Rhythmicity;
243 vTurbulenceNED(eDown)+= delta;
244 if (in.DistanceAGL/in.wingspan < 3.0)
245 vTurbulenceNED(eDown) *= in.DistanceAGL/in.wingspan * 0.3333;
248 vTurbulenceNED(eNorth) = sin( delta * 3.0 );
249 vTurbulenceNED(eEast) = cos( delta * 3.0 );
252 vTurbPQR(eP) += delta * 0.04;
262 if (probability_of_exceedence_index == 0 || in.V == 0) {
263 vTurbulenceNED(eNorth) = vTurbulenceNED(eEast) = vTurbulenceNED(eDown) = 0.0;
264 vTurbPQR(eP) = vTurbPQR(eQ) = vTurbPQR(eR) = 0.0;
269 double b_w = in.wingspan, L_u, L_w, sig_u, sig_w;
271 if (b_w == 0.) b_w = 30.;
274 if (h <= 10.) h = 10;
278 L_u = h/pow(0.177 + 0.000823*h, 1.2);
280 sig_w = 0.1*windspeed_at_20ft;
281 sig_u = sig_w/pow(0.177 + 0.000823*h, 0.4);
282 }
else if (h <= 2000) {
284 L_u = L_w = 1000 + (h-1000.)/1000.*750.;
285 sig_u = sig_w = 0.1*windspeed_at_20ft
286 + (h-1000.)/1000.*(POE_Table->
GetValue(probability_of_exceedence_index, h) - 0.1*windspeed_at_20ft);
289 sig_u = sig_w = POE_Table->
GetValue(probability_of_exceedence_index, h);
293 T_V = in.totalDeltaT,
294 sig_p = 1.9/sqrt(L_w*b_w)*sig_w,
297 L_p = sqrt(L_w*b_w)/2.6,
301 tau_q = 4*b_w/M_PI/in.V,
302 tau_r =3*b_w/M_PI/in.V,
303 nu_u = generator->GetNormalRandomNumber(),
304 nu_v = generator->GetNormalRandomNumber(),
305 nu_w = generator->GetNormalRandomNumber(),
306 nu_p = generator->GetNormalRandomNumber(),
307 xi_u=0, xi_v=0, xi_w=0, xi_p=0, xi_q=0, xi_r=0;
311 if (turbType == ttTustin) {
316 C_BL = 1/tau_u/tan(T_V/2/tau_u),
317 C_BLp = 1/tau_p/tan(T_V/2/tau_p),
318 C_BLq = 1/tau_q/tan(T_V/2/tau_q),
319 C_BLr = 1/tau_r/tan(T_V/2/tau_r);
325 xi_u = -(1 - C_BL*tau_u)/(1 + C_BL*tau_u)*xi_u_km1
326 + sig_u*sqrt(2*tau_u/T_V)/(1 + C_BL*tau_u)*(nu_u + nu_u_km1);
327 xi_v = -2*(
sqr(omega_v) -
sqr(C_BL))/
sqr(omega_v + C_BL)*xi_v_km1
328 -
sqr(omega_v - C_BL)/
sqr(omega_v + C_BL) * xi_v_km2
329 + sig_u*sqrt(3*omega_v/T_V)/
sqr(omega_v + C_BL)*(
330 (C_BL + omega_v/sqrt(3.))*nu_v
331 + 2/sqrt(3.)*omega_v*nu_v_km1
332 + (omega_v/sqrt(3.) - C_BL)*nu_v_km2);
333 xi_w = -2*(
sqr(omega_w) -
sqr(C_BL))/
sqr(omega_w + C_BL)*xi_w_km1
334 -
sqr(omega_w - C_BL)/
sqr(omega_w + C_BL) * xi_w_km2
335 + sig_w*sqrt(3*omega_w/T_V)/
sqr(omega_w + C_BL)*(
336 (C_BL + omega_w/sqrt(3.))*nu_w
337 + 2/sqrt(3.)*omega_w*nu_w_km1
338 + (omega_w/sqrt(3.) - C_BL)*nu_w_km2);
339 xi_p = -(1 - C_BLp*tau_p)/(1 + C_BLp*tau_p)*xi_p_km1
340 + sig_p*sqrt(2*tau_p/T_V)/(1 + C_BLp*tau_p) * (nu_p + nu_p_km1);
341 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
342 + C_BLq/in.V/(1 + 4*b_w*C_BLq/M_PI/in.V) * (xi_w - xi_w_km1);
343 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
344 + C_BLr/in.V/(1 + 3*b_w*C_BLr/M_PI/in.V) * (xi_v - xi_v_km1);
346 }
else if (turbType == ttMilspec) {
349 xi_u = (1 - T_V/tau_u) *xi_u_km1 + sig_u*sqrt(2*T_V/tau_u)*nu_u;
350 xi_v = (1 - 2*T_V/tau_u)*xi_v_km1 + sig_u*sqrt(4*T_V/tau_u)*nu_v;
351 xi_w = (1 - 2*T_V/tau_w)*xi_w_km1 + sig_w*sqrt(4*T_V/tau_w)*nu_w;
352 xi_p = (1 - T_V/tau_p) *xi_p_km1 + sig_p*sqrt(2*T_V/tau_p)*nu_p;
353 xi_q = (1 - T_V/tau_q) *xi_q_km1 + M_PI/4/b_w*(xi_w - xi_w_km1);
354 xi_r = (1 - T_V/tau_r) *xi_r_km1 + M_PI/3/b_w*(xi_v - xi_v_km1);
358 double cospsi = cos(psiw), sinpsi = sin(psiw);
359 vTurbulenceNED(eNorth) = cospsi*xi_u + sinpsi*xi_v;
360 vTurbulenceNED(eEast) = -sinpsi*xi_u + cospsi*xi_v;
361 vTurbulenceNED(eDown) = xi_w;
363 vTurbPQR(eP) = cospsi*xi_p + sinpsi*xi_q;
364 vTurbPQR(eQ) = -sinpsi*xi_p + cospsi*xi_q;
368 vTurbPQR = in.Tl2b*vTurbPQR;
371 xi_u_km1 = xi_u; nu_u_km1 = nu_u;
372 xi_v_km2 = xi_v_km1; xi_v_km1 = xi_v; nu_v_km2 = nu_v_km1; nu_v_km1 = nu_v;
373 xi_w_km2 = xi_w_km1; xi_w_km1 = xi_w; nu_w_km2 = nu_w_km1; nu_w_km1 = nu_w;
374 xi_p_km1 = xi_p; nu_p_km1 = nu_p;
383 TurbDirection = atan2( vTurbulenceNED(eEast), vTurbulenceNED(eNorth))*radtodeg;
389double FGWinds::CosineGustProfile(
double startDuration,
double steadyDuration,
double endDuration,
double elapsedTime)
392 if (elapsedTime >= 0 && elapsedTime <= startDuration) {
393 factor = (1.0 - cos(M_PI*elapsedTime/startDuration))/2.0;
394 }
else if (elapsedTime > startDuration && (elapsedTime <= (startDuration + steadyDuration))) {
396 }
else if (elapsedTime > (startDuration + steadyDuration) && elapsedTime <= (startDuration + steadyDuration + endDuration)) {
397 factor = (1-cos(M_PI*(1-(elapsedTime-(startDuration + steadyDuration))/endDuration)))/2.0;
407void FGWinds::CosineGust()
409 struct OneMinusCosineProfile& profile = oneMinusCosineGust.
gustProfile;
411 double factor = CosineGustProfile( profile.startupDuration,
412 profile.steadyDuration,
414 profile.elapsedTime);
437 profile.elapsedTime += in.totalDeltaT;
439 if (profile.elapsedTime > (profile.startupDuration + profile.steadyDuration + profile.endDuration)) {
440 profile.Running =
false;
441 profile.elapsedTime = 0.0;
443 vCosineGust.InitMatrix(0);
449void FGWinds::NumberOfUpDownburstCells(
int num)
451 for (
unsigned int i=0; i<UpDownBurstCells.size();i++)
delete UpDownBurstCells[i];
452 UpDownBurstCells.clear();
454 for (
int i=0; i<num; i++) UpDownBurstCells.push_back(
new struct UpDownBurst);
463double FGWinds::DistanceFromRingCenter(
double lat,
double lon)
465 double deltaLat = in.latitude - lat;
466 double deltaLong = in.longitude - lon;
467 double dLat2 = deltaLat/2.0;
468 double dLong2 = deltaLong/2.0;
469 double a = sin(dLat2)*sin(dLat2)
470 + cos(lat)*cos(in.latitude)*sin(dLong2)*sin(dLong2);
471 double c = 2.0*atan2(sqrt(a), sqrt(1.0 - a));
472 double d = in.planetRadius*c;
478void FGWinds::UpDownBurst()
481 for (
unsigned int i=0; i<UpDownBurstCells.size(); i++) {
482 DistanceFromRingCenter(UpDownBurstCells[i]->ringLatitude, UpDownBurstCells[i]->ringLongitude);
492void FGWinds::SetRandomSeed(
int sr)
495 generator = std::make_shared<RandomNumberGenerator>(*RandomSeed);
498int FGWinds::GetRandomSeed(
void)
const {
502 return FDMExec->SRand();
507void FGWinds::bind(
void)
511 PropertyManager->Tie(
"atmosphere/wind-north-fps",
this, eNorth,
513 PropertyManager->Tie(
"atmosphere/wind-east-fps",
this, eEast,
515 PropertyManager->Tie(
"atmosphere/wind-down-fps",
this, eDown,
517 PropertyManager->Tie(
"atmosphere/wind-mag-fps",
this, &FGWinds::GetWindspeed,
518 &FGWinds::SetWindspeed);
521 PropertyManager->Tie(
"atmosphere/gust-north-fps",
this, eNorth,
523 PropertyManager->Tie(
"atmosphere/gust-east-fps",
this, eEast,
525 PropertyManager->Tie(
"atmosphere/gust-down-fps",
this, eDown,
529 PropertyManager->Tie<
FGWinds,
double>(
"atmosphere/cosine-gust/startup-duration-sec",
531 PropertyManager->Tie<
FGWinds,
double>(
"atmosphere/cosine-gust/steady-duration-sec",
533 PropertyManager->Tie<
FGWinds,
double>(
"atmosphere/cosine-gust/end-duration-sec",
535 PropertyManager->Tie<
FGWinds,
double>(
"atmosphere/cosine-gust/magnitude-ft_sec",
537 PropertyManager->Tie<
FGWinds, eGustFrame>(
"atmosphere/cosine-gust/frame",
this,
539 PropertyManager->Tie<
FGWinds,
double>(
"atmosphere/cosine-gust/X-velocity-ft_sec",
541 PropertyManager->Tie<
FGWinds,
double>(
"atmosphere/cosine-gust/Y-velocity-ft_sec",
543 PropertyManager->Tie<
FGWinds,
double>(
"atmosphere/cosine-gust/Z-velocity-ft_sec",
545 PropertyManager->Tie<
FGWinds,
bool>(
"atmosphere/cosine-gust/start",
this,
549 PropertyManager->Tie<
FGWinds,
int>(
"atmosphere/updownburst/number-of-cells",
550 this,
nullptr, &FGWinds::NumberOfUpDownburstCells);
560 PropertyManager->Tie(
"atmosphere/turb-north-fps",
this, eNorth,
562 PropertyManager->Tie(
"atmosphere/turb-east-fps",
this, eEast,
564 PropertyManager->Tie(
"atmosphere/turb-down-fps",
this, eDown,
567 PropertyManager->Tie(
"atmosphere/p-turb-rad_sec",
this, eP, &FGWinds::GetTurbPQR);
568 PropertyManager->Tie(
"atmosphere/q-turb-rad_sec",
this, eQ, &FGWinds::GetTurbPQR);
569 PropertyManager->Tie(
"atmosphere/r-turb-rad_sec",
this, eR, &FGWinds::GetTurbPQR);
570 PropertyManager->Tie(
"atmosphere/turb-type",
this, &FGWinds::GetTurbType, &
FGWinds::SetTurbType);
571 PropertyManager->Tie(
"atmosphere/turb-rate",
this, &FGWinds::GetTurbRate, &FGWinds::SetTurbRate);
572 PropertyManager->Tie(
"atmosphere/turb-gain",
this, &FGWinds::GetTurbGain, &FGWinds::SetTurbGain);
573 PropertyManager->Tie(
"atmosphere/turb-rhythmicity",
this, &FGWinds::GetRhythmicity,
574 &FGWinds::SetRhythmicity);
577 PropertyManager->Tie(
"atmosphere/turbulence/milspec/windspeed_at_20ft_AGL-fps",
578 this, &FGWinds::GetWindspeed20ft,
579 &FGWinds::SetWindspeed20ft);
580 PropertyManager->Tie(
"atmosphere/turbulence/milspec/severity",
581 this, &FGWinds::GetProbabilityOfExceedence,
590 PropertyManager->Tie(
"atmosphere/randomseed",
this, &FGWinds::GetRandomSeed, &FGWinds::SetRandomSeed);
612void FGWinds::Debug(
int from)
614 if (debug_lvl <= 0)
return;
620 if (debug_lvl & 2 ) {
621 FGLogging log(LogLevel::DEBUG);
622 if (from == 0) log <<
"Instantiated: FGWinds" << endl;
623 if (from == 1) log <<
"Destroyed: FGWinds" << endl;
625 if (debug_lvl & 4 ) {
627 if (debug_lvl & 8 ) {
629 if (debug_lvl & 16) {
631 if (debug_lvl & 128) {
633 if (debug_lvl & 64) {
double Magnitude(void) const
Length of the vector.
FGColumnVector3 & Normalize(void)
Normalize.
Encapsulates the JSBSim simulation executive.
double GetSimTime(void) const
Returns the cumulative simulation time in seconds.
FGMatrix33 Inverse(void) const
Return the inverse of the matrix.
Base class for all scheduled JSBSim models.
virtual bool Run(bool Holding)
Runs the model; called by the Executive.
double GetValue(void) const
Get the current table value.
virtual const FGColumnVector3 & GetWindNED(void) const
Retrieves the wind components in NED frame.
virtual void GustMagnitude(double mag)
Specifies the magnitude of the gust in feet/second.
virtual const FGColumnVector3 & GetTotalWindNED(void) const
Retrieves the total wind components in NED frame.
virtual void GustXComponent(double x)
Specifies the X component of velocity in the specified gust frame (ft/sec).
virtual void SetTurbNED(int idx, double turb)
Sets a turbulence component in NED frame.
virtual void SetWindPsi(double dir)
Sets the direction that the wind is coming from.
virtual void StartupGustDuration(double dur)
Specifies the duration of the startup portion of the gust.
bool Run(bool Holding) override
Runs the winds model; called by the Executive Can pass in a value indicating if the executive is dire...
virtual void SetTurbType(tType tt)
Turbulence models available: ttNone, ttStandard, ttBerndt, ttCulp, ttMilspec, ttTustin.
virtual double GetTurbNED(int idx) const
Retrieves a turbulence component in NED frame.
virtual void SteadyGustDuration(double dur)
Specifies the length of time that the gust is at a steady, full strength.
virtual void EndGustDuration(double dur)
Specifies the length of time it takes for the gust to return to zero velocity.
virtual void SetWindNED(double wN, double wE, double wD)
Sets the wind components in NED frame.
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
virtual void GustFrame(eGustFrame gFrame)
Specifies the frame that the gust direction vector components are specified in.
FGWinds(FGFDMExec *)
Constructor.
virtual void GustYComponent(double y)
Specifies the Y component of velocity in the specified gust frame (ft/sec).
virtual void StartGust(bool running)
Initiates the execution of the gust.
virtual const FGColumnVector3 & GetGustNED(void) const
Retrieves the gust components in NED frame.
virtual double GetWindPsi(void) const
Retrieves the direction that the wind is coming from.
virtual void SetGustNED(int idx, double gust)
Sets a gust component in NED frame.
Main namespace for the JSBSim Flight Dynamics Model.
constexpr double sqr(double x)
simply square a value
FGColumnVector3 vWindTransformed
struct OneMinusCosineProfile gustProfile