Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add rebinning to a constant resolution for spectra and extinction data #100

Merged
merged 18 commits into from
Mar 2, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
22 changes: 10 additions & 12 deletions docs/calc_extinction.rst
Original file line number Diff line number Diff line change
Expand Up @@ -244,19 +244,17 @@ exists.
.. plot::

import pkg_resources

import numpy as np
import matplotlib.pyplot as plt

from measure_extinction.stardata import StarData
from measure_extinction.extdata import ExtData

# get the location of the data files
data_path = pkg_resources.resource_filename('measure_extinction',
'data/')
data_path = pkg_resources.resource_filename("measure_extinction", "data/")

# read in the observed data on the star
redstar = StarData('hd229238.dat', path=data_path)
compstar = StarData('hd204172.dat', path=data_path)
redstar = StarData("hd229238.dat", path=data_path)
compstar = StarData("hd204172.dat", path=data_path)

# calculate the extinction curve
extdata = ExtData()
Expand All @@ -272,12 +270,12 @@ exists.
extdata.plot(ax)

# finish configuring the plot
ax.set_title('HD229238/HD204172 extinction')
ax.set_xscale('log')
ax.set_xlabel('$\lambda$ [$\mu m$]')
ax.set_ylabel('$E(\lambda - V)/E(B-V)$')
ax.tick_params('both', length=10, width=2, which='major')
ax.tick_params('both', length=5, width=1, which='minor')
ax.set_title("HD229238/HD204172 extinction")
ax.set_xscale("log")
ax.set_xlabel(r"$\lambda$ [$\mu m$]")
ax.set_ylabel(r"$E(\lambda - V)/E(B-V)$")
ax.tick_params("both", length=10, width=2, which="major")
ax.tick_params("both", length=5, width=1, which="minor")

# use the whitespace better
fig.tight_layout()
Expand Down
131 changes: 131 additions & 0 deletions docs/code_capabilities.rst
Original file line number Diff line number Diff line change
@@ -0,0 +1,131 @@

=================
Code Capabilities
=================

The code is built around two classes.
StarData is for photometry and spectroscopy of a single star.
ExtData is for calculating and storing the extinction curve for a single
sightline.

For both classes, the data is stored by the source of origin of the data. Possible
origins are band photometry ("BAND") or spectroscopy from the
International Ultraviolet Explorer ("IUE"),
Far-Ultraviolet Spectroscopic Explorer ("FUSE"),
Hubble Space Telescope Imaging Spectrograph ("STIS"),
NASA's Infrared Telescope SpeX Instrument ("SpeX_SXD", "SpeX_LXD"),
and Spitzer Infrared Spectrograph ("IRS")
The strings in "" give the dictionary key.

This package assumes the spectral data from a specific source is on the
same wavelength grid for all stars/sightlines.
This simplifies the extinction calculations.
Code to do this for some data sources can be found in
`measure_extinction.merge_obsspec`.

StarData
========

The data stored as part of `measure_extinction.stardata.StarData`
is the photometry and spectroscopy for a
single star, each dataset stored in a class itself.
The photometry is stored in the dedicated BandData class to
provide specific capabilities needed for handling photometric data.
The spectroscopy is stored in SpecData class.

The details of the data stored include:

* file: name of DAT file
* path: path of DAT file
* data: dictionary of BandData and SpecData objects containing the stellar data
* photonly: only photometry used (all spectroscopic data ignored)
* sptype: spectral type of star from DAT file
* use_corfac: boolean for determining if corfacs should be used
* corfac: dictionary of correction factors for spectral data
* dereddened: boolean if the data has been dereddened
* dereddenParams: dictionary of FM90 and CCM89 dereddening parameters
* model_parameters: stellar atmosphere model parameters (for model based data)

karllark marked this conversation as resolved.
Show resolved Hide resolved
Member Functions
----------------

The member functions of StarData include:

* read: read in the data on a star based on a DAT formatted file
* deredden: deredden the data based on FM90 (UV) and CCM89 (optical/NIR)
* get_flat_data_arrays: get waves, fluxes, uncs as vectors combining all the data
* plot: plot the stellar data given a matplotlib plot object

