CamiFITS.jl

FITS stands for Flexible Image Transport System. This is an open standard origionally developed for the astronomy community to store telescope images together with tables of spectral information. Over the years it has developed into a scientific standard [W. D. Pence et al., A&A, 524 (2010) A42]. The standard is maintained by the FITS Support Office at NASA/GSFC [FITS standard - Version 4.0]. The website also offers a FITS Conformance Verifier.

CamiFITS offers the basic FITS functionality for scientific users not requiring celestal coordinates. Optional Conforming Extentions are under development. The user can create, read and extend .fits files as well as create, edit and delete user-defined metainformation.

Disclaimer 2023-03-13: The author is currently writing the manual. In this process the code is tested, both with regard to FITS conformance and runtest coverage. Known issues remain to be solved and the package certainly did not reach a stable form.

Table of contents

Installation

The package is installed using the Julia package manager

julia> using Pkg; Pkg.add("CamiFITS")

julia> using CamiFITS

Manual

Introduction

A FITS file consists of a sequence of one or more header-data-units (FITS_HDUs), each containing a FITS_data block preceeded by FITS_header records of metainformation.

We distinguish between IMAGE and TABLE HDU data types. The first HDU in a .fits file is called the PRIMARY HDU.

By the command f = fits_read("filnam.fits") we asign a collection of FITS_HDU objects from the file "filnam.fits" to the variable f.

The elements of the HDU collection f are f[1], f[2], ..., with f[1] representing the PRIMARY HDU. The structure of HDU f[i] can be printed in formated form using the command FITS_info(f[i]).

FITS files are created using the command fits_create.

The simplest FITS file

The simplest file conforming to the FITS standard consists of a single HDU containing an empty data field of the type Any[].

julia> filename = "empty.fits";

julia> fits_create(filename; protect=false)

julia> f = fits_read(filename);

julia> fits_info(f[1])

File: minimal.fits
hdu: 1
hdutype: PRIMARY
DataType: Any
Datasize: (0,)

Metainformation:
SIMPLE  =                    T / file does conform to FITS standard
NAXIS   =                    0 / number of data axes
EXTEND  =                    T / FITS dataset may contain extensions
COMMENT    Primary FITS HDU    / http://fits.gsfc.nasa.gov/iaufwg
END

Any[]

julia> rm(filename); f = nothing

The FITS file for a single matrix

We first create the data field in the form of a 3x3 matrix:

julia> filename = "matrix.fits";

julia> data = [11,21,31,12,22,23,13,23,33];

julia> data = reshape(data,(3,3,1))
3×3×1 Array{Int64, 3}:
[:, :, 1] =
 11  12  13
 21  22  23
 31  23  33

We next create and inspact the FITS file for the matrix data

julia> fits_create(filename, data; protect=false)

julia> f = fits_read(filename);

julia> fits_info(f[1])

File: matrix.fits
hdu: 1
hdutype: PRIMARY
DataType: Int64
Datasize: (3, 3, 1)

Metainformation:
SIMPLE  =                    T / file does conform to FITS standard
BITPIX  =                   64 / number of bits per data pixel
NAXIS   =                    3 / number of data axes
NAXIS1  =                    3 / length of data axis 1
NAXIS2  =                    3 / length of data axis 2
NAXIS3  =                    1 / length of data axis 3
BZERO   =                  0.0 / offset data range to that of unsigned integer
BSCALE  =                  1.0 / default scaling factor
EXTEND  =                    T / FITS dataset may contain extensions
COMMENT    Primary FITS HDU    / http://fits.gsfc.nasa.gov/iaufwg
END

3×3×1 Array{Int64, 3}:
[:, :, 1] =
 11  12  13
 21  22  23
 31  23  33

julia> rm(filename); f = nothing

The keywords NAXIS1, NAXIS2 and NAXIS3 represent the dimensions of the data matrix in $x$, $y$ and $z$ direction.

