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