BandData
========

Member Functions
----------------

SpecData
========

A spectrum is stored as part of `measure_extinction.stardata.SpecData`.
The details of the data stored are:

* waves: wavelengths with units
* n_waves: number of wavelengths
* wave_range: min/max of waves with units
* fluxes: fluxes versus wavelengths with units
* uncs: flux uncertainties versus wavelength with units
* npts: number of measurements versus wavelength

Member Functions
----------------

The member functions of SpecData include:

* read_spectra: generic read for a spectrum in a specifically formatted FITS file
* read_fuse: read a FUSE spectrum
* read_iue: read an IUE spectrum (includes cutting data > 3200 A)
* read_stis: read a Hubble/STIS spectrum
* read_spex: read a IRTF/SpeX spectrum (includes scaling by corfacs)
* read_irs: read a Spitzer/IRS spectrum (includes scaling by corfacs, cutting above some wavelength)
* rebin_constres: rebin spectrum to a constant input resolution

ExtData
=======

The data that is stored as part of `measure_extinction.extdata.ExtData`
is the extinction curve from each source and associated details.
Specific details stored include (where src = data source),

* type: type of extinction measurement (e.g., elx, elxebv, alax)
* type_rel_band: photometric band that is used for the relative extinction measurement (usually V band)
* waves[src]: wavelengths
* exts[src]: extinction versus wavelength
* uncs[src]: uncertainty in the extinction versus wavelength
* npts[src]: number of measurements at each wavelength
* names[src]: names of the photometric bands (only present for names["BAND"])
* columns: dictionary of gas or dust column measurements (e.g., A(V), R(V), N(HI))
* red_file: filename of the reddened star
* comp_file: filename of the comparison star

Member Functions
----------------

The member functions of ExtData include:

* calc_elx: calculate E(lambda-X) for all data sources present in both reddened and comparison StarData objects
* calc_elx_bands: calculate E(lambda-X) for the "BAND" data
* calc_elx_spectra: calculate E(lambda-X) for a single spectral data source
* calc_EBV: determine E(B-V) based on E(lambda-V) "BAND" data
* calc_AV: determine A(V) based on E(lambda-V) "BAND" data
* calc_RV: determine R(V) using calc_EBV and calc_AV as needed
* trans_elv_elvebv: calculate E(lambda-V)/E(B-V) based on E(lambda-V) and E(B-V)
* trans_elv_alav: calculate A(lambda)/A(V) based on E(lambda-V) and A(V)
* rebin_constres: rebin one data source extinction data to a constant input resolution
* get_fitdata: gets a tuple of vectors with (wavelengths, exts, uncs) that is useful for fitting
* save: save the extinction curve to a FITS file
* read: read an extinction curve from a FITS file
* plot: plot the extinction curve given a matplotlib plot object

.. note::
ExtData partially supports extinction curves relative to an arbitrary band.
Some member functions only support extinction curve relative to V band.
Work continues to update the code to allow arbitrary bands for all functions.
3 changes: 2 additions & 1 deletion docs/index.rst
Original file line number Diff line number Diff line change
Expand Up @@ -29,8 +29,9 @@ User Documentation
:maxdepth: 2

Extinction Calculation <calc_extinction.rst>
Observed Data Formats <data_formats.rst>
Capabilities <code_capabilities.rst>
Plotting <plotting.rst>
Observed Data Formats <data_formats.rst>
Stellar Models as Standards <model_standards.rst>

Reporting Issues
Expand Down
88 changes: 76 additions & 12 deletions measure_extinction/extdata.py
Original file line number Diff line number Diff line change
Expand Up @@ -13,6 +13,8 @@

from dust_extinction.conversions import AxAvToExv

from measure_extinction.merge_obsspec import _wavegrid

__all__ = ["ExtData", "AverageExtData"]


Expand Down Expand Up @@ -266,7 +268,7 @@ class ExtData:
comp_file : string
comparison star filename

columns : list of tuples of column measurements
columns : dict of tuples of column measurements
measurements are A(V), R(V), N(HI), etc.
tuples are measurement, uncertainty