The matrix elements are referred to as pixels and their bit size is represented by the keyword BITPIX. In the above example the pixel value is given by the matrix indices.

API

FITS - Types

CamiFITS.FITS_HDUType
FITS_HDU{T,V}

Object to hold a single "Header-Data Unit" (HDU).

The fields are

  • .filename: name of the corresponding FITS file (::String)
  • .hduindex:: identifier (a file may contain more than one HDU) (:Int)
  • .header: the header object where T=FITS_header (::T)
  • .dataobject: the data object where V=FITS_data (::V)
CamiFITS.FITS_headerType
FITS_header

Object to hold the header information of a FITS_HDU.

The fields are:

  • .hduindex: identifier (a file may contain more than one HDU) (::Int)
  • .records: the header formated as an array of strings of 80 ASCII characters (::Array{String,1})
  • .keys: keys[i] - key corresponding to records[i] (record of index i) (::Array{String,1})
  • .values: value[i] - corresponding to records[i] (::Array{Any,1})
  • .comments: comments[i] - comment corresponding to records[i] (::String)
  • .dict: Dictionary key[i] => value[i] (::Dict{String,Any})
  • .maps: Dictionary key[i] => i (::Dict{String,Int})
CamiFITS.FITS_dataType
FITS_data

Object to hold the data of the FITS_HDU of given hduindex and hdutype.

The fields are:

  • .hduindex: identifier (a file may contain more than one HDU) (::Int)
  • .hdutype: accepted types are 'PRIMARY', 'IMAGE' and 'TABLE' (::String)
  • .data: in the from appropriate for the hdutype (::Any)
CamiFITS.FITS_tableType
FITS_table

Object to hold the data of a TABLE HDU (a FITS_HDU for ASCII tables). It contains the data in the form of records (rows) of ASCII strings.

The fields are:

  • .hduindex: identifier (a file may contain more than one HDU) (::Int)
  • .rows: the table formated as an array of rows of ASCII strings (::Array{String,1})
CamiFITS.FITS_nameType
FITS_name

FITS object to hold the decomposed name of a .fits file.

The fields are:

  • .name: for 'p#.fits' this is 'p#.fits' (::String)
  • .prefix: for 'p#.fits' this is 'p' (::String)
  • .numerator: for 'p#.fits' this is '#', a serial number (e.g., '3') or a range (e.g., '3-7') (::String)
  • .extension: for 'p#.fits' this is '.fits' (::String)

FITS - HDU Methods

CamiFITS.fits_infoMethod
fits_info(hdu)

Print metafinformation and data of given FITS_HDU

Example:

strExample = "remove.fits"
data = [11,21,31,12,22,23,13,23,33]
data = reshape(data,(3,3,1))
fits_create(strExample, data; protect=false)

f = fits_read(strExample)
fits_info(f[1])

  File: remove.fits
  hdu: 1
  hdutype: PRIMARY
  DataType: Int64
  Datasize: (3, 3, 1)

  Metainformation:
  SIMPLE  =                    T / file does conform to FITS standard
  BITPIX  =                   64 / number of bits per data pixel
  NAXIS   =                    3 / number of data axes
  NAXIS1  =                    3 / length of data axis 1
  NAXIS2  =                    3 / length of data axis 2
  NAXIS3  =                    1 / length of data axis 3
  BZERO   =                  0.0 / offset data range to that of unsigned integer
  BSCALE  =                  1.0 / default scaling factor
  EXTEND  =                    T / FITS dataset may contain extensions
  COMMENT    Primary FITS HDU    / http://fits.gsfc.nasa.gov/iaufwg
  END

  3×3×1 Array{Int64, 3}:
  [:, :, 1] =
   11  12  13
   21  22  23
   31  23  33
CamiFITS.parse_FITS_TABLEMethod
parse_FITS_TABLE(hdu)

