hum_subs.py 20.5 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26
import numpy as np
import sys
from util_subs import (CtoK)

def VaporPressure(temp, P, phase, meth):
    """
    Calculate the saturation vapor pressure

    For temperatures above 0 deg C the vapor pressure over liquid water
    is calculated.

    The optional parameter 'liquid' changes the calculation to vapor pressure
    over liquid water over the entire temperature range.

    The current default fomulas are Hyland and Wexler for liquid and
    Goff Gratch for ice.

    Parameters
    ----------
    temp : float
        Temperature in [deg C]
    phase : str
        'liquid' : Calculate vapor pressure over liqiud water or
        'ice' : Calculate vapor pressure over ice
    meth : str
        formula to be used
27
        Hardy               : vaporpressure formula from Hardy (1998)
28 29
        MagnusTetens        : vaporpressure formula from Magnus Tetens
        GoffGratch          : vaporpressure formula from Goff Gratch
30 31 32 33
        Buck                : vaporpressure formula from Buck (1981)
        Buck2               : vaporpressure formula from the Buck (2012)
        WMO                 : vaporpressure formula from WMO (1988)
        WMO2018             : vaporpressure formula from WMO (2018)
34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50
        Wexler              : vaporpressure formula from Wexler (1976)
        Sonntag             : vaporpressure formula from Sonntag (1994)
        Bolton              : vaporpressure formula from Bolton (1980)
        HylandWexler        : vaporpressure formula from Hyland and Wexler (1983)
        IAPWS               : vaporpressure formula from IAPWS (2002)
        Preining            : vaporpressure formula from Preining (2002)
        MurphyKoop          : vaporpressure formula from Murphy and Koop (2005)

    Returns
    -------
    P : float
        Saturation vapor pressure in [hPa]


    Author
    ------
    Ported to Python and modified by S. Biri from Holger Voemel's original
51
    """
52 53 54 55 56 57 58 59 60

    Psat = np.zeros(temp.size)*np.nan
    if (np.nanmin(temp) > 200):  # if Ta in Kelvin convert to Celsius
        temp = temp-273.16
    T = np.copy(temp)+273.16    #  Most formulas use T in [K]
                                #  Formulas using [C] use the variable temp
    #  Calculate saturation pressure over liquid water
    if (phase == 'liquid'):
        if (meth == 'HylandWexler' or meth == ''):
61 62 63 64
            """
            Source Hyland, R. W. and A. Wexler, Formulations for the
            Thermodynamic Properties of the saturated Phases of H2O from
            173.15K to 473.15K, ASHRAE Trans, 89(2A), 500-519, 1983."""
65 66 67 68 69
            Psat = np.exp(-0.58002206e4/T+0.13914993e1-0.48640239e-1*T +
                          0.41764768e-4*np.power(T, 2) -
                          0.14452093e-7*np.power(T, 3) +
                          0.65459673e1*np.log(T))/100
        if (meth == 'Hardy'):
70 71 72 73 74
            """
            Source Hardy, B., 1998, ITS-90 Formulations for Vapor Pressure,
            Frostpoint Temperature, Dewpoint Temperature, and Enhancement
            Factors in the Range -100 to +100°C, The Proceedings of the Third
            International Symposium on Humidity & Moisture, London, England"""
75 76 77 78 79 80 81
            Psat = np.exp(-2.8365744e3/np.power(T, 2)-6.028076559e3/T +
                          1.954263612e1-2.737830188e-2*T +
                          1.6261698e-5*np.power(T, 2) +
                          7.0229056e-10*np.power(T, 3) -
                          1.8680009e-13*np.power(T, 4) +
                          2.7150305*np.log(T))/100
        if (meth == 'Preining'):
82 83 84
            """Source : Vehkamaeki, H., M. Kulmala, I. Napari, K. E. J.
            Lehtinen, C.Timmreck, M. Noppel, and A. Laaksonen (2002),
            J. Geophys. Res., 107, doi:10.1029/2002JD002184."""
