JSBSim Flight Dynamics Model  1.0 (02 March 2017)
An Open Source Flight Dynamics and Control Software Library in C++
FGStandardAtmosphere Class Reference

Models the 1976 U.S. More...

#include <FGStandardAtmosphere.h>

+ Inheritance diagram for FGStandardAtmosphere:
+ Collaboration diagram for FGStandardAtmosphere:

Public Member Functions

 FGStandardAtmosphere (FGFDMExec *)
 Constructor.
 
virtual ~FGStandardAtmosphere ()
 Destructor.
 
bool InitModel (void)
 
virtual void PrintStandardAtmosphereTable ()
 Prints the U.S. Standard Atmosphere table.
 
Temperature access functions.

There are several ways to get the temperature, and several modeled temperature values that can be retrieved.

The U.S. Standard Atmosphere temperature either at a specified altitude, or at sea level can be retrieved. These two temperatures do NOT include the effects of any bias or delta gradient that may have been supplied by the user. The modeled temperature and the modeled temperature at sea level can also be retrieved. These two temperatures DO include the effects of an optionally user-supplied bias or delta gradient.

virtual double GetTemperature (double altitude) const
 Returns the actual modeled temperature in degrees Rankine at a specified altitude. More...
 
virtual double GetStdTemperature (double altitude) const
 Returns the standard temperature in degrees Rankine at a specified altitude. More...
 
virtual double GetStdTemperatureSL () const
 Returns the standard sea level temperature in degrees Rankine. More...
 
virtual double GetStdTemperatureRatio (double h) const
 Returns the ratio of the standard temperature at the supplied altitude over the standard sea level temperature. More...
 
virtual double GetTemperatureBias (eTemperature to) const
 Returns the temperature bias over the sea level value in degrees Rankine.
 
virtual double GetTemperatureDeltaGradient (eTemperature to)
 Returns the temperature gradient to be applied on top of the standard temperature gradient. More...
 
virtual void SetTemperatureSL (double t, eTemperature unit=eFahrenheit)
 Sets the Sea Level temperature, if it is to be different than the standard. More...
 
virtual void SetTemperature (double t, double h, eTemperature unit=eFahrenheit)
 Sets the temperature at the supplied altitude, if it is to be different than the standard temperature. More...
 
virtual void SetTemperatureBias (eTemperature unit, double t)
 Sets the temperature bias to be added to the standard temperature at all altitudes. More...
 
virtual void SetSLTemperatureGradedDelta (eTemperature unit, double t)
 Sets a Sea Level temperature delta that is ramped out by 86 km. More...
 
virtual void SetTemperatureGradedDelta (double t, double h, eTemperature unit=eFahrenheit)
 Sets the temperature delta value at the supplied altitude/elevation above sea level, to be added to the standard temperature and ramped out by 86 km. More...
 
virtual void ResetSLTemperature ()
 This function resets the model to apply no bias or delta gradient to the temperature. More...
 
Pressure access functions.
virtual double GetPressure (double altitude) const
 Returns the pressure at a specified altitude in psf.
 
virtual double GetStdPressure100K (double altitude) const
 Returns the standard pressure at a specified altitude in psf.
 
virtual double GetStdPressure (double altitude) const
 Returns the standard pressure at the specified altitude.
 
virtual void SetPressureSL (ePressure unit, double pressure)
 Sets the sea level pressure for modeling an off-standard pressure profile. More...
 
virtual void ResetSLPressure ()
 Resets the sea level to the Standard sea level pressure, and recalculates dependent parameters so that the pressure calculations are standard. More...
 
Density access functions.
virtual double GetStdDensity (double altitude) const
 Returns the standard density at a specified altitude.
 
- Public Member Functions inherited from FGAtmosphere
 FGAtmosphere (FGFDMExec *)
 Constructor.
 
virtual ~FGAtmosphere ()
 Destructor.
 
virtual double GetDensityAltitude () const
 
virtual double GetPressureAltitude () const
 
