wcslib (8.2.2)

(root)/
include/
wcslib-8.2.2/
wcs.h
       1  /*============================================================================
       2    WCSLIB 8.2 - an implementation of the FITS WCS standard.
       3    Copyright (C) 1995-2023, Mark Calabretta
       4  
       5    This file is part of WCSLIB.
       6  
       7    WCSLIB is free software: you can redistribute it and/or modify it under the
       8    terms of the GNU Lesser General Public License as published by the Free
       9    Software Foundation, either version 3 of the License, or (at your option)
      10    any later version.
      11  
      12    WCSLIB is distributed in the hope that it will be useful, but WITHOUT ANY
      13    WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
      14    FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public License for
      15    more details.
      16  
      17    You should have received a copy of the GNU Lesser General Public License
      18    along with WCSLIB.  If not, see http://www.gnu.org/licenses.
      19  
      20    Author: Mark Calabretta, Australia Telescope National Facility, CSIRO.
      21    http://www.atnf.csiro.au/people/Mark.Calabretta
      22    $Id: wcs.h,v 8.2.1.1 2023/11/16 10:05:57 mcalabre Exp mcalabre $
      23  *=============================================================================
      24  *
      25  * WCSLIB 8.2 - C routines that implement the FITS World Coordinate System
      26  * (WCS) standard.  Refer to the README file provided with WCSLIB for an
      27  * overview of the library.
      28  *
      29  *
      30  * Summary of the wcs routines
      31  * ---------------------------
      32  * Routines in this suite implement the FITS World Coordinate System (WCS)
      33  * standard which defines methods to be used for computing world coordinates
      34  * from image pixel coordinates, and vice versa.  The standard, and proposed
      35  * extensions for handling distortions, are described in
      36  *
      37  =   "Representations of world coordinates in FITS",
      38  =   Greisen, E.W., & Calabretta, M.R. 2002, A&A, 395, 1061 (WCS Paper I)
      39  =
      40  =   "Representations of celestial coordinates in FITS",
      41  =   Calabretta, M.R., & Greisen, E.W. 2002, A&A, 395, 1077 (WCS Paper II)
      42  =
      43  =   "Representations of spectral coordinates in FITS",
      44  =   Greisen, E.W., Calabretta, M.R., Valdes, F.G., & Allen, S.L.
      45  =   2006, A&A, 446, 747 (WCS Paper III)
      46  =
      47  =   "Representations of distortions in FITS world coordinate systems",
      48  =   Calabretta, M.R. et al. (WCS Paper IV, draft dated 2004/04/22),
      49  =   available from http://www.atnf.csiro.au/people/Mark.Calabretta
      50  =
      51  =   "Mapping on the HEALPix grid",
      52  =   Calabretta, M.R., & Roukema, B.F. 2007, MNRAS, 381, 865 (WCS Paper V)
      53  =
      54  =   "Representing the 'Butterfly' Projection in FITS -- Projection Code XPH",
      55  =   Calabretta, M.R., & Lowe, S.R. 2013, PASA, 30, e050 (WCS Paper VI)
      56  =
      57  =   "Representations of time coordinates in FITS -
      58  =    Time and relative dimension in space",
      59  =   Rots, A.H., Bunclark, P.S., Calabretta, M.R., Allen, S.L.,
      60  =   Manchester, R.N., & Thompson, W.T. 2015, A&A, 574, A36 (WCS Paper VII)
      61  *
      62  * These routines are based on the wcsprm struct which contains all information
      63  * needed for the computations.  The struct contains some members that must be
      64  * set by the user, and others that are maintained by these routines, somewhat
      65  * like a C++ class but with no encapsulation.
      66  *
      67  * wcsnpv(), wcsnps(), wcsini(), wcsinit(), wcssub(), wcsfree(), and wcstrim(),
      68  * are provided to manage the wcsprm struct, wcssize() computes its total size
      69  * including allocated memory, and wcsprt() prints its contents.  Refer to the
      70  * description of the wcsprm struct for an explanation of the anticipated usage
      71  * of these routines.  wcscopy(), which does a deep copy of one wcsprm struct
      72  * to another, is defined as a preprocessor macro function that invokes
      73  * wcssub().
      74  *
      75  * wcsperr() prints the error message(s) (if any) stored in a wcsprm struct,
      76  * and the linprm, celprm, prjprm, spcprm, and tabprm structs that it contains.
      77  *
      78  * A setup routine, wcsset(), computes intermediate values in the wcsprm struct
      79  * from parameters in it that were supplied by the user.  The struct always
      80  * needs to be set up by wcsset() but this need not be called explicitly -
      81  * refer to the explanation of wcsprm::flag.
      82  *
      83  * wcsp2s() and wcss2p() implement the WCS world coordinate transformations.
      84  * In fact, they are high level driver routines for the WCS linear,
      85  * logarithmic, celestial, spectral and tabular transformation routines
      86  * described in lin.h, log.h, cel.h, spc.h and tab.h.
      87  *
      88  * Given either the celestial longitude or latitude plus an element of the
      89  * pixel coordinate a hybrid routine, wcsmix(), iteratively solves for the
      90  * unknown elements.
      91  *
      92  * wcsccs() changes the celestial coordinate system of a wcsprm struct, for
      93  * example, from equatorial to galactic, and wcssptr() translates the spectral
      94  * axis.  For example, a 'FREQ' axis may be translated into 'ZOPT-F2W' and vice
      95  * versa.
      96  *
      97  * wcslib_version() returns the WCSLIB version number.
      98  *
      99  * Quadcube projections:
     100  * ---------------------
     101  *   The quadcube projections (TSC, CSC, QSC) may be represented in FITS in
     102  *   either of two ways:
     103  *
     104  *     a: The six faces may be laid out in one plane and numbered as follows:
     105  *
     106  =                                 0
     107  =
     108  =                        4  3  2  1  4  3  2
     109  =
     110  =                                 5
     111  *
     112  *        Faces 2, 3 and 4 may appear on one side or the other (or both).  The
     113  *        world-to-pixel routines map faces 2, 3 and 4 to the left but the
     114  *        pixel-to-world routines accept them on either side.
     115  *
     116  *     b: The "COBE" convention in which the six faces are stored in a
     117  *        three-dimensional structure using a CUBEFACE axis indexed from
     118  *        0 to 5 as above.
     119  *
     120  *   These routines support both methods; wcsset() determines which is being
     121  *   used by the presence or absence of a CUBEFACE axis in ctype[].  wcsp2s()
     122  *   and wcss2p() translate the CUBEFACE axis representation to the single
     123  *   plane representation understood by the lower-level WCSLIB projection
     124  *   routines.
     125  *
     126  *
     127  * wcsnpv() - Memory allocation for PVi_ma
     128  * ---------------------------------------
     129  * wcsnpv() sets or gets the value of NPVMAX (default 64).  This global
     130  * variable controls the number of pvcard structs, for holding PVi_ma
     131  * keyvalues, that wcsini() should allocate space for.  It is also used by
     132  * wcsinit() as the default value of npvmax.
     133  *
     134  * PLEASE NOTE: This function is not thread-safe.
     135  *
     136  * Given:
     137  *   n         int       Value of NPVMAX; ignored if < 0.  Use a value less
     138  *                       than zero to get the current value.
     139  *
     140  * Function return value:
     141  *             int       Current value of NPVMAX.
     142  *
     143  *
     144  * wcsnps() - Memory allocation for PSi_ma
     145  * ---------------------------------------
     146  * wcsnps() sets or gets the value of NPSMAX (default 8).  This global variable
     147  * controls the number of pscard structs, for holding PSi_ma keyvalues, that
     148  * wcsini() should allocate space for.  It is also used by wcsinit() as the
     149  * default value of npsmax.
     150  *
     151  * PLEASE NOTE: This function is not thread-safe.
     152  *
     153  * Given:
     154  *   n         int       Value of NPSMAX; ignored if < 0.  Use a value less
     155  *                       than zero to get the current value.
     156  *
     157  * Function return value:
     158  *             int       Current value of NPSMAX.
     159  *
     160  *
     161  * wcsini() - Default constructor for the wcsprm struct
     162  * ----------------------------------------------------
     163  * wcsini() is a thin wrapper on wcsinit().  It invokes it with npvmax,
     164  * npsmax, and ndpmax set to -1 which causes it to use the values of the
     165  * global variables NDPMAX, NPSMAX, and NDPMAX.  It is thereby potentially
     166  * thread-unsafe if these variables are altered dynamically via wcsnpv(),
     167  * wcsnps(), and disndp().  Use wcsinit() for a thread-safe alternative in
     168  * this case.
     169  *
     170  *
     171  * wcsinit() - Default constructor for the wcsprm struct
     172  * -----------------------------------------------------
     173  * wcsinit() optionally allocates memory for arrays in a wcsprm struct and sets
     174  * all members of the struct to default values.
     175  *
     176  * PLEASE NOTE: every wcsprm struct should be initialized by wcsinit(),
     177  * possibly repeatedly.  On the first invokation, and only the first
     178  * invokation, wcsprm::flag must be set to -1 to initialize memory management,
     179  * regardless of whether wcsinit() will actually be used to allocate memory.
     180  *
     181  * Given:
     182  *   alloc     int       If true, allocate memory unconditionally for the
     183  *                       crpix, etc. arrays.  Please note that memory is never
     184  *                       allocated by wcsinit() for the auxprm, tabprm, nor
     185  *                       wtbarr structs.
     186  *
     187  *                       If false, it is assumed that pointers to these arrays
     188  *                       have been set by the user except if they are null
     189  *                       pointers in which case memory will be allocated for
     190  *                       them regardless.  (In other words, setting alloc true
     191  *                       saves having to initalize these pointers to zero.)
     192  *
     193  *   naxis     int       The number of world coordinate axes.  This is used to
     194  *                       determine the length of the various wcsprm vectors and
     195  *                       matrices and therefore the amount of memory to
     196  *                       allocate for them.
     197  *
     198  * Given and returned:
     199  *   wcs       struct wcsprm*
     200  *                       Coordinate transformation parameters.
     201  *
     202  *                       Note that, in order to initialize memory management,
     203  *                       wcsprm::flag should be set to -1 when wcs is
     204  *                       initialized for the first time (memory leaks may
     205  *                       result if it had already been initialized).
     206  *
     207  * Given:
     208  *   npvmax    int       The number of PVi_ma keywords to allocate space for.
     209  *                       If set to -1, the value of the global variable NPVMAX
     210  *                       will be used.  This is potentially thread-unsafe if
     211  *                       wcsnpv() is being used dynamically to alter its value.
     212  *
     213  *   npsmax    int       The number of PSi_ma keywords to allocate space for.
     214  *                       If set to -1, the value of the global variable NPSMAX
     215  *                       will be used.  This is potentially thread-unsafe if
     216  *                       wcsnps() is being used dynamically to alter its value.
     217  *
     218  *   ndpmax    int       The number of DPja or DQia keywords to allocate space
     219  *                       for.  If set to -1, the value of the global variable
     220  *                       NDPMAX will be used.  This is potentially
     221  *                       thread-unsafe if disndp() is being used dynamically to
     222  *                       alter its value.
     223  *
     224  * Function return value:
     225  *             int       Status return value:
     226  *                         0: Success.
     227  *                         1: Null wcsprm pointer passed.
     228  *                         2: Memory allocation failed.
     229  *
     230  *                       For returns > 1, a detailed error message is set in
     231  *                       wcsprm::err if enabled, see wcserr_enable().
     232  *
     233  *
     234  * wcsauxi() - Default constructor for the auxprm struct
     235  * -----------------------------------------------------
     236  * wcsauxi() optionally allocates memory for an auxprm struct, attaches it to
     237  * wcsprm, and sets all members of the struct to default values.
     238  *
     239  * Given:
     240  *   alloc     int       If true, allocate memory unconditionally for the
     241  *                       auxprm struct.
     242  *
     243  *                       If false, it is assumed that wcsprm::aux has already
     244  *                       been set to point to an auxprm struct, in which case
     245  *                       the user is responsible for managing that memory.
     246  *                       However, if wcsprm::aux is a null pointer, memory will
     247  *                       be allocated regardless.  (In other words, setting
     248  *                       alloc true saves having to initalize the pointer to
     249  *                       zero.)
     250  *
     251  * Given and returned:
     252  *   wcs       struct wcsprm*
     253  *                       Coordinate transformation parameters.
     254  *
     255  * Function return value:
     256  *             int       Status return value:
     257  *                         0: Success.
     258  *                         1: Null wcsprm pointer passed.
     259  *                         2: Memory allocation failed.
     260  *
     261  *
     262  * wcssub() - Subimage extraction routine for the wcsprm struct
     263  * ------------------------------------------------------------
     264  * wcssub() extracts the coordinate description for a subimage from a wcsprm
     265  * struct.  It does a deep copy, using wcsinit() to allocate memory for its
     266  * arrays if required.  Only the "information to be provided" part of the
     267  * struct is extracted.  Consequently, wcsset() need not have been, and won't
     268  * be invoked on the struct from which the subimage is extracted.  A call to
     269  * wcsset() is required to set up the subimage struct.
     270  *
     271  * The world coordinate system of the subimage must be separable in the sense
     272  * that the world coordinates at any point in the subimage must depend only on
     273  * the pixel coordinates of the axes extracted.  In practice, this means that
     274  * the linear transformation matrix of the original image must not contain
     275  * non-zero off-diagonal terms that associate any of the subimage axes with any
     276  * of the non-subimage axes.  Likewise, if any distortions are associated with
     277  * the subimage axes, they must not depend on any of the axes that are not
     278  * being extracted.
     279  *
     280  * Note that while the required elements of the tabprm array are extracted, the
     281  * wtbarr array is not.  (Thus it is not appropriate to call wcssub() after
     282  * wcstab() but before filling the tabprm structs - refer to wcshdr.h.)
     283  *
     284  * wcssub() can also add axes to a wcsprm struct.  The new axes will be created
     285  * using the defaults set by wcsinit() which produce a simple, unnamed, linear
     286  * axis with world coordinate equal to the pixel coordinate.  These default
     287  * values can be changed afterwards, before invoking wcsset().
     288  *
     289  * Given:
     290  *   alloc     int       If true, allocate memory for the crpix, etc. arrays in
     291  *                       the destination.  Otherwise, it is assumed that
     292  *                       pointers to these arrays have been set by the user
     293  *                       except if they are null pointers in which case memory
     294  *                       will be allocated for them regardless.
     295  *
     296  *   wcssrc    const struct wcsprm*
     297  *                       Struct to extract from.
     298  *
     299  * Given and returned:
     300  *   nsub      int*
     301  *   axes      int[]     Vector of length *nsub containing the image axis
     302  *                       numbers (1-relative) to extract.  Order is
     303  *                       significant; axes[0] is the axis number of the input
     304  *                       image that corresponds to the first axis in the
     305  *                       subimage, etc.
     306  *
     307  *                       Use an axis number of 0 to create a new axis using
     308  *                       the defaults set by wcsinit().  They can be changed
     309  *                       later.
     310  *
     311  *                       nsub (the pointer) may be set to zero, and so also may
     312  *                       *nsub, which is interpreted to mean all axes in the
     313  *                       input image; the number of axes will be returned if
     314  *                       nsub != 0x0.  axes itself (the pointer) may be set to
     315  *                       zero to indicate the first *nsub axes in their
     316  *                       original order.
     317  *
     318  *                       Set both nsub (or *nsub) and axes to zero to do a deep
     319  *                       copy of one wcsprm struct to another.
     320  *
     321  *                       Subimage extraction by coordinate axis type may be
     322  *                       done by setting the elements of axes[] to the
     323  *                       following special preprocessor macro values:
     324  *
     325  *                         WCSSUB_LONGITUDE: Celestial longitude.
     326  *                         WCSSUB_LATITUDE:  Celestial latitude.
     327  *                         WCSSUB_CUBEFACE:  Quadcube CUBEFACE axis.
     328  *                         WCSSUB_SPECTRAL:  Spectral axis.
     329  *                         WCSSUB_STOKES:    Stokes axis.
     330  *                         WCSSUB_TIME:      Time axis.
     331  *
     332  *                       Refer to the notes (below) for further usage examples.
     333  *
     334  *                       On return, *nsub will be set to the number of axes in
     335  *                       the subimage; this may be zero if there were no axes
     336  *                       of the required type(s) (in which case no memory will
     337  *                       be allocated).  axes[] will contain the axis numbers
     338  *                       that were extracted, or 0 for newly created axes.  The
     339  *                       vector length must be sufficient to contain all axis
     340  *                       numbers.  No checks are performed to verify that the
     341  *                       coordinate axes are consistent, this is done by
     342  *                       wcsset().
     343  *
     344  *   wcsdst    struct wcsprm*
     345  *                       Struct describing the subimage.  wcsprm::flag should
     346  *                       be set to -1 if wcsdst was not previously initialized
     347  *                       (memory leaks may result if it was previously
     348  *                       initialized).
     349  *
     350  * Function return value:
     351  *             int       Status return value:
     352  *                         0: Success.
     353  *                         1: Null wcsprm pointer passed.
     354  *                         2: Memory allocation failed.
     355  *                        12: Invalid subimage specification.
     356  *                        13: Non-separable subimage coordinate system.
     357  *
     358  *                       For returns > 1, a detailed error message is set in
     359  *                       wcsprm::err if enabled, see wcserr_enable().
     360  *
     361  * Notes:
     362  *   1: Combinations of subimage axes of particular types may be extracted in
     363  *      the same order as they occur in the input image by combining
     364  *      preprocessor codes, for example
     365  *
     366  =        *nsub = 1;
     367  =        axes[0] = WCSSUB_LONGITUDE | WCSSUB_LATITUDE | WCSSUB_SPECTRAL;
     368  *
     369  *      would extract the longitude, latitude, and spectral axes in the same
     370  *      order as the input image.  If one of each were present, *nsub = 3 would
     371  *      be returned.
     372  *
     373  *      For convenience, WCSSUB_CELESTIAL is defined as the combination
     374  *      WCSSUB_LONGITUDE | WCSSUB_LATITUDE | WCSSUB_CUBEFACE.
     375  *
     376  *      The codes may also be negated to extract all but the types specified,
     377  *      for example
     378  *
     379  =        *nsub = 4;
     380  =        axes[0] = WCSSUB_LONGITUDE;
     381  =        axes[1] = WCSSUB_LATITUDE;
     382  =        axes[2] = WCSSUB_CUBEFACE;
     383  =        axes[3] = -(WCSSUB_SPECTRAL | WCSSUB_STOKES);
     384  *
     385  *      The last of these specifies all axis types other than spectral or
     386  *      Stokes.  Extraction is done in the order specified by axes[] a
     387  *      longitude axis (if present) would be extracted first (via axes[0]) and
     388  *      not subsequently (via axes[3]).  Likewise for the latitude and cubeface
     389  *      axes in this example.
     390  *
     391  *      From the foregoing, it is apparent that the value of *nsub returned may
     392  *      be less than or greater than that given.  However, it will never exceed
     393  *      the number of axes in the input image (plus the number of newly-created
     394  *      axes if any were specified on input).
     395  *
     396  *
     397  * wcscompare() - Compare two wcsprm structs for equality
     398  * ------------------------------------------------------
     399  * wcscompare() compares two wcsprm structs for equality.
     400  *
     401  * Given:
     402  *   cmp       int       A bit field controlling the strictness of the
     403  *                       comparison.  When 0, all fields must be identical.
     404  *
     405  *                       The following constants may be or'ed together to
     406  *                       relax the comparison:
     407  *                         WCSCOMPARE_ANCILLARY: Ignore ancillary keywords
     408  *                           that don't change the WCS transformation, such
     409  *                           as DATE-OBS or EQUINOX.
     410  *                         WCSCOMPARE_TILING: Ignore integral differences in
     411  *                           CRPIXja.  This is the 'tiling' condition, where
     412  *                           two WCSes cover different regions of the same
     413  *                           map projection and align on the same map grid.
     414  *                         WCSCOMPARE_CRPIX: Ignore any differences at all in
     415  *                           CRPIXja.  The two WCSes cover different regions
     416  *                           of the same map projection but may not align on
     417  *                           the same map grid.  Overrides WCSCOMPARE_TILING.
     418  *
     419  *   tol       double    Tolerance for comparison of floating-point values.
     420  *                       For example, for tol == 1e-6, all floating-point
     421  *                       values in the structs must be equal to the first 6
     422  *                       decimal places.  A value of 0 implies exact equality.
     423  *
     424  *   wcs1      const struct wcsprm*
     425  *                       The first wcsprm struct to compare.
     426  *
     427  *   wcs2      const struct wcsprm*
     428  *                       The second wcsprm struct to compare.
     429  *
     430  * Returned:
     431  *   equal     int*      Non-zero when the given structs are equal.
     432  *
     433  * Function return value:
     434  *             int       Status return value:
     435  *                         0: Success.
     436  *                         1: Null pointer passed.
     437  *
     438  *
     439  * wcscopy() macro - Copy routine for the wcsprm struct
     440  * ----------------------------------------------------
     441  * wcscopy() does a deep copy of one wcsprm struct to another.  As of
     442  * WCSLIB 3.6, it is implemented as a preprocessor macro that invokes
     443  * wcssub() with the nsub and axes pointers both set to zero.
     444  *
     445  *
     446  * wcsfree() - Destructor for the wcsprm struct
     447  * --------------------------------------------
     448  * wcsfree() frees memory allocated for the wcsprm arrays by wcsinit() and/or
     449  * wcsset().  wcsinit() records the memory it allocates and wcsfree() will only
     450  * attempt to free this.
     451  *
     452  * PLEASE NOTE: wcsfree() must not be invoked on a wcsprm struct that was not
     453  * initialized by wcsinit().
     454  *
     455  * Given and returned:
     456  *   wcs       struct wcsprm*
     457  *                       Coordinate transformation parameters.
     458  *
     459  * Function return value:
     460  *             int       Status return value:
     461  *                         0: Success.
     462  *                         1: Null wcsprm pointer passed.
     463  *
     464  *
     465  * wcstrim() - Free unused arrays in the wcsprm struct
     466  * ---------------------------------------------------
     467  * wcstrim() frees memory allocated by wcsinit() for arrays in the wcsprm
     468  * struct that remains unused after it has been set up by wcsset().
     469  *
     470  * The free'd array members are associated with FITS WCS keyrecords that are
     471  * rarely used and usually just bloat the struct: wcsprm::crota, wcsprm::colax,
     472  * wcsprm::cname, wcsprm::crder, wcsprm::csyer, wcsprm::czphs, and
     473  * wcsprm::cperi.  If unused, wcsprm::pv, wcsprm::ps, and wcsprm::cd are also
     474  * freed.
     475  *
     476  * Once these arrays have been freed, a test such as
     477  =
     478  =        if (!undefined(wcs->cname[i])) {...}
     479  =
     480  * must be protected as follows
     481  =
     482  =        if (wcs->cname && !undefined(wcs->cname[i])) {...}
     483  =
     484  * In addition, if wcsprm::npv is non-zero but less than wcsprm::npvmax, then
     485  * the unused space in wcsprm::pv will be recovered (using realloc()).
     486  * Likewise for wcsprm::ps.
     487  *
     488  * Given and returned:
     489  *   wcs       struct wcsprm*
     490  *                       Coordinate transformation parameters.
     491  *
     492  * Function return value:
     493  *             int       Status return value:
     494  *                         0: Success.
     495  *                         1: Null wcsprm pointer passed.
     496  *                        14: wcsprm struct is unset.
     497  *
     498  *
     499  * wcssize() - Compute the size of a wcsprm struct
     500  * -----------------------------------------------
     501  * wcssize() computes the full size of a wcsprm struct, including allocated
     502  * memory.
     503  *
     504  * Given:
     505  *   wcs       const struct wcsprm*
     506  *                       Coordinate transformation parameters.
     507  *
     508  *                       If NULL, the base size of the struct and the allocated
     509  *                       size are both set to zero.
     510  *
     511  * Returned:
     512  *   sizes     int[2]    The first element is the base size of the struct as
     513  *                       returned by sizeof(struct wcsprm).  The second element
     514  *                       is the total allocated size, in bytes, assuming that
     515  *                       the allocation was done by wcsini().  This figure
     516  *                       includes memory allocated for members of constituent
     517  *                       structs, such as wcsprm::lin.
     518  *
     519  *                       It is not an error for the struct not to have been set
     520  *                       up via wcsset(), which normally results in additional
     521  *                       memory allocation. 
     522  *
     523  * Function return value:
     524  *             int       Status return value:
     525  *                         0: Success.
     526  *
     527  *
     528  * auxsize() - Compute the size of a auxprm struct
     529  * -----------------------------------------------
     530  * auxsize() computes the full size of an auxprm struct, including allocated
     531  * memory.
     532  *
     533  * Given:
     534  *   aux       const struct auxprm*
     535  *                       Auxiliary coordinate information.
     536  *
     537  *                       If NULL, the base size of the struct and the allocated
     538  *                       size are both set to zero.
     539  *
     540  * Returned:
     541  *   sizes     int[2]    The first element is the base size of the struct as
     542  *                       returned by sizeof(struct auxprm).  The second element
     543  *                       is the total allocated size, in bytes, currently zero.
     544  *
     545  * Function return value:
     546  *             int       Status return value:
     547  *                         0: Success.
     548  *
     549  *
     550  * wcsprt() - Print routine for the wcsprm struct
     551  * ----------------------------------------------
     552  * wcsprt() prints the contents of a wcsprm struct using wcsprintf().  Mainly
     553  * intended for diagnostic purposes.
     554  *
     555  * Given:
     556  *   wcs       const struct wcsprm*
     557  *                       Coordinate transformation parameters.
     558  *
     559  * Function return value:
     560  *             int       Status return value:
     561  *                         0: Success.
     562  *                         1: Null wcsprm pointer passed.
     563  *
     564  *
     565  * wcsperr() - Print error messages from a wcsprm struct
     566  * -----------------------------------------------------
     567  * wcsperr() prints the error message(s), if any, stored in a wcsprm struct,
     568  * and the linprm, celprm, prjprm, spcprm, and tabprm structs that it contains.
     569  * If there are no errors then nothing is printed.  It uses wcserr_prt(), q.v.
     570  *
     571  * Given:
     572  *   wcs       const struct wcsprm*
     573  *                       Coordinate transformation parameters.
     574  *
     575  *   prefix    const char *
     576  *                       If non-NULL, each output line will be prefixed with
     577  *                       this string.
     578  *
     579  * Function return value:
     580  *             int       Status return value:
     581  *                         0: Success.
     582  *                         1: Null wcsprm pointer passed.
     583  *
     584  *
     585  * wcsbchk() - Enable/disable bounds checking
     586  * ------------------------------------------
     587  * wcsbchk() is used to control bounds checking in the projection routines.
     588  * Note that wcsset() always enables bounds checking.  wcsbchk() will invoke
     589  * wcsset() on the wcsprm struct beforehand if necessary.
     590  *
     591  * Given and returned:
     592  *   wcs       struct wcsprm*
     593  *                       Coordinate transformation parameters.
     594  *
     595  * Given:
     596  *   bounds    int       If bounds&1 then enable strict bounds checking for the
     597  *                       spherical-to-Cartesian (s2x) transformation for the
     598  *                       AZP, SZP, TAN, SIN, ZPN, and COP projections.
     599  *
     600  *                       If bounds&2 then enable strict bounds checking for the
     601  *                       Cartesian-to-spherical (x2s) transformation for the
     602  *                       HPX and XPH projections.
     603  *
     604  *                       If bounds&4 then enable bounds checking on the native
     605  *                       coordinates returned by the Cartesian-to-spherical
     606  *                       (x2s) transformations using prjchk().
     607  *
     608  *                       Zero it to disable all checking.
     609  *
     610  * Function return value:
     611  *             int       Status return value:
     612  *                         0: Success.
     613  *                         1: Null wcsprm pointer passed.
     614  *
     615  *
     616  * wcsset() - Setup routine for the wcsprm struct
     617  * ----------------------------------------------
     618  * wcsset() sets up a wcsprm struct according to information supplied within
     619  * it (refer to the description of the wcsprm struct).
     620  *
     621  * wcsset() recognizes the NCP projection and converts it to the equivalent SIN
     622  * projection and likewise translates GLS into SFL.  It also translates the
     623  * AIPS spectral types ('FREQ-LSR', 'FELO-HEL', etc.), possibly changing the
     624  * input header keywords wcsprm::ctype and/or wcsprm::specsys if necessary.
     625  *
     626  * Note that this routine need not be called directly; it will be invoked by
     627  * wcsp2s() and wcss2p() if the wcsprm::flag is anything other than a
     628  * predefined magic value.
     629  *
     630  * Given and returned:
     631  *   wcs       struct wcsprm*
     632  *                       Coordinate transformation parameters.
     633  *
     634  * Function return value:
     635  *             int       Status return value:
     636  *                         0: Success.
     637  *                         1: Null wcsprm pointer passed.
     638  *                         2: Memory allocation failed.
     639  *                         3: Linear transformation matrix is singular.
     640  *                         4: Inconsistent or unrecognized coordinate axis
     641  *                            types.
     642  *                         5: Invalid parameter value.
     643  *                         6: Invalid coordinate transformation parameters.
     644  *                         7: Ill-conditioned coordinate transformation
     645  *                            parameters.
     646  *
     647  *                       For returns > 1, a detailed error message is set in
     648  *                       wcsprm::err if enabled, see wcserr_enable().
     649  *
     650  * Notes:
     651  *   1: wcsset() always enables strict bounds checking in the projection
     652  *      routines (via a call to prjini()).  Use wcsbchk() to modify
     653  *      bounds-checking after wcsset() is invoked.
     654  *
     655  *
     656  * wcsp2s() - Pixel-to-world transformation
     657  * ----------------------------------------
     658  * wcsp2s() transforms pixel coordinates to world coordinates.
     659  *
     660  * Given and returned:
     661  *   wcs       struct wcsprm*
     662  *                       Coordinate transformation parameters.
     663  *
     664  * Given:
     665  *   ncoord,
     666  *   nelem     int       The number of coordinates, each of vector length
     667  *                       nelem but containing wcs.naxis coordinate elements.
     668  *                       Thus nelem must equal or exceed the value of the
     669  *                       NAXIS keyword unless ncoord == 1, in which case nelem
     670  *                       is not used.
     671  *
     672  *   pixcrd    const double[ncoord][nelem]
     673  *                       Array of pixel coordinates.
     674  *
     675  * Returned:
     676  *   imgcrd    double[ncoord][nelem]
     677  *                       Array of intermediate world coordinates.  For
     678  *                       celestial axes, imgcrd[][wcs.lng] and
     679  *                       imgcrd[][wcs.lat] are the projected x-, and
     680  *                       y-coordinates in pseudo "degrees".  For spectral
     681  *                       axes, imgcrd[][wcs.spec] is the intermediate spectral
     682  *                       coordinate, in SI units.  For time axes,
     683  *                       imgcrd[][wcs.time] is the intermediate time
     684  *                       coordinate.
     685  *
     686  *   phi,theta double[ncoord]
     687  *                       Longitude and latitude in the native coordinate system
     688  *                       of the projection [deg].
     689  *
     690  *   world     double[ncoord][nelem]
     691  *                       Array of world coordinates.  For celestial axes,
     692  *                       world[][wcs.lng] and world[][wcs.lat] are the
     693  *                       celestial longitude and latitude [deg].  For spectral
     694  *                       axes, world[][wcs.spec] is the spectral coordinate, in
     695  *                       SI units.  For time axes, world[][wcs.time] is the
     696  *                       time coordinate.
     697  *
     698  *   stat      int[ncoord]
     699  *                       Status return value for each coordinate:
     700  *                         0: Success.
     701  *                        1+: A bit mask indicating invalid pixel coordinate
     702  *                            element(s).
     703  *
     704  * Function return value:
     705  *             int       Status return value:
     706  *                         0: Success.
     707  *                         1: Null wcsprm pointer passed.
     708  *                         2: Memory allocation failed.
     709  *                         3: Linear transformation matrix is singular.
     710  *                         4: Inconsistent or unrecognized coordinate axis
     711  *                            types.
     712  *                         5: Invalid parameter value.
     713  *                         6: Invalid coordinate transformation parameters.
     714  *                         7: Ill-conditioned coordinate transformation
     715  *                            parameters.
     716  *                         8: One or more of the pixel coordinates were
     717  *                            invalid, as indicated by the stat vector.
     718  *
     719  *                       For returns > 1, a detailed error message is set in
     720  *                       wcsprm::err if enabled, see wcserr_enable().
     721  *
     722  *
     723  * wcss2p() - World-to-pixel transformation
     724  * ----------------------------------------
     725  * wcss2p() transforms world coordinates to pixel coordinates.
     726  *
     727  * Given and returned:
     728  *   wcs       struct wcsprm*
     729  *                       Coordinate transformation parameters.
     730  *
     731  * Given:
     732  *   ncoord,
     733  *   nelem     int       The number of coordinates, each of vector length nelem
     734  *                       but containing wcs.naxis coordinate elements.  Thus
     735  *                       nelem must equal or exceed the value of the NAXIS
     736  *                       keyword unless ncoord == 1, in which case nelem is not
     737  *                       used.
     738  *
     739  *   world     const double[ncoord][nelem]
     740  *                       Array of world coordinates.  For celestial axes,
     741  *                       world[][wcs.lng] and world[][wcs.lat] are the
     742  *                       celestial longitude and latitude [deg]. For spectral
     743  *                       axes, world[][wcs.spec] is the spectral coordinate, in
     744  *                       SI units.  For time axes, world[][wcs.time] is the
     745  *                       time coordinate.
     746  *
     747  * Returned:
     748  *   phi,theta double[ncoord]
     749  *                       Longitude and latitude in the native coordinate
     750  *                       system of the projection [deg].
     751  *
     752  *   imgcrd    double[ncoord][nelem]
     753  *                       Array of intermediate world coordinates.  For
     754  *                       celestial axes, imgcrd[][wcs.lng] and
     755  *                       imgcrd[][wcs.lat] are the projected x-, and
     756  *                       y-coordinates in pseudo "degrees".  For quadcube
     757  *                       projections with a CUBEFACE axis the face number is
     758  *                       also returned in imgcrd[][wcs.cubeface].  For
     759  *                       spectral axes, imgcrd[][wcs.spec] is the intermediate
     760  *                       spectral coordinate, in SI units.  For time axes,
     761  *                       imgcrd[][wcs.time] is the intermediate time
     762  *                       coordinate.
     763  *
     764  *   pixcrd    double[ncoord][nelem]
     765  *                       Array of pixel coordinates.
     766  *
     767  *   stat      int[ncoord]
     768  *                       Status return value for each coordinate:
     769  *                         0: Success.
     770  *                        1+: A bit mask indicating invalid world coordinate
     771  *                            element(s).
     772  *
     773  * Function return value:
     774  *             int       Status return value:
     775  *                         0: Success.
     776  *                         1: Null wcsprm pointer passed.
     777  *                         2: Memory allocation failed.
     778  *                         3: Linear transformation matrix is singular.
     779  *                         4: Inconsistent or unrecognized coordinate axis
     780  *                            types.
     781  *                         5: Invalid parameter value.
     782  *                         6: Invalid coordinate transformation parameters.
     783  *                         7: Ill-conditioned coordinate transformation
     784  *                            parameters.
     785  *                         9: One or more of the world coordinates were
     786  *                            invalid, as indicated by the stat vector.
     787  *
     788  *                       For returns > 1, a detailed error message is set in
     789  *                       wcsprm::err if enabled, see wcserr_enable().
     790  *
     791  *
     792  * wcsmix() - Hybrid coordinate transformation
     793  * -------------------------------------------
     794  * wcsmix(), given either the celestial longitude or latitude plus an element
     795  * of the pixel coordinate, solves for the remaining elements by iterating on
     796  * the unknown celestial coordinate element using wcss2p().  Refer also to the
     797  * notes below.
     798  *
     799  * Given and returned:
     800  *   wcs       struct wcsprm*
     801  *                       Indices for the celestial coordinates obtained
     802  *                       by parsing the wcsprm::ctype[].
     803  *
     804  * Given:
     805  *   mixpix    int       Which element of the pixel coordinate is given.
     806  *
     807  *   mixcel    int       Which element of the celestial coordinate is given:
     808  *                         1: Celestial longitude is given in
     809  *                            world[wcs.lng], latitude returned in
     810  *                            world[wcs.lat].
     811  *                         2: Celestial latitude is given in
     812  *                            world[wcs.lat], longitude returned in
     813  *                            world[wcs.lng].
     814  *
     815  *   vspan     const double[2]
     816  *                       Solution interval for the celestial coordinate [deg].
     817  *                       The ordering of the two limits is irrelevant.
     818  *                       Longitude ranges may be specified with any convenient
     819  *                       normalization, for example [-120,+120] is the same as
     820  *                       [240,480], except that the solution will be returned
     821  *                       with the same normalization, i.e. lie within the
     822  *                       interval specified.
     823  *
     824  *   vstep     const double
     825  *                       Step size for solution search [deg].  If zero, a
     826  *                       sensible, although perhaps non-optimal default will be
     827  *                       used.
     828  *
     829  *   viter     int       If a solution is not found then the step size will be
     830  *                       halved and the search recommenced.  viter controls how
     831  *                       many times the step size is halved.  The allowed range
     832  *                       is 5 - 10.
     833  *
     834  * Given and returned:
     835  *   world     double[naxis]
     836  *                       World coordinate elements.  world[wcs.lng] and
     837  *                       world[wcs.lat] are the celestial longitude and
     838  *                       latitude [deg].  Which is given and which returned
     839  *                       depends on the value of mixcel.  All other elements
     840  *                       are given.
     841  *
     842  * Returned:
     843  *   phi,theta double[naxis]
     844  *                       Longitude and latitude in the native coordinate
     845  *                       system of the projection [deg].
     846  *
     847  *   imgcrd    double[naxis]
     848  *                       Image coordinate elements.  imgcrd[wcs.lng] and
     849  *                       imgcrd[wcs.lat] are the projected x-, and
     850  *                       y-coordinates in pseudo "degrees".
     851  *
     852  * Given and returned:
     853  *   pixcrd    double[naxis]
     854  *                       Pixel coordinate.  The element indicated by mixpix is
     855  *                       given and the remaining elements are returned.
     856  *
     857  * Function return value:
     858  *             int       Status return value:
     859  *                         0: Success.
     860  *                         1: Null wcsprm pointer passed.
     861  *                         2: Memory allocation failed.
     862  *                         3: Linear transformation matrix is singular.
     863  *                         4: Inconsistent or unrecognized coordinate axis
     864  *                            types.
     865  *                         5: Invalid parameter value.
     866  *                         6: Invalid coordinate transformation parameters.
     867  *                         7: Ill-conditioned coordinate transformation
     868  *                            parameters.
     869  *                        10: Invalid world coordinate.
     870  *                        11: No solution found in the specified interval.
     871  *
     872  *                       For returns > 1, a detailed error message is set in
     873  *                       wcsprm::err if enabled, see wcserr_enable().
     874  *
     875  * Notes:
     876  *   1: Initially the specified solution interval is checked to see if it's a
     877  *      "crossing" interval.  If it isn't, a search is made for a crossing
     878  *      solution by iterating on the unknown celestial coordinate starting at
     879  *      the upper limit of the solution interval and decrementing by the
     880  *      specified step size.  A crossing is indicated if the trial value of the
     881  *      pixel coordinate steps through the value specified.  If a crossing
     882  *      interval is found then the solution is determined by a modified form of
     883  *      "regula falsi" division of the crossing interval.  If no crossing
     884  *      interval was found within the specified solution interval then a search
     885  *      is made for a "non-crossing" solution as may arise from a point of
     886  *      tangency.  The process is complicated by having to make allowance for
     887  *      the discontinuities that occur in all map projections.
     888  *
     889  *      Once one solution has been determined others may be found by subsequent
     890  *      invokations of wcsmix() with suitably restricted solution intervals.
     891  *
     892  *      Note the circumstance that arises when the solution point lies at a
     893  *      native pole of a projection in which the pole is represented as a
     894  *      finite curve, for example the zenithals and conics.  In such cases two
     895  *      or more valid solutions may exist but wcsmix() only ever returns one.
     896  *
     897  *      Because of its generality wcsmix() is very compute-intensive.  For
     898  *      compute-limited applications more efficient special-case solvers could
     899  *      be written for simple projections, for example non-oblique cylindrical
     900  *      projections.
     901  *
     902  *
     903  * wcsccs() - Change celestial coordinate system
     904  * ---------------------------------------------
     905  * wcsccs() changes the celestial coordinate system of a wcsprm struct.  For
     906  * example, from equatorial to galactic coordinates.
     907  *
     908  * Parameters that define the spherical coordinate transformation, essentially
     909  * being three Euler angles, must be provided.  Thereby wcsccs() does not need
     910  * prior knowledge of specific celestial coordinate systems.  It also has the
     911  * advantage of making it completely general.
     912  *
     913  * Auxiliary members of the wcsprm struct relating to equatorial celestial
     914  * coordinate systems may also be changed.
     915  *
     916  * Only orthodox spherical coordinate systems are supported.  That is, they
     917  * must be right-handed, with latitude increasing from zero at the equator to
     918  * +90 degrees at the pole.  This precludes systems such as aziumuth and zenith
     919  * distance, which, however, could be handled as negative azimuth and
     920  * elevation.
     921  *
     922  * PLEASE NOTE: Information in the wcsprm struct relating to the original
     923  * coordinate system will be overwritten and therefore lost.  If this is
     924  * undesirable, invoke wcsccs() on a copy of the struct made with wcssub().
     925  * The wcsprm struct is reset on return with an explicit call to wcsset().
     926  *
     927  * Given and returned:
     928  *   wcs       struct wcsprm*
     929  *                       Coordinate transformation parameters.  Particular
     930  *                       "values to be given" elements of the wcsprm struct
     931  *                       are modified.
     932  *
     933  * Given:
     934  *   lng2p1,
     935  *   lat2p1    double    Longitude and latitude in the new celestial coordinate
     936  *                       system of the pole (i.e. latitude +90) of the original
     937  *                       system [deg].  See notes 1 and 2 below.
     938  *
     939  *   lng1p2    double    Longitude in the original celestial coordinate system
     940  *                       of the pole (i.e. latitude +90) of the new system
     941  *                       [deg].  See note 1 below.
     942  *
     943  *   clng,clat const char*
     944  *                       Longitude and latitude identifiers of the new CTYPEia
     945  *                       celestial axis codes, without trailing dashes.  For
     946  *                       example, "RA" and "DEC" or "GLON" and "GLAT".  Up to
     947  *                       four characters are used, longer strings need not be
     948  *                       null-terminated.
     949  *
     950  *   radesys   const char*
     951  *                       Used when transforming to equatorial coordinates,
     952  *                       identified by clng == "RA" and clat = "DEC".  May be
     953  *                       set to the null pointer to preserve the current value.
     954  *                       Up to 71 characters are used, longer strings need not
     955  *                       be null-terminated.
     956  *
     957  *                       If the new coordinate system is anything other than
     958  *                       equatorial, then wcsprm::radesys will be cleared.
     959  *
     960  *   equinox   double    Used when transforming to equatorial coordinates.  May
     961  *                       be set to zero to preserve the current value.
     962  *
     963  *                       If the new coordinate system is not equatorial, then
     964  *                       wcsprm::equinox will be marked as undefined.
     965  *
     966  *   alt       const char*
     967  *                       Character code for alternate coordinate descriptions
     968  *                       (i.e. the 'a' in keyword names such as CTYPEia).  This
     969  *                       is blank for the primary coordinate description, or
     970  *                       one of the 26 upper-case letters, A-Z.  May be set to
     971  *                       the null pointer, or null string if no change is
     972  *                       required.
     973  *
     974  * Function return value:
     975  *             int       Status return value:
     976  *                         0: Success.
     977  *                         1: Null wcsprm pointer passed.
     978  *                        12: Invalid subimage specification (no celestial
     979  *                            axes).
     980  *
     981  * Notes:
     982  *   1: Follows the prescription given in WCS Paper II, Sect. 2.7 for changing
     983  *      celestial coordinates.
     984  *
     985  *      The implementation takes account of indeterminacies that arise in that
     986  *      prescription in the particular cases where one of the poles of the new
     987  *      system is at the fiducial point, or one of them is at the native pole.
     988  *
     989  *   2: If lat2p1 == +90, i.e. where the poles of the two coordinate systems
     990  *      coincide, then the spherical coordinate transformation becomes a simple
     991  *      change in origin of longitude given by
     992  *      lng2 = lng1 + (lng2p1 - lng1p2 - 180), and lat2 = lat1, where
     993  *      (lng2,lat2) are coordinates in the new system, and (lng1,lat1) are
     994  *      coordinates in the original system.
     995  *
     996  *      Likewise, if lat2p1 == -90, then lng2 = -lng1 + (lng2p1 + lng1p2), and
     997  *      lat2 = -lat1.
     998  *
     999  *   3: For example, if the original coordinate system is B1950 equatorial and
    1000  *      the desired new coordinate system is galactic, then
    1001  *
    1002  *      - (lng2p1,lat2p1) are the galactic coordinates of the B1950 celestial
    1003  *        pole, defined by the IAU to be (123.0,+27.4), and lng1p2 is the B1950
    1004  *        right ascension of the galactic pole, defined as 192.25.  Clearly
    1005  *        these coordinates are fixed for a particular coordinate
    1006  *        transformation.
    1007  *
    1008  *      - (clng,clat) would be 'GLON' and 'GLAT', these being the FITS standard
    1009  *        identifiers for galactic coordinates.
    1010  *
    1011  *      - Since the new coordinate system is not equatorial, wcsprm::radesys
    1012  *        and wcsprm::equinox will be cleared.
    1013  *
    1014  *   4. The coordinates required for some common transformations (obtained from
    1015  *      https://ned.ipac.caltech.edu/coordinate_calculator) are as follows:
    1016  *
    1017  =      (123.0000,+27.4000) galactic coordinates of B1950 celestial pole,
    1018  =      (192.2500,+27.4000) B1950 equatorial coordinates of galactic pole.
    1019  *
    1020  =      (122.9319,+27.1283) galactic coordinates of J2000 celestial pole,
    1021  =      (192.8595,+27.1283) J2000 equatorial coordinates of galactic pole.
    1022  *
    1023  =      (359.6774,+89.7217) B1950 equatorial coordinates of J2000 pole,
    1024  =      (180.3162,+89.7217) J2000 equatorial coordinates of B1950 pole.
    1025  *
    1026  =      (270.0000,+66.5542) B1950 equatorial coordinates of B1950 ecliptic pole,
    1027  =      ( 90.0000,+66.5542) B1950 ecliptic coordinates of B1950 celestial pole.
    1028  *
    1029  =      (270.0000,+66.5607) J2000 equatorial coordinates of J2000 ecliptic pole,
    1030  =      ( 90.0000,+66.5607) J2000 ecliptic coordinates of J2000 celestial pole.
    1031  *
    1032  =      ( 26.7315,+15.6441) supergalactic coordinates of B1950 celestial pole,
    1033  =      (283.1894,+15.6441) B1950 equatorial coordinates of supergalactic pole.
    1034  *
    1035  =      ( 26.4505,+15.7089) supergalactic coordinates of J2000 celestial pole,
    1036  =      (283.7542,+15.7089) J2000 equatorial coordinates of supergalactic pole.
    1037  *
    1038  *
    1039  * wcssptr() - Spectral axis translation
    1040  * -------------------------------------
    1041  * wcssptr() translates the spectral axis in a wcsprm struct.  For example, a
    1042  * 'FREQ' axis may be translated into 'ZOPT-F2W' and vice versa.
    1043  *
    1044  * PLEASE NOTE: Information in the wcsprm struct relating to the original
    1045  * coordinate system will be overwritten and therefore lost.  If this is
    1046  * undesirable, invoke wcssptr() on a copy of the struct made with wcssub().
    1047  * The wcsprm struct is reset on return with an explicit call to wcsset().
    1048  *
    1049  * Given and returned:
    1050  *   wcs       struct wcsprm*
    1051  *                       Coordinate transformation parameters.
    1052  *
    1053  *   i         int*      Index of the spectral axis (0-relative).  If given < 0
    1054  *                       it will be set to the first spectral axis identified
    1055  *                       from the ctype[] keyvalues in the wcsprm struct.
    1056  *
    1057  *   ctype     char[9]   Desired spectral CTYPEia.  Wildcarding may be used as
    1058  *                       for the ctypeS2 argument to spctrn() as described in
    1059  *                       the prologue of spc.h, i.e. if the final three
    1060  *                       characters are specified as "???", or if just the
    1061  *                       eighth character is specified as '?', the correct
    1062  *                       algorithm code will be substituted and returned.
    1063  *
    1064  * Function return value:
    1065  *             int       Status return value:
    1066  *                         0: Success.
    1067  *                         1: Null wcsprm pointer passed.
    1068  *                         2: Memory allocation failed.
    1069  *                         3: Linear transformation matrix is singular.
    1070  *                         4: Inconsistent or unrecognized coordinate axis
    1071  *                            types.
    1072  *                         5: Invalid parameter value.
    1073  *                         6: Invalid coordinate transformation parameters.
    1074  *                         7: Ill-conditioned coordinate transformation
    1075  *                            parameters.
    1076  *                        12: Invalid subimage specification (no spectral
    1077  *                            axis).
    1078  *
    1079  *                       For returns > 1, a detailed error message is set in
    1080  *                       wcsprm::err if enabled, see wcserr_enable().
    1081  *
    1082  *
    1083  * wcslib_version() - WCSLIB version number
    1084  * ----------------------------------------
    1085  * wcslib_version() returns the WCSLIB version number.
    1086  *
    1087  * The major version number changes when the ABI changes or when the license
    1088  * conditions change.  ABI changes typically result from a change to the
    1089  * contents of one of the structs.  The major version number is used to
    1090  * distinguish between incompatible versions of the sharable library.
    1091  *
    1092  * The minor version number changes with new functionality or bug fixes that do
    1093  * not involve a change in the ABI.
    1094  *
    1095  * The auxiliary version number (which is often absent) signals changes to the
    1096  * documentation, test suite, build procedures, or any other change that does
    1097  * not affect the compiled library.
    1098  *
    1099  * Returned:
    1100  *   vers[3]   int[3]    The broken-down version number:
    1101  *                         0: Major version number.
    1102  *                         1: Minor version number.
    1103  *                         2: Auxiliary version number (zero if absent).
    1104  *                       May be given as a null pointer if not required.
    1105  *
    1106  * Function return value:
    1107  *             char*     A null-terminated, statically allocated string
    1108  *                       containing the version number in the usual form, i.e.
    1109  *                       "<major>.<minor>.<auxiliary>".
    1110  *
    1111  *
    1112  * wcsprm struct - Coordinate transformation parameters
    1113  * ----------------------------------------------------
    1114  * The wcsprm struct contains information required to transform world
    1115  * coordinates.  It consists of certain members that must be set by the user
    1116  * ("given") and others that are set by the WCSLIB routines ("returned").
    1117  * While the addresses of the arrays themselves may be set by wcsinit() if it
    1118  * (optionally) allocates memory, their contents must be set by the user.
    1119  *
    1120  * Some parameters that are given are not actually required for transforming
    1121  * coordinates.  These are described as "auxiliary"; the struct simply provides
    1122  * a place to store them, though they may be used by wcshdo() in constructing a
    1123  * FITS header from a wcsprm struct.  Some of the returned values are supplied
    1124  * for informational purposes and others are for internal use only as
    1125  * indicated.
    1126  *
    1127  * In practice, it is expected that a WCS parser would scan the FITS header to
    1128  * determine the number of coordinate axes.  It would then use wcsinit() to
    1129  * allocate memory for arrays in the wcsprm struct and set default values.
    1130  * Then as it reread the header and identified each WCS keyrecord it would load
    1131  * the value into the relevant wcsprm array element.  This is essentially what
    1132  * wcspih() does - refer to the prologue of wcshdr.h.  As the final step,
    1133  * wcsset() is invoked, either directly or indirectly, to set the derived
    1134  * members of the wcsprm struct.  wcsset() strips off trailing blanks in all
    1135  * string members and null-fills the character array.
    1136  *
    1137  *   int flag
    1138  *     (Given and returned) This flag must be set to zero whenever any of the
    1139  *     following wcsprm struct members are set or changed:
    1140  *
    1141  *       - wcsprm::naxis (q.v., not normally set by the user),
    1142  *       - wcsprm::crpix,
    1143  *       - wcsprm::pc,
    1144  *       - wcsprm::cdelt,
    1145  *       - wcsprm::crval,
    1146  *       - wcsprm::cunit,
    1147  *       - wcsprm::ctype,
    1148  *       - wcsprm::lonpole,
    1149  *       - wcsprm::latpole,
    1150  *       - wcsprm::restfrq,
    1151  *       - wcsprm::restwav,
    1152  *       - wcsprm::npv,
    1153  *       - wcsprm::pv,
    1154  *       - wcsprm::nps,
    1155  *       - wcsprm::ps,
    1156  *       - wcsprm::cd,
    1157  *       - wcsprm::crota,
    1158  *       - wcsprm::altlin,
    1159  *       - wcsprm::ntab,
    1160  *       - wcsprm::nwtb,
    1161  *       - wcsprm::tab,
    1162  *       - wcsprm::wtb.
    1163  *
    1164  *     This signals the initialization routine, wcsset(), to recompute the
    1165  *     returned members of the linprm, celprm, spcprm, and tabprm structs.
    1166  *     wcsset() will reset flag to indicate that this has been done.
    1167  *
    1168  *     PLEASE NOTE: flag should be set to -1 when wcsinit() is called for the
    1169  *     first time for a particular wcsprm struct in order to initialize memory
    1170  *     management.  It must ONLY be used on the first initialization otherwise
    1171  *     memory leaks may result.
    1172  *
    1173  *   int naxis
    1174  *     (Given or returned) Number of pixel and world coordinate elements.
    1175  *
    1176  *     If wcsinit() is used to initialize the linprm struct (as would normally
    1177  *     be the case) then it will set naxis from the value passed to it as a
    1178  *     function argument.  The user should not subsequently modify it.
    1179  *
    1180  *   double *crpix
    1181  *     (Given) Address of the first element of an array of double containing
    1182  *     the coordinate reference pixel, CRPIXja.
    1183  *
    1184  *   double *pc
    1185  *     (Given) Address of the first element of the PCi_ja (pixel coordinate)
    1186  *     transformation matrix.  The expected order is
    1187  *
    1188  =       struct wcsprm wcs;
    1189  =       wcs.pc = {PC1_1, PC1_2, PC2_1, PC2_2};
    1190  *
    1191  *     This may be constructed conveniently from a 2-D array via
    1192  *
    1193  =       double m[2][2] = {{PC1_1, PC1_2},
    1194  =                         {PC2_1, PC2_2}};
    1195  *
    1196  *     which is equivalent to
    1197  *
    1198  =       double m[2][2];
    1199  =       m[0][0] = PC1_1;
    1200  =       m[0][1] = PC1_2;
    1201  =       m[1][0] = PC2_1;
    1202  =       m[1][1] = PC2_2;
    1203  *
    1204  *     The storage order for this 2-D array is the same as for the 1-D array,
    1205  *     whence
    1206  *
    1207  =       wcs.pc = *m;
    1208  *
    1209  *     would be legitimate.
    1210  *
    1211  *   double *cdelt
    1212  *     (Given) Address of the first element of an array of double containing
    1213  *     the coordinate increments, CDELTia.
    1214  *
    1215  *   double *crval
    1216  *     (Given) Address of the first element of an array of double containing
    1217  *     the coordinate reference values, CRVALia.
    1218  *
    1219  *   char (*cunit)[72]
    1220  *     (Given) Address of the first element of an array of char[72] containing
    1221  *     the CUNITia keyvalues which define the units of measurement of the
    1222  *     CRVALia, CDELTia, and CDi_ja keywords.
    1223  *
    1224  *     As CUNITia is an optional header keyword, cunit[][72] may be left blank
    1225  *     but otherwise is expected to contain a standard units specification as
    1226  *     defined by WCS Paper I.  Utility function wcsutrn(), described in
    1227  *     wcsunits.h, is available to translate commonly used non-standard units
    1228  *     specifications but this must be done as a separate step before invoking
    1229  *     wcsset().
    1230  *
    1231  *     For celestial axes, if cunit[][72] is not blank, wcsset() uses
    1232  *     wcsunits() to parse it and scale cdelt[], crval[], and cd[][*] to
    1233  *     degrees.  It then resets cunit[][72] to "deg".
    1234  *
    1235  *     For spectral axes, if cunit[][72] is not blank, wcsset() uses wcsunits()
    1236  *     to parse it and scale cdelt[], crval[], and cd[][*] to SI units.  It
    1237  *     then resets cunit[][72] accordingly.
    1238  *
    1239  *     wcsset() ignores cunit[][72] for other coordinate types; cunit[][72] may
    1240  *     be used to label coordinate values.
    1241  *
    1242  *     These variables accomodate the longest allowed string-valued FITS
    1243  *     keyword, being limited to 68 characters, plus the null-terminating
    1244  *     character.
    1245  *
    1246  *   char (*ctype)[72]
    1247  *     (Given) Address of the first element of an array of char[72] containing
    1248  *     the coordinate axis types, CTYPEia.
    1249  *
    1250  *     The ctype[][72] keyword values must be in upper case and there must be
    1251  *     zero or one pair of matched celestial axis types, and zero or one
    1252  *     spectral axis.  The ctype[][72] strings should be padded with blanks on
    1253  *     the right and null-terminated so that they are at least eight characters
    1254  *     in length.
    1255  *
    1256  *     These variables accomodate the longest allowed string-valued FITS
    1257  *     keyword, being limited to 68 characters, plus the null-terminating
    1258  *     character.
    1259  *
    1260  *   double lonpole
    1261  *     (Given and returned) The native longitude of the celestial pole, phi_p,
    1262  *     given by LONPOLEa [deg] or by PVi_2a [deg] attached to the longitude
    1263  *     axis which takes precedence if defined, and ...
    1264  *   double latpole
    1265  *     (Given and returned) ... the native latitude of the celestial pole,
    1266  *     theta_p, given by LATPOLEa [deg] or by PVi_3a [deg] attached to the
    1267  *     longitude axis which takes precedence if defined.
    1268  *
    1269  *     lonpole and latpole may be left to default to values set by wcsinit()
    1270  *     (see celprm::ref), but in any case they will be reset by wcsset() to
    1271  *     the values actually used.  Note therefore that if the wcsprm struct is
    1272  *     reused without resetting them, whether directly or via wcsinit(), they
    1273  *     will no longer have their default values.
    1274  *
    1275  *   double restfrq
    1276  *     (Given) The rest frequency [Hz], and/or ...
    1277  *   double restwav
    1278  *     (Given) ... the rest wavelength in vacuo [m], only one of which need be
    1279  *     given, the other should be set to zero.
    1280  *
    1281  *   int npv
    1282  *     (Given) The number of entries in the wcsprm::pv[] array.
    1283  *
    1284  *   int npvmax
    1285  *     (Given or returned) The length of the wcsprm::pv[] array.
    1286  *
    1287  *     npvmax will be set by wcsinit() if it allocates memory for wcsprm::pv[],
    1288  *     otherwise it must be set by the user.  See also wcsnpv().
    1289  *
    1290  *   struct pvcard *pv
    1291  *     (Given) Address of the first element of an array of length npvmax of
    1292  *     pvcard structs.
    1293  *
    1294  *     As a FITS header parser encounters each PVi_ma keyword it should load it
    1295  *     into a pvcard struct in the array and increment npv.  wcsset()
    1296  *     interprets these as required.
    1297  *
    1298  *     Note that, if they were not given, wcsset() resets the entries for
    1299  *     PVi_1a, PVi_2a, PVi_3a, and PVi_4a for longitude axis i to match
    1300  *     phi_0 and theta_0 (the native longitude and latitude of the reference
    1301  *     point), LONPOLEa and LATPOLEa respectively.
    1302  *
    1303  *   int nps
    1304  *     (Given) The number of entries in the wcsprm::ps[] array.
    1305  *
    1306  *   int npsmax
    1307  *     (Given or returned) The length of the wcsprm::ps[] array.
    1308  *
    1309  *     npsmax will be set by wcsinit() if it allocates memory for wcsprm::ps[],
    1310  *     otherwise it must be set by the user.  See also wcsnps().
    1311  *
    1312  *   struct pscard *ps
    1313  *     (Given) Address of the first element of an array of length npsmax of
    1314  *     pscard structs.
    1315  *
    1316  *     As a FITS header parser encounters each PSi_ma keyword it should load it
    1317  *     into a pscard struct in the array and increment nps.  wcsset()
    1318  *     interprets these as required (currently no PSi_ma keyvalues are
    1319  *     recognized).
    1320  *
    1321  *   double *cd
    1322  *     (Given) For historical compatibility, the wcsprm struct supports two
    1323  *     alternate specifications of the linear transformation matrix, those
    1324  *     associated with the CDi_ja keywords, and ...
    1325  *   double *crota
    1326  *     (Given) ... those associated with the CROTAi keywords.  Although these
    1327  *     may not formally co-exist with PCi_ja, the approach taken here is simply
    1328  *     to ignore them if given in conjunction with PCi_ja.
    1329  *
    1330  *   int altlin
    1331  *     (Given) altlin is a bit flag that denotes which of the PCi_ja, CDi_ja
    1332  *     and CROTAi keywords are present in the header:
    1333  *
    1334  *     - Bit 0: PCi_ja is present.
    1335  *
    1336  *     - Bit 1: CDi_ja is present.
    1337  *
    1338  *       Matrix elements in the IRAF convention are equivalent to the product
    1339  *       CDi_ja = CDELTia * PCi_ja, but the defaults differ from that of the
    1340  *       PCi_ja matrix.  If one or more CDi_ja keywords are present then all
    1341  *       unspecified CDi_ja default to zero.  If no CDi_ja (or CROTAi) keywords
    1342  *       are present, then the header is assumed to be in PCi_ja form whether
    1343  *       or not any PCi_ja keywords are present since this results in an
    1344  *       interpretation of CDELTia consistent with the original FITS
    1345  *       specification.
    1346  *
    1347  *       While CDi_ja may not formally co-exist with PCi_ja, it may co-exist
    1348  *       with CDELTia and CROTAi which are to be ignored.
    1349  *
    1350  *     - Bit 2: CROTAi is present.
    1351  *
    1352  *       In the AIPS convention, CROTAi may only be associated with the
    1353  *       latitude axis of a celestial axis pair.  It specifies a rotation in
    1354  *       the image plane that is applied AFTER the CDELTia; any other CROTAi
    1355  *       keywords are ignored.
    1356  *
    1357  *       CROTAi may not formally co-exist with PCi_ja.
    1358  *
    1359  *       CROTAi and CDELTia may formally co-exist with CDi_ja but if so are to
    1360  *       be ignored.
    1361  *
    1362  *     - Bit 3: PCi_ja + CDELTia was derived from CDi_ja by wcspcx().
    1363  *
    1364  *       This bit is set by wcspcx() when it derives PCi_ja and CDELTia from
    1365  *       CDi_ja via an orthonormal decomposition.  In particular, it signals
    1366  *       wcsset() not to replace PCi_ja by a copy of CDi_ja with CDELTia set
    1367  *       to unity.
    1368  *
    1369  *     CDi_ja and CROTAi keywords, if found, are to be stored in the wcsprm::cd
    1370  *     and wcsprm::crota arrays which are dimensioned similarly to wcsprm::pc
    1371  *     and wcsprm::cdelt.  FITS header parsers should use the following
    1372  *     procedure:
    1373  *
    1374  *     - Whenever a PCi_ja keyword is encountered: altlin |= 1;
    1375  *
    1376  *     - Whenever a CDi_ja keyword is encountered: altlin |= 2;
    1377  *
    1378  *     - Whenever a CROTAi keyword is encountered: altlin |= 4;
    1379  *
    1380  *     If none of these bits are set the PCi_ja representation results, i.e.
    1381  *     wcsprm::pc and wcsprm::cdelt will be used as given.
    1382  *
    1383  *     These alternate specifications of the linear transformation matrix are
    1384  *     translated immediately to PCi_ja by wcsset() and are invisible to the
    1385  *     lower-level WCSLIB routines.  In particular, unless bit 3 is also set,
    1386  *     wcsset() resets wcsprm::cdelt to unity if CDi_ja is present (and no
    1387  *     PCi_ja).
    1388  *
    1389  *     If CROTAi are present but none is associated with the latitude axis
    1390  *     (and no PCi_ja or CDi_ja), then wcsset() reverts to a unity PCi_ja
    1391  *     matrix.
    1392  *
    1393  *   int velref
    1394  *     (Given) AIPS velocity code VELREF, refer to spcaips().
    1395  *
    1396  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1397  *     wcsprm::velref is changed.
    1398  *
    1399  *   char alt[4]
    1400  *     (Given, auxiliary) Character code for alternate coordinate descriptions
    1401  *     (i.e. the 'a' in keyword names such as CTYPEia).  This is blank for the
    1402  *     primary coordinate description, or one of the 26 upper-case letters,
    1403  *     A-Z.
    1404  *
    1405  *     An array of four characters is provided for alignment purposes, only the
    1406  *     first is used.
    1407  *
    1408  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1409  *     wcsprm::alt is changed.
    1410  *
    1411  *   int colnum
    1412  *     (Given, auxiliary) Where the coordinate representation is associated
    1413  *     with an image-array column in a FITS binary table, this variable may be
    1414  *     used to record the relevant column number.
    1415  *
    1416  *     It should be set to zero for an image header or pixel list.
    1417  *
    1418  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1419  *     wcsprm::colnum is changed.
    1420  *
    1421  *   int *colax
    1422  *     (Given, auxiliary) Address of the first element of an array of int
    1423  *     recording the column numbers for each axis in a pixel list.
    1424  *
    1425  *     The array elements should be set to zero for an image header or image
    1426  *     array in a binary table.
    1427  *
    1428  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1429  *     wcsprm::colax is changed.
    1430  *
    1431  *   char (*cname)[72]
    1432  *     (Given, auxiliary) The address of the first element of an array of
    1433  *     char[72] containing the coordinate axis names, CNAMEia.
    1434  *
    1435  *     These variables accomodate the longest allowed string-valued FITS
    1436  *     keyword, being limited to 68 characters, plus the null-terminating
    1437  *     character.
    1438  *
    1439  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1440  *     wcsprm::cname is changed.
    1441  *
    1442  *   double *crder
    1443  *     (Given, auxiliary) Address of the first element of an array of double
    1444  *     recording the random error in the coordinate value, CRDERia.
    1445  *
    1446  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1447  *     wcsprm::crder is changed.
    1448  *
    1449  *   double *csyer
    1450  *     (Given, auxiliary) Address of the first element of an array of double
    1451  *     recording the systematic error in the coordinate value, CSYERia.
    1452  *
    1453  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1454  *     wcsprm::csyer is changed.
    1455  *
    1456  *   double *czphs
    1457  *     (Given, auxiliary) Address of the first element of an array of double
    1458  *     recording the time at the zero point of a phase axis, CZPHSia.
    1459  *
    1460  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1461  *     wcsprm::czphs is changed.
    1462  *
    1463  *   double *cperi
    1464  *     (Given, auxiliary) Address of the first element of an array of double
    1465  *     recording the period of a phase axis, CPERIia.
    1466  *
    1467  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1468  *     wcsprm::cperi is changed.
    1469  *
    1470  *   char wcsname[72]
    1471  *     (Given, auxiliary) The name given to the coordinate representation,
    1472  *     WCSNAMEa.  This variable accomodates the longest allowed string-valued
    1473  *     FITS keyword, being limited to 68 characters, plus the null-terminating
    1474  *     character.
    1475  *
    1476  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1477  *     wcsprm::wcsname is changed.
    1478  *
    1479  *   char timesys[72]
    1480  *     (Given, auxiliary) TIMESYS keyvalue, being the time scale (UTC, TAI,
    1481  *     etc.) in which all other time-related auxiliary header values are
    1482  *     recorded.  Also defines the time scale for an image axis with CTYPEia
    1483  *     set to 'TIME'.
    1484  *
    1485  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1486  *     wcsprm::timesys is changed.
    1487  *
    1488  *   char trefpos[72]
    1489  *     (Given, auxiliary) TREFPOS keyvalue, being the location in space where
    1490  *     the recorded time is valid.
    1491  *
    1492  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1493  *     wcsprm::trefpos is changed.
    1494  *
    1495  *   char trefdir[72]
    1496  *     (Given, auxiliary) TREFDIR keyvalue, being the reference direction used
    1497  *     in calculating a pathlength delay.
    1498  *
    1499  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1500  *     wcsprm::trefdir is changed.
    1501  *
    1502  *   char plephem[72]
    1503  *     (Given, auxiliary) PLEPHEM keyvalue, being the Solar System ephemeris
    1504  *     used for calculating a pathlength delay.
    1505  *
    1506  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1507  *     wcsprm::plephem is changed.
    1508  *
    1509  *   char timeunit[72]
    1510  *     (Given, auxiliary) TIMEUNIT keyvalue, being the time units in which
    1511  *     the following header values are expressed: TSTART, TSTOP, TIMEOFFS,
    1512  *     TIMSYER, TIMRDER, TIMEDEL.  It also provides the default value for
    1513  *     CUNITia for time axes.
    1514  *
    1515  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1516  *     wcsprm::timeunit is changed.
    1517  *
    1518  *   char dateref[72]
    1519  *     (Given, auxiliary) DATEREF keyvalue, being the date of a reference epoch
    1520  *     relative to which other time measurements refer.
    1521  *
    1522  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1523  *     wcsprm::dateref is changed.
    1524  *
    1525  *   double mjdref[2]
    1526  *     (Given, auxiliary) MJDREF keyvalue, equivalent to DATEREF expressed as
    1527  *     a Modified Julian Date (MJD = JD - 2400000.5).  The value is given as
    1528  *     the sum of the two-element vector, allowing increased precision.
    1529  *
    1530  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1531  *     wcsprm::mjdref is changed.
    1532  *
    1533  *   double timeoffs
    1534  *     (Given, auxiliary) TIMEOFFS keyvalue, being a time offset, which may be
    1535  *     used, for example, to provide a uniform clock correction for times
    1536  *     referenced to DATEREF.
    1537  *
    1538  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1539  *     wcsprm::timeoffs is changed.
    1540  *
    1541  *   char dateobs[72]
    1542  *     (Given, auxiliary) DATE-OBS keyvalue, being the date at the start of the
    1543  *     observation unless otherwise explained in the DATE-OBS keycomment, in
    1544  *     ISO format, yyyy-mm-ddThh:mm:ss.
    1545  *
    1546  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1547  *     wcsprm::dateobs is changed.
    1548  *
    1549  *   char datebeg[72]
    1550  *     (Given, auxiliary) DATE-BEG keyvalue, being the date at the start of the
    1551  *     observation in ISO format, yyyy-mm-ddThh:mm:ss.
    1552  *
    1553  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1554  *     wcsprm::datebeg is changed.
    1555  *
    1556  *   char dateavg[72]
    1557  *     (Given, auxiliary) DATE-AVG keyvalue, being the date at a representative
    1558  *     mid-point of the observation in ISO format, yyyy-mm-ddThh:mm:ss.
    1559  *
    1560  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1561  *     wcsprm::dateavg is changed.
    1562  *
    1563  *   char dateend[72]
    1564  *     (Given, auxiliary) DATE-END keyvalue, baing the date at the end of the
    1565  *     observation in ISO format, yyyy-mm-ddThh:mm:ss.
    1566  *
    1567  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1568  *     wcsprm::dateend is changed.
    1569  *
    1570  *   double mjdobs
    1571  *     (Given, auxiliary) MJD-OBS keyvalue, equivalent to DATE-OBS expressed
    1572  *     as a Modified Julian Date (MJD = JD - 2400000.5).
    1573  *
    1574  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1575  *     wcsprm::mjdobs is changed.
    1576  *
    1577  *   double mjdbeg
    1578  *     (Given, auxiliary) MJD-BEG keyvalue, equivalent to DATE-BEG expressed
    1579  *     as a Modified Julian Date (MJD = JD - 2400000.5).
    1580  *
    1581  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1582  *     wcsprm::mjdbeg is changed.
    1583  *
    1584  *   double mjdavg
    1585  *     (Given, auxiliary) MJD-AVG keyvalue, equivalent to DATE-AVG expressed
    1586  *     as a Modified Julian Date (MJD = JD - 2400000.5).
    1587  *
    1588  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1589  *     wcsprm::mjdavg is changed.
    1590  *
    1591  *   double mjdend
    1592  *     (Given, auxiliary) MJD-END keyvalue, equivalent to DATE-END expressed
    1593  *     as a Modified Julian Date (MJD = JD - 2400000.5).
    1594  *
    1595  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1596  *     wcsprm::mjdend is changed.
    1597  *
    1598  *   double jepoch
    1599  *     (Given, auxiliary) JEPOCH keyvalue, equivalent to DATE-OBS expressed
    1600  *     as a Julian epoch.
    1601  *
    1602  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1603  *     wcsprm::jepoch is changed.
    1604  *
    1605  *   double bepoch
    1606  *     (Given, auxiliary) BEPOCH keyvalue, equivalent to DATE-OBS expressed
    1607  *     as a Besselian epoch
    1608  *
    1609  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1610  *     wcsprm::bepoch is changed.
    1611  *
    1612  *   double tstart
    1613  *     (Given, auxiliary) TSTART keyvalue, equivalent to DATE-BEG expressed
    1614  *     as a time in units of TIMEUNIT relative to DATEREF+TIMEOFFS.
    1615  *
    1616  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1617  *     wcsprm::tstart is changed.
    1618  *
    1619  *   double tstop
    1620  *     (Given, auxiliary) TSTOP keyvalue, equivalent to DATE-END expressed
    1621  *     as a time in units of TIMEUNIT relative to DATEREF+TIMEOFFS.
    1622  *
    1623  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1624  *     wcsprm::tstop is changed.
    1625  *
    1626  *   double xposure
    1627  *     (Given, auxiliary) XPOSURE keyvalue, being the effective exposure time
    1628  *     in units of TIMEUNIT.
    1629  *
    1630  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1631  *     wcsprm::xposure is changed.
    1632  *
    1633  *   double telapse
    1634  *     (Given, auxiliary) TELAPSE keyvalue, equivalent to the elapsed time
    1635  *     between DATE-BEG and DATE-END, in units of TIMEUNIT.
    1636  *
    1637  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1638  *     wcsprm::telapse is changed.
    1639  *
    1640  *   double timsyer
    1641  *     (Given, auxiliary) TIMSYER keyvalue, being the absolute error of the
    1642  *     time values, in units of TIMEUNIT.
    1643  *
    1644  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1645  *     wcsprm::timsyer is changed.
    1646  *
    1647  *   double timrder
    1648  *     (Given, auxiliary) TIMRDER keyvalue, being the accuracy of time stamps
    1649  *     relative to each other, in units of TIMEUNIT.
    1650  *
    1651  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1652  *     wcsprm::timrder is changed.
    1653  *
    1654  *   double timedel
    1655  *     (Given, auxiliary) TIMEDEL keyvalue, being the resolution of the time
    1656  *     stamps.
    1657  *
    1658  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1659  *     wcsprm::timedel is changed.
    1660  *
    1661  *   double timepixr
    1662  *     (Given, auxiliary) TIMEPIXR keyvalue, being the relative position of the
    1663  *     time stamps in binned time intervals, a value between 0.0 and 1.0.
    1664  *
    1665  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1666  *     wcsprm::timepixr is changed.
    1667  *
    1668  *   double obsgeo[6]
    1669  *     (Given, auxiliary) Location of the observer in a standard terrestrial
    1670  *     reference frame.  The first three give ITRS Cartesian coordinates
    1671  *     OBSGEO-X [m],   OBSGEO-Y [m],   OBSGEO-Z [m], and the second three give
    1672  *     OBSGEO-L [deg], OBSGEO-B [deg], OBSGEO-H [m], which are related through
    1673  *     a standard transformation.
    1674  *
    1675  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1676  *     wcsprm::obsgeo is changed.
    1677  *
    1678  *   char obsorbit[72]
    1679  *     (Given, auxiliary) OBSORBIT keyvalue, being the URI, URL, or name of an
    1680  *     orbit ephemeris file giving spacecraft coordinates relating to TREFPOS.
    1681  *
    1682  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1683  *     wcsprm::obsorbit is changed.
    1684  *
    1685  *   char radesys[72]
    1686  *     (Given, auxiliary) The equatorial or ecliptic coordinate system type,
    1687  *     RADESYSa.
    1688  *
    1689  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1690  *     wcsprm::radesys is changed.
    1691  *
    1692  *   double equinox
    1693  *     (Given, auxiliary) The equinox associated with dynamical equatorial or
    1694  *     ecliptic coordinate systems, EQUINOXa (or EPOCH in older headers).  Not
    1695  *     applicable to ICRS equatorial or ecliptic coordinates.
    1696  *
    1697  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1698  *     wcsprm::equinox is changed.
    1699  *
    1700  *   char specsys[72]
    1701  *     (Given, auxiliary) Spectral reference frame (standard of rest),
    1702  *     SPECSYSa.
    1703  *
    1704  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1705  *     wcsprm::specsys is changed.
    1706  *
    1707  *   char ssysobs[72]
    1708  *     (Given, auxiliary) The spectral reference frame in which there is no
    1709  *     differential variation in the spectral coordinate across the
    1710  *     field-of-view, SSYSOBSa.
    1711  *
    1712  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1713  *     wcsprm::ssysobs is changed.
    1714  *
    1715  *   double velosys
    1716  *     (Given, auxiliary) The relative radial velocity [m/s] between the
    1717  *     observer and the selected standard of rest in the direction of the
    1718  *     celestial reference coordinate, VELOSYSa.
    1719  *
    1720  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1721  *     wcsprm::velosys is changed.
    1722  *
    1723  *   double zsource
    1724  *     (Given, auxiliary) The redshift, ZSOURCEa, of the source.
    1725  *
    1726  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1727  *     wcsprm::zsource is changed.
    1728  *
    1729  *   char ssyssrc[72]
    1730  *     (Given, auxiliary) The spectral reference frame (standard of rest),
    1731  *     SSYSSRCa, in which wcsprm::zsource was measured.
    1732  *
    1733  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1734  *     wcsprm::ssyssrc is changed.
    1735  *
    1736  *   double velangl
    1737  *     (Given, auxiliary) The angle [deg] that should be used to decompose an
    1738  *     observed velocity into radial and transverse components.
    1739  *
    1740  *     It is not necessary to reset the wcsprm struct (via wcsset()) when
    1741  *     wcsprm::velangl is changed.
    1742  *
    1743  *   struct auxprm *aux
    1744  *     (Given, auxiliary) This struct holds auxiliary coordinate system
    1745  *     information of a specialist nature.  While these parameters may be
    1746  *     widely recognized within particular fields of astronomy, they differ
    1747  *     from the above auxiliary parameters in not being defined by any of the
    1748  *     FITS WCS standards.  Collecting them together in a separate struct that
    1749  *     is allocated only when required helps to control bloat in the size of
    1750  *     the wcsprm struct.
    1751  *
    1752  *   int ntab
    1753  *     (Given) See wcsprm::tab.
    1754  *
    1755  *   int nwtb
    1756  *     (Given) See wcsprm::wtb.
    1757  *
    1758  *   struct tabprm *tab
    1759  *     (Given) Address of the first element of an array of ntab tabprm structs
    1760  *     for which memory has been allocated.  These are used to store tabular
    1761  *     transformation parameters.
    1762  *
    1763  *     Although technically wcsprm::ntab and tab are "given", they will
    1764  *     normally be set by invoking wcstab(), whether directly or indirectly.
    1765  *
    1766  *     The tabprm structs contain some members that must be supplied and others
    1767  *     that are derived.  The information to be supplied comes primarily from
    1768  *     arrays stored in one or more FITS binary table extensions.  These
    1769  *     arrays, referred to here as "wcstab arrays", are themselves located by
    1770  *     parameters stored in the FITS image header.
    1771  *
    1772  *   struct wtbarr *wtb
    1773  *     (Given) Address of the first element of an array of nwtb wtbarr structs
    1774  *     for which memory has been allocated.  These are used in extracting
    1775  *     wcstab arrays from a FITS binary table.
    1776  *
    1777  *     Although technically wcsprm::nwtb and wtb are "given", they will
    1778  *     normally be set by invoking wcstab(), whether directly or indirectly.
    1779  *
    1780  *   char lngtyp[8]
    1781  *     (Returned) Four-character WCS celestial longitude and ...
    1782  *   char lattyp[8]
    1783  *     (Returned) ... latitude axis types. e.g. "RA", "DEC", "GLON", "GLAT",
    1784  *     etc. extracted from 'RA--', 'DEC-', 'GLON', 'GLAT', etc. in the first
    1785  *     four characters of CTYPEia but with trailing dashes removed.  (Declared
    1786  *     as char[8] for alignment reasons.)
    1787  *
    1788  *   int lng
    1789  *     (Returned) Index for the longitude coordinate, and ...
    1790  *   int lat
    1791  *     (Returned) ... index for the latitude coordinate, and ...
    1792  *   int spec
    1793  *     (Returned) ... index for the spectral coordinate, and ...
    1794  *   int time
    1795  *     (Returned) ... index for the time coordinate in the imgcrd[][] and
    1796  *     world[][] arrays in the API of wcsp2s(), wcss2p() and wcsmix().
    1797  *
    1798  *     These may also serve as indices into the pixcrd[][] array provided that
    1799  *     the PCi_ja matrix does not transpose axes.
    1800  *
    1801  *   int cubeface
    1802  *     (Returned) Index into the pixcrd[][] array for the CUBEFACE axis.  This
    1803  *     is used for quadcube projections where the cube faces are stored on a
    1804  *     separate axis (see wcs.h).
    1805  *
    1806  *   int *types
    1807  *     (Returned) Address of the first element of an array of int containing a
    1808  *     four-digit type code for each axis.
    1809  *
    1810  *     - First digit (i.e. 1000s):
    1811  *       - 0: Non-specific coordinate type.
    1812  *       - 1: Stokes coordinate.
    1813  *       - 2: Celestial coordinate (including CUBEFACE).
    1814  *       - 3: Spectral coordinate.
    1815  *       - 4: Time coordinate.
    1816  *
    1817  *     - Second digit (i.e. 100s):
    1818  *       - 0: Linear axis.
    1819  *       - 1: Quantized axis (STOKES, CUBEFACE).
    1820  *       - 2: Non-linear celestial axis.
    1821  *       - 3: Non-linear spectral axis.
    1822  *       - 4: Logarithmic axis.
    1823  *       - 5: Tabular axis.
    1824  *
    1825  *     - Third digit (i.e. 10s):
    1826  *       - 0: Group number, e.g. lookup table number, being an index into the
    1827  *            tabprm array (see above).
    1828  *
    1829  *     - The fourth digit is used as a qualifier depending on the axis type.
    1830  *
    1831  *       - For celestial axes:
    1832  *         - 0: Longitude coordinate.
    1833  *         - 1: Latitude coordinate.
    1834  *         - 2: CUBEFACE number.
    1835  *
    1836  *       - For lookup tables: the axis number in a multidimensional table.
    1837  *
    1838  *     CTYPEia in "4-3" form with unrecognized algorithm code will have its
    1839  *     type set to -1 and generate an error.
    1840  *
    1841  *   struct linprm lin
    1842  *     (Returned) Linear transformation parameters (usage is described in the
    1843  *     prologue to lin.h).
    1844  *
    1845  *   struct celprm cel
    1846  *     (Returned) Celestial transformation parameters (usage is described in
    1847  *     the prologue to cel.h).
    1848  *
    1849  *   struct spcprm spc
    1850  *     (Returned) Spectral transformation parameters (usage is described in the
    1851  *     prologue to spc.h).
    1852  *
    1853  *   struct wcserr *err
    1854  *     (Returned) If enabled, when an error status is returned, this struct
    1855  *     contains detailed information about the error, see wcserr_enable().
    1856  *
    1857  *   int m_flag
    1858  *     (For internal use only.)
    1859  *   int m_naxis
    1860  *     (For internal use only.)
    1861  *   double *m_crpix
    1862  *     (For internal use only.)
    1863  *   double *m_pc
    1864  *     (For internal use only.)
    1865  *   double *m_cdelt
    1866  *     (For internal use only.)
    1867  *   double *m_crval
    1868  *     (For internal use only.)
    1869  *   char (*m_cunit)[72]
    1870  *     (For internal use only.)
    1871  *   char (*m_ctype)[72]
    1872  *     (For internal use only.)
    1873  *   struct pvcard *m_pv
    1874  *     (For internal use only.)
    1875  *   struct pscard *m_ps
    1876  *     (For internal use only.)
    1877  *   double *m_cd
    1878  *     (For internal use only.)
    1879  *   double *m_crota
    1880  *     (For internal use only.)
    1881  *   int *m_colax
    1882  *     (For internal use only.)
    1883  *   char (*m_cname)[72]
    1884  *     (For internal use only.)
    1885  *   double *m_crder
    1886  *     (For internal use only.)
    1887  *   double *m_csyer
    1888  *     (For internal use only.)
    1889  *   double *m_czphs
    1890  *     (For internal use only.)
    1891  *   double *m_cperi
    1892  *     (For internal use only.)
    1893  *   struct tabprm *m_tab
    1894  *     (For internal use only.)
    1895  *   struct wtbarr *m_wtb
    1896  *     (For internal use only.)
    1897  *
    1898  *
    1899  * pvcard struct - Store for PVi_ma keyrecords
    1900  * -------------------------------------------
    1901  * The pvcard struct is used to pass the parsed contents of PVi_ma keyrecords
    1902  * to wcsset() via the wcsprm struct.
    1903  *
    1904  * All members of this struct are to be set by the user.
    1905  *
    1906  *   int i
    1907  *     (Given) Axis number (1-relative), as in the FITS PVi_ma keyword.  If
    1908  *     i == 0, wcsset() will replace it with the latitude axis number.
    1909  *
    1910  *   int m
    1911  *     (Given) Parameter number (non-negative), as in the FITS PVi_ma keyword.
    1912  *
    1913  *   double value
    1914  *     (Given) Parameter value.
    1915  *
    1916  *
    1917  * pscard struct - Store for PSi_ma keyrecords
    1918  * -------------------------------------------
    1919  * The pscard struct is used to pass the parsed contents of PSi_ma keyrecords
    1920  * to wcsset() via the wcsprm struct.
    1921  *
    1922  * All members of this struct are to be set by the user.
    1923  *
    1924  *   int i
    1925  *     (Given) Axis number (1-relative), as in the FITS PSi_ma keyword.
    1926  *
    1927  *   int m
    1928  *     (Given) Parameter number (non-negative), as in the FITS PSi_ma keyword.
    1929  *
    1930  *   char value[72]
    1931  *     (Given) Parameter value.
    1932  *
    1933  *
    1934  * auxprm struct - Additional auxiliary parameters
    1935  * -----------------------------------------------
    1936  * The auxprm struct holds auxiliary coordinate system information of a
    1937  * specialist nature.  It is anticipated that this struct will expand in future
    1938  * to accomodate additional parameters.
    1939  *
    1940  * All members of this struct are to be set by the user.
    1941  *
    1942  *   double rsun_ref
    1943  *     (Given, auxiliary) Reference radius of the Sun used in coordinate
    1944  *     calculations (m).
    1945  *
    1946  *   double dsun_obs
    1947  *     (Given, auxiliary) Distance between the centre of the Sun and the
    1948  *     observer (m).
    1949  *
    1950  *   double crln_obs
    1951  *     (Given, auxiliary) Carrington heliographic longitude of the observer
    1952  *     (deg).
    1953  *
    1954  *   double hgln_obs
    1955  *     (Given, auxiliary) Stonyhurst heliographic longitude of the observer
    1956  *     (deg).
    1957  *
    1958  *   double hglt_obs
    1959  *     (Given, auxiliary) Heliographic latitude (Carrington or Stonyhurst) of
    1960  *     the observer (deg).
    1961  *
    1962  *   double a_radius
    1963  *     Length of the semi-major axis of a triaxial ellipsoid approximating the
    1964  *     shape of a body (e.g. planet) in the solar system (m).
    1965  *
    1966  *   double b_radius
    1967  *     Length of the intermediate axis, normal to the semi-major and semi-minor
    1968  *     axes, of a triaxial ellipsoid approximating the shape of a body (m).
    1969  *
    1970  *   double c_radius
    1971  *     Length of the semi-minor axis, normal to the semi-major axis, of a
    1972  *     triaxial ellipsoid approximating the shape of a body (m).
    1973  *
    1974  *   double blon_obs
    1975  *     Bodycentric longitude of the observer in the coordinate system fixed to
    1976  *     the planet or other solar system body (deg, in range 0 to 360).
    1977  *
    1978  *   double blat_obs
    1979  *     Bodycentric latitude of the observer in the coordinate system fixed to
    1980  *     the planet or other solar system body (deg).
    1981  *
    1982  *   double bdis_obs
    1983  *     Bodycentric distance of the observer (m).
    1984  *
    1985  * Global variable: const char *wcs_errmsg[] - Status return messages
    1986  * ------------------------------------------------------------------
    1987  * Error messages to match the status value returned from each function.
    1988  *
    1989  *===========================================================================*/
    1990  
    1991  #ifndef WCSLIB_WCS
    1992  #define WCSLIB_WCS
    1993  
    1994  #include "lin.h"
    1995  #include "cel.h"
    1996  #include "spc.h"
    1997  
    1998  #ifdef __cplusplus
    1999  extern "C" {
    2000  #define wtbarr wtbarr_s		// See prologue of wtbarr.h.
    2001  #endif
    2002  
    2003  #define WCSSUB_LONGITUDE 0x1001
    2004  #define WCSSUB_LATITUDE  0x1002
    2005  #define WCSSUB_CUBEFACE  0x1004
    2006  #define WCSSUB_CELESTIAL 0x1007
    2007  #define WCSSUB_SPECTRAL  0x1008
    2008  #define WCSSUB_STOKES    0x1010
    2009  #define WCSSUB_TIME      0x1020
    2010  
    2011  
    2012  #define WCSCOMPARE_ANCILLARY 0x0001
    2013  #define WCSCOMPARE_TILING    0x0002
    2014  #define WCSCOMPARE_CRPIX     0x0004
    2015  
    2016  
    2017  extern const char *wcs_errmsg[];
    2018  
    2019  enum wcs_errmsg_enum {
    2020    WCSERR_SUCCESS         =  0,	// Success.
    2021    WCSERR_NULL_POINTER    =  1,	// Null wcsprm pointer passed.
    2022    WCSERR_MEMORY          =  2,	// Memory allocation failed.
    2023    WCSERR_SINGULAR_MTX    =  3,	// Linear transformation matrix is singular.
    2024    WCSERR_BAD_CTYPE       =  4,	// Inconsistent or unrecognized coordinate
    2025  				// axis type.
    2026    WCSERR_BAD_PARAM       =  5,	// Invalid parameter value.
    2027    WCSERR_BAD_COORD_TRANS =  6,	// Unrecognized coordinate transformation
    2028  				// parameter.
    2029    WCSERR_ILL_COORD_TRANS =  7,	// Ill-conditioned coordinate transformation
    2030  				// parameter.
    2031    WCSERR_BAD_PIX         =  8,	// One or more of the pixel coordinates were
    2032  				// invalid.
    2033    WCSERR_BAD_WORLD       =  9,	// One or more of the world coordinates were
    2034  				// invalid.
    2035    WCSERR_BAD_WORLD_COORD = 10,	// Invalid world coordinate.
    2036    WCSERR_NO_SOLUTION     = 11,	// No solution found in the specified
    2037  				// interval.
    2038    WCSERR_BAD_SUBIMAGE    = 12,	// Invalid subimage specification.
    2039    WCSERR_NON_SEPARABLE   = 13,	// Non-separable subimage coordinate system.
    2040    WCSERR_UNSET           = 14 	// wcsprm struct is unset.
    2041  };
    2042  
    2043  
    2044  // Struct used for storing PVi_ma keywords.
    2045  struct pvcard {
    2046    int i;			// Axis number, as in PVi_ma (1-relative).
    2047    int m;			// Parameter number, ditto  (0-relative).
    2048    double value;			// Parameter value.
    2049  };
    2050  
    2051  // Size of the pvcard struct in int units, used by the Fortran wrappers.
    2052  #define PVLEN (sizeof(struct pvcard)/sizeof(int))
    2053  
    2054  // Struct used for storing PSi_ma keywords.
    2055  struct pscard {
    2056    int i;			// Axis number, as in PSi_ma (1-relative).
    2057    int m;			// Parameter number, ditto  (0-relative).
    2058    char value[72];		// Parameter value.
    2059  };
    2060  
    2061  // Size of the pscard struct in int units, used by the Fortran wrappers.
    2062  #define PSLEN (sizeof(struct pscard)/sizeof(int))
    2063  
    2064  // Struct used to hold additional auxiliary parameters.
    2065  struct auxprm {
    2066    double rsun_ref;              // Solar radius.
    2067    double dsun_obs;              // Distance from Sun centre to observer.
    2068    double crln_obs;              // Carrington heliographic lng of observer.
    2069    double hgln_obs;              // Stonyhurst heliographic lng of observer.
    2070    double hglt_obs;              // Heliographic latitude of observer.
    2071  
    2072    double a_radius;              // Semi-major axis of solar system body.
    2073    double b_radius;              // Semi-intermediate axis of solar system body.
    2074    double c_radius;              // Semi-minor axis of solar system body.
    2075    double blon_obs;              // Bodycentric longitude of observer.
    2076    double blat_obs;              // Bodycentric latitude of observer.
    2077    double bdis_obs;              // Bodycentric distance of observer.
    2078    double dummy[2];              // Reserved for future use.
    2079  };
    2080  
    2081  // Size of the auxprm struct in int units, used by the Fortran wrappers.
    2082  #define AUXLEN (sizeof(struct auxprm)/sizeof(int))
    2083  
    2084  
    2085  struct wcsprm {
    2086    // Initialization flag (see the prologue above).
    2087    //--------------------------------------------------------------------------
    2088    int    flag;			// Set to zero to force initialization.
    2089  
    2090    // FITS header keyvalues to be provided (see the prologue above).
    2091    //--------------------------------------------------------------------------
    2092    int    naxis;			// Number of axes (pixel and coordinate).
    2093    double *crpix;		// CRPIXja keyvalues for each pixel axis.
    2094    double *pc;			// PCi_ja  linear transformation matrix.
    2095    double *cdelt;		// CDELTia keyvalues for each coord axis.
    2096    double *crval;		// CRVALia keyvalues for each coord axis.
    2097  
    2098    char   (*cunit)[72];		// CUNITia keyvalues for each coord axis.
    2099    char   (*ctype)[72];		// CTYPEia keyvalues for each coord axis.
    2100  
    2101    double lonpole;		// LONPOLEa keyvalue.
    2102    double latpole;		// LATPOLEa keyvalue.
    2103  
    2104    double restfrq;		// RESTFRQa keyvalue.
    2105    double restwav;		// RESTWAVa keyvalue.
    2106  
    2107    int    npv;			// Number of PVi_ma keywords, and the
    2108    int    npvmax;		// number for which space was allocated.
    2109    struct pvcard *pv;		// PVi_ma keywords for each i and m.
    2110  
    2111    int    nps;			// Number of PSi_ma keywords, and the
    2112    int    npsmax;		// number for which space was allocated.
    2113    struct pscard *ps;		// PSi_ma keywords for each i and m.
    2114  
    2115    // Alternative header keyvalues (see the prologue above).
    2116    //--------------------------------------------------------------------------
    2117    double *cd;			// CDi_ja linear transformation matrix.
    2118    double *crota;		// CROTAi keyvalues for each coord axis.
    2119    int    altlin;		// Alternative representations
    2120  				//   Bit 0: PCi_ja is present,
    2121  				//   Bit 1: CDi_ja is present,
    2122  				//   Bit 2: CROTAi is present.
    2123    int    velref;		// AIPS velocity code, VELREF.
    2124  
    2125    // Auxiliary coordinate system information of a general nature.  Not
    2126    // used by WCSLIB.  Refer to the prologue comments above for a brief
    2127    // explanation of these values.
    2128    char   alt[4];
    2129    int    colnum;
    2130    int    *colax;
    2131  				// Auxiliary coordinate axis information.
    2132    char   (*cname)[72];
    2133    double *crder;
    2134    double *csyer;
    2135    double *czphs;
    2136    double *cperi;
    2137  
    2138    char   wcsname[72];
    2139  				// Time reference system and measurement.
    2140    char   timesys[72], trefpos[72], trefdir[72], plephem[72];
    2141    char   timeunit[72];
    2142    char   dateref[72];
    2143    double mjdref[2];
    2144    double timeoffs;
    2145  				// Data timestamps and durations.
    2146    char   dateobs[72], datebeg[72], dateavg[72], dateend[72];
    2147    double mjdobs, mjdbeg, mjdavg, mjdend;
    2148    double jepoch, bepoch;
    2149    double tstart, tstop;
    2150    double xposure, telapse;
    2151  				// Timing accuracy.
    2152    double timsyer, timrder;
    2153    double timedel, timepixr;
    2154  				// Spatial & celestial reference frame.
    2155    double obsgeo[6];
    2156    char   obsorbit[72];
    2157    char   radesys[72];
    2158    double equinox;
    2159    char   specsys[72];
    2160    char   ssysobs[72];
    2161    double velosys;
    2162    double zsource;
    2163    char   ssyssrc[72];
    2164    double velangl;
    2165  
    2166    // Additional auxiliary coordinate system information of a specialist
    2167    // nature.  Not used by WCSLIB.  Refer to the prologue comments above.
    2168    struct auxprm *aux;
    2169  
    2170    // Coordinate lookup tables (see the prologue above).
    2171    //--------------------------------------------------------------------------
    2172    int    ntab;			// Number of separate tables.
    2173    int    nwtb;			// Number of wtbarr structs.
    2174    struct tabprm *tab;		// Tabular transformation parameters.
    2175    struct wtbarr *wtb;		// Array of wtbarr structs.
    2176  
    2177    //--------------------------------------------------------------------------
    2178    // Information derived from the FITS header keyvalues by wcsset().
    2179    //--------------------------------------------------------------------------
    2180    char   lngtyp[8], lattyp[8];	// Celestial axis types, e.g. RA, DEC.
    2181    int    lng, lat, spec, time;	// Longitude, latitude, spectral, and time
    2182  				// axis indices (0-relative).
    2183    int    cubeface;		// True if there is a CUBEFACE axis.
    2184    int    dummy;			// Dummy for alignment purposes.
    2185    int    *types;		// Coordinate type codes for each axis.
    2186  
    2187    struct linprm lin;		//    Linear transformation parameters.
    2188    struct celprm cel;		// Celestial transformation parameters.
    2189    struct spcprm spc;		//  Spectral transformation parameters.
    2190  
    2191    //--------------------------------------------------------------------------
    2192    //             THE REMAINDER OF THE WCSPRM STRUCT IS PRIVATE.
    2193    //--------------------------------------------------------------------------
    2194  
    2195    // Error handling, if enabled.
    2196    //--------------------------------------------------------------------------
    2197    struct wcserr *err;
    2198  
    2199    // Memory management.
    2200    //--------------------------------------------------------------------------
    2201    int    m_flag, m_naxis;
    2202    double *m_crpix, *m_pc, *m_cdelt, *m_crval;
    2203    char  (*m_cunit)[72], (*m_ctype)[72];
    2204    struct pvcard *m_pv;
    2205    struct pscard *m_ps;
    2206    double *m_cd, *m_crota;
    2207    int    *m_colax;
    2208    char  (*m_cname)[72];
    2209    double *m_crder, *m_csyer, *m_czphs, *m_cperi;
    2210    struct auxprm *m_aux;
    2211    struct tabprm *m_tab;
    2212    struct wtbarr *m_wtb;
    2213  };
    2214  
    2215  // Size of the wcsprm struct in int units, used by the Fortran wrappers.
    2216  #define WCSLEN (sizeof(struct wcsprm)/sizeof(int))
    2217  
    2218  
    2219  int wcsnpv(int n);
    2220  
    2221  int wcsnps(int n);
    2222  
    2223  int wcsini(int alloc, int naxis, struct wcsprm *wcs);
    2224  
    2225  int wcsinit(int alloc, int naxis, struct wcsprm *wcs, int npvmax, int npsmax,
    2226              int ndpmax);
    2227  
    2228  int wcsauxi(int alloc, struct wcsprm *wcs);
    2229  
    2230  int wcssub(int alloc, const struct wcsprm *wcssrc, int *nsub, int axes[],
    2231             struct wcsprm *wcsdst);
    2232  
    2233  int wcscompare(int cmp, double tol, const struct wcsprm *wcs1,
    2234                 const struct wcsprm *wcs2, int *equal);
    2235  
    2236  int wcsfree(struct wcsprm *wcs);
    2237  
    2238  int wcstrim(struct wcsprm *wcs);
    2239  
    2240  int wcssize(const struct wcsprm *wcs, int sizes[2]);
    2241  
    2242  int auxsize(const struct auxprm *aux, int sizes[2]);
    2243  
    2244  int wcsprt(const struct wcsprm *wcs);
    2245  
    2246  int wcsperr(const struct wcsprm *wcs, const char *prefix);
    2247  
    2248  int wcsbchk(struct wcsprm *wcs, int bounds);
    2249  
    2250  int wcsset(struct wcsprm *wcs);
    2251  
    2252  int wcsp2s(struct wcsprm *wcs, int ncoord, int nelem, const double pixcrd[],
    2253             double imgcrd[], double phi[], double theta[], double world[],
    2254             int stat[]);
    2255  
    2256  int wcss2p(struct wcsprm *wcs, int ncoord, int nelem, const double world[],
    2257             double phi[], double theta[], double imgcrd[], double pixcrd[],
    2258             int stat[]);
    2259  
    2260  int wcsmix(struct wcsprm *wcs, int mixpix, int mixcel, const double vspan[2],
    2261             double vstep, int viter, double world[], double phi[],
    2262             double theta[], double imgcrd[], double pixcrd[]);
    2263  
    2264  int wcsccs(struct wcsprm *wcs, double lng2p1, double lat2p1, double lng1p2,
    2265             const char *clng, const char *clat, const char *radesys,
    2266             double equinox, const char *alt);
    2267  
    2268  int wcssptr(struct wcsprm *wcs, int *i, char ctype[9]);
    2269  
    2270  const char* wcslib_version(int vers[3]);
    2271  
    2272  // Defined mainly for backwards compatibility, use wcssub() instead.
    2273  #define wcscopy(alloc, wcssrc, wcsdst) wcssub(alloc, wcssrc, 0x0, 0x0, wcsdst)
    2274  
    2275  
    2276  // Deprecated.
    2277  #define wcsini_errmsg wcs_errmsg
    2278  #define wcssub_errmsg wcs_errmsg
    2279  #define wcscopy_errmsg wcs_errmsg
    2280  #define wcsfree_errmsg wcs_errmsg
    2281  #define wcsprt_errmsg wcs_errmsg
    2282  #define wcsset_errmsg wcs_errmsg
    2283  #define wcsp2s_errmsg wcs_errmsg
    2284  #define wcss2p_errmsg wcs_errmsg
    2285  #define wcsmix_errmsg wcs_errmsg
    2286  
    2287  #ifdef __cplusplus
    2288  #undef wtbarr
    2289  }
    2290  #endif
    2291  
    2292  #endif // WCSLIB_WCS