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