Parse FITS_TABLE (ASCII table) into a Vector of its columns for further processing by the user. Default formatting in ISO 2004 FORTRAN data format specified by keys "TFORMS1" - "TFORMSn"). Display formatting in ISO 2004 FORTRAN data format ("TDISP1" - "TDISPn") prepared for user editing.

Example:

strExample = "example.fits"
data = [10, 20, 30]
fits_create(strExample, data; protect=false)

t1 = Float16[1.01E-6,2.0E-6,3.0E-6,4.0E-6,5.0E-6]
t2 = [0x0000043e, 0x0000040c, 0x0000041f, 0x0000042e, 0x0000042f]
t3 = [1.23,2.12,3.,4.,5.]
t4 = ['a','b','c','d','e']
t5 = ["a","bb","ccc","dddd","ABCeeaeeEEEEEEEEEEEE"]
data = [t1,t2,t3,t4,t5]
fits_extend(strExample, data, "TABLE")

f = fits_read(strExample)
d = f[2].header.dict
d = [get(d,"TFORM$i",0) for i=1:5]; println(strip.(d))
  SubString{String}["'E6.1    '", "'I4      '", "'F4.2    '", "'A1      '", "'A20     '"]

f[2].dataobject.data                            # this is the table hdu
  5-element Vector{String}:
   "1.0e-6 1086 1.23 a a                    "
   "2.0e-6 1036 2.12 b bb                   "
   "3.0e-6 1055 3.0  c ccc                  "
   "4.0e-6 1070 4.0  d dddd                 "
   "5.0e-6 1071 5.0  e ABCeeaeeEEEEEEEEEEEE "

parse_FITS_TABLE(f[2])
  5-element Vector{Vector{T} where T}:
   [1.0e-6, 2.0e-6, 3.0e-6, 4.0e-6, 5.0e-6]
   [1086, 1036, 1055, 1070, 1071]
   [1.23, 2.12, 3.0, 4.0, 5.0]
   ["a", "b", "c", "d", "e"]
   ["a                   ", "bb                  ", "ccc                 ", "dddd                ", "ABCeeaeeEEEEEEEEEEEE"]

FITS - File Methods

CamiFITS.cast_FITS_nameMethod
cast_FITS_name(str::String)

Decompose the FITS filename 'filnam.fits' into its name, prefix, numerator and extension.

Examples:

strExample = "T23.01.fits"
f = cast_FITS_name(strExample)
FITS_name("T23.01", "T23.", "01", ".fits")

f.name, f.prefix, f.numerator, f.extension
("T23.01", "T23.", "01", ".fits")
CamiFITS.fits_combineMethod
fits_combine(strFirst, strLast [; protect=true])

Copy "filenameFirst" to "filenameLast" (with mandatory ".fits" extension)

Key:

  • protect::Bool: overwrite protection

Example:

fits_combine("T01.fits", "T22.fits")
  'T01-T22.fits': file created
CamiFITS.fits_copyFunction
fits_copy(filenameA [, filenameB="" [; protect=true]])

Copy "filenameA" to "filenameB" (with mandatory ".fits" extension) Key:

  • protect::Bool: overwrite protection

Examples:

fits_copy("T01.fits")
  'T01.fits' was saved as 'T01 - Copy.fits'

fits_copy("T01.fits", "T01a.fits")
  FitsError: 'T01a.fits' in use (set ';protect=false' to lift overwrite protection)

fits_copy("T01.fits", "T01a.fits"; protect=false)
  'T01.fits' was saved as 'T01a.fits'
CamiFITS.fits_createFunction
fits_create(filename [, data [; protect=true]])

Create FITS file of given filename [, optional data block [, default overwrite protection]] and return Array of HDUs. Key:

  • protect::Bool: overwrite protection

Examples:

strExample = "minimal.fits"
fits_create(strExample; protect=false)

f = fits_read(strExample)
a = f[1].dataobject.data
b = f[1].header.keys
println(a);println(b)
  Any[]
  ["SIMPLE", "NAXIS", "EXTEND", "COMMENT", "END"]

