Interfacing with FITS Files

As many astronomical images are stored as FITS files, GalSim includes functionality for reading and writing these files with GalSim Image instances.

We include routines for reading and writing an individual Image to/from FITS files, and also routines for handling multiple Image instances in a single FITS file.

We also have a wrapper around the FITS header information to make it work more like a Python dict, called FitsHeader.

Note

These routines are largely wrappers of the astropy.io.fits package. They are still fairly useful for connecting GalSim objects with the AstroPy API. However, they used to be critically important for providing a stable API across different PyFITS and then AstroPy versions. For instance, now the astropy.io.fits.Header API is very similar to our own FitsHeader, but we used to have many checks for different PyFITS and AstroPy versions to call things in different ways while maintaining an intuitive front-end user interface.

Reading FITS Files

galsim.fits.read(file_name=None, dir=None, hdu_list=None, hdu=None, compression='auto', read_header=False, suppress_warning=True)[source]

Construct an Image from a FITS file or HDUList.

The normal usage for this function is to read a fits file and return the image contained therein, automatically decompressing it if necessary. However, you may also pass it an HDUList, in which case it will select the indicated hdu (with the hdu parameter) from that.

Not all FITS pixel types are supported – only short, int, unsigned short, unsigned int, float, and double.

If the FITS header has keywords that start with GS_, these will be used to initialize the bounding box and WCS. If these are absent, the code will try to read whatever WCS is given in the FITS header. cf. galsim.wcs.readFromFitsHeader. The default bounding box will have (xmin,ymin) at (1,1). The default WCS, if there is no WCS information in the FITS file, will be PixelScale(1.0).

This function is called as im = galsim.fits.read(...)

Parameters:
  • file_name – The name of the file to read in. [Either file_name or hdu_list is required.]

  • dir – Optionally a directory name can be provided if file_name does not already include it. [default: None]

  • hdu_list – Either an astropy.io.fits.HDUList, astropy.io.fits.PrimaryHDU, or astropy.io.fits..ImageHDU. In the former case, the hdu in the list will be selected. In the latter two cases, the hdu parameter is ignored. [Either file_name or hdu_list is required.]

  • hdu – The number of the HDU to return. [default: None, which means to return either the primary or first extension as appropriate for the given compression. (e.g. for ‘rice’, the first extension is the one you normally want.)]

  • compression

    Which decompression scheme to use (if any). Options are:

    • None or ‘none’ = no decompression

    • ’rice’ = use rice decompression in tiles

    • ’gzip’ = use gzip to decompress the full file

    • ’bzip2’ = use bzip2 to decompress the full file

    • ’gzip_tile’ = use gzip decompression in tiles

    • ’hcompress’ = use hcompress decompression in tiles

    • ’plio’ = use plio decompression in tiles

    • ’auto’ = determine the decompression from the extension of the file name (requires file_name to be given).

      • ’.fz’ => ‘rice’

      • ’.gz’ => ‘gzip’

      • ’.bz2’ => ‘bzip2’

      • otherwise None

    [default: ‘auto’]

  • read_header – Whether to read the header and store it as image.header.[default: False]

  • suppress_warning – Whether to suppress a warning that the WCS could not be read from the FITS header, so the WCS defaulted to either a PixelScale or AffineTransform. [default: True]

Returns:

the image as an Image instance.

galsim.fits.readMulti(file_name=None, dir=None, hdu_list=None, compression='auto', read_headers=False, suppress_warning=True)[source]

Construct a list of Image instances from a FITS file or HDUList.

The normal usage for this function is to read a fits file and return a list of all the images contained therein, automatically decompressing them if necessary. However, you may also pass it an HDUList, in which case it will build the images from these directly.

Not all FITS pixel types are supported – only short, int, unsigned short, unsigned int, float, and double.

If the FITS header has keywords that start with GS_, these will be used to initialize the bounding box and WCS. If these are absent, the code will try to read whatever WCS is given in the FITS header. cf. galsim.wcs.readFromFitsHeader. The default bounding box will have (xmin,ymin) at (1,1). The default WCS, if there is no WCS information in the FITS file, will be PixelScale(1.0).

