source: lmdz_wrf/trunk/tools/diag_tools.py @ 2178

Last change on this file since 2178 was 2140, checked in by lfita, 6 years ago

Adding:

  • `timemax': When a given variable [varname] got its maximum
File size: 88.7 KB
Line 
1# Tools for the compute of diagnostics
2# L. Fita, CIMA. CONICET-UBA, CNRS UMI-IFAECI, Buenos Aires, Argentina
3
4# Available general pupose diagnostics (model independent) providing (varv1, varv2, ..., dimns, dimvns)
5# compute_accum: Function to compute the accumulation of a variable
6# compute_cllmh: Function to compute cllmh: low/medium/hight cloud fraction following
7#   newmicro.F90 from LMDZ compute_clt(cldfra, pres, dimns, dimvns)
8# compute_clt: Function to compute the total cloud fraction following 'newmicro.F90' from LMDZ
9# compute_clivi: Function to compute cloud-ice water path (clivi)
10# compute_clwvl: Function to compute condensed water path (clwvl)
11# compute_deaccum: Function to compute the deaccumulation of a variable
12# compute_mslp: Function to compute mslp: mean sea level pressure following p_interp.F90 from WRF
13# compute_OMEGAw: Function to transform OMEGA [Pas-1] to velocities [ms-1]
14# compute_prw: Function to compute water vapour path (prw)
15# compute_rh: Function to compute relative humidity following 'Tetens' equation (T,P) ...'
16# compute_td: Function to compute the dew point temperature
17# compute_turbulence: Function to compute the rubulence term of the Taylor's decomposition ...'
18# C_diagnostic: Class to compute generic variables
19# compute_wd: Function to compute the wind direction 3D
20# compute_wds: Function to compute the wind direction
21# compute_wss: Function to compute the wind speed
22# compute_WRFhur: Function to compute WRF relative humidity following Teten's equation
23# compute_WRFta: Function to compute WRF air temperature
24# compute_WRFtd: Function to compute WRF dew-point air temperature
25# compute_WRFua: Function to compute geographical rotated WRF x-wind
26# compute_WRFva: Function to compute geographical rotated WRF y-wind
27# compute_WRFuava: Function to compute geographical rotated WRF 3D winds
28# compute_WRFuas: Function to compute geographical rotated WRF 2-meter x-wind
29# compute_WRFvas: Function to compute geographical rotated WRF 2-meter y-wind
30# compute_WRFuasvas: Fucntion to compute geographical rotated WRF 2-meter winds
31# derivate_centered: Function to compute the centered derivate of a given field
32# Forcompute_cllmh: Function to compute cllmh: low/medium/hight cloud fraction
33#   following newmicro.F90 from LMDZ via Fortran subroutine
34# Forcompute_clt: Function to compute the total cloud fraction following
35#   'newmicro.F90' from LMDZ via a Fortran module
36# Forcompute_fog_K84: Computation of fog and visibility following Kunkel, (1984)
37# Forcompute_fog_RUC: Computation of fog and visibility following RUC method Smirnova, (2000)
38# Forcompute_fog_FRAML50: fog and visibility following Gultepe and Milbrandt, (2010)
39# Forcompute_potevap_orPM: Function to compute potential evapotranspiration following
40#   Penman-Monteith formulation implemented in ORCHIDEE
41# Forcompute_psl_ptarget: Function to compute the sea-level pressure following
42#   target_pressure value found in `p_interp.F'
43# Forcompute_zmla_gen: Function to compute the boundary layer height following a
44#   generic method with Fortran
45# Forcompute_zwind: Function to compute the wind at a given height following the
46#   power law method
47# Forcompute_zwind_log: Function to compute the wind at a given height following the
48#   logarithmic law method
49# Forcompute_zwindMO: Function to compute the wind at a given height following the
50#   Monin-Obukhov theory
51# W_diagnostic: Class to compute WRF diagnostics variables
52
53# Others just providing variable values
54# var_cllmh: Fcuntion to compute cllmh on a 1D column
55# var_clt: Function to compute the total cloud fraction following 'newmicro.F90' from
56#   LMDZ using 1D vertical column values
57# var_convini: Function returns convective initialization (pr(t) > 0.0001) in time units
58# var_hur: Function to compute relative humidity following 'August - Roche - Magnus'
59#   formula
60# var_hur_Uhus: Function to compute relative humidity following 'August-Roche-Magnus'
61#   formula using hus
62# var_mslp: Fcuntion to compute mean sea-level pressure
63# var_td: Function to compute dew-point air temperature from temperature and pressure
64#   values
65# var_td_Uhus: Function to compute dew-point air temperature from temperature and
66#   pressure values using hus
67# var_timemax: This function returns the time at which variable reaches its maximum in time
68#   units
69# var_timeoverthres: This function returns the time at which (varv(t) > thres) in time units
70# var_virtualTemp: This function returns virtual temperature in K,
71# var_WRFtime: Function to copmute CFtimes from WRFtime variable
72# var_wd: Function to compute the wind direction
73# var_wd: Function to compute the wind speed
74# rotational_z: z-component of the rotatinoal of horizontal vectorial field
75# turbulence_var: Function to compute the Taylor's decomposition turbulence term from
76#   a a given variable
77
78import numpy as np
79from netCDF4 import Dataset as NetCDFFile
80import os
81import re
82import nc_var_tools as ncvar
83import generic_tools as gen
84import datetime as dtime
85import module_ForDiag as fdin
86import module_ForDef as fdef
87
88main = 'diag_tools.py'
89errormsg = 'ERROR -- error -- ERROR -- error'
90warnmsg = 'WARNING -- warning -- WARNING -- warning'
91
92# Constants
93grav = fdef.module_definitions.grav
94
95# Available WRFiag
96Wavailablediags = ['hur', 'p', 'ta', 'td', 'ua', 'va', 'uas', 'vas', 'wd', 'ws', 'zg']
97
98# Available General diagnostics
99Cavailablediags = ['hur', 'hur_Uhus', 'td', 'td_Uhus', 'wd', 'ws']
100
101# Gneral information
102##
103def reduce_spaces(string):
104    """ Function to give words of a line of text removing any extra space
105    """
106    values = string.replace('\n','').split(' ')
107    vals = []
108    for val in values:
109         if len(val) > 0:
110             vals.append(val)
111
112    return vals
113
114def variable_combo(varn,combofile):
115    """ Function to provide variables combination from a given variable name
116      varn= name of the variable
117      combofile= ASCII file with the combination of variables
118        [varn] [combo]
119          [combo]: '@' separated list of variables to use to generate [varn]
120            [WRFdt] to get WRF time-step (from general attributes)
121    >>> variable_combo('WRFprls','/home/lluis/PY/diagnostics.inf')
122    deaccum@RAINNC@XTIME@prnc
123    """
124    fname = 'variable_combo'
125
126    if varn == 'h':
127        print fname + '_____________________________________________________________'
128        print variable_combo.__doc__
129        quit()
130
131    if not os.path.isfile(combofile):
132        print errormsg
133        print '  ' + fname + ": file with combinations '" + combofile +              \
134          "' does not exist!!"
135        quit(-1)
136
137    objf = open(combofile, 'r')
138
139    found = False
140    for line in objf:
141        linevals = reduce_spaces(line)
142        varnf = linevals[0]
143        combo = linevals[1].replace('\n','')
144        if varn == varnf: 
145            found = True
146            break
147
148    if not found:
149        print errormsg
150        print '  ' + fname + ": variable '" + varn + "' not found in '" + combofile +\
151          "' !!"
152        combo='ERROR'
153
154    objf.close()
155
156    return combo
157
158# Mathematical operators
159##
160def compute_accum(varv, dimns, dimvns):
161    """ Function to compute the accumulation of a variable
162    compute_accum(varv, dimnames, dimvns)
163      [varv]= values to accum (assuming [t,])
164      [dimns]= list of the name of the dimensions of the [varv]
165      [dimvns]= list of the name of the variables with the values of the
166        dimensions of [varv]
167    """
168    fname = 'compute_accum'
169
170    deacdims = dimns[:]
171    deacvdims = dimvns[:]
172
173    slicei = []
174    slicee = []
175
176    Ndims = len(varv.shape)
177    for iid in range(0,Ndims):
178        slicei.append(slice(0,varv.shape[iid]))
179        slicee.append(slice(0,varv.shape[iid]))
180
181    slicee[0] = np.arange(varv.shape[0])
182    slicei[0] = np.arange(varv.shape[0])
183    slicei[0][1:varv.shape[0]] = np.arange(varv.shape[0]-1)
184
185    vari = varv[tuple(slicei)]
186    vare = varv[tuple(slicee)]
187
188    ac = vari*0.
189    for it in range(1,varv.shape[0]):
190        ac[it,] = ac[it-1,] + vare[it,]
191
192    return ac, deacdims, deacvdims
193
194def compute_deaccum(varv, dimns, dimvns):
195    """ Function to compute the deaccumulation of a variable
196    compute_deaccum(varv, dimnames, dimvns)
197      [varv]= values to deaccum (assuming [t,])
198      [dimns]= list of the name of the dimensions of the [varv]
199      [dimvns]= list of the name of the variables with the values of the
200        dimensions of [varv]
201    """
202    fname = 'compute_deaccum'
203
204    deacdims = dimns[:]
205    deacvdims = dimvns[:]
206
207    slicei = []
208    slicee = []
209
210    Ndims = len(varv.shape)
211    for iid in range(0,Ndims):
212        slicei.append(slice(0,varv.shape[iid]))
213        slicee.append(slice(0,varv.shape[iid]))
214
215    slicee[0] = np.arange(varv.shape[0])
216    slicei[0] = np.arange(varv.shape[0])
217    slicei[0][1:varv.shape[0]] = np.arange(varv.shape[0]-1)
218
219    vari = varv[tuple(slicei)]
220    vare = varv[tuple(slicee)]
221
222    deac = vare - vari
223
224    return deac, deacdims, deacvdims
225
226def derivate_centered(var,dim,dimv):
227    """ Function to compute the centered derivate of a given field
228      centered derivate(n) = (var(n-1) + var(n+1))/(2*dn).
229    [var]= variable
230    [dim]= which dimension to compute the derivate
231    [dimv]= dimension values (can be of different dimension of [var])
232    >>> derivate_centered(np.arange(16).reshape(4,4)*1.,1,1.)
233    [[  0.   1.   2.   0.]
234     [  0.   5.   6.   0.]
235     [  0.   9.  10.   0.]
236     [  0.  13.  14.   0.]]
237    """
238
239    fname = 'derivate_centered'
240   
241    vark = var.dtype
242
243    if hasattr(dimv, "__len__"):
244# Assuming that the last dimensions of var [..., N, M] are the same of dimv [N, M]
245        if len(var.shape) != len(dimv.shape):
246            dimvals = np.zeros((var.shape), dtype=vark)
247            if len(var.shape) - len(dimv.shape) == 1:
248                for iz in range(var.shape[0]):
249                    dimvals[iz,] = dimv
250            elif len(var.shape) - len(dimv.shape) == 2:
251                for it in range(var.shape[0]):
252                    for iz in range(var.shape[1]):
253                        dimvals[it,iz,] = dimv
254            else:
255                print errormsg
256                print '  ' + fname + ': dimension difference between variable',      \
257                  var.shape,'and variable with dimension values',dimv.shape,         \
258                  ' not ready !!!'
259                quit(-1)
260        else:
261            dimvals = dimv
262    else:
263# dimension values are identical everywhere!
264# from: http://stackoverflow.com/questions/16807011/python-how-to-identify-if-a-variable-is-an-array-or-a-scalar   
265        dimvals = np.ones((var.shape), dtype=vark)*dimv
266
267    derivate = np.zeros((var.shape), dtype=vark)
268    if dim > len(var.shape) - 1:
269        print errormsg
270        print '  ' + fname + ': dimension',dim,' too big for given variable of ' +   \
271          'shape:', var.shape,'!!!'
272        quit(-1)
273
274    slicebef = []
275    sliceaft = []
276    sliceder = []
277
278    for id in range(len(var.shape)):
279        if id == dim:
280            slicebef.append(slice(0,var.shape[id]-2))
281            sliceaft.append(slice(2,var.shape[id]))
282            sliceder.append(slice(1,var.shape[id]-1))
283        else:
284            slicebef.append(slice(0,var.shape[id]))
285            sliceaft.append(slice(0,var.shape[id]))
286            sliceder.append(slice(0,var.shape[id]))
287
288    if hasattr(dimv, "__len__"):
289        derivate[tuple(sliceder)] = (var[tuple(slicebef)] + var[tuple(sliceaft)])/   \
290          ((dimvals[tuple(sliceaft)] - dimvals[tuple(slicebef)]))
291        print (dimvals[tuple(sliceaft)] - dimvals[tuple(slicebef)])
292    else:
293        derivate[tuple(sliceder)] = (var[tuple(slicebef)] + var[tuple(sliceaft)])/   \
294          (2.*dimv)
295
296#    print 'before________'
297#    print var[tuple(slicebef)]
298
299#    print 'after________'
300#    print var[tuple(sliceaft)]
301
302    return derivate
303
304def rotational_z(Vx,Vy,pos):
305    """ z-component of the rotatinoal of horizontal vectorial field
306    \/ x (Vx,Vy,Vz) = \/xVy - \/yVx
307    [Vx]= Variable component x
308    [Vy]=  Variable component y
309    [pos]= poisition of the grid points
310    >>> rotational_z(np.arange(16).reshape(4,4)*1., np.arange(16).reshape(4,4)*1., 1.)
311    [[  0.   1.   2.   0.]
312     [ -4.   0.   0.  -7.]
313     [ -8.   0.   0. -11.]
314     [  0.  13.  14.   0.]]
315    """
316
317    fname =  'rotational_z'
318
319    ndims = len(Vx.shape)
320    rot1 = derivate_centered(Vy,ndims-1,pos)
321    rot2 = derivate_centered(Vx,ndims-2,pos)
322
323    rot = rot1 - rot2
324
325    return rot
326
327# Diagnostics
328##
329def Forcompute_cape_afwa(ta, hur, pa, zg, hgt, parcelm, dimns, dimvns):
330    """ Function to compute the CAPE, CIN, ZLFC, PLFC, LI following WRF
331      'phys/module_diaf_afwa.F' methodology
332    Forcompute_cape_afwa(ta, hur, pa, hgt, zsfc, parcelm, dimns, dimvns)
333      [ta]= air-temperature values (assuming [[t],z,y,x]) [K]
334      [pa]= pressure values (assuming [[t],z,y,x]) [Pa]
335      [zg]= gopotential height (assuming [[t],z,y,x]) [gpm]
336      [hgt]= topographical height (assuming [m]
337      [parcelm]= method of air-parcel to use
338      [dimns]= list of the name of the dimensions of [pa]
339      [dimvns]= list of the name of the variables with the values of the
340        dimensions of [pa]
341    """
342    fname = 'Forcompute_cape_afwa'
343
344    psldims = dimns[:]
345    pslvdims = dimvns[:]
346
347    if len(pa.shape) == 4:
348        cape = np.zeros((pa.shape[0],pa.shape[2],pa.shape[3]), dtype=np.float)
349        cin = np.zeros((pa.shape[0],pa.shape[2],pa.shape[3]), dtype=np.float)
350        zlfc = np.zeros((pa.shape[0],pa.shape[2],pa.shape[3]), dtype=np.float)
351        plfc = np.zeros((pa.shape[0],pa.shape[2],pa.shape[3]), dtype=np.float)
352        li = np.zeros((pa.shape[0],pa.shape[2],pa.shape[3]), dtype=np.float)
353
354        dx = pa.shape[3]
355        dy = pa.shape[2]
356        dz = pa.shape[1]
357        dt = pa.shape[0]
358        psldims.pop(1)
359        pslvdims.pop(1)
360
361        pcape,pcin,pzlfc,pplfc,pli= fdin.module_fordiagnostics.compute_cape_afwa4d(  \
362          ta=ta[:].transpose(), hur=hur[:].transpose(), press=pa[:].transpose(),     \
363          zg=zg[:].transpose(), hgt=hgt.transpose(), parcelmethod=parcelm, d1=dx,    \
364          d2=dy, d3=dz, d4=dt)
365        cape = pcape.transpose()
366        cin = pcin.transpose()
367        zlfc = pzlfc.transpose()
368        plfc = pplfc.transpose()
369        li = pli.transpose()
370    else:
371        print errormsg
372        print '  ' + fname + ': rank', len(pa.shape), 'not ready !!'
373        print '  it only computes 4D [t,z,y,x] rank values'
374        quit(-1)
375
376    return cape, cin, zlfc, plfc, li, psldims, pslvdims
377
378
379
380def var_clt(cfra):
381    """ Function to compute the total cloud fraction following 'newmicro.F90' from
382      LMDZ using 1D vertical column values
383      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
384    """
385    ZEPSEC=1.0E-12
386
387    fname = 'var_clt'
388
389    zclear = 1.
390    zcloud = 0.
391
392    dz = cfra.shape[0]
393    for iz in range(dz):
394        zclear =zclear*(1.-np.max([cfra[iz],zcloud]))/(1.-np.min([zcloud,1.-ZEPSEC]))
395        clt = 1. - zclear
396        zcloud = cfra[iz]
397
398    return clt
399
400def compute_clt(cldfra, dimns, dimvns):
401    """ Function to compute the total cloud fraction following 'newmicro.F90' from
402      LMDZ
403    compute_clt(cldfra, dimnames)
404      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
405      [dimns]= list of the name of the dimensions of [cldfra]
406      [dimvns]= list of the name of the variables with the values of the
407        dimensions of [cldfra]
408    """
409    fname = 'compute_clt'
410
411    cltdims = dimns[:]
412    cltvdims = dimvns[:]
413
414    if len(cldfra.shape) == 4:
415        clt = np.zeros((cldfra.shape[0],cldfra.shape[2],cldfra.shape[3]),            \
416          dtype=np.float)
417        dx = cldfra.shape[3]
418        dy = cldfra.shape[2]
419        dz = cldfra.shape[1]
420        dt = cldfra.shape[0]
421        cltdims.pop(1)
422        cltvdims.pop(1)
423
424        for it in range(dt):
425            for ix in range(dx):
426                for iy in range(dy):
427                    zclear = 1.
428                    zcloud = 0.
429                    gen.percendone(it*dx*dy + ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
430                    clt[it,iy,ix] = var_clt(cldfra[it,:,iy,ix])
431
432    else:
433        clt = np.zeros((cldfra.shape[1],cldfra.shape[2]), dtype=np.float)
434        dx = cldfra.shape[2]
435        dy = cldfra.shape[1]
436        dy = cldfra.shape[0]
437        cltdims.pop(0)
438        cltvdims.pop(0)
439        for ix in range(dx):
440            for iy in range(dy):
441                zclear = 1.
442                zcloud = 0.
443                gen.percendone(ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
444                clt[iy,ix] = var_clt(cldfra[:,iy,ix])
445
446    return clt, cltdims, cltvdims
447
448def Forcompute_clt(cldfra, dimns, dimvns):
449    """ Function to compute the total cloud fraction following 'newmicro.F90' from
450      LMDZ via a Fortran module
451    compute_clt(cldfra, dimnames)
452      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
453      [dimns]= list of the name of the dimensions of [cldfra]
454      [dimvns]= list of the name of the variables with the values of the
455        dimensions of [cldfra]
456    """
457    fname = 'Forcompute_clt'
458
459    cltdims = dimns[:]
460    cltvdims = dimvns[:]
461
462
463    if len(cldfra.shape) == 4:
464        clt = np.zeros((cldfra.shape[0],cldfra.shape[2],cldfra.shape[3]),            \
465          dtype=np.float)
466        dx = cldfra.shape[3]
467        dy = cldfra.shape[2]
468        dz = cldfra.shape[1]
469        dt = cldfra.shape[0]
470        cltdims.pop(1)
471        cltvdims.pop(1)
472
473        clt = fdin.module_fordiagnostics.compute_clt4d2(cldfra[:])
474
475    else:
476        clt = np.zeros((cldfra.shape[1],cldfra.shape[2]), dtype=np.float)
477        dx = cldfra.shape[2]
478        dy = cldfra.shape[1]
479        dy = cldfra.shape[0]
480        cltdims.pop(0)
481        cltvdims.pop(0)
482
483        clt = fdin.module_fordiagnostics.compute_clt3d1(cldfra[:])
484
485    return clt, cltdims, cltvdims
486
487def var_cllmh(cfra, p):
488    """ Fcuntion to compute cllmh on a 1D column
489    """
490
491    fname = 'var_cllmh'
492
493    ZEPSEC =1.0E-12
494    prmhc = 440.*100.
495    prmlc = 680.*100.
496
497    zclearl = 1.
498    zcloudl = 0.
499    zclearm = 1.
500    zcloudm = 0.
501    zclearh = 1.
502    zcloudh = 0.
503
504    dvz = cfra.shape[0]
505
506    cllmh = np.ones((3), dtype=np.float)
507
508    for iz in range(dvz):
509        if p[iz] < prmhc:
510            cllmh[2] = cllmh[2]*(1.-np.max([cfra[iz], zcloudh]))/(1.-                \
511              np.min([zcloudh,1.-ZEPSEC]))
512            zcloudh = cfra[iz]
513        elif p[iz] >= prmhc and p[iz] < prmlc:
514            cllmh[1] = cllmh[1]*(1.-np.max([cfra[iz], zcloudm]))/(1.-                \
515              np.min([zcloudm,1.-ZEPSEC]))
516            zcloudm = cfra[iz]
517        elif p[iz] >= prmlc:
518            cllmh[0] = cllmh[0]*(1.-np.max([cfra[iz], zcloudl]))/(1.-                \
519              np.min([zcloudl,1.-ZEPSEC]))
520            zcloudl = cfra[iz]
521
522    cllmh = 1.- cllmh
523
524    return cllmh
525
526def Forcompute_cllmh(cldfra, pres, dimns, dimvns):
527    """ Function to compute cllmh: low/medium/hight cloud fraction following newmicro.F90 from LMDZ via Fortran subroutine
528    compute_clt(cldfra, pres, dimns, dimvns)
529      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
530      [pres] = pressure field
531      [dimns]= list of the name of the dimensions of [cldfra]
532      [dimvns]= list of the name of the variables with the values of the
533        dimensions of [cldfra]
534    """
535    fname = 'Forcompute_cllmh'
536
537    cllmhdims = dimns[:]
538    cllmhvdims = dimvns[:]
539
540    if len(cldfra.shape) == 4:
541        dx = cldfra.shape[3]
542        dy = cldfra.shape[2]
543        dz = cldfra.shape[1]
544        dt = cldfra.shape[0]
545        cllmhdims.pop(1)
546        cllmhvdims.pop(1)
547
548        cllmh = fdin.module_fordiagnostics.compute_cllmh4d2(cldfra[:], pres[:])
549       
550    else:
551        dx = cldfra.shape[2]
552        dy = cldfra.shape[1]
553        dz = cldfra.shape[0]
554        cllmhdims.pop(0)
555        cllmhvdims.pop(0)
556
557        cllmh = fdin.module_fordiagnostics.compute_cllmh3d1(cldfra[:], pres[:])
558
559    return cllmh, cllmhdims, cllmhvdims
560
561def compute_cllmh(cldfra, pres, dimns, dimvns):
562    """ Function to compute cllmh: low/medium/hight cloud fraction following newmicro.F90 from LMDZ
563    compute_clt(cldfra, pres, dimns, dimvns)
564      [cldfra]= cloud fraction values (assuming [[t],z,y,x])
565      [pres] = pressure field
566      [dimns]= list of the name of the dimensions of [cldfra]
567      [dimvns]= list of the name of the variables with the values of the
568        dimensions of [cldfra]
569    """
570    fname = 'compute_cllmh'
571
572    cllmhdims = dimns[:]
573    cllmhvdims = dimvns[:]
574
575    if len(cldfra.shape) == 4:
576        dx = cldfra.shape[3]
577        dy = cldfra.shape[2]
578        dz = cldfra.shape[1]
579        dt = cldfra.shape[0]
580        cllmhdims.pop(1)
581        cllmhvdims.pop(1)
582
583        cllmh = np.ones(tuple([3, dt, dy, dx]), dtype=np.float)
584
585        for it in range(dt):
586            for ix in range(dx):
587                for iy in range(dy):
588                    gen.percendone(it*dx*dy + ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
589                    cllmh[:,it,iy,ix] = var_cllmh(cldfra[it,:,iy,ix], pres[it,:,iy,ix])
590       
591    else:
592        dx = cldfra.shape[2]
593        dy = cldfra.shape[1]
594        dz = cldfra.shape[0]
595        cllmhdims.pop(0)
596        cllmhvdims.pop(0)
597
598        cllmh = np.ones(tuple([3, dy, dx]), dtype=np.float)
599
600        for ix in range(dx):
601            for iy in range(dy):
602                gen.percendone(ix*dy + iy,dx*dy, 5, 'diagnosted')
603                cllmh[:,iy,ix] = var_cllmh(cldfra[:,iy,ix], pres[:,iy,ix])
604
605    return cllmh, cllmhdims, cllmhvdims
606
607def compute_clivi(dens, qtot, dimns, dimvns):
608    """ Function to compute cloud-ice water path (clivi)
609      [dens] = density [in kgkg-1] (assuming [t],z,y,x)
610      [qtot] = added mixing ratio of all cloud-ice species in [kgkg-1] (assuming [t],z,y,x)
611      [dimns]= list of the name of the dimensions of [q]
612      [dimvns]= list of the name of the variables with the values of the
613        dimensions of [q]
614    """
615    fname = 'compute_clivi'
616
617    clividims = dimns[:]
618    clivivdims = dimvns[:]
619
620    if len(qtot.shape) == 4:
621        clividims.pop(1)
622        clivivdims.pop(1)
623    else:
624        clividims.pop(0)
625        clivivdims.pop(0)
626
627    data1 = dens*qtot
628    clivi = np.sum(data1, axis=1)
629
630    return clivi, clividims, clivivdims
631
632
633def compute_clwvl(dens, qtot, dimns, dimvns):
634    """ Function to compute condensed water path (clwvl)
635      [dens] = density [in kgkg-1] (assuming [t],z,y,x)
636      [qtot] = added mixing ratio of all cloud-water species in [kgkg-1] (assuming [t],z,y,x)
637      [dimns]= list of the name of the dimensions of [q]
638      [dimvns]= list of the name of the variables with the values of the
639        dimensions of [q]
640    """
641    fname = 'compute_clwvl'
642
643    clwvldims = dimns[:]
644    clwvlvdims = dimvns[:]
645
646    if len(qtot.shape) == 4:
647        clwvldims.pop(1)
648        clwvlvdims.pop(1)
649    else:
650        clwvldims.pop(0)
651        clwvlvdims.pop(0)
652
653    data1 = dens*qtot
654    clwvl = np.sum(data1, axis=1)
655
656    return clwvl, clwvldims, clwvlvdims
657
658def var_virtualTemp (temp,rmix):
659    """ This function returns virtual temperature in K,
660      temp: temperature [K]
661      rmix: mixing ratio in [kgkg-1]
662    """
663
664    fname = 'var_virtualTemp'
665
666    virtual=temp*(0.622+rmix)/(0.622*(1.+rmix))
667
668    return virtual
669
670def var_convini(pr, time, dimns, dimvns):
671    """ This function returns convective initialization (pr(t) > 0.0001) in time units
672      pr: precipitation fux [kgm-2s-1]
673      time: time in CF coordinates
674    """
675    fname = 'var_convini'
676
677    dt = pr.shape[0]
678    dy = pr.shape[1]
679    dx = pr.shape[2]
680
681    vardims = dimns[:]
682    varvdims = dimvns[:]
683
684    vardims.pop(0)
685    varvdims.pop(0)
686
687    prmin = 0.0001
688    convini = np.ones((dy, dx), dtype=np.float)*gen.fillValueF
689    for it in range(dt):
690# NOT working ?
691#        convini = np.where(convini == gen.fillValueF and pr[it,:,:] >= prmin,        \
692#          time[it], fillValueF)
693        for j in range(dy):
694            for i in range(dx):
695                if convini[j,i] == gen.fillValueF and pr[it,j,i] >= prmin:
696                    convini[j,i] = time[it]
697                    break
698
699    return convini, vardims, varvdims
700
701def var_timemax(varv, time, dimns, dimvns):
702    """ This function returns the time at which variable reaches its maximum in time
703        units
704      varv: values of the variable to use
705      time: time in CF coordinates
706    """
707    fname = 'var_timemax'
708
709    dt = varv.shape[0]
710    dy = varv.shape[1]
711    dx = varv.shape[2]
712
713    vardims = dimns[:]
714    varvdims = dimvns[:]
715
716    vardims.pop(0)
717    varvdims.pop(0)
718
719    timemax = np.ones((dy, dx), dtype=np.float)*gen.fillValueF
720    varmax = np.max(varv, axis=0)
721    for j in range(dy):
722        for i in range(dx):
723            it = gen.index_vec(varv[:,j,i], varmax[j,i])
724            timemax[j,i] = time[it]
725
726    return timemax, vardims, varvdims
727
728def var_timeoverthres(varv, time, thres, dimns, dimvns):
729    """ This function returns the time at which (varv(t) > thres) in time units
730      varv: values of the variable to use
731      time: time in CF coordinates
732      thres: threshold to overpass
733    """
734    fname = 'var_timeoverthres'
735
736    dt = varv.shape[0]
737    dy = varv.shape[1]
738    dx = varv.shape[2]
739
740    vardims = dimns[:]
741    varvdims = dimvns[:]
742
743    vardims.pop(0)
744    varvdims.pop(0)
745
746    timeoverthres = np.ones((dy, dx), dtype=np.float)*gen.fillValueF
747    for it in range(dt):
748        for j in range(dy):
749            for i in range(dx):
750                if timeoverthres[j,i] == gen.fillValueF and varv[it,j,i] >= thres:
751                    timeoverthres[j,i] = time[it]
752                    break
753
754    return timeoverthres, vardims, varvdims
755
756def Forcompute_zint(var, zinterlev, zweights, dimns, dimvns):
757    """ Function to compute a vertical integration of volumetric quantities
758    Forcompute_mrso(smois, dsoil, dimns, dimvns)
759      [var]= values (assuming [[t],z,y,x]) [volumetric units]
760      [zinterlev]= depth of each layer (assuming [z]) [same z units as var]
761      [zweights]= weights to apply to each level (just in case...)
762      [dimns]= list of the name of the dimensions of [smois]
763      [dimvns]= list of the name of the variables with the values of the
764        dimensions of [smois]
765    """
766    fname = 'Forcompute_zint'
767
768    vardims = dimns[:]
769    varvdims = dimvns[:]
770
771    if len(var.shape) == 4:
772        zint = np.zeros((var.shape[0],var.shape[2],var.shape[3]), dtype=np.float)
773        dx = var.shape[3]
774        dy = var.shape[2]
775        dz = var.shape[1]
776        dt = var.shape[0]
777        vardims.pop(1)
778        varvdims.pop(1)
779
780        zintvart=fdin.module_fordiagnostics.compute_zint4d(var4d=var[:].transpose(), \
781          dlev=zinterlev[:].transpose(), zweight=zweights[:].transpose(), d1=dx,     \
782          d2=dy, d3=dz, d4=dt)
783        zintvar = zintvart.transpose()
784    else:
785        print errormsg
786        print '  ' + fname + ': rank', len(var.shape), 'not ready !!'
787        print '  it only computes 4D [t,z,y,x] rank values'
788        quit(-1)
789
790    return zintvar, vardims, varvdims
791
792def var_mslp(pres, psfc, ter, tk, qv):
793    """ Function to compute mslp on a 1D column
794    """
795
796    fname = 'var_mslp'
797
798    N = 1.0
799    expon=287.04*.0065/9.81
800    pref = 40000.
801
802# First find where about 400 hPa is located
803    dz=len(pres) 
804
805    kref = -1
806    pinc = pres[0] - pres[dz-1]
807
808    if pinc < 0.:
809        for iz in range(1,dz):
810            if pres[iz-1] >= pref and pres[iz] < pref: 
811                kref = iz
812                break
813    else:
814        for iz in range(dz-1):
815            if pres[iz] >= pref and pres[iz+1] < pref: 
816                kref = iz
817                break
818
819    if kref == -1:
820        print errormsg
821        print '  ' + fname + ': no reference pressure:',pref,'found!!'
822        print '    values:',pres[:]
823        quit(-1)
824
825    mslp = 0.
826
827# We are below both the ground and the lowest data level.
828
829# First, find the model level that is closest to a "target" pressure
830# level, where the "target" pressure is delta-p less that the local
831# value of a horizontally smoothed surface pressure field.  We use
832# delta-p = 150 hPa here. A standard lapse rate temperature profile
833# passing through the temperature at this model level will be used
834# to define the temperature profile below ground.  This is similar
835# to the Benjamin and Miller (1990) method, using 
836# 700 hPa everywhere for the "target" pressure.
837
838# ptarget = psfc - 15000.
839    ptarget = 70000.
840    dpmin=1.e4
841    kupper = 0
842    if pinc > 0.:
843        for iz in range(dz-1,0,-1):
844            kupper = iz
845            dp=np.abs( pres[iz] - ptarget )
846            if dp < dpmin: exit
847            dpmin = np.min([dpmin, dp])
848    else:
849        for iz in range(dz):
850            kupper = iz
851            dp=np.abs( pres[iz] - ptarget )
852            if dp < dpmin: exit
853            dpmin = np.min([dpmin, dp])
854
855    pbot=np.max([pres[0], psfc])
856#    zbot=0.
857
858#    tbotextrap=tk(i,j,kupper,itt)*(pbot/pres_field(i,j,kupper,itt))**expon
859#    tvbotextrap=virtual(tbotextrap,qv(i,j,1,itt))
860
861#    data_out(i,j,itt,1) = (zbot+tvbotextrap/.0065*(1.-(interp_levels(1)/pbot)**expon))
862    tbotextrap = tk[kupper]*(psfc/ptarget)**expon
863    tvbotextrap = var_virtualTemp(tbotextrap, qv[kupper])
864    mslp = psfc*( (tvbotextrap+0.0065*ter)/tvbotextrap)**(1./expon)
865
866    return mslp
867
868def compute_mslp(pressure, psurface, terrain, temperature, qvapor, dimns, dimvns):
869    """ Function to compute mslp: mean sea level pressure following p_interp.F90 from WRF
870    var_mslp(pres, ter, tk, qv, dimns, dimvns)
871      [pressure]= pressure field [Pa] (assuming [[t],z,y,x])
872      [psurface]= surface pressure field [Pa]
873      [terrain]= topography [m]
874      [temperature]= temperature [K]
875      [qvapor]= water vapour mixing ratio [kgkg-1]
876      [dimns]= list of the name of the dimensions of [cldfra]
877      [dimvns]= list of the name of the variables with the values of the
878        dimensions of [pres]
879    """
880
881    fname = 'compute_mslp'
882
883    mslpdims = list(dimns[:])
884    mslpvdims = list(dimvns[:])
885
886    if len(pressure.shape) == 4:
887        mslpdims.pop(1)
888        mslpvdims.pop(1)
889    else:
890        mslpdims.pop(0)
891        mslpvdims.pop(0)
892
893    if len(pressure.shape) == 4:
894        dx = pressure.shape[3]
895        dy = pressure.shape[2]
896        dz = pressure.shape[1]
897        dt = pressure.shape[0]
898
899        mslpv = np.zeros(tuple([dt, dy, dx]), dtype=np.float)
900
901# Terrain... to 2D !
902        terval = np.zeros(tuple([dy, dx]), dtype=np.float)
903        if len(terrain.shape) == 3:
904            terval = terrain[0,:,:]
905        else:
906            terval = terrain
907
908        for ix in range(dx):
909            for iy in range(dy):
910                if terval[iy,ix] > 0.:
911                    for it in range(dt):
912                        mslpv[it,iy,ix] = var_mslp(pressure[it,:,iy,ix],             \
913                          psurface[it,iy,ix], terval[iy,ix], temperature[it,:,iy,ix],\
914                          qvapor[it,:,iy,ix])
915
916                        gen.percendone(it*dx*dy + ix*dy + iy, dx*dy*dt, 5, 'diagnosted')
917                else:
918                    mslpv[:,iy,ix] = psurface[:,iy,ix]
919
920    else:
921        dx = pressure.shape[2]
922        dy = pressure.shape[1]
923        dz = pressure.shape[0]
924
925        mslpv = np.zeros(tuple([dy, dx]), dtype=np.float)
926
927# Terrain... to 2D !
928        terval = np.zeros(tuple([dy, dx]), dtype=np.float)
929        if len(terrain.shape) == 3:
930            terval = terrain[0,:,:]
931        else:
932            terval = terrain
933
934        for ix in range(dx):
935            for iy in range(dy):
936                gen.percendone(ix*dy + iy,dx*dy, 5, 'diagnosted')
937                if terval[iy,ix] > 0.:
938                    mslpv[iy,ix] = var_mslp(pressure[:,iy,ix], psurface[iy,ix],      \
939                      terval[iy,ix], temperature[:,iy,ix], qvapor[:,iy,ix])
940                else:
941                    mslpv[iy,ix] = psfc[iy,ix]
942
943    return mslpv, mslpdims, mslpvdims
944
945def Forcompute_psl_ecmwf(ps, hgt, ta1, pa2, unpa1, dimns, dimvns):
946    """ Function to compute the sea-level pressure following Mats Hamrud and Philippe Courtier [Pa]
947    Forcompute_psl_ptarget(ps, hgt, ta1, pa2, unpa1, dimns, dimvns)
948      [ps]= surface pressure values (assuming [[t],y,x]) [Pa]
949      [hgt]= opography (assuming [y,x]) [m]
950      [ta1]= air-temperature values at first half-mass level (assuming [[t],y,x]) [K]
951      [pa2]= pressure values at second full-mass levels (assuming [[t],y,x]) [Pa]
952      [unpa1]= pressure values at first half-mass levels (assuming [[t],y,x]) [Pa]
953      [dimns]= list of the name of the dimensions of [pa]
954      [dimvns]= list of the name of the variables with the values of the
955        dimensions of [pa]
956    """
957    fname = 'Forcompute_psl_ecmwf'
958
959    vardims = dimns[:]
960    varvdims = dimvns[:]
961
962    if len(pa2.shape) == 3:
963        psl = np.zeros((pa2.shape[0],pa2.shape[1],pa2.shape[2]), dtype=np.float)
964        dx = pa2.shape[2]
965        dy = pa2.shape[1]
966        dt = pa2.shape[0]
967        pslt= fdin.module_fordiagnostics.compute_psl_ecmwf( ps=ps[:].transpose(),    \
968          hgt=hgt[:].transpose(), t=ta1[:].transpose(), press=pa2[:].transpose(),    \
969          unpress=unpa1[:].transpose(), d1=dx, d2=dy, d4=dt)
970        psl = pslt.transpose()
971    else:
972        print errormsg
973        print '  ' + fname + ': rank', len(pa2.shape), 'not ready !!'
974        print '  it only computes 3D [t,y,x] rank values'
975        quit(-1)
976
977    return psl, vardims, varvdims
978
979def Forcompute_psl_ptarget(pa, ps, ta, hgt, qv, target_pressure, dimns, dimvns):
980    """ Function to compute the sea-level pressure following target_pressure value
981      found in `p_interp.F'
982    Forcompute_psl_ptarget(pa, ps, ta, hgt, qv, dimns, dimvns)
983      [pa]= pressure values (assuming [[t],z,y,x]) [Pa]
984      [ps]= surface pressure values (assuming [[t],y,x]) [Pa]
985      [ta]= air-temperature values (assuming [[t],z,y,x]) [K]
986      [hgt]= opography (assuming [y,x]) [m]
987      [qv]= water vapour mixing ratio (assuming [[t],z,y,x]) [kgkg-1]
988      [dimns]= list of the name of the dimensions of [pa]
989      [dimvns]= list of the name of the variables with the values of the
990        dimensions of [pa]
991    """
992    fname = 'Forcompute_psl_ptarget'
993
994    psldims = dimns[:]
995    pslvdims = dimvns[:]
996
997    if len(pa.shape) == 4:
998        psl = np.zeros((pa.shape[0],pa.shape[2],pa.shape[3]), dtype=np.float)
999        dx = pa.shape[3]
1000        dy = pa.shape[2]
1001        dz = pa.shape[1]
1002        dt = pa.shape[0]
1003        psldims.pop(1)
1004        pslvdims.pop(1)
1005
1006        pslt= fdin.module_fordiagnostics.compute_psl_ptarget4d2(                     \
1007          press=pa[:].transpose(), ps=ps[:].transpose(), hgt=hgt[:].transpose(),     \
1008          ta=ta[:].transpose(), qv=qv[:].transpose(), ptarget=target_pressure,       \
1009          d1=dx, d2=dy, d3=dz, d4=dt)
1010        psl = pslt.transpose()
1011    else:
1012        print errormsg
1013        print '  ' + fname + ': rank', len(pa.shape), 'not ready !!'
1014        print '  it only computes 4D [t,z,y,x] rank values'
1015        quit(-1)
1016
1017    return psl, psldims, pslvdims
1018
1019def Forcompute_zmla_gen(theta, qratio, zpl, hgt, dimns, dimvns):
1020    """ Function to compute the boundary layer height following a generic method
1021      with Fortran
1022    Forcompute_zmla_gen(theta, qratio, zpl, hgt, zmla, dimns, dimvns)
1023      [theta]= potential air-temperature values (assuming [[t],z,y,x]) [K]
1024      [qratio]= water mixing ratio (assuming [[t],z,y,x]) [kgkg-1]
1025      [zpl]= height from sea level (assuming [[t],z,y,x]) [m]
1026      [hgt]= topographical height (assuming [m]
1027      [zmla]= boundary layer height [m]
1028      [dimns]= list of the name of the dimensions of [theta]
1029      [dimvns]= list of the name of the variables with the values of the
1030        dimensions of [theta]
1031    """
1032    fname = 'Forcompute_zmla_gen'
1033
1034    zmladims = dimns[:]
1035    zmlavdims = dimvns[:]
1036
1037    if len(theta.shape) == 4:
1038        zmla= np.zeros((theta.shape[0],theta.shape[2],theta.shape[3]), dtype=np.float)
1039
1040        dx = theta.shape[3]
1041        dy = theta.shape[2]
1042        dz = theta.shape[1]
1043        dt = theta.shape[0]
1044        zmladims.pop(1)
1045        zmlavdims.pop(1)
1046
1047        pzmla= fdin.module_fordiagnostics.compute_zmla_generic4d(                    \
1048          tpot=theta[:].transpose(), qratio=qratio[:].transpose(),                   \
1049          z=zpl[:].transpose(), hgt=hgt.transpose(), d1=dx, d2=dy, d3=dz, d4=dt)
1050        zmla = pzmla.transpose()
1051    else:
1052        print errormsg
1053        print '  ' + fname + ': rank', len(theta.shape), 'not ready !!'
1054        print '  it only computes 4D [t,z,y,x] rank values'
1055        quit(-1)
1056
1057    return zmla, zmladims, zmlavdims
1058
1059def Forcompute_zwind(ua, va, z, uas, vas, sina, cosa, zval, dimns, dimvns):
1060    """ Function to compute the wind at a given height following the power law method
1061    Forcompute_zwind(ua, va, zsl, uas, vas, hgt, sina, cosa, zval, dimns, dimvns)
1062      [ua]= x-component of unstaggered 3D wind (assuming [[t],z,y,x]) [ms-1]
1063      [va]= y-component of unstaggered 3D wind (assuming [[t],z,y,x]) [ms-1]
1064      [z]= height above surface [m]
1065      [uas]= x-component of unstaggered 10 m wind (assuming [[t],z,y,x]) [ms-1]
1066      [vas]= y-component of unstaggered 10 m wind (assuming [[t],z,y,x]) [ms-1]
1067      [sina]= local sine of map rotation [1.]
1068      [cosa]= local cosine of map rotation [1.]
1069      [zval]= desired height for winds [m]
1070      [dimns]= list of the name of the dimensions of [ua]
1071      [dimvns]= list of the name of the variables with the values of the
1072        dimensions of [ua]
1073    """
1074    fname = 'Forcompute_zwind'
1075
1076    vardims = dimns[:]
1077    varvdims = dimvns[:]
1078
1079    if len(ua.shape) == 4:
1080        var1= np.zeros((ua.shape[0],ua.shape[2],ua.shape[3]), dtype=np.float)
1081        var2= np.zeros((ua.shape[0],ua.shape[2],ua.shape[3]), dtype=np.float)
1082
1083        dx = ua.shape[3]
1084        dy = ua.shape[2]
1085        dz = ua.shape[1]
1086        dt = ua.shape[0]
1087        vardims.pop(1)
1088        varvdims.pop(1)
1089
1090        pvar1, pvar2= fdin.module_fordiagnostics.compute_zwind4d(ua=ua.transpose(),  \
1091          va=va[:].transpose(), z=z[:].transpose(), uas=uas.transpose(),             \
1092          vas=vas.transpose(), sina=sina.transpose(), cosa=cosa.transpose(),         \
1093          zextrap=zval, d1=dx, d2=dy, d3=dz, d4=dt)
1094        var1 = pvar1.transpose()
1095        var2 = pvar2.transpose()
1096    else:
1097        print errormsg
1098        print '  ' + fname + ': rank', len(ua.shape), 'not ready !!'
1099        print '  it only computes 4D [t,z,y,x] rank values'
1100        quit(-1)
1101
1102    return var1, var2, vardims, varvdims
1103
1104def Forcompute_zwind_log(ua, va, z, uas, vas, sina, cosa, zval, dimns, dimvns):
1105    """ Function to compute the wind at a given height following the logarithmic law method
1106    Forcompute_zwind(ua, va, zsl, uas, vas, hgt, sina, cosa, zval, dimns, dimvns)
1107      [ua]= x-component of unstaggered 3D wind (assuming [[t],z,y,x]) [ms-1]
1108      [va]= y-component of unstaggered 3D wind (assuming [[t],z,y,x]) [ms-1]
1109      [z]= height above surface [m]
1110      [uas]= x-component of unstaggered 10 m wind (assuming [[t],z,y,x]) [ms-1]
1111      [vas]= y-component of unstaggered 10 m wind (assuming [[t],z,y,x]) [ms-1]
1112      [sina]= local sine of map rotation [1.]
1113      [cosa]= local cosine of map rotation [1.]
1114      [zval]= desired height for winds [m]
1115      [dimns]= list of the name of the dimensions of [ua]
1116      [dimvns]= list of the name of the variables with the values of the
1117        dimensions of [ua]
1118    """
1119    fname = 'Forcompute_zwind_log'
1120
1121    vardims = dimns[:]
1122    varvdims = dimvns[:]
1123
1124    if len(ua.shape) == 4:
1125        var1= np.zeros((ua.shape[0],ua.shape[2],ua.shape[3]), dtype=np.float)
1126        var2= np.zeros((ua.shape[0],ua.shape[2],ua.shape[3]), dtype=np.float)
1127
1128        dx = ua.shape[3]
1129        dy = ua.shape[2]
1130        dz = ua.shape[1]
1131        dt = ua.shape[0]
1132        vardims.pop(1)
1133        varvdims.pop(1)
1134
1135        pvar1, pvar2= fdin.module_fordiagnostics.compute_zwind_log4d(                \
1136          ua=ua.transpose(), va=va[:].transpose(), z=z[:].transpose(),               \
1137          uas=uas.transpose(), vas=vas.transpose(), sina=sina.transpose(),           \
1138          cosa=cosa.transpose(), zextrap=zval, d1=dx, d2=dy, d3=dz, d4=dt)
1139        var1 = pvar1.transpose()
1140        var2 = pvar2.transpose()
1141    else:
1142        print errormsg
1143        print '  ' + fname + ': rank', len(ua.shape), 'not ready !!'
1144        print '  it only computes 4D [t,z,y,x] rank values'
1145        quit(-1)
1146
1147    return var1, var2, vardims, varvdims
1148
1149def Forcompute_zwindMO(ust, znt, rmol, uas, vas, sina, cosa, zval, dimns, dimvns):
1150    """ Function to compute the wind at a given height following the Monin-Obukhov theory
1151    Forcompute_zwind(ust, znt, rmol, uas, vas, sina, cosa, zval, dimns, dimvns)
1152      [ust]: u* in similarity theory (assuming [[t],y,x]) [ms-1]
1153      [znt]: thermal time-varying roughness length (assuming [[t],y,x]) [m]
1154      [rmol]: inverse of Obukhov length (assuming [[t],y,x]) [m-1]
1155      [uas]= x-component of unstaggered 10 m wind (assuming [[t],y,x]) [ms-1]
1156      [vas]= y-component of unstaggered 10 m wind (assuming [[t],y,x]) [ms-1]
1157      [sina]= local sine of map rotation [1.]
1158      [cosa]= local cosine of map rotation [1.]
1159      [zval]= desired height for winds [m]
1160      [dimns]= list of the name of the dimensions of [uas]
1161      [dimvns]= list of the name of the variables with the values of the
1162        dimensions of [uas]
1163    """
1164    fname = 'Forcompute_zwindMO'
1165
1166    vardims = dimns[:]
1167    varvdims = dimvns[:]
1168
1169    if len(uas.shape) == 3:
1170        var1= np.zeros((uas.shape[0],uas.shape[1],uas.shape[2]), dtype=np.float)
1171        var2= np.zeros((uas.shape[0],uas.shape[1],uas.shape[2]), dtype=np.float)
1172
1173        dx = uas.shape[2]
1174        dy = uas.shape[1]
1175        dt = uas.shape[0]
1176
1177        pvar1, pvar2 = fdin.module_fordiagnostics.compute_zwindmo3d(                 \
1178          ust=ust.transpose(), znt=znt[:].transpose(), rmol=rmol[:].transpose(),     \
1179          uas=uas.transpose(), vas=vas.transpose(), sina=sina.transpose(),           \
1180          cosa=cosa.transpose(), newz=zval, d1=dx, d2=dy, d3=dt)
1181        var1 = pvar1.transpose()
1182        var2 = pvar2.transpose()
1183    else:
1184        print errormsg
1185        print '  ' + fname + ': rank', len(uas.shape), 'not ready !!'
1186        print '  it only computes 3D [t,y,x] rank values'
1187        quit(-1)
1188
1189    return var1, var2, vardims, varvdims
1190
1191def Forcompute_potevap_orPM(rho1, ust, uas, vas, tas, ps, qv1, dimns, dimvns):
1192    """ Function to compute potential evapotranspiration following Penman-Monteith
1193      formulation implemented in ORCHIDEE in src_sechiba/enerbil.f90
1194    Forcompute_potevap_orPM(rho1, uas, vas, tas, ps, qv2, qv1, dimns, dimvns)
1195      [rho1]= air-density at the first layer (assuming [[t],y,m]) [kgm-3]
1196      [ust]= u* in similarity theory (assuming [[t],y,x]) [ms-1]
1197      [uas]= x-component of unstaggered 10 m wind (assuming [[t],y,x]) [ms-1]
1198      [vas]= y-component of unstaggered 10 m wind (assuming [[t],y,x]) [ms-1]
1199      [tas]= 2m air temperature [K]
1200      [ps]= surface pressure [Pa]
1201      [qv1]= mixing ratio at the first atmospheric layer [kgkg-1]
1202      [dimns]= list of the name of the dimensions of [uas]
1203      [dimvns]= list of the name of the variables with the values of the
1204        dimensions of [uas]
1205    """
1206    fname = 'Forcompute_potevap_orPM'
1207
1208    vardims = dimns[:]
1209    varvdims = dimvns[:]
1210
1211    if len(uas.shape) == 3:
1212        var1= np.zeros((uas.shape[0],uas.shape[1],uas.shape[2]), dtype=np.float)
1213        var2= np.zeros((uas.shape[0],uas.shape[1],uas.shape[2]), dtype=np.float)
1214
1215        dx = uas.shape[2]
1216        dy = uas.shape[1]
1217        dt = uas.shape[0]
1218
1219        pvar = fdin.module_fordiagnostics.compute_potevap_orpm3d(                    \
1220          rho1=rho1.transpose(), ust=ust.transpose(), uas=uas.transpose(),           \
1221          vas=vas.transpose(), tas=tas.transpose(), ps=ps.transpose(),               \
1222          qv1=qv1.transpose(), d1=dx, d2=dy, d3=dt)
1223        var = pvar.transpose()
1224    else:
1225        print errormsg
1226        print '  ' + fname + ': rank', len(uas.shape), 'not ready !!'
1227        print '  it only computes 3D [t,y,x] rank values'
1228        quit(-1)
1229
1230    return var, vardims, varvdims
1231
1232def Forcompute_fog_K84(qcloud, qice, dimns, dimvns):
1233    """ Function to compute fog and visibility following Kunkel, (1984)
1234    Forcompute_fog_K84(qcloud, qice, dimns, dimvns)
1235      [qcloud]= cloud mixing ratio [kgk-1]
1236      [qice]= ice mixing ratio [kgk-1]
1237      [dimns]= list of the name of the dimensions of [uas]
1238      [dimvns]= list of the name of the variables with the values of the
1239        dimensions of [qcloud]
1240    """
1241    fname = 'Forcompute_fog_K84'
1242
1243    vardims = dimns[:]
1244    varvdims = dimvns[:]
1245
1246    if len(qcloud.shape) == 4:
1247        var= np.zeros((qcloud.shape[0],qcloud.shape[2],qcloud.shape[3]), dtype=np.float)
1248
1249        dx = qcloud.shape[3]
1250        dy = qcloud.shape[2]
1251        dz = qcloud.shape[1]
1252        dt = qcloud.shape[0]
1253        vardims.pop(1)
1254        varvdims.pop(1)
1255
1256        pvar1, pvar2 = fdin.module_fordiagnostics.compute_fog_k84(                   \
1257          qc=qcloud[:,0,:,:].transpose(), qi=qice[:,0,:,:].transpose(), d1=dx, d2=dy,\
1258          d3=dt)
1259        var1 = pvar1.transpose()
1260        var2 = pvar2.transpose()
1261    else:
1262        print errormsg
1263        print '  ' + fname + ': rank', len(qcloud.shape), 'not ready !!'
1264        print '  it only computes 4D [t,z,y,x] rank values'
1265        quit(-1)
1266
1267    return var1, var2, vardims, varvdims
1268
1269def Forcompute_fog_RUC(qvapor, temp, pres, dimns, dimvns):
1270    """ Function to compute fog and visibility following RUC method Smirnova, (2000)
1271    Forcompute_fog_RUC(qcloud, qice, dimns, dimvns)
1272      [qvapor]= water vapor mixing ratio [kgk-1]
1273      [temp]= temperature [K]
1274      [pres]= pressure [Pa]
1275      [dimns]= list of the name of the dimensions of [uas]
1276      [dimvns]= list of the name of the variables with the values of the
1277        dimensions of [qcloud]
1278    """
1279    fname = 'Forcompute_fog_RUC'
1280
1281    vardims = dimns[:]
1282    varvdims = dimvns[:]
1283
1284    if len(qvapor.shape) == 4:
1285        var= np.zeros((qvapor.shape[0],qvapor.shape[2],qvapor.shape[3]), dtype=np.float)
1286
1287        dx = qvapor.shape[3]
1288        dy = qvapor.shape[2]
1289        dz = qvapor.shape[1]
1290        dt = qvapor.shape[0]
1291        vardims.pop(1)
1292        varvdims.pop(1)
1293
1294        pvar1, pvar2 = fdin.module_fordiagnostics.compute_fog_ruc(                   \
1295          qv=qvapor[:,0,:,:].transpose(), ta=temp[:,0,:,:].transpose(),              \
1296          pres=pres[:,0,:,:].transpose(), d1=dx, d2=dy, d3=dt)
1297        var1 = pvar1.transpose()
1298        var2 = pvar2.transpose()
1299    elif len(qvapor.shape) == 3:
1300        var= np.zeros((qvapor.shape[0],qvapor.shape[1],qvapor.shape[2]), dtype=np.float)
1301
1302        dx = qvapor.shape[2]
1303        dy = qvapor.shape[1]
1304        dt = qvapor.shape[0]
1305
1306        pvar1, pvar2 = fdin.module_fordiagnostics.compute_fog_ruc(                   \
1307          qv=qvapor[:].transpose(), ta=temp[:].transpose(), pres=pres[:].transpose(),\
1308          d1=dx, d2=dy, d3=dt)
1309        var1 = pvar1.transpose()
1310        var2 = pvar2.transpose()
1311    else:
1312        print errormsg
1313        print '  ' + fname + ': rank', len(qcloud.shape), 'not ready !!'
1314        print '  it only computes 4D [t,z,y,x] or 3D [t,z,y,x] rank values'
1315        quit(-1)
1316
1317    return var1, var2, vardims, varvdims
1318
1319def Forcompute_fog_FRAML50(qvapor, temp, pres, dimns, dimvns):
1320    """ Function to compute fog (vis < 1km) and visibility following FRAM-L 50 % prob
1321         Gultepe, and Milbrandt, (2010), J. Appl. Meteor. Climatol.
1322    Forcompute_fog_FRAML50(qvapor, temp, pres, dimns, dimvns)
1323      [qvapor]= vapor mixing ratio [kgk-1]
1324      [temp]= temperature [K]
1325      [pres]= pressure [Pa]
1326      [dimns]= list of the name of the dimensions of [uas]
1327      [dimvns]= list of the name of the variables with the values of the
1328        dimensions of [qvapor]
1329    """
1330    fname = 'Forcompute_fog_FRAML50'
1331
1332    vardims = dimns[:]
1333    varvdims = dimvns[:]
1334
1335    if len(qvapor.shape) == 4:
1336        var= np.zeros((qvapor.shape[0],qvapor.shape[2],qvapor.shape[3]), dtype=np.float)
1337
1338        dx = qvapor.shape[3]
1339        dy = qvapor.shape[2]
1340        dz = qvapor.shape[1]
1341        dt = qvapor.shape[0]
1342        vardims.pop(1)
1343        varvdims.pop(1)
1344
1345        pvar1, pvar2 = fdin.module_fordiagnostics.compute_fog_framl50(               \
1346          qv=qvapor[:,0,:,:].transpose(), ta=temp[:,0,:,:].transpose(),              \
1347          pres=pres[:,0,:,:].transpose(), d1=dx, d2=dy, d3=dt)
1348        var1 = pvar1.transpose()
1349        var2 = pvar2.transpose()
1350    elif len(qvapor.shape) == 3:
1351        var= np.zeros((qvapor.shape[0],qvapor.shape[1],qvapor.shape[2]), dtype=np.float)
1352
1353        dx = qvapor.shape[2]
1354        dy = qvapor.shape[1]
1355        dt = qvapor.shape[0]
1356
1357        pvar1, pvar2 = fdin.module_fordiagnostics.compute_fog_framl50(               \
1358          qv=qvapor[:].transpose(), ta=temp[:].transpose(), pres=pres[:].transpose(),\
1359          d1=dx, d2=dy, d3=dt)
1360        var1 = pvar1.transpose()
1361        var2 = pvar2.transpose()
1362    else:
1363        print errormsg
1364        print '  ' + fname + ': rank', len(qvapor.shape), 'not ready !!'
1365        print '  it only computes 4D [t,z,y,x] or 3D [t,y,x] rank values'
1366        quit(-1)
1367
1368    return var1, var2, vardims, varvdims
1369
1370####### ###### ##### #### ### ## # END Fortran diagnostics
1371
1372def compute_OMEGAw(omega, p, t, dimns, dimvns):
1373    """ Function to transform OMEGA [Pas-1] to velocities [ms-1]
1374      tacking: https://www.ncl.ucar.edu/Document/Functions/Contributed/omega_to_w.shtml
1375      [omega] = vertical velocity [in ms-1] (assuming [t],z,y,x)
1376      [p] = pressure in [Pa] (assuming [t],z,y,x)
1377      [t] = temperature in [K] (assuming [t],z,y,x)
1378      [dimns]= list of the name of the dimensions of [q]
1379      [dimvns]= list of the name of the variables with the values of the
1380        dimensions of [q]
1381    """
1382    fname = 'compute_OMEGAw'
1383
1384    rgas = 287.058            # J/(kg-K) => m2/(s2 K)
1385    g    = 9.80665            # m/s2
1386
1387    wdims = dimns[:]
1388    wvdims = dimvns[:]
1389
1390    rho  = p/(rgas*t)         # density => kg/m3
1391    w    = -omega/(rho*g)     
1392
1393    return w, wdims, wvdims
1394
1395def compute_prw(dens, q, dimns, dimvns):
1396    """ Function to compute water vapour path (prw)
1397      [dens] = density [in kgkg-1] (assuming [t],z,y,x)
1398      [q] = mixing ratio in [kgkg-1] (assuming [t],z,y,x)
1399      [dimns]= list of the name of the dimensions of [q]
1400      [dimvns]= list of the name of the variables with the values of the
1401        dimensions of [q]
1402    """
1403    fname = 'compute_prw'
1404
1405    prwdims = dimns[:]
1406    prwvdims = dimvns[:]
1407
1408    if len(q.shape) == 4:
1409        prwdims.pop(1)
1410        prwvdims.pop(1)
1411    else:
1412        prwdims.pop(0)
1413        prwvdims.pop(0)
1414
1415    data1 = dens*q
1416    prw = np.sum(data1, axis=1)
1417
1418    return prw, prwdims, prwvdims
1419
1420def compute_rh(p, t, q, dimns, dimvns):
1421    """ Function to compute relative humidity following 'Tetens' equation (T,P) ...'
1422      [t]= temperature (assuming [[t],z,y,x] in [K])
1423      [p] = pressure field (assuming in [hPa])
1424      [q] = mixing ratio in [kgkg-1]
1425      [dimns]= list of the name of the dimensions of [t]
1426      [dimvns]= list of the name of the variables with the values of the
1427        dimensions of [t]
1428    """
1429    fname = 'compute_rh'
1430
1431    rhdims = dimns[:]
1432    rhvdims = dimvns[:]
1433
1434    data1 = 10.*0.6112*np.exp(17.67*(t-273.16)/(t-29.65))
1435    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
1436
1437    rh = q/data2
1438
1439    return rh, rhdims, rhvdims
1440
1441def compute_td(p, temp, qv, dimns, dimvns):
1442    """ Function to compute the dew point temperature
1443      [p]= pressure [Pa]
1444      [temp]= temperature [C]
1445      [qv]= mixing ratio [kgkg-1]
1446      [dimns]= list of the name of the dimensions of [p]
1447      [dimvns]= list of the name of the variables with the values of the
1448        dimensions of [p]
1449    """
1450    fname = 'compute_td'
1451
1452#    print '    ' + fname + ': computing dew-point temperature from TS as t and Tetens...'
1453# tacking from: http://en.wikipedia.org/wiki/Dew_point
1454    tk = temp
1455    data1 = 10.*0.6112*np.exp(17.67*(tk-273.16)/(tk-29.65))
1456    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
1457
1458    rh = qv/data2
1459               
1460    pa = rh * data1
1461    td = 257.44*np.log(pa/6.1121)/(18.678-np.log(pa/6.1121))
1462
1463    tddims = dimns[:]
1464    tdvdims = dimvns[:]
1465
1466    return td, tddims, tdvdims
1467
1468def var_WRFtime(timewrfv, refdate='19491201000000', tunitsval='minutes'):
1469    """ Function to copmute CFtimes from WRFtime variable
1470      refdate= [YYYYMMDDMIHHSS] format of reference date
1471      tunitsval= CF time units
1472      timewrfv= matrix string values of WRF 'Times' variable
1473    """
1474    fname = 'var_WRFtime'
1475
1476    yrref=refdate[0:4]
1477    monref=refdate[4:6]
1478    dayref=refdate[6:8]
1479    horref=refdate[8:10]
1480    minref=refdate[10:12]
1481    secref=refdate[12:14]
1482
1483    refdateS = yrref + '-' + monref + '-' + dayref + ' ' + horref + ':' + minref +   \
1484      ':' + secref
1485
1486    dt = timewrfv.shape[0]
1487    WRFtime = np.zeros((dt), dtype=np.float)
1488
1489    for it in range(dt):
1490        wrfdates = gen.datetimeStr_conversion(timewrfv[it,:],'WRFdatetime', 'matYmdHMS')
1491        WRFtime[it] = gen.realdatetime1_CFcompilant(wrfdates, refdate, tunitsval)
1492
1493    tunits = tunitsval + ' since ' + refdateS
1494
1495    return WRFtime, tunits
1496
1497def turbulence_var(varv, dimvn, dimn):
1498    """ Function to compute the Taylor's decomposition turbulence term from a a given variable
1499      x*=<x^2>_t-(<X>_t)^2
1500    turbulence_var(varv,dimn)
1501      varv= values of the variable
1502      dimvn= names of the dimension of the variable
1503      dimn= names of the dimensions (as a dictionary with 'X', 'Y', 'Z', 'T')
1504    >>> turbulence_var(np.arange((27)).reshape(3,3,3),['time','y','x'],{'T':'time', 'Y':'y', 'X':'x'})
1505    [[ 54.  54.  54.]
1506     [ 54.  54.  54.]
1507     [ 54.  54.  54.]]
1508    """
1509    fname = 'turbulence_varv'
1510
1511    timedimid = dimvn.index(dimn['T'])
1512
1513    varv2 = varv*varv
1514
1515    vartmean = np.mean(varv, axis=timedimid)
1516    var2tmean = np.mean(varv2, axis=timedimid)
1517
1518    varvturb = var2tmean - (vartmean*vartmean)
1519
1520    return varvturb
1521
1522def compute_turbulence(v, dimns, dimvns):
1523    """ Function to compute the rubulence term of the Taylor's decomposition ...'
1524      x*=<x^2>_t-(<X>_t)^2
1525      [v]= variable (assuming [[t],z,y,x])
1526      [dimns]= list of the name of the dimensions of [v]
1527      [dimvns]= list of the name of the variables with the values of the
1528        dimensions of [v]
1529    """
1530    fname = 'compute_turbulence'
1531
1532    turbdims = dimns[:]
1533    turbvdims = dimvns[:]
1534
1535    turbdims.pop(0)
1536    turbvdims.pop(0)
1537
1538    v2 = v*v
1539
1540    vartmean = np.mean(v, axis=0)
1541    var2tmean = np.mean(v2, axis=0)
1542
1543    turb = var2tmean - (vartmean*vartmean)
1544
1545    return turb, turbdims, turbvdims
1546
1547def compute_wd(u, v, dimns, dimvns):
1548    """ Function to compute the wind direction
1549      [u]= W-E wind direction [ms-1, knot, ...]
1550      [v]= N-S wind direction [ms-1, knot, ...]
1551      [dimns]= list of the name of the dimensions of [u]
1552      [dimvns]= list of the name of the variables with the values of the
1553        dimensions of [u]
1554    """
1555    fname = 'compute_wds'
1556
1557#    print '    ' + fname + ': computing wind direction as ATAN2(v,u) ...'
1558    theta = np.arctan2(v,u)
1559    theta = np.where(theta < 0., theta + 2.*np.pi, theta)
1560
1561    var = 360.*theta/(2.*np.pi)
1562
1563    vardims = dimns[:]
1564    varvdims = dimvns[:]
1565
1566    return var, vardims, varvdims
1567
1568def compute_wds(u, v, dimns, dimvns):
1569    """ Function to compute the wind direction
1570      [u]= W-E wind direction [ms-1, knot, ...]
1571      [v]= N-S wind direction [ms-1, knot, ...]
1572      [dimns]= list of the name of the dimensions of [u]
1573      [dimvns]= list of the name of the variables with the values of the
1574        dimensions of [u]
1575    """
1576    fname = 'compute_wds'
1577
1578#    print '    ' + fname + ': computing wind direction as ATAN2(v,u) ...'
1579    theta = np.arctan2(v,u)
1580    theta = np.where(theta < 0., theta + 2.*np.pi, theta)
1581
1582    wds = 360.*theta/(2.*np.pi)
1583
1584    wdsdims = dimns[:]
1585    wdsvdims = dimvns[:]
1586
1587    return wds, wdsdims, wdsvdims
1588
1589def compute_wss(u, v, dimns, dimvns):
1590    """ Function to compute the wind speed
1591      [u]= W-E wind direction [ms-1, knot, ...]
1592      [v]= N-S wind direction [ms-1, knot, ...]
1593      [dimns]= list of the name of the dimensions of [u]
1594      [dimvns]= list of the name of the variables with the values of the
1595        dimensions of [u]
1596    """
1597    fname = 'compute_wss'
1598
1599#    print '    ' + fname + ': computing wind speed as SQRT(v**2 + u**2) ...'
1600    wss = np.sqrt(u*u + v*v)
1601
1602    wssdims = dimns[:]
1603    wssvdims = dimvns[:]
1604
1605    return wss, wssdims, wssvdims
1606
1607def timeunits_seconds(dtu):
1608    """ Function to transform a time units to seconds
1609    timeunits_seconds(timeuv)
1610      [dtu]= time units value to transform in seconds
1611    """
1612    fname='timunits_seconds'
1613
1614    if dtu == 'years':
1615        times = 365.*24.*3600.
1616    elif dtu == 'weeks':
1617        times = 7.*24.*3600.
1618    elif dtu == 'days':
1619        times = 24.*3600.
1620    elif dtu == 'hours':
1621        times = 3600.
1622    elif dtu == 'minutes':
1623        times = 60.
1624    elif dtu == 'seconds':
1625        times = 1.
1626    elif dtu == 'miliseconds':
1627        times = 1./1000.
1628    else:
1629        print errormsg
1630        print '  ' + fname  + ": time units '" + dtu + "' not ready !!"
1631        quit(-1)
1632
1633    return times
1634
1635def compute_WRFhur(t, p, qv, dimns, dimvns):
1636    """ Function to compute WRF relative humidity following Teten's equation
1637      t= orginal WRF temperature
1638      p= original WRF pressure (P + PB)
1639      formula:
1640          temp = theta*(p/p0)**(R/Cp)
1641
1642    """
1643    fname = 'compute_WRFtd'
1644
1645    tk = (t+300.)*(p/fdef.module_definitions.p0ref)**(fdef.module_definitions.rcp)
1646    data1 = 10.*0.6112*np.exp(17.67*(tk-273.16)/(tk-29.65))
1647    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
1648
1649    rh = qv/data2
1650               
1651    dnamesvar = ['Time','bottom_top','south_north','west_east']
1652    dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dimns,dimvns)
1653
1654    return rh, dnamesvar, dvnamesvar
1655
1656def compute_WRFua(u, v, sina, cosa, dimns, dimvns):
1657    """ Function to compute geographical rotated WRF 3D winds
1658      u= orginal WRF x-wind
1659      v= orginal WRF y-wind
1660      sina= original WRF local sinus of map rotation
1661      cosa= original WRF local cosinus of map rotation
1662      formula:
1663        ua = u*cosa-va*sina
1664        va = u*sina+va*cosa
1665    """
1666    fname = 'compute_WRFua'
1667
1668    var0 = u
1669    var1 = v
1670    var2 = sina
1671    var3 = cosa
1672
1673    # un-staggering variables
1674    unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1675    ua = np.zeros(tuple(unstgdims), dtype=np.float)
1676    unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1677    unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1678    unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1679    unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1680
1681    for iz in range(var0.shape[1]):
1682        ua[:,iz,:,:] = unstgvar0[:,iz,:,:]*var3 - unstgvar1[:,iz,:,:]*var2
1683
1684    dnamesvar = ['Time','bottom_top','south_north','west_east']
1685    dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dimns,dimvns)
1686
1687    return ua, dnamesvar, dvnamesvar
1688
1689def compute_WRFva(u, v, sina, cosa, dimns, dimvns):
1690    """ Function to compute geographical rotated WRF 3D winds
1691      u= orginal WRF x-wind
1692      v= orginal WRF y-wind
1693      sina= original WRF local sinus of map rotation
1694      cosa= original WRF local cosinus of map rotation
1695      formula:
1696        ua = u*cosa-va*sina
1697        va = u*sina+va*cosa
1698    """
1699    fname = 'compute_WRFva'
1700
1701    var0 = u
1702    var1 = v
1703    var2 = sina
1704    var3 = cosa
1705
1706    # un-staggering variables
1707    unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1708    va = np.zeros(tuple(unstgdims), dtype=np.float)
1709    unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1710    unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1711    unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1712    unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1713
1714    for iz in range(var0.shape[1]):
1715        va[:,iz,:,:] = unstgvar0[:,iz,:,:]*var2 + unstgvar1[:,iz,:,:]*var3
1716
1717    dnamesvar = ['Time','bottom_top','south_north','west_east']
1718    dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dimns,dimvns)
1719
1720    return va, dnamesvar, dvnamesvar
1721
1722def compute_WRFuava(u, v, sina, cosa, dimns, dimvns):
1723    """ Function to compute geographical rotated WRF 3D winds
1724      u= orginal WRF x-wind
1725      v= orginal WRF y-wind
1726      sina= original WRF local sinus of map rotation
1727      cosa= original WRF local cosinus of map rotation
1728      formula:
1729        ua = u*cosa-va*sina
1730        va = u*sina+va*cosa
1731    """
1732    fname = 'compute_WRFuava'
1733
1734    var0 = u
1735    var1 = v
1736    var2 = sina
1737    var3 = cosa
1738
1739    # un-staggering variables
1740    unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1741    ua = np.zeros(tuple(unstgdims), dtype=np.float)
1742    va = np.zeros(tuple(unstgdims), dtype=np.float)
1743    unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1744    unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1745    unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1746    unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1747
1748    for iz in range(var0.shape[1]):
1749        ua[:,iz,:,:] = unstgvar0[:,iz,:,:]*var3 - unstgvar1[:,iz,:,:]*var2
1750        va[:,iz,:,:] = unstgvar0[:,iz,:,:]*var2 + unstgvar1[:,iz,:,:]*var3
1751
1752    dnamesvar = ['Time','bottom_top','south_north','west_east']
1753    dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dimns,dimvns)
1754
1755    return ua, va, dnamesvar, dvnamesvar
1756
1757def compute_WRFuas(u10, v10, sina, cosa, dimns, dimvns):
1758    """ Function to compute geographical rotated WRF 2-meter x-wind
1759      u10= orginal WRF 10m x-wind
1760      v10= orginal WRF 10m y-wind
1761      sina= original WRF local sinus of map rotation
1762      cosa= original WRF local cosinus of map rotation
1763      formula:
1764        uas = u10*cosa-va10*sina
1765        vas = u10*sina+va10*cosa
1766    """
1767    fname = 'compute_WRFuas'
1768
1769    var0 = u10
1770    var1 = v10
1771    var2 = sina
1772    var3 = cosa
1773
1774    uas = np.zeros(var0.shape, dtype=np.float)
1775    vas = np.zeros(var0.shape, dtype=np.float)
1776
1777    uas = var0*var3 - var1*var2
1778
1779    dnamesvar = ['Time','south_north','west_east']
1780    dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dimns,dimvns)
1781
1782    return uas, dnamesvar, dvnamesvar
1783
1784def compute_WRFvas(u10, v10, sina, cosa, dimns, dimvns):
1785    """ Function to compute geographical rotated WRF 2-meter y-wind
1786      u10= orginal WRF 10m x-wind
1787      v10= orginal WRF 10m y-wind
1788      sina= original WRF local sinus of map rotation
1789      cosa= original WRF local cosinus of map rotation
1790      formula:
1791        uas = u10*cosa-va10*sina
1792        vas = u10*sina+va10*cosa
1793    """
1794    fname = 'compute_WRFvas'
1795
1796    var0 = u10
1797    var1 = v10
1798    var2 = sina
1799    var3 = cosa
1800
1801    uas = np.zeros(var0.shape, dtype=np.float)
1802    vas = np.zeros(var0.shape, dtype=np.float)
1803
1804    vas = var0*var2 + var1*var3
1805
1806    dnamesvar = ['Time','south_north','west_east']
1807    dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dimns,dimvns)
1808
1809    return vas, dnamesvar, dvnamesvar
1810
1811def compute_WRFuasvas(u10, v10, sina, cosa, dimns, dimvns):
1812    """ Function to compute geographical rotated WRF 2-meter winds
1813      u10= orginal WRF 10m x-wind
1814      v10= orginal WRF 10m y-wind
1815      sina= original WRF local sinus of map rotation
1816      cosa= original WRF local cosinus of map rotation
1817      formula:
1818        uas = u10*cosa-va10*sina
1819        vas = u10*sina+va10*cosa
1820    """
1821    fname = 'compute_WRFuasvas'
1822
1823    var0 = u10
1824    var1 = v10
1825    var2 = sina
1826    var3 = cosa
1827
1828    uas = np.zeros(var0.shape, dtype=np.float)
1829    vas = np.zeros(var0.shape, dtype=np.float)
1830
1831    uas = var0*var3 - var1*var2
1832    vas = var0*var2 + var1*var3
1833
1834    dnamesvar = ['Time','south_north','west_east']
1835    dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dimns,dimvns)
1836
1837    return uas, vas, dnamesvar, dvnamesvar
1838
1839def compute_WRFta(t, p, dimns, dimvns):
1840    """ Function to compute WRF air temperature
1841      t= orginal WRF temperature
1842      p= original WRF pressure (P + PB)
1843      formula:
1844          temp = theta*(p/p0)**(R/Cp)
1845
1846    """
1847    fname = 'compute_WRFta'
1848
1849    ta = (t+300.)*(p/fdef.module_definitions.p0ref)**(fdef.module_definitions.rcp)
1850
1851    dnamesvar = ['Time','bottom_top','south_north','west_east']
1852    dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dimns,dimvns)
1853
1854    return ta, dnamesvar, dvnamesvar
1855
1856def compute_WRFtd(t, p, qv, dimns, dimvns):
1857    """ Function to compute WRF dew-point air temperature
1858      t= orginal WRF temperature
1859      p= original WRF pressure (P + PB)
1860      formula:
1861          temp = theta*(p/p0)**(R/Cp)
1862
1863    """
1864    fname = 'compute_WRFtd'
1865
1866    tk = (t+300.)*(p/fdef.module_definitions.p0ref)**(fdef.module_definitions.rcp)
1867    data1 = 10.*0.6112*np.exp(17.67*(tk-273.16)/(tk-29.65))
1868    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
1869
1870    rh = qv/data2
1871               
1872    pa = rh * data1
1873    td = 257.44*np.log(pa/6.1121)/(18.678-np.log(pa/6.1121))
1874
1875    dnamesvar = ['Time','bottom_top','south_north','west_east']
1876    dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dimns,dimvns)
1877
1878    return td, dnamesvar, dvnamesvar
1879
1880def compute_WRFwd(u, v, sina, cosa, dimns, dimvns):
1881    """ Function to compute the wind direction
1882      u= W-E wind direction [ms-1]
1883      v= N-S wind direction [ms-1]
1884      sina= original WRF local sinus of map rotation
1885      cosa= original WRF local cosinus of map rotation
1886    """
1887    fname = 'compute_WRFwd'
1888    var0 = u
1889    var1 = v
1890    var2 = sina
1891    var3 = cosa
1892
1893    # un-staggering variables
1894    unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1895    ua = np.zeros(tuple(unstgdims), dtype=np.float)
1896    va = np.zeros(tuple(unstgdims), dtype=np.float)
1897    unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1898    unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1899    unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1900    unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1901
1902    for iz in range(var0.shape[1]):
1903        ua[:,iz,:,:] = unstgvar0[:,iz,:,:]*var3 - unstgvar1[:,iz,:,:]*var2
1904        va[:,iz,:,:] = unstgvar0[:,iz,:,:]*var2 + unstgvar1[:,iz,:,:]*var3
1905
1906    theta = np.arctan2(va,ua)
1907    theta = np.where(theta < 0., theta + 2.*np.pi, theta)
1908
1909    wd = 360.*theta/(2.*np.pi)
1910
1911    dnamesvar = ['Time','bottom_top','south_north','west_east']
1912    dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dimns,dimvns)
1913
1914    return wd
1915
1916####### Variables (as they arrive without dimensions staff)
1917
1918def var_hur(p, t, q):
1919    """ Function to compute relative humidity following 'August - Roche - Magnus' formula
1920      [t]= temperature (assuming [[t],z,y,x] in [K])
1921      [p] = pressure field (assuming in [Pa])
1922      [q] = mixing ratio in [kgkg-1]
1923      ref.: M. G. Lawrence, BAMS, 2005, 225
1924    >>> print var_rh(101300., 290., 3.)
1925    0.250250256174
1926    """
1927    fname = 'var_hur'
1928
1929    # Enthalpy of vaporization [Jkg-1]
1930    L = 2.453*10.**6
1931    # Gas constant for water vapor [JK-1Kg-1]
1932    Rv = 461.5 
1933
1934    # Computing saturated pressure
1935    #es = 10.*0.6112*np.exp(17.67*(t-273.16)/(t-29.65))
1936    #es = 6.11*np.exp(L/Rv*(1.-273./t)/273.)
1937
1938    # August - Roche - Magnus formula
1939    es = 6.1094*np.exp(17.625*(t-273.15)/((t-273.15)+243.04))
1940
1941    # Saturated mixing ratio [g/kg]
1942    ws = 0.622*es/(0.01*p-es)
1943
1944    hur = q/(ws*1000.)
1945
1946    #print 'q:', q[5,5,5,5], 't:', t[5,5,5,5], 'p:', p[5,5,5,5]
1947    #print 'es:', es[5,5,5,5], 'ws:', ws[5,5,5,5], 'hur:', hur[5,5,5,5]
1948
1949    return hur
1950
1951def var_hur_Uhus(p, t, hus):
1952    """ Function to compute relative humidity following 'August - Roche - Magnus' formula using hus
1953      [t]= temperature (assuming [[t],z,y,x] in [K])
1954      [p] = pressure field (assuming in [Pa])
1955      [hus] = specific humidty [1]
1956      ref.: M. G. Lawrence, BAMS, 2005, 225
1957    >>> print var_rh(101300., 290., 3.)
1958    0.250250256174
1959    """
1960    fname = 'var_hur_Uhus'
1961
1962    # Enthalpy of vaporization [Jkg-1]
1963    L = 2.453*10.**6
1964    # Gas constant for water vapor [JK-1Kg-1]
1965    Rv = 461.5 
1966
1967    # Computing saturated pressure
1968    #es = 10.*0.6112*np.exp(17.67*(t-273.16)/(t-29.65))
1969    #es = 6.11*np.exp(L/Rv*(1.-273./t)/273.)
1970
1971    # August - Roche - Magnus formula
1972    es = 6.1094*np.exp(17.625*(t-273.15)/((t-273.15)+243.04))
1973
1974    # Saturated mixing ratio [g/kg]
1975    ws = 0.622*es/(0.01*p-es)
1976
1977    # Mixing ratio
1978    q = hus/(1.+hus)
1979
1980    hur = q/ws
1981
1982    #print 'q:', q[5,5,5,5], 't:', t[5,5,5,5], 'p:', p[5,5,5,5]
1983    #print 'es:', es[5,5,5,5], 'ws:', ws[5,5,5,5], 'hur:', hur[5,5,5,5]
1984
1985    return hur
1986
1987def var_td(t, p, qv):
1988    """ Function to compute dew-point air temperature from temperature and pressure values
1989      t= temperature [K]
1990      p= pressure (Pa)
1991      formula:
1992          temp = theta*(p/p0)**(R/Cp)
1993
1994    """
1995    fname = 'compute_td'
1996
1997    tk = (t)*(p/fdef.module_definitions.p0ref)**(fdef.module_definitions.rcp)
1998    data1 = 10.*0.6112*np.exp(17.67*(tk-273.16)/(tk-29.65))
1999    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
2000
2001    rh = qv/data2
2002               
2003    pa = rh * data1
2004    td = 257.44*np.log(pa/6.1121)/(18.678-np.log(pa/6.1121))
2005
2006    return td
2007
2008def var_td_Uhus(t, p, hus):
2009    """ Function to compute dew-point air temperature from temperature and pressure values using hus
2010      t= temperature [K]
2011      hus= specific humidity [1]
2012      formula:
2013          temp = theta*(p/p0)**(R/Cp)
2014
2015    """
2016    fname = 'compute_td'
2017
2018    tk = (t)*(p/fdef.module_definitions.p0ref)**(fdef.module_definitions.rcp)
2019    data1 = 10.*0.6112*np.exp(17.67*(tk-273.16)/(tk-29.65))
2020    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
2021
2022    qv = hus/(1.+hus)
2023
2024    rh = qv/data2
2025               
2026    pa = rh * data1
2027    td = 257.44*np.log(pa/6.1121)/(18.678-np.log(pa/6.1121))
2028
2029    return td
2030
2031def var_wd(u, v):
2032    """ Function to compute the wind direction
2033      [u]= W-E wind direction [ms-1, knot, ...]
2034      [v]= N-S wind direction [ms-1, knot, ...]
2035    """
2036    fname = 'var_wd'
2037
2038    theta = np.arctan2(v,u)
2039    theta = np.where(theta < 0., theta + 2.*np.pi, theta)
2040
2041    wd = 360.*theta/(2.*np.pi)
2042
2043    return wd
2044
2045def var_ws(u, v):
2046    """ Function to compute the wind speed
2047      [u]= W-E wind direction [ms-1, knot, ...]
2048      [v]= N-S wind direction [ms-1, knot, ...]
2049    """
2050    fname = 'var_ws'
2051
2052    ws = np.sqrt(u*u + v*v)
2053
2054    return ws
2055
2056class C_diagnostic(object):
2057    """ Class to compute generic variables
2058      Cdiag: name of the diagnostic to compute
2059      ncobj: netcdf object with data
2060      sfcvars: dictionary with CF equivalencies of surface variables inside file
2061      vars3D: dictionary with CF equivalencies of 3D variables inside file
2062      dictdims: dictionary with CF equivalencies of dimensions inside file
2063        self.values = Values of the diagnostic
2064        self.dims = Dimensions of the diagnostic
2065        self.units = units of the diagnostic
2066        self.incvars = list of variables from the input netCDF object
2067    """ 
2068    def __init__(self, Cdiag, ncobj, sfcvars, vars3D, dictdims):
2069        fname = 'C_diagnostic'
2070        self.values = None
2071        self.dims = None
2072        self.incvars = ncobj.variables
2073        self.units = None
2074
2075        if Cdiag == 'hur':
2076            """ Computing relative humidity
2077            """
2078            vn = 'hur'
2079            CF3Dvars = ['ta', 'plev', 'q']
2080            for v3D in CF3Dvars:
2081                if not vars3D.has_key(v3D):
2082                    print gen.errormsg
2083                    print '  ' + fname + ": missing variable '" + v3D +              \
2084                      "' attribution to compute '" + vn + "' !!"
2085                    print '  Equivalence of 3D variables provided _______'
2086                    gen.printing_dictionary(vars3D)
2087                    quit(-1)
2088                if not self.incvars.has_key(vars3D[v3D]):
2089                    print gen.errormsg
2090                    print '  ' + fname + ": missing variable '" + vars3D[v3D] +      \
2091                      "' in input file to compute '" + vn + "' !!"
2092                    print '  available variables:', self.incvars.keys()
2093                    print '  looking for variables _______' 
2094                    gen.printing_dictionary(vars3D)
2095                    quit(-1)
2096
2097            ta = ncobj.variables[vars3D['ta']][:]
2098            p = ncobj.variables[vars3D['plev']][:]
2099            q = ncobj.variables[vars3D['q']][:]
2100
2101            if len(ta.shape) != len(p.shape):
2102                p = gen.fill_Narray(p, ta*0., filldim=[0,2,3])
2103
2104            self.values = var_hur(p, ta, q)
2105            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2106              dictdims['lon']]
2107            self.units = '1'
2108
2109        if Cdiag == 'hur_Uhus':
2110            """ Computing relative humidity using hus
2111            """
2112            vn = 'hur'
2113            CF3Dvars = ['ta', 'plev', 'hus']
2114            for v3D in CF3Dvars:
2115                if not vars3D.has_key(v3D):
2116                    print gen.errormsg
2117                    print '  ' + fname + ": missing variable '" + v3D +              \
2118                      "' attribution to compute '" + vn + "' !!"
2119                    print '  Equivalence of 3D variables provided _______'
2120                    gen.printing_dictionary(vars3D)
2121                    quit(-1)
2122                if not self.incvars.has_key(vars3D[v3D]):
2123                    print gen.errormsg
2124                    print '  ' + fname + ": missing variable '" + vars3D[v3D] +      \
2125                      "' in input file to compute '" + vn + "' !!"
2126                    print '  available variables:', self.incvars.keys()
2127                    print '  looking for variables _______' 
2128                    gen.printing_dictionary(vars3D)
2129                    quit(-1)
2130
2131            ta = ncobj.variables[vars3D['ta']][:]
2132            p = ncobj.variables[vars3D['plev']][:]
2133            hus = ncobj.variables[vars3D['hus']][:]
2134
2135            if len(ta.shape) != len(p.shape):
2136                p = gen.fill_Narray(p, ta*0., filldim=[0,2,3])
2137
2138            self.values = var_hur_Uhus(p, ta, hus)
2139            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2140              dictdims['lon']]
2141            self.units = '1'
2142
2143        elif Cdiag == 'td':
2144            """ Computing dew-point temperature
2145            """
2146            vn = 'td'
2147            CF3Dvars = ['ta', 'plev', 'q']
2148            for v3D in CF3Dvars:
2149                if not vars3D.has_key(v3D):
2150                    print gen.errormsg
2151                    print '  ' + fname + ": missing variable '" + v3D +              \
2152                      "' attribution to compute '" + vn + "' !!"
2153                    print '  Equivalence of 3D variables provided _______'
2154                    gen.printing_dictionary(vars3D)
2155                    quit(-1)
2156                if not self.incvars.has_key(vars3D[v3D]):
2157                    print gen.errormsg
2158                    print '  ' + fname + ": missing variable '" + vars3D[v3D] +      \
2159                      "' in input file to compute '" + vn + "' !!"
2160                    print '  available variables:', self.incvars.keys()
2161                    print '  looking for variables _______' 
2162                    gen.printing_dictionary(vars3D)
2163                    quit(-1)
2164
2165            ta = ncobj.variables[vars3D['ta']][:]
2166            p = ncobj.variables[vars3D['plev']][:]
2167            q = ncobj.variables[vars3D['q']][:]
2168
2169            if len(ta.shape) != len(p.shape):
2170                p = gen.fill_Narray(p, ta*0., filldim=[0,2,3])
2171
2172            self.values = var_td(ta, p, q)
2173            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2174              dictdims['lon']]
2175            self.units = 'K'
2176
2177        elif Cdiag == 'td_Uhus':
2178            """ Computing dew-point temperature
2179            """
2180            vn = 'td'
2181            CF3Dvars = ['ta', 'plev', 'hus']
2182            for v3D in CF3Dvars:
2183                if not vars3D.has_key(v3D):
2184                    print gen.errormsg
2185                    print '  ' + fname + ": missing variable '" + v3D +              \
2186                      "' attribution to compute '" + vn + "' !!"
2187                    print '  Equivalence of 3D variables provided _______'
2188                    gen.printing_dictionary(vars3D)
2189                    quit(-1)
2190                if not self.incvars.has_key(vars3D[v3D]):
2191                    print gen.errormsg
2192                    print '  ' + fname + ": missing variable '" + vars3D[v3D] +      \
2193                      "' in input file to compute '" + vn + "' !!"
2194                    print '  available variables:', self.incvars.keys()
2195                    print '  looking for variables _______' 
2196                    gen.printing_dictionary(vars3D)
2197                    quit(-1)
2198
2199            ta = ncobj.variables[vars3D['ta']][:]
2200            p = ncobj.variables[vars3D['plev']][:]
2201            hus = ncobj.variables[vars3D['hus']][:]
2202
2203            if len(ta.shape) != len(p.shape):
2204                p = gen.fill_Narray(p, ta*0., filldim=[0,2,3])
2205
2206            self.values = var_td_Uhus(ta, p, hus)
2207            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2208              dictdims['lon']]
2209            self.units = 'K'
2210
2211        elif Cdiag == 'wd':
2212            """ Computing wind direction
2213            """
2214            vn = 'wd'
2215            CF3Dvars = ['ua', 'va']
2216            for v3D in CF3Dvars:
2217                if not vars3D.has_key(v3D):
2218                    print gen.errormsg
2219                    print '  ' + fname + ": missing variable '" + v3D +              \
2220                      "self.' attribution to compute '" + vn + "' !!"
2221                    print '  Equivalence of 3D variables provided _______'
2222                    gen.printing_dictionary(vars3D)
2223                    quit(-1)
2224                if not self.incvars.has_key(vars3D[v3D]):
2225                    print gen.errormsg
2226                    print '  ' + fname + ": missing variable '" + vars3D[v3D] +      \
2227                      "' in input file to compute '" + vn + "' !!"
2228                    print '  available variables:', self.incvars.keys()
2229                    print '  looking for variables _______' 
2230                    gen.printing_dictionary(vars3D)
2231                    quit(-1)
2232   
2233            ua = ncobj.variables[vars3D['ua']][:]
2234            va = ncobj.variables[vars3D['va']][:]
2235   
2236            self.values = var_wd(ua, va)
2237            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2238              dictdims['lon']]
2239            self.units = 'degree'
2240
2241        elif Cdiag == 'ws':
2242            """ Computing wind speed
2243            """
2244            vn = 'ws'
2245            CF3Dvars = ['ua', 'va']
2246            for v3D in CF3Dvars:
2247                if not vars3D.has_key(v3D):
2248                    print gen.errormsg
2249                    print '  ' + fname + ": missing variable '" + v3D +              \
2250                      "' attribution to compute '" + vn + "' !!"
2251                    print '  Equivalence of 3D variables provided _______'
2252                    gen.printing_dictionary(vars3D)
2253                    quit(-1)
2254                if not self.incvars.has_key(vars3D[v3D]):
2255                    print gen.errormsg
2256                    print '  ' + fname + ": missing variable '" + vars3D[v3D] +      \
2257                      "' in input file to compute '" + vn + "' !!"
2258                    print '  available variables:', self.incvars.keys()
2259                    print '  looking for variables _______' 
2260                    gen.printing_dictionary(vars3D)
2261                    quit(-1)
2262
2263            ua = ncobj.variables[vars3D['ua']][:]
2264            va = ncobj.variables[vars3D['va']][:]
2265
2266            self.values = var_ws(ua, va)
2267            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2268              dictdims['lon']]
2269            self.units = ncobj.variables[vars3D['ua']].units
2270
2271        else:
2272            print gen.errormsg
2273            print '  ' + fname + ": variable '" + Wdiag + "' not ready !!"
2274            print '  available ones:', Cavailablediags
2275            quit(-1)
2276
2277class W_diagnostic(object):
2278    """ Class to compute WRF diagnostics variables
2279      Wdiag: name of the diagnostic to compute
2280      ncobj: netcdf object with data
2281      sfcvars: dictionary with CF equivalencies of surface variables inside file
2282      vars3D: dictionary with CF equivalencies of 3D variables inside file
2283      indims: list of dimensions inside file
2284      invardims: list of dimension-variables inside file
2285      dictdims: dictionary with CF equivalencies of dimensions inside file
2286        self.values = Values of the diagnostic
2287        self.dims = Dimensions of the diagnostic
2288        self.units = units of the diagnostic
2289        self.incvars = list of variables from the input netCDF object
2290    """   
2291    def __init__(self, Wdiag, ncobj, sfcvars, vars3D, indims, invardims, dictdims):
2292        fname = 'W_diagnostic'
2293
2294        self.values = None
2295        self.dims = None
2296        self.incvars = ncobj.variables
2297        self.units = None
2298
2299        if Wdiag == 'hur':
2300            """ Computing relative humidity
2301            """
2302            vn = 'hur'
2303            CF3Dvars = ['ta', 'hus']
2304            for v3D in CF3Dvars:
2305                if not vars3D.has_key(v3D):
2306                    print gen.errormsg
2307                    print '  ' + fname + ": missing variable '" + v3D +              \
2308                      "' attribution to compute '" + vn + "' !!"
2309                    print '  Equivalence of 3D variables provided _______'
2310                    gen.printing_dictionary(vars3D)
2311                    quit(-1)
2312
2313            ta = ncobj.variables['T'][:]
2314            p = ncobj.variables['P'][:] + ncobj.variables['PB'][:]
2315            hur = ncobj.variables['QVAPOR'][:]
2316   
2317            vals, dims, vdims = compute_WRFhur(ta, p, hur, indims, invardims)
2318            self.values = vals
2319            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2320              dictdims['lon']]
2321            self.units = '1'
2322
2323        elif Wdiag == 'p':
2324            """ Computing air pressure
2325            """
2326            vn = 'p'
2327
2328            self.values = ncobj.variables['PB'][:] + ncobj.variables['P'][:]
2329            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2330              dictdims['lon']]
2331            self.units = ncobj.variables['PB'].units
2332
2333        elif Wdiag == 'ta':
2334            """ Computing air temperature
2335            """
2336            vn = 'ta'
2337            CF3Dvars = ['ta']
2338            for v3D in CF3Dvars:
2339                if not vars3D.has_key(v3D):
2340                    print gen.errormsg
2341                    print '  ' + fname + ": missing variable '" + v3D +              \
2342                      "' attribution to compute '" + vn + "' !!"
2343                    print '  Equivalence of 3D variables provided _______'
2344                    gen.printing_dictionary(vars3D)
2345                    quit(-1)
2346
2347            ta = ncobj.variables['T'][:]
2348            p = ncobj.variables['P'][:] + ncobj.variables['PB'][:]
2349   
2350            vals, dims, vdims = compute_WRFta(ta, p, indims, invardims)
2351            self.values = vals
2352            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2353              dictdims['lon']]
2354            self.units = 'K'
2355
2356        elif Wdiag == 'td':
2357            """ Computing dew-point temperature
2358            """
2359            vn = 'td'
2360            CF3Dvars = ['ta', 'hus']
2361            for v3D in CF3Dvars:
2362                if not vars3D.has_key(v3D):
2363                    print gen.errormsg
2364                    print '  ' + fname + ": missing variable '" + v3D +              \
2365                      "' attribution to compute '" + vn + "' !!"
2366                    print '  Equivalence of 3D variables provided _______'
2367                    gen.printing_dictionary(vars3D)
2368                    quit(-1)
2369
2370            ta = ncobj.variables['T'][:]
2371            p = ncobj.variables['P'][:] + ncobj.variables['PB'][:]
2372            hus = ncobj.variables['QVAPOR'][:]
2373   
2374            vals, dims, vdims = compute_WRFtd(ta, p, hus, indims, invardims)
2375            self.values = vals
2376            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2377              dictdims['lon']]
2378            self.units = 'K'
2379   
2380        elif Wdiag == 'ua':
2381            """ Computing x-wind
2382            """
2383            vn = 'ua'
2384            CF3Dvars = ['ua', 'va']
2385            for v3D in CF3Dvars:
2386                if not vars3D.has_key(v3D):
2387                    print gen.errormsg
2388                    print '  ' + fname + ": missing variable '" + v3D +              \
2389                      "' attribution to compute '" + vn + "' !!"
2390                    print '  Equivalence of 3D variables provided _______'
2391                    gen.printing_dictionary(vars3D)
2392                    quit(-1)
2393
2394            ua = ncobj.variables['U'][:]
2395            va = ncobj.variables['V'][:]
2396            sina = ncobj.variables['SINALPHA'][:]
2397            cosa = ncobj.variables['COSALPHA'][:]
2398   
2399            vals, dims, vdims = compute_WRFua(ua, va, sina, cosa, indims, invardims)
2400            self.values = vals
2401            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2402              dictdims['lon']]
2403            self.units = ncobj.variables['U'].units
2404   
2405        elif Wdiag == 'uas':
2406            """ Computing 10m x-wind
2407            """
2408            vn = 'uas'
2409            CFsfcvars = ['uas', 'vas']
2410            for vsf in CFsfcvars:
2411                if not sfcvars.has_key(vsf):
2412                    print gen.errormsg
2413                    print '  ' + fname + ": missing variable '" + vsf +              \
2414                      "' attribution to compute '" + vn + "' !!"
2415                    print '  Equivalence of sfc variables provided _______'
2416                    gen.printing_dictionary(sfcvars)
2417                    quit(-1)
2418   
2419            uas = ncobj.variables['U10'][:]
2420            vas = ncobj.variables['V10'][:]
2421            sina = ncobj.variables['SINALPHA'][:]
2422            cosa = ncobj.variables['COSALPHA'][:]
2423   
2424            vals,dims,vdims = compute_WRFuas(uas, vas, sina, cosa, indims, invardims)
2425            self.values = vals
2426            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2427              dictdims['lon']]
2428            self.units = ncobj.variables['U10'].units
2429   
2430        elif Wdiag == 'va':
2431            """ Computing y-wind
2432            """
2433            vn = 'ua'
2434            CF3Dvars = ['ua', 'va']
2435            for v3D in CF3Dvars:
2436                if not vars3D.has_key(v3D):
2437                    print gen.errormsg
2438                    print '  ' + fname + ": missing variable '" + v3D +              \
2439                      "' attribution to compute '" + vn + "' !!"
2440                    print '  Equivalence of 3D variables provided _______'
2441                    gen.printing_dictionary(vars3D)
2442                    quit(-1)
2443   
2444            ua = ncobj.variables['U'][:]
2445            va = ncobj.variables['V'][:]
2446            sina = ncobj.variables['SINALPHA'][:]
2447            cosa = ncobj.variables['COSALPHA'][:]
2448   
2449            vals, dims, vdims = compute_WRFva(ua, va, sina, cosa, indims, invardims)
2450            self.values = vals
2451            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2452              dictdims['lon']]
2453            self.units = ncobj.variables['U'].units
2454   
2455        elif Wdiag == 'vas':
2456            """ Computing 10m y-wind
2457            """
2458            vn = 'uas'
2459            CFsfcvars = ['uas', 'vas']
2460            for vsf in CFsfcvars:
2461                if not sfcvars.has_key(vsf):
2462                    print gen.errormsg
2463                    print '  ' + fname + ": missing variable '" + vsf +              \
2464                      "' attribution to compute '" + vn + "' !!"
2465                    print '  Equivalence of sfc variables provided _______'
2466                    gen.printing_dictionary(sfcvars)
2467                    quit(-1)
2468   
2469            uas = ncobj.variables['U10'][:]
2470            vas = ncobj.variables['V10'][:]
2471            sina = ncobj.variables['SINALPHA'][:]
2472            cosa = ncobj.variables['COSALPHA'][:]
2473   
2474            vals,dims,vdims = compute_WRFvas(uas, vas, sina, cosa, indims, invardims)
2475            self.values = vals
2476            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2477              dictdims['lon']]
2478            self.units = ncobj.variables['U10'].units
2479
2480        elif Wdiag == 'wd':
2481            """ Computing wind direction
2482            """
2483            vn = 'wd'
2484            CF3Dvars = ['ua', 'va']
2485            for v3D in CF3Dvars:
2486                if not vars3D.has_key(v3D):
2487                    print gen.errormsg
2488                    print '  ' + fname + ": missing variable '" + v3D +              \
2489                      "' attribution to compute '" + vn + "' !!"
2490                    print '  Equivalence of 3D variables provided _______'
2491                    gen.printing_dictionary(vars3D)
2492                    quit(-1)
2493
2494            ua = ncobj.variables['U10'][:]
2495            va = ncobj.variables['V10'][:]
2496            sina = ncobj.variables['SINALPHA'][:]
2497            cosa = ncobj.variables['COSALPHA'][:]
2498   
2499            vals, dims, vdims = compute_WRFwd(ua, va, sina, cosa, indims, invardims)
2500            self.values = vals
2501            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2502              dictdims['lon']]
2503            self.units = 'degree'
2504
2505        elif Wdiag == 'ws':
2506            """ Computing wind speed
2507            """
2508            vn = 'ws'
2509            CF3Dvars = ['ua', 'va']
2510            for v3D in CF3Dvars:
2511                if not vars3D.has_key(v3D):
2512                    print gen.errormsg
2513                    print '  ' + fname + ": missing variable '" + v3D +              \
2514                      "' attribution to compute '" + vn + "' !!"
2515                    print '  Equivalence of 3D variables provided _______'
2516                    gen.printing_dictionary(vars3D)
2517                    quit(-1)
2518   
2519            ua = ncobj.variables['U10'][:]
2520            va = ncobj.variables['V10'][:]
2521   
2522            self.values = var_ws(ua, va)
2523            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2524              dictdims['lon']]
2525            self.units = ncobj.variables['U10'].units
2526
2527        elif Wdiag == 'zg':
2528            """ Computing geopotential
2529            """
2530            vn = 'zg'
2531
2532            self.values = ncobj.variables['PHB'][:] + ncobj.variables['PH'][:]
2533            self.dims = [dictdims['time'], dictdims['plev'], dictdims['lat'],        \
2534              dictdims['lon']]
2535            self.units = ncobj.variables['PHB'].units
2536
2537        else:
2538            print gen.errormsg
2539            print '  ' + fname + ": variable '" + Wdiag + "' not ready !!"
2540            print '  available ones:', Wavailablediags
2541            quit(-1)
Note: See TracBrowser for help on using the repository browser.