1/* The FITSIO software was written by William Pence at the High Energy */
2/* Astrophysic Science Archive Research Center (HEASARC) at the NASA */
3/* Goddard Space Flight Center. */
4/*
5
6Copyright (Unpublished--all rights reserved under the copyright laws of
7the United States), U.S. Government as represented by the Administrator
8of the National Aeronautics and Space Administration. No copyright is
9claimed in the United States under Title 17, U.S. Code.
10
11Permission to freely use, copy, modify, and distribute this software
12and its documentation without fee is hereby granted, provided that this
13copyright notice and disclaimer of warranty appears in all copies.
14
15DISCLAIMER:
16
17THE SOFTWARE IS PROVIDED 'AS IS' WITHOUT ANY WARRANTY OF ANY KIND,
18EITHER EXPRESSED, IMPLIED, OR STATUTORY, INCLUDING, BUT NOT LIMITED TO,
19ANY WARRANTY THAT THE SOFTWARE WILL CONFORM TO SPECIFICATIONS, ANY
20IMPLIED WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR
21PURPOSE, AND FREEDOM FROM INFRINGEMENT, AND ANY WARRANTY THAT THE
22DOCUMENTATION WILL CONFORM TO THE SOFTWARE, OR ANY WARRANTY THAT THE
23SOFTWARE WILL BE ERROR FREE. IN NO EVENT SHALL NASA BE LIABLE FOR ANY
24DAMAGES, INCLUDING, BUT NOT LIMITED TO, DIRECT, INDIRECT, SPECIAL OR
25CONSEQUENTIAL DAMAGES, ARISING OUT OF, RESULTING FROM, OR IN ANY WAY
26CONNECTED WITH THIS SOFTWARE, WHETHER OR NOT BASED UPON WARRANTY,
27CONTRACT, TORT , OR OTHERWISE, WHETHER OR NOT INJURY WAS SUSTAINED BY
28PERSONS OR PROPERTY OR OTHERWISE, AND WHETHER OR NOT LOSS WAS SUSTAINED
29FROM, OR AROSE OUT OF THE RESULTS OF, OR USE OF, THE SOFTWARE OR
30SERVICES PROVIDED HEREUNDER."
31
32*/
33
34#ifndef _FITSIO_H
35#define _FITSIO_H
36
37#define CFITSIO_VERSION 3.36
38#define CFITSIO_MINOR 36
39#define CFITSIO_MAJOR 3
40#define CFITSIO_SONAME 2
41
42/* the SONAME is incremented in a new release if the binary shared */
43/* library (on linux and Mac systems) is not backward compatible */
44/* with the previous release of CFITSIO */
45
46
47/* CFITS_API is defined below for use on Windows systems. */
48/* It is used to identify the public functions which should be exported. */
49/* This has no effect on non-windows platforms where "WIN32" is not defined */
50
51#if defined (WIN32)
52 #if defined(cfitsio_EXPORTS)
53 #define CFITS_API __declspec(dllexport)
54 #else
55 #define CFITS_API /* __declspec(dllimport) */
56 #endif /* CFITS_API */
57#else /* defined (WIN32) */
58 #define CFITS_API
59#endif
60
61#include <stdio.h>
62
63/* the following was provided by Michael Greason (GSFC) to fix a */
64/* C/Fortran compatibility problem on an SGI Altix system running */
65/* SGI ProPack 4 [this is a Novell SuSE Enterprise 9 derivative] */
66/* and using the Intel C++ and Fortran compilers (version 9.1) */
67#if defined(__INTEL_COMPILER) && defined(__itanium__)
68# define mipsFortran 1
69# define _MIPS_SZLONG 64
70#endif
71
72#if defined(linux) || defined(__APPLE__) || defined(__sgi)
73# include <sys/types.h> /* apparently needed on debian linux systems */
74#endif /* to define off_t */
75
76#include <stdlib.h> /* apparently needed to define size_t with gcc 2.8.1 */
77#include <limits.h> /* needed for LLONG_MAX and INT64_MAX definitions */
78
79/* Define the datatype for variables which store file offset values. */
80/* The newer 'off_t' datatype should be used for this purpose, but some */
81/* older compilers do not recognize this type, in which case we use 'long' */
82/* instead. Note that _OFF_T is defined (or not) in stdio.h depending */
83/* on whether _LARGEFILE_SOURCE is defined in sys/feature_tests.h */
84/* (at least on Solaris platforms using cc) */
85
86/* Debian systems require: "(defined(linux) && defined(__off_t_defined))" */
87/* the mingw-w64 compiler requires: "(defined(__MINGW32__) && defined(_OFF_T_DEFINED))" */
88#if defined(_OFF_T) \
89 || (defined(linux) && defined(__off_t_defined)) \
90 || (defined(__MINGW32__) && defined(_OFF_T_DEFINED)) \
91 || defined(_MIPS_SZLONG) || defined(__APPLE__) || defined(_AIX)
92# define OFF_T off_t
93#elif defined(_MSC_VER) && (_MSC_VER>= 1400)
94# define OFF_T long long
95#else
96# define OFF_T long
97#endif
98
99/* this block determines if the the string function name is
100 strtol or strtoll, and whether to use %ld or %lld in printf statements */
101
102/*
103 The following 2 cases for that Athon64 were removed on 4 Jan 2006;
104 they appear to be incorrect now that LONGLONG is always typedef'ed
105 to 'long long'
106 || defined(__ia64__) \
107 || defined(__x86_64__) \
108*/
109#if (defined(__alpha) && ( defined(__unix__) || defined(__NetBSD__) )) \
110 || defined(__sparcv9) || (defined(__sparc__) && defined(__arch64__)) \
111 || defined(__powerpc64__) || defined(__64BIT__) \
112 || (defined(_MIPS_SZLONG) && _MIPS_SZLONG == 64) \
113 || defined( _MSC_VER)|| defined(__BORLANDC__)
114
115# define USE_LL_SUFFIX 0
116#else
117# define USE_LL_SUFFIX 1
118#endif
119
120/*
121 Determine what 8-byte integer data type is available.
122 'long long' is now supported by most compilers, but
123 older MS Visual C++ compilers before V7.0 use '__int64' instead.
124*/
125
126#ifndef LONGLONG_TYPE /* this may have been previously defined */
127#if defined(_MSC_VER) /* Microsoft Visual C++ */
128
129#if (_MSC_VER < 1300) /* versions earlier than V7.0 do not have 'long long' */
130 typedef __int64 LONGLONG;
131#else /* newer versions do support 'long long' */
132 typedef long long LONGLONG;
133#endif
134
135#elif defined( __BORLANDC__) /* for the Borland 5.5 compiler, in particular */
136 typedef __int64 LONGLONG;
137#else
138 typedef long long LONGLONG;
139#endif
140
141#define LONGLONG_TYPE
142#endif
143
144#ifndef LONGLONG_MAX
145
146#ifdef LLONG_MAX
147/* Linux and Solaris definition */
148#define LONGLONG_MAX LLONG_MAX
149#define LONGLONG_MIN LLONG_MIN
150
151#elif defined(LONG_LONG_MAX)
152#define LONGLONG_MAX LONG_LONG_MAX
153#define LONGLONG_MIN LONG_LONG_MIN
154
155#elif defined(__LONG_LONG_MAX__)
156/* Mac OS X & CYGWIN defintion */
157#define LONGLONG_MAX __LONG_LONG_MAX__
158#define LONGLONG_MIN (-LONGLONG_MAX -1LL)
159
160#elif defined(INT64_MAX)
161/* windows definition */
162#define LONGLONG_MAX INT64_MAX
163#define LONGLONG_MIN INT64_MIN
164
165#elif defined(_I64_MAX)
166/* windows definition */
167#define LONGLONG_MAX _I64_MAX
168#define LONGLONG_MIN _I64_MIN
169
170#elif (defined(__alpha) && ( defined(__unix__) || defined(__NetBSD__) )) \
171 || defined(__sparcv9) \
172 || defined(__ia64__) \
173 || defined(__x86_64__) \
174 || defined(_SX) \
175 || defined(__powerpc64__) || defined(__64BIT__) \
176 || (defined(_MIPS_SZLONG) && _MIPS_SZLONG == 64)
177/* sizeof(long) = 64 */
178#define LONGLONG_MAX 9223372036854775807L /* max 64-bit integer */
179#define LONGLONG_MIN (-LONGLONG_MAX -1L) /* min 64-bit integer */
180
181#else
182/* define a default value, even if it is never used */
183#define LONGLONG_MAX 9223372036854775807LL /* max 64-bit integer */
184#define LONGLONG_MIN (-LONGLONG_MAX -1LL) /* min 64-bit integer */
185
186#endif
187#endif /* end of ndef LONGLONG_MAX section */
188
189
190/* ================================================================= */
191
192
193/* The following exclusion if __CINT__ is defined is needed for ROOT */
194#ifndef __CINT__
195#include "longnam.h"
196#endif
197
198#define NIOBUF 40 /* number of IO buffers to create (default = 40) */
199 /* !! Significantly increasing NIOBUF may degrade performance !! */
200
201#define IOBUFLEN 2880 /* size in bytes of each IO buffer (DONT CHANGE!) */
202
203/* global variables */
204
205#define FLEN_FILENAME 1025 /* max length of a filename */
206#define FLEN_KEYWORD 72 /* max length of a keyword (HIERARCH convention) */
207#define FLEN_CARD 81 /* length of a FITS header card */
208#define FLEN_VALUE 71 /* max length of a keyword value string */
209#define FLEN_COMMENT 73 /* max length of a keyword comment string */
210#define FLEN_ERRMSG 81 /* max length of a FITSIO error message */
211#define FLEN_STATUS 31 /* max length of a FITSIO status text string */
212
213#define TBIT 1 /* codes for FITS table data types */
214#define TBYTE 11
215#define TSBYTE 12
216#define TLOGICAL 14
217#define TSTRING 16
218#define TUSHORT 20
219#define TSHORT 21
220#define TUINT 30
221#define TINT 31
222#define TULONG 40
223#define TLONG 41
224#define TINT32BIT 41 /* used when returning datatype of a column */
225#define TFLOAT 42
226#define TLONGLONG 81
227#define TDOUBLE 82
228#define TCOMPLEX 83
229#define TDBLCOMPLEX 163
230
231#define TYP_STRUC_KEY 10
232#define TYP_CMPRS_KEY 20
233#define TYP_SCAL_KEY 30
234#define TYP_NULL_KEY 40
235#define TYP_DIM_KEY 50
236#define TYP_RANG_KEY 60
237#define TYP_UNIT_KEY 70
238#define TYP_DISP_KEY 80
239#define TYP_HDUID_KEY 90
240#define TYP_CKSUM_KEY 100
241#define TYP_WCS_KEY 110
242#define TYP_REFSYS_KEY 120
243#define TYP_COMM_KEY 130
244#define TYP_CONT_KEY 140
245#define TYP_USER_KEY 150
246
247
248#define INT32BIT int /* 32-bit integer datatype. Currently this */
249 /* datatype is an 'int' on all useful platforms */
250 /* however, it is possible that that are cases */
251 /* where 'int' is a 2-byte integer, in which case */
252 /* INT32BIT would need to be defined as 'long'. */
253
254#define BYTE_IMG 8 /* BITPIX code values for FITS image types */
255#define SHORT_IMG 16
256#define LONG_IMG 32
257#define LONGLONG_IMG 64
258#define FLOAT_IMG -32
259#define DOUBLE_IMG -64
260 /* The following 2 codes are not true FITS */
261 /* datatypes; these codes are only used internally */
262 /* within cfitsio to make it easier for users */
263 /* to deal with unsigned integers. */
264#define SBYTE_IMG 10
265#define USHORT_IMG 20
266#define ULONG_IMG 40
267
268#define IMAGE_HDU 0 /* Primary Array or IMAGE HDU */
269#define ASCII_TBL 1 /* ASCII table HDU */
270#define BINARY_TBL 2 /* Binary table HDU */
271#define ANY_HDU -1 /* matches any HDU type */
272
273#define READONLY 0 /* options when opening a file */
274#define READWRITE 1
275
276/* adopt a hopefully obscure number to use as a null value flag */
277/* could be problems if the FITS files contain data with these values */
278#define FLOATNULLVALUE -9.11912E-36F
279#define DOUBLENULLVALUE -9.1191291391491E-36
280
281/* compression algorithm codes */
282#define NO_DITHER -1
283#define SUBTRACTIVE_DITHER_1 1
284#define SUBTRACTIVE_DITHER_2 2
285#define MAX_COMPRESS_DIM 6
286#define RICE_1 11
287#define GZIP_1 21
288#define GZIP_2 22
289#define PLIO_1 31
290#define HCOMPRESS_1 41
291#define BZIP2_1 51 /* not publicly supported; only for test purposes */
292#define NOCOMPRESS -1
293
294#ifndef TRUE
295#define TRUE 1
296#endif
297
298#ifndef FALSE
299#define FALSE 0
300#endif
301
302#define CASESEN 1 /* do case-sensitive string match */
303#define CASEINSEN 0 /* do case-insensitive string match */
304
305#define GT_ID_ALL_URI 0 /* hierarchical grouping parameters */
306#define GT_ID_REF 1
307#define GT_ID_POS 2
308#define GT_ID_ALL 3
309#define GT_ID_REF_URI 11
310#define GT_ID_POS_URI 12
311
312#define OPT_RM_GPT 0
313#define OPT_RM_ENTRY 1
314#define OPT_RM_MBR 2
315#define OPT_RM_ALL 3
316
317#define OPT_GCP_GPT 0
318#define OPT_GCP_MBR 1
319#define OPT_GCP_ALL 2
320
321#define OPT_MCP_ADD 0
322#define OPT_MCP_NADD 1
323#define OPT_MCP_REPL 2
324#define OPT_MCP_MOV 3
325
326#define OPT_MRG_COPY 0
327#define OPT_MRG_MOV 1
328
329#define OPT_CMT_MBR 1
330#define OPT_CMT_MBR_DEL 11
331
332typedef struct /* structure used to store table column information */
333{
334 char ttype[70]; /* column name = FITS TTYPEn keyword; */
335 LONGLONG tbcol; /* offset in row to first byte of each column */
336 int tdatatype; /* datatype code of each column */
337 LONGLONG trepeat; /* repeat count of column; number of elements */
338 double tscale; /* FITS TSCALn linear scaling factor */
339 double tzero; /* FITS TZEROn linear scaling zero point */
340 LONGLONG tnull; /* FITS null value for int image or binary table cols */
341 char strnull[20]; /* FITS null value string for ASCII table columns */
342 char tform[10]; /* FITS tform keyword value */
343 long twidth; /* width of each ASCII table column */
344}tcolumn;
345
346#define VALIDSTRUC 555 /* magic value used to identify if structure is valid */
347
348typedef struct /* structure used to store basic FITS file information */
349{
350 int filehandle; /* handle returned by the file open function */
351 int driver; /* defines which set of I/O drivers should be used */
352 int open_count; /* number of opened 'fitsfiles' using this structure */
353 char *filename; /* file name */
354 int validcode; /* magic value used to verify that structure is valid */
355 int only_one; /* flag meaning only copy the specified extension */
356 LONGLONG filesize; /* current size of the physical disk file in bytes */
357 LONGLONG logfilesize; /* logical size of file, including unflushed buffers */
358 int lasthdu; /* is this the last HDU in the file? 0 = no, else yes */
359 LONGLONG bytepos; /* current logical I/O pointer position in file */
360 LONGLONG io_pos; /* current I/O pointer position in the physical file */
361 int curbuf; /* number of I/O buffer currently in use */
362 int curhdu; /* current HDU number; 0 = primary array */
363 int hdutype; /* 0 = primary array, 1 = ASCII table, 2 = binary table */
364 int writemode; /* 0 = readonly, 1 = readwrite */
365 int maxhdu; /* highest numbered HDU known to exist in the file */
366 int MAXHDU; /* dynamically allocated dimension of headstart array */
367 LONGLONG *headstart; /* byte offset in file to start of each HDU */
368 LONGLONG headend; /* byte offest in file to end of the current HDU header */
369 LONGLONG ENDpos; /* byte offest to where the END keyword was last written */
370 LONGLONG nextkey; /* byte offset in file to beginning of next keyword */
371 LONGLONG datastart; /* byte offset in file to start of the current data unit */
372 int imgdim; /* dimension of image; cached for fast access */
373 LONGLONG imgnaxis[99]; /* length of each axis; cached for fast access */
374 int tfield; /* number of fields in the table (primary array has 2 */
375 int startcol; /* used by ffgcnn to record starting column number */
376 LONGLONG origrows; /* original number of rows (value of NAXIS2 keyword) */
377 LONGLONG numrows; /* number of rows in the table (dynamically updated) */
378 LONGLONG rowlength; /* length of a table row or image size (bytes) */
379 tcolumn *tableptr; /* pointer to the table structure */
380 LONGLONG heapstart; /* heap start byte relative to start of data unit */
381 LONGLONG heapsize; /* size of the heap, in bytes */
382
383 /* the following elements are related to compressed images */
384
385 /* these record the 'requested' options to be used when the image is compressed */
386 int request_compress_type; /* requested image compression algorithm */
387 long request_tilesize[MAX_COMPRESS_DIM]; /* requested tiling size */
388 float request_quantize_level; /* requested quantize level */
389 int request_quantize_method ; /* requested quantizing method */
390 int request_dither_seed; /* starting offset into the array of random dithering */
391 int request_lossy_int_compress; /* lossy compress integer image as if float image? */
392 int request_huge_hdu; /* use '1Q' rather then '1P' variable length arrays */
393 float request_hcomp_scale; /* requested HCOMPRESS scale factor */
394 int request_hcomp_smooth; /* requested HCOMPRESS smooth parameter */
395
396 /* these record the actual options that were used when the image was compressed */
397 int compress_type; /* type of compression algorithm */
398 long tilesize[MAX_COMPRESS_DIM]; /* size of compression tiles */
399 float quantize_level; /* floating point quantization level */
400 int quantize_method; /* floating point pixel quantization algorithm */
401 int dither_seed; /* starting offset into the array of random dithering */
402
403 /* other compression parameters */
404 int compressimg; /* 1 if HDU contains a compressed image, else 0 */
405 char zcmptype[12]; /* compression type string */
406 int zbitpix; /* FITS data type of image (BITPIX) */
407 int zndim; /* dimension of image */
408 long znaxis[MAX_COMPRESS_DIM]; /* length of each axis */
409 long maxtilelen; /* max number of pixels in each image tile */
410 long maxelem; /* maximum byte length of tile compressed arrays */
411
412 int cn_compressed; /* column number for COMPRESSED_DATA column */
413 int cn_uncompressed; /* column number for UNCOMPRESSED_DATA column */
414 int cn_gzip_data; /* column number for GZIP2 lossless compressed data */
415 int cn_zscale; /* column number for ZSCALE column */
416 int cn_zzero; /* column number for ZZERO column */
417 int cn_zblank; /* column number for the ZBLANK column */
418
419 double zscale; /* scaling value, if same for all tiles */
420 double zzero; /* zero pt, if same for all tiles */
421 double cn_bscale; /* value of the BSCALE keyword in header */
422 double cn_bzero; /* value of the BZERO keyword (may be reset) */
423 double cn_actual_bzero; /* actual value of the BZERO keyword */
424 int zblank; /* value for null pixels, if not a column */
425
426 int rice_blocksize; /* first compression parameter: Rice pixels/block */
427 int rice_bytepix; /* 2nd compression parameter: Rice bytes/pixel */
428 float hcomp_scale; /* 1st hcompress compression parameter */
429 int hcomp_smooth; /* 2nd hcompress compression parameter */
430
431 int *tilerow; /* row number of the array of uncompressed tiledata */
432 long *tiledatasize; /* length of the array of tile data in bytes */
433 int *tiletype; /* datatype of the array of tile (TINT, TSHORT, etc) */
434 void **tiledata; /* array of uncompressed tile of data, for row *tilerow */
435 void **tilenullarray; /* array of optional array of null value flags */
436 int *tileanynull; /* anynulls in the array of tile? */
437
438 char *iobuffer; /* pointer to FITS file I/O buffers */
439 long bufrecnum[NIOBUF]; /* file record number of each of the buffers */
440 int dirty[NIOBUF]; /* has the corresponding buffer been modified? */
441 int ageindex[NIOBUF]; /* relative age of each buffer */
442} FITSfile;
443
444typedef struct /* structure used to store basic HDU information */
445{
446 int HDUposition; /* HDU position in file; 0 = first HDU */
447 FITSfile *Fptr; /* pointer to FITS file structure */
448}fitsfile;
449
450typedef struct /* structure for the iterator function column information */
451{
452 /* elements required as input to fits_iterate_data: */
453
454 fitsfile *fptr; /* pointer to the HDU containing the column */
455 int colnum; /* column number in the table (use name if < 1) */
456 char colname[70]; /* name (= TTYPEn value) of the column (optional) */
457 int datatype; /* output datatype (converted if necessary */
458 int iotype; /* = InputCol, InputOutputCol, or OutputCol */
459
460 /* output elements that may be useful for the work function: */
461
462 void *array; /* pointer to the array (and the null value) */
463 long repeat; /* binary table vector repeat value */
464 long tlmin; /* legal minimum data value */
465 long tlmax; /* legal maximum data value */
466 char tunit[70]; /* physical unit string */
467 char tdisp[70]; /* suggested display format */
468
469} iteratorCol;
470
471#define InputCol 0 /* flag for input only iterator column */
472#define InputOutputCol 1 /* flag for input and output iterator column */
473#define OutputCol 2 /* flag for output only iterator column */
474
475/*=============================================================================
476*
477* The following wtbarr typedef is used in the fits_read_wcstab() routine,
478* which is intended for use with the WCSLIB library written by Mark
479* Calabretta, http://www.atnf.csiro.au/~mcalabre/index.html
480*
481* In order to maintain WCSLIB and CFITSIO as independent libraries it
482* was not permissible for any CFITSIO library code to include WCSLIB
483* header files, or vice versa. However, the CFITSIO function
484* fits_read_wcstab() accepts an array of structs defined by wcs.h within
485* WCSLIB. The problem then was to define this struct within fitsio.h
486* without including wcs.h, especially noting that wcs.h will often (but
487* not always) be included together with fitsio.h in an applications
488* program that uses fits_read_wcstab().
489*
490* Of the various possibilities, the solution adopted was for WCSLIB to
491* define "struct wtbarr" while fitsio.h defines "typedef wtbarr", a
492* untagged struct with identical members. This allows both wcs.h and
493* fitsio.h to define a wtbarr data type without conflict by virtue of
494* the fact that structure tags and typedef names share different
495* namespaces in C. Therefore, declarations within WCSLIB look like
496*
497* struct wtbarr *w;
498*
499* while within CFITSIO they are simply
500*
501* wtbarr *w;
502*
503* but as suggested by the commonality of the names, these are really the
504* same aggregate data type. However, in passing a (struct wtbarr *) to
505* fits_read_wcstab() a cast to (wtbarr *) is formally required.
506*===========================================================================*/
507
508#ifndef WCSLIB_GETWCSTAB
509#define WCSLIB_GETWCSTAB
510
511typedef struct {
512 int i; /* Image axis number. */
513 int m; /* Array axis number for index vectors. */
514 int kind; /* Array type, 'c' (coord) or 'i' (index). */
515 char extnam[72]; /* EXTNAME of binary table extension. */
516 int extver; /* EXTVER of binary table extension. */
517 int extlev; /* EXTLEV of binary table extension. */
518 char ttype[72]; /* TTYPEn of column containing the array. */
519 long row; /* Table row number. */
520 int ndim; /* Expected array dimensionality. */
521 int *dimlen; /* Where to write the array axis lengths. */
522 double **arrayp; /* Where to write the address of the array */
523 /* allocated to store the array. */
524} wtbarr;
525
526int CFITS_API fits_read_wcstab(fitsfile *fptr, int nwtb, wtbarr *wtb, int *status);
527
528#endif /* WCSLIB_GETWCSTAB */
529
530/* error status codes */
531
532#define CREATE_DISK_FILE -106 /* create disk file, without extended filename syntax */
533#define OPEN_DISK_FILE -105 /* open disk file, without extended filename syntax */
534#define SKIP_TABLE -104 /* move to 1st image when opening file */
535#define SKIP_IMAGE -103 /* move to 1st table when opening file */
536#define SKIP_NULL_PRIMARY -102 /* skip null primary array when opening file */
537#define USE_MEM_BUFF -101 /* use memory buffer when opening file */
538#define OVERFLOW_ERR -11 /* overflow during datatype conversion */
539#define PREPEND_PRIMARY -9 /* used in ffiimg to insert new primary array */
540#define SAME_FILE 101 /* input and output files are the same */
541#define TOO_MANY_FILES 103 /* tried to open too many FITS files */
542#define FILE_NOT_OPENED 104 /* could not open the named file */
543#define FILE_NOT_CREATED 105 /* could not create the named file */
544#define WRITE_ERROR 106 /* error writing to FITS file */
545#define END_OF_FILE 107 /* tried to move past end of file */
546#define READ_ERROR 108 /* error reading from FITS file */
547#define FILE_NOT_CLOSED 110 /* could not close the file */
548#define ARRAY_TOO_BIG 111 /* array dimensions exceed internal limit */
549#define READONLY_FILE 112 /* Cannot write to readonly file */
550#define MEMORY_ALLOCATION 113 /* Could not allocate memory */
551#define BAD_FILEPTR 114 /* invalid fitsfile pointer */
552#define NULL_INPUT_PTR 115 /* NULL input pointer to routine */
553#define SEEK_ERROR 116 /* error seeking position in file */
554
555#define BAD_URL_PREFIX 121 /* invalid URL prefix on file name */
556#define TOO_MANY_DRIVERS 122 /* tried to register too many IO drivers */
557#define DRIVER_INIT_FAILED 123 /* driver initialization failed */
558#define NO_MATCHING_DRIVER 124 /* matching driver is not registered */
559#define URL_PARSE_ERROR 125 /* failed to parse input file URL */
560#define RANGE_PARSE_ERROR 126 /* failed to parse input file URL */
561
562#define SHARED_ERRBASE (150)
563#define SHARED_BADARG (SHARED_ERRBASE + 1)
564#define SHARED_NULPTR (SHARED_ERRBASE + 2)
565#define SHARED_TABFULL (SHARED_ERRBASE + 3)
566#define SHARED_NOTINIT (SHARED_ERRBASE + 4)
567#define SHARED_IPCERR (SHARED_ERRBASE + 5)
568#define SHARED_NOMEM (SHARED_ERRBASE + 6)
569#define SHARED_AGAIN (SHARED_ERRBASE + 7)
570#define SHARED_NOFILE (SHARED_ERRBASE + 8)
571#define SHARED_NORESIZE (SHARED_ERRBASE + 9)
572
573#define HEADER_NOT_EMPTY 201 /* header already contains keywords */
574#define KEY_NO_EXIST 202 /* keyword not found in header */
575#define KEY_OUT_BOUNDS 203 /* keyword record number is out of bounds */
576#define VALUE_UNDEFINED 204 /* keyword value field is blank */
577#define NO_QUOTE 205 /* string is missing the closing quote */
578#define BAD_INDEX_KEY 206 /* illegal indexed keyword name */
579#define BAD_KEYCHAR 207 /* illegal character in keyword name or card */
580#define BAD_ORDER 208 /* required keywords out of order */
581#define NOT_POS_INT 209 /* keyword value is not a positive integer */
582#define NO_END 210 /* couldn't find END keyword */
583#define BAD_BITPIX 211 /* illegal BITPIX keyword value*/
584#define BAD_NAXIS 212 /* illegal NAXIS keyword value */
585#define BAD_NAXES 213 /* illegal NAXISn keyword value */
586#define BAD_PCOUNT 214 /* illegal PCOUNT keyword value */
587#define BAD_GCOUNT 215 /* illegal GCOUNT keyword value */
588#define BAD_TFIELDS 216 /* illegal TFIELDS keyword value */
589#define NEG_WIDTH 217 /* negative table row size */
590#define NEG_ROWS 218 /* negative number of rows in table */
591#define COL_NOT_FOUND 219 /* column with this name not found in table */
592#define BAD_SIMPLE 220 /* illegal value of SIMPLE keyword */
593#define NO_SIMPLE 221 /* Primary array doesn't start with SIMPLE */
594#define NO_BITPIX 222 /* Second keyword not BITPIX */
595#define NO_NAXIS 223 /* Third keyword not NAXIS */
596#define NO_NAXES 224 /* Couldn't find all the NAXISn keywords */
597#define NO_XTENSION 225 /* HDU doesn't start with XTENSION keyword */
598#define NOT_ATABLE 226 /* the CHDU is not an ASCII table extension */
599#define NOT_BTABLE 227 /* the CHDU is not a binary table extension */
600#define NO_PCOUNT 228 /* couldn't find PCOUNT keyword */
601#define NO_GCOUNT 229 /* couldn't find GCOUNT keyword */
602#define NO_TFIELDS 230 /* couldn't find TFIELDS keyword */
603#define NO_TBCOL 231 /* couldn't find TBCOLn keyword */
604#define NO_TFORM 232 /* couldn't find TFORMn keyword */
605#define NOT_IMAGE 233 /* the CHDU is not an IMAGE extension */
606#define BAD_TBCOL 234 /* TBCOLn keyword value < 0 or > rowlength */
607#define NOT_TABLE 235 /* the CHDU is not a table */
608#define COL_TOO_WIDE 236 /* column is too wide to fit in table */
609#define COL_NOT_UNIQUE 237 /* more than 1 column name matches template */
610#define BAD_ROW_WIDTH 241 /* sum of column widths not = NAXIS1 */
611#define UNKNOWN_EXT 251 /* unrecognizable FITS extension type */
612#define UNKNOWN_REC 252 /* unrecognizable FITS record */
613#define END_JUNK 253 /* END keyword is not blank */
614#define BAD_HEADER_FILL 254 /* Header fill area not blank */
615#define BAD_DATA_FILL 255 /* Data fill area not blank or zero */
616#define BAD_TFORM 261 /* illegal TFORM format code */
617#define BAD_TFORM_DTYPE 262 /* unrecognizable TFORM datatype code */
618#define BAD_TDIM 263 /* illegal TDIMn keyword value */
619#define BAD_HEAP_PTR 264 /* invalid BINTABLE heap address */
620
621#define BAD_HDU_NUM 301 /* HDU number < 1 or > MAXHDU */
622#define BAD_COL_NUM 302 /* column number < 1 or > tfields */
623#define NEG_FILE_POS 304 /* tried to move before beginning of file */
624#define NEG_BYTES 306 /* tried to read or write negative bytes */
625#define BAD_ROW_NUM 307 /* illegal starting row number in table */
626#define BAD_ELEM_NUM 308 /* illegal starting element number in vector */
627#define NOT_ASCII_COL 309 /* this is not an ASCII string column */
628#define NOT_LOGICAL_COL 310 /* this is not a logical datatype column */
629#define BAD_ATABLE_FORMAT 311 /* ASCII table column has wrong format */
630#define BAD_BTABLE_FORMAT 312 /* Binary table column has wrong format */
631#define NO_NULL 314 /* null value has not been defined */
632#define NOT_VARI_LEN 317 /* this is not a variable length column */
633#define BAD_DIMEN 320 /* illegal number of dimensions in array */
634#define BAD_PIX_NUM 321 /* first pixel number greater than last pixel */
635#define ZERO_SCALE 322 /* illegal BSCALE or TSCALn keyword = 0 */
636#define NEG_AXIS 323 /* illegal axis length < 1 */
637
638#define NOT_GROUP_TABLE 340
639#define HDU_ALREADY_MEMBER 341
640#define MEMBER_NOT_FOUND 342
641#define GROUP_NOT_FOUND 343
642#define BAD_GROUP_ID 344
643#define TOO_MANY_HDUS_TRACKED 345
644#define HDU_ALREADY_TRACKED 346
645#define BAD_OPTION 347
646#define IDENTICAL_POINTERS 348
647#define BAD_GROUP_ATTACH 349
648#define BAD_GROUP_DETACH 350
649
650#define BAD_I2C 401 /* bad int to formatted string conversion */
651#define BAD_F2C 402 /* bad float to formatted string conversion */
652#define BAD_INTKEY 403 /* can't interprete keyword value as integer */
653#define BAD_LOGICALKEY 404 /* can't interprete keyword value as logical */
654#define BAD_FLOATKEY 405 /* can't interprete keyword value as float */
655#define BAD_DOUBLEKEY 406 /* can't interprete keyword value as double */
656#define BAD_C2I 407 /* bad formatted string to int conversion */
657#define BAD_C2F 408 /* bad formatted string to float conversion */
658#define BAD_C2D 409 /* bad formatted string to double conversion */
659#define BAD_DATATYPE 410 /* bad keyword datatype code */
660#define BAD_DECIM 411 /* bad number of decimal places specified */
661#define NUM_OVERFLOW 412 /* overflow during datatype conversion */
662
663# define DATA_COMPRESSION_ERR 413 /* error in imcompress routines */
664# define DATA_DECOMPRESSION_ERR 414 /* error in imcompress routines */
665# define NO_COMPRESSED_TILE 415 /* compressed tile doesn't exist */
666
667#define BAD_DATE 420 /* error in date or time conversion */
668
669#define PARSE_SYNTAX_ERR 431 /* syntax error in parser expression */
670#define PARSE_BAD_TYPE 432 /* expression did not evaluate to desired type */
671#define PARSE_LRG_VECTOR 433 /* vector result too large to return in array */
672#define PARSE_NO_OUTPUT 434 /* data parser failed not sent an out column */
673#define PARSE_BAD_COL 435 /* bad data encounter while parsing column */
674#define PARSE_BAD_OUTPUT 436 /* Output file not of proper type */
675
676#define ANGLE_TOO_BIG 501 /* celestial angle too large for projection */
677#define BAD_WCS_VAL 502 /* bad celestial coordinate or pixel value */
678#define WCS_ERROR 503 /* error in celestial coordinate calculation */
679#define BAD_WCS_PROJ 504 /* unsupported type of celestial projection */
680#define NO_WCS_KEY 505 /* celestial coordinate keywords not found */
681#define APPROX_WCS_KEY 506 /* approximate WCS keywords were calculated */
682
683#define NO_CLOSE_ERROR 999 /* special value used internally to switch off */
684 /* the error message from ffclos and ffchdu */
685
686/*------- following error codes are used in the grparser.c file -----------*/
687#define NGP_ERRBASE (360) /* base chosen so not to interfere with CFITSIO */
688#define NGP_OK (0)
689#define NGP_NO_MEMORY (NGP_ERRBASE + 0) /* malloc failed */
690#define NGP_READ_ERR (NGP_ERRBASE + 1) /* read error from file */
691#define NGP_NUL_PTR (NGP_ERRBASE + 2) /* null pointer passed as argument */
692#define NGP_EMPTY_CURLINE (NGP_ERRBASE + 3) /* line read seems to be empty */
693#define NGP_UNREAD_QUEUE_FULL (NGP_ERRBASE + 4) /* cannot unread more then 1 line (or single line twice) */
694#define NGP_INC_NESTING (NGP_ERRBASE + 5) /* too deep include file nesting (inf. loop ?) */
695#define NGP_ERR_FOPEN (NGP_ERRBASE + 6) /* fopen() failed, cannot open file */
696#define NGP_EOF (NGP_ERRBASE + 7) /* end of file encountered */
697#define NGP_BAD_ARG (NGP_ERRBASE + 8) /* bad arguments passed */
698#define NGP_TOKEN_NOT_EXPECT (NGP_ERRBASE + 9) /* token not expected here */
699
700/* The following exclusion if __CINT__ is defined is needed for ROOT */
701#ifndef __CINT__
702/* the following 3 lines are needed to support C++ compilers */
703#ifdef __cplusplus
704extern "C" {
705#endif
706#endif
707
708int CFITS2Unit( fitsfile *fptr );
709CFITS_API fitsfile* CUnit2FITS(int unit);
710
711/*---------------- FITS file URL parsing routines -------------*/
712int CFITS_API fits_get_token (char **ptr, char *delimiter, char *token, int *isanumber);
713int CFITS_API fits_get_token2(char **ptr, char *delimiter, char **token, int *isanumber, int *status);
714char CFITS_API *fits_split_names(char *list);
715int CFITS_API ffiurl( char *url, char *urltype, char *infile,
716 char *outfile, char *extspec, char *rowfilter,
717 char *binspec, char *colspec, int *status);
718int CFITS_API ffifile (char *url, char *urltype, char *infile,
719 char *outfile, char *extspec, char *rowfilter,
720 char *binspec, char *colspec, char *pixfilter, int *status);
721int CFITS_API ffifile2 (char *url, char *urltype, char *infile,
722 char *outfile, char *extspec, char *rowfilter,
723 char *binspec, char *colspec, char *pixfilter, char *compspec, int *status);
724int CFITS_API ffrtnm(char *url, char *rootname, int *status);
725int CFITS_API ffexist(const char *infile, int *exists, int *status);
726int CFITS_API ffexts(char *extspec, int *extnum, char *extname, int *extvers,
727 int *hdutype, char *colname, char *rowexpress, int *status);
728int CFITS_API ffextn(char *url, int *extension_num, int *status);
729int CFITS_API ffurlt(fitsfile *fptr, char *urlType, int *status);
730int CFITS_API ffbins(char *binspec, int *imagetype, int *haxis,
731 char colname[4][FLEN_VALUE], double *minin,
732 double *maxin, double *binsizein,
733 char minname[4][FLEN_VALUE], char maxname[4][FLEN_VALUE],
734 char binname[4][FLEN_VALUE], double *weight, char *wtname,
735 int *recip, int *status);
736int CFITS_API ffbinr(char **binspec, char *colname, double *minin,
737 double *maxin, double *binsizein, char *minname,
738 char *maxname, char *binname, int *status);
739int CFITS_API fits_copy_cell2image(fitsfile *fptr, fitsfile *newptr, char *colname,
740 long rownum, int *status);
741int CFITS_API fits_copy_image2cell(fitsfile *fptr, fitsfile *newptr, char *colname,
742 long rownum, int copykeyflag, int *status);
743int CFITS_API fits_copy_pixlist2image(fitsfile *infptr, fitsfile *outfptr, int firstkey, /* I - first HDU record number to start with */
744 int naxis, int *colnum, int *status);
745int CFITS_API ffimport_file( char *filename, char **contents, int *status );
746int CFITS_API ffrwrg( char *rowlist, LONGLONG maxrows, int maxranges, int *numranges,
747 long *minrow, long *maxrow, int *status);
748int CFITS_API ffrwrgll( char *rowlist, LONGLONG maxrows, int maxranges, int *numranges,
749 LONGLONG *minrow, LONGLONG *maxrow, int *status);
750/*---------------- FITS file I/O routines -------------*/
751int CFITS_API fits_init_cfitsio(void);
752int CFITS_API ffomem(fitsfile **fptr, const char *name, int mode, void **buffptr,
753 size_t *buffsize, size_t deltasize,
754 void *(*mem_realloc)(void *p, size_t newsize),
755 int *status);
756int CFITS_API ffopen(fitsfile **fptr, const char *filename, int iomode, int *status);
757int CFITS_API ffopentest(int soname, fitsfile **fptr, const char *filename, int iomode, int *status);
758
759int CFITS_API ffdopn(fitsfile **fptr, const char *filename, int iomode, int *status);
760int CFITS_API fftopn(fitsfile **fptr, const char *filename, int iomode, int *status);
761int CFITS_API ffiopn(fitsfile **fptr, const char *filename, int iomode, int *status);
762int CFITS_API ffdkopn(fitsfile **fptr, const char *filename, int iomode, int *status);
763int CFITS_API ffreopen(fitsfile *openfptr, fitsfile **newfptr, int *status);
764int CFITS_API ffinit( fitsfile **fptr, const char *filename, int *status);
765int CFITS_API ffdkinit(fitsfile **fptr, const char *filename, int *status);
766int CFITS_API ffimem(fitsfile **fptr, void **buffptr,
767 size_t *buffsize, size_t deltasize,
768 void *(*mem_realloc)(void *p, size_t newsize),
769 int *status);
770int CFITS_API fftplt(fitsfile **fptr, const char *filename, const char *tempname,
771 int *status);
772int CFITS_API ffflus(fitsfile *fptr, int *status);
773int CFITS_API ffflsh(fitsfile *fptr, int clearbuf, int *status);
774int CFITS_API ffclos(fitsfile *fptr, int *status);
775int CFITS_API ffdelt(fitsfile *fptr, int *status);
776int CFITS_API ffflnm(fitsfile *fptr, char *filename, int *status);
777int CFITS_API ffflmd(fitsfile *fptr, int *filemode, int *status);
778int CFITS_API fits_delete_iraf_file(const char *filename, int *status);
779
780/*---------------- utility routines -------------*/
781
782float CFITS_API ffvers(float *version);
783void CFITS_API ffupch(char *string);
784void CFITS_API ffgerr(int status, char *errtext);
785void CFITS_API ffpmsg(const char *err_message);
786void CFITS_API ffpmrk(void);
787int CFITS_API ffgmsg(char *err_message);
788void CFITS_API ffcmsg(void);
789void CFITS_API ffcmrk(void);
790void CFITS_API ffrprt(FILE *stream, int status);
791void CFITS_API ffcmps(char *templt, char *colname, int casesen, int *match,
792 int *exact);
793int CFITS_API fftkey(const char *keyword, int *status);
794int CFITS_API fftrec(char *card, int *status);
795int CFITS_API ffnchk(fitsfile *fptr, int *status);
796int CFITS_API ffkeyn(const char *keyroot, int value, char *keyname, int *status);
797int CFITS_API ffnkey(int value, const char *keyroot, char *keyname, int *status);
798int CFITS_API ffgkcl(char *card);
799int CFITS_API ffdtyp(const char *cval, char *dtype, int *status);
800int CFITS_API ffinttyp(char *cval, int *datatype, int *negative, int *status);
801int CFITS_API ffpsvc(char *card, char *value, char *comm, int *status);
802int CFITS_API ffgknm(char *card, char *name, int *length, int *status);
803int CFITS_API ffgthd(char *tmplt, char *card, int *hdtype, int *status);
804int CFITS_API ffmkky(const char *keyname, char *keyval, const char *comm, char *card, int *status);
805int CFITS_API fits_translate_keyword(char *inrec, char *outrec, char *patterns[][2],
806 int npat, int n_value, int n_offset, int n_range, int *pat_num,
807 int *i, int *j, int *m, int *n, int *status);
808int CFITS_API fits_translate_keywords(fitsfile *infptr, fitsfile *outfptr,
809 int firstkey, char *patterns[][2],
810 int npat, int n_value, int n_offset, int n_range, int *status);
811int CFITS_API ffasfm(char *tform, int *datacode, long *width, int *decim, int *status);
812int CFITS_API ffbnfm(char *tform, int *datacode, long *repeat, long *width, int *status);
813int CFITS_API ffbnfmll(char *tform, int *datacode, LONGLONG *repeat, long *width, int *status);
814int CFITS_API ffgabc(int tfields, char **tform, int space, long *rowlen, long *tbcol,
815 int *status);
816int CFITS_API fits_get_section_range(char **ptr,long *secmin,long *secmax,long *incre,
817 int *status);
818/* ffmbyt should not normally be used in application programs, but it is
819 defined here as a publicly available routine because there are a few
820 rare cases where it is needed
821*/
822int CFITS_API ffmbyt(fitsfile *fptr, LONGLONG bytpos, int ignore_err, int *status);
823/*----------------- write single keywords --------------*/
824int CFITS_API ffpky(fitsfile *fptr, int datatype, const char *keyname, void *value,
825 const char *comm, int *status);
826int CFITS_API ffprec(fitsfile *fptr, const char *card, int *status);
827int CFITS_API ffpcom(fitsfile *fptr, const char *comm, int *status);
828int CFITS_API ffpunt(fitsfile *fptr, const char *keyname, const char *unit, int *status);
829int CFITS_API ffphis(fitsfile *fptr, const char *history, int *status);
830int CFITS_API ffpdat(fitsfile *fptr, int *status);
831int CFITS_API ffverifydate(int year, int month, int day, int *status);
832int CFITS_API ffgstm(char *timestr, int *timeref, int *status);
833int CFITS_API ffgsdt(int *day, int *month, int *year, int *status);
834int CFITS_API ffdt2s(int year, int month, int day, char *datestr, int *status);
835int CFITS_API fftm2s(int year, int month, int day, int hour, int minute, double second,
836 int decimals, char *datestr, int *status);
837int CFITS_API ffs2dt(char *datestr, int *year, int *month, int *day, int *status);
838int CFITS_API ffs2tm(char *datestr, int *year, int *month, int *day, int *hour,
839 int *minute, double *second, int *status);
840int CFITS_API ffpkyu(fitsfile *fptr, const char *keyname, const char *comm, int *status);
841int CFITS_API ffpkys(fitsfile *fptr, const char *keyname, const char *value, const char *comm,int *status);
842int CFITS_API ffpkls(fitsfile *fptr, const char *keyname, const char *value, const char *comm,int *status);
843int CFITS_API ffplsw(fitsfile *fptr, int *status);
844int CFITS_API ffpkyl(fitsfile *fptr, const char *keyname, int value, const char *comm, int *status);
845int CFITS_API ffpkyj(fitsfile *fptr, const char *keyname, LONGLONG value, const char *comm, int *status);
846int CFITS_API ffpkyf(fitsfile *fptr, const char *keyname, float value, int decim, const char *comm,
847 int *status);
848int CFITS_API ffpkye(fitsfile *fptr, const char *keyname, float value, int decim, const char *comm,
849 int *status);
850int CFITS_API ffpkyg(fitsfile *fptr, const char *keyname, double value, int decim, const char *comm,
851 int *status);
852int CFITS_API ffpkyd(fitsfile *fptr, const char *keyname, double value, int decim, const char *comm,
853 int *status);
854int CFITS_API ffpkyc(fitsfile *fptr, const char *keyname, float *value, int decim, const char *comm,
855 int *status);
856int CFITS_API ffpkym(fitsfile *fptr, const char *keyname, double *value, int decim, const char *comm,
857 int *status);
858int CFITS_API ffpkfc(fitsfile *fptr, const char *keyname, float *value, int decim, const char *comm,
859 int *status);
860int CFITS_API ffpkfm(fitsfile *fptr, const char *keyname, double *value, int decim, const char *comm,
861 int *status);
862int CFITS_API ffpkyt(fitsfile *fptr, const char *keyname, long intval, double frac, const char *comm,
863 int *status);
864int CFITS_API ffptdm( fitsfile *fptr, int colnum, int naxis, long naxes[], int *status);
865int CFITS_API ffptdmll( fitsfile *fptr, int colnum, int naxis, LONGLONG naxes[], int *status);
866
867/*----------------- write array of keywords --------------*/
868int CFITS_API ffpkns(fitsfile *fptr, const char *keyroot, int nstart, int nkey, char *value[],
869 char *comm[], int *status);
870int CFITS_API ffpknl(fitsfile *fptr, const char *keyroot, int nstart, int nkey, int *value,
871 char *comm[], int *status);
872int CFITS_API ffpknj(fitsfile *fptr, const char *keyroot, int nstart, int nkey, long *value,
873 char *comm[], int *status);
874int CFITS_API ffpknjj(fitsfile *fptr, const char *keyroot, int nstart, int nkey, LONGLONG *value,
875 char *comm[], int *status);
876int CFITS_API ffpknf(fitsfile *fptr, const char *keyroot, int nstart, int nkey, float *value,
877 int decim, char *comm[], int *status);
878int CFITS_API ffpkne(fitsfile *fptr, const char *keyroot, int nstart, int nkey, float *value,
879 int decim, char *comm[], int *status);
880int CFITS_API ffpkng(fitsfile *fptr, const char *keyroot, int nstart, int nkey, double *value,
881 int decim, char *comm[], int *status);
882int CFITS_API ffpknd(fitsfile *fptr, const char *keyroot, int nstart, int nkey, double *value,
883 int decim, char *comm[], int *status);
884int CFITS_API ffcpky(fitsfile *infptr,fitsfile *outfptr,int incol,int outcol,
885 char *rootname, int *status);
886
887/*----------------- write required header keywords --------------*/
888int CFITS_API ffphps( fitsfile *fptr, int bitpix, int naxis, long naxes[], int *status);
889int CFITS_API ffphpsll( fitsfile *fptr, int bitpix, int naxis, LONGLONG naxes[], int *status);
890int CFITS_API ffphpr( fitsfile *fptr, int simple, int bitpix, int naxis, long naxes[],
891 LONGLONG pcount, LONGLONG gcount, int extend, int *status);
892int CFITS_API ffphprll( fitsfile *fptr, int simple, int bitpix, int naxis, LONGLONG naxes[],
893 LONGLONG pcount, LONGLONG gcount, int extend, int *status);
894int CFITS_API ffphtb(fitsfile *fptr, LONGLONG naxis1, LONGLONG naxis2, int tfields, char **ttype,
895 long *tbcol, char **tform, char **tunit, const char *extname, int *status);
896int CFITS_API ffphbn(fitsfile *fptr, LONGLONG naxis2, int tfields, char **ttype,
897 char **tform, char **tunit, const char *extname, LONGLONG pcount, int *status);
898int CFITS_API ffphext( fitsfile *fptr, const char *xtension, int bitpix, int naxis, long naxes[],
899 LONGLONG pcount, LONGLONG gcount, int *status);
900/*----------------- write template keywords --------------*/
901int CFITS_API ffpktp(fitsfile *fptr, const char *filename, int *status);
902
903/*------------------ get header information --------------*/
904int CFITS_API ffghsp(fitsfile *fptr, int *nexist, int *nmore, int *status);
905int CFITS_API ffghps(fitsfile *fptr, int *nexist, int *position, int *status);
906
907/*------------------ move position in header -------------*/
908int CFITS_API ffmaky(fitsfile *fptr, int nrec, int *status);
909int CFITS_API ffmrky(fitsfile *fptr, int nrec, int *status);
910
911/*------------------ read single keywords -----------------*/
912int CFITS_API ffgnxk(fitsfile *fptr, char **inclist, int ninc, char **exclist,
913 int nexc, char *card, int *status);
914int CFITS_API ffgrec(fitsfile *fptr, int nrec, char *card, int *status);
915int CFITS_API ffgcrd(fitsfile *fptr, const char *keyname, char *card, int *status);
916int CFITS_API ffgstr(fitsfile *fptr, const char *string, char *card, int *status);
917int CFITS_API ffgunt(fitsfile *fptr, const char *keyname, char *unit, int *status);
918int CFITS_API ffgkyn(fitsfile *fptr, int nkey, char *keyname, char *keyval, char *comm,
919 int *status);
920int CFITS_API ffgkey(fitsfile *fptr, const char *keyname, char *keyval, char *comm,
921 int *status);
922
923int CFITS_API ffgky( fitsfile *fptr, int datatype, const char *keyname, void *value,
924 char *comm, int *status);
925int CFITS_API ffgkys(fitsfile *fptr, const char *keyname, char *value, char *comm, int *status);
926int CFITS_API ffgkls(fitsfile *fptr, const char *keyname, char **value, char *comm, int *status);
927int CFITS_API fffree(void *value, int *status);
928int CFITS_API fffkls(char *value, int *status);
929int CFITS_API ffgkyl(fitsfile *fptr, const char *keyname, int *value, char *comm, int *status);
930int CFITS_API ffgkyj(fitsfile *fptr, const char *keyname, long *value, char *comm, int *status);
931int CFITS_API ffgkyjj(fitsfile *fptr, const char *keyname, LONGLONG *value, char *comm, int *status);
932int CFITS_API ffgkye(fitsfile *fptr, const char *keyname, float *value, char *comm,int *status);
933int CFITS_API ffgkyd(fitsfile *fptr, const char *keyname, double *value,char *comm,int *status);
934int CFITS_API ffgkyc(fitsfile *fptr, const char *keyname, float *value, char *comm,int *status);
935int CFITS_API ffgkym(fitsfile *fptr, const char *keyname, double *value,char *comm,int *status);
936int CFITS_API ffgkyt(fitsfile *fptr, const char *keyname, long *ivalue, double *dvalue,
937 char *comm, int *status);
938int CFITS_API ffgtdm(fitsfile *fptr, int colnum, int maxdim, int *naxis, long naxes[],
939 int *status);
940int CFITS_API ffgtdmll(fitsfile *fptr, int colnum, int maxdim, int *naxis, LONGLONG naxes[],
941 int *status);
942int CFITS_API ffdtdm(fitsfile *fptr, char *tdimstr, int colnum, int maxdim,
943 int *naxis, long naxes[], int *status);
944int CFITS_API ffdtdmll(fitsfile *fptr, char *tdimstr, int colnum, int maxdim,
945 int *naxis, LONGLONG naxes[], int *status);
946
947/*------------------ read array of keywords -----------------*/
948int CFITS_API ffgkns(fitsfile *fptr, const char *keyname, int nstart, int nmax, char *value[],
949 int *nfound, int *status);
950int CFITS_API ffgknl(fitsfile *fptr, const char *keyname, int nstart, int nmax, int *value,
951 int *nfound, int *status);
952int CFITS_API ffgknj(fitsfile *fptr, const char *keyname, int nstart, int nmax, long *value,
953 int *nfound, int *status);
954int CFITS_API ffgknjj(fitsfile *fptr, const char *keyname, int nstart, int nmax, LONGLONG *value,
955 int *nfound, int *status);
956int CFITS_API ffgkne(fitsfile *fptr, const char *keyname, int nstart, int nmax, float *value,
957 int *nfound, int *status);
958int CFITS_API ffgknd(fitsfile *fptr, const char *keyname, int nstart, int nmax, double *value,
959 int *nfound, int *status);
960int CFITS_API ffh2st(fitsfile *fptr, char **header, int *status);
961int CFITS_API ffhdr2str( fitsfile *fptr, int exclude_comm, char **exclist,
962 int nexc, char **header, int *nkeys, int *status);
963int CFITS_API ffcnvthdr2str( fitsfile *fptr, int exclude_comm, char **exclist,
964 int nexc, char **header, int *nkeys, int *status);
965
966/*----------------- read required header keywords --------------*/
967int CFITS_API ffghpr(fitsfile *fptr, int maxdim, int *simple, int *bitpix, int *naxis,
968 long naxes[], long *pcount, long *gcount, int *extend, int *status);
969
970int CFITS_API ffghprll(fitsfile *fptr, int maxdim, int *simple, int *bitpix, int *naxis,
971 LONGLONG naxes[], long *pcount, long *gcount, int *extend, int *status);
972
973int CFITS_API ffghtb(fitsfile *fptr,int maxfield, long *naxis1, long *naxis2,
974 int *tfields, char **ttype, long *tbcol, char **tform, char **tunit,
975 char *extname, int *status);
976
977int CFITS_API ffghtbll(fitsfile *fptr,int maxfield, LONGLONG *naxis1, LONGLONG *naxis2,
978 int *tfields, char **ttype, LONGLONG *tbcol, char **tform, char **tunit,
979 char *extname, int *status);
980
981
982int CFITS_API ffghbn(fitsfile *fptr, int maxfield, long *naxis2, int *tfields,
983 char **ttype, char **tform, char **tunit, char *extname,
984 long *pcount, int *status);
985
986int CFITS_API ffghbnll(fitsfile *fptr, int maxfield, LONGLONG *naxis2, int *tfields,
987 char **ttype, char **tform, char **tunit, char *extname,
988 LONGLONG *pcount, int *status);
989
990/*--------------------- update keywords ---------------*/
991int CFITS_API ffuky(fitsfile *fptr, int datatype, const char *keyname, void *value,
992 const char *comm, int *status);
993int CFITS_API ffucrd(fitsfile *fptr, const char *keyname, const char *card, int *status);
994int CFITS_API ffukyu(fitsfile *fptr, const char *keyname, const char *comm, int *status);
995int CFITS_API ffukys(fitsfile *fptr, const char *keyname, const char *value, const char *comm, int *status);
996int CFITS_API ffukls(fitsfile *fptr, const char *keyname, const char *value, const char *comm, int *status);
997int CFITS_API ffukyl(fitsfile *fptr, const char *keyname, int value, const char *comm, int *status);
998int CFITS_API ffukyj(fitsfile *fptr, const char *keyname, LONGLONG value, const char *comm, int *status);
999int CFITS_API ffukyf(fitsfile *fptr, const char *keyname, float value, int decim, const char *comm,
1000 int *status);
1001int CFITS_API ffukye(fitsfile *fptr, const char *keyname, float value, int decim, const char *comm,
1002 int *status);
1003int CFITS_API ffukyg(fitsfile *fptr, const char *keyname, double value, int decim, const char *comm,
1004 int *status);
1005int CFITS_API ffukyd(fitsfile *fptr, const char *keyname, double value, int decim, const char *comm,
1006 int *status);
1007int CFITS_API ffukyc(fitsfile *fptr, const char *keyname, float *value, int decim, const char *comm,
1008 int *status);
1009int CFITS_API ffukym(fitsfile *fptr, const char *keyname, double *value, int decim, const char *comm,
1010 int *status);
1011int CFITS_API ffukfc(fitsfile *fptr, const char *keyname, float *value, int decim, const char *comm,
1012 int *status);
1013int CFITS_API ffukfm(fitsfile *fptr, const char *keyname, double *value, int decim, const char *comm,
1014 int *status);
1015
1016/*--------------------- modify keywords ---------------*/
1017int CFITS_API ffmrec(fitsfile *fptr, int nkey, const char *card, int *status);
1018int CFITS_API ffmcrd(fitsfile *fptr, const char *keyname, const char *card, int *status);
1019int CFITS_API ffmnam(fitsfile *fptr, const char *oldname, const char *newname, int *status);
1020int CFITS_API ffmcom(fitsfile *fptr, const char *keyname, const char *comm, int *status);
1021int CFITS_API ffmkyu(fitsfile *fptr, const char *keyname, const char *comm, int *status);
1022int CFITS_API ffmkys(fitsfile *fptr, const char *keyname, const char *value, const char *comm,int *status);
1023int CFITS_API ffmkls(fitsfile *fptr, const char *keyname, const char *value, const char *comm,int *status);
1024int CFITS_API ffmkyl(fitsfile *fptr, const char *keyname, int value, const char *comm, int *status);
1025int CFITS_API ffmkyj(fitsfile *fptr, const char *keyname, LONGLONG value, const char *comm, int *status);
1026int CFITS_API ffmkyf(fitsfile *fptr, const char *keyname, float value, int decim, const char *comm,
1027 int *status);
1028int CFITS_API ffmkye(fitsfile *fptr, const char *keyname, float value, int decim, const char *comm,
1029 int *status);
1030int CFITS_API ffmkyg(fitsfile *fptr, const char *keyname, double value, int decim, const char *comm,
1031 int *status);
1032int CFITS_API ffmkyd(fitsfile *fptr, const char *keyname, double value, int decim, const char *comm,
1033 int *status);
1034int CFITS_API ffmkyc(fitsfile *fptr, const char *keyname, float *value, int decim, const char *comm,
1035 int *status);
1036int CFITS_API ffmkym(fitsfile *fptr, const char *keyname, double *value, int decim, const char *comm,
1037 int *status);
1038int CFITS_API ffmkfc(fitsfile *fptr, const char *keyname, float *value, int decim, const char *comm,
1039 int *status);
1040int CFITS_API ffmkfm(fitsfile *fptr, const char *keyname, double *value, int decim, const char *comm,
1041 int *status);
1042
1043/*--------------------- insert keywords ---------------*/
1044int CFITS_API ffirec(fitsfile *fptr, int nkey, const char *card, int *status);
1045int CFITS_API ffikey(fitsfile *fptr, const char *card, int *status);
1046int CFITS_API ffikyu(fitsfile *fptr, const char *keyname, const char *comm, int *status);
1047int CFITS_API ffikys(fitsfile *fptr, const char *keyname, const char *value, const char *comm,int *status);
1048int CFITS_API ffikls(fitsfile *fptr, const char *keyname, const char *value, const char *comm,int *status);
1049int CFITS_API ffikyl(fitsfile *fptr, const char *keyname, int value, const char *comm, int *status);
1050int CFITS_API ffikyj(fitsfile *fptr, const char *keyname, LONGLONG value, const char *comm, int *status);
1051int CFITS_API ffikyf(fitsfile *fptr, const char *keyname, float value, int decim, const char *comm,
1052 int *status);
1053int CFITS_API ffikye(fitsfile *fptr, const char *keyname, float value, int decim, const char *comm,
1054 int *status);
1055int CFITS_API ffikyg(fitsfile *fptr, const char *keyname, double value, int decim, const char *comm,
1056 int *status);
1057int CFITS_API ffikyd(fitsfile *fptr, const char *keyname, double value, int decim, const char *comm,
1058 int *status);
1059int CFITS_API ffikyc(fitsfile *fptr, const char *keyname, float *value, int decim, const char *comm,
1060 int *status);
1061int CFITS_API ffikym(fitsfile *fptr, const char *keyname, double *value, int decim, const char *comm,
1062 int *status);
1063int CFITS_API ffikfc(fitsfile *fptr, const char *keyname, float *value, int decim, const char *comm,
1064 int *status);
1065int CFITS_API ffikfm(fitsfile *fptr, const char *keyname, double *value, int decim, const char *comm,
1066 int *status);
1067
1068/*--------------------- delete keywords ---------------*/
1069int CFITS_API ffdkey(fitsfile *fptr, const char *keyname, int *status);
1070int CFITS_API ffdstr(fitsfile *fptr, const char *string, int *status);
1071int CFITS_API ffdrec(fitsfile *fptr, int keypos, int *status);
1072
1073/*--------------------- get HDU information -------------*/
1074int CFITS_API ffghdn(fitsfile *fptr, int *chdunum);
1075int CFITS_API ffghdt(fitsfile *fptr, int *exttype, int *status);
1076int CFITS_API ffghad(fitsfile *fptr, long *headstart, long *datastart, long *dataend,
1077 int *status);
1078int CFITS_API ffghadll(fitsfile *fptr, LONGLONG *headstart, LONGLONG *datastart,
1079 LONGLONG *dataend, int *status);
1080int CFITS_API ffghof(fitsfile *fptr, OFF_T *headstart, OFF_T *datastart, OFF_T *dataend,
1081 int *status);
1082int CFITS_API ffgipr(fitsfile *fptr, int maxaxis, int *imgtype, int *naxis,
1083 long *naxes, int *status);
1084int CFITS_API ffgiprll(fitsfile *fptr, int maxaxis, int *imgtype, int *naxis,
1085 LONGLONG *naxes, int *status);
1086int CFITS_API ffgidt(fitsfile *fptr, int *imgtype, int *status);
1087int CFITS_API ffgiet(fitsfile *fptr, int *imgtype, int *status);
1088int CFITS_API ffgidm(fitsfile *fptr, int *naxis, int *status);
1089int CFITS_API ffgisz(fitsfile *fptr, int nlen, long *naxes, int *status);
1090int CFITS_API ffgiszll(fitsfile *fptr, int nlen, LONGLONG *naxes, int *status);
1091
1092/*--------------------- HDU operations -------------*/
1093int CFITS_API ffmahd(fitsfile *fptr, int hdunum, int *exttype, int *status);
1094int CFITS_API ffmrhd(fitsfile *fptr, int hdumov, int *exttype, int *status);
1095int CFITS_API ffmnhd(fitsfile *fptr, int exttype, char *hduname, int hduvers,
1096 int *status);
1097int CFITS_API ffthdu(fitsfile *fptr, int *nhdu, int *status);
1098int CFITS_API ffcrhd(fitsfile *fptr, int *status);
1099int CFITS_API ffcrim(fitsfile *fptr, int bitpix, int naxis, long *naxes, int *status);
1100int CFITS_API ffcrimll(fitsfile *fptr, int bitpix, int naxis, LONGLONG *naxes, int *status);
1101int CFITS_API ffcrtb(fitsfile *fptr, int tbltype, LONGLONG naxis2, int tfields, char **ttype,
1102 char **tform, char **tunit, const char *extname, int *status);
1103int CFITS_API ffiimg(fitsfile *fptr, int bitpix, int naxis, long *naxes, int *status);
1104int CFITS_API ffiimgll(fitsfile *fptr, int bitpix, int naxis, LONGLONG *naxes, int *status);
1105int CFITS_API ffitab(fitsfile *fptr, LONGLONG naxis1, LONGLONG naxis2, int tfields, char **ttype,
1106 long *tbcol, char **tform, char **tunit, const char *extname, int *status);
1107int CFITS_API ffibin(fitsfile *fptr, LONGLONG naxis2, int tfields, char **ttype, char **tform,
1108 char **tunit, const char *extname, LONGLONG pcount, int *status);
1109int CFITS_API ffrsim(fitsfile *fptr, int bitpix, int naxis, long *naxes, int *status);
1110int CFITS_API ffrsimll(fitsfile *fptr, int bitpix, int naxis, LONGLONG *naxes, int *status);
1111int CFITS_API ffdhdu(fitsfile *fptr, int *hdutype, int *status);
1112int CFITS_API ffcopy(fitsfile *infptr, fitsfile *outfptr, int morekeys, int *status);
1113int CFITS_API ffcpfl(fitsfile *infptr, fitsfile *outfptr, int prev, int cur, int follow,
1114 int *status);
1115int CFITS_API ffcphd(fitsfile *infptr, fitsfile *outfptr, int *status);
1116int CFITS_API ffcpdt(fitsfile *infptr, fitsfile *outfptr, int *status);
1117int CFITS_API ffchfl(fitsfile *fptr, int *status);
1118int CFITS_API ffcdfl(fitsfile *fptr, int *status);
1119int CFITS_API ffwrhdu(fitsfile *fptr, FILE *outstream, int *status);
1120
1121int CFITS_API ffrdef(fitsfile *fptr, int *status);
1122int CFITS_API ffhdef(fitsfile *fptr, int morekeys, int *status);
1123int CFITS_API ffpthp(fitsfile *fptr, long theap, int *status);
1124
1125int CFITS_API ffcsum(fitsfile *fptr, long nrec, unsigned long *sum, int *status);
1126void CFITS_API ffesum(unsigned long sum, int complm, char *ascii);
1127unsigned long CFITS_API ffdsum(char *ascii, int complm, unsigned long *sum);
1128int CFITS_API ffpcks(fitsfile *fptr, int *status);
1129int CFITS_API ffupck(fitsfile *fptr, int *status);
1130int CFITS_API ffvcks(fitsfile *fptr, int *datastatus, int *hdustatus, int *status);
1131int CFITS_API ffgcks(fitsfile *fptr, unsigned long *datasum, unsigned long *hdusum,
1132 int *status);
1133
1134/*--------------------- define scaling or null values -------------*/
1135int CFITS_API ffpscl(fitsfile *fptr, double scale, double zero, int *status);
1136int CFITS_API ffpnul(fitsfile *fptr, LONGLONG nulvalue, int *status);
1137int CFITS_API fftscl(fitsfile *fptr, int colnum, double scale, double zero, int *status);
1138int CFITS_API fftnul(fitsfile *fptr, int colnum, LONGLONG nulvalue, int *status);
1139int CFITS_API ffsnul(fitsfile *fptr, int colnum, char *nulstring, int *status);
1140
1141/*--------------------- get column information -------------*/
1142int CFITS_API ffgcno(fitsfile *fptr, int casesen, char *templt, int *colnum,
1143 int *status);
1144int CFITS_API ffgcnn(fitsfile *fptr, int casesen, char *templt, char *colname,
1145 int *colnum, int *status);
1146
1147int CFITS_API ffgtcl(fitsfile *fptr, int colnum, int *typecode, long *repeat,
1148 long *width, int *status);
1149int CFITS_API ffgtclll(fitsfile *fptr, int colnum, int *typecode, LONGLONG *repeat,
1150 LONGLONG *width, int *status);
1151int CFITS_API ffeqty(fitsfile *fptr, int colnum, int *typecode, long *repeat,
1152 long *width, int *status);
1153int CFITS_API ffeqtyll(fitsfile *fptr, int colnum, int *typecode, LONGLONG *repeat,
1154 LONGLONG *width, int *status);
1155int CFITS_API ffgncl(fitsfile *fptr, int *ncols, int *status);
1156int CFITS_API ffgnrw(fitsfile *fptr, long *nrows, int *status);
1157int CFITS_API ffgnrwll(fitsfile *fptr, LONGLONG *nrows, int *status);
1158int CFITS_API ffgacl(fitsfile *fptr, int colnum, char *ttype, long *tbcol,
1159 char *tunit, char *tform, double *tscal, double *tzero,
1160 char *tnull, char *tdisp, int *status);
1161int CFITS_API ffgbcl(fitsfile *fptr, int colnum, char *ttype, char *tunit,
1162 char *dtype, long *repeat, double *tscal, double *tzero,
1163 long *tnull, char *tdisp, int *status);
1164int CFITS_API ffgbclll(fitsfile *fptr, int colnum, char *ttype, char *tunit,
1165 char *dtype, LONGLONG *repeat, double *tscal, double *tzero,
1166 LONGLONG *tnull, char *tdisp, int *status);
1167int CFITS_API ffgrsz(fitsfile *fptr, long *nrows, int *status);
1168int CFITS_API ffgcdw(fitsfile *fptr, int colnum, int *width, int *status);
1169
1170/*--------------------- read primary array or image elements -------------*/
1171int CFITS_API ffgpxv(fitsfile *fptr, int datatype, long *firstpix, LONGLONG nelem,
1172 void *nulval, void *array, int *anynul, int *status);
1173int CFITS_API ffgpxvll(fitsfile *fptr, int datatype, LONGLONG *firstpix, LONGLONG nelem,
1174 void *nulval, void *array, int *anynul, int *status);
1175int CFITS_API ffgpxf(fitsfile *fptr, int datatype, long *firstpix, LONGLONG nelem,
1176 void *array, char *nullarray, int *anynul, int *status);
1177int CFITS_API ffgpxfll(fitsfile *fptr, int datatype, LONGLONG *firstpix, LONGLONG nelem,
1178 void *array, char *nullarray, int *anynul, int *status);
1179int CFITS_API ffgsv(fitsfile *fptr, int datatype, long *blc, long *trc, long *inc,
1180 void *nulval, void *array, int *anynul, int *status);
1181
1182int CFITS_API ffgpv(fitsfile *fptr, int datatype, LONGLONG firstelem, LONGLONG nelem,
1183 void *nulval, void *array, int *anynul, int *status);
1184int CFITS_API ffgpf(fitsfile *fptr, int datatype, LONGLONG firstelem, LONGLONG nelem,
1185 void *array, char *nullarray, int *anynul, int *status);
1186int CFITS_API ffgpvb(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem, unsigned
1187 char nulval, unsigned char *array, int *anynul, int *status);
1188int CFITS_API ffgpvsb(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem, signed
1189 char nulval, signed char *array, int *anynul, int *status);
1190int CFITS_API ffgpvui(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1191 unsigned short nulval, unsigned short *array, int *anynul,
1192 int *status);
1193int CFITS_API ffgpvi(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1194 short nulval, short *array, int *anynul, int *status);
1195int CFITS_API ffgpvuj(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1196 unsigned long nulval, unsigned long *array, int *anynul,
1197 int *status);
1198int CFITS_API ffgpvj(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1199 long nulval, long *array, int *anynul, int *status);
1200int CFITS_API ffgpvjj(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1201 LONGLONG nulval, LONGLONG *array, int *anynul, int *status);
1202int CFITS_API ffgpvuk(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1203 unsigned int nulval, unsigned int *array, int *anynul, int *status);
1204int CFITS_API ffgpvk(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1205 int nulval, int *array, int *anynul, int *status);
1206int CFITS_API ffgpve(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1207 float nulval, float *array, int *anynul, int *status);
1208int CFITS_API ffgpvd(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1209 double nulval, double *array, int *anynul, int *status);
1210
1211int CFITS_API ffgpfb(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1212 unsigned char *array, char *nularray, int *anynul, int *status);
1213int CFITS_API ffgpfsb(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1214 signed char *array, char *nularray, int *anynul, int *status);
1215int CFITS_API ffgpfui(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1216 unsigned short *array, char *nularray, int *anynul, int *status);
1217int CFITS_API ffgpfi(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1218 short *array, char *nularray, int *anynul, int *status);
1219int CFITS_API ffgpfuj(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1220 unsigned long *array, char *nularray, int *anynul, int *status);
1221int CFITS_API ffgpfj(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1222 long *array, char *nularray, int *anynul, int *status);
1223int CFITS_API ffgpfjj(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1224 LONGLONG *array, char *nularray, int *anynul, int *status);
1225int CFITS_API ffgpfuk(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1226 unsigned int *array, char *nularray, int *anynul, int *status);
1227int CFITS_API ffgpfk(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1228 int *array, char *nularray, int *anynul, int *status);
1229int CFITS_API ffgpfe(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1230 float *array, char *nularray, int *anynul, int *status);
1231int CFITS_API ffgpfd(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1232 double *array, char *nularray, int *anynul, int *status);
1233
1234int CFITS_API ffg2db(fitsfile *fptr, long group, unsigned char nulval, LONGLONG ncols,
1235 LONGLONG naxis1, LONGLONG naxis2, unsigned char *array,
1236 int *anynul, int *status);
1237int CFITS_API ffg2dsb(fitsfile *fptr, long group, signed char nulval, LONGLONG ncols,
1238 LONGLONG naxis1, LONGLONG naxis2, signed char *array,
1239 int *anynul, int *status);
1240int CFITS_API ffg2dui(fitsfile *fptr, long group, unsigned short nulval, LONGLONG ncols,
1241 LONGLONG naxis1, LONGLONG naxis2, unsigned short *array,
1242 int *anynul, int *status);
1243int CFITS_API ffg2di(fitsfile *fptr, long group, short nulval, LONGLONG ncols,
1244 LONGLONG naxis1, LONGLONG naxis2, short *array,
1245 int *anynul, int *status);
1246int CFITS_API ffg2duj(fitsfile *fptr, long group, unsigned long nulval, LONGLONG ncols,
1247 LONGLONG naxis1, LONGLONG naxis2, unsigned long *array,
1248 int *anynul, int *status);
1249int CFITS_API ffg2dj(fitsfile *fptr, long group, long nulval, LONGLONG ncols,
1250 LONGLONG naxis1, LONGLONG naxis2, long *array,
1251 int *anynul, int *status);
1252int CFITS_API ffg2djj(fitsfile *fptr, long group, LONGLONG nulval, LONGLONG ncols,
1253 LONGLONG naxis1, LONGLONG naxis2, LONGLONG *array,
1254 int *anynul, int *status);
1255int CFITS_API ffg2duk(fitsfile *fptr, long group, unsigned int nulval, LONGLONG ncols,
1256 LONGLONG naxis1, LONGLONG naxis2, unsigned int *array,
1257 int *anynul, int *status);
1258int CFITS_API ffg2dk(fitsfile *fptr, long group, int nulval, LONGLONG ncols,
1259 LONGLONG naxis1, LONGLONG naxis2, int *array,
1260 int *anynul, int *status);
1261int CFITS_API ffg2de(fitsfile *fptr, long group, float nulval, LONGLONG ncols,
1262 LONGLONG naxis1, LONGLONG naxis2, float *array,
1263 int *anynul, int *status);
1264int CFITS_API ffg2dd(fitsfile *fptr, long group, double nulval, LONGLONG ncols,
1265 LONGLONG naxis1, LONGLONG naxis2, double *array,
1266 int *anynul, int *status);
1267
1268int CFITS_API ffg3db(fitsfile *fptr, long group, unsigned char nulval, LONGLONG ncols,
1269 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1270 unsigned char *array, int *anynul, int *status);
1271int CFITS_API ffg3dsb(fitsfile *fptr, long group, signed char nulval, LONGLONG ncols,
1272 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1273 signed char *array, int *anynul, int *status);
1274int CFITS_API ffg3dui(fitsfile *fptr, long group, unsigned short nulval, LONGLONG ncols,
1275 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1276 unsigned short *array, int *anynul, int *status);
1277int CFITS_API ffg3di(fitsfile *fptr, long group, short nulval, LONGLONG ncols,
1278 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1279 short *array, int *anynul, int *status);
1280int CFITS_API ffg3duj(fitsfile *fptr, long group, unsigned long nulval, LONGLONG ncols,
1281 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1282 unsigned long *array, int *anynul, int *status);
1283int CFITS_API ffg3dj(fitsfile *fptr, long group, long nulval, LONGLONG ncols,
1284 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1285 long *array, int *anynul, int *status);
1286int CFITS_API ffg3djj(fitsfile *fptr, long group, LONGLONG nulval, LONGLONG ncols,
1287 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1288 LONGLONG *array, int *anynul, int *status);
1289int CFITS_API ffg3duk(fitsfile *fptr, long group, unsigned int nulval, LONGLONG ncols,
1290 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1291 unsigned int *array, int *anynul, int *status);
1292int CFITS_API ffg3dk(fitsfile *fptr, long group, int nulval, LONGLONG ncols,
1293 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1294 int *array, int *anynul, int *status);
1295int CFITS_API ffg3de(fitsfile *fptr, long group, float nulval, LONGLONG ncols,
1296 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1297 float *array, int *anynul, int *status);
1298int CFITS_API ffg3dd(fitsfile *fptr, long group, double nulval, LONGLONG ncols,
1299 LONGLONG nrows, LONGLONG naxis1, LONGLONG naxis2, LONGLONG naxis3,
1300 double *array, int *anynul, int *status);
1301
1302int CFITS_API ffgsvb(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1303 long *trc, long *inc, unsigned char nulval, unsigned char *array,
1304 int *anynul, int *status);
1305int CFITS_API ffgsvsb(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1306 long *trc, long *inc, signed char nulval, signed char *array,
1307 int *anynul, int *status);
1308int CFITS_API ffgsvui(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1309 long *trc, long *inc, unsigned short nulval, unsigned short *array,
1310 int *anynul, int *status);
1311int CFITS_API ffgsvi(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1312 long *trc, long *inc, short nulval, short *array, int *anynul, int *status);
1313int CFITS_API ffgsvuj(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1314 long *trc, long *inc, unsigned long nulval, unsigned long *array,
1315 int *anynul, int *status);
1316int CFITS_API ffgsvj(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1317 long *trc, long *inc, long nulval, long *array, int *anynul, int *status);
1318int CFITS_API ffgsvjj(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1319 long *trc, long *inc, LONGLONG nulval, LONGLONG *array, int *anynul,
1320 int *status);
1321int CFITS_API ffgsvuk(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1322 long *trc, long *inc, unsigned int nulval, unsigned int *array,
1323 int *anynul, int *status);
1324int CFITS_API ffgsvk(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1325 long *trc, long *inc, int nulval, int *array, int *anynul, int *status);
1326int CFITS_API ffgsve(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1327 long *trc, long *inc, float nulval, float *array, int *anynul, int *status);
1328int CFITS_API ffgsvd(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1329 long *trc, long *inc, double nulval, double *array, int *anynul,
1330 int *status);
1331
1332int CFITS_API ffgsfb(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1333 long *trc, long *inc, unsigned char *array, char *flagval,
1334 int *anynul, int *status);
1335int CFITS_API ffgsfsb(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1336 long *trc, long *inc, signed char *array, char *flagval,
1337 int *anynul, int *status);
1338int CFITS_API ffgsfui(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1339 long *trc, long *inc, unsigned short *array, char *flagval, int *anynul,
1340 int *status);
1341int CFITS_API ffgsfi(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1342 long *trc, long *inc, short *array, char *flagval, int *anynul, int *status);
1343int CFITS_API ffgsfuj(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1344 long *trc, long *inc, unsigned long *array, char *flagval, int *anynul,
1345 int *status);
1346int CFITS_API ffgsfj(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1347 long *trc, long *inc, long *array, char *flagval, int *anynul, int *status);
1348int CFITS_API ffgsfjj(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1349 long *trc, long *inc, LONGLONG *array, char *flagval, int *anynul,
1350 int *status);
1351int CFITS_API ffgsfuk(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1352 long *trc, long *inc, unsigned int *array, char *flagval, int *anynul,
1353 int *status);
1354int CFITS_API ffgsfk(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1355 long *trc, long *inc, int *array, char *flagval, int *anynul, int *status);
1356int CFITS_API ffgsfe(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1357 long *trc, long *inc, float *array, char *flagval, int *anynul, int *status);
1358int CFITS_API ffgsfd(fitsfile *fptr, int colnum, int naxis, long *naxes, long *blc,
1359 long *trc, long *inc, double *array, char *flagval, int *anynul,
1360 int *status);
1361
1362int CFITS_API ffggpb(fitsfile *fptr, long group, long firstelem, long nelem,
1363 unsigned char *array, int *status);
1364int CFITS_API ffggpsb(fitsfile *fptr, long group, long firstelem, long nelem,
1365 signed char *array, int *status);
1366int CFITS_API ffggpui(fitsfile *fptr, long group, long firstelem, long nelem,
1367 unsigned short *array, int *status);
1368int CFITS_API ffggpi(fitsfile *fptr, long group, long firstelem, long nelem,
1369 short *array, int *status);
1370int CFITS_API ffggpuj(fitsfile *fptr, long group, long firstelem, long nelem,
1371 unsigned long *array, int *status);
1372int CFITS_API ffggpj(fitsfile *fptr, long group, long firstelem, long nelem,
1373 long *array, int *status);
1374int CFITS_API ffggpjj(fitsfile *fptr, long group, long firstelem, long nelem,
1375 LONGLONG *array, int *status);
1376int CFITS_API ffggpuk(fitsfile *fptr, long group, long firstelem, long nelem,
1377 unsigned int *array, int *status);
1378int CFITS_API ffggpk(fitsfile *fptr, long group, long firstelem, long nelem,
1379 int *array, int *status);
1380int CFITS_API ffggpe(fitsfile *fptr, long group, long firstelem, long nelem,
1381 float *array, int *status);
1382int CFITS_API ffggpd(fitsfile *fptr, long group, long firstelem, long nelem,
1383 double *array, int *status);
1384
1385/*--------------------- read column elements -------------*/
1386int CFITS_API ffgcv( fitsfile *fptr, int datatype, int colnum, LONGLONG firstrow,
1387 LONGLONG firstelem, LONGLONG nelem, void *nulval, void *array, int *anynul,
1388 int *status);
1389int CFITS_API ffgcf( fitsfile *fptr, int datatype, int colnum, LONGLONG firstrow,
1390 LONGLONG firstelem, LONGLONG nelem, void *array, char *nullarray,
1391 int *anynul, int *status);
1392int CFITS_API ffgcvs(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1393 LONGLONG nelem, char *nulval, char **array, int *anynul, int *status);
1394int CFITS_API ffgcl (fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1395 LONGLONG nelem, char *array, int *status);
1396int CFITS_API ffgcvl (fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1397 LONGLONG nelem, char nulval, char *array, int *anynul, int *status);
1398int CFITS_API ffgcvb(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1399 LONGLONG nelem, unsigned char nulval, unsigned char *array,
1400 int *anynul, int *status);
1401int CFITS_API ffgcvsb(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1402 LONGLONG nelem, signed char nulval, signed char *array,
1403 int *anynul, int *status);
1404int CFITS_API ffgcvui(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1405 LONGLONG nelem, unsigned short nulval, unsigned short *array,
1406 int *anynul, int *status);
1407int CFITS_API ffgcvi(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1408 LONGLONG nelem, short nulval, short *array, int *anynul, int *status);
1409int CFITS_API ffgcvuj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1410 LONGLONG nelem, unsigned long nulval, unsigned long *array, int *anynul,
1411 int *status);
1412int CFITS_API ffgcvj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1413 LONGLONG nelem, long nulval, long *array, int *anynul, int *status);
1414int CFITS_API ffgcvjj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1415 LONGLONG nelem, LONGLONG nulval, LONGLONG *array, int *anynul,
1416 int *status);
1417int CFITS_API ffgcvuk(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1418 LONGLONG nelem, unsigned int nulval, unsigned int *array, int *anynul,
1419 int *status);
1420int CFITS_API ffgcvk(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1421 LONGLONG nelem, int nulval, int *array, int *anynul, int *status);
1422int CFITS_API ffgcve(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1423 LONGLONG nelem, float nulval, float *array, int *anynul, int *status);
1424int CFITS_API ffgcvd(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1425 LONGLONG nelem, double nulval, double *array, int *anynul, int *status);
1426int CFITS_API ffgcvc(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1427 LONGLONG nelem, float nulval, float *array, int *anynul, int *status);
1428int CFITS_API ffgcvm(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1429 LONGLONG nelem, double nulval, double *array, int *anynul, int *status);
1430
1431int CFITS_API ffgcx(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstbit,
1432 LONGLONG nbits, char *larray, int *status);
1433int CFITS_API ffgcxui(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG nrows,
1434 long firstbit, int nbits, unsigned short *array, int *status);
1435int CFITS_API ffgcxuk(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG nrows,
1436 long firstbit, int nbits, unsigned int *array, int *status);
1437
1438int CFITS_API ffgcfs(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1439 LONGLONG nelem, char **array, char *nularray, int *anynul, int *status);
1440int CFITS_API ffgcfl(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1441 LONGLONG nelem, char *array, char *nularray, int *anynul, int *status);
1442int CFITS_API ffgcfb(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1443 LONGLONG nelem, unsigned char *array, char *nularray, int *anynul, int *status);
1444int CFITS_API ffgcfsb(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1445 LONGLONG nelem, signed char *array, char *nularray, int *anynul, int *status);
1446int CFITS_API ffgcfui(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1447 LONGLONG nelem, unsigned short *array, char *nularray, int *anynul,
1448 int *status);
1449int CFITS_API ffgcfi(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1450 LONGLONG nelem, short *array, char *nularray, int *anynul, int *status);
1451int CFITS_API ffgcfuj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1452 LONGLONG nelem, unsigned long *array, char *nularray, int *anynul,
1453 int *status);
1454int CFITS_API ffgcfj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1455 LONGLONG nelem, long *array, char *nularray, int *anynul, int *status);
1456int CFITS_API ffgcfjj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1457 LONGLONG nelem, LONGLONG *array, char *nularray, int *anynul, int *status);
1458int CFITS_API ffgcfuk(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1459 LONGLONG nelem, unsigned int *array, char *nularray, int *anynul,
1460 int *status);
1461int CFITS_API ffgcfk(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1462 LONGLONG nelem, int *array, char *nularray, int *anynul, int *status);
1463int CFITS_API ffgcfe(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1464 LONGLONG nelem, float *array, char *nularray, int *anynul, int *status);
1465int CFITS_API ffgcfd(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1466 LONGLONG nelem, double *array, char *nularray, int *anynul, int *status);
1467int CFITS_API ffgcfc(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1468 LONGLONG nelem, float *array, char *nularray, int *anynul, int *status);
1469int CFITS_API ffgcfm(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1470 LONGLONG nelem, double *array, char *nularray, int *anynul, int *status);
1471
1472int CFITS_API ffgdes(fitsfile *fptr, int colnum, LONGLONG rownum, long *length,
1473 long *heapaddr, int *status);
1474int CFITS_API ffgdesll(fitsfile *fptr, int colnum, LONGLONG rownum, LONGLONG *length,
1475 LONGLONG *heapaddr, int *status);
1476int CFITS_API ffgdess(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG nrows, long *length,
1477 long *heapaddr, int *status);
1478int CFITS_API ffgdessll(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG nrows, LONGLONG *length,
1479 LONGLONG *heapaddr, int *status);
1480int CFITS_API ffpdes(fitsfile *fptr, int colnum, LONGLONG rownum, LONGLONG length,
1481 LONGLONG heapaddr, int *status);
1482int CFITS_API fftheap(fitsfile *fptr, LONGLONG *heapsize, LONGLONG *unused, LONGLONG *overlap,
1483 int *valid, int *status);
1484int CFITS_API ffcmph(fitsfile *fptr, int *status);
1485
1486int CFITS_API ffgtbb(fitsfile *fptr, LONGLONG firstrow, LONGLONG firstchar, LONGLONG nchars,
1487 unsigned char *values, int *status);
1488
1489int CFITS_API ffgextn(fitsfile *fptr, LONGLONG offset, LONGLONG nelem, void *array, int *status);
1490int CFITS_API ffpextn(fitsfile *fptr, LONGLONG offset, LONGLONG nelem, void *array, int *status);
1491
1492/*------------ write primary array or image elements -------------*/
1493int CFITS_API ffppx(fitsfile *fptr, int datatype, long *firstpix, LONGLONG nelem,
1494 void *array, int *status);
1495int CFITS_API ffppxll(fitsfile *fptr, int datatype, LONGLONG *firstpix, LONGLONG nelem,
1496 void *array, int *status);
1497int CFITS_API ffppxn(fitsfile *fptr, int datatype, long *firstpix, LONGLONG nelem,
1498 void *array, void *nulval, int *status);
1499int CFITS_API ffppxnll(fitsfile *fptr, int datatype, LONGLONG *firstpix, LONGLONG nelem,
1500 void *array, void *nulval, int *status);
1501int CFITS_API ffppr(fitsfile *fptr, int datatype, LONGLONG firstelem,
1502 LONGLONG nelem, void *array, int *status);
1503int CFITS_API ffpprb(fitsfile *fptr, long group, LONGLONG firstelem,
1504 LONGLONG nelem, unsigned char *array, int *status);
1505int CFITS_API ffpprsb(fitsfile *fptr, long group, LONGLONG firstelem,
1506 LONGLONG nelem, signed char *array, int *status);
1507int CFITS_API ffpprui(fitsfile *fptr, long group, LONGLONG firstelem,
1508 LONGLONG nelem, unsigned short *array, int *status);
1509int CFITS_API ffppri(fitsfile *fptr, long group, LONGLONG firstelem,
1510 LONGLONG nelem, short *array, int *status);
1511int CFITS_API ffppruj(fitsfile *fptr, long group, LONGLONG firstelem,
1512 LONGLONG nelem, unsigned long *array, int *status);
1513int CFITS_API ffpprj(fitsfile *fptr, long group, LONGLONG firstelem,
1514 LONGLONG nelem, long *array, int *status);
1515int CFITS_API ffppruk(fitsfile *fptr, long group, LONGLONG firstelem,
1516 LONGLONG nelem, unsigned int *array, int *status);
1517int CFITS_API ffpprk(fitsfile *fptr, long group, LONGLONG firstelem,
1518 LONGLONG nelem, int *array, int *status);
1519int CFITS_API ffppre(fitsfile *fptr, long group, LONGLONG firstelem,
1520 LONGLONG nelem, float *array, int *status);
1521int CFITS_API ffpprd(fitsfile *fptr, long group, LONGLONG firstelem,
1522 LONGLONG nelem, double *array, int *status);
1523int CFITS_API ffpprjj(fitsfile *fptr, long group, LONGLONG firstelem,
1524 LONGLONG nelem, LONGLONG *array, int *status);
1525
1526int CFITS_API ffppru(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1527 int *status);
1528int CFITS_API ffpprn(fitsfile *fptr, LONGLONG firstelem, LONGLONG nelem, int *status);
1529
1530int CFITS_API ffppn(fitsfile *fptr, int datatype, LONGLONG firstelem, LONGLONG nelem,
1531 void *array, void *nulval, int *status);
1532int CFITS_API ffppnb(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1533 unsigned char *array, unsigned char nulval, int *status);
1534int CFITS_API ffppnsb(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1535 signed char *array, signed char nulval, int *status);
1536int CFITS_API ffppnui(fitsfile *fptr, long group, LONGLONG firstelem,
1537 LONGLONG nelem, unsigned short *array, unsigned short nulval,
1538 int *status);
1539int CFITS_API ffppni(fitsfile *fptr, long group, LONGLONG firstelem,
1540 LONGLONG nelem, short *array, short nulval, int *status);
1541int CFITS_API ffppnj(fitsfile *fptr, long group, LONGLONG firstelem,
1542 LONGLONG nelem, long *array, long nulval, int *status);
1543int CFITS_API ffppnuj(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1544 unsigned long *array, unsigned long nulval, int *status);
1545int CFITS_API ffppnuk(fitsfile *fptr, long group, LONGLONG firstelem, LONGLONG nelem,
1546 unsigned int *array, unsigned int nulval, int *status);
1547int CFITS_API ffppnk(fitsfile *fptr, long group, LONGLONG firstelem,
1548 LONGLONG nelem, int *array, int nulval, int *status);
1549int CFITS_API ffppne(fitsfile *fptr, long group, LONGLONG firstelem,
1550 LONGLONG nelem, float *array, float nulval, int *status);
1551int CFITS_API ffppnd(fitsfile *fptr, long group, LONGLONG firstelem,
1552 LONGLONG nelem, double *array, double nulval, int *status);
1553int CFITS_API ffppnjj(fitsfile *fptr, long group, LONGLONG firstelem,
1554 LONGLONG nelem, LONGLONG *array, LONGLONG nulval, int *status);
1555
1556int CFITS_API ffp2db(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1557 LONGLONG naxis2, unsigned char *array, int *status);
1558int CFITS_API ffp2dsb(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1559 LONGLONG naxis2, signed char *array, int *status);
1560int CFITS_API ffp2dui(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1561 LONGLONG naxis2, unsigned short *array, int *status);
1562int CFITS_API ffp2di(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1563 LONGLONG naxis2, short *array, int *status);
1564int CFITS_API ffp2duj(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1565 LONGLONG naxis2, unsigned long *array, int *status);
1566int CFITS_API ffp2dj(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1567 LONGLONG naxis2, long *array, int *status);
1568int CFITS_API ffp2duk(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1569 LONGLONG naxis2, unsigned int *array, int *status);
1570int CFITS_API ffp2dk(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1571 LONGLONG naxis2, int *array, int *status);
1572int CFITS_API ffp2de(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1573 LONGLONG naxis2, float *array, int *status);
1574int CFITS_API ffp2dd(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1575 LONGLONG naxis2, double *array, int *status);
1576int CFITS_API ffp2djj(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG naxis1,
1577 LONGLONG naxis2, LONGLONG *array, int *status);
1578
1579int CFITS_API ffp3db(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1580 LONGLONG naxis2, LONGLONG naxis3, unsigned char *array, int *status);
1581int CFITS_API ffp3dsb(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1582 LONGLONG naxis2, LONGLONG naxis3, signed char *array, int *status);
1583int CFITS_API ffp3dui(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1584 LONGLONG naxis2, LONGLONG naxis3, unsigned short *array, int *status);
1585int CFITS_API ffp3di(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1586 LONGLONG naxis2, LONGLONG naxis3, short *array, int *status);
1587int CFITS_API ffp3duj(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1588 LONGLONG naxis2, LONGLONG naxis3, unsigned long *array, int *status);
1589int CFITS_API ffp3dj(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1590 LONGLONG naxis2, LONGLONG naxis3, long *array, int *status);
1591int CFITS_API ffp3duk(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1592 LONGLONG naxis2, LONGLONG naxis3, unsigned int *array, int *status);
1593int CFITS_API ffp3dk(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1594 LONGLONG naxis2, LONGLONG naxis3, int *array, int *status);
1595int CFITS_API ffp3de(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1596 LONGLONG naxis2, LONGLONG naxis3, float *array, int *status);
1597int CFITS_API ffp3dd(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1598 LONGLONG naxis2, LONGLONG naxis3, double *array, int *status);
1599int CFITS_API ffp3djj(fitsfile *fptr, long group, LONGLONG ncols, LONGLONG nrows, LONGLONG naxis1,
1600 LONGLONG naxis2, LONGLONG naxis3, LONGLONG *array, int *status);
1601
1602int CFITS_API ffpss(fitsfile *fptr, int datatype,
1603 long *fpixel, long *lpixel, void *array, int *status);
1604int CFITS_API ffpssb(fitsfile *fptr, long group, long naxis, long *naxes,
1605 long *fpixel, long *lpixel, unsigned char *array, int *status);
1606int CFITS_API ffpsssb(fitsfile *fptr, long group, long naxis, long *naxes,
1607 long *fpixel, long *lpixel, signed char *array, int *status);
1608int CFITS_API ffpssui(fitsfile *fptr, long group, long naxis, long *naxes,
1609 long *fpixel, long *lpixel, unsigned short *array, int *status);
1610int CFITS_API ffpssi(fitsfile *fptr, long group, long naxis, long *naxes,
1611 long *fpixel, long *lpixel, short *array, int *status);
1612int CFITS_API ffpssuj(fitsfile *fptr, long group, long naxis, long *naxes,
1613 long *fpixel, long *lpixel, unsigned long *array, int *status);
1614int CFITS_API ffpssj(fitsfile *fptr, long group, long naxis, long *naxes,
1615 long *fpixel, long *lpixel, long *array, int *status);
1616int CFITS_API ffpssuk(fitsfile *fptr, long group, long naxis, long *naxes,
1617 long *fpixel, long *lpixel, unsigned int *array, int *status);
1618int CFITS_API ffpssk(fitsfile *fptr, long group, long naxis, long *naxes,
1619 long *fpixel, long *lpixel, int *array, int *status);
1620int CFITS_API ffpsse(fitsfile *fptr, long group, long naxis, long *naxes,
1621 long *fpixel, long *lpixel, float *array, int *status);
1622int CFITS_API ffpssd(fitsfile *fptr, long group, long naxis, long *naxes,
1623 long *fpixel, long *lpixel, double *array, int *status);
1624int CFITS_API ffpssjj(fitsfile *fptr, long group, long naxis, long *naxes,
1625 long *fpixel, long *lpixel, LONGLONG *array, int *status);
1626
1627int CFITS_API ffpgpb(fitsfile *fptr, long group, long firstelem,
1628 long nelem, unsigned char *array, int *status);
1629int CFITS_API ffpgpsb(fitsfile *fptr, long group, long firstelem,
1630 long nelem, signed char *array, int *status);
1631int CFITS_API ffpgpui(fitsfile *fptr, long group, long firstelem,
1632 long nelem, unsigned short *array, int *status);
1633int CFITS_API ffpgpi(fitsfile *fptr, long group, long firstelem,
1634 long nelem, short *array, int *status);
1635int CFITS_API ffpgpuj(fitsfile *fptr, long group, long firstelem,
1636 long nelem, unsigned long *array, int *status);
1637int CFITS_API ffpgpj(fitsfile *fptr, long group, long firstelem,
1638 long nelem, long *array, int *status);
1639int CFITS_API ffpgpuk(fitsfile *fptr, long group, long firstelem,
1640 long nelem, unsigned int *array, int *status);
1641int CFITS_API ffpgpk(fitsfile *fptr, long group, long firstelem,
1642 long nelem, int *array, int *status);
1643int CFITS_API ffpgpe(fitsfile *fptr, long group, long firstelem,
1644 long nelem, float *array, int *status);
1645int CFITS_API ffpgpd(fitsfile *fptr, long group, long firstelem,
1646 long nelem, double *array, int *status);
1647int CFITS_API ffpgpjj(fitsfile *fptr, long group, long firstelem,
1648 long nelem, LONGLONG *array, int *status);
1649
1650/*--------------------- iterator functions -------------*/
1651int CFITS_API fits_iter_set_by_name(iteratorCol *col, fitsfile *fptr, char *colname,
1652 int datatype, int iotype);
1653int CFITS_API fits_iter_set_by_num(iteratorCol *col, fitsfile *fptr, int colnum,
1654 int datatype, int iotype);
1655int CFITS_API fits_iter_set_file(iteratorCol *col, fitsfile *fptr);
1656int CFITS_API fits_iter_set_colname(iteratorCol *col, char *colname);
1657int CFITS_API fits_iter_set_colnum(iteratorCol *col, int colnum);
1658int CFITS_API fits_iter_set_datatype(iteratorCol *col, int datatype);
1659int CFITS_API fits_iter_set_iotype(iteratorCol *col, int iotype);
1660
1661CFITS_API fitsfile * fits_iter_get_file(iteratorCol *col);
1662char CFITS_API * fits_iter_get_colname(iteratorCol *col);
1663int CFITS_API fits_iter_get_colnum(iteratorCol *col);
1664int CFITS_API fits_iter_get_datatype(iteratorCol *col);
1665int CFITS_API fits_iter_get_iotype(iteratorCol *col);
1666void CFITS_API *fits_iter_get_array(iteratorCol *col);
1667long CFITS_API fits_iter_get_tlmin(iteratorCol *col);
1668long CFITS_API fits_iter_get_tlmax(iteratorCol *col);
1669long CFITS_API fits_iter_get_repeat(iteratorCol *col);
1670char CFITS_API *fits_iter_get_tunit(iteratorCol *col);
1671char CFITS_API *fits_iter_get_tdisp(iteratorCol *col);
1672
1673int CFITS_API ffiter(int ncols, iteratorCol *data, long offset, long nPerLoop,
1674 int (*workFn)( long totaln, long offset, long firstn,
1675 long nvalues, int narrays, iteratorCol *data, void *userPointer),
1676 void *userPointer, int *status);
1677
1678/*--------------------- write column elements -------------*/
1679int CFITS_API ffpcl(fitsfile *fptr, int datatype, int colnum, LONGLONG firstrow,
1680 LONGLONG firstelem, LONGLONG nelem, void *array, int *status);
1681int CFITS_API ffpcls(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1682 LONGLONG nelem, char **array, int *status);
1683int CFITS_API ffpcll(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1684 LONGLONG nelem, char *array, int *status);
1685int CFITS_API ffpclb(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1686 LONGLONG nelem, unsigned char *array, int *status);
1687int CFITS_API ffpclsb(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1688 LONGLONG nelem, signed char *array, int *status);
1689int CFITS_API ffpclui(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1690 LONGLONG nelem, unsigned short *array, int *status);
1691int CFITS_API ffpcli(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1692 LONGLONG nelem, short *array, int *status);
1693int CFITS_API ffpcluj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1694 LONGLONG nelem, unsigned long *array, int *status);
1695int CFITS_API ffpclj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1696 LONGLONG nelem, long *array, int *status);
1697int CFITS_API ffpcluk(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1698 LONGLONG nelem, unsigned int *array, int *status);
1699int CFITS_API ffpclk(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1700 LONGLONG nelem, int *array, int *status);
1701int CFITS_API ffpcle(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1702 LONGLONG nelem, float *array, int *status);
1703int CFITS_API ffpcld(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1704 LONGLONG nelem, double *array, int *status);
1705int CFITS_API ffpclc(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1706 LONGLONG nelem, float *array, int *status);
1707int CFITS_API ffpclm(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1708 LONGLONG nelem, double *array, int *status);
1709int CFITS_API ffpclu(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1710 LONGLONG nelem, int *status);
1711int CFITS_API ffprwu(fitsfile *fptr, LONGLONG firstrow, LONGLONG nrows, int *status);
1712int CFITS_API ffpcljj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1713 LONGLONG nelem, LONGLONG *array, int *status);
1714int CFITS_API ffpclx(fitsfile *fptr, int colnum, LONGLONG frow, long fbit, long nbit,
1715 char *larray, int *status);
1716
1717int CFITS_API ffpcn(fitsfile *fptr, int datatype, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1718 LONGLONG nelem, void *array, void *nulval, int *status);
1719int CFITS_API ffpcns( fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1720 LONGLONG nelem, char **array, char *nulvalue, int *status);
1721int CFITS_API ffpcnl( fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1722 LONGLONG nelem, char *array, char nulvalue, int *status);
1723int CFITS_API ffpcnb(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1724 LONGLONG nelem, unsigned char *array, unsigned char nulvalue,
1725 int *status);
1726int CFITS_API ffpcnsb(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1727 LONGLONG nelem, signed char *array, signed char nulvalue,
1728 int *status);
1729int CFITS_API ffpcnui(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1730 LONGLONG nelem, unsigned short *array, unsigned short nulvalue,
1731 int *status);
1732int CFITS_API ffpcni(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1733 LONGLONG nelem, short *array, short nulvalue, int *status);
1734int CFITS_API ffpcnuj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1735 LONGLONG nelem, unsigned long *array, unsigned long nulvalue,
1736 int *status);
1737int CFITS_API ffpcnj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1738 LONGLONG nelem, long *array, long nulvalue, int *status);
1739int CFITS_API ffpcnuk(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1740 LONGLONG nelem, unsigned int *array, unsigned int nulvalue,
1741 int *status);
1742int CFITS_API ffpcnk(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1743 LONGLONG nelem, int *array, int nulvalue, int *status);
1744int CFITS_API ffpcne(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1745 LONGLONG nelem, float *array, float nulvalue, int *status);
1746int CFITS_API ffpcnd(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1747 LONGLONG nelem, double *array, double nulvalue, int *status);
1748int CFITS_API ffpcnjj(fitsfile *fptr, int colnum, LONGLONG firstrow, LONGLONG firstelem,
1749 LONGLONG nelem, LONGLONG *array, LONGLONG nulvalue, int *status);
1750int CFITS_API ffptbb(fitsfile *fptr, LONGLONG firstrow, LONGLONG firstchar, LONGLONG nchars,
1751 unsigned char *values, int *status);
1752
1753int CFITS_API ffirow(fitsfile *fptr, LONGLONG firstrow, LONGLONG nrows, int *status);
1754int CFITS_API ffdrow(fitsfile *fptr, LONGLONG firstrow, LONGLONG nrows, int *status);
1755int CFITS_API ffdrrg(fitsfile *fptr, char *ranges, int *status);
1756int CFITS_API ffdrws(fitsfile *fptr, long *rownum, long nrows, int *status);
1757int CFITS_API ffdrwsll(fitsfile *fptr, LONGLONG *rownum, LONGLONG nrows, int *status);
1758int CFITS_API fficol(fitsfile *fptr, int numcol, char *ttype, char *tform, int *status);
1759int CFITS_API fficls(fitsfile *fptr, int firstcol, int ncols, char **ttype,
1760 char **tform, int *status);
1761int CFITS_API ffmvec(fitsfile *fptr, int colnum, LONGLONG newveclen, int *status);
1762int CFITS_API ffdcol(fitsfile *fptr, int numcol, int *status);
1763int CFITS_API ffcpcl(fitsfile *infptr, fitsfile *outfptr, int incol, int outcol,
1764 int create_col, int *status);
1765int CFITS_API ffcprw(fitsfile *infptr, fitsfile *outfptr, LONGLONG firstrow,
1766 LONGLONG nrows, int *status);
1767
1768/*--------------------- WCS Utilities ------------------*/
1769int CFITS_API ffgics(fitsfile *fptr, double *xrval, double *yrval, double *xrpix,
1770 double *yrpix, double *xinc, double *yinc, double *rot,
1771 char *type, int *status);
1772int CFITS_API ffgicsa(fitsfile *fptr, char version, double *xrval, double *yrval, double *xrpix,
1773 double *yrpix, double *xinc, double *yinc, double *rot,
1774 char *type, int *status);
1775int CFITS_API ffgtcs(fitsfile *fptr, int xcol, int ycol, double *xrval,
1776 double *yrval, double *xrpix, double *yrpix, double *xinc,
1777 double *yinc, double *rot, char *type, int *status);
1778int CFITS_API ffwldp(double xpix, double ypix, double xref, double yref,
1779 double xrefpix, double yrefpix, double xinc, double yinc,
1780 double rot, char *type, double *xpos, double *ypos, int *status);
1781int CFITS_API ffxypx(double xpos, double ypos, double xref, double yref,
1782 double xrefpix, double yrefpix, double xinc, double yinc,
1783 double rot, char *type, double *xpix, double *ypix, int *status);
1784
1785/* WCS support routines (provide interface to Doug Mink's WCS library */
1786int CFITS_API ffgiwcs(fitsfile *fptr, char **header, int *status);
1787int CFITS_API ffgtwcs(fitsfile *fptr, int xcol, int ycol, char **header, int *status);
1788
1789/*--------------------- lexical parsing routines ------------------*/
1790int CFITS_API fftexp( fitsfile *fptr, char *expr, int maxdim,
1791 int *datatype, long *nelem, int *naxis,
1792 long *naxes, int *status );
1793
1794int CFITS_API fffrow( fitsfile *infptr, char *expr,
1795 long firstrow, long nrows,
1796 long *n_good_rows, char *row_status, int *status);
1797
1798int CFITS_API ffffrw( fitsfile *fptr, char *expr, long *rownum, int *status);
1799
1800int CFITS_API fffrwc( fitsfile *fptr, char *expr, char *timeCol,
1801 char *parCol, char *valCol, long ntimes,
1802 double *times, char *time_status, int *status );
1803
1804int CFITS_API ffsrow( fitsfile *infptr, fitsfile *outfptr, char *expr,
1805 int *status);
1806
1807int CFITS_API ffcrow( fitsfile *fptr, int datatype, char *expr,
1808 long firstrow, long nelements, void *nulval,
1809 void *array, int *anynul, int *status );
1810
1811int CFITS_API ffcalc_rng( fitsfile *infptr, char *expr, fitsfile *outfptr,
1812 char *parName, char *parInfo, int nRngs,
1813 long *start, long *end, int *status );
1814
1815int CFITS_API ffcalc( fitsfile *infptr, char *expr, fitsfile *outfptr,
1816 char *parName, char *parInfo, int *status );
1817
1818 /* ffhist is not really intended as a user-callable routine */
1819 /* but it may be useful for some specialized applications */
1820 /* ffhist2 is a newer version which is strongly recommended instead of ffhist */
1821
1822int CFITS_API ffhist(fitsfile **fptr, char *outfile, int imagetype, int naxis,
1823 char colname[4][FLEN_VALUE],
1824 double *minin, double *maxin, double *binsizein,
1825 char minname[4][FLEN_VALUE], char maxname[4][FLEN_VALUE],
1826 char binname[4][FLEN_VALUE],
1827 double weightin, char wtcol[FLEN_VALUE],
1828 int recip, char *rowselect, int *status);
1829int CFITS_API ffhist2(fitsfile **fptr, char *outfile, int imagetype, int naxis,
1830 char colname[4][FLEN_VALUE],
1831 double *minin, double *maxin, double *binsizein,
1832 char minname[4][FLEN_VALUE], char maxname[4][FLEN_VALUE],
1833 char binname[4][FLEN_VALUE],
1834 double weightin, char wtcol[FLEN_VALUE],
1835 int recip, char *rowselect, int *status);
1836
1837int CFITS_API fits_select_image_section(fitsfile **fptr, char *outfile,
1838 char *imagesection, int *status);
1839int CFITS_API fits_copy_image_section(fitsfile *infptr, fitsfile *outfile,
1840 char *imagesection, int *status);
1841
1842int CFITS_API fits_calc_binning(fitsfile *fptr, int naxis, char colname[4][FLEN_VALUE],
1843 double *minin, double *maxin, double *binsizein,
1844 char minname[4][FLEN_VALUE], char maxname[4][FLEN_VALUE],
1845 char binname[4][FLEN_VALUE], int *colnum, long *haxes, float *amin,
1846 float *amax, float *binsize, int *status);
1847
1848int CFITS_API fits_write_keys_histo(fitsfile *fptr, fitsfile *histptr,
1849 int naxis, int *colnum, int *status);
1850int CFITS_API fits_rebin_wcs( fitsfile *fptr, int naxis, float *amin, float *binsize,
1851 int *status);
1852int CFITS_API fits_make_hist(fitsfile *fptr, fitsfile *histptr, int bitpix,int naxis,
1853 long *naxes, int *colnum, float *amin, float *amax, float *binsize,
1854 float weight, int wtcolnum, int recip, char *selectrow, int *status);
1855
1856typedef struct
1857{
1858 /* input(s) */
1859 int count;
1860 char ** path;
1861 char ** tag;
1862 fitsfile ** ifptr;
1863
1864 char * expression;
1865
1866 /* output control */
1867 int bitpix;
1868 long blank;
1869 fitsfile * ofptr;
1870 char keyword[FLEN_KEYWORD];
1871 char comment[FLEN_COMMENT];
1872} PixelFilter;
1873
1874
1875int CFITS_API fits_pixel_filter (PixelFilter * filter, int * status);
1876
1877
1878/*--------------------- grouping routines ------------------*/
1879
1880int CFITS_API ffgtcr(fitsfile *fptr, char *grpname, int grouptype, int *status);
1881int CFITS_API ffgtis(fitsfile *fptr, char *grpname, int grouptype, int *status);
1882int CFITS_API ffgtch(fitsfile *gfptr, int grouptype, int *status);
1883int CFITS_API ffgtrm(fitsfile *gfptr, int rmopt, int *status);
1884int CFITS_API ffgtcp(fitsfile *infptr, fitsfile *outfptr, int cpopt, int *status);
1885int CFITS_API ffgtmg(fitsfile *infptr, fitsfile *outfptr, int mgopt, int *status);
1886int CFITS_API ffgtcm(fitsfile *gfptr, int cmopt, int *status);
1887int CFITS_API ffgtvf(fitsfile *gfptr, long *firstfailed, int *status);
1888int CFITS_API ffgtop(fitsfile *mfptr,int group,fitsfile **gfptr,int *status);
1889int CFITS_API ffgtam(fitsfile *gfptr, fitsfile *mfptr, int hdupos, int *status);
1890int CFITS_API ffgtnm(fitsfile *gfptr, long *nmembers, int *status);
1891int CFITS_API ffgmng(fitsfile *mfptr, long *nmembers, int *status);
1892int CFITS_API ffgmop(fitsfile *gfptr, long member, fitsfile **mfptr, int *status);
1893int CFITS_API ffgmcp(fitsfile *gfptr, fitsfile *mfptr, long member, int cpopt,
1894 int *status);
1895int CFITS_API ffgmtf(fitsfile *infptr, fitsfile *outfptr, long member, int tfopt,
1896 int *status);
1897int CFITS_API ffgmrm(fitsfile *fptr, long member, int rmopt, int *status);
1898
1899/*--------------------- group template parser routines ------------------*/
1900
1901int CFITS_API fits_execute_template(fitsfile *ff, char *ngp_template, int *status);
1902
1903int CFITS_API fits_img_stats_short(short *array,long nx, long ny, int nullcheck,
1904 short nullvalue,long *ngoodpix, short *minvalue, short *maxvalue, double *mean,
1905 double *sigma, double *noise1, double *noise2, double *noise3, double *noise5, int *status);
1906int CFITS_API fits_img_stats_int(int *array,long nx, long ny, int nullcheck,
1907 int nullvalue,long *ngoodpix, int *minvalue, int *maxvalue, double *mean,
1908 double *sigma, double *noise1, double *noise2, double *noise3, double *noise5, int *status);
1909int CFITS_API fits_img_stats_float(float *array, long nx, long ny, int nullcheck,
1910 float nullvalue,long *ngoodpix, float *minvalue, float *maxvalue, double *mean,
1911 double *sigma, double *noise1, double *noise2, double *noise3, double *noise5, int *status);
1912
1913/*--------------------- image compression routines ------------------*/
1914
1915int CFITS_API fits_set_compression_type(fitsfile *fptr, int ctype, int *status);
1916int CFITS_API fits_set_tile_dim(fitsfile *fptr, int ndim, long *dims, int *status);
1917int CFITS_API fits_set_noise_bits(fitsfile *fptr, int noisebits, int *status);
1918int CFITS_API fits_set_quantize_level(fitsfile *fptr, float qlevel, int *status);
1919int CFITS_API fits_set_hcomp_scale(fitsfile *fptr, float scale, int *status);
1920int CFITS_API fits_set_hcomp_smooth(fitsfile *fptr, int smooth, int *status);
1921int CFITS_API fits_set_quantize_method(fitsfile *fptr, int method, int *status);
1922int CFITS_API fits_set_quantize_dither(fitsfile *fptr, int dither, int *status);
1923int CFITS_API fits_set_dither_seed(fitsfile *fptr, int seed, int *status);
1924int CFITS_API fits_set_dither_offset(fitsfile *fptr, int offset, int *status);
1925int CFITS_API fits_set_lossy_int(fitsfile *fptr, int lossy_int, int *status);
1926int CFITS_API fits_set_huge_hdu(fitsfile *fptr, int huge, int *status);
1927int CFITS_API fits_set_compression_pref(fitsfile *infptr, fitsfile *outfptr, int *status);
1928
1929int CFITS_API fits_get_compression_type(fitsfile *fptr, int *ctype, int *status);
1930int CFITS_API fits_get_tile_dim(fitsfile *fptr, int ndim, long *dims, int *status);
1931int CFITS_API fits_get_quantize_level(fitsfile *fptr, float *qlevel, int *status);
1932int CFITS_API fits_get_noise_bits(fitsfile *fptr, int *noisebits, int *status);
1933int CFITS_API fits_get_hcomp_scale(fitsfile *fptr, float *scale, int *status);
1934int CFITS_API fits_get_hcomp_smooth(fitsfile *fptr, int *smooth, int *status);
1935int CFITS_API fits_get_dither_seed(fitsfile *fptr, int *seed, int *status);
1936
1937int CFITS_API fits_img_compress(fitsfile *infptr, fitsfile *outfptr, int *status);
1938int CFITS_API fits_compress_img(fitsfile *infptr, fitsfile *outfptr, int compress_type,
1939 long *tilesize, int parm1, int parm2, int *status);
1940int CFITS_API fits_is_compressed_image(fitsfile *fptr, int *status);
1941int CFITS_API fits_is_reentrant(void);
1942int CFITS_API fits_decompress_img (fitsfile *infptr, fitsfile *outfptr, int *status);
1943int CFITS_API fits_img_decompress_header(fitsfile *infptr, fitsfile *outfptr, int *status);
1944int CFITS_API fits_img_decompress (fitsfile *infptr, fitsfile *outfptr, int *status);
1945
1946/* H-compress routines */
1947int CFITS_API fits_hcompress(int *a, int nx, int ny, int scale, char *output,
1948 long *nbytes, int *status);
1949int CFITS_API fits_hcompress64(LONGLONG *a, int nx, int ny, int scale, char *output,
1950 long *nbytes, int *status);
1951int CFITS_API fits_hdecompress(unsigned char *input, int smooth, int *a, int *nx,
1952 int *ny, int *scale, int *status);
1953int CFITS_API fits_hdecompress64(unsigned char *input, int smooth, LONGLONG *a, int *nx,
1954 int *ny, int *scale, int *status);
1955
1956int CFITS_API fits_compress_table (fitsfile *infptr, fitsfile *outfptr, int *status);
1957int CFITS_API fits_uncompress_table(fitsfile *infptr, fitsfile *outfptr, int *status);
1958
1959/* The following exclusion if __CINT__ is defined is needed for ROOT */
1960#ifndef __CINT__
1961#ifdef __cplusplus
1962}
1963#endif
1964#endif
1965
1966#endif
1967
1968