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