85 86 87
            Psat = np.exp(-7235.424651/T+77.34491296+5.7113e-3*T -
                          8.2*np.log(T))/100
        if (meth == 'Wexler'):
88 89 90 91
            """Wexler, A., Vapor Pressure Formulation for Water in Range 0 to
            100 C. A Revision, Journal of Research of the National Bureau of
            Standards - A. Physics and Chemistry, September - December 1976,
            Vol. 80A, Nos.5 and 6, 775-785"""
92 93 94 95 96 97 98 99
            Psat = np.exp(-0.29912729e4*np.power(T, -2) -
                          0.60170128e4*np.power(T, -1) +
                          0.1887643854e2-0.28354721e-1*T +
                          0.17838301e-4*np.power(T, 2) -
                          0.84150417e-9*np.power(T, 3) +
                          0.44412543e-12*np.power(T, 4) +  # This line was corrected from '-' to '+' following the original citation. (HV 20140819). The change makes only negligible difference
                          2.858487*np.log(T))/100
        if ((meth == 'GoffGratch') or (meth == 'MartiMauersberger')):
100 101 102 103 104
            """Marti and Mauersberger don't have a vapor pressure curve over
            liquid. Using Goff Gratch instead; Goff Gratch formulation Source :
            Smithsonian Meteorological Tables, 5th edition, p. 350, 1984
            From original source: Goff & Gratch (1946), p. 107)
            """
105 106 107 108 109 110 111
            Ts = 373.16  # steam point temperature in K
            ews = 1013.246  # saturation pressure at steam point temperature
            Psat = np.power(10, -7.90298*(Ts/T-1)+5.02808*np.log10(Ts/T) -
                            1.3816e-7*(np.power(10, 11.344*(1-T/Ts))-1) +
                            8.1328e-3*(np.power(10, -3.49149*(Ts/T-1))-1) +
                            np.log10(ews))
        if (meth == 'MagnusTetens'):
112
            """Source: Murray, F. W., On the computation of \
113
                         saturation vapor pressure, J. Appl. Meteorol., \
114
                         6, 203-204, 1967."""
115 116 117 118 119
            Psat = np.power(10, 7.5*(temp)/(temp+237.5)+0.7858)
            # Murray quotes this as the original formula and
            Psat = 6.1078*np.exp(17.269388*temp/(temp+237.3))
            # this as the mathematical aquivalent in the form of base e.
        if (meth == 'Buck'):
120 121 122
            """Bucks vapor pressure formulation based on Tetens formula.
            Source: Buck, A. L., New equations for computing vapor pressure and
            enhancement factor, J. Appl. Meteorol., 20, 1527-1532, 1981."""
123 124 125
            Psat = (6.1121*np.exp(17.502*temp/(240.97+temp)) *
                    (1.0007+(3.46e-6*P)))
        if (meth == 'Buck2'):
126 127
            """Bucks vapor pressure formulation based on Tetens formula.
            Source: Buck Research, Model CR-1A Hygrometer Operating Manual,
128
            May 2012"""
129 130 131
            Psat = (6.1121*np.exp((18.678-(temp)/234.5)*(temp)/(257.14+temp)) *
                    (1+1e-4*(7.2+P*(0.0320)+5.9e-6*np.power(T, 2))))
        if (meth == 'WMO'):
132 133 134 135 136 137
            """Intended WMO formulation, originally published by Goff (1957)
            incorrectly referenced by WMO technical regulations, WMO-NO 49,
            Vol I, General Meteorological Standards and Recommended Practices,
            App. A, Corrigendum Aug 2000. and incorrectly referenced by WMO
            technical regulations, WMO-NO 49, Vol I, General Meteorological
            Standards and Recommended Practices, App. A, 1988."""
138 139 140
            Ts = 273.16  # triple point temperature in K
            Psat = np.power(10, 10.79574*(1-Ts/T)-5.028*np.log10(T/Ts) +
                            1.50475e-4*(1-np.power(10, -8.2969*(T/Ts-1))) +
141
                            0.42873e-3*(np.power(10, 4.76955*(1-Ts/T))-1) +  # in eq. 13 is -4.76955; in aerobulk is like this
142
                            0.78614)