bool Run (bool Holding)
 Runs the atmosphere forces model; called by the Executive. More...
 
virtual double GetTemperature () const
 Returns the actual, modeled temperature at the current altitude in degrees Rankine. More...
 
virtual double GetTemperatureSL () const
 Returns the actual, modeled sea level temperature in degrees Rankine. More...
 
virtual double GetTemperatureRatio () const
 Returns the ratio of the at-current-altitude temperature as modeled over the sea level value. More...
 
virtual double GetTemperatureRatio (double h) const
 Returns the ratio of the temperature as modeled at the supplied altitude over the sea level value. More...
 
virtual double GetPressure (void) const
 Returns the pressure in psf.
 
virtual double GetPressureSL (ePressure to=ePSF) const
 
virtual double GetPressureRatio (void) const
 Returns the ratio of at-altitude pressure over the sea level value.
 
virtual double GetDensity (void) const
 Returns the density in slugs/ft^3. More...
 
virtual double GetDensity (double altitude) const
 Returns the density in slugs/ft^3 at a given altitude in ft. More...
 
virtual double GetDensitySL (void) const
 Returns the sea level density in slugs/ft^3.
 
virtual double GetDensityRatio (void) const
 Returns the ratio of at-altitude density over the sea level value.
 
virtual double GetSoundSpeed (void) const
 Returns the speed of sound in ft/sec.
 
virtual double GetSoundSpeed (double altitude) const
 Returns the speed of sound in ft/sec at a given altitude in ft.
 
virtual double GetSoundSpeedSL (void) const
 Returns the sea level speed of sound in ft/sec.
 
virtual double GetSoundSpeedRatio (void) const
 Returns the ratio of at-altitude sound speed over the sea level value.
 
virtual double GetAbsoluteViscosity (void) const
 Returns the absolute viscosity.
 
virtual double GetKinematicViscosity (void) const
 Returns the kinematic viscosity.
 
- Public Member Functions inherited from FGModel
 FGModel (FGFDMExec *)
 Constructor.
 
virtual ~FGModel ()
 Destructor.
 
virtual SGPath FindFullPathName (const SGPath &path) const
 
FGFDMExecGetExec (void)
 
unsigned int GetRate (void)
 Get the output rate for the model in frames.
 
void SetPropertyManager (FGPropertyManager *fgpm)
 
void SetRate (unsigned int tt)
 Set the ouput rate for the model in frames.
 
- Public Member Functions inherited from FGModelFunctions
std::string GetFunctionStrings (const std::string &delimeter) const
 Gets the strings for the current set of functions. More...
 
std::string GetFunctionValues (const std::string &delimeter) const
 Gets the function values. More...
 
FGFunctionGetPreFunction (const std::string &name)
 Get one of the "pre" function. More...
 
bool Load (Element *el, FGPropertyManager *PropertyManager, std::string prefix="")
 
void PostLoad (Element *el, FGPropertyManager *PropertyManager, std::string prefix="")
 
void PreLoad (Element *el, FGPropertyManager *PropertyManager, std::string prefix="")
 
void RunPostFunctions (void)
 
void RunPreFunctions (void)
 
- Public Member Functions inherited from FGJSBBase
 FGJSBBase ()
 Constructor for FGJSBBase.
 
virtual ~FGJSBBase ()
 Destructor for FGJSBBase.
 
void disableHighLighting (void)
 Disables highlighting in the console output.
 
std::string GetVersion (void)
 Returns the version number of JSBSim. More...
 
void PutMessage (const Message &msg)
 Places a Message structure on the Message queue. More...
 
void PutMessage (const std::string &text)
 Creates a message with the given text and places it on the queue. More...
 
void PutMessage (const std::string &text, bool bVal)
 Creates a message with the given text and boolean value and places it on the queue. More...
 
void PutMessage (const std::string &text, int iVal)
 Creates a message with the given text and integer value and places it on the queue. More...
 
