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

Last change on this file since 781 was 781, checked in by aslmd, 12 years ago

UTIL PYTHON. added a gfortran script to compile mcd fortran stuff in python. added comments here and there.

  • Property svn:executable set to *
File size: 59.4 KB
RevLine 
[349]1#######################
2##### PLANETOPLOT #####
3#######################
[345]4
[349]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
[392]7### A. Colaitis  -- LMD --    11/2011 -- Mostly minor improvements and inter-plot operation capabilities + zrecast interpolation for gcm
[448]8### A. Spiga     -- LMD -- 11~12/2011 -- Extended multivar subplot capabilities + cosmetic changes + general cleaning and tests
[435]9### A. Colaitis  -- LMD --    12/2011 -- Added movie capability [mencoder must be installed]
[475]10### A. Spiga     -- LMD --    12/2011 -- Added HTML animated page capability + general tests of consistency [winds, etc...] + consistent generic movie loop
[525]11### J. Leconte   -- LMD --    02/2012 -- Added area weighted averaging. Compatibility with terrestrial gcm.
[579]12### A. Spiga     -- LMD --    03/2012 -- Cleaning and improved comments
[725]13### T. Navarro   -- LMD --    04/2012 -- Added capabilities (e.g. histograms for difference maps)
14### A. Colaitis  -- LMD -- 05~06/2012 -- Added capabilities for analysis of mesoscale files (wind speed, etc...)
15### A. Spiga     -- LMD -- 04~07/2012 -- Added larger support of files + planets. Corrected a few bugs. Cleaning and improved comments
[781]16### A. Colaitis  -- LMD --    08/2012 -- Added functionalities for further analysis: spectra, hodographs, etc... plus improved flexibility (e.g. grid)
[725]17
[345]18def planetoplot (namefiles,\
[351]19           level=0,\
[350]20           vertmode=0,\
[345]21           proj=None,\
22           back=None,\
23           target=None,
24           stride=3,\
25           var=None,\
[760]26           clb=None,\
[399]27           winds=False,\
[345]28           addchar=None,\
29           vmin=None,\
30           vmax=None,\
31           tile=False,\
32           zoom=None,\
33           display=True,\
34           hole=False,\
35           save="gui",\
36           anomaly=False,\
37           var2=None,\
38           ndiv=10,\
39           mult=1.,\
[578]40           zetitle=["fill"],\
[345]41           slon=None,\
42           slat=None,\
43           svert=None,\
[359]44           stime=None,\
45           outputname=None,\
46           resolution=200,\
47           ope=None,\
48           fileref=None,\
49           minop=0.,\
50           maxop=0.,\
[363]51           titleref="fill",\
[369]52           invert_y=False,\
53           xaxis=[None,None],\
[372]54           yaxis=[None,None],\
[382]55           ylog=False,\
[763]56           xlog=False,\
[385]57           yintegral=False,\
[388]58           blat=None,\
[451]59           blon=None,\
[418]60           tsat=False,\
[430]61           flagnolow=False,\
[444]62           mrate=None,\
[453]63           mquality=False,\
[457]64           trans=1,\
[468]65           zarea=None,\
[483]66           axtime=None,\
[502]67           redope=None,\
[507]68           seevar=False,\
[502]69           xlab=None,\
[569]70           ylab=None,\
[612]71           lbls=None,\
72           lstyle=None,\
[647]73           cross=None,\
[721]74           facwind=1.,\
[760]75           trycol=False,\
[763]76           streamflag=False,\
77           analysis=None):
[345]78
79    ####################################################################################################################
80    ### Colorbars http://www.scipy.org/Cookbook/Matplotlib/Show_colormaps?action=AttachFile&do=get&target=colormaps3.png
81
82    #################################
83    ### Load librairies and functions
84    from netCDF4 import Dataset
85    from myplot import getcoord2d,define_proj,makeplotres,simplinterv,vectorfield,ptitle,latinterv,getproj,wrfinterv,dumpbdy,\
86                       fmtvar,definecolorvec,defcolorb,getprefix,putpoints,calculate_bounds,errormess,definesubplot,\
87                       zoomset,getcoorddef,getwinddef,whatkindfile,reducefield,bounds,getstralt,getfield,smooth,nolow,\
[612]88                       getname,localtime,check_localtime,polarinterv,getsindex,define_axis,determineplot,readslices,bidimfind,getlschar,hole_bounds,\
[760]89                       getdimfromvar,select_getfield
[754]90    from mymath import deg,max,min,mean,writeascii,fig2data,fig2img
[363]91    import matplotlib as mpl
[725]92    from matplotlib.pyplot import contour,contourf,hist, text,subplot, figure, rcParams, savefig, colorbar, \
93                                  pcolor, show, plot, clabel, title, close, legend, xlabel, axis, ylabel, \
94                                  subplots_adjust, axes, clabel
[345]95    from matplotlib.cm import get_cmap
[490]96    #from mpl_toolkits.basemap import cm
[345]97    import numpy as np
98    from numpy.core.defchararray import find
[763]99    from scipy.stats import gaussian_kde,describe
[430]100    from videosink import VideoSink
[673]101    from times import sol2ls
[444]102    import subprocess
[490]103    #from singlet import singlet
[587]104    from itertools import cycle
[647]105    import os
[359]106
[579]107#########################
108### PRELIMINARY STUFF ###
109#########################
110### we say hello
[392]111    print "********************************************"
112    print "********** WELCOME TO PLANETOPLOT **********"
113    print "********************************************"
[579]114### we ensure namefiles and var are arrays
[392]115    if not isinstance(namefiles, np.ndarray): namefiles = [namefiles]
116    if not isinstance(var, np.ndarray):       var = [var]
[579]117### we initialize a few variables
118    initime=-1 ; sslon = None ; sslat = None
[760]119    k = 0 ; firstfile = True ; count = 0
120### we perform sanity checks and correct for insufficient information
[582]121    if slon is not None: sslon = np.zeros([1,2])
122    if slat is not None: sslat = np.zeros([1,2])
[760]123    if clb is None:            clb = ["def"]*len(var)
124    elif len(clb) < len(var):  clb = [clb[0]]*len(var) ; print "WARNING: less color than vars! setting all to 1st value."
125### we set option trycol i.e. the user wants to try a set of colorbars
126    if trycol: clb = ["Greys","Blues","YlOrRd","jet","spectral","hot","RdBu","RdYlBu","Paired"] ; zetitle = clb ; var = [var[0]]*9
[579]127### we avoid specific cases not yet implemented
128    if mrate is not None and len(var) > 1: errormess("multivar not allowed in movies. should be fixed soon!")
[760]129### we prepare number of plot fields [zelen] and number of plot instances [numplot] according to user choices
[579]130### --> we support multifile and multivar plots : files + vars separated by commas are on the same figure
[763]131    nlon, nlat, nvert, ntime, mapmode, nslices = determineplot(slon, slat, svert, stime, redope)
[392]132    zelen = len(namefiles)*len(var)
[763]133### we correct number of plot fields for possible operation (substract, etc...)
134    if ope is not None:
135        if fileref is not None:       zelen = 3*len(var)*len(namefiles)
136        elif "var" in ope:            zelen = zelen + 1
[424]137    numplot = zelen*nslices
[392]138    print "********** FILES, SLICES, VARS, TOTAL PLOTS: ", len(namefiles), nslices, len(var), numplot
[579]139    print "********** MAPMODE: ", mapmode
140### we define the arrays for plot fields -- which will be placed within multiplot loops
[763]141    all_var  = [[]]*zelen ; all_var2  = [[]]*zelen
142    all_title = [[]]*zelen ; all_varname = [[]]*zelen ; all_namefile = [[]]*zelen ; all_colorb = [[]]*zelen
143    all_time = [[]]*zelen ; all_vert = [[]]*zelen ; all_lat = [[]]*zelen ; all_lon = [[]]*zelen
[579]144    all_windu = [[]]*zelen ; all_windv = [[]]*zelen
[763]145    plot_x = [[]]*zelen ; plot_y = [[]]*zelen ; multiplot = [[]]*zelen
146### tool (should be move to mymath)
147    getVar = lambda searchList, ind: [searchList[i] for i in ind]
[579]148#############################
149### LOOP OVER PLOT FIELDS ###
150#############################
[763]151   
[392]152    for nnn in range(len(namefiles)):
153     for vvv in range(len(var)): 
154
[579]155    ### we load each NETCDF objects in namefiles
[468]156      namefile = namefiles[nnn] 
[345]157      nc  = Dataset(namefile)
[579]158    ### we explore the variables in the file
[507]159      varinfile = nc.variables.keys()
[748]160      if seevar: print varinfile ; exit() 
[579]161    ### we define the type of file we have (gcm, meso, etc...)
[507]162      typefile = whatkindfile(nc)
[579]163      if firstfile:                 typefile0 = typefile
164      elif typefile != typefile0:   errormess("Not the same kind of files !", [typefile0, typefile])
[718]165    ### we care for input file being 1D simulations
166      is1d=999 
[637]167      if "longitude" in nc.dimensions and "latitude" in nc.dimensions: is1d = len(nc.variables["longitude"][:])*len(nc.variables["latitude"][:])
168      elif "lon" in nc.dimensions and "lat" in nc.dimensions: is1d = len(nc.variables["lon"][:])*len(nc.variables["lat"][:])
169      if typefile in ['gcm','earthgcm'] and is1d == 1:       mapmode=0 ; winds=False
[579]170    ### we create default vert and time prescriptions if not here in case mapping mode is on (lat/lon maps)
[548]171      if redope is None and mapmode == 1:
[477]172          if svert is None:  svert = readslices(str(level)) ; nvert=1
173          if stime is None and mrate is None:
174             stime = readslices(str(0)) ; ntime=1 ## this is a default choice
175             print "WELL... nothing about time axis. I took default: first time reference stored in file."
[579]176    ### we get the names of variables to be read. in case only one available, we choose this one.
[725]177    ### (we take out of the test specific names e.g. UV is not in the file but used to ask a wind speed computation)
[760]178      varname = select_getfield(zvarname=var[vvv],znc=nc,ztypefile=typefile,mode='check')
[579]179    ### we get the names of wind variables to be read (if any)
[349]180      if winds:                                                   
[345]181         [uchar,vchar,metwind] = getwinddef(nc)             
182         if uchar == 'not found': winds = False
[579]183    ### we tell the user that either no var or no wind is not acceptable
[392]184      if not varname and not winds: errormess("please set at least winds or var",printvar=nc.variables)
[579]185    ### we get the coordinates lat/lon to be used
[399]186      [lon2d,lat2d] = getcoorddef(nc)
[579]187    ### we get an adapted map projection if none is provided by the user
[548]188      if proj == None:   proj = getproj(nc)                 
[579]189    ### we define plot boundaries given projection or user choices
[489]190      if firstfile:
191         if proj in ["npstere","spstere"]: [wlon,wlat] = polarinterv(lon2d,lat2d)
192         elif proj in ["lcc","laea"]:      [wlon,wlat] = wrfinterv(lon2d,lat2d)
193         else:                             [wlon,wlat] = simplinterv(lon2d,lat2d)
194         if zoom:                          [wlon,wlat] = zoomset(wlon,wlat,zoom)
[558]195         elif zarea is not None:           [wlon,wlat] = latinterv(area=zarea)
[489]196
[725]197#############################################################
198############ WE LOAD 4D DIMENSIONS : x, y, z, t #############
199#############################################################
200
[781]201      ###should be available for GCM or simple files ?
202      ##if ope in ["cat"] and nnn > 0:    count = time[-1] + 1  ## so that a cat is possible with simple subscripts
203
[725]204    ### TYPE 1 : GCM files or simple files
[558]205      if typefile in ["gcm","earthgcm","ecmwf"]:
[725]206      ### this is needed for continuity
[610]207          if slon is not None: sslon = slon 
208          if slat is not None: sslat = slat 
[725]209      ### we define lat/lon vectors. we get what was done in getcoorddef.
[735]210          lat = lat2d[:,0] ; lon = lon2d[0,:]
[725]211      ### we define areas. this is needed for calculate means and weight with area. this is not compulsory (see reduce_field).
212          if "aire" in nc.variables:      area = nc.variables["aire"][:,:]
213          ### --> add a line here if your reference is not present
214          else:                           area = None
215      ### we define altitude vector. either it is referenced or it is guessed based on last variable's dimensions.
[721]216          if "altitude" in nc.variables:   vert = nc.variables["altitude"][:]
217          elif "Alt" in nc.variables:      vert = nc.variables["Alt"][:]
218          elif "lev" in nc.variables:      vert = nc.variables["lev"][:]
219          elif "presnivs" in nc.variables: vert = nc.variables["presnivs"][:]
[725]220          ### --> add a line here if your reference is not present
[724]221          else: 
[725]222              print "No altitude found. Try to build a simple axis." ; dadim = getdimfromvar(nc) 
[724]223              if   len(dadim) == 4:  print "-- 4D field. Assume z is dim 2." ; vert = np.arange(dadim[-3])
224              elif len(dadim) == 3:  print "-- 3D field. Assume z is dim 1." ; vert = [0.]
[725]225              else:                  vert = [0.]
226      ### we define time vector. either it is referenced or it is guessed based on last variable's dimensions.
[525]227          if "Time" in nc.variables:            time = nc.variables["Time"][:]
[748]228          elif "time_counter" in nc.variables:  time = nc.variables["time_counter"][:]/86400. #### convert from s to days
[525]229          elif "time" in nc.variables:          time = nc.variables["time"][:]
[725]230          ### --> add a line here if your reference is not present
[724]231          else: 
232              print "No time found. Try to build a simple axis. Assume t is dim 1." 
233              dadim = getdimfromvar(nc)
234              if   len(dadim) == 4:  time = np.arange(dadim[-4])
235              elif len(dadim) == 3:  time = np.arange(dadim[-3])
[725]236              else:                  time = [0.] #errormess("no time axis found.")
237          ### (SPECIFIC. convert to Ls for Martian GCM files.)
[640]238          if axtime in ["ls"]:
239              print "converting to Ls ..."
240              for iii in range(len(time)):
[673]241                time[iii] = sol2ls(time[iii])
[640]242                if iii > 0:
[725]243                  while abs(time[iii]-time[iii-1]) > 300: time[iii] = time[iii]+360
244          ### (a case where the user would like to set local times. e.g. : 1D plot with no longitude reference)
[494]245          if axtime in ["lt"]:
246              if initime == -1: initime=input("Please type initial local time:")
[725]247              time = (initime+time*24)%24 ; print "LOCAL TIMES.... ", time
248
249    ### TYPE 2 : MESOSCALE FILES
[548]250      elif typefile in ['meso','geo']:
[725]251          ### area not active with mesoscale files
252          area = None 
253          ### HACK TO GET THE CORRECT LAT/LON FROM MESOSCALE FILES WITH 2D LAT/LON ARRAYS
254          ### principle: calculate correct indices then repopulate slon and slat
[490]255          if slon is not None or slat is not None:
[725]256              show_topo_map_user = ((save == 'png') and (typefile == 'meso') and ("HGT" in varinfile) and display)
257              if firstfile and show_topo_map_user:  iwantawhereplot = nc     #show a topo map with a cross on the chosen point
258              else:                                 iwantawhereplot = None   #do not show anything, just select indices
[490]259              numlon = 1 ; numlat = 1 
[610]260              if slon is not None:   numlon = slon.shape[1]   
261              if slat is not None:   numlat = slat.shape[1]
[490]262              indices = np.ones([numlon,numlat,2]) ; vlon = None ; vlat = None
263              for iii in range(numlon): 
264               for jjj in range(numlat):
[610]265                 if slon is not None:  vlon = slon[0][iii]  ### note: slon[:][0] does not work
266                 if slat is not None:  vlat = slat[0][jjj]  ### note: slon[:][0] does not work
[490]267                 indices[iii,jjj,:] = bidimfind(lon2d,lat2d,vlon,vlat,file=iwantawhereplot) 
268                 lonp,latp = ( lon2d[indices[iii,jjj,0],indices[iii,jjj,1]] , lat2d[indices[iii,jjj,0],indices[iii,jjj,1]] )
[587]269              ### possible bug here if several --lat
[490]270              for iii in range(numlon):
271               for jjj in range(numlat):
[610]272                 if slon is not None: sslon[0][iii] = indices[iii,0,1] #...this is idx
273                 if slat is not None: sslat[0][jjj] = indices[0,jjj,0] #...this is idy
[490]274              lonp,latp = ( lon2d[indices[0,0,0],indices[0,0,1]] , lat2d[indices[0,0,0],indices[0,0,1]] )
[725]275          ### we get rid of boundary relaxation zone for plots. important to do that now and not before.
276          if typefile in ['meso'] and mapmode == 1: lon2d = dumpbdy(lon2d,6) ; lat2d = dumpbdy(lat2d,6) 
277          ### we read the keyword for vertical dimension. we take care for vertical staggering.
[393]278          if varname in ['PHTOT','W']:    vertdim='BOTTOM-TOP_PATCH_END_STAG'
279          else:                           vertdim='BOTTOM-TOP_PATCH_END_UNSTAG'
[490]280          if (var2 is not None and var2 not in ['PHTOT','W']): dumped_vert_stag=True ; vertdim='BOTTOM-TOP_PATCH_END_UNSTAG'
281          else:                                                dumped_vert_stag=False
[725]282          ### we read the keyword for horizontal dimensions. we take care for horizontal staggering.
[393]283          if varname in ['V']:  latdim='SOUTH-NORTH_PATCH_END_STAG'
284          else:                 latdim='SOUTH-NORTH_PATCH_END_UNSTAG'
285          if varname in ['U']:  londim='WEST-EAST_PATCH_END_STAG'
286          else:                 londim='WEST-EAST_PATCH_END_UNSTAG'
[402]287          lon = np.arange(0,getattr(nc,londim),1) ; lat = np.arange(0,getattr(nc,latdim),1)
[725]288          ### we define the time axis and take care of various specificities (lt, ls, sol) or issues (concatenation)
[468]289          if axtime in ["ls","sol"]:
290              lstab, soltab, lttab = getlschar ( namefile, getaxis = True )
291              if axtime == "ls":      time = lstab
292              elif axtime == "sol":   time = soltab
293          else:
[665]294              if ope in ["cat"] and nnn > 0:    count = time[-1] + 1  ## so that a cat is possible with simple subscripts
295              else:                             count = 0
[613]296              if "Times" in nc.dimensions:   time = count + np.arange(0,len(nc.dimensions["Times"]),1)
297              elif "Time" in nc.dimensions:  time = count + np.arange(0,len(nc.dimensions["Time"]),1)
298              else:                          time = count + np.arange(0,1,1)
[489]299          if axtime in ["lt"]:
[569]300              ftime = np.zeros(len(time))
301              for i in range(len(time)): 
[763]302                 ftime[i] = localtime ( time[i], slon , nc )
[569]303              time=ftime
304              time=check_localtime(time)
[489]305              print "LOCAL TIMES.... ", time
[725]306          ### we define the vertical axis (or lack thereof) and cover possibilities for it to be altitude, pressure, geopotential. quite SPECIFIC.
[426]307          if typefile in ['geo']:   vert = [0.] ; stime = readslices(str(0))
308          else:
309              if vertmode is None:  vertmode=0
[687]310              if vertmode == 0:     
311                  if "vert" in nc.variables: vert = nc.variables["vert"][:]/1000. ; vertmode = 1
312                  else:                      vert = np.arange(0,getattr(nc,vertdim),1)
[562]313              elif vertmode == -1:  vert = nc.variables["PHTOT"][0,:,0,0]/3.72 ; vert = np.array(vert[0:len(vert)-1]) #; print vert
[610]314              elif vertmode == 1 or vertmode == 2:  vert = nc.variables["vert"][:]        ## pressure in Pa
315              else:                                 vert = nc.variables["vert"][:]/1000.  ## altitude in km
[725]316####################################################################
317############ END of WE LOAD 4D DIMENSIONS : x, y, z, t #############
318####################################################################
[345]319
[760]320    ### we fill the arrays of varname, namefile, time, colorbar at the current step considered (NB: why use both k and nnn ?)
[402]321      all_varname[k] = varname
322      all_namefile[k] = namefile
[406]323      all_time[k] = time
[763]324      all_vert[k] = vert
325      all_lat[k] = lat
326      all_lon[k] =  lon
[760]327      all_colorb[k] = clb[vvv]
[763]328      if var2:  all_var2[k], plot_x[k], plot_y[k] = select_getfield(zvarname=var2,znc=nc,ztypefile=typefile,mode='getvar',ztsat=tsat,ylon=all_lon[k],ylat=all_lat[k],yalt=all_vert[k],ytime=all_time[k],analysis=analysis)
[451]329      if winds: all_windu[k] = getfield(nc,uchar) ; all_windv[k] = getfield(nc,vchar)
[725]330    ### we fill the arrays of fields to be plotted at the current step considered
[763]331      all_var[k], plot_x[k], plot_y[k] = select_getfield(zvarname=all_varname[k],znc=nc,ztypefile=typefile,mode='getvar',ztsat=tsat,ylon=all_lon[k],ylat=all_lat[k],yalt=all_vert[k],ytime=all_time[k],analysis=analysis)
332
333      # last line of for namefile in namefiles
[468]334      print "**** GOT SUBDATA:",k," NAMEFILE:",namefile," VAR:",varname, var2 ; k += 1 ; firstfile = False
[345]335
[725]336    ### note to self : stopped comment rewriting here.
337
[359]338    ##################################
[763]339    ### Operation on files (I) with _var
340    if ope is not None and "var" in ope:
341         print "********** OPERATION: ",ope
342         if len(namefiles) > 1: errormess("for this operation... please set only one file !")
343         if len(var) > 2:       errormess("not sure this works for more than 2 vars... please check.")
344         if   "div_var" in ope: all_var[k] = all_var[k-2] / all_var[k-1] ; insert = '_div_'
345         elif "mul_var" in ope: all_var[k] = all_var[k-2] * all_var[k-1] ; insert = '_mul_'
346         elif "add_var" in ope: all_var[k] = all_var[k-2] + all_var[k-1] ; insert = '_add_'
347         elif "sub_var" in ope: all_var[k] = all_var[k-2] - all_var[k-1] ; insert = '_sub_'
348         else:                    errormess(ope+" : non-implemented operation. Check pp.py --help")
349         plot_x[k] = None ; plot_y[k] = None
350         all_time[k] = all_time[k-1] ; all_vert[k] = all_vert[k-1] ; all_lat[k] = all_lat[k-1] ; all_lon[k] = all_lon[k-1] ; all_namefile[k] = all_namefile[k-1]
351         all_varname[k] = all_varname[k-2] + insert + all_varname[k-1]
352         if len(clb) >= zelen: all_colorb[k] = clb[-1]   # last additional user-defined color is for operation plot
353         else:                 all_colorb[k] = all_colorb[k-1]  # if no additional user-defined color... set same as var
354         ### only the operation plot. do not mention colorb so that it is user-defined?
355         if "only" in ope:
356             numplot = 1 ; all_var[0] = all_var[k]
357             all_time[0] = all_time[k] ; all_vert[0] = all_vert[k] ; all_lat[0] = all_lat[k] ; all_lon[0] = all_lon[k] ; all_namefile[0] = all_namefile[k]
358             all_varname[0] = all_varname[k-2] + insert + all_varname[k-1]
[721]359
360
[345]361    ##################################
[763]362    ### Operation on files (II) without _var
363    # we re-iterate on the plots to set operation subplots to make it compatible with multifile (using the same ref file)
364    # (k+1)%3==0 is the index of operation plots
365    # (k+2)%3==0 is the index of reference plots
366    # (k+3)%3==0 is the index of first plots
367    opefirstpass=True
368    if ope is not None and "var" not in ope:
369       print "********** OPERATION: ",ope
370       for k in np.arange(zelen):
371               if len(var) > 1: errormess("for this operation... please set only one var !")
372               if ope in ["-","+","-%","-_only","+_only","-%_only","-_histo"]:
373                  if fileref is None: errormess("fileref is missing!")
374                  else:ncref = Dataset(fileref)
375
376                  if opefirstpass: ## first plots
377                     for ll in np.arange(len(namefiles)):
378                        print "SETTING FIRST PLOT"
379                        all_varname[3*ll] = all_varname[ll] ; all_time[3*ll] = all_time[ll] ; all_vert[3*ll] = all_vert[ll] ; all_lat[3*ll] = all_lat[ll] ; all_lon[3*ll] = all_lon[ll] ; all_namefile[3*ll] = all_namefile[ll] ; all_var2[3*ll] = all_var2[ll] ; all_colorb[3*ll] = all_colorb[ll] ; all_var[3*ll] = all_var[ll]
380                        if plot_y[ll] is not None: plot_y[3*ll] = plot_y[ll] ; plot_x[3*ll] = plot_x[ll]
381                        else: plot_y[3*ll] = None ; plot_x[3*ll] = None
382                        opefirstpass=False
383
384                  if (k+2)%3==0: ## reference plots
385                        print "SETTING REFERENCE PLOT"
386                        all_varname[k] = all_varname[k-1] ; all_time[k] = all_time[k-1] ; all_vert[k] = all_vert[k-1] ; all_lat[k] = all_lat[k-1] ; all_lon[k] = all_lon[k-1] ; all_namefile[k] = all_namefile[k-1] ; all_var2[k] = all_var2[k-1] ; all_colorb[k] = all_colorb[k-1]
387                        all_var[k], plot_x[k], plot_y[k] = select_getfield(zvarname=all_varname[k-1],znc=ncref,ztypefile=typefile,mode='getvar',ztsat=tsat,ylon=all_lon[k],ylat=all_lat[k],yalt=all_vert[k],ytime=all_time[k],analysis=analysis)
388                        if winds: all_windu[k] = getfield(ncref,uchar) ; all_windv[k] = getfield(ncref,vchar)
389
390                  if (k+1)%3==0: ## operation plots
391                     print "SETTING OPERATION PLOT"
392                     all_varname[k] = all_varname[k-1] ; all_time[k] = all_time[k-1] ; all_vert[k] = all_vert[k-1] ; all_lat[k] = all_lat[k-1] ; all_lon[k] = all_lon[k-1] ; all_namefile[k] = all_namefile[k-1] ; all_var2[k] = all_var2[k-1]
393                     if ope in ["-","-_only","-_histo"]:
394                         all_var[k]= all_var[k-2] - all_var[k-1]
395                         if plot_y[k-1] is not None and plot_y[k-2] is not None: plot_y[k] = plot_y[k-2] - plot_y[k-1]
396                         if plot_y[k-2] is None: plot_y[k] = None; plot_x[k] = None
397                     elif ope in ["+","+_only"]:   
398                         all_var[k]= all_var[k-2] + all_var[k-1]
399                         if plot_y[k-1] is not None and plot_y[k-2] is not None: plot_y[k] = plot_y[k-2] + plot_y[k-1]
400                         if plot_y[k-2] is None: plot_y[k] = None; plot_x[k] = None
401                     elif ope in ["-%","-%_only"]:
402                         masked = np.ma.masked_where(all_var[k-1] == 0,all_var[k-1])
403                         masked.set_fill_value([np.NaN])
404                         all_var[k]= 100.*(all_var[k-2] - masked)/masked
405                         if plot_y[k-1] is not None and plot_y[k-2] is not None: 
406                            masked = np.ma.masked_where(plot_y[k-1] == 0,plot_y[k-1])
407                            masked.set_fill_value([np.NaN])
408                            plot_y[k]= 100.*(plot_y[k-2] - masked)/masked
409                         if plot_y[k-2] is None: plot_y[k] = None; plot_x[k] = None
410                     if len(clb) >= zelen: all_colorb[k] = clb[-1]
411                     else: all_colorb[k] = "RdBu_r" # if no additional user-defined color... set a good default one
412                     if winds: all_windu[k] = all_windu[k-2]-all_windu[k-1] ; all_windv[k] = all_windv[k-2] - all_windv[k-1]
413
414               elif ope in ["cat"]:
415                  tabtime = all_time[0];tab = all_var[0];k = 1
416                  if var2: tab2 = all_var2[0]
417                  while k != len(namefiles) and len(all_time[k]) != 0:
418                      if var2: tab2 = np.append(tab2,all_var2[k],axis=0) 
419                      tabtime = np.append(tabtime,all_time[k]) ; tab = np.append(tab,all_var[k],axis=0) ; k += 1
420                  all_time[0] = np.array(tabtime) ; all_var[0] = np.array(tab) ; numplot = 1
421                  if var2: all_var2[0] = np.array(tab2)
422               else: errormess(ope+" : non-implemented operation. Check pp.py --help")
423       if "only" in ope:
424           numplot = 1 ; all_var[0] = all_var[k]
425           all_time[0] = all_time[k] ; all_vert[0] = all_vert[k] ; all_lat[0] = all_lat[k] ; all_lon[0] = all_lon[k] ; all_namefile[0] = all_namefile[k] ; plot_x[0]=plot_x[k] ; plot_y[0]=plot_y[k]
426           all_varname[0] = all_varname[k]
427
428    ##################################
[345]429    ### Open a figure and set subplots
430    fig = figure()
431    subv,subh = definesubplot( numplot, fig ) 
[763]432    if ope in ['-','-%','-_histo'] and len(namefiles) ==1 : subv,subh = 2,2
433    elif ope in ['-','-%'] and len(namefiles)>1 : subv, subh = len(namefiles), 3
[345]434 
435    #################################
436    ### Time loop for plotting device
[763]437    nplot = 1 ; error = False ; firstpass = True 
[612]438    if lstyle is not None: linecycler = cycle(lstyle)
439    else: linecycler = cycle(["-","--","-.",":"])
[392]440    print "********************************************"
[345]441    while error is False:
[458]442     
[763]443       print "********** PLOT", nplot, " OF ",numplot
[424]444       if nplot > numplot: break
[345]445
[424]446       ####################################################################
[392]447       ## get all indexes to be taken into account for this subplot and then reduce field
448       ## 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
[423]449       if ope is not None:
[763]450           if fileref is not None:      index_f = ((nplot-1)//(nlon*nlat*nvert*ntime))%(3*len(namefiles))  ## OK only 1 var,  see test in the beginning
[423]451           elif "var" in ope:           index_f = ((nplot-1)//(nlon*nlat*nvert*ntime))%(len(var)+1)        ## OK only 1 file, see test in the beginning
[451]452           elif "cat" in ope:           index_f = 0
[763]453       elif not firstpass:
454          if len(namefiles) > 1 and len(var) == 1 and which == "unidim": pass
455          else: yeah = len(namefiles)*len(var) ; index_f = ((nplot-1)//(nlon*nlat*nvert*ntime))%yeah
456       else: yeah = len(namefiles)*len(var) ; index_f = ((nplot-1)//(nlon*nlat*nvert*ntime))%yeah
457       time = all_time[index_f] ; vert = all_vert[index_f] ; lat = all_lat[index_f] ; lon = all_lon[index_f]
458       indexlon  = getsindex(sslon,(nplot-1)%nlon,lon)
459       indexlat  = getsindex(sslat,((nplot-1)//nlon)%nlat,lat)
460       indexvert = getsindex(svert,((nplot-1)//(nlon*nlat))%nvert,vert)
461       plotx=plot_x[index_f] ; ploty=plot_y[index_f]
[435]462       if mrate is not None:                 indextime = None 
463       else:                                 indextime = getsindex(stime,((nplot-1)//(nlon*nlat*nvert))%ntime,time)
[426]464       ltst = None 
[763]465       if typefile in ['meso'] and indextime is not None and len(indextime) < 2: ltst = localtime ( indextime, slon , nc)
[468]466       print "********** INDEX LON:",indexlon," LAT:",indexlat," VERT:",indexvert," TIME:",indextime
[490]467       ##var = nc.variables["phisinit"][:,:]
468       ##contourf(np.transpose(var),30,cmap = get_cmap(name="Greys_r") ) ; axis('off') ; plot(indexlat,indexlon,'mx',mew=4.0,ms=20.0)
469       ##show()
470       ##exit()
471       #truc = True
472       #truc = False
473       #if truc: indexvert = None
[424]474       ####################################################################
[475]475       ########## REDUCE FIELDS
476       ####################################################################
[448]477       error = False
[392]478       varname = all_varname[index_f]
[763]479       which=''
[448]480       if varname:   ### what is shaded.
[436]481           what_I_plot, error = reducefield( all_var[index_f], d4=indextime, d1=indexlon, d2=indexlat, d3=indexvert, \
[717]482                                             yint=yintegral, alt=vert, anomaly=anomaly, redope=redope, mesharea=area, unidim=is1d)
[516]483           if mult != 2718.:  what_I_plot = what_I_plot*mult
484           else:              what_I_plot = np.log10(what_I_plot) ; print "log plot"
[760]485                 
[451]486       if var2:      ### what is contoured.
[448]487           what_I_plot_contour, error = reducefield( all_var2[index_f], d4=indextime, d1=indexlon, d2=indexlat , d3=indexvert, \
[760]488                                                     yint=yintegral, alt=vert )
[451]489       if winds:     ### what is plot as vectors.
490           vecx, error = reducefield( all_windu[index_f], d4=indextime, d3=indexvert, yint=yintegral, alt=vert)
491           vecy, error = reducefield( all_windv[index_f], d4=indextime, d3=indexvert, yint=yintegral, alt=vert)
[490]492           if varname in [uchar,vchar]: what_I_plot = np.sqrt( np.square(vecx) + np.square(vecy) ) ; varname = "wind"
[763]493 
494       if plotx is not None:
495          plotx, error = reducefield( plotx, d4=indextime, d1=indexlon, d2=indexlat, d3=indexvert, \
496                                             yint=yintegral, alt=vert, anomaly=anomaly, redope=redope, mesharea=area, unidim=is1d)
497          ploty, error = reducefield( ploty, d4=indextime, d1=indexlon, d2=indexlat, d3=indexvert, \
498                                             yint=yintegral, alt=vert, anomaly=anomaly, redope=redope, mesharea=area, unidim=is1d)
499          which='xy'
[490]500       #####################################################################
501       #if truc:
502       #   nx = what_I_plot.shape[2] ; ny = what_I_plot.shape[1] ; nz = what_I_plot.shape[0]
503       #   for k in range(nz): print k,' over ',nz ; what_I_plot[k,:,:] = what_I_plot[k,:,:] / smooth(what_I_plot[k,:,:],12)
504       #   for iii in range(nx):
505       #    for jjj in range(ny):
506       #     deviation = what_I_plot[:,jjj,iii] ; mx = max(deviation) ; mn = min(deviation)
507       #     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
508       #     else:                                                   what_I_plot[0,jjj,iii]     = 0.
509       #     if np.abs(what_I_plot[0,jjj,iii]) > 1.5:
510       #         print iii,jjj,what_I_plot[0,jjj,iii],int(abs(1.-mx)*100.),int(abs(1.-mn)*100.)
511       #         plot(rel)
512       #   show()
513       #   anomaly = True ### pour avoir les bons reglages plots
514       #   what_I_plot = what_I_plot[0,:,:] 
515       #####################################################################
516
[448]517       ####################################################################
[458]518       ### General plot settings
[763]519       changesubplot = (numplot > 1) and (len(what_I_plot.shape) != 1) and (which != "xy")  ## default for 1D plots: superimposed. to be reworked for better flexibility.
[518]520       if changesubplot: subplot(subv,subh,nplot) #; subplots_adjust(wspace=0,hspace=0)
[760]521       colorb = all_colorb[index_f]
[458]522       ####################################################################
[475]523       if error:
524               errormess("There is an error in reducing field !")
525       else:
[448]526               ticks = ndiv + 1 
[405]527               fvar = varname
528               if anomaly: fvar = 'anomaly'
[475]529               ###
530               if mapmode == 0:    ### could this be moved inside imov loop ?
[430]531                   itime=indextime
[587]532                   if len(what_I_plot.shape) == 3: itime=[0]
[475]533                   m = None ; x = None ; y = None
[608]534                   latyeah = lat ; lonyeah = lon
535                   if typefile in ['meso']:
536                       # now that the section is determined we can set the real lat
537                       # ... or for now, a temptative one.
538                       if slon is not None: latyeah = lat2d[:,0]
539                       if slat is not None: lonyeah = lon2d[0,:]
540                   what_I_plot, x, y = define_axis(lonyeah,latyeah,vert,time,indexlon,indexlat,indexvert,\
[763]541                         itime,what_I_plot, len(all_var[index_f].shape),vertmode,redope)
[475]542               ###
[763]543               if (fileref is not None) and ((index_f+1)%3 == 0):    zevmin, zevmax = calculate_bounds(what_I_plot,vmin=minop,vmax=maxop)
[405]544               else:                                                   zevmin, zevmax = calculate_bounds(what_I_plot,vmin=vmin,vmax=vmax)
[760]545               #if (fileref is not None) and (index_f == numplot-1):    colorb = "RdBu_r"
[518]546               if colorb in ["def","nobar","onebar"]:                  palette = get_cmap(name=defcolorb(fvar.upper()))
[436]547               else:                                                   palette = get_cmap(name=colorb)
[490]548               #palette = cm.GMT_split
[489]549               ##### 1. ELIMINATE 0D or >3D CASES
550               if len(what_I_plot.shape) == 0:   
551                 print "VALUE VALUE VALUE VALUE ::: ",what_I_plot
552                 save = 'donothing'
553               elif len(what_I_plot.shape) >= 4:
[475]554                 print "WARNING!!! ",len(what_I_plot.shape),"-D PLOT NOT SUPPORTED !!! dimensions: ",what_I_plot.shape
555                 errormess("Are you sure you did not forget to prescribe a dimension ?")
556               ##### 2. HANDLE simple 1D/2D field and movies of 1D/2D fields
557               else:
[434]558                 if mrate is not None: iend=len(time)-1
[448]559                 else:                 iend=0
[763]560                 imov = 0
561                 if analysis in ['density','histo','fft','histodensity']: which="xy"
562                 elif len(what_I_plot.shape) == 3:
563                    if var2 and which == '':               which = "contour" ## have to start with contours rather than shading
564                    elif which == '':                      which = "regular"
[475]565                    if mrate is None:      errormess("3D field. Use --rate RATE for movie or specify --time TIME. Exit.")
566                 elif len(what_I_plot.shape) == 2:
[763]567                    if var2 and which == '':               which = "contour" ## have to start with contours rather than shading
568                    elif which == '':                      which = "regular"
569                    if mrate is not None and which == '':  which = "unidim"
570                 elif len(what_I_plot.shape) == 1 and which == '' :
[475]571                    which = "unidim"
[562]572                    if what_I_plot.shape[-1] == 1:      print "VALUE VALUE VALUE VALUE ::: ", what_I_plot[0] ; save = 'donothing'
[763]573##                 if which == "unidim" and len(namefiles) > 1: numplot = 1 # this case is similar to several vars from one file
[475]574                 ##### IMOV LOOP #### IMOV LOOP
[430]575                 while imov <= iend:
[448]576                    print "-> frame ",imov+1, which
577                    if which == "regular":   
578                        if mrate is None:                                   what_I_plot_frame = what_I_plot
579                        else:                                               what_I_plot_frame = what_I_plot[imov,:,:]
[451]580                        if winds:
581                            if mrate is None:                                   vecx_frame = vecx ; vecy_frame = vecy
582                            else:                                               vecx_frame = vecx[imov,:,:] ; vecy_frame = vecy[imov,:,:]
[448]583                    elif which == "contour": 
584                        if mrate is None or what_I_plot_contour.ndim < 3:   what_I_plot_frame = what_I_plot_contour
585                        else:                                               what_I_plot_frame = what_I_plot_contour[imov,:,:]
[475]586                    elif which == "unidim":
587                        if mrate is None:                                   what_I_plot_frame = what_I_plot
588                        else:                                               what_I_plot_frame = what_I_plot[:,imov]  ## because swapaxes
589                    #if mrate is not None:     
590                    if mapmode == 1: 
[763]591                        m = define_proj(proj,wlon,wlat,back=back,blat=blat,blon=blon)  ## this is dirty, defined above but out of imov loop
592                        x, y = m(lon2d, lat2d)                                         ## this is dirty, defined above but out of imov loop
[548]593                    if typefile in ['meso'] and mapmode == 1:   what_I_plot_frame = dumpbdy(what_I_plot_frame,6,condition=True)
[464]594#                   if typefile in ['mesoideal']:    what_I_plot_frame = dumpbdy(what_I_plot_frame,0,stag='W',condition=dumped_vert_stag)
[444]595
[475]596                    if which == "unidim":
[763]597#                        if lbls is not None: lbl=lbls[nplot-1]
598                        if lbls is not None: lbl=lbls[index_f]
[612]599                        else:
600                           lbl = ""
601                           if indexlat is not None:  lbl = lbl + " ix" + str(indexlat[0])
602                           if indexlon is not None:  lbl = lbl + " iy" + str(indexlon[0])
603                           if indexvert is not None: lbl = lbl + " iz" + str(indexvert[0])
604                           if indextime is not None: lbl = lbl + " it" + str(indextime[0])
[647]605                           if lbl == "": lbl = all_namefile[index_f]
[569]606
[475]607                        if mrate is not None: x = y  ## because swapaxes...
[490]608                        #what_I_plot_frame = np.diff(what_I_plot_frame, n=1) ; x = x[1:]
[587]609                     
610                        zeline = next(linecycler) ## "-" for simple lines
611                        if tile:      zemarker = 'x'
612                        else:         zemarker = None 
[579]613                        this_is_a_regular_plot = (indexvert is not None) or (indextime is None) or (indexlat is None) or (indexlon is None)
[587]614                        if this_is_a_regular_plot:   plot(x,what_I_plot_frame,zeline,label=lbl,marker=zemarker)  ## vertical profile
615                        else:                        plot(what_I_plot_frame,x,zeline,label=lbl,marker=zemarker)  ## regular plot
[475]616                        if nplot > 1: legend(loc='best')
[502]617                        if indextime is None and axtime is not None and xlab is None:    xlabel(axtime.upper()) ## define the right label
[642]618                        if save == 'txt':  writeascii(np.transpose([x,np.transpose(what_I_plot)]),'profile'+str(nplot*1000+imov)+'.txt')
[475]619
[763]620                    elif which == "xy":
621                        if lbls is not None: lbl=lbls[index_f]
622                        else: lbl=None
623                        if analysis is not None:
624                           if analysis == 'histo': 
625                               if zelen == 1:
626                                  mpl.pyplot.hist(ploty.flatten(),bins=ndiv,normed=True, alpha=0.67, facecolor = 'green', label = lbls)
627                                  mpl.pyplot.legend()
628                                  mpl.pyplot.grid(True)
629                                  mpl.pyplot.title(zetitle)
630                               else:
631                                  multiplot[index_f]=ploty.flatten()
632                                  if index_f == zelen-1: 
633                                     if ope is not None: multiplot = getVar(multiplot,3*np.arange((index_f+1)/3)+2) ## we only compute histograms for the operation plots
634                                     mpl.pyplot.hist(multiplot,bins=ndiv,normed=True, alpha=0.75, label = lbls) 
635                                     mpl.pyplot.legend()
636                                     mpl.pyplot.grid(True)
637                                     mpl.pyplot.title(zetitle)
638                           elif analysis in ['density','histodensity']:
639                                  if ope is not None and (index_f+1)%3 !=0: pass
640                                  else:
641                                     plotx = np.linspace(min(ploty.flatten()),max(ploty.flatten()),1000)
642                                     density = gaussian_kde(ploty.flatten())
643   #                                  density.covariance_factor = lambda : .25  # adjust the covariance factor to change the bandwidth if needed
644   #                                  density._compute_covariance()
645                                        # display the mean and variance of the kde:
646                                     sample = density.resample(size=20000)
647                                     n, (smin, smax), sm, sv, ss, sk = describe(sample[0])
648                                     mpl.pyplot.plot(plotx,density(plotx), label = lbl+'\nmean: '+str(sm)[0:5]+'   std: '+str(np.sqrt(sv))[0:5]+'\nskewness: '+str(ss)[0:5]+'   kurtosis: '+str(sk)[0:5])
649                                     if analysis == 'histodensity':  # plot both histo and density (to assess the rightness of the kernel density estimate for exemple) and display the estimated variance
650                                        mpl.pyplot.hist(ploty.flatten(),bins=ndiv,normed=True, alpha=0.30, label = lbl)
651                                     if index_f == zelen-1: mpl.pyplot.legend() ; mpl.pyplot.title(zetitle)
652                           else:
653                              plot(plotx,ploty,label = lbl)
654                              if index_f == zelen-1: mpl.pyplot.legend() ; mpl.pyplot.title(zetitle)
655                              mpl.pyplot.grid(True)
656                        else:
657                           plot(plotx,ploty,label = lbl)
658                           if index_f == zelen-1: mpl.pyplot.legend() ; mpl.pyplot.title(zetitle)
659                        if varname == 'hodograph':
660                            a=0
661                            for ii in np.arange(len(time)): 
662                               if a%6 == 0: mpl.pyplot.text(plotx[ii],ploty[ii],time[ii]) 
663                               a=a+1
664                            mpl.pyplot.grid(True)
665
[569]666                    elif which == "regular":
[647]667                   
668                        # plot stream lines if there is a stream file and a vert/lat slice. Might not work with movies ??
669                        if streamflag and sslat is None and svert is None:
670                             streamfile = all_namefile[index_f].replace('.nc','_stream.nc')
671                             teststream = os.path.exists(streamfile)
672                             if teststream:
673                                print 'INFO: Using stream file',streamfile, 'for stream lines'
674                                ncstream = Dataset(streamfile)
675                                psi = getfield(ncstream,'psi')
676                                psi = psi[0,:,:,0] # time and longitude are dummy dimensions
677                                if psi.shape[1] != len(x) or psi.shape[0] != len(y):
678                                    errormess('stream file does not fit! Dimensions: '+str(psi.shape)+' '+str(x.shape)+' '+str(y.shape))
679                                zelevels = np.arange(-1.e10,1.e10,1.e9)
680                                zemin = np.min(abs(zelevels))
681                                zemax = np.max(abs(zelevels))
682                                zewidth  =  (abs(zelevels)-zemin)*(5.- 0.5)/(zemax - zemin) + 0.5 # linewidth ranges from 5 to 0.5
683                                cs = contour( x,y,psi, zelevels, colors='k', linewidths = zewidth)
684                                clabel(cs, inline=True, fontsize = 4.*rcParams['font.size']/5., fmt="%1.1e")
685                             else:
686                                print 'WARNING: STREAM FILE',streamfile, 'DOES NOT EXIST !'
687                             
[448]688                        if hole:         what_I_plot_frame = hole_bounds(what_I_plot_frame,zevmin,zevmax)
689                        else:            what_I_plot_frame = bounds(what_I_plot_frame,zevmin,zevmax)
690                        if flagnolow:    what_I_plot_frame = nolow(what_I_plot_frame)
691                        if not tile:
692                            #zelevels = np.linspace(zevmin*(1. + 1.e-7),zevmax*(1. - 1.e-7)) #,num=20)
693                            zelevels = np.linspace(zevmin,zevmax,num=ticks)
[548]694                            #what_I_plot_frame = smooth(what_I_plot_frame,100)
[453]695                            if mapmode == 1:       m.contourf( x, y, what_I_plot_frame, zelevels, cmap = palette, alpha=trans)
696                            elif mapmode == 0:     contourf( x, y, what_I_plot_frame, zelevels, cmap = palette, alpha=trans)
[612]697                            if cross is not None and mapmode == 1:
698                               idx,idy=m(cross[0][0],cross[0][1])
699                               mpl.pyplot.plot([idx],[idy],'+k',mew=0.5,ms=18)
[448]700                        else:
[458]701                            if mapmode == 1:       m.pcolor( x, y, what_I_plot_frame, cmap = palette, vmin=zevmin, vmax=zevmax, alpha=trans)
702                            elif mapmode == 0:     pcolor( x, y, what_I_plot_frame, cmap = palette, vmin=zevmin, vmax=zevmax, alpha=trans)
[647]703                       
[475]704
[518]705                        if colorb not in ['nobar','onebar']:       
[763]706                            if (fileref is not None) and ((index_f+1)%3 == 0):   daformat = "%.3f" 
[468]707                            elif mult != 1:                                        daformat = "%.1f"
[448]708                            else:                                                  daformat = fmtvar(fvar.upper())
[603]709                            #if proj in ['moll','cyl']:  zeorientation="horizontal" ; zepad = 0.05
710                            if proj in ['moll']:        zeorientation="horizontal" ; zepad = 0.05
[599]711                            else:                       zeorientation="vertical" ; zepad = 0.03
712                            zecb = colorbar( fraction=0.05,pad=zepad,format=daformat,orientation=zeorientation,\
[468]713                                      ticks=np.linspace(zevmin,zevmax,num=min([ticks/2+1,21])),extend='neither',spacing='proportional' ) 
[578]714                            if zeorientation == "horizontal" and zetitle[0] != "fill": zecb.ax.set_xlabel(zetitle[index_f]) ; zetitle[index_f]=""
[451]715                        if winds:
[548]716                            if typefile in ['meso']:
[451]717                                [vecx_frame,vecy_frame] = [dumpbdy(vecx_frame,6,stag=uchar,condition=True), dumpbdy(vecy_frame,6,stag=vchar,condition=True)]
718                                key = True
[637]719                                if fvar in ['UV','uv','uvmet']: key = False
[451]720                            elif typefile in ['gcm']:
721                                key = False
[548]722                            if metwind and mapmode == 1:   [vecx_frame,vecy_frame] = m.rotate_vector(vecx_frame, vecy_frame, lon2d, lat2d)
[721]723                            if trans != 0.0:   colorvec = definecolorvec(colorb) 
724                            else:              colorvec = definecolorvec(back) 
[451]725                            vectorfield(vecx_frame, vecy_frame, x, y, stride=stride, csmooth=2,\
726                                             #scale=15., factor=300., color=colorvec, key=key)
[721]727                                             scale=20., factor=250./facwind, color=colorvec, key=key)
[451]728                                                              #200.         ## or csmooth=stride
[721]729                        ### THIS IS A QUITE SPECIFIC PIECE (does not work for mesoscale files)
[763]730                        if ope == '-_histo' and nplot == numplot: # this should work as long as ope is '-' guarantees 3 plots for 4 panels without contour
[640]731                            subplot(subv,subh,nplot+1)
[642]732                            rcParams["legend.fontsize"] = 'xx-large'
[641]733                            if indexlat is None:
734                                latmin = -50.; latmax = 50. # latitude range for histogram of difference
[640]735                                zeindexlat = (lat<latmax)*(lat>latmin)
[721]736                                if typefile in ['meso']: zeindexlat = 10
[642]737                                # this follows the define_axis logic in myplot.py:
738                                if indextime is None or indexlon is None: what_I_plot_frame = what_I_plot_frame[zeindexlat,:]
739                                else: what_I_plot_frame = what_I_plot_frame[:,zeindexlat]
[640]740                            toplot = np.ravel(what_I_plot_frame[np.isnan(what_I_plot_frame)==False])
[641]741                            zebins = np.linspace(zevmin,zevmax,num=30)
[640]742                            hist(toplot,bins=zebins,histtype='step',linewidth=2,color='k',normed=True)
[647]743                            zestd = np.std(toplot);zemean = np.mean(toplot)
[641]744                            zebins = np.linspace(zevmin,zevmax,num=300)
[640]745                            zegauss = (1./(zestd * np.sqrt(2 * np.pi)) * np.exp( - (zebins - zemean)**2 / (2 * zestd**2) ) )
[642]746                            plot(zebins, zegauss, linewidth=1, color='r',label="mean: "+str(zemean)[0:5]+"\nstd: "+str(zestd)[0:5])
747                            legend(loc=0,frameon=False)
[640]748                            subplot(subv,subh,nplot) # go back to last plot for title of contour difference
[763]749                        if ope is not None and "only" not in ope: title("fig(1) "+ope+" fig(2)")
750                        elif ope is not None and "only" in ope: title("fig(1) "+ope[0]+" fig(2)")
[642]751                           
[448]752                    elif which == "contour":
[562]753                        zevminc, zevmaxc = calculate_bounds(what_I_plot_frame, vmin=min(what_I_plot_frame), vmax=max(what_I_plot_frame))
[448]754                        zelevels = np.linspace(zevminc,zevmaxc,ticks/2) #20)
[721]755                        ### another dirty specific stuff in the wall
[748]756                        if var2 == 'HGT':        zelevels = np.arange(-10000.,30000.,500.) #1000.)
[721]757                        elif var2 == 'tpot':     zelevels = np.arange(270,370,5)
758                        elif var2 == 'tk':       zelevels = np.arange(150,250,5)
759                        elif var2 == 'wstar':    zelevels = np.arange(0,10,1.0)
760                        elif var2 == 'zmax_th':  zelevels = np.arange(0,10,2.0) ; what_I_plot_frame = what_I_plot_frame / 1000.
761                        ###
[608]762                        if mapmode == 0:   
[637]763                            what_I_plot_frame, x, y = define_axis( lonyeah,latyeah,vert,time,indexlon,indexlat,indexvert,\
[763]764                                                              itime,what_I_plot_frame, len(all_var2[index_f].shape),vertmode,redope)
[637]765                            ## is this needed? only if len(all_var2[index_f].shape) != len(all_var[index_f].shape)
766                            cs = contour( x,y,what_I_plot_frame, zelevels, colors='k', linewidths = 0.33)#, alpha=0.5, linestyles='solid')
[721]767                            #cs = contour( x,y,what_I_plot_frame, zelevels, colors='w', linewidths = 0.5)#, alpha=0.5, linestyles='solid')
[637]768                            clabel(cs, inline=1, fontsize = 4.*rcParams['font.size']/5., fmt=fmtvar(var2.upper()))
[721]769                        elif mapmode == 1: 
770                            cs = m.contour( x,y,what_I_plot_frame, zelevels, colors='k', linewidths = 0.33)#, alpha=0.5, linestyles='solid')
771                            #clabel(cs, inline=0, fontsize = rcParams['font.size'], fmt="%.0f") #fmtvar(var2.upper()))
[763]772                    if which in ["regular","unidim","xy"]:
[448]773
[763]774                        if nplot > 1 and which in ["unidim","xy"]:
[475]775                           pass  ## because we superimpose nplot instances
776                        else:
777                           # Axis directives for movie frames [including the first one).
778                           zxmin, zxmax = xaxis ; zymin, zymax = yaxis
779                           if zxmin is not None: mpl.pyplot.xlim(xmin=zxmin)
780                           if zxmax is not None: mpl.pyplot.xlim(xmax=zxmax)
781                           if zymin is not None: mpl.pyplot.ylim(ymin=zymin)
782                           if zymax is not None: mpl.pyplot.ylim(ymax=zymax)
[763]783                           if ylog and not xlog:      mpl.pyplot.semilogy()
784                           if xlog and not ylog:      mpl.pyplot.semilogx()
785                           if xlog and ylog: 
786                                mpl.pyplot.xscale('log') 
787                                mpl.pyplot.yscale('log')
[475]788                           if invert_y:  ax = mpl.pyplot.gca() ; ax.set_ylim(ax.get_ylim()[::-1])
[502]789                           if xlab is not None: xlabel(xlab)
790                           if ylab is not None: ylabel(ylab)
[475]791
[448]792                        if mrate is not None:
[451]793                           ### THIS IS A MENCODER MOVIE
794                           if mrate > 0:
795                             figframe=mpl.pyplot.gcf()
796                             if mquality:   figframe.set_dpi(600.)
797                             else:          figframe.set_dpi(200.)
798                             mframe=fig2img(figframe)
799                             if imov == 0:
800                                moviename='movie' ;W,H = figframe.canvas.get_width_height()
801                                video = VideoSink((H,W), moviename, rate=mrate, byteorder="rgba")
802                             video.run(mframe) ; close()
803                             if imov == iend: video.close()                           
804                           ### THIS IS A WEBPAGE MOVIE
805                           else:
806                             nameframe = "image"+str(1000+imov)
807                             makeplotres(nameframe,res=100.,disp=False) ; close()
808                             if imov == 0: myfile = open("zepics", 'w')
809                             myfile.write("modImages["+str(imov)+"] = '"+nameframe+"_100.png';"+ '\n')
810                             if imov == iend:
811                                 myfile.write("first_image = 0;"+ '\n')
812                                 myfile.write("last_image = "+str(iend)+";"+ '\n')
813                                 myfile.close()
[475]814                        if var2 and which == "regular":  which = "contour"
[448]815                        imov = imov+1
816                    elif which == "contour":
817                        which = "regular"
818
[345]819       ### Next subplot
[483]820       zevarname = varname
821       if redope is not None: zevarname = zevarname + "_" + redope
822       basename = getname(var=zevarname,var2=var2,winds=winds,anomaly=anomaly)
[462]823       if len(what_I_plot.shape) > 3:
824           basename = basename + getstralt(nc,level) 
[451]825       if mrate is not None: basename = "movie_" + basename
[548]826       if typefile in ['meso']:
[569]827            if sslon is not None: basename = basename + "_lon_" + str(int(round(lonp)))
828            if sslat is not None: basename = basename + "_lat_" + str(int(round(latp)))
[349]829            plottitle = basename
[468]830            ### dans le nouveau systeme time=ls,sol,lt cette ligne pourrait ne servir a rien (ou deplacer au dessus)
831            if addchar and indextime is not None:   [addchar,gogol,gogol2] = getlschar ( all_namefile[index_f] )  ;  plottitle = plottitle + addchar
[485]832            ### en fait redope is None doit etre remplace par : n'est ni maxt ni mint
833            if redope is None and ltst is not None and ( (mapmode == 0) or (proj in ["lcc","laea","merc","nsper"]) ):  plottitle = plottitle + "_LT" + str(ltst)
[349]834       else:
[359]835            if fileref is not None:
[763]836               if (index_f+1)%3 == 0:     plottitle = basename+' '+"fig(1) "+ope+" fig(2)"
837               elif (index_f+2)%3 == 0:   plottitle = basename+' '+fileref
838               else:                        plottitle = basename+' '+all_namefile[index_f]
[647]839            else:                            plottitle = basename+' '+all_namefile[index_f]
[468]840       if mult != 1:                         plottitle = '{:.0e}'.format(mult) + "*" + plottitle
[578]841       if zetitle[0] != "fill":                 
[580]842          if index_f < len(zetitle): plottitle = zetitle[index_f]
843          else: plottitle = zetitle[0]
[578]844          if titleref is "fill":             titleref=zetitle[index_f]
[763]845          if fileref is not None and which != "xy":
846             if (index_f+2)%3 == 0:        plottitle = titleref
847             if (index_f+1)%3 == 0:        plottitle = "fig(1) "+ope+" fig(2)"
[402]848#       if indexlon is not None:      plottitle = plottitle + " lon: " + str(min(lon[indexlon])) +" "+ str(max(lon[indexlon]))
849#       if indexlat is not None:      plottitle = plottitle + " lat: " + str(min(lat[indexlat])) +" "+ str(max(lat[indexlat]))
850#       if indexvert is not None:     plottitle = plottitle + " vert: " + str(min(vert[indexvert])) +" "+ str(max(vert[indexvert]))
851#       if indextime is not None:     plottitle = plottitle + " time: " + str(min(time[indextime])) +" "+ str(max(time[indextime]))
[518]852       if colorb != "onebar": title( plottitle )
[402]853       if nplot >= numplot: error = True
[345]854       nplot += 1
855
[763]856       if len(namefiles) > 1 and len(var) == 1 and which == "unidim": index_f=index_f+1
857       firstpass=False
858
[518]859    if colorb == "onebar":
860        cax = axes([0.1, 0.2, 0.8, 0.03]) # a ameliorer
861        zecb = colorbar(cax=cax, orientation="horizontal", format=fmtvar(fvar.upper()),\
862                 ticks=np.linspace(zevmin,zevmax,num=min([ticks/2+1,21])),extend='neither',spacing='proportional')
[578]863        if zetitle[0] != "fill": zecb.ax.set_xlabel(zetitle[index_f]) ; zetitle[index_f]=""
[345]864
[612]865
[345]866    ##########################################################################
867    ### Save the figure in a file in the data folder or an user-defined folder
[359]868    if outputname is None:
[548]869       if typefile in ['meso']:   prefix = getprefix(nc)
[359]870       elif typefile in ['gcm']:            prefix = 'LMD_GCM_'
871       else:                                prefix = ''
[345]872    ###
[359]873       zeplot = prefix + basename
[721]874       if zoom:            zeplot = zeplot + "zoom"+str(abs(zoom))
[359]875       if addchar:         zeplot = zeplot + addchar
876       if numplot <= 0:    zeplot = zeplot + "_LT"+str(abs(numplot))
[345]877    ###
[359]878       if not target:      zeplot = namefile[0:find(namefile,'wrfout')] + zeplot
879       else:               zeplot = target + "/" + zeplot 
[345]880    ###
[359]881    else:
882       zeplot=outputname
883
[436]884    if mrate is None:
885        pad_inches_value = 0.35
[760]886        if wlon[1]-wlon[0] < 2.: pad_inches_value = 0.5  # LOCAL MODE (small values)
[436]887        print "********** SAVE ", save
888        if save == 'png': 
889            if display: makeplotres(zeplot,res=100.,pad_inches_value=pad_inches_value) #,erase=True)  ## a miniature
890            makeplotres(zeplot,res=resolution,pad_inches_value=pad_inches_value,disp=False)
891        elif save in ['eps','svg','pdf']:     makeplotres(zeplot,pad_inches_value=pad_inches_value,disp=False,ext=save)
892        elif save == 'gui':                   show()
[489]893        elif save == 'donothing':             pass
[436]894        elif save == 'txt':                   print "Saved results in txt file." 
895        else: 
896            print "INFO: save mode not supported. using gui instead."
897            show()
[345]898
[447]899    ###################################
900    #### Getting more out of this video -- PROBLEMS WITH CREATED VIDEOS
901    #
902    #if mrate is not None:
903    #    print "Re-encoding movie.. first pass"
904    #    video.first_pass(filename=moviename,quality=mquality,rate=mrate)
905    #    print "Re-encoding movie.. second pass"
906    #    video.second_pass(filename=moviename,quality=mquality,rate=mrate)   
[444]907
[345]908    ###############
909    ### Now the end
910    return zeplot
Note: See TracBrowser for help on using the repository browser.