[ncl-talk] converting pressure levels to hybrid height coordiantes ( output of HadGEM2-ES).

Dennis Shea shea at ucar.edu
Tue Jun 23 14:08:35 MDT 2015


Hello Chao Tang,

I'm sure you can get what you want  .....
It may be relatively simple or it *may* involve multiple steps (multiple
files).

[1]
I do not work with the CMIP5 archive. Certainly not the "HadGEM2-ES
<http://view.es-doc.org/?renderMethod=name&type=cim.1.software.ModelComponent&name=HadGEM2-ES&project=CMIP5>
: Hadley Global Environment Model 2 - Earth System" 6-hrly output.

If there is a file corresponding to

    hus_6hrLev_HadGEM2-ES_historical_r1i1p1_198412010600.nc


        double hus(time, lev, lat, lon) ;

                hus:standard_name = "specific_humidity" ;

                hus:long_name = "Specific Humidity" ;



that contains pressure values at each grid for all levels, for example,

    pressure_6hrLev_HadGEM2-ES_historical_r1i1p1_198412010600.nc

        double pressure(time, lev, lat, lon) ;

                pressure:standard_name = "pressure" ;

                pressure:long_name = "Pa" ;             <=== or "hPa"


then the process is fairly simple:

       level = (/1, 2, 3, 5, 7, 10, 20, 30, 50, 70, 100, 125, 150, 175,
200, 225,                \

                  250, 300, 350, 400, 450, 500, 550, 600, 650, 700, 750,
775, 800, 825,  \

                  850, 875, 900, 925, 950, 975, 1000 /)*100

      level at units = Pa"


       fhus =
addfile("hus_6hrLev_HadGEM2-ES_historical_r1i1p1_198412010600.nc","r")

       fprs =
addfile("pressure_6hrLev_HadGEM2-ES_historical_r1i1p1_198412010600.nc","r")

       q     = fhus->hus                      ; (time,lev,lat,lon)

       p     = fprs->pressure              ; (time,lev,lat,lon)

  ;;
https://www.ncl.ucar.edu/Document/Functions/Contributed/int2p_n_Wrap.shtml

       linlog = -1   ; any option

   q_isobaric = *int2p_n_Wrap*(p,q,level,linlog,0) ; (time,lev,lat,lon)

       printVarSummary(q_isobaric)


[2]
If [1] is not available, the user (you) must calculate the pressure values
yourself.
Again, I do not work with HadGEM2-ES data.

By definition, the hybrid height coordinate is

             z(t,k,j,i) = a(k) + b(k)*orog(t,j,i)

<-- or, if t is degenerate ('t' is of size 1)

             z(k,j,i) = a(k) + b(k)*orog(j,i)

'b_bnds' and 'orog' are on the 'hus' file.

=====

 On the file you provided, using NCL syntax

      b = 0.5*( b_bnds(:,0) + b_bnds(:,1) )   ; b[*]

So, 'we' have the b(k)  ..


What are the 'a(k)' ? Can we use the following?

           a = lev      ; a(k) = lev(k)

In this case heights at each level and grid point can readily be computed
using the above formula.

         z(time, lev, lat, lon)    or    z(lev, lat, lon)

====
Still, there is no 'pressure' information available. Values must be
computed. This is the key.

How can the *pressures* at each hybrid height coordinate be obtained?
Certainly, the pressure levels could be *estimated* from the standard
atmosphere

http://www.ncl.ucar.edu/Document/Functions/Built-in/stdatmus_z2tdp.shtml
However, it is very unlikely that this is what you want.
====

The hydrostatic equation could be used if surface pressure values are
available. This requires orography, temperature, humidity and surface
pressure.

Use the hydrostatic formula:

     z2-z1 = (T*R/g) ln (p1/p2)