void PutMessage (const std::string &text, double dVal)
 Creates a message with the given text and double value and places it on the queue. More...
 
int SomeMessages (void)
 Reads the message on the queue (but does not delete it). More...
 
void ProcessMessage (void)
 Reads the message on the queue and removes it from the queue. More...
 
MessageProcessNextMessage (void)
 Reads the next message on the queue and removes it from the queue. More...
 

Protected Member Functions

virtual void bind (void)
 
void CalculateLapseRates ()
 Recalculate the lapse rate vectors when the temperature profile is altered in a way that would change the lapse rates, such as when a gradient is applied. More...
 
void CalculatePressureBreakpoints ()
 Calculate (or recalculate) the atmospheric pressure breakpoints at the altitudes in the standard temperature table. More...
 
void Debug (int from)
 
- Protected Member Functions inherited from FGAtmosphere
void Calculate (double altitude)
 Calculate the atmosphere for the given altitude.
 
virtual double ConvertFromPSF (double t, ePressure unit=ePSF) const
 
virtual double ConvertToPSF (double t, ePressure unit=ePSF) const
 
virtual double ConvertToRankine (double t, eTemperature unit) const
 
- Protected Member Functions inherited from FGModel
virtual bool Load (Element *el)
 Loads this model. More...
 
- Protected Member Functions inherited from FGJSBBase
void Debug (int)
 

Protected Attributes

double GradientFadeoutAltitude
 
std::vector< double > LapseRateVector
 
std::vector< double > PressureBreakpointVector
 
FGTableStdAtmosTemperatureTable
 
double StdSLdensity
 
double StdSLpressure
 
double StdSLsoundspeed
 
double StdSLtemperature
 
double TemperatureBias
 
double TemperatureDeltaGradient
 
- Protected Attributes inherited from FGAtmosphere
const double Beta
 
double Density
 
double DensityAltitude
 
double KinematicViscosity
 
double Pressure
 
double PressureAltitude
 
double rSLdensity
 
double rSLpressure
 
double rSLsoundspeed
 
double rSLtemperature
 
double SLdensity
 
double SLpressure
 
double SLsoundspeed
 
double SLtemperature
 
double Soundspeed
 
const double SutherlandConstant
 
double Temperature
 
double Viscosity
 
- Protected Attributes inherited from FGModel
unsigned int exe_ctr
 
FGFDMExecFDMExec
 
FGPropertyManagerPropertyManager
 
unsigned int rate
 
- Protected Attributes inherited from FGModelFunctions
FGPropertyReader LocalProperties
 
std::vector< FGFunction * > PostFunctions
 
std::vector< FGFunction * > PreFunctions
 

Additional Inherited Members

- Public Types inherited from FGAtmosphere
enum  ePressure {
  eNoPressUnit =0, ePSF, eMillibars, ePascals,
  eInchesHg
}
 Enums for specifying pressure units.
 
enum  eTemperature {
  eNoTempUnit =0, eFahrenheit, eCelsius, eRankine,
  eKelvin
}
 Enums for specifying temperature units.
 
- Public Types inherited from FGJSBBase
enum  { eL = 1, eM, eN }
 Moments L, M, N.
 
enum  { eP = 1, eQ, eR }
 Rates P, Q, R.
 
enum  { eU = 1, eV, eW }
 Velocities U, V, W.
 
enum  { eX = 1, eY, eZ }
 Positions X, Y, Z.
 
enum  { ePhi = 1, eTht, ePsi }
 Euler angles Phi, Theta, Psi.
 
enum  { eDrag = 1, eSide, eLift }
 Stability axis forces, Drag, Side force, Lift.
 
enum  { eRoll = 1, ePitch, eYaw }
 Local frame orientation Roll, Pitch, Yaw.
 
enum  { eNorth = 1, eEast, eDown }
 Local frame position North, East, Down.
 
enum  { eLat = 1, eLong, eRad }
 Locations Radius, Latitude, Longitude.
 
enum  {
  inNone = 0, inDegrees, inRadians, inMeters,
  inFeet
}
 Conversion specifiers.
 
