These routines read or write data values in the primary data array (i.e., the first HDU in the FITS file) or an IMAGE extension. Automatic data type conversion is performed for if the data type of the FITS array (as defined by the BITPIX keyword) differs from the data type of the array in the calling routine. The data values are automatically scaled by the BSCALE and BZERO header values as they are being written or read from the FITS array. Unlike the basic routines described in the previous chapter, most of these routines specifically support the FITS random groups format. See Appendix B for the definition of the parameters used in these routines.
The more primitive reading and writing routines (i. e., ffppr_, ffppn_, ffppn, ffgpv_, or ffgpf_) simply treat the primary array as a long 1-dimensional array of pixels, ignoring the intrinsic dimensionality of the array. When dealing with a 2D image, for example, the application program must calculate the pixel offset in the 1-D array that corresponds to any particular X, Y coordinate in the image. C programmers should note that the ordering of arrays in FITS files, and hence in all the CFITSIO calls, is more similar to the dimensionality of arrays in Fortran rather than C. For instance if a FITS image has NAXIS1 = 100 and NAXIS2 = 50, then a 2-D array just large enough to hold the image should be declared as array[50][100] and not as array[100][50].
For convenience, higher-level routines are also provided to specifically deal with 2D images (ffp2d_ and ffg2d_) and 3D data cubes (ffp3d_ and ffg3d_). The dimensionality of the FITS image is passed by the naxis1, naxis2, and naxis3 parameters and the declared dimensions of the program array are passed in the dim1 and dim2 parameters. Note that the dimensions of the program array may be larger than the dimensions of the FITS array. For example if a FITS image with NAXIS1 = NAXIS2 = 400 is read into a program array which is dimensioned as 512 x 512 pixels, then the image will just fill the lower left corner of the array with pixels in the range 1 - 400 in the X an Y directions. This has the effect of taking a contiguous set of pixel value in the FITS array and writing them to a non-contiguous array in program memory (i.e., there are now some blank pixels around the edge of the image in the program array).
The most general set of routines (ffpss_, ffgsv_, and ffgsf_) may be used to transfer a rectangular subset of the pixels in a FITS N-dimensional image to or from an array which has been declared in the calling program. The fpixel and lpixel parameters are integer arrays which specify the starting and ending pixel coordinate in each dimension (starting with 1, not 0) of the FITS image that is to be read or written. It is important to note that these are the starting and ending pixels in the FITS image, not in the declared array in the program. The array parameter in these routines is treated simply as a large one-dimensional array of the appropriate data type containing the pixel values; The pixel values in the FITS array are read/written from/to this program array in strict sequence without any gaps; it is up to the calling routine to correctly interpret the dimensionality of this array. The two FITS reading routines (ffgsv_ and ffgsf_ ) also have an `inc' parameter which defines the data sampling interval in each dimension of the FITS array. For example, if inc[0]=2 and inc[1]=3 when reading a 2-dimensional FITS image, then only every other pixel in the first dimension and every 3rd pixel in the second dimension will be returned to the 'array' parameter.
Two types of routines are provided to read the data array which differ in the way undefined pixels are handled. The first type of routines (e.g., ffgpv_) simply return an array of data elements in which undefined pixels are set equal to a value specified by the user in the `nulval' parameter. An additional feature of these routines is that if the user sets nulval = 0, then no checks for undefined pixels will be performed, thus reducing the amount of CPU processing. The second type of routines (e.g., ffgpf_) returns the data element array and, in addition, a char array that indicates whether the value of the corresponding data pixel is undefined (= 1) or defined (= 0). The latter type of routines may be more convenient to use in some circumstances, however, it requires an additional array of logical values which can be unwieldy when working with large data arrays.
int fits_write_img / ffppr (fitsfile *fptr, int datatype, LONGLONG firstelem, LONGLONG nelements, DTYPE *array, int *status); int fits_write_img_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffppr[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelements, DTYPE *array, > int *status); int fits_write_imgnull / ffppn (fitsfile *fptr, int datatype, LONGLONG firstelem, LONGLONG nelements, DTYPE *array, DTYPE *nulval, > int *status); int fits_write_imgnull_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffppn[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelements, DTYPE *array, DTYPE nulval, > int *status);
int fits_write_img_null / ffppru (fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelements, > int *status)
int fits_write_grppar_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffpgp[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, long firstelem, long nelements, > DTYPE *array, int *status)
int fits_write_2d_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffp2d[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, LONGLONG dim1, LONGLONG naxis1, LONGLONG naxis2, DTYPE *array, > int *status) int fits_write_3d_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffp3d[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, LONGLONG dim1, LONGLONG dim2, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3, DTYPE *array, > int *status)
int fits_write_subset_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffpss[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, long naxis, long *naxes, long *fpixel, long *lpixel, DTYPE *array, > int *status)
int fits_read_img / ffgpv (fitsfile *fptr, int datatype, long firstelem, long nelements, DTYPE *nulval, > DTYPE *array, int *anynul, int *status) int fits_read_img_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffgpv[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, long firstelem, long nelements, DTYPE nulval, > DTYPE *array, int *anynul, int *status) int fits_read_imgnull / ffgpf (fitsfile *fptr, int datatype, long firstelem, long nelements, > DTYPE *array, char *nullarray, int *anynul, int *status) int fits_read_imgnull_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffgpf[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, long firstelem, long nelements, > DTYPE *array, char *nullarray, int *anynul, int *status)
int fits_read_grppar_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffggp[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, long firstelem, long nelements, > DTYPE *array, int *status)
int fits_read_2d_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffg2d[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, DTYPE nulval, LONGLONG dim1, LONGLONG naxis1, LONGLONG naxis2, > DTYPE *array, int *anynul, int *status) int fits_read_3d_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffg3d[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, long group, DTYPE nulval, LONGLONG dim1, LONGLONG dim2, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3, > DTYPE *array, int *anynul, int *status)
int fits_read_subset_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffgsv[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, int group, int naxis, long *naxes, long *fpixel, long *lpixel, long *inc, DTYPE nulval, > DTYPE *array, int *anynul, int *status) int fits_read_subsetnull_[byt, sht, usht, int, uint, lng, ulng, lnglng, ulnglng, flt, dbl] / ffgsf[b,i,ui,k,uk,j,uj,jj,ujj,e,d] (fitsfile *fptr, int group, int naxis, long *naxes, long *fpixel, long *lpixel, long *inc, > DTYPE *array, char *nullarray, int *anynul, int *status)