48#include "math/FGTable.h"
72constexpr double sqr(
double x) {
return x*x; }
75 :
FGModel(fdmex), generator(fdmex->GetRandomGenerator())
79 MagnitudedAccelDt = MagnitudeAccel = Magnitude = TurbDirection = 0.0;
84 spike = target_time = strength = 0.0;
85 wind_from_clockwise = 0.0;
88 vGustNED.InitMatrix();
89 vTurbulenceNED.InitMatrix();
90 vCosineGust.InitMatrix();
93 windspeed_at_20ft = 0.;
94 probability_of_exceedence_index = 0;
99 << 500.0 << 1750.0 << 3750.0 << 7500.0 << 15000.0 << 25000.0 << 35000.0 << 45000.0 << 55000.0 << 65000.0 << 75000.0 << 80000.0
100 << 1 << 3.2 << 2.2 << 1.5 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0
101 << 2 << 4.2 << 3.6 << 3.3 << 1.6 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0
102 << 3 << 6.6 << 6.9 << 7.4 << 6.7 << 4.6 << 2.7 << 0.4 << 0.0 << 0.0 << 0.0 << 0.0 << 0.0
103 << 4 << 8.6 << 9.6 << 10.6 << 10.1 << 8.0 << 6.6 << 5.0 << 4.2 << 2.7 << 0.0 << 0.0 << 0.0
104 << 5 << 11.8 << 13.0 << 16.0 << 15.1 << 11.6 << 9.7 << 8.1 << 8.2 << 7.9 << 4.9 << 3.2 << 2.1
105 << 6 << 15.6 << 17.6 << 23.0 << 23.6 << 22.1 << 20.0 << 16.0 << 15.1 << 12.1 << 7.9 << 6.2 << 5.1
106 << 7 << 18.7 << 21.5 << 28.4 << 30.2 << 30.7 << 31.0 << 25.2 << 23.1 << 17.5 << 10.7 << 8.4 << 7.2;
122bool FGWinds::InitModel(
void)
124 if (!FGModel::InitModel())
return false;
128 vGustNED.InitMatrix();
129 vTurbulenceNED.InitMatrix();
130 vCosineGust.InitMatrix();
135 xi_u_km1 = nu_u_km1 = 0;
136 xi_v_km1 = xi_v_km2 = nu_v_km1 = nu_v_km2 = 0;
137 xi_w_km1 = xi_w_km2 = nu_w_km1 = nu_w_km2 = 0;
138 xi_p_km1 = nu_p_km1 = 0;
139 xi_q_km1 = xi_r_km1 = 0;
149 if (Holding)
return false;
151 if (turbType != ttNone)
152 Turbulence(in.AltitudeASL);
154 vTurbulenceNED.InitMatrix();
157 vTotalWindNED = vWindNED + vGustNED + vCosineGust + vTurbulenceNED;
160 if (vWindNED(eX) != 0.0) psiw = atan2( vWindNED(eY), vWindNED(eX) );
161 if (psiw < 0) psiw += 2*M_PI;
171void FGWinds::SetWindspeed(
double speed)
175 vWindNED(eNorth) = speed;
177 vWindNED(eNorth) = speed * cos(psiw);
178 vWindNED(eEast) = speed * sin(psiw);
179 vWindNED(eDown) = 0.0;
185double FGWinds::GetWindspeed(
void)
const
196 double mag = GetWindspeed();
203void FGWinds::Turbulence(
double h)
209 vTurbPQR(eP) = wind_from_clockwise;
210 if (TurbGain == 0.0)
return;
213 if (TurbGain < 0.0) TurbGain = 0.0;
214 if (TurbGain > 1.0) TurbGain = 1.0;
215 if (TurbRate < 0.0) TurbRate = 0.0;
216 if (TurbRate > 30.0) TurbRate = 30.0;
217 if (Rhythmicity < 0.0) Rhythmicity = 0.0;
218 if (Rhythmicity > 1.0) Rhythmicity = 1.0;
222 double sinewave = sin( time * TurbRate * 6.283185307 );
225 if (target_time == 0.0) {
226 strength = random = generator->GetUniformRandomNumber();
227 target_time = time + 0.71 + (random * 0.5);
229 if (time > target_time) {
237 vTurbulenceNED.InitMatrix();
238 double delta = strength * max_vs * TurbGain * (1-Rhythmicity) * spike;
241 vTurbulenceNED(eDown) = sinewave * max_vs * TurbGain * Rhythmicity;
242 vTurbulenceNED(eDown)+= delta;
243 if (in.DistanceAGL/in.wingspan < 3.0)
244 vTurbulenceNED(eDown) *= in.DistanceAGL/in.wingspan * 0.3333;
247 vTurbulenceNED(eNorth) = sin( delta * 3.0 );
248 vTurbulenceNED(eEast) = cos( delta * 3.0 );
251 vTurbPQR(eP) += delta * 0.04;
261 if (probability_of_exceedence_index == 0 || in.V == 0) {
262 vTurbulenceNED(eNorth) = vTurbulenceNED(eEast) = vTurbulenceNED(eDown) = 0.0;
263 vTurbPQR(eP) = vTurbPQR(eQ) = vTurbPQR(eR) = 0.0;
268 double b_w = in.wingspan, L_u, L_w, sig_u, sig_w;
270 if (b_w == 0.) b_w = 30.;
273 if (h <= 10.) h = 10;
277 L_u = h/pow(0.177 + 0.000823*h, 1.2);
279 sig_w = 0.1*windspeed_at_20ft;
280 sig_u = sig_w/pow(0.177 + 0.000823*h, 0.4);
281 }
else if (h <= 2000) {
283 L_u = L_w = 1000 + (h-1000.)/1000.*750.;
284 sig_u = sig_w = 0.1*windspeed_at_20ft
285 + (h-1000.)/1000.*(POE_Table->
GetValue(probability_of_exceedence_index, h) - 0.1*windspeed_at_20ft);
288 sig_u = sig_w = POE_Table->
GetValue(probability_of_exceedence_index, h);
292 T_V = in.totalDeltaT,
293 sig_p = 1.9/sqrt(L_w*b_w)*sig_w,
296 L_p = sqrt(L_w*b_w)/2.6,
300 tau_q = 4*b_w/M_PI/in.V,
301 tau_r =3*b_w/M_PI/in.V,
302 nu_u = generator->GetNormalRandomNumber(),
303 nu_v = generator->GetNormalRandomNumber(),
304 nu_w = generator->GetNormalRandomNumber(),
305 nu_p = generator->GetNormalRandomNumber(),
306 xi_u=0, xi_v=0, xi_w=0, xi_p=0, xi_q=0, xi_r=0;
310 if (turbType == ttTustin) {
315 C_BL = 1/tau_u/tan(T_V/2/tau_u),
316 C_BLp = 1/tau_p/tan(T_V/2/tau_p),
317 C_BLq = 1/tau_q/tan(T_V/2/tau_q),
318 C_BLr = 1/tau_r/tan(T_V/2/tau_r);
324 xi_u = -(1 - C_BL*tau_u)/(1 + C_BL*tau_u)*xi_u_km1
325 + sig_u*sqrt(2*tau_u/T_V)/(1 + C_BL*tau_u)*(nu_u + nu_u_km1);
326 xi_v = -2*(sqr(omega_v) - sqr(C_BL))/sqr(omega_v + C_BL)*xi_v_km1
327 - sqr(omega_v - C_BL)/sqr(omega_v + C_BL) * xi_v_km2
328 + sig_u*sqrt(3*omega_v/T_V)/sqr(omega_v + C_BL)*(
329 (C_BL + omega_v/sqrt(3.))*nu_v
330 + 2/sqrt(3.)*omega_v*nu_v_km1
331 + (omega_v/sqrt(3.) - C_BL)*nu_v_km2);
332 xi_w = -2*(sqr(omega_w) - sqr(C_BL))/sqr(omega_w + C_BL)*xi_w_km1
333 - sqr(omega_w - C_BL)/sqr(omega_w + C_BL) * xi_w_km2
334 + sig_w*sqrt(3*omega_w/T_V)/sqr(omega_w + C_BL)*(
335 (C_BL + omega_w/sqrt(3.))*nu_w
336 + 2/sqrt(3.)*omega_w*nu_w_km1
337 + (omega_w/sqrt(3.) - C_BL)*nu_w_km2);
338 xi_p = -(1 - C_BLp*tau_p)/(1 + C_BLp*tau_p)*xi_p_km1
339 + sig_p*sqrt(2*tau_p/T_V)/(1 + C_BLp*tau_p) * (nu_p + nu_p_km1);
340 xi_q = -(1 - 4*b_w*C_BLq/M_PI/in.V)/(1 + 4*b_w*C_BLq/M_PI/in.V) * xi_q_km1
341 + C_BLq/in.V/(1 + 4*b_w*C_BLq/M_PI/in.V) * (xi_w - xi_w_km1);
342 xi_r = - (1 - 3*b_w*C_BLr/M_PI/in.V)/(1 + 3*b_w*C_BLr/M_PI/in.V) * xi_r_km1
343 + C_BLr/in.V/(1 + 3*b_w*C_BLr/M_PI/in.V) * (xi_v - xi_v_km1);
345 }
else if (turbType == ttMilspec) {
348 xi_u = (1 - T_V/tau_u) *xi_u_km1 + sig_u*sqrt(2*T_V/tau_u)*nu_u;
349 xi_v = (1 - 2*T_V/tau_u)*xi_v_km1 + sig_u*sqrt(4*T_V/tau_u)*nu_v;
350 xi_w = (1 - 2*T_V/tau_w)*xi_w_km1 + sig_w*sqrt(4*T_V/tau_w)*nu_w;
351 xi_p = (1 - T_V/tau_p) *xi_p_km1 + sig_p*sqrt(2*T_V/tau_p)*nu_p;
352 xi_q = (1 - T_V/tau_q) *xi_q_km1 + M_PI/4/b_w*(xi_w - xi_w_km1);
353 xi_r = (1 - T_V/tau_r) *xi_r_km1 + M_PI/3/b_w*(xi_v - xi_v_km1);
357 double cospsi = cos(psiw), sinpsi = sin(psiw);
358 vTurbulenceNED(eNorth) = cospsi*xi_u + sinpsi*xi_v;
359 vTurbulenceNED(eEast) = -sinpsi*xi_u + cospsi*xi_v;
360 vTurbulenceNED(eDown) = xi_w;
362 vTurbPQR(eP) = cospsi*xi_p + sinpsi*xi_q;
363 vTurbPQR(eQ) = -sinpsi*xi_p + cospsi*xi_q;
367 vTurbPQR = in.Tl2b*vTurbPQR;
370 xi_u_km1 = xi_u; nu_u_km1 = nu_u;
371 xi_v_km2 = xi_v_km1; xi_v_km1 = xi_v; nu_v_km2 = nu_v_km1; nu_v_km1 = nu_v;
372 xi_w_km2 = xi_w_km1; xi_w_km1 = xi_w; nu_w_km2 = nu_w_km1; nu_w_km1 = nu_w;
373 xi_p_km1 = xi_p; nu_p_km1 = nu_p;
382 TurbDirection = atan2( vTurbulenceNED(eEast), vTurbulenceNED(eNorth))*radtodeg;
388double FGWinds::CosineGustProfile(
double startDuration,
double steadyDuration,
double endDuration,
double elapsedTime)
391 if (elapsedTime >= 0 && elapsedTime <= startDuration) {
392 factor = (1.0 - cos(M_PI*elapsedTime/startDuration))/2.0;
393 }
else if (elapsedTime > startDuration && (elapsedTime <= (startDuration + steadyDuration))) {
395 }
else if (elapsedTime > (startDuration + steadyDuration) && elapsedTime <= (startDuration + steadyDuration + endDuration)) {
396 factor = (1-cos(M_PI*(1-(elapsedTime-(startDuration + steadyDuration))/endDuration)))/2.0;
406void FGWinds::CosineGust()
408 struct OneMinusCosineProfile& profile = oneMinusCosineGust.
gustProfile;
410 double factor = CosineGustProfile( profile.startupDuration,
411 profile.steadyDuration,
413 profile.elapsedTime);
436 profile.elapsedTime += in.totalDeltaT;
438 if (profile.elapsedTime > (profile.startupDuration + profile.steadyDuration + profile.endDuration)) {
439 profile.Running =
false;
440 profile.elapsedTime = 0.0;
442 vCosineGust.InitMatrix(0);
448void FGWinds::NumberOfUpDownburstCells(
int num)
450 for (
unsigned int i=0; i<UpDownBurstCells.size();i++)
delete UpDownBurstCells[i];
451 UpDownBurstCells.clear();
453 for (
int i=0; i<num; i++) UpDownBurstCells.push_back(
new struct UpDownBurst);
462double FGWinds::DistanceFromRingCenter(
double lat,
double lon)
464 double deltaLat = in.latitude - lat;
465 double deltaLong = in.longitude - lon;
466 double dLat2 = deltaLat/2.0;
467 double dLong2 = deltaLong/2.0;
468 double a = sin(dLat2)*sin(dLat2)
469 + cos(lat)*cos(in.latitude)*sin(dLong2)*sin(dLong2);
470 double c = 2.0*atan2(sqrt(a), sqrt(1.0 - a));
471 double d = in.planetRadius*c;
477void FGWinds::UpDownBurst()
480 for (
unsigned int i=0; i<UpDownBurstCells.size(); i++) {
481 DistanceFromRingCenter(UpDownBurstCells[i]->ringLatitude, UpDownBurstCells[i]->ringLongitude);
491void FGWinds::SetRandomSeed(
int sr)
494 generator = std::make_shared<RandomNumberGenerator>(*RandomSeed);
497int FGWinds::GetRandomSeed(
void)
const {
501 return FDMExec->SRand();
506void FGWinds::bind(
void)
508 typedef double (FGWinds::*PMF)(int) const;
509 typedef int (FGWinds::*PMFt)(void) const;
510 typedef void (FGWinds::*PMFd)(int,double);
511 typedef void (FGWinds::*PMFi)(int);
512 typedef double (FGWinds::*Ptr)(void) const;
522 PropertyManager->Tie(
"atmosphere/wind-mag-fps",
this, &FGWinds::GetWindspeed,
523 &FGWinds::SetWindspeed);
538 PropertyManager->Tie(
"atmosphere/cosine-gust/frame",
this, (PMFt)0L, (PMFi)&
FGWinds::GustFrame);
542 PropertyManager->Tie(
"atmosphere/cosine-gust/start",
this,
static_cast<bool (
FGWinds::*)(
void) const
>(
nullptr), &
FGWinds::StartGust);
545 PropertyManager->Tie(
"atmosphere/updownburst/number-of-cells",
this, (PMFt)0L, &FGWinds::NumberOfUpDownburstCells);
562 PropertyManager->Tie(
"atmosphere/p-turb-rad_sec",
this,1, (PMF)&FGWinds::GetTurbPQR);
563 PropertyManager->Tie(
"atmosphere/q-turb-rad_sec",
this,2, (PMF)&FGWinds::GetTurbPQR);
564 PropertyManager->Tie(
"atmosphere/r-turb-rad_sec",
this,3, (PMF)&FGWinds::GetTurbPQR);
565 PropertyManager->Tie(
"atmosphere/turb-type",
this, (PMFt)&FGWinds::GetTurbType, (PMFi)&
FGWinds::SetTurbType);
566 PropertyManager->Tie(
"atmosphere/turb-rate",
this, &FGWinds::GetTurbRate, &FGWinds::SetTurbRate);
567 PropertyManager->Tie(
"atmosphere/turb-gain",
this, &FGWinds::GetTurbGain, &FGWinds::SetTurbGain);
568 PropertyManager->Tie(
"atmosphere/turb-rhythmicity",
this, &FGWinds::GetRhythmicity,
569 &FGWinds::SetRhythmicity);
572 PropertyManager->Tie(
"atmosphere/turbulence/milspec/windspeed_at_20ft_AGL-fps",
573 this, &FGWinds::GetWindspeed20ft,
574 &FGWinds::SetWindspeed20ft);
575 PropertyManager->Tie(
"atmosphere/turbulence/milspec/severity",
576 this, &FGWinds::GetProbabilityOfExceedence,
585 PropertyManager->Tie(
"atmosphere/randomseed",
this, (PMFt)&FGWinds::GetRandomSeed, &FGWinds::SetRandomSeed);
607void FGWinds::Debug(
int from)
609 if (debug_lvl <= 0)
return;
615 if (debug_lvl & 2 ) {
616 if (from == 0) cout <<
"Instantiated: FGWinds" << endl;
617 if (from == 1) cout <<
"Destroyed: FGWinds" << endl;
619 if (debug_lvl & 4 ) {
621 if (debug_lvl & 8 ) {
623 if (debug_lvl & 16) {
625 if (debug_lvl & 128) {
627 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.
FGColumnVector3 vWindTransformed
struct OneMinusCosineProfile gustProfile