143 144 145 146
        if (meth == 'WMO2018'):
            """WMO 2018 edition. Annex 4.B, eq. 4.B.1, 4.B.2, 4.B.5 """
            Psat = 6.112*np.exp(17.62*temp/(243.12+temp))*(1.0016+3.15e-6*P -
                                                           0.074/P)
147
        if (meth == 'Sonntag'):
148 149
            """Source: Sonntag, D., Advancements in the field of hygrometry,
            Meteorol. Z., N. F., 3, 51-66, 1994."""
150 151 152 153
            Psat = np.exp(-6096.9385*np.power(T, -1)+16.635794 -
                          2.711193e-2*T+1.673952e-5*np.power(T, 2) +
                          2.433502*np.log(T))#*(1.0016+P*3.15e-6-0.074/P)
        if (meth == 'Bolton'):
154 155 156
            """Source: Bolton, D., The computation of equivalent potential
            temperature, Monthly Weather Report, 108, 1046-1053, 1980.
            equation (10)"""
157 158
            Psat = 6.112*np.exp(17.67*temp/(temp+243.5))
        if (meth == 'IAPWS'):
159 160 161 162 163 164 165
            """Source: Wagner W. and A. Pruss (2002), The IAPWS formulation
            1995 for the thermodynamic properties of ordinary water substance
            for general and scientific use, J. Phys. Chem. Ref. Data, 31(2),
            387-535. This is the 'official' formulation from the International
            Association for the Properties of  Water and Steam The valid range
            of this formulation is 273.16 <= T <= 647.096 K and is based on the
            ITS90 temperature scale."""
166 167 168 169 170 171 172 173 174
            Tc = 647.096   # K   : Temperature at the critical point
            Pc = 22.064e4  #  hPa : Vapor pressure at the critical point
            nu = (1-T/Tc)
            a1, a2, a3 = -7.85951783, 1.84408259, -11.7866497
            a4, a5, a6 = 22.6807411, -15.9618719, 1.80122502
            Psat = (Pc*np.exp(Tc/T*(a1*nu+a2*np.power(nu, 1.5) +
                    a3*np.power(nu, 3)+a4*np.power(nu, 3.5) +
                    a5*np.power(nu, 4)+ a6*np.power(nu, 7.5))))
        if (meth == 'MurphyKoop'):
175 176 177
            """Source : Murphy and Koop, Review of the vapour pressure of ice
            and supercooled water for atmospheric applications, Q. J. R.
            Meteorol. Soc (2005), 131, pp. 1539-1565."""
178 179 180 181 182
            Psat = np.exp(54.842763-6763.22/T-4.210*np.log(T)+0.000367*T +
                          np.tanh(0.0415*(T-218.8))*(53.878-1331.22/T -
                          9.44523*np.log(T)+0.014025*T))/100
    # Calculate saturation pressure over ice ----------------------------------
    elif (phase == 'ice'):
183 184
        """Default uses Goff Gratch over ice. There is little ambiguity in the
        ice saturation curve. Goff Gratch is widely used."""
185
        if (meth == 'MartiMauersberger'):
186 187 188
            """Source : Marti, J. and K Mauersberger, A survey and new
            measurements of ice vapor pressure at temperatures between 170 and
            250 K, GRL 20, 363-366, 1993."""
189 190
            Psat = np.power(10, -2663.5/T+12.537)/100
        if (meth == 'HylandWexler'):
191 192 193
            """Source Hyland, R. W. and A. Wexler, Formulations for the
            Thermodynamic Properties of the saturated Phases of H2O from
            173.15K to 473.15K, ASHRAE Trans, 89(2A), 500-519, 1983."""
194 195 196 197 198 199
            Psat = np.exp(-0.56745359e4/T+0.63925247e1-0.96778430e-2*T +
                          0.62215701e-6*np.power(T, 2) +
                          0.20747825e-8*np.power(T, 3) -
                          0.9484024e-12*np.power(T, 4) +
                          0.41635019e1*np.log(T))/100
        if (meth == 'Wexler'):
