nrlmsise00

Python/NumPy interface

msise_model(time, alt, lat, lon, f107a, f107, ap)

Interface to gtd7() [1] and gtd7d() [2]

msise_flat(*args, **kwargs)

Flattened interface to gtd7() [1] and gtd7d() [2]

gtd7_flat(*args, **kwargs)

Flattened variant of the MSIS gtd7() function

gtd7d_flat(*args, **kwargs)

Flattened variant of the MSIS gtd7d() function

scale_height(alt, lat, molw, temp)

Atmospheric scale height

Python interface based on the NRLMSISE-00 C version [0]

nrlmsise00.gtd7_flat(*args, **kwargs)

Flattened variant of the MSIS gtd7() function

Returns a single 11-element numpy.ndarray instead of the two lists. All arguments except the keywords flags and ap_a can be numpy.ndarray to facilitate calculations at many locations/times.

gtd7(year, doy, sec, alt, g_lat, g_long, lst, f107A, f107, ap, ap_a=None, flags=None)

MSIS Neutral Atmosphere Empircial Model from the surface to lower exosphere.

Parameters:
  • year (int) – Year, but has no real effect, more important is doy.

  • doy (int) – Day of the year.

  • sec (float) – Seconds into the day (UT).

  • alt (float) – Altitude in [km].

  • g_lat (float) – Geodetic latitude in [degrees N].

  • g_long (float) – Geodetic longitude in [degrees E].

  • lst (float) – Apparent local solar timei [h].

  • f107A (float) – 81 day average of 10.7 cm radio flux (centered on doy) at the actual distance of the Earth from the Sun rather than the radio flux at 1 AU.

  • f107 (float) – Daily F10.7 flux for previous day at position of Earth. Like f107A at the actual distance of the Earth from the Sun rather than the radio flux at 1 AU.

  • ap (float) – Daily geomagnetic ap index.

  • ap_a (list of 7 floats, optional) –

    Array containing the following magnetic values:

    1. daily AP

    2. 3 hr AP index for current time

    3. 3 hr AP index for 3 hrs before current time

    4. 3 hr AP index for 6 hrs before current time

    5. 3 hr AP index for 9 hrs before current time

    6. Average of eight 3 hr AP indicies from 12 to 33 hrs

      prior to current time

    7. Average of eight 3 hr AP indicies from 36 to 57 hrs

      prior to current time

  • flags (list of 24 int, optional) –

    Sets the model’s internal switches array. Quote from the NRLMSISE-00 source code: Switches: to turn on and off particular variations use these switches. 0 is off, 1 is on, and 2 is main effects off but cross terms on.

    Standard values are 0 for switch 0 and 1 for switches 1 to 23. The array ‘switches’ needs to be set accordingly by the calling program. The arrays sw and swc are set internally.

    switches[i]:

    1. output in meters and kilograms instead of centimetres and grams

    2. F10.7 effect on mean

    3. time independent

    4. symmetrical annual

    5. symmetrical semiannual

    6. asymmetrical annual

    7. asymmetrical semiannual

    8. diurnal

    9. semidiurnal

    10. daily ap

      [when this is set to -1 (!) the pointer ap_a in struct nrlmsise_input must point to a struct ap_array]

    11. all UT/long effects

    12. longitudinal

    13. UT and mixed UT/long

    14. mixed AP/UT/LONG

    15. terdiurnal

    16. departures from diffusive equilibrium

    17. all TINF var

    18. all TLB var

    19. all TN1 var

    20. all S var

    21. all TN2 var

    22. all NLB var

    23. all TN3 var

    24. turbo scale height var

Returns:

  • densities (list) – the NRLMSISE-00 densities:

    • d[0] - HE NUMBER DENSITY(CM-3)

    • d[1] - O NUMBER DENSITY(CM-3)

    • d[2] - N2 NUMBER DENSITY(CM-3)

    • d[3] - O2 NUMBER DENSITY(CM-3)

    • d[4] - AR NUMBER DENSITY(CM-3)

    • d[5] - TOTAL MASS DENSITY(GM/CM3) [includes d[8] in td7d]

    • d[6] - H NUMBER DENSITY(CM-3)

    • d[7] - N NUMBER DENSITY(CM-3)

    • d[8] - Anomalous oxygen NUMBER DENSITY(CM-3)

    O, H, and N are set to zero below 72.5 km

    d[5], TOTAL MASS DENSITY, is NOT the same for subroutines GTD7 and GTD7D SUBROUTINE GTD7 – d[5] is the sum of the mass densities of the species labeled by indices 0-4 and 6-7 in output variable d. This includes He, O, N2, O2, Ar, H, and N but does NOT include anomalous oxygen (species index 8).

  • temperatures (list) – the NRLMSISE-00 temperatures:

    • t[0] - EXOSPHERIC TEMPERATURE

    • t[1] - TEMPERATURE AT ALT

    t[0], Exospheric temperature, is set to global average for altitudes below 120 km. The 120 km gradient is left at global average value for altitudes below 72 km.

nrlmsise00.gtd7d_flat(*args, **kwargs)

Flattened variant of the MSIS gtd7d() function

Returns a single 11-element numpy.ndarray instead of the two lists. All arguments except the keywords flags and ap_a can be numpy.ndarray to facilitate calculations at many locations/times.

gtd7d(*args, **kwargs)

MSIS Neutral Atmosphere Empircial Model from the surface to lower exosphere.

