1 | #################################################### |
---|
2 | ### A Python Class for the Mars Climate Database ### |
---|
3 | ### ---------------------------------------------### |
---|
4 | ### Aymeric SPIGA 17-21/04/2012 ### |
---|
5 | ### ---------------------------------------------### |
---|
6 | ### (see mcdtest.py for examples of use) ### |
---|
7 | #################################################### |
---|
8 | |
---|
9 | import numpy as np |
---|
10 | import fmcd |
---|
11 | import matplotlib.pyplot as mpl |
---|
12 | import myplot |
---|
13 | |
---|
14 | |
---|
15 | class mcd(): |
---|
16 | |
---|
17 | def __repr__(self): |
---|
18 | # print out a help string when help is invoked on the object |
---|
19 | whatprint = 'MCD object. \"help(mcd)\" for more information\n' |
---|
20 | return whatprint |
---|
21 | |
---|
22 | ######################## |
---|
23 | ### Default settings ### |
---|
24 | ######################## |
---|
25 | |
---|
26 | def __init__(self): |
---|
27 | # default settings |
---|
28 | ## 0. general stuff |
---|
29 | self.name = "MCD v4.3 output" |
---|
30 | #self.dset = '/home/aymeric/Science/MCD_v4.3/data/' |
---|
31 | self.dset = '/home/marshttp/MCD_v4.3/data/' |
---|
32 | ## 1. spatio-temporal coordinates |
---|
33 | self.lat = 0. |
---|
34 | self.lon = 0. |
---|
35 | self.loct = 0. |
---|
36 | self.xdate = 0. # see datekey |
---|
37 | self.xz = 10. # see zkey |
---|
38 | ## 1bis. related settings |
---|
39 | self.zkey = 3 # specify that xz is the altitude above surface (m) |
---|
40 | self.datekey = 1 # 0 = "Earth time": xdate is given in Julian days (localtime must be set to zero) |
---|
41 | # 1 = "Mars date": xdate is the value of Ls |
---|
42 | ## 2. climatological options |
---|
43 | self.dust = 2 #our best guess MY24 scenario, with solar average conditions |
---|
44 | self.hrkey = 1 #set high resolution mode on (hrkey=0 to set high resolution off) |
---|
45 | ## 3. additional settings for advanced use |
---|
46 | self.extvarkey = 1 #extra output variables (1: yes, 0: no) |
---|
47 | self.perturkey = 0 #integer perturkey ! perturbation type (0: none) |
---|
48 | self.seedin = 1 #random number generator seed (unused if perturkey=0) |
---|
49 | self.gwlength = 0. #gravity Wave wavelength (unused if perturkey=0) |
---|
50 | ## outputs. just to define attributes. |
---|
51 | ## --> in update |
---|
52 | self.pres = None ; self.dens = None ; self.temp = None ; self.zonwind = None ; self.merwind = None ; self.meanvar = None ; self.extvar = None |
---|
53 | self.seedout = None ; self.ierr = None |
---|
54 | ## --> in prepare |
---|
55 | self.xcoord = None ; self.ycoord = None |
---|
56 | self.prestab = None ; self.denstab = None ; self.temptab = None |
---|
57 | self.zonwindtab = None ; self.merwindtab = None ; self.meanvartab = None ; self.extvartab = None |
---|
58 | |
---|
59 | def viking1(self): self.name = "Viking 1 site. MCD v4.3 output" ; self.lat = 22.48 ; self.lon = -49.97 ; self.xdate = 97. |
---|
60 | def viking2(self): self.name = "Viking 2 site. MCD v4.3 output" ; self.lat = 47.97 ; self.lon = -225.74 ; self.xdate = 117.6 |
---|
61 | |
---|
62 | def getextvarlab(self,num): |
---|
63 | whichfield = { \ |
---|
64 | 91: "Pressure (Pa)", \ |
---|
65 | 92: "Density (kg/m3)", \ |
---|
66 | 93: "Temperature (K)", \ |
---|
67 | 94: "W-E wind component (m/s)", \ |
---|
68 | 95: "S-N wind component (m/s)", \ |
---|
69 | 1: "Radial distance from planet center (m)",\ |
---|
70 | 2: "Altitude above areoid (Mars geoid) (m)",\ |
---|
71 | 3: "Altitude above local surface (m)",\ |
---|
72 | 4: "orographic height (m) (surface altitude above areoid)",\ |
---|
73 | 5: "Ls, solar longitude of Mars (deg)",\ |
---|
74 | 6: "LST local true solar time (hrs)",\ |
---|
75 | 7: "Universal solar time (LST at lon=0) (hrs)",\ |
---|
76 | 8: "Air heat capacity Cp (J kg-1 K-1)",\ |
---|
77 | 9: "gamma=Cp/Cv Ratio of specific heats",\ |
---|
78 | 10: "density RMS day to day variations (kg/m^3)",\ |
---|
79 | 11: "[not defined]",\ |
---|
80 | 12: "[not defined]",\ |
---|
81 | 13: "scale height H(p) (m)",\ |
---|
82 | 14: "GCM orography (m)",\ |
---|
83 | 15: "surface temperature (K)",\ |
---|
84 | 16: "daily maximum mean surface temperature (K)",\ |
---|
85 | 17: "daily minimum mean surface temperature (K)",\ |
---|
86 | 18: "surf. temperature RMS day to day variations (K)",\ |
---|
87 | 19: "surface pressure (high resolution if hireskey=1)",\ |
---|
88 | 20: "GCM surface pressure (Pa)",\ |
---|
89 | 21: "atmospheric pressure RMS day to day variations (Pa)",\ |
---|
90 | 22: "surface pressure RMS day to day variations (Pa)",\ |
---|
91 | 23: "temperature RMS day to day variations (K)",\ |
---|
92 | 24: "zonal wind RMS day to day variations (m/s)",\ |
---|
93 | 25: "meridional wind RMS day to day variations (m/s)",\ |
---|
94 | 26: "vertical wind component (m/s) >0 when downwards!",\ |
---|
95 | 27: "vertical wind RMS day to day variations (m/s)",\ |
---|
96 | 28: "small scale perturbation (gravity wave) (kg/m^3)",\ |
---|
97 | 29: "q2: turbulent kinetic energy (m2/s2)",\ |
---|
98 | 30: "[not defined]",\ |
---|
99 | 31: "thermal IR flux to surface (W/m2)",\ |
---|
100 | 32: "solar flux to surface (W/m2)",\ |
---|
101 | 33: "thermal IR flux to space (W/m2)",\ |
---|
102 | 34: "solar flux reflected to space (W/m2)",\ |
---|
103 | 35: "surface CO2 ice layer (kg/m2)",\ |
---|
104 | 36: "DOD: Dust column visible optical depth",\ |
---|
105 | 37: "Dust mass mixing ratio (kg/kg)",\ |
---|
106 | 38: "DOD RMS day to day variations",\ |
---|
107 | 39: "DOD total standard deviation over season",\ |
---|
108 | 40: "Water vapor column (kg/m2)",\ |
---|
109 | 41: "Water vapor vol. mixing ratio (mol/mol)",\ |
---|
110 | 42: "Water ice column (kg/m2)",\ |
---|
111 | 43: "Water ice mixing ratio (mol/mol)",\ |
---|
112 | 44: "O3 ozone vol. mixing ratio (mol/mol)",\ |
---|
113 | 45: "[CO2] vol. mixing ratio (mol/mol)",\ |
---|
114 | 46: "[O] vol. mixing ratio (mol/mol)",\ |
---|
115 | 47: "[N2] vol. mixing ratio (mol/mol)",\ |
---|
116 | 48: "[CO] vol. mixing ratio (mol/mol)",\ |
---|
117 | 49: "R: Molecular gas constant (J K-1 kg-1)",\ |
---|
118 | 50: "Air viscosity estimation (N s m-2)" |
---|
119 | } |
---|
120 | if num not in whichfield: myplot.errormess("Incorrect subscript in extvar.") |
---|
121 | return whichfield[num] |
---|
122 | |
---|
123 | def convertlab(self,num): |
---|
124 | ## a conversion from text inquiries to extvar numbers. to be completed. |
---|
125 | if num == "p": num = 91 |
---|
126 | elif num == "rho": num = 92 |
---|
127 | elif num == "t": num = 93 |
---|
128 | elif num == "u": num = 94 |
---|
129 | elif num == "v": num = 95 |
---|
130 | elif num == "tsurf": num = 15 |
---|
131 | elif num == "topo": num = 4 |
---|
132 | elif num == "h": num = 13 |
---|
133 | elif num == "ps": num = 19 |
---|
134 | elif num == "tau": num = 36 |
---|
135 | elif num == "mtot": num = 40 |
---|
136 | elif num == "icetot": num = 42 |
---|
137 | elif num == "ps_ddv": num = 22 |
---|
138 | elif num == "h2ovap": num = 41 |
---|
139 | elif num == "h2oice": num = 43 |
---|
140 | elif num == "cp": num = 8 |
---|
141 | elif num == "rho_ddv": num = 10 |
---|
142 | elif num == "tsurfmx": num = 16 |
---|
143 | elif num == "tsurfmn": num = 17 |
---|
144 | elif num == "lwdown": num = 31 |
---|
145 | elif num == "swdown": num = 32 |
---|
146 | elif num == "lwup": num = 33 |
---|
147 | elif num == "swup": num = 34 |
---|
148 | elif num == "o3": num = 44 |
---|
149 | elif num == "o": num = 46 |
---|
150 | elif num == "co": num = 48 |
---|
151 | elif num == "visc": num = 50 |
---|
152 | elif num == "co2ice": num = 35 |
---|
153 | elif not isinstance(num, np.int): myplot.errormess("field reference not found.") |
---|
154 | return num |
---|
155 | |
---|
156 | ################### |
---|
157 | ### One request ### |
---|
158 | ################### |
---|
159 | |
---|
160 | def update(self): |
---|
161 | # retrieve fields from MCD (call_mcd). more info in fmcd.call_mcd.__doc__ |
---|
162 | (self.pres, self.dens, self.temp, self.zonwind, self.merwind, \ |
---|
163 | self.meanvar, self.extvar, self.seedout, self.ierr) \ |
---|
164 | = \ |
---|
165 | fmcd.call_mcd(self.zkey,self.xz,self.lon,self.lat,self.hrkey, \ |
---|
166 | self.datekey,self.xdate,self.loct,self.dset,self.dust, \ |
---|
167 | self.perturkey,self.seedin,self.gwlength,self.extvarkey ) |
---|
168 | ## we use the end of extvar (unused) to store meanvar. this is convenient for getextvar(lab) |
---|
169 | self.extvar[90] = self.pres ; self.extvar[91] = self.dens |
---|
170 | self.extvar[92] = self.temp ; self.extvar[93] = self.zonwind ; self.extvar[94] = self.merwind |
---|
171 | |
---|
172 | def printset(self): |
---|
173 | # print main settings |
---|
174 | print "zkey",self.zkey,"xz",self.xz,"lon",self.lon,"lat",self.lat,"hrkey",self.hrkey, \ |
---|
175 | "xdate",self.xdate,"loct",self.loct,"dust",self.dust |
---|
176 | |
---|
177 | def getnameset(self): |
---|
178 | # set a name referring to settings [convenient for databases] |
---|
179 | name = str(self.zkey)+str(self.xz)+str(self.lon)+str(self.lat)+str(self.hrkey)+str(self.datekey)+str(self.xdate)+str(self.loct)+str(self.dust) |
---|
180 | return name |
---|
181 | |
---|
182 | def printcoord(self): |
---|
183 | # print requested space-time coordinates |
---|
184 | print "LAT",self.lat,"LON",self.lon,"LOCT",self.loct,"XDATE",self.xdate |
---|
185 | |
---|
186 | def printmeanvar(self): |
---|
187 | # print mean MCD variables |
---|
188 | print "Pressure = %5.3f pascals. " % (self.pres) |
---|
189 | print "Density = %5.3f kilograms per cubic meter. " % (self.dens) |
---|
190 | print "Temperature = %3.0f kelvins (%4.0f degrees celsius)." % (self.temp,self.temp-273.15) |
---|
191 | print "Zonal wind = %5.3f meters per second." % (self.zonwind) |
---|
192 | print "Meridional wind = %5.3f meters per second." % (self.merwind) |
---|
193 | |
---|
194 | def printextvar(self,num): |
---|
195 | # print extra MCD variables |
---|
196 | num = self.convertlab(num) |
---|
197 | print self.getextvarlab(num) + " ..... " + str(self.extvar[num-1]) |
---|
198 | |
---|
199 | def printallextvar(self): |
---|
200 | # print all extra MCD variables |
---|
201 | for i in range(50): self.printextvar(i+1) |
---|
202 | |
---|
203 | def htmlprinttabextvar(self,tabtodo): |
---|
204 | print "Results from the Mars Climate Database" |
---|
205 | print "<ul>" |
---|
206 | for i in range(len(tabtodo)): print "<li>" ; self.printextvar(tabtodo[i]) ; print "</li>" |
---|
207 | print "</ul>" |
---|
208 | print "<hr>" |
---|
209 | print "SETTINGS<br />" |
---|
210 | self.printcoord() |
---|
211 | self.printset() |
---|
212 | |
---|
213 | def printmcd(self): |
---|
214 | # 1. call MCD 2. print settings 3. print mean vars |
---|
215 | self.update() |
---|
216 | self.printcoord() |
---|
217 | print "-------------------------------------------" |
---|
218 | self.printmeanvar() |
---|
219 | |
---|
220 | ######################## |
---|
221 | ### Several requests ### |
---|
222 | ######################## |
---|
223 | |
---|
224 | def prepare(self,ndx=None,ndy=None): |
---|
225 | ### prepare I/O arrays for 1d slices |
---|
226 | if ndx is None: print "No dimension in prepare. Exit. Set at least ndx." ; exit() |
---|
227 | else: self.xcoord = np.ones(ndx) |
---|
228 | if ndy is None: dashape = (ndx) ; dashapemean = (ndx,6) ; dashapeext = (ndx,101) ; self.ycoord = None |
---|
229 | else: dashape = (ndx,ndy) ; dashapemean = (ndx,ndy,6) ; dashapeext = (ndx,ndy,101) ; self.ycoord = np.ones(ndy) |
---|
230 | self.prestab = np.ones(dashape) ; self.denstab = np.ones(dashape) ; self.temptab = np.ones(dashape) |
---|
231 | self.zonwindtab = np.ones(dashape) ; self.merwindtab = np.ones(dashape) |
---|
232 | self.meanvartab = np.ones(dashapemean) ; self.extvartab = np.ones(dashapeext) |
---|
233 | |
---|
234 | def getextvar(self,num): |
---|
235 | ### get a given var in extvartab |
---|
236 | try: field=self.extvartab[:,:,num] |
---|
237 | except: field=self.extvartab[:,num] |
---|
238 | return field |
---|
239 | |
---|
240 | def definefield(self,choice): |
---|
241 | ### for analysis or plot purposes, set field and field label from user-defined choice |
---|
242 | choice = self.convertlab(choice) |
---|
243 | field = self.getextvar(choice); fieldlab = self.getextvarlab(choice) |
---|
244 | return field,fieldlab |
---|
245 | |
---|
246 | ################### |
---|
247 | ### 1D analysis ### |
---|
248 | ################### |
---|
249 | |
---|
250 | def put1d(self,i): |
---|
251 | ## fill in subscript i in output arrays |
---|
252 | ## (arrays must have been correctly defined through prepare) |
---|
253 | if self.prestab is None: myplot.errormess("arrays must be prepared first through self.prepare") |
---|
254 | self.prestab[i] = self.pres ; self.denstab[i] = self.dens ; self.temptab[i] = self.temp |
---|
255 | self.zonwindtab[i] = self.zonwind ; self.merwindtab[i] = self.merwind |
---|
256 | self.meanvartab[i,1:5] = self.meanvar[0:4] ## note: var numbering according to MCD manual is kept |
---|
257 | self.extvartab[i,1:100] = self.extvar[0:99] ## note: var numbering according to MCD manual is kept |
---|
258 | |
---|
259 | def diurnal(self,nd=13,start=0.,end=24.): |
---|
260 | ### retrieve a local time slice |
---|
261 | self.xlabel = "Local time (Martian hour)" |
---|
262 | self.prepare(ndx=nd) ; self.xcoord = np.linspace(start,end,nd) |
---|
263 | for i in range(nd): self.loct = self.xcoord[i] ; self.update() ; self.put1d(i) |
---|
264 | |
---|
265 | def zonal(self,nd=37,start=-180.,end=180.): |
---|
266 | ### retrieve a longitude slice |
---|
267 | self.xlabel = "East longitude (degrees)" |
---|
268 | self.prepare(ndx=nd) ; self.xcoord = np.linspace(start,end,nd) |
---|
269 | for i in range(nd): self.lon = self.xcoord[i] ; self.update() ; self.put1d(i) |
---|
270 | |
---|
271 | def meridional(self,nd=19,start=-90.,end=90.): |
---|
272 | ### retrieve a latitude slice |
---|
273 | self.xlabel = "North latitude (degrees)" |
---|
274 | self.prepare(ndx=nd) ; self.xcoord = np.linspace(start,end,nd) |
---|
275 | for i in range(nd): self.lat = self.xcoord[i] ; self.update() ; self.put1d(i) |
---|
276 | |
---|
277 | def profile(self,nd=20,start=0.,end=120000.,tabperso=None): |
---|
278 | ### retrieve an altitude slice (profile) |
---|
279 | self.xlabel = "Altitude (m)" |
---|
280 | if tabperso is not None: nd = len(tabperso) |
---|
281 | self.prepare(ndx=nd) |
---|
282 | if tabperso is None: self.xcoord = np.linspace(start,end,nd) |
---|
283 | else: self.xcoord = tabperso |
---|
284 | for i in range(nd): self.xz = self.xcoord[i] ; self.update() ; self.put1d(i) |
---|
285 | |
---|
286 | def seasonal(self,nd=12,start=0.,end=360.): |
---|
287 | ### retrieve a seasonal slice |
---|
288 | self.xlabel = "Areocentric longitude (degrees)" |
---|
289 | self.prepare(ndx=nd) ; self.xcoord = np.linspace(start,end,nd) |
---|
290 | for i in range(nd): self.xdate = self.xcoord[i] ; self.update() ; self.put1d(i) |
---|
291 | |
---|
292 | def makeplot1d(self,choice,vertplot=0): |
---|
293 | ### one 1D plot is created for the user-defined variable in choice. |
---|
294 | (field, fieldlab) = self.definefield(choice) |
---|
295 | if vertplot != 1: absc = self.xcoord ; ordo = field ; ordolab = fieldlab ; absclab = self.xlabel |
---|
296 | else: ordo = self.xcoord ; absc = field ; absclab = fieldlab ; ordolab = self.xlabel |
---|
297 | mpl.plot(absc,ordo,'-bo') ; mpl.ylabel(ordolab) ; mpl.xlabel(absclab) #; mpl.xticks(query.xcoord) |
---|
298 | mpl.figtext(0.5, 0.01, "Mars Climate Database (c) LMD/OU/IAA/ESA/CNES", ha='center') |
---|
299 | |
---|
300 | def plot1d(self,tabtodo,vertplot=0): |
---|
301 | ### complete 1D figure with possible multiplots |
---|
302 | if isinstance(tabtodo,np.str): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
303 | if isinstance(tabtodo,np.int): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
304 | fig = mpl.figure() ; subv,subh = myplot.definesubplot( len(tabtodo) , fig ) |
---|
305 | for i in range(len(tabtodo)): mpl.subplot(subv,subh,i+1).grid(True, linestyle=':', color='grey') ; self.makeplot1d(tabtodo[i],vertplot) |
---|
306 | |
---|
307 | def htmlplot1d(self,tabtodo,vertplot=0,figname="temp.png"): |
---|
308 | ### complete 1D figure with possible multiplots |
---|
309 | ### added in 09/2012 for online MCD |
---|
310 | ### see http://www.dalkescientific.com/writings/diary/archive/2005/04/23/matplotlib_without_gui.html |
---|
311 | from matplotlib.figure import Figure |
---|
312 | from matplotlib.backends.backend_agg import FigureCanvasAgg |
---|
313 | if isinstance(tabtodo,np.str): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
314 | if isinstance(tabtodo,np.int): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
315 | fig = Figure(figsize=(8,8)) ; subv,subh = myplot.definesubplot( len(tabtodo) , fig ) |
---|
316 | for i in range(len(tabtodo)): |
---|
317 | yeah = fig.add_subplot(subv,subh,i+1) #.grid(True, linestyle=':', color='grey') |
---|
318 | choice = tabtodo[i] |
---|
319 | (field, fieldlab) = self.definefield(choice) |
---|
320 | if vertplot != 1: absc = self.xcoord ; ordo = field ; ordolab = fieldlab ; absclab = self.xlabel |
---|
321 | else: ordo = self.xcoord ; absc = field ; absclab = fieldlab ; ordolab = self.xlabel |
---|
322 | yeah.plot(absc,ordo,'-bo') #; mpl.xticks(query.xcoord) |
---|
323 | ax = fig.gca() ; ax.set_ylabel(ordolab) ; ax.set_xlabel(absclab) |
---|
324 | fig.text(0.5, 0.01, "Mars Climate Database (c) LMD/OU/IAA/ESA/CNES", ha='center') |
---|
325 | canvas = FigureCanvasAgg(fig) |
---|
326 | # The size * the dpi gives the final image size |
---|
327 | # a4"x4" image * 80 dpi ==> 320x320 pixel image |
---|
328 | canvas.print_figure(figname, dpi=80) |
---|
329 | |
---|
330 | ################### |
---|
331 | ### 2D analysis ### |
---|
332 | ################### |
---|
333 | |
---|
334 | def latlon(self,ndx=37,startx=-180.,endx=180.,ndy=19,starty=-90.,endy=90.,fixedlt=False): |
---|
335 | ### retrieve a latitude/longitude slice |
---|
336 | ### default is: local time is not fixed. user-defined local time is at longitude 0. |
---|
337 | self.xlabel = "East longitude (degrees)" ; self.ylabel = "North latitude (degrees)" |
---|
338 | self.prepare(ndx=ndx,ndy=ndy) |
---|
339 | self.xcoord = np.linspace(startx,endx,ndx) ; self.ycoord = np.linspace(starty,endy,ndy) |
---|
340 | if not fixedlt: umst = self.loct |
---|
341 | for i in range(ndx): |
---|
342 | for j in range(ndy): |
---|
343 | self.lon = self.xcoord[i] ; self.lat = self.ycoord[j] |
---|
344 | if not fixedlt: self.loct = (umst + self.lon/15.) % 24 |
---|
345 | self.update() ; self.put2d(i,j) |
---|
346 | if not fixedlt: self.loct = umst |
---|
347 | |
---|
348 | def put2d(self,i,j): |
---|
349 | ## fill in subscript i,j in output arrays |
---|
350 | ## (arrays must have been correctly defined through prepare) |
---|
351 | if self.prestab is None: myplot.errormess("arrays must be prepared first through self.prepare") |
---|
352 | self.prestab[i,j] = self.pres ; self.denstab[i,j] = self.dens ; self.temptab[i,j] = self.temp |
---|
353 | self.zonwindtab[i,j] = self.zonwind ; self.merwindtab[i,j] = self.merwind |
---|
354 | self.meanvartab[i,j,1:5] = self.meanvar[0:4] ## note: var numbering according to MCD manual is kept |
---|
355 | self.extvartab[i,j,1:100] = self.extvar[0:99] ## note: var numbering according to MCD manual is kept |
---|
356 | |
---|
357 | def makemap2d(self,choice,incwind=False,fixedlt=False): |
---|
358 | ### one 2D map is created for the user-defined variable in choice. |
---|
359 | self.latlon(fixedlt=fixedlt) ## a map is implicitely a lat-lon plot. otherwise it is a plot (cf. makeplot2d) |
---|
360 | (field, fieldlab) = self.definefield(choice) |
---|
361 | if incwind: |
---|
362 | (windx, fieldlabwx) = self.definefield("u") |
---|
363 | (windy, fieldlabwy) = self.definefield("v") |
---|
364 | myplot.maplatlon(self.xcoord,self.ycoord,field,title=fieldlab,proj="cyl",vecx=windx,vecy=windy) |
---|
365 | else: |
---|
366 | myplot.maplatlon(self.xcoord,self.ycoord,field,title=fieldlab,proj="moll") |
---|
367 | mpl.figtext(0.5, 0.0, "Mars Climate Database (c) LMD/OU/IAA/ESA/CNES", ha='center') |
---|
368 | |
---|
369 | def map2d(self,tabtodo,incwind=False,fixedlt=False): |
---|
370 | ### complete 2D figure with possible multiplots |
---|
371 | if isinstance(tabtodo,np.str): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
372 | if isinstance(tabtodo,np.int): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
373 | fig = mpl.figure() |
---|
374 | subv,subh = myplot.definesubplot( len(tabtodo) , fig ) |
---|
375 | for i in range(len(tabtodo)): mpl.subplot(subv,subh,i+1) ; self.makemap2d(tabtodo[i],incwind=incwind,fixedlt=fixedlt) |
---|
376 | |
---|
377 | def htmlmap2d(self,tabtodo,incwind=False,fixedlt=False,figname="temp.png"): |
---|
378 | ### complete 2D figure with possible multiplots |
---|
379 | ### added in 09/2012 for online MCD |
---|
380 | ### see http://www.dalkescientific.com/writings/diary/archive/2005/04/23/matplotlib_without_gui.html |
---|
381 | from matplotlib.figure import Figure |
---|
382 | from matplotlib.backends.backend_agg import FigureCanvasAgg |
---|
383 | from matplotlib.cm import get_cmap |
---|
384 | if isinstance(tabtodo,np.str): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
385 | if isinstance(tabtodo,np.int): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
386 | fig = Figure(figsize=(8,8)) ; subv,subh = myplot.definesubplot( len(tabtodo) , fig ) |
---|
387 | |
---|
388 | ### topocontours |
---|
389 | fieldc = self.getextvar(self.convertlab("topo")) |
---|
390 | |
---|
391 | for i in range(len(tabtodo)): |
---|
392 | yeah = fig.add_subplot(subv,subh,i+1) |
---|
393 | choice = tabtodo[i] |
---|
394 | self.latlon(fixedlt=fixedlt) ## a map is implicitely a lat-lon plot. otherwise it is a plot (cf. makeplot2d) |
---|
395 | (field, fieldlab) = self.definefield(choice) |
---|
396 | if incwind: (windx, fieldlabwx) = self.definefield("u") ; (windy, fieldlabwy) = self.definefield("v") |
---|
397 | |
---|
398 | proj="cyl" ; colorb="jet" ; ndiv=20 ; zeback="molabw" ; trans=1.0 #0.6 |
---|
399 | title="" ; vecx=None ; vecy=None ; stride=2 |
---|
400 | lon = self.xcoord |
---|
401 | lat = self.ycoord |
---|
402 | |
---|
403 | ### get lon and lat in 2D version. get lat/lon intervals |
---|
404 | #numdim = len(np.array(lon).shape) |
---|
405 | #if numdim == 2: [lon2d,lat2d] = [lon,lat] |
---|
406 | #elif numdim == 1: [lon2d,lat2d] = np.meshgrid(lon,lat) |
---|
407 | #else: errormess("lon and lat arrays must be 1D or 2D") |
---|
408 | #[wlon,wlat] = myplot.latinterv() |
---|
409 | ### define projection and background. define x and y given the projection |
---|
410 | #m = basemap.Basemap(projection='moll') marche pas |
---|
411 | #m = myplot.define_proj(proj,wlon,wlat,back=zeback,blat=None,blon=None) |
---|
412 | #x, y = m(lon2d, lat2d) |
---|
413 | ### TEMP |
---|
414 | x = lon ; y = lat |
---|
415 | ## define field. bound field. |
---|
416 | what_I_plot = np.transpose(field) |
---|
417 | zevmin, zevmax = myplot.calculate_bounds(what_I_plot) ## vmin=min(what_I_plot_frame), vmax=max(what_I_plot_frame)) |
---|
418 | what_I_plot = myplot.bounds(what_I_plot,zevmin,zevmax) |
---|
419 | ## define contour field levels. define color palette |
---|
420 | ticks = ndiv + 1 |
---|
421 | zelevels = np.linspace(zevmin,zevmax,ticks) |
---|
422 | palette = get_cmap(name=colorb) |
---|
423 | ## contours topo |
---|
424 | zelevc = np.linspace(-8000.,20000.,20) |
---|
425 | yeah.contour( x, y, np.transpose(fieldc), zelevc, colors='black',linewidths = 0.4) |
---|
426 | # contour field |
---|
427 | c = yeah.contourf( x, y, what_I_plot, zelevels, cmap = palette, alpha = trans ) |
---|
428 | Figure.colorbar(fig,c,orientation='vertical',format="%.1e") |
---|
429 | ax = fig.gca() ; ax.set_title(fieldlab) ; ax.set_ylabel("Latitude") ; ax.set_xlabel("Longitude") |
---|
430 | ax.set_xticks(np.arange(-180,181,45)) ; ax.set_xbound(lower=-180, upper=180) |
---|
431 | ax.set_yticks(np.arange(-90,91,30)) ; ax.set_ybound(lower=-90, upper=90) |
---|
432 | if incwind: |
---|
433 | [x2d,y2d] = np.meshgrid(x,y) |
---|
434 | yeah.quiver(x2d,y2d,np.transpose(windx),np.transpose(windy)) |
---|
435 | fig.text(0.5, 0.01, "Mars Climate Database (c) LMD/OU/IAA/ESA/CNES", ha='center') |
---|
436 | canvas = FigureCanvasAgg(fig) |
---|
437 | # The size * the dpi gives the final image size |
---|
438 | # a4"x4" image * 80 dpi ==> 320x320 pixel image |
---|
439 | canvas.print_figure(figname, dpi=80) |
---|
440 | |
---|
441 | |
---|
442 | ### TODO: makeplot2d, plot2d, passer plot settings |
---|
443 | |
---|