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" |
---|
30 | self.ack = "Mars Climate Database (c) LMD/OU/IAA/ESA/CNES" |
---|
31 | #self.dset = '/home/aymeric/Science/MCD_v4.3/data/' |
---|
32 | self.dset = '/home/marshttp/MCD_v4.3/data/' |
---|
33 | ## 1. spatio-temporal coordinates |
---|
34 | self.lat = 0. |
---|
35 | self.lats = None |
---|
36 | self.late = None |
---|
37 | self.lon = 0. |
---|
38 | self.lons = None |
---|
39 | self.lone = None |
---|
40 | self.loct = 0. |
---|
41 | self.locts = None |
---|
42 | self.locte = None |
---|
43 | self.xdate = 0. # see datekey |
---|
44 | self.xdates = None |
---|
45 | self.xdatee = None |
---|
46 | self.xz = 10. # see zkey |
---|
47 | self.xzs = None |
---|
48 | self.xze = None |
---|
49 | ## 1bis. related settings |
---|
50 | self.zkey = 3 # specify that xz is the altitude above surface (m) |
---|
51 | # zkey : <integer> type of vertical coordinate xz |
---|
52 | # 1 = radius from centre of planet (m) |
---|
53 | # 2 = height above areoid (m) (MOLA zero datum) |
---|
54 | # 3 = height above surface (m) |
---|
55 | # 4 = pressure level (Pa) |
---|
56 | # 5 = altitude above mean Mars Radius(=3396000m) (m) |
---|
57 | self.datekey = 1 # 0 = "Earth time": xdate is given in Julian days (localtime must be set to zero) |
---|
58 | # 1 = "Mars date": xdate is the value of Ls |
---|
59 | ## 2. climatological options |
---|
60 | self.dust = 2 #our best guess MY24 scenario, with solar average conditions |
---|
61 | self.hrkey = 1 #set high resolution mode on (hrkey=0 to set high resolution off) |
---|
62 | ## 3. additional settings for advanced use |
---|
63 | self.extvarkey = 1 #extra output variables (1: yes, 0: no) |
---|
64 | self.perturkey = 0 #integer perturkey ! perturbation type (0: none) |
---|
65 | self.seedin = 1 #random number generator seed (unused if perturkey=0) |
---|
66 | self.gwlength = 0. #gravity Wave wavelength (unused if perturkey=0) |
---|
67 | ## outputs. just to define attributes. |
---|
68 | ## --> in update |
---|
69 | self.pres = None ; self.dens = None ; self.temp = None ; self.zonwind = None ; self.merwind = None ; self.meanvar = None ; self.extvar = None |
---|
70 | self.seedout = None ; self.ierr = None |
---|
71 | ## --> in prepare |
---|
72 | self.xcoord = None ; self.ycoord = None |
---|
73 | self.prestab = None ; self.denstab = None ; self.temptab = None |
---|
74 | self.zonwindtab = None ; self.merwindtab = None ; self.meanvartab = None ; self.extvartab = None |
---|
75 | ## plot stuff |
---|
76 | self.xlabel = None ; self.ylabel = None |
---|
77 | self.vertplot = False |
---|
78 | self.fmt = "%.2e" |
---|
79 | |
---|
80 | def viking1(self): self.name = "Viking 1 site. MCD v4.3 output" ; self.lat = 22.48 ; self.lon = -49.97 ; self.xdate = 97. |
---|
81 | def viking2(self): self.name = "Viking 2 site. MCD v4.3 output" ; self.lat = 47.97 ; self.lon = -225.74 ; self.xdate = 117.6 |
---|
82 | |
---|
83 | def getdustlabel(self): |
---|
84 | if self.dust == 1: self.dustlabel = "MY24 minimum solar scenario" |
---|
85 | elif self.dust == 2: self.dustlabel = "MY24 average solar scenario" |
---|
86 | elif self.dust == 3: self.dustlabel = "MY24 maximum solar scenario" |
---|
87 | elif self.dust == 4: self.dustlabel = "dust storm minimum solar scenario" |
---|
88 | elif self.dust == 5: self.dustlabel = "dust storm average solar scenario" |
---|
89 | elif self.dust == 6: self.dustlabel = "dust storm maximum solar scenario" |
---|
90 | elif self.dust == 7: self.dustlabel = "warm scenario (dusty, maximum solar)" |
---|
91 | elif self.dust == 8: self.dustlabel = "cold scenario (low dust, minimum solar)" |
---|
92 | |
---|
93 | def gettitle(self,oneline=False): |
---|
94 | self.getdustlabel() |
---|
95 | self.title = self.name + " with " + self.dustlabel + "." |
---|
96 | if self.datekey == 1: self.title = self.title + " Ls " + str(self.xdate) + "deg." |
---|
97 | elif self.datekey == 0: self.title = self.title + " JD " + str(self.xdate) + "." |
---|
98 | if not oneline: self.title = self.title + "\n" |
---|
99 | if self.lats is None: self.title = self.title + " Latitude " + str(self.lat) + "E" |
---|
100 | if self.lons is None: self.title = self.title + " Longitude " + str(self.lon) + "N" |
---|
101 | if self.xzs is None: |
---|
102 | self.vertunits() |
---|
103 | self.title = self.title + " Altitude " + str(self.xz) + " " + self.vunits |
---|
104 | if self.locts is None: self.title = self.title + " Local time " + str(self.loct) + "h" |
---|
105 | |
---|
106 | def getextvarlab(self,num): |
---|
107 | whichfield = { \ |
---|
108 | 91: "Pressure (Pa)", \ |
---|
109 | 92: "Density (kg/m3)", \ |
---|
110 | 93: "Temperature (K)", \ |
---|
111 | 94: "W-E wind component (m/s)", \ |
---|
112 | 95: "S-N wind component (m/s)", \ |
---|
113 | 1: "Radial distance from planet center (m)",\ |
---|
114 | 2: "Altitude above areoid (Mars geoid) (m)",\ |
---|
115 | 3: "Altitude above local surface (m)",\ |
---|
116 | 4: "orographic height (m) (surf alt above areoid)",\ |
---|
117 | 5: "Ls, solar longitude of Mars (deg)",\ |
---|
118 | 6: "LST local true solar time (hrs)",\ |
---|
119 | 7: "Universal solar time (LST at lon=0) (hrs)",\ |
---|
120 | 8: "Air heat capacity Cp (J kg-1 K-1)",\ |
---|
121 | 9: "gamma=Cp/Cv Ratio of specific heats",\ |
---|
122 | 10: "density RMS day to day variations (kg/m^3)",\ |
---|
123 | 11: "[not defined]",\ |
---|
124 | 12: "[not defined]",\ |
---|
125 | 13: "scale height H(p) (m)",\ |
---|
126 | 14: "GCM orography (m)",\ |
---|
127 | 15: "surface temperature (K)",\ |
---|
128 | 16: "daily max mean surface temperature (K)",\ |
---|
129 | 17: "daily min mean surface temperature (K)",\ |
---|
130 | 18: "surf. temperature RMS day to day variations (K)",\ |
---|
131 | 19: "surface pressure (high resolution if hireskey=1)",\ |
---|
132 | 20: "GCM surface pressure (Pa)",\ |
---|
133 | 21: "atmospheric pressure RMS day to day variations (Pa)",\ |
---|
134 | 22: "surface pressure RMS day to day variations (Pa)",\ |
---|
135 | 23: "temperature RMS day to day variations (K)",\ |
---|
136 | 24: "zonal wind RMS day to day variations (m/s)",\ |
---|
137 | 25: "meridional wind RMS day to day variations (m/s)",\ |
---|
138 | 26: "vertical wind component (m/s) >0 when downwards!",\ |
---|
139 | 27: "vertical wind RMS day to day variations (m/s)",\ |
---|
140 | 28: "small scale perturbation (gravity wave) (kg/m^3)",\ |
---|
141 | 29: "q2: turbulent kinetic energy (m2/s2)",\ |
---|
142 | 30: "[not defined]",\ |
---|
143 | 31: "thermal IR flux to surface (W/m2)",\ |
---|
144 | 32: "solar flux to surface (W/m2)",\ |
---|
145 | 33: "thermal IR flux to space (W/m2)",\ |
---|
146 | 34: "solar flux reflected to space (W/m2)",\ |
---|
147 | 35: "surface CO2 ice layer (kg/m2)",\ |
---|
148 | 36: "DOD: Dust column visible optical depth",\ |
---|
149 | 37: "Dust mass mixing ratio (kg/kg)",\ |
---|
150 | 38: "DOD RMS day to day variations",\ |
---|
151 | 39: "DOD total standard deviation over season",\ |
---|
152 | 40: "Water vapor column (kg/m2)",\ |
---|
153 | 41: "Water vapor vol. mixing ratio (mol/mol)",\ |
---|
154 | 42: "Water ice column (kg/m2)",\ |
---|
155 | 43: "Water ice mixing ratio (mol/mol)",\ |
---|
156 | 44: "O3 ozone vol. mixing ratio (mol/mol)",\ |
---|
157 | 45: "[CO2] vol. mixing ratio (mol/mol)",\ |
---|
158 | 46: "[O] vol. mixing ratio (mol/mol)",\ |
---|
159 | 47: "[N2] vol. mixing ratio (mol/mol)",\ |
---|
160 | 48: "[CO] vol. mixing ratio (mol/mol)",\ |
---|
161 | 49: "R: Molecular gas constant (J K-1 kg-1)",\ |
---|
162 | 50: "Air viscosity estimation (N s m-2)" |
---|
163 | } |
---|
164 | if num not in whichfield: myplot.errormess("Incorrect subscript in extvar.") |
---|
165 | |
---|
166 | dastuff = whichfield[num] |
---|
167 | |
---|
168 | if "(K)" in dastuff: self.fmt="%.0f" |
---|
169 | elif "(Pa)" in dastuff: self.fmt="%.1f" |
---|
170 | elif "(W/m2)" in dastuff: self.fmt="%.0f" |
---|
171 | elif "(m/s)" in dastuff: self.fmt="%.1f" |
---|
172 | else: self.fmt="%.2e" |
---|
173 | |
---|
174 | return dastuff |
---|
175 | |
---|
176 | def convertlab(self,num): |
---|
177 | ## a conversion from text inquiries to extvar numbers. to be completed. |
---|
178 | if num == "p": num = 91 |
---|
179 | elif num == "rho": num = 92 |
---|
180 | elif num == "t": num = 93 |
---|
181 | elif num == "u": num = 94 |
---|
182 | elif num == "v": num = 95 |
---|
183 | elif num == "tsurf": num = 15 |
---|
184 | elif num == "topo": num = 4 |
---|
185 | elif num == "h": num = 13 |
---|
186 | elif num == "ps": num = 19 |
---|
187 | elif num == "tau": num = 36 |
---|
188 | elif num == "mtot": num = 40 |
---|
189 | elif num == "icetot": num = 42 |
---|
190 | elif num == "ps_ddv": num = 22 |
---|
191 | elif num == "h2ovap": num = 41 |
---|
192 | elif num == "h2oice": num = 43 |
---|
193 | elif num == "cp": num = 8 |
---|
194 | elif num == "rho_ddv": num = 10 |
---|
195 | elif num == "tsurfmx": num = 16 |
---|
196 | elif num == "tsurfmn": num = 17 |
---|
197 | elif num == "lwdown": num = 31 |
---|
198 | elif num == "swdown": num = 32 |
---|
199 | elif num == "lwup": num = 33 |
---|
200 | elif num == "swup": num = 34 |
---|
201 | elif num == "o3": num = 44 |
---|
202 | elif num == "o": num = 46 |
---|
203 | elif num == "co": num = 48 |
---|
204 | elif num == "visc": num = 50 |
---|
205 | elif num == "co2ice": num = 35 |
---|
206 | elif not isinstance(num, np.int): myplot.errormess("field reference not found.") |
---|
207 | return num |
---|
208 | |
---|
209 | ################### |
---|
210 | ### One request ### |
---|
211 | ################### |
---|
212 | |
---|
213 | def update(self): |
---|
214 | # retrieve fields from MCD (call_mcd). more info in fmcd.call_mcd.__doc__ |
---|
215 | ## sanity first |
---|
216 | self.loct = abs(self.loct)%24 |
---|
217 | if self.locts is not None and self.locte is not None: |
---|
218 | self.locts = abs(self.locts)%24 |
---|
219 | self.locte = abs(self.locte)%24 |
---|
220 | if self.locts == self.locte: self.locte = self.locts + 24 |
---|
221 | if self.lat > 90.: self.lat = 90. |
---|
222 | if self.lat < -90.: self.lat = -90. |
---|
223 | if self.lats is not None and self.late is not None: |
---|
224 | if abs(self.lats) > 90.: self.lats = 90. |
---|
225 | if abs(self.late) > 90.: self.late = 90. |
---|
226 | if abs(self.lats) < -90.: self.lats = -90. |
---|
227 | if abs(self.late) < -90.: self.late = -90. |
---|
228 | ## now MCD request |
---|
229 | (self.pres, self.dens, self.temp, self.zonwind, self.merwind, \ |
---|
230 | self.meanvar, self.extvar, self.seedout, self.ierr) \ |
---|
231 | = \ |
---|
232 | fmcd.call_mcd(self.zkey,self.xz,self.lon,self.lat,self.hrkey, \ |
---|
233 | self.datekey,self.xdate,self.loct,self.dset,self.dust, \ |
---|
234 | self.perturkey,self.seedin,self.gwlength,self.extvarkey ) |
---|
235 | ## we use the end of extvar (unused) to store meanvar. this is convenient for getextvar(lab) |
---|
236 | self.extvar[90] = self.pres ; self.extvar[91] = self.dens |
---|
237 | self.extvar[92] = self.temp ; self.extvar[93] = self.zonwind ; self.extvar[94] = self.merwind |
---|
238 | ## treat missing values |
---|
239 | if self.temp == -999: self.extvar[:] = np.NaN ; self.meanvar[:] = np.NaN |
---|
240 | |
---|
241 | def printset(self): |
---|
242 | # print main settings |
---|
243 | print "zkey",self.zkey,"xz",self.xz,"lon",self.lon,"lat",self.lat,"hrkey",self.hrkey, \ |
---|
244 | "xdate",self.xdate,"loct",self.loct,"dust",self.dust |
---|
245 | |
---|
246 | def getnameset(self): |
---|
247 | # set a name referring to settings [convenient for databases] |
---|
248 | strlat = str(self.lat)+str(self.lats)+str(self.late) |
---|
249 | strlon = str(self.lon)+str(self.lons)+str(self.lone) |
---|
250 | strxz = str(self.xz)+str(self.xzs)+str(self.xze) |
---|
251 | strloct = str(self.loct)+str(self.locts)+str(self.locte) |
---|
252 | name = str(self.zkey)+strxz+strlon+strlat+str(self.hrkey)+str(self.datekey)+str(self.xdate)+strloct+str(self.dust) |
---|
253 | return name |
---|
254 | |
---|
255 | def printcoord(self): |
---|
256 | # print requested space-time coordinates |
---|
257 | print "LAT",self.lat,"LON",self.lon,"LOCT",self.loct,"XDATE",self.xdate |
---|
258 | |
---|
259 | def printmeanvar(self): |
---|
260 | # print mean MCD variables |
---|
261 | print "Pressure = %5.3f pascals. " % (self.pres) |
---|
262 | print "Density = %5.3f kilograms per cubic meter. " % (self.dens) |
---|
263 | print "Temperature = %3.0f kelvins (%4.0f degrees celsius)." % (self.temp,self.temp-273.15) |
---|
264 | print "Zonal wind = %5.3f meters per second." % (self.zonwind) |
---|
265 | print "Meridional wind = %5.3f meters per second." % (self.merwind) |
---|
266 | print "Total horizontal wind = %5.3f meters per second." % ( np.sqrt(self.zonwind**2 + self.merwind**2) ) |
---|
267 | |
---|
268 | def printextvar(self,num): |
---|
269 | # print extra MCD variables |
---|
270 | num = self.convertlab(num) |
---|
271 | print self.getextvarlab(num) + " ..... " + str(self.extvar[num-1]) |
---|
272 | |
---|
273 | def printallextvar(self): |
---|
274 | # print all extra MCD variables |
---|
275 | for i in range(50): self.printextvar(i+1) |
---|
276 | |
---|
277 | def htmlprinttabextvar(self,tabtodo): |
---|
278 | self.gettitle() |
---|
279 | print "<hr>" |
---|
280 | print self.title |
---|
281 | print "<hr>" |
---|
282 | print "<ul>" |
---|
283 | for i in range(len(tabtodo)): print "<li>" ; self.printextvar(tabtodo[i]) ; print "</li>" |
---|
284 | print "</ul>" |
---|
285 | print "<hr>" |
---|
286 | print self.ack |
---|
287 | print "<hr>" |
---|
288 | #print "SETTINGS<br />" |
---|
289 | #self.printcoord() |
---|
290 | #self.printset() |
---|
291 | |
---|
292 | def printmcd(self): |
---|
293 | # 1. call MCD 2. print settings 3. print mean vars |
---|
294 | self.update() |
---|
295 | self.printcoord() |
---|
296 | print "-------------------------------------------" |
---|
297 | self.printmeanvar() |
---|
298 | |
---|
299 | ######################## |
---|
300 | ### Several requests ### |
---|
301 | ######################## |
---|
302 | |
---|
303 | def prepare(self,ndx=None,ndy=None): |
---|
304 | ### prepare I/O arrays for 1d slices |
---|
305 | if ndx is None: print "No dimension in prepare. Exit. Set at least ndx." ; exit() |
---|
306 | else: self.xcoord = np.ones(ndx) |
---|
307 | if ndy is None: dashape = (ndx) ; dashapemean = (ndx,6) ; dashapeext = (ndx,101) ; self.ycoord = None |
---|
308 | else: dashape = (ndx,ndy) ; dashapemean = (ndx,ndy,6) ; dashapeext = (ndx,ndy,101) ; self.ycoord = np.ones(ndy) |
---|
309 | self.prestab = np.ones(dashape) ; self.denstab = np.ones(dashape) ; self.temptab = np.ones(dashape) |
---|
310 | self.zonwindtab = np.ones(dashape) ; self.merwindtab = np.ones(dashape) |
---|
311 | self.meanvartab = np.ones(dashapemean) ; self.extvartab = np.ones(dashapeext) |
---|
312 | |
---|
313 | def getextvar(self,num): |
---|
314 | ### get a given var in extvartab |
---|
315 | try: field=self.extvartab[:,:,num] |
---|
316 | except: field=self.extvartab[:,num] |
---|
317 | return field |
---|
318 | |
---|
319 | def definefield(self,choice): |
---|
320 | ### for analysis or plot purposes, set field and field label from user-defined choice |
---|
321 | choice = self.convertlab(choice) |
---|
322 | field = self.getextvar(choice); fieldlab = self.getextvarlab(choice) |
---|
323 | ## fix for possibly slightly negative tracers |
---|
324 | if "(mol/mol)" in fieldlab or "(kg/kg)" in fieldlab or "(kg/m2)" in fieldlab or "(W/m2)" in fieldlab: |
---|
325 | ind = np.where(field < 1.e-30) |
---|
326 | if ind != -1: field[ind] = 1.e-30 ## 0 does not work everywhere. |
---|
327 | return field,fieldlab |
---|
328 | |
---|
329 | def ininterv(self,dstart,dend,nd,start=None,end=None,yaxis=False,vertcoord=False): |
---|
330 | ### user-defined start and end are used to create xcoord (or ycoord) vector |
---|
331 | if start is not None and end is not None: first, second = self.correctbounds(start,end,vertcoord) |
---|
332 | else: first, second = self.correctbounds(dstart,dend,vertcoord) |
---|
333 | if self.zkey != 4 or not vertcoord: tabtab = np.linspace(first,second,nd) |
---|
334 | else: tabtab = np.logspace(first,second,nd) |
---|
335 | if not yaxis: self.xcoord = tabtab |
---|
336 | else: self.ycoord = tabtab |
---|
337 | |
---|
338 | def correctbounds(self,start,end,vertcoord): |
---|
339 | if self.zkey != 4 or not vertcoord: |
---|
340 | # regular altitudes |
---|
341 | if start > end: first = end ; second = start |
---|
342 | else: first = start ; second = end |
---|
343 | else: |
---|
344 | # pressure: reversed avis |
---|
345 | if start < end: first = np.log10(end) ; second = np.log10(start) |
---|
346 | else: first = np.log10(start) ; second = np.log10(end) |
---|
347 | return first, second |
---|
348 | |
---|
349 | def vertlabel(self): |
---|
350 | if self.zkey == 1: self.xlabel = "radius from centre of planet (m)" |
---|
351 | elif self.zkey == 2: self.xlabel = "height above areoid (m) (MOLA zero datum)" |
---|
352 | elif self.zkey == 3: self.xlabel = "height above surface (m)" |
---|
353 | elif self.zkey == 4: self.xlabel = "pressure level (Pa)" |
---|
354 | elif self.zkey == 5: self.xlabel = "altitude above mean Mars Radius(=3396000m) (m)" |
---|
355 | |
---|
356 | def vertunits(self): |
---|
357 | if self.zkey == 1: self.vunits = "m CP" |
---|
358 | elif self.zkey == 2: self.vunits = "m AMR" |
---|
359 | elif self.zkey == 3: self.vunits = "m ALS" |
---|
360 | elif self.zkey == 4: self.vunits = "Pa" |
---|
361 | elif self.zkey == 5: self.vunits = "m AMMRad" |
---|
362 | |
---|
363 | def vertaxis(self,number,yaxis=False): |
---|
364 | if self.zkey == 2: self.ininterv(-5000.,100000.,number,start=self.xzs,end=self.xze,yaxis=yaxis,vertcoord=True) |
---|
365 | elif self.zkey == 3: self.ininterv(0.,120000.,number,start=self.xzs,end=self.xze,yaxis=yaxis,vertcoord=True) |
---|
366 | elif self.zkey == 5: self.ininterv(-5000.,100000.,number,start=self.xzs,end=self.xze,yaxis=yaxis,vertcoord=True) |
---|
367 | elif self.zkey == 4: self.ininterv(1000.,0.001,number,start=self.xzs,end=self.xze,yaxis=yaxis,vertcoord=True) |
---|
368 | elif self.zkey == 1: self.ininterv(3396000,3596000,number,start=self.xzs,end=self.xze,yaxis=yaxis,vertcoord=True) |
---|
369 | |
---|
370 | ################### |
---|
371 | ### 1D analysis ### |
---|
372 | ################### |
---|
373 | |
---|
374 | def put1d(self,i): |
---|
375 | ## fill in subscript i in output arrays |
---|
376 | ## (arrays must have been correctly defined through prepare) |
---|
377 | if self.prestab is None: myplot.errormess("arrays must be prepared first through self.prepare") |
---|
378 | self.prestab[i] = self.pres ; self.denstab[i] = self.dens ; self.temptab[i] = self.temp |
---|
379 | self.zonwindtab[i] = self.zonwind ; self.merwindtab[i] = self.merwind |
---|
380 | self.meanvartab[i,1:5] = self.meanvar[0:4] ## note: var numbering according to MCD manual is kept |
---|
381 | self.extvartab[i,1:100] = self.extvar[0:99] ## note: var numbering according to MCD manual is kept |
---|
382 | |
---|
383 | def diurnal(self,nd=13): |
---|
384 | ### retrieve a local time slice |
---|
385 | save = self.loct |
---|
386 | self.xlabel = "Local time (Martian hour)" |
---|
387 | self.prepare(ndx=nd) ; self.ininterv(0.,24.,nd,start=self.locts,end=self.locte) |
---|
388 | for i in range(nd): self.loct = self.xcoord[i] ; self.update() ; self.put1d(i) |
---|
389 | self.loct = save |
---|
390 | |
---|
391 | def zonal(self,nd=37): |
---|
392 | ### retrieve a longitude slice |
---|
393 | save = self.lon |
---|
394 | self.xlabel = "East longitude (degrees)" |
---|
395 | self.prepare(ndx=nd) ; self.ininterv(-180.,180.,nd,start=self.lons,end=self.lone) |
---|
396 | for i in range(nd): self.lon = self.xcoord[i] ; self.update() ; self.put1d(i) |
---|
397 | self.lon = save |
---|
398 | |
---|
399 | def meridional(self,nd=19): |
---|
400 | ### retrieve a latitude slice |
---|
401 | save = self.lat |
---|
402 | self.xlabel = "North latitude (degrees)" |
---|
403 | self.prepare(ndx=nd) ; self.ininterv(-90.,90.,nd,start=self.lats,end=self.late) |
---|
404 | for i in range(nd): self.lat = self.xcoord[i] ; self.update() ; self.put1d(i) |
---|
405 | self.lat = save |
---|
406 | |
---|
407 | def profile(self,nd=20,tabperso=None): |
---|
408 | ### retrieve an altitude slice (profile) |
---|
409 | save = self.xz |
---|
410 | self.vertlabel() |
---|
411 | self.vertplot = True |
---|
412 | if tabperso is not None: nd = len(tabperso) |
---|
413 | correct = False |
---|
414 | self.prepare(ndx=nd) ; self.vertaxis(nd) |
---|
415 | if tabperso is not None: self.xcoord = tabperso |
---|
416 | for i in range(nd): self.xz = self.xcoord[i] ; self.update() ; self.put1d(i) |
---|
417 | self.xz = save |
---|
418 | |
---|
419 | def seasonal(self,nd=12): |
---|
420 | ### retrieve a seasonal slice |
---|
421 | save = self.xdate |
---|
422 | self.xlabel = "Areocentric longitude (degrees)" |
---|
423 | self.prepare(ndx=nd) ; self.ininterv(0.,360.,nd,start=self.xdates,end=self.xdatee) |
---|
424 | for i in range(nd): self.xdate = self.xcoord[i] ; self.update() ; self.put1d(i) |
---|
425 | self.xdate = save |
---|
426 | |
---|
427 | def getascii(self,tabtodo,filename="output.txt"): |
---|
428 | ### print out values in an ascii file |
---|
429 | if isinstance(tabtodo,np.str): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
430 | if isinstance(tabtodo,np.int): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
431 | asciifile = open(filename, "w") |
---|
432 | for i in range(len(tabtodo)): |
---|
433 | (field, fieldlab) = self.definefield(tabtodo[i]) |
---|
434 | self.gettitle(oneline=True) |
---|
435 | asciifile.write("### " + self.title + "\n") |
---|
436 | asciifile.write("### " + self.ack + "\n") |
---|
437 | asciifile.write("### Column 1 is " + self.xlabel + "\n") |
---|
438 | asciifile.write("### Column 2 is " + fieldlab + "\n") |
---|
439 | for ix in range(len(self.xcoord)): |
---|
440 | asciifile.write("%15.5e%15.5e\n" % ( self.xcoord[ix], field[ix] ) ) |
---|
441 | asciifile.close() |
---|
442 | return |
---|
443 | |
---|
444 | def makeplot1d(self,choice): |
---|
445 | ### one 1D plot is created for the user-defined variable in choice. |
---|
446 | (field, fieldlab) = self.definefield(choice) |
---|
447 | if not self.vertplot: absc = self.xcoord ; ordo = field ; ordolab = fieldlab ; absclab = self.xlabel |
---|
448 | else: ordo = self.xcoord ; absc = field ; absclab = fieldlab ; ordolab = self.xlabel |
---|
449 | mpl.plot(absc,ordo,'-bo') ; mpl.ylabel(ordolab) ; mpl.xlabel(absclab) #; mpl.xticks(query.xcoord) |
---|
450 | if self.zkey == 4: mpl.semilogy() ; ax = mpl.gca() ; ax.set_ylim(ax.get_ylim()[::-1]) |
---|
451 | mpl.figtext(0.5, 0.01, self.ack, ha='center') |
---|
452 | |
---|
453 | def plot1d(self,tabtodo): |
---|
454 | ### complete 1D figure with possible multiplots |
---|
455 | if isinstance(tabtodo,np.str): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
456 | if isinstance(tabtodo,np.int): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
457 | fig = mpl.figure() ; subv,subh = myplot.definesubplot( len(tabtodo) , fig ) |
---|
458 | for i in range(len(tabtodo)): mpl.subplot(subv,subh,i+1).grid(True, linestyle=':', color='grey') ; self.makeplot1d(tabtodo[i]) |
---|
459 | |
---|
460 | def htmlplot1d(self,tabtodo,figname="temp.png",title=""): |
---|
461 | ### complete 1D figure with possible multiplots |
---|
462 | ### added in 09/2012 for online MCD |
---|
463 | ### see http://www.dalkescientific.com/writings/diary/archive/2005/04/23/matplotlib_without_gui.html |
---|
464 | from matplotlib.figure import Figure |
---|
465 | from matplotlib.backends.backend_agg import FigureCanvasAgg |
---|
466 | if isinstance(tabtodo,np.str): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
467 | if isinstance(tabtodo,np.int): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
468 | |
---|
469 | howmanyplots = len(tabtodo) |
---|
470 | if howmanyplots == 1: fig = Figure(figsize=(16,8)) |
---|
471 | elif howmanyplots == 2: fig = Figure(figsize=(8,8)) |
---|
472 | elif howmanyplots == 3: fig = Figure(figsize=(8,16)) |
---|
473 | elif howmanyplots == 4: fig = Figure(figsize=(16,8)) |
---|
474 | |
---|
475 | subv,subh = myplot.definesubplot( len(tabtodo) , fig ) |
---|
476 | for i in range(len(tabtodo)): |
---|
477 | yeah = fig.add_subplot(subv,subh,i+1) #.grid(True, linestyle=':', color='grey') |
---|
478 | choice = tabtodo[i] |
---|
479 | (field, fieldlab) = self.definefield(choice) |
---|
480 | if not self.vertplot: absc = self.xcoord ; ordo = field ; ordolab = fieldlab ; absclab = self.xlabel |
---|
481 | else: ordo = self.xcoord ; absc = field ; absclab = fieldlab ; ordolab = self.xlabel |
---|
482 | yeah.plot(absc,ordo,'-bo') #; mpl.xticks(query.xcoord) |
---|
483 | ax = fig.gca() ; ax.set_ylabel(ordolab) ; ax.set_xlabel(absclab) |
---|
484 | if self.zkey == 4: ax.set_yscale('log') ; ax.set_ylim(ax.get_ylim()[::-1]) |
---|
485 | |
---|
486 | if self.lats is not None: ax.set_xticks(np.arange(-90,91,15)) ; ax.set_xbound(lower=self.lats, upper=self.late) |
---|
487 | elif self.lons is not None: ax.set_xticks(np.arange(-360,361,30)) ; ax.set_xbound(lower=self.lons, upper=self.lone) |
---|
488 | elif self.locts is not None: ax.set_xticks(np.arange(0,26,2)) ; ax.set_xbound(lower=self.locts, upper=self.locte) |
---|
489 | |
---|
490 | self.gettitle() |
---|
491 | fig.text(0.5, 0.95, self.title, ha='center') |
---|
492 | fig.text(0.5, 0.01, self.ack, ha='center') |
---|
493 | canvas = FigureCanvasAgg(fig) |
---|
494 | # The size * the dpi gives the final image size |
---|
495 | # a4"x4" image * 80 dpi ==> 320x320 pixel image |
---|
496 | canvas.print_figure(figname, dpi=80) |
---|
497 | |
---|
498 | ################### |
---|
499 | ### 2D analysis ### |
---|
500 | ################### |
---|
501 | |
---|
502 | def latlon(self,ndx=37,ndy=19,fixedlt=False): |
---|
503 | ### retrieve a latitude/longitude slice |
---|
504 | ### default is: local time is not fixed. user-defined local time is at longitude 0. |
---|
505 | save1 = self.lon ; save2 = self.lat ; save3 = self.loct |
---|
506 | self.xlabel = "East longitude (degrees)" ; self.ylabel = "North latitude (degrees)" |
---|
507 | self.prepare(ndx=ndx,ndy=ndy) |
---|
508 | self.ininterv(-180.,180.,ndx,start=self.lons,end=self.lone) |
---|
509 | self.ininterv(-90., 90.,ndy,start=self.lats,end=self.late,yaxis=True) |
---|
510 | if not fixedlt: umst = self.loct |
---|
511 | for i in range(ndx): |
---|
512 | for j in range(ndy): |
---|
513 | self.lon = self.xcoord[i] ; self.lat = self.ycoord[j] |
---|
514 | if not fixedlt: self.loct = (umst + self.lon/15.) % 24 |
---|
515 | self.update() ; self.put2d(i,j) |
---|
516 | if not fixedlt: self.loct = umst |
---|
517 | self.lon = save1 ; self.lat = save2 ; self.loct = save3 |
---|
518 | |
---|
519 | def lonalt(self,ndx=37,ndy=20,fixedlt=False): |
---|
520 | ### retrieve a longitude/altitude slice |
---|
521 | save1 = self.lon ; save2 = self.xz ; save3 = self.loct |
---|
522 | self.vertlabel() ; self.ylabel = self.xlabel |
---|
523 | self.xlabel = "East longitude (degrees)" |
---|
524 | self.prepare(ndx=ndx,ndy=ndy) |
---|
525 | self.ininterv(-180.,180.,ndx,start=self.lons,end=self.lone) |
---|
526 | self.vertaxis(ndy,yaxis=True) |
---|
527 | if not fixedlt: umst = self.loct |
---|
528 | for i in range(ndx): |
---|
529 | for j in range(ndy): |
---|
530 | self.lon = self.xcoord[i] ; self.xz = self.ycoord[j] |
---|
531 | if not fixedlt: self.loct = (umst + self.lon/15.) % 24 |
---|
532 | self.update() ; self.put2d(i,j) |
---|
533 | if not fixedlt: self.loct = umst |
---|
534 | self.lon = save1 ; self.xz = save2 ; self.loct = save3 |
---|
535 | |
---|
536 | def latalt(self,ndx=19,ndy=20,fixedlt=False): |
---|
537 | ### retrieve a latitude/altitude slice |
---|
538 | save1 = self.lat ; save2 = self.xz ; save3 = self.loct |
---|
539 | self.vertlabel() ; self.ylabel = self.xlabel |
---|
540 | self.xlabel = "North latitude (degrees)" |
---|
541 | self.prepare(ndx=ndx,ndy=ndy) |
---|
542 | self.ininterv(-90.,90.,ndx,start=self.lats,end=self.late) |
---|
543 | self.vertaxis(ndy,yaxis=True) |
---|
544 | if not fixedlt: umst = self.loct |
---|
545 | for i in range(ndx): |
---|
546 | for j in range(ndy): |
---|
547 | self.lat = self.xcoord[i] ; self.xz = self.ycoord[j] |
---|
548 | if not fixedlt: self.loct = (umst + self.lon/15.) % 24 |
---|
549 | self.update() ; self.put2d(i,j) |
---|
550 | if not fixedlt: self.loct = umst |
---|
551 | self.lat = save1 ; self.xz = save2 ; self.loct = save3 |
---|
552 | |
---|
553 | def put2d(self,i,j): |
---|
554 | ## fill in subscript i,j in output arrays |
---|
555 | ## (arrays must have been correctly defined through prepare) |
---|
556 | if self.prestab is None: myplot.errormess("arrays must be prepared first through self.prepare") |
---|
557 | self.prestab[i,j] = self.pres ; self.denstab[i,j] = self.dens ; self.temptab[i,j] = self.temp |
---|
558 | self.zonwindtab[i,j] = self.zonwind ; self.merwindtab[i,j] = self.merwind |
---|
559 | self.meanvartab[i,j,1:5] = self.meanvar[0:4] ## note: var numbering according to MCD manual is kept |
---|
560 | self.extvartab[i,j,1:100] = self.extvar[0:99] ## note: var numbering according to MCD manual is kept |
---|
561 | |
---|
562 | def makemap2d(self,choice,incwind=False,fixedlt=False,proj="cyl"): |
---|
563 | ### one 2D map is created for the user-defined variable in choice. |
---|
564 | self.latlon(fixedlt=fixedlt) ## a map is implicitely a lat-lon plot. otherwise it is a plot (cf. makeplot2d) |
---|
565 | if choice == "wind" or incwind: |
---|
566 | (windx, fieldlabwx) = self.definefield("u") |
---|
567 | (windy, fieldlabwy) = self.definefield("v") |
---|
568 | if choice == "wind": |
---|
569 | field = np.sqrt(windx*windx + windy*windy) |
---|
570 | fieldlab = "Horizontal wind speed (m/s)" |
---|
571 | else: |
---|
572 | (field, fieldlab) = self.definefield(choice) |
---|
573 | if incwind: myplot.maplatlon(self.xcoord,self.ycoord,field,title=fieldlab,proj=proj,vecx=windx,vecy=windy) #,stride=1) |
---|
574 | else: myplot.maplatlon(self.xcoord,self.ycoord,field,title=fieldlab,proj=proj) |
---|
575 | mpl.figtext(0.5, 0.0, self.ack, ha='center') |
---|
576 | |
---|
577 | def map2d(self,tabtodo,incwind=False,fixedlt=False,proj="cyl"): |
---|
578 | ### complete 2D figure with possible multiplots |
---|
579 | if isinstance(tabtodo,np.str): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
580 | if isinstance(tabtodo,np.int): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
581 | fig = mpl.figure() |
---|
582 | subv,subh = myplot.definesubplot( len(tabtodo) , fig ) |
---|
583 | for i in range(len(tabtodo)): mpl.subplot(subv,subh,i+1) ; self.makemap2d(tabtodo[i],incwind=incwind,fixedlt=fixedlt,proj=proj) |
---|
584 | |
---|
585 | def htmlmap2d(self,tabtodo,incwind=False,fixedlt=False,figname="temp.png",title="",back="zMOL"): |
---|
586 | ### complete 2D figure with possible multiplots |
---|
587 | ### added in 09/2012 for online MCD |
---|
588 | ### see http://www.dalkescientific.com/writings/diary/archive/2005/04/23/matplotlib_without_gui.html |
---|
589 | from matplotlib.figure import Figure |
---|
590 | from matplotlib.backends.backend_agg import FigureCanvasAgg |
---|
591 | from matplotlib.cm import get_cmap |
---|
592 | from matplotlib import rcParams |
---|
593 | |
---|
594 | #from mpl_toolkits.basemap import Basemap |
---|
595 | |
---|
596 | from Scientific.IO import NetCDF |
---|
597 | filename = "/home/marshttp/surface.nc" |
---|
598 | zefile = NetCDF.NetCDFFile(filename, 'r') |
---|
599 | fieldc = zefile.variables[back] |
---|
600 | yc = zefile.variables['latitude'] |
---|
601 | xc = zefile.variables['longitude'] |
---|
602 | ## plutot que fieldc = self.getextvar(self.convertlab("topo")) |
---|
603 | |
---|
604 | if isinstance(tabtodo,np.str): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
605 | if isinstance(tabtodo,np.int): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
606 | |
---|
607 | howmanyplots = len(tabtodo) |
---|
608 | if howmanyplots == 1: fig = Figure(figsize=(16,8)) |
---|
609 | elif howmanyplots == 2: fig = Figure(figsize=(8,8)) |
---|
610 | elif howmanyplots == 3: fig = Figure(figsize=(8,16)) |
---|
611 | elif howmanyplots == 4: fig = Figure(figsize=(16,8)) |
---|
612 | |
---|
613 | subv,subh = myplot.definesubplot( len(tabtodo) , fig ) |
---|
614 | |
---|
615 | for i in range(len(tabtodo)): |
---|
616 | yeah = fig.add_subplot(subv,subh,i+1) |
---|
617 | choice = tabtodo[i] |
---|
618 | self.latlon(fixedlt=fixedlt,ndx=64,ndy=48) |
---|
619 | ## a map is implicitely a lat-lon plot. otherwise it is a plot (cf. makeplot2d) |
---|
620 | (field, fieldlab) = self.definefield(choice) |
---|
621 | if incwind: (windx, fieldlabwx) = self.definefield("u") ; (windy, fieldlabwy) = self.definefield("v") |
---|
622 | |
---|
623 | proj="moll" ; colorb="jet" ; ndiv=20 ; zeback="molabw" ; trans=1.0 #0.6 |
---|
624 | title="" ; vecx=None ; vecy=None ; stride=2 |
---|
625 | lon = self.xcoord |
---|
626 | lat = self.ycoord |
---|
627 | |
---|
628 | #[lon2d,lat2d] = np.meshgrid(lon,lat) |
---|
629 | ##### define projection and background. define x and y given the projection |
---|
630 | ##[wlon,wlat] = myplot.latinterv() |
---|
631 | ##yeahm = myplot.define_proj(proj,wlon,wlat,back=zeback,blat=None,blon=None) |
---|
632 | ##x, y = yeahm(lon2d, lat2d) |
---|
633 | #map = Basemap(projection='ortho',lat_0=45,lon_0=-100) |
---|
634 | #x, y = map(lon2d, lat2d) |
---|
635 | |
---|
636 | #### TEMP |
---|
637 | x = lon ; y = lat |
---|
638 | |
---|
639 | ## define field. bound field. |
---|
640 | what_I_plot = np.transpose(field) |
---|
641 | zevmin, zevmax = myplot.calculate_bounds(what_I_plot) ## vmin=min(what_I_plot_frame), vmax=max(what_I_plot_frame)) |
---|
642 | what_I_plot = myplot.bounds(what_I_plot,zevmin,zevmax) |
---|
643 | ## define contour field levels. define color palette |
---|
644 | ticks = ndiv + 1 |
---|
645 | zelevels = np.linspace(zevmin,zevmax,ticks) |
---|
646 | palette = get_cmap(name=colorb) |
---|
647 | |
---|
648 | # You can set negative contours to be solid instead of dashed: |
---|
649 | rcParams['contour.negative_linestyle'] = 'solid' |
---|
650 | ## contours topo |
---|
651 | zelevc = np.linspace(-9.,20.,11) |
---|
652 | yeah.contour( xc, yc, fieldc, zelevc, colors='black',linewidths = 0.4) |
---|
653 | yeah.contour( np.array(xc) + 360., yc, fieldc, zelevc, colors='black',linewidths = 0.4) |
---|
654 | yeah.contour( np.array(xc) - 360., yc, fieldc, zelevc, colors='black',linewidths = 0.4) |
---|
655 | # contour field |
---|
656 | c = yeah.contourf( x, y, what_I_plot, zelevels, cmap = palette, alpha = trans ) |
---|
657 | clb = Figure.colorbar(fig,c,orientation='vertical',format=self.fmt,ticks=np.linspace(zevmin,zevmax,num=min([ticks/2+1,21]))) |
---|
658 | clb.set_label(fieldlab) |
---|
659 | if incwind: |
---|
660 | [x2d,y2d] = np.meshgrid(x,y) |
---|
661 | yeah.quiver(x2d,y2d,np.transpose(windx),np.transpose(windy)) |
---|
662 | ax = fig.gca() ; ax.set_ylabel("Latitude") ; ax.set_xlabel("Longitude") |
---|
663 | ax.set_xticks(np.arange(-360,361,45)) ; ax.set_xbound(lower=self.lons, upper=self.lone) |
---|
664 | ax.set_yticks(np.arange(-90,91,30)) ; ax.set_ybound(lower=self.lats, upper=self.late) |
---|
665 | self.gettitle() |
---|
666 | fig.text(0.5, 0.95, self.title, ha='center') |
---|
667 | fig.text(0.5, 0.01, self.ack, ha='center') |
---|
668 | canvas = FigureCanvasAgg(fig) |
---|
669 | # The size * the dpi gives the final image size |
---|
670 | # a4"x4" image * 80 dpi ==> 320x320 pixel image |
---|
671 | canvas.print_figure(figname, dpi=80) |
---|
672 | |
---|
673 | def htmlplot2d(self,tabtodo,fixedlt=False,figname="temp.png",title=""): |
---|
674 | ### complete 2D figure with possible multiplots |
---|
675 | ### added in 10/2012 for online MCD |
---|
676 | ### see http://www.dalkescientific.com/writings/diary/archive/2005/04/23/matplotlib_without_gui.html |
---|
677 | from matplotlib.figure import Figure |
---|
678 | from matplotlib.backends.backend_agg import FigureCanvasAgg |
---|
679 | from matplotlib.cm import get_cmap |
---|
680 | if isinstance(tabtodo,np.str): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
681 | if isinstance(tabtodo,np.int): tabtodo=[tabtodo] ## so that asking one element without [] is possible. |
---|
682 | |
---|
683 | howmanyplots = len(tabtodo) |
---|
684 | if howmanyplots == 1: fig = Figure(figsize=(16,8)) |
---|
685 | elif howmanyplots == 2: fig = Figure(figsize=(8,8)) |
---|
686 | elif howmanyplots == 3: fig = Figure(figsize=(8,16)) |
---|
687 | elif howmanyplots == 4: fig = Figure(figsize=(16,8)) |
---|
688 | |
---|
689 | subv,subh = myplot.definesubplot( len(tabtodo) , fig ) |
---|
690 | |
---|
691 | for i in range(len(tabtodo)): |
---|
692 | yeah = fig.add_subplot(subv,subh,i+1) |
---|
693 | choice = tabtodo[i] |
---|
694 | |
---|
695 | if self.lons is not None: self.lonalt(fixedlt=fixedlt,ndx=64,ndy=35) |
---|
696 | elif self.lats is not None: self.latalt(fixedlt=fixedlt,ndx=48,ndy=35) |
---|
697 | |
---|
698 | (field, fieldlab) = self.definefield(choice) |
---|
699 | |
---|
700 | colorb="jet" ; ndiv=20 ; title="" |
---|
701 | |
---|
702 | ## define field. bound field. |
---|
703 | what_I_plot = np.transpose(field) |
---|
704 | zevmin, zevmax = myplot.calculate_bounds(what_I_plot) ## vmin=min(what_I_plot_frame), vmax=max(what_I_plot_frame)) |
---|
705 | what_I_plot = myplot.bounds(what_I_plot,zevmin,zevmax) |
---|
706 | ## define contour field levels. define color palette |
---|
707 | ticks = ndiv + 1 |
---|
708 | zelevels = np.linspace(zevmin,zevmax,ticks) |
---|
709 | palette = get_cmap(name=colorb) |
---|
710 | # contour field |
---|
711 | c = yeah.contourf( self.xcoord, self.ycoord, what_I_plot, zelevels, cmap = palette ) |
---|
712 | clb = Figure.colorbar(fig,c,orientation='vertical',format=self.fmt,ticks=np.linspace(zevmin,zevmax,num=min([ticks/2+1,21]))) |
---|
713 | clb.set_label(fieldlab) |
---|
714 | ax = fig.gca() ; ax.set_ylabel(self.ylabel) ; ax.set_xlabel(self.xlabel) |
---|
715 | |
---|
716 | if self.zkey == 4: |
---|
717 | ax.set_yscale('log') ; ax.set_ylim(ax.get_ylim()[::-1]) |
---|
718 | else: |
---|
719 | #ax.set_yticks(np.arange(self.xzs,self.xze,10000.)) ; |
---|
720 | ax.set_ybound(lower=self.xzs, upper=self.xze) |
---|
721 | |
---|
722 | if self.lons is not None: ax.set_xticks(np.arange(-360,361,45)) ; ax.set_xbound(lower=self.lons, upper=self.lone) |
---|
723 | elif self.lats is not None: ax.set_xticks(np.arange(-90,91,30)) ; ax.set_xbound(lower=self.lats, upper=self.late) |
---|
724 | |
---|
725 | self.gettitle() |
---|
726 | fig.text(0.5, 0.95, self.title, ha='center') |
---|
727 | fig.text(0.5, 0.01, self.ack, ha='center') |
---|
728 | canvas = FigureCanvasAgg(fig) |
---|
729 | # The size * the dpi gives the final image size |
---|
730 | # a4"x4" image * 80 dpi ==> 320x320 pixel image |
---|
731 | canvas.print_figure(figname, dpi=80) |
---|
732 | |
---|
733 | |
---|
734 | ### TODO: makeplot2d, plot2d, passer plot settings |
---|
735 | |
---|