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

Last change on this file since 613 was 613, checked in by aslmd, 13 years ago

UTIL PYTHON: an example script for pp.py (whole ExoMars? report). minor corrections. MESOSCALE : dust storm def files

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