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

Last change on this file since 1647 was 1647, checked in by lfita, 7 years ago

Adding removal of not desired variables from 'ACRAINTOT'

File size: 64.3 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#   To be used with module_ForDiagnostics.F90, module_ForDiagnosticsVars.F90, module_generic.F90
5#      foudre: f2py -m module_ForDiagnostics --f90exec=/usr/bin/gfortran-4.7 -c module_generic.F90 module_ForDiagnosticsVars.F90 module_ForDiagnostics.F90 >& run_f2py.log
6#      ciclad: f2py --f90flags="-fPIC" --f90exec=/usr/bin/gfortran -L/opt/canopy-1.3.0/Canopy_64bit/System/lib/ -L/usr/lib64/ -L/opt/canopy-1.3.0/Canopy_64bit/System/lib/ -m module_ForDiagnostics -c module_generic.F90 module_ForDiagnosticsVars.F90 module_ForDiagnostics.F90 >& run_f2py.log
7
8## 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
9## e.g. # diagnostics.py -f /home/lluis/PY/diagnostics.inf -d variable_combo -v WRFprc
10
11from optparse import OptionParser
12import numpy as np
13from netCDF4 import Dataset as NetCDFFile
14import os
15import re
16import nc_var_tools as ncvar
17import generic_tools as gen
18import datetime as dtime
19import module_ForDiag as fdin
20
21main = 'diagnostics.py'
22errormsg = 'ERROR -- error -- ERROR -- error'
23warnmsg = 'WARNING -- warning -- WARNING -- warning'
24
25# Constants
26grav = 9.81
27
28# Gneral information
29##
30def reduce_spaces(string):
31    """ Function to give words of a line of text removing any extra space
32    """
33    values = string.replace('\n','').split(' ')
34    vals = []
35    for val in values:
36         if len(val) > 0:
37             vals.append(val)
38
39    return vals
40
41def variable_combo(varn,combofile):
42    """ Function to provide variables combination from a given variable name
43      varn= name of the variable
44      combofile= ASCII file with the combination of variables
45        [varn] [combo]
46          [combo]: '@' separated list of variables to use to generate [varn]
47            [WRFdt] to get WRF time-step (from general attributes)
48    >>> variable_combo('WRFprls','/home/lluis/PY/diagnostics.inf')
49    deaccum@RAINNC@XTIME@prnc
50    """
51    fname = 'variable_combo'
52
53    if varn == 'h':
54        print fname + '_____________________________________________________________'
55        print variable_combo.__doc__
56        quit()
57
58    if not os.path.isfile(combofile):
59        print errormsg
60        print '  ' + fname + ": file with combinations '" + combofile +              \
61          "' does not exist!!"
62        quit(-1)
63
64    objf = open(combofile, 'r')
65
66    found = False
67    for line in objf:
68        linevals = reduce_spaces(line)
69        varnf = linevals[0]
70        combo = linevals[1].replace('\n','')
71        if varn == varnf: 
72            found = True
73            break
74
75    if not found:
76        print errormsg
77        print '  ' + fname + ": variable '" + varn + "' not found in '" + combofile +\
78          "' !!"
79        combo='ERROR'
80
81    objf.close()
82
83    return combo
84
85# Mathematical operators
86##
87def compute_accum(varv, dimns, dimvns):
88    """ Function to compute the accumulation of a variable
89    compute_accum(varv, dimnames, dimvns)
90      [varv]= values to accum (assuming [t,])
91      [dimns]= list of the name of the dimensions of the [varv]
92      [dimvns]= list of the name of the variables with the values of the
93        dimensions of [varv]
94    """
95    fname = 'compute_accum'
96
97    deacdims = dimns[:]
98    deacvdims = dimvns[:]
99
100    slicei = []
101    slicee = []
102
103    Ndims = len(varv.shape)
104    for iid in range(0,Ndims):
105        slicei.append(slice(0,varv.shape[iid]))
106        slicee.append(slice(0,varv.shape[iid]))
107
108    slicee[0] = np.arange(varv.shape[0])
109    slicei[0] = np.arange(varv.shape[0])
110    slicei[0][1:varv.shape[0]] = np.arange(varv.shape[0]-1)
111
112    vari = varv[tuple(slicei)]
113    vare = varv[tuple(slicee)]
114
115    ac = vari*0.
116    for it in range(1,varv.shape[0]):
117        ac[it,] = ac[it-1,] + vare[it,]
118
119    return ac, deacdims, deacvdims
120
121def compute_deaccum(varv, dimns, dimvns):
122    """ Function to compute the deaccumulation of a variable
123    compute_deaccum(varv, dimnames, dimvns)
124      [varv]= values to deaccum (assuming [t,])
125      [dimns]= list of the name of the dimensions of the [varv]
126      [dimvns]= list of the name of the variables with the values of the
127        dimensions of [varv]
128    """
129    fname = 'compute_deaccum'
130
131    deacdims = dimns[:]
132    deacvdims = dimvns[:]
133
134    slicei = []
135    slicee = []
136
137    Ndims = len(varv.shape)
138    for iid in range(0,Ndims):
139        slicei.append(slice(0,varv.shape[iid]))
140        slicee.append(slice(0,varv.shape[iid]))
141
142    slicee[0] = np.arange(varv.shape[0])
143    slicei[0] = np.arange(varv.shape[0])
144    slicei[0][1:varv.shape[0]] = np.arange(varv.shape[0]-1)
145
146    vari = varv[tuple(slicei)]
147    vare = varv[tuple(slicee)]
148
149    deac = vare - vari
150
151    return deac, deacdims, deacvdims
152
153def derivate_centered(var,dim,dimv):
154    """ Function to compute the centered derivate of a given field
155      centered derivate(n) = (var(n-1) + var(n+1))/(2*dn).
156    [var]= variable
157    [dim]= which dimension to compute the derivate
158    [dimv]= dimension values (can be of different dimension of [var])
159    >>> derivate_centered(np.arange(16).reshape(4,4)*1.,1,1.)
160    [[  0.   1.   2.   0.]
161     [  0.   5.   6.   0.]
162     [  0.   9.  10.   0.]
163     [  0.  13.  14.   0.]]
164    """
165
166    fname = 'derivate_centered'
167   
168    vark = var.dtype
169
170    if hasattr(dimv, "__len__"):
171# Assuming that the last dimensions of var [..., N, M] are the same of dimv [N, M]
172        if len(var.shape) != len(dimv.shape):
173            dimvals = np.zeros((var.shape), dtype=vark)
174            if len(var.shape) - len(dimv.shape) == 1:
175                for iz in range(var.shape[0]):
176                    dimvals[iz,] = dimv
177            elif len(var.shape) - len(dimv.shape) == 2:
178                for it in range(var.shape[0]):
179                    for iz in range(var.shape[1]):
180                        dimvals[it,iz,] = dimv
181            else:
182                print errormsg
183                print '  ' + fname + ': dimension difference between variable',      \
184                  var.shape,'and variable with dimension values',dimv.shape,         \
185                  ' not ready !!!'
186                quit(-1)
187        else:
188            dimvals = dimv
189    else:
190# dimension values are identical everywhere!
191# from: http://stackoverflow.com/questions/16807011/python-how-to-identify-if-a-variable-is-an-array-or-a-scalar   
192        dimvals = np.ones((var.shape), dtype=vark)*dimv
193
194    derivate = np.zeros((var.shape), dtype=vark)
195    if dim > len(var.shape) - 1:
196        print errormsg
197        print '  ' + fname + ': dimension',dim,' too big for given variable of ' +   \
198          'shape:', var.shape,'!!!'
199        quit(-1)
200
201    slicebef = []
202    sliceaft = []
203    sliceder = []
204
205    for id in range(len(var.shape)):
206        if id == dim:
207            slicebef.append(slice(0,var.shape[id]-2))
208            sliceaft.append(slice(2,var.shape[id]))
209            sliceder.append(slice(1,var.shape[id]-1))
210        else:
211            slicebef.append(slice(0,var.shape[id]))
212            sliceaft.append(slice(0,var.shape[id]))
213            sliceder.append(slice(0,var.shape[id]))
214
215    if hasattr(dimv, "__len__"):
216        derivate[tuple(sliceder)] = (var[tuple(slicebef)] + var[tuple(sliceaft)])/   \
217          ((dimvals[tuple(sliceaft)] - dimvals[tuple(slicebef)]))
218        print (dimvals[tuple(sliceaft)] - dimvals[tuple(slicebef)])
219    else:
220        derivate[tuple(sliceder)] = (var[tuple(slicebef)] + var[tuple(sliceaft)])/   \
221          (2.*dimv)
222
223#    print 'before________'
224#    print var[tuple(slicebef)]
225
226#    print 'after________'
227#    print var[tuple(sliceaft)]
228
229    return derivate
230
231def rotational_z(Vx,Vy,pos):
232    """ z-component of the rotatinoal of horizontal vectorial field
233    \/ x (Vx,Vy,Vz) = \/xVy - \/yVx
234    [Vx]= Variable component x
235    [Vy]=  Variable component y
236    [pos]= poisition of the grid points
237    >>> rotational_z(np.arange(16).reshape(4,4)*1., np.arange(16).reshape(4,4)*1., 1.)
238    [[  0.   1.   2.   0.]
239     [ -4.   0.   0.  -7.]
240     [ -8.   0.   0. -11.]
241     [  0.  13.  14.   0.]]
242    """
243
244    fname =  'rotational_z'
245
246    ndims = len(Vx.shape)
247    rot1 = derivate_centered(Vy,ndims-1,pos)
248    rot2 = derivate_centered(Vx,ndims-2,pos)
249
250    rot = rot1 - rot2
251
252    return rot
253
254# Diagnostics
255##
256
257def var_clt(cfra):
258    """ Function to compute the total cloud fraction following 'newmicro.F90' from
259      LMDZ using 1D vertical column values
260      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
261    """
262    ZEPSEC=1.0E-12
263
264    fname = 'var_clt'
265
266    zclear = 1.
267    zcloud = 0.
268
269    dz = cfra.shape[0]
270    for iz in range(dz):
271        zclear =zclear*(1.-np.max([cfra[iz],zcloud]))/(1.-np.min([zcloud,1.-ZEPSEC]))
272        clt = 1. - zclear
273        zcloud = cfra[iz]
274
275    return clt
276
277def compute_clt(cldfra, dimns, dimvns):
278    """ Function to compute the total cloud fraction following 'newmicro.F90' from
279      LMDZ
280    compute_clt(cldfra, dimnames)
281      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
282      [dimns]= list of the name of the dimensions of [cldfra]
283      [dimvns]= list of the name of the variables with the values of the
284        dimensions of [cldfra]
285    """
286    fname = 'compute_clt'
287
288    cltdims = dimns[:]
289    cltvdims = dimvns[:]
290
291    if len(cldfra.shape) == 4:
292        clt = np.zeros((cldfra.shape[0],cldfra.shape[2],cldfra.shape[3]),            \
293          dtype=np.float)
294        dx = cldfra.shape[3]
295        dy = cldfra.shape[2]
296        dz = cldfra.shape[1]
297        dt = cldfra.shape[0]
298        cltdims.pop(1)
299        cltvdims.pop(1)
300
301        for it in range(dt):
302            for ix in range(dx):
303                for iy in range(dy):
304                    zclear = 1.
305                    zcloud = 0.
306                    gen.percendone(it*dx*dy + ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
307                    clt[it,iy,ix] = var_clt(cldfra[it,:,iy,ix])
308
309    else:
310        clt = np.zeros((cldfra.shape[1],cldfra.shape[2]), dtype=np.float)
311        dx = cldfra.shape[2]
312        dy = cldfra.shape[1]
313        dy = cldfra.shape[0]
314        cltdims.pop(0)
315        cltvdims.pop(0)
316        for ix in range(dx):
317            for iy in range(dy):
318                zclear = 1.
319                zcloud = 0.
320                gen.percendone(ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
321                clt[iy,ix] = var_clt(cldfra[:,iy,ix])
322
323    return clt, cltdims, cltvdims
324
325def Forcompute_clt(cldfra, dimns, dimvns):
326    """ Function to compute the total cloud fraction following 'newmicro.F90' from
327      LMDZ via a Fortran module
328    compute_clt(cldfra, dimnames)
329      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
330      [dimns]= list of the name of the dimensions of [cldfra]
331      [dimvns]= list of the name of the variables with the values of the
332        dimensions of [cldfra]
333    """
334    fname = 'Forcompute_clt'
335
336    cltdims = dimns[:]
337    cltvdims = dimvns[:]
338
339
340    if len(cldfra.shape) == 4:
341        clt = np.zeros((cldfra.shape[0],cldfra.shape[2],cldfra.shape[3]),            \
342          dtype=np.float)
343        dx = cldfra.shape[3]
344        dy = cldfra.shape[2]
345        dz = cldfra.shape[1]
346        dt = cldfra.shape[0]
347        cltdims.pop(1)
348        cltvdims.pop(1)
349
350        clt = fdin.module_fordiagnostics.compute_clt4d2(cldfra[:])
351
352    else:
353        clt = np.zeros((cldfra.shape[1],cldfra.shape[2]), dtype=np.float)
354        dx = cldfra.shape[2]
355        dy = cldfra.shape[1]
356        dy = cldfra.shape[0]
357        cltdims.pop(0)
358        cltvdims.pop(0)
359
360        clt = fdin.module_fordiagnostics.compute_clt3d1(cldfra[:])
361
362    return clt, cltdims, cltvdims
363
364def var_cllmh(cfra, p):
365    """ Fcuntion to compute cllmh on a 1D column
366    """
367
368    fname = 'var_cllmh'
369
370    ZEPSEC =1.0E-12
371    prmhc = 440.*100.
372    prmlc = 680.*100.
373
374    zclearl = 1.
375    zcloudl = 0.
376    zclearm = 1.
377    zcloudm = 0.
378    zclearh = 1.
379    zcloudh = 0.
380
381    dvz = cfra.shape[0]
382
383    cllmh = np.ones((3), dtype=np.float)
384
385    for iz in range(dvz):
386        if p[iz] < prmhc:
387            cllmh[2] = cllmh[2]*(1.-np.max([cfra[iz], zcloudh]))/(1.-                \
388              np.min([zcloudh,1.-ZEPSEC]))
389            zcloudh = cfra[iz]
390        elif p[iz] >= prmhc and p[iz] < prmlc:
391            cllmh[1] = cllmh[1]*(1.-np.max([cfra[iz], zcloudm]))/(1.-                \
392              np.min([zcloudm,1.-ZEPSEC]))
393            zcloudm = cfra[iz]
394        elif p[iz] >= prmlc:
395            cllmh[0] = cllmh[0]*(1.-np.max([cfra[iz], zcloudl]))/(1.-                \
396              np.min([zcloudl,1.-ZEPSEC]))
397            zcloudl = cfra[iz]
398
399    cllmh = 1.- cllmh
400
401    return cllmh
402
403def Forcompute_cllmh(cldfra, pres, dimns, dimvns):
404    """ Function to compute cllmh: low/medium/hight cloud fraction following newmicro.F90 from LMDZ via Fortran subroutine
405    compute_clt(cldfra, pres, dimns, dimvns)
406      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
407      [pres] = pressure field
408      [dimns]= list of the name of the dimensions of [cldfra]
409      [dimvns]= list of the name of the variables with the values of the
410        dimensions of [cldfra]
411    """
412    fname = 'Forcompute_cllmh'
413
414    cllmhdims = dimns[:]
415    cllmhvdims = dimvns[:]
416
417    if len(cldfra.shape) == 4:
418        dx = cldfra.shape[3]
419        dy = cldfra.shape[2]
420        dz = cldfra.shape[1]
421        dt = cldfra.shape[0]
422        cllmhdims.pop(1)
423        cllmhvdims.pop(1)
424
425        cllmh = fdin.module_fordiagnostics.compute_cllmh4d2(cldfra[:], pres[:])
426       
427    else:
428        dx = cldfra.shape[2]
429        dy = cldfra.shape[1]
430        dz = cldfra.shape[0]
431        cllmhdims.pop(0)
432        cllmhvdims.pop(0)
433
434        cllmh = fdin.module_fordiagnostics.compute_cllmh3d1(cldfra[:], pres[:])
435
436    return cllmh, cllmhdims, cllmhvdims
437
438def compute_cllmh(cldfra, pres, dimns, dimvns):
439    """ Function to compute cllmh: low/medium/hight cloud fraction following newmicro.F90 from LMDZ
440    compute_clt(cldfra, pres, dimns, dimvns)
441      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
442      [pres] = pressure field
443      [dimns]= list of the name of the dimensions of [cldfra]
444      [dimvns]= list of the name of the variables with the values of the
445        dimensions of [cldfra]
446    """
447    fname = 'compute_cllmh'
448
449    cllmhdims = dimns[:]
450    cllmhvdims = dimvns[:]
451
452    if len(cldfra.shape) == 4:
453        dx = cldfra.shape[3]
454        dy = cldfra.shape[2]
455        dz = cldfra.shape[1]
456        dt = cldfra.shape[0]
457        cllmhdims.pop(1)
458        cllmhvdims.pop(1)
459
460        cllmh = np.ones(tuple([3, dt, dy, dx]), dtype=np.float)
461
462        for it in range(dt):
463            for ix in range(dx):
464                for iy in range(dy):
465                    gen.percendone(it*dx*dy + ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
466                    cllmh[:,it,iy,ix] = var_cllmh(cldfra[it,:,iy,ix], pres[it,:,iy,ix])
467       
468    else:
469        dx = cldfra.shape[2]
470        dy = cldfra.shape[1]
471        dz = cldfra.shape[0]
472        cllmhdims.pop(0)
473        cllmhvdims.pop(0)
474
475        cllmh = np.ones(tuple([3, dy, dx]), dtype=np.float)
476
477        for ix in range(dx):
478            for iy in range(dy):
479                gen.percendone(ix*dy + iy,dx*dy, 5, 'diagnosted')
480                cllmh[:,iy,ix] = var_cllmh(cldfra[:,iy,ix], pres[:,iy,ix])
481
482    return cllmh, cllmhdims, cllmhvdims
483
484def compute_clivi(dens, qtot, dimns, dimvns):
485    """ Function to compute cloud-ice water path (clivi)
486      [dens] = density [in kgkg-1] (assuming [t],z,y,x)
487      [qtot] = added mixing ratio of all cloud-ice species in [kgkg-1] (assuming [t],z,y,x)
488      [dimns]= list of the name of the dimensions of [q]
489      [dimvns]= list of the name of the variables with the values of the
490        dimensions of [q]
491    """
492    fname = 'compute_clivi'
493
494    clividims = dimns[:]
495    clivivdims = dimvns[:]
496
497    if len(qtot.shape) == 4:
498        clividims.pop(1)
499        clivivdims.pop(1)
500    else:
501        clividims.pop(0)
502        clivivdims.pop(0)
503
504    data1 = dens*qtot
505    clivi = np.sum(data1, axis=1)
506
507    return clivi, clividims, clivivdims
508
509
510def compute_clwvl(dens, qtot, dimns, dimvns):
511    """ Function to compute condensed water path (clwvl)
512      [dens] = density [in kgkg-1] (assuming [t],z,y,x)
513      [qtot] = added mixing ratio of all cloud-water species in [kgkg-1] (assuming [t],z,y,x)
514      [dimns]= list of the name of the dimensions of [q]
515      [dimvns]= list of the name of the variables with the values of the
516        dimensions of [q]
517    """
518    fname = 'compute_clwvl'
519
520    clwvldims = dimns[:]
521    clwvlvdims = dimvns[:]
522
523    if len(qtot.shape) == 4:
524        clwvldims.pop(1)
525        clwvlvdims.pop(1)
526    else:
527        clwvldims.pop(0)
528        clwvlvdims.pop(0)
529
530    data1 = dens*qtot
531    clwvl = np.sum(data1, axis=1)
532
533    return clwvl, clwvldims, clwvlvdims
534
535def var_virtualTemp (temp,rmix):
536    """ This function returns virtual temperature in K,
537      temp: temperature [K]
538      rmix: mixing ratio in [kgkg-1]
539    """
540
541    fname = 'var_virtualTemp'
542
543    virtual=temp*(0.622+rmix)/(0.622*(1.+rmix))
544
545    return virtual
546
547
548def var_mslp(pres, psfc, ter, tk, qv):
549    """ Function to compute mslp on a 1D column
550    """
551
552    fname = 'var_mslp'
553
554    N = 1.0
555    expon=287.04*.0065/9.81
556    pref = 40000.
557
558# First find where about 400 hPa is located
559    dz=len(pres) 
560
561    kref = -1
562    pinc = pres[0] - pres[dz-1]
563
564    if pinc < 0.:
565        for iz in range(1,dz):
566            if pres[iz-1] >= pref and pres[iz] < pref: 
567                kref = iz
568                break
569    else:
570        for iz in range(dz-1):
571            if pres[iz] >= pref and pres[iz+1] < pref: 
572                kref = iz
573                break
574
575    if kref == -1:
576        print errormsg
577        print '  ' + fname + ': no reference pressure:',pref,'found!!'
578        print '    values:',pres[:]
579        quit(-1)
580
581    mslp = 0.
582
583# We are below both the ground and the lowest data level.
584
585# First, find the model level that is closest to a "target" pressure
586# level, where the "target" pressure is delta-p less that the local
587# value of a horizontally smoothed surface pressure field.  We use
588# delta-p = 150 hPa here. A standard lapse rate temperature profile
589# passing through the temperature at this model level will be used
590# to define the temperature profile below ground.  This is similar
591# to the Benjamin and Miller (1990) method, using 
592# 700 hPa everywhere for the "target" pressure.
593
594# ptarget = psfc - 15000.
595    ptarget = 70000.
596    dpmin=1.e4
597    kupper = 0
598    if pinc > 0.:
599        for iz in range(dz-1,0,-1):
600            kupper = iz
601            dp=np.abs( pres[iz] - ptarget )
602            if dp < dpmin: exit
603            dpmin = np.min([dpmin, dp])
604    else:
605        for iz in range(dz):
606            kupper = iz
607            dp=np.abs( pres[iz] - ptarget )
608            if dp < dpmin: exit
609            dpmin = np.min([dpmin, dp])
610
611    pbot=np.max([pres[0], psfc])
612#    zbot=0.
613
614#    tbotextrap=tk(i,j,kupper,itt)*(pbot/pres_field(i,j,kupper,itt))**expon
615#    tvbotextrap=virtual(tbotextrap,qv(i,j,1,itt))
616
617#    data_out(i,j,itt,1) = (zbot+tvbotextrap/.0065*(1.-(interp_levels(1)/pbot)**expon))
618    tbotextrap = tk[kupper]*(psfc/ptarget)**expon
619    tvbotextrap = var_virtualTemp(tbotextrap, qv[kupper])
620    mslp = psfc*( (tvbotextrap+0.0065*ter)/tvbotextrap)**(1./expon)
621
622    return mslp
623
624def compute_mslp(pressure, psurface, terrain, temperature, qvapor, dimns, dimvns):
625    """ Function to compute mslp: mean sea level pressure following p_interp.F90 from WRF
626    var_mslp(pres, ter, tk, qv, dimns, dimvns)
627      [pressure]= pressure field [Pa] (assuming [[t],z,y,x])
628      [psurface]= surface pressure field [Pa]
629      [terrain]= topography [m]
630      [temperature]= temperature [K]
631      [qvapor]= water vapour mixing ratio [kgkg-1]
632      [dimns]= list of the name of the dimensions of [cldfra]
633      [dimvns]= list of the name of the variables with the values of the
634        dimensions of [pres]
635    """
636
637    fname = 'compute_mslp'
638
639    mslpdims = list(dimns[:])
640    mslpvdims = list(dimvns[:])
641
642    if len(pressure.shape) == 4:
643        mslpdims.pop(1)
644        mslpvdims.pop(1)
645    else:
646        mslpdims.pop(0)
647        mslpvdims.pop(0)
648
649    if len(pressure.shape) == 4:
650        dx = pressure.shape[3]
651        dy = pressure.shape[2]
652        dz = pressure.shape[1]
653        dt = pressure.shape[0]
654
655        mslpv = np.zeros(tuple([dt, dy, dx]), dtype=np.float)
656
657# Terrain... to 2D !
658        terval = np.zeros(tuple([dy, dx]), dtype=np.float)
659        if len(terrain.shape) == 3:
660            terval = terrain[0,:,:]
661        else:
662            terval = terrain
663
664        for ix in range(dx):
665            for iy in range(dy):
666                if terval[iy,ix] > 0.:
667                    for it in range(dt):
668                        mslpv[it,iy,ix] = var_mslp(pressure[it,:,iy,ix],             \
669                          psurface[it,iy,ix], terval[iy,ix], temperature[it,:,iy,ix],\
670                          qvapor[it,:,iy,ix])
671
672                        gen.percendone(it*dx*dy + ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
673                else:
674                    mslpv[:,iy,ix] = psurface[:,iy,ix]
675
676    else:
677        dx = pressure.shape[2]
678        dy = pressure.shape[1]
679        dz = pressure.shape[0]
680
681        mslpv = np.zeros(tuple([dy, dx]), dtype=np.float)
682
683# Terrain... to 2D !
684        terval = np.zeros(tuple([dy, dx]), dtype=np.float)
685        if len(terrain.shape) == 3:
686            terval = terrain[0,:,:]
687        else:
688            terval = terrain
689
690        for ix in range(dx):
691            for iy in range(dy):
692                gen.percendone(ix*dy + iy,dx*dy, 5, 'diagnosted')
693                if terval[iy,ix] > 0.:
694                    mslpv[iy,ix] = var_mslp(pressure[:,iy,ix], psurface[iy,ix],          \
695                      terval[iy,ix], temperature[:,iy,ix], qvapor[:,iy,ix])
696                else:
697                    mslpv[iy,ix] = psfc[iy,ix]
698
699    return mslpv, mslpdims, mslpvdims
700
701def compute_OMEGAw(omega, p, t, dimns, dimvns):
702    """ Function to transform OMEGA [Pas-1] to velocities [ms-1]
703      tacking: https://www.ncl.ucar.edu/Document/Functions/Contributed/omega_to_w.shtml
704      [omega] = vertical velocity [in ms-1] (assuming [t],z,y,x)
705      [p] = pressure in [Pa] (assuming [t],z,y,x)
706      [t] = temperature in [K] (assuming [t],z,y,x)
707      [dimns]= list of the name of the dimensions of [q]
708      [dimvns]= list of the name of the variables with the values of the
709        dimensions of [q]
710    """
711    fname = 'compute_OMEGAw'
712
713    rgas = 287.058            # J/(kg-K) => m2/(s2 K)
714    g    = 9.80665            # m/s2
715
716    wdims = dimns[:]
717    wvdims = dimvns[:]
718
719    rho  = p/(rgas*t)         # density => kg/m3
720    w    = -omega/(rho*g)     
721
722    return w, wdims, wvdims
723
724def compute_prw(dens, q, dimns, dimvns):
725    """ Function to compute water vapour path (prw)
726      [dens] = density [in kgkg-1] (assuming [t],z,y,x)
727      [q] = mixing ratio in [kgkg-1] (assuming [t],z,y,x)
728      [dimns]= list of the name of the dimensions of [q]
729      [dimvns]= list of the name of the variables with the values of the
730        dimensions of [q]
731    """
732    fname = 'compute_prw'
733
734    prwdims = dimns[:]
735    prwvdims = dimvns[:]
736
737    if len(q.shape) == 4:
738        prwdims.pop(1)
739        prwvdims.pop(1)
740    else:
741        prwdims.pop(0)
742        prwvdims.pop(0)
743
744    data1 = dens*q
745    prw = np.sum(data1, axis=1)
746
747    return prw, prwdims, prwvdims
748
749def compute_rh(p, t, q, dimns, dimvns):
750    """ Function to compute relative humidity following 'Tetens' equation (T,P) ...'
751      [t]= temperature (assuming [[t],z,y,x] in [K])
752      [p] = pressure field (assuming in [hPa])
753      [q] = mixing ratio in [kgkg-1]
754      [dimns]= list of the name of the dimensions of [t]
755      [dimvns]= list of the name of the variables with the values of the
756        dimensions of [t]
757    """
758    fname = 'compute_rh'
759
760    rhdims = dimns[:]
761    rhvdims = dimvns[:]
762
763    data1 = 10.*0.6112*np.exp(17.67*(t-273.16)/(t-29.65))
764    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
765
766    rh = q/data2
767
768    return rh, rhdims, rhvdims
769
770def compute_td(p, temp, qv, dimns, dimvns):
771    """ Function to compute the dew point temperature
772      [p]= pressure [Pa]
773      [temp]= temperature [C]
774      [qv]= mixing ratio [kgkg-1]
775      [dimns]= list of the name of the dimensions of [p]
776      [dimvns]= list of the name of the variables with the values of the
777        dimensions of [p]
778    """
779    fname = 'compute_td'
780
781#    print '    ' + fname + ': computing dew-point temperature from TS as t and Tetens...'
782# tacking from: http://en.wikipedia.org/wiki/Dew_point
783    tk = temp
784    data1 = 10.*0.6112*np.exp(17.67*(tk-273.16)/(tk-29.65))
785    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
786
787    rh = qv/data2
788               
789    pa = rh * data1
790    td = 257.44*np.log(pa/6.1121)/(18.678-np.log(pa/6.1121))
791
792    tddims = dimns[:]
793    tdvdims = dimvns[:]
794
795    return td, tddims, tdvdims
796
797def turbulence_var(varv, dimvn, dimn):
798    """ Function to compute the Taylor's decomposition turbulence term from a a given variable
799      x*=<x^2>_t-(<X>_t)^2
800    turbulence_var(varv,dimn)
801      varv= values of the variable
802      dimvn= names of the dimension of the variable
803      dimn= names of the dimensions (as a dictionary with 'X', 'Y', 'Z', 'T')
804    >>> turbulence_var(np.arange((27)).reshape(3,3,3),['time','y','x'],{'T':'time', 'Y':'y', 'X':'x'})
805    [[ 54.  54.  54.]
806     [ 54.  54.  54.]
807     [ 54.  54.  54.]]
808    """
809    fname = 'turbulence_varv'
810
811    timedimid = dimvn.index(dimn['T'])
812
813    varv2 = varv*varv
814
815    vartmean = np.mean(varv, axis=timedimid)
816    var2tmean = np.mean(varv2, axis=timedimid)
817
818    varvturb = var2tmean - (vartmean*vartmean)
819
820    return varvturb
821
822def compute_turbulence(v, dimns, dimvns):
823    """ Function to compute the rubulence term of the Taylor's decomposition ...'
824      x*=<x^2>_t-(<X>_t)^2
825      [v]= variable (assuming [[t],z,y,x])
826      [dimns]= list of the name of the dimensions of [v]
827      [dimvns]= list of the name of the variables with the values of the
828        dimensions of [v]
829    """
830    fname = 'compute_turbulence'
831
832    turbdims = dimns[:]
833    turbvdims = dimvns[:]
834
835    turbdims.pop(0)
836    turbvdims.pop(0)
837
838    v2 = v*v
839
840    vartmean = np.mean(v, axis=0)
841    var2tmean = np.mean(v2, axis=0)
842
843    turb = var2tmean - (vartmean*vartmean)
844
845    return turb, turbdims, turbvdims
846
847def compute_wds(u, v, dimns, dimvns):
848    """ Function to compute the wind direction
849      [u]= W-E wind direction [ms-1, knot, ...]
850      [v]= N-S wind direction [ms-1, knot, ...]
851      [dimns]= list of the name of the dimensions of [u]
852      [dimvns]= list of the name of the variables with the values of the
853        dimensions of [u]
854    """
855    fname = 'compute_wds'
856
857#    print '    ' + fname + ': computing wind direction as ATAN2(v,u) ...'
858    theta = np.arctan2(v,u)
859    theta = np.where(theta < 0., theta + 2.*np.pi, theta)
860
861    wds = 360.*theta/(2.*np.pi)
862
863    wdsdims = dimns[:]
864    wdsvdims = dimvns[:]
865
866    return wds, wdsdims, wdsvdims
867
868def compute_wss(u, v, dimns, dimvns):
869    """ Function to compute the wind speed
870      [u]= W-E wind direction [ms-1, knot, ...]
871      [v]= N-S wind direction [ms-1, knot, ...]
872      [dimns]= list of the name of the dimensions of [u]
873      [dimvns]= list of the name of the variables with the values of the
874        dimensions of [u]
875    """
876    fname = 'compute_wss'
877
878#    print '    ' + fname + ': computing wind speed as SQRT(v**2 + u**2) ...'
879    wss = np.sqrt(u*u + v*v)
880
881    wssdims = dimns[:]
882    wssvdims = dimvns[:]
883
884    return wss, wssdims, wssvdims
885
886def timeunits_seconds(dtu):
887    """ Function to transform a time units to seconds
888    timeunits_seconds(timeuv)
889      [dtu]= time units value to transform in seconds
890    """
891    fname='timunits_seconds'
892
893    if dtu == 'years':
894        times = 365.*24.*3600.
895    elif dtu == 'weeks':
896        times = 7.*24.*3600.
897    elif dtu == 'days':
898        times = 24.*3600.
899    elif dtu == 'hours':
900        times = 3600.
901    elif dtu == 'minutes':
902        times = 60.
903    elif dtu == 'seconds':
904        times = 1.
905    elif dtu == 'miliseconds':
906        times = 1./1000.
907    else:
908        print errormsg
909        print '  ' + fname  + ": time units '" + dtu + "' not ready !!"
910        quit(-1)
911
912    return times
913
914####### ###### ##### #### ### ## #
915comboinf="\nIF -d 'variable_combo', provides information of the combination to obtain -v [varn] with the ASCII file with the combinations as -f [combofile]"
916
917parser = OptionParser()
918parser.add_option("-f", "--netCDF_file", dest="ncfile", help="file to use", metavar="FILE")
919parser.add_option("-d", "--dimensions", dest="dimns", 
920  help="[dimtn]@[dtvn],[dimzn]@[dzvn],[...,[dimxn]@[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 ('WRFtime', for WRF time copmutation)" + comboinf, 
921  metavar="LABELS")
922parser.add_option("-v", "--variables", dest="varns", 
923  help=" [varn1]|[var11]@[...[varN1]],[...,[varnM]|[var1M]@[...[varLM]]] ',' list of variables to compute [varnK] and its necessary ones [var1K]...[varPK]", metavar="VALUES")
924
925(opts, args) = parser.parse_args()
926
927#######    #######
928## MAIN
929    #######
930availdiags = ['ACRAINTOT', 'accum', 'clt', 'cllmh', 'deaccum', 'LMDZrh', 'mslp',     \
931  'OMEGAw', 'RAINTOT',                                                               \
932  'rvors', 'td', 'turbulence', 'WRFclivi', 'WRFclwvl', 'WRFgeop', 'WRFp',            \
933  'WRFrvors', 'ws', 'wds', 'wss', 'WRFheight', 'WRFheightrel', 'WRFua', 'WRFva']
934
935methods = ['accum', 'deaccum']
936
937# Variables not to check
938NONcheckingvars = ['cllmh', 'deaccum', 'TSrhs', 'TStd', 'TSwds', 'TSwss', 'WRFbils', \
939  'WRFclivi', 'WRFclwvl', 'WRFdens', 'WRFgeop',                                      \
940  'WRFp', 'WRFtd',                                                                   \
941  'WRFpos', 'WRFprc', 'WRFprls', 'WRFrh', 'LMDZrh', 'LMDZrhs', 'WRFrhs', 'WRFrvors', \
942  'WRFt', 'WRFtime', 'WRFua', 'WRFva', 'WRFwds', 'WRFwss', 'WRFheight']
943
944NONchkvardims = ['WRFtime']
945
946ofile = 'diagnostics.nc'
947
948dimns = opts.dimns
949varns = opts.varns
950
951# Special method. knowing variable combination
952##
953if opts.dimns == 'variable_combo':
954    print warnmsg
955    print '  ' + main + ': knowing variable combination !!!'
956    combination = variable_combo(opts.varns,opts.ncfile)
957    print '     COMBO: ' + combination
958    quit(-1)
959
960if not os.path.isfile(opts.ncfile):
961    print errormsg
962    print '   ' + main + ": file '" + opts.ncfile + "' does not exist !!"
963    quit(-1)
964
965ncobj = NetCDFFile(opts.ncfile, 'r')
966
967# Looking for specific variables that might be use in more than one diagnostic
968WRFgeop_compute = False
969WRFp_compute = False
970WRFt_compute = False
971WRFrh_compute = False
972WRFght_compute = False
973WRFdens_compute = False
974WRFpos_compute = False
975WRFtime_compute = False
976
977# File creation
978newnc = NetCDFFile(ofile,'w')
979
980# dimensions
981dimvalues = dimns.split(',')
982dnames = []
983dvnames = []
984
985for dimval in dimvalues:
986    dn = dimval.split('@')[0]
987    dnv = dimval.split('@')[1]
988    dnames.append(dn)
989    dvnames.append(dnv)
990    # Is there any dimension-variable which should be computed?
991    if dnv == 'WRFgeop':WRFgeop_compute = True
992    if dnv == 'WRFp': WRFp_compute = True
993    if dnv == 'WRFt': WRFt_compute = True
994    if dnv == 'WRFrh': WRFrh_compute = True
995    if dnv == 'WRFght': WRFght_compute = True
996    if dnv == 'WRFdens': WRFdens_compute = True
997    if dnv == 'WRFpos': WRFpos_compute = True
998    if dnv == 'WRFtime': WRFtime_compute = True
999
1000# diagnostics to compute
1001diags = varns.split(',')
1002Ndiags = len(diags)
1003
1004for idiag in range(Ndiags):
1005    if diags[idiag].split('|')[1].find('@') == -1:
1006        depvars = diags[idiag].split('|')[1]
1007        if depvars == 'WRFgeop':WRFgeop_compute = True
1008        if depvars == 'WRFp': WRFp_compute = True
1009        if depvars == 'WRFt': WRFt_compute = True
1010        if depvars == 'WRFrh': WRFrh_compute = True
1011        if depvars == 'WRFght': WRFght_compute = True
1012        if depvars == 'WRFdens': WRFdens_compute = True
1013        if depvars == 'WRFpos': WRFpos_compute = True
1014        if depvars == 'WRFtime': WRFtime_compute = True
1015    else:
1016        depvars = diags[idiag].split('|')[1].split('@')
1017        if gen.searchInlist(depvars, 'WRFgeop'): WRFgeop_compute = True
1018        if gen.searchInlist(depvars, 'WRFp'): WRFp_compute = True
1019        if gen.searchInlist(depvars, 'WRFt'): WRFt_compute = True
1020        if gen.searchInlist(depvars, 'WRFrh'): WRFrh_compute = True
1021        if gen.searchInlist(depvars, 'WRFght'): WRFght_compute = True
1022        if gen.searchInlist(depvars, 'WRFdens'): WRFdens_compute = True
1023        if gen.searchInlist(depvars, 'WRFpos'): WRFpos_compute = True
1024        if gen.searchInlist(depvars, 'WRFtime'): WRFtime_compute = True
1025
1026# Dictionary with the new computed variables to be able to add them
1027dictcompvars = {}
1028if WRFgeop_compute:
1029    print '  ' + main + ': Retrieving geopotential value from WRF as PH + PHB'
1030    dimv = ncobj.variables['PH'].shape
1031    WRFgeop = ncobj.variables['PH'][:] + ncobj.variables['PHB'][:]
1032
1033    # Attributes of the variable
1034    Vvals = gen.variables_values('WRFgeop')
1035    dictcompvars['WRFgeop'] = {'name': Vvals[0], 'standard_name': Vvals[1],          \
1036      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
1037
1038if WRFp_compute:
1039    print '  ' + main + ': Retrieving pressure value from WRF as P + PB'
1040    dimv = ncobj.variables['P'].shape
1041    WRFp = ncobj.variables['P'][:] + ncobj.variables['PB'][:]
1042
1043    # Attributes of the variable
1044    Vvals = gen.variables_values('WRFp')
1045    dictcompvars['WRFgeop'] = {'name': Vvals[0], 'standard_name': Vvals[1],          \
1046      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
1047
1048if WRFght_compute:
1049    print '    ' + main + ': computing geopotential height from WRF as PH + PHB ...' 
1050    WRFght = ncobj.variables['PH'][:] + ncobj.variables['PHB'][:]
1051
1052    # Attributes of the variable
1053    Vvals = gen.variables_values('WRFght')
1054    dictcompvars['WRFgeop'] = {'name': Vvals[0], 'standard_name': Vvals[1],          \
1055      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
1056
1057if WRFrh_compute:
1058    print '    ' + main + ": computing relative humidity from WRF as 'Tetens'" +     \
1059      ' equation (T,P) ...'
1060    p0=100000.
1061    p=ncobj.variables['P'][:] + ncobj.variables['PB'][:]
1062    tk = (ncobj.variables['T'][:] + 300.)*(p/p0)**(2./7.)
1063    qv = ncobj.variables['QVAPOR'][:]
1064
1065    data1 = 10.*0.6112*np.exp(17.67*(tk-273.16)/(tk-29.65))
1066    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
1067
1068    WRFrh = qv/data2
1069
1070    # Attributes of the variable
1071    Vvals = gen.variables_values('WRFrh')
1072    dictcompvars['WRFrh'] = {'name': Vvals[0], 'standard_name': Vvals[1],            \
1073      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
1074
1075if WRFt_compute:
1076    print '    ' + main + ': computing temperature from WRF as inv_potT(T + 300) ...'
1077    p0=100000.
1078    p=ncobj.variables['P'][:] + ncobj.variables['PB'][:]
1079
1080    WRFt = (ncobj.variables['T'][:] + 300.)*(p/p0)**(2./7.)
1081
1082    # Attributes of the variable
1083    Vvals = gen.variables_values('WRFt')
1084    dictcompvars['WRFt'] = {'name': Vvals[0], 'standard_name': Vvals[1],             \
1085      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
1086
1087if WRFdens_compute:
1088    print '    ' + main + ': computing air density from WRF as ((MU + MUB) * ' +     \
1089      'DNW)/g ...'
1090
1091# Just we need in in absolute values: Size of the central grid cell
1092##    dxval = ncobj.getncattr('DX')
1093##    dyval = ncobj.getncattr('DY')
1094##    mapfac = ncobj.variables['MAPFAC_M'][:]
1095##    area = dxval*dyval*mapfac
1096
1097    mu = (ncobj.variables['MU'][:] + ncobj.variables['MUB'][:])
1098    dnw = ncobj.variables['DNW'][:]
1099
1100    WRFdens = np.zeros((mu.shape[0], dnw.shape[1], mu.shape[1], mu.shape[2]),        \
1101      dtype=np.float)
1102    levval = np.zeros((mu.shape[1], mu.shape[2]), dtype=np.float)
1103
1104    for it in range(mu.shape[0]):
1105        for iz in range(dnw.shape[1]):
1106            levval.fill(np.abs(dnw[it,iz]))
1107            WRFdens[it,iz,:,:] = levval
1108            WRFdens[it,iz,:,:] = mu[it,:,:]*WRFdens[it,iz,:,:]/grav
1109
1110    # Attributes of the variable
1111    Vvals = gen.variables_values('WRFdens')
1112    dictcompvars['WRFdens'] = {'name': Vvals[0], 'standard_name': Vvals[1],          \
1113      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
1114
1115if WRFpos_compute:
1116# WRF positions from the lowest-leftest corner of the matrix
1117    print '    ' + main + ': computing position from MAPFAC_M as sqrt(DY*j**2 + ' +  \
1118      'DX*x**2)*MAPFAC_M ...'
1119
1120    mapfac = ncobj.variables['MAPFAC_M'][:]
1121
1122    distx = np.float(ncobj.getncattr('DX'))
1123    disty = np.float(ncobj.getncattr('DY'))
1124
1125    print 'distx:',distx,'disty:',disty
1126
1127    dx = mapfac.shape[2]
1128    dy = mapfac.shape[1]
1129    dt = mapfac.shape[0]
1130
1131    WRFpos = np.zeros((dt, dy, dx), dtype=np.float)
1132
1133    for i in range(1,dx):
1134        WRFpos[0,0,i] = distx*i/mapfac[0,0,i]
1135    for j in range(1,dy):
1136        i=0
1137        WRFpos[0,j,i] = WRFpos[0,j-1,i] + disty/mapfac[0,j,i]
1138        for i in range(1,dx):
1139#            WRFpos[0,j,i] = np.sqrt((disty*j)**2. + (distx*i)**2.)/mapfac[0,j,i]
1140#            WRFpos[0,j,i] = np.sqrt((disty*j)**2. + (distx*i)**2.)
1141             WRFpos[0,j,i] = WRFpos[0,j,i-1] + distx/mapfac[0,j,i]
1142
1143    for it in range(1,dt):
1144        WRFpos[it,:,:] = WRFpos[0,:,:]
1145
1146if WRFtime_compute:
1147    print '    ' + main + ': computing time from WRF as CFtime(Times) ...'
1148
1149    refdate='19491201000000'
1150    tunitsval='minutes'
1151
1152    timeobj = ncobj.variables['Times']
1153    timewrfv = timeobj[:]
1154
1155    yrref=refdate[0:4]
1156    monref=refdate[4:6]
1157    dayref=refdate[6:8]
1158    horref=refdate[8:10]
1159    minref=refdate[10:12]
1160    secref=refdate[12:14]
1161
1162    refdateS = yrref + '-' + monref + '-' + dayref + ' ' + horref + ':' + minref +   \
1163      ':' + secref
1164
1165    dt = timeobj.shape[0]
1166    WRFtime = np.zeros((dt), dtype=np.float)
1167
1168    for it in range(dt):
1169        wrfdates = gen.datetimeStr_conversion(timewrfv[it,:],'WRFdatetime', 'matYmdHMS')
1170        WRFtime[it] = gen.realdatetime1_CFcompilant(wrfdates, refdate, tunitsval)
1171
1172    tunits = tunitsval + ' since ' + refdateS
1173
1174    # Attributes of the variable
1175    dictcompvars['WRFtime'] = {'name': 'time', 'standard_name': 'time',              \
1176      'long_name': 'time', 'units': tunits, 'calendar': 'gregorian'}
1177
1178### ## #
1179# Going for the diagnostics
1180### ## #
1181print '  ' + main + ' ...'
1182varsadd = []
1183
1184for idiag in range(Ndiags):
1185    print '    diagnostic:',diags[idiag]
1186    diag = diags[idiag].split('|')[0]
1187    depvars = diags[idiag].split('|')[1].split('@')
1188    if diags[idiag].split('|')[1].find('@') != -1:
1189        depvars = diags[idiag].split('|')[1].split('@')
1190        if depvars[0] == 'deaccum': diag='deaccum'
1191        if depvars[0] == 'accum': diag='accum'
1192        for depv in depvars:
1193            if not ncobj.variables.has_key(depv) and not                             \
1194              gen.searchInlist(NONcheckingvars, depv) and                            \
1195              not gen.searchInlist(methods, depv) and not depvars[0] == 'deaccum'    \
1196              and not depvars[0] == 'accum':
1197                print errormsg
1198                print '  ' + main + ": file '" + opts.ncfile +                       \
1199                  "' does not have variable '" + depv + "' !!"
1200                quit(-1)
1201    else:
1202        depvars = diags[idiag].split('|')[1]
1203        if not ncobj.variables.has_key(depvars) and not                              \
1204          gen.searchInlist(NONcheckingvars, depvars) and                             \
1205          not gen.searchInlist(methods, depvars):
1206            print errormsg
1207            print '  ' + main + ": file '" + opts.ncfile +                           \
1208              "' does not have variable '" + depvars + "' !!"
1209            quit(-1)
1210
1211    print "\n    Computing '" + diag + "' from: ", depvars, '...'
1212
1213# acraintot: accumulated total precipitation from WRF RAINC, RAINNC
1214    if diag == 'ACRAINTOT':
1215           
1216        var0 = ncobj.variables[depvars[0]]
1217        var1 = ncobj.variables[depvars[1]]
1218        diagout = var0[:] + var1[:]
1219
1220        dnamesvar = var0.dimensions
1221        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1222
1223        # Removing the nonChecking variable-dimensions from the initial list
1224        varsadd = []
1225        for nonvd in NONchkvardims:
1226            if gen.searchInlist(dvnamesvar,nonvd): dvnamesvar.remove(nonvd)
1227            varsadd.append(nonvd)
1228
1229        ncvar.insert_variable(ncobj, 'pracc', diagout, dnamesvar, dvnamesvar, newnc)
1230
1231# accum: acumulation of any variable as (Variable, time [as [tunits]
1232#   from/since ....], newvarname)
1233    elif diag == 'accum':
1234
1235        var0 = ncobj.variables[depvars[0]]
1236        var1 = ncobj.variables[depvars[1]]
1237
1238        dnamesvar = var0.dimensions
1239        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1240
1241        diagout, diagoutd, diagoutvd = compute_accum(var0,dnamesvar,dvnamesvar)
1242
1243        CFvarn = ncvar.variables_values(depvars[0])[0]
1244
1245# Removing the flux
1246        if depvars[1] == 'XTIME':
1247            dtimeunits = var1.getncattr('description')
1248            tunits = dtimeunits.split(' ')[0]
1249        else:
1250            dtimeunits = var1.getncattr('units')
1251            tunits = dtimeunits.split(' ')[0]
1252
1253        dtime = (var1[1] - var1[0])*timeunits_seconds(tunits)
1254
1255        ncvar.insert_variable(ncobj, CFvarn + 'acc', diagout*dtime, diagoutd, diagoutvd, newnc)
1256
1257# cllmh with cldfra, pres
1258    elif diag == 'cllmh':
1259           
1260        var0 = ncobj.variables[depvars[0]]
1261        if depvars[1] == 'WRFp':
1262            var1 = WRFp
1263        else:
1264            var01 = ncobj.variables[depvars[1]]
1265            if len(size(var1.shape)) < len(size(var0.shape)):
1266                var1 = np.brodcast_arrays(var01,var0)[0]
1267            else:
1268                var1 = var01
1269
1270        diagout, diagoutd, diagoutvd = Forcompute_cllmh(var0,var1,dnames,dvnames)
1271
1272        # Removing the nonChecking variable-dimensions from the initial list
1273        varsadd = []
1274        for nonvd in NONchkvardims:
1275            if gen.searchInlist(diagoutvd,nonvd): diagoutvd.remove(nonvd)
1276            varsadd.append(nonvd)
1277
1278        ncvar.insert_variable(ncobj, 'cll', diagout[0,:], diagoutd, diagoutvd, newnc)
1279        ncvar.insert_variable(ncobj, 'clm', diagout[1,:], diagoutd, diagoutvd, newnc)
1280        ncvar.insert_variable(ncobj, 'clh', diagout[2,:], diagoutd, diagoutvd, newnc)
1281
1282# clt with cldfra
1283    elif diag == 'clt':
1284           
1285        var0 = ncobj.variables[depvars]
1286        diagout, diagoutd, diagoutvd = Forcompute_clt(var0,dnames,dvnames)
1287
1288        # Removing the nonChecking variable-dimensions from the initial list
1289        varsadd = []
1290        for nonvd in NONchkvardims:
1291            if gen.searchInlist(diagoutvd,nonvd): diagoutvd.remove(nonvd)
1292            varsadd.append(nonvd)
1293           
1294        ncvar.insert_variable(ncobj, 'clt', diagout, diagoutd, diagoutvd, newnc)
1295
1296# deaccum: deacumulation of any variable as (Variable, time [as [tunits]
1297#   from/since ....], newvarname)
1298    elif diag == 'deaccum':
1299
1300        var0 = ncobj.variables[depvars[1]]
1301        var1 = ncobj.variables[depvars[2]]
1302
1303        dnamesvar = var0.dimensions
1304        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1305
1306        diagout, diagoutd, diagoutvd = compute_deaccum(var0,dnamesvar,dvnamesvar)
1307
1308# Transforming to a flux
1309        if depvars[2] == 'XTIME':
1310            dtimeunits = var1.getncattr('description')
1311            tunits = dtimeunits.split(' ')[0]
1312        else:
1313            dtimeunits = var1.getncattr('units')
1314            tunits = dtimeunits.split(' ')[0]
1315
1316        dtime = (var1[1] - var1[0])*timeunits_seconds(tunits)
1317        ncvar.insert_variable(ncobj, depvars[3], diagout/dtime, diagoutd, diagoutvd, newnc)
1318
1319# LMDZrh (pres, t, r)
1320    elif diag == 'LMDZrh':
1321           
1322        var0 = ncobj.variables[depvars[0]][:]
1323        var1 = ncobj.variables[depvars[1]][:]
1324        var2 = ncobj.variables[depvars[2]][:]
1325
1326        diagout, diagoutd, diagoutvd = compute_rh(var0,var1,var2,dnames,dvnames)
1327        ncvar.insert_variable(ncobj, 'hur', diagout, diagoutd, diagoutvd, newnc)
1328
1329# LMDZrhs (psol, t2m, q2m)
1330    elif diag == 'LMDZrhs':
1331           
1332        var0 = ncobj.variables[depvars[0]][:]
1333        var1 = ncobj.variables[depvars[1]][:]
1334        var2 = ncobj.variables[depvars[2]][:]
1335
1336        dnamesvar = ncobj.variables[depvars[0]].dimensions
1337        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1338
1339        diagout, diagoutd, diagoutvd = compute_rh(var0,var1,var2,dnamesvar,dvnamesvar)
1340
1341        ncvar.insert_variable(ncobj, 'hurs', diagout, diagoutd, diagoutvd, newnc)
1342
1343# mslp: mean sea level pressure (pres, psfc, terrain, temp, qv)
1344    elif diag == 'mslp' or diag == 'WRFmslp':
1345           
1346        var1 = ncobj.variables[depvars[1]][:]
1347        var2 = ncobj.variables[depvars[2]][:]
1348        var4 = ncobj.variables[depvars[4]][:]
1349
1350        if diag == 'WRFmslp':
1351            var0 = WRFp
1352            var3 = WRFt
1353            dnamesvar = ncobj.variables['P'].dimensions
1354        else:
1355            var0 = ncobj.variables[depvars[0]][:]
1356            var3 = ncobj.variables[depvars[3]][:]
1357            dnamesvar = ncobj.variables[depvars[0]].dimensions
1358
1359        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1360
1361        diagout, diagoutd, diagoutvd = compute_mslp(var0, var1, var2, var3, var4,    \
1362          dnamesvar, dvnamesvar)
1363
1364        # Removing the nonChecking variable-dimensions from the initial list
1365        varsadd = []
1366        diagoutvd = list(dvnames)
1367        for nonvd in NONchkvardims:
1368            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1369            varsadd.append(nonvd)
1370        ncvar.insert_variable(ncobj, 'psl', diagout, diagoutd, diagoutvd, newnc)
1371
1372# OMEGAw (omega, p, t) from NCL formulation (https://www.ncl.ucar.edu/Document/Functions/Contributed/omega_to_w.shtml)
1373    elif diag == 'OMEGAw':
1374           
1375        var0 = ncobj.variables[depvars[0]][:]
1376        var1 = ncobj.variables[depvars[1]][:]
1377        var2 = ncobj.variables[depvars[2]][:]
1378
1379        dnamesvar = ncobj.variables[depvars[0]].dimensions
1380        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1381
1382        diagout, diagoutd, diagoutvd = compute_OMEGAw(var0,var1,var2,dnamesvar,dvnamesvar)
1383
1384        ncvar.insert_variable(ncobj, 'wa', diagout, diagoutd, diagoutvd, newnc)
1385
1386# raintot: instantaneous total precipitation from WRF as (RAINC + RAINC) / dTime
1387    elif diag == 'RAINTOT':
1388
1389        var0 = ncobj.variables[depvars[0]]
1390        var1 = ncobj.variables[depvars[1]]
1391        if depvars[2] != 'WRFtime':
1392            var2 = ncobj.variables[depvars[2]]
1393        else:
1394            var2 = np.arange(var0.shape[0], dtype=int)
1395
1396        var = var0[:] + var1[:]
1397
1398        dnamesvar = var0.dimensions
1399        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1400
1401        diagout, diagoutd, diagoutvd = compute_deaccum(var,dnamesvar,dvnamesvar)
1402
1403# Transforming to a flux
1404        if var2.shape[0] > 1:
1405            if depvars[2] != 'WRFtime':
1406                dtimeunits = var2.getncattr('units')
1407                tunits = dtimeunits.split(' ')[0]
1408   
1409                dtime = (var2[1] - var2[0])*timeunits_seconds(tunits)
1410            else:
1411                var2 = ncobj.variables['Times']
1412                time1 = var2[0,:]
1413                time2 = var2[1,:]
1414                tmf1 = ''
1415                tmf2 = ''
1416                for ic in range(len(time1)):
1417                    tmf1 = tmf1 + time1[ic]
1418                    tmf2 = tmf2 + time2[ic]
1419                dtdate1 = dtime.datetime.strptime(tmf1,"%Y-%m-%d_%H:%M:%S")
1420                dtdate2 = dtime.datetime.strptime(tmf2,"%Y-%m-%d_%H:%M:%S")
1421                diffdate12 = dtdate2 - dtdate1
1422                dtime = diffdate12.total_seconds()
1423                print 'dtime:',dtime
1424        else:
1425            print warnmsg
1426            print '  ' + main + ": only 1 time-step for '" + diag + "' !!"
1427            print '    leaving a zero value!'
1428            diagout = var0[:]*0.
1429            dtime=1.
1430
1431        # Removing the nonChecking variable-dimensions from the initial list
1432        varsadd = []
1433        for nonvd in NONchkvardims:
1434            if gen.searchInlist(diagoutvd,nonvd): diagoutvd.remove(nonvd)
1435            varsadd.append(nonvd)
1436           
1437        ncvar.insert_variable(ncobj, 'pr', diagout/dtime, diagoutd, diagoutvd, newnc)
1438
1439# rhs (psfc, t, q) from TimeSeries files
1440    elif diag == 'TSrhs':
1441           
1442        p0=100000.
1443        var0 = ncobj.variables[depvars[0]][:]
1444        var1 = (ncobj.variables[depvars[1]][:])*(var0/p0)**(2./7.)
1445        var2 = ncobj.variables[depvars[2]][:]
1446
1447        dnamesvar = ncobj.variables[depvars[0]].dimensions
1448        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1449
1450        diagout, diagoutd, diagoutvd = compute_rh(var0,var1,var2,dnamesvar,dvnamesvar)
1451
1452        ncvar.insert_variable(ncobj, 'hurs', diagout, diagoutd, diagoutvd, newnc)
1453
1454# td (psfc, t, q) from TimeSeries files
1455    elif diag == 'TStd' or diag == 'td':
1456           
1457        var0 = ncobj.variables[depvars[0]][:]
1458        var1 = ncobj.variables[depvars[1]][:] - 273.15
1459        var2 = ncobj.variables[depvars[2]][:]
1460
1461        dnamesvar = ncobj.variables[depvars[0]].dimensions
1462        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1463
1464        diagout, diagoutd, diagoutvd = compute_td(var0,var1,var2,dnamesvar,dvnamesvar)
1465
1466        ncvar.insert_variable(ncobj, 'tds', diagout, diagoutd, diagoutvd, newnc)
1467
1468# td (psfc, t, q) from TimeSeries files
1469    elif diag == 'TStdC' or diag == 'tdC':
1470           
1471        var0 = ncobj.variables[depvars[0]][:]
1472# Temperature is already in degrees Celsius
1473        var1 = ncobj.variables[depvars[1]][:]
1474        var2 = ncobj.variables[depvars[2]][:]
1475
1476        dnamesvar = ncobj.variables[depvars[0]].dimensions
1477        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1478
1479        diagout, diagoutd, diagoutvd = compute_td(var0,var1,var2,dnamesvar,dvnamesvar)
1480
1481        ncvar.insert_variable(ncobj, 'tds', diagout, diagoutd, diagoutvd, newnc)
1482
1483# wds (u, v)
1484    elif diag == 'TSwds' or diag == 'wds' :
1485 
1486        var0 = ncobj.variables[depvars[0]][:]
1487        var1 = ncobj.variables[depvars[1]][:]
1488
1489        dnamesvar = ncobj.variables[depvars[0]].dimensions
1490        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1491
1492        diagout, diagoutd, diagoutvd = compute_wds(var0,var1,dnamesvar,dvnamesvar)
1493
1494        ncvar.insert_variable(ncobj, 'wds', diagout, diagoutd, diagoutvd, newnc)
1495
1496# wss (u, v)
1497    elif diag == 'TSwss' or diag == 'wss':
1498           
1499        var0 = ncobj.variables[depvars[0]][:]
1500        var1 = ncobj.variables[depvars[1]][:]
1501
1502        dnamesvar = ncobj.variables[depvars[0]].dimensions
1503        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1504
1505        diagout, diagoutd, diagoutvd = compute_wss(var0,var1,dnamesvar,dvnamesvar)
1506
1507        ncvar.insert_variable(ncobj, 'wss', diagout, diagoutd, diagoutvd, newnc)
1508
1509# turbulence (var)
1510    elif diag == 'turbulence':
1511
1512        var0 = ncobj.variables[depvars][:]
1513
1514        dnamesvar = list(ncobj.variables[depvars].dimensions)
1515        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1516
1517        diagout, diagoutd, diagoutvd = compute_turbulence(var0,dnamesvar,dvnamesvar)
1518        valsvar = gen.variables_values(depvars)
1519
1520        newvarn = depvars + 'turb'
1521        print main + '; Lluis newvarn:', newvarn
1522        ncvar.insert_variable(ncobj, newvarn, diagout, diagoutd, 
1523          diagoutvd, newnc)
1524        print main + '; Lluis variables:', newnc.variables.keys()
1525        varobj = newnc.variables[newvarn]
1526        attrv = varobj.long_name
1527        attr = varobj.delncattr('long_name')
1528        newattr = ncvar.set_attribute(varobj, 'long_name', attrv +                   \
1529          " Taylor decomposition turbulence term")
1530
1531# WRFbils fom WRF as HFX + LH
1532    elif diag == 'WRFbils':
1533           
1534        var0 = ncobj.variables[depvars[0]][:]
1535        var1 = ncobj.variables[depvars[1]][:]
1536
1537        diagout = var0 + var1
1538        dnamesvar = list(ncobj.variables[depvars[0]].dimensions)
1539        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1540
1541        ncvar.insert_variable(ncobj, 'bils', diagout, dnamesvar, dvnamesvar, newnc)
1542
1543# WRFclivi WRF water vapour path WRFdens, QICE, QGRAUPEL, QHAIL
1544    elif diag == 'WRFclivi':
1545           
1546        var0 = WRFdens
1547        qtot = ncobj.variables[depvars[1]]
1548        qtotv = qtot[:]
1549        Nspecies = len(depvars) - 2
1550        for iv in range(Nspecies):
1551            if ncobj.variables.has_key(depvars[iv+2]):
1552                var1 = ncobj.variables[depvars[iv+2]][:]
1553                qtotv = qtotv + var1
1554
1555        dnamesvar = list(qtot.dimensions)
1556        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1557
1558        diagout, diagoutd, diagoutvd = compute_clivi(var0, qtotv, dnamesvar,dvnamesvar)
1559
1560        # Removing the nonChecking variable-dimensions from the initial list
1561        varsadd = []
1562        diagoutvd = list(dvnames)
1563        for nonvd in NONchkvardims:
1564            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1565            varsadd.append(nonvd)
1566        ncvar.insert_variable(ncobj, 'clivi', diagout, diagoutd, diagoutvd, newnc)
1567
1568# WRFclwvl WRF water cloud-condensed path WRFdens, QCLOUD, QICE, QGRAUPEL, QHAIL
1569    elif diag == 'WRFclwvl':
1570           
1571        var0 = WRFdens
1572        qtot = ncobj.variables[depvars[1]]
1573        qtotv = ncobj.variables[depvars[1]]
1574        Nspecies = len(depvars) - 2
1575        for iv in range(Nspecies):
1576            if ncobj.variables.has_key(depvars[iv+2]):
1577                var1 = ncobj.variables[depvars[iv+2]]
1578                qtotv = qtotv + var1[:]
1579
1580        dnamesvar = list(qtot.dimensions)
1581        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1582
1583        diagout, diagoutd, diagoutvd = compute_clwvl(var0, qtotv, dnamesvar,dvnamesvar)
1584
1585        # Removing the nonChecking variable-dimensions from the initial list
1586        varsadd = []
1587        diagoutvd = list(dvnames)
1588        for nonvd in NONchkvardims:
1589            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1590            varsadd.append(nonvd)
1591        ncvar.insert_variable(ncobj, 'clwvl', diagout, diagoutd, diagoutvd, newnc)
1592
1593# WRFgeop geopotential from WRF as PH + PHB
1594    elif diag == 'WRFgeop':
1595        var0 = ncobj.variables[depvars[0]][:]
1596        var1 = ncobj.variables[depvars[1]][:]
1597
1598        # de-staggering geopotential
1599        diagout0 = var0 + var1
1600        dt = diagout0.shape[0]
1601        dz = diagout0.shape[1]
1602        dy = diagout0.shape[2]
1603        dx = diagout0.shape[3]
1604
1605        diagout = np.zeros((dt,dz-1,dy,dx), dtype=np.float)
1606        diagout = 0.5*(diagout0[:,1:dz,:,:]+diagout0[:,0:dz-1,:,:])
1607
1608        # Removing the nonChecking variable-dimensions from the initial list
1609        varsadd = []
1610        diagoutvd = list(dvnames)
1611        for nonvd in NONchkvardims:
1612            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1613            varsadd.append(nonvd)
1614
1615        ncvar.insert_variable(ncobj, 'zg', diagout, dnames, diagoutvd, newnc)
1616
1617# WRFp pressure from WRF as P + PB
1618    elif diag == 'WRFp':
1619           
1620        diagout = WRFp
1621
1622        ncvar.insert_variable(ncobj, 'pres', diagout, dnames, dvnames, newnc)
1623
1624# WRFpos
1625    elif diag == 'WRFpos':
1626           
1627        dnamesvar = ncobj.variables['MAPFAC_M'].dimensions
1628        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1629
1630        ncvar.insert_variable(ncobj, 'WRFpos', WRFpos, dnamesvar, dvnamesvar, newnc)
1631
1632# WRFprw WRF water vapour path WRFdens, QVAPOR
1633    elif diag == 'WRFprw':
1634           
1635        var0 = WRFdens
1636        var1 = ncobj.variables[depvars[1]]
1637
1638        dnamesvar = list(var1.dimensions)
1639        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1640
1641        diagout, diagoutd, diagoutvd = compute_prw(var0, var1, dnamesvar,dvnamesvar)
1642
1643        # Removing the nonChecking variable-dimensions from the initial list
1644        varsadd = []
1645        diagoutvd = list(dvnames)
1646        for nonvd in NONchkvardims:
1647            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1648            varsadd.append(nonvd)
1649        ncvar.insert_variable(ncobj, 'prw', diagout, diagoutd, diagoutvd, newnc)
1650
1651# WRFrh (P, T, QVAPOR)
1652    elif diag == 'WRFrh':
1653           
1654        dnamesvar = list(ncobj.variables[depvars[2]].dimensions)
1655        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1656
1657        ncvar.insert_variable(ncobj, 'hur', WRFrh, dnames, dvnames, newnc)
1658
1659# WRFrhs (PSFC, T2, Q2)
1660    elif diag == 'WRFrhs':
1661           
1662        var0 = ncobj.variables[depvars[0]][:]
1663        var1 = ncobj.variables[depvars[1]][:]
1664        var2 = ncobj.variables[depvars[2]][:]
1665
1666        dnamesvar = list(ncobj.variables[depvars[2]].dimensions)
1667        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1668
1669        diagout, diagoutd, diagoutvd = compute_rh(var0,var1,var2,dnamesvar,dvnamesvar)
1670        ncvar.insert_variable(ncobj, 'hurs', diagout, diagoutd, diagoutvd, newnc)
1671
1672# rvors (u10, v10, WRFpos)
1673    elif diag == 'WRFrvors':
1674           
1675        var0 = ncobj.variables[depvars[0]]
1676        var1 = ncobj.variables[depvars[1]]
1677
1678        diagout = rotational_z(var0, var1, distx)
1679
1680        dnamesvar = ncobj.variables[depvars[0]].dimensions
1681        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1682
1683        ncvar.insert_variable(ncobj, 'rvors', diagout, dnamesvar, dvnamesvar, newnc)
1684
1685# WRFt (T, P, PB)
1686    elif diag == 'WRFt':
1687        var0 = ncobj.variables[depvars[0]][:]
1688        var1 = ncobj.variables[depvars[1]][:]
1689        var2 = ncobj.variables[depvars[2]][:]
1690
1691        p0=100000.
1692        p=var1 + var2
1693
1694        WRFt = (var0 + 300.)*(p/p0)**(2./7.)
1695
1696        dnamesvar = list(ncobj.variables[depvars[0]].dimensions)
1697        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1698
1699        # Removing the nonChecking variable-dimensions from the initial list
1700        varsadd = []
1701        diagoutvd = list(dvnames)
1702        for nonvd in NONchkvardims:
1703            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1704            varsadd.append(nonvd)
1705
1706        ncvar.insert_variable(ncobj, 'ta', WRFt, dnames, diagoutvd, newnc)
1707
1708# WRFua (U, V, SINALPHA, COSALPHA) to be rotated !!
1709    elif diag == 'WRFua':
1710        var0 = ncobj.variables[depvars[0]][:]
1711        var1 = ncobj.variables[depvars[1]][:]
1712        var2 = ncobj.variables[depvars[2]][:]
1713        var3 = ncobj.variables[depvars[3]][:]
1714
1715        # un-staggering variables
1716        unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1717        ua = np.zeros(tuple(unstgdims), dtype=np.float)
1718        unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1719        unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1720        unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1721        unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1722
1723        for iz in range(var0.shape[1]):
1724            ua[:,iz,:,:] = unstgvar0[:,iz,:,:]*var3 - unstgvar1[:,iz,:,:]*var2
1725
1726        dnamesvar = ['Time','bottom_top','south_north','west_east']
1727        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1728
1729        # Removing the nonChecking variable-dimensions from the initial list
1730        varsadd = []
1731        diagoutvd = list(dvnames)
1732        for nonvd in NONchkvardims:
1733            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1734            varsadd.append(nonvd)
1735
1736        ncvar.insert_variable(ncobj, 'ua', ua, dnames, diagoutvd, newnc)
1737
1738# WRFua (U, V, SINALPHA, COSALPHA) to be rotated !!
1739    elif diag == 'WRFva':
1740        var0 = ncobj.variables[depvars[0]][:]
1741        var1 = ncobj.variables[depvars[1]][:]
1742        var2 = ncobj.variables[depvars[2]][:]
1743        var3 = ncobj.variables[depvars[3]][:]
1744
1745        # un-staggering variables
1746        unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1747        va = np.zeros(tuple(unstgdims), dtype=np.float)
1748        unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1749        unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1750        unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1751        unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1752        for iz in range(var0.shape[1]):
1753            va[:,iz,:,:] = unstgvar0[:,iz,:,:]*var2 + unstgvar1[:,iz,:,:]*var3
1754
1755        dnamesvar = ['Time','bottom_top','south_north','west_east']
1756        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1757
1758        # Removing the nonChecking variable-dimensions from the initial list
1759        varsadd = []
1760        diagoutvd = list(dvnames)
1761        for nonvd in NONchkvardims:
1762            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1763            varsadd.append(nonvd)
1764        ncvar.insert_variable(ncobj, 'va', va, dnames, diagoutvd, newnc)
1765
1766# WRFtime
1767    elif diag == 'WRFtime':
1768           
1769        diagout = WRFtime
1770
1771        dnamesvar = ['Time']
1772        dvnamesvar = ['Times']
1773
1774        ncvar.insert_variable(ncobj, 'time', diagout, dnamesvar, dvnamesvar, newnc)
1775
1776# ws (U, V)
1777    elif diag == 'ws':
1778           
1779        var0 = ncobj.variables[depvars[0]][:]
1780        var1 = ncobj.variables[depvars[1]][:]
1781        # un-staggering variables
1782        unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1783        va = np.zeros(tuple(unstgdims), dtype=np.float)
1784        unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1785        unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1786        unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1787        unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1788
1789        dnamesvar = ['Time','bottom_top','south_north','west_east']
1790        diagout = np.sqrt(unstgvar0*unstgvar0 + unstgvar1*unstgvar1)
1791
1792#        dnamesvar = ncobj.variables[depvars[0]].dimensions
1793        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1794
1795        # Removing the nonChecking variable-dimensions from the initial list
1796        varsadd = []
1797        diagoutvd = list(dvnamesvar)
1798        for nonvd in NONchkvardims:
1799            if gen.searchInlist(dvnamesvar,nonvd): diagoutvd.remove(nonvd)
1800            varsadd.append(nonvd)
1801        ncvar.insert_variable(ncobj, 'ws', diagout, dnamesvar, diagoutvd, newnc)
1802
1803# wss (u10, v10)
1804    elif diag == 'wss':
1805           
1806        var0 = ncobj.variables[depvars[0]][:]
1807        var1 = ncobj.variables[depvars[1]][:]
1808
1809        diagout = np.sqrt(var0*var0 + var1*var1)
1810
1811        dnamesvar = ncobj.variables[depvars[0]].dimensions
1812        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1813
1814        ncvar.insert_variable(ncobj, 'wss', diagout, dnamesvar, dvnamesvar, newnc)
1815
1816# WRFheight height from WRF geopotential as WRFGeop/g
1817    elif diag == 'WRFheight':
1818           
1819        diagout = WRFgeop/grav
1820
1821        # Removing the nonChecking variable-dimensions from the initial list
1822        varsadd = []
1823        diagoutvd = list(dvnames)
1824        for nonvd in NONchkvardims:
1825            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1826            varsadd.append(nonvd)
1827
1828        ncvar.insert_variable(ncobj, 'zhgt', diagout, dnames, diagoutvd, newnc)
1829
1830# WRFheightrel relative-height from WRF geopotential as WRFgeop(PH + PHB)/g-HGT 'WRFheightrel|PH@PHB@HGT
1831    elif diag == 'WRFheightrel':
1832        var0 = ncobj.variables[depvars[0]][:]
1833        var1 = ncobj.variables[depvars[1]][:]
1834        var2 = ncobj.variables[depvars[2]][:]
1835
1836        dimz = var0.shape[1]
1837        diagout = np.zeros(tuple(var0.shape), dtype=np.float)
1838        for iz in range(dimz):
1839            diagout[:,iz,:,:] = (var0[:,iz,:,:]+ var1[:,iz,:,:])/grav - var2
1840
1841        # Removing the nonChecking variable-dimensions from the initial list
1842        varsadd = []
1843        diagoutvd = list(dvnames)
1844        for nonvd in NONchkvardims:
1845            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1846            varsadd.append(nonvd)
1847
1848        ncvar.insert_variable(ncobj, 'zhgtrel', diagout, dnames, diagoutvd, newnc)
1849
1850    else:
1851        print errormsg
1852        print '  ' + main + ": diagnostic '" + diag + "' not ready!!!"
1853        print '    available diagnostics: ', availdiags
1854        quit(-1)
1855
1856    newnc.sync()
1857    # Adding that additional variables required to compute some diagnostics which
1858    #   where not in the original file
1859    for vadd in varsadd:
1860        if not gen.searchInlist(newnc.variables.keys(),vadd):
1861            attrs = dictcompvars[vadd]
1862            vvn = attrs['name']
1863            if not gen.searchInlist(newnc.variables.keys(), vvn):
1864                iidvn = dvnames.index(vadd)
1865                dnn = dnames[iidvn]
1866                if vadd == 'WRFtime':
1867                    dvarvals = WRFtime[:]
1868                newvar = newnc.createVariable(vvn, 'f8', (dnn))
1869                newvar[:] = dvarvals
1870                for attn in attrs.keys():
1871                    if attn != 'name':
1872                        attv = attrs[attn]
1873                        ncvar.set_attribute(newvar, attn, attv)
1874
1875#   end of diagnostics
1876
1877# Global attributes
1878##
1879atvar = ncvar.set_attribute(newnc, 'program', 'diagnostics.py')
1880atvar = ncvar.set_attribute(newnc, 'version', '1.0')
1881atvar = ncvar.set_attribute(newnc, 'author', 'Fita Borrell, Lluis')
1882atvar = ncvar.set_attribute(newnc, 'institution', 'Laboratoire Meteorologie ' +      \
1883  'Dynamique')
1884atvar = ncvar.set_attribute(newnc, 'university', 'Universite Pierre et Marie ' +     \
1885  'Curie -- Jussieu')
1886atvar = ncvar.set_attribute(newnc, 'centre', 'Centre national de la recherche ' +    \
1887  'scientifique')
1888atvar = ncvar.set_attribute(newnc, 'city', 'Paris')
1889atvar = ncvar.set_attribute(newnc, 'original_file', opts.ncfile)
1890
1891gorigattrs = ncobj.ncattrs()
1892
1893for attr in gorigattrs:
1894    attrv = ncobj.getncattr(attr)
1895    atvar = ncvar.set_attribute(newnc, attr, attrv)
1896
1897ncobj.close()
1898newnc.close()
1899
1900print '\n' + main + ': successfull writting of diagnostics file "' + ofile + '" !!!'
Note: See TracBrowser for help on using the repository browser.