This function is called as im = galsim.fits.readMulti(...)

Note

This function along with writeMulti can be used to effect the equivalent of a simple version of fpack or funpack. To Rice compress a fits file, you can call:

fname = 'some_image_file.fits'
galsim.fits.writeMulti(galsim.fits.readMulti(fname, read_headers=True), fname+'.fz')

To uncompress:

fname = 'some_image_file.fits.fz'
galsim.fits.writeMulti(galsim.fits.readMulti(fname, read_headers=True), fname[:-3])
Parameters:
  • file_name – The name of the file to read in. [Either file_name or hdu_list is required.]

  • dir – Optionally a directory name can be provided if file_name does not already include it. [default: None]

  • hdu_list – An astropy.io.fits.HDUList from which to read the images. [Either file_name or hdu_list is required.]

  • compression

    Which decompression scheme to use (if any). Options are:

    • None or ‘none’ = no decompression

    • ’rice’ = use rice decompression in tiles

    • ’gzip’ = use gzip to decompress the full file

    • ’bzip2’ = use bzip2 to decompress the full file

    • ’gzip_tile’ = use gzip decompression in tiles

    • ’hcompress’ = use hcompress decompression in tiles

    • ’plio’ = use plio decompression in tiles

    • ’auto’ = determine the decompression from the extension of the file name (requires file_name to be given).

      • ’.fz’ => ‘rice’

      • ’.gz’ => ‘gzip’

      • ’.bz2’ => ‘bzip2’

      • otherwise None

    [default: ‘auto’]

  • read_headers – Whether to read the headers and store them as image.header.[default: False]

  • suppress_warning – Whether to suppress a warning that the WCS could not be read from the FITS header, so the WCS defaulted to either a PixelScale or AffineTransform. [default: True]

Returns:

a Python list of Image instances.

galsim.fits.readCube(file_name=None, dir=None, hdu_list=None, hdu=None, compression='auto', suppress_warning=True)[source]

Construct a Python list of Image instances from a FITS data cube.

Not all FITS pixel types are supported – only short, int, unsigned short, unsigned int, float, and double.

If the FITS header has keywords that start with GS_, these will be used to initialize the bounding box and WCS. If these are absent, the code will try to read whatever WCS is given in the FITS header. cf. galsim.wcs.readFromFitsHeader. The default bounding box will have (xmin,ymin) at (1,1). The default WCS, if there is no WCS information in the FITS file, will be PixelScale(1.0).

This function is called as image_list = galsim.fits.readCube(...)

Parameters:
  • file_name – The name of the file to read in. [Either file_name or hdu_list is required.]

  • dir – Optionally a directory name can be provided if file_name does not already include it. [default: None]

  • hdu_list – Either an astropy.io.fits.HDUList, an astropy.io.fits.PrimaryHDU, or astropy.io.fits.ImageHDU. In the former case, the hdu in the list will be selected. In the latter two cases, the hdu parameter is ignored. [Either file_name or hdu_list is required.]

  • hdu – The number of the HDU to return. [default: None, which means to return either the primary or first extension as appropriate for the given compression. (e.g. for rice, the first extension is the one you normally want.)]

  • compression

    Which decompression scheme to use (if any). Options are:

    • None or ‘none’ = no decompression

    • ’rice’ = use rice decompression in tiles

    • ’gzip’ = use gzip to decompress the full file

    • ’bzip2’ = use bzip2 to decompress the full file

    • ’gzip_tile’ = use gzip decompression in tiles

    • ’hcompress’ = use hcompress decompression in tiles

    • ’plio’ = use plio decompression in tiles

    • ’auto’ = determine the decompression from the extension of the file name (requires file_name to be given).

      • ’.fz’ => ‘rice’

      • ’.gz’ => ‘gzip’

      • ’.bz2’ => ‘bzip2’

      • otherwise None

    [default: ‘auto’]

  • suppress_warning – Whether to suppress a warning that the WCS could not be read from the FITS header, so the WCS defaulted to either a PixelScale or AffineTransform. [default: True]

