COESA_density Function

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


Calls

proc~~coesa_density~~CallsGraph proc~coesa_density COESA_density proc~pressure_upper pressure_upper proc~coesa_density->proc~pressure_upper proc~mean_molecular_weight_upper mean_molecular_weight_upper proc~coesa_density->proc~mean_molecular_weight_upper proc~geopotential_altitude geopotential_altitude proc~coesa_density->proc~geopotential_altitude proc~mean_molecular_weight_lower mean_molecular_weight_lower proc~coesa_density->proc~mean_molecular_weight_lower proc~temperature_upper temperature_upper proc~coesa_density->proc~temperature_upper proc~pressure_lower pressure_lower proc~coesa_density->proc~pressure_lower proc~temperature_lower temperature_lower proc~coesa_density->proc~temperature_lower proc~interpolation_index interpolation_index proc~pressure_upper->proc~interpolation_index proc~interpolation_scale_factors interpolation_scale_factors proc~pressure_upper->proc~interpolation_scale_factors proc~mean_molecular_weight_upper->proc~interpolation_index proc~mean_molecular_weight_upper->proc~interpolation_scale_factors proc~mean_molecular_weight_ratio_lower mean_molecular_weight_ratio_lower proc~mean_molecular_weight_lower->proc~mean_molecular_weight_ratio_lower proc~find find proc~pressure_lower->proc~find proc~tm Tm proc~temperature_lower->proc~tm proc~interpolation_index->proc~find proc~tm->proc~find proc~mean_molecular_weight_ratio_lower->proc~find

Contents

Source Code


Source Code

function COESA_density(Z) result(density)
    !! a simple version that just returns density
    real(wp),intent(in) :: Z !! altitude in meters
    real(wp) :: density !! kg/m^3
    real(wp) :: H,M,T,P
    if (Z < -5000.0_wp) then
        error stop "altitude below lower bound of -5000 m"
    else if (Z > 1000000.0_wp) then
        error stop "altitude above upper bound of 1000000 m"
    else if (Z < 86000.0_wp) then
        H = geopotential_altitude(Z)
        M = mean_molecular_weight_lower(Z)
        T = temperature_lower(H, M)
        P = pressure_lower(H)
    else
        T = temperature_upper(Z)
        P = pressure_upper(Z)
        M = mean_molecular_weight_upper(Z)
    end if
    density = P * M / (Rstar * T)
end function COESA_density