source: lmdz_wrf/trunk/tools/diagnostics.py @ 463

Last change on this file since 463 was 449, checked in by lfita, 10 years ago

Fixing minor typos

File size: 36.6 KB
Line 
1# Pthon script to comput diagnostics
2# L. Fita, LMD. CNR, UPMC-Jussieu, Paris, France
3# File diagnostics.inf provides the combination of variables to get the desired diagnostic
4#
5## e.g. # diagnostics.py -d 'Time@time,bottom_top@ZNU,south_north@XLAT,west_east@XLONG' -v 'clt|CLDFRA,cllmh|CLDFRA@WRFp,RAINTOT|RAINC@RAINNC@XTIME' -f WRF_LMDZ/NPv31/wrfout_d01_1980-03-01_00:00:00
6## e.g. # diagnostics.py -f /home/lluis/PY/diagnostics.inf -d variable_combo -v WRFprc
7
8from optparse import OptionParser
9import numpy as np
10from netCDF4 import Dataset as NetCDFFile
11import os
12import re
13import nc_var_tools as ncvar
14import datetime as dt
15
16main = 'diagnostics.py'
17errormsg = 'ERROR -- error -- ERROR -- error'
18warnmsg = 'WARNING -- warning -- WARNING -- warning'
19
20# Gneral information
21##
22def reduce_spaces(string):
23    """ Function to give words of a line of text removing any extra space
24    """
25    values = string.replace('\n','').split(' ')
26    vals = []
27    for val in values:
28         if len(val) > 0:
29             vals.append(val)
30
31    return vals
32
33def variable_combo(varn,combofile):
34    """ Function to provide variables combination from a given variable name
35      varn= name of the variable
36      combofile= ASCII file with the combination of variables
37        [varn] [combo]
38          [combo]: '@' separated list of variables to use to generate [varn]
39            [WRFdt] to get WRF time-step (from general attributes)
40    >>> variable_combo('WRFprls','/home/lluis/PY/diagnostics.inf')
41    deaccum@RAINNC@XTIME@prnc
42    """
43    fname = 'variable_combo'
44
45    if varn == 'h':
46        print fname + '_____________________________________________________________'
47        print variable_combo.__doc__
48        quit()
49
50    if not os.path.isfile(combofile):
51        print errormsg
52        print '  ' + fname + ": file with combinations '" + combofile +              \
53          "' does not exist!!"
54        quit(-1)
55
56    objf = open(combofile, 'r')
57
58    found = False
59    for line in objf:
60        linevals = reduce_spaces(line)
61        varnf = linevals[0]
62        combo = linevals[1].replace('\n','')
63        if varn == varnf: 
64            found = True
65            break
66
67    if not found:
68        print errormsg
69        print '  ' + fname + ": variable '" + varn + "' not found in '" + combofile +\
70          "' !!"
71        combo='ERROR'
72
73    objf.close()
74
75    return combo
76
77# Mathematical operators
78##
79def compute_deaccum(varv, dimns, dimvns):
80    """ Function to compute the deaccumulation of a variable
81    compute_deaccum(varv, dimnames, dimvns)
82      [varv]= values to deaccum (assuming [t,])
83      [dimns]= list of the name of the dimensions of the [varv]
84      [dimvns]= list of the name of the variables with the values of the
85        dimensions of [varv]
86    """
87    fname = 'compute_deaccum'
88
89    deacdims = dimns[:]
90    deacvdims = dimvns[:]
91
92    slicei = []
93    slicee = []
94
95    Ndims = len(varv.shape)
96    for iid in range(0,Ndims):
97        slicei.append(slice(0,varv.shape[iid]))
98        slicee.append(slice(0,varv.shape[iid]))
99
100    slicee[0] = np.arange(varv.shape[0])
101    slicei[0] = np.arange(varv.shape[0])
102    slicei[0][1:varv.shape[0]] = np.arange(varv.shape[0]-1)
103
104    vari = varv[tuple(slicei)]
105    vare = varv[tuple(slicee)]
106
107    deac = vare - vari
108
109    return deac, deacdims, deacvdims
110
111def derivate_centered(var,dim,dimv):
112    """ Function to compute the centered derivate of a given field
113      centered derivate(n) = (var(n-1) + var(n+1))/(2*dn).
114    [var]= variable
115    [dim]= which dimension to compute the derivate
116    [dimv]= dimension values (can be of different dimension of [var])
117    >>> derivate_centered(np.arange(16).reshape(4,4)*1.,1,1.)
118    [[  0.   1.   2.   0.]
119     [  0.   5.   6.   0.]
120     [  0.   9.  10.   0.]
121     [  0.  13.  14.   0.]]
122    """
123
124    fname = 'derivate_centered'
125   
126    vark = var.dtype
127
128    if hasattr(dimv, "__len__"):
129# Assuming that the last dimensions of var [..., N, M] are the same of dimv [N, M]
130        if len(var.shape) != len(dimv.shape):
131            dimvals = np.zeros((var.shape), dtype=vark)
132            if len(var.shape) - len(dimv.shape) == 1:
133                for iz in range(var.shape[0]):
134                    dimvals[iz,] = dimv
135            elif len(var.shape) - len(dimv.shape) == 2:
136                for it in range(var.shape[0]):
137                    for iz in range(var.shape[1]):
138                        dimvals[it,iz,] = dimv
139            else:
140                print errormsg
141                print '  ' + fname + ': dimension difference between variable',      \
142                  var.shape,'and variable with dimension values',dimv.shape,         \
143                  ' not ready !!!'
144                quit(-1)
145        else:
146            dimvals = dimv
147    else:
148# dimension values are identical everywhere!
149# from: http://stackoverflow.com/questions/16807011/python-how-to-identify-if-a-variable-is-an-array-or-a-scalar   
150        dimvals = np.ones((var.shape), dtype=vark)*dimv
151
152    derivate = np.zeros((var.shape), dtype=vark)
153    if dim > len(var.shape) - 1:
154        print errormsg
155        print '  ' + fname + ': dimension',dim,' too big for given variable of ' +   \
156          'shape:', var.shape,'!!!'
157        quit(-1)
158
159    slicebef = []
160    sliceaft = []
161    sliceder = []
162
163    for id in range(len(var.shape)):
164        if id == dim:
165            slicebef.append(slice(0,var.shape[id]-2))
166            sliceaft.append(slice(2,var.shape[id]))
167            sliceder.append(slice(1,var.shape[id]-1))
168        else:
169            slicebef.append(slice(0,var.shape[id]))
170            sliceaft.append(slice(0,var.shape[id]))
171            sliceder.append(slice(0,var.shape[id]))
172
173    if hasattr(dimv, "__len__"):
174        derivate[tuple(sliceder)] = (var[tuple(slicebef)] + var[tuple(sliceaft)])/   \
175          ((dimvals[tuple(sliceaft)] - dimvals[tuple(slicebef)]))
176        print (dimvals[tuple(sliceaft)] - dimvals[tuple(slicebef)])
177    else:
178        derivate[tuple(sliceder)] = (var[tuple(slicebef)] + var[tuple(sliceaft)])/   \
179          (2.*dimv)
180
181#    print 'before________'
182#    print var[tuple(slicebef)]
183
184#    print 'after________'
185#    print var[tuple(sliceaft)]
186
187    return derivate
188
189def rotational_z(Vx,Vy,pos):
190    """ z-component of the rotatinoal of horizontal vectorial field
191    \/ x (Vx,Vy,Vz) = \/xVy - \/yVx
192    [Vx]= Variable component x
193    [Vy]=  Variable component y
194    [pos]= poisition of the grid points
195    >>> rotational_z(np.arange(16).reshape(4,4)*1., np.arange(16).reshape(4,4)*1., 1.)
196    [[  0.   1.   2.   0.]
197     [ -4.   0.   0.  -7.]
198     [ -8.   0.   0. -11.]
199     [  0.  13.  14.   0.]]
200    """
201
202    fname =  'rotational_z'
203
204    ndims = len(Vx.shape)
205    rot1 = derivate_centered(Vy,ndims-1,pos)
206    rot2 = derivate_centered(Vx,ndims-2,pos)
207
208    rot = rot1 - rot2
209
210    return rot
211
212# Diagnostics
213##
214
215def var_clt(cfra):
216    """ Function to compute the total cloud fraction following 'newmicro.F90' from
217      LMDZ using 1D vertical column values
218      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
219    """
220    ZEPSEC=1.0E-12
221
222    fname = 'var_clt'
223
224    zclear = 1.
225    zcloud = 0.
226
227    dz = cfra.shape[0]
228    for iz in range(dz):
229        zclear =zclear*(1.-np.max([cfra[iz],zcloud]))/(1.-np.min([zcloud,1.-ZEPSEC]))
230        clt = 1. - zclear
231        zcloud = cfra[iz]
232
233    return clt
234
235def compute_clt(cldfra, dimns, dimvns):
236    """ Function to compute the total cloud fraction following 'newmicro.F90' from
237      LMDZ
238    compute_clt(cldfra, dimnames)
239      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
240      [dimns]= list of the name of the dimensions of [cldfra]
241      [dimvns]= list of the name of the variables with the values of the
242        dimensions of [cldfra]
243    """
244    fname = 'compute_clt'
245
246    cltdims = dimns[:]
247    cltvdims = dimvns[:]
248
249    if len(cldfra.shape) == 4:
250        clt = np.zeros((cldfra.shape[0],cldfra.shape[2],cldfra.shape[3]),            \
251          dtype=np.float)
252        dx = cldfra.shape[3]
253        dy = cldfra.shape[2]
254        dz = cldfra.shape[1]
255        dt = cldfra.shape[0]
256        cltdims.pop(1)
257        cltvdims.pop(1)
258
259        for it in range(dt):
260            for ix in range(dx):
261                for iy in range(dy):
262                    zclear = 1.
263                    zcloud = 0.
264                    ncvar.percendone(it*dx*dy + ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
265                    clt[it,iy,ix] = var_clt(cldfra[it,:,iy,ix])
266
267    else:
268        clt = np.zeros((cldfra.shape[1],cldfra.shape[2]), dtype=np.float)
269        dx = cldfra.shape[2]
270        dy = cldfra.shape[1]
271        dy = cldfra.shape[0]
272        cltdims.pop(0)
273        cltvdims.pop(0)
274        for ix in range(dx):
275            for iy in range(dy):
276                zclear = 1.
277                zcloud = 0.
278                ncvar.percendone(ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
279                clt[iy,ix] = var_clt(cldfra[:,iy,ix])
280
281    return clt, cltdims, cltvdims
282
283def var_cllmh(cfra, p):
284    """ Fcuntion to compute cllmh on a 1D column
285    """
286
287    fname = 'var_cllmh'
288
289    ZEPSEC =1.0E-12
290    prmhc = 440.*100.
291    prmlc = 680.*100.
292
293    zclearl = 1.
294    zcloudl = 0.
295    zclearm = 1.
296    zcloudm = 0.
297    zclearh = 1.
298    zcloudh = 0.
299
300    dvz = cfra.shape[0]
301
302    cllmh = np.ones((3), dtype=np.float)
303
304    for iz in range(dvz):
305        if p[iz] < prmhc:
306            cllmh[2] = cllmh[2]*(1.-np.max([cfra[iz], zcloudh]))/(1.-                \
307              np.min([zcloudh,1.-ZEPSEC]))
308            zcloudh = cfra[iz]
309        elif p[iz] >= prmhc and p[iz] < prmlc:
310            cllmh[1] = cllmh[1]*(1.-np.max([cfra[iz], zcloudm]))/(1.-                \
311              np.min([zcloudm,1.-ZEPSEC]))
312            zcloudm = cfra[iz]
313        elif p[iz] >= prmlc:
314            cllmh[0] = cllmh[0]*(1.-np.max([cfra[iz], zcloudl]))/(1.-                \
315              np.min([zcloudl,1.-ZEPSEC]))
316            zcloudl = cfra[iz]
317
318    cllmh = 1.- cllmh
319
320    return cllmh
321
322def compute_cllmh(cldfra, pres, dimns, dimvns):
323    """ Function to compute cllmh: low/medium/hight cloud fraction following newmicro.F90 from LMDZ
324    compute_clt(cldfra, pres, dimns, dimvns)
325      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
326      [pres] = pressure field
327      [dimns]= list of the name of the dimensions of [cldfra]
328      [dimvns]= list of the name of the variables with the values of the
329        dimensions of [cldfra]
330    """
331    fname = 'compute_cllmh'
332
333    cllmhdims = dimns[:]
334    cllmhvdims = dimvns[:]
335
336    if len(cldfra.shape) == 4:
337        dx = cldfra.shape[3]
338        dy = cldfra.shape[2]
339        dz = cldfra.shape[1]
340        dt = cldfra.shape[0]
341        cllmhdims.pop(1)
342        cllmhvdims.pop(1)
343
344        cllmh = np.ones(tuple([3, dt, dy, dx]), dtype=np.float)
345
346        for it in range(dt):
347            for ix in range(dx):
348                for iy in range(dy):
349                    ncvar.percendone(it*dx*dy + ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
350                    cllmh[:,it,iy,ix] = var_cllmh(cldfra[it,:,iy,ix], pres[it,:,iy,ix])
351       
352    else:
353        dx = cldfra.shape[2]
354        dy = cldfra.shape[1]
355        dz = cldfra.shape[0]
356        cllmhdims.pop(0)
357        cllmhvdims.pop(0)
358
359        cllmh = np.ones(tuple([3, dy, dx]), dtype=np.float)
360
361        for ix in range(dx):
362            for iy in range(dy):
363                ncvar.percendone(ix*dy + iy,dx*dy, 5, 'diagnosted')
364                cllmh[:,iy,ix] = var_cllmh(cldfra[:,iy,ix], pres[:,iy,ix])
365
366    return cllmh, cllmhdims, cllmhvdims
367
368def var_virtualTemp (temp,rmix):
369    """ This function returns virtual temperature in K,
370      temp: temperature [K]
371      rmix: mixing ratio in [kgkg-1]
372    """
373
374    fname = 'var_virtualTemp'
375
376    virtual=temp*(0.622+rmix)/(0.622*(1.+rmix))
377
378    return virtual
379
380
381def var_mslp(pres, psfc, ter, tk, qv):
382    """ Function to compute mslp on a 1D column
383    """
384
385    fname = 'var_mslp'
386
387    N = 1.0
388    expon=287.04*.0065/9.81
389    pref = 40000.
390
391# First find where about 400 hPa is located
392    dz=len(pres) 
393
394    kref = -1
395    pinc = pres[0] - pres[dz-1]
396
397    if pinc < 0.:
398        for iz in range(1,dz):
399            if pres[iz-1] >= pref and pres[iz] < pref: 
400                kref = iz
401                break
402    else:
403        for iz in range(dz-1):
404            if pres[iz] >= pref and pres[iz+1] < pref: 
405                kref = iz
406                break
407
408    if kref == -1:
409        print errormsg
410        print '  ' + fname + ': no reference pressure:',pref,'found!!'
411        print '    values:',pres[:]
412        quit(-1)
413
414    mslp = 0.
415
416# We are below both the ground and the lowest data level.
417
418# First, find the model level that is closest to a "target" pressure
419# level, where the "target" pressure is delta-p less that the local
420# value of a horizontally smoothed surface pressure field.  We use
421# delta-p = 150 hPa here. A standard lapse rate temperature profile
422# passing through the temperature at this model level will be used
423# to define the temperature profile below ground.  This is similar
424# to the Benjamin and Miller (1990) method, using 
425# 700 hPa everywhere for the "target" pressure.
426
427# ptarget = psfc - 15000.
428    ptarget = 70000.
429    dpmin=1.e4
430    kupper = 0
431    if pinc > 0.:
432        for iz in range(dz-1,0,-1):
433            kupper = iz
434            dp=np.abs( pres[iz] - ptarget )
435            if dp < dpmin: exit
436            dpmin = np.min([dpmin, dp])
437    else:
438        for iz in range(dz):
439            kupper = iz
440            dp=np.abs( pres[iz] - ptarget )
441            if dp < dpmin: exit
442            dpmin = np.min([dpmin, dp])
443
444    pbot=np.max([pres[0], psfc])
445#    zbot=0.
446
447#    tbotextrap=tk(i,j,kupper,itt)*(pbot/pres_field(i,j,kupper,itt))**expon
448#    tvbotextrap=virtual(tbotextrap,qv(i,j,1,itt))
449
450#    data_out(i,j,itt,1) = (zbot+tvbotextrap/.0065*(1.-(interp_levels(1)/pbot)**expon))
451    tbotextrap = tk[kupper]*(psfc/ptarget)**expon
452    tvbotextrap = var_virtualTemp(tbotextrap, qv[kupper])
453    mslp = psfc*( (tvbotextrap+0.0065*ter)/tvbotextrap)**(1./expon)
454
455    return mslp
456
457def compute_mslp(pressure, psurface, terrain, temperature, qvapor, dimns, dimvns):
458    """ Function to compute mslp: mean sea level pressure following p_interp.F90 from WRF
459    var_mslp(pres, ter, tk, qv, dimns, dimvns)
460      [pressure]= pressure field [Pa] (assuming [[t],z,y,x])
461      [psurface]= surface pressure field [Pa]
462      [terrain]= topography [m]
463      [temperature]= temperature [K]
464      [qvapor]= water vapour mixing ratio [kgkg-1]
465      [dimns]= list of the name of the dimensions of [cldfra]
466      [dimvns]= list of the name of the variables with the values of the
467        dimensions of [pres]
468    """
469
470    fname = 'compute_mslp'
471
472    mslpdims = list(dimns[:])
473    mslpvdims = list(dimvns[:])
474
475    if len(pressure.shape) == 4:
476        mslpdims.pop(1)
477        mslpvdims.pop(1)
478    else:
479        mslpdims.pop(0)
480        mslpvdims.pop(0)
481
482    if len(pressure.shape) == 4:
483        dx = pressure.shape[3]
484        dy = pressure.shape[2]
485        dz = pressure.shape[1]
486        dt = pressure.shape[0]
487
488        mslpv = np.zeros(tuple([dt, dy, dx]), dtype=np.float)
489
490# Terrain... to 2D !
491        terval = np.zeros(tuple([dy, dx]), dtype=np.float)
492        if len(terrain.shape) == 3:
493            terval = terrain[0,:,:]
494        else:
495            terval = terrain
496
497        for ix in range(dx):
498            for iy in range(dy):
499                if terval[iy,ix] > 0.:
500                    for it in range(dt):
501                        mslpv[it,iy,ix] = var_mslp(pressure[it,:,iy,ix],             \
502                          psurface[it,iy,ix], terval[iy,ix], temperature[it,:,iy,ix],\
503                          qvapor[it,:,iy,ix])
504
505                        ncvar.percendone(it*dx*dy + ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
506                else:
507                    mslpv[:,iy,ix] = psurface[:,iy,ix]
508
509    else:
510        dx = pressure.shape[2]
511        dy = pressure.shape[1]
512        dz = pressure.shape[0]
513
514        mslpv = np.zeros(tuple([dy, dx]), dtype=np.float)
515
516# Terrain... to 2D !
517        terval = np.zeros(tuple([dy, dx]), dtype=np.float)
518        if len(terrain.shape) == 3:
519            terval = terrain[0,:,:]
520        else:
521            terval = terrain
522
523        for ix in range(dx):
524            for iy in range(dy):
525                ncvar.percendone(ix*dy + iy,dx*dy, 5, 'diagnosted')
526                if terval[iy,ix] > 0.:
527                    mslpv[iy,ix] = var_mslp(pressure[:,iy,ix], psurface[iy,ix],          \
528                      terval[iy,ix], temperature[:,iy,ix], qvapor[:,iy,ix])
529                else:
530                    mslpv[iy,ix] = psfc[iy,ix]
531
532    return mslpv, mslpdims, mslpvdims
533
534def compute_prw(dens, q, dimns, dimvns):
535    """ Function to compute water vapour path (prw)
536      [dens] = density [in kgkg-1] (assuming [t],z,y,x)
537      [q] = mixing ratio in [kgkg-1] (assuming [t],z,y,x)
538      [dimns]= list of the name of the dimensions of [q]
539      [dimvns]= list of the name of the variables with the values of the
540        dimensions of [q]
541    """
542    fname = 'compute_prw'
543
544    prwdims = dimns[:]
545    prwvdims = dimvns[:]
546
547    if len(q.shape) == 4:
548        prwdims.pop(1)
549        prwvdims.pop(1)
550    else:
551        prwdims.pop(0)
552        prwvdims.pop(0)
553
554    data1 = dens*q
555    prw = np.sum(data1, axis=1)
556
557    return prw, prwdims, prwvdims
558
559def compute_rh(p, t, q, dimns, dimvns):
560    """ Function to compute relative humidity following 'Tetens' equation (T,P) ...'
561      [t]= temperature (assuming [[t],z,y,x] in [K])
562      [p] = pressure field (assuming in [hPa])
563      [q] = mixing ratio in [kgkg-1]
564      [dimns]= list of the name of the dimensions of [t]
565      [dimvns]= list of the name of the variables with the values of the
566        dimensions of [t]
567    """
568    fname = 'compute_rh'
569
570    rhdims = dimns[:]
571    rhvdims = dimvns[:]
572
573    data1 = 10.*0.6112*np.exp(17.67*(t-273.16)/(t-29.65))
574    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
575
576    rh = q/data2
577
578    return rh, rhdims, rhvdims
579
580def turbulence_var(varv, dimvn, dimn):
581    """ Function to compute the Taylor's decomposition turbulence term from a a given variable
582      x*=<x^2>_t-(<X>_t)^2
583    turbulence_var(varv,dimn)
584      varv= values of the variable
585      dimvn= names of the dimension of the variable
586      dimn= names of the dimensions (as a dictionary with 'X', 'Y', 'Z', 'T')
587    >>> turbulence_var(np.arange((27)).reshape(3,3,3),['time','y','x'],{'T':'time', 'Y':'y', 'X':'x'})
588    [[ 54.  54.  54.]
589     [ 54.  54.  54.]
590     [ 54.  54.  54.]]
591    """
592    fname = 'turbulence_varv'
593
594    timedimid = dimvn.index(dimn['T'])
595
596    varv2 = varv*varv
597
598    vartmean = np.mean(varv, axis=timedimid)
599    var2tmean = np.mean(varv2, axis=timedimid)
600
601    varvturb = var2tmean - (vartmean*vartmean)
602
603    return varvturb
604
605def compute_turbulence(v, dimns, dimvns):
606    """ Function to compute the rubulence term of the Taylor's decomposition ...'
607      x*=<x^2>_t-(<X>_t)^2
608      [v]= variable (assuming [[t],z,y,x])
609      [dimns]= list of the name of the dimensions of [v]
610      [dimvns]= list of the name of the variables with the values of the
611        dimensions of [v]
612    """
613    fname = 'compute_turbulence'
614
615    turbdims = dimns[:]
616    turbvdims = dimvns[:]
617
618    turbdims.pop(0)
619    turbvdims.pop(0)
620
621    v2 = v*v
622
623    vartmean = np.mean(v, axis=0)
624    var2tmean = np.mean(v2, axis=0)
625
626    turb = var2tmean - (vartmean*vartmean)
627
628    return turb, turbdims, turbvdims
629
630def timeunits_seconds(dtu):
631    """ Function to transform a time units to seconds
632    timeunits_seconds(timeuv)
633      [dtu]= time units value to transform in seconds
634    """
635    fname='timunits_seconds'
636
637    if dtu == 'years':
638        times = 365.*24.*3600.
639    elif dtu == 'weeks':
640        times = 7.*24.*3600.
641    elif dtu == 'days':
642        times = 24.*3600.
643    elif dtu == 'hours':
644        times = 3600.
645    elif dtu == 'minutes':
646        times = 60.
647    elif dtu == 'seconds':
648        times = 1.
649    elif dtu == 'miliseconds':
650        times = 1./1000.
651    else:
652        print errormsg
653        print '  ' + fname  + ": time units '" + dtu + "' not ready !!"
654        quit(-1)
655
656    return times
657
658####### ###### ##### #### ### ## #
659comboinf="\nIF -d 'variable_combo', provides information of the combination to obtain -v [varn] with the ASCII file with the combinations as -f [combofile]"
660
661parser = OptionParser()
662parser.add_option("-f", "--netCDF_file", dest="ncfile", help="file to use", metavar="FILE")
663parser.add_option("-d", "--dimensions", dest="dimns", 
664  help="[dimxn]@[dxvn],[dimyn]@[dxvn],[...,[dimtn]@[dxvn]], ',' list with the couples [dimDn]@[dDvn], [dimDn], name of the dimension D and name of the variable [dDvn] with the values of the dimension" + comboinf, 
665  metavar="LABELS")
666parser.add_option("-v", "--variables", dest="varns", 
667  help=" [varn1]|[var11]@[...[varN1]],[...,[varnM]|[var1M]@[...[varLM]]] ',' list of variables to compute [varnK] and its necessary ones [var1K]...[varPK]", metavar="VALUES")
668
669(opts, args) = parser.parse_args()
670
671#######    #######
672## MAIN
673    #######
674availdiags = ['ACRAINTOT', 'clt', 'cllmh', 'deaccum', 'LMDZrh', 'mslp', 'RAINTOT',   \
675  'rvors', 'turbulence', 'WRFrvors']
676
677# Variables not to check
678NONcheckingvars = ['cllmh', 'deaccum', 'WRFbils', 'WRFdens', 'WRFgeop', 'WRFp',      \
679  'WRFpos', 'WRFprc', 'WRFprls', 'WRFrh', 'LMDZrh', 'LMDZrhs', 'WRFrhs', 'WRFrvors', \
680  'WRFt', 'WRFtime']
681
682ofile = 'diagnostics.nc'
683
684dimns = opts.dimns
685varns = opts.varns
686
687# Special method. knowing variable combination
688##
689if opts.dimns == 'variable_combo':
690    print warnmsg
691    print '  ' + main + ': knowing variable combination !!!'
692    combination = variable_combo(opts.varns,opts.ncfile)
693    print '     COMBO: ' + combination
694    quit(-1)
695
696if not os.path.isfile(opts.ncfile):
697    print errormsg
698    print '   ' + main + ": file '" + opts.ncfile + "' does not exist !!"
699    quit(-1)
700
701ncobj = NetCDFFile(opts.ncfile, 'r')
702
703# File creation
704newnc = NetCDFFile(ofile,'w')
705
706# dimensions
707dimvalues = dimns.split(',')
708dnames = []
709dvnames = []
710
711for dimval in dimvalues:
712    dnames.append(dimval.split('@')[0])
713    dvnames.append(dimval.split('@')[1])
714
715# diagnostics to compute
716diags = varns.split(',')
717Ndiags = len(diags)
718
719# Looking for specific variables that might be use in more than one diagnostic
720WRFp_compute = False
721WRFt_compute = False
722WRFrh_compute = False
723WRFght_compute = False
724WRFdens_compute = False
725WRFpos_compute = False
726
727for idiag in range(Ndiags):
728    if diags[idiag].split('|')[1].find('@') == -1:
729        depvars = diags[idiag].split('|')[1]
730        if depvars == 'WRFp': WRFp_compute = True
731        if depvars == 'WRFt': WRFt_compute = True
732        if depvars == 'WRFrh': WRFrh_compute = True
733        if depvars == 'WRFght': WRFght_compute = True
734        if depvars == 'WRFdens': WRFdens_compute = True
735        if depvars == 'WRFpos': WRFpos_compute = True
736
737    else:
738        depvars = diags[idiag].split('|')[1].split('@')
739        if ncvar.searchInlist(depvars, 'WRFp'): WRFp_compute = True
740        if ncvar.searchInlist(depvars, 'WRFt'): WRFt_compute = True
741        if ncvar.searchInlist(depvars, 'WRFrh'): WRFrh_compute = True
742        if ncvar.searchInlist(depvars, 'WRFght'): WRFght_compute = True
743        if ncvar.searchInlist(depvars, 'WRFdens'): WRFdens_compute = True
744        if ncvar.searchInlist(depvars, 'WRFpos'): WRFpos_compute = True
745
746if WRFp_compute:
747    print '  ' + main + ': Retrieving pressure value from WRF as P + PB'
748    dimv = ncobj.variables['P'].shape
749    WRFp = ncobj.variables['P'][:] + ncobj.variables['PB'][:]
750
751if WRFght_compute:
752    print '    ' + main + ': computing geopotential height from WRF as PH + PHB ...' 
753    WRFght = ncobj.variables['PH'][:] + ncobj.variables['PHB'][:]
754
755if WRFrh_compute:
756    print '    ' + main + ": computing relative humidity from WRF as 'Tetens'" +     \
757      ' equation (T,P) ...'
758    p0=100000.
759    p=ncobj.variables['P'][:] + ncobj.variables['PB'][:]
760    tk = (ncobj.variables['T'][:] + 300.)*(p/p0)**(2./7.)
761    qv = ncobj.variables['QVAPOR'][:]
762
763    data1 = 10.*0.6112*np.exp(17.67*(tk-273.16)/(tk-29.65))
764    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
765
766    WRFrh = qv/data2
767
768if WRFt_compute:
769    print '    ' + main + ': computing temperature from WRF as inv_potT(T + 300) ...'
770    p0=100000.
771    p=ncobj.variables['P'][:] + ncobj.variables['PB'][:]
772
773    WRFt = (ncobj.variables['T'][:] + 300.)*(p/p0)**(2./7.)
774
775if WRFdens_compute:
776    print '    ' + main + ': computing air density from WRF as ((MU + MUB) * ' +     \
777      'DNW)/g ...'
778    grav = 9.81
779
780# Just we need in in absolute values: Size of the central grid cell
781##    dxval = ncobj.getncattr('DX')
782##    dyval = ncobj.getncattr('DY')
783##    mapfac = ncobj.variables['MAPFAC_M'][:]
784##    area = dxval*dyval*mapfac
785
786    mu = (ncobj.variables['MU'][:] + ncobj.variables['MUB'][:])
787    dnw = ncobj.variables['DNW'][:]
788
789    WRFdens = np.zeros((mu.shape[0], dnw.shape[1], mu.shape[1], mu.shape[2]),        \
790      dtype=np.float)
791    levval = np.zeros((mu.shape[1], mu.shape[2]), dtype=np.float)
792
793    for it in range(mu.shape[0]):
794        for iz in range(dnw.shape[1]):
795            levval.fill(np.abs(dnw[it,iz]))
796            WRFdens[it,iz,:,:] = levval
797            WRFdens[it,iz,:,:] = mu[it,:,:]*WRFdens[it,iz,:,:]/grav
798
799if WRFpos_compute:
800# WRF positions from the lowest-leftest corner of the matrix
801    print '    ' + main + ': computing position from MAPFAC_M as sqrt(DY*j**2 + ' +  \
802      'DX*x**2)*MAPFAC_M ...'
803
804    mapfac = ncobj.variables['MAPFAC_M'][:]
805
806    distx = np.float(ncobj.getncattr('DX'))
807    disty = np.float(ncobj.getncattr('DY'))
808
809    print 'distx:',distx,'disty:',disty
810
811    dx = mapfac.shape[2]
812    dy = mapfac.shape[1]
813    dt = mapfac.shape[0]
814
815    WRFpos = np.zeros((dt, dy, dx), dtype=np.float)
816
817    for i in range(1,dx):
818        WRFpos[0,0,i] = distx*i/mapfac[0,0,i]
819    for j in range(1,dy):
820        i=0
821        WRFpos[0,j,i] = WRFpos[0,j-1,i] + disty/mapfac[0,j,i]
822        for i in range(1,dx):
823#            WRFpos[0,j,i] = np.sqrt((disty*j)**2. + (distx*i)**2.)/mapfac[0,j,i]
824#            WRFpos[0,j,i] = np.sqrt((disty*j)**2. + (distx*i)**2.)
825             WRFpos[0,j,i] = WRFpos[0,j,i-1] + distx/mapfac[0,j,i]
826
827    for it in range(1,dt):
828        WRFpos[it,:,:] = WRFpos[0,:,:]
829
830### ## #
831# Going for the diagnostics
832### ## #
833print '  ' + main + ' ...'
834
835for idiag in range(Ndiags):
836    print '    diagnostic:',diags[idiag]
837    diag = diags[idiag].split('|')[0]
838    depvars = diags[idiag].split('|')[1].split('@')
839    if diags[idiag].split('|')[1].find('@') != -1:
840        depvars = diags[idiag].split('|')[1].split('@')
841        if depvars[0] == 'deaccum': diag='deaccum'
842        for depv in depvars:
843            if not ncobj.variables.has_key(depv) and not                             \
844              ncvar.searchInlist(NONcheckingvars, depv) and depvars[0] != 'deaccum':
845                print errormsg
846                print '  ' + main + ": file '" + opts.ncfile +                       \
847                  "' does not have variable '" + depv + "' !!"
848                quit(-1)
849    else:
850        depvars = diags[idiag].split('|')[1]
851        if not ncobj.variables.has_key(depvars) and not                              \
852          ncvar.searchInlist(NONcheckingvars, depvars) and depvars[0] != 'deaccum':
853            print errormsg
854            print '  ' + main + ": file '" + opts.ncfile +                           \
855              "' does not have variable '" + depvars + "' !!"
856            quit(-1)
857
858    print "\n    Computing '" + diag + "' from: ", depvars, '...'
859
860# acraintot: accumulated total precipitation from WRF RAINC, RAINNC
861    if diag == 'ACRAINTOT':
862           
863        var0 = ncobj.variables[depvars[0]]
864        var1 = ncobj.variables[depvars[1]]
865        diagout = var0[:] + var1[:]
866
867        dnamesvar = var0.dimensions
868        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
869
870        ncvar.insert_variable(ncobj, 'acpr', diagout, dnamesvar, dvnamesvar, newnc)
871
872# cllmh with cldfra, pres
873    elif diag == 'cllmh':
874           
875        var0 = ncobj.variables[depvars[0]]
876        if depvars[1] == 'WRFp':
877            var1 = WRFp
878        else:
879            var01 = ncobj.variables[depvars[1]]
880            if len(size(var1.shape)) < len(size(var0.shape)):
881                var1 = np.brodcast_arrays(var01,var0)[0]
882            else:
883                var1 = var01
884
885        diagout, diagoutd, diagoutvd = compute_cllmh(var0,var1,dnames,dvnames)
886        ncvar.insert_variable(ncobj, 'cll', diagout[0,:], diagoutd, diagoutvd, newnc)
887        ncvar.insert_variable(ncobj, 'clm', diagout[1,:], diagoutd, diagoutvd, newnc)
888        ncvar.insert_variable(ncobj, 'clh', diagout[2,:], diagoutd, diagoutvd, newnc)
889
890# clt with cldfra
891    elif diag == 'clt':
892           
893        var0 = ncobj.variables[depvars]
894        diagout, diagoutd, diagoutvd = compute_clt(var0,dnames,dvnames)
895        ncvar.insert_variable(ncobj, 'clt', diagout, diagoutd, diagoutvd, newnc)
896
897# deaccum: deacumulation of any variable as (Variable, time [as [tunits]
898#   from/since ....], newvarname)
899    elif diag == 'deaccum':
900
901        var0 = ncobj.variables[depvars[1]]
902        var1 = ncobj.variables[depvars[2]]
903
904        dnamesvar = var0.dimensions
905        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
906
907        diagout, diagoutd, diagoutvd = compute_deaccum(var0,dnamesvar,dvnamesvar)
908
909# Transforming to a flux
910        if depvars[2] == 'XTIME':
911            dtimeunits = var1.getncattr('description')
912            tunits = dtimeunits.split(' ')[0]
913        else:
914            dtimeunits = var1.getncattr('units')
915            tunits = dtimeunits.split(' ')[0]
916
917        dtime = (var1[1] - var1[0])*timeunits_seconds(tunits)
918        ncvar.insert_variable(ncobj, depvars[3], diagout/dtime, diagoutd, diagoutvd, newnc)
919
920# LMDZrh (pres, t, r)
921    elif diag == 'LMDZrh':
922           
923        var0 = ncobj.variables[depvars[0]][:]
924        var1 = ncobj.variables[depvars[1]][:]
925        var2 = ncobj.variables[depvars[2]][:]
926
927        diagout, diagoutd, diagoutvd = compute_rh(var0,var1,var2,dnames,dvnames)
928        ncvar.insert_variable(ncobj, 'hus', diagout, diagoutd, diagoutvd, newnc)
929
930# LMDZrhs (psol, t2m, q2m)
931    elif diag == 'LMDZrhs':
932           
933        var0 = ncobj.variables[depvars[0]][:]
934        var1 = ncobj.variables[depvars[1]][:]
935        var2 = ncobj.variables[depvars[2]][:]
936
937        dnamesvar = ncobj.variables[depvars[0]].dimensions
938        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
939
940        diagout, diagoutd, diagoutvd = compute_rh(var0,var1,var2,dnamesvar,dvnamesvar)
941
942        ncvar.insert_variable(ncobj, 'huss', diagout, diagoutd, diagoutvd, newnc)
943
944# mslp: mean sea level pressure (pres, psfc, terrain, temp, qv)
945    elif diag == 'mslp' or diag == 'WRFmslp':
946           
947        var1 = ncobj.variables[depvars[1]][:]
948        var2 = ncobj.variables[depvars[2]][:]
949        var4 = ncobj.variables[depvars[4]][:]
950
951        if diag == 'WRFmslp':
952            var0 = WRFp
953            var3 = WRFt
954            dnamesvar = ncobj.variables['P'].dimensions
955        else:
956            var0 = ncobj.variables[depvars[0]][:]
957            var3 = ncobj.variables[depvars[3]][:]
958            dnamesvar = ncobj.variables[depvars[0]].dimensions
959
960        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
961
962        diagout, diagoutd, diagoutvd = compute_mslp(var0, var1, var2, var3, var4,    \
963          dnamesvar, dvnamesvar)
964
965        ncvar.insert_variable(ncobj, 'psl', diagout, diagoutd, diagoutvd, newnc)
966
967# raintot: instantaneous total precipitation from WRF as (RAINC + RAINC) / dTime
968    elif diag == 'RAINTOT':
969
970        var0 = ncobj.variables[depvars[0]]
971        var1 = ncobj.variables[depvars[1]]
972        if depvars[2] != 'WRFtime':
973            var2 = ncobj.variables[depvars[2]]
974
975        var = var0[:] + var1[:]
976
977        dnamesvar = var0.dimensions
978        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
979
980        diagout, diagoutd, diagoutvd = compute_deaccum(var,dnamesvar,dvnamesvar)
981
982# Transforming to a flux
983        if depvars[2] != 'WRFtime':
984            dtimeunits = var2.getncattr('units')
985            tunits = dtimeunits.split(' ')[0]
986
987            dtime = (var2[1] - var2[0])*timeunits_seconds(tunits)
988        else:
989            var2 = ncobj.variables['Times']
990            time1 = var2[0,:]
991            time2 = var2[1,:]
992            tmf1 = ''
993            tmf2 = ''
994            for ic in range(len(time1)):
995                tmf1 = tmf1 + time1[ic]
996                tmf2 = tmf2 + time2[ic]
997            dtdate1 = dt.datetime.strptime(tmf1,"%Y-%m-%d_%H:%M:%S")
998            dtdate2 = dt.datetime.strptime(tmf2,"%Y-%m-%d_%H:%M:%S")
999            diffdate12 = dtdate2 - dtdate1
1000            dtime = diffdate12.total_seconds()
1001            print 'dtime:',dtime
1002
1003        ncvar.insert_variable(ncobj, 'pr', diagout/dtime, diagoutd, diagoutvd, newnc)
1004
1005# turbulence (var)
1006    elif diag == 'turbulence':
1007
1008        var0 = ncobj.variables[depvars][:]
1009
1010        dnamesvar = list(ncobj.variables[depvars].dimensions)
1011        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1012
1013        diagout, diagoutd, diagoutvd = compute_turbulence(var0,dnamesvar,dvnamesvar)
1014        valsvar = ncvar.variables_values(depvars)
1015
1016        ncvar.insert_variable(ncobj, valsvar[0] + 'turb', diagout, diagoutd, 
1017          diagoutvd, newnc)
1018        varobj = newnc.variables[valsvar[0] + 'turb']
1019        attrv = varobj.long_name
1020        attr = varobj.delncattr('long_name')
1021        newattr = ncvar.set_attribute(varobj, 'long_name', attrv +                   \
1022          " Taylor decomposition turbulence term")
1023
1024# WRFbils fom WRF as HFX + LH
1025    elif diag == 'WRFbils':
1026           
1027        var0 = ncobj.variables[depvars[0]][:]
1028        var1 = ncobj.variables[depvars[1]][:]
1029
1030        diagout = var0 + var1
1031
1032        ncvar.insert_variable(ncobj, 'bils', diagout, dnames, dvnames, newnc)
1033
1034# WRFp pressure from WRF as P + PB
1035    elif diag == 'WRFp':
1036           
1037        diagout = WRFp
1038
1039        ncvar.insert_variable(ncobj, 'pres', diagout, dnames, dvnames, newnc)
1040
1041# WRFpos
1042    elif diag == 'WRFpos':
1043           
1044        dnamesvar = ncobj.variables['MAPFAC_M'].dimensions
1045        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1046
1047        ncvar.insert_variable(ncobj, 'WRFpos', WRFpos, dnamesvar, dvnamesvar, newnc)
1048
1049# WRFprw WRF water vapour path WRFdens, QVAPOR
1050    elif diag == 'WRFprw':
1051           
1052        var0 = WRFdens
1053        var1 = ncobj.variables[depvars[1]]
1054
1055        dnamesvar = list(var1.dimensions)
1056        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1057
1058        diagout, diagoutd, diagoutvd = compute_prw(var0, var1, dnamesvar,dvnamesvar)
1059
1060        ncvar.insert_variable(ncobj, 'prw', diagout, diagoutd, diagoutvd, newnc)
1061
1062# WRFrh (P, T, QVAPOR)
1063    elif diag == 'WRFrh':
1064           
1065        dnamesvar = list(ncobj.variables[depvars[2]].dimensions)
1066        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1067
1068        ncvar.insert_variable(ncobj, 'hus', WRFrh, dnames, dvnames, newnc)
1069
1070# WRFrhs (PSFC, T2, Q2)
1071    elif diag == 'WRFrhs':
1072           
1073        var0 = ncobj.variables[depvars[0]][:]
1074        var1 = ncobj.variables[depvars[1]][:]
1075        var2 = ncobj.variables[depvars[2]][:]
1076
1077        dnamesvar = list(ncobj.variables[depvars[2]].dimensions)
1078        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1079
1080        diagout, diagoutd, diagoutvd = compute_rh(var0,var1,var2,dnamesvar,dvnamesvar)
1081        ncvar.insert_variable(ncobj, 'huss', diagout, diagoutd, diagoutvd, newnc)
1082
1083# rvors (u10, v10, WRFpos)
1084    elif diag == 'WRFrvors':
1085           
1086        var0 = ncobj.variables[depvars[0]]
1087        var1 = ncobj.variables[depvars[1]]
1088
1089        diagout = rotational_z(var0, var1, distx)
1090
1091        dnamesvar = ncobj.variables[depvars[0]].dimensions
1092        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1093
1094        ncvar.insert_variable(ncobj, 'rvors', diagout, dnamesvar, dvnamesvar, newnc)
1095
1096# wss (u10, v10)
1097    elif diag == 'wss':
1098           
1099        var0 = ncobj.variables[depvars[0]][:]
1100        var1 = ncobj.variables[depvars[1]][:]
1101
1102        diagout = np.sqrt(var0*var0 + var1*var1)
1103
1104        dnamesvar = ncobj.variables[depvars[0]].dimensions
1105        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1106
1107        print 'dnamesvar',dnamesvar
1108        print 'dnames',dnames
1109        print 'dvnames',dvnames
1110        print 'dvnamesvar',dvnamesvar
1111
1112        ncvar.insert_variable(ncobj, 'wss', diagout, dnamesvar, dvnamesvar, newnc)
1113
1114    else:
1115        print errormsg
1116        print '  ' + main + ": diagnostic '" + diag + "' not ready!!!"
1117        print '    available diagnostics: ', availdiags
1118        quit(-1)
1119
1120    newnc.sync()
1121
1122#   end of diagnostics
1123
1124# Global attributes
1125##
1126atvar = ncvar.set_attribute(newnc, 'program', 'diagnostics.py')
1127atvar = ncvar.set_attribute(newnc, 'version', '1.0')
1128atvar = ncvar.set_attribute(newnc, 'author', 'Fita Borrell, Lluis')
1129atvar = ncvar.set_attribute(newnc, 'institution', 'Laboratoire Meteorologie ' +      \
1130  'Dynamique')
1131atvar = ncvar.set_attribute(newnc, 'university', 'Universite Pierre et Marie ' +     \
1132  'Curie -- Jussieu')
1133atvar = ncvar.set_attribute(newnc, 'centre', 'Centre national de la recherche ' +    \
1134  'scientifique')
1135atvar = ncvar.set_attribute(newnc, 'city', 'Paris')
1136atvar = ncvar.set_attribute(newnc, 'original_file', opts.ncfile)
1137
1138gorigattrs = ncobj.ncattrs()
1139
1140for attr in gorigattrs:
1141    attrv = ncobj.getncattr(attr)
1142    atvar = ncvar.set_attribute(newnc, attr, attrv)
1143
1144ncobj.close()
1145newnc.close()
1146
1147print '\n' + main + ': successfull writting of diagnostics file "' + ofile + '" !!!'
Note: See TracBrowser for help on using the repository browser.