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