strExample = "remove.fits"
data = [11,21,31,12,22,23,13,23,33]
data = reshape(data,(3,3,1))
fits_create(strExample, data; protect=false)

f = fits_read(strExample)
fits_info(f[1])

  File: remove.fits
  hdu: 1
  hdutype: PRIMARY
  DataType: Int64
  Datasize: (3, 3, 1)

  Metainformation:
  SIMPLE  =                    T / file does conform to FITS standard
  BITPIX  =                   64 / number of bits per data pixel
  NAXIS   =                    3 / number of data axes
  NAXIS1  =                    3 / length of data axis 1
  NAXIS2  =                    3 / length of data axis 2
  NAXIS3  =                    1 / length of data axis 3
  BZERO   =                  0.0 / offset data range to that of unsigned integer
  BSCALE  =                  1.0 / default scaling factor
  EXTEND  =                    T / FITS dataset may contain extensions
  COMMENT    Primary FITS HDU    / http://fits.gsfc.nasa.gov/iaufwg
  END

  3×3×1 Array{Int64, 3}:
  [:, :, 1] =
   11  12  13
   21  22  23
   31  23  33
CamiFITS.fits_extendFunction
fits_extend(filename, data_extend, hdutype="IMAGE")

Extend the FITS file of given filename with the data of hdutype from data_extend and return Array of HDUs.

Examples:

strExample = "test_example.fits"
data = [0x0000043e, 0x0000040c, 0x0000041f]
fits_create(strExample, data; protect=false)

f = fits_read(strExample)
a = Float16[1.01E-6,2.0E-6,3.0E-6,4.0E-6,5.0E-6]
b = [0x0000043e, 0x0000040c, 0x0000041f, 0x0000042e, 0x0000042f]
c = [1.23,2.12,3.,4.,5.]
d = ['a','b','c','d','e']
e = ["a","bb","ccc","dddd","ABCeeaeeEEEEEEEEEEEE"]
data = [a,b,c,d,e]
fits_extend(strExample, data, "TABLE")

f = fits_read(strExample)
f[2].dataobject.data
  5-element Vector{String}:
   "1.0e-6 1086 1.23 a a                    "
   "2.0e-6 1036 2.12 b bb                   "
   "3.0e-6 1055 3.0  c ccc                  "
   "4.0e-6 1070 4.0  d dddd                 "
   "5.0e-6 1071 5.0  e ABCeeaeeEEEEEEEEEEEE "

rm(strExample); f = data = a = b = c = d = e = nothing
CamiFITS.fits_readMethod
fits_read(filename)

Read FITS file and return Array of FITS_HDUs

Example:

strExample = "minimal.fits"
fits_create(strExample; protect=false)

f = fits_read(strExample)
f[1].dataobject.data
  Any[]

rm(strExample); f = nothing

FITS - Key Methods

CamiFITS.fits_add_keyMethod
fits_add_key(filename, hduindex, key, value, comment)

Add a header record of given 'key, value and comment' to 'HDU[hduindex]' of file with name 'filename'

Example:

strExample="minimal.fits"
fits_create(strExample; protect=false)
fits_add_key(strExample, 1, "KEYNEW1", true, "FITS dataset may contain extension")

f = fits_read(strExample)
fits_info(f[1])

  File: minimal.fits
  hdu: 1
  hdutype: PRIMARY
  DataType: Any
  Datasize: (0,)

  Metainformation:
  SIMPLE  =                    T / file does conform to FITS standard
  NAXIS   =                    0 / number of data axes
  EXTEND  =                    T / FITS dataset may contain extensions
  COMMENT    Primary FITS HDU    / http://fits.gsfc.nasa.gov/iaufwg
  KEYNEW1 =                    T / FITS dataset may contain extension
  END

  Any[]
CamiFITS.fits_delete_keyMethod
fits_delete_key(filename, hduindex, key)

