JSBSim Flight Dynamics Model 1.2.2 (22 Mar 2025)
An Open Source Flight Dynamics and Control Software Library in C++
Loading...
Searching...
No Matches
FGMSIS Class Reference

Detailed Description

Models the MSIS-00 atmosphere.

This is a wrapper for the NRL-MSIS-00 model 2001:

This C++ format model uses the NRLMSISE-00 C source code package - release 20020503

The NRLMSISE-00 model was developed by Mike Picone, Alan Hedin, and Doug Drob. They also wrote a NRLMSISE-00 distribution package in FORTRAN which is available at http://uap-www.nrl.navy.mil/models_web/msis/msis_home.htm

Dominik Brodowski implemented and maintains this C version. You can reach him at devel.nosp@m.@bro.nosp@m.do.de. See the file "DOCUMENTATION" for details, and check http://www.brodo.de/english/pub/nrlmsise/index.html for updated releases of this package.

Author
David Culp

Definition at line 83 of file FGMSIS.h.

#include <FGMSIS.h>

+ Inheritance diagram for FGMSIS:
+ Collaboration diagram for FGMSIS:

Public Member Functions

 FGMSIS (FGFDMExec *)
 Constructor.
 
 ~FGMSIS ()
 Destructor.
 
virtual double GetDensity (double altitude) const
 Returns the density in slugs/ft^3 at a given altitude in ft.
 
double GetDensity (double altitude) const override
 Returns the density in slugs/ft^3 at a given altitude in ft.
 
virtual double GetDensity (void) const
 Returns the density in slugs/ft^3.
 
double GetPressure (double altitude) const override
 Returns the pressure at a specified altitude in psf.
 
virtual double GetPressure (double altitude) const=0
 Returns the pressure at a specified altitude in psf.
 
virtual double GetPressure (void) const
 Returns the pressure in psf.
 
virtual double GetSoundSpeed (double altitude) const
 Returns the speed of sound in ft/sec at a given altitude in ft.
 
double GetSoundSpeed (double altitude) const override
 Returns the speed of sound in ft/sec at a given altitude in ft.
 
virtual double GetSoundSpeed (void) const
 Returns the speed of sound in ft/sec.
 
virtual double GetTemperature () const
 Returns the actual, modeled temperature at the current altitude in degrees Rankine.
 
double GetTemperature (double altitude) const override
 Returns the actual modeled temperature in degrees Rankine at a specified altitude.
 
virtual double GetTemperature (double altitude) const=0
 Returns the actual modeled temperature in degrees Rankine at a specified altitude.
 
bool InitModel (void) override
 Runs the MSIS-00 atmosphere model; called by the Executive Can pass in a value indicating if the executive is directing the simulation to Hold.
 
bool Load (Element *el) override
 
- Public Member Functions inherited from FGStandardAtmosphere
 FGStandardAtmosphere (FGFDMExec *)
 Constructor.
 
virtual ~FGStandardAtmosphere ()
 Destructor.
 
bool InitModel (void) override
 
double GetTemperature (double altitude) const override
 Returns the actual modeled temperature in degrees Rankine at a specified altitude.
 
virtual double GetStdTemperature (double altitude) const
 Returns the standard temperature in degrees Rankine at a specified altitude.
 
virtual double GetStdTemperatureSL () const
 Returns the standard sea level temperature in degrees Rankine.
 
virtual double GetStdTemperatureRatio (double h) const
 Returns the ratio of the standard temperature at the supplied altitude over the standard sea level temperature.
 
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.
 
void SetTemperatureSL (double t, eTemperature unit=eFahrenheit) override
 Sets the Sea Level temperature, if it is to be different than the standard.
 
void SetTemperature (double t, double h, eTemperature unit=eFahrenheit) override
 Sets the temperature at the supplied altitude, if it is to be different than the standard temperature.
 
virtual void SetTemperatureBias (eTemperature unit, double t)
 Sets the temperature bias to be added to the standard temperature at all altitudes.
 
virtual void SetSLTemperatureGradedDelta (eTemperature unit, double t)
 Sets a Sea Level temperature delta that is ramped out by 86 km.
 
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.
 