Returns:

a Python list of Image instances.

galsim.fits.readFile(file_name, dir=None, hdu=None, compression='auto')[source]

Read in a Pyfits hdu_list from a FITS file, taking care of the GalSim compression options.

If you want to do something different with an hdu or hdu_list than one of our other read functions, you can use this function. It handles the compression options in the standard GalSim way and just returns the hdu (and hdu_list) for you to use as you see fit.

This function is called as:

>>> hdu, hdu_list, fin = galsim.fits.readFile(...)

The first item in the returned tuple is the specified hdu (or the primary if none was specifically requested). The other two are returned so you can properly close them. They are the full HDUList and possibly a file handle. The appropriate cleanup can be done with:

>>> galsim.fits.closeHDUList(hdu_list, fin)
Parameters:
  • file_name – The name of the file to read in.

  • dir – Optionally a directory name can be provided if file_name does not already include it. [default: None]

  • hdu – The number of the HDU to return. [default: None, which means to return either the primary or first extension as appropriate for the given compression. (e.g. for rice, the first extension is the one you normally want.)]

  • compression

    Which decompression scheme to use (if any). Options are:

    • None or ‘none’ = no decompression

    • ’rice’ = use rice decompression in tiles

    • ’gzip’ = use gzip to decompress the full file

    • ’bzip2’ = use bzip2 to decompress the full file

    • ’gzip_tile’ = use gzip decompression in tiles

    • ’hcompress’ = use hcompress decompression in tiles

    • ’plio’ = use plio decompression in tiles

    • ’auto’ = determine the decompression from the extension of the file name (requires file_name to be given).

      • ’.fz’ => ‘rice’

      • ’.gz’ => ‘gzip’

      • ’.bz2’ => ‘bzip2’

      • otherwise None

    [default: ‘auto’]

Returns:

(hdu, hdu_list, fin).

Return type:

a tuple with three items

galsim.fits.closeHDUList(hdu_list, fin)[source]

If necessary, close the file handle that was opened to read in the hdu_list

Writing FITS Files

galsim.fits.write(image, file_name=None, dir=None, hdu_list=None, clobber=True, compression='auto')[source]

Write a single image to a FITS file.

Write the Image instance image to a FITS file, with details depending on the arguments. This function can be called directly as galsim.fits.write(image, ...), with the image as the first argument, or as an Image method: image.write(...).

Parameters:
  • image – The Image to write to file. Per the description of this method, it may be given explicitly via galsim.fits.write(image, ...) or the method may be called directly as an image method, image.write(...). Note that if the image has a ‘header’ attribute containing a FitsHeader, then the FitsHeader is written to the header in the PrimaryHDU, followed by the WCS as usual.

  • file_name – The name of the file to write to. [Either file_name or hdu_list is required.]

  • dir – Optionally a directory name can be provided if file_name does not already include it. [default: None]

  • hdu_list – An astropy.io.fits.HDUList. If this is provided instead of file_name, then the Image is appended to the end of the HDUList as a new HDU. In that case, the user is responsible for calling either hdu_list.writeto(...) or galsim.fits.writeFile(...) afterwards. [Either file_name or hdu_list is required.]

  • clobber – Setting clobber=True will silently overwrite existing files. [default: True]

  • compression

    Which compression scheme to use (if any). Options are:

    • None or ‘none’ = no compression

    • ’rice’ = use rice compression in tiles (preserves header readability)

    • ’gzip’ = use gzip to compress the full file

    • ’bzip2’ = use bzip2 to compress the full file

    • ’gzip_tile’ = use gzip in tiles (preserves header readability)

    • ’hcompress’ = use hcompress in tiles (only valid for 2-d images)

    • ’plio’ = use plio compression in tiles (only valid for pos integer data)

    • ’auto’ = determine the compression from the extension of the file name (requires file_name to be given):

      • ’.fz’ => ‘rice’

      • ’.gz’ => ‘gzip’

      • ’.bz2’ => ‘bzip2’

      • otherwise None

    [default: ‘auto’]