where z2 and z1 are the calculated height values, g=gravity, R= dry gas
constant

    p2 = p1*exp[(-(g/RT)*(z2-z1)]

Technically, 'T' should be mean *layer virtual temperature* (which requites
'q')

For example, if z1=0 (mean sea level) and 'pmsl' is the sea level pressure,
the pressure at an arbitrary height 'z' is

     p = pmsl*exp(-g*z/(R*T))

--
In your case, initially, z1 could initially be orography [ 'orog' ]  and
'pmsl' could be surface pressure ('ps').
Then, use the calculated 'z' and 'layer T' to create a 'pressure' array at
each grid point and level.
Then us [1] for calculation.

Good luck













On Mon, Jun 22, 2015 at 12:10 AM, chao.tang <chao.tang at univ-reunion.fr>
wrote:

> Hi Dennis et all,
>
>
> I send you these two file, I want to convert the vertical levels from
> hybrid height coordinate to pressure level.
>
> As exampled bellow I attached two files with just one time step.
>
>
>
>
> the first file (HadGEM2 one ) is in hybrid height coordinate, the second
> (EIN75) is in pressure levels.
> Please note that they cover different domains, the EIN75 one was cutted to
> reduce the size for easier transfer by email.
>
> thanks very much.
>
>
> *Best regards,*
> *Chao TANG*
> *--------------------------------------*
> *--------------------------------------**---------------------*
> *Laboratoire d'Energétique, d'Electronique et Procédés*
> *Unité de Recherche EA 4079*
> Université de La Réunion, LE2P
> 15 avenue René Cassin CS 92003
> 97744 ST DENIS CEDEX 9
> Mobile: 0693 91 57 38 <secretariat.le2p at univ-reunion.fr>
> Courriel: chao.tang at univ-reunion.fr  <secretariat.le2p at univ-reunion.fr>
> *--------------------------------------*
> *--------------------------------------**---------------------*
>
> On 15 Jun 2015, at 19:25, Dennis Shea <shea at ucar.edu> wrote:
>
> From the file and dump you sent offline:
>
> double lev(lev) ;
> lev:bounds = "lev_bnds" ;
> lev:units = "m" ;
> lev:axis = "Z" ;
> lev:positive = "up" ;
> lev:long_name = "hybrid height coordinate" ;
> lev:standard_name = "atmosphere_hybrid_height_coordinate" ;
> lev:formula = "z(k,j,i) = a(k) + b(k)*orog(j,i)" ;
> lev:formula_terms = "a: lev b: b orog: orog" ;
>
> The "b" coefficient variable  is present
>
> double b(lev) ;
> b:long_name = "vertical coordinate formula term: b(k)" ;
>
> b = 0.99771648645401, 0.990881502628326, 0.979542553424835,
> 0.9637770652771,
>     0.943695485591888, 0.919438362121582, 0.891178011894226,
>     0.859118342399597, 0.823493480682373, 0.784570515155792,
>     0.742646217346191, 0.698050200939178, 0.651142716407776,
>     0.602314412593842, 0.55198872089386, 0.500619947910309,
> 0.44869339466095,
>     0.39672577381134, 0.34526526927948, 0.294891387224197,
> 0.24621507525444,
>     0.199878215789795, 0.156554222106934, 0.116947874426842,
>     0.0817952379584312, 0.0518637150526047, 0.0279368180781603,
>     0.0107164792716503, 0.00130179093685001, 0, 0, 0, 0, 0, 0, 0, 0, 0 ;
>
>
> However, there is no "a" coefficient variable present (nor is there an
> "a_bnds" coefficient variable).
>
> If you send another file, please only send one time-step containing the
> 'a' and 'be' coefficients.
> ... that is all that is needed
>
> Also, what might a sample pressure level variable look like"
>
> Is it isobaric levels p(lev) or a 4D pressure variable: p(time,lev,lat,lon)
>
> If the former [ p(lev] ] just include in an email; if 4D please send as a
> netCDF file.
>
> THX
>
>
> On Mon, Jun 15, 2015 at 6:16 AM, chao.tang <chao.tang at univ-reunion.fr>
> wrote:
>
>> Hello everyone,
>>
>>
>> I was searching for the solution to converting pressure levels to hybrid
>> height coordinates (unit: m, output of HadGEM2-ES).
>>
>> Any idea would be really appreciated.
>>
>> thanks in advance, best wishes.
>>
>>
>>
>> *Best regards,*
>> *Chao TANG*
>> *--------------------------------------*
>> *--------------------------------------**---------------------*
>> *Laboratoire d'Energétique, d'Electronique et Procédés*
>> *Unité de Recherche EA 4079*
>> Université de La Réunion, LE2P
>> 15 avenue René Cassin CS 92003
>> 97744 ST DENIS CEDEX 9
>> Mobile: 0693 91 57 38 <secretariat.le2p at univ-reunion.fr>
>> Courriel: chao.tang at univ-reunion.fr  <secretariat.le2p at univ-reunion.fr>
>> *--------------------------------------*
>> *--------------------------------------**---------------------*
>>
>>
>> _______________________________________________
>> ncl-talk mailing list
>> ncl-talk at ucar.edu
>> List instructions, subscriber options, unsubscribe:
>> http://mailman.ucar.edu/mailman/listinfo/ncl-talk
>>
>>
>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: http://mailman.ucar.edu/pipermail/ncl-talk/attachments/20150623/02faf71b/attachment.html 


More information about the ncl-talk mailing list