- Static Public Member Functions inherited from FGJSBBase
static double CelsiusToFahrenheit (double celsius)
 Converts from degrees Celsius to degrees Fahrenheit. More...
 
static double CelsiusToKelvin (double celsius)
 Converts from degrees Celsius to degrees Kelvin. More...
 
static double CelsiusToRankine (double celsius)
 Converts from degrees Celsius to degrees Rankine. More...
 
static double Constrain (double min, double value, double max)
 Constrain a value between a minimum and a maximum value.
 
static bool EqualToRoundoff (double a, double b)
 Finite precision comparison. More...
 
static bool EqualToRoundoff (float a, float b)
 Finite precision comparison. More...
 
static bool EqualToRoundoff (float a, double b)
 Finite precision comparison. More...
 
static bool EqualToRoundoff (double a, float b)
 Finite precision comparison. More...
 
static double FahrenheitToCelsius (double fahrenheit)
 Converts from degrees Fahrenheit to degrees Celsius. More...
 
static double FeetToMeters (double measure)
 Converts from feet to meters. More...
 
static double GaussianRandomNumber (void)
 
static double KelvinToCelsius (double kelvin)
 Converts from degrees Kelvin to degrees Celsius. More...
 
static double KelvinToFahrenheit (double kelvin)
 Converts from degrees Kelvin to degrees Fahrenheit. More...
 
static double KelvinToRankine (double kelvin)
 Converts from degrees Kelvin to degrees Rankine. More...
 
static double MachFromVcalibrated (double vcas, double p, double psl, double rhosl)
 Calculate the Mach number from the calibrated airspeed. More...
 
static double PitotTotalPressure (double mach, double p)
 Compute the total pressure in front of the Pitot tube. More...
 
static double RankineToCelsius (double rankine)
 Converts from degrees Rankine to degrees Celsius. More...
 
static double RankineToKelvin (double rankine)
 Converts from degrees Rankine to degrees Kelvin. More...
 
static double sign (double num)
 
static double VcalibratedFromMach (double mach, double p, double psl, double rhosl)
 Calculate the calibrated airspeed from the Mach number. More...
 
- Public Attributes inherited from FGAtmosphere
struct JSBSim::FGAtmosphere::Inputs in
 
- Public Attributes inherited from FGModel
std::string Name
 
- Static Public Attributes inherited from FGJSBBase
static short debug_lvl = 1
 
static char highint [5] = {27, '[', '1', 'm', '\0' }
 highlights text
 
static char halfint [5] = {27, '[', '2', 'm', '\0' }
 low intensity text
 
static char normint [6] = {27, '[', '2', '2', 'm', '\0' }
 normal intensity text
 
static char reset [5] = {27, '[', '0', 'm', '\0' }
 resets text properties
 
static char underon [5] = {27, '[', '4', 'm', '\0' }
 underlines text
 
static char underoff [6] = {27, '[', '2', '4', 'm', '\0' }
 underline off
 
static char fgblue [6] = {27, '[', '3', '4', 'm', '\0' }
 blue text
 
static char fgcyan [6] = {27, '[', '3', '6', 'm', '\0' }
 cyan text
 
static char fgred [6] = {27, '[', '3', '1', 'm', '\0' }
 red text
 
static char fggreen [6] = {27, '[', '3', '2', 'm', '\0' }
 green text
 
static char fgdef [6] = {27, '[', '3', '9', 'm', '\0' }
 default text
 
- Static Protected Member Functions inherited from FGJSBBase
static std::string CreateIndexedPropertyName (const std::string &Property, int index)
 
- Static Protected Attributes inherited from FGJSBBase
static const double degtorad = 0.017453292519943295769236907684886
 
static const double fpstokts = 1.0/ktstofps
 
static const double fttom = 0.3048
 
static int gaussian_random_number_phase = 0
 
static const double hptoftlbssec = 550.0
 
static const double in3tom3 = 1.638706E-5
 