200 201
            """Wexler, A., Vapor pressure formulation for ice, Journal of
            Research of the National Bureau of Standards-A. 81A, 5-20, 1977."""
202 203 204 205 206
            Psat = np.exp(-0.58653696e4*np.power(T, -1)+0.22241033e2 +
                          0.13749042e-1*T-0.34031775e-4*np.power(T, 2) +
                          0.26967687e-7*np.power(T, 3) +
                          0.6918651*np.log(T))/100
        if (meth == 'Hardy'):
207 208 209 210 211 212 213
            """Source Hardy, B., 1998, ITS-90 Formulations for Vapor Pressure,
            Frostpoint Temperature, Dewpoint Temperature, and Enhancement
            Factors in the Range -100 to +100°C, The Proceedings of the Third
            International Symposium on Humidity & Moisture, London, England.
            These coefficients are updated to ITS90 based on the work by
            Bob Hardy at Thunder Scientific: http://www.thunderscientific.com/
            tech_info/reflibrary/its90formulas.pdf """
214 215 216 217 218
            Psat = np.exp(-0.58666426e4*np.power(T, -1)+0.2232870244e2 +
                          0.139387003e-1*T-0.34262402e-4*np.power(T, 2) +
                          0.27040955e-7*np.power(T, 3) +
                          0.67063522e-1*np.log(T))/100
        if (meth == 'GoffGratch' or meth == '' or meth == 'IAPWS'):
219 220 221 222
            """IAPWS does not provide a vapor pressure formulation over ice use
            Goff Gratch instead.
            Source : Smithsonian Meteorological Tables, 5th edition, p. 350,
            1984"""
223 224 225 226 227
            ei0 = 6.1071  # mbar
            T0 = 273.16   # triple point in K
            Psat = np.power(10, -9.09718*(T0/T-1)-3.56654*np.log10(T0/T) +
                            0.876793*(1-T/T0)+np.log10(ei0))
        if (meth == 'MagnusTetens'):
228 229
            """Source: Murray, F. W., On the computation of saturation vapor
            pressure, J. Appl. Meteorol., 6, 203-204, 1967."""
230 231 232 233 234
            Psat = np.power(10, 9.5*temp/(265.5+temp)+0.7858)
            #  Murray quotes this as the original formula and
            Psat = 6.1078*np.exp(21.8745584*(T-273.16)/(T-7.66))
            # this as the mathematical aquivalent in the form of base e.
        if (meth == 'Buck'):
235 236 237
            """Bucks vapor pressure formulation based on Tetens formula.
            Source: Buck, A. L., New equations for computing vapor pressure and
            enhancement factor, J. Appl. Meteorol., 20, 1527-1532, 1981."""
238 239 240
            Psat = (6.1115*np.exp(22.452*temp/(272.55+temp)) *
                    (1.0003+(4.18e-6*P)))
        if (meth == 'Buck2'):
241 242 243
            """Bucks vapor pressure formulation based on Tetens formula.
            Source: Buck Research, Model CR-1A Hygrometer Operating Manual,
            Sep 2001"""
244 245 246
            Psat = (6.1115*np.exp((23.036-temp/333.7)*temp/(279.82+temp)) *
                    (1+1e-4*(2.2+P*(0.0383+6.4e-6*np.power(T, 2)))))
        if (meth == 'CIMO'):
247 248 249
            """Source: Annex 4B, Guide to Meteorological Instruments and
            Methods of Observation, WMO Publication No 8, 7th edition, Geneva,
            2008. (CIMO Guide)"""
250 251 252
            Psat = (6.112*np.exp(22.46*temp/(272.62+temp)) *
                    (1.0016+3.15e-6*P-0.074/P))
        if (meth == 'WMO' or meth == 'WMO2000'):
