write_fits_cut4

This IDL facility writes out a cut sky HEALPix map into a FITS file according to the HEALPix convention. The format used for the FITS file follows the one used for Boomerang98 and is adapted from COBE/DMR. This routine can be used to store polarized maps, where the information relative to the Stokes parameters I, Q and U are placed in extension 0, 1 and 2 respectively by successive invocation of the routine.

Location in HEALPix directory tree: src/idl/fits/write_fits_cut4.pro 


FORMAT

IDL> WRITE_FITS_CUT4 , File, Pixel, Signal[, N_Obs, Serror, COORDSYS=, EXTENSION=, HDR=, /NESTED, NSIDE=, ORDERING=, /POLARISATION, /RING, UNITS=, XHDR=, HELP=]


QUALIFIERS

File
name of a FITS file in which the map is to be written

Pixel
(LONG or LONG64 vector),
index of observed (or valid) pixels

Signal
(FLOAT or DOUBLE vector, same size as Pixel),
value of signal in each observed pixel

N_Obs
(LONG or INT or LONG64 vector, Optional, same size as Pixel),
number of observation per pixel.
If absent, the field N_OBS will take a value of 1 in the output file. If set to a scalar constant, N_OBS will take this value in the output file

Serror
(FLOAT or DOUBLE vector, Optional, same size as Pixel)
rms of signal in pixel, for white noise, this is $\propto 1/\sqrt{{\rm n\_obs}}$
If absent, the field SERROR will take a value of 0.0 in the output file. If set to a scalar constant, SERROR will take this value in the output file


KEYWORDS

COORDSYS=
(optional),
if set to either 'C', 'E' or 'G', specifies that the Healpix coordinate system is respectively Celestial=equatorial, Ecliptic or Galactic. (The relevant keyword is then added/updated in the extension header, but the map is NOT rotated)

EXTENSION=
(optional),
(0 based) extension number in which to write data. default:0. If set to 0 (or not set) a new file is written from scratch. If set to a value larger than 1, the corresponding extension is added or updated, as long as all previous extensions already exist. All extensions of the same file should use the same ORDERING, NSIDE and COORDSYS.

HDR=
(optional),
String array containing the information to be put in the primary header.

/NESTED
(optional) if set, specifies that the map is in the NESTED ordering scheme
see also:Ordering and Ring

NSIDE=
(optional),
scalar integer, HEALPix resolution parameter of the data set. The resolution parameter should be made available to the FITS file, either thru this qualifier, or via the header (see XHDR).

ORDERING=
(optional),
if set to either 'ring' or 'nested' (case un-sensitive), specifies that the map is respectively in RING or NESTED ordering scheme
see also:Nested and Ring
The ordering information should be made available to the FITS file, either thru a combination of Ordering/Ring/Nested, or via the header (see XHDR).

/POLARISATION
specifies that file will contain the I, Q and U polarisation Stokes parameter in extensions 0, 1 and 2 respectively, and sets the FITS header keywords accordingly

/RING
if set, specifies that the map is in the RING ordering scheme
see also:Ordering and Nested

UNITS=
(optional),
string describing the physical units of the data set (only applies to Signal and Serror)

XHDR=
(optional),
String array containing the information to be put in the extension header.

HELP=
(optional),
if set, an extensive help is displayed, and no file is written


DESCRIPTION

For more information on the FITS file format supported in HEALPix, including the one implemented in write_fits_cut4 , see https://healpix.sourceforge.io/data/examples/healpix_fits_specs.pdf.


RELATED ROUTINES

This section lists the routines related to write_fits_cut4

idl
version 6.4 or more is necessary to run write_fits_cut4
read_fits_cut4
This HEALPix IDL facility can be used to read in maps written by write_fits_cut4 .
write_fits_cut4, write_fits_partial, write_fits_map
write_tqu, write_fits_sb
HEALPix IDL routines to write cut-sky and partial maps, full-sky maps, polarized full-sky maps and arbitrary data sets into FITS files
sxaddpar
This IDL routine (included in HEALPix package) can be used to update or add FITS keywords to the header in HDR and XHDR


EXAMPLE # 1:

write_fits_cut4 , 'map_cut.fits', pixel, temperature, /ring,nside=32, /pol  
writes in 'map_cut.fits' a FITS file containing the temperature measured in a set of HEALPix pixel.

EXAMPLE # 2:

write_fits_cut4 , 'tqu_cut.fits', pixel, temperature, n_t, s_t, $
/ring, nside=32, /pol
write_fits_cut4 , 'tqu_cut.fits', pixel, qstokes, n_q, s_q, $
/ring, nside=32, /pol, ext=1
write_fits_cut4 , 'tqu_cut.fits', pixel, ustokes, n_u, s_u, $
/ring, nside=32, /pol, ext=2
writes in 'tqu_cut.fits' a FITS file with three extensions, each of them containing information on the observed pixel, the measured signal, the number of observations and noise per pixel, for the three Stokes parameters I, Q and U respectively. The HEALPix ring ordered scheme and the resolution Nside=32 is assumed.

Version 3.82, 2022-07-28