Skip to content

Commit c8e3a4b

Browse files
committed
tests
1 parent 0cd7697 commit c8e3a4b

File tree

4 files changed

+56
-59
lines changed

4 files changed

+56
-59
lines changed

specutils/io/default_loaders/tabular_fits.py

Lines changed: 13 additions & 16 deletions
Original file line numberDiff line numberDiff line change
@@ -47,11 +47,11 @@ def tabular_fits_loader(file_obj, column_mapping=None, hdu=1, store_data_header=
4747
4848
Parameters
4949
----------
50-
file_obj : str, file-like, or HDUList
51-
FITS file name, object (provided from name by Astropy I/O Registry),
52-
or HDUList (as resulting from astropy.io.fits.open()).
50+
file_obj : str, file-like, or `~astropy.io.fits.HDUList`
51+
FITS file name, object (provided from name by Astropy I/O Registry), or
52+
`~astropy.io.fits.HDUList` (as resulting from `astropy.io.fits.open`).
5353
hdu : int
54-
The HDU of the fits file (default: 1st extension) to read from
54+
The HDU of the fits file (default: 1st extension) to read
5555
store_data_header : bool
5656
Defaults to ``False``, which stores the primary header in ``Spectrum.meta['header']``.
5757
Set to ``True`` to instead store the header from the specified data HDU.
@@ -86,11 +86,8 @@ def tabular_fits_loader(file_obj, column_mapping=None, hdu=1, store_data_header=
8686
else:
8787
tab.meta = hdulist[0].header
8888

89-
# Determine if there is a correlation matrix
90-
correl = None
91-
if 'CORREL' in [h.name for h in hdulist]:
92-
correl = Table.read(hdulist['CORREL'])
93-
correl.meta = hdulist['CORREL'].header
89+
# Determine if there is a covariance matrix
90+
covar = Table.read(hdulist['COVAR']) if 'COVAR' in [h.name for h in hdulist] else None
9491

9592
# Minimal checks for wcs consistency with table data -
9693
# assume 1D spectral axis (having shape (0, NAXIS1),
@@ -102,9 +99,9 @@ def tabular_fits_loader(file_obj, column_mapping=None, hdu=1, store_data_header=
10299
# If no column mapping is given, attempt to parse the file using
103100
# unit information
104101
if column_mapping is None:
105-
return generic_spectrum_from_table(tab, wcs=wcs, correl=correl, **kwargs)
102+
return generic_spectrum_from_table(tab, wcs=wcs, covar=covar, **kwargs)
106103

107-
return spectrum_from_column_mapping(tab, column_mapping, wcs=wcs, correl=correl)
104+
return spectrum_from_column_mapping(tab, column_mapping, wcs=wcs, covar=covar)
108105

109106

110107
@custom_writer("tabular-fits")
@@ -175,11 +172,11 @@ def tabular_fits_writer(spectrum, file_name, hdu=1, update_header=False, store_d
175172
colnames = [dispname, "flux"]
176173

177174
# Include uncertainty - units to be inferred from spectrum.flux
178-
correl = None
175+
covar = None
179176
if spectrum.uncertainty is not None:
180177
if isinstance(spectrum.uncertainty, Covariance):
181-
var, correl = spectrum.uncertainty.to_tables()
182-
columns.append(np.sqrt(var) * funit)
178+
covar = spectrum.uncertainty.to_table()
179+
columns.append(np.sqrt(spectrum.uncertainty.variance) * funit)
183180
colnames.append("uncertainty")
184181
else:
185182
try:
@@ -224,8 +221,8 @@ def tabular_fits_writer(spectrum, file_name, hdu=1, update_header=False, store_d
224221
hdu1 = fits.BinTableHDU(data=tab, name='DATA')
225222

226223
hdulist = fits.HDUList([hdu0, hdu1])
227-
if correl is not None:
228-
hdulist.append(fits.BinTableHDU(data=correl, name='CORREL'))
224+
if covar is not None:
225+
hdulist.append(fits.BinTableHDU(data=covar, name='COVAR'))
229226

230227
# TODO: Use output_verify options to check for valid FITS
231228
hdulist.writeto(file_name, **kwargs)

specutils/io/parsing_utils.py

Lines changed: 33 additions & 35 deletions
Original file line numberDiff line numberDiff line change
@@ -52,7 +52,7 @@ def read_fileobj_or_hdulist(*args, **kwargs):
5252
hdulist.close()
5353

5454

55-
def spectrum_from_column_mapping(table, column_mapping, wcs=None, correl=None, verbose=False):
55+
def spectrum_from_column_mapping(table, column_mapping, wcs=None, covar=None, verbose=False):
5656
"""
5757
Given a table and a mapping of the table column names to attributes
5858
on the Spectrum object, parse the information into a Spectrum.
@@ -76,10 +76,9 @@ def spectrum_from_column_mapping(table, column_mapping, wcs=None, correl=None, v
7676
wcs : :class:`~astropy.wcs.WCS` or :class:`gwcs.WCS`
7777
WCS object passed to the Spectrum initializer.
7878
79-
correl : `astropy.table.Table`, optional
80-
Table with correlation matrix for the uncertainties in coordinate
81-
format; see `~astropy.nddata.Covariance`. If None, uncertainties are
82-
assumed to be uncorrelated.
79+
covar : `astropy.table.Table`, optional
80+
Table providing a covariance matrix for the uncertainties in coordinate
81+
format; see `~astropy.nddata.Covariance`.
8382
8483
verbose : bool
8584
Print extra info.
@@ -136,28 +135,27 @@ def spectrum_from_column_mapping(table, column_mapping, wcs=None, correl=None, v
136135
spec_kwargs.setdefault(kwarg_name, kwarg_val)
137136

138137
# Ensure that the uncertainties are a subclass of NDUncertainty
139-
if spec_kwargs.get('uncertainty') is None and correl is not None:
138+
if spec_kwargs.get('uncertainty') is None and covar is not None:
140139
warnings.warn('Unable to parse uncertainty from provided table. Ignoring provided '
141-
'correlation matrix data.')
142-
correl = None
140+
'covariance matrix data.')
141+
covar = None
143142
if spec_kwargs.get('uncertainty') is not None:
144-
if correl is not None:
145-
err = spec_kwargs.get('uncertainty')
143+
if covar is not None:
146144
try:
147-
spec_kwargs['uncertainty'] = Covariance.from_tables(err**2, correl)
145+
spec_kwargs['uncertainty'] = Covariance.from_table(covar)
148146
except (ValueError, TypeError):
149-
warnings.warn('Unable to parse correlation table into a Covariance object. '
150-
'Ignoring correlation matrix data.')
151-
correl = None
147+
warnings.warn('Unable to parse table with covariance data into into a Covariance '
148+
'object. Ignoring covariance matrix data.')
149+
covar = None
152150
# NOTE: This is not an `else` or `elif` block in order to catch the
153-
# change to correl=None when handling the exception above.
154-
if correl is None:
151+
# change to covar=None when handling the exception above.
152+
if covar is None:
155153
spec_kwargs['uncertainty'] = StdDevUncertainty(spec_kwargs.get('uncertainty'))
156154

157155
return Spectrum(**spec_kwargs, wcs=wcs, meta={'header': table.meta})
158156

159157

160-
def generic_spectrum_from_table(table, wcs=None, correl=None, **kwargs):
158+
def generic_spectrum_from_table(table, wcs=None, covar=None, **kwargs):
161159
"""
162160
Load spectrum from an Astropy table into a Spectrum object.
163161
Uses the following logic to figure out which column is which:
@@ -174,15 +172,15 @@ def generic_spectrum_from_table(table, wcs=None, correl=None, **kwargs):
174172
175173
Parameters
176174
----------
177-
table : `astropy.table.QTable`
178-
Tabulated data with units
175+
table : :class:`~astropy.table.Table`
176+
Table containing a column of ``flux``, and optionally ``spectral_axis``
177+
and ``uncertainty`` as defined above.
179178
wcs : :class:`~astropy.wcs.WCS`, optional
180179
A FITS WCS object. If this is present, the machinery will fall back
181-
to using the wcs to find the dispersion information.
182-
correl : `astropy.table.Table`, optional
183-
Table with correlation matrix for the uncertainties in coordinate
184-
format; see `~astropy.nddata.Covariance`. If None, uncertainties are
185-
assumed to be uncorrelated.
180+
to using the ``wcs`` to find the dispersion information.
181+
covar : `astropy.table.Table`, optional
182+
Table providing a covariance matrix for the uncertainties in coordinate
183+
format; see `~astropy.nddata.Covariance`.
186184
187185
Returns
188186
-------
@@ -297,24 +295,24 @@ def _find_spectral_column(table, columns_to_search, spectral_axis):
297295
if table[err_column].ndim > 1:
298296
err = table[err_column].T
299297
elif flux.ndim > 1: # Repeat uncertainties over all flux columns
300-
if correl is not None:
301-
warnings.warn('When applying correlated errors, the dimensionality of the error '
298+
if covar is not None:
299+
warnings.warn('When applying covariance, the dimensionality of the error '
302300
'array must match the dimensionality of the flux array. Ignoring '
303-
'correlated errors.')
304-
correl = None
301+
'covariance.')
302+
covar = None
305303
err = np.tile(table[err_column], flux.shape[0], 1)
306304
else:
307305
err = table[err_column]
308-
if correl is not None:
306+
if covar is not None:
309307
try:
310-
err = Covariance.from_tables(err**2, correl)
308+
err = Covariance.from_table(covar)
311309
except (ValueError, TypeError):
312-
warnings.warn('Unable to parse correlation table into a Covariance object. '
313-
'Ignoring correlation matrix data.')
314-
correl = None
310+
warnings.warn('Unable to parse covariance table into a Covariance object. '
311+
'Ignoring covariance matrix data.')
312+
covar = None
315313
# NOTE: This is not an `else` or `elif` block in order to catch the
316-
# change to correl=None when handling the exception above.
317-
if correl is None:
314+
# change to covar=None when handling the exception above.
315+
if covar is None:
318316
err = StdDevUncertainty(err.to(err.unit))
319317
if np.min(table[err_column]) <= 0.:
320318
warnings.warn("Standard Deviation has values of 0 or less", AstropyUserWarning)

specutils/spectra/spectrum.py

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -506,7 +506,7 @@ def __getitem__(self, item):
506506
new_meta = deepcopy(self.meta)
507507

508508
if isinstance(self.uncertainty, Covariance):
509-
new_unc = self.uncertainty.sub_matrix(item)
509+
new_unc = self.uncertainty.match_to_data_slice(item)
510510
elif self.uncertainty is not None:
511511
new_unc = self.uncertainty[item]
512512
else:
@@ -933,7 +933,7 @@ def __str__(self):
933933
# Add information about uncertainties if available
934934
if self.uncertainty:
935935
_arr = (
936-
self.uncertainty.toarray()
936+
self.uncertainty.to_dense()
937937
if isinstance(self.uncertainty, Covariance)
938938
else self.uncertainty.array
939939
)

specutils/tests/test_loaders.py

Lines changed: 8 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -610,25 +610,25 @@ def test_tabular_fits_cov_io(tmp_path):
610610
np.full(1000-1, 0.5, dtype=float),
611611
np.full(1000-2, 0.2, dtype=float),
612612
]
613-
cov = Covariance(array=sparse.diags(cov_diags, [0, 1, 2]), unit=u.Jy**2)
613+
cov = Covariance(array=sparse.diags(cov_diags, [0, 1, 2]), unit=u.Jy**2, assume_symmetric=True)
614614

615615
# Create the Spectrum1D
616616
spectrum = Spectrum1D(flux=flux, spectral_axis=wave, uncertainty=cov)
617617
# Simple test of ingestion
618-
assert np.array_equal(spectrum.uncertainty.toarray(), cov.toarray()), \
618+
assert np.array_equal(spectrum.uncertainty.to_dense(), cov.to_dense()), \
619619
'Covariance not included correctly'
620620
# Write it
621621
tmpfile = str(tmp_path / '_tst.fits')
622622
spectrum.write(tmpfile, format='tabular-fits')
623623

624624
# Check the output file
625625
with fits.open(tmpfile) as hdu:
626-
assert len(hdu) == 3, 'Should contain 3 extensions, primary, DATA, CORREL'
627-
assert np.array_equal([h.name for h in hdu], ['PRIMARY', 'DATA', 'CORREL']), \
626+
assert len(hdu) == 3, 'Should contain 3 extensions, primary, DATA, COVAR'
627+
assert np.array_equal([h.name for h in hdu], ['PRIMARY', 'DATA', 'COVAR']), \
628628
'Extension names are wrong'
629629
assert all([h.__class__.__name__ == 'BinTableHDU' for h in hdu[1:]]), \
630630
'Data extensions should both be BinTableHDU'
631-
assert len(hdu['CORREL'].data) == cov.stored_nnz, \
631+
assert len(hdu['COVAR'].data) == cov.stored_nnz, \
632632
'Number of non-zero cov elements mismatch'
633633

634634
# Read it
@@ -637,7 +637,9 @@ def test_tabular_fits_cov_io(tmp_path):
637637
assert spectrum.spectral_axis.unit == _spectrum.spectral_axis.unit
638638
assert quantity_allclose(spectrum.spectral_axis, _spectrum.spectral_axis)
639639
assert quantity_allclose(spectrum.flux, _spectrum.flux)
640-
assert np.array_equal(spectrum.uncertainty.toarray(), _spectrum.uncertainty.toarray())
640+
assert np.array_equal(spectrum.uncertainty.to_dense(), _spectrum.uncertainty.to_dense())
641+
642+
# TODO: Add test for multidim
641643

642644

643645
@pytest.mark.parametrize("ndim", range(1, 4))

0 commit comments

Comments
 (0)