diff --git a/datamodel.txt b/datamodel.txt index f61271a..5284d2b 100644 --- a/datamodel.txt +++ b/datamodel.txt @@ -673,6 +673,64 @@ The format will be identical to that of the pfsObject files. -=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=- +Physical parameter measurements for Galactic Archeology targets including flux-calibrated +combined spectra and synthetic stellar template fits. + + "pfsGAObject-%05d-%05d-%s-%016x-%03d-0x%016x.fits" + % (catId, tract, patch, objId, nVisit % 1000, pfsVisitHash) + +The format is similar to pfsObject files, but with additional HDUs and a different fluxTable format. + +HDU #j VELCORR Velocity corrections used at each visit [BINARY FITS TABLE] +HDU #i STELLARPARAMS Fundamental stellar parameter measurements [BINARY FITS TABLE] +HDU #l STELLARCOVAR Covariance matrix of stellar parameters [BINARY FITS TABLE] n*n +HDU #k ABUND Single element abundances [BINARY FITS TABLE] +HDU #m ABUNDCOVAR Covarinace matrix of single element abundances [32-bit FLOAT] n*n + +In the data tables outlined below, we allow for the possibility of multiple measurements of the same +physical parameters indicated by the `method` field but we provide the full covariance matrix for +the primary method only. The column `covarId` indicates the position of the parameter within the +corresponding covariance matrix. + +The VELCORR table lists the velocity corrections used at each visit: + + visit visit identifier 32-bit int + JD Julian Date of the visit 32-bit float + helio Heliocentric correction at the visit [km s-1] 32-bit float + bary Barycentric correction at the visit [km s-2] 32-bit float + +The STELLARPARAMS table lists the measured fundamental parameters + + method Method of measuring the parameters string + frame Reference frame for measuring the velocity string + param Stellar parameter string + covarId Index of parameter within the covariance matrix 8-bit uint + unit Physical unit of the parameter string + value Measured value of the parameter 32-bit float + valueErr Uncertainty of the measured value 32-bit float + flag Measurement flag (true means bad) bool + status Flags describing the quality of the measurement string + +The ABUND table lists the measured single element abundances + + method Method of measuring the parameters string + element Element name string + covarId Index of parameter within the covariance matrix 8-bit uint + value Measured abundance 32-bit float + valueErr Uncertainty of the measured abundance 32-bit float + flag Measurement flag (true means bad) bool + status Flags describing the quality of the measurement string + +The FLUXTABLE of pfsGAObject files is extended and includes the following additional columns: + + model Best fit fluxed template model 32-bit float + cont Continuum model 32-bit float + norm_flux Continuum-normalized flux 32-bit float + norm_error Error of continuum-normalized flux 32-bit float + norm_model Best fit continuum-normalized model 32-bit float + +-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=- + 2D spectral line measurements The 2D pipeline measures the position and flux of sky and arc lines as part of the reduction process. diff --git a/python/pfs/datamodel/__init__.py b/python/pfs/datamodel/__init__.py index 1f2b6c1..8d7f19f 100644 --- a/python/pfs/datamodel/__init__.py +++ b/python/pfs/datamodel/__init__.py @@ -15,3 +15,4 @@ from .pfsTable import * from .pfsFocalPlaneFunction import * from .pfsFiberNorms import * +from .ga import * \ No newline at end of file diff --git a/python/pfs/datamodel/fluxTable.py b/python/pfs/datamodel/fluxTable.py index c60c1cd..192827f 100644 --- a/python/pfs/datamodel/fluxTable.py +++ b/python/pfs/datamodel/fluxTable.py @@ -38,11 +38,11 @@ class FluxTable: _hduName = "FLUX_TABLE" # HDU name to use def __init__(self, wavelength, flux, error, mask, flags): - dims = np.array([len(wavelength.shape), len(flux.shape), len(error.shape), len(mask.shape)]) - lengths = set([wavelength.shape, flux.shape, error.shape, mask.shape]) - if np.any(dims != 1) or len(lengths) > 1: - raise RuntimeError("Bad shapes for wavelength,flux,error,mask: %s,%s,%s,%s" % - (wavelength.shape, flux.shape, error.shape, mask.shape)) + self.checkShapes(wavelength=wavelength, + flux=flux, + error=error, + mask=mask) + self.wavelength = wavelength self.flux = flux self.error = error @@ -53,6 +53,17 @@ def __len__(self): """Return number of elements""" return len(self.wavelength) + def checkShapes(self, **kwargs): + keys = list(sorted(kwargs.keys())) + dims = np.array([len(kwargs[k].shape) for k in keys]) + lengths = set([kwargs[k].shape for k in keys]) + + if np.any(dims != 1) or len(lengths) > 1: + names = ','.join(keys) + shapes = ','.join([str(kwargs[k].shape) for k in keys]) + raise RuntimeError("Bad shapes for %s: %s" % + (names, shapes)) + def toFits(self, fits): """Write to a FITS file diff --git a/python/pfs/datamodel/ga.py b/python/pfs/datamodel/ga.py new file mode 100644 index 0000000..69ad3a7 --- /dev/null +++ b/python/pfs/datamodel/ga.py @@ -0,0 +1,321 @@ +import numpy as np + +from .notes import makeNotesClass, Notes +from .pfsFiberArray import PfsFiberArray +from .fluxTable import FluxTable +from .pfsTable import PfsTable, Column +from .utils import inheritDocstrings +from .utils import astropyHeaderToDict, astropyHeaderFromDict +from .masks import MaskHelper + +__all__ = [ + "VelocityCorrections", + "StellarParams", + "Abundances", + "GAFluxTable", + "PfsGAObjectNotes", + "PfsGAObject", +] + + +class VelocityCorrections(PfsTable): + """A table of velocity corrections applied to the individual visits.""" + + damdVer = 2 + schema = [ + Column("visit", np.int32, "ID of the visit these corrections apply for", -1), + Column("JD", np.float32, "Julian date of the visit", -1), + Column("helio", np.float32, "Heliocentric correction", np.nan), + Column("bary", np.float32, "Barycentric correction", np.nan), + ] + fitsExtName = 'VELCORR' + + +class StellarParams(PfsTable): + """List of measured stellar parameters for a target.""" + + damdVer = 2 + schema = [ + Column("method", str, "Line-of-sight velocity measurement method", ""), + Column("frame", str, "Reference frame of velocity: helio, bary", ""), + Column("param", str, "Stellar parameter: v_los, M_H, T_eff, log_g, a_M", ""), + Column("covarId", np.uint8, "Param position within covariance matrix", -1), + Column("unit", str, "Physical unit of parameter", ""), + Column("value", np.float32, "Stellar parameter value", np.nan), + Column("valueErr", np.float32, "Stellar parameter error", np.nan), + # TODO: add quantiles or similar for MCMC results + Column("flag", bool, "Measurement flag (true means bad)", False), + Column("status", str, "Measurement flags", ""), + ] + fitsExtName = 'STELLARPARAM' + + +class Abundances(PfsTable): + """List of measured abundance parameters for stellar targets.""" + + damdVer = 2 + schema = [ + Column("method", str, "Abundance measurement method", ""), + Column("element", str, "Chemical element the abundance is measured for", ""), + Column("covarId", np.uint8, "Param position within covariance matrix", -1), + Column("value", np.float32, "Abundance value", np.nan), + Column("valueErr", np.float32, "Abundance error", np.nan), + # TODO: will we have systematic errors? + Column("flag", bool, "Measurement flag (true means bad)", False), + Column("status", str, "Measurement flags", ""), + ] + fitsExtName = 'ABUND' + + +class GAFluxTable(FluxTable): + """Table of coadded fluxes at near-original sampling and model fits + + Merged and coadded spectra have been resampled to a standard wavelength + sampling. This representation provides coadded fluxes at approximately the + native wavelength sampling, for those that want the data with a minimum of + resampling. This is mostly of use for single exposures and coadds made from + back-to-back exposures with the same top-end configuration. For coadds made + from exposures with different top-end configurations, the different + wavelength samplings obtained from the different fibers means there's no + single native wavelength sampling, and so this is less useful. + + This is like a `pfs.datamodel.PfsSimpleSpectrum`, except that it includes a + variance array, and is written to a FITS HDU rather than a file (so it can + be incorporated within a `pfs.datamodel.PfsSpectrum`). + + Parameters + ---------- + wavelength : `numpy.ndarray` of `float` + Array of wavelengths. + flux : `numpy.ndarray` of `float` + Array of fluxes. + error : `numpy.ndarray` of `float` + Array of flux errors. + model : `numpy.ndarray` of `float` + Array of best-fit model flux. + cont : `numpy.ndarray` of `float` + Array of continuum model. + norm_flux : `numpy.ndarray` of `float` + Array of continuum-normalized flux. + norm_error : `numpy.ndarray` of `float` + Array of continuum-normalized flux error. + norm_model : `numpy.ndarray` of `float` + Array of continuum-normalized model. + mask : `numpy.ndarray` of `int` + Array of mask pixels. + flags : `pfs.datamodel.MaskHelper` + Helper for dealing with symbolic names for mask values. + """ + _hduName = "FLUX_TABLE" # HDU name to use + + def __init__(self, wavelength, flux, error, model, cont, norm_flux, norm_error, norm_model, mask, flags): + self.checkShapes(wavelength=wavelength, + flux=flux, + error=error, + model=model, + cont=cont, + norm_flux=norm_flux, + norm_error=norm_error, + norm_model=norm_model, + mask=mask) + + self.wavelength = wavelength + self.flux = flux + self.error = error + self.model = model + self.cont = cont + self.norm_flux = norm_flux + self.norm_error = norm_error + self.norm_model = norm_model + self.mask = mask + self.flags = flags + + def toFits(self, fits): + """Write to a FITS file + + Parameters + ---------- + fits : `astropy.io.fits.HDUList` + Opened FITS file. + """ + # NOTE: When making any changes to this method that modify the output + # format, increment the DAMD_VER header value and record the change in + # the versions.txt file. + from astropy.io.fits import BinTableHDU, Column + header = astropyHeaderFromDict(self.flags.toFitsHeader()) + header['DAMD_VER'] = (1, "GAFluxTable datamodel version") + hdu = BinTableHDU.from_columns([ + Column("wavelength", "D", array=self.wavelength), + Column("flux", "E", array=self.flux), + Column("error", "E", array=self.error), + Column("model", "E", array=self.model), + Column("cont", "E", array=self.cont), + Column("norm_flux", "E", array=self.norm_flux), + Column("norm_error", "E", array=self.norm_error), + Column("norm_model", "E", array=self.norm_model), + Column("mask", "K", array=self.mask), + ], header=header, name=self._hduName) + fits.append(hdu) + + @classmethod + def fromFits(cls, fits): + """Construct from a FITS file + + Parameters + ---------- + fits : `astropy.io.fits.HDUList` + Opened FITS file. + + Returns + ------- + self : `FluxTable` + Constructed `FluxTable`. + """ + hdu = fits[cls._hduName] + header = astropyHeaderToDict(hdu.header) + flags = MaskHelper.fromFitsHeader(header) + return cls(hdu.data["wavelength"].astype(float), + hdu.data["flux"].astype(np.float32), + hdu.data["error"].astype(np.float32), + hdu.data["model"].astype(np.float32), + hdu.data["cont"].astype(np.float32), + hdu.data["norm_flux"].astype(np.float32), + hdu.data["norm_error"].astype(np.float32), + hdu.data["norm_model"].astype(np.float32), + hdu.data["mask"].astype(np.int32), + flags) + + +PfsGAObjectNotes = makeNotesClass( + "PfsGAObjectNotes", + [] +) + + +@inheritDocstrings +class PfsGAObject(PfsFiberArray): + """Coadded spectrum of a GA target with derived quantities. + + Produced by ˙˙gapipe`` + + Parameters + ---------- + target : `pfs.datamodel.Target` + Target information. + observations : `pfs.datamodel.Observations` + Observations of the target. + wavelength : `numpy.ndarray` of `float` + Array of wavelengths. + flux : `numpy.ndarray` of `float` + Array of fluxes. + mask : `numpy.ndarray` of `int` + Array of mask pixels. + sky : `numpy.ndarray` of `float` + Array of sky values. + covar : `numpy.ndarray` of `float` + Near-diagonal (diagonal and either side) part of the covariance matrix. + covar2 : `numpy.ndarray` of `float` + Low-resolution non-sparse covariance estimate. + flags : `MaskHelper` + Helper for dealing with symbolic names for mask values. + metadata : `dict` (`str`: POD), optional + Keyword-value pairs for the header. + fluxTable : `pfs.datamodel.GAFluxTable`, optional + Table of coadded fluxes and continuum-normalized flux from contributing observations. + stellarParams: `pfs.datamodel.StellarParams`, optional + Table of measured stellar parameters. + velocityCorrections: `pfs.datamodel.VelocityCorrections`, optional + Table of velocity corrections applied to the individual visits. + abundances: `pfs.datamodel.Abundances`, optional + Table of measured abundance parameters. + paramsCovar: `numpy.ndarray` of `float`, optional + Covariance matrix for stellar parameters. + abundCovar: `numpy.ndarray` of `float`, optional + Covariance matrix for abundance parameters. + notes : `Notes`, optional + Reduction notes. + """ + + filenameFormat = ("pfsGAObject-%(catId)05d-%(tract)05d-%(patch)s-%(objId)016x" + "-%(nVisit)03d-0x%(pfsVisitHash)016x.fits") + filenameRegex = r"^pfsGAObject-(\d{5})-(\d{5})-(.*)-([0-9a-f]{16})-(\d{3})-0x([0-9a-f]{16})\.fits.*$" + filenameKeys = [("catId", int), ("tract", int), ("patch", str), ("objId", int), + ("nVisit", int), ("pfsVisitHash", int)] + NotesClass = PfsGAObjectNotes + FluxTableClass = GAFluxTable + + StellarParamsFitsExtName = "STELLARCOVAR" + AbundancesFitsExtName = "ABUNDCOVAR" + + def __init__( + self, + target, + observations, + wavelength, + flux, + mask, + sky, + covar, + covar2, + flags, + metadata=None, + fluxTable=None, + stellarParams=None, + velocityCorrections=None, + abundances=None, + paramsCovar=None, + abundCovar=None, + notes: Notes = None, + ): + super().__init__(target, observations, wavelength, flux, mask, sky, + covar, covar2, flags, metadata=metadata, fluxTable=fluxTable, notes=notes) + + self.stellarParams = stellarParams + self.velocityCorrections = velocityCorrections + self.abundances = abundances + self.paramsCovar = paramsCovar + self.abundCovar = abundCovar + + def validate(self): + """Validate that all the arrays are of the expected shape""" + super().validate() + + # TODO: write any validation code + + @classmethod + def _readImpl(cls, fits): + data = super()._readImpl(fits) + + # TODO: handle missing extensions + + data["velocityCorrections"] = VelocityCorrections.readHdu(fits) + data["stellarParams"] = StellarParams.readHdu(fits) + if cls.StellarParamsFitsExtName in fits: + data["paramsCovar"] = fits[cls.StellarParamsFitsExtName].data.astype(np.float32) + data["abundances"] = Abundances.readHdu(fits) + if cls.AbundancesFitsExtName in fits: + data["abundCovar"] = fits[cls.AbundancesFitsExtName].data.astype(np.float32) + + return data + + def _writeImpl(self, fits): + from astropy.io.fits import ImageHDU + + header = super()._writeImpl(fits) + + if self.velocityCorrections is not None: + self.velocityCorrections.writeHdu(fits) + if self.stellarParams is not None: + self.stellarParams.writeHdu(fits) + if self.paramsCovar is not None: + fits.append(ImageHDU(self.paramsCovar.astype(np.float32), + header=header, + name=self.StellarParamsFitsExtName)) + if self.abundances is not None: + self.abundances.writeHdu(fits) + if self.abundCovar is not None: + fits.append(ImageHDU(self.abundCovar.astype(np.float32), + header=header, + name=self.AbundancesFitsExtName)) + + return header diff --git a/python/pfs/datamodel/pfsFiberArray.py b/python/pfs/datamodel/pfsFiberArray.py index 56105f2..6486efb 100644 --- a/python/pfs/datamodel/pfsFiberArray.py +++ b/python/pfs/datamodel/pfsFiberArray.py @@ -44,6 +44,7 @@ class PfsFiberArray(PfsSimpleSpectrum): """ filenameFormat = None # Subclasses should override NotesClass: Type[Notes] # Subclasses should override + FluxTableClass: Type[FluxTable] = FluxTable # Subclasses may override def __init__( self, @@ -114,7 +115,7 @@ def _readImpl(cls, fits): data["covar"] = fits["COVAR"].data.astype(np.float32) data["covar2"] = fits["COVAR2"].data.astype(np.float32) try: - fluxTable = FluxTable.fromFits(fits) + fluxTable = cls.FluxTableClass.fromFits(fits) except KeyError as exc: # Only want to catch "Extension XXX not found." if not exc.args[0].startswith("Extension"): @@ -151,3 +152,4 @@ def _writeImpl(self, fits): self.notes.writeFits(fits) if self.fluxTable is not None: self.fluxTable.toFits(fits) + return header diff --git a/python/pfs/datamodel/pfsFluxReference.py b/python/pfs/datamodel/pfsFluxReference.py index 97262a0..bc8c159 100644 --- a/python/pfs/datamodel/pfsFluxReference.py +++ b/python/pfs/datamodel/pfsFluxReference.py @@ -210,7 +210,7 @@ def fits_getdata(hdulist, name, dtype=None, needHeader=False): data["identity"] = Identity.fromFits(fd) data["metadata"] = astropyHeaderToDict(fd[0].header) data["wavelength"] = WavelengthArray.fromFitsHeader( - wcsHeader, data["flux"].shape[1], dtype=np.float + wcsHeader, data["flux"].shape[1], dtype=float ) data["fitFlagNames"] = MaskHelper.fromFitsHeader(flagHeader) diff --git a/tests/test_docstrings.py b/tests/test_docstrings.py index 4e590fd..8b1553e 100644 --- a/tests/test_docstrings.py +++ b/tests/test_docstrings.py @@ -5,18 +5,21 @@ import lsst.utils.tests from pfs.datamodel.drp import PfsArm, PfsMerged, PfsReference, PfsSingle, PfsObject +from pfs.datamodel.ga import PfsGAObject class DocstringsTestCase(unittest.TestCase): def testDocstrings(self): - for cls in (PfsArm, PfsMerged, PfsReference, PfsSingle, PfsObject): + for cls in (PfsArm, PfsMerged, PfsReference, PfsSingle, PfsObject, PfsGAObject): for name, attr in inspect.getmembers(cls): if not hasattr(attr, "__doc__") or not attr.__doc__: continue docstring = attr.__doc__ for base in cls.__mro__[1:-1]: - self.assertNotIn(base.__name__, docstring, - f"{cls.__name__}.{name}.__doc__ contains {base.__name__}: {docstring}") + if not name.endswith("Class"): + self.assertNotIn( + base.__name__, docstring, + f"{cls.__name__}.{name}.__doc__ contains {base.__name__}: {docstring}") class TestMemory(lsst.utils.tests.MemoryTestCase): diff --git a/tests/test_pfsGAObject.py b/tests/test_pfsGAObject.py new file mode 100644 index 0000000..2353c64 --- /dev/null +++ b/tests/test_pfsGAObject.py @@ -0,0 +1,186 @@ +import os +import re +import numpy as np +from unittest import TestCase + +from pfs.datamodel import Target, TargetType +from pfs.datamodel import Observations +from pfs.datamodel import MaskHelper +from pfs.datamodel import GAFluxTable +from pfs.datamodel import PfsGAObject, StellarParams, Abundances, VelocityCorrections + + +class PfsGAObjectTestCase(TestCase): + """ Check the format of example datamodel files are + consistent with that specified in the corresponding + datamodel classes. + """ + + def makePfsGAObject(self): + """Construct a PfsGAObject with dummy values for testing.""" + + catId = 12345 + tract = 1 + patch = '1,1' + objId = 123456789 + ra = -100.63654 + dec = -68.591576 + targetType = TargetType.SCIENCE + + target = Target(catId, tract, patch, objId, ra, dec, targetType) + + visit = np.array([83219, 83220]) + arm = np.array(['b', 'm']) + spectrograph = np.array([1, 1]) + pfsDesignId = np.array([8854764194165386399, 8854764194165386400]) + fiberId = np.array([476, 476]) + pfiNominal = np.array([[ra, dec], [ra, dec]]) + pfiCenter = np.array([[ra, dec], [ra, dec]]) + + observations = Observations(visit, arm, spectrograph, pfsDesignId, fiberId, pfiNominal, pfiCenter) + + npix = 4096 + wavelength = np.concatenate([ + np.linspace(380, 650, npix, dtype=np.float32), + np.linspace(710, 885, npix, dtype=np.float32) + ]) + flux = np.zeros_like(wavelength) + error = np.zeros_like(wavelength) + model = np.zeros_like(wavelength) + cont = np.zeros_like(wavelength) + norm_flux = np.zeros_like(wavelength) + norm_error = np.zeros_like(wavelength) + norm_model = np.zeros_like(wavelength) + mask = np.zeros_like(wavelength, dtype=np.int32) + sky = np.zeros_like(wavelength) + covar = np.zeros((3, wavelength.size), dtype=np.float32) # Tridiagonal covariance matrix of flux + covar2 = np.zeros((1, 1), dtype=np.float32) # ? + + flags = MaskHelper() # + metadata = {} # Key-value pairs to put in the header + fluxTable = GAFluxTable(wavelength, flux, error, model, cont, + norm_flux, norm_error, norm_model, + mask, flags) + + stellarParams = StellarParams( + method=np.array(['rvfit', 'rvfit', 'rvfit', 'rvfit', 'rvfit']), + frame=np.array(['helio', '', '', '', '']), + param=np.array(['v_los', 'Fe_H', 'T_eff', 'log_g', 'a_Fe']), + covarId=np.array([0, 1, 2, 3, 4]), + unit=np.array(['km s-1', 'dex', 'K', '', 'dex']), + value=np.array([0.0, 0.0, 0.0, 0.0, 0.0]), + valueErr=np.array([0.0, 0.0, 0.0, 0.0, 0.0]), + flag=np.array([False, False, False, False, False]), + status=np.array(['', '', '', '', '']), + ) + + velocityCorrections = VelocityCorrections( + visit=visit, + JD=np.zeros_like(visit, dtype=np.float32), + helio=np.zeros_like(visit, dtype=np.float32), + bary=np.zeros_like(visit, dtype=np.float32), + ) + + abundances = Abundances( + method=np.array(['chemfit', 'chemfit', 'chemfit']), + element=np.array(['Mg', 'Ti', 'Si']), + covarId=np.array([0, 1, 2]), + value=np.array([0.0, 0.0, 0.0]), + valueErr=np.array([0.0, 0.0, 0.0]), + flag=np.array([False, False, False]), + status=np.array(['', '', '']), + ) + + paramsCovar = np.eye(3, dtype=np.float32) + abundCovar = np.eye(4, dtype=np.float32) + notes = None + + return PfsGAObject(target, observations, + wavelength, flux, mask, sky, covar, covar2, + flags, metadata, + fluxTable, + stellarParams, + velocityCorrections, + abundances, + paramsCovar, + abundCovar, + notes) + + def assertPfsGAObject(self, lhs, rhs): + np.testing.assert_array_equal(lhs.observations.visit, rhs.observations.visit) + + # TODO: add more tests here + + def extractAttributes(self, cls, fileName): + matches = re.search(cls.filenameRegex, fileName) + if not matches: + self.fail( + "Unable to parse filename: {} using regex {}" + .format(fileName, cls.filenameRegex)) + + # Cannot use algorithm in PfsSpectra._parseFilename(), + # specifically cls.filenameKeys, due to ambiguity in parsing + # integers in hex format (eg objId). Need to parse cls.filenameFormat + ff = re.search(r'^[a-zA-Z]+(.*)\.fits', cls.filenameFormat)[1] + cmps = re.findall(r'-{0,1}(0x){0,1}\%\((\w+)\)\d*(\w)', ff) + fmts = [(kk, tt) for ox, kk, tt in cmps] + + d = {} + for (kk, tt), vv in zip(fmts, matches.groups()): + if tt == 'd': + ii = int(vv) + elif tt == 'x': + ii = int(vv, 16) + elif tt == 's': + ii = vv + d[kk] = ii + return d + + def test_filenameRegex(self): + d = self.extractAttributes( + PfsGAObject, + 'pfsGAObject-07621-01234-2,2-02468ace1234abcd-003-0x0123456789abcdef.fits') + self.assertEqual(d['catId'], 7621) + self.assertEqual(d['tract'], 1234) + self.assertEqual(d['patch'], '2,2') + self.assertEqual(d['objId'], 163971054118939597) + self.assertEqual(d['nVisit'], 3) + self.assertEqual(d['pfsVisitHash'], 81985529216486895) + + def test_getIdentity(self): + """Construct a PfsGAObject and get its identity.""" + + pfsGAObject = self.makePfsGAObject() + identity = pfsGAObject.getIdentity() + filename = pfsGAObject.filenameFormat % identity + + self.assertEqual('pfsGAObject-12345-00001-1,1-00000000075bcd15-002-0x05a95bc24d8ce16f.fits', filename) + + def test_validate(self): + """Construct a PfsGAObject and run validation.""" + + pfsGAObject = self.makePfsGAObject() + pfsGAObject.validate() + + def test_writeFits_fromFits(self): + """Construct a PfsGAObject and save it to a FITS file.""" + + pfsGAObject = self.makePfsGAObject() + + dirName = os.path.splitext(__file__)[0] + if not os.path.exists(dirName): + os.makedirs(dirName) + + id = pfsGAObject.getIdentity() + filename = os.path.join(dirName, pfsGAObject.filenameFormat % id) + if os.path.exists(filename): + os.unlink(filename) + + try: + pfsGAObject.writeFits(filename) + other = PfsGAObject.readFits(filename) + self.assertPfsGAObject(pfsGAObject, other) + except Exception as e: + raise e + finally: + os.unlink(filename)