253 254 255 256 257
            """There is no typo issue in the WMO formulations for ice.
            WMO formulation, which is very similar to Goff & Gratch.
            Source : WMO technical regulations, WMO-NO 49, Vol I, General
            Meteorological Standards and Recommended Practices, Aug 2000,
            App. A."""
258 259 260 261
            T0 = 273.16  # triple point temperature in K
            Psat = np.power(10, -9.09685*(T0/T-1)-3.56654*np.log10(T0/T) +
                            0.87682*(1-T/T0)+0.78614)
        if (meth == 'Sonntag'):
262 263
            """Source: Sonntag, D., Advancements in the field of hygrometry,
            Meteorol. Z., N. F., 3, 51-66, 1994."""
264 265 266
            Psat = np.exp(-6024.5282*np.power(T, -1)+24.721994+1.0613868e-2*T -
                          1.3198825e-5*np.power(T, 2)-0.49382577*np.log(T))
        if (meth == 'MurphyKoop'):
267 268 269
            """Source : Murphy and Koop, Review of the vapour pressure of ice
            and supercooled water for atmospheric applications, Q. J. R.
            Meteorol. Soc (2005), 131, pp. 1539-1565."""
270 271 272 273 274
            Psat = np.exp(9.550426-5723.265/T+3.53068*np.log(T) -
                          0.00728332*T)/100

        s = np.where(temp > 0)
        if (s.size[0] >= 1):
275 276 277 278 279
            """Independent of the formula used for ice, use Hyland Wexler
            (water) for temperatures above freezing (see above).
            Source Hyland, R. W. and A. Wexler, Formulations for the
            Thermodynamic Properties of the saturated Phases of H2O from
            173.15K to 473.15K, ASHRAE Trans, 89(2A), 500-519, 1983."""
280 281 282 283 284 285 286 287 288 289 290
            Psat_w = np.exp(-0.58002206e4/T+0.13914993e1-0.48640239e-1*T +
                            0.41764768e-4*np.power(T, 2) -
                            0.14452093e-7*np.power(T, 3) +
                            0.65459673e1*np.log(T))/100
            Psat[s] = Psat_w[s]

    return Psat
#-----------------------------------------------------------------------------


def qsat_sea(T, P, qmeth):
sbiri's avatar
sbiri committed
291 292
    """
    computes surface saturation specific humidity [g/kg]
293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317

    Parameters
    ----------
    T : float
        temperature ($^\\circ$\\,C)
    P : float
        pressure (mb)
    qmeth : str
        method to calculate vapor pressure

    Returns
    -------
    qs : float
    """
    T = np.asarray(T)
    if (np.nanmin(T) > 200):  # if Ta in Kelvin convert to Celsius
        T = T-CtoK
    ex = VaporPressure(T, P, 'liquid', qmeth)
    es = 0.98*ex  # reduction at sea surface
    qs = 622*es/(P-0.378*es)
    return qs
#-----------------------------------------------------------------------------


def qsat_air(T, P, rh, qmeth):
sbiri's avatar
sbiri committed
318 319
    """
    computes saturation specific humidity [g/kg]
320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355

    Parameters
    ----------
    T : float
        temperature ($^\circ$\,C)
    P : float
        pressure (mb)
    rh : float
       relative humidity (%)
    qmeth : str
        method to calculate vapor pressure

    Returns
    -------
    q : float
    em : float
    """
    T = np.asarray(T)
    if (np.nanmin(T) > 200):  # if Ta in Kelvin convert to Celsius
        T = T-CtoK
    es = VaporPressure(T, P, 'liquid', qmeth)
    em = 0.01*rh*es
    q = 622*em/(P-0.378*em)
    return q
#------------------------------------------------------------------------------


def get_hum(hum, T, sst, P, qmeth):
    """
    Get specific humidity output

    Parameters
    ----------
    hum : array
        humidity input switch 2x1 [x, values] default is relative humidity
            x='rh' : relative humidity in %
sbiri's avatar
sbiri committed
356
            x='q' : specific humidity (kg/kg)
357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374
            x='Td' : dew point temperature (K)
    T : float
        air temperature in K
    sst : float
        sea surface temperature in K
    P : float
        air pressure at sea level in hPa
    qmeth : str
        method to calculate specific humidity from vapor pressure

    Returns
    -------
    qair : float
        specific humidity of air
    qsea : float
        specific humidity over sea surface

    """
