COESA_module Module

1976 Standard Atmosphere Model

History

  • Original Julia code from: https://github.com/danielmatz/COESA.jl [MIT License]
  • Jacob Williams, converted to Fortran 4/30/2020

Uses

  • module~~coesa_module~~UsesGraph module~coesa_module COESA_module iso_fortran_env iso_fortran_env module~coesa_module->iso_fortran_env

Contents


Variables

TypeVisibilityAttributesNameInitial
real(kind=wp), private, parameter:: r0 =6356766.0_wp
real(kind=wp), private, parameter:: g0 =9.80665_wp
real(kind=wp), private, parameter:: M0 =28.9644_wp
real(kind=wp), private, parameter:: Rstar =8.31432e3_wp
real(kind=wp), private, parameter:: gamma =1.4_wp
real(kind=wp), private, parameter, dimension(*):: Hb =[0, 11, 20, 32, 47, 51, 71]*1000.0_wp
real(kind=wp), private, parameter, dimension(*):: Lmb =[-6.5_wp, 0.0_wp, 1.0_wp, 2.8_wp, 0.0_wp, -2.8_wp, -2.0_wp]/1000.0_wp
real(kind=wp), private, parameter, dimension(*):: Ztable =[80.0_wp, 80.5_wp, 81.0_wp, 81.5_wp, 82.0_wp, 82.5_wp, 83.0_wp, 83.5_wp, 84.0_wp, 84.5_wp, 85.0_wp, 85.5_wp, 86.0_wp]*1000
real(kind=wp), private, parameter, dimension(*):: Mratiotable =[1.0_wp, 0.999996_wp, 0.999989_wp, 0.999971_wp, 0.999941_wp, 0.999909_wp, 0.999870_wp, 0.999829_wp, 0.999786_wp, 0.999741_wp, 0.999694_wp, 0.999641_wp, 0.999579_wp]
real(kind=wp), private, parameter, dimension(*):: Ztableupper =[86000.0_wp, 87000.0_wp, 88000.0_wp, 89000.0_wp, 90000.0_wp, 91000.0_wp, 93000.0_wp, 95000.0_wp, 97000.0_wp, 99000.0_wp, 101000.0_wp, 103000.0_wp, 105000.0_wp, 107000.0_wp, 109000.0_wp, 110000.0_wp, 111000.0_wp, 112000.0_wp, 113000.0_wp, 114000.0_wp, 115000.0_wp, 116000.0_wp, 117000.0_wp, 118000.0_wp, 119000.0_wp, 120000.0_wp, 125000.0_wp, 130000.0_wp, 135000.0_wp, 140000.0_wp, 145000.0_wp, 150000.0_wp, 160000.0_wp, 170000.0_wp, 180000.0_wp, 190000.0_wp, 200000.0_wp, 210000.0_wp, 220000.0_wp, 230000.0_wp, 240000.0_wp, 250000.0_wp, 260000.0_wp, 270000.0_wp, 280000.0_wp, 290000.0_wp, 300000.0_wp, 310000.0_wp, 320000.0_wp, 330000.0_wp, 340000.0_wp, 350000.0_wp, 360000.0_wp, 370000.0_wp, 380000.0_wp, 390000.0_wp, 400000.0_wp, 410000.0_wp, 420000.0_wp, 430000.0_wp, 440000.0_wp, 450000.0_wp, 460000.0_wp, 470000.0_wp, 480000.0_wp, 490000.0_wp, 500000.0_wp, 525000.0_wp, 550000.0_wp, 575000.0_wp, 600000.0_wp, 625000.0_wp, 650000.0_wp, 675000.0_wp, 700000.0_wp, 725000.0_wp, 750000.0_wp, 775000.0_wp, 800000.0_wp, 825000.0_wp, 850000.0_wp, 875000.0_wp, 900000.0_wp, 925000.0_wp, 950000.0_wp, 975000.0_wp, 1000000.0_wp]
real(kind=wp), private, parameter, dimension(*):: Ptableupper =[3.7338E-1_wp, 3.1259E-1_wp, 2.6173E-1_wp, 2.1919E-1_wp, 1.8359E-1_wp, 1.5381E-1_wp, 1.0801E-1_wp, 7.5966E-2_wp, 5.3571E-2_wp, 3.7948E-2_wp, 2.7192E-2_wp, 1.9742E-2_wp, 1.4477E-2_wp, 1.0751E-2_wp, 8.1142E-3_wp, 7.1042E-3_wp, 6.2614E-3_wp, 5.5547E-3_wp, 4.9570E-3_wp, 4.4473E-3_wp, 4.0096E-3_wp, 3.6312E-3_wp, 3.3022E-3_wp, 3.0144E-3_wp, 2.7615E-3_wp, 2.5382E-3_wp, 1.7354E-3_wp, 1.2505E-3_wp, 9.3568E-4_wp, 7.2028E-4_wp, 5.6691E-4_wp, 4.5422E-4_wp, 3.0395E-4_wp, 2.1210E-4_wp, 1.5271E-4_wp, 1.1266E-4_wp, 8.4736E-5_wp, 6.4756E-5_wp, 5.0149E-5_wp, 3.9276E-5_wp, 3.1059E-5_wp, 2.4767E-5_wp, 1.9894E-5_wp, 1.6083E-5_wp, 1.3076E-5_wp, 1.0683E-5_wp, 8.7704E-6_wp, 7.2285E-6_wp, 5.9796E-6_wp, 4.9630E-6_wp, 4.1320E-6_wp, 3.4498E-6_wp, 2.8878E-6_wp, 2.4234E-6_wp, 2.0384E-6_wp, 1.7184E-6_wp, 1.4518E-6_wp, 1.2291E-6_wp, 1.0427E-6_wp, 8.8645E-7_wp, 7.5517E-7_wp, 6.4468E-7_wp, 5.5155E-7_wp, 4.7292E-7_wp, 4.0642E-7_wp, 3.5011E-7_wp, 3.0236E-7_wp, 2.1200E-7_wp, 1.5137E-7_wp, 1.1028E-7_wp, 8.2130E-8_wp, 6.2601E-8_wp, 4.8865E-8_wp, 3.9048E-8_wp, 3.1908E-8_wp, 2.6611E-8_wp, 2.2599E-8_wp, 1.9493E-8_wp, 1.7036E-8_wp, 1.5051E-8_wp, 1.3415E-8_wp, 1.2043E-8_wp, 1.0873E-8_wp, 9.8635E-9_wp, 8.9816E-9_wp, 8.2043E-9_wp, 7.5138E-9_wp]
real(kind=wp), private, parameter, dimension(*):: logPtableupper =log(Ptableupper)
real(kind=wp), private, parameter, dimension(*):: Mtableupper =[28.95_wp, 28.95_wp, 28.94_wp, 28.93_wp, 28.91_wp, 28.89_wp, 28.82_wp, 28.73_wp, 28.62_wp, 28.48_wp, 28.30_wp, 28.10_wp, 27.88_wp, 27.64_wp, 27.39_wp, 27.27_wp, 27.14_wp, 27.02_wp, 26.90_wp, 26.79_wp, 26.68_wp, 26.58_wp, 26.48_wp, 26.38_wp, 26.29_wp, 26.20_wp, 25.80_wp, 25.44_wp, 25.09_wp, 24.75_wp, 24.42_wp, 24.10_wp, 23.49_wp, 22.90_wp, 22.34_wp, 21.81_wp, 21.30_wp, 20.83_wp, 20.37_wp, 19.95_wp, 19.56_wp, 19.19_wp, 18.85_wp, 18.53_wp, 18.24_wp, 17.97_wp, 17.73_wp, 17.50_wp, 17.29_wp, 17.09_wp, 16.91_wp, 16.74_wp, 16.57_wp, 16.42_wp, 16.27_wp, 16.13_wp, 15.98_wp, 15.84_wp, 15.70_wp, 15.55_wp, 15.40_wp, 15.25_wp, 15.08_wp, 14.91_wp, 14.73_wp, 14.54_wp, 14.33_wp, 13.76_wp, 13.09_wp, 12.34_wp, 11.51_wp, 10.62_wp, 9.72_wp, 8.83_wp, 8.00_wp, 7.24_wp, 6.58_wp, 6.01_wp, 5.54_wp, 5.16_wp, 4.85_wp, 4.60_wp, 4.40_wp, 4.25_wp, 4.12_wp, 4.02_wp, 3.94_wp]
real(kind=wp), private, parameter, dimension(*):: Tmb =[0.288150000000000000E+003, 0.216650000000000000E+003, 0.216650000000000000E+003, 0.228650000000000000E+003, 0.270650000000000000E+003, 0.270650000000000000E+003, 0.214650000000000000E+003]
real(kind=wp), private, parameter, dimension(*):: Pb =[0.101325000000000000E+006, 0.226320639734629302E+005, 0.547488866967777956E+004, 0.868018684755227330E+003, 0.110906305554965877E+003, 0.669388731186872660E+002, 0.395642042804072865E+001]

