wcslib (8.2.2)
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "https://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml" lang="en-US">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=11"/>
<meta name="generator" content="Doxygen 1.9.8"/>
<meta name="viewport" content="width=device-width, initial-scale=1"/>
<title>WCSLIB: wcsfix.h Source File</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="search/search.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="search/searchdata.js"></script>
<script type="text/javascript" src="search/search.js"></script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
<tbody>
<tr id="projectrow">
<td id="projectalign">
<div id="projectname">WCSLIB<span id="projectnumber"> 8.2.2</span>
</div>
</td>
</tr>
</tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.9.8 -->
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
var searchBox = new SearchBox("searchBox", "search/",'.html');
/* @license-end */
</script>
<script type="text/javascript" src="menudata.js"></script>
<script type="text/javascript" src="menu.js"></script>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
$(function() {
initMenu('',true,false,'search.php','Search');
$(document).ready(function() { init_search(); });
});
/* @license-end */
</script>
<div id="main-nav"></div>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
$(document).ready(function() { init_codefold(0); });
/* @license-end */
</script>
<!-- window showing the filter options -->
<div id="MSearchSelectWindow"
onmouseover="return searchBox.OnSearchSelectShow()"
onmouseout="return searchBox.OnSearchSelectHide()"
onkeydown="return searchBox.OnSearchSelectKey(event)">
</div>
<!-- iframe showing the search results (closed by default) -->
<div id="MSearchResultsWindow">
<div id="MSearchResults">
<div class="SRPage">
<div id="SRIndex">
<div id="SRResults"></div>
<div class="SRStatus" id="Loading">Loading...</div>
<div class="SRStatus" id="Searching">Searching...</div>
<div class="SRStatus" id="NoMatches">No Matches</div>
</div>
</div>
</div>
</div>
<div id="nav-path" class="navpath">
<ul>
<li class="navelem"><a class="el" href="dir_af99bb3c152a306abd27951285ad1127.html">C</a></li> </ul>
</div>
</div><!-- top -->
<div class="header">
<div class="headertitle"><div class="title">wcsfix.h</div></div>
</div><!--header-->
<div class="contents">
<a href="wcsfix_8h.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a id="l00001" name="l00001"></a><span class="lineno"> 1</span><span class="comment">/*============================================================================</span></div>
<div class="line"><a id="l00002" name="l00002"></a><span class="lineno"> 2</span><span class="comment"> WCSLIB 8.2 - an implementation of the FITS WCS standard.</span></div>
<div class="line"><a id="l00003" name="l00003"></a><span class="lineno"> 3</span><span class="comment"> Copyright (C) 1995-2023, Mark Calabretta</span></div>
<div class="line"><a id="l00004" name="l00004"></a><span class="lineno"> 4</span><span class="comment"></span> </div>
<div class="line"><a id="l00005" name="l00005"></a><span class="lineno"> 5</span><span class="comment"> This file is part of WCSLIB.</span></div>
<div class="line"><a id="l00006" name="l00006"></a><span class="lineno"> 6</span><span class="comment"></span> </div>
<div class="line"><a id="l00007" name="l00007"></a><span class="lineno"> 7</span><span class="comment"> WCSLIB is free software: you can redistribute it and/or modify it under the</span></div>
<div class="line"><a id="l00008" name="l00008"></a><span class="lineno"> 8</span><span class="comment"> terms of the GNU Lesser General Public License as published by the Free</span></div>
<div class="line"><a id="l00009" name="l00009"></a><span class="lineno"> 9</span><span class="comment"> Software Foundation, either version 3 of the License, or (at your option)</span></div>
<div class="line"><a id="l00010" name="l00010"></a><span class="lineno"> 10</span><span class="comment"> any later version.</span></div>
<div class="line"><a id="l00011" name="l00011"></a><span class="lineno"> 11</span><span class="comment"></span> </div>
<div class="line"><a id="l00012" name="l00012"></a><span class="lineno"> 12</span><span class="comment"> WCSLIB is distributed in the hope that it will be useful, but WITHOUT ANY</span></div>
<div class="line"><a id="l00013" name="l00013"></a><span class="lineno"> 13</span><span class="comment"> WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS</span></div>
<div class="line"><a id="l00014" name="l00014"></a><span class="lineno"> 14</span><span class="comment"> FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for</span></div>
<div class="line"><a id="l00015" name="l00015"></a><span class="lineno"> 15</span><span class="comment"> more details.</span></div>
<div class="line"><a id="l00016" name="l00016"></a><span class="lineno"> 16</span><span class="comment"></span> </div>
<div class="line"><a id="l00017" name="l00017"></a><span class="lineno"> 17</span><span class="comment"> You should have received a copy of the GNU Lesser General Public License</span></div>
<div class="line"><a id="l00018" name="l00018"></a><span class="lineno"> 18</span><span class="comment"> along with WCSLIB. If not, see http://www.gnu.org/licenses.</span></div>
<div class="line"><a id="l00019" name="l00019"></a><span class="lineno"> 19</span><span class="comment"></span> </div>
<div class="line"><a id="l00020" name="l00020"></a><span class="lineno"> 20</span><span class="comment"> Author: Mark Calabretta, Australia Telescope National Facility, CSIRO.</span></div>
<div class="line"><a id="l00021" name="l00021"></a><span class="lineno"> 21</span><span class="comment"> http://www.atnf.csiro.au/people/Mark.Calabretta</span></div>
<div class="line"><a id="l00022" name="l00022"></a><span class="lineno"> 22</span><span class="comment"> $Id: wcsfix.h,v 8.2.1.1 2023/11/16 10:05:57 mcalabre Exp mcalabre $</span></div>
<div class="line"><a id="l00023" name="l00023"></a><span class="lineno"> 23</span><span class="comment">*=============================================================================</span></div>
<div class="line"><a id="l00024" name="l00024"></a><span class="lineno"> 24</span><span class="comment">*</span></div>
<div class="line"><a id="l00025" name="l00025"></a><span class="lineno"> 25</span><span class="comment">* WCSLIB 8.2 - C routines that implement the FITS World Coordinate System</span></div>
<div class="line"><a id="l00026" name="l00026"></a><span class="lineno"> 26</span><span class="comment">* (WCS) standard. Refer to the README file provided with WCSLIB for an</span></div>
<div class="line"><a id="l00027" name="l00027"></a><span class="lineno"> 27</span><span class="comment">* overview of the library.</span></div>
<div class="line"><a id="l00028" name="l00028"></a><span class="lineno"> 28</span><span class="comment">*</span></div>
<div class="line"><a id="l00029" name="l00029"></a><span class="lineno"> 29</span><span class="comment">*</span></div>
<div class="line"><a id="l00030" name="l00030"></a><span class="lineno"> 30</span><span class="comment">* Summary of the wcsfix routines</span></div>
<div class="line"><a id="l00031" name="l00031"></a><span class="lineno"> 31</span><span class="comment">* ------------------------------</span></div>
<div class="line"><a id="l00032" name="l00032"></a><span class="lineno"> 32</span><span class="comment">* Routines in this suite identify and translate various forms of construct</span></div>
<div class="line"><a id="l00033" name="l00033"></a><span class="lineno"> 33</span><span class="comment">* known to occur in FITS headers that violate the FITS World Coordinate System</span></div>
<div class="line"><a id="l00034" name="l00034"></a><span class="lineno"> 34</span><span class="comment">* (WCS) standard described in</span></div>
<div class="line"><a id="l00035" name="l00035"></a><span class="lineno"> 35</span><span class="comment">*</span></div>
<div class="line"><a id="l00036" name="l00036"></a><span class="lineno"> 36</span><span class="comment">= "Representations of world coordinates in FITS",</span></div>
<div class="line"><a id="l00037" name="l00037"></a><span class="lineno"> 37</span><span class="comment">= Greisen, E.W., & Calabretta, M.R. 2002, A&A, 395, 1061 (WCS Paper I)</span></div>
<div class="line"><a id="l00038" name="l00038"></a><span class="lineno"> 38</span><span class="comment">=</span></div>
<div class="line"><a id="l00039" name="l00039"></a><span class="lineno"> 39</span><span class="comment">= "Representations of celestial coordinates in FITS",</span></div>
<div class="line"><a id="l00040" name="l00040"></a><span class="lineno"> 40</span><span class="comment">= Calabretta, M.R., & Greisen, E.W. 2002, A&A, 395, 1077 (WCS Paper II)</span></div>
<div class="line"><a id="l00041" name="l00041"></a><span class="lineno"> 41</span><span class="comment">=</span></div>
<div class="line"><a id="l00042" name="l00042"></a><span class="lineno"> 42</span><span class="comment">= "Representations of spectral coordinates in FITS",</span></div>
<div class="line"><a id="l00043" name="l00043"></a><span class="lineno"> 43</span><span class="comment">= Greisen, E.W., Calabretta, M.R., Valdes, F.G., & Allen, S.L.</span></div>
<div class="line"><a id="l00044" name="l00044"></a><span class="lineno"> 44</span><span class="comment">= 2006, A&A, 446, 747 (WCS Paper III)</span></div>
<div class="line"><a id="l00045" name="l00045"></a><span class="lineno"> 45</span><span class="comment">=</span></div>
<div class="line"><a id="l00046" name="l00046"></a><span class="lineno"> 46</span><span class="comment">= "Representations of time coordinates in FITS -</span></div>
<div class="line"><a id="l00047" name="l00047"></a><span class="lineno"> 47</span><span class="comment">= Time and relative dimension in space",</span></div>
<div class="line"><a id="l00048" name="l00048"></a><span class="lineno"> 48</span><span class="comment">= Rots, A.H., Bunclark, P.S., Calabretta, M.R., Allen, S.L.,</span></div>
<div class="line"><a id="l00049" name="l00049"></a><span class="lineno"> 49</span><span class="comment">= Manchester, R.N., & Thompson, W.T. 2015, A&A, 574, A36 (WCS Paper VII)</span></div>
<div class="line"><a id="l00050" name="l00050"></a><span class="lineno"> 50</span><span class="comment">*</span></div>
<div class="line"><a id="l00051" name="l00051"></a><span class="lineno"> 51</span><span class="comment">* Repairs effected by these routines range from the translation of</span></div>
<div class="line"><a id="l00052" name="l00052"></a><span class="lineno"> 52</span><span class="comment">* non-standard values for standard WCS keywords, to the repair of malformed</span></div>
<div class="line"><a id="l00053" name="l00053"></a><span class="lineno"> 53</span><span class="comment">* coordinate representations. Some routines are also provided to check the</span></div>
<div class="line"><a id="l00054" name="l00054"></a><span class="lineno"> 54</span><span class="comment">* consistency of pairs of keyvalues that define the same measure in two</span></div>
<div class="line"><a id="l00055" name="l00055"></a><span class="lineno"> 55</span><span class="comment">* different ways, for example, as a date and an MJD.</span></div>
<div class="line"><a id="l00056" name="l00056"></a><span class="lineno"> 56</span><span class="comment">*</span></div>
<div class="line"><a id="l00057" name="l00057"></a><span class="lineno"> 57</span><span class="comment">* A separate routine, wcspcx(), "regularizes" the linear transformation matrix</span></div>
<div class="line"><a id="l00058" name="l00058"></a><span class="lineno"> 58</span><span class="comment">* component (PCi_j) of the coordinate transformation to make it more human-</span></div>
<div class="line"><a id="l00059" name="l00059"></a><span class="lineno"> 59</span><span class="comment">* readable. Where a coordinate description was constructed from CDi_j, it</span></div>
<div class="line"><a id="l00060" name="l00060"></a><span class="lineno"> 60</span><span class="comment">* decomposes it into PCi_j + CDELTi in a meaningful way. Optionally, it can</span></div>
<div class="line"><a id="l00061" name="l00061"></a><span class="lineno"> 61</span><span class="comment">* also diagonalize the PCi_j matrix (as far as possible), i.e. undo a</span></div>
<div class="line"><a id="l00062" name="l00062"></a><span class="lineno"> 62</span><span class="comment">* transposition of axes in the intermediate pixel coordinate system.</span></div>
<div class="line"><a id="l00063" name="l00063"></a><span class="lineno"> 63</span><span class="comment">*</span></div>
<div class="line"><a id="l00064" name="l00064"></a><span class="lineno"> 64</span><span class="comment">* Non-standard keyvalues:</span></div>
<div class="line"><a id="l00065" name="l00065"></a><span class="lineno"> 65</span><span class="comment">* -----------------------</span></div>
<div class="line"><a id="l00066" name="l00066"></a><span class="lineno"> 66</span><span class="comment">* AIPS-convention celestial projection types, NCP and GLS, and spectral</span></div>
<div class="line"><a id="l00067" name="l00067"></a><span class="lineno"> 67</span><span class="comment">* types, 'FREQ-LSR', 'FELO-HEL', etc., set in CTYPEia are translated</span></div>
<div class="line"><a id="l00068" name="l00068"></a><span class="lineno"> 68</span><span class="comment">* on-the-fly by wcsset() but without modifying the relevant ctype[], pv[] or</span></div>
<div class="line"><a id="l00069" name="l00069"></a><span class="lineno"> 69</span><span class="comment">* specsys members of the wcsprm struct. That is, only the information</span></div>
<div class="line"><a id="l00070" name="l00070"></a><span class="lineno"> 70</span><span class="comment">* extracted from ctype[] is translated when wcsset() fills in wcsprm::cel</span></div>
<div class="line"><a id="l00071" name="l00071"></a><span class="lineno"> 71</span><span class="comment">* (celprm struct) or wcsprm::spc (spcprm struct).</span></div>
<div class="line"><a id="l00072" name="l00072"></a><span class="lineno"> 72</span><span class="comment">*</span></div>
<div class="line"><a id="l00073" name="l00073"></a><span class="lineno"> 73</span><span class="comment">* On the other hand, these routines do change the values of wcsprm::ctype[],</span></div>
<div class="line"><a id="l00074" name="l00074"></a><span class="lineno"> 74</span><span class="comment">* wcsprm::pv[], wcsprm::specsys and other wcsprm struct members as</span></div>
<div class="line"><a id="l00075" name="l00075"></a><span class="lineno"> 75</span><span class="comment">* appropriate to produce the same result as if the FITS header itself had</span></div>
<div class="line"><a id="l00076" name="l00076"></a><span class="lineno"> 76</span><span class="comment">* been translated.</span></div>
<div class="line"><a id="l00077" name="l00077"></a><span class="lineno"> 77</span><span class="comment">*</span></div>
<div class="line"><a id="l00078" name="l00078"></a><span class="lineno"> 78</span><span class="comment">* Auxiliary WCS header information not used directly by WCSLIB may also be</span></div>
<div class="line"><a id="l00079" name="l00079"></a><span class="lineno"> 79</span><span class="comment">* translated. For example, the older DATE-OBS date format (wcsprm::dateobs)</span></div>
<div class="line"><a id="l00080" name="l00080"></a><span class="lineno"> 80</span><span class="comment">* is recast to year-2000 standard form, and MJD-OBS (wcsprm::mjdobs) will be</span></div>
<div class="line"><a id="l00081" name="l00081"></a><span class="lineno"> 81</span><span class="comment">* deduced from it if not already set.</span></div>
<div class="line"><a id="l00082" name="l00082"></a><span class="lineno"> 82</span><span class="comment">*</span></div>
<div class="line"><a id="l00083" name="l00083"></a><span class="lineno"> 83</span><span class="comment">* Certain combinations of keyvalues that result in malformed coordinate</span></div>
<div class="line"><a id="l00084" name="l00084"></a><span class="lineno"> 84</span><span class="comment">* systems, as described in Sect. 7.3.4 of Paper I, may also be repaired.</span></div>
<div class="line"><a id="l00085" name="l00085"></a><span class="lineno"> 85</span><span class="comment">* These are handled by cylfix().</span></div>
<div class="line"><a id="l00086" name="l00086"></a><span class="lineno"> 86</span><span class="comment">*</span></div>
<div class="line"><a id="l00087" name="l00087"></a><span class="lineno"> 87</span><span class="comment">* Non-standard keywords:</span></div>
<div class="line"><a id="l00088" name="l00088"></a><span class="lineno"> 88</span><span class="comment">* ----------------------</span></div>
<div class="line"><a id="l00089" name="l00089"></a><span class="lineno"> 89</span><span class="comment">* The AIPS-convention CROTAn keywords are recognized as quasi-standard</span></div>
<div class="line"><a id="l00090" name="l00090"></a><span class="lineno"> 90</span><span class="comment">* and as such are accomodated by wcsprm::crota[] and translated to</span></div>
<div class="line"><a id="l00091" name="l00091"></a><span class="lineno"> 91</span><span class="comment">* wcsprm::pc[][] by wcsset(). These are not dealt with here, nor are any</span></div>
<div class="line"><a id="l00092" name="l00092"></a><span class="lineno"> 92</span><span class="comment">* other non-standard keywords since these routines work only on the contents</span></div>
<div class="line"><a id="l00093" name="l00093"></a><span class="lineno"> 93</span><span class="comment">* of a wcsprm struct and do not deal with FITS headers per se. In</span></div>
<div class="line"><a id="l00094" name="l00094"></a><span class="lineno"> 94</span><span class="comment">* particular, they do not identify or translate CD00i00j, PC00i00j, PROJPn,</span></div>
<div class="line"><a id="l00095" name="l00095"></a><span class="lineno"> 95</span><span class="comment">* EPOCH, VELREF or VSOURCEa keywords; this may be done by the FITS WCS</span></div>
<div class="line"><a id="l00096" name="l00096"></a><span class="lineno"> 96</span><span class="comment">* header parser supplied with WCSLIB, refer to wcshdr.h.</span></div>
<div class="line"><a id="l00097" name="l00097"></a><span class="lineno"> 97</span><span class="comment">*</span></div>
<div class="line"><a id="l00098" name="l00098"></a><span class="lineno"> 98</span><span class="comment">* wcsfix() and wcsfixi() apply all of the corrections handled by the following</span></div>
<div class="line"><a id="l00099" name="l00099"></a><span class="lineno"> 99</span><span class="comment">* specific functions, which may also be invoked separately:</span></div>
<div class="line"><a id="l00100" name="l00100"></a><span class="lineno"> 100</span><span class="comment">*</span></div>
<div class="line"><a id="l00101" name="l00101"></a><span class="lineno"> 101</span><span class="comment">* - cdfix(): Sets the diagonal element of the CDi_ja matrix to 1.0 if all</span></div>
<div class="line"><a id="l00102" name="l00102"></a><span class="lineno"> 102</span><span class="comment">* CDi_ja keywords associated with a particular axis are omitted.</span></div>
<div class="line"><a id="l00103" name="l00103"></a><span class="lineno"> 103</span><span class="comment">*</span></div>
<div class="line"><a id="l00104" name="l00104"></a><span class="lineno"> 104</span><span class="comment">* - datfix(): recast an older DATE-OBS date format in dateobs to year-2000</span></div>
<div class="line"><a id="l00105" name="l00105"></a><span class="lineno"> 105</span><span class="comment">* standard form. Derive dateref from mjdref if not already set.</span></div>
<div class="line"><a id="l00106" name="l00106"></a><span class="lineno"> 106</span><span class="comment">* Alternatively, if dateref is set and mjdref isn't, then derive mjdref</span></div>
<div class="line"><a id="l00107" name="l00107"></a><span class="lineno"> 107</span><span class="comment">* from it. If both are set, then check consistency. Likewise for dateobs</span></div>
<div class="line"><a id="l00108" name="l00108"></a><span class="lineno"> 108</span><span class="comment">* and mjdobs; datebeg and mjdbeg; dateavg and mjdavg; and dateend and</span></div>
<div class="line"><a id="l00109" name="l00109"></a><span class="lineno"> 109</span><span class="comment">* mjdend.</span></div>
<div class="line"><a id="l00110" name="l00110"></a><span class="lineno"> 110</span><span class="comment">*</span></div>
<div class="line"><a id="l00111" name="l00111"></a><span class="lineno"> 111</span><span class="comment">* - obsfix(): if only one half of obsgeo[] is set, then derive the other</span></div>
<div class="line"><a id="l00112" name="l00112"></a><span class="lineno"> 112</span><span class="comment">* half from it. If both halves are set, then check consistency.</span></div>
<div class="line"><a id="l00113" name="l00113"></a><span class="lineno"> 113</span><span class="comment">*</span></div>
<div class="line"><a id="l00114" name="l00114"></a><span class="lineno"> 114</span><span class="comment">* - unitfix(): translate some commonly used but non-standard unit strings in</span></div>
<div class="line"><a id="l00115" name="l00115"></a><span class="lineno"> 115</span><span class="comment">* the CUNITia keyvalues, e.g. 'DEG' -> 'deg'.</span></div>
<div class="line"><a id="l00116" name="l00116"></a><span class="lineno"> 116</span><span class="comment">*</span></div>
<div class="line"><a id="l00117" name="l00117"></a><span class="lineno"> 117</span><span class="comment">* - spcfix(): translate AIPS-convention spectral types, 'FREQ-LSR',</span></div>
<div class="line"><a id="l00118" name="l00118"></a><span class="lineno"> 118</span><span class="comment">* 'FELO-HEL', etc., in ctype[] as set from CTYPEia.</span></div>
<div class="line"><a id="l00119" name="l00119"></a><span class="lineno"> 119</span><span class="comment">*</span></div>
<div class="line"><a id="l00120" name="l00120"></a><span class="lineno"> 120</span><span class="comment">* - celfix(): translate AIPS-convention celestial projection types, NCP and</span></div>
<div class="line"><a id="l00121" name="l00121"></a><span class="lineno"> 121</span><span class="comment">* GLS, in ctype[] as set from CTYPEia.</span></div>
<div class="line"><a id="l00122" name="l00122"></a><span class="lineno"> 122</span><span class="comment">*</span></div>
<div class="line"><a id="l00123" name="l00123"></a><span class="lineno"> 123</span><span class="comment">* - cylfix(): fixes WCS keyvalues for malformed cylindrical projections that</span></div>
<div class="line"><a id="l00124" name="l00124"></a><span class="lineno"> 124</span><span class="comment">* suffer from the problem described in Sect. 7.3.4 of Paper I.</span></div>
<div class="line"><a id="l00125" name="l00125"></a><span class="lineno"> 125</span><span class="comment">*</span></div>
<div class="line"><a id="l00126" name="l00126"></a><span class="lineno"> 126</span><span class="comment">*</span></div>
<div class="line"><a id="l00127" name="l00127"></a><span class="lineno"> 127</span><span class="comment">* wcsfix() - Translate a non-standard WCS struct</span></div>
<div class="line"><a id="l00128" name="l00128"></a><span class="lineno"> 128</span><span class="comment">* ----------------------------------------------</span></div>
<div class="line"><a id="l00129" name="l00129"></a><span class="lineno"> 129</span><span class="comment">* wcsfix() is identical to wcsfixi(), but lacks the info argument.</span></div>
<div class="line"><a id="l00130" name="l00130"></a><span class="lineno"> 130</span><span class="comment">*</span></div>
<div class="line"><a id="l00131" name="l00131"></a><span class="lineno"> 131</span><span class="comment">*</span></div>
<div class="line"><a id="l00132" name="l00132"></a><span class="lineno"> 132</span><span class="comment">* wcsfixi() - Translate a non-standard WCS struct</span></div>
<div class="line"><a id="l00133" name="l00133"></a><span class="lineno"> 133</span><span class="comment">* -----------------------------------------------</span></div>
<div class="line"><a id="l00134" name="l00134"></a><span class="lineno"> 134</span><span class="comment">* wcsfixi() applies all of the corrections handled separately by cdfix(),</span></div>
<div class="line"><a id="l00135" name="l00135"></a><span class="lineno"> 135</span><span class="comment">* datfix(), obsfix(), unitfix(), spcfix(), celfix(), and cylfix().</span></div>
<div class="line"><a id="l00136" name="l00136"></a><span class="lineno"> 136</span><span class="comment">*</span></div>
<div class="line"><a id="l00137" name="l00137"></a><span class="lineno"> 137</span><span class="comment">* Given:</span></div>
<div class="line"><a id="l00138" name="l00138"></a><span class="lineno"> 138</span><span class="comment">* ctrl int Do potentially unsafe translations of non-standard</span></div>
<div class="line"><a id="l00139" name="l00139"></a><span class="lineno"> 139</span><span class="comment">* unit strings as described in the usage notes to</span></div>
<div class="line"><a id="l00140" name="l00140"></a><span class="lineno"> 140</span><span class="comment">* wcsutrn().</span></div>
<div class="line"><a id="l00141" name="l00141"></a><span class="lineno"> 141</span><span class="comment">*</span></div>
<div class="line"><a id="l00142" name="l00142"></a><span class="lineno"> 142</span><span class="comment">* naxis const int []</span></div>
<div class="line"><a id="l00143" name="l00143"></a><span class="lineno"> 143</span><span class="comment">* Image axis lengths. If this array pointer is set to</span></div>
<div class="line"><a id="l00144" name="l00144"></a><span class="lineno"> 144</span><span class="comment">* zero then cylfix() will not be invoked.</span></div>
<div class="line"><a id="l00145" name="l00145"></a><span class="lineno"> 145</span><span class="comment">*</span></div>
<div class="line"><a id="l00146" name="l00146"></a><span class="lineno"> 146</span><span class="comment">* Given and returned:</span></div>
<div class="line"><a id="l00147" name="l00147"></a><span class="lineno"> 147</span><span class="comment">* wcs struct wcsprm*</span></div>
<div class="line"><a id="l00148" name="l00148"></a><span class="lineno"> 148</span><span class="comment">* Coordinate transformation parameters.</span></div>
<div class="line"><a id="l00149" name="l00149"></a><span class="lineno"> 149</span><span class="comment">*</span></div>
<div class="line"><a id="l00150" name="l00150"></a><span class="lineno"> 150</span><span class="comment">* Returned:</span></div>
<div class="line"><a id="l00151" name="l00151"></a><span class="lineno"> 151</span><span class="comment">* stat int [NWCSFIX]</span></div>
<div class="line"><a id="l00152" name="l00152"></a><span class="lineno"> 152</span><span class="comment">* Status returns from each of the functions. Use the</span></div>
<div class="line"><a id="l00153" name="l00153"></a><span class="lineno"> 153</span><span class="comment">* preprocessor macros NWCSFIX to dimension this vector</span></div>
<div class="line"><a id="l00154" name="l00154"></a><span class="lineno"> 154</span><span class="comment">* and CDFIX, DATFIX, OBSFIX, UNITFIX, SPCFIX, CELFIX,</span></div>
<div class="line"><a id="l00155" name="l00155"></a><span class="lineno"> 155</span><span class="comment">* and CYLFIX to access its elements. A status value</span></div>
<div class="line"><a id="l00156" name="l00156"></a><span class="lineno"> 156</span><span class="comment">* of -2 is set for functions that were not invoked.</span></div>
<div class="line"><a id="l00157" name="l00157"></a><span class="lineno"> 157</span><span class="comment">*</span></div>
<div class="line"><a id="l00158" name="l00158"></a><span class="lineno"> 158</span><span class="comment">* info struct wcserr [NWCSFIX]</span></div>
<div class="line"><a id="l00159" name="l00159"></a><span class="lineno"> 159</span><span class="comment">* Status messages from each of the functions. Use the</span></div>
<div class="line"><a id="l00160" name="l00160"></a><span class="lineno"> 160</span><span class="comment">* preprocessor macros NWCSFIX to dimension this vector</span></div>
<div class="line"><a id="l00161" name="l00161"></a><span class="lineno"> 161</span><span class="comment">* and CDFIX, DATFIX, OBSFIX, UNITFIX, SPCFIX, CELFIX,</span></div>
<div class="line"><a id="l00162" name="l00162"></a><span class="lineno"> 162</span><span class="comment">* and CYLFIX to access its elements.</span></div>
<div class="line"><a id="l00163" name="l00163"></a><span class="lineno"> 163</span><span class="comment">*</span></div>
<div class="line"><a id="l00164" name="l00164"></a><span class="lineno"> 164</span><span class="comment">* Note that the memory allocated by wcsfixi() for the</span></div>
<div class="line"><a id="l00165" name="l00165"></a><span class="lineno"> 165</span><span class="comment">* message in each wcserr struct (wcserr::msg, if</span></div>
<div class="line"><a id="l00166" name="l00166"></a><span class="lineno"> 166</span><span class="comment">* non-zero) must be freed by the user. See</span></div>
<div class="line"><a id="l00167" name="l00167"></a><span class="lineno"> 167</span><span class="comment">* wcsdealloc().</span></div>
<div class="line"><a id="l00168" name="l00168"></a><span class="lineno"> 168</span><span class="comment">*</span></div>
<div class="line"><a id="l00169" name="l00169"></a><span class="lineno"> 169</span><span class="comment">* Function return value:</span></div>
<div class="line"><a id="l00170" name="l00170"></a><span class="lineno"> 170</span><span class="comment">* int Status return value:</span></div>
<div class="line"><a id="l00171" name="l00171"></a><span class="lineno"> 171</span><span class="comment">* 0: Success.</span></div>
<div class="line"><a id="l00172" name="l00172"></a><span class="lineno"> 172</span><span class="comment">* 1: One or more of the translation functions</span></div>
<div class="line"><a id="l00173" name="l00173"></a><span class="lineno"> 173</span><span class="comment">* returned an error.</span></div>
<div class="line"><a id="l00174" name="l00174"></a><span class="lineno"> 174</span><span class="comment">*</span></div>
<div class="line"><a id="l00175" name="l00175"></a><span class="lineno"> 175</span><span class="comment">*</span></div>
<div class="line"><a id="l00176" name="l00176"></a><span class="lineno"> 176</span><span class="comment">* cdfix() - Fix erroneously omitted CDi_ja keywords</span></div>
<div class="line"><a id="l00177" name="l00177"></a><span class="lineno"> 177</span><span class="comment">* -------------------------------------------------</span></div>
<div class="line"><a id="l00178" name="l00178"></a><span class="lineno"> 178</span><span class="comment">* cdfix() sets the diagonal element of the CDi_ja matrix to unity if all</span></div>
<div class="line"><a id="l00179" name="l00179"></a><span class="lineno"> 179</span><span class="comment">* CDi_ja keywords associated with a given axis were omitted. According to WCS</span></div>
<div class="line"><a id="l00180" name="l00180"></a><span class="lineno"> 180</span><span class="comment">* Paper I, if any CDi_ja keywords at all are given in a FITS header then those</span></div>
<div class="line"><a id="l00181" name="l00181"></a><span class="lineno"> 181</span><span class="comment">* not given default to zero. This results in a singular matrix with an</span></div>
<div class="line"><a id="l00182" name="l00182"></a><span class="lineno"> 182</span><span class="comment">* intersecting row and column of zeros.</span></div>
<div class="line"><a id="l00183" name="l00183"></a><span class="lineno"> 183</span><span class="comment">*</span></div>
<div class="line"><a id="l00184" name="l00184"></a><span class="lineno"> 184</span><span class="comment">* cdfix() is expected to be invoked before wcsset(), which will fail if these</span></div>
<div class="line"><a id="l00185" name="l00185"></a><span class="lineno"> 185</span><span class="comment">* errors have not been corrected.</span></div>
<div class="line"><a id="l00186" name="l00186"></a><span class="lineno"> 186</span><span class="comment">*</span></div>
<div class="line"><a id="l00187" name="l00187"></a><span class="lineno"> 187</span><span class="comment">* Given and returned:</span></div>
<div class="line"><a id="l00188" name="l00188"></a><span class="lineno"> 188</span><span class="comment">* wcs struct wcsprm*</span></div>
<div class="line"><a id="l00189" name="l00189"></a><span class="lineno"> 189</span><span class="comment">* Coordinate transformation parameters.</span></div>
<div class="line"><a id="l00190" name="l00190"></a><span class="lineno"> 190</span><span class="comment">*</span></div>
<div class="line"><a id="l00191" name="l00191"></a><span class="lineno"> 191</span><span class="comment">* Function return value:</span></div>
<div class="line"><a id="l00192" name="l00192"></a><span class="lineno"> 192</span><span class="comment">* int Status return value:</span></div>
<div class="line"><a id="l00193" name="l00193"></a><span class="lineno"> 193</span><span class="comment">* -1: No change required (not an error).</span></div>
<div class="line"><a id="l00194" name="l00194"></a><span class="lineno"> 194</span><span class="comment">* 0: Success.</span></div>
<div class="line"><a id="l00195" name="l00195"></a><span class="lineno"> 195</span><span class="comment">* 1: Null wcsprm pointer passed.</span></div>
<div class="line"><a id="l00196" name="l00196"></a><span class="lineno"> 196</span><span class="comment">*</span></div>
<div class="line"><a id="l00197" name="l00197"></a><span class="lineno"> 197</span><span class="comment">*</span></div>
<div class="line"><a id="l00198" name="l00198"></a><span class="lineno"> 198</span><span class="comment">* datfix() - Translate DATE-OBS and derive MJD-OBS or vice versa</span></div>
<div class="line"><a id="l00199" name="l00199"></a><span class="lineno"> 199</span><span class="comment">* --------------------------------------------------------------</span></div>
<div class="line"><a id="l00200" name="l00200"></a><span class="lineno"> 200</span><span class="comment">* datfix() translates the old DATE-OBS date format set in wcsprm::dateobs to</span></div>
<div class="line"><a id="l00201" name="l00201"></a><span class="lineno"> 201</span><span class="comment">* year-2000 standard form (yyyy-mm-ddThh:mm:ss). It derives wcsprm::dateref</span></div>
<div class="line"><a id="l00202" name="l00202"></a><span class="lineno"> 202</span><span class="comment">* from wcsprm::mjdref if not already set. Alternatively, if dateref is set</span></div>
<div class="line"><a id="l00203" name="l00203"></a><span class="lineno"> 203</span><span class="comment">* and mjdref isn't, then it derives mjdref from it. If both are set but</span></div>
<div class="line"><a id="l00204" name="l00204"></a><span class="lineno"> 204</span><span class="comment">* disagree by more than 0.001 day (86.4 seconds) then an error status is</span></div>
<div class="line"><a id="l00205" name="l00205"></a><span class="lineno"> 205</span><span class="comment">* returned. Likewise for wcsprm::dateobs and wcsprm::mjdobs; wcsprm::datebeg</span></div>
<div class="line"><a id="l00206" name="l00206"></a><span class="lineno"> 206</span><span class="comment">* and wcsprm::mjdbeg; wcsprm::dateavg and wcsprm::mjdavg; and wcsprm::dateend</span></div>
<div class="line"><a id="l00207" name="l00207"></a><span class="lineno"> 207</span><span class="comment">* and wcsprm::mjdend.</span></div>
<div class="line"><a id="l00208" name="l00208"></a><span class="lineno"> 208</span><span class="comment">*</span></div>
<div class="line"><a id="l00209" name="l00209"></a><span class="lineno"> 209</span><span class="comment">* If neither dateobs nor mjdobs are set, but wcsprm::jepoch (primarily) or</span></div>
<div class="line"><a id="l00210" name="l00210"></a><span class="lineno"> 210</span><span class="comment">* wcsprm::bepoch is, then both are derived from it. If jepoch and/or bepoch</span></div>
<div class="line"><a id="l00211" name="l00211"></a><span class="lineno"> 211</span><span class="comment">* are set but disagree with dateobs or mjdobs by more than 0.000002 year</span></div>
<div class="line"><a id="l00212" name="l00212"></a><span class="lineno"> 212</span><span class="comment">* (63.2 seconds), an informative message is produced.</span></div>
<div class="line"><a id="l00213" name="l00213"></a><span class="lineno"> 213</span><span class="comment">*</span></div>
<div class="line"><a id="l00214" name="l00214"></a><span class="lineno"> 214</span><span class="comment">* The translations done by datfix() do not affect and are not affected by</span></div>
<div class="line"><a id="l00215" name="l00215"></a><span class="lineno"> 215</span><span class="comment">* wcsset().</span></div>
<div class="line"><a id="l00216" name="l00216"></a><span class="lineno"> 216</span><span class="comment">*</span></div>
<div class="line"><a id="l00217" name="l00217"></a><span class="lineno"> 217</span><span class="comment">* Given and returned:</span></div>
<div class="line"><a id="l00218" name="l00218"></a><span class="lineno"> 218</span><span class="comment">* wcs struct wcsprm*</span></div>
<div class="line"><a id="l00219" name="l00219"></a><span class="lineno"> 219</span><span class="comment">* Coordinate transformation parameters.</span></div>
<div class="line"><a id="l00220" name="l00220"></a><span class="lineno"> 220</span><span class="comment">* wcsprm::dateref and/or wcsprm::mjdref may be changed.</span></div>
<div class="line"><a id="l00221" name="l00221"></a><span class="lineno"> 221</span><span class="comment">* wcsprm::dateobs and/or wcsprm::mjdobs may be changed.</span></div>
<div class="line"><a id="l00222" name="l00222"></a><span class="lineno"> 222</span><span class="comment">* wcsprm::datebeg and/or wcsprm::mjdbeg may be changed.</span></div>
<div class="line"><a id="l00223" name="l00223"></a><span class="lineno"> 223</span><span class="comment">* wcsprm::dateavg and/or wcsprm::mjdavg may be changed.</span></div>
<div class="line"><a id="l00224" name="l00224"></a><span class="lineno"> 224</span><span class="comment">* wcsprm::dateend and/or wcsprm::mjdend may be changed.</span></div>
<div class="line"><a id="l00225" name="l00225"></a><span class="lineno"> 225</span><span class="comment">*</span></div>
<div class="line"><a id="l00226" name="l00226"></a><span class="lineno"> 226</span><span class="comment">* Function return value:</span></div>
<div class="line"><a id="l00227" name="l00227"></a><span class="lineno"> 227</span><span class="comment">* int Status return value:</span></div>
<div class="line"><a id="l00228" name="l00228"></a><span class="lineno"> 228</span><span class="comment">* -1: No change required (not an error).</span></div>
<div class="line"><a id="l00229" name="l00229"></a><span class="lineno"> 229</span><span class="comment">* 0: Success.</span></div>
<div class="line"><a id="l00230" name="l00230"></a><span class="lineno"> 230</span><span class="comment">* 1: Null wcsprm pointer passed.</span></div>
<div class="line"><a id="l00231" name="l00231"></a><span class="lineno"> 231</span><span class="comment">* 5: Invalid parameter value.</span></div>
<div class="line"><a id="l00232" name="l00232"></a><span class="lineno"> 232</span><span class="comment">*</span></div>
<div class="line"><a id="l00233" name="l00233"></a><span class="lineno"> 233</span><span class="comment">* For returns >= 0, a detailed message, whether</span></div>
<div class="line"><a id="l00234" name="l00234"></a><span class="lineno"> 234</span><span class="comment">* informative or an error message, may be set in</span></div>
<div class="line"><a id="l00235" name="l00235"></a><span class="lineno"> 235</span><span class="comment">* wcsprm::err if enabled, see wcserr_enable(), with</span></div>
<div class="line"><a id="l00236" name="l00236"></a><span class="lineno"> 236</span><span class="comment">* wcsprm::err.status set to FIXERR_DATE_FIX.</span></div>
<div class="line"><a id="l00237" name="l00237"></a><span class="lineno"> 237</span><span class="comment">*</span></div>
<div class="line"><a id="l00238" name="l00238"></a><span class="lineno"> 238</span><span class="comment">* Notes:</span></div>
<div class="line"><a id="l00239" name="l00239"></a><span class="lineno"> 239</span><span class="comment">* 1: The MJD algorithms used by datfix() are from D.A. Hatcher, 1984, QJRAS,</span></div>
<div class="line"><a id="l00240" name="l00240"></a><span class="lineno"> 240</span><span class="comment">* 25, 53-55, as modified by P.T. Wallace for use in SLALIB subroutines</span></div>
<div class="line"><a id="l00241" name="l00241"></a><span class="lineno"> 241</span><span class="comment">* CLDJ and DJCL.</span></div>
<div class="line"><a id="l00242" name="l00242"></a><span class="lineno"> 242</span><span class="comment">*</span></div>
<div class="line"><a id="l00243" name="l00243"></a><span class="lineno"> 243</span><span class="comment">*</span></div>
<div class="line"><a id="l00244" name="l00244"></a><span class="lineno"> 244</span><span class="comment">* obsfix() - complete the OBSGEO-[XYZLBH] vector of observatory coordinates</span></div>
<div class="line"><a id="l00245" name="l00245"></a><span class="lineno"> 245</span><span class="comment">* -------------------------------------------------------------------------</span></div>
<div class="line"><a id="l00246" name="l00246"></a><span class="lineno"> 246</span><span class="comment">* obsfix() completes the wcsprm::obsgeo vector of observatory coordinates.</span></div>
<div class="line"><a id="l00247" name="l00247"></a><span class="lineno"> 247</span><span class="comment">* That is, if only the (x,y,z) Cartesian coordinate triplet or the (l,b,h)</span></div>
<div class="line"><a id="l00248" name="l00248"></a><span class="lineno"> 248</span><span class="comment">* geodetic coordinate triplet are set, then it derives the other triplet from</span></div>
<div class="line"><a id="l00249" name="l00249"></a><span class="lineno"> 249</span><span class="comment">* it. If both triplets are set, then it checks for consistency at the level</span></div>
<div class="line"><a id="l00250" name="l00250"></a><span class="lineno"> 250</span><span class="comment">* of 1 metre.</span></div>
<div class="line"><a id="l00251" name="l00251"></a><span class="lineno"> 251</span><span class="comment">*</span></div>
<div class="line"><a id="l00252" name="l00252"></a><span class="lineno"> 252</span><span class="comment">* The operations done by obsfix() do not affect and are not affected by</span></div>
<div class="line"><a id="l00253" name="l00253"></a><span class="lineno"> 253</span><span class="comment">* wcsset().</span></div>
<div class="line"><a id="l00254" name="l00254"></a><span class="lineno"> 254</span><span class="comment">*</span></div>
<div class="line"><a id="l00255" name="l00255"></a><span class="lineno"> 255</span><span class="comment">* Given:</span></div>
<div class="line"><a id="l00256" name="l00256"></a><span class="lineno"> 256</span><span class="comment">* ctrl int Flag that controls behaviour if one triplet is</span></div>
<div class="line"><a id="l00257" name="l00257"></a><span class="lineno"> 257</span><span class="comment">* defined and the other is only partially defined:</span></div>
<div class="line"><a id="l00258" name="l00258"></a><span class="lineno"> 258</span><span class="comment">* 0: Reset only the undefined elements of an</span></div>
<div class="line"><a id="l00259" name="l00259"></a><span class="lineno"> 259</span><span class="comment">* incomplete coordinate triplet.</span></div>
<div class="line"><a id="l00260" name="l00260"></a><span class="lineno"> 260</span><span class="comment">* 1: Reset all elements of an incomplete triplet.</span></div>
<div class="line"><a id="l00261" name="l00261"></a><span class="lineno"> 261</span><span class="comment">* 2: Don't make any changes, check for consistency</span></div>
<div class="line"><a id="l00262" name="l00262"></a><span class="lineno"> 262</span><span class="comment">* only. Returns an error if either of the two</span></div>
<div class="line"><a id="l00263" name="l00263"></a><span class="lineno"> 263</span><span class="comment">* triplets is incomplete.</span></div>
<div class="line"><a id="l00264" name="l00264"></a><span class="lineno"> 264</span><span class="comment">*</span></div>
<div class="line"><a id="l00265" name="l00265"></a><span class="lineno"> 265</span><span class="comment">* Given and returned:</span></div>
<div class="line"><a id="l00266" name="l00266"></a><span class="lineno"> 266</span><span class="comment">* wcs struct wcsprm*</span></div>
<div class="line"><a id="l00267" name="l00267"></a><span class="lineno"> 267</span><span class="comment">* Coordinate transformation parameters.</span></div>
<div class="line"><a id="l00268" name="l00268"></a><span class="lineno"> 268</span><span class="comment">* wcsprm::obsgeo may be changed.</span></div>
<div class="line"><a id="l00269" name="l00269"></a><span class="lineno"> 269</span><span class="comment">*</span></div>
<div class="line"><a id="l00270" name="l00270"></a><span class="lineno"> 270</span><span class="comment">* Function return value:</span></div>
<div class="line"><a id="l00271" name="l00271"></a><span class="lineno"> 271</span><span class="comment">* int Status return value:</span></div>
<div class="line"><a id="l00272" name="l00272"></a><span class="lineno"> 272</span><span class="comment">* -1: No change required (not an error).</span></div>
<div class="line"><a id="l00273" name="l00273"></a><span class="lineno"> 273</span><span class="comment">* 0: Success.</span></div>
<div class="line"><a id="l00274" name="l00274"></a><span class="lineno"> 274</span><span class="comment">* 1: Null wcsprm pointer passed.</span></div>
<div class="line"><a id="l00275" name="l00275"></a><span class="lineno"> 275</span><span class="comment">* 5: Invalid parameter value.</span></div>
<div class="line"><a id="l00276" name="l00276"></a><span class="lineno"> 276</span><span class="comment">*</span></div>
<div class="line"><a id="l00277" name="l00277"></a><span class="lineno"> 277</span><span class="comment">* For returns >= 0, a detailed message, whether</span></div>
<div class="line"><a id="l00278" name="l00278"></a><span class="lineno"> 278</span><span class="comment">* informative or an error message, may be set in</span></div>
<div class="line"><a id="l00279" name="l00279"></a><span class="lineno"> 279</span><span class="comment">* wcsprm::err if enabled, see wcserr_enable(), with</span></div>
<div class="line"><a id="l00280" name="l00280"></a><span class="lineno"> 280</span><span class="comment">* wcsprm::err.status set to FIXERR_OBS_FIX.</span></div>
<div class="line"><a id="l00281" name="l00281"></a><span class="lineno"> 281</span><span class="comment">*</span></div>
<div class="line"><a id="l00282" name="l00282"></a><span class="lineno"> 282</span><span class="comment">* Notes:</span></div>
<div class="line"><a id="l00283" name="l00283"></a><span class="lineno"> 283</span><span class="comment">* 1: While the International Terrestrial Reference System (ITRS) is based</span></div>
<div class="line"><a id="l00284" name="l00284"></a><span class="lineno"> 284</span><span class="comment">* solely on Cartesian coordinates, it recommends the use of the GRS80</span></div>
<div class="line"><a id="l00285" name="l00285"></a><span class="lineno"> 285</span><span class="comment">* ellipsoid in converting to geodetic coordinates. However, while WCS</span></div>
<div class="line"><a id="l00286" name="l00286"></a><span class="lineno"> 286</span><span class="comment">* Paper III recommends ITRS Cartesian coordinates, Paper VII prescribes</span></div>
<div class="line"><a id="l00287" name="l00287"></a><span class="lineno"> 287</span><span class="comment">* the use of the IAU(1976) ellipsoid for geodetic coordinates, and</span></div>
<div class="line"><a id="l00288" name="l00288"></a><span class="lineno"> 288</span><span class="comment">* consequently that is what is used here.</span></div>
<div class="line"><a id="l00289" name="l00289"></a><span class="lineno"> 289</span><span class="comment">*</span></div>
<div class="line"><a id="l00290" name="l00290"></a><span class="lineno"> 290</span><span class="comment">* 2: For reference, parameters of commonly used global reference ellipsoids:</span></div>
<div class="line"><a id="l00291" name="l00291"></a><span class="lineno"> 291</span><span class="comment">*</span></div>
<div class="line"><a id="l00292" name="l00292"></a><span class="lineno"> 292</span><span class="comment">= a (m) 1/f Standard</span></div>
<div class="line"><a id="l00293" name="l00293"></a><span class="lineno"> 293</span><span class="comment">= --------- ------------- --------------------------------</span></div>
<div class="line"><a id="l00294" name="l00294"></a><span class="lineno"> 294</span><span class="comment">= 6378140 298.2577 IAU(1976)</span></div>
<div class="line"><a id="l00295" name="l00295"></a><span class="lineno"> 295</span><span class="comment">= 6378137 298.257222101 GRS80</span></div>
<div class="line"><a id="l00296" name="l00296"></a><span class="lineno"> 296</span><span class="comment">= 6378137 298.257223563 WGS84</span></div>
<div class="line"><a id="l00297" name="l00297"></a><span class="lineno"> 297</span><span class="comment">= 6378136 298.257 IERS(1989)</span></div>
<div class="line"><a id="l00298" name="l00298"></a><span class="lineno"> 298</span><span class="comment">= 6378136.6 298.25642 IERS(2003,2010), IAU(2009/2012)</span></div>
<div class="line"><a id="l00299" name="l00299"></a><span class="lineno"> 299</span><span class="comment">*</span></div>
<div class="line"><a id="l00300" name="l00300"></a><span class="lineno"> 300</span><span class="comment">* where f = (a - b) / a is the flattening, and a and b are the semi-major</span></div>
<div class="line"><a id="l00301" name="l00301"></a><span class="lineno"> 301</span><span class="comment">* and semi-minor radii in metres.</span></div>
<div class="line"><a id="l00302" name="l00302"></a><span class="lineno"> 302</span><span class="comment">*</span></div>
<div class="line"><a id="l00303" name="l00303"></a><span class="lineno"> 303</span><span class="comment">* 3: The transformation from geodetic (lng,lat,hgt) to Cartesian (x,y,z) is</span></div>
<div class="line"><a id="l00304" name="l00304"></a><span class="lineno"> 304</span><span class="comment">*</span></div>
<div class="line"><a id="l00305" name="l00305"></a><span class="lineno"> 305</span><span class="comment">= x = (n + hgt)*coslng*coslat,</span></div>
<div class="line"><a id="l00306" name="l00306"></a><span class="lineno"> 306</span><span class="comment">= y = (n + hgt)*sinlng*coslat,</span></div>
<div class="line"><a id="l00307" name="l00307"></a><span class="lineno"> 307</span><span class="comment">= z = (n*(1.0 - e^2) + hgt)*sinlat,</span></div>
<div class="line"><a id="l00308" name="l00308"></a><span class="lineno"> 308</span><span class="comment">*</span></div>
<div class="line"><a id="l00309" name="l00309"></a><span class="lineno"> 309</span><span class="comment">* where the "prime vertical radius", n, is a function of latitude</span></div>
<div class="line"><a id="l00310" name="l00310"></a><span class="lineno"> 310</span><span class="comment">*</span></div>
<div class="line"><a id="l00311" name="l00311"></a><span class="lineno"> 311</span><span class="comment">= n = a / sqrt(1 - (e*sinlat)^2),</span></div>
<div class="line"><a id="l00312" name="l00312"></a><span class="lineno"> 312</span><span class="comment">*</span></div>
<div class="line"><a id="l00313" name="l00313"></a><span class="lineno"> 313</span><span class="comment">* and a, the equatorial radius, and e^2 = (2 - f)*f, the (first)</span></div>
<div class="line"><a id="l00314" name="l00314"></a><span class="lineno"> 314</span><span class="comment">* eccentricity of the ellipsoid, are constants. obsfix() inverts these</span></div>
<div class="line"><a id="l00315" name="l00315"></a><span class="lineno"> 315</span><span class="comment">* iteratively by writing</span></div>
<div class="line"><a id="l00316" name="l00316"></a><span class="lineno"> 316</span><span class="comment">*</span></div>
<div class="line"><a id="l00317" name="l00317"></a><span class="lineno"> 317</span><span class="comment">= x = rho*coslng*coslat,</span></div>
<div class="line"><a id="l00318" name="l00318"></a><span class="lineno"> 318</span><span class="comment">= y = rho*sinlng*coslat,</span></div>
<div class="line"><a id="l00319" name="l00319"></a><span class="lineno"> 319</span><span class="comment">= zeta = rho*sinlat,</span></div>
<div class="line"><a id="l00320" name="l00320"></a><span class="lineno"> 320</span><span class="comment">*</span></div>
<div class="line"><a id="l00321" name="l00321"></a><span class="lineno"> 321</span><span class="comment">* where</span></div>
<div class="line"><a id="l00322" name="l00322"></a><span class="lineno"> 322</span><span class="comment">*</span></div>
<div class="line"><a id="l00323" name="l00323"></a><span class="lineno"> 323</span><span class="comment">= rho = n + hgt,</span></div>
<div class="line"><a id="l00324" name="l00324"></a><span class="lineno"> 324</span><span class="comment">= = sqrt(x^2 + y^2 + zeta^2),</span></div>
<div class="line"><a id="l00325" name="l00325"></a><span class="lineno"> 325</span><span class="comment">= zeta = z / (1 - n*e^2/rho),</span></div>
<div class="line"><a id="l00326" name="l00326"></a><span class="lineno"> 326</span><span class="comment">*</span></div>
<div class="line"><a id="l00327" name="l00327"></a><span class="lineno"> 327</span><span class="comment">* and iterating over the value of zeta. Since e is small, a good first</span></div>
<div class="line"><a id="l00328" name="l00328"></a><span class="lineno"> 328</span><span class="comment">* approximation is given by zeta = z.</span></div>
<div class="line"><a id="l00329" name="l00329"></a><span class="lineno"> 329</span><span class="comment">*</span></div>
<div class="line"><a id="l00330" name="l00330"></a><span class="lineno"> 330</span><span class="comment">*</span></div>
<div class="line"><a id="l00331" name="l00331"></a><span class="lineno"> 331</span><span class="comment">* unitfix() - Correct aberrant CUNITia keyvalues</span></div>
<div class="line"><a id="l00332" name="l00332"></a><span class="lineno"> 332</span><span class="comment">* ----------------------------------------------</span></div>
<div class="line"><a id="l00333" name="l00333"></a><span class="lineno"> 333</span><span class="comment">* unitfix() applies wcsutrn() to translate non-standard CUNITia keyvalues,</span></div>
<div class="line"><a id="l00334" name="l00334"></a><span class="lineno"> 334</span><span class="comment">* e.g. 'DEG' -> 'deg', also stripping off unnecessary whitespace.</span></div>
<div class="line"><a id="l00335" name="l00335"></a><span class="lineno"> 335</span><span class="comment">*</span></div>
<div class="line"><a id="l00336" name="l00336"></a><span class="lineno"> 336</span><span class="comment">* unitfix() is expected to be invoked before wcsset(), which will fail if</span></div>
<div class="line"><a id="l00337" name="l00337"></a><span class="lineno"> 337</span><span class="comment">* non-standard CUNITia keyvalues have not been translated.</span></div>
<div class="line"><a id="l00338" name="l00338"></a><span class="lineno"> 338</span><span class="comment">*</span></div>
<div class="line"><a id="l00339" name="l00339"></a><span class="lineno"> 339</span><span class="comment">* Given:</span></div>
<div class="line"><a id="l00340" name="l00340"></a><span class="lineno"> 340</span><span class="comment">* ctrl int Do potentially unsafe translations described in the</span></div>
<div class="line"><a id="l00341" name="l00341"></a><span class="lineno"> 341</span><span class="comment">* usage notes to wcsutrn().</span></div>
<div class="line"><a id="l00342" name="l00342"></a><span class="lineno"> 342</span><span class="comment">*</span></div>
<div class="line"><a id="l00343" name="l00343"></a><span class="lineno"> 343</span><span class="comment">* Given and returned:</span></div>
<div class="line"><a id="l00344" name="l00344"></a><span class="lineno"> 344</span><span class="comment">* wcs struct wcsprm*</span></div>
<div class="line"><a id="l00345" name="l00345"></a><span class="lineno"> 345</span><span class="comment">* Coordinate transformation parameters.</span></div>
<div class="line"><a id="l00346" name="l00346"></a><span class="lineno"> 346</span><span class="comment">*</span></div>
<div class="line"><a id="l00347" name="l00347"></a><span class="lineno"> 347</span><span class="comment">* Function return value:</span></div>
<div class="line"><a id="l00348" name="l00348"></a><span class="lineno"> 348</span><span class="comment">* int Status return value:</span></div>
<div class="line"><a id="l00349" name="l00349"></a><span class="lineno"> 349</span><span class="comment">* -1: No change required (not an error).</span></div>
<div class="line"><a id="l00350" name="l00350"></a><span class="lineno"> 350</span><span class="comment">* 0: Success (an alias was applied).</span></div>
<div class="line"><a id="l00351" name="l00351"></a><span class="lineno"> 351</span><span class="comment">* 1: Null wcsprm pointer passed.</span></div>
<div class="line"><a id="l00352" name="l00352"></a><span class="lineno"> 352</span><span class="comment">*</span></div>
<div class="line"><a id="l00353" name="l00353"></a><span class="lineno"> 353</span><span class="comment">* When units are translated (i.e. 0 is returned), an</span></div>
<div class="line"><a id="l00354" name="l00354"></a><span class="lineno"> 354</span><span class="comment">* informative message is set in wcsprm::err if enabled,</span></div>
<div class="line"><a id="l00355" name="l00355"></a><span class="lineno"> 355</span><span class="comment">* see wcserr_enable(), with wcsprm::err.status set to</span></div>
<div class="line"><a id="l00356" name="l00356"></a><span class="lineno"> 356</span><span class="comment">* FIXERR_UNITS_ALIAS.</span></div>
<div class="line"><a id="l00357" name="l00357"></a><span class="lineno"> 357</span><span class="comment">*</span></div>
<div class="line"><a id="l00358" name="l00358"></a><span class="lineno"> 358</span><span class="comment">*</span></div>
<div class="line"><a id="l00359" name="l00359"></a><span class="lineno"> 359</span><span class="comment">* spcfix() - Translate AIPS-convention spectral types</span></div>
<div class="line"><a id="l00360" name="l00360"></a><span class="lineno"> 360</span><span class="comment">* ---------------------------------------------------</span></div>
<div class="line"><a id="l00361" name="l00361"></a><span class="lineno"> 361</span><span class="comment">* spcfix() translates AIPS-convention spectral coordinate types,</span></div>
<div class="line"><a id="l00362" name="l00362"></a><span class="lineno"> 362</span><span class="comment">* '{FREQ,FELO,VELO}-{LSR,HEL,OBS}' (e.g. 'FREQ-OBS', 'FELO-HEL', 'VELO-LSR')</span></div>
<div class="line"><a id="l00363" name="l00363"></a><span class="lineno"> 363</span><span class="comment">* set in wcsprm::ctype[], subject to VELREF set in wcsprm::velref.</span></div>
<div class="line"><a id="l00364" name="l00364"></a><span class="lineno"> 364</span><span class="comment">*</span></div>
<div class="line"><a id="l00365" name="l00365"></a><span class="lineno"> 365</span><span class="comment">* Note that if wcs::specsys is already set then it will not be overridden.</span></div>
<div class="line"><a id="l00366" name="l00366"></a><span class="lineno"> 366</span><span class="comment">*</span></div>
<div class="line"><a id="l00367" name="l00367"></a><span class="lineno"> 367</span><span class="comment">* AIPS-convention spectral types set in CTYPEia are translated on-the-fly by</span></div>
<div class="line"><a id="l00368" name="l00368"></a><span class="lineno"> 368</span><span class="comment">* wcsset() but without modifying wcsprm::ctype[] or wcsprm::specsys. That is,</span></div>
<div class="line"><a id="l00369" name="l00369"></a><span class="lineno"> 369</span><span class="comment">* only the information extracted from wcsprm::ctype[] is translated when</span></div>
<div class="line"><a id="l00370" name="l00370"></a><span class="lineno"> 370</span><span class="comment">* wcsset() fills in wcsprm::spc (spcprm struct). spcfix() modifies</span></div>
<div class="line"><a id="l00371" name="l00371"></a><span class="lineno"> 371</span><span class="comment">* wcsprm::ctype[] so that if the header is subsequently written out, e.g. by</span></div>
<div class="line"><a id="l00372" name="l00372"></a><span class="lineno"> 372</span><span class="comment">* wcshdo(), then it will contain translated CTYPEia keyvalues.</span></div>
<div class="line"><a id="l00373" name="l00373"></a><span class="lineno"> 373</span><span class="comment">*</span></div>
<div class="line"><a id="l00374" name="l00374"></a><span class="lineno"> 374</span><span class="comment">* The operations done by spcfix() do not affect and are not affected by</span></div>
<div class="line"><a id="l00375" name="l00375"></a><span class="lineno"> 375</span><span class="comment">* wcsset().</span></div>
<div class="line"><a id="l00376" name="l00376"></a><span class="lineno"> 376</span><span class="comment">*</span></div>
<div class="line"><a id="l00377" name="l00377"></a><span class="lineno"> 377</span><span class="comment">* Given and returned:</span></div>
<div class="line"><a id="l00378" name="l00378"></a><span class="lineno"> 378</span><span class="comment">* wcs struct wcsprm*</span></div>
<div class="line"><a id="l00379" name="l00379"></a><span class="lineno"> 379</span><span class="comment">* Coordinate transformation parameters. wcsprm::ctype[]</span></div>
<div class="line"><a id="l00380" name="l00380"></a><span class="lineno"> 380</span><span class="comment">* and/or wcsprm::specsys may be changed.</span></div>
<div class="line"><a id="l00381" name="l00381"></a><span class="lineno"> 381</span><span class="comment">*</span></div>
<div class="line"><a id="l00382" name="l00382"></a><span class="lineno"> 382</span><span class="comment">* Function return value:</span></div>
<div class="line"><a id="l00383" name="l00383"></a><span class="lineno"> 383</span><span class="comment">* int Status return value:</span></div>
<div class="line"><a id="l00384" name="l00384"></a><span class="lineno"> 384</span><span class="comment">* -1: No change required (not an error).</span></div>
<div class="line"><a id="l00385" name="l00385"></a><span class="lineno"> 385</span><span class="comment">* 0: Success.</span></div>
<div class="line"><a id="l00386" name="l00386"></a><span class="lineno"> 386</span><span class="comment">* 1: Null wcsprm pointer passed.</span></div>
<div class="line"><a id="l00387" name="l00387"></a><span class="lineno"> 387</span><span class="comment">* 2: Memory allocation failed.</span></div>
<div class="line"><a id="l00388" name="l00388"></a><span class="lineno"> 388</span><span class="comment">* 3: Linear transformation matrix is singular.</span></div>
<div class="line"><a id="l00389" name="l00389"></a><span class="lineno"> 389</span><span class="comment">* 4: Inconsistent or unrecognized coordinate axis</span></div>
<div class="line"><a id="l00390" name="l00390"></a><span class="lineno"> 390</span><span class="comment">* types.</span></div>
<div class="line"><a id="l00391" name="l00391"></a><span class="lineno"> 391</span><span class="comment">* 5: Invalid parameter value.</span></div>
<div class="line"><a id="l00392" name="l00392"></a><span class="lineno"> 392</span><span class="comment">* 6: Invalid coordinate transformation parameters.</span></div>
<div class="line"><a id="l00393" name="l00393"></a><span class="lineno"> 393</span><span class="comment">* 7: Ill-conditioned coordinate transformation</span></div>
<div class="line"><a id="l00394" name="l00394"></a><span class="lineno"> 394</span><span class="comment">* parameters.</span></div>
<div class="line"><a id="l00395" name="l00395"></a><span class="lineno"> 395</span><span class="comment">*</span></div>
<div class="line"><a id="l00396" name="l00396"></a><span class="lineno"> 396</span><span class="comment">* For returns >= 0, a detailed message, whether</span></div>
<div class="line"><a id="l00397" name="l00397"></a><span class="lineno"> 397</span><span class="comment">* informative or an error message, may be set in</span></div>
<div class="line"><a id="l00398" name="l00398"></a><span class="lineno"> 398</span><span class="comment">* wcsprm::err if enabled, see wcserr_enable(), with</span></div>
<div class="line"><a id="l00399" name="l00399"></a><span class="lineno"> 399</span><span class="comment">* wcsprm::err.status set to FIXERR_SPC_UPDTE.</span></div>
<div class="line"><a id="l00400" name="l00400"></a><span class="lineno"> 400</span><span class="comment">*</span></div>
<div class="line"><a id="l00401" name="l00401"></a><span class="lineno"> 401</span><span class="comment">*</span></div>
<div class="line"><a id="l00402" name="l00402"></a><span class="lineno"> 402</span><span class="comment">* celfix() - Translate AIPS-convention celestial projection types</span></div>
<div class="line"><a id="l00403" name="l00403"></a><span class="lineno"> 403</span><span class="comment">* ---------------------------------------------------------------</span></div>
<div class="line"><a id="l00404" name="l00404"></a><span class="lineno"> 404</span><span class="comment">* celfix() translates AIPS-convention celestial projection types, NCP and</span></div>
<div class="line"><a id="l00405" name="l00405"></a><span class="lineno"> 405</span><span class="comment">* GLS, set in the ctype[] member of the wcsprm struct.</span></div>
<div class="line"><a id="l00406" name="l00406"></a><span class="lineno"> 406</span><span class="comment">*</span></div>
<div class="line"><a id="l00407" name="l00407"></a><span class="lineno"> 407</span><span class="comment">* Two additional pv[] keyvalues are created when translating NCP, and three</span></div>
<div class="line"><a id="l00408" name="l00408"></a><span class="lineno"> 408</span><span class="comment">* are created when translating GLS with non-zero reference point. If the pv[]</span></div>
<div class="line"><a id="l00409" name="l00409"></a><span class="lineno"> 409</span><span class="comment">* array was initially allocated by wcsini() then the array will be expanded if</span></div>
<div class="line"><a id="l00410" name="l00410"></a><span class="lineno"> 410</span><span class="comment">* necessary. Otherwise, error 2 will be returned if sufficient empty slots</span></div>
<div class="line"><a id="l00411" name="l00411"></a><span class="lineno"> 411</span><span class="comment">* are not already available for use.</span></div>
<div class="line"><a id="l00412" name="l00412"></a><span class="lineno"> 412</span><span class="comment">*</span></div>
<div class="line"><a id="l00413" name="l00413"></a><span class="lineno"> 413</span><span class="comment">* AIPS-convention celestial projection types set in CTYPEia are translated</span></div>
<div class="line"><a id="l00414" name="l00414"></a><span class="lineno"> 414</span><span class="comment">* on-the-fly by wcsset() but without modifying wcsprm::ctype[], wcsprm::pv[],</span></div>
<div class="line"><a id="l00415" name="l00415"></a><span class="lineno"> 415</span><span class="comment">* or wcsprm::npv. That is, only the information extracted from</span></div>
<div class="line"><a id="l00416" name="l00416"></a><span class="lineno"> 416</span><span class="comment">* wcsprm::ctype[] is translated when wcsset() fills in wcsprm::cel (celprm</span></div>
<div class="line"><a id="l00417" name="l00417"></a><span class="lineno"> 417</span><span class="comment">* struct). celfix() modifies wcsprm::ctype[], wcsprm::pv[], and wcsprm::npv</span></div>
<div class="line"><a id="l00418" name="l00418"></a><span class="lineno"> 418</span><span class="comment">* so that if the header is subsequently written out, e.g. by wcshdo(), then it</span></div>
<div class="line"><a id="l00419" name="l00419"></a><span class="lineno"> 419</span><span class="comment">* will contain translated CTYPEia keyvalues and the relevant PVi_ma.</span></div>
<div class="line"><a id="l00420" name="l00420"></a><span class="lineno"> 420</span><span class="comment">*</span></div>
<div class="line"><a id="l00421" name="l00421"></a><span class="lineno"> 421</span><span class="comment">* The operations done by celfix() do not affect and are not affected by</span></div>
<div class="line"><a id="l00422" name="l00422"></a><span class="lineno"> 422</span><span class="comment">* wcsset(). However, it uses information in the wcsprm struct provided by</span></div>
<div class="line"><a id="l00423" name="l00423"></a><span class="lineno"> 423</span><span class="comment">* wcsset(), and will invoke it if necessary.</span></div>
<div class="line"><a id="l00424" name="l00424"></a><span class="lineno"> 424</span><span class="comment">*</span></div>
<div class="line"><a id="l00425" name="l00425"></a><span class="lineno"> 425</span><span class="comment">* Given and returned:</span></div>
<div class="line"><a id="l00426" name="l00426"></a><span class="lineno"> 426</span><span class="comment">* wcs struct wcsprm*</span></div>
<div class="line"><a id="l00427" name="l00427"></a><span class="lineno"> 427</span><span class="comment">* Coordinate transformation parameters. wcsprm::ctype[]</span></div>
<div class="line"><a id="l00428" name="l00428"></a><span class="lineno"> 428</span><span class="comment">* and/or wcsprm::pv[] may be changed.</span></div>
<div class="line"><a id="l00429" name="l00429"></a><span class="lineno"> 429</span><span class="comment">*</span></div>
<div class="line"><a id="l00430" name="l00430"></a><span class="lineno"> 430</span><span class="comment">* Function return value:</span></div>
<div class="line"><a id="l00431" name="l00431"></a><span class="lineno"> 431</span><span class="comment">* int Status return value:</span></div>
<div class="line"><a id="l00432" name="l00432"></a><span class="lineno"> 432</span><span class="comment">* -1: No change required (not an error).</span></div>
<div class="line"><a id="l00433" name="l00433"></a><span class="lineno"> 433</span><span class="comment">* 0: Success.</span></div>
<div class="line"><a id="l00434" name="l00434"></a><span class="lineno"> 434</span><span class="comment">* 1: Null wcsprm pointer passed.</span></div>
<div class="line"><a id="l00435" name="l00435"></a><span class="lineno"> 435</span><span class="comment">* 2: Memory allocation failed.</span></div>
<div class="line"><a id="l00436" name="l00436"></a><span class="lineno"> 436</span><span class="comment">* 3: Linear transformation matrix is singular.</span></div>
<div class="line"><a id="l00437" name="l00437"></a><span class="lineno"> 437</span><span class="comment">* 4: Inconsistent or unrecognized coordinate axis</span></div>
<div class="line"><a id="l00438" name="l00438"></a><span class="lineno"> 438</span><span class="comment">* types.</span></div>
<div class="line"><a id="l00439" name="l00439"></a><span class="lineno"> 439</span><span class="comment">* 5: Invalid parameter value.</span></div>
<div class="line"><a id="l00440" name="l00440"></a><span class="lineno"> 440</span><span class="comment">* 6: Invalid coordinate transformation parameters.</span></div>
<div class="line"><a id="l00441" name="l00441"></a><span class="lineno"> 441</span><span class="comment">* 7: Ill-conditioned coordinate transformation</span></div>
<div class="line"><a id="l00442" name="l00442"></a><span class="lineno"> 442</span><span class="comment">* parameters.</span></div>
<div class="line"><a id="l00443" name="l00443"></a><span class="lineno"> 443</span><span class="comment">*</span></div>
<div class="line"><a id="l00444" name="l00444"></a><span class="lineno"> 444</span><span class="comment">* For returns > 1, a detailed error message is set in</span></div>
<div class="line"><a id="l00445" name="l00445"></a><span class="lineno"> 445</span><span class="comment">* wcsprm::err if enabled, see wcserr_enable().</span></div>
<div class="line"><a id="l00446" name="l00446"></a><span class="lineno"> 446</span><span class="comment">*</span></div>
<div class="line"><a id="l00447" name="l00447"></a><span class="lineno"> 447</span><span class="comment">*</span></div>
<div class="line"><a id="l00448" name="l00448"></a><span class="lineno"> 448</span><span class="comment">* cylfix() - Fix malformed cylindrical projections</span></div>
<div class="line"><a id="l00449" name="l00449"></a><span class="lineno"> 449</span><span class="comment">* ------------------------------------------------</span></div>
<div class="line"><a id="l00450" name="l00450"></a><span class="lineno"> 450</span><span class="comment">* cylfix() fixes WCS keyvalues for malformed cylindrical projections that</span></div>
<div class="line"><a id="l00451" name="l00451"></a><span class="lineno"> 451</span><span class="comment">* suffer from the problem described in Sect. 7.3.4 of Paper I.</span></div>
<div class="line"><a id="l00452" name="l00452"></a><span class="lineno"> 452</span><span class="comment">*</span></div>
<div class="line"><a id="l00453" name="l00453"></a><span class="lineno"> 453</span><span class="comment">* cylfix() requires the wcsprm struct to have been set up by wcsset(), and</span></div>
<div class="line"><a id="l00454" name="l00454"></a><span class="lineno"> 454</span><span class="comment">* will invoke it if necessary. After modification, the struct is reset on</span></div>
<div class="line"><a id="l00455" name="l00455"></a><span class="lineno"> 455</span><span class="comment">* return with an explicit call to wcsset().</span></div>
<div class="line"><a id="l00456" name="l00456"></a><span class="lineno"> 456</span><span class="comment">*</span></div>
<div class="line"><a id="l00457" name="l00457"></a><span class="lineno"> 457</span><span class="comment">* Given:</span></div>
<div class="line"><a id="l00458" name="l00458"></a><span class="lineno"> 458</span><span class="comment">* naxis const int []</span></div>
<div class="line"><a id="l00459" name="l00459"></a><span class="lineno"> 459</span><span class="comment">* Image axis lengths.</span></div>
<div class="line"><a id="l00460" name="l00460"></a><span class="lineno"> 460</span><span class="comment">*</span></div>
<div class="line"><a id="l00461" name="l00461"></a><span class="lineno"> 461</span><span class="comment">* Given and returned:</span></div>
<div class="line"><a id="l00462" name="l00462"></a><span class="lineno"> 462</span><span class="comment">* wcs struct wcsprm*</span></div>
<div class="line"><a id="l00463" name="l00463"></a><span class="lineno"> 463</span><span class="comment">* Coordinate transformation parameters.</span></div>
<div class="line"><a id="l00464" name="l00464"></a><span class="lineno"> 464</span><span class="comment">*</span></div>
<div class="line"><a id="l00465" name="l00465"></a><span class="lineno"> 465</span><span class="comment">* Function return value:</span></div>
<div class="line"><a id="l00466" name="l00466"></a><span class="lineno"> 466</span><span class="comment">* int Status return value:</span></div>
<div class="line"><a id="l00467" name="l00467"></a><span class="lineno"> 467</span><span class="comment">* -1: No change required (not an error).</span></div>
<div class="line"><a id="l00468" name="l00468"></a><span class="lineno"> 468</span><span class="comment">* 0: Success.</span></div>
<div class="line"><a id="l00469" name="l00469"></a><span class="lineno"> 469</span><span class="comment">* 1: Null wcsprm pointer passed.</span></div>
<div class="line"><a id="l00470" name="l00470"></a><span class="lineno"> 470</span><span class="comment">* 2: Memory allocation failed.</span></div>
<div class="line"><a id="l00471" name="l00471"></a><span class="lineno"> 471</span><span class="comment">* 3: Linear transformation matrix is singular.</span></div>
<div class="line"><a id="l00472" name="l00472"></a><span class="lineno"> 472</span><span class="comment">* 4: Inconsistent or unrecognized coordinate axis</span></div>
<div class="line"><a id="l00473" name="l00473"></a><span class="lineno"> 473</span><span class="comment">* types.</span></div>
<div class="line"><a id="l00474" name="l00474"></a><span class="lineno"> 474</span><span class="comment">* 5: Invalid parameter value.</span></div>
<div class="line"><a id="l00475" name="l00475"></a><span class="lineno"> 475</span><span class="comment">* 6: Invalid coordinate transformation parameters.</span></div>
<div class="line"><a id="l00476" name="l00476"></a><span class="lineno"> 476</span><span class="comment">* 7: Ill-conditioned coordinate transformation</span></div>
<div class="line"><a id="l00477" name="l00477"></a><span class="lineno"> 477</span><span class="comment">* parameters.</span></div>
<div class="line"><a id="l00478" name="l00478"></a><span class="lineno"> 478</span><span class="comment">* 8: All of the corner pixel coordinates are invalid.</span></div>
<div class="line"><a id="l00479" name="l00479"></a><span class="lineno"> 479</span><span class="comment">* 9: Could not determine reference pixel coordinate.</span></div>
<div class="line"><a id="l00480" name="l00480"></a><span class="lineno"> 480</span><span class="comment">* 10: Could not determine reference pixel value.</span></div>
<div class="line"><a id="l00481" name="l00481"></a><span class="lineno"> 481</span><span class="comment">*</span></div>
<div class="line"><a id="l00482" name="l00482"></a><span class="lineno"> 482</span><span class="comment">* For returns > 1, a detailed error message is set in</span></div>
<div class="line"><a id="l00483" name="l00483"></a><span class="lineno"> 483</span><span class="comment">* wcsprm::err if enabled, see wcserr_enable().</span></div>
<div class="line"><a id="l00484" name="l00484"></a><span class="lineno"> 484</span><span class="comment">*</span></div>
<div class="line"><a id="l00485" name="l00485"></a><span class="lineno"> 485</span><span class="comment">*</span></div>
<div class="line"><a id="l00486" name="l00486"></a><span class="lineno"> 486</span><span class="comment">* wcspcx() - regularize PCi_j</span></div>
<div class="line"><a id="l00487" name="l00487"></a><span class="lineno"> 487</span><span class="comment">* ---------------------------</span></div>
<div class="line"><a id="l00488" name="l00488"></a><span class="lineno"> 488</span><span class="comment">* wcspcx() "regularizes" the linear transformation matrix component of the</span></div>
<div class="line"><a id="l00489" name="l00489"></a><span class="lineno"> 489</span><span class="comment">* coordinate transformation (PCi_ja) to make it more human-readable.</span></div>
<div class="line"><a id="l00490" name="l00490"></a><span class="lineno"> 490</span><span class="comment">*</span></div>
<div class="line"><a id="l00491" name="l00491"></a><span class="lineno"> 491</span><span class="comment">* Normally, upon encountering a FITS header containing a CDi_ja matrix,</span></div>
<div class="line"><a id="l00492" name="l00492"></a><span class="lineno"> 492</span><span class="comment">* wcsset() simply treats it as PCi_ja and sets CDELTia to unity. However,</span></div>
<div class="line"><a id="l00493" name="l00493"></a><span class="lineno"> 493</span><span class="comment">* wcspcx() decomposes CDi_ja into PCi_ja and CDELTia in such a way that</span></div>
<div class="line"><a id="l00494" name="l00494"></a><span class="lineno"> 494</span><span class="comment">* CDELTia form meaningful scaling parameters. In practice, the residual</span></div>
<div class="line"><a id="l00495" name="l00495"></a><span class="lineno"> 495</span><span class="comment">* PCi_ja matrix will often then be orthogonal, i.e. unity, or describing a</span></div>
<div class="line"><a id="l00496" name="l00496"></a><span class="lineno"> 496</span><span class="comment">* pure rotation, axis permutation, or reflection, or a combination thereof.</span></div>
<div class="line"><a id="l00497" name="l00497"></a><span class="lineno"> 497</span><span class="comment">*</span></div>
<div class="line"><a id="l00498" name="l00498"></a><span class="lineno"> 498</span><span class="comment">* The decomposition is based on normalizing the length in the transformed</span></div>
<div class="line"><a id="l00499" name="l00499"></a><span class="lineno"> 499</span><span class="comment">* system (i.e. intermediate pixel coordinates) of the orthonormal basis</span></div>
<div class="line"><a id="l00500" name="l00500"></a><span class="lineno"> 500</span><span class="comment">* vectors of the pixel coordinate system. This deviates slightly from the</span></div>
<div class="line"><a id="l00501" name="l00501"></a><span class="lineno"> 501</span><span class="comment">* prescription given by Eq. (4) of WCS Paper I, namely Sum(j=1,N)(PCi_ja)² = 1,</span></div>
<div class="line"><a id="l00502" name="l00502"></a><span class="lineno"> 502</span><span class="comment">* in replacing the sum over j with the sum over i. Consequently, the columns</span></div>
<div class="line"><a id="l00503" name="l00503"></a><span class="lineno"> 503</span><span class="comment">* of PCi_ja will consist of unit vectors. In practice, especially in cubes</span></div>
<div class="line"><a id="l00504" name="l00504"></a><span class="lineno"> 504</span><span class="comment">* and higher dimensional images, at least some pairs of these unit vectors, if</span></div>
<div class="line"><a id="l00505" name="l00505"></a><span class="lineno"> 505</span><span class="comment">* not all, will often be orthogonal or close to orthogonal.</span></div>
<div class="line"><a id="l00506" name="l00506"></a><span class="lineno"> 506</span><span class="comment">*</span></div>
<div class="line"><a id="l00507" name="l00507"></a><span class="lineno"> 507</span><span class="comment">* The sign of CDELTia is chosen to make the PCi_ja matrix as close to the,</span></div>
<div class="line"><a id="l00508" name="l00508"></a><span class="lineno"> 508</span><span class="comment">* possibly permuted, unit matrix as possible, except that where the coordinate</span></div>
<div class="line"><a id="l00509" name="l00509"></a><span class="lineno"> 509</span><span class="comment">* description contains a pair of celestial axes, the sign of CDELTia is set</span></div>
<div class="line"><a id="l00510" name="l00510"></a><span class="lineno"> 510</span><span class="comment">* negative for the longitude axis and positive for the latitude axis.</span></div>
<div class="line"><a id="l00511" name="l00511"></a><span class="lineno"> 511</span><span class="comment">*</span></div>
<div class="line"><a id="l00512" name="l00512"></a><span class="lineno"> 512</span><span class="comment">* Optionally, rows of the PCi_ja matrix may also be permuted to diagonalize</span></div>
<div class="line"><a id="l00513" name="l00513"></a><span class="lineno"> 513</span><span class="comment">* it as far as possible, thus undoing any transposition of axes in the</span></div>
<div class="line"><a id="l00514" name="l00514"></a><span class="lineno"> 514</span><span class="comment">* intermediate pixel coordinate system.</span></div>
<div class="line"><a id="l00515" name="l00515"></a><span class="lineno"> 515</span><span class="comment">*</span></div>
<div class="line"><a id="l00516" name="l00516"></a><span class="lineno"> 516</span><span class="comment">* If the coordinate description contains a celestial plane, then the angle of</span></div>
<div class="line"><a id="l00517" name="l00517"></a><span class="lineno"> 517</span><span class="comment">* rotation of each of the basis vectors associated with the celestial axes is</span></div>
<div class="line"><a id="l00518" name="l00518"></a><span class="lineno"> 518</span><span class="comment">* returned. For a pure rotation the two angles should be identical. Any</span></div>
<div class="line"><a id="l00519" name="l00519"></a><span class="lineno"> 519</span><span class="comment">* difference between them is a measure of axis skewness.</span></div>
<div class="line"><a id="l00520" name="l00520"></a><span class="lineno"> 520</span><span class="comment">*</span></div>
<div class="line"><a id="l00521" name="l00521"></a><span class="lineno"> 521</span><span class="comment">* The decomposition is not performed for axes involving a sequent distortion</span></div>
<div class="line"><a id="l00522" name="l00522"></a><span class="lineno"> 522</span><span class="comment">* function that is defined in terms of CDi_ja, such as TPV, TNX, or ZPX, which</span></div>
<div class="line"><a id="l00523" name="l00523"></a><span class="lineno"> 523</span><span class="comment">* always are. The independent variables of the polynomial are therefore</span></div>
<div class="line"><a id="l00524" name="l00524"></a><span class="lineno"> 524</span><span class="comment">* intermediate world coordinates rather than intermediate pixel coordinates.</span></div>
<div class="line"><a id="l00525" name="l00525"></a><span class="lineno"> 525</span><span class="comment">* Because sequent distortions are always applied before CDELTia, if CDi_ja was</span></div>
<div class="line"><a id="l00526" name="l00526"></a><span class="lineno"> 526</span><span class="comment">* translated to PCi_ja plus CDELTia, then the distortion would be altered</span></div>
<div class="line"><a id="l00527" name="l00527"></a><span class="lineno"> 527</span><span class="comment">* unless the polynomial coefficients were also adjusted to account for the</span></div>
<div class="line"><a id="l00528" name="l00528"></a><span class="lineno"> 528</span><span class="comment">* change of scale.</span></div>
<div class="line"><a id="l00529" name="l00529"></a><span class="lineno"> 529</span><span class="comment">*</span></div>
<div class="line"><a id="l00530" name="l00530"></a><span class="lineno"> 530</span><span class="comment">* wcspcx() requires the wcsprm struct to have been set up by wcsset(), and</span></div>
<div class="line"><a id="l00531" name="l00531"></a><span class="lineno"> 531</span><span class="comment">* will invoke it if necessary. The wcsprm struct is reset on return with an</span></div>
<div class="line"><a id="l00532" name="l00532"></a><span class="lineno"> 532</span><span class="comment">* explicit call to wcsset().</span></div>
<div class="line"><a id="l00533" name="l00533"></a><span class="lineno"> 533</span><span class="comment">*</span></div>
<div class="line"><a id="l00534" name="l00534"></a><span class="lineno"> 534</span><span class="comment">* Given and returned:</span></div>
<div class="line"><a id="l00535" name="l00535"></a><span class="lineno"> 535</span><span class="comment">* wcs struct wcsprm*</span></div>
<div class="line"><a id="l00536" name="l00536"></a><span class="lineno"> 536</span><span class="comment">* Coordinate transformation parameters.</span></div>
<div class="line"><a id="l00537" name="l00537"></a><span class="lineno"> 537</span><span class="comment">*</span></div>
<div class="line"><a id="l00538" name="l00538"></a><span class="lineno"> 538</span><span class="comment">* Given:</span></div>
<div class="line"><a id="l00539" name="l00539"></a><span class="lineno"> 539</span><span class="comment">* dopc int If 1, then PCi_ja and CDELTia, as given, will be</span></div>
<div class="line"><a id="l00540" name="l00540"></a><span class="lineno"> 540</span><span class="comment">* recomposed according to the above prescription. If 0,</span></div>
<div class="line"><a id="l00541" name="l00541"></a><span class="lineno"> 541</span><span class="comment">* the operation is restricted to decomposing CDi_ja.</span></div>
<div class="line"><a id="l00542" name="l00542"></a><span class="lineno"> 542</span><span class="comment">*</span></div>
<div class="line"><a id="l00543" name="l00543"></a><span class="lineno"> 543</span><span class="comment">* permute int If 1, then after decomposition (or recomposition),</span></div>
<div class="line"><a id="l00544" name="l00544"></a><span class="lineno"> 544</span><span class="comment">* permute rows of PCi_ja to make the axes of the</span></div>
<div class="line"><a id="l00545" name="l00545"></a><span class="lineno"> 545</span><span class="comment">* intermediate pixel coordinate system match as closely</span></div>
<div class="line"><a id="l00546" name="l00546"></a><span class="lineno"> 546</span><span class="comment">* as possible those of the pixel coordinates. That is,</span></div>
<div class="line"><a id="l00547" name="l00547"></a><span class="lineno"> 547</span><span class="comment">* make it as close to a diagonal matrix as possible.</span></div>
<div class="line"><a id="l00548" name="l00548"></a><span class="lineno"> 548</span><span class="comment">* However, celestial axes are special in always being</span></div>
<div class="line"><a id="l00549" name="l00549"></a><span class="lineno"> 549</span><span class="comment">* paired, with the longitude axis preceding the latitude</span></div>
<div class="line"><a id="l00550" name="l00550"></a><span class="lineno"> 550</span><span class="comment">* axis.</span></div>
<div class="line"><a id="l00551" name="l00551"></a><span class="lineno"> 551</span><span class="comment">*</span></div>
<div class="line"><a id="l00552" name="l00552"></a><span class="lineno"> 552</span><span class="comment">* All WCS entities indexed by i, such as CTYPEia,</span></div>
<div class="line"><a id="l00553" name="l00553"></a><span class="lineno"> 553</span><span class="comment">* CRVALia, CDELTia, etc., including coordinate lookup</span></div>
<div class="line"><a id="l00554" name="l00554"></a><span class="lineno"> 554</span><span class="comment">* tables, will also be permuted as necessary to account</span></div>
<div class="line"><a id="l00555" name="l00555"></a><span class="lineno"> 555</span><span class="comment">* for the change to PCi_ja. This does not apply to</span></div>
<div class="line"><a id="l00556" name="l00556"></a><span class="lineno"> 556</span><span class="comment">* CRPIXja, nor prior distortion functions. These</span></div>
<div class="line"><a id="l00557" name="l00557"></a><span class="lineno"> 557</span><span class="comment">* operate on pixel coordinates, which are not affected</span></div>
<div class="line"><a id="l00558" name="l00558"></a><span class="lineno"> 558</span><span class="comment">* by the permutation.</span></div>
<div class="line"><a id="l00559" name="l00559"></a><span class="lineno"> 559</span><span class="comment">*</span></div>
<div class="line"><a id="l00560" name="l00560"></a><span class="lineno"> 560</span><span class="comment">* Returned:</span></div>
<div class="line"><a id="l00561" name="l00561"></a><span class="lineno"> 561</span><span class="comment">* rotn double[2] Rotation angle [deg] of each basis vector associated</span></div>
<div class="line"><a id="l00562" name="l00562"></a><span class="lineno"> 562</span><span class="comment">* with the celestial axes. For a pure rotation the two</span></div>
<div class="line"><a id="l00563" name="l00563"></a><span class="lineno"> 563</span><span class="comment">* angles should be identical. Any difference between</span></div>
<div class="line"><a id="l00564" name="l00564"></a><span class="lineno"> 564</span><span class="comment">* them is a measure of axis skewness.</span></div>
<div class="line"><a id="l00565" name="l00565"></a><span class="lineno"> 565</span><span class="comment">*</span></div>
<div class="line"><a id="l00566" name="l00566"></a><span class="lineno"> 566</span><span class="comment">* May be set to the NULL pointer if this information is</span></div>
<div class="line"><a id="l00567" name="l00567"></a><span class="lineno"> 567</span><span class="comment">* not required.</span></div>
<div class="line"><a id="l00568" name="l00568"></a><span class="lineno"> 568</span><span class="comment">*</span></div>
<div class="line"><a id="l00569" name="l00569"></a><span class="lineno"> 569</span><span class="comment">* Function return value:</span></div>
<div class="line"><a id="l00570" name="l00570"></a><span class="lineno"> 570</span><span class="comment">* int Status return value:</span></div>
<div class="line"><a id="l00571" name="l00571"></a><span class="lineno"> 571</span><span class="comment">* 0: Success.</span></div>
<div class="line"><a id="l00572" name="l00572"></a><span class="lineno"> 572</span><span class="comment">* 1: Null wcsprm pointer passed.</span></div>
<div class="line"><a id="l00573" name="l00573"></a><span class="lineno"> 573</span><span class="comment">* 2: Memory allocation failed.</span></div>
<div class="line"><a id="l00574" name="l00574"></a><span class="lineno"> 574</span><span class="comment">* 5: CDi_j matrix not used.</span></div>
<div class="line"><a id="l00575" name="l00575"></a><span class="lineno"> 575</span><span class="comment">* 6: Sequent distortion function present.</span></div>
<div class="line"><a id="l00576" name="l00576"></a><span class="lineno"> 576</span><span class="comment">*</span></div>
<div class="line"><a id="l00577" name="l00577"></a><span class="lineno"> 577</span><span class="comment">*</span></div>
<div class="line"><a id="l00578" name="l00578"></a><span class="lineno"> 578</span><span class="comment">* Global variable: const char *wcsfix_errmsg[] - Status return messages</span></div>
<div class="line"><a id="l00579" name="l00579"></a><span class="lineno"> 579</span><span class="comment">* ---------------------------------------------------------------------</span></div>
<div class="line"><a id="l00580" name="l00580"></a><span class="lineno"> 580</span><span class="comment">* Error messages to match the status value returned from each function.</span></div>
<div class="line"><a id="l00581" name="l00581"></a><span class="lineno"> 581</span><span class="comment">*</span></div>
<div class="line"><a id="l00582" name="l00582"></a><span class="lineno"> 582</span><span class="comment">*===========================================================================*/</span></div>
<div class="line"><a id="l00583" name="l00583"></a><span class="lineno"> 583</span> </div>
<div class="line"><a id="l00584" name="l00584"></a><span class="lineno"> 584</span><span class="preprocessor">#ifndef WCSLIB_WCSFIX</span></div>
<div class="line"><a id="l00585" name="l00585"></a><span class="lineno"> 585</span><span class="preprocessor">#define WCSLIB_WCSFIX</span></div>
<div class="line"><a id="l00586" name="l00586"></a><span class="lineno"> 586</span> </div>
<div class="line"><a id="l00587" name="l00587"></a><span class="lineno"> 587</span><span class="preprocessor">#include "<a class="code" href="wcs_8h.html">wcs.h</a>"</span></div>
<div class="line"><a id="l00588" name="l00588"></a><span class="lineno"> 588</span><span class="preprocessor">#include "<a class="code" href="wcserr_8h.html">wcserr.h</a>"</span></div>
<div class="line"><a id="l00589" name="l00589"></a><span class="lineno"> 589</span> </div>
<div class="line"><a id="l00590" name="l00590"></a><span class="lineno"> 590</span><span class="preprocessor">#ifdef __cplusplus</span></div>
<div class="line"><a id="l00591" name="l00591"></a><span class="lineno"> 591</span><span class="keyword">extern</span> <span class="stringliteral">"C"</span> {</div>
<div class="line"><a id="l00592" name="l00592"></a><span class="lineno"> 592</span><span class="preprocessor">#endif</span></div>
<div class="line"><a id="l00593" name="l00593"></a><span class="lineno"> 593</span> </div>
<div class="line"><a id="l00594" name="l00594"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#af23e7b02522c40fa5dfbf3d569348844"> 594</a></span><span class="preprocessor">#define CDFIX 0</span></div>
<div class="line"><a id="l00595" name="l00595"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a7181ebe5e9f0a4058642c56dc848bd5c"> 595</a></span><span class="preprocessor">#define DATFIX 1</span></div>
<div class="line"><a id="l00596" name="l00596"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#ac293192284a2143faaf5b3b1d0be9262"> 596</a></span><span class="preprocessor">#define OBSFIX 2</span></div>
<div class="line"><a id="l00597" name="l00597"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a8f4a947e2605b35ffa92f08b113d60b2"> 597</a></span><span class="preprocessor">#define UNITFIX 3</span></div>
<div class="line"><a id="l00598" name="l00598"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0816c5f2354ee6c0044e11867d7558ea"> 598</a></span><span class="preprocessor">#define SPCFIX 4</span></div>
<div class="line"><a id="l00599" name="l00599"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#af1b99efe520fbd2d4bd0e5a35f87e186"> 599</a></span><span class="preprocessor">#define CELFIX 5</span></div>
<div class="line"><a id="l00600" name="l00600"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a4d37e0274dff84649cba075b8761b3fa"> 600</a></span><span class="preprocessor">#define CYLFIX 6</span></div>
<div class="line"><a id="l00601" name="l00601"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0ed13e54c3eacb9325afbae78ef33b61"> 601</a></span><span class="preprocessor">#define NWCSFIX 7</span></div>
<div class="line"><a id="l00602" name="l00602"></a><span class="lineno"> 602</span> </div>
<div class="line"><a id="l00603" name="l00603"></a><span class="lineno"> 603</span><span class="keyword">extern</span> <span class="keyword">const</span> <span class="keywordtype">char</span> *<a class="code hl_variable" href="wcsfix_8h.html#a256ce6281894f65dd15396cc0994e875">wcsfix_errmsg</a>[];</div>
<div class="line"><a id="l00604" name="l00604"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a3229b126ed844da0a2d4f7abff1de7d0"> 604</a></span><span class="preprocessor">#define cylfix_errmsg wcsfix_errmsg</span></div>
<div class="line"><a id="l00605" name="l00605"></a><span class="lineno"> 605</span> </div>
<div class="foldopen" id="foldopen00606" data-start="{" data-end="};">
<div class="line"><a id="l00606" name="l00606"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183"> 606</a></span><span class="keyword">enum</span> <a class="code hl_enumeration" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183">wcsfix_errmsg_enum</a> {</div>
<div class="line"><a id="l00607" name="l00607"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1e29d170017673feccf26343a7fd9528"> 607</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1e29d170017673feccf26343a7fd9528">FIXERR_OBSGEO_FIX</a> = -5, <span class="comment">// Observatory coordinates amended.</span></div>
<div class="line"><a id="l00608" name="l00608"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1edfc4a839295befd132ce460d253311"> 608</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1edfc4a839295befd132ce460d253311">FIXERR_DATE_FIX</a> = -4, <span class="comment">// Date string reformatted.</span></div>
<div class="line"><a id="l00609" name="l00609"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a7906f35912dd0fb39954bfd5140672a7"> 609</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a7906f35912dd0fb39954bfd5140672a7">FIXERR_SPC_UPDATE</a> = -3, <span class="comment">// Spectral axis type modified.</span></div>
<div class="line"><a id="l00610" name="l00610"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1e255b3aa269ff0c251ada8a6c3f9602"> 610</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1e255b3aa269ff0c251ada8a6c3f9602">FIXERR_UNITS_ALIAS</a> = -2, <span class="comment">// Units alias translation.</span></div>
<div class="line"><a id="l00611" name="l00611"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183aec3fdc50ed9f4ca8d80d7ce7751ef0e3"> 611</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183aec3fdc50ed9f4ca8d80d7ce7751ef0e3">FIXERR_NO_CHANGE</a> = -1, <span class="comment">// No change.</span></div>
<div class="line"><a id="l00612" name="l00612"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183aee9fbc64e56bb6d307d06d8ef8e8b244"> 612</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183aee9fbc64e56bb6d307d06d8ef8e8b244">FIXERR_SUCCESS</a> = 0, <span class="comment">// Success.</span></div>
<div class="line"><a id="l00613" name="l00613"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183af574a836e251e8a0257da97580bb9354"> 613</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183af574a836e251e8a0257da97580bb9354">FIXERR_NULL_POINTER</a> = 1, <span class="comment">// Null wcsprm pointer passed.</span></div>
<div class="line"><a id="l00614" name="l00614"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1e4cf4eeb3cd2f4d8c2c1f040aa62f6c"> 614</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1e4cf4eeb3cd2f4d8c2c1f040aa62f6c">FIXERR_MEMORY</a> = 2, <span class="comment">// Memory allocation failed.</span></div>
<div class="line"><a id="l00615" name="l00615"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a8553bf40509263e3c3a198810f83d26e"> 615</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a8553bf40509263e3c3a198810f83d26e">FIXERR_SINGULAR_MTX</a> = 3, <span class="comment">// Linear transformation matrix is singular.</span></div>
<div class="line"><a id="l00616" name="l00616"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a421fc9b9a2aac54bc832b3c1180f8f07"> 616</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a421fc9b9a2aac54bc832b3c1180f8f07">FIXERR_BAD_CTYPE</a> = 4, <span class="comment">// Inconsistent or unrecognized coordinate</span></div>
<div class="line"><a id="l00617" name="l00617"></a><span class="lineno"> 617</span> <span class="comment">// axis types.</span></div>
<div class="line"><a id="l00618" name="l00618"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a26d787caed068586fbef3d3c0fbce41f"> 618</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a26d787caed068586fbef3d3c0fbce41f">FIXERR_BAD_PARAM</a> = 5, <span class="comment">// Invalid parameter value.</span></div>
<div class="line"><a id="l00619" name="l00619"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a5dd410d6f1a55543c4f7d0f82435eb40"> 619</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a5dd410d6f1a55543c4f7d0f82435eb40">FIXERR_BAD_COORD_TRANS</a> = 6, <span class="comment">// Invalid coordinate transformation</span></div>
<div class="line"><a id="l00620" name="l00620"></a><span class="lineno"> 620</span> <span class="comment">// parameters.</span></div>
<div class="line"><a id="l00621" name="l00621"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a81b5390b4f770515ae950d9e382b2885"> 621</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a81b5390b4f770515ae950d9e382b2885">FIXERR_ILL_COORD_TRANS</a> = 7, <span class="comment">// Ill-conditioned coordinate transformation</span></div>
<div class="line"><a id="l00622" name="l00622"></a><span class="lineno"> 622</span> <span class="comment">// parameters.</span></div>
<div class="line"><a id="l00623" name="l00623"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a3f4b7a9a303943f6c12ea51cce2240cf"> 623</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a3f4b7a9a303943f6c12ea51cce2240cf">FIXERR_BAD_CORNER_PIX</a> = 8, <span class="comment">// All of the corner pixel coordinates are</span></div>
<div class="line"><a id="l00624" name="l00624"></a><span class="lineno"> 624</span> <span class="comment">// invalid.</span></div>
<div class="line"><a id="l00625" name="l00625"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183ad6bf7801d043f41f67c54677d6cfcb75"> 625</a></span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183ad6bf7801d043f41f67c54677d6cfcb75">FIXERR_NO_REF_PIX_COORD</a> = 9, <span class="comment">// Could not determine reference pixel</span></div>
<div class="line"><a id="l00626" name="l00626"></a><span class="lineno"> 626</span> <span class="comment">// coordinate.</span></div>
<div class="line"><a id="l00627" name="l00627"></a><span class="lineno"> 627</span> <a class="code hl_enumvalue" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a15a9e5f9cbb559ef53018e9aade43e88">FIXERR_NO_REF_PIX_VAL</a> = 10 <span class="comment">// Could not determine reference pixel value.</span></div>
<div class="line"><a id="l00628" name="l00628"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a15a9e5f9cbb559ef53018e9aade43e88"> 628</a></span>};</div>
</div>
<div class="line"><a id="l00629" name="l00629"></a><span class="lineno"> 629</span> </div>
<div class="line"><a id="l00630" name="l00630"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a89e1b5b4d2fa89af03f5d1143352b05f"> 630</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="wcsfix_8h.html#a89e1b5b4d2fa89af03f5d1143352b05f">wcsfix</a>(<span class="keywordtype">int</span> ctrl, <span class="keyword">const</span> <span class="keywordtype">int</span> naxis[], <span class="keyword">struct</span> <a class="code hl_struct" href="structwcsprm.html">wcsprm</a> *wcs, <span class="keywordtype">int</span> stat[]);</div>
<div class="line"><a id="l00631" name="l00631"></a><span class="lineno"> 631</span> </div>
<div class="line"><a id="l00632" name="l00632"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a62298e0fb06332a282d9daab718a1286"> 632</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="wcsfix_8h.html#a62298e0fb06332a282d9daab718a1286">wcsfixi</a>(<span class="keywordtype">int</span> ctrl, <span class="keyword">const</span> <span class="keywordtype">int</span> naxis[], <span class="keyword">struct</span> <a class="code hl_struct" href="structwcsprm.html">wcsprm</a> *wcs, <span class="keywordtype">int</span> stat[],</div>
<div class="line"><a id="l00633" name="l00633"></a><span class="lineno"> 633</span> <span class="keyword">struct</span> <a class="code hl_struct" href="structwcserr.html">wcserr</a> info[]);</div>
<div class="line"><a id="l00634" name="l00634"></a><span class="lineno"> 634</span> </div>
<div class="line"><a id="l00635" name="l00635"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a25714f1558ecbee6c1b1fef0abf8ea7f"> 635</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="wcsfix_8h.html#a25714f1558ecbee6c1b1fef0abf8ea7f">cdfix</a>(<span class="keyword">struct</span> <a class="code hl_struct" href="structwcsprm.html">wcsprm</a> *wcs);</div>
<div class="line"><a id="l00636" name="l00636"></a><span class="lineno"> 636</span> </div>
<div class="line"><a id="l00637" name="l00637"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a77b614a15de67b42040c2be46cbfca1a"> 637</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="wcsfix_8h.html#a77b614a15de67b42040c2be46cbfca1a">datfix</a>(<span class="keyword">struct</span> <a class="code hl_struct" href="structwcsprm.html">wcsprm</a> *wcs);</div>
<div class="line"><a id="l00638" name="l00638"></a><span class="lineno"> 638</span> </div>
<div class="line"><a id="l00639" name="l00639"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a52691a17e224d8f4f1d1f897798efa49"> 639</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="wcsfix_8h.html#a52691a17e224d8f4f1d1f897798efa49">obsfix</a>(<span class="keywordtype">int</span> ctrl, <span class="keyword">struct</span> <a class="code hl_struct" href="structwcsprm.html">wcsprm</a> *wcs);</div>
<div class="line"><a id="l00640" name="l00640"></a><span class="lineno"> 640</span> </div>
<div class="line"><a id="l00641" name="l00641"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a883167275c4d3855ba453364db3d8d66"> 641</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="wcsfix_8h.html#a883167275c4d3855ba453364db3d8d66">unitfix</a>(<span class="keywordtype">int</span> ctrl, <span class="keyword">struct</span> <a class="code hl_struct" href="structwcsprm.html">wcsprm</a> *wcs);</div>
<div class="line"><a id="l00642" name="l00642"></a><span class="lineno"> 642</span> </div>
<div class="line"><a id="l00643" name="l00643"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#af011e4065b6179e19d2964bc9646b6af"> 643</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="wcsfix_8h.html#af011e4065b6179e19d2964bc9646b6af">spcfix</a>(<span class="keyword">struct</span> <a class="code hl_struct" href="structwcsprm.html">wcsprm</a> *wcs);</div>
<div class="line"><a id="l00644" name="l00644"></a><span class="lineno"> 644</span> </div>
<div class="line"><a id="l00645" name="l00645"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#ac1df72303f64e50d5e3cb320c126443b"> 645</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="wcsfix_8h.html#ac1df72303f64e50d5e3cb320c126443b">celfix</a>(<span class="keyword">struct</span> <a class="code hl_struct" href="structwcsprm.html">wcsprm</a> *wcs);</div>
<div class="line"><a id="l00646" name="l00646"></a><span class="lineno"> 646</span> </div>
<div class="line"><a id="l00647" name="l00647"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a07281faacbec1df800a417bf157751d7"> 647</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="wcsfix_8h.html#a07281faacbec1df800a417bf157751d7">cylfix</a>(<span class="keyword">const</span> <span class="keywordtype">int</span> naxis[], <span class="keyword">struct</span> <a class="code hl_struct" href="structwcsprm.html">wcsprm</a> *wcs);</div>
<div class="line"><a id="l00648" name="l00648"></a><span class="lineno"> 648</span> </div>
<div class="line"><a id="l00649" name="l00649"></a><span class="lineno"><a class="line" href="wcsfix_8h.html#a7b70e8127424672db38ce3e653e46467"> 649</a></span><span class="keywordtype">int</span> <a class="code hl_function" href="wcsfix_8h.html#a7b70e8127424672db38ce3e653e46467">wcspcx</a>(<span class="keyword">struct</span> <a class="code hl_struct" href="structwcsprm.html">wcsprm</a> *wcs, <span class="keywordtype">int</span> dopc, <span class="keywordtype">int</span> permute, <span class="keywordtype">double</span> rotn[2]);</div>
<div class="line"><a id="l00650" name="l00650"></a><span class="lineno"> 650</span> </div>
<div class="line"><a id="l00651" name="l00651"></a><span class="lineno"> 651</span> </div>
<div class="line"><a id="l00652" name="l00652"></a><span class="lineno"> 652</span><span class="preprocessor">#ifdef __cplusplus</span></div>
<div class="line"><a id="l00653" name="l00653"></a><span class="lineno"> 653</span>}</div>
<div class="line"><a id="l00654" name="l00654"></a><span class="lineno"> 654</span><span class="preprocessor">#endif</span></div>
<div class="line"><a id="l00655" name="l00655"></a><span class="lineno"> 655</span> </div>
<div class="line"><a id="l00656" name="l00656"></a><span class="lineno"> 656</span><span class="preprocessor">#endif </span><span class="comment">// WCSLIB_WCSFIX</span></div>
<div class="ttc" id="astructwcserr_html"><div class="ttname"><a href="structwcserr.html">wcserr</a></div><div class="ttdoc">Error message handling.</div><div class="ttdef"><b>Definition</b> wcserr.h:243</div></div>
<div class="ttc" id="astructwcsprm_html"><div class="ttname"><a href="structwcsprm.html">wcsprm</a></div><div class="ttdoc">Coordinate transformation parameters.</div><div class="ttdef"><b>Definition</b> wcs.h:2085</div></div>
<div class="ttc" id="awcs_8h_html"><div class="ttname"><a href="wcs_8h.html">wcs.h</a></div></div>
<div class="ttc" id="awcserr_8h_html"><div class="ttname"><a href="wcserr_8h.html">wcserr.h</a></div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183">wcsfix_errmsg_enum</a></div><div class="ttdeci">wcsfix_errmsg_enum</div><div class="ttdef"><b>Definition</b> wcsfix.h:606</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a15a9e5f9cbb559ef53018e9aade43e88"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a15a9e5f9cbb559ef53018e9aade43e88">FIXERR_NO_REF_PIX_VAL</a></div><div class="ttdeci">@ FIXERR_NO_REF_PIX_VAL</div><div class="ttdef"><b>Definition</b> wcsfix.h:627</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a1e255b3aa269ff0c251ada8a6c3f9602"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1e255b3aa269ff0c251ada8a6c3f9602">FIXERR_UNITS_ALIAS</a></div><div class="ttdeci">@ FIXERR_UNITS_ALIAS</div><div class="ttdef"><b>Definition</b> wcsfix.h:610</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a1e29d170017673feccf26343a7fd9528"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1e29d170017673feccf26343a7fd9528">FIXERR_OBSGEO_FIX</a></div><div class="ttdeci">@ FIXERR_OBSGEO_FIX</div><div class="ttdef"><b>Definition</b> wcsfix.h:607</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a1e4cf4eeb3cd2f4d8c2c1f040aa62f6c"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1e4cf4eeb3cd2f4d8c2c1f040aa62f6c">FIXERR_MEMORY</a></div><div class="ttdeci">@ FIXERR_MEMORY</div><div class="ttdef"><b>Definition</b> wcsfix.h:614</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a1edfc4a839295befd132ce460d253311"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a1edfc4a839295befd132ce460d253311">FIXERR_DATE_FIX</a></div><div class="ttdeci">@ FIXERR_DATE_FIX</div><div class="ttdef"><b>Definition</b> wcsfix.h:608</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a26d787caed068586fbef3d3c0fbce41f"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a26d787caed068586fbef3d3c0fbce41f">FIXERR_BAD_PARAM</a></div><div class="ttdeci">@ FIXERR_BAD_PARAM</div><div class="ttdef"><b>Definition</b> wcsfix.h:618</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a3f4b7a9a303943f6c12ea51cce2240cf"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a3f4b7a9a303943f6c12ea51cce2240cf">FIXERR_BAD_CORNER_PIX</a></div><div class="ttdeci">@ FIXERR_BAD_CORNER_PIX</div><div class="ttdef"><b>Definition</b> wcsfix.h:623</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a421fc9b9a2aac54bc832b3c1180f8f07"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a421fc9b9a2aac54bc832b3c1180f8f07">FIXERR_BAD_CTYPE</a></div><div class="ttdeci">@ FIXERR_BAD_CTYPE</div><div class="ttdef"><b>Definition</b> wcsfix.h:616</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a5dd410d6f1a55543c4f7d0f82435eb40"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a5dd410d6f1a55543c4f7d0f82435eb40">FIXERR_BAD_COORD_TRANS</a></div><div class="ttdeci">@ FIXERR_BAD_COORD_TRANS</div><div class="ttdef"><b>Definition</b> wcsfix.h:619</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a7906f35912dd0fb39954bfd5140672a7"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a7906f35912dd0fb39954bfd5140672a7">FIXERR_SPC_UPDATE</a></div><div class="ttdeci">@ FIXERR_SPC_UPDATE</div><div class="ttdef"><b>Definition</b> wcsfix.h:609</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a81b5390b4f770515ae950d9e382b2885"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a81b5390b4f770515ae950d9e382b2885">FIXERR_ILL_COORD_TRANS</a></div><div class="ttdeci">@ FIXERR_ILL_COORD_TRANS</div><div class="ttdef"><b>Definition</b> wcsfix.h:621</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183a8553bf40509263e3c3a198810f83d26e"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183a8553bf40509263e3c3a198810f83d26e">FIXERR_SINGULAR_MTX</a></div><div class="ttdeci">@ FIXERR_SINGULAR_MTX</div><div class="ttdef"><b>Definition</b> wcsfix.h:615</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183ad6bf7801d043f41f67c54677d6cfcb75"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183ad6bf7801d043f41f67c54677d6cfcb75">FIXERR_NO_REF_PIX_COORD</a></div><div class="ttdeci">@ FIXERR_NO_REF_PIX_COORD</div><div class="ttdef"><b>Definition</b> wcsfix.h:625</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183aec3fdc50ed9f4ca8d80d7ce7751ef0e3"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183aec3fdc50ed9f4ca8d80d7ce7751ef0e3">FIXERR_NO_CHANGE</a></div><div class="ttdeci">@ FIXERR_NO_CHANGE</div><div class="ttdef"><b>Definition</b> wcsfix.h:611</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183aee9fbc64e56bb6d307d06d8ef8e8b244"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183aee9fbc64e56bb6d307d06d8ef8e8b244">FIXERR_SUCCESS</a></div><div class="ttdeci">@ FIXERR_SUCCESS</div><div class="ttdef"><b>Definition</b> wcsfix.h:612</div></div>
<div class="ttc" id="awcsfix_8h_html_a0399bbea1e28abad3259a8ea05b25183af574a836e251e8a0257da97580bb9354"><div class="ttname"><a href="wcsfix_8h.html#a0399bbea1e28abad3259a8ea05b25183af574a836e251e8a0257da97580bb9354">FIXERR_NULL_POINTER</a></div><div class="ttdeci">@ FIXERR_NULL_POINTER</div><div class="ttdef"><b>Definition</b> wcsfix.h:613</div></div>
<div class="ttc" id="awcsfix_8h_html_a07281faacbec1df800a417bf157751d7"><div class="ttname"><a href="wcsfix_8h.html#a07281faacbec1df800a417bf157751d7">cylfix</a></div><div class="ttdeci">int cylfix(const int naxis[], struct wcsprm *wcs)</div><div class="ttdoc">Fix malformed cylindrical projections.</div></div>
<div class="ttc" id="awcsfix_8h_html_a256ce6281894f65dd15396cc0994e875"><div class="ttname"><a href="wcsfix_8h.html#a256ce6281894f65dd15396cc0994e875">wcsfix_errmsg</a></div><div class="ttdeci">const char * wcsfix_errmsg[]</div><div class="ttdoc">Status return messages.</div></div>
<div class="ttc" id="awcsfix_8h_html_a25714f1558ecbee6c1b1fef0abf8ea7f"><div class="ttname"><a href="wcsfix_8h.html#a25714f1558ecbee6c1b1fef0abf8ea7f">cdfix</a></div><div class="ttdeci">int cdfix(struct wcsprm *wcs)</div><div class="ttdoc">Fix erroneously omitted CDi_ja keywords.</div></div>
<div class="ttc" id="awcsfix_8h_html_a52691a17e224d8f4f1d1f897798efa49"><div class="ttname"><a href="wcsfix_8h.html#a52691a17e224d8f4f1d1f897798efa49">obsfix</a></div><div class="ttdeci">int obsfix(int ctrl, struct wcsprm *wcs)</div><div class="ttdoc">complete the OBSGEO-[XYZLBH] vector of observatory coordinates.</div></div>
<div class="ttc" id="awcsfix_8h_html_a62298e0fb06332a282d9daab718a1286"><div class="ttname"><a href="wcsfix_8h.html#a62298e0fb06332a282d9daab718a1286">wcsfixi</a></div><div class="ttdeci">int wcsfixi(int ctrl, const int naxis[], struct wcsprm *wcs, int stat[], struct wcserr info[])</div><div class="ttdoc">Translate a non-standard WCS struct.</div></div>
<div class="ttc" id="awcsfix_8h_html_a77b614a15de67b42040c2be46cbfca1a"><div class="ttname"><a href="wcsfix_8h.html#a77b614a15de67b42040c2be46cbfca1a">datfix</a></div><div class="ttdeci">int datfix(struct wcsprm *wcs)</div><div class="ttdoc">Translate DATE-OBS and derive MJD-OBS or vice versa.</div></div>
<div class="ttc" id="awcsfix_8h_html_a7b70e8127424672db38ce3e653e46467"><div class="ttname"><a href="wcsfix_8h.html#a7b70e8127424672db38ce3e653e46467">wcspcx</a></div><div class="ttdeci">int wcspcx(struct wcsprm *wcs, int dopc, int permute, double rotn[2])</div><div class="ttdoc">regularize PCi_j.</div></div>
<div class="ttc" id="awcsfix_8h_html_a883167275c4d3855ba453364db3d8d66"><div class="ttname"><a href="wcsfix_8h.html#a883167275c4d3855ba453364db3d8d66">unitfix</a></div><div class="ttdeci">int unitfix(int ctrl, struct wcsprm *wcs)</div><div class="ttdoc">Correct aberrant CUNITia keyvalues.</div></div>
<div class="ttc" id="awcsfix_8h_html_a89e1b5b4d2fa89af03f5d1143352b05f"><div class="ttname"><a href="wcsfix_8h.html#a89e1b5b4d2fa89af03f5d1143352b05f">wcsfix</a></div><div class="ttdeci">int wcsfix(int ctrl, const int naxis[], struct wcsprm *wcs, int stat[])</div><div class="ttdoc">Translate a non-standard WCS struct.</div></div>
<div class="ttc" id="awcsfix_8h_html_ac1df72303f64e50d5e3cb320c126443b"><div class="ttname"><a href="wcsfix_8h.html#ac1df72303f64e50d5e3cb320c126443b">celfix</a></div><div class="ttdeci">int celfix(struct wcsprm *wcs)</div><div class="ttdoc">Translate AIPS-convention celestial projection types.</div></div>
<div class="ttc" id="awcsfix_8h_html_af011e4065b6179e19d2964bc9646b6af"><div class="ttname"><a href="wcsfix_8h.html#af011e4065b6179e19d2964bc9646b6af">spcfix</a></div><div class="ttdeci">int spcfix(struct wcsprm *wcs)</div><div class="ttdoc">Translate AIPS-convention spectral types.</div></div>
</div><!-- fragment --></div><!-- contents -->
<!-- start footer part -->
<hr class="footer"/><address class="footer"><small>
Generated on Wed Nov 29 2023 19:09:57 for WCSLIB by <a href="https://www.doxygen.org/index.html"><img class="footer" src="doxygen.svg" width="104" height="31" alt="doxygen"/></a> 1.9.8
</small></address>
</body>
</html>