375
    if (hum[0] not in ['rh', 'q', 'Td']):
376 377 378 379 380 381 382
        sys.exit("unknown humidity input")
        qair, qsea = np.nan, np.nan
    elif (hum[0] == 'rh'):
        RH = hum[1]
        if (np.all(RH < 1)):
            sys.exit("input relative humidity units should be \%")
            qair, qsea = np.nan, np.nan
383
        # RH = np.where(RH > 100, np.nan, RH)  # ensure RH <=100
sbiri's avatar
sbiri committed
384 385
        qsea = qsat_sea(sst, P, qmeth)/1000    # surface water q (kg/kg)
        qair = qsat_air(T, P, RH, qmeth)/1000  # q of air (kg/kg)
386 387
    elif (hum[0] == 'q'):
        qair = hum[1]
sbiri's avatar
sbiri committed
388
        qsea = qsat_sea(sst, P, qmeth)/1000  # surface water q (kg/kg)
389 390 391 392 393 394 395
    elif (hum[0] == 'Td'):
        Td = hum[1] # dew point temperature (K)
        Td = np.where(Td < 200, np.copy(Td)+CtoK, np.copy(Td))
        T = np.where(T < 200, np.copy(T)+CtoK, np.copy(T))
        esd = 611.21*np.exp(17.502*((Td-273.16)/(Td-32.19)))
        es = 611.21*np.exp(17.502*((T-273.16)/(T-32.19)))
        RH = 100*esd/es
396
        # RH = np.where(RH > 100, np.nan, RH)  # ensure RH <=100
sbiri's avatar
sbiri committed
397 398
        qair = qsat_air(T, P, RH, qmeth)/1000  # q of air (kg/kg)
        qsea = qsat_sea(sst, P, qmeth)/1000    # surface water q (kg/kg)
399
    return qair, qsea
400
#------------------------------------------------------------------------------
401 402


403
def gamma(opt, sst, t, q, cp):
404
    """
405
    Computes the adiabatic lapse-rate
406 407 408

    Parameters
    ----------
409 410 411 412
    opt : str
        type of adiabatic lapse rate dry or "moist"
        dry has options to be constant "dry_c", for dry air "dry", or
        for unsaturated air with water vapor "dry_v"
413 414 415 416 417 418
    sst : float
        sea surface temperature [K]
    t : float
        air temperature [K]
    q : float
        specific humidity of air [kg/kg]
419 420
    cp : float
        specific capacity of air at constant Pressure [kJ/(kg*K)]
421 422 423 424 425 426 427 428 429 430
    Returns
    -------
    gamma : float
        lapse rate [K/m]

    """
    if (np.nanmin(sst) < 200):  # if sst in Celsius convert to Kelvin
        sst = sst+CtoK
    if (np.nanmin(t) < 200):  # if sst in Celsius convert to Kelvin
        t = t+CtoK
431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447
    if (opt == "moist"):
        t = np.maximum(t, 180)
        q = np.maximum(q,  1e-6)
        w = q/(1-q) # mixing ratio w = q/(1-q)
        iRT = 1/(287.05*t)
        # latent heat of vaporization of water as a function of temperature
        lv = (2.501-0.00237*(sst-CtoK))*1e6
        gamma = 9.8*(1+lv*w*iRT)/(1005+np.power(lv, 2)*w*(287.05/461.495) *
                                  iRT/t)
    elif (opt == "dry_c"):
        gamma = 0.0098*np.ones(t.shape)
    elif (opt == "dry"):
        gamma = 9.81/cp
    elif (opt == "dry_v"):
        w = q/(1-q) # mixing ratio
        f_v = 1-0.85*w #(1+w)/(1+w*)
        gamma = f_v*9.81/cp
448 449
    return gamma
#------------------------------------------------------------------------------