virtual void ResetSLTemperature ()
 This function resets the model to apply no bias or delta gradient to the temperature.
 
double GetPressure (double altitude) const override
 Returns the pressure at a specified altitude in psf.
 
virtual double GetStdPressure (double altitude) const
 Returns the standard pressure at the specified altitude.
 
void SetPressureSL (ePressure unit, double pressure) override
 Sets the sea level pressure for modeling an off-standard pressure profile.
 
virtual void ResetSLPressure ()
 Resets the sea level to the Standard sea level pressure, and recalculates dependent parameters so that the pressure calculations are standard.
 
virtual double GetStdDensity (double altitude) const
 Returns the standard density at a specified altitude.
 
void SetDewPoint (eTemperature unit, double dewpoint)
 Sets the dew point.
 
double GetDewPoint (eTemperature to) const
 Returns the dew point.
 
void SetVaporPressure (ePressure unit, double Pv)
 Sets the partial pressure of water vapor.
 
double GetVaporPressure (ePressure to) const
 Returns the partial pressure of water vapor.
 
double GetSaturatedVaporPressure (ePressure to) const
 Returns the saturated pressure of water vapor.
 
void SetRelativeHumidity (double RH)
 Sets the relative humidity.
 
double GetRelativeHumidity (void) const
 Returns the relative humidity in percent.
 
void SetVaporMassFractionPPM (double frac)
 Sets the vapor mass per million of dry air mass units.
 
double GetVaporMassFractionPPM (void) const
 Returns the vapor mass per million of dry air mass units (ppm).
 
virtual void PrintStandardAtmosphereTable ()
 Prints the U.S. Standard Atmosphere table.
 
- Public Member Functions inherited from FGAtmosphere
 FGAtmosphere (FGFDMExec *)
 Constructor.
 
virtual ~FGAtmosphere ()
 Destructor.
 
bool Run (bool Holding) override
 Runs the atmosphere forces model; called by the Executive.
 
virtual double GetTemperatureSL () const
 Returns the actual, modeled sea level temperature in degrees Rankine.
 
virtual double GetTemperatureRatio () const
 Returns the ratio of the at-current-altitude temperature as modeled over the sea level value.
 
virtual double GetTemperatureRatio (double h) const
 Returns the ratio of the temperature as modeled at the supplied altitude over the sea level value.
 
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 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 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.
 
virtual double GetDensityAltitude () const
 
virtual double GetPressureAltitude () const
 
- Public Member Functions inherited from FGModel
 FGModel (FGFDMExec *)
 Constructor.
 
 ~FGModel () override
 Destructor.
 
virtual SGPath FindFullPathName (const SGPath &path) const
 
FGFDMExecGetExec (void)
 
const std::string & GetName (void)
 
unsigned int GetRate (void)
 Get the output rate for the model in frames.
 
