Package astLib :: Module astSED
[hide private]
[frames] | no frames]

Source Code for Module astLib.astSED

  1  # -*- coding: utf-8 -*- 
  2  """module for performing calculations on Spectral Energy Distributions (SEDs) 
  3   
  4  (c) 2007-2009 Matt Hilton  
  5   
  6  U{http://astlib.sourceforge.net} 
  7   
  8  This module provides classes for manipulating SEDs, in particular the Bruzual & Charlot 2003 and Maraston 
  9  2005 stellar population synthesis models are currently supported. Functions are provided for calculating  
 10  the evolution of broadband colours and magnitudes in these models with redshift etc.. It is not especially 
 11  fast at present. 
 12   
 13  @var VEGA: The SED of Vega, used for calculation of magnitudes on the Vega system. 
 14  @type VEGA: astSED.SED object 
 15  @var AB: Flat spectrum SED, used for calculation of magnitudes on the AB system. 
 16  @type AB: astSED.SED object 
 17   
 18  """ 
 19   
 20  #------------------------------------------------------------------------------------------------------------ 
 21  import sys 
 22  import numpy 
 23  import math 
 24  try: 
 25      from scipy import interpolate 
 26      from scipy import ndimage 
 27  except: 
 28      print "WARNING: astSED: failed to import scipy modules - some functions will not work." 
 29  import astLib 
 30  from astLib import astCalc 
 31  import os 
 32  try: 
 33      import matplotlib 
 34      from matplotlib import pylab 
 35      matplotlib.interactive(False) 
 36  except: 
 37      print "WARNING: astSED: failed to import matplotlib - some functions will not work." 
 38       
 39  #------------------------------------------------------------------------------------------------------------ 
40 -class Passband:
41 """This class describes a filter transmission curve. Passband objects are created by loading data from 42 from text files containing wavelength in angstroms in the first column, relative transmission efficiency 43 in the second column (whitespace delimited). For example, to create a Passband object for the 2MASS J 44 filter: 45 46 passband=astSED.Passband("J_2MASS.res") 47 48 where "J_2MASS.res" is a file in the current working directory that describes the filter. 49 50 Wavelength units can be specified as 'angstroms', 'nanometres' or 'microns'; if either of the latter, 51 they will be converted to angstroms. 52 53 """
54 - def __init__(self, fileName, renormalise = True, inputUnits = 'angstroms'):
55 56 inFile=file(fileName, "rb") 57 lines=inFile.readlines() 58 59 wavelength=[] 60 transmission=[] 61 for line in lines: 62 63 if line[0] != "#" and len(line) > 3: 64 65 bits=line.split() 66 transmission.append(float(bits[1])) 67 wavelength.append(float(bits[0])) 68 69 self.wavelength=numpy.array(wavelength) 70 self.transmission=numpy.array(transmission) 71 72 if inputUnits == 'nanometres': 73 self.wavelength=self.wavelength*10.0 74 elif inputUnits == 'microns': 75 self.wavelength=self.wavelength*10000.0 76 77 if renormalise == True: 78 self.transmission=self.transmission/numpy.trapz(self.transmission, self.wavelength) 79 80 # Store a ready-to-go interpolation object to speed calculation of fluxes up 81 self.interpolator=interpolate.interp1d(self.wavelength, self.transmission, kind='linear')
82
83 - def asList(self):
84 """Returns a two dimensional list of [wavelength, transmission], suitable for plotting by gnuplot. 85 86 @rtype: list 87 @return: list in format [wavelength, transmission] 88 89 """ 90 91 listData=[] 92 for l, f in zip(self.wavelength, self.transmission): 93 listData.append([l, f]) 94 95 return listData
96
97 - def rescale(self, maxTransmission):
98 """Rescales the passband so that maximum value of the transmission is equal to maxTransmission. 99 Useful for plotting. 100 101 @type maxTransmission: float 102 @param maxTransmission: maximum value of rescaled transmission curve 103 104 """ 105 106 self.transmission=self.transmission*(maxTransmission/self.transmission.max())
107
108 - def plot(self, xmin = 'min', xmax = 'max', maxTransmission = None):
109 """Plots the passband, rescaling the maximum of the tranmission curve to maxTransmission if 110 required. 111 112 @type xmin: float or 'min' 113 @param xmin: minimum of the wavelength range of the plot 114 @type xmax: float or 'max' 115 @param xmax: maximum of the wavelength range of the plot 116 @type maxTransmission: float 117 @param maxTransmission: maximum value of rescaled transmission curve 118 119 """ 120 121 if maxTransmission != None: 122 self.rescale(maxTransmission) 123 124 pylab.matplotlib.interactive(True) 125 pylab.plot(self.wavelength, self.transmission) 126 127 if xmin == 'min': 128 xmin=self.wavelength.min() 129 if xmax == 'max': 130 xmax=self.wavelength.max() 131 132 pylab.xlim(xmin, xmax) 133 pylab.xlabel("Wavelength") 134 pylab.ylabel("Relative Flux")
135 136 #------------------------------------------------------------------------------------------------------------
137 -class SED:
138 """This class describes a Spectral Energy Distribution (SED). 139 140 To create a SED object, lists (or numpy arrays) of wavelength and relative flux must be provided. The SED 141 can optionally be redshifted. Typically, the units of SEDs are Angstroms - flux calculations using 142 Passband and SED objects specified with different wavelength units will be incorrect. 143 144 The L{StellarPopulation} class (and derivatives) can be used to extract SEDs for specified ages from e.g. 145 the Bruzual & Charlot 2003 or Maraston 2005 models. 146 147 """ 148
149 - def __init__(self, wavelength = [], flux = [], z = 0.0000001, ageGyr = None, renormalise = False, label = None):
150 151 # We keep a copy of the wavelength, flux at z = 0, as it's more robust to copy that 152 # to self.wavelength, flux and redshift it, rather than repeatedly redshifting the same 153 # arrays back and forth 154 self.z0wavelength=numpy.array(wavelength) 155 self.z0flux=numpy.array(flux) 156 self.wavelength=numpy.array(wavelength) 157 self.flux=numpy.array(flux) 158 self.z=0.0000001 159 self.label=label # plain text label, handy for using in photo-z codes 160 161 if renormalise == True: 162 self.renormalise() 163 164 if z != 0.0: 165 self.redshift(z) 166 167 self.ageGyr=ageGyr
168
169 - def loadFromFile(self, fileName):
170 """Loads SED from a white space delimited file in the format wavelength, flux. Lines beginning with 171 # are ignored. 172 173 @type fileName: string 174 @param fileName: path to file containing wavelength, flux data 175 176 """ 177 178 inFile=file(fileName, "rb") 179 lines=inFile.readlines() 180 inFile.close() 181 wavelength=[] 182 flux=[] 183 wholeLines=[] 184 for line in lines: 185 if line[0] != "#" and len(line) > 3: 186 bits=line.split() 187 wavelength.append(float(bits[0])) 188 flux.append(float(bits[1])) 189 190 # Sort SED so wavelength is in ascending order 191 if wavelength[0] > wavelength[-1]: 192 wavelength.reverse() 193 flux.reverse() 194 195 self.z0wavelength=numpy.array(wavelength) 196 self.z0flux=numpy.array(flux) 197 self.wavelength=numpy.array(wavelength) 198 self.flux=numpy.array(flux)
199
200 - def writeToFile(self, fileName):
201 """Writes SED to a white space delimited file in the format wavelength, flux. 202 203 @type fileName: string 204 @param fileName: path to file 205 206 """ 207 208 outFile=file(fileName, "wb") 209 for l, f in zip(self.wavelength, self.flux): 210 outFile.write(str(l)+" "+str(f)+"\n") 211 outFile.close()
212
213 - def asList(self):
214 """Returns a two dimensional list of [wavelength, flux], suitable for plotting by gnuplot. 215 216 @rtype: list 217 @return: list in format [wavelength, flux] 218 219 """ 220 221 listData=[] 222 for l, f in zip(self.wavelength, self.flux): 223 listData.append([l, f]) 224 225 return listData
226
227 - def plot(self, xmin = 'min', xmax = 'max'):
228 """Produces a simple (wavelength, flux) plot of the SED. 229 230 @type xmin: float or 'min' 231 @param xmin: minimum of the wavelength range of the plot 232 @type xmax: float or 'max' 233 @param xmax: maximum of the wavelength range of the plot 234 235 """ 236 237 pylab.matplotlib.interactive(True) 238 pylab.plot(self.wavelength, self.flux) 239 240 if xmin == 'min': 241 xmin=self.wavelength.min() 242 if xmax == 'max': 243 xmax=self.wavelength.max() 244 245 pylab.xlim(xmin, xmax) 246 pylab.xlabel("Wavelength") 247 pylab.ylabel("Relative Flux")
248
249 - def smooth(self, smoothPix):
250 """Smooths SED.flux with a uniform (boxcar) filter of width smoothPix. Cannot be undone. 251 252 @type smoothPix: int 253 254 """ 255 smoothed=ndimage.uniform_filter1d(self.flux, smoothPix) 256 self.flux=smoothed
257
258 - def redshift(self, z):
259 """Redshifts the SED to redshift z. 260 261 @type z: float 262 @param z: redshift 263 264 """ 265 266 # We have to conserve energy so the area under the redshifted SED has to be equal to 267 # the area under the unredshifted SED, otherwise magnitude calculations will be wrong 268 # when comparing SEDs at different zs 269 self.wavelength=self.z0wavelength 270 self.flux=self.z0flux 271 272 z0TotalFlux=numpy.trapz(self.z0wavelength, self.z0flux) 273 self.wavelength=self.wavelength*(1.0+z) 274 zTotalFlux=numpy.trapz(self.wavelength, self.flux) 275 self.flux=self.flux*(z0TotalFlux/zTotalFlux) 276 self.z=z
277
278 - def renormalise(self, minWavelength = 'min', maxWavelength = 'max'):
279 """Normalises the SED such that the area under the specified wavelength range is equal to 1. 280 281 @type minWavelength: float or 'min' 282 @param minWavelength: minimum wavelength of range over which to normalise SED 283 @type maxWavelength: float or 'max' 284 @param maxWavelength: maximum wavelength of range over which to normalise SED 285 286 """ 287 if minWavelength == 'min': 288 minWavelength=self.wavelength.min() 289 if maxWavelength == 'max': 290 maxWavelength=self.wavelength.max() 291 292 lowCut=numpy.greater(self.wavelength, minWavelength) 293 highCut=numpy.less(self.wavelength, maxWavelength) 294 totalCut=numpy.logical_and(lowCut, highCut) 295 sedFluxSlice=self.flux[totalCut] 296 sedWavelengthSlice=self.wavelength[totalCut] 297 298 self.flux=self.flux/numpy.trapz(abs(sedFluxSlice), sedWavelengthSlice)#self.wavelength)
299
300 - def matchFlux(self, matchSED, minWavelength, maxWavelength):
301 """Matches the flux in the wavelength range given by minWavelength, maxWavelength to the 302 flux in the same region in matchSED. Useful for plotting purposes. 303 304 @type matchSED: astSED.SED object 305 @param matchSED: SED to match flux to 306 @type minWavelength: float 307 @param minWavelength: minimum of range in which to match flux of current SED to matchSED 308 @type maxWavelength: float 309 @param maxWavelength: maximum of range in which to match flux of current SED to matchSED 310 311 """ 312 313 interpMatch=interpolate.interp1d(matchSED.wavelength, matchSED.flux, kind='linear') 314 interpSelf=interpolate.interp1d(self.wavelength, self.flux, kind='linear') 315 316 wavelengthRange=numpy.arange(minWavelength, maxWavelength, 5.0) 317 318 matchFlux=numpy.trapz(interpMatch(wavelengthRange), wavelengthRange) 319 selfFlux=numpy.trapz(interpSelf(wavelengthRange), wavelengthRange) 320 321 self.flux=self.flux*(matchFlux/selfFlux)
322 323
324 - def calcFlux(self, passband):
325 """Calculates relative flux in the given passband. 326 327 @type passband: astSED.Passband object 328 @param passband: filter passband through which to calculate the flux from the SED 329 330 """ 331 lowCut=numpy.greater(self.wavelength, passband.wavelength.min()) 332 highCut=numpy.less(self.wavelength, passband.wavelength.max()) 333 totalCut=numpy.logical_and(lowCut, highCut) 334 sedFluxSlice=self.flux[totalCut] 335 sedWavelengthSlice=self.wavelength[totalCut] 336 337 # Use linear interpolation to rebin the passband to the same dimensions as the 338 # part of the SED we're interested in 339 sedInBand=passband.interpolator(sedWavelengthSlice)*sedFluxSlice 340 341 totalFlux=numpy.trapz(sedInBand, sedWavelengthSlice) # sum using trapezoid rule 342 343 return totalFlux
344
345 - def calcAbsMag(self, passband, magType = "Vega"):
346 """Calculates the relative absolute magnitude in the given passband. By 'relative', we mean 347 that we calculate -2.5*log10*flux, relative to either the Vega or AB magnitude systems; i.e. 348 we do not assume the flux is coming from an object 10 pc away (as would be the case for a 349 true absolute magnitude), as calculated magnitudes need to be normalised to observations 350 in any case for them to be meaningful. 351 352 @type passband: astSED.Passband object 353 @param passband: filter passband through which to calculate the magnitude from the SED 354 @type magType: string 355 @param magType: either "Vega" or "AB" 356 @rtype: float 357 @return: relative absolute magnitude on the specified magnitude system 358 359 """ 360 f1=self.calcFlux(passband) 361 if magType == "Vega": 362 f2=VEGA.calcFlux(passband) 363 elif magType == "AB": 364 f2=AB.calcFlux(passband) 365 366 mag=-2.5*math.log10(f1/f2) 367 368 return mag
369
370 - def calcAppMag(self, passband, magType = "Vega"):
371 """Calculates the relative apparent magnitude in the given passband. By 'relative', we 372 mean that we calculate the apparent mag as 5.0*log10*(dl*1e5)-2.5*log10*flux, where 373 flux is relative to either the Vega or AB systems, and dl is the luminosity distance 374 in Mpc at the redshift of the SED. Calculated magnitudes need to be normalised to 375 observations to make them physically meaningful. 376 377 @type passband: astSED.Passband object 378 @param passband: filter passband through which to calculate the magnitude from the SED 379 @type magType: string 380 @param magType: either "Vega" or "AB" 381 @rtype: float 382 @return: relative absolute magnitude on the specified magnitude system 383 384 """ 385 f1=self.calcFlux(passband) 386 if magType == "Vega": 387 f2=VEGA.calcFlux(passband) 388 elif magType == "AB": 389 f2=AB.calcFlux(passband) 390 391 mag=-2.5*math.log10(f1/f2) 392 393 appMag=5.0*math.log10(astCalc.dl(self.z)*1e5)+mag 394 395 return appMag
396
397 - def calcColour(self, passband1, passband2, magType = "Vega"):
398 """Calculates the colour passband1-passband2. 399 400 @type passband1: astSED.Passband object 401 @param passband1: filter passband through which to calculate the first magnitude 402 @type passband2: astSED.Passband object 403 @param passband1: filter passband through which to calculate the second magnitude 404 @type magType: string 405 @param magType: either "Vega" or "AB" 406 @rtype: float 407 @return: colour defined by passband1 - passband2 on the specified magnitude system 408 409 """ 410 mag1=self.calcAppMag(passband1, magType = magType) 411 mag2=self.calcAppMag(passband2, magType = magType) 412 413 colour=mag1-mag2 414 return colour
415 416 #------------------------------------------------------------------------------------------------------------
417 -class VegaSED(SED):
418 """This class stores the SED of Vega, used for calculation of magnitudes on the Vega system. 419 420 The Vega SED used is taken from Bohlin 2007 (http://adsabs.harvard.edu/abs/2007ASPC..364..315B), and is 421 available from the STScI CALSPEC library (http://www.stsci.edu/hst/observatory/cdbs/calspec.html). 422 423 """ 424 425 # Vega SED has to be renormalised if other SEDs are not, 426 # otherwise can get overflow problems where dividing a big (10e29) number by a small one (10e-10) 427 # such as when calculating Vega mags?
428 - def __init__(self, renormalise = True):
429 430 VEGA_SED_PATH=astLib.__path__[0]+os.path.sep+"data"+os.path.sep+"bohlin2006_Vega.sed" # from HST CALSPEC 431 432 inFile=file(VEGA_SED_PATH, "rb") 433 lines=inFile.readlines() 434 435 wavelength=[] 436 flux=[] 437 for line in lines: 438 439 if line[0] != "#" and len(line) > 3: 440 441 bits=line.split() 442 flux.append(float(bits[1])) 443 wavelength.append(float(bits[0])) 444 445 self.wavelength=numpy.array(wavelength) 446 self.flux=numpy.array(flux) 447 448 # We may want to redshift reference SEDs to calculate rest-frame colors from SEDs at different zs 449 self.z0wavelength=numpy.array(wavelength) 450 self.z0flux=numpy.array(flux) 451 self.z=0.0000001 452 453 if renormalise == True: 454 self.flux=self.flux/numpy.trapz(self.flux, self.wavelength) 455 self.z0flux=self.z0flux/numpy.trapz(self.z0flux, self.z0wavelength)
456 457 #------------------------------------------------------------------------------------------------------------
458 -class StellarPopulation:
459 """This class describes a stellar population model, either a Simple Stellar Population (SSP) or a 460 Composite Stellar Population (CSP), such as the models of Bruzual & Charlot 2003 or Maraston 2005. 461 462 The constructor for this class can be used for generic SSPs or CSPs stored in white space delimited text 463 files, containing columns for age, wavelength, and flux. Columns are counted from 0 ... n. Lines starting 464 with # are ignored. 465 466 The classes L{M05Model} (for Maraston 2005 models), L{BC03Model} (for Bruzual & Charlot 2003 models) are 467 derived from this class. The only difference between them is the code used to load in the model data. 468 469 """
470 - def __init__(self, fileName, ageColumn = 0, wavelengthColumn = 1, fluxColumn = 2):
471 472 inFile=file(fileName, "rb") 473 lines=inFile.readlines() 474 inFile.close() 475 476 # Extract a list of model ages and valid wavelengths from the file 477 self.ages=[] 478 self.wavelengths=[] 479 for line in lines: 480 if line[0] !="#" and len(line) > 3: 481 bits=line.split() 482 age=float(bits[ageColumn]) 483 wavelength=float(bits[wavelengthColumn]) 484 if age not in self.ages: 485 self.ages.append(age) 486 if wavelength not in self.wavelengths: 487 self.wavelengths.append(wavelength) 488 489 # Construct a grid of flux - rows correspond to each wavelength, columns to age 490 self.fluxGrid=numpy.zeros([len(self.ages), len(self.wavelengths)]) 491 for line in lines: 492 if line[0] !="#" and len(line) > 3: 493 bits=line.split() 494 sedAge=float(bits[ageColumn]) 495 sedWavelength=float(bits[wavelengthColumn]) 496 sedFlux=float(bits[fluxColumn]) 497 498 row=self.ages.index(sedAge) 499 column=self.wavelengths.index(sedWavelength) 500 501 self.fluxGrid[row][column]=sedFlux
502
503 - def getSED(self, ageGyr, z = 0.0, renormalise = False, label = None):
504 """Extract a SED for given age. Do linear interpolation between models if necessary. 505 506 @type ageGyr: float 507 @param ageGyr: age of the SED in Gyr 508 @type z: float 509 @param z: redshift the SED from z = 0 to z = z 510 @type renormalise: bool 511 @param renormalise: normalise the SED to have area 1 512 @rtype: astSED.SED object 513 @return: SED 514 515 """ 516 517 if ageGyr in self.ages: 518 519 flux=self.fluxGrid[self.ages.index(ageGyr)] 520 sed=SED(self.wavelengths, flux, z = z, renormalise = renormalise, label = label) 521 return sed 522 523 else: 524 525 # Use interpolation, iterating over each wavelength column 526 flux=[] 527 for i in range(len(self.wavelengths)): 528 interpolator=interpolate.interp1d(self.ages, self.fluxGrid[:,i], kind='linear') 529 sedFlux=interpolator(ageGyr)[0] 530 flux.append(sedFlux) 531 sed=SED(self.wavelengths, flux, z = z, renormalise = renormalise, label = label) 532 return sed
533
534 - def getColourEvolution(self, passband1, passband2, zFormation, zStepSize = 0.05, magType = "Vega"):
535 """Calculates the evolution of the colour observed through passband1 - passband2 for the 536 StellarPopulation with redshift, from z = 0 to z = zFormation. 537 538 @type passband1: astSED.Passband object 539 @param passband1: filter passband through which to calculate the first magnitude 540 @type passband2: astSED.Passband object 541 @param passband2: filter passband through which to calculate the second magnitude 542 @type zFormation: float 543 @param zFormation: formation redshift of the StellarPopulation 544 @type zStepSize: float 545 @param zStepSize: size of interval in z at which to calculate model colours 546 @type magType: string 547 @param magType: either "Vega" or "AB" 548 @rtype: dictionary 549 @return: dictionary of numpy.arrays in format {'z', 'colour'} 550 551 """ 552 553 zSteps=int(math.ceil(zFormation/zStepSize)) 554 zData=[] 555 colourData=[] 556 for i in range(1, zSteps): 557 zc=i*zStepSize 558 age=astCalc.tl(zFormation)-astCalc.tl(zc) 559 sed=self.getSED(age, z = zc) 560 colour=sed.calcColour(passband1, passband2, magType = magType) 561 zData.append(zc) 562 colourData.append(colour) 563 564 zData=numpy.array(zData) 565 colourData=numpy.array(colourData) 566 567 return {'z': zData, 'colour': colourData}
568
569 - def getMagEvolution(self, passband, magNormalisation, zNormalisation, zFormation, zStepSize = 0.05, 570 onePlusZSteps = False, magType = "Vega"):
571 """Calculates the evolution with redshift (from z = 0 to z = zFormation) of apparent magnitude 572 in the observed frame through the passband for the StellarPopulation, normalised to magNormalisation 573 (apparent) at z = zNormalisation. 574 575 @type passband: astSED.Passband object 576 @param passband: filter passband through which to calculate the magnitude 577 @type magNormalisation: float 578 @param magNormalisation: sets the apparent magnitude of the SED at zNormalisation 579 @type zNormalisation: float 580 @param zNormalisation: the redshift at which the magnitude normalisation is carried out 581 @type zFormation: float 582 @param zFormation: formation redshift of the StellarPopulation 583 @type zStepSize: float 584 @param zStepSize: size of interval in z at which to calculate model magnitudes 585 @type onePlusZSteps: bool 586 @param onePlusZSteps: if True, zSteps are (1+z)*zStepSize, otherwise zSteps are linear 587 @type magType: string 588 @param magType: either "Vega" or "AB" 589 @rtype: dictionary 590 @return: dictionary of numpy.arrays in format {'z', 'mag'} 591 592 """ 593 594 # Count upwards in z steps as interpolation doesn't work if array ordered z decreasing 595 zSteps=int(math.ceil(zFormation/zStepSize)) 596 zData=[] 597 magData=[] 598 absMagData=[] 599 zc0=0.0 600 for i in range(1, zSteps): 601 if onePlusZSteps == False: 602 zc=i*zStepSize 603 else: 604 zc=zc0+(1+zc0)*zStepSize 605 zc0=zc 606 if zc >= zFormation: 607 break 608 age=astCalc.tl(zFormation)-astCalc.tl(zc) 609 sed=self.getSED(age, z = zc) 610 mag=sed.calcAppMag(passband, magType = magType) 611 zData.append(zc) 612 magData.append(mag) 613 absMagData.append(sed.calcAbsMag(passband)) 614 615 zData=numpy.array(zData) 616 magData=numpy.array(magData) 617 618 # Do the normalisation 619 interpolator=interpolate.interp1d(zData, magData, kind='linear') 620 modelNormMag=interpolator(zNormalisation)[0] 621 normConstant=magNormalisation-modelNormMag 622 magData=magData+normConstant 623 624 return {'z': zData, 'mag': magData}
625
626 - def calcEvolutionCorrection(self, zFrom, zTo, zFormation, passband, magType = "Vega"):
627 """Calculates the evolution correction in magnitudes in the rest frame through the passband 628 from redshift zFrom to redshift zTo, where the stellarPopulation is assumed to be formed 629 at redshift zFormation. 630 631 @type zFrom: float 632 @param zFormation: redshift to evolution correct from 633 @type zTo: float 634 @param zTo: redshift to evolution correct to 635 @type zFormation: float 636 @param zFormation: formation redshift of the StellarPopulation 637 @type passband: astSED.Passband object 638 @param passband: filter passband through which to calculate magnitude 639 @type magType: string 640 @param magType: either "Vega" or "AB" 641 @rtype: float 642 @return: evolution correction in magnitudes in the rest frame 643 644 """ 645 646 ageFrom=astCalc.tl(zFormation)-astCalc.tl(zFrom) 647 ageTo=astCalc.tl(zFormation)-astCalc.tl(zTo) 648 649 fromSED=self.getSED(ageFrom) 650 toSED=self.getSED(ageTo) 651 652 fromMag=fromSED.calcAbsMag(passband, magType = magType) 653 toMag=toSED.calcAbsMag(passband, magType = magType) 654 655 return fromMag-toMag
656 657 #------------------------------------------------------------------------------------------------------------
658 -class M05Model(StellarPopulation):
659 """This class describes a Maraston 2005 stellar population model. To load a composite stellar population 660 model (CSP) for a tau = 0.1 Gyr burst of star formation, solar metallicity, Salpeter IMF: 661 662 m05csp = astSED.M05Model(M05_DIR+"/csp_e_0.10_z02_salp.sed_agb") 663 664 where M05_DIR is set to point to the directory where the Maraston 2005 models are stored on your system. 665 666 The file format of the Maraston 2005 simple stellar poulation (SSP) models is different to the file 667 format used for the CSPs, and this needs to be specified using the fileType parameter. To load a SSP with 668 solar metallicity, red horizontal branch morphology: 669 670 m05ssp = astSED.M05Model(M05_DIR+"/sed.ssz002.rhb", fileType = "ssp") 671 672 """
673 - def __init__(self, fileName, fileType = "csp"):
674 675 inFile=file(fileName, "rb") 676 lines=inFile.readlines() 677 inFile.close() 678 679 if fileType == "csp": 680 ageColumn=0 681 wavelengthColumn=1 682 fluxColumn=2 683 elif fileType == "ssp": 684 ageColumn=0 685 wavelengthColumn=2 686 fluxColumn=3 687 else: 688 raise Exception, "fileType must be 'ssp' or 'csp'" 689 690 # Extract a list of model ages and valid wavelengths from the file 691 self.ages=[] 692 self.wavelengths=[] 693 for line in lines: 694 if line[0] !="#" and len(line) > 3: 695 bits=line.split() 696 age=float(bits[ageColumn]) 697 wavelength=float(bits[wavelengthColumn]) 698 if age not in self.ages: 699 self.ages.append(age) 700 if wavelength not in self.wavelengths: 701 self.wavelengths.append(wavelength) 702 703 # Construct a grid of flux - rows correspond to each wavelength, columns to age 704 self.fluxGrid=numpy.zeros([len(self.ages), len(self.wavelengths)]) 705 for line in lines: 706 if line[0] !="#" and len(line) > 3: 707 bits=line.split() 708 sedAge=float(bits[ageColumn]) 709 sedWavelength=float(bits[wavelengthColumn]) 710 sedFlux=float(bits[fluxColumn]) 711 712 row=self.ages.index(sedAge) 713 column=self.wavelengths.index(sedWavelength) 714 715 self.fluxGrid[row][column]=sedFlux
716 717 #------------------------------------------------------------------------------------------------------------
718 -class BC03Model(StellarPopulation):
719 """This class describes a Bruzual & Charlot 2003 stellar population model, extracted from a GALAXEV .ised 720 file using the galaxevpl program that is included in GALAXEV. The file format is white space delimited, 721 with wavelength in the first column. Subsequent columns contain the model fluxes for SEDs of different 722 ages, as specified when running galaxevpl. The age corresponding to each flux column is taken from the 723 comment line beginning "# Age (yr)", and is converted to Gyr. 724 725 For example, to load a tau = 0.1 Gyr burst of star formation, solar metallicity, Salpeter IMF model 726 stored in a file (created by galaxevpl) called "csp_lr_solar_0p1Gyr.136": 727 728 bc03model = BC03Model("csp_lr_solar_0p1Gyr.136") 729 730 """ 731
732 - def __init__(self, fileName):
733 734 inFile=file(fileName, "rb") 735 lines=inFile.readlines() 736 inFile.close() 737 738 # Extract a list of model ages - BC03 ages are in years, so convert to Gyr 739 self.ages=[] 740 for line in lines: 741 if line.find("# Age (yr)") != -1: 742 rawAges=line[line.find("# Age (yr)")+10:].split() 743 for age in rawAges: 744 self.ages.append(float(age)/1e9) 745 break 746 747 # Extract a list of valid wavelengths from the file 748 self.wavelengths=[] 749 for line in lines: 750 if line[0] !="#" and len(line) > 3: 751 bits=line.split() 752 self.wavelengths.append(float(bits[0])) 753 754 # Construct a grid of flux - rows correspond to each wavelength, columns to age 755 self.fluxGrid=numpy.zeros([len(self.ages), len(self.wavelengths)]) 756 for line in lines: 757 if line[0] !="#" and len(line) > 3: 758 bits=line.split() 759 ageFluxes=bits[1:] 760 sedWavelength=float(bits[0]) 761 column=self.wavelengths.index(sedWavelength) 762 763 for row in range(len(ageFluxes)): 764 sedFlux=float(ageFluxes[row]) 765 self.fluxGrid[row][column]=sedFlux
766 767 #------------------------------------------------------------------------------------------------------------ 768 # Data 769 VEGA=VegaSED() 770 771 # Normalisation factor for Vega and AB SEDs, so can mix them and work out mag differences 772 _norm=3652.935 # for the Bohlin 2006 Vega SED, mean over several bands, uncertainty ~1 per cent 773 AB=VegaSED() 774 AB.flux=numpy.ones(len(AB.wavelength)) 775 AB.flux=(AB.flux/(AB.wavelength*AB.wavelength)) 776 AB.flux=AB.flux*_norm 777 AB.z0flux=AB.flux[:] 778