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