void SetPropertyManager (std::shared_ptr< 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.
 
std::string GetFunctionValues (const std::string &delimeter) const
 Gets the function values.
 
std::shared_ptr< FGFunctionGetPreFunction (const std::string &name)
 Get one of the "pre" function.
 
bool Load (Element *el, FGFDMExec *fdmex, std::string prefix="")
 
void PostLoad (Element *el, FGFDMExec *fdmex, std::string prefix="")
 
void PreLoad (Element *el, FGFDMExec *fdmex, 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.
 

Protected Member Functions

void Calculate (double altitude) override
 Calculate the atmosphere for the given altitude.
 
void Compute (double altitude, double &pression, double &temperature, double &density, double &Rair) const
 
- Protected Member Functions inherited from FGStandardAtmosphere
void Calculate (double altitude) override
 Calculate the atmosphere for the given altitude.
 
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.
 
void CalculatePressureBreakpoints (double SLpress)
 Calculate (or recalculate) the atmospheric pressure breakpoints at the altitudes in the standard temperature table.
 
void CalculateStdDensityBreakpoints ()
 Calculate the atmospheric density breakpoints at the altitudes in the standard temperature table.
 
double GeopotentialAltitude (double geometalt) const
 Convert a geometric altitude to a geopotential altitude.
 
double GeometricAltitude (double geopotalt) const
 Convert a geopotential altitude to a geometric altitude.
 
double CalculateDensityAltitude (double density, double geometricAlt) override
 Calculates the density altitude given any temperature or pressure bias.
 
double CalculatePressureAltitude (double pressure, double geometricAlt) override
 Calculates the pressure altitude given any temperature or pressure bias.
 
double CalculateVaporPressure (double temperature)
 Calculate the pressure of water vapor with the Magnus formula.
 
void ValidateVaporMassFraction (double geometricAlt)
 Validate the value of the vapor mass fraction.
 
void CalculateSLDensity (void)
 Calculate the SL density.
 
void CalculateSLSoundSpeedAndDensity (void)
 Calculate the SL density and sound speed.
 
void bind (void) override
 
void Debug (int from) override
 
- Protected Member Functions inherited from FGAtmosphere
double ConvertToRankine (double t, eTemperature unit) const
 Converts to Rankine from one of several unit systems.
 
double ConvertFromRankine (double t, eTemperature unit) const
 Converts from Rankine to one of several unit systems.
 
double ConvertToPSF (double t, ePressure unit=ePSF) const
 Converts to PSF (pounds per square foot) from one of several unit systems.
 
double ConvertFromPSF (double t, ePressure unit=ePSF) const
 Converts from PSF (pounds per square foot) to one of several unit systems.
 
double ValidatePressure (double p, const std::string &msg, bool quiet=false) const
 Check that the pressure is within plausible boundaries.
 
double ValidateTemperature (double t, const std::string &msg, bool quiet=false) const
 Check that the temperature is within plausible boundaries.
 
- Protected Member Functions inherited from FGModel
bool Upload (Element *el, bool preLoad)
 Uploads this model in memory.
 

Protected Attributes

double day_of_year = 1.0
 
struct nrlmsise_flags flags
 
struct nrlmsise_input input
 
double seconds_in_day = 0.0
 
- Protected Attributes inherited from FGStandardAtmosphere
double StdSLtemperature
 Standard sea level conditions.
 
double StdSLdensity
 
double StdSLpressure
 
double StdSLsoundspeed
 
double TemperatureBias
 
double TemperatureDeltaGradient
 
double GradientFadeoutAltitude
 
double VaporMassFraction
 
double SaturatedVaporPressure
 
FGTable StdAtmosTemperatureTable
 
FGTable MaxVaporMassFraction
 
std::vector< double > LapseRates
 
std::vector< double > PressureBreakpoints
 
std::vector< double > StdPressureBreakpoints
 
std::vector< double > StdDensityBreakpoints
 
std::vector< double > StdLapseRates
 
- Protected Attributes inherited from FGAtmosphere
double SLtemperature = 1.8
 
double SLdensity = 1.0
 
double SLpressure = 1.0
 
double SLsoundspeed = 1.0
 
double Temperature = 1.8
 
double Density = 0.0
 
double Pressure = 0.0
 
double Soundspeed = 0.0
 
double PressureAltitude = 0.0
 
double DensityAltitude = 0.0
 
double Viscosity = 0.0
 
double KinematicViscosity = 0.0
 
double Reng = Reng0
 
- Protected Attributes inherited from FGModel
unsigned int exe_ctr
 
FGFDMExecFDMExec
 
std::string Name
 
std::shared_ptr< FGPropertyManagerPropertyManager
 
unsigned int rate
 
- Protected Attributes inherited from FGModelFunctions
FGPropertyReader LocalProperties
 
std::vector< std::shared_ptr< FGFunction > > PostFunctions
 
std::vector< std::shared_ptr< FGFunction > > PreFunctions
 

Additional Inherited Members

- Public Types inherited from FGAtmosphere
enum  ePressure {
  eNoPressUnit =0 , ePSF , eMillibars , ePascals ,
  eInchesHg
}
 Enums for specifying pressure units. More...
 
enum  eTemperature {
  eNoTempUnit =0 , eFahrenheit , eCelsius , eRankine ,
  eKelvin
}
 Enums for specifying temperature units. More...
 
- Public Types inherited from FGJSBBase
enum  { eL = 1 , eM , eN }
 Moments L, M, N. More...
 
enum  { eP = 1 , eQ , eR }
 Rates P, Q, R. More...
 
enum  { eU = 1 , eV , eW }
 Velocities U, V, W. More...
 
enum  { eX = 1 , eY , eZ }
 Positions X, Y, Z. More...
 
enum  { ePhi = 1 , eTht , ePsi }
 Euler angles Phi, Theta, Psi. More...
 
enum  { eDrag = 1 , eSide , eLift }
 Stability axis forces, Drag, Side force, Lift. More...
 
enum  { eRoll = 1 , ePitch , eYaw }
 Local frame orientation Roll, Pitch, Yaw. More...
 
enum  { eNorth = 1 , eEast , eDown }
 Local frame position North, East, Down. More...
 
enum  { eLat = 1 , eLong , eRad }
 Locations Radius, Latitude, Longitude. More...
 
enum  {
  inNone = 0 , inDegrees , inRadians , inMeters ,
  inFeet
}
 Conversion specifiers. More...
 
- Static Public Member Functions inherited from FGJSBBase
static const std::string & GetVersion (void)
 Returns the version number of JSBSim.
 
static constexpr double KelvinToFahrenheit (double kelvin)
 Converts from degrees Kelvin to degrees Fahrenheit.
 
static constexpr double CelsiusToRankine (double celsius)
 Converts from degrees Celsius to degrees Rankine.
 
static constexpr double RankineToCelsius (double rankine)
 Converts from degrees Rankine to degrees Celsius.
 
static constexpr double KelvinToRankine (double kelvin)
 Converts from degrees Kelvin to degrees Rankine.
 
static constexpr double RankineToKelvin (double rankine)
 Converts from degrees Rankine to degrees Kelvin.
 
static constexpr double FahrenheitToCelsius (double fahrenheit)
 Converts from degrees Fahrenheit to degrees Celsius.
 
static constexpr double CelsiusToFahrenheit (double celsius)
 Converts from degrees Celsius to degrees Fahrenheit.
 
static constexpr double CelsiusToKelvin (double celsius)
 Converts from degrees Celsius to degrees Kelvin.
 
static constexpr double KelvinToCelsius (double kelvin)
 Converts from degrees Kelvin to degrees Celsius.
 
static constexpr double FeetToMeters (double measure)
 Converts from feet to meters.
 
static bool EqualToRoundoff (double a, double b)
 Finite precision comparison.
 
static bool EqualToRoundoff (float a, float b)
 Finite precision comparison.
 
static bool EqualToRoundoff (float a, double b)
 Finite precision comparison.
 
static bool EqualToRoundoff (double a, float b)
 Finite precision comparison.
 
static constexpr double Constrain (double min, double value, double max)
 Constrain a value between a minimum and a maximum value.
 
static constexpr double sign (double num)
 
- Public Attributes inherited from FGAtmosphere
struct JSBSim::FGAtmosphere::Inputs in
 
const double StdDaySLsoundspeed
 
- Static Public Attributes inherited from FGAtmosphere
static constexpr double StdDaySLtemperature = 518.67
 
static constexpr double StdDaySLpressure = 2116.228
 
static constexpr double SHRatio = 1.4
 
static constexpr double StdDaySLdensity = StdDaySLpressure / (Reng0 * StdDaySLtemperature)
 
- Static Public Attributes inherited from FGJSBBase
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 short debug_lvl = 1
 
- Static Protected Member Functions inherited from FGJSBBase
static std::string CreateIndexedPropertyName (const std::string &Property, int index)
 
- Static Protected Attributes inherited from FGStandardAtmosphere
static constexpr double EarthRadius = 6356766.0 / fttom
 Earth radius in ft as defined for ISA 1976.
 
static constexpr double a = 611.2/psftopa
 Sonntag constants based on ref [2].
 
static constexpr double b = 17.62
 
static constexpr double c = 243.12
 
static constexpr double Mwater = 18.016 * kgtoslug / 1000.0
 Mean molecular weight for water - slug/mol.
 
static constexpr double Rdry = Rstar / Mair
 
static constexpr double Rwater = Rstar / Mwater
 
- Static Protected Attributes inherited from FGAtmosphere
static constexpr double SutherlandConstant = 198.72
 
static constexpr double Beta = 2.269690E-08
 
static constexpr double Rstar = 8.31432 * kgtoslug / KelvinToRankine(fttom * fttom)
 Universal gas constant - ft*lbf/R/mol.
 
static constexpr double Mair = 28.9645 * kgtoslug / 1000.0
 Mean molecular weight for air - slug/mol.
 
static constexpr double g0 = 9.80665 / fttom
 Sea-level acceleration of gravity - ft/s^2.
 
static constexpr double Reng0 = Rstar / Mair
 Specific gas constant for air - ft*lbf/slug/R.
 
- Static Protected Attributes inherited from FGJSBBase
static constexpr double radtodeg = 180. / M_PI
 
static constexpr double degtorad = M_PI / 180.
 
static constexpr double hptoftlbssec = 550.0
 
static constexpr double psftoinhg = 0.014138
 
static constexpr double psftopa = 47.88
 
static constexpr double fttom = 0.3048
 
static constexpr double ktstofps = 1852./(3600*fttom)
 
static constexpr double fpstokts = 1.0 / ktstofps
 
static constexpr double inchtoft = 1.0/12.0
 
static constexpr double m3toft3 = 1.0/(fttom*fttom*fttom)
 
static constexpr double in3tom3 = inchtoft*inchtoft*inchtoft/m3toft3
 
static constexpr double inhgtopa = 3386.38
 
static constexpr double slugtolb = 32.174049
 Note that definition of lbtoslug by the inverse of slugtolb and not to a different constant you can also get from some tables will make lbtoslug*slugtolb == 1 up to the magnitude of roundoff.
 
static constexpr double lbtoslug = 1.0/slugtolb
 
static constexpr double kgtolb = 2.20462
 
static constexpr double kgtoslug = 0.06852168
 
static const std::string needed_cfg_version = "2.0"
 
static const std::string JSBSim_version = JSBSIM_VERSION " " __DATE__ " " __TIME__
 

Constructor & Destructor Documentation

◆ FGMSIS()

FGMSIS ( FGFDMExec fdmex)

Constructor.

Definition at line 73 of file FGMSIS.cpp.

74{
75 Name = "MSIS";
76
77 for(unsigned int i=0; i<24; ++i)
78 flags.switches[i] = 1;
79 input.year = 0; // Ignored by NRLMSIS
80 input.f107A = 150.;
81 input.f107 = 150.;
82 input.ap = 4.;
83 input.ap_a = nullptr;
84
85 Debug(0);
86}
FGStandardAtmosphere(FGFDMExec *)
Constructor.

◆ ~FGMSIS()

~FGMSIS ( )

Destructor.

Definition at line 90 of file FGMSIS.cpp.

91{
92 Debug(1);
93}

Member Function Documentation

◆ Calculate()

void Calculate ( double  altitude)
overrideprotectedvirtual

Calculate the atmosphere for the given altitude.

Reimplemented from FGAtmosphere.

Definition at line 124 of file FGMSIS.cpp.

125{
126 double SLRair = 0.0;
127
128 Compute(0.0, SLpressure, SLtemperature, SLdensity, SLRair);
129 Compute(altitude, Pressure, Temperature, Density, Reng);
130
131 SLsoundspeed = sqrt(SHRatio*SLRair*SLtemperature);
132 Soundspeed = sqrt(SHRatio*Reng*Temperature);
133 PressureAltitude = CalculatePressureAltitude(Pressure, altitude);
134 DensityAltitude = CalculateDensityAltitude(Density, altitude);
135
136 Viscosity = Beta * pow(Temperature, 1.5) / (SutherlandConstant + Temperature);
137 KinematicViscosity = Viscosity / Density;
138 // SaturatedVaporPressure = CalculateVaporPressure(Temperature);
139 // ValidateVaporMassFraction(altitude);
140}
double CalculateDensityAltitude(double density, double geometricAlt) override
Calculates the density altitude given any temperature or pressure bias.
double CalculatePressureAltitude(double pressure, double geometricAlt) override
Calculates the pressure altitude given any temperature or pressure bias.
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ Compute()

void Compute ( double  altitude,
double &  pression,
double &  temperature,
double &  density,
double &  Rair 
) const
protected

Definition at line 144 of file FGMSIS.cpp.

146{
147 constexpr double fttokm = fttom / 1000.;
148 constexpr double kgm3_to_slugft3 = kgtoslug / m3toft3;
149 constexpr double gtoslug = kgtoslug / 1000.;
150 // Molecular weight (g/mol)
151 // N2 O2 O He H Ar N OA
152 const double species_mmol[8] {28.0134, 31.9988, 31.9988/2.0, 4.0, 1.0, 39.948,
153 28.0134/2.0, 31.9988/2.0};
154
155 double dn[10] {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0};
156 double h = altitude*fttokm;
157 double lat = in.GeodLatitudeDeg;
158 double lon = in.LongitudeDeg;
159
160 // Compute epoch
161 double utc_seconds = seconds_in_day + FDMExec->GetSimTime();
162 unsigned int days = utc_seconds / 86400.;
163 utc_seconds -= days * 86400.;
164 double today = day_of_year + days;
165 unsigned int year = today / 365.;
166 today -= year * 365.;
167
168 struct nrlmsise_output output;
169
170 input.doy = today;
171 input.sec = utc_seconds;
172 input.alt = h;
173 input.g_lat = lat;
174 input.g_long = lon;
175 input.lst = utc_seconds/3600 + lon/15; // Local Solar Time (hours)
176 assert(flags.switches[9] != -1); // Make sure that input.ap is used.
177
178 gtd7(&input, &flags, &output);
179
180 temperature = KelvinToRankine(output.t[1]);
181 density = output.d[5] * kgm3_to_slugft3;
182 dn[1] = output.d[2]; // N2
183 dn[2] = output.d[3]; // O2
184 dn[3] = output.d[1]; // O
185 dn[4] = output.d[0]; // He
186 dn[5] = output.d[6]; // H
187 dn[6] = output.d[4]; // Ar
188 dn[7] = output.d[7]; // N
189 // SUBROUTINE GTD7 does NOT include anomalous oxygen so we drop it from
190 // the molar mass computation as well for consistency.
191 dn[8] = 0.0; // OA
192
193 // Compute specific gas constant for air
194 double mmol = 0.0;
195 double qty_mol = 0.0;
196
197 for (unsigned int i=1; i<9; ++i) {
198 mmol += dn[i]*species_mmol[i-1];
199 qty_mol += dn[i];
200 }
201 double mair = mmol * gtoslug / qty_mol;
202 Rair = Rstar / mair;
203
204 pressure = density * Rair * temperature;
205}
static constexpr double Rstar
Universal gas constant - ft*lbf/R/mol.
double GetSimTime(void) const
Returns the cumulative simulation time in seconds.
Definition FGFDMExec.h:549
static constexpr double KelvinToRankine(double kelvin)
Converts from degrees Kelvin to degrees Rankine.
Definition FGJSBBase.h:206

◆ GetDensity() [1/3]

double GetDensity ( double  altitude) const
virtual

Returns the density in slugs/ft^3 at a given altitude in ft.

Reimplemented from FGAtmosphere.

Definition at line 170 of file FGAtmosphere.cpp.

183{
184 return GetPressure(altitude)/(Reng * GetTemperature(altitude));
185}
virtual double GetPressure(void) const
Returns the pressure in psf.
virtual double GetTemperature() const
Returns the actual, modeled temperature at the current altitude in degrees Rankine.

◆ GetDensity() [2/3]

double GetDensity ( double  altitude) const
inlineoverridevirtual

Returns the density in slugs/ft^3 at a given altitude in ft.

Reimplemented from FGAtmosphere.

Definition at line 117 of file FGMSIS.h.

117 {
118 double t, p, rho, R;
119 Compute(altitude, p, t, rho, R);
120 return rho;
121 }

◆ GetDensity() [3/3]

virtual double GetDensity ( void  ) const
inlinevirtual

Returns the density in slugs/ft^3.

This function may only be used if Run() is called first.

Reimplemented from FGAtmosphere.

Definition at line 167 of file FGAtmosphere.h.

167{return Density;}

◆ GetPressure() [1/3]

double GetPressure ( double  altitude) const
inlineoverridevirtual

Returns the pressure at a specified altitude in psf.

Implements FGAtmosphere.

Definition at line 110 of file FGMSIS.h.

110 {
111 double t, p, rho, R;
112 Compute(altitude, p, t, rho, R);
113 return p;
114 }

◆ GetPressure() [2/3]

virtual double GetPressure ( double  altitude) const
virtual

Returns the pressure at a specified altitude in psf.

Implements FGAtmosphere.

◆ GetPressure() [3/3]

virtual double GetPressure ( void  ) const
inlinevirtual

Returns the pressure in psf.

Reimplemented from FGAtmosphere.

Definition at line 144 of file FGAtmosphere.h.

144{return Pressure;}

◆ GetSoundSpeed() [1/3]

double GetSoundSpeed ( double  altitude) const
virtual

Returns the speed of sound in ft/sec at a given altitude in ft.

Reimplemented from FGAtmosphere.

Definition at line 186 of file FGAtmosphere.cpp.

190{
191 return sqrt(SHRatio * Reng * GetTemperature(altitude));
192}

◆ GetSoundSpeed() [2/3]

double GetSoundSpeed ( double  altitude) const
inlineoverridevirtual

Returns the speed of sound in ft/sec at a given altitude in ft.

Reimplemented from FGAtmosphere.

Definition at line 124 of file FGMSIS.h.

124 {
125 double t, p, rho, R;
126 Compute(altitude, p, t, rho, R);
127 return sqrt(FGAtmosphere::SHRatio*R*t);
128 }

◆ GetSoundSpeed() [3/3]

virtual double GetSoundSpeed ( void  ) const
inlinevirtual

Returns the speed of sound in ft/sec.

Reimplemented from FGAtmosphere.

Definition at line 183 of file FGAtmosphere.h.

183{return Soundspeed;}

◆ GetTemperature() [1/3]

virtual double GetTemperature ( ) const
inlinevirtual

Returns the actual, modeled temperature at the current altitude in degrees Rankine.

Returns
Modeled temperature in degrees Rankine.

Reimplemented from FGAtmosphere.

Definition at line 109 of file FGAtmosphere.h.

109{return Temperature;}

◆ GetTemperature() [2/3]

double GetTemperature ( double  altitude) const
inlineoverridevirtual

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 103 of file FGMSIS.h.

103 {
104 double t, p, rho, R;
105 Compute(altitude, p, t, rho, R);
106 return t;
107 }

◆ GetTemperature() [3/3]

virtual 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.

◆ InitModel()

bool InitModel ( void  )
overridevirtual

Runs the MSIS-00 atmosphere model; called by the Executive Can pass in a value indicating if the executive is directing the simulation to Hold.

Parameters
Holdingif true, the executive has been directed to hold the sim from advancing time. Some models may ignore this flag, such as the Input model, which may need to be active to listen on a socket for the "Resume" command to be given.
Returns
false if no error

Reimplemented from FGAtmosphere.

Definition at line 97 of file FGMSIS.cpp.

98{
99 FGStandardAtmosphere::InitModel();
100
101 Calculate(0.0);
102
103 return true;
104}
void Calculate(double altitude) override
Calculate the atmosphere for the given altitude.
Definition FGMSIS.cpp:124
+ Here is the call graph for this function:

◆ Load()

bool Load ( Element el)
overridevirtual

Reimplemented from FGModel.

Definition at line 108 of file FGMSIS.cpp.

109{
110 if (!Upload(el, true)) return false;
111
112 if (el->FindElement("day"))
113 day_of_year = el->FindElementValueAsNumber("day");
114 if (el->FindElement("utc"))
115 seconds_in_day = el->FindElementValueAsNumber("utc");
116
117 Debug(3);
118
119 return true;
120}
bool Upload(Element *el, bool preLoad)
Uploads this model in memory.
Definition FGModel.cpp:110

Member Data Documentation

◆ day_of_year

double day_of_year = 1.0
protected

Definition at line 135 of file FGMSIS.h.

◆ flags

struct nrlmsise_flags flags
protected

Definition at line 138 of file FGMSIS.h.

◆ input

struct nrlmsise_input input
protected

Definition at line 139 of file FGMSIS.h.

◆ seconds_in_day

double seconds_in_day = 0.0
protected

Definition at line 136 of file FGMSIS.h.


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