static const double inchtoft = 0.08333333
 
static const double inhgtopa = 3386.38
 
static const std::string JSBSim_version = "1.0 " __DATE__ " " __TIME__
 
static const double kgtolb = 2.20462
 
static const double kgtoslug = 0.06852168
 
static const double ktstofps = 1.68781
 
static const double lbtoslug = 1.0/slugtolb
 
static Message localMsg
 
static const double m3toft3 = 1.0/(fttom*fttom*fttom)
 
static double Mair = 28.9645
 
static unsigned int messageId = 0
 
static std::queue< MessageMessages
 
static const std::string needed_cfg_version = "2.0"
 
static const double psftoinhg = 0.014138
 
static const double psftopa = 47.88
 
static const double radtodeg = 57.295779513082320876798154814105
 
static double Reng = 1716.56
 
static double Rstar = 1545.348
 
static const double SHRatio = 1.40
 
static const double slugtolb = 32.174049
 

Detailed Description

Models the 1976 U.S.

Standard Atmosphere, with the ability to modify the temperature and pressure. A base feature of the model is the temperature profile that is stored as an FGTable object with this data:

GeoMet Alt Temp GeoPot Alt GeoMet Alt
(ft) (deg R) (km) (km)
--------- -------- ---------- ----------
0.0 518.67 // 0.000 0.000
36151.6 390.0 // 11.000 11.019
65823.5 390.0 // 20.000 20.063
105518.4 411.6 // 32.000 32.162
155347.8 487.2 // 47.000 47.350
168677.8 487.2 // 51.000 51.413
235570.9 386.4 // 71.000 71.802
282152.2 336.5; // 84.852 86.000