galsim.fits.writeMulti(image_list, file_name=None, dir=None, hdu_list=None, clobber=True, compression='auto')[source]

Write a Python list of images to a multi-extension FITS file.

The details of how the images are written to file depends on the arguments.

Note

This function along with readMulti can be used to effect the equivalent of a simple version of fpack or funpack. To Rice compress a fits file, you can call:

fname = 'some_image_file.fits'
galsim.fits.writeMulti(galsim.fits.readMulti(fname, read_headers=True), fname+'.fz')

To uncompress:

fname = 'some_image_file.fits.fz'
galsim.fits.writeMulti(galsim.fits.readMulti(fname, read_headers=True), fname[:-3])
Parameters:
  • image_list – A Python list of Image instances. (For convenience, some items in this list may be HDUs already. Any Image will be converted into an astropy.io.fits.HDU.)

  • file_name – The name of the file to write to. [Either file_name or hdu_list is required.]

  • dir – Optionally a directory name can be provided if file_name does not already include it. [default: None]

  • hdu_list – An astropy.io.fits.HDUList. If this is provided instead of file_name, then the Image is appended to the end of the HDUList as a new HDU. In that case, the user is responsible for calling either hdu_list.writeto(...) or galsim.fits.writeFile(...) afterwards. [Either file_name or hdu_list is required.]

  • clobber – Setting clobber=True will silently overwrite existing files. [default: True]

  • compression

    Which compression scheme to use (if any). Options are:

    • None or ‘none’ = no compression

    • ’rice’ = use rice compression in tiles (preserves header readability)

    • ’gzip’ = use gzip to compress the full file

    • ’bzip2’ = use bzip2 to compress the full file

    • ’gzip_tile’ = use gzip in tiles (preserves header readability)

    • ’hcompress’ = use hcompress in tiles (only valid for 2-d images)

    • ’plio’ = use plio compression in tiles (only valid for pos integer data)

    • ’auto’ = determine the compression from the extension of the file name (requires file_name to be given):

      • ’.fz’ => ‘rice’

      • ’.gz’ => ‘gzip’

      • ’.bz2’ => ‘bzip2’

      • otherwise None

    [default: ‘auto’]

galsim.fits.writeCube(image_list, file_name=None, dir=None, hdu_list=None, clobber=True, compression='auto')[source]

Write a Python list of images to a FITS file as a data cube.

The details of how the images are written to file depends on the arguments. Unlike for writeMulti, when writing a data cube it is necessary that each Image in image_list has the same size (nx, ny). No check is made to confirm that all images have the same origin and pixel scale (or WCS).

In fact, the WCS of the first image is the one that gets put into the FITS header (since only one WCS can be put into a FITS header). Thus, if the images have different WCS functions, only the first one will be rendered correctly by plotting programs such as ds9. The FITS standard does not support any way to have the various images in a data cube to have different WCS solutions.

Parameters:
  • image_list – The image_list can also be either an array of NumPy arrays or a 3d NumPy array, in which case this is written to the fits file directly. In the former case, no explicit check is made that the NumPy arrays are all the same shape, but a NumPy exception will be raised which we let pass upstream unmolested.

  • file_name – The name of the file to write to. [Either file_name or hdu_list is required.]

  • dir – Optionally a directory name can be provided if file_name does not already include it. [default: None]

  • hdu_list – An astropy.io.fits.HDUList. If this is provided instead of file_name, then the cube is appended to the end of the HDUList as a new HDU. In that case, the user is responsible for calling either hdu_list.writeto(...) or galsim.fits.writeFile(...) afterwards. [Either file_name or hdu_list is required.]

  • clobber – Setting clobber=True will silently overwrite existing files. [default: True]

  • compression

    Which compression scheme to use (if any). Options are:

    • None or ‘none’ = no compression

    • ’rice’ = use rice compression in tiles (preserves header readability)

    • ’gzip’ = use gzip to compress the full file

    • ’bzip2’ = use bzip2 to compress the full file

    • ’gzip_tile’ = use gzip in tiles (preserves header readability)

    • ’hcompress’ = use hcompress in tiles (only valid for 2-d images)

    • ’plio’ = use plio compression in tiles (only valid for pos integer data)

    • ’auto’ = determine the compression from the extension of the file name (requires file_name to be given):

      • ’.fz’ => ‘rice’

      • ’.gz’ => ‘gzip’

      • ’.bz2’ => ‘bzip2’

      • otherwise None

    [default: ‘auto’]