Derived Types

type, public :: State

Components

TypeVisibilityAttributesNameInitial
real(kind=wp), public :: mean_molecular_weight =0.0_wp
real(kind=wp), public :: temperature =0.0_wp
real(kind=wp), public :: pressure =0.0_wp
real(kind=wp), public :: speed_of_sound =0.0_wp

Type-Bound Procedures

procedure, public :: density

Functions

private pure function find(x, xvec) result(i)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: x
real(kind=wp), intent(in), dimension(:):: xvec

Return Value integer

private pure function density(s)

Arguments

TypeIntentOptionalAttributesName
class(State), intent(in) :: s

Return Value real(kind=wp)

private pure function geopotential_altitude(Z)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: Z

Return Value real(kind=wp)

private pure function Tm(H)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: H

Return Value real(kind=wp)

private pure function temperature_lower(H, M)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: H
real(kind=wp), intent(in) :: M

Return Value real(kind=wp)

private pure function pressure_lower(H)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: H

Return Value real(kind=wp)

private pure function speed_of_sound_lower(T, M)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: T
real(kind=wp), intent(in) :: M

Return Value real(kind=wp)

private pure function mean_molecular_weight_ratio_lower(Z)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: Z

Return Value real(kind=wp)

private pure function mean_molecular_weight_lower(Z)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: Z

Return Value real(kind=wp)

private pure function temperature_upper(Z)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: Z

Return Value real(kind=wp)

private pure function interpolation_index(Z) result(i)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: Z

Return Value integer

private pure function pressure_upper(Z)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: Z

Return Value real(kind=wp)

private pure function mean_molecular_weight_upper(Z)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: Z

Return Value real(kind=wp)

private pure function speed_of_sound_86km()

Arguments

None

Return Value real(kind=wp)

private pure function dynamic_viscosity(Z, T)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: Z

altitude (m)

real(kind=wp), intent(in) :: T

temperature (K)

Return Value real(kind=wp)

public function COESA_atmosphere(Z) result(s)

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: Z

Return Value type(State)

public function COESA_density(Z) result(density)

a simple version that just returns density

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: Z

altitude in meters

Return Value real(kind=wp)

kg/m^3


Subroutines

private subroutine initialize()

Arguments

None

private pure subroutine interpolation_scale_factors(i, Z, scale0, scale1, scale2)

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: i
real(kind=wp), intent(in) :: Z
real(kind=wp), intent(out) :: scale0
real(kind=wp), intent(out) :: scale1
real(kind=wp), intent(out) :: scale2