The pressure is calculated at lower altitudes through the use of two equations that are presented in the U.S. Standard Atmosphere document (see references). Density, kinematic viscosity, speed of sound, etc., are all calculated based on various constants and temperature and pressure. At higher altitudes (above 86 km (282152 ft) a different and more complicated method of calculating pressure is used. The temperature may be modified through the use of several methods. Ultimately, these access methods allow the user to modify the sea level standard temperature, and/or the sea level standard pressure, so that the entire profile will be consistently and accurately calculated.

Properties

  • atmosphere/delta-T
  • atmosphere/T-sl-dev-F
Author
Jon Berndt
See also
"U.S. Standard Atmosphere, 1976", NASA TM-X-74335
Version
Id
FGStandardAtmosphere.h,v 1.17 2012/04/13 13:18:27 jberndt Exp

Definition at line 103 of file FGStandardAtmosphere.h.

Member Function Documentation

◆ CalculateLapseRates()

void CalculateLapseRates ( )
protected

Recalculate the lapse rate vectors when the temperature profile is altered in a way that would change the lapse rates, such as when a gradient is applied.

This function is also called to initialize the lapse rate vector.

Definition at line 377 of file FGStandardAtmosphere.cpp.

378 {
379  for (unsigned int bh=0; bh<LapseRateVector.size(); bh++)
380  {
381  double t0 = (*StdAtmosTemperatureTable)(bh+1,1);
382  double t1 = (*StdAtmosTemperatureTable)(bh+2,1);
383  double h0 = (*StdAtmosTemperatureTable)(bh+1,0);
384  double h1 = (*StdAtmosTemperatureTable)(bh+2,0);
385  LapseRateVector[bh] = (t1 - t0) / (h1 - h0) + TemperatureDeltaGradient;
386  }
387 }
+ Here is the caller graph for this function:

◆ CalculatePressureBreakpoints()

void CalculatePressureBreakpoints ( )
protected

Calculate (or recalculate) the atmospheric pressure breakpoints at the altitudes in the standard temperature table.

Definition at line 391 of file FGStandardAtmosphere.cpp.

392 {
393  for (unsigned int b=0; b<PressureBreakpointVector.size()-1; b++) {
394  double BaseTemp = (*StdAtmosTemperatureTable)(b+1,1);
395  double BaseAlt = (*StdAtmosTemperatureTable)(b+1,0);
396  double UpperAlt = (*StdAtmosTemperatureTable)(b+2,0);
397  double deltaH = UpperAlt - BaseAlt;
398  double Tmb = BaseTemp
399  + TemperatureBias
400  + (GradientFadeoutAltitude - BaseAlt)*TemperatureDeltaGradient;
401  if (LapseRateVector[b] != 0.00) {
402  double Lmb = LapseRateVector[b];
403  double Exp = Mair/(Rstar*Lmb);
404  double factor = Tmb/(Tmb + Lmb*deltaH);
405  PressureBreakpointVector[b+1] = PressureBreakpointVector[b]*pow(factor, Exp);
406  } else {
407  PressureBreakpointVector[b+1] = PressureBreakpointVector[b]*exp(-Mair*deltaH/(Rstar*Tmb));
408  }
409  }
410 }
+ Here is the caller graph for this function:

◆ GetStdTemperature()

double GetStdTemperature ( double  altitude) const
virtual

Returns the standard temperature in degrees Rankine at a specified altitude.

Parameters
altitudeThe altitude in feet above sea level (ASL) to get the temperature at.
Returns
The STANDARD temperature in degrees Rankine at the specified altitude.

Definition at line 213 of file FGStandardAtmosphere.cpp.

214 {
215  double Lk9 = 0.00658368; // deg R per foot
216  double Tinf = 1800.0; // Same as 1000 Kelvin
217  double temp = Tinf;
218 
219  if (altitude < 298556.4) { // 91 km - station 8
220 
221  double GeoPotAlt = (altitude*20855531.5)/(20855531.5+altitude);
222  temp = StdAtmosTemperatureTable->GetValue(GeoPotAlt);
223 
224  } else if (altitude < 360892.4) { // 110 km - station 9
225 
226  temp = 473.7429 - 137.38176 * sqrt(1.0 - pow((altitude - 298556.4)/65429.462, 2.0));
227 
228  } else if (altitude < 393700.8) { // 120 km - station 10
229 
230  temp = 432 + Lk9 * (altitude - 360892.4);
231 
232  } else if (altitude < 3280839.9) { // 1000 km station 12
233 
234  double lambda = 0.00001870364;
235  double eps = (altitude - 393700.8) * (20855531.5 + 393700.8) / (20855531.5 + altitude);
236  temp = Tinf - (Tinf - 648.0) * exp(-lambda*eps);
237 
238  }
239 
240  return temp;
241 }
+ Here is the caller graph for this function:

◆ GetStdTemperatureRatio()

virtual double GetStdTemperatureRatio ( double  h) const
inlinevirtual

Returns the ratio of the standard temperature at the supplied altitude over the standard sea level temperature.

Definition at line 139 of file FGStandardAtmosphere.h.

139 { return GetStdTemperature(h)*rSLtemperature; }
virtual double GetStdTemperature(double altitude) const
Returns the standard temperature in degrees Rankine at a specified altitude.
+ Here is the call graph for this function:

◆ GetStdTemperatureSL()

virtual double GetStdTemperatureSL ( ) const
inlinevirtual

Returns the standard sea level temperature in degrees Rankine.

Returns
The STANDARD temperature at sea level in degrees Rankine.

Definition at line 135 of file FGStandardAtmosphere.h.

135 { return GetStdTemperature(0.0); }
virtual double GetStdTemperature(double altitude) const
Returns the standard temperature in degrees Rankine at a specified altitude.
+ Here is the call graph for this function:

◆ GetTemperature()

double GetTemperature ( double  altitude) const
virtual

Returns the actual modeled temperature in degrees Rankine at a specified altitude.

Parameters
altitudeThe altitude above sea level (ASL) in feet.
Returns
Modeled temperature in degrees Rankine at the specified altitude.

Implements FGAtmosphere.

Definition at line 199 of file FGStandardAtmosphere.cpp.

200 {
201  double GeoPotAlt = (altitude*20855531.5)/(20855531.5+altitude);
202 
203  double T = StdAtmosTemperatureTable->GetValue(GeoPotAlt) + TemperatureBias;
204  if (altitude <= GradientFadeoutAltitude)
205  T += TemperatureDeltaGradient * (GradientFadeoutAltitude - altitude);
206 
207  return T;
208 }

◆ GetTemperatureDeltaGradient()

virtual double GetTemperatureDeltaGradient ( eTemperature  to)
inlinevirtual

Returns the temperature gradient to be applied on top of the standard temperature gradient.

Definition at line 147 of file FGStandardAtmosphere.h.

148  { if (to == eCelsius || to == eKelvin) return TemperatureDeltaGradient/1.80; else return TemperatureDeltaGradient; }
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ResetSLPressure()

void ResetSLPressure ( )
virtual

Resets the sea level to the Standard sea level pressure, and recalculates dependent parameters so that the pressure calculations are standard.

Definition at line 423 of file FGStandardAtmosphere.cpp.

424 {
425  PressureBreakpointVector[0] = StdSLpressure; // psf
427 }
void CalculatePressureBreakpoints()
Calculate (or recalculate) the atmospheric pressure breakpoints at the altitudes in the standard temp...
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ResetSLTemperature()

void ResetSLTemperature ( )
virtual

This function resets the model to apply no bias or delta gradient to the temperature.

The delta gradient and bias values are reset to 0.0, and the standard temperature is used for the entire temperature profile at all altitudes.

Definition at line 414 of file FGStandardAtmosphere.cpp.

415 {
416  TemperatureBias = TemperatureDeltaGradient = 0.0;
419 }
void CalculateLapseRates()
Recalculate the lapse rate vectors when the temperature profile is altered in a way that would change...
void CalculatePressureBreakpoints()
Calculate (or recalculate) the atmospheric pressure breakpoints at the altitudes in the standard temp...
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SetPressureSL()

void SetPressureSL ( ePressure  unit,
double  pressure 
)
virtual

Sets the sea level pressure for modeling an off-standard pressure profile.

This could be useful in the case where the pressure at an airfield is known or set for a particular simulation run.

Parameters
pressureThe pressure in the units specified.
unitthe unit of measure that the specified pressure is supplied in.

Reimplemented from FGAtmosphere.

Definition at line 187 of file FGStandardAtmosphere.cpp.

188 {
189  double press = ConvertToPSF(pressure, unit);
190 
191  PressureBreakpointVector[0] = press;
193 }
void CalculatePressureBreakpoints()
Calculate (or recalculate) the atmospheric pressure breakpoints at the altitudes in the standard temp...
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SetSLTemperatureGradedDelta()

void SetSLTemperatureGradedDelta ( eTemperature  unit,
double  t 
)
virtual

Sets a Sea Level temperature delta that is ramped out by 86 km.

The value of the delta is used to calculate a delta gradient that is applied to the temperature at all altitudes below 86 km (282152 ft). For instance, if a temperature of 20 degrees F is supplied, the delta gradient would be 20/282152 - or, about 7.09E-5 degrees/ft. At sea level, the full 20 degrees would be added to the standard temperature, but that 20 degree delta would be reduced by 7.09E-5 degrees for every foot of altitude above sea level, so that by 86 km, there would be no further delta added to the standard temperature. The graded delta can be used along with the a bias to tailor the temperature profile as desired.

Parameters
tthe sea level temperature delta value in the unit provided.
unitthe unit of the temperature.

Definition at line 329 of file FGStandardAtmosphere.cpp.

330 {
331  SetTemperatureGradedDelta(deltemp, 0.0, unit);
332 }
virtual void SetTemperatureGradedDelta(double t, double h, eTemperature unit=eFahrenheit)
Sets the temperature delta value at the supplied altitude/elevation above sea level, to be added to the standard temperature and ramped out by 86 km.
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SetTemperature()

void SetTemperature ( double  t,
double  h,
eTemperature  unit = eFahrenheit 
)
virtual

Sets the temperature at the supplied altitude, if it is to be different than the standard temperature.

This function will calculate a bias - a difference - from the standard atmosphere temperature at the supplied altitude and will apply that calculated bias to the entire temperature profile. If a graded delta is present, that will be included in the calculation - that is, regardless of any graded delta present, a temperature bias will be determined so that the temperature requested in this function call will be reached.

Parameters
tThe temperature value in the unit provided.
hThe altitude in feet above sea level.
unitThe unit of the temperature.

Implements FGAtmosphere.

Definition at line 295 of file FGStandardAtmosphere.cpp.

296 {
297  double targetSLtemp = ConvertToRankine(t, unit);
298 
299  TemperatureBias = 0.0;
300  TemperatureBias = targetSLtemp - GetTemperature(h);
302 }
void CalculatePressureBreakpoints()
Calculate (or recalculate) the atmospheric pressure breakpoints at the altitudes in the standard temp...
virtual double GetTemperature() const
Returns the actual, modeled temperature at the current altitude in degrees Rankine.
Definition: FGAtmosphere.h:117
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SetTemperatureBias()

void SetTemperatureBias ( eTemperature  unit,
double  t 
)
virtual

Sets the temperature bias to be added to the standard temperature at all altitudes.

This function sets the bias - the difference - from the standard atmosphere temperature. This bias applies to the entire temperature profile. Another way to set the temperature bias is to use the SetSLTemperature function, which replaces the value calculated by this function with a calculated bias.

Parameters
tthe temperature value in the unit provided.
unitthe unit of the temperature.

Definition at line 306 of file FGStandardAtmosphere.cpp.

307 {
308  if (unit == eCelsius || unit == eKelvin)
309  t *= 1.80; // If temp delta "t" is given in metric, scale up to English
310 
311  TemperatureBias = t;
313 }
void CalculatePressureBreakpoints()
Calculate (or recalculate) the atmospheric pressure breakpoints at the altitudes in the standard temp...
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SetTemperatureGradedDelta()

void SetTemperatureGradedDelta ( double  t,
double  h,
eTemperature  unit = eFahrenheit 
)
virtual

Sets the temperature delta value at the supplied altitude/elevation above sea level, to be added to the standard temperature and ramped out by 86 km.

This function computes the sea level delta from the standard atmosphere temperature at sea level.

Parameters
tthe temperature skew value in the unit provided.
unitthe unit of the temperature.

Definition at line 341 of file FGStandardAtmosphere.cpp.

342 {
343  if (unit == eCelsius || unit == eKelvin)
344  deltemp *= 1.80; // If temp delta "t" is given in metric, scale up to English
345 
346  TemperatureDeltaGradient = deltemp/(GradientFadeoutAltitude - h);
349 }
void CalculateLapseRates()
Recalculate the lapse rate vectors when the temperature profile is altered in a way that would change...
void CalculatePressureBreakpoints()
Calculate (or recalculate) the atmospheric pressure breakpoints at the altitudes in the standard temp...
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SetTemperatureSL()

void SetTemperatureSL ( double  t,
eTemperature  unit = eFahrenheit 
)
virtual

Sets the Sea Level temperature, if it is to be different than the standard.

This function will calculate a bias - a difference - from the standard atmosphere temperature and will apply that bias to the entire temperature profile. This is one way to set the temperature bias. Using the SetTemperatureBias function will replace the value calculated by this function.

Parameters
tthe temperature value in the unit provided.
unitthe unit of the temperature.

Reimplemented from FGAtmosphere.

Definition at line 321 of file FGStandardAtmosphere.cpp.

322 {
323  SetTemperature(t, 0.0, unit);
324 }
virtual void SetTemperature(double t, double h, eTemperature unit=eFahrenheit)
Sets the temperature at the supplied altitude, if it is to be different than the standard temperature...
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

The documentation for this class was generated from the following files: