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-2011 Matt Hilton  
   5   
   6  U{http://astlib.sourceforge.net} 
   7   
   8  This module provides classes for manipulating SEDs, in particular the Bruzual & Charlot 2003, Maraston 
   9  2005, and Percival et al 2009 stellar population synthesis models are currently supported. Functions are  
  10  provided for calculating the evolution of colours and magnitudes in these models with redshift etc., and  
  11  for fitting broadband photometry using these models. 
  12   
  13  @var VEGA: The SED of Vega, used for calculation of magnitudes on the Vega system. 
  14  @type VEGA: L{SED} object 
  15  @var AB: Flat spectrum SED, used for calculation of magnitudes on the AB system. 
  16  @type AB: L{SED} object 
  17  @var SOL: The SED of the Sun. 
  18  @type SOL: L{SED} object 
  19   
  20  """ 
  21   
  22  #------------------------------------------------------------------------------------------------------------ 
  23  import sys 
  24  import numpy 
  25  import math 
  26  import operator 
  27  try: 
  28      from scipy import interpolate 
  29      from scipy import ndimage 
  30      from scipy import optimize 
  31  except: 
  32      print "WARNING: astSED: failed to import scipy modules - some functions will not work." 
  33  import astLib 
  34  from astLib import astCalc 
  35  import os 
  36  try: 
  37      import matplotlib 
  38      from matplotlib import pylab 
  39      matplotlib.interactive(False) 
  40  except: 
  41      print "WARNING: astSED: failed to import matplotlib - some functions will not work." 
  42  import glob 
  43   
  44  #------------------------------------------------------------------------------------------------------------ 
45 -class Passband:
46 """This class describes a filter transmission curve. Passband objects are created by loading data from 47 from text files containing wavelength in angstroms in the first column, relative transmission efficiency 48 in the second column (whitespace delimited). For example, to create a Passband object for the 2MASS J 49 filter: 50 51 passband=astSED.Passband("J_2MASS.res") 52 53 where "J_2MASS.res" is a file in the current working directory that describes the filter. 54 55 Wavelength units can be specified as 'angstroms', 'nanometres' or 'microns'; if either of the latter, 56 they will be converted to angstroms. 57 58 """
59 - def __init__(self, fileName, normalise = True, inputUnits = 'angstroms'):
60 61 inFile=file(fileName, "rb") 62 lines=inFile.readlines() 63 64 wavelength=[] 65 transmission=[] 66 for line in lines: 67 68 if line[0] != "#" and len(line) > 3: 69 70 bits=line.split() 71 transmission.append(float(bits[1])) 72 wavelength.append(float(bits[0])) 73 74 self.wavelength=numpy.array(wavelength) 75 self.transmission=numpy.array(transmission) 76 77 if inputUnits == 'angstroms': 78 pass 79 elif inputUnits == 'nanometres': 80 self.wavelength=self.wavelength*10.0 81 elif inputUnits == 'microns': 82 self.wavelength=self.wavelength*10000.0 83 elif inputUnits == 'mm': 84 self.wavelength=self.wavelength*1e7 85 elif inputUnits == 'GHz': 86 self.wavelength=3e8/(self.wavelength*1e9) 87 self.wavelength=self.wavelength*1e10 88 else: 89 raise Exception, "didn't understand passband input units" 90 91 # Sort into ascending order of wavelength otherwise normalisation will be wrong 92 merged=numpy.array([self.wavelength, self.transmission]).transpose() 93 sortedMerged=numpy.array(sorted(merged, key=operator.itemgetter(0))) 94 self.wavelength=sortedMerged[:, 0] 95 self.transmission=sortedMerged[:, 1] 96 97 if normalise == True: 98 self.transmission=self.transmission/numpy.trapz(self.transmission, self.wavelength) 99 100 # Store a ready-to-go interpolation object to speed calculation of fluxes up 101 self.interpolator=interpolate.interp1d(self.wavelength, self.transmission, kind='linear')
102
103 - def asList(self):
104 """Returns a two dimensional list of [wavelength, transmission], suitable for plotting by gnuplot. 105 106 @rtype: list 107 @return: list in format [wavelength, transmission] 108 109 """ 110 111 listData=[] 112 for l, f in zip(self.wavelength, self.transmission): 113 listData.append([l, f]) 114 115 return listData
116
117 - def rescale(self, maxTransmission):
118 """Rescales the passband so that maximum value of the transmission is equal to maxTransmission. 119 Useful for plotting. 120 121 @type maxTransmission: float 122 @param maxTransmission: maximum value of rescaled transmission curve 123 124 """ 125 126 self.transmission=self.transmission*(maxTransmission/self.transmission.max())
127
128 - def plot(self, xmin = 'min', xmax = 'max', maxTransmission = None):
129 """Plots the passband, rescaling the maximum of the tranmission curve to maxTransmission if 130 required. 131 132 @type xmin: float or 'min' 133 @param xmin: minimum of the wavelength range of the plot 134 @type xmax: float or 'max' 135 @param xmax: maximum of the wavelength range of the plot 136 @type maxTransmission: float 137 @param maxTransmission: maximum value of rescaled transmission curve 138 139 """ 140 141 if maxTransmission != None: 142 self.rescale(maxTransmission) 143 144 pylab.matplotlib.interactive(True) 145 pylab.plot(self.wavelength, self.transmission) 146 147 if xmin == 'min': 148 xmin=self.wavelength.min() 149 if xmax == 'max': 150 xmax=self.wavelength.max() 151 152 pylab.xlim(xmin, xmax) 153 pylab.xlabel("Wavelength") 154 pylab.ylabel("Relative Flux")
155
156 - def effectiveWavelength(self):
157 """Calculates effective wavelength for the passband. This is the same as equation (3) of 158 Carter et al. 2009. 159 160 @rtype: float 161 @return: effective wavelength of the passband, in Angstroms 162 163 """ 164 165 a=numpy.trapz(self.transmission*self.wavelength) 166 b=numpy.trapz(self.transmission/self.wavelength) 167 effWavelength=numpy.sqrt(a/b) 168 169 return effWavelength
170 171 #------------------------------------------------------------------------------------------------------------
172 -class TopHatPassband(Passband):
173 """This class generates a passband with a top hat response between the given wavelengths. 174 175 """ 176
177 - def __init__(self, wavelengthMin, wavelengthMax, normalise = True):
178 """Generates a passband object with top hat response between wavelengthMin, wavelengthMax. 179 Units are assumed to be Angstroms. 180 181 @type wavelengthMin: float 182 @param wavelengthMin: minimum of the wavelength range of the passband 183 @type wavelengthMax: float 184 @param wavelengthMax: maximum of the wavelength range of the passband 185 @type normalise: bool 186 @param normalise: if True, scale such that total area under the passband over the wavelength 187 range is 1. 188 189 """ 190 191 self.wavelength=numpy.arange(wavelengthMin, wavelengthMax+10, 10, dtype = float) 192 self.transmission=numpy.ones(self.wavelength.shape, dtype = float) 193 194 if normalise == True: 195 self.transmission=self.transmission/numpy.trapz(self.transmission, self.wavelength) 196 197 # Store a ready-to-go interpolation object to speed calculation of fluxes up 198 self.interpolator=interpolate.interp1d(self.wavelength, self.transmission, kind='linear')
199 200 201 #------------------------------------------------------------------------------------------------------------
202 -class SED:
203 """This class describes a Spectral Energy Distribution (SED). 204 205 To create a SED object, lists (or numpy arrays) of wavelength and relative flux must be provided. The SED 206 can optionally be redshifted. The wavelength units of SEDs are assumed to be Angstroms - flux 207 calculations using Passband and SED objects specified with different wavelength units will be incorrect. 208 209 The L{StellarPopulation} class (and derivatives) can be used to extract SEDs for specified ages from e.g. 210 the Bruzual & Charlot 2003 or Maraston 2005 models. 211 212 """ 213
214 - def __init__(self, wavelength = [], flux = [], z = 0.0, ageGyr = None, normalise = False, label = None):
215 216 # We keep a copy of the wavelength, flux at z = 0, as it's more robust to copy that 217 # to self.wavelength, flux and redshift it, rather than repeatedly redshifting the same 218 # arrays back and forth 219 self.z0wavelength=numpy.array(wavelength) 220 self.z0flux=numpy.array(flux) 221 self.wavelength=numpy.array(wavelength) 222 self.flux=numpy.array(flux) 223 self.z=z 224 self.label=label # plain text label, handy for using in photo-z codes 225 226 # Store the intrinsic (i.e. unextincted) flux in case we change extinction 227 self.EBMinusV=0.0 228 self.intrinsic_z0flux=numpy.array(flux) 229 230 if normalise == True: 231 self.normalise() 232 233 if z != 0.0: 234 self.redshift(z) 235 236 self.ageGyr=ageGyr
237 238
239 - def copy(self):
240 """Copies the SED, returning a new SED object 241 242 @rtype: L{SED} object 243 @return: SED 244 245 """ 246 247 newSED=SED(wavelength = self.z0wavelength, flux = self.z0flux, z = self.z, ageGyr = self.ageGyr, 248 normalise = False, label = self.label) 249 250 return newSED
251 252
253 - def loadFromFile(self, fileName):
254 """Loads SED from a white space delimited file in the format wavelength, flux. Lines beginning with 255 # are ignored. 256 257 @type fileName: string 258 @param fileName: path to file containing wavelength, flux data 259 260 """ 261 262 inFile=file(fileName, "rb") 263 lines=inFile.readlines() 264 inFile.close() 265 wavelength=[] 266 flux=[] 267 wholeLines=[] 268 for line in lines: 269 if line[0] != "#" and len(line) > 3: 270 bits=line.split() 271 wavelength.append(float(bits[0])) 272 flux.append(float(bits[1])) 273 274 # Sort SED so wavelength is in ascending order 275 if wavelength[0] > wavelength[-1]: 276 wavelength.reverse() 277 flux.reverse() 278 279 self.z0wavelength=numpy.array(wavelength) 280 self.z0flux=numpy.array(flux) 281 self.wavelength=numpy.array(wavelength) 282 self.flux=numpy.array(flux)
283
284 - def writeToFile(self, fileName):
285 """Writes SED to a white space delimited file in the format wavelength, flux. 286 287 @type fileName: string 288 @param fileName: path to file 289 290 """ 291 292 outFile=file(fileName, "wb") 293 for l, f in zip(self.wavelength, self.flux): 294 outFile.write(str(l)+" "+str(f)+"\n") 295 outFile.close()
296
297 - def asList(self):
298 """Returns a two dimensional list of [wavelength, flux], suitable for plotting by gnuplot. 299 300 @rtype: list 301 @return: list in format [wavelength, flux] 302 303 """ 304 305 listData=[] 306 for l, f in zip(self.wavelength, self.flux): 307 listData.append([l, f]) 308 309 return listData
310
311 - def plot(self, xmin = 'min', xmax = 'max'):
312 """Produces a simple (wavelength, flux) plot of the SED. 313 314 @type xmin: float or 'min' 315 @param xmin: minimum of the wavelength range of the plot 316 @type xmax: float or 'max' 317 @param xmax: maximum of the wavelength range of the plot 318 319 """ 320 321 pylab.matplotlib.interactive(True) 322 pylab.plot(self.wavelength, self.flux) 323 324 if xmin == 'min': 325 xmin=self.wavelength.min() 326 if xmax == 'max': 327 xmax=self.wavelength.max() 328 329 # Sensible y scale 330 plotMask=numpy.logical_and(numpy.greater(self.wavelength, xmin), numpy.less(self.wavelength, xmax)) 331 plotMax=self.flux[plotMask].max() 332 pylab.ylim(0, plotMax*1.1) 333 pylab.xlim(xmin, xmax) 334 pylab.xlabel("Wavelength") 335 pylab.ylabel("Relative Flux")
336
337 - def integrate(self, wavelengthMin = 'min', wavelengthMax = 'max'):
338 """Calculates flux in SED within given wavelength range. 339 340 @type wavelengthMin: float or 'min' 341 @param wavelengthMin: minimum of the wavelength range 342 @type wavelengthMax: float or 'max' 343 @param wavelengthMax: maximum of the wavelength range 344 @rtype: float 345 @return: relative flux 346 347 """ 348 349 if wavelengthMin == 'min': 350 wavelengthMin=self.wavelength.min() 351 if wavelengthMax == 'max': 352 wavelengthMax=self.wavelength.max() 353 354 mask=numpy.logical_and(numpy.greater(self.wavelength, wavelengthMin), \ 355 numpy.less(self.wavelength, wavelengthMax)) 356 flux=numpy.trapz(self.flux[mask], self.wavelength[mask]) 357 358 return flux
359
360 - def smooth(self, smoothPix):
361 """Smooths SED.flux with a uniform (boxcar) filter of width smoothPix. Cannot be undone. 362 363 @type smoothPix: int 364 @param smoothPix: size of uniform filter applied to SED, in pixels 365 366 """ 367 smoothed=ndimage.uniform_filter1d(self.flux, smoothPix) 368 self.flux=smoothed
369
370 - def redshift(self, z):
371 """Redshifts the SED to redshift z. 372 373 @type z: float 374 @param z: redshift 375 376 """ 377 378 # We have to conserve energy so the area under the redshifted SED has to be equal to 379 # the area under the unredshifted SED, otherwise magnitude calculations will be wrong 380 # when comparing SEDs at different zs 381 self.wavelength=numpy.zeros(self.z0wavelength.shape[0]) 382 self.flux=numpy.zeros(self.z0flux.shape[0]) 383 self.wavelength=self.wavelength+self.z0wavelength 384 self.flux=self.flux+self.z0flux 385 386 z0TotalFlux=numpy.trapz(self.z0wavelength, self.z0flux) 387 self.wavelength=self.wavelength*(1.0+z) 388 zTotalFlux=numpy.trapz(self.wavelength, self.flux) 389 self.flux=self.flux*(z0TotalFlux/zTotalFlux) 390 self.z=z
391
392 - def normalise(self, minWavelength = 'min', maxWavelength = 'max'):
393 """Normalises the SED such that the area under the specified wavelength range is equal to 1. 394 395 @type minWavelength: float or 'min' 396 @param minWavelength: minimum wavelength of range over which to normalise SED 397 @type maxWavelength: float or 'max' 398 @param maxWavelength: maximum wavelength of range over which to normalise SED 399 400 """ 401 if minWavelength == 'min': 402 minWavelength=self.wavelength.min() 403 if maxWavelength == 'max': 404 maxWavelength=self.wavelength.max() 405 406 lowCut=numpy.greater(self.wavelength, minWavelength) 407 highCut=numpy.less(self.wavelength, maxWavelength) 408 totalCut=numpy.logical_and(lowCut, highCut) 409 sedFluxSlice=self.flux[totalCut] 410 sedWavelengthSlice=self.wavelength[totalCut] 411 412 self.flux=self.flux/numpy.trapz(abs(sedFluxSlice), sedWavelengthSlice)#self.wavelength)
413
414 - def normaliseToMag(self, ABMag, passband):
415 """Normalises the SED to match the flux equivalent to the given AB magnitude in the given passband. 416 417 @type ABMag: float 418 @param ABMag: AB magnitude to which the SED is to be normalised at the given passband 419 @type passband: an L{Passband} object 420 @param passband: passband at which normalisation to AB magnitude is calculated 421 422 """ 423 424 magFlux=mag2Flux(ABMag, 0.0, passband) 425 sedFlux=self.calcFlux(passband) 426 norm=magFlux[0]/sedFlux 427 self.flux=self.flux*norm 428 self.z0flux=self.z0flux*norm
429
430 - def matchFlux(self, matchSED, minWavelength, maxWavelength):
431 """Matches the flux in the wavelength range given by minWavelength, maxWavelength to the 432 flux in the same region in matchSED. Useful for plotting purposes. 433 434 @type matchSED: L{SED} object 435 @param matchSED: SED to match flux to 436 @type minWavelength: float 437 @param minWavelength: minimum of range in which to match flux of current SED to matchSED 438 @type maxWavelength: float 439 @param maxWavelength: maximum of range in which to match flux of current SED to matchSED 440 441 """ 442 443 interpMatch=interpolate.interp1d(matchSED.wavelength, matchSED.flux, kind='linear') 444 interpSelf=interpolate.interp1d(self.wavelength, self.flux, kind='linear') 445 446 wavelengthRange=numpy.arange(minWavelength, maxWavelength, 5.0) 447 448 matchFlux=numpy.trapz(interpMatch(wavelengthRange), wavelengthRange) 449 selfFlux=numpy.trapz(interpSelf(wavelengthRange), wavelengthRange) 450 451 self.flux=self.flux*(matchFlux/selfFlux)
452 453
454 - def calcFlux(self, passband):
455 """Calculates flux in the given passband. 456 457 @type passband: L{Passband} object 458 @param passband: filter passband through which to calculate the flux from the SED 459 @rtype: float 460 @return: flux 461 462 """ 463 lowCut=numpy.greater(self.wavelength, passband.wavelength.min()) 464 highCut=numpy.less(self.wavelength, passband.wavelength.max()) 465 totalCut=numpy.logical_and(lowCut, highCut) 466 sedFluxSlice=self.flux[totalCut] 467 sedWavelengthSlice=self.wavelength[totalCut] 468 469 # Use linear interpolation to rebin the passband to the same dimensions as the 470 # part of the SED we're interested in 471 sedInBand=passband.interpolator(sedWavelengthSlice)*sedFluxSlice 472 totalFlux=numpy.trapz(sedInBand*sedWavelengthSlice, sedWavelengthSlice) 473 totalFlux=totalFlux/numpy.trapz(passband.interpolator(sedWavelengthSlice)\ 474 *sedWavelengthSlice, sedWavelengthSlice) 475 476 return totalFlux
477
478 - def calcMag(self, passband, addDistanceModulus = True, magType = "Vega"):
479 """Calculates magnitude in the given passband. If addDistanceModulus == True, 480 then the distance modulus (5.0*log10*(dl*1e5), where dl is the luminosity distance 481 in Mpc at the redshift of the L{SED}) is added. 482 483 @type passband: L{Passband} object 484 @param passband: filter passband through which to calculate the magnitude from the SED 485 @type addDistanceModulus: bool 486 @param addDistanceModulus: if True, adds 5.0*log10*(dl*1e5) to the mag returned, where 487 dl is the luminosity distance (Mpc) corresponding to the SED z 488 @type magType: string 489 @param magType: either "Vega" or "AB" 490 @rtype: float 491 @return: magnitude through the given passband on the specified magnitude system 492 493 """ 494 f1=self.calcFlux(passband) 495 if magType == "Vega": 496 f2=VEGA.calcFlux(passband) 497 elif magType == "AB": 498 f2=AB.calcFlux(passband) 499 500 mag=-2.5*math.log10(f1/f2) 501 if magType == "Vega": 502 mag=mag+0.026 # Add 0.026 because Vega has V=0.026 (e.g. Bohlin & Gilliland 2004) 503 504 if self.z > 0.0 and addDistanceModulus == True: 505 appMag=5.0*math.log10(astCalc.dl(self.z)*1e5)+mag 506 else: 507 appMag=mag 508 509 return appMag
510
511 - def calcColour(self, passband1, passband2, magType = "Vega"):
512 """Calculates the colour passband1-passband2. 513 514 @type passband1: L{Passband} object 515 @param passband1: filter passband through which to calculate the first magnitude 516 @type passband2: L{Passband} object 517 @param passband1: filter passband through which to calculate the second magnitude 518 @type magType: string 519 @param magType: either "Vega" or "AB" 520 @rtype: float 521 @return: colour defined by passband1 - passband2 on the specified magnitude system 522 523 """ 524 mag1=self.calcMag(passband1, magType = magType, addDistanceModulus = True) 525 mag2=self.calcMag(passband2, magType = magType, addDistanceModulus = True) 526 527 colour=mag1-mag2 528 return colour
529
530 - def getSEDDict(self, passbands):
531 """This is a convenience function for pulling out fluxes from a SED for a given set of passbands 532 in the same format as made by L{mags2SEDDict} - designed to make fitting code simpler. 533 534 @type passbands: list of L{Passband} objects 535 @param passbands: list of passbands through which fluxes will be calculated 536 537 """ 538 539 flux=[] 540 wavelength=[] 541 for p in passbands: 542 flux.append(self.calcFlux(p)) 543 wavelength.append(p.effectiveWavelength()) 544 545 SEDDict={} 546 SEDDict['flux']=numpy.array(flux) 547 SEDDict['wavelength']=numpy.array(wavelength) 548 549 return SEDDict
550
551 - def extinctionCalzetti(self, EBMinusV):
552 """Applies the Calzetti et al. 2000 (ApJ, 533, 682) extinction law to the SED with the given 553 E(B-V) amount of extinction. R_v' = 4.05 is assumed (see equation (5) of Calzetti et al.). 554 555 @type EBMinusV: float 556 @param EBMinusV: extinction E(B-V), in magnitudes 557 558 """ 559 560 self.EBMinusV=EBMinusV 561 562 # All done in rest frame 563 self.z0flux=self.intrinsic_z0flux 564 565 # Allow us to set EBMinusV == 0 to turn extinction off 566 if EBMinusV > 0: 567 # Note that EBMinusV is assumed to be Es as in equations (2) - (5) 568 # Note here wavelength units have to be microns for constants to make sense 569 RvPrime=4.05 # equation (5) of Calzetti et al. 2000 570 shortWavelengthMask=numpy.logical_and(numpy.greater_equal(self.z0wavelength, 1200), \ 571 numpy.less(self.z0wavelength, 6300)) 572 longWavelengthMask=numpy.logical_and(numpy.greater_equal(self.z0wavelength, 6300), \ 573 numpy.less_equal(self.z0wavelength, 22000)) 574 wavelengthMicrons=numpy.array(self.z0wavelength/10000.0, dtype=numpy.float64) 575 kPrime=numpy.zeros(self.z0wavelength.shape[0], dtype=numpy.float64) 576 kPrimeLong=(2.659*(-1.857 \ 577 +1.040/wavelengthMicrons \ 578 ))+RvPrime 579 kPrimeShort=(2.659*(-2.156 \ 580 +1.509/wavelengthMicrons \ 581 -0.198/wavelengthMicrons**2 \ 582 +0.011/wavelengthMicrons**3 \ 583 ))+RvPrime 584 kPrime[longWavelengthMask]=kPrimeLong[longWavelengthMask] 585 kPrime[shortWavelengthMask]=kPrimeShort[shortWavelengthMask] 586 587 # Here we extrapolate kPrime in similar way to what HYPERZ does 588 # Short wavelengths 589 try: 590 interpolator=interpolate.interp1d(self.z0wavelength, kPrimeShort, kind='linear') 591 slope=(interpolator(1100.0)-interpolator(1200.0))/(1100.0-1200.0) 592 intercept=interpolator(1200.0)-(slope*1200.0) 593 mask=numpy.less(self.z0wavelength, 1200.0) 594 kPrime[mask]=slope*self.z0wavelength[mask]+intercept 595 596 # Long wavelengths 597 interpolator=interpolate.interp1d(self.z0wavelength, kPrimeLong, kind='linear') 598 slope=(interpolator(21900.0)-interpolator(22000.0))/(21900.0-22000.0) 599 intercept=interpolator(21900.0)-(slope*21900.0) 600 mask=numpy.greater(self.z0wavelength, 22000.0) 601 kPrime[mask]=slope*self.z0wavelength[mask]+intercept 602 except: 603 raise Exception, "This SED has a wavelength range that doesn't cover ~1200-22000 Angstroms" 604 605 # Never let go negative 606 kPrime[numpy.less_equal(kPrime, 0.0)]=1e-6 607 608 reddening=numpy.power(10, 0.4*EBMinusV*kPrime) 609 self.z0flux=self.z0flux/reddening 610 611 self.redshift(self.z)
612 613 #------------------------------------------------------------------------------------------------------------
614 -class VegaSED(SED):
615 """This class stores the SED of Vega, used for calculation of magnitudes on the Vega system. 616 617 The Vega SED used is taken from Bohlin 2007 (http://adsabs.harvard.edu/abs/2007ASPC..364..315B), and is 618 available from the STScI CALSPEC library (http://www.stsci.edu/hst/observatory/cdbs/calspec.html). 619 620 """ 621
622 - def __init__(self, normalise = False):
623 624 VEGA_SED_PATH=astLib.__path__[0]+os.path.sep+"data"+os.path.sep+"bohlin2006_Vega.sed" # from HST CALSPEC 625 626 inFile=file(VEGA_SED_PATH, "rb") 627 lines=inFile.readlines() 628 629 wavelength=[] 630 flux=[] 631 for line in lines: 632 633 if line[0] != "#" and len(line) > 3: 634 635 bits=line.split() 636 flux.append(float(bits[1])) 637 wavelength.append(float(bits[0])) 638 639 self.wavelength=numpy.array(wavelength) 640 self.flux=numpy.array(flux, dtype=numpy.float64) 641 642 # We may want to redshift reference SEDs to calculate rest-frame colors from SEDs at different zs 643 self.z0wavelength=numpy.array(wavelength) 644 self.z0flux=numpy.array(flux, dtype=numpy.float64) 645 self.z=0.0
646 647 #if normalise == True: 648 #self.flux=self.flux/numpy.trapz(self.flux, self.wavelength) 649 #self.z0flux=self.z0flux/numpy.trapz(self.z0flux, self.z0wavelength) 650 651 #------------------------------------------------------------------------------------------------------------
652 -class StellarPopulation:
653 """This class describes a stellar population model, either a Simple Stellar Population (SSP) or a 654 Composite Stellar Population (CSP), such as the models of Bruzual & Charlot 2003 or Maraston 2005. 655 656 The constructor for this class can be used for generic SSPs or CSPs stored in white space delimited text 657 files, containing columns for age, wavelength, and flux. Columns are counted from 0 ... n. Lines starting 658 with # are ignored. 659 660 The classes L{M05Model} (for Maraston 2005 models), L{BC03Model} (for Bruzual & Charlot 2003 models), and 661 L{P09Model} (for Percival et al. 2009 models) are derived from this class. The only difference between 662 them is the code used to load in the model data. 663 664 """
665 - def __init__(self, fileName, ageColumn = 0, wavelengthColumn = 1, fluxColumn = 2):
666 667 inFile=file(fileName, "rb") 668 lines=inFile.readlines() 669 inFile.close() 670 671 self.fileName=fileName 672 673 # Extract a list of model ages and valid wavelengths from the file 674 self.ages=[] 675 self.wavelengths=[] 676 for line in lines: 677 if line[0] !="#" and len(line) > 3: 678 bits=line.split() 679 age=float(bits[ageColumn]) 680 wavelength=float(bits[wavelengthColumn]) 681 if age not in self.ages: 682 self.ages.append(age) 683 if wavelength not in self.wavelengths: 684 self.wavelengths.append(wavelength) 685 686 # Construct a grid of flux - rows correspond to each wavelength, columns to age 687 self.fluxGrid=numpy.zeros([len(self.ages), len(self.wavelengths)]) 688 for line in lines: 689 if line[0] !="#" and len(line) > 3: 690 bits=line.split() 691 sedAge=float(bits[ageColumn]) 692 sedWavelength=float(bits[wavelengthColumn]) 693 sedFlux=float(bits[fluxColumn]) 694 695 row=self.ages.index(sedAge) 696 column=self.wavelengths.index(sedWavelength) 697 698 self.fluxGrid[row][column]=sedFlux
699
700 - def getSED(self, ageGyr, z = 0.0, normalise = False, label = None):
701 """Extract a SED for given age. Do linear interpolation between models if necessary. 702 703 @type ageGyr: float 704 @param ageGyr: age of the SED in Gyr 705 @type z: float 706 @param z: redshift the SED from z = 0 to z = z 707 @type normalise: bool 708 @param normalise: normalise the SED to have area 1 709 @rtype: L{SED} object 710 @return: SED 711 712 """ 713 714 if ageGyr in self.ages: 715 716 flux=self.fluxGrid[self.ages.index(ageGyr)] 717 sed=SED(self.wavelengths, flux, z = z, normalise = normalise, label = label) 718 return sed 719 720 else: 721 722 # Use interpolation, iterating over each wavelength column 723 flux=[] 724 for i in range(len(self.wavelengths)): 725 interpolator=interpolate.interp1d(self.ages, self.fluxGrid[:,i], kind='linear') 726 sedFlux=interpolator(ageGyr) 727 flux.append(sedFlux) 728 sed=SED(self.wavelengths, flux, z = z, normalise = normalise, label = label) 729 return sed
730
731 - def getColourEvolution(self, passband1, passband2, zFormation, zStepSize = 0.05, magType = "Vega"):
732 """Calculates the evolution of the colour observed through passband1 - passband2 for the 733 StellarPopulation with redshift, from z = 0 to z = zFormation. 734 735 @type passband1: L{Passband} object 736 @param passband1: filter passband through which to calculate the first magnitude 737 @type passband2: L{Passband} object 738 @param passband2: filter passband through which to calculate the second magnitude 739 @type zFormation: float 740 @param zFormation: formation redshift of the StellarPopulation 741 @type zStepSize: float 742 @param zStepSize: size of interval in z at which to calculate model colours 743 @type magType: string 744 @param magType: either "Vega" or "AB" 745 @rtype: dictionary 746 @return: dictionary of numpy.arrays in format {'z', 'colour'} 747 748 """ 749 750 zSteps=int(math.ceil(zFormation/zStepSize)) 751 zData=[] 752 colourData=[] 753 for i in range(1, zSteps): 754 zc=i*zStepSize 755 age=astCalc.tl(zFormation)-astCalc.tl(zc) 756 sed=self.getSED(age, z = zc) 757 colour=sed.calcColour(passband1, passband2, magType = magType) 758 zData.append(zc) 759 colourData.append(colour) 760 761 zData=numpy.array(zData) 762 colourData=numpy.array(colourData) 763 764 return {'z': zData, 'colour': colourData}
765
766 - def getMagEvolution(self, passband, magNormalisation, zNormalisation, zFormation, zStepSize = 0.05, 767 onePlusZSteps = False, magType = "Vega"):
768 """Calculates the evolution with redshift (from z = 0 to z = zFormation) of apparent magnitude 769 in the observed frame through the passband for the StellarPopulation, normalised to magNormalisation 770 (apparent) at z = zNormalisation. 771 772 @type passband: L{Passband} object 773 @param passband: filter passband through which to calculate the magnitude 774 @type magNormalisation: float 775 @param magNormalisation: sets the apparent magnitude of the SED at zNormalisation 776 @type zNormalisation: float 777 @param zNormalisation: the redshift at which the magnitude normalisation is carried out 778 @type zFormation: float 779 @param zFormation: formation redshift of the StellarPopulation 780 @type zStepSize: float 781 @param zStepSize: size of interval in z at which to calculate model magnitudes 782 @type onePlusZSteps: bool 783 @param onePlusZSteps: if True, zSteps are (1+z)*zStepSize, otherwise zSteps are linear 784 @type magType: string 785 @param magType: either "Vega" or "AB" 786 @rtype: dictionary 787 @return: dictionary of numpy.arrays in format {'z', 'mag'} 788 789 """ 790 791 # Count upwards in z steps as interpolation doesn't work if array ordered z decreasing 792 zSteps=int(math.ceil(zFormation/zStepSize)) 793 zData=[] 794 magData=[] 795 absMagData=[] 796 zc0=0.0 797 for i in range(1, zSteps): 798 if onePlusZSteps == False: 799 zc=i*zStepSize 800 else: 801 zc=zc0+(1+zc0)*zStepSize 802 zc0=zc 803 if zc >= zFormation: 804 break 805 age=astCalc.tl(zFormation)-astCalc.tl(zc) 806 sed=self.getSED(age, z = zc) 807 mag=sed.calcMag(passband, magType = magType, addDistanceModulus = True) 808 zData.append(zc) 809 magData.append(mag) 810 absMagData.append(sed.calcMag(passband, addDistanceModulus = False)) 811 812 zData=numpy.array(zData) 813 magData=numpy.array(magData) 814 815 # Do the normalisation 816 interpolator=interpolate.interp1d(zData, magData, kind='linear') 817 modelNormMag=interpolator(zNormalisation) 818 normConstant=magNormalisation-modelNormMag 819 magData=magData+normConstant 820 821 return {'z': zData, 'mag': magData}
822
823 - def calcEvolutionCorrection(self, zFrom, zTo, zFormation, passband, magType = "Vega"):
824 """Calculates the evolution correction in magnitudes in the rest frame through the passband 825 from redshift zFrom to redshift zTo, where the stellarPopulation is assumed to be formed 826 at redshift zFormation. 827 828 @type zFrom: float 829 @param zFormation: redshift to evolution correct from 830 @type zTo: float 831 @param zTo: redshift to evolution correct to 832 @type zFormation: float 833 @param zFormation: formation redshift of the StellarPopulation 834 @type passband: L{Passband} object 835 @param passband: filter passband through which to calculate magnitude 836 @type magType: string 837 @param magType: either "Vega" or "AB" 838 @rtype: float 839 @return: evolution correction in magnitudes in the rest frame 840 841 """ 842 843 ageFrom=astCalc.tl(zFormation)-astCalc.tl(zFrom) 844 ageTo=astCalc.tl(zFormation)-astCalc.tl(zTo) 845 846 fromSED=self.getSED(ageFrom) 847 toSED=self.getSED(ageTo) 848 849 fromMag=fromSED.calcMag(passband, magType = magType, addDistanceModulus = False) 850 toMag=toSED.calcMag(passband, magType = magType, addDistanceModulus = False) 851 852 return fromMag-toMag
853 854 #------------------------------------------------------------------------------------------------------------
855 -class M05Model(StellarPopulation):
856 """This class describes a Maraston 2005 stellar population model. To load a composite stellar population 857 model (CSP) for a tau = 0.1 Gyr burst of star formation, solar metallicity, Salpeter IMF: 858 859 m05csp = astSED.M05Model(M05_DIR+"/csp_e_0.10_z02_salp.sed_agb") 860 861 where M05_DIR is set to point to the directory where the Maraston 2005 models are stored on your system. 862 863 The file format of the Maraston 2005 simple stellar poulation (SSP) models is different to the file 864 format used for the CSPs, and this needs to be specified using the fileType parameter. To load a SSP with 865 solar metallicity, red horizontal branch morphology: 866 867 m05ssp = astSED.M05Model(M05_DIR+"/sed.ssz002.rhb", fileType = "ssp") 868 869 The wavelength units of SEDs from M05 models are Angstroms, with flux in units of erg/s/Angstrom. 870 871 """
872 - def __init__(self, fileName, fileType = "csp"):
873 874 self.modelFamily="M05" 875 876 inFile=file(fileName, "rb") 877 lines=inFile.readlines() 878 inFile.close() 879 880 self.fileName=fileName 881 882 if fileType == "csp": 883 ageColumn=0 884 wavelengthColumn=1 885 fluxColumn=2 886 elif fileType == "ssp": 887 ageColumn=0 888 wavelengthColumn=2 889 fluxColumn=3 890 else: 891 raise Exception, "fileType must be 'ssp' or 'csp'" 892 893 # Extract a list of model ages and valid wavelengths from the file 894 self.ages=[] 895 self.wavelengths=[] 896 for line in lines: 897 if line[0] !="#" and len(line) > 3: 898 bits=line.split() 899 age=float(bits[ageColumn]) 900 wavelength=float(bits[wavelengthColumn]) 901 if age not in self.ages: 902 self.ages.append(age) 903 if wavelength not in self.wavelengths: 904 self.wavelengths.append(wavelength) 905 906 # Construct a grid of flux - rows correspond to each wavelength, columns to age 907 self.fluxGrid=numpy.zeros([len(self.ages), len(self.wavelengths)]) 908 for line in lines: 909 if line[0] !="#" and len(line) > 3: 910 bits=line.split() 911 sedAge=float(bits[ageColumn]) 912 sedWavelength=float(bits[wavelengthColumn]) 913 sedFlux=float(bits[fluxColumn]) 914 915 row=self.ages.index(sedAge) 916 column=self.wavelengths.index(sedWavelength) 917 918 self.fluxGrid[row][column]=sedFlux
919 920 #------------------------------------------------------------------------------------------------------------
921 -class BC03Model(StellarPopulation):
922 """This class describes a Bruzual & Charlot 2003 stellar population model, extracted from a GALAXEV .ised 923 file using the galaxevpl program that is included in GALAXEV. The file format is white space delimited, 924 with wavelength in the first column. Subsequent columns contain the model fluxes for SEDs of different 925 ages, as specified when running galaxevpl. The age corresponding to each flux column is taken from the 926 comment line beginning "# Age (yr)", and is converted to Gyr. 927 928 For example, to load a tau = 0.1 Gyr burst of star formation, solar metallicity, Salpeter IMF model 929 stored in a file (created by galaxevpl) called "csp_lr_solar_0p1Gyr.136": 930 931 bc03model = BC03Model("csp_lr_solar_0p1Gyr.136") 932 933 The wavelength units of SEDs from BC03 models are Angstroms. Flux is converted into units of 934 erg/s/Angstrom (the units in the files output by galaxevpl are LSun/Angstrom). 935 936 """ 937
938 - def __init__(self, fileName):
939 940 self.modelFamily="BC03" 941 self.fileName=fileName 942 943 inFile=file(fileName, "rb") 944 lines=inFile.readlines() 945 inFile.close() 946 947 # Extract a list of model ages - BC03 ages are in years, so convert to Gyr 948 self.ages=[] 949 for line in lines: 950 if line.find("# Age (yr)") != -1: 951 rawAges=line[line.find("# Age (yr)")+10:].split() 952 for age in rawAges: 953 self.ages.append(float(age)/1e9) 954 955 # Extract a list of valid wavelengths from the file 956 # If we have many ages in the file, this is more complicated... 957 lambdaLinesCount=0 958 startFluxDataLine=None 959 for i in range(len(lines)): 960 line=lines[i] 961 if "# Lambda(A)" in line: 962 lambdaLinesCount=lambdaLinesCount+1 963 if line[0] != "#" and len(line) > 3 and startFluxDataLine == None: 964 startFluxDataLine=i 965 self.wavelengths=[] 966 for i in range(startFluxDataLine, len(lines), lambdaLinesCount): 967 line=lines[i] 968 bits=line.split() 969 self.wavelengths.append(float(bits[0])) 970 971 # Construct a grid of flux - rows correspond to each wavelength, columns to age 972 self.fluxGrid=numpy.zeros([len(self.ages), len(self.wavelengths)]) 973 for i in range(startFluxDataLine, len(lines), lambdaLinesCount): 974 line=lines[i] 975 bits=[] 976 for k in range(i, i+lambdaLinesCount): 977 bits=bits+lines[k].split() 978 ageFluxes=bits[1:] 979 sedWavelength=float(bits[0]) 980 column=self.wavelengths.index(sedWavelength) 981 for row in range(len(ageFluxes)): 982 sedFlux=float(ageFluxes[row]) 983 self.fluxGrid[row][column]=sedFlux 984 985 # Convert flux into erg/s/Angstrom - native units of galaxevpl files are LSun/Angstrom 986 self.fluxGrid=self.fluxGrid*3.826e33
987 988 #------------------------------------------------------------------------------------------------------------
989 -class P09Model(StellarPopulation):
990 """This class describes a Percival et al 2009 (BaSTI; http://albione.oa-teramo.inaf.it) stellar 991 population model. We assume that the synthetic spectra for each model are unpacked under the directory 992 pointed to by fileName. 993 994 The wavelength units of SEDs from P09 models are converted to Angstroms. Flux is converted into units of 995 erg/s/Angstrom (the units in the BaSTI low-res spectra are 4.3607e-33 erg/s/m). 996 997 """ 998
999 - def __init__(self, fileName):
1000 1001 self.modelFamily="P09" 1002 1003 files=glob.glob(fileName+os.path.sep+"*.t??????") 1004 1005 self.fileName=fileName 1006 1007 # Map end of filenames to ages in Gyr 1008 extensionAgeMap={} 1009 self.ages=[] 1010 for f in files: 1011 ext=f.split(".")[-1] 1012 ageGyr=float(f[-5:])/1e3 1013 self.ages.append(ageGyr) 1014 extensionAgeMap[ext]=ageGyr 1015 self.ages.sort() 1016 1017 # Construct a grid of flux - rows correspond to each wavelength, columns to age 1018 self.wavelengths=None 1019 self.fluxGrid=None 1020 for i in range(len(self.ages)): 1021 for e in extensionAgeMap.keys(): 1022 if extensionAgeMap[e] == self.ages[i]: 1023 inFileName=glob.glob(fileName+os.path.sep+"*."+e)[0] 1024 inFile=file(inFileName, "rb") 1025 lines=inFile.readlines() 1026 inFile.close() 1027 wavelength=[] 1028 flux=[] 1029 for line in lines: 1030 bits=line.split() 1031 wavelength.append(float(bits[0])*10.0) # units in file are nm, not angstroms 1032 flux.append(float(bits[1])) 1033 if self.wavelengths == None: 1034 self.wavelengths=wavelength 1035 if self.fluxGrid == None: 1036 self.fluxGrid=numpy.zeros([len(self.ages), len(self.wavelengths)]) 1037 self.fluxGrid[i]=flux 1038 1039 # Convert flux into erg/s/Angstrom - native units in BaSTI files are 4.3607e-33 erg/s/m 1040 self.fluxGrid=self.fluxGrid/4.3607e-33/1e10
1041 1042 #------------------------------------------------------------------------------------------------------------
1043 -def makeModelSEDDictList(modelList, z, passbandsList, labelsList = [], EBMinusVList = [0.0], forceYoungerThanUniverse = True):
1044 """This routine makes a list of SEDDict dictionaries (see L{mags2SEDDict}) for fitting using 1045 L{fitSEDDict}. This speeds up the fitting as this allows us to calculate model SED magnitudes only once, 1046 if all objects to be fitted are at the same redshift. We add some meta data to the modelSEDDicts (e.g. 1047 the model file names). 1048 1049 The effect of extinction by dust (assuming the Calzetti et al. 2000 law) can be included by giving a list 1050 of E(B-V) values. 1051 1052 If forceYoungerThanUniverse == True, ages which are older than the universe at the given z will not be 1053 included. 1054 1055 @type modelList: list of L{StellarPopulation} model objects 1056 @param modelList: list of StellarPopulation models to include 1057 @type z: float 1058 @param z: redshift to apply to all stellar population models in modelList 1059 @type EBMinusVList: list 1060 @param EBMinusVList: list of E(B-V) extinction values to apply to all models, in magnitudes 1061 @type labelsList: list 1062 @param labelsList: optional list used for labelling passbands in output SEDDicts 1063 @type forceYoungerThanUniverse: bool 1064 @param forceYoungerThanUniverse: if True, do not allow models that exceed the age of the universe at z 1065 @rtype: list 1066 @return: list of dictionaries containing model fluxes, to be used as input to L{fitSEDDict}. 1067 1068 """ 1069 1070 # Otherwise if this is the case we won't actually make any model SEDDicts ... 1071 if EBMinusVList == []: 1072 EBMinusVList=[0.0] 1073 1074 modelSEDDictList=[] 1075 for m in range(len(modelList)): 1076 testAges=numpy.array(modelList[m].ages) 1077 if forceYoungerThanUniverse == True: 1078 testAges=testAges[numpy.logical_and(numpy.less(testAges, astCalc.tz(z)), numpy.greater(testAges, 0))] 1079 for t in testAges: 1080 s=modelList[m].getSED(t, z = z, label=modelList[m].fileName+" - age="+str(t)+" Gyr") 1081 for EBMinusV in EBMinusVList: 1082 try: 1083 s.extinctionCalzetti(EBMinusV) 1084 except: 1085 raise Exception, "Model %s has a wavelength range that doesn't cover ~1200-22000 Angstroms" % (modelList[m].fileName) 1086 modelSEDDict=s.getSEDDict(passbandsList) 1087 modelSEDDict['labels']=labelsList 1088 modelSEDDict['E(B-V)']=EBMinusV 1089 modelSEDDict['ageGyr']=t 1090 modelSEDDict['z']=z 1091 modelSEDDict['fileName']=modelList[m].fileName 1092 modelSEDDict['modelListIndex']=m 1093 modelSEDDictList.append(modelSEDDict) 1094 1095 return modelSEDDictList
1096 1097 #------------------------------------------------------------------------------------------------------------
1098 -def fitSEDDict(SEDDict, modelSEDDictList):
1099 """Fits the given SED dictionary (made using L{mags2SEDDict}) with the given list of model SED 1100 dictionaries. The latter should be made using L{makeModelSEDDictList}, and entries for fluxes should 1101 correspond directly between the model and SEDDict. 1102 1103 Returns a dictionary with best fit values. 1104 1105 @type SEDDict: dictionary, in format of L{mags2SEDDict} 1106 @param SEDDict: dictionary of observed fluxes and uncertainties, in format of L{mags2SEDDict} 1107 @type modelSEDDictList: list of dictionaries, in format of L{makeModelSEDDictList} 1108 @param modelSEDDictList: list of dictionaries containing fluxes of models to be fitted to the observed 1109 fluxes listed in the SEDDict. This should be made using L{makeModelSEDDictList}. 1110 @rtype: dictionary 1111 @return: results of the fitting - keys: 1112 - 'minChiSq': minimum chi squared value of best fit 1113 - 'chiSqContrib': corresponding contribution at each passband to the minimum chi squared value 1114 - 'ageGyr': the age in Gyr of the best fitting model 1115 - 'modelFileName': the file name of the stellar population model corresponding to the best fit 1116 - 'modelListIndex': the index of the best fitting model in the input modelSEDDictList 1117 - 'norm': the normalisation that the best fit model should be multiplied by to match the SEDDict 1118 - 'z': the redshift of the best fit model 1119 - 'E(B-V)': the extinction, E(B-V), in magnitudes, of the best fit model 1120 1121 """ 1122 1123 modelFlux=[] 1124 for modelSEDDict in modelSEDDictList: 1125 modelFlux.append(modelSEDDict['flux']) 1126 modelFlux=numpy.array(modelFlux) 1127 sedFlux=numpy.array([SEDDict['flux']]*len(modelSEDDictList)) 1128 sedFluxErr=numpy.array([SEDDict['fluxErr']]*len(modelSEDDictList)) 1129 1130 # Analytic expression below is for normalisation at minimum chi squared (see note book) 1131 norm=numpy.sum((modelFlux*sedFlux)/(sedFluxErr**2), axis=1)/numpy.sum(modelFlux**2/sedFluxErr**2, axis=1) 1132 norms=numpy.array([norm]*modelFlux.shape[1]).transpose() 1133 chiSq=numpy.sum(((sedFlux-norms*modelFlux)**2)/sedFluxErr**2, axis=1) 1134 chiSq[numpy.isnan(chiSq)]=1e6 # throw these out, should check this out and handle more gracefully 1135 minChiSq=chiSq.min() 1136 bestMatchIndex=numpy.equal(chiSq, minChiSq).nonzero()[0][0] 1137 bestNorm=norm[bestMatchIndex] 1138 bestChiSq=minChiSq 1139 bestChiSqContrib=((sedFlux[bestMatchIndex]-norms[bestMatchIndex]*modelFlux[bestMatchIndex])**2)\ 1140 /sedFluxErr[bestMatchIndex]**2 1141 1142 resultsDict={'minChiSq': bestChiSq, 1143 'chiSqContrib': bestChiSqContrib, 1144 'allChiSqValues': chiSq, 1145 'ageGyr': modelSEDDictList[bestMatchIndex]['ageGyr'], 1146 'modelFileName': modelSEDDictList[bestMatchIndex]['fileName'], 1147 'modelListIndex': modelSEDDictList[bestMatchIndex]['modelListIndex'], 1148 'norm': bestNorm, 1149 'z': modelSEDDictList[bestMatchIndex]['z'], 1150 'E(B-V)': modelSEDDictList[bestMatchIndex]['E(B-V)']} 1151 1152 return resultsDict
1153 1154 #------------------------------------------------------------------------------------------------------------
1155 -def mags2SEDDict(ABMags, ABMagErrs, passbands):
1156 """Takes a set of corresponding AB magnitudes, uncertainties, and passbands, and 1157 returns a dictionary with keys 'flux', 'fluxErr', 'wavelength'. Fluxes are in units of 1158 erg/s/cm^2/Angstrom, wavelength in Angstroms. These dictionaries are the staple diet of the 1159 L{fitSEDDict} routine. 1160 1161 @type ABMags: list or numpy array 1162 @param ABMags: AB magnitudes, specified in corresponding order to passbands and ABMagErrs 1163 @type ABMagErrs: list or numpy array 1164 @param ABMagErrs: AB magnitude errors, specified in corresponding order to passbands and ABMags 1165 @type passbands: list of L{Passband} objects 1166 @param passbands: passband objects, specified in corresponding order to ABMags and ABMagErrs 1167 @rtype: dictionary 1168 @return: dictionary with keys {'flux', 'fluxErr', 'wavelength'}, suitable for input to L{fitSEDDict} 1169 1170 """ 1171 1172 flux=[] 1173 fluxErr=[] 1174 wavelength=[] 1175 for m, e, p in zip(ABMags, ABMagErrs, passbands): 1176 f, err=mag2Flux(m, e, p) 1177 flux.append(f) 1178 fluxErr.append(err) 1179 wavelength.append(p.effectiveWavelength()) 1180 1181 SEDDict={} 1182 SEDDict['flux']=numpy.array(flux) 1183 SEDDict['fluxErr']=numpy.array(fluxErr) 1184 SEDDict['wavelength']=numpy.array(wavelength) 1185 1186 return SEDDict
1187 1188 #------------------------------------------------------------------------------------------------------------
1189 -def mag2Flux(ABMag, ABMagErr, passband):
1190 """Converts given AB magnitude and uncertainty into flux, in erg/s/cm^2/Angstrom. 1191 1192 @type ABMag: float 1193 @param ABMag: magnitude on AB system in passband 1194 @type ABMagErr: float 1195 @param ABMagErr: uncertainty in AB magnitude in passband 1196 @type passband: L{Passband} object 1197 @param passband: L{Passband} object at which ABMag was measured 1198 @rtype: list 1199 @return: [flux, fluxError], in units of erg/s/cm^2/Angstrom 1200 1201 """ 1202 1203 fluxJy=(10**23.0)*10**(-(ABMag+48.6)/2.5) # AB mag 1204 aLambda=3e-13 # for conversion to erg s-1 cm-2 angstrom-1 with lambda in microns 1205 effLMicron=passband.effectiveWavelength()*(1e-10/1e-6) 1206 fluxWLUnits=aLambda*fluxJy/effLMicron**2 1207 1208 fluxJyErr=(10**23.0)*10**(-(ABMag-ABMagErr+48.6)/2.5) # AB mag 1209 fluxWLUnitsErr=aLambda*fluxJyErr/effLMicron**2 1210 fluxWLUnitsErr=fluxWLUnitsErr-fluxWLUnits 1211 1212 return [fluxWLUnits, fluxWLUnitsErr]
1213 1214 #------------------------------------------------------------------------------------------------------------
1215 -def flux2Mag(flux, fluxErr, passband):
1216 """Converts given flux and uncertainty in erg/s/cm^2/Angstrom into AB magnitudes. 1217 1218 @type flux: float 1219 @param flux: flux in erg/s/cm^2/Angstrom in passband 1220 @type fluxErr: float 1221 @param fluxErr: uncertainty in flux in passband, in erg/s/cm^2/Angstrom 1222 @type passband: L{Passband} object 1223 @param passband: L{Passband} object at which ABMag was measured 1224 @rtype: list 1225 @return: [ABMag, ABMagError], in AB magnitudes 1226 1227 """ 1228 1229 # aLambda = 3x10-5 for effective wavelength in angstroms 1230 aLambda=3e-13 # for conversion to erg s-1 cm-2 angstrom-1 with lambda in microns 1231 effLMicron=passband.effectiveWavelength()*(1e-10/1e-6) 1232 1233 fluxJy=(flux*effLMicron**2)/aLambda 1234 mag=-2.5*numpy.log10(fluxJy/10**23)-48.6 1235 1236 fluxErrJy=(fluxErr*effLMicron**2)/aLambda 1237 magErr=mag-(-2.5*numpy.log10((fluxJy+fluxErrJy)/10**23)-48.6) 1238 1239 return [mag, magErr]
1240 1241 #------------------------------------------------------------------------------------------------------------
1242 -def mag2Jy(ABMag):
1243 """Converts an AB magnitude into flux density in Jy 1244 1245 @type ABMag: float 1246 @param ABMag: AB magnitude 1247 @rtype: float 1248 @return: flux density in Jy 1249 1250 """ 1251 1252 fluxJy=((10**23)*10**(-(float(ABMag)+48.6)/2.5)) 1253 1254 return fluxJy
1255 1256 1257 #------------------------------------------------------------------------------------------------------------
1258 -def Jy2Mag(fluxJy):
1259 """Converts flux density in Jy into AB magnitude 1260 1261 @type fluxJy: float 1262 @param fluxJy: flux density in Jy 1263 @rtype: float 1264 @return: AB magnitude 1265 1266 """ 1267 1268 ABMag=-2.5*(numpy.log10(fluxJy)-23.0)-48.6 1269 1270 return ABMag
1271 1272 #------------------------------------------------------------------------------------------------------------ 1273 # Data 1274 VEGA=VegaSED() 1275 1276 # AB SED has constant flux density 3631 Jy 1277 AB=SED(wavelength = numpy.logspace(1, 8, 1e5), flux = numpy.ones(1e6)) 1278 AB.flux=(3e-5*3631)/(AB.wavelength**2) 1279 AB.z0flux=AB.flux[:] 1280 1281 # Solar SED from HST CALSPEC (http://www.stsci.edu/hst/observatory/cdbs/calspec.html) 1282 SOL=SED() 1283 SOL.loadFromFile(astLib.__path__[0]+os.path.sep+"data"+os.path.sep+"sun_reference_stis_001.ascii") 1284