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

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

Adding:

  • `WRFbnds': grid cell bounds from XLONG_U, XLAT_U, XLONG_V, XLAT_V as intersection of their related parallels and meridians
File size: 82.4 KB
Line 
1# Python script to comput diagnostics
2# From L. Fita work in different places: CCRC (Australia), LMD (France)
3# More information at: http://www.xn--llusfb-5va.cat/python/PyNCplot
4#
5# pyNCplot and its component nc_var.py comes with ABSOLUTELY NO WARRANTY.
6# This work is licendes under a Creative Commons
7#   Attribution-ShareAlike 4.0 International License (http://creativecommons.org/licenses/by-sa/4.0)
8#
9# L. Fita, CIMA. CONICET-UBA, CNRS UMI-IFAECI, C.A. Buenos Aires, Argentina
10# File diagnostics.inf provides the combination of variables to get the desired diagnostic
11#   To be used with module_ForDiagnostics.F90, module_ForDiagnosticsVars.F90, module_generic.F90
12#      foudre: f2py -m module_ForDiagnostics --f90exec=/usr/bin/gfortran-4.7 -c module_generic.F90 module_ForDiagnosticsVars.F90 module_ForDiagnostics.F90 >& run_f2py.log
13#      ciclad: f2py --f90flags="-fPIC" --f90exec=/usr/bin/gfortran -L/opt/canopy-1.3.0/Canopy_64bit/System/lib/ -L/usr/lib64/ -L/opt/canopy-1.3.0/Canopy_64bit/System/lib/ -m module_ForDiagnostics -c module_generic.F90 module_ForDiagnosticsVars.F90 module_ForDiagnostics.F90 >& run_f2py.log
14
15## e.g. # diagnostics.py -d 'Time@WRFtime,bottom_top@ZNU,south_north@XLAT,west_east@XLONG' -v 'clt|CLDFRA,cllmh|CLDFRA@WRFp,RAINTOT|RAINC@RAINNC@RAINSH@XTIME' -f WRF_LMDZ/NPv31/wrfout_d01_1980-03-01_00:00:00
16## e.g. # diagnostics.py -f /home/lluis/PY/diagnostics.inf -d variable_combo -v WRFprc
17
18# Available general pupose diagnostics (model independent) providing (varv1, varv2, ..., dimns, dimvns)
19# compute_accum: Function to compute the accumulation of a variable
20# compute_cllmh: Function to compute cllmh: low/medium/hight cloud fraction following
21#   newmicro.F90 from LMDZ compute_clt(cldfra, pres, dimns, dimvns)
22# compute_clt: Function to compute the total cloud fraction following 'newmicro.F90' from LMDZ
23# compute_clivi: Function to compute cloud-ice water path (clivi)
24# compute_clwvl: Function to compute condensed water path (clwvl)
25# compute_deaccum: Function to compute the deaccumulation of a variable
26# compute_mslp: Function to compute mslp: mean sea level pressure following p_interp.F90 from WRF
27# compute_OMEGAw: Function to transform OMEGA [Pas-1] to velocities [ms-1]
28# compute_prw: Function to compute water vapour path (prw)
29# compute_range_faces: Function to compute faces [uphill, valley, downhill] of sections of a mountain
30#   range, along a given face
31# compute_rh: Function to compute relative humidity following 'Tetens' equation (T,P) ...'
32# compute_td: Function to compute the dew point temperature
33# compute_turbulence: Function to compute the rubulence term of the Taylor's decomposition ...'
34# compute_wds: Function to compute the wind direction
35# compute_wss: Function to compute the wind speed
36# compute_WRFuava: Function to compute geographical rotated WRF 3D winds
37# compute_WRFuasvas: Fucntion to compute geographical rotated WRF 2-meter winds
38# derivate_centered: Function to compute the centered derivate of a given field
39# def Forcompute_cllmh: Function to compute cllmh: low/medium/hight cloud fraction following newmicro.F90 from LMDZ via Fortran subroutine
40# Forcompute_clt: Function to compute the total cloud fraction following 'newmicro.F90' from LMDZ via a Fortran module
41# Forcompute_psl_ptarget: Function to compute the sea-level pressure following target_pressure value found in `p_interp.F'
42
43# Others just providing variable values
44# var_cllmh: Fcuntion to compute cllmh on a 1D column
45# var_clt: Function to compute the total cloud fraction following 'newmicro.F90' from LMDZ using 1D vertical column values
46# var_mslp: Fcuntion to compute mean sea-level pressure
47# var_virtualTemp: This function returns virtual temperature in K,
48# var_WRFtime: Function to copmute CFtimes from WRFtime variable
49# rotational_z: z-component of the rotatinoal of horizontal vectorial field
50# turbulence_var: Function to compute the Taylor's decomposition turbulence term from a a given variable
51# timeoverthres: When a given variable [varname] overpass a given [value]. Being [CFvarn] the name of the diagnostics in
52#   variables_values.dat
53# timemax ([varname], time). When a given variable [varname] got its maximum
54
55from optparse import OptionParser
56import numpy as np
57from netCDF4 import Dataset as NetCDFFile
58import os
59import re
60import nc_var_tools as ncvar
61import generic_tools as gen
62import datetime as dtime
63import module_ForDiag as fdin
64import module_ForDef as fdef
65import diag_tools as diag
66
67main = 'diagnostics.py'
68errormsg = 'ERROR -- error -- ERROR -- error'
69warnmsg = 'WARNING -- warning -- WARNING -- warning'
70
71# Constants
72grav = 9.81
73
74
75####### ###### ##### #### ### ## #
76comboinf="\nIF -d 'variable_combo', provides information of the combination to obtain -v [varn] with the ASCII file with the combinations as -f [combofile]"
77
78parser = OptionParser()
79parser.add_option("-f", "--netCDF_file", dest="ncfile", help="file to use", metavar="FILE")
80parser.add_option("-d", "--dimensions", dest="dimns", 
81  help="[dimtn]@[dtvn],[dimzn]@[dzvn],[...,[dimxn]@[dxvn]], ',' list with the couples [dimDn]@[dDvn], [dimDn], name of the dimension D and name of the variable [dDvn] with the values of the dimension ('WRFtime', for WRF time copmutation). NOTE: same order as in file!!!!" + comboinf, 
82  metavar="LABELS")
83parser.add_option("-v", "--variables", dest="varns", 
84  help=" [varn1]|[var11]@[...[varN1]],[...,[varnM]|[var1M]@[...[varLM]]] ',' list of variables to compute [varnK] and its necessary ones [var1K]...[varPK]", metavar="VALUES")
85
86(opts, args) = parser.parse_args()
87
88#######    #######
89## MAIN
90    #######
91availdiags = ['ACRAINTOT', 'accum', 'clt', 'cllmh', 'convini', 'deaccum', 'fog_K84', \
92  'fog_RUC', 'LMDZrh', 'mslp', 'OMEGAw', 'RAINTOT', 'range_faces',                   \
93  'rvors', 'td', 'timemax', 'timeoverthres', 'turbulence', 'uavaFROMwswd',           \
94  'WRFbnds', 'WRFcape_afwa', 'WRFclivi', 'WRFclwvi', 'WRF_denszint', 'WRFgeop',      \
95  'WRFmrso', 'WRFmrsos', 'WRFpotevap_orPM', 'WRFp', 'WRFpsl_ecmwf',                  \
96  'WRFpsl_ptarget', 'WRFrvors', 'WRFslw', 'ws', 'wds', 'wss', 'WRFheight',           \
97  'WRFheightrel', 'WRFtda', 'WRFtdas', 'WRFua', 'WRFva', 'WRFzwind', 'WRFzwind_log', \
98  'WRFzwindMO']
99
100methods = ['accum', 'deaccum']
101
102# Variables not to check
103NONcheckingvars = ['accum', 'cllmh', 'deaccum', 'face', 'LONLATdxdy', 'TSrhs',       \
104  'TStd', 'TSwds', 'TSwss',                                                          \
105  'WRFbils',  'WRFbnds',                                                             \
106  'WRFclivi', 'WRFclwvi', 'WRFdens', 'WRFdx', 'WRFdxdy', 'WRFdxdywps', 'WRFdy',      \
107  'WRFgeop', 'WRFp', 'WRFtd',                                                        \
108  'WRFpos', 'WRFprc', 'WRFprls', 'WRFrh', 'LMDZrh', 'LMDZrhs',                       \
109  'WRFrhs', 'WRFrvors',                                                              \
110  'WRFt', 'WRFtime', 'WRFua', 'WRFva', 'WRFwds', 'WRFwss', 'WRFheight', 'WRFz',      \
111  'WRFzg']
112
113# diagnostics not to check their dependeny
114NONcheckdepvars = ['accum', 'deaccum', 'timeoverthres', 'WRF_denszint',              \
115  'WRFzwind_log', 'WRFzwind', 'WRFzwindMO']
116
117NONchkvardims = ['WRFtime']
118
119ofile = 'diagnostics.nc'
120
121dimns = opts.dimns
122varns = opts.varns
123
124# Special method. knowing variable combination
125##
126if opts.dimns == 'variable_combo':
127    print warnmsg
128    print '  ' + main + ': knowing variable combination !!!'
129    combination = variable_combo(opts.varns,opts.ncfile)
130    print '     COMBO: ' + combination
131    quit(-1)
132
133if opts.ncfile is None:
134    print errormsg
135    print '   ' + main + ": No file provided !!"
136    print '     is mandatory to provide a file -f [filename]'
137    quit(-1)
138
139if opts.dimns is None:
140    print errormsg
141    print '   ' + main + ": No description of dimensions are provided !!"
142    print '     is mandatory to provide description of dimensions as -d [dimn]@[vardimname],... '
143    quit(-1)
144
145if opts.varns is None:
146    print errormsg
147    print '   ' + main + ": No variable to diagnose is provided !!"
148    print '     is mandatory to provide a variable to diagnose as -v [diagn]|[varn1]@... '
149    quit(-1)
150
151if not os.path.isfile(opts.ncfile):
152    print errormsg
153    print '   ' + main + ": file '" + opts.ncfile + "' does not exist !!"
154    quit(-1)
155
156ncobj = NetCDFFile(opts.ncfile, 'r')
157
158# Looking for specific variables that might be use in more than one diagnostic
159WRFgeop_compute = False
160WRFp_compute = False
161WRFt_compute = False
162WRFrh_compute = False
163WRFght_compute = False
164WRFdens_compute = False
165WRFpos_compute = False
166WRFtime_compute = False
167WRFz_compute = False
168WRFdxdy_compute = False
169WRFdxdywps_compute = False
170LONLATdxdy_compute = False
171
172# File creation
173newnc = NetCDFFile(ofile,'w')
174
175# dimensions
176dimvalues = dimns.split(',')
177dnames = []
178dvnames = []
179
180for dimval in dimvalues:
181    dn = dimval.split('@')[0]
182    dnv = dimval.split('@')[1]
183    dnames.append(dn)
184    dvnames.append(dnv)
185    # Is there any dimension-variable which should be computed?
186    if dnv == 'WRFgeop':WRFgeop_compute = True
187    if dnv == 'WRFp': WRFp_compute = True
188    if dnv == 'WRFt': WRFt_compute = True
189    if dnv == 'WRFrh': WRFrh_compute = True
190    if dnv == 'WRFght': WRFght_compute = True
191    if dnv == 'WRFdens': WRFdens_compute = True
192    if dnv == 'WRFpos': WRFpos_compute = True
193    if dnv == 'WRFtime': WRFtime_compute = True
194    if dnv == 'WRFz':WRFz_compute = True
195    if dnv == 'WRFdxdy':WRFdxdy_compute = True
196    if dnv == 'WRFdxdywps':WRFdxdywps_compute = True
197    if dnv == 'LONLATdxdy':LONLATdxdy_compute = True
198
199# diagnostics to compute
200diags = varns.split(',')
201Ndiags = len(diags)
202
203for idiag in range(Ndiags):
204    if diags[idiag].split('|')[1].find('@') == -1:
205        depvars = diags[idiag].split('|')[1]
206        if depvars == 'WRFgeop':WRFgeop_compute = True
207        if depvars == 'WRFp': WRFp_compute = True
208        if depvars == 'WRFt': WRFt_compute = True
209        if depvars == 'WRFrh': WRFrh_compute = True
210        if depvars == 'WRFght': WRFght_compute = True
211        if depvars == 'WRFdens': WRFdens_compute = True
212        if depvars == 'WRFpos': WRFpos_compute = True
213        if depvars == 'WRFtime': WRFtime_compute = True
214        if depvars == 'WRFz': WRFz_compute = True
215    else:
216        depvars = diags[idiag].split('|')[1].split('@')
217        if gen.searchInlist(depvars, 'WRFgeop'): WRFgeop_compute = True
218        if gen.searchInlist(depvars, 'WRFp'): WRFp_compute = True
219        if gen.searchInlist(depvars, 'WRFt'): WRFt_compute = True
220        if gen.searchInlist(depvars, 'WRFrh'): WRFrh_compute = True
221        if gen.searchInlist(depvars, 'WRFght'): WRFght_compute = True
222        if gen.searchInlist(depvars, 'WRFdens'): WRFdens_compute = True
223        if gen.searchInlist(depvars, 'WRFpos'): WRFpos_compute = True
224        if gen.searchInlist(depvars, 'WRFtime'): WRFtime_compute = True
225        if gen.searchInlist(depvars, 'WRFz'): WRFz_compute = True
226        if gen.searchInlist(depvars, 'WRFdxdy'): WRFdxdy_compute = True
227        if gen.searchInlist(depvars, 'WRFdxdywps'): WRFdxdywps_compute = True
228        if gen.searchInlist(depvars, 'LONLATdxdy'): LONLATdxdy_compute = True
229
230# Dictionary with the new computed variables to be able to add them
231dictcompvars = {}
232if WRFgeop_compute:
233    print '  ' + main + ': Retrieving geopotential value from WRF as PH + PHB'
234    dimv = ncobj.variables['PH'].shape
235    WRFgeop = ncobj.variables['PH'][:] + ncobj.variables['PHB'][:]
236
237    # Attributes of the variable
238    Vvals = gen.variables_values('WRFgeop')
239    dictcompvars['WRFgeop'] = {'name': Vvals[0], 'standard_name': Vvals[1],          \
240      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
241
242if WRFp_compute:
243    print '  ' + main + ': Retrieving pressure value from WRF as P + PB'
244    dimv = ncobj.variables['P'].shape
245    WRFp = ncobj.variables['P'][:] + ncobj.variables['PB'][:]
246
247    # Attributes of the variable
248    Vvals = gen.variables_values('WRFp')
249    dictcompvars['WRFgeop'] = {'name': Vvals[0], 'standard_name': Vvals[1],          \
250      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
251
252if WRFght_compute:
253    print '    ' + main + ': computing geopotential height from WRF as PH + PHB ...' 
254    WRFght = ncobj.variables['PH'][:] + ncobj.variables['PHB'][:]
255
256    # Attributes of the variable
257    Vvals = gen.variables_values('WRFght')
258    dictcompvars['WRFgeop'] = {'name': Vvals[0], 'standard_name': Vvals[1],          \
259      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
260
261if WRFrh_compute:
262    print '    ' + main + ": computing relative humidity from WRF as 'Tetens'" +     \
263      ' equation (T,P) ...'
264    p0=100000.
265    p=ncobj.variables['P'][:] + ncobj.variables['PB'][:]
266    tk = (ncobj.variables['T'][:] + 300.)*(p/p0)**(2./7.)
267    qv = ncobj.variables['QVAPOR'][:]
268
269    data1 = 10.*0.6112*np.exp(17.67*(tk-273.16)/(tk-29.65))
270    data2 = 0.622*data1/(0.01*p-(1.-0.622)*data1)
271
272    WRFrh = qv/data2
273
274    # Attributes of the variable
275    Vvals = gen.variables_values('WRFrh')
276    dictcompvars['WRFrh'] = {'name': Vvals[0], 'standard_name': Vvals[1],            \
277      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
278
279if WRFt_compute:
280    print '    ' + main + ': computing temperature from WRF as inv_potT(T + 300) ...'
281    p0=100000.
282    p=ncobj.variables['P'][:] + ncobj.variables['PB'][:]
283
284    WRFt = (ncobj.variables['T'][:] + 300.)*(p/p0)**(2./7.)
285
286    # Attributes of the variable
287    Vvals = gen.variables_values('WRFt')
288    dictcompvars['WRFt'] = {'name': Vvals[0], 'standard_name': Vvals[1],             \
289      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
290
291if WRFdens_compute:
292    print '    ' + main + ': computing air density from WRF as ((MU + MUB) * ' +     \
293      'DNW)/g ...'
294
295# Just we need in in absolute values: Size of the central grid cell
296##    dxval = ncobj.getncattr('DX')
297##    dyval = ncobj.getncattr('DY')
298##    mapfac = ncobj.variables['MAPFAC_M'][:]
299##    area = dxval*dyval*mapfac
300
301    mu = (ncobj.variables['MU'][:] + ncobj.variables['MUB'][:])
302    dnw = ncobj.variables['DNW'][:]
303
304    WRFdens = np.zeros((mu.shape[0], dnw.shape[1], mu.shape[1], mu.shape[2]),        \
305      dtype=np.float)
306    levval = np.zeros((mu.shape[1], mu.shape[2]), dtype=np.float)
307
308    for it in range(mu.shape[0]):
309        for iz in range(dnw.shape[1]):
310            levval.fill(np.abs(dnw[it,iz]))
311            WRFdens[it,iz,:,:] = levval
312            WRFdens[it,iz,:,:] = mu[it,:,:]*WRFdens[it,iz,:,:]/grav
313
314    # Attributes of the variable
315    Vvals = gen.variables_values('WRFdens')
316    dictcompvars['WRFdens'] = {'name': Vvals[0], 'standard_name': Vvals[1],          \
317      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
318
319if WRFpos_compute:
320# WRF positions from the lowest-leftest corner of the matrix
321    print '    ' + main + ': computing position from MAPFAC_M as sqrt(DY*j**2 + ' +  \
322      'DX*x**2)*MAPFAC_M ...'
323
324    mapfac = ncobj.variables['MAPFAC_M'][:]
325
326    distx = np.float(ncobj.getncattr('DX'))
327    disty = np.float(ncobj.getncattr('DY'))
328
329    print 'distx:',distx,'disty:',disty
330
331    dx = mapfac.shape[2]
332    dy = mapfac.shape[1]
333    dt = mapfac.shape[0]
334
335    WRFpos = np.zeros((dt, dy, dx), dtype=np.float)
336
337    for i in range(1,dx):
338        WRFpos[0,0,i] = distx*i/mapfac[0,0,i]
339    for j in range(1,dy):
340        i=0
341        WRFpos[0,j,i] = WRFpos[0,j-1,i] + disty/mapfac[0,j,i]
342        for i in range(1,dx):
343#            WRFpos[0,j,i] = np.sqrt((disty*j)**2. + (distx*i)**2.)/mapfac[0,j,i]
344#            WRFpos[0,j,i] = np.sqrt((disty*j)**2. + (distx*i)**2.)
345             WRFpos[0,j,i] = WRFpos[0,j,i-1] + distx/mapfac[0,j,i]
346
347    for it in range(1,dt):
348        WRFpos[it,:,:] = WRFpos[0,:,:]
349
350if WRFtime_compute:
351    print '    ' + main + ': computing time from WRF as CFtime(Times) ...'
352
353    refdate='19491201000000'
354    tunitsval='minutes'
355
356    timeobj = ncobj.variables['Times']
357    timewrfv = timeobj[:]
358
359    yrref=refdate[0:4]
360    monref=refdate[4:6]
361    dayref=refdate[6:8]
362    horref=refdate[8:10]
363    minref=refdate[10:12]
364    secref=refdate[12:14]
365
366    refdateS = yrref + '-' + monref + '-' + dayref + ' ' + horref + ':' + minref +   \
367      ':' + secref
368
369
370    if len(timeobj.shape) == 2:
371        dt = timeobj.shape[0]
372    else:
373        dt = 1
374    WRFtime = np.zeros((dt), dtype=np.float)
375
376    if len(timeobj.shape) == 2:
377        for it in range(dt):
378            wrfdates = gen.datetimeStr_conversion(timewrfv[it,:],'WRFdatetime',      \
379              'matYmdHMS')
380            WRFtime[it] = gen.realdatetime1_CFcompilant(wrfdates, refdate, tunitsval)
381    else:
382        wrfdates = gen.datetimeStr_conversion(timewrfv[:],'WRFdatetime',             \
383          'matYmdHMS')
384        WRFtime[0] = gen.realdatetime1_CFcompilant(wrfdates, refdate, tunitsval)
385
386    tunits = tunitsval + ' since ' + refdateS
387
388    # Attributes of the variable
389    dictcompvars['WRFtime'] = {'name': 'time', 'standard_name': 'time',              \
390      'long_name': 'time', 'units': tunits, 'calendar': 'gregorian'}
391
392if WRFz_compute:
393    print '  ' + main + ': Retrieving z: height above surface value from WRF as ' +  \
394      'unstagger(PH + PHB)/9.8-hgt'
395    dimv = ncobj.variables['PH'].shape
396    WRFzg = (ncobj.variables['PH'][:] + ncobj.variables['PHB'][:])/9.8
397
398    unzgd = (dimv[0], dimv[1]-1, dimv[2], dimv[3])
399    unzg = np.zeros(unzgd, dtype=np.float)
400    unzg = 0.5*(WRFzg[:,0:dimv[1]-1,:,:] + WRFzg[:,1:dimv[1],:,:])
401
402    WRFz = np.zeros(unzgd, dtype=np.float)
403    for iz in range(dimv[1]-1):
404        WRFz[:,iz,:,:] = unzg[:,iz,:,:] - ncobj.variables['HGT'][:]
405
406    # Attributes of the variable
407    Vvals = gen.variables_values('WRFz')
408    dictcompvars['WRFz'] = {'name': Vvals[0], 'standard_name': Vvals[1],             \
409      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
410
411if WRFdxdy_compute:
412    print '  ' + main + ': Retrieving dxdy: real distance between grid points ' +    \
413      'from WRF as dx=(XLONG(i+1)-XLONG(i))*DX/MAPFAC_M, dy=(XLAT(j+1)-XLAT(i))*DY/'+\
414      'MAPFAC_M, ds=sqrt(dx**2+dy**2)'
415    dimv = ncobj.variables['XLONG'].shape
416    WRFlon = ncobj.variables['XLONG'][0,:,:]
417    WRFlat = ncobj.variables['XLAT'][0,:,:]
418    WRFmapfac_m = ncobj.variables['MAPFAC_M'][0,:,:]
419    DX = ncobj.DX
420    DY = ncobj.DY
421
422    dimx = dimv[2]
423    dimy = dimv[1]
424
425    WRFdx = np.zeros((dimy,dimx), dtype=np.float)
426    WRFdy = np.zeros((dimy,dimx), dtype=np.float)
427
428    WRFdx[:,0:dimx-1]=(WRFlon[:,1:dimx]-WRFlon[:,0:dimx-1])*DX/WRFmapfac_m[:,0:dimx-1]
429    WRFdy[0:dimy-1,:]=(WRFlat[1:dimy,:]-WRFlat[0:dimy-1,:])*DY/WRFmapfac_m[0:dimy-1,:]
430    WRFds = np.sqrt(WRFdx**2 + WRFdy**2)
431
432    # Attributes of the variable
433    Vvals = gen.variables_values('WRFdx')
434    dictcompvars['WRFdx'] = {'name': Vvals[0], 'standard_name': Vvals[1],             \
435      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
436    Vvals = gen.variables_values('WRFdy')
437    dictcompvars['WRFdy'] = {'name': Vvals[0], 'standard_name': Vvals[1],             \
438      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
439    Vvals = gen.variables_values('WRFds')
440    dictcompvars['WRFds'] = {'name': Vvals[0], 'standard_name': Vvals[1],            \
441      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
442
443if WRFdxdywps_compute:
444    print '  ' + main + ': Retrieving dxdy: real distance between grid points ' +    \
445      'from wpsWRF as dx=(XLONG_M(i+1)-XLONG_M(i))*DX/MAPFAC_M, ' +                  \
446      'dy=(XLAT_M(j+1)-XLAT_M(i))*DY/MAPFAC_M, ds=sqrt(dx**2+dy**2)'
447    dimv = ncobj.variables['XLONG_M'].shape
448    WRFlon = ncobj.variables['XLONG_M'][0,:,:]
449    WRFlat = ncobj.variables['XLAT_M'][0,:,:]
450    WRFmapfac_m = ncobj.variables['MAPFAC_M'][0,:,:]
451    DX = ncobj.DX
452    DY = ncobj.DY
453
454    dimx = dimv[2]
455    dimy = dimv[1]
456
457    WRFdx = np.zeros((dimy,dimx), dtype=np.float)
458    WRFdy = np.zeros((dimy,dimx), dtype=np.float)
459
460    WRFdx[:,0:dimx-1]=(WRFlon[:,1:dimx]-WRFlon[:,0:dimx-1])*DX/WRFmapfac_m[:,0:dimx-1]
461    WRFdy[0:dimy-1,:]=(WRFlat[1:dimy,:]-WRFlat[0:dimy-1,:])*DY/WRFmapfac_m[0:dimy-1,:]
462    WRFds = np.sqrt(WRFdx**2 + WRFdy**2)
463
464    # Attributes of the variable
465    Vvals = gen.variables_values('WRFdx')
466    dictcompvars['WRFdx'] = {'name': Vvals[0], 'standard_name': Vvals[1],             \
467      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
468    Vvals = gen.variables_values('WRFdy')
469    dictcompvars['WRFdy'] = {'name': Vvals[0], 'standard_name': Vvals[1],             \
470      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
471    Vvals = gen.variables_values('WRFds')
472    dictcompvars['WRFds'] = {'name': Vvals[0], 'standard_name': Vvals[1],            \
473      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
474
475if LONLATdxdy_compute:
476    print '  ' + main + ': Retrieving dxdy: real distance between grid points ' +    \
477      'from a regular lonlat projection as dx=(lon[i+1]-lon[i])*raddeg*Rearth*' +    \
478      'cos(abs(lat[i])); dy=(lat[j+1]-lat[i])*raddeg*Rearth; ds=sqrt(dx**2+dy**2); '+\
479      'raddeg = pi/180; Rearth=6370.0e03'
480    dimv = ncobj.variables['lon'].shape
481    lon = ncobj.variables['lon'][:]
482    lat = ncobj.variables['lat'][:]
483
484    WRFlon, WRFlat = gen.lonlat2D(lon,lat)
485
486    dimx = WRFlon.shape[1]
487    dimy = WRFlon.shape[0]
488
489    WRFdx = np.zeros((dimy,dimx), dtype=np.float)
490    WRFdy = np.zeros((dimy,dimx), dtype=np.float)
491
492    raddeg = np.pi/180.
493
494    Rearth = fdef.module_definitions.earthradii
495
496    WRFdx[:,0:dimx-1]=(WRFlon[:,1:dimx]-WRFlon[:,0:dimx-1])*raddeg*Rearth*           \
497      np.cos(np.abs(WRFlat[:,0:dimx-1]*raddeg))
498    WRFdy[0:dimy-1,:]=(WRFlat[1:dimy,:]-WRFlat[0:dimy-1,:])*raddeg*Rearth
499    WRFds = np.sqrt(WRFdx**2 + WRFdy**2)
500
501    # Attributes of the variable
502    Vvals = gen.variables_values('WRFdx')
503    dictcompvars['WRFdx'] = {'name': Vvals[0], 'standard_name': Vvals[1],             \
504      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
505    Vvals = gen.variables_values('WRFdy')
506    dictcompvars['WRFdy'] = {'name': Vvals[0], 'standard_name': Vvals[1],             \
507      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
508    Vvals = gen.variables_values('WRFds')
509    dictcompvars['WRFds'] = {'name': Vvals[0], 'standard_name': Vvals[1],             \
510      'long_name': Vvals[4].replace('|',' '), 'units': Vvals[5]}
511
512### ## #
513# Going for the diagnostics
514### ## #
515print '  ' + main + ' ...'
516varsadd = []
517
518for idiag in range(Ndiags):
519    print '    diagnostic:',diags[idiag]
520    diagn = diags[idiag].split('|')[0]
521    depvars = diags[idiag].split('|')[1].split('@')
522    if not gen.searchInlist(NONcheckdepvars, diagn):
523        if diags[idiag].split('|')[1].find('@') != -1:
524            depvars = diags[idiag].split('|')[1].split('@')
525            if depvars[0] == 'deaccum': diagn='deaccum'
526            if depvars[0] == 'accum': diagn='accum'
527            for depv in depvars:
528                # Checking without extra arguments of a depending variable (':', separated)
529                if depv.find(':') != -1: depv=depv.split(':')[0]
530                if not ncobj.variables.has_key(depv) and not                         \
531                  gen.searchInlist(NONcheckingvars, depv) and                        \
532                  not gen.searchInlist(methods, depv) and not depvars[0] == 'deaccum'\
533                  and not depvars[0] == 'accum' and not depv[0:2] == 'z=':
534                    print errormsg
535                    print '  ' + main + ": file '" + opts.ncfile +                   \
536                      "' does not have variable '" + depv + "' !!"
537                    quit(-1)
538        else:
539            depvars = diags[idiag].split('|')[1]
540            if not ncobj.variables.has_key(depvars) and not                          \
541              gen.searchInlist(NONcheckingvars, depvars) and                         \
542              not gen.searchInlist(methods, depvars):
543                print errormsg
544                print '  ' + main + ": file '" + opts.ncfile +                       \
545                  "' does not have variable '" + depvars + "' !!"
546                quit(-1)
547
548    print "\n    Computing '" + diagn + "' from: ", depvars, '...'
549
550# acraintot: accumulated total precipitation from WRF RAINC, RAINNC, RAINSH
551    if diagn == 'ACRAINTOT':
552           
553        var0 = ncobj.variables[depvars[0]]
554        var1 = ncobj.variables[depvars[1]]
555        var2 = ncobj.variables[depvars[2]]
556
557        diagout = var0[:] + var1[:] + var2[:]
558
559        dnamesvar = var0.dimensions
560        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
561
562        # Removing the nonChecking variable-dimensions from the initial list
563        varsadd = []
564        for nonvd in NONchkvardims:
565            if gen.searchInlist(dvnamesvar,nonvd): dvnamesvar.remove(nonvd)
566            varsadd.append(nonvd)
567
568        ncvar.insert_variable(ncobj, 'pracc', diagout, dnamesvar, dvnamesvar, newnc)
569
570# accum: acumulation of any variable as (Variable, time [as [tunits]
571#   from/since ....], newvarname)
572    elif diagn == 'accum':
573
574        var0 = ncobj.variables[depvars[0]]
575        var1 = ncobj.variables[depvars[1]]
576
577        dnamesvar = var0.dimensions
578        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
579
580        diagout, diagoutd, diagoutvd = diag.compute_accum(var0,dnamesvar,dvnamesvar)
581        # Removing the nonChecking variable-dimensions from the initial list
582        varsadd = []
583        diagoutvd = list(dvnames)
584        for nonvd in NONchkvardims:
585            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
586            varsadd.append(nonvd)
587
588        CFvarn = ncvar.variables_values(depvars[0])[0]
589
590# Removing the flux
591        if depvars[1] == 'XTIME':
592            dtimeunits = var1.getncattr('description')
593            tunits = dtimeunits.split(' ')[0]
594        else:
595            dtimeunits = var1.getncattr('units')
596            tunits = dtimeunits.split(' ')[0]
597
598        dtime = (var1[1] - var1[0])*diag.timeunits_seconds(tunits)
599
600        ncvar.insert_variable(ncobj, CFvarn + 'acc', diagout*dtime, diagoutd, diagoutvd, newnc)
601
602# cllmh with cldfra, pres
603    elif diagn == 'cllmh':
604           
605        var0 = ncobj.variables[depvars[0]]
606        if depvars[1] == 'WRFp':
607            var1 = WRFp
608        else:
609            var01 = ncobj.variables[depvars[1]]
610            if len(size(var1.shape)) < len(size(var0.shape)):
611                var1 = np.brodcast_arrays(var01,var0)[0]
612            else:
613                var1 = var01
614
615        diagout, diagoutd, diagoutvd = diag.Forcompute_cllmh(var0,var1,dnames,dvnames)
616
617        # Removing the nonChecking variable-dimensions from the initial list
618        varsadd = []
619        for nonvd in NONchkvardims:
620            if gen.searchInlist(diagoutvd,nonvd): diagoutvd.remove(nonvd)
621            varsadd.append(nonvd)
622
623        ncvar.insert_variable(ncobj, 'cll', diagout[0,:], diagoutd, diagoutvd, newnc)
624        ncvar.insert_variable(ncobj, 'clm', diagout[1,:], diagoutd, diagoutvd, newnc)
625        ncvar.insert_variable(ncobj, 'clh', diagout[2,:], diagoutd, diagoutvd, newnc)
626
627# clt with cldfra
628    elif diagn == 'clt':
629           
630        var0 = ncobj.variables[depvars]
631        diagout, diagoutd, diagoutvd = diag.Forcompute_clt(var0,dnames,dvnames)
632
633        # Removing the nonChecking variable-dimensions from the initial list
634        varsadd = []
635        for nonvd in NONchkvardims:
636            if gen.searchInlist(diagoutvd,nonvd): diagoutvd.remove(nonvd)
637            varsadd.append(nonvd)
638           
639        ncvar.insert_variable(ncobj, 'clt', diagout, diagoutd, diagoutvd, newnc)
640
641# convini (pr, time)
642    elif diagn == 'convini':
643           
644        var0 = ncobj.variables[depvars[0]][:]
645        var1 = ncobj.variables[depvars[1]][:]
646        otime = ncobj.variables[depvars[1]]
647
648        dnamesvar = ncobj.variables[depvars[0]].dimensions
649        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
650
651        diagout, diagoutd, diagoutvd  = diag.var_convini(var0, var1, dnames, dvnames)
652
653        ncvar.insert_variable(ncobj, 'convini', diagout, diagoutd, diagoutvd, newnc, \
654          fill=gen.fillValueF)
655        # Getting the right units
656        ovar = newnc.variables['convini']
657        if gen.searchInlist(otime.ncattrs(), 'units'): 
658            tunits = otime.getncattr('units')
659            ncvar.set_attribute(ovar, 'units', tunits)
660            newnc.sync()
661
662# deaccum: deacumulation of any variable as (Variable, time [as [tunits]
663#   from/since ....], newvarname)
664    elif diagn == 'deaccum':
665
666        var0 = ncobj.variables[depvars[0]]
667        var1 = ncobj.variables[depvars[1]]
668
669        dnamesvar = var0.dimensions
670        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
671
672        diagout, diagoutd, diagoutvd = diag.compute_deaccum(var0,dnamesvar,dvnamesvar)
673        # Removing the nonChecking variable-dimensions from the initial list
674        varsadd = []
675        diagoutvd = list(dvnames)
676        for nonvd in NONchkvardims:
677            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
678            varsadd.append(nonvd)
679
680# Transforming to a flux
681        if depvars[1] == 'XTIME':
682            dtimeunits = var1.getncattr('description')
683            tunits = dtimeunits.split(' ')[0]
684        else:
685            dtimeunits = var1.getncattr('units')
686            tunits = dtimeunits.split(' ')[0]
687
688        dtime = (var1[1] - var1[0])*diag.timeunits_seconds(tunits)
689        ncvar.insert_variable(ncobj, depvars[2], diagout/dtime, diagoutd, diagoutvd, \
690          newnc)
691
692# fog_K84: Computation of fog and visibility following Kunkel, (1984) as QCLOUD, QICE
693    elif diagn == 'fog_K84':
694
695        var0 = ncobj.variables[depvars[0]]
696        var1 = ncobj.variables[depvars[1]]
697
698        dnamesvar = list(var0.dimensions)
699        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
700
701        diag1, diag2, diagoutd, diagoutvd = diag.Forcompute_fog_K84(var0, var1,      \
702          dnamesvar, dvnamesvar)
703        # Removing the nonChecking variable-dimensions from the initial list
704        varsadd = []
705        diagoutvd = list(dvnames)
706        for nonvd in NONchkvardims:
707            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
708            varsadd.append(nonvd)
709        ncvar.insert_variable(ncobj, 'fog', diag1, diagoutd, diagoutvd, newnc)
710        ncvar.insert_variable(ncobj, 'fogvisblty', diag2, diagoutd, diagoutvd, newnc)
711
712# fog_RUC: Computation of fog and visibility following Kunkel, (1984) as QVAPOR,
713#    WRFt, WRFp or Q2, T2, PSFC
714    elif diagn == 'fog_RUC':
715
716        var0 = ncobj.variables[depvars[0]]
717        print gen.infmsg
718        if depvars[1] == 'WRFt':
719            print '  ' + main + ": computing '" + diagn + "' using 3D variables !!"
720            var1 = WRFt
721            var2 = WRFp
722        else:
723            print '  ' + main + ": computing '" + diagn + "' using 2D variables !!"
724            var1 = ncobj.variables[depvars[1]]
725            var2 = ncobj.variables[depvars[2]]
726
727        dnamesvar = list(var0.dimensions)
728        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
729
730        diag1, diag2, diagoutd, diagoutvd = diag.Forcompute_fog_RUC(var0, var1, var2,\
731          dnamesvar, dvnamesvar)
732        # Removing the nonChecking variable-dimensions from the initial list
733        varsadd = []
734        diagoutvd = list(dvnames)
735        for nonvd in NONchkvardims:
736            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
737            varsadd.append(nonvd)
738        ncvar.insert_variable(ncobj, 'fog', diag1, diagoutd, diagoutvd, newnc)
739        ncvar.insert_variable(ncobj, 'fogvisblty', diag2, diagoutd, diagoutvd, newnc)
740
741# fog_FRAML50: Computation of fog and visibility following Gultepe, I. and
742#   J.A. Milbrandt, 2010 as QVAPOR, WRFt, WRFp or Q2, T2, PSFC
743    elif diagn == 'fog_FRAML50':
744
745        var0 = ncobj.variables[depvars[0]]
746        print gen.infmsg
747        if depvars[1] == 'WRFt':
748            print '  ' + main + ": computing '" + diagn + "' using 3D variables !!"
749            var1 = WRFt
750            var2 = WRFp
751        else:
752            print '  ' + main + ": computing '" + diagn + "' using 2D variables !!"
753            var1 = ncobj.variables[depvars[1]]
754            var2 = ncobj.variables[depvars[2]]
755
756        dnamesvar = list(var0.dimensions)
757        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
758
759        diag1, diag2, diagoutd, diagoutvd = diag.Forcompute_fog_FRAML50(var0, var1,  \
760          var2, dnamesvar, dvnamesvar)
761        # Removing the nonChecking variable-dimensions from the initial list
762        varsadd = []
763        diagoutvd = list(dvnames)
764        for nonvd in NONchkvardims:
765            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
766            varsadd.append(nonvd)
767        ncvar.insert_variable(ncobj, 'fog', diag1, diagoutd, diagoutvd, newnc)
768        ncvar.insert_variable(ncobj, 'fogvisblty', diag2, diagoutd, diagoutvd, newnc)
769
770# LMDZrh (pres, t, r)
771    elif diagn == 'LMDZrh':
772           
773        var0 = ncobj.variables[depvars[0]][:]
774        var1 = ncobj.variables[depvars[1]][:]
775        var2 = ncobj.variables[depvars[2]][:]
776
777        diagout, diagoutd, diagoutvd = diag.compute_rh(var0,var1,var2,dnames,dvnames)
778        ncvar.insert_variable(ncobj, 'hur', diagout, diagoutd, diagoutvd, newnc)
779
780# LMDZrhs (psol, t2m, q2m)
781    elif diagn == 'LMDZrhs':
782           
783        var0 = ncobj.variables[depvars[0]][:]
784        var1 = ncobj.variables[depvars[1]][:]
785        var2 = ncobj.variables[depvars[2]][:]
786
787        dnamesvar = ncobj.variables[depvars[0]].dimensions
788        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
789
790        diagout, diagoutd, diagoutvd = diag.compute_rh(var0,var1,var2,dnamesvar,dvnamesvar)
791
792        ncvar.insert_variable(ncobj, 'hurs', diagout, diagoutd, diagoutvd, newnc)
793
794# range_faces: LON, LAT, HGT, WRFdxdy, 'face:['WE'/'SN']:[dsfilt]:[dsnewrange]:[hvalleyrange]'
795    elif diagn == 'range_faces':
796           
797        var0 = ncobj.variables[depvars[0]][:]
798        var1 = ncobj.variables[depvars[1]][:]
799        var2 = ncobj.variables[depvars[2]][:]
800        face = depvars[4].split(':')[1]
801        dsfilt = np.float(depvars[4].split(':')[2])
802        dsnewrange = np.float(depvars[4].split(':')[3])
803        hvalleyrange = np.float(depvars[4].split(':')[4])
804
805        dnamesvar = list(ncobj.variables[depvars[2]].dimensions)
806        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
807        lon, lat = gen.lonlat2D(var0, var1)
808        if len(var2.shape) == 3:
809            print warnmsg
810            print '  ' + diagn + ": shapping to 2D variable '" + depvars[2] + "' !!"
811            hgt = var2[0,:,:]
812            dnamesvar.pop(0)
813            dvnamesvar.pop(0)           
814        else:
815            hgt = var2[:]
816
817        orogmax, ptorogmax, dhgt, peaks, valleys, origfaces, diagout, diagoutd,      \
818          diagoutvd, rng, rngorogmax, ptrngorogmax= diag.Forcompute_range_faces(lon, \
819           lat, hgt, WRFdx, WRFdy, WRFds, face, dsfilt, dsnewrange, hvalleyrange,    \
820           dnamesvar, dvnamesvar)
821
822        # Removing the nonChecking variable-dimensions from the initial list
823        varsadd = []
824        diagoutvd = list(dvnames)
825        for nonvd in NONchkvardims:
826            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
827            varsadd.append(nonvd)
828
829        ncvar.insert_variable(ncobj, 'dx', WRFdx, diagoutd, diagoutvd, newnc)
830        ncvar.insert_variable(ncobj, 'dy', WRFdy, diagoutd, diagoutvd, newnc)
831        ncvar.insert_variable(ncobj, 'ds', WRFds, diagoutd, diagoutvd, newnc)
832
833        # adding variables to output file
834        if face == 'WE': axis = 'lon'
835        elif face == 'SN': axis = 'lat'
836
837        ncvar.insert_variable(ncobj, 'range', rng, diagoutd, diagoutvd, newnc,       \
838          fill=gen.fillValueI)
839        ovar = newnc.variables['range']
840        ncvar.set_attribute(ovar, 'deriv', axis)
841
842        ncvar.insert_variable(ncobj, 'orogmax', rngorogmax, diagoutd, diagoutvd,     \
843          newnc, fill=gen.fillValueF)
844        newnc.renameVariable('orogmax', 'rangeorogmax')
845        ovar = newnc.variables['rangeorogmax']
846        ncvar.set_attribute(ovar, 'deriv', axis)
847        stdn = ovar.standard_name
848        ncvar.set_attribute(ovar, 'standard_name', 'range_' + stdn)
849        Ln = ovar.long_name
850        ncvar.set_attribute(ovar, 'long_name', 'range ' + stdn)
851
852        ncvar.insert_variable(ncobj, 'ptorogmax', ptrngorogmax, diagoutd, diagoutvd, \
853          newnc)
854        newnc.renameVariable('ptorogmax', 'rangeptorogmax')
855        ovar = newnc.variables['rangeptorogmax']
856        ncvar.set_attribute(ovar, 'deriv', axis)
857        stdn = ovar.standard_name
858        ncvar.set_attribute(ovar, 'standard_name', 'range_' + stdn)
859        Ln = ovar.long_name
860        ncvar.set_attribute(ovar, 'long_name', 'range ' + stdn)
861
862        ncvar.insert_variable(ncobj, 'orogmax', orogmax, diagoutd, diagoutvd,        \
863          newnc)
864        ovar = newnc.variables['orogmax']
865        ncvar.set_attribute(ovar, 'deriv', axis)
866
867        ncvar.insert_variable(ncobj, 'ptorogmax', ptorogmax, diagoutd, diagoutvd,    \
868          newnc)
869        ovar = newnc.variables['ptorogmax']
870        ncvar.set_attribute(ovar, 'deriv', axis)
871
872        ncvar.insert_variable(ncobj, 'orogderiv', dhgt, diagoutd, diagoutvd, newnc)
873        ovar = newnc.variables['orogderiv']
874        ncvar.set_attribute(ovar, 'deriv', axis)
875
876        ncvar.insert_variable(ncobj, 'peak', peaks, diagoutd, diagoutvd, newnc)
877        ncvar.insert_variable(ncobj, 'valley', valleys, diagoutd, diagoutvd, newnc)
878
879        ncvar.insert_variable(ncobj, 'rangefaces', diagout, diagoutd, diagoutvd,    \
880          newnc)
881        newnc.renameVariable('rangefaces', 'rangefacesfilt')
882        ovar = newnc.variables['rangefacesfilt']
883        ncvar.set_attribute(ovar, 'face', face)
884        ncvar.set_attributek(ovar, 'dist_filter', dsfilt, 'F')
885
886        ncvar.insert_variable(ncobj, 'rangefaces', origfaces, diagoutd, diagoutvd,  \
887          newnc, fill=gen.fillValueI)
888        ovar = newnc.variables['rangefaces']
889        ncvar.set_attribute(ovar, 'face', face)
890        ncvar.set_attributek(ovar, 'dist_newrange', dsnewrange, 'F')
891        ncvar.set_attributek(ovar, 'h_valley_newrange', hvalleyrange, 'F')
892
893# cell_bnds: grid cell bounds from XLONG_U, XLAT_U, XLONG_V, XLAT_V as intersection
894#   of their related parallels and meridians
895    elif diagn == 'WRFbnds':
896           
897        var0 = ncobj.variables[depvars[0]][0,:,:]
898        var1 = ncobj.variables[depvars[1]][0,:,:]
899        var2 = ncobj.variables[depvars[2]][0,:,:]
900        var3 = ncobj.variables[depvars[3]][0,:,:]
901
902        dnamesvar = []
903        dnamesvar.append(ncobj.variables[depvars[2]].dimensions[2])
904        dnamesvar.append(ncobj.variables[depvars[0]].dimensions[1])
905        dnamesvar.append('bnds')
906        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
907
908        cellbndsx, cellbndsy, diagoutd, diagoutvd = diag.Forcompute_cellbnds(var0,   \
909          var1, var2, var3, dnamesvar, dvnamesvar)
910
911        # Removing the nonChecking variable-dimensions from the initial list
912        varsadd = []
913        diagoutvd = list(dvnames)
914        for nonvd in NONchkvardims:
915            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
916            varsadd.append(nonvd)
917        # creation of bounds dimension
918        newdim = newnc.createDimension('bnds', 4)
919
920        ncvar.insert_variable(ncobj, 'lon_bnds', cellbndsx, diagoutd, diagoutvd, newnc)
921        ncvar.insert_variable(ncobj, 'lat_bnds', cellbndsy, diagoutd, diagoutvd, newnc)
922
923# mrso: total soil moisture SMOIS, DZS
924    elif diagn == 'WRFmrso':
925           
926        var0 = ncobj.variables[depvars[0]][:]
927        var10 = ncobj.variables[depvars[1]][:]
928        dnamesvar = list(ncobj.variables[depvars[0]].dimensions)
929        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
930
931        var1 = var0.copy()*0.
932        var2 = var0.copy()*0.+1.
933        # Must be a better way....
934        for j in range(var0.shape[2]):
935          for i in range(var0.shape[3]):
936              var1[:,:,j,i] = var10
937
938        diagout, diagoutd, diagoutvd = diag.Forcompute_zint(var0, var1, var2,        \
939          dnamesvar, dvnamesvar)
940
941        # Removing the nonChecking variable-dimensions from the initial list
942        varsadd = []
943        diagoutvd = list(dvnames)
944        for nonvd in NONchkvardims:
945            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
946            varsadd.append(nonvd)
947        ncvar.insert_variable(ncobj, 'mrso', diagout, diagoutd, diagoutvd, newnc)
948
949# mrsos: First layer soil moisture SMOIS, DZS
950    elif diagn == 'WRFmrsos':
951           
952        var0 = ncobj.variables[depvars[0]][:]
953        var1 = ncobj.variables[depvars[1]][:]
954        diagoutd = list(ncobj.variables[depvars[0]].dimensions)
955        diagoutvd = ncvar.var_dim_dimv(diagoutd,dnames,dvnames)
956
957        diagoutd.pop(1)
958        diagoutvd.pop(1)
959
960        diagout= np.zeros((var0.shape[0],var0.shape[2],var0.shape[3]), dtype=np.float)
961
962        # Must be a better way....
963        for j in range(var0.shape[2]):
964          for i in range(var0.shape[3]):
965              diagout[:,j,i] = var0[:,0,j,i]*var1[:,0]
966
967        # Removing the nonChecking variable-dimensions from the initial list
968        varsadd = []
969        diagoutvd = list(dvnames)
970        for nonvd in NONchkvardims:
971            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
972            varsadd.append(nonvd)
973        ncvar.insert_variable(ncobj, 'mrsos', diagout, diagoutd, diagoutvd, newnc)
974
975# mslp: mean sea level pressure (pres, psfc, terrain, temp, qv)
976    elif diagn == 'mslp' or diagn == 'WRFmslp':
977           
978        var1 = ncobj.variables[depvars[1]][:]
979        var2 = ncobj.variables[depvars[2]][:]
980        var4 = ncobj.variables[depvars[4]][:]
981
982        if diagn == 'WRFmslp':
983            var0 = WRFp
984            var3 = WRFt
985            dnamesvar = ncobj.variables['P'].dimensions
986        else:
987            var0 = ncobj.variables[depvars[0]][:]
988            var3 = ncobj.variables[depvars[3]][:]
989            dnamesvar = ncobj.variables[depvars[0]].dimensions
990
991        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
992
993        diagout, diagoutd, diagoutvd = diag.compute_mslp(var0, var1, var2, var3, var4,    \
994          dnamesvar, dvnamesvar)
995
996        # Removing the nonChecking variable-dimensions from the initial list
997        varsadd = []
998        diagoutvd = list(dvnames)
999        for nonvd in NONchkvardims:
1000            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1001            varsadd.append(nonvd)
1002        ncvar.insert_variable(ncobj, 'psl', diagout, diagoutd, diagoutvd, newnc)
1003
1004# OMEGAw (omega, p, t) from NCL formulation (https://www.ncl.ucar.edu/Document/Functions/Contributed/omega_to_w.shtml)
1005    elif diagn == 'OMEGAw':
1006           
1007        var0 = ncobj.variables[depvars[0]][:]
1008        var1 = ncobj.variables[depvars[1]][:]
1009        var2 = ncobj.variables[depvars[2]][:]
1010
1011        dnamesvar = ncobj.variables[depvars[0]].dimensions
1012        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1013
1014        diagout, diagoutd, diagoutvd = diag.compute_OMEGAw(var0,var1,var2,dnamesvar,dvnamesvar)
1015
1016        ncvar.insert_variable(ncobj, 'wa', diagout, diagoutd, diagoutvd, newnc)
1017
1018# raintot: instantaneous total precipitation from WRF as (RAINC + RAINC + RAINSH) / dTime
1019    elif diagn == 'RAINTOT':
1020
1021        var0 = ncobj.variables[depvars[0]]
1022        var1 = ncobj.variables[depvars[1]]
1023        var2 = ncobj.variables[depvars[2]]
1024
1025        if depvars[3] != 'WRFtime':
1026            var3 = ncobj.variables[depvars[3]]
1027        else:
1028            var3 = np.arange(var0.shape[0], dtype=int)
1029
1030        var = var0[:] + var1[:] + var2[:]
1031
1032        dnamesvar = var0.dimensions
1033        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1034
1035        diagout, diagoutd, diagoutvd = diag.compute_deaccum(var,dnamesvar,dvnamesvar)
1036
1037# Transforming to a flux
1038        if var3.shape[0] > 1:
1039            if depvars[3] != 'WRFtime':
1040                dtimeunits = var3.getncattr('units')
1041                tunits = dtimeunits.split(' ')[0]
1042   
1043                dtime = (var3[1] - var3[0])*diag.timeunits_seconds(tunits)
1044            else:
1045                var3 = ncobj.variables['Times']
1046                time1 = var3[0,:]
1047                time2 = var3[1,:]
1048                tmf1 = ''
1049                tmf2 = ''
1050                for ic in range(len(time1)):
1051                    tmf1 = tmf1 + time1[ic]
1052                    tmf2 = tmf2 + time2[ic]
1053                dtdate1 = dtime.datetime.strptime(tmf1,"%Y-%m-%d_%H:%M:%S")
1054                dtdate2 = dtime.datetime.strptime(tmf2,"%Y-%m-%d_%H:%M:%S")
1055                diffdate12 = dtdate2 - dtdate1
1056                dtime = diffdate12.total_seconds()
1057                print 'dtime:',dtime
1058        else:
1059            print warnmsg
1060            print '  ' + main + ": only 1 time-step for '" + diag + "' !!"
1061            print '    leaving a zero value!'
1062            diagout = var0[:]*0.
1063            dtime=1.
1064
1065        # Removing the nonChecking variable-dimensions from the initial list
1066        varsadd = []
1067        for nonvd in NONchkvardims:
1068            if gen.searchInlist(diagoutvd,nonvd): diagoutvd.remove(nonvd)
1069            varsadd.append(nonvd)
1070           
1071        ncvar.insert_variable(ncobj, 'pr', diagout/dtime, diagoutd, diagoutvd, newnc)
1072
1073# timemax ([varname], time). When a given variable [varname] got its maximum
1074    elif diagn == 'timemax':
1075           
1076        var0 = ncobj.variables[depvars[0]][:]
1077        var1 = ncobj.variables[depvars[1]][:]
1078
1079        otime = ncobj.variables[depvars[1]]
1080
1081        dnamesvar = ncobj.variables[depvars[0]].dimensions
1082        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1083
1084        diagout, diagoutd, diagoutvd  = diag.var_timemax(var0, var1, dnames,         \
1085          dvnames)
1086
1087        ncvar.insert_variable(ncobj, 'timemax', diagout, diagoutd, diagoutvd, newnc, \
1088          fill=gen.fillValueF)
1089        # Getting the right units
1090        ovar = newnc.variables['timemax']
1091        if gen.searchInlist(otime.ncattrs(), 'units'): 
1092            tunits = otime.getncattr('units')
1093            ncvar.set_attribute(ovar, 'units', tunits)
1094            newnc.sync()
1095        ncvar.set_attribute(ovar, 'variable', depvars[0])
1096
1097# timeoverthres ([varname], time, [value], [CFvarn]). When a given variable [varname]   
1098#   overpass a given [value]. Being [CFvarn] the name of the diagnostics in
1099#   variables_values.dat
1100    elif diagn == 'timeoverthres':
1101           
1102        var0 = ncobj.variables[depvars[0]][:]
1103        var1 = ncobj.variables[depvars[1]][:]
1104        var2 = np.float(depvars[2])
1105        var3 = depvars[3]
1106
1107        otime = ncobj.variables[depvars[1]]
1108
1109        dnamesvar = ncobj.variables[depvars[0]].dimensions
1110        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1111
1112        diagout, diagoutd, diagoutvd  = diag.var_timeoverthres(var0, var1, var2,     \
1113          dnames, dvnames)
1114
1115        ncvar.insert_variable(ncobj, var3, diagout, diagoutd, diagoutvd, newnc, \
1116          fill=gen.fillValueF)
1117        # Getting the right units
1118        ovar = newnc.variables[var3]
1119        if gen.searchInlist(otime.ncattrs(), 'units'): 
1120            tunits = otime.getncattr('units')
1121            ncvar.set_attribute(ovar, 'units', tunits)
1122            newnc.sync()
1123        ncvar.set_attribute(ovar, 'overpassed_threshold', var2)
1124
1125# rhs (psfc, t, q) from TimeSeries files
1126    elif diagn == 'TSrhs':
1127           
1128        p0=100000.
1129        var0 = ncobj.variables[depvars[0]][:]
1130        var1 = (ncobj.variables[depvars[1]][:])*(var0/p0)**(2./7.)
1131        var2 = ncobj.variables[depvars[2]][:]
1132
1133        dnamesvar = ncobj.variables[depvars[0]].dimensions
1134        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1135
1136        diagout, diagoutd, diagoutvd = diag.compute_rh(var0,var1,var2,dnamesvar,dvnamesvar)
1137
1138        ncvar.insert_variable(ncobj, 'hurs', diagout, diagoutd, diagoutvd, newnc)
1139
1140# slw: total soil liquid water SH2O, DZS
1141    elif diagn == 'WRFslw':
1142           
1143        var0 = ncobj.variables[depvars[0]][:]
1144        var10 = ncobj.variables[depvars[1]][:]
1145        dnamesvar = list(ncobj.variables[depvars[0]].dimensions)
1146        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1147
1148        var1 = var0.copy()*0.
1149        var2 = var0.copy()*0.+1.
1150        # Must be a better way....
1151        for j in range(var0.shape[2]):
1152          for i in range(var0.shape[3]):
1153              var1[:,:,j,i] = var10
1154
1155        diagout, diagoutd, diagoutvd = diag.Forcompute_zint(var0, var1, var2,        \
1156          dnamesvar, dvnamesvar)
1157
1158        # Removing the nonChecking variable-dimensions from the initial list
1159        varsadd = []
1160        diagoutvd = list(dvnames)
1161        for nonvd in NONchkvardims:
1162            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1163            varsadd.append(nonvd)
1164        ncvar.insert_variable(ncobj, 'slw', diagout, diagoutd, diagoutvd, newnc)
1165
1166# td (psfc, t, q) from TimeSeries files
1167    elif diagn == 'TStd' or diagn == 'td':
1168           
1169        var0 = ncobj.variables[depvars[0]][:]
1170        var1 = ncobj.variables[depvars[1]][:] - 273.15
1171        var2 = ncobj.variables[depvars[2]][:]
1172
1173        dnamesvar = ncobj.variables[depvars[0]].dimensions
1174        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1175
1176        diagout, diagoutd, diagoutvd = diag.compute_td(var0,var1,var2,dnamesvar,dvnamesvar)
1177
1178        ncvar.insert_variable(ncobj, 'tdas', diagout, diagoutd, diagoutvd, newnc)
1179
1180# td (psfc, t, q) from TimeSeries files
1181    elif diagn == 'TStdC' or diagn == 'tdC':
1182           
1183        var0 = ncobj.variables[depvars[0]][:]
1184# Temperature is already in degrees Celsius
1185        var1 = ncobj.variables[depvars[1]][:]
1186        var2 = ncobj.variables[depvars[2]][:]
1187
1188        dnamesvar = ncobj.variables[depvars[0]].dimensions
1189        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1190
1191        diagout, diagoutd, diagoutvd = diag.compute_td(var0,var1,var2,dnamesvar,dvnamesvar)
1192
1193        # Removing the nonChecking variable-dimensions from the initial list
1194        varsadd = []
1195        for nonvd in NONchkvardims:
1196            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1197            varsadd.append(nonvd)
1198
1199        ncvar.insert_variable(ncobj, 'tdas', diagout, diagoutd, diagoutvd, newnc)
1200
1201# wds (u, v)
1202    elif diagn == 'TSwds' or diagn == 'wds' :
1203 
1204        var0 = ncobj.variables[depvars[0]][:]
1205        var1 = ncobj.variables[depvars[1]][:]
1206
1207        dnamesvar = ncobj.variables[depvars[0]].dimensions
1208        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1209
1210        diagout, diagoutd, diagoutvd = diag.compute_wds(var0,var1,dnamesvar,dvnamesvar)
1211
1212        # Removing the nonChecking variable-dimensions from the initial list
1213        varsadd = []
1214        for nonvd in NONchkvardims:
1215            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1216            varsadd.append(nonvd)
1217
1218        ncvar.insert_variable(ncobj, 'wds', diagout, diagoutd, diagoutvd, newnc)
1219
1220# wss (u, v)
1221    elif diagn == 'TSwss' or diagn == 'wss':
1222           
1223        var0 = ncobj.variables[depvars[0]][:]
1224        var1 = ncobj.variables[depvars[1]][:]
1225
1226        dnamesvar = ncobj.variables[depvars[0]].dimensions
1227        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1228
1229        diagout, diagoutd, diagoutvd = diag.compute_wss(var0,var1,dnamesvar,dvnamesvar)
1230
1231        # Removing the nonChecking variable-dimensions from the initial list
1232        varsadd = []
1233        for nonvd in NONchkvardims:
1234            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1235            varsadd.append(nonvd)
1236
1237        ncvar.insert_variable(ncobj, 'wss', diagout, diagoutd, diagoutvd, newnc)
1238
1239# turbulence (var)
1240    elif diagn == 'turbulence':
1241
1242        var0 = ncobj.variables[depvars][:]
1243
1244        dnamesvar = list(ncobj.variables[depvars].dimensions)
1245        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1246
1247        diagout, diagoutd, diagoutvd = diag.compute_turbulence(var0,dnamesvar,dvnamesvar)
1248        valsvar = gen.variables_values(depvars)
1249
1250        newvarn = depvars + 'turb'
1251        ncvar.insert_variable(ncobj, newvarn, diagout, diagoutd, 
1252          diagoutvd, newnc)
1253        varobj = newnc.variables[newvarn]
1254        attrv = varobj.long_name
1255        attr = varobj.delncattr('long_name')
1256        newattr = ncvar.set_attribute(varobj, 'long_name', attrv +                   \
1257          " Taylor decomposition turbulence term")
1258
1259# ua va from ws wd (deg)
1260    elif diagn == 'uavaFROMwswd':
1261           
1262        var0 = ncobj.variables[depvars[0]][:]
1263        var1 = ncobj.variables[depvars[1]][:]
1264
1265        ua = var0*np.cos(var1*np.pi/180.)
1266        va = var0*np.sin(var1*np.pi/180.)
1267
1268        dnamesvar = ncobj.variables[depvars[0]].dimensions
1269        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1270
1271        ncvar.insert_variable(ncobj, 'ua', ua, dnamesvar, dvnamesvar, newnc)
1272        ncvar.insert_variable(ncobj, 'va', va, dnamesvar, dvnamesvar, newnc)
1273
1274# ua va from obs ws wd (deg)
1275    elif diagn == 'uavaFROMobswswd':
1276           
1277        var0 = ncobj.variables[depvars[0]][:]
1278        var1 = ncobj.variables[depvars[1]][:]
1279
1280        ua = var0*np.cos((var1+180.)*np.pi/180.)
1281        va = var0*np.sin((var1+180.)*np.pi/180.)
1282
1283        dnamesvar = ncobj.variables[depvars[0]].dimensions
1284        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1285
1286        ncvar.insert_variable(ncobj, 'ua', ua, dnamesvar, dvnamesvar, newnc)
1287        ncvar.insert_variable(ncobj, 'va', va, dnamesvar, dvnamesvar, newnc)
1288
1289# WRFbils fom WRF as HFX + LH
1290    elif diagn == 'WRFbils':
1291           
1292        var0 = ncobj.variables[depvars[0]][:]
1293        var1 = ncobj.variables[depvars[1]][:]
1294
1295        diagout = var0 + var1
1296        dnamesvar = list(ncobj.variables[depvars[0]].dimensions)
1297        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1298
1299        ncvar.insert_variable(ncobj, 'bils', diagout, dnamesvar, dvnamesvar, newnc)
1300
1301# WRFcape_afwa CAPE, CIN, ZLFC, PLFC, LI following WRF 'phys/module_diaf_afwa.F'
1302#   methodology as WRFt, WRFrh, WRFp, WRFgeop, HGT
1303    elif diagn == 'WRFcape_afwa':
1304        var0 = WRFt
1305        var1 = WRFrh
1306        var2 = WRFp
1307        dz = WRFgeop.shape[1]
1308        # de-staggering
1309        var3 = 0.5*(WRFgeop[:,0:dz-1,:,:]+WRFgeop[:,1:dz,:,:])/9.8
1310        var4 = ncobj.variables[depvars[4]][0,:,:]
1311
1312        dnamesvar = list(ncobj.variables['T'].dimensions)
1313        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1314
1315        diagout = np.zeros(var0.shape, dtype=np.float)
1316        diagout1, diagout2, diagout3, diagout4, diagout5, diagoutd, diagoutvd =      \
1317          diag.Forcompute_cape_afwa(var0, var1, var2, var3, var4, 3, dnamesvar,      \
1318          dvnamesvar)
1319
1320        # Removing the nonChecking variable-dimensions from the initial list
1321        varsadd = []
1322        for nonvd in NONchkvardims:
1323            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1324            varsadd.append(nonvd)
1325
1326        ncvar.insert_variable(ncobj, 'cape', diagout1, diagoutd, diagoutvd, newnc)
1327        ncvar.insert_variable(ncobj, 'cin', diagout2, diagoutd, diagoutvd, newnc)
1328        ncvar.insert_variable(ncobj, 'zlfc', diagout3, diagoutd, diagoutvd, newnc)
1329        ncvar.insert_variable(ncobj, 'plfc', diagout4, diagoutd, diagoutvd, newnc)
1330        ncvar.insert_variable(ncobj, 'li', diagout5, diagoutd, diagoutvd, newnc)
1331
1332# WRFclivi WRF water vapour path WRFdens, QICE, QGRAUPEL, QHAIL
1333    elif diagn == 'WRFclivi':
1334           
1335        var0 = WRFdens
1336        qtot = ncobj.variables[depvars[1]]
1337        qtotv = qtot[:]
1338        Nspecies = len(depvars) - 2
1339        for iv in range(Nspecies):
1340            if ncobj.variables.has_key(depvars[iv+2]):
1341                var1 = ncobj.variables[depvars[iv+2]][:]
1342                qtotv = qtotv + var1
1343
1344        dnamesvar = list(qtot.dimensions)
1345        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1346
1347        diagout, diagoutd, diagoutvd = diag.compute_clivi(var0, qtotv, dnamesvar,dvnamesvar)
1348
1349        # Removing the nonChecking variable-dimensions from the initial list
1350        varsadd = []
1351        diagoutvd = list(dvnames)
1352        for nonvd in NONchkvardims:
1353            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1354            varsadd.append(nonvd)
1355        ncvar.insert_variable(ncobj, 'clivi', diagout, diagoutd, diagoutvd, newnc)
1356
1357# WRFclwvi WRF water cloud-condensed path WRFdens, QCLOUD, QICE, QGRAUPEL, QHAIL
1358    elif diagn == 'WRFclwvi':
1359           
1360        var0 = WRFdens
1361        qtot = ncobj.variables[depvars[1]]
1362        qtotv = ncobj.variables[depvars[1]]
1363        Nspecies = len(depvars) - 2
1364        for iv in range(Nspecies):
1365            if ncobj.variables.has_key(depvars[iv+2]):
1366                var1 = ncobj.variables[depvars[iv+2]]
1367                qtotv = qtotv + var1[:]
1368
1369        dnamesvar = list(qtot.dimensions)
1370        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1371
1372        diagout, diagoutd, diagoutvd = diag.compute_clwvl(var0, qtotv, dnamesvar,dvnamesvar)
1373
1374        # Removing the nonChecking variable-dimensions from the initial list
1375        varsadd = []
1376        diagoutvd = list(dvnames)
1377        for nonvd in NONchkvardims:
1378            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1379            varsadd.append(nonvd)
1380        ncvar.insert_variable(ncobj, 'clwvi', diagout, diagoutd, diagoutvd, newnc)
1381
1382# WRF_denszint WRF vertical integration as WRFdens, sum(Q[water species1], ..., Q[water speciesN]), varn=[varN]
1383    elif diagn == 'WRF_denszint':
1384           
1385        var0 = WRFdens
1386        varn = depvars[1].split('=')[1]
1387        qtot = ncobj.variables[depvars[2]]
1388        qtotv = ncobj.variables[depvars[2]]
1389        Nspecies = len(depvars) - 2
1390        for iv in range(Nspecies):
1391            if ncobj.variables.has_key(depvars[iv+2]):
1392                var1 = ncobj.variables[depvars[iv+2]]
1393                qtotv = qtotv + var1[:]
1394
1395        dnamesvar = list(qtot.dimensions)
1396        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1397
1398        diagout, diagoutd, diagoutvd = diag.compute_clwvl(var0, qtotv, dnamesvar,dvnamesvar)
1399
1400        # Removing the nonChecking variable-dimensions from the initial list
1401        varsadd = []
1402        diagoutvd = list(dvnames)
1403        for nonvd in NONchkvardims:
1404            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1405            varsadd.append(nonvd)
1406        ncvar.insert_variable(ncobj, varn, diagout, diagoutd, diagoutvd, newnc)
1407
1408# WRFgeop geopotential from WRF as PH + PHB
1409    elif diagn == 'WRFgeop':
1410        var0 = ncobj.variables[depvars[0]][:]
1411        var1 = ncobj.variables[depvars[1]][:]
1412
1413        # de-staggering geopotential
1414        diagout0 = var0 + var1
1415        dt = diagout0.shape[0]
1416        dz = diagout0.shape[1]
1417        dy = diagout0.shape[2]
1418        dx = diagout0.shape[3]
1419
1420        diagout = np.zeros((dt,dz-1,dy,dx), dtype=np.float)
1421        diagout = 0.5*(diagout0[:,1:dz,:,:]+diagout0[:,0:dz-1,:,:])
1422
1423        # Removing the nonChecking variable-dimensions from the initial list
1424        varsadd = []
1425        diagoutvd = list(dvnames)
1426        for nonvd in NONchkvardims:
1427            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1428            varsadd.append(nonvd)
1429
1430        ncvar.insert_variable(ncobj, 'zg', diagout, dnames, diagoutvd, newnc)
1431
1432# WRFpotevap_orPM potential evapotranspiration following Penman-Monteith formulation
1433#   implemented in ORCHIDEE (in src_sechiba/enerbil.f90) as: WRFdens, UST, U10, V10, T2, PSFC, QVAPOR
1434    elif diagn == 'WRFpotevap_orPM':
1435        var0 = WRFdens[:,0,:,:]
1436        var1 = ncobj.variables[depvars[1]][:]
1437        var2 = ncobj.variables[depvars[2]][:]
1438        var3 = ncobj.variables[depvars[3]][:]
1439        var4 = ncobj.variables[depvars[4]][:]
1440        var5 = ncobj.variables[depvars[5]][:]
1441        var6 = ncobj.variables[depvars[6]][:,0,:,:]
1442
1443        dnamesvar = list(ncobj.variables[depvars[1]].dimensions)
1444        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1445
1446        diagout = np.zeros(var1.shape, dtype=np.float)
1447        diagout, diagoutd, diagoutvd = diag.Forcompute_potevap_orPM(var0, var1, var2,\
1448          var3, var4, var5, var6, dnamesvar, dvnamesvar)
1449
1450        # Removing the nonChecking variable-dimensions from the initial list
1451        varsadd = []
1452        for nonvd in NONchkvardims:
1453            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1454            varsadd.append(nonvd)
1455
1456        ncvar.insert_variable(ncobj, 'evspsblpot', diagout, diagoutd, diagoutvd, newnc)
1457
1458# WRFmslp_ptarget sea-level pressure following ECMWF method as PSFC, HGT, WRFt, WRFp, ZNU, ZNW
1459    elif diagn == 'WRFpsl_ecmwf':
1460        var0 = ncobj.variables[depvars[0]][:]
1461        var1 = ncobj.variables[depvars[1]][0,:,:]
1462        var2 = WRFt[:,0,:,:]
1463        var4 = WRFp[:,0,:,:]
1464        var5 = ncobj.variables[depvars[4]][0,:]
1465        var6 = ncobj.variables[depvars[5]][0,:]
1466
1467        # This is quite too appriximate!! passing pressure at half-levels to 2nd full
1468        #   level, using eta values at full (ZNW) and half (ZNU) mass levels
1469        var3 = WRFp[:,0,:,:] + (var6[1] - var5[0])*(WRFp[:,1,:,:] - WRFp[:,0,:,:])/  \
1470          (var5[1]-var5[0])
1471
1472        dnamesvar = list(ncobj.variables[depvars[0]].dimensions)
1473        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1474
1475        diagout = np.zeros(var0.shape, dtype=np.float)
1476        diagout, diagoutd, diagoutvd = diag.Forcompute_psl_ecmwf(var0, var1, var2,   \
1477          var3, var4, dnamesvar, dvnamesvar)
1478
1479        # Removing the nonChecking variable-dimensions from the initial list
1480        varsadd = []
1481        for nonvd in NONchkvardims:
1482            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1483            varsadd.append(nonvd)
1484
1485        ncvar.insert_variable(ncobj, 'psl', diagout, diagoutd, diagoutvd, newnc)
1486
1487# WRFmslp_ptarget sea-level pressure following ptarget method as WRFp, PSFC, WRFt, HGT, QVAPOR
1488    elif diagn == 'WRFpsl_ptarget':
1489        var0 = WRFp
1490        var1 = ncobj.variables[depvars[1]][:]
1491        var2 = WRFt
1492        var3 = ncobj.variables[depvars[3]][0,:,:]
1493        var4 = ncobj.variables[depvars[4]][:]
1494
1495        dnamesvar = list(ncobj.variables[depvars[4]].dimensions)
1496        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1497
1498        diagout = np.zeros(var0.shape, dtype=np.float)
1499        diagout, diagoutd, diagoutvd = diag.Forcompute_psl_ptarget(var0, var1, var2, \
1500          var3, var4, 700000., dnamesvar, dvnamesvar)
1501
1502        # Removing the nonChecking variable-dimensions from the initial list
1503        varsadd = []
1504        for nonvd in NONchkvardims:
1505            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1506            varsadd.append(nonvd)
1507
1508        ncvar.insert_variable(ncobj, 'psl', diagout, diagoutd, diagoutvd, newnc)
1509
1510# WRFp pressure from WRF as P + PB
1511    elif diagn == 'WRFp':
1512        var0 = ncobj.variables[depvars[0]][:]
1513        var1 = ncobj.variables[depvars[1]][:]
1514           
1515        diagout = var0 + var1
1516        diagoutd = list(ncobj.variables[depvars[0]].dimensions)
1517        diagoutvd = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1518
1519        # Removing the nonChecking variable-dimensions from the initial list
1520        varsadd = []
1521        for nonvd in NONchkvardims:
1522            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1523            varsadd.append(nonvd)
1524
1525        ncvar.insert_variable(ncobj, 'pres', diagout, diagoutd, diagoutvd, newnc)
1526
1527# WRFpos
1528    elif diagn == 'WRFpos':
1529           
1530        dnamesvar = ncobj.variables['MAPFAC_M'].dimensions
1531        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1532
1533        ncvar.insert_variable(ncobj, 'WRFpos', WRFpos, dnamesvar, dvnamesvar, newnc)
1534
1535# WRFprw WRF water vapour path WRFdens, QVAPOR
1536    elif diagn == 'WRFprw':
1537           
1538        var0 = WRFdens
1539        var1 = ncobj.variables[depvars[1]]
1540
1541        dnamesvar = list(var1.dimensions)
1542        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1543
1544        diagout, diagoutd, diagoutvd = diag.compute_prw(var0, var1, dnamesvar,dvnamesvar)
1545
1546        # Removing the nonChecking variable-dimensions from the initial list
1547        varsadd = []
1548        diagoutvd = list(dvnames)
1549        for nonvd in NONchkvardims:
1550            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1551            varsadd.append(nonvd)
1552        ncvar.insert_variable(ncobj, 'prw', diagout, diagoutd, diagoutvd, newnc)
1553
1554# WRFrh (P, T, QVAPOR)
1555    elif diagn == 'WRFrh':
1556           
1557        dnamesvar = list(ncobj.variables[depvars[2]].dimensions)
1558        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1559
1560        ncvar.insert_variable(ncobj, 'hur', WRFrh, dnames, dvnames, newnc)
1561
1562# WRFrhs (PSFC, T2, Q2)
1563    elif diagn == 'WRFrhs':
1564           
1565        var0 = ncobj.variables[depvars[0]][:]
1566        var1 = ncobj.variables[depvars[1]][:]
1567        var2 = ncobj.variables[depvars[2]][:]
1568
1569        dnamesvar = list(ncobj.variables[depvars[2]].dimensions)
1570        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1571
1572        diagout, diagoutd, diagoutvd = diag.compute_rh(var0,var1,var2,dnamesvar,dvnamesvar)
1573        ncvar.insert_variable(ncobj, 'hurs', diagout, diagoutd, diagoutvd, newnc)
1574
1575# rvors (u10, v10, WRFpos)
1576    elif diagn == 'WRFrvors':
1577           
1578        var0 = ncobj.variables[depvars[0]]
1579        var1 = ncobj.variables[depvars[1]]
1580
1581        diagout = rotational_z(var0, var1, distx)
1582
1583        dnamesvar = ncobj.variables[depvars[0]].dimensions
1584        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1585
1586        ncvar.insert_variable(ncobj, 'rvors', diagout, dnamesvar, dvnamesvar, newnc)
1587
1588# WRFt (T, P, PB)
1589    elif diagn == 'WRFt':
1590        var0 = ncobj.variables[depvars[0]][:]
1591        var1 = ncobj.variables[depvars[1]][:]
1592        var2 = ncobj.variables[depvars[2]][:]
1593
1594        p0=100000.
1595        p=var1 + var2
1596
1597        WRFt = (var0 + 300.)*(p/p0)**(2./7.)
1598
1599        dnamesvar = list(ncobj.variables[depvars[0]].dimensions)
1600        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1601
1602        # Removing the nonChecking variable-dimensions from the initial list
1603        varsadd = []
1604        diagoutvd = list(dvnames)
1605        for nonvd in NONchkvardims:
1606            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1607            varsadd.append(nonvd)
1608
1609        ncvar.insert_variable(ncobj, 'ta', WRFt, dnames, diagoutvd, newnc)
1610
1611# WRFtda (WRFrh, WRFt)
1612    elif diagn == 'WRFtda':
1613        ARM2 = fdef.module_definitions.arm2
1614        ARM3 = fdef.module_definitions.arm3
1615
1616        gammatarh = np.log(WRFrh) + ARM2*(WRFt-273.15)/((WRFt-273.15)+ARM3)
1617        td = ARM3*gammatarh/(ARM2-gammatarh)
1618
1619        dnamesvar = list(ncobj.variables['T'].dimensions)
1620        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1621
1622        # Removing the nonChecking variable-dimensions from the initial list
1623        varsadd = []
1624        diagoutvd = list(dvnames)
1625        for nonvd in NONchkvardims:
1626            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1627            varsadd.append(nonvd)
1628
1629        ncvar.insert_variable(ncobj, 'tda', td, dnames, diagoutvd, newnc)
1630
1631# WRFtdas (PSFC, T2, Q2)
1632    elif diagn == 'WRFtdas':
1633        ARM2 = fdef.module_definitions.arm2
1634        ARM3 = fdef.module_definitions.arm3
1635
1636        var0 = ncobj.variables[depvars[0]][:]
1637        var1 = ncobj.variables[depvars[1]][:]
1638        var2 = ncobj.variables[depvars[2]][:]
1639
1640        dnamesvar = list(ncobj.variables[depvars[1]].dimensions)
1641        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1642
1643        rhs, diagoutd, diagoutvd = diag.compute_rh(var0,var1,var2,dnamesvar,dvnamesvar)
1644
1645        gammatarhs = np.log(rhs) + ARM2*(var1-273.15)/((var1-273.15)+ARM3)
1646        tdas = ARM3*gammatarhs/(ARM2-gammatarhs) + 273.15
1647
1648        # Removing the nonChecking variable-dimensions from the initial list
1649        varsadd = []
1650        diagoutvd = list(dvnames)
1651        for nonvd in NONchkvardims:
1652            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1653            varsadd.append(nonvd)
1654
1655        ncvar.insert_variable(ncobj, 'tdas', tdas, dnames, diagoutvd, newnc)
1656
1657# WRFua (U, V, SINALPHA, COSALPHA) to be rotated !!
1658    elif diagn == 'WRFua':
1659        var0 = ncobj.variables[depvars[0]][:]
1660        var1 = ncobj.variables[depvars[1]][:]
1661        var2 = ncobj.variables[depvars[2]][:]
1662        var3 = ncobj.variables[depvars[3]][:]
1663
1664        # un-staggering variables
1665        if len(var0.shape) == 4:
1666            unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1667        elif len(var0.shape) == 3:
1668            # Asuming sunding point (dimt, dimz, dimstgx)
1669            unstgdims = [var0.shape[0], var0.shape[1]]
1670
1671        ua = np.zeros(tuple(unstgdims), dtype=np.float)
1672        unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1673        unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1674
1675        if len(var0.shape) == 4:
1676            unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1677            unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1678
1679            for iz in range(var0.shape[1]):
1680                ua[:,iz,:,:] = unstgvar0[:,iz,:,:]*var3 - unstgvar1[:,iz,:,:]*var2
1681
1682            dnamesvar = ['Time','bottom_top','south_north','west_east']
1683
1684        elif len(var0.shape) == 3:
1685            unstgvar0 = 0.5*(var0[:,:,0] + var0[:,:,1])
1686            unstgvar1 = 0.5*(var1[:,:,0] + var1[:,:,1])
1687            for iz in range(var0.shape[1]):
1688                ua[:,iz] = unstgvar0[:,iz]*var3 - unstgvar1[:,iz]*var2
1689
1690            dnamesvar = ['Time','bottom_top']
1691
1692        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1693
1694        # Removing the nonChecking variable-dimensions from the initial list
1695        varsadd = []
1696        diagoutvd = list(dvnames)
1697        for nonvd in NONchkvardims:
1698            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1699            varsadd.append(nonvd)
1700
1701        ncvar.insert_variable(ncobj, 'ua', ua, dnames, diagoutvd, newnc)
1702
1703# WRFua (U, V, SINALPHA, COSALPHA) to be rotated !!
1704    elif diagn == 'WRFva':
1705        var0 = ncobj.variables[depvars[0]][:]
1706        var1 = ncobj.variables[depvars[1]][:]
1707        var2 = ncobj.variables[depvars[2]][:]
1708        var3 = ncobj.variables[depvars[3]][:]
1709
1710        # un-staggering variables
1711        if len(var0.shape) == 4:
1712            unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1713        elif len(var0.shape) == 3:
1714            # Asuming sunding point (dimt, dimz, dimstgx)
1715            unstgdims = [var0.shape[0], var0.shape[1]]
1716
1717        va = np.zeros(tuple(unstgdims), dtype=np.float)
1718        unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1719        unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1720
1721        if len(var0.shape) == 4:
1722            unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1723            unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1724
1725            for iz in range(var0.shape[1]):
1726                va[:,iz,:,:] = unstgvar0[:,iz,:,:]*var2 + unstgvar1[:,iz,:,:]*var3
1727
1728            dnamesvar = ['Time','bottom_top','south_north','west_east']
1729
1730        elif len(var0.shape) == 3:
1731            unstgvar0 = 0.5*(var0[:,:,0] + var0[:,:,1])
1732            unstgvar1 = 0.5*(var1[:,:,0] + var1[:,:,1])
1733            for iz in range(var0.shape[1]):
1734                va[:,iz] = unstgvar0[:,iz]*var2 + unstgvar1[:,iz]*var3
1735
1736            dnamesvar = ['Time','bottom_top']
1737
1738        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1739
1740        # Removing the nonChecking variable-dimensions from the initial list
1741        varsadd = []
1742        diagoutvd = list(dvnames)
1743        for nonvd in NONchkvardims:
1744            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1745            varsadd.append(nonvd)
1746        ncvar.insert_variable(ncobj, 'va', va, dnames, diagoutvd, newnc)
1747
1748
1749# WRFwd (U, V, SINALPHA, COSALPHA) to be rotated !!
1750    elif diagn == 'WRFwd':
1751        var0 = ncobj.variables[depvars[0]][:]
1752        var1 = ncobj.variables[depvars[1]][:]
1753        var2 = ncobj.variables[depvars[2]][:]
1754        var3 = ncobj.variables[depvars[3]][:]
1755
1756        # un-staggering variables
1757        if len(var0.shape) == 4:
1758            unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1759        elif len(var0.shape) == 3:
1760            # Asuming sunding point (dimt, dimz, dimstgx)
1761            unstgdims = [var0.shape[0], var0.shape[1]]
1762
1763        ua = np.zeros(tuple(unstgdims), dtype=np.float)
1764        va = np.zeros(tuple(unstgdims), dtype=np.float)
1765        unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1766        unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1767
1768        if len(var0.shape) == 4:
1769            unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1770            unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1771
1772            for iz in range(var0.shape[1]):
1773                ua[:,iz,:,:] = unstgvar0[:,iz,:,:]*var3 - unstgvar1[:,iz,:,:]*var2
1774                va[:,iz,:,:] = unstgvar0[:,iz,:,:]*var2 + unstgvar1[:,iz,:,:]*var3
1775
1776            dnamesvar = ['Time','bottom_top','south_north','west_east']
1777
1778        elif len(var0.shape) == 3:
1779            unstgvar0 = 0.5*(var0[:,:,0] + var0[:,:,1])
1780            unstgvar1 = 0.5*(var1[:,:,0] + var1[:,:,1])
1781            for iz in range(var0.shape[1]):
1782                ua[:,iz] = unstgvar0[:,iz]*var3 - unstgvar1[:,iz]*var2
1783                va[:,iz] = unstgvar0[:,iz]*var2 + unstgvar1[:,iz]*var3
1784
1785            dnamesvar = ['Time','bottom_top']
1786
1787        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1788
1789        wd, dnames, dvnames = diag.compute_wd(ua, va, dnamesvar, dvnamesvar)
1790
1791        # Removing the nonChecking variable-dimensions from the initial list
1792        varsadd = []
1793        diagoutvd = list(dvnames)
1794        for nonvd in NONchkvardims:
1795            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1796            varsadd.append(nonvd)
1797
1798        ncvar.insert_variable(ncobj, 'wd', wd, dnames, diagoutvd, newnc)
1799
1800# WRFtime
1801    elif diagn == 'WRFtime':
1802           
1803        diagout = WRFtime
1804
1805        dnamesvar = ['Time']
1806        dvnamesvar = ['Times']
1807
1808        ncvar.insert_variable(ncobj, 'time', diagout, dnamesvar, dvnamesvar, newnc)
1809
1810# ws (U, V)
1811    elif diagn == 'ws':
1812           
1813        # un-staggering variables
1814        if len(var0.shape) == 4:
1815            unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1816        elif len(var0.shape) == 3:
1817            # Asuming sunding point (dimt, dimz, dimstgx)
1818            unstgdims = [var0.shape[0], var0.shape[1]]
1819
1820        ua = np.zeros(tuple(unstgdims), dtype=np.float)
1821        va = np.zeros(tuple(unstgdims), dtype=np.float)
1822        unstgvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1823        unstgvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1824
1825        if len(var0.shape) == 4:
1826            unstgvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1827            unstgvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1828
1829            for iz in range(var0.shape[1]):
1830                ua[:,iz,:,:] = unstgvar0[:,iz,:,:]*var3 - unstgvar1[:,iz,:,:]*var2
1831                va[:,iz,:,:] = unstgvar0[:,iz,:,:]*var2 + unstgvar1[:,iz,:,:]*var3
1832
1833            dnamesvar = ['Time','bottom_top','south_north','west_east']
1834
1835        elif len(var0.shape) == 3:
1836            unstgvar0 = 0.5*(var0[:,:,0] + var0[:,:,1])
1837            unstgvar1 = 0.5*(var1[:,:,0] + var1[:,:,1])
1838            for iz in range(var0.shape[1]):
1839                ua[:,iz] = unstgvar0[:,iz]*var3 - unstgvar1[:,iz]*var2
1840                va[:,iz] = unstgvar0[:,iz]*var2 + unstgvar1[:,iz]*var3
1841
1842            dnamesvar = ['Time','bottom_top']
1843
1844        diagout = np.sqrt(unstgvar0*unstgvar0 + unstgvar1*unstgvar1)
1845
1846#        dnamesvar = ncobj.variables[depvars[0]].dimensions
1847        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1848
1849        # Removing the nonChecking variable-dimensions from the initial list
1850        varsadd = []
1851        diagoutvd = list(dvnamesvar)
1852        for nonvd in NONchkvardims:
1853            if gen.searchInlist(dvnamesvar,nonvd): diagoutvd.remove(nonvd)
1854            varsadd.append(nonvd)
1855        ncvar.insert_variable(ncobj, 'ws', diagout, dnamesvar, diagoutvd, newnc)
1856
1857# wss (u10, v10)
1858    elif diagn == 'wss':
1859           
1860        var0 = ncobj.variables[depvars[0]][:]
1861        var1 = ncobj.variables[depvars[1]][:]
1862
1863        diagout = np.sqrt(var0*var0 + var1*var1)
1864
1865        dnamesvar = ncobj.variables[depvars[0]].dimensions
1866        dvnamesvar = ncvar.var_dim_dimv(dnamesvar,dnames,dvnames)
1867
1868        ncvar.insert_variable(ncobj, 'wss', diagout, dnamesvar, dvnamesvar, newnc)
1869
1870# WRFheight height from WRF geopotential as WRFGeop/g
1871    elif diagn == 'WRFheight':
1872           
1873        diagout = WRFgeop/grav
1874
1875        # Removing the nonChecking variable-dimensions from the initial list
1876        varsadd = []
1877        diagoutvd = list(dvnames)
1878        for nonvd in NONchkvardims:
1879            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1880            varsadd.append(nonvd)
1881
1882        ncvar.insert_variable(ncobj, 'zhgt', diagout, dnames, diagoutvd, newnc)
1883
1884# WRFheightrel relative-height from WRF geopotential as WRFgeop(PH + PHB)/g-HGT 'WRFheightrel|PH@PHB@HGT
1885    elif diagn == 'WRFheightrel':
1886        var0 = ncobj.variables[depvars[0]][:]
1887        var1 = ncobj.variables[depvars[1]][:]
1888        var2 = ncobj.variables[depvars[2]][:]
1889
1890        dimz = var0.shape[1]
1891        diagout = np.zeros(tuple(var0.shape), dtype=np.float)
1892        for iz in range(dimz):
1893            diagout[:,iz,:,:] = (var0[:,iz,:,:]+ var1[:,iz,:,:])/grav - var2
1894
1895        # Removing the nonChecking variable-dimensions from the initial list
1896        varsadd = []
1897        diagoutvd = list(dvnames)
1898        for nonvd in NONchkvardims:
1899            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1900            varsadd.append(nonvd)
1901
1902        ncvar.insert_variable(ncobj, 'zhgtrel', diagout, dnames, diagoutvd, newnc)
1903
1904# WRFzmla_gen generic boundary layer hieght computation from WRF theta, QVAPOR, WRFgeop, HGT,
1905    elif diagn == 'WRFzmlagen':
1906        var0 = ncobj.variables[depvars[0]][:]+300.
1907        var1 = ncobj.variables[depvars[1]][:]
1908        dimz = var0.shape[1]
1909        var2 = WRFgeop[:,1:dimz+1,:,:]/9.8
1910        var3 = ncobj.variables[depvars[3]][0,:,:]
1911
1912        diagout, diagoutd, diagoutvd = diag.Forcompute_zmla_gen(var0,var1,var2,var3, \
1913          dnames,dvnames)
1914
1915        # Removing the nonChecking variable-dimensions from the initial list
1916        varsadd = []
1917        for nonvd in NONchkvardims:
1918            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1919            varsadd.append(nonvd)
1920
1921        ncvar.insert_variable(ncobj, 'zmla', diagout, diagoutd, diagoutvd, newnc)
1922
1923# WRFzwind wind extrapolation at a given height using power law computation from WRF
1924#   U, V, WRFz, U10, V10, SINALPHA, COSALPHA, z=[zval]
1925    elif diagn == 'WRFzwind':
1926        var0 = ncobj.variables[depvars[0]][:]
1927        var1 = ncobj.variables[depvars[1]][:]
1928        var2 = WRFz
1929        var3 = ncobj.variables[depvars[3]][:]
1930        var4 = ncobj.variables[depvars[4]][:]
1931        var5 = ncobj.variables[depvars[5]][0,:,:]
1932        var6 = ncobj.variables[depvars[6]][0,:,:]
1933        var7 = np.float(depvars[7].split('=')[1])
1934
1935        # un-staggering 3D winds
1936        unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1937        va = np.zeros(tuple(unstgdims), dtype=np.float)
1938        unvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1939        unvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1940        unvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1941        unvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1942
1943        diagout1, diagout2, diagoutd, diagoutvd = diag.Forcompute_zwind(unvar0,      \
1944          unvar1, var2, var3, var4, var5, var6, var7, dnames, dvnames)
1945
1946        # Removing the nonChecking variable-dimensions from the initial list
1947        varsadd = []
1948        for nonvd in NONchkvardims:
1949            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1950            varsadd.append(nonvd)
1951
1952        ncvar.insert_variable(ncobj, 'uaz', diagout1, diagoutd, diagoutvd, newnc)
1953        ncvar.insert_variable(ncobj, 'vaz', diagout2, diagoutd, diagoutvd, newnc)
1954
1955# WRFzwind wind extrapolation at a given hieght using logarithmic law computation
1956#   from WRF U, V, WRFz, U10, V10, SINALPHA, COSALPHA, z=[zval]
1957    elif diagn == 'WRFzwind_log':
1958        var0 = ncobj.variables[depvars[0]][:]
1959        var1 = ncobj.variables[depvars[1]][:]
1960        var2 = WRFz
1961        var3 = ncobj.variables[depvars[3]][:]
1962        var4 = ncobj.variables[depvars[4]][:]
1963        var5 = ncobj.variables[depvars[5]][0,:,:]
1964        var6 = ncobj.variables[depvars[6]][0,:,:]
1965        var7 = np.float(depvars[7].split('=')[1])
1966
1967        # un-staggering 3D winds
1968        unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1]
1969        va = np.zeros(tuple(unstgdims), dtype=np.float)
1970        unvar0 = np.zeros(tuple(unstgdims), dtype=np.float)
1971        unvar1 = np.zeros(tuple(unstgdims), dtype=np.float)
1972        unvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]])
1973        unvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:])
1974
1975        diagout1, diagout2, diagoutd, diagoutvd = diag.Forcompute_zwind_log(unvar0,  \
1976          unvar1, var2, var3, var4, var5, var6, var7, dnames, dvnames)
1977
1978        # Removing the nonChecking variable-dimensions from the initial list
1979        varsadd = []
1980        for nonvd in NONchkvardims:
1981            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
1982            varsadd.append(nonvd)
1983
1984        ncvar.insert_variable(ncobj, 'uaz', diagout1, diagoutd, diagoutvd, newnc)
1985        ncvar.insert_variable(ncobj, 'vaz', diagout2, diagoutd, diagoutvd, newnc)
1986
1987# WRFzwindMO wind extrapolation at a given height computation using Monin-Obukhow
1988#   theory from WRF UST, ZNT, RMOL, U10, V10, SINALPHA, COSALPHA, z=[zval]
1989#   NOTE: only useful for [zval] < 80. m
1990    elif diagn == 'WRFzwindMO':
1991        var0 = ncobj.variables[depvars[0]][:]
1992        var1 = ncobj.variables[depvars[1]][:]
1993        var2 = ncobj.variables[depvars[2]][:]
1994        var3 = ncobj.variables[depvars[3]][:]
1995        var4 = ncobj.variables[depvars[4]][:]
1996        var5 = ncobj.variables[depvars[5]][0,:,:]
1997        var6 = ncobj.variables[depvars[6]][0,:,:]
1998        var7 = np.float(depvars[7].split('=')[1])
1999
2000        diagout1, diagout2, diagoutd, diagoutvd = diag.Forcompute_zwindMO(var0, var1,\
2001          var2, var3, var4, var5, var6, var7, dnames, dvnames)
2002
2003        # Removing the nonChecking variable-dimensions from the initial list
2004        varsadd = []
2005        for nonvd in NONchkvardims:
2006            if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd)
2007            varsadd.append(nonvd)
2008
2009        ncvar.insert_variable(ncobj, 'uaz', diagout1, diagoutd, diagoutvd, newnc)
2010        ncvar.insert_variable(ncobj, 'vaz', diagout2, diagoutd, diagoutvd, newnc)
2011
2012    else:
2013        print errormsg
2014        print '  ' + main + ": diagnostic '" + diagn + "' not ready!!!"
2015        print '    available diagnostics: ', availdiags
2016        quit(-1)
2017
2018    newnc.sync()
2019    # Adding that additional variables required to compute some diagnostics which
2020    #   where not in the original file
2021    print '  adding additional variables...'
2022    for vadd in varsadd:
2023        if not gen.searchInlist(newnc.variables.keys(),vadd) and                     \
2024          dictcompvars.has_key(vadd):
2025            attrs = dictcompvars[vadd]
2026            vvn = attrs['name']
2027            if not gen.searchInlist(newnc.variables.keys(), vvn):
2028                iidvn = dvnames.index(vadd)
2029                dnn = dnames[iidvn]
2030                if vadd == 'WRFtime':
2031                    dvarvals = WRFtime[:]
2032                newvar = newnc.createVariable(vvn, 'f8', (dnn))
2033                newvar[:] = dvarvals
2034                for attn in attrs.keys():
2035                    if attn != 'name':
2036                        attv = attrs[attn]
2037                        ncvar.set_attribute(newvar, attn, attv)
2038
2039#   end of diagnostics
2040
2041# Global attributes
2042##
2043ncvar.add_global_PyNCplot(newnc, main, None, '2.0')
2044
2045gorigattrs = ncobj.ncattrs()
2046for attr in gorigattrs:
2047    attrv = ncobj.getncattr(attr)
2048    atvar = ncvar.set_attribute(newnc, attr, attrv)
2049
2050ncobj.close()
2051newnc.close()
2052
2053print '\n' + main + ': successfull writting of diagnostics file "' + ofile + '" !!!'
Note: See TracBrowser for help on using the repository browser.