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