galsim.fits.writeFile(file_name, hdu_list, dir=None, clobber=True, compression='auto')[source]

Write a Pyfits hdu_list to a FITS file, taking care of the GalSim compression options.

If you have used the write(), writeMulti() or writeCube() functions with the hdu_list option rather than writing directly to a file, you may subsequently use the command hdu_list.writeto(...). However, it may be more convenient to use this function, writeFile() instead, since it treats the compression option consistently with how that option is handled in the above functions.

Parameters:
  • file_name – The name of the file to write to.

  • hdu_list – An astropy.io.fits.HDUList.

  • dir – Optionally a directory name can be provided if file_name does not already include it. [default: None]

  • clobber – Setting clobber=True will silently overwrite existing files. [default: True]

  • compression

    Which compression scheme to use (if any). Options are:

    • None or ‘none’ = no compression

    • ’gzip’ = use gzip to compress the full file

    • ’bzip2’ = use bzip2 to compress the full file

    • ’auto’ = determine the compression from the extension of the file name (requires file_name to be given):

      • ’.gz’ => ‘gzip’

      • ’.bz2’ => ‘bzip2’

      • otherwise None

    Note that the other options, such as ‘rice’, that operate on the image directly are not available at this point. If you want to use one of them, it must be applied when writing each hdu. [default: ‘auto’]

FITS Headers

class galsim.fits.FitsHeader(header=None, file_name=None, dir=None, hdu_list=None, hdu=None, compression='auto', text_file=False)[source]

A class storing key/value pairs from a FITS Header

This class works a lot like the regular read() function, but rather than returning the image part of the FITS file, it gives you access to the header information.

After construction, you can access a header value by:

>>> value = fits_header[key]

or write to it with:

>>> fits_header[key] = value                # If you just want to set a value.
>>> fits_header[key] = (value, comment)     # If you want to include a comment field.

In fact, most of the normal functions available for a dict are available::

>>> keys = fits_header.keys()
>>> items = fits_header.items()
>>> for key in fits_header:
>>>     value = fits_header[key]
>>> value = fits_header.get(key, default)
>>> del fits_header[key]
>>> etc.

Note

This used to be a particularly useful abstraction, since PyFITS and then AstroPy used to keep changing their syntax for how to write to a fits header rather often, so this class had numerous checks for which version of PPyFITS or AstroPy was installed and call things the right way depending on the version. Thus, it was able to maintain a stable, intuitive API that would work with any version on the backend. We no longer support PyFITS or older versions of AstroPy, so now much of the syntax of this class is very similar in interface to the current version of astropy.io.fits.Header. Indeed it is now a rather light wrapper around their Header class with just a few convenience features to make it easier to work with GalSim objects.

The underlying Header object is available as a .header attribute:

>>> apy_header = fits_header.header

A FitsHeader may be constructed from a file name, an open PyFits (or astropy.io.fits) HDUList object, or a PyFits (or astropy.io.fits) Header object. It can also be constructed with no parameters, in which case a blank Header will be constructed with no keywords yet if you want to add the keywords you want by hand.:

>>> h1 = galsim.FitsHeader(file_name = file_name)
>>> h2 = galsim.FitsHeader(header = header)
>>> h3 = galsim.FitsHeader(hdu_list = hdu_list)
>>> h4 = galsim.FitsHeader()

For convenience, the first parameter may be unnamed as either a header or a file_name:

