source: trunk/UTIL/PYTHON/planetoplot.py @ 584

Last change on this file since 584 was 582, checked in by acolaitis, 13 years ago

Python. Bug correction (another one, yay \! )

  • Property svn:executable set to *
File size: 41.1 KB
Line 
1#######################
2##### PLANETOPLOT #####
3#######################
4
5### A. Spiga     -- LMD -- 06~09/2011 -- General building and mapping capabilities
6### T. Navarro   -- LMD -- 10~11/2011 -- Improved use for GCM and added sections + 1Dplot capabilities
7### A. Colaitis  -- LMD --    11/2011 -- Mostly minor improvements and inter-plot operation capabilities + zrecast interpolation for gcm
8### A. Spiga     -- LMD -- 11~12/2011 -- Extended multivar subplot capabilities + cosmetic changes + general cleaning and tests
9### A. Colaitis  -- LMD --    12/2011 -- Added movie capability [mencoder must be installed]
10### A. Spiga     -- LMD --    12/2011 -- Added HTML animated page capability + general tests of consistency [winds, etc...] + consistent generic movie loop
11### J. Leconte   -- LMD --    02/2012 -- Added area weighted averaging. Compatibility with terrestrial gcm.
12### A. Spiga     -- LMD --    03/2012 -- Cleaning and improved comments
13def planetoplot (namefiles,\
14           level=0,\
15           vertmode=0,\
16           proj=None,\
17           back=None,\
18           target=None,
19           stride=3,\
20           var=None,\
21           colorb="def",\
22           winds=False,\
23           addchar=None,\
24           interv=[0,1],\
25           vmin=None,\
26           vmax=None,\
27           tile=False,\
28           zoom=None,\
29           display=True,\
30           hole=False,\
31           save="gui",\
32           anomaly=False,\
33           var2=None,\
34           ndiv=10,\
35           mult=1.,\
36           zetitle=["fill"],\
37           slon=None,\
38           slat=None,\
39           svert=None,\
40           stime=None,\
41           outputname=None,\
42           resolution=200,\
43           ope=None,\
44           fileref=None,\
45           minop=0.,\
46           maxop=0.,\
47           titleref="fill",\
48           invert_y=False,\
49           xaxis=[None,None],\
50           yaxis=[None,None],\
51           ylog=False,\
52           yintegral=False,\
53           blat=None,\
54           blon=None,\
55           tsat=False,\
56           flagnolow=False,\
57           mrate=None,\
58           mquality=False,\
59           trans=1,\
60           zarea=None,\
61           axtime=None,\
62           redope=None,\
63           seevar=False,\
64           xlab=None,\
65           ylab=None,\
66           lbls=None):
67
68    ####################################################################################################################
69    ### Colorbars http://www.scipy.org/Cookbook/Matplotlib/Show_colormaps?action=AttachFile&do=get&target=colormaps3.png
70
71    #################################
72    ### Load librairies and functions
73    from netCDF4 import Dataset
74    from myplot import getcoord2d,define_proj,makeplotres,simplinterv,vectorfield,ptitle,latinterv,getproj,wrfinterv,dumpbdy,\
75                       fmtvar,definecolorvec,defcolorb,getprefix,putpoints,calculate_bounds,errormess,definesubplot,\
76                       zoomset,getcoorddef,getwinddef,whatkindfile,reducefield,bounds,getstralt,getfield,smooth,nolow,\
77                       getname,localtime,check_localtime,polarinterv,getsindex,define_axis,determineplot,readslices,bidimfind,getlschar,hole_bounds,windamplitude
78    from mymath import deg,max,min,mean,get_tsat,writeascii,fig2data,fig2img
79    import matplotlib as mpl
80    from matplotlib.pyplot import contour,contourf, subplot, figure, rcParams, savefig, colorbar, pcolor, show, plot, clabel, title, close, legend, xlabel, axis, ylabel, subplots_adjust, axes
81    from matplotlib.cm import get_cmap
82    #from mpl_toolkits.basemap import cm
83    import numpy as np
84    from numpy.core.defchararray import find
85    from videosink import VideoSink
86    import subprocess
87    #from singlet import singlet
88
89#########################
90### PRELIMINARY STUFF ###
91#########################
92### we say hello
93    print "********************************************"
94    print "********** WELCOME TO PLANETOPLOT **********"
95    print "********************************************"
96### we ensure namefiles and var are arrays
97    if not isinstance(namefiles, np.ndarray): namefiles = [namefiles]
98    if not isinstance(var, np.ndarray):       var = [var]
99### we initialize a few variables
100    initime=-1 ; sslon = None ; sslat = None
101    if slon is not None: sslon = np.zeros([1,2])
102    if slat is not None: sslat = np.zeros([1,2])
103    k = 0 ; firstfile = True ; count = 0
104### we avoid specific cases not yet implemented
105    if mrate is not None and len(var) > 1: errormess("multivar not allowed in movies. should be fixed soon!")
106### we prepare number of plot fields [zelen] and plot instances [numplot] according to user choices
107### --> we support multifile and multivar plots : files + vars separated by commas are on the same figure
108    nlon, nlat, nvert, ntime, mapmode, nslices = determineplot(slon, slat, svert, stime)
109    zelen = len(namefiles)*len(var)
110    numplot = zelen*nslices
111    print "********** FILES, SLICES, VARS, TOTAL PLOTS: ", len(namefiles), nslices, len(var), numplot
112    print "********** MAPMODE: ", mapmode
113### we correct number of plot fields for possible operation (substract, etc...)
114    if ope is not None:
115        if fileref is not None:       zelen = zelen + 2
116        elif "var" in ope:            zelen = zelen + 1
117### we define the arrays for plot fields -- which will be placed within multiplot loops
118    all_var  = [[]]*zelen ; all_var2  = [[]]*zelen
119    all_title = [[]]*zelen ; all_varname = [[]]*zelen ; all_namefile = [[]]*zelen ; all_time = [[]]*zelen
120    all_windu = [[]]*zelen ; all_windv = [[]]*zelen
121 
122#############################
123### LOOP OVER PLOT FIELDS ###
124#############################
125    for nnn in range(len(namefiles)):
126     for vvv in range(len(var)): 
127
128    ### we load each NETCDF objects in namefiles
129      namefile = namefiles[nnn] 
130      nc  = Dataset(namefile)
131    ### we explore the variables in the file
132      varinfile = nc.variables.keys()
133      if seevar: print varinfile ; exit()
134    ### we define the type of file we have (gcm, meso, etc...)
135      typefile = whatkindfile(nc)
136      if firstfile:                 typefile0 = typefile
137      elif typefile != typefile0:   errormess("Not the same kind of files !", [typefile0, typefile])
138    ### we care for input file being 1D simulations
139      if typefile in ['gcm'] and len(nc.variables["longitude"][:]) is 1 and len(nc.variables["latitude"][:]) is 1:       mapmode=0 ; winds=False
140      elif typefile in ['earthgcm'] and len(nc.variables["lon"][:]) is 1 and len(nc.variables["lat"][:]) is 1:           mapmode=0 ; winds=False
141    ### we create default vert and time prescriptions if not here in case mapping mode is on (lat/lon maps)
142      if redope is None and mapmode == 1:
143          if svert is None:  svert = readslices(str(level)) ; nvert=1
144          if stime is None and mrate is None:
145             stime = readslices(str(0)) ; ntime=1 ## this is a default choice
146             print "WELL... nothing about time axis. I took default: first time reference stored in file."
147    ### we get the names of variables to be read. in case only one available, we choose this one.
148      varname=var[vvv] 
149      if ((varname not in nc.variables) and not ((typefile in ['meso']) and (varname in ['UV','uv','uvmet']))): 
150          if len(varinfile) == 1:   varname = varinfile[0] 
151          else:                     varname = False
152    ### we get the names of wind variables to be read (if any)
153      if winds:                                                   
154         [uchar,vchar,metwind] = getwinddef(nc)             
155         if uchar == 'not found': winds = False
156    ### we tell the user that either no var or no wind is not acceptable
157      if not varname and not winds: errormess("please set at least winds or var",printvar=nc.variables)
158    ### we get the coordinates lat/lon to be used
159      [lon2d,lat2d] = getcoorddef(nc)
160    ### we get an adapted map projection if none is provided by the user
161      if proj == None:   proj = getproj(nc)                 
162    ### we define plot boundaries given projection or user choices
163      if firstfile:
164         if proj in ["npstere","spstere"]: [wlon,wlat] = polarinterv(lon2d,lat2d)
165         elif proj in ["lcc","laea"]:      [wlon,wlat] = wrfinterv(lon2d,lat2d)
166         else:                             [wlon,wlat] = simplinterv(lon2d,lat2d)
167         if zoom:                          [wlon,wlat] = zoomset(wlon,wlat,zoom)
168         elif zarea is not None:           [wlon,wlat] = latinterv(area=zarea)
169
170##########################################################
171############ LOAD 4D DIMENSIONS : x, y, z, t #############
172##########################################################
173      if typefile in ["gcm","earthgcm","ecmwf"]:
174          ### LAT
175          if "latitude" in nc.variables:  lat = nc.variables["latitude"][:]
176          elif "lat" in nc.variables:     lat = nc.variables["lat"][:]
177          ### LON
178          if "longitude" in nc.variables: lon = nc.variables["longitude"][:]
179          elif "lon" in nc.variables:     lon = nc.variables["lon"][:]
180          ### ALT
181          if "altitude" in nc.variables:  vert = nc.variables["altitude"][:]
182          elif "Alt" in nc.variables:     vert = nc.variables["Alt"][:]
183          elif "lev" in nc.variables:     vert = nc.variables["lev"][:]
184          else:                           vert = [0.]
185          ### AIRE (to weight means with the area)
186          if "aire" in nc.variables:      area = nc.variables["aire"][:,:] 
187          else:                           area = None
188          ### TIME
189          if "Time" in nc.variables:            time = nc.variables["Time"][:]
190          elif "time" in nc.variables:          time = nc.variables["time"][:]
191          elif "time_counter" in nc.variables:  time = nc.variables["time_counter"][:]/86400. #### time counter cinverstion from s-> days
192          else:                                 errormess("no time axis found.")
193          if axtime in ["ls","sol"]:   errormess("not supported. should not be too difficult though.")
194          # for 1D plots (no need for longitude computation):
195          if axtime in ["lt"]:
196              if initime == -1: initime=input("Please type initial local time:")
197              time = (initime+time*24)%24
198              print "LOCAL TIMES.... ", time
199      elif typefile in ['meso','geo']:
200          area = None ## not active for the moment
201          ###### STUFF TO GET THE CORRECT LAT/LON FROM MESOSCALE FILES WITH 2D LAT/LON ARRAYS
202          ###### principle: calculate correct indices then repopulate slon and slat
203          if slon is not None or slat is not None:
204              if firstfile and save == 'png' and typefile == 'meso' and "HGT" in varinfile:   iwantawhereplot = nc     #show a topo map with a cross on the chosen point
205              else:                                                                           iwantawhereplot = None   #do not show anything, just select indices
206              numlon = 1 ; numlat = 1 
207              if slon is not None:   numlon = slon.shape[0]   
208              if slat is not None:   numlat = slat.shape[0]
209              indices = np.ones([numlon,numlat,2]) ; vlon = None ; vlat = None
210              for iii in range(numlon): 
211               for jjj in range(numlat):
212                 if slon is not None:  vlon = slon[iii][0]  ### note: slon[:][0] does not work
213                 if slat is not None:  vlat = slat[jjj][0]  ### note: slon[:][0] does not work
214                 indices[iii,jjj,:] = bidimfind(lon2d,lat2d,vlon,vlat,file=iwantawhereplot) 
215                 lonp,latp = ( lon2d[indices[iii,jjj,0],indices[iii,jjj,1]] , lat2d[indices[iii,jjj,0],indices[iii,jjj,1]] )
216              for iii in range(numlon):
217               for jjj in range(numlat):
218                 if slon is not None: sslon[iii][0] = indices[iii,0,1] ; sslon[iii][1] = indices[iii,0,1]  #...this is idx
219                 if slat is not None: sslat[jjj][0] = indices[0,jjj,0] ; sslat[jjj][1] = indices[0,jjj,0]  #...this is idy
220              lonp,latp = ( lon2d[indices[0,0,0],indices[0,0,1]] , lat2d[indices[0,0,0],indices[0,0,1]] )
221          ######
222          if typefile in ['meso'] and mapmode == 1: lon2d = dumpbdy(lon2d,6) ; lat2d = dumpbdy(lat2d,6)  ### important to do that now and not before
223          ######
224          if varname in ['PHTOT','W']:    vertdim='BOTTOM-TOP_PATCH_END_STAG'
225          else:                           vertdim='BOTTOM-TOP_PATCH_END_UNSTAG'
226          if (var2 is not None and var2 not in ['PHTOT','W']): dumped_vert_stag=True ; vertdim='BOTTOM-TOP_PATCH_END_UNSTAG'
227          else:                                                dumped_vert_stag=False
228          if varname in ['V']:  latdim='SOUTH-NORTH_PATCH_END_STAG'
229          else:                 latdim='SOUTH-NORTH_PATCH_END_UNSTAG'
230          if varname in ['U']:  londim='WEST-EAST_PATCH_END_STAG'
231          else:                 londim='WEST-EAST_PATCH_END_UNSTAG'
232          lon = np.arange(0,getattr(nc,londim),1) ; lat = np.arange(0,getattr(nc,latdim),1)
233          ###
234          if axtime in ["ls","sol"]:
235              lstab, soltab, lttab = getlschar ( namefile, getaxis = True )
236              if axtime == "ls":      time = lstab
237              elif axtime == "sol":   time = soltab
238          else:
239              if "Times" in nc.variables:   time = count + np.arange(0,len(nc.variables["Times"]),1)
240              elif "Time" in nc.variables:  time = count + np.arange(0,len(nc.variables["Time"]),1)
241              else:                         time = count + np.arange(0,1,1)
242              if ope in ["cat"]:
243                 if nnn > 0:  count = time[-1] + 1  ## so that a cat is possible with simple subscripts
244                 else:        count = 0
245              else: count = 0
246          if axtime in ["lt"]:
247              ftime = np.zeros(len(time))
248              for i in range(len(time)): 
249                 ftime[i] = localtime ( interv[0]+time[i]*interv[1], 0.5*(wlon[0]+wlon[1]) )
250              time=ftime
251              time=check_localtime(time)
252              print "LOCAL TIMES.... ", time
253
254          ###
255          if typefile in ['geo']:   vert = [0.] ; stime = readslices(str(0))
256          else:
257              if vertmode is None:  vertmode=0
258              if vertmode == 0:     vert = np.arange(0,getattr(nc,vertdim),1)
259              elif vertmode == -1:  vert = nc.variables["PHTOT"][0,:,0,0]/3.72 ; vert = np.array(vert[0:len(vert)-1]) #; print vert
260              else:                 vert = nc.variables["vert"][:]
261       #if firstfile:
262       #   lat0 = lat
263       #elif len(lat0) != len(lat):
264       #   errormess("Not the same latitude lengths !", [len(lat0), len(lat)])
265       #elif sum((lat == lat0) == False) != 0:
266       #   errormess("Not the same latitudes !", [lat,lat0])
267       ## Faire d'autre checks sur les compatibilites entre fichiers!!
268##########################################################
269##########################################################
270##########################################################
271
272      all_varname[k] = varname
273      all_namefile[k] = namefile
274      all_time[k] = time
275      if var2: all_var2[k] = getfield(nc,var2)
276      if winds: all_windu[k] = getfield(nc,uchar) ; all_windv[k] = getfield(nc,vchar)
277      #####################
278      ##### GETFIELDS #####
279      #####################
280      ##### SPECIFIC CASES
281      ##### 1. saturation temperature
282      if varname in ["temp","t","T_nadir_nit","T_nadir_day","temp_day","temp_night"] and tsat:
283          tt=getfield(nc,varname) ; print "computing Tsat-T, I ASSUME Z-AXIS IS PRESSURE"
284          if type(tt).__name__=='MaskedArray':  tt.set_fill_value([np.NaN]) ; tinput=tt.filled()
285          else:                                 tinput=tt
286          all_var[k]=get_tsat(vert,tinput,zlon=lon,zlat=lat,zalt=vert,ztime=time)
287      ##### 2. wind amplitude
288      elif ((varname in ['UV','uv','uvmet']) and (typefile in ['meso']) and (varname not in nc.variables)):
289          all_var[k]=windamplitude(nc)
290      else:
291      ##### GENERAL STUFF HERE
292          all_var[k] = getfield(nc,varname)
293      print "**** GOT SUBDATA:",k," NAMEFILE:",namefile," VAR:",varname, var2 ; k += 1 ; firstfile = False
294      #### End of for namefile in namefiles
295
296    ##################################
297    ### Operation on files
298    if ope is not None:
299        print "********** OPERATION: ",ope
300        if "var" not in ope:
301             if len(var) > 1: errormess("for this operation... please set only one var !")
302             if ope in ["-","+","-%"]:
303                if fileref is not None:   
304                   ncref = Dataset(fileref)
305                   if ((all_varname[k-1] in ['UV','uv','uvmet']) and (typefile in ['meso']) and (all_varname[k-1] not in ncref.variables)):
306                      all_var[k] = windamplitude(ncref)
307                   else:
308                      all_var[k] = getfield(ncref,all_varname[k-1])
309                   all_varname[k] = all_varname[k-1] ; all_time[k] = all_time[k-1] ; all_namefile[k] = all_namefile[k-1] ; all_var2[k] = all_var2[k-1]
310                else: errormess("fileref is missing!") 
311                if ope == "-":     all_var[k+1]= all_var[k-1] - all_var[k]
312                elif ope == "+":   all_var[k+1]= all_var[k-1] + all_var[k]
313                elif ope == "-%":
314                    masked = np.ma.masked_where(all_var[k] == 0,all_var[k])
315                    masked.set_fill_value([np.NaN])
316                    all_var[k+1]= 100.*(all_var[k-1] - masked)/masked
317                all_varname[k+1] = all_varname[k] ; all_time[k+1] = all_time[k] ; all_namefile[k+1] = all_namefile[k] ; all_var2[k+1] = all_var2[k] ; numplot = numplot+2
318             elif ope in ["cat"]:
319                tabtime = all_time[0];tab = all_var[0];k = 1
320                if var2: tab2 = all_var2[0]
321                while k != len(namefiles) and len(all_time[k]) != 0:
322                    if var2: tab2 = np.append(tab2,all_var2[k],axis=0) 
323                    tabtime = np.append(tabtime,all_time[k]) ; tab = np.append(tab,all_var[k],axis=0) ; k += 1
324                all_time[0] = np.array(tabtime) ; all_var[0] = np.array(tab) ; numplot = 1
325                if var2: all_var2[0] = np.array(tab2)
326             else: errormess(ope+" : non-implemented operation. Check pp.py --help")
327        else:
328             if len(namefiles) > 1: errormess("for this operation... please set only one file !") 
329             if len(var) > 2:       errormess("not sure this works for more than 2 vars... please check.")
330             if   ope in ["div_var"]: all_var[k] = all_var[k-2] / all_var[k-1] ; insert = '_div_'
331             elif ope in ["mul_var"]: all_var[k] = all_var[k-2] * all_var[k-1] ; insert = '_mul_'
332             elif ope in ["add_var"]: all_var[k] = all_var[k-2] + all_var[k-1] ; insert = '_add_'
333             elif ope in ["sub_var"]: all_var[k] = all_var[k-2] - all_var[k-1] ; insert = '_sub_'
334             else:                    errormess(ope+" : non-implemented operation. Check pp.py --help")
335             numplot = numplot + 1 ; all_time[k] = all_time[k-1] ; all_namefile[k] = all_namefile[k-1]
336             all_varname[k] = all_varname[k-2] + insert + all_varname[k-1] 
337    ##################################
338    ### Open a figure and set subplots
339    fig = figure()
340    subv,subh = definesubplot( numplot, fig ) 
341    if ope in ['-','-%']: subv,subh = 2,2
342 
343    #################################
344    ### Time loop for plotting device
345    nplot = 1 ; error = False 
346    print "********************************************"
347    while error is False:
348     
349       print "********** NPLOT", nplot
350       if nplot > numplot: break
351
352       ####################################################################
353       ## get all indexes to be taken into account for this subplot and then reduce field
354       ## We plot 1) all lon slices 2) all lat slices 3) all vert slices 4) all time slices and then go to the next slice
355       indexlon  = getsindex(sslon,(nplot-1)%nlon,lon)
356       indexlat  = getsindex(sslat,((nplot-1)//nlon)%nlat,lat)
357       indexvert = getsindex(svert,((nplot-1)//(nlon*nlat))%nvert,vert) 
358       if ope is not None:
359           if fileref is not None:      index_f = ((nplot-1)//(nlon*nlat*nvert*ntime))%(len(namefiles)+2)  ## OK only 1 var,  see test in the beginning
360           elif "var" in ope:           index_f = ((nplot-1)//(nlon*nlat*nvert*ntime))%(len(var)+1)        ## OK only 1 file, see test in the beginning
361           elif "cat" in ope:           index_f = 0
362       else:                            yeah = len(namefiles)*len(var) ; index_f = ((nplot-1)//(nlon*nlat*nvert*ntime))%yeah
363       time = all_time[index_f]
364       if mrate is not None:                 indextime = None 
365       else:                                 indextime = getsindex(stime,((nplot-1)//(nlon*nlat*nvert))%ntime,time)
366       ltst = None 
367       if typefile in ['meso'] and indextime is not None:  ltst = localtime ( interv[0]+indextime*interv[1], 0.5*(wlon[0]+wlon[1]) ) 
368       print "********** INDEX LON:",indexlon," LAT:",indexlat," VERT:",indexvert," TIME:",indextime
369       ##var = nc.variables["phisinit"][:,:]
370       ##contourf(np.transpose(var),30,cmap = get_cmap(name="Greys_r") ) ; axis('off') ; plot(indexlat,indexlon,'mx',mew=4.0,ms=20.0)
371       ##show()
372       ##exit()
373       #truc = True
374       #truc = False
375       #if truc: indexvert = None
376       ####################################################################
377       ########## REDUCE FIELDS
378       ####################################################################
379       error = False
380       varname = all_varname[index_f]
381       if varname:   ### what is shaded.
382           what_I_plot, error = reducefield( all_var[index_f], d4=indextime, d1=indexlon, d2=indexlat, d3=indexvert, \
383                                             yint=yintegral, alt=vert, anomaly=anomaly, redope=redope, mesharea=area )
384           if mult != 2718.:  what_I_plot = what_I_plot*mult
385           else:              what_I_plot = np.log10(what_I_plot) ; print "log plot"
386       if var2:      ### what is contoured.
387           what_I_plot_contour, error = reducefield( all_var2[index_f], d4=indextime, d1=indexlon, d2=indexlat , d3=indexvert, \
388                                             yint=yintegral, alt=vert )
389       if winds:     ### what is plot as vectors.
390           vecx, error = reducefield( all_windu[index_f], d4=indextime, d3=indexvert, yint=yintegral, alt=vert)
391           vecy, error = reducefield( all_windv[index_f], d4=indextime, d3=indexvert, yint=yintegral, alt=vert)
392           if varname in [uchar,vchar]: what_I_plot = np.sqrt( np.square(vecx) + np.square(vecy) ) ; varname = "wind"
393
394       #####################################################################
395       #if truc:
396       #   nx = what_I_plot.shape[2] ; ny = what_I_plot.shape[1] ; nz = what_I_plot.shape[0]
397       #   for k in range(nz): print k,' over ',nz ; what_I_plot[k,:,:] = what_I_plot[k,:,:] / smooth(what_I_plot[k,:,:],12)
398       #   for iii in range(nx):
399       #    for jjj in range(ny):
400       #     deviation = what_I_plot[:,jjj,iii] ; mx = max(deviation) ; mn = min(deviation)
401       #     if iii > 6 and iii < nx-6 and jjj > 6 and jjj < ny-6:   what_I_plot[0,jjj,iii],rel = singlet(deviation,vert/1000.)  ### z must be in km
402       #     else:                                                   what_I_plot[0,jjj,iii]     = 0.
403       #     if np.abs(what_I_plot[0,jjj,iii]) > 1.5:
404       #         print iii,jjj,what_I_plot[0,jjj,iii],int(abs(1.-mx)*100.),int(abs(1.-mn)*100.)
405       #         plot(rel)
406       #   show()
407       #   anomaly = True ### pour avoir les bons reglages plots
408       #   what_I_plot = what_I_plot[0,:,:] 
409       #####################################################################
410
411       ####################################################################
412       ### General plot settings
413       changesubplot = (numplot > 1) and (len(what_I_plot.shape) != 1)  ## default for 1D plots: superimposed. to be reworked for better flexibility.
414       if changesubplot: subplot(subv,subh,nplot) #; subplots_adjust(wspace=0,hspace=0)
415       ####################################################################
416       if error:
417               errormess("There is an error in reducing field !")
418       else:
419               ticks = ndiv + 1 
420               fvar = varname
421               if anomaly: fvar = 'anomaly'
422               ###
423               if mapmode == 0:    ### could this be moved inside imov loop ?
424                   itime=indextime
425                   if len(what_I_plot.shape) is 3: itime=[0]
426                   m = None ; x = None ; y = None
427                   what_I_plot, x, y = define_axis(lon,lat,vert,time,indexlon,indexlat,indexvert,\
428                         itime,what_I_plot, len(all_var[index_f].shape),vertmode)
429               ###
430               if (fileref is not None) and (index_f == numplot-1):    zevmin, zevmax = calculate_bounds(what_I_plot,vmin=minop,vmax=maxop)
431               else:                                                   zevmin, zevmax = calculate_bounds(what_I_plot,vmin=vmin,vmax=vmax)
432               if (fileref is not None) and (index_f == numplot-1):    colorb = "RdBu_r"
433               if colorb in ["def","nobar","onebar"]:                  palette = get_cmap(name=defcolorb(fvar.upper()))
434               else:                                                   palette = get_cmap(name=colorb)
435               #palette = cm.GMT_split
436               ##### 1. ELIMINATE 0D or >3D CASES
437               if len(what_I_plot.shape) == 0:   
438                 print "VALUE VALUE VALUE VALUE ::: ",what_I_plot
439                 save = 'donothing'
440               elif len(what_I_plot.shape) >= 4:
441                 print "WARNING!!! ",len(what_I_plot.shape),"-D PLOT NOT SUPPORTED !!! dimensions: ",what_I_plot.shape
442                 errormess("Are you sure you did not forget to prescribe a dimension ?")
443               ##### 2. HANDLE simple 1D/2D field and movies of 1D/2D fields
444               else:
445                 if mrate is not None: iend=len(time)-1
446                 else:                 iend=0
447                 imov = 0 
448                 if len(what_I_plot.shape) == 3:
449                    if var2:               which = "contour" ## have to start with contours rather than shading
450                    else:                  which = "regular"
451                    if mrate is None:      errormess("3D field. Use --rate RATE for movie or specify --time TIME. Exit.")
452                 elif len(what_I_plot.shape) == 2:
453                    if var2:               which = "contour" ## have to start with contours rather than shading
454                    else:                  which = "regular"
455                    if mrate is not None:  which = "unidim"
456                 elif len(what_I_plot.shape) == 1:
457                    which = "unidim"
458                    if what_I_plot.shape[-1] == 1:      print "VALUE VALUE VALUE VALUE ::: ", what_I_plot[0] ; save = 'donothing'
459                 ##### IMOV LOOP #### IMOV LOOP
460                 while imov <= iend:
461                    print "-> frame ",imov+1, which
462                    if which == "regular":   
463                        if mrate is None:                                   what_I_plot_frame = what_I_plot
464                        else:                                               what_I_plot_frame = what_I_plot[imov,:,:]
465                        if winds:
466                            if mrate is None:                                   vecx_frame = vecx ; vecy_frame = vecy
467                            else:                                               vecx_frame = vecx[imov,:,:] ; vecy_frame = vecy[imov,:,:]
468                    elif which == "contour": 
469                        if mrate is None or what_I_plot_contour.ndim < 3:   what_I_plot_frame = what_I_plot_contour
470                        else:                                               what_I_plot_frame = what_I_plot_contour[imov,:,:]
471                    elif which == "unidim":
472                        if mrate is None:                                   what_I_plot_frame = what_I_plot
473                        else:                                               what_I_plot_frame = what_I_plot[:,imov]  ## because swapaxes
474                    #if mrate is not None:     
475                    if mapmode == 1: 
476                            m = define_proj(proj,wlon,wlat,back=back,blat=blat,blon=blon)  ## this is dirty, defined above but out of imov loop
477                            x, y = m(lon2d, lat2d)                                         ## this is dirty, defined above but out of imov loop
478                    if typefile in ['meso'] and mapmode == 1:   what_I_plot_frame = dumpbdy(what_I_plot_frame,6,condition=True)
479#                   if typefile in ['mesoideal']:    what_I_plot_frame = dumpbdy(what_I_plot_frame,0,stag='W',condition=dumped_vert_stag)
480
481                    if which == "unidim":
482                        lbl = ""
483                        if indexlat is not None:  lbl = lbl + " ix" + str(indexlat[0])
484                        if indexlon is not None:  lbl = lbl + " iy" + str(indexlon[0])
485                        if indexvert is not None: lbl = lbl + " iz" + str(indexvert[0])
486                        if indextime is not None: lbl = lbl + " it" + str(indextime[0])
487                        if lbl == "": lbl = namefiles[index_f]
488
489                        if lbls is not None: lbl=lbls[index_f]
490
491                        if mrate is not None: x = y  ## because swapaxes...
492                        #what_I_plot_frame = np.diff(what_I_plot_frame, n=1) ; x = x[1:]
493                       
494                        if not tile:  zeline='-'
495                        else:         zeline=','
496                        this_is_a_regular_plot = (indexvert is not None) or (indextime is None) or (indexlat is None) or (indexlon is None)
497                        if this_is_a_regular_plot:   plot(x,what_I_plot_frame,zeline,label=lbl)  ## vertical profile
498                        else:                        plot(what_I_plot_frame,x,zeline,label=lbl)  ## regular plot
499                        if nplot > 1: legend(loc='best')
500                        if indextime is None and axtime is not None and xlab is None:    xlabel(axtime.upper()) ## define the right label
501                        if save == 'txt':  writeascii(np.transpose(what_I_plot),'profile'+str(nplot*1000+imov)+'.txt')
502
503                    elif which == "regular":
504                        if hole:         what_I_plot_frame = hole_bounds(what_I_plot_frame,zevmin,zevmax)
505                        else:            what_I_plot_frame = bounds(what_I_plot_frame,zevmin,zevmax)
506                        if flagnolow:    what_I_plot_frame = nolow(what_I_plot_frame)
507                        if not tile:
508                            #zelevels = np.linspace(zevmin*(1. + 1.e-7),zevmax*(1. - 1.e-7)) #,num=20)
509                            zelevels = np.linspace(zevmin,zevmax,num=ticks)
510                            #what_I_plot_frame = smooth(what_I_plot_frame,100)
511                            if mapmode == 1:       m.contourf( x, y, what_I_plot_frame, zelevels, cmap = palette, alpha=trans)
512                            elif mapmode == 0:     contourf( x, y, what_I_plot_frame, zelevels, cmap = palette, alpha=trans)
513                        else:
514                            if mapmode == 1:       m.pcolor( x, y, what_I_plot_frame, cmap = palette, vmin=zevmin, vmax=zevmax, alpha=trans)
515                            elif mapmode == 0:     pcolor( x, y, what_I_plot_frame, cmap = palette, vmin=zevmin, vmax=zevmax, alpha=trans)
516
517                        if colorb not in ['nobar','onebar']:       
518                            if (fileref is not None) and (index_f == numplot-1):   daformat = "%.3f" 
519                            elif mult != 1:                                        daformat = "%.1f"
520                            else:                                                  daformat = fmtvar(fvar.upper())
521                            if proj in ['moll']:  zeorientation="horizontal"
522                            else:                 zeorientation="vertical"
523                            zecb = colorbar( fraction=0.05,pad=0.03,format=daformat,orientation=zeorientation,\
524                                      ticks=np.linspace(zevmin,zevmax,num=min([ticks/2+1,21])),extend='neither',spacing='proportional' ) 
525                            if zeorientation == "horizontal" and zetitle[0] != "fill": zecb.ax.set_xlabel(zetitle[index_f]) ; zetitle[index_f]=""
526                        if winds:
527                            if typefile in ['meso']:
528                                [vecx_frame,vecy_frame] = [dumpbdy(vecx_frame,6,stag=uchar,condition=True), dumpbdy(vecy_frame,6,stag=vchar,condition=True)]
529                                key = True
530                            elif typefile in ['gcm']:
531                                key = False
532                            if metwind and mapmode == 1:   [vecx_frame,vecy_frame] = m.rotate_vector(vecx_frame, vecy_frame, lon2d, lat2d)
533                            if var:       colorvec = definecolorvec(back)
534                            else:         colorvec = definecolorvec(colorb)
535                            vectorfield(vecx_frame, vecy_frame, x, y, stride=stride, csmooth=2,\
536                                             #scale=15., factor=300., color=colorvec, key=key)
537                                             scale=20., factor=250., color=colorvec, key=key)
538                                                              #200.         ## or csmooth=stride
539                    elif which == "contour":
540                        zevminc, zevmaxc = calculate_bounds(what_I_plot_frame, vmin=min(what_I_plot_frame), vmax=max(what_I_plot_frame))
541                        zelevels = np.linspace(zevminc,zevmaxc,ticks/2) #20)
542                        if var2 == 'HGT': zelevels = np.arange(-10000.,30000.,2000.)
543                        if mapmode == 0:   
544                            what_I_plot_frame, x, y = define_axis( lon,lat,vert,time,indexlon,indexlat,indexvert,\
545                                                              itime,what_I_plot_frame, len(all_var2[index_f].shape),vertmode )
546                            cs = contour( x,y,what_I_plot_frame, zelevels, colors='k', linewidths = 1 ) #0.33 colors='w' )# , alpha=0.5)#, linestyles='solid')
547                        elif mapmode == 1:  cs = m.contour( x,y,what_I_plot_frame, zelevels, colors='k', linewidths = 1 ) #0.33 colors='w' )# , alpha=0.5)# , linestyles='solid')
548                    if which in ["regular","unidim"]:
549
550                        if nplot > 1 and which == "unidim":
551                           pass  ## because we superimpose nplot instances
552                        else:
553                           # Axis directives for movie frames [including the first one).
554                           zxmin, zxmax = xaxis ; zymin, zymax = yaxis
555                           if zxmin is not None: mpl.pyplot.xlim(xmin=zxmin)
556                           if zxmax is not None: mpl.pyplot.xlim(xmax=zxmax)
557                           if zymin is not None: mpl.pyplot.ylim(ymin=zymin)
558                           if zymax is not None: mpl.pyplot.ylim(ymax=zymax)
559                           if ylog:      mpl.pyplot.semilogy()
560                           if invert_y:  ax = mpl.pyplot.gca() ; ax.set_ylim(ax.get_ylim()[::-1])
561                           if xlab is not None: xlabel(xlab)
562                           if ylab is not None: ylabel(ylab)
563
564                        if mrate is not None:
565                           ### THIS IS A MENCODER MOVIE
566                           if mrate > 0:
567                             figframe=mpl.pyplot.gcf()
568                             if mquality:   figframe.set_dpi(600.)
569                             else:          figframe.set_dpi(200.)
570                             mframe=fig2img(figframe)
571                             if imov == 0:
572                                moviename='movie' ;W,H = figframe.canvas.get_width_height()
573                                video = VideoSink((H,W), moviename, rate=mrate, byteorder="rgba")
574                             video.run(mframe) ; close()
575                             if imov == iend: video.close()                           
576                           ### THIS IS A WEBPAGE MOVIE
577                           else:
578                             nameframe = "image"+str(1000+imov)
579                             makeplotres(nameframe,res=100.,disp=False) ; close()
580                             if imov == 0: myfile = open("zepics", 'w')
581                             myfile.write("modImages["+str(imov)+"] = '"+nameframe+"_100.png';"+ '\n')
582                             if imov == iend:
583                                 myfile.write("first_image = 0;"+ '\n')
584                                 myfile.write("last_image = "+str(iend)+";"+ '\n')
585                                 myfile.close()
586                        if var2 and which == "regular":  which = "contour"
587                        imov = imov+1
588                    elif which == "contour":
589                        which = "regular"
590
591       ### Next subplot
592       zevarname = varname
593       if redope is not None: zevarname = zevarname + "_" + redope
594       basename = getname(var=zevarname,var2=var2,winds=winds,anomaly=anomaly)
595       if len(what_I_plot.shape) > 3:
596           basename = basename + getstralt(nc,level) 
597       if mrate is not None: basename = "movie_" + basename
598       if typefile in ['meso']:
599            if sslon is not None: basename = basename + "_lon_" + str(int(round(lonp)))
600            if sslat is not None: basename = basename + "_lat_" + str(int(round(latp)))
601            plottitle = basename
602            ### dans le nouveau systeme time=ls,sol,lt cette ligne pourrait ne servir a rien (ou deplacer au dessus)
603            if addchar and indextime is not None:   [addchar,gogol,gogol2] = getlschar ( all_namefile[index_f] )  ;  plottitle = plottitle + addchar
604            ### en fait redope is None doit etre remplace par : n'est ni maxt ni mint
605            if redope is None and ltst is not None and ( (mapmode == 0) or (proj in ["lcc","laea","merc","nsper"]) ):  plottitle = plottitle + "_LT" + str(ltst)
606       else:
607            if fileref is not None:
608                if index_f is numplot-1:     plottitle = basename+' '+"fig(1) "+ope+" fig(2)"
609                elif index_f is numplot-2:   plottitle = basename+' '+fileref
610                else:                        plottitle = basename+' '+namefiles[0]#index_f]
611            else:                            plottitle = basename+' '+namefiles[0]#index_f]
612       if mult != 1:                         plottitle = '{:.0e}'.format(mult) + "*" + plottitle
613       if zetitle[0] != "fill":                 
614          if index_f < len(zetitle): plottitle = zetitle[index_f]
615          else: plottitle = zetitle[0]
616          if titleref is "fill":             titleref=zetitle[index_f]
617          if fileref is not None:
618             if index_f is numplot-2:        plottitle = titleref
619             if index_f is numplot-1:        plottitle = "fig(1) "+ope+" fig(2)"
620#       if indexlon is not None:      plottitle = plottitle + " lon: " + str(min(lon[indexlon])) +" "+ str(max(lon[indexlon]))
621#       if indexlat is not None:      plottitle = plottitle + " lat: " + str(min(lat[indexlat])) +" "+ str(max(lat[indexlat]))
622#       if indexvert is not None:     plottitle = plottitle + " vert: " + str(min(vert[indexvert])) +" "+ str(max(vert[indexvert]))
623#       if indextime is not None:     plottitle = plottitle + " time: " + str(min(time[indextime])) +" "+ str(max(time[indextime]))
624       if colorb != "onebar": title( plottitle )
625       if nplot >= numplot: error = True
626       nplot += 1
627
628    if colorb == "onebar":
629        cax = axes([0.1, 0.2, 0.8, 0.03]) # a ameliorer
630        zecb = colorbar(cax=cax, orientation="horizontal", format=fmtvar(fvar.upper()),\
631                 ticks=np.linspace(zevmin,zevmax,num=min([ticks/2+1,21])),extend='neither',spacing='proportional')
632        if zetitle[0] != "fill": zecb.ax.set_xlabel(zetitle[index_f]) ; zetitle[index_f]=""
633
634     
635    ##########################################################################
636    ### Save the figure in a file in the data folder or an user-defined folder
637    if outputname is None:
638       if typefile in ['meso']:   prefix = getprefix(nc)
639       elif typefile in ['gcm']:            prefix = 'LMD_GCM_'
640       else:                                prefix = ''
641    ###
642       zeplot = prefix + basename
643       if addchar:         zeplot = zeplot + addchar
644       if numplot <= 0:    zeplot = zeplot + "_LT"+str(abs(numplot))
645    ###
646       if not target:      zeplot = namefile[0:find(namefile,'wrfout')] + zeplot
647       else:               zeplot = target + "/" + zeplot 
648    ###
649    else:
650       zeplot=outputname
651
652    if mrate is None:
653        pad_inches_value = 0.35
654        print "********** SAVE ", save
655        if save == 'png': 
656            if display: makeplotres(zeplot,res=100.,pad_inches_value=pad_inches_value) #,erase=True)  ## a miniature
657            makeplotres(zeplot,res=resolution,pad_inches_value=pad_inches_value,disp=False)
658        elif save in ['eps','svg','pdf']:     makeplotres(zeplot,pad_inches_value=pad_inches_value,disp=False,ext=save)
659        elif save == 'gui':                   show()
660        elif save == 'donothing':             pass
661        elif save == 'txt':                   print "Saved results in txt file." 
662        else: 
663            print "INFO: save mode not supported. using gui instead."
664            show()
665
666    ###################################
667    #### Getting more out of this video -- PROBLEMS WITH CREATED VIDEOS
668    #
669    #if mrate is not None:
670    #    print "Re-encoding movie.. first pass"
671    #    video.first_pass(filename=moviename,quality=mquality,rate=mrate)
672    #    print "Re-encoding movie.. second pass"
673    #    video.second_pass(filename=moviename,quality=mquality,rate=mrate)   
674
675    ###############
676    ### Now the end
677    return zeplot
Note: See TracBrowser for help on using the repository browser.