This subroutine provides Effective Total Mass Density for output d[5] which includes contributions from ‘anomalous oxygen’ which can affect satellite drag above 500 km. See ‘returns’ for additional details.

Parameters:
  • *args – Same as for gtd7().

  • **kwargs – Same as for gtd7().

Returns:

densities, temperatures – See documentation for gtd7(), except for d[5]:

SUBROUTINE GTD7D – d[5] is the ‘effective total mass density for drag’ and is the sum of the mass densities of all species in this model, INCLUDING anomalous oxygen.

Return type:

(list, list)

nrlmsise00.msise_flat(*args, **kwargs)[source]

Flattened interface to gtd7() [1] and gtd7d() [2]

Calls the C model function using a datetime.datetime instance to calculate the day of year, seconds and so on. The method keyword decides which version to call, the difference being that gtd7d() includes anomalous oxygen in the total mass density (d[5]), gtd7() does not.

Parameters:
  • time (datetime.datetime) – Date and time as a datetime.dateime.

  • alt (float) – Altitude in km.

  • lat (float) – Latitude in degrees north.

  • lon (float) – Longitude in degrees east.

  • f107a (float) – The observed f107a (81-day running mean of f107) centred at date.

  • f107 (float) – The observed f107 value on the previous day.

  • ap (float) – The ap value at date.

  • lst (float, optional) – The local solar time, can be different from the calculated one.

  • ap_a (list, optional) – List of length 7 containing ap values to be used when flags[9] is set to -1, otherwise no effect.

  • flags (list, optional) – List of length 24 setting the NRLMSIS switches explicitly.

  • method (string, optional) – Set to “gtd7d” to use gtd7d() (which includes anomalous oxygen in the total mass density) instead of the “standard” gtd7() function without it.

Returns:

  • densities (list of floats) –

    1. He number density [cm^-3]

    2. O number density [cm^-3]

    3. N2 number density [cm^-3]

    4. O2 number density [cm^-3]

    5. AR number density [cm^-3]

    6. total mass density [g cm^-3] (includes d[8] in gtd7d)

    7. H number density [cm^-3]

    8. N number density [cm^-3]

    9. Anomalous oxygen number density [cm^-3]

  • temperatures (list of floats) –

    1. Exospheric temperature [K]

    2. Temperature at alt [K]

Note

No date and time conversion will be attempted, the input needs to be converted before, dateutil or pandas provide convenient functions.

The local solar time is calculated from time and longitude, except when lst is set, then that value is used.

The solar and geomagnetic indices have to be provided, so far the values are not included in the module.

nrlmsise00.msise_model(time, alt, lat, lon, f107a, f107, ap, lst=None, ap_a=None, flags=None, method='gtd7')[source]

Interface to gtd7() [1] and gtd7d() [2]

Calls the C model function using a datetime.datetime instance to calculate the day of year, seconds and so on. The method keyword decides which version to call, the difference being that gtd7d() includes anomalous oxygen in the total mass density (d[5]), gtd7() does not.

Parameters:
  • time (datetime.datetime) – Date and time as a datetime.dateime.

  • alt (float) – Altitude in km.

  • lat (float) – Latitude in degrees north.

  • lon (float) – Longitude in degrees east.

  • f107a (float) – The observed f107a (81-day running mean of f107) centred at date.

  • f107 (float) – The observed f107 value on the previous day.

  • ap (float) – The ap value at date.

  • lst (float, optional) – The local solar time, can be different from the calculated one.

  • ap_a (list, optional) – List of length 7 containing ap values to be used when flags[9] is set to -1, otherwise no effect.

  • flags (list, optional) – List of length 24 setting the NRLMSIS switches explicitly.

  • method (string, optional) – Set to “gtd7d” to use gtd7d() (which includes anomalous oxygen in the total mass density) instead of the “standard” gtd7() function without it.

Returns:

  • densities (list of floats) –

    1. He number density [cm^-3]

    2. O number density [cm^-3]

    3. N2 number density [cm^-3]

    4. O2 number density [cm^-3]

    5. AR number density [cm^-3]

    6. total mass density [g cm^-3] (includes d[8] in gtd7d)

    7. H number density [cm^-3]

    8. N number density [cm^-3]

    9. Anomalous oxygen number density [cm^-3]

  • temperatures (list of floats) –

    1. Exospheric temperature [K]

    2. Temperature at alt [K]

Note

No date and time conversion will be attempted, the input needs to be converted before, dateutil or pandas provide convenient functions.

The local solar time is calculated from time and longitude, except when lst is set, then that value is used.

The solar and geomagnetic indices have to be provided, so far the values are not included in the module.

nrlmsise00.scale_height(alt, lat, molw, temp)[source]

Atmospheric scale height

Extracted from the C-code for easy access, with constants updated to the standard SI values. It is reasonably fast and numpy-broadcasting compatible.

Parameters:
  • alt (float or array_like) – Altitude in [km].

  • lat (float or array_like) – Geodetic latitude in [degrees N].

  • molw (float or array_like) – Molecular mass at alt in [kg / mol]. Can be derived by dividing the total mass density by the sum of number densities, e.g. from gtd7(), multiplied by Avogadro’s constant (6.0221407e23 / mol).

  • temp (float or array_like) – Temperature at alt in [K].

Returns:

scale_height – Scale height in [m].

Return type:

float or array_like

4-D Xarray interface