Delete a header record of given key, value and comment to FITS_HDU[hduindex] of file with name 'filename'

Examples:

strExample="minimal.fits"
fits_create(strExample; protect=false)
fits_add_key(strExample, 1, "KEYNEW1", true, "this is record 5")

f = fits_read(strExample)
get(f[1].header.maps,"KEYNEW1",0)
  5

fits_delete_key(strExample, 1, "KEYNEW1")

f = fits_read(strExample)
get(f[1].header.maps,"KEYNEW1",0)
  0

fits_delete_key(filnam, 1, "NAXIS")
 'NAXIS': cannot be deleted (key protected under FITS standard)
CamiFITS.fits_edit_keyMethod
fits_edit_key(filename, hduindex, key, value, comment)

Edit a header record of given 'key, value and comment' to 'HDU[hduindex]' of file with name 'filename'

Example:

data = DateTime("2020-01-01T00:00:00.000")
strExample="minimal.fits"
fits_create(strExample; protect=false)
fits_add_key(strExample, 1, "KEYNEW1", true, "this is record 5")
fits_edit_key(strExample, 1, "KEYNEW1", data, "record 5 changed to a DateTime type")

f = fits_read(strExample)
fits_info(f[1])

  File: minimal.fits
  hdu: 1
  hdutype: PRIMARY
  DataType: Any
  Datasize: (0,)

  Metainformation:
  SIMPLE  =                    T / file does conform to FITS standard
  NAXIS   =                    0 / number of data axes
  EXTEND  =                    T / FITS dataset may contain extensions
  COMMENT    Primary FITS HDU    / http://fits.gsfc.nasa.gov/iaufwg
  KEYNEW1 = '2020-01-01T00:00:00' / record 5 changed to a DateTime type
  END

  Any[]
CamiFITS.fits_rename_keyMethod
fits_rename_key(filename, hduindex, keyold, kewnew)

Rename the key of a header record of file with name 'filename'

Example:

strExample="minimal.fits"
fits_create(strExample; protect=false)
fits_add_key(strExample, 1, "KEYNEW1", true, "this is record 5")
fits_rename_key(strExample, 1, "KEYNEW1",  "KEYNEW2")

f = fits_read(strExample)
fits_info(f[1])

  File: minimal.fits
  hdu: 1
  hdutype: PRIMARY
  DataType: Any
  Datasize: (0,)

  Metainformation:
  SIMPLE  =                    T / file does conform to FITS standard
  NAXIS   =                    0 / number of data axes
  EXTEND  =                    T / FITS dataset may contain extensions
  COMMENT    Primary FITS HDU    / http://fits.gsfc.nasa.gov/iaufwg
  KEYNEW2 =                    T / this is record 5
  END

  Any[]

FORTRAN

CamiFITS.FORTRAN_formatType
FORTRAN_format

Object to hold a FORTRAN format specifier decomposed in its fields.

Accepted datatype specifiers are: Aw, Iw, Fw.d, Ew.d, Dw.d

Accepted output formating specifiers are: Aw, Iw.m, Bw.m, Ow.m, Zw.m, Fw.d, Ew.dEe, ENw.d, ESw.d, Gw.dEe, Dw.dEe. Notation: w - width, m (optional) - minimum number of digits, d - number of digits to right of decimal, e - number of digits in exponent N/S (optional) indicates engineering/scientific formating of the E type.

The fields are:

  • .Type: primary FORTRAN datatype (::String)
  • .TypeChar: primary FORTRAN datatype character (::Char)
  • .EngSci: secundary datatype character - N for engineering/ S for scientific (::Union{Char,Nothing})
  • .width: width of numeric field (::Int)
  • .nmin: minimum number of digits displayed (::Int)
  • .ndec: number of digits to right of decimal (::Int)
  • .nexp: number of digits in exponent (::Int)
CamiFITS.cast_FORTRAN_formatMethod
cast_FORTRAN_format(format::String)

