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