Expand Down Expand Up @@ -534,9 +536,6 @@ def calc_AV_JHK(self):
- extrapolate from J, H, & K photometry
- assumes functional form from Rieke, Rieke, & Paul (1989)

Parameters
----------

Returns
-------
Updates self.columns["AV"]
Expand Down Expand Up @@ -573,9 +572,6 @@ def calc_RV(self):
"""
Calculate R(V) from the observed extinction curve

Parameters
----------

Returns
-------
Updates self.columns["RV"]
Expand Down Expand Up @@ -626,8 +622,10 @@ def trans_elv_elvebv(self, ebv=None):
fullebv = _get_column_plus_unc(ebv)

for curname in self.exts.keys():
self.uncs[curname] = (self.exts[curname] / fullebv[0]) * np.sqrt(
np.square(self.uncs[curname] / self.exts[curname])
# only compute where there is data and exts is not zero
gvals = (self.exts[curname] != 0) & (self.npts[curname] > 0)
self.uncs[curname][gvals] = (self.exts[curname][gvals] / fullebv[0]) * np.sqrt(
np.square(self.uncs[curname][gvals] / self.exts[curname][gvals])
karllark marked this conversation as resolved.
Show resolved Hide resolved
+ np.square(fullebv[1] / fullebv[0])
)
self.exts[curname] /= fullebv[0]
Expand Down Expand Up @@ -698,6 +696,75 @@ def trans_elv_alav(self, av=None, akav=0.112):

self.type = "alax"

def rebin_constres(self, source, waverange, resolution):
"""
Rebin the source extinction curve to a fixed spectral resolution
and min/max wavelength range.

Parameters
----------
source : str
source of extinction (i.e. "IUE", "IRS")
waverange : 2 element array of astropy Quantities
Min/max of wavelength range with units
resolution : float
Spectral resolution of rebinned extinction curve

Returns
-------
measure_extinction ExtData
Object with source extinction curve rebinned

"""
if source == "BAND":
raise ValueError("BAND extinction cannot be rebinned")

if source not in self.exts.keys():
warnings.warn(f"{source} extinction not present")
else:
# setup new wavelength grid
full_wave, full_wave_min, full_wave_max = _wavegrid(
resolution, waverange.to(u.micron).value
karllark marked this conversation as resolved.
Show resolved Hide resolved
)
n_waves = len(full_wave)

# setup the new rebinned vectors
new_waves = full_wave * u.micron
new_exts = np.zeros((n_waves), dtype=float)
new_uncs = np.zeros((n_waves), dtype=float)
new_npts = np.zeros((n_waves), dtype=int)

# check if uncertainties defined and set temporarily to 1 if not
# needed to avoid infinite weights
nouncs = False
if np.sum(self.uncs[source] > 0.0) == 0:
nouncs = True
self.uncs[source] = np.full((len(self.waves[source])), 1.0)

# rebin using a weighted average
owaves = self.waves[source].to(u.micron).value
for k in range(n_waves):
(indxs,) = np.where(
(owaves >= full_wave_min[k])
& (owaves < full_wave_max[k])
& (self.npts[source] > 0.0)
)
if len(indxs) > 0:
weights = 1.0 / np.square(self.uncs[source][indxs])
sweights = np.sum(weights)
new_exts[k] = np.sum(weights * self.exts[source][indxs]) / sweights
new_uncs[k] = 1.0 / np.sqrt(sweights)
new_npts[k] = np.sum(self.npts[source][indxs])

if nouncs:
new_uncs = np.full((n_waves), 0.0)

# update source values
self.waves[source] = new_waves
self.exts[source] = new_exts
self.uncs[source] = new_uncs
self.npts[source] = new_npts

def get_fitdata(
self,
req_datasources,
Expand Down Expand Up @@ -1374,9 +1441,6 @@ def fit_band_ext(self):
"""
Fit the observed NIR extinction curve with a powerlaw model, based on the band data between 1 and 40 micron

Parameters
----------

Returns
-------
Updates self.model["waves", "exts", "residuals", "params"] and self.columns["AV"] with the fitting results:
Expand Down
Loading