Decompose the format specifier format into its fields and cast this into the FORTRAN_format object. Allowed format specifiers are of the types: Aw, Iw.m, Bw.m, Ow.m, Zw.m, Fw.d, Ew.dEe, ENw.d, ESw.d, Gw.dEe, Dw.dEe, with: w - width, m(optional) - minimum number of digits, d - number of digits to right of decimal, e - number of digits in exponent; N/S (optional) indicates engineering/scientific formating of the E type.

Examples:

f = cast_FORTRAN_format("I10")
  FORTRAN_format("Iw", 'I', nothing, 10, 0, 0, 0)

f = cast_FORTRAN_format("I10.12")
  FORTRAN_format("Iw.m", 'I', nothing, 10, 12, 0, 0)

f = cast_FORTRAN_format("E10.5E3")
  FORTRAN_format("Ew.dEe", 'E', nothing, 10, 0, 5, 3)

f.Type, f.TypeChar, f.EngSci, f.width, f.nmin, f.ndec, f.nexp
  ("Ew.dEe", 'E', nothing, 10, 0, 5, 3)
CamiFITS.cast_FORTRAN_datatypeMethod
cast_FORTRAN_datatype(format::String)

Decompose the format specifier format into its fields and cast this into the FORTRAN_format object. Allowed format specifiers are of the types: Aw, Iw, Fw.d, Ew.d, Dw.d, where: w - width, d - number of digits to right of decimal point.

Examples:

f = cast_FORTRAN_datatype("I10")
  FORTRAN_format("Iw", 'I', nothing, 10, 0, 0, 0)

f = cast_FORTRAN_datatypet("F10.4")
  FORTRAN_format("Fw.d", 'F', nothing, 10, 0, 4, 0)

f = cast_FORTRAN_datatype("E10.5")
  FORTRAN_format("Ew.d", 'E', nothing, 10, 0, 5, 0)

f.Type, f.TypeChar, f.EngSci, f.width, f.nmin, f.ndec, f.nexp
  ("Ew.d", 'E', nothing, 10, 0, 5, 0)

Plotting

CamiFITS.step125Method
step125(x)

Step used for deviding the number x in steps according to 1-2-5 scheme

Examples:

step125.([5,10,21.3,50,100.1])
5-element Vector{Int64}:
  1
  2
  5
 10
 20
CamiFITS.select125Method
select125(x)

Select elements of the collection x by index according to 1-2-5 scheme

Examples:

x = [1,2,4,6,8,10,13,16,18,20,40,60,80,100]
select125(x)
 [2, 6, 10, 16, 20, 60, 100]

x = string.(x)
select125(x)
 ["2", "6", "10", "16", "20", "60", "100"]

x = 1:100
select125(x)
 [20, 40, 60, 80, 100]
CamiFITS.stepsMethod
steps(x)

Heatmap range transformation for steplength specification vector x

Examples:

x = [4,2,6]
steps(x)
 [0, 4, 6, 12]
CamiFITS.stepcentersMethod
stepcenters(x)

Stepcenter positions for steplength specification vector x

Examples:

x = [4,2,6]
stepcenters(x)
 [2.0, 5.0, 9.0]
CamiFITS.stepedgesMethod
stepedges(x)

Stepedges for steplength specification vector x

Examples:

x = [4,2,6]
stepedges(x)
 [0, 4, 6, 12]
CamiFITS.edgesFunction
edges(px [, Δx[, x0]])

Heatmap range transformation from pixel coordinates to physical coordinates, with pixelsize Δx and offset x0, both in physical units.

Examples:

px = 1:5
Δx = 2.5
x0 = 2.5
edges(px)
 [0.5, 1.5, 2.5, 3.5, 4.5]

edges(px, Δx)
 [1.25, 3.75, 6.25, 8.75, 11.25]

edges(px, Δx, x0)
 [-1.25, 1.25, 3.75, 6.25, 8.75]

Index