>>> h1 = galsim.FitsHeader(file_name)
>>> h2 = galsim.FitsHeader(header)
Parameters:
  • header – An astropy.io.fits.Header object or in fact any dict-like object or list of (key,value) pairs. [default: None]

  • file_name – The name of the file to read in. [default: None]

  • dir – Optionally a directory name can be provided if file_name does not already include it. [default: None]

  • hdu_list – Either an astropy.io.fits.HDUList, an astropy.io.fits.PrimaryHDU, or astropy.io.fits.ImageHDU. In the former case, the hdu in the list will be selected. In the latter two cases, the hdu parameter is ignored. [default: None]

  • hdu – The number of the HDU to return. [default: None, which means to return either the primary or first extension as appropriate for the given compression. (e.g. for rice, the first extension is the one you normally want.)]

  • compression

    Which decompression scheme to use (if any). Options are:

    • None or ‘none’ = no decompression

    • ’rice’ = use rice decompression in tiles

    • ’gzip’ = use gzip to decompress the full file

    • ’bzip2’ = use bzip2 to decompress the full file

    • ’gzip_tile’ = use gzip decompression in tiles

    • ’hcompress’ = use hcompress decompression in tiles

    • ’plio’ = use plio decompression in tiles

    • ’auto’ = determine the decompression from the extension of the file name (requires file_name to be given).

      • ’.fz’ => ‘rice’

      • ’.gz’ => ‘gzip’

      • ’.bz2’ => ‘bzip2’

      • otherwise None

    [default: ‘auto’]

  • text_file – Normally a file is taken to be a fits file, but you can also give it a text file with the header information (like the .head file output from SCamp). In this case you should set text_file = True to tell GalSim to parse the file this way. [default: False]

append(key, value='', comment=None, useblanks=True)[source]

Append an item to the end of the header.

This breaks convention a bit by treating the header more like a list than a dict, but sometimes that is necessary to get the header structured the way you want it.

Parameters:
  • key – The key of the entry to append

  • value – The value of the entry to append [default: ‘’]

  • comment – A comment field if desired [default: None]

  • useblanks – If there are blank entries currently at the end, should they be overwritten with the new entry? [default: True]

clear()[source]

Clear all values in the header. Works like dict.clear.

comment(key)[source]

Get the comment field for the given key.

Parameter:

key: The header key for which to get the comment field.

Returns:

the comment field.

extend(other, replace=False, useblanks=True)[source]

Extend this FitsHeader with items from another FitsHeader.

Equivalent to appending all the other’s items to the end of this one with the exception that it ignores items that are already in self. If you want to replace existing values rather than ignore duplicates, use replace=True.

Parameters:
  • other – Another FitsHeader object.

  • replace – Replace duplicate entries rather than ignore them. [default: False]

  • useblanks – If there are blank entries currently at the end, should they be overwritten with the new entry? [default: True]

get(key, default=None)[source]

Get the value of a given key. Works like dict.get.

Parameters:
  • key – The header key for which to get the value

  • default – Optionally, A value to use if the key is not present. [default: None]

Returns:

the value of the given key

items()[source]

Get all header items. Works like dict.items.

Returns:

A list of (key, value) tuples.

iteritems()[source]

Synonym for self.items()

iterkeys()[source]

Synonym for self.keys()

itervalues()[source]

Synonym for self.values()

keys()[source]

Get all header keys. Works like dict.keys

Returns:

A list of keys.

pop(key, default=None)[source]

Pop off a value from the header. Works like dict.pop.

Parameters:
  • key – The header key for which to get the value

  • default – Optionally, A value to use if the key is not present. [default: None]

Returns:

the value of the given key

update(dict2)[source]

Update the header with a dict-like object. Works like dict.update.

If there are any items in dict2 that are duplicates of items already in the header, the current items will ber overwritten.

Parameters:

dict2 – Another header or dict-like object with keys and values to update.

values()[source]

Get all header values. Works like dict.values.

Returns:

A list of values.