1 /*
2 *class++
3 *  Name:
4 *     SpecFrame
5 
6 *  Purpose:
7 *     Spectral coordinate system description.
8 
9 *  Constructor Function:
10 c     astSpecFrame
11 f     AST_SPECFRAME
12 
13 *  Description:
14 *     A SpecFrame is a specialised form of one-dimensional Frame which
15 *     represents various coordinate systems used to describe positions within
16 *     an electro-magnetic spectrum. The particular coordinate system to be
17 *     used is specified by setting the SpecFrame's System attribute (the
18 *     default is wavelength) qualified, as necessary, by other attributes
19 *     such as the rest frequency, the standard of rest, the epoch of
20 *     observation, units, etc (see the description of the System attribute
21 *     for details).
22 *
23 *     By setting a value for thr SpecOrigin attribute, a SpecFrame can be made
24 *     to represent offsets from a given spectral position, rather than absolute
25 *     spectral values.
26 
27 *  Inheritance:
28 *     The SpecFrame class inherits from the Frame class.
29 
30 *  Attributes:
31 *     In addition to those attributes common to all Frames, every
32 *     SpecFrame also has the following attributes:
33 *
34 *     - AlignSpecOffset: Align SpecFrames using the offset coordinate system?
35 *     - AlignStdOfRest: Standard of rest in which to align SpecFrames
36 *     - RefDec: Declination of the source (FK5 J2000)
37 *     - RefRA: Right ascension of the source (FK5 J2000)
38 *     - RestFreq: Rest frequency
39 *     - SourceSys: Source velocity spectral system
40 *     - SourceVel: Source velocity
41 *     - SourceVRF: Source velocity rest frame
42 *     - SpecOrigin: The zero point for SpecFrame axis values
43 *     - StdOfRest: Standard of rest
44 *
45 *     Several of the Frame attributes inherited by the SpecFrame class
46 *     refer to a specific axis of the Frame (for instance Unit(axis),
47 *     Label(axis), etc). Since a SpecFrame is strictly one-dimensional,
48 *     it allows these attributes to be specified without an axis index.
49 *     So for instance, "Unit" is allowed in place of "Unit(1)".
50 
51 *  Functions:
52 c     In addition to those functions applicable to all Frames, the
53 c     following functions may also be applied to all SpecFrames:
54 f     In addition to those routines applicable to all Frames, the
55 f     following routines may also be applied to all SpecFrames:
56 *
57 c     - astSetRefPos: Set reference position in any celestial system
58 f     - AST_SETREFPOS: Set reference position in any celestial system
59 c     - astGetRefPos: Get reference position in any celestial system
60 f     - AST_GETREFPOS: Get reference position in any celestial system
61 
62 *  Copyright:
63 *     Copyright (C) 1997-2006 Council for the Central Laboratory of the
64 *     Research Councils
65 
66 *  Licence:
67 *     This program is free software: you can redistribute it and/or
68 *     modify it under the terms of the GNU Lesser General Public
69 *     License as published by the Free Software Foundation, either
70 *     version 3 of the License, or (at your option) any later
71 *     version.
72 *
73 *     This program is distributed in the hope that it will be useful,
74 *     but WITHOUT ANY WARRANTY; without even the implied warranty of
75 *     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
76 *     GNU Lesser General Public License for more details.
77 *
78 *     You should have received a copy of the GNU Lesser General
79 *     License along with this program.  If not, see
80 *     <http://www.gnu.org/licenses/>.
81 
82 *  Authors:
83 *     DSB: David S. Berry (Starlink)
84 
85 *  History:
86 *     4-NOV-2002 (DSB):
87 *        Original version.
88 *     2-FEB-2005 (DSB):
89 *        - Avoid using astStore to allocate more storage than is supplied
90 *        in the "data" pointer. This can cause access violations since
91 *        astStore will then read beyond the end of the "data" area.
92 *     22-MAR-2005 (DSB):
93 *        - Re-structure MakeSpecMapping in order to avoid unnecessary
94 *        access to SpecFrame attributes which may not be set, and to
95 *        check that all required attributes have been set if UseDefs is
96 *        zero.
97 *     23-MAR-2005 (DSB):
98 *        - Added missing rest frames to SorEqual.
99 *     12-AUG-2005 (DSB):
100 *        - Remove GeoLon and GeoLat attributes. Use the new ObsLon and
101 *        ObsLat attributes in the parent Frame class instead. Note, for
102 *        backward compatibility the public attribute accessors and the
103 *        astLoadSpecFrame functions still recogonise GeoLon and GeoLat,
104 *        but use the ObsLat/ObsLon attributes internally.
105 *     14-FEB-2006 (DSB):
106 *        Override astGetObjSize.
107 *     1-MAR-2006 (DSB):
108 *        Replace astSetPermMap within DEBUG blocks by astBeginPM/astEndPM.
109 *     6-OCT-2006 (DSB):
110 *        Guard against annulling null pointers in subFrame.
111 *     18-OCT-2006 (DSB):
112 *        Added SpecOrigin and AlignSpecOffset attributes.
113 *     23-OCT-2006 (DSB):
114 *        Fix memory leak caused by addition of SpecOrigin and AlignSpecOffset
115 *        attributes.
116 *     15-NOV-2006 (DSB):
117 *        Only write out SpecOrigin if it is not bad.
118 *     8-JAN-2006 (DSB):
119 *        - SubFrame: Copy the SourceSystem and SourceStdOfRest attributes
120 *        to the System and StdOfRest attributes of the "align_frm"
121 *        SpecFrame before calling MakeSpecMapping. Previously, the
122 *        values assigned to  SourceSystem and SourceStdOfRest were
123 *        ignored, and alignment was always performed in the templates System
124 *        and StdOfRest.
125 *        - MakeSpecMapping: Correct logic used to decide if steps 2 and 7
126 *        can be cancelled.
127 *        - OriginSystem: Clear the AlignSpecOffset attributes before
128 *        finding the Mapping between the old and new Systems.
129 *     16-JAN-2006 (DSB):
130 *        Fix bug in Dump that caused SrcVRF not to be written out.
131 *     31-JAN-2007 (DSB):
132 *        Modified so that a SpecFrame can be used as a template to find a
133 *        SpecFrame contained within a CmpFrame. This involves changes in
134 *        Match and the removal of the local versions of SetMaxAxes and
135 *        SetMinAxes.
136 *     8-AUG-2007 (DSB):
137 *        Changed Overlay to avoid the possibility of making permanent
138 *        changes to the supplied template Frame.
139 *     3-SEP-2007 (DSB):
140 *        In SubFrame, since AlignSystem is extended by the SpecFrame class
141 *        it needs to be cleared before invoking the parent SubFrame
142 *        method in cases where the result Frame is not a SkyFrame.
143 *     2-OCT-2007 (DSB):
144 *        In Overlay, clear AlignSystem as well as System before calling
145 *        the parent overlay method.
146 *     4-SEP-2009 (DSB):
147 *        In MakeSpecMapping, in order to produce alignment that is not
148 *        affected by the epoch or reference position, make the alignment
149 *        frame adapt to the epoch and reference position of the target
150 *        and result Frames.
151 *     14-SEP-2009 (DSB):
152 *        In MakeSpecMapping, extend the 4-SEP-2009 fix to cover other
153 *        attributes that define the available rest frames (e.g.
154 *        SourceVRF, SourceVel, ObsLat, ObsLon, ObsAlt).
155 *     16-SEP-2009 (DSB):
156 *        In MakeSpecMapping, retain the original alignment frame attribute
157 *        values if we are restoring the integrity of a FrameSet.
158 *     29-APR-2011 (DSB):
159 *        Prevent astFindFrame from matching a subclass template against a
160 *        superclass target.
161 *class--
162 */
163 
164 /* Module Macros. */
165 /* ============== */
166 /* Set the name of the class we are implementing. This indicates to
167    the header files that define class interfaces that they should make
168    "protected" symbols available. */
169 #define astCLASS SpecFrame
170 
171 /* Define the first and last acceptable System values. */
172 #define FIRST_SYSTEM AST__FREQ
173 #define LAST_SYSTEM AST__VREL
174 
175 /* Define the first and last acceptable StdOfRest values. */
176 #define FIRST_SOR AST__TPSOR
177 #define LAST_SOR AST__SCSOR
178 
179 /* The supported spectral coordinate systems fall into two groups;
180    "relative", and "absolute". The relative systems define each axis
181    value with respect to the rest frequency, whereas the absolute systems
182    have axis values which do not depend on the rest frequency. Define a
183    macro which returns one if the specified system is absolute, and zero
184    otherwise. */
185 #define ABS_SYSTEM(sys) \
186       ( ( sys == AST__ENERGY || \
187           sys == AST__WAVENUM || \
188           sys == AST__WAVELEN || \
189           sys == AST__AIRWAVE || \
190           sys == AST__FREQ ) ? 1 : 0 )
191 
192 /* Macros which return the maximum and minimum of two values. */
193 #define MAX(aa,bb) ((aa)>(bb)?(aa):(bb))
194 #define MIN(aa,bb) ((aa)<(bb)?(aa):(bb))
195 
196 /* Macro to check for equality of floating point values. We cannot
197    compare bad values directory because of the danger of floating point
198    exceptions, so bad values are dealt with explicitly. */
199 #define EQUAL(aa,bb) (((aa)==AST__BAD)?(((bb)==AST__BAD)?1:0):(((bb)==AST__BAD)?0:(fabs((aa)-(bb))<=1.0E5*MAX((fabs(aa)+fabs(bb))*DBL_EPSILON,DBL_MIN))))
200 
201 /* Define other numerical constants for use in this module. */
202 #define GETATTRIB_BUFF_LEN 50
203 #define GETLABEL_BUFF_LEN 200
204 #define GETSYMBOL_BUFF_LEN 20
205 #define GETTITLE_BUFF_LEN 200
206 
207 /* Header files. */
208 /* ============= */
209 /* Interface definitions. */
210 /* ---------------------- */
211 
212 #include "globals.h"             /* Thread-safe global data access */
213 #include "error.h"               /* Error reporting facilities */
214 #include "memory.h"              /* Memory allocation facilities */
215 #include "unit.h"                /* Units management facilities */
216 #include "globals.h"             /* Thread-safe global data access */
217 #include "object.h"              /* Base Object class */
218 #include "specmap.h"             /* Spectral coordinate Mappings */
219 #include "frame.h"               /* Parent Frame class */
220 #include "skyframe.h"            /* Celestial coordinate frames */
221 #include "specframe.h"           /* Interface definition for this class */
222 #include "mapping.h"             /* Coordinate Mappings */
223 #include "cmpmap.h"              /* Compound Mappings */
224 #include "unitmap.h"             /* Unit Mappings */
225 #include "pal.h"                 /* SlaLib interface */
226 #include "shiftmap.h"            /* Change of origin */
227 
228 /* Error code definitions. */
229 /* ----------------------- */
230 #include "ast_err.h"             /* AST error codes */
231 
232 /* C header files. */
233 /* --------------- */
234 #include <stdarg.h>
235 #include <stdio.h>
236 #include <string.h>
237 #include <ctype.h>
238 #include <stddef.h>
239 #include <math.h>
240 #include <limits.h>
241 
242 /* Module Variables. */
243 /* ================= */
244 
245 /* Address of this static variable is used as a unique identifier for
246    member of this class. */
247 static int class_check;
248 
249 /* Pointers to parent class methods which are used or extended by this
250    class. */
251 static int (* parent_getobjsize)( AstObject *, int * );
252 static AstSystemType (* parent_getalignsystem)( AstFrame *, int * );
253 static AstSystemType (* parent_getsystem)( AstFrame *, int * );
254 static const char *(* parent_getattrib)( AstObject *, const char *, int * );
255 static const char *(* parent_getdomain)( AstFrame *, int * );
256 static const char *(* parent_getlabel)( AstFrame *, int, int * );
257 static const char *(* parent_getsymbol)( AstFrame *, int, int * );
258 static const char *(* parent_gettitle)( AstFrame *, int * );
259 static const char *(* parent_getunit)( AstFrame *, int, int * );
260 static int (* parent_match)( AstFrame *, AstFrame *, int, int **, int **, AstMapping **, AstFrame **, int * );
261 static int (* parent_subframe)( AstFrame *, AstFrame *, int, const int *, const int *, AstMapping **, AstFrame **, int * );
262 static int (* parent_testattrib)( AstObject *, const char *, int * );
263 static void (* parent_setunit)( AstFrame *, int, const char *, int * );
264 static void (* parent_clearattrib)( AstObject *, const char *, int * );
265 static void (* parent_overlay)( AstFrame *, const int *, AstFrame *, int * );
266 static void (* parent_setattrib)( AstObject *, const char *, int * );
267 static void (* parent_setsystem)( AstFrame *, AstSystemType, int * );
268 static void (* parent_clearsystem)( AstFrame *, int * );
269 static void (* parent_clearunit)( AstFrame *, int, int * );
270 
271 /* Define a variable to hold a SkyFrame which will be used for formatting
272    and unformatting sky positions, etc. */
273 static AstSkyFrame *skyframe;
274 
275 /* Define macros for accessing each item of thread specific global data. */
276 #ifdef THREAD_SAFE
277 
278 /* Define how to initialise thread-specific globals. */
279 #define GLOBAL_inits \
280    globals->Class_Init = 0; \
281    globals->GetAttrib_Buff[ 0 ] = 0; \
282    globals->GetLabel_Buff[ 0 ] = 0; \
283    globals->GetSymbol_Buff[ 0 ] = 0; \
284    globals->GetTitle_Buff[ 0 ] = 0; \
285 
286 /* Create the function that initialises global data for this module. */
287 astMAKE_INITGLOBALS(SpecFrame)
288 
289 /* Define macros for accessing each item of thread specific global data. */
290 #define class_init astGLOBAL(SpecFrame,Class_Init)
291 #define class_vtab astGLOBAL(SpecFrame,Class_Vtab)
292 #define getattrib_buff astGLOBAL(SpecFrame,GetAttrib_Buff)
293 #define getlabel_buff astGLOBAL(SpecFrame,GetLabel_Buff)
294 #define getsymbol_buff astGLOBAL(SpecFrame,GetSymbol_Buff)
295 #define gettitle_buff astGLOBAL(SpecFrame,GetTitle_Buff)
296 
297 
298 
299 static pthread_mutex_t mutex2 = PTHREAD_MUTEX_INITIALIZER;
300 #define LOCK_MUTEX2 pthread_mutex_lock( &mutex2 );
301 #define UNLOCK_MUTEX2 pthread_mutex_unlock( &mutex2 );
302 
303 /* If thread safety is not needed, declare and initialise globals at static
304    variables. */
305 #else
306 
307 /* Buffer returned by GetAttrib. */
308 static char getattrib_buff[ 51 ];
309 
310 /* Default GetLabel string buffer */
311 static char getlabel_buff[ 201 ];
312 
313 /* Default GetSymbol buffer */
314 static char getsymbol_buff[ 21 ];
315 
316 /* Default Title string buffer */
317 static char gettitle_buff[ 201 ];
318 
319 
320 /* Define the class virtual function table and its initialisation flag
321    as static variables. */
322 static AstSpecFrameVtab class_vtab;   /* Virtual function table */
323 static int class_init = 0;       /* Virtual function table initialised? */
324 
325 #define LOCK_MUTEX2
326 #define UNLOCK_MUTEX2
327 
328 #endif
329 
330 
331 /* Prototypes for Private Member Functions. */
332 /* ======================================== */
333 static AstStdOfRestType StdOfRestCode( const char *, int * );
334 static int GetObjSize( AstObject *, int * );
335 static AstSystemType GetAlignSystem( AstFrame *, int * );
336 static AstSystemType SystemCode( AstFrame *, const char *, int * );
337 static AstSystemType ValidateSystem( AstFrame *, AstSystemType, const char *, int * );
338 static const char *DefUnit( AstSystemType, const char *, const char *, int * );
339 static const char *GetDomain( AstFrame *, int * );
340 static const char *GetLabel( AstFrame *, int, int * );
341 static const char *GetSymbol( AstFrame *, int, int * );
342 static const char *GetTitle( AstFrame *, int * );
343 static const char *GetUnit( AstFrame *, int, int * );
344 static const char *SpecMapUnit( AstSystemType, const char *, const char *, int * );
345 static const char *StdOfRestString( AstStdOfRestType, int * );
346 static const char *SystemLabel( AstSystemType, int * );
347 static const char *SystemString( AstFrame *, AstSystemType, int * );
348 static double ConvertSourceVel( AstSpecFrame *, AstStdOfRestType, AstSystemType, int * );
349 static int EqualSor( AstSpecFrame *, AstSpecFrame *, int * );
350 static int GetActiveUnit( AstFrame *, int * );
351 static int MakeSpecMapping( AstSpecFrame *, AstSpecFrame *, AstSpecFrame *, int, AstMapping **, int * );
352 static int Match( AstFrame *, AstFrame *, int, int **, int **, AstMapping **, AstFrame **, int * );
353 static int SorConvert( AstSpecFrame *, AstSpecFrame *, AstSpecMap *, int * );
354 static int SubFrame( AstFrame *, AstFrame *, int, const int *, const int *, AstMapping **, AstFrame **, int * );
355 static int TestActiveUnit( AstFrame *, int * );
356 static void ClearUnit( AstFrame *, int, int * );
357 static void Copy( const AstObject *, AstObject *, int * );
358 static void Delete( AstObject *, int * );
359 static void Dump( AstObject *, AstChannel *, int * );
360 static void GetRefPos( AstSpecFrame *, AstSkyFrame *, double *, double *, int * );
361 static void Overlay( AstFrame *, const int *, AstFrame *, int * );
362 static void SetRefPos( AstSpecFrame *, AstSkyFrame *, double, double, int * );
363 static void SetUnit( AstFrame *, int, const char *, int * );
364 static void VerifyAttrs( AstSpecFrame *, const char *, const char *, const char *, int * );
365 static double ToUnits( AstSpecFrame *, const char *, double, const char *, int * );
366 static void OriginStdOfRest( AstSpecFrame *, AstStdOfRestType, const char *, int * );
367 static void OriginSystem( AstSpecFrame *, AstSystemType, const char *, int * );
368 
369 static AstSystemType GetSystem( AstFrame *, int * );
370 static void SetSystem( AstFrame *, AstSystemType, int * );
371 static void ClearSystem( AstFrame *, int * );
372 
373 static const char *GetAttrib( AstObject *, const char *, int * );
374 static int TestAttrib( AstObject *, const char *, int * );
375 static void ClearAttrib( AstObject *, const char *, int * );
376 static void SetAttrib( AstObject *, const char *, int * );
377 
378 static AstStdOfRestType GetAlignStdOfRest( AstSpecFrame *, int * );
379 static int TestAlignStdOfRest( AstSpecFrame *, int * );
380 static void ClearAlignStdOfRest( AstSpecFrame *, int * );
381 static void SetAlignStdOfRest( AstSpecFrame *, AstStdOfRestType, int * );
382 
383 static AstStdOfRestType GetStdOfRest( AstSpecFrame *, int * );
384 static int TestStdOfRest( AstSpecFrame *, int * );
385 static void ClearStdOfRest( AstSpecFrame *, int * );
386 static void SetStdOfRest( AstSpecFrame *, AstStdOfRestType, int * );
387 
388 static double GetRestFreq( AstSpecFrame *, int * );
389 static int TestRestFreq( AstSpecFrame *, int * );
390 static void ClearRestFreq( AstSpecFrame *, int * );
391 static void SetRestFreq( AstSpecFrame *, double, int * );
392 
393 static double GetSourceVel( AstSpecFrame *, int * );
394 static int TestSourceVel( AstSpecFrame *, int * );
395 static void ClearSourceVel( AstSpecFrame *, int * );
396 static void SetSourceVel( AstSpecFrame *, double, int * );
397 
398 static double GetRefRA( AstSpecFrame *, int * );
399 static int TestRefRA( AstSpecFrame *, int * );
400 static void ClearRefRA( AstSpecFrame *, int * );
401 static void SetRefRA( AstSpecFrame *, double, int * );
402 
403 static double GetRefDec( AstSpecFrame *, int * );
404 static int TestRefDec( AstSpecFrame *, int * );
405 static void ClearRefDec( AstSpecFrame *, int * );
406 static void SetRefDec( AstSpecFrame *, double, int * );
407 
408 static AstStdOfRestType GetSourceVRF( AstSpecFrame *, int * );
409 static int TestSourceVRF( AstSpecFrame *, int * );
410 static void ClearSourceVRF( AstSpecFrame *, int * );
411 static void SetSourceVRF( AstSpecFrame *, AstStdOfRestType, int * );
412 
413 static AstSystemType GetSourceSys( AstSpecFrame *, int * );
414 static int TestSourceSys( AstSpecFrame *, int * );
415 static void ClearSourceSys( AstSpecFrame *, int * );
416 static void SetSourceSys( AstSpecFrame *, AstSystemType, int * );
417 
418 static double GetSpecOrigin( AstSpecFrame *, int * );
419 static int TestSpecOrigin( AstSpecFrame *, int * );
420 static void ClearSpecOrigin( AstSpecFrame *, int * );
421 static void SetSpecOrigin( AstSpecFrame *, double, int * );
422 static double GetSpecOriginCur( AstSpecFrame *, int * );
423 
424 static int GetAlignSpecOffset( AstSpecFrame *, int * );
425 static int TestAlignSpecOffset( AstSpecFrame *, int * );
426 static void SetAlignSpecOffset( AstSpecFrame *, int, int * );
427 static void ClearAlignSpecOffset( AstSpecFrame *, int * );
428 
429 /* Member functions. */
430 /* ================= */
431 
ClearAttrib(AstObject * this_object,const char * attrib,int * status)432 static void ClearAttrib( AstObject *this_object, const char *attrib, int *status ) {
433 /*
434 *  Name:
435 *     ClearAttrib
436 
437 *  Purpose:
438 *     Clear an attribute value for a SpecFrame.
439 
440 *  Type:
441 *     Private function.
442 
443 *  Synopsis:
444 *     #include "specframe.h"
445 *     void ClearAttrib( AstObject *this, const char *attrib, int *status )
446 
447 *  Class Membership:
448 *     SpecFrame member function (over-rides the astClearAttrib protected
449 *     method inherited from the Frame class).
450 
451 *  Description:
452 *     This function clears the value of a specified attribute for a
453 *     SpecFrame, so that the default value will subsequently be used.
454 
455 *  Parameters:
456 *     this
457 *        Pointer to the SpecFrame.
458 *     attrib
459 *        Pointer to a null terminated string specifying the attribute
460 *        name.  This should be in lower case with no surrounding white
461 *        space.
462 *     status
463 *        Pointer to the inherited status variable.
464 
465 *  Notes:
466 *     - This function uses one-based axis numbering so that it is
467 *     suitable for external (public) use.
468 */
469 
470 /* Local Variables: */
471    AstSpecFrame *this;           /* Pointer to the SpecFrame structure */
472    char *new_attrib;             /* Pointer value to new attribute name */
473    int len;                      /* Length of attrib string */
474 
475 /* Check the global error status. */
476    if ( !astOK ) return;
477 
478 /* Obtain a pointer to the SpecFrame structure. */
479    this = (AstSpecFrame *) this_object;
480 
481 /* Obtain the length of the "attrib" string. */
482    len = strlen( attrib );
483 
484 /* Check the attribute name and clear the appropriate attribute. */
485 
486 /* First look for axis attributes defined by the Frame class. Since a
487    SpecFrame has only 1 axis, we allow these attributes to be specified
488    without a trailing "(axis)" string. */
489    if ( !strcmp( attrib, "direction" ) ||
490         !strcmp( attrib, "bottom" ) ||
491         !strcmp( attrib, "top" ) ||
492         !strcmp( attrib, "format" ) ||
493         !strcmp( attrib, "label" ) ||
494         !strcmp( attrib, "symbol" ) ||
495         !strcmp( attrib, "unit" ) ) {
496 
497 /* Create a new attribute name from the original by appending the string
498    "(1)" and then use the parent ClearAttrib method. */
499       new_attrib = astMalloc( len + 4 );
500       if( new_attrib ) {
501          memcpy( new_attrib, attrib, len );
502          memcpy( new_attrib + len, "(1)", 4 );
503          (*parent_clearattrib)( this_object, new_attrib, status );
504          new_attrib = astFree( new_attrib );
505       }
506 
507 /* AlignStdOfRest. */
508 /* --------------- */
509    } else if ( !strcmp( attrib, "alignstdofrest" ) ) {
510       astClearAlignStdOfRest( this );
511 
512 /* GeoLat. */
513 /* ------- */
514 /* Retained for backward compatibility with older versions of AST in which
515    SpecFrame had GeoLon/Lat attributes (now ObsLon/Lat are used instead). */
516    } else if ( !strcmp( attrib, "geolat" ) ) {
517       astClearAttrib( this, "obslat" );
518 
519 /* GeoLon. */
520 /* ------- */
521    } else if ( !strcmp( attrib, "geolon" ) ) {
522       astClearAttrib( this, "obslon" );
523 
524 /* RefDec. */
525 /* ---------- */
526    } else if ( !strcmp( attrib, "refdec" ) ) {
527       astClearRefDec( this );
528 
529 /* RefRA. */
530 /* --------- */
531    } else if ( !strcmp( attrib, "refra" ) ) {
532       astClearRefRA( this );
533 
534 /* RestFreq. */
535 /* --------- */
536    } else if ( !strcmp( attrib, "restfreq" ) ) {
537       astClearRestFreq( this );
538 
539 /* SourceVel. */
540 /* ---------- */
541    } else if ( !strcmp( attrib, "sourcevel" ) ) {
542       astClearSourceVel( this );
543 
544 /* SpecOrigin. */
545 /* ---------- */
546    } else if ( !strcmp( attrib, "specorigin" ) ) {
547       astClearSpecOrigin( this );
548 
549 /* AlignSpecOffset. */
550 /* ---------------- */
551    } else if ( !strcmp( attrib, "alignspecoffset" ) ) {
552       astClearAlignSpecOffset( this );
553 
554 /* SourceVRF */
555 /* --------- */
556    } else if ( !strcmp( attrib, "sourcevrf" ) ) {
557       astClearSourceVRF( this );
558 
559 /* SourceSys */
560 /* --------- */
561    } else if ( !strcmp( attrib, "sourcesys" ) ) {
562       astClearSourceSys( this );
563 
564 /* StdOfRest. */
565 /* ---------- */
566    } else if ( !strcmp( attrib, "stdofrest" ) ) {
567       astClearStdOfRest( this );
568 
569 /* If the attribute is not recognised, pass it on to the parent method
570    for further interpretation. */
571    } else {
572       (*parent_clearattrib)( this_object, attrib, status );
573    }
574 }
575 
ClearSystem(AstFrame * this_frame,int * status)576 static void ClearSystem( AstFrame *this_frame, int *status ) {
577 /*
578 *  Name:
579 *     ClearSystem
580 
581 *  Purpose:
582 *     Clear the System attribute for a SpecFrame.
583 
584 *  Type:
585 *     Private function.
586 
587 *  Synopsis:
588 *     #include "specframe.h"
589 *     void ClearSystem( AstFrame *this_frame, int *status )
590 
591 *  Class Membership:
592 *     SpecFrame member function (over-rides the astClearSystem protected
593 *     method inherited from the Frame class).
594 
595 *  Description:
596 *     This function clears the System attribute for a SpecFrame.
597 
598 *  Parameters:
599 *     this
600 *        Pointer to the SpecFrame.
601 *     status
602 *        Pointer to the inherited status variable.
603 
604 */
605 
606 /* Local Variables: */
607    AstSpecFrame *this;           /* Pointer to SpecFrame structure */
608    AstSystemType newsys;         /* System after clearing */
609    AstSystemType oldsys;         /* System before clearing */
610 
611 /* Check the global error status. */
612    if ( !astOK ) return;
613 
614 /* Obtain a pointer to the SpecFrame structure. */
615    this = (AstSpecFrame *) this_frame;
616 
617 /* Save the original system */
618    oldsys = astGetSystem( this_frame );
619 
620 /* Use the parent ClearSystem method to clear the System value. */
621    (*parent_clearsystem)( this_frame, status );
622 
623 /* Get the default System. */
624    newsys = astGetSystem( this_frame );
625 
626 /* If the system has actually changed. */
627    if( newsys != oldsys ) {
628 
629 /* Changing the System value will in general require the Units to change
630    as well. If the used has previously specified the units to be used with
631    the new system, then re-instate them (they are stored in the "usedunits"
632    array in the SpecFrame structure). Otherwise, clear the units so that
633    the default units will eb used with the new System. */
634       if( (int) newsys < this->nuunits && this->usedunits &&
635           this->usedunits[ (int) newsys ] ) {
636          astSetUnit( this, 0, this->usedunits[ (int) newsys ] );
637       } else {
638          astClearUnit( this, 0 );
639       }
640 
641 /* Also, clear all attributes which have system-specific defaults. */
642       astClearLabel( this_frame, 0 );
643       astClearSymbol( this_frame, 0 );
644       astClearTitle( this_frame );
645 
646 /* Modify the SpecOrigin value to use the new System */
647       OriginSystem( this, oldsys, "astClearSystem", status );
648 
649    }
650 
651 }
652 
ClearStdOfRest(AstSpecFrame * this,int * status)653 static void ClearStdOfRest( AstSpecFrame *this, int *status ) {
654 /*
655 *+
656 *  Name:
657 *     astClearStdOfRest
658 
659 *  Purpose:
660 *     Clear the StdOfRest attribute for a SpecFrame.
661 
662 *  Type:
663 *     Protected function.
664 
665 *  Synopsis:
666 *     #include "timeframe.h"
667 *     void astClearStdOfRest( AstSpecFrame *this )
668 
669 *  Class Membership:
670 *     SpecFrame virtual function
671 
672 *  Description:
673 *     This function clears the StdOfRest attribute for a SpecFrame.
674 
675 *  Parameters:
676 *     this
677 *        Pointer to the SpecFrame.
678 
679 *-
680 */
681 
682 /* Check the global error status. */
683    if ( !astOK ) return;
684 
685 /* Modify the SpecOrigin value stored in the SpecFrame structure to refer to the
686    default rest frame (heliocentric). */
687    OriginStdOfRest( this, AST__HLSOR, "astClearStdOfRest", status );
688 
689 /* Store a bad value for the standard of rest in the SpecFrame structure. */
690    this->stdofrest = AST__BADSOR;
691 }
692 
693 
ClearUnit(AstFrame * this_frame,int axis,int * status)694 static void ClearUnit( AstFrame *this_frame, int axis, int *status ) {
695 /*
696 *  Name:
697 *     ClearUnit
698 
699 *  Purpose:
700 *     Clear the value of the Unit string for a SpecFrame's axis.
701 
702 *  Type:
703 *     Private function.
704 
705 *  Synopsis:
706 *     #include "specframe.h"
707 *     void ClearUnit( AstFrame *this_frame, int axis )
708 
709 *  Class Membership:
710 *     SpecFrame member function (over-rides the astClearUnit method inherited
711 *     from the Frame class).
712 
713 *  Description:
714 *     This function clears the Unit string for a specified axis of a
715 *     SpecFrame. It also clears the UsedUnit item in the SpecFrame
716 *     structure corresponding to the current System.
717 
718 *  Parameters:
719 *     this
720 *        Pointer to the SpecFrame.
721 *     axis
722 *        The number of the axis (zero-based).
723 */
724 
725 /* Local Variables: */
726    AstSpecFrame *this;           /* Pointer to the SpecFrame structure */
727    int system;                   /* The SpecFrame's System value */
728 
729 /* Check the global error status. */
730    if ( !astOK ) return;
731 
732 /* Obtain a pointer to the SpecFrame structure. */
733    this = (AstSpecFrame *) this_frame;
734 
735 /* Validate the axis index. */
736    astValidateAxis( this, axis, 1, "astClearUnit" );
737 
738 /* Clear the UsedUnit item for the current System, if current set. */
739    system = (int) astGetSystem( this );
740    if( system < this->nuunits && this->usedunits ) {
741       this->usedunits[ system ] = astFree( this->usedunits[ system ] );
742    }
743 
744 /* Use the parent method to clear the Unit attribute of the axis. */
745    (*parent_clearunit)( this_frame, axis, status );
746 }
747 
ConvertSourceVel(AstSpecFrame * this,AstStdOfRestType newsor,AstSystemType newsys,int * status)748 static double ConvertSourceVel( AstSpecFrame *this, AstStdOfRestType newsor,
749                                 AstSystemType newsys, int *status ) {
750 /*
751 *  Name:
752 *     ConvertSourceVel
753 
754 *  Purpose:
755 *     Convert the SourceVel value to a specified rest frame and spectral
756 *     system.
757 
758 *  Type:
759 *     Private function.
760 
761 *  Synopsis:
762 *     #include "specframe.h"
763 *     double ConvertSourceVel( AstSpecFrame *this, AstStdOfRestType newsor,
764 *                              AstSystemType newsys, int *status )
765 
766 *  Class Membership:
767 *     SpecFrame member function
768 
769 *  Description:
770 *     This function convert the SourceVel value to a specified rest frame
771 *     and spectral system, and returns the new value.
772 
773 *  Parameters:
774 *     this
775 *        Pointer to the SpecFrame.
776 *     newsor
777 *        The rest frame in which the source velocity is required.
778 *     newsys
779 *        The spectral system (AST__VREL or AST__REDSHIFT) in which the
780 *        source velocity is required.
781 *     status
782 *        Pointer to the inherited status variable.
783 
784 *  Returned Value:
785 *     The converted source velocity (m/s), or redshift.
786 
787 *  Notes:
788 *     - This function returns zero if an error occurs.
789 */
790 
791 /* Local Variables: */
792    AstSpecFrame *from;     /* Pointer to a source SpecFrame */
793    AstSpecFrame *to;       /* Pointer to a destination SpecFrame */
794    AstSpecMap *specmap;    /* Pointer to a SpecMap */
795    AstStdOfRestType sor;   /* Standard of rest in which SourceVel is defined */
796    AstSystemType sys;      /* Spectral system in which SourceVel is defined */
797    double ret;             /* The returned value */
798    double rf;              /* Rest frequency (Hz) */
799    double temp;            /* Temporary storage */
800 
801 /* Initialise */
802    ret = 0.0;
803 
804 /* Check the global error status. */
805    if ( !astOK ) return ret;
806 
807 /* Get the value of the SourceVel attribute. This will be a velocity in m/s
808    (relativistic, radio or optical), or unitless redshift or beta factor,
809    depending on the current value of SourceSys. */
810    ret = astGetSourceVel( this );
811 
812 /* Check it can be used (depends on whether a value has been set and
813    whether the UseDefs attribute is zero). */
814    VerifyAttrs( this, "convert source velocity to a new standard of rest",
815                 "SourceVel", "astMatch", status );
816 
817 /* Get the rest frame and spectral system to which value refers. */
818    sor = astGetSourceVRF( this );
819    sys = astGetSourceSys( this );
820 
821 /* If necessary, convert to the requested rest frame and spectral system. */
822    if( sor != newsor || sys != newsys ) {
823 
824 /* Verify that usable value is available for the RestFreq attribute. An
825    error is reported if not. */
826       VerifyAttrs( this, "convert source velocity to a new standard of rest",
827                    "RestFreq", "astMatch", status );
828 
829 /* Take two copies of the supplied SpecFrame and set their StdOfRest
830    attributes to the required values. */
831       from = astCopy( this );
832       astSetStdOfRest( from, sor );
833 
834       to = astCopy( this );
835       astSetStdOfRest( to, newsor );
836 
837 /* Initialise a new SpecMap to describe the conversion. The new SpecMap
838    initially represents a UnitMap. */
839       specmap = astSpecMap( 1, 0, "", status );
840 
841 /* Add a conversion from the spectral system in which the SourceVEl value
842    is stored, to relativistic velocity. */
843       if( sys == AST__VRADIO ) {
844          astSpecAdd( specmap, "VRTOVL", NULL );
845 
846       } else if( sys == AST__VOPTICAL ) {
847          astSpecAdd( specmap, "VOTOVL", NULL );
848 
849       } else if( sys == AST__REDSHIFT ) {
850          astSpecAdd( specmap, "ZOTOVL", NULL );
851 
852       } else if( sys == AST__BETA ) {
853          astSpecAdd( specmap, "BTTOVL", NULL );
854       }
855 
856 /* Add a conversion from velocity to frequency since SorConvert converts
857    frequencies. */
858       rf = astGetRestFreq( this );
859       astSpecAdd( specmap, "VLTOFR", &rf );
860 
861 /* Now add a conversion from frequency in the SourveVRF standard of rest to
862    frequency in the required rest frame. */
863       SorConvert( from, to, specmap, status );
864 
865 /* Add a conversion from frequency back to velocity. Note, the value of the
866    rest frequency does not affect the overall conversion. */
867       astSpecAdd( specmap, "FRTOVL", &rf );
868 
869 /* Add a conversion from relativistic velocity to the required spectral
870    system, if needed. */
871       if( newsys == AST__VRADIO ) {
872          astSpecAdd( specmap, "VLTOVR", NULL );
873 
874       } else if( newsys == AST__VOPTICAL ) {
875          astSpecAdd( specmap, "VLTOVO", NULL );
876 
877       } else if( newsys == AST__REDSHIFT ) {
878          astSpecAdd( specmap, "VLTOZO", NULL );
879 
880       } else if( newsys == AST__BETA ) {
881          astSpecAdd( specmap, "VLTOBT", NULL );
882       }
883 
884 /* Use the SpecMap to convert the source velocity in the SourceVRF
885    standard of rest and SourceSys spectral system to the required rest
886    frame and spectral system. */
887       temp = ret;
888       astTran1( specmap, 1, &temp, 1, &ret );
889 
890 /* Free resources */
891       specmap = astAnnul( specmap );
892       to = astAnnul( to );
893       from = astAnnul( from );
894    }
895 
896 /* Return zero if an error has occurred. */
897    if( !astOK ) ret = 0.0;
898 
899 /* Return the answer. */
900    return ret;
901 
902 }
903 
DefUnit(AstSystemType system,const char * method,const char * class,int * status)904 static const char *DefUnit( AstSystemType system, const char *method,
905                             const char *class, int *status ){
906 /*
907 *  Name:
908 *     DefUnit
909 
910 *  Purpose:
911 *     Return the default units for a spectral coordinate system type.
912 
913 *  Type:
914 *     Private function.
915 
916 *  Synopsis:
917 *     #include "specframe.h"
918 *     const char *DefUnit( AstSystemType system, const char *method,
919 *                          const char *class, int *status )
920 
921 *  Class Membership:
922 *     SpecFrame member function.
923 
924 *  Description:
925 *     This function returns a textual representation of the default
926 *     units associated with the specified spectral coordinate system.
927 
928 *  Parameters:
929 *     system
930 *        The spectral coordinate system.
931 *     method
932 *        Pointer to a string holding the name of the calling method.
933 *        This is only for use in constructing error messages.
934 *     class
935 *        Pointer to a string holding the name of the supplied object class.
936 *        This is only for use in constructing error messages.
937 *     status
938 *        Pointer to the inherited status variable.
939 
940 *  Returned Value:
941 *     As tring describing the default units. This string follows the
942 *     units syntax described in FITS WCS paper I "Representations of world
943 *     coordinates in FITS" (Greisen & Calabretta).
944 
945 *  Notes:
946 *     - A NULL pointer is returned if this function is invoked with
947 *     the global error status set or if it should fail for any reason.
948 */
949 
950 /* Local Variables: */
951    const char *result;           /* Value to return */
952 
953 /* Initialize */
954    result = NULL;
955 
956 /* Check the global error status. */
957    if ( !astOK ) return result;
958 
959 /* Get an identifier for the default units. */
960    if( system == AST__FREQ ) {
961       result = "GHz";
962    } else if( system == AST__ENERGY ) {
963       result = "J";
964    } else if( system == AST__WAVENUM ) {
965       result = "1/m";
966    } else if( system == AST__WAVELEN ) {
967       result = "Angstrom";
968    } else if( system == AST__AIRWAVE ) {
969       result = "Angstrom";
970    } else if( system == AST__VRADIO ) {
971       result = "km/s";
972    } else if( system == AST__VOPTICAL ) {
973       result = "km/s";
974    } else if( system == AST__REDSHIFT ) {
975       result = "";
976    } else if( system == AST__BETA ) {
977       result = "";
978    } else if( system == AST__VREL ) {
979       result = "km/s";
980 
981 /* Report an error if the coordinate system was not recognised. */
982    } else {
983       astError( AST__SCSIN, "%s(%s): Corrupt %s contains illegal System "
984                 "identification code (%d).", status, method, class, class,
985                 (int) system );
986    }
987 
988 /* Return the result. */
989    return result;
990 }
991 
EqualSor(AstSpecFrame * this,AstSpecFrame * that,int * status)992 static int EqualSor( AstSpecFrame *this, AstSpecFrame *that, int *status ) {
993 /*
994 *  Name:
995 *     EqualSor
996 
997 *  Purpose:
998 *     Do two SpecFrames use the same standard of rest?
999 
1000 *  Type:
1001 *     Private function.
1002 
1003 *  Synopsis:
1004 *     #include "specframe.h"
1005 *     int EqualSor( AstSpecFrame *this, AstSpecFrame *that, int *status )
1006 
1007 *  Class Membership:
1008 *     SpecFrame member function
1009 
1010 *  Description:
1011 *    This function returns non-zero if the two supplied SpecFrames use
1012 *    the same standard of rest.
1013 
1014 *  Parameters:
1015 *     this
1016 *        Pointer to the first SpecFrame.
1017 *     that
1018 *        Pointer to the second SpecFrame.
1019 *     status
1020 *        Pointer to the inherited status variable.
1021 
1022 *  Returned Value:
1023 *     Non-zero if the two SpecFrames use the same standard of rest. Zero
1024 *     otherwise.
1025 
1026 */
1027 
1028 /* Local Variables: */
1029    AstStdOfRestType sor;             /* Standard of rest */
1030    int result;                       /* Value to return */
1031 
1032 /* Check the global error status. */
1033    if ( !astOK ) return 0;
1034 
1035 /* Initialise. */
1036    result = 1;
1037 
1038 /* Compare StdOfRest attributes. */
1039    sor = astGetStdOfRest( this );
1040    if( astGetStdOfRest( that ) != sor ) {
1041       result = 0;
1042 
1043 /* If the standards of rest are equal we need to check the the attributes
1044    which specify the precise rest frame. */
1045    } else {
1046 
1047 /* The reference RA and Dec need to be equal */
1048       if( !EQUAL( astGetRefRA( this ), astGetRefRA( that ) ) ||
1049           !EQUAL( astGetRefDec( this ), astGetRefDec( that ) ) ) {
1050          result = 0;
1051 
1052 /* For source rest frame, the source velocities, rest frames and systems must
1053    be equal */
1054       } else if( sor == AST__SCSOR ){
1055          if( !EQUAL( astGetSourceVel( this ), astGetSourceVel( that ) ) ||
1056                      astGetSourceVRF( this ) != astGetSourceVRF( that ) ||
1057                      astGetSourceSys( this ) != astGetSourceSys( that ) ) {
1058             result = 0;
1059          }
1060 
1061 /* For geocentric, barycentric and heliocentric rest frames, the epochs must
1062    be the same */
1063       } else if( sor == AST__GESOR || sor == AST__BYSOR || sor == AST__HLSOR ){
1064          if( !EQUAL( astGetEpoch( this ), astGetEpoch( that ) ) ) result = 0;
1065 
1066 /* For topocentric rest frame, the epoch and position of the observer must be
1067    the same */
1068       } else if( sor == AST__TPSOR ){
1069          if( !EQUAL( astGetEpoch( this ), astGetEpoch( that ) ) ||
1070              !EQUAL( astGetObsAlt( this ), astGetObsAlt( that ) ) ||
1071              !EQUAL( astGetObsLon( this ), astGetObsLon( that ) ) ||
1072              !EQUAL( astGetObsLat( this ), astGetObsLat( that ) ) ) result = 0;
1073 
1074       } else if( sor != AST__LKSOR && sor != AST__LDSOR &&
1075                  sor != AST__GLSOR && sor != AST__LGSOR && astOK ) {
1076          astError( AST__INTER, "SorEqual(SpecFrame): Function SorEqual "
1077                    "does not yet support rest frame %d (AST internal "
1078                    "programming error)", status, sor );
1079       }
1080    }
1081 
1082 /* Return the result */
1083    return result;
1084 }
1085 
GetObjSize(AstObject * this_object,int * status)1086 static int GetObjSize( AstObject *this_object, int *status ) {
1087 /*
1088 *  Name:
1089 *     GetObjSize
1090 
1091 *  Purpose:
1092 *     Return the in-memory size of an Object.
1093 
1094 *  Type:
1095 *     Private function.
1096 
1097 *  Synopsis:
1098 *     #include "specframe.h"
1099 *     int GetObjSize( AstObject *this, int *status )
1100 
1101 *  Class Membership:
1102 *     SpecFrame member function (over-rides the astGetObjSize protected
1103 *     method inherited from the parent class).
1104 
1105 *  Description:
1106 *     This function returns the in-memory size of the supplied SpecFrame,
1107 *     in bytes.
1108 
1109 *  Parameters:
1110 *     this
1111 *        Pointer to the SpecFrame.
1112 *     status
1113 *        Pointer to the inherited status variable.
1114 
1115 *  Returned Value:
1116 *     The Object size, in bytes.
1117 
1118 *  Notes:
1119 *     - A value of zero will be returned if this function is invoked
1120 *     with the global status set, or if it should fail for any reason.
1121 */
1122 
1123 /* Local Variables: */
1124    AstSpecFrame *this;         /* Pointer to SpecFrame structure */
1125    int result;                /* Result value to return */
1126    int i;
1127 
1128 /* Initialise. */
1129    result = 0;
1130 
1131 /* Check the global error status. */
1132    if ( !astOK ) return result;
1133 
1134 /* Obtain a pointers to the SpecFrame structure. */
1135    this = (AstSpecFrame *) this_object;
1136 
1137 /* Invoke the GetObjSize method inherited from the parent class, and then
1138    add on any components of the class structure defined by thsi class
1139    which are stored in dynamically allocated memory. */
1140    result = (*parent_getobjsize)( this_object, status );
1141    if( this->usedunits ) {
1142       for( i = 0; i < this->nuunits; i++ ) {
1143          result += astTSizeOf( this->usedunits[ i ] );
1144       }
1145       result += astTSizeOf( this->usedunits );
1146    }
1147 
1148 /* If an error occurred, clear the result value. */
1149    if ( !astOK ) result = 0;
1150 
1151 /* Return the result, */
1152    return result;
1153 }
1154 
GetActiveUnit(AstFrame * this_frame,int * status)1155 static int GetActiveUnit( AstFrame *this_frame, int *status ) {
1156 /*
1157 *  Name:
1158 *     GetActiveUnit
1159 
1160 *  Purpose:
1161 *     Obtain the value of the ActiveUnit flag for a SpecFrame.
1162 
1163 *  Type:
1164 *     Private function.
1165 
1166 *  Synopsis:
1167 *     #include "specframe.h"
1168 *     int GetActiveUnit( AstFrame *this_frame, int *status )
1169 
1170 *  Class Membership:
1171 *     SpecFrame member function (over-rides the astGetActiveUnit protected
1172 *     method inherited from the Frame class).
1173 
1174 *  Description:
1175 *    This function returns the value of the ActiveUnit flag for a
1176 *    SpecFrame, which is always 1.
1177 
1178 *  Parameters:
1179 *     this
1180 *        Pointer to the SpecFrame.
1181 *     status
1182 *        Pointer to the inherited status variable.
1183 
1184 *  Returned Value:
1185 *     The value to use for the ActiveUnit flag (1).
1186 
1187 */
1188    return 1;
1189 }
1190 
GetAttrib(AstObject * this_object,const char * attrib,int * status)1191 static const char *GetAttrib( AstObject *this_object, const char *attrib, int *status ) {
1192 /*
1193 *  Name:
1194 *     GetAttrib
1195 
1196 *  Purpose:
1197 *     Get the value of a specified attribute for a SpecFrame.
1198 
1199 *  Type:
1200 *     Private function.
1201 
1202 *  Synopsis:
1203 *     #include "specframe.h"
1204 *     const char *GetAttrib( AstObject *this, const char *attrib, int *status )
1205 
1206 *  Class Membership:
1207 *     SpecFrame member function (over-rides the protected astGetAttrib
1208 *     method inherited from the Frame class).
1209 
1210 *  Description:
1211 *     This function returns a pointer to the value of a specified
1212 *     attribute for a SpecFrame, formatted as a character string.
1213 
1214 *  Parameters:
1215 *     this
1216 *        Pointer to the SpecFrame.
1217 *     attrib
1218 *        Pointer to a null-terminated string containing the name of
1219 *        the attribute whose value is required. This name should be in
1220 *        lower case, with all white space removed.
1221 *     status
1222 *        Pointer to the inherited status variable.
1223 
1224 *  Returned Value:
1225 *     - Pointer to a null-terminated string containing the attribute
1226 *     value.
1227 
1228 *  Notes:
1229 *     - This function uses one-based axis numbering so that it is
1230 *     suitable for external (public) use.
1231 *     - The returned string pointer may point at memory allocated
1232 *     within the SpecFrame, or at static memory. The contents of the
1233 *     string may be over-written or the pointer may become invalid
1234 *     following a further invocation of the same function or any
1235 *     modification of the SpecFrame. A copy of the string should
1236 *     therefore be made if necessary.
1237 *     - A NULL pointer will be returned if this function is invoked
1238 *     with the global error status set, or if it should fail for any
1239 *     reason.
1240 */
1241 
1242 /* Local Variables: */
1243    astDECLARE_GLOBALS            /* Declare the thread specific global data */
1244    AstSpecFrame *this;           /* Pointer to the SpecFrame structure */
1245    AstStdOfRestType sor;         /* Standard of rest */
1246    AstSystemType sys;            /* Spectral system */
1247    char *new_attrib;             /* Pointer value to new attribute name */
1248    const char *result;           /* Pointer value to return */
1249    double dval;                  /* Attribute value */
1250    int ival;                     /* Attribute value */
1251    int len;                      /* Length of attrib string */
1252 
1253 /* Initialise. */
1254    result = NULL;
1255 
1256 /* Check the global error status. */
1257    if ( !astOK ) return result;
1258 
1259 /* Get a pointer to the structure holding thread-specific global data. */
1260    astGET_GLOBALS(this_object);
1261 
1262 /* Obtain a pointer to the SpecFrame structure. */
1263    this = (AstSpecFrame *) this_object;
1264 
1265 /* Create an FK5 J2000 SkyFrame which will be used for formatting and
1266    unformatting sky positions, etc. */
1267    LOCK_MUTEX2
1268    if( !skyframe ) {
1269       astBeginPM;
1270       skyframe = astSkyFrame( "system=FK5,equinox=J2000", status );
1271       astEndPM;
1272    }
1273    UNLOCK_MUTEX2
1274 
1275 /* Obtain the length of the attrib string. */
1276    len = strlen( attrib );
1277 
1278 /* Compare "attrib" with each recognised attribute name in turn,
1279    obtaining the value of the required attribute. If necessary, write
1280    the value into "getattrib_buff" as a null-terminated string in an appropriate
1281    format.  Set "result" to point at the result string. */
1282 
1283 /* First look for axis attributes defined by the Frame class. Since a
1284    SpecFrame has only 1 axis, we allow these attributes to be specified
1285    without a trailing "(axis)" string. */
1286    if ( !strcmp( attrib, "direction" ) ||
1287         !strcmp( attrib, "bottom" ) ||
1288         !strcmp( attrib, "top" ) ||
1289         !strcmp( attrib, "format" ) ||
1290         !strcmp( attrib, "label" ) ||
1291         !strcmp( attrib, "symbol" ) ||
1292         !strcmp( attrib, "unit" ) ) {
1293 
1294 /* Create a new attribute name from the original by appending the string
1295    "(1)" and then use the parent GetAttrib method. */
1296       new_attrib = astMalloc( len + 4 );
1297       if( new_attrib ) {
1298          memcpy( new_attrib, attrib, len );
1299          memcpy( new_attrib + len, "(1)", 4 );
1300          result = (*parent_getattrib)( this_object, new_attrib, status );
1301          new_attrib = astFree( new_attrib );
1302       }
1303 
1304 /* AlignStdOfRest. */
1305 /* --------------- */
1306 /* Obtain the AlignStdOfRest code and convert to a string. */
1307    } else if ( !strcmp( attrib, "alignstdofrest" ) ) {
1308       sor = astGetAlignStdOfRest( this );
1309       if ( astOK ) {
1310          result = StdOfRestString( sor, status );
1311 
1312 /* Report an error if the value was not recognised. */
1313          if ( !result ) {
1314             astError( AST__SCSIN,
1315                      "astGetAttrib(%s): Corrupt %s contains invalid AlignStdOfRest "
1316                      "identification code (%d).", status, astGetClass( this ),
1317                      astGetClass( this ), (int) sor );
1318          }
1319       }
1320 
1321 /* AlignSpecOffset */
1322 /* --------------- */
1323    } else if ( !strcmp( attrib, "alignspecoffset" ) ) {
1324       ival = astGetAlignSpecOffset( this );
1325       if ( astOK ) {
1326          (void) sprintf( getattrib_buff, "%d", ival );
1327          result = getattrib_buff;
1328       }
1329 
1330 /* GeoLat. */
1331 /* ------- */
1332 /* Retained for backward compatibility with older versions of AST in which
1333    SpecFrame had GeoLon/Lat attributes (now ObsLon/Lat are used instead). */
1334    } else if ( !strcmp( attrib, "geolat" ) ) {
1335       result = astGetAttrib( this, "obslat" );
1336 
1337 /* GeoLon. */
1338 /* ------- */
1339    } else if ( !strcmp( attrib, "geolon" ) ) {
1340       result = astGetAttrib( this, "obslon" );
1341 
1342 /* RefDec. */
1343 /* ------- */
1344 /* Convert to a string using the SkyFrame Format method. */
1345    } else if ( !strcmp( attrib, "refdec" ) ) {
1346       dval = astGetRefDec( this );
1347       if ( astOK ) {
1348          result = astFormat( skyframe, 1, dval );
1349       }
1350 
1351 /* RefRA. */
1352 /* ------ */
1353 /* Convert to a string using the SkyFrame Format method. */
1354    } else if ( !strcmp( attrib, "refra" ) ) {
1355       dval = astGetRefRA( this );
1356       if ( astOK ) {
1357          result = astFormat( skyframe, 0, dval );
1358       }
1359 
1360 /* RestFreq. */
1361 /* --------- */
1362    } else if ( !strcmp( attrib, "restfreq" ) ) {
1363       dval = astGetRestFreq( this );
1364       if ( astOK ) {
1365          (void) sprintf( getattrib_buff, "%.*g", DBL_DIG, dval*1.0E-9 );
1366          result = getattrib_buff;
1367       }
1368 
1369 /* SourceVel */
1370 /* --------- */
1371    } else if ( !strcmp( attrib, "sourcevel" ) ) {
1372       dval = astGetSourceVel( this );
1373       if ( astOK ) {
1374 
1375 /* Convert from m/s to km/s if the SourceVel value is a velocity. . */
1376          if( astGetSourceSys( this ) == AST__VREL ||
1377              astGetSourceSys( this ) == AST__VRADIO ||
1378              astGetSourceSys( this ) == AST__VOPTICAL ) dval *= 1.0E-3;
1379 
1380 /* Format */
1381          (void) sprintf( getattrib_buff, "%.*g", DBL_DIG, dval );
1382          result = getattrib_buff;
1383 
1384       }
1385 
1386 /* SpecOrigin. */
1387 /* ----------- */
1388    } else if ( !strcmp( attrib, "specorigin" ) ) {
1389       dval = GetSpecOriginCur( this, status );
1390       if( astOK ) {
1391          (void) sprintf( getattrib_buff, "%.*g", DBL_DIG, dval );
1392          result = getattrib_buff;
1393       }
1394 
1395 
1396 /* SourceVRF */
1397 /* ----------*/
1398    } else if ( !strcmp( attrib, "sourcevrf" ) ) {
1399       sor = astGetSourceVRF( this );
1400       if ( astOK ) {
1401          result = StdOfRestString( sor, status );
1402 
1403 /* Report an error if the value was not recognised. */
1404          if ( !result ) {
1405             astError( AST__SCSIN,
1406                      "astGetAttrib(%s): Corrupt %s contains invalid SourceVRF "
1407                      "identification code (%d).", status, astGetClass( this ),
1408                      astGetClass( this ), (int) sor );
1409          }
1410       }
1411 
1412 /* SourceSys */
1413 /* ----------*/
1414    } else if ( !strcmp( attrib, "sourcesys" ) ) {
1415       sys = astGetSourceSys( this );
1416       if ( astOK ) {
1417          result = SystemString( (AstFrame *) this, sys, status );
1418 
1419 /* Report an error if the value was not recognised. */
1420          if ( !result ) {
1421             astError( AST__SCSIN,
1422                      "astGetAttrib(%s): Corrupt %s contains invalid SourceSys "
1423                      "identification code (%d).", status, astGetClass( this ),
1424                      astGetClass( this ), (int) sys );
1425          }
1426       }
1427 
1428 /* StdOfRest. */
1429 /* ---------- */
1430 /* Obtain the StdOfRest code and convert to a string. */
1431    } else if ( !strcmp( attrib, "stdofrest" ) ) {
1432       sor = astGetStdOfRest( this );
1433       if ( astOK ) {
1434          result = StdOfRestString( sor, status );
1435 
1436 /* Report an error if the value was not recognised. */
1437          if ( !result ) {
1438             astError( AST__SCSIN,
1439                      "astGetAttrib(%s): Corrupt %s contains invalid StdOfRest "
1440                      "identification code (%d).", status, astGetClass( this ),
1441                      astGetClass( this ), (int) sor );
1442          }
1443       }
1444 
1445 /* If the attribute name was not recognised, pass it on to the parent
1446    method for further interpretation. */
1447    } else {
1448       result = (*parent_getattrib)( this_object, attrib, status );
1449    }
1450 
1451 /* Return the result. */
1452    return result;
1453 }
1454 
GetDomain(AstFrame * this_frame,int * status)1455 static const char *GetDomain( AstFrame *this_frame, int *status ) {
1456 /*
1457 *  Name:
1458 *     GetDomain
1459 
1460 *  Purpose:
1461 *     Obtain a pointer to the Domain attribute string for a SpecFrame.
1462 
1463 *  Type:
1464 *     Private function.
1465 
1466 *  Synopsis:
1467 *     #include "specframe.h"
1468 *     const char *GetDomain( AstFrame *this, int *status )
1469 
1470 *  Class Membership:
1471 *     SpecFrame member function (over-rides the astGetDomain protected
1472 *     method inherited from the Frame class).
1473 
1474 *  Description:
1475 *    This function returns a pointer to the Domain attribute string
1476 *    for a SpecFrame.
1477 
1478 *  Parameters:
1479 *     this
1480 *        Pointer to the SpecFrame.
1481 *     status
1482 *        Pointer to the inherited status variable.
1483 
1484 *  Returned Value:
1485 *     A pointer to a constant null-terminated string containing the
1486 *     Domain value.
1487 
1488 *  Notes:
1489 *     - The returned pointer or the string it refers to may become
1490 *     invalid following further invocation of this function or
1491 *     modification of the SpecFrame.
1492 *     - A NULL pointer is returned if this function is invoked with
1493 *     the global error status set or if it should fail for any reason.
1494 */
1495 
1496 /* Local Variables: */
1497    AstSpecFrame *this;           /* Pointer to SpecFrame structure */
1498    const char *result;           /* Pointer value to return */
1499 
1500 /* Initialise. */
1501    result = NULL;
1502 
1503 /* Check the global error status. */
1504    if ( !astOK ) return result;
1505 
1506 /* Obtain a pointer to the SpecFrame structure. */
1507    this = (AstSpecFrame *) this_frame;
1508 
1509 /* If a Domain attribute string has been set, invoke the parent method
1510    to obtain a pointer to it. */
1511    if ( astTestDomain( this ) ) {
1512       result = (*parent_getdomain)( this_frame, status );
1513 
1514 /* Otherwise, provide a pointer to a suitable default string. */
1515    } else {
1516       result = "SPECTRUM";
1517    }
1518 
1519 /* Return the result. */
1520    return result;
1521 }
1522 
GetLabel(AstFrame * this,int axis,int * status)1523 static const char *GetLabel( AstFrame *this, int axis, int *status ) {
1524 /*
1525 *  Name:
1526 *     GetLabel
1527 
1528 *  Purpose:
1529 *     Access the Label string for a SpecFrame axis.
1530 
1531 *  Type:
1532 *     Private function.
1533 
1534 *  Synopsis:
1535 *     #include "specframe.h"
1536 *     const char *GetLabel( AstFrame *this, int axis, int *status )
1537 
1538 *  Class Membership:
1539 *     SpecFrame member function (over-rides the astGetLabel method inherited
1540 *     from the Frame class).
1541 
1542 *  Description:
1543 *     This function returns a pointer to the Label string for a specified axis
1544 *     of a SpecFrame.
1545 
1546 *  Parameters:
1547 *     this
1548 *        Pointer to the SpecFrame.
1549 *     axis
1550 *        Axis index (zero-based) identifying the axis for which information is
1551 *        required.
1552 *     status
1553 *        Pointer to the inherited status variable.
1554 
1555 *  Returned Value:
1556 *     Pointer to a constant null-terminated character string containing the
1557 *     requested information.
1558 
1559 *  Notes:
1560 *     -  A NULL pointer will be returned if this function is invoked with the
1561 *     global error status set, or if it should fail for any reason.
1562 */
1563 
1564 /* Local Variables: */
1565    astDECLARE_GLOBALS            /* Declare the thread specific global data */
1566    AstMapping *map;              /* Mapping between units */
1567    AstSystemType system;         /* Code identifying type of spectral coordinates */
1568    char *new_lab;                /* Modified label string */
1569    const char *result;           /* Pointer to label string */
1570    double orig;                  /* Spec origin */
1571 
1572 /* Check the global error status. */
1573    if ( !astOK ) return NULL;
1574 
1575 /* Get a pointer to the structure holding thread-specific global data. */
1576    astGET_GLOBALS(this);
1577 
1578 /* Initialise. */
1579    result = NULL;
1580 
1581 /* Validate the axis index. */
1582    astValidateAxis( this, axis, 1, "astGetLabel" );
1583 
1584 /* Check if a value has been set for the required axis label string. If so,
1585    invoke the parent astGetLabel method to obtain a pointer to it. */
1586    if ( astTestLabel( this, axis ) ) {
1587       result = (*parent_getlabel)( this, axis, status );
1588 
1589 /* Otherwise, identify the spectral coordinate system described by the
1590    SpecFrame. */
1591    } else {
1592       system = astGetSystem( this );
1593 
1594 /* If OK, supply a pointer to a suitable default label string. */
1595       if ( astOK ) {
1596          result = strcpy( getlabel_buff, SystemLabel( system, status ) );
1597          getlabel_buff[ 0 ] = toupper( getlabel_buff[ 0 ] );
1598 
1599 /* If a non-zero SpecOrigin has been specified, include the offset now. */
1600          orig = GetSpecOriginCur( (AstSpecFrame *) this, status );
1601          if( orig != 0.0 ) {
1602             sprintf( getlabel_buff + strlen( getlabel_buff ), " offset from %s",
1603                      astFormat( this, 0, orig ) );
1604          }
1605 
1606 /* Modify this default to take account of the current value of the Unit
1607    attribute, if set. */
1608          if( astTestUnit( this, axis ) ) {
1609 
1610 /* Find a Mapping from the default Units for the current System, to the
1611    units indicated by the Unit attribute. This Mapping is used to modify
1612    the existing default label appropriately. For instance, if the default
1613    units is "Hz" and the actual units is "log(Hz)", then the default label
1614    of "Frequency" is changed to "log( frequency )". */
1615             map = astUnitMapper( DefUnit( system, "astGetLabel",
1616                                           astGetClass( this ), status ),
1617                                  astGetUnit( this, axis ), result,
1618                                  &new_lab );
1619             if( new_lab ) {
1620                result = strcpy( getlabel_buff, new_lab );
1621                new_lab = astFree( new_lab );
1622             }
1623 
1624 /* Annul the unused Mapping. */
1625             if( map ) map = astAnnul( map );
1626 
1627          }
1628       }
1629    }
1630 
1631 /* Return the result. */
1632    return result;
1633 }
1634 
GetRefPos(AstSpecFrame * this,AstSkyFrame * frm,double * lon,double * lat,int * status)1635 static void GetRefPos( AstSpecFrame *this, AstSkyFrame *frm, double *lon,
1636                        double *lat, int *status ){
1637 /*
1638 *++
1639 *  Name:
1640 c     astGetRefPos
1641 f     AST_GETREFPOS
1642 
1643 *  Purpose:
1644 *     Return the reference position in a specified celestial coordinate system.
1645 
1646 *  Type:
1647 *     Public virtual function.
1648 
1649 *  Synopsis:
1650 c     #include "specframe.h"
1651 c     void astGetRefPos( AstSpecFrame *this, AstSkyFrame *frm, double *lon,
1652 c                        double *lat )
1653 f     CALL AST_GETREFPOS( THIS, FRM, LON, LAT, STATUS )
1654 
1655 *  Class Membership:
1656 *     Frame method.
1657 
1658 *  Description:
1659 c     This function
1660 f     This routine
1661 *     returns the reference position (specified by attributes RefRA and
1662 *     RefDec) converted to the celestial coordinate system represented by
1663 *     a supplied SkyFrame. The celestial longitude and latitude values
1664 *     are returned in radians.
1665 
1666 *  Parameters:
1667 c     this
1668 f     THIS = INTEGER (Given)
1669 *        Pointer to the SpecFrame.
1670 c     frm
1671 f     FRM = INTEGER (Given)
1672 *        Pointer to the SkyFrame which defines the required celestial
1673 *        coordinate system.
1674 c        If NULL
1675 f        If AST__NULL
1676 *        is supplied, then the longitude and latitude values are returned
1677 *        as FK5 J2000 RA and Dec values.
1678 c     lon
1679 f     LON = DOUBLE PRECISION (Returned)
1680 c        A pointer to a double in which to store the
1681 f        The
1682 *        longitude of the reference point, in the coordinate system
1683 *        represented by the supplied SkyFrame (radians).
1684 c     lat
1685 f     LAT = DOUBLE PRECISION (Returned)
1686 c        A pointer to a double in which to store the
1687 f        The
1688 *        latitude of the reference point, in the coordinate system
1689 *        represented by the supplied SkyFrame (radians).
1690 f     STATUS = INTEGER (Given and Returned)
1691 f        The global status.
1692 
1693 *  Notes:
1694 *     - Values of AST__BAD will be returned if this function is
1695 c     invoked with the AST error status set, or if it should fail for
1696 f     invoked with STATUS set to an error value, or if it should fail for
1697 *     any reason.
1698 *--
1699 */
1700 
1701 /* Local Variables: */
1702    AstFrameSet *fs;            /* Conversion FrameSet */
1703    AstFrame *fb;               /* Base Frame */
1704    AstFrame *fc;               /* Current Frame */
1705    double xin[ 1 ];            /* Axis 1 values */
1706    double yin[ 1 ];            /* Axis 2 values */
1707    double xout[ 1 ];           /* Axis 1 values */
1708    double yout[ 1 ];           /* Axis 2 values */
1709 
1710 /* Initialise. */
1711    if( lon ) *lon = AST__BAD;
1712    if( lat ) *lat = AST__BAD;
1713 
1714 /* Check the global error status. */
1715    if ( !astOK ) return;
1716 
1717 /* If no SkyFrame was supplied, just return the stored RefRA and RefDec
1718    values. */
1719    if( !frm ) {
1720       if( lon ) *lon = astGetRefRA( this );
1721       if( lat ) *lat = astGetRefDec( this );
1722 
1723 /* Otherwise, convert the stored values to the requested system. */
1724    } else {
1725 
1726 /* Create an FK5 J2000 SkyFrame which will be used for formatting and
1727    unformatting sky positions, etc. */
1728       LOCK_MUTEX2
1729       if( !skyframe ) {
1730          astBeginPM;
1731          skyframe = astSkyFrame( "system=FK5,equinox=J2000", status );
1732          astEndPM;
1733       }
1734       UNLOCK_MUTEX2
1735 
1736 /* Find the Mapping from the SkyFrame which describes the internal format
1737    in which the RefRA and RefDec attribute values are stored, to the
1738    supplied Frame. */
1739       fs = astFindFrame( skyframe, frm, "" );
1740 
1741 /* If alignment was possible, use the Mapping to transform the internal
1742    RefRA and RefDec values. Check for axis permutatuion. */
1743       if( fs ) {
1744          fb = astGetFrame( fs, AST__BASE );
1745          if( astGetLonAxis( fb ) == 0 ) {
1746             xin[ 0 ] = astGetRefRA( this );
1747             yin[ 0 ] = astGetRefDec( this );
1748          } else {
1749             yin[ 0 ] = astGetRefRA( this );
1750             xin[ 0 ] = astGetRefDec( this );
1751          }
1752          astTran2( fs, 1, xin, yin, 1, xout, yout );
1753 
1754 /* Store the returned values, checking to see if the axes of the supplied
1755    SkyFrame have been permuted. */
1756          fc = astGetFrame( fs, AST__CURRENT );
1757          if( astGetLonAxis( fc ) == 0 ) {
1758             if( lon ) *lon = xout[ 0 ];
1759             if( lat ) *lat = yout[ 0 ];
1760          } else {
1761             if( lon ) *lon = yout[ 0 ];
1762             if( lat ) *lat = xout[ 0 ];
1763          }
1764 
1765 /* Annul object references. */
1766          fc = astAnnul( fc );
1767          fb = astAnnul( fb );
1768          fs = astAnnul( fs );
1769       }
1770    }
1771 }
1772 
GetSymbol(AstFrame * this,int axis,int * status)1773 static const char *GetSymbol( AstFrame *this, int axis, int *status ) {
1774 /*
1775 *  Name:
1776 *     GetSymbol
1777 
1778 *  Purpose:
1779 *     Obtain a pointer to the Symbol string for a SpecFrame axis.
1780 
1781 *  Type:
1782 *     Private function.
1783 
1784 *  Synopsis:
1785 *     #include "specframe.h"
1786 *     const char *GetSymbol( AstFrame *this, int axis, int *status )
1787 
1788 *  Class Membership:
1789 *     SpecFrame member function (over-rides the astGetSymbol method inherited
1790 *     from the Frame class).
1791 
1792 *  Description:
1793 *     This function returns a pointer to the Symbol string for a specified axis
1794 *     of a SpecFrame.
1795 
1796 *  Parameters:
1797 *     this
1798 *        Pointer to the SpecFrame.
1799 *     axis
1800 *        Axis index (zero-based) identifying the axis for which information is
1801 *        required.
1802 *     status
1803 *        Pointer to the inherited status variable.
1804 
1805 *  Returned Value:
1806 *     Pointer to a constant null-terminated character string containing the
1807 *     requested information.
1808 
1809 *  Notes:
1810 *     -  A NULL pointer will be returned if this function is invoked with the
1811 *     global error status set, or if it should fail for any reason.
1812 */
1813 
1814 /* Local Variables: */
1815    astDECLARE_GLOBALS            /* Declare the thread specific global data */
1816    AstMapping *map;              /* Mapping between units */
1817    AstSystemType system;         /* Code identifying type of sky coordinates */
1818    char *new_sym;                /* Modified symbol string */
1819    const char *result;           /* Pointer to symbol string */
1820 
1821 /* Check the global error status. */
1822    if ( !astOK ) return NULL;
1823 
1824 /* Get a pointer to the structure holding thread-specific global data. */
1825    astGET_GLOBALS(this);
1826 
1827 /* Initialise. */
1828    result = NULL;
1829 
1830 /* Validate the axis index. */
1831    astValidateAxis( this, axis, 1, "astGetSymbol" );
1832 
1833 /* Check if a value has been set for the required axis symbol string. If so,
1834    invoke the parent astGetSymbol method to obtain a pointer to it. */
1835    if ( astTestSymbol( this, axis ) ) {
1836       result = (*parent_getsymbol)( this, axis, status );
1837 
1838 /* Otherwise, identify the sky coordinate system described by the SpecFrame. */
1839    } else {
1840       system = astGetSystem( this );
1841 
1842 /* If OK, supply a pointer to a suitable default Symbol string. */
1843       if ( astOK ) {
1844 
1845          if( system == AST__FREQ ) {
1846 	    result = "FREQ";
1847          } else if( system == AST__ENERGY ) {
1848 	    result = "ENER";
1849          } else if( system == AST__WAVENUM ) {
1850 	    result = "WAVN";
1851          } else if( system == AST__WAVELEN ) {
1852 	    result = "WAVE";
1853          } else if( system == AST__AIRWAVE ) {
1854 	    result = "AWAV";
1855          } else if( system == AST__VRADIO ) {
1856 	    result = "VRAD";
1857          } else if( system == AST__VOPTICAL ) {
1858 	    result = "VOPT";
1859          } else if( system == AST__REDSHIFT ) {
1860 	    result = "ZOPT";
1861          } else if( system == AST__BETA ) {
1862 	    result = "BETA";
1863          } else if( system == AST__VREL ) {
1864 	    result = "VELO";
1865 
1866 /* Report an error if the coordinate system was not recognised. */
1867          } else {
1868 	    astError( AST__SCSIN, "astGetSymbol(%s): Corrupt %s contains "
1869 		      "invalid System identification code (%d).", status,
1870                       astGetClass( this ), astGetClass( this ), (int) system );
1871          }
1872 
1873 /* Modify this default to take account of the current value of the Unit
1874    attribute, if set. */
1875          if( astTestUnit( this, axis ) ) {
1876 
1877 /* Find a Mapping from the default Units for the current System, to the
1878    units indicated by the Unit attribute. This Mapping is used to modify
1879    the existing default symbol appropriately. For instance, if the default
1880    units is "Hz" and the actual units is "log(Hz)", then the default symbol
1881    of "nu" is changed to "log( nu )". */
1882             map = astUnitMapper( DefUnit( system, "astGetSymbol",
1883                                           astGetClass( this ), status ),
1884                                  astGetUnit( this, axis ), result,
1885                                  &new_sym );
1886             if( new_sym ) {
1887                result = strcpy( getsymbol_buff, new_sym );
1888                new_sym = astFree( new_sym );
1889             }
1890 
1891 /* Annul the unused Mapping. */
1892             if( map ) map = astAnnul( map );
1893 
1894          }
1895       }
1896    }
1897 
1898 /* Return the result. */
1899    return result;
1900 }
1901 
GetAlignSystem(AstFrame * this_frame,int * status)1902 static AstSystemType GetAlignSystem( AstFrame *this_frame, int *status ) {
1903 /*
1904 *  Name:
1905 *     GetAlignSystem
1906 
1907 *  Purpose:
1908 *     Obtain the AlignSystem attribute for a SpecFrame.
1909 
1910 *  Type:
1911 *     Private function.
1912 
1913 *  Synopsis:
1914 *     #include "Specframe.h"
1915 *     AstSystemType GetAlignSystem( AstFrame *this_frame, int *status )
1916 
1917 *  Class Membership:
1918 *     SpecFrame member function (over-rides the astGetAlignSystem protected
1919 *     method inherited from the Frame class).
1920 
1921 *  Description:
1922 *     This function returns the AlignSystem attribute for a SpecFrame.
1923 
1924 *  Parameters:
1925 *     this
1926 *        Pointer to the SpecFrame.
1927 *     status
1928 *        Pointer to the inherited status variable.
1929 
1930 *  Returned Value:
1931 *     The AlignSystem value.
1932 
1933 */
1934 
1935 /* Local Variables: */
1936    AstSpecFrame *this;           /* Pointer to SpecFrame structure */
1937    AstSystemType result;         /* Value to return */
1938 
1939 /* Initialise. */
1940    result = AST__BADSYSTEM;
1941 
1942 /* Check the global error status. */
1943    if ( !astOK ) return result;
1944 
1945 /* Obtain a pointer to the SpecFrame structure. */
1946    this = (AstSpecFrame *) this_frame;
1947 
1948 /* If a AlignSystem attribute has been set, invoke the parent method to obtain
1949    it. */
1950    if ( astTestAlignSystem( this ) ) {
1951       result = (*parent_getalignsystem)( this_frame, status );
1952 
1953 /* Otherwise, provide a suitable default. */
1954    } else {
1955       result = AST__WAVELEN;
1956    }
1957 
1958 /* Return the result. */
1959    return result;
1960 }
1961 
GetSystem(AstFrame * this_frame,int * status)1962 static AstSystemType GetSystem( AstFrame *this_frame, int *status ) {
1963 /*
1964 *  Name:
1965 *     GetSystem
1966 
1967 *  Purpose:
1968 *     Obtain the System attribute for a SpecFrame.
1969 
1970 *  Type:
1971 *     Private function.
1972 
1973 *  Synopsis:
1974 *     #include "specframe.h"
1975 *     AstSystemType GetSystem( AstFrame *this_frame, int *status )
1976 
1977 *  Class Membership:
1978 *     SpecFrame member function (over-rides the astGetSystem protected
1979 *     method inherited from the Frame class).
1980 
1981 *  Description:
1982 *     This function returns the System attribute for a SpecFrame.
1983 
1984 *  Parameters:
1985 *     this
1986 *        Pointer to the SpecFrame.
1987 *     status
1988 *        Pointer to the inherited status variable.
1989 
1990 *  Returned Value:
1991 *     The System value.
1992 
1993 *  Notes:
1994 *     - AST__BADSYSTEM is returned if this function is invoked with
1995 *     the global error status set or if it should fail for any reason.
1996 */
1997 
1998 /* Local Variables: */
1999    AstSpecFrame *this;           /* Pointer to SpecFrame structure */
2000    AstSystemType result;         /* Value to return */
2001 
2002 /* Initialise. */
2003    result = AST__BADSYSTEM;
2004 
2005 /* Check the global error status. */
2006    if ( !astOK ) return result;
2007 
2008 /* Obtain a pointer to the SpecFrame structure. */
2009    this = (AstSpecFrame *) this_frame;
2010 
2011 /* If a System attribute has been set, invoke the parent method to obtain
2012    it. */
2013    if ( astTestSystem( this ) ) {
2014       result = (*parent_getsystem)( this_frame, status );
2015 
2016 /* Otherwise, provide a suitable default. */
2017    } else {
2018       result = AST__WAVELEN;
2019    }
2020 
2021 /* Return the result. */
2022    return result;
2023 }
2024 
GetTitle(AstFrame * this_frame,int * status)2025 static const char *GetTitle( AstFrame *this_frame, int *status ) {
2026 /*
2027 *  Name:
2028 *     GetTitle
2029 
2030 *  Purpose:
2031 *     Obtain a pointer to the Title string for a SpecFrame.
2032 
2033 *  Type:
2034 *     Private function.
2035 
2036 *  Synopsis:
2037 *     #include "specframe.h"
2038 *     const char *GetTitle( AstFrame *this_frame, int *status )
2039 
2040 *  Class Membership:
2041 *     SpecFrame member function (over-rides the astGetTitle method inherited
2042 *     from the Frame class).
2043 
2044 *  Description:
2045 *     This function returns a pointer to the Title string for a SpecFrame.
2046 *     A pointer to a suitable default string is returned if no Title value has
2047 *     previously been set.
2048 
2049 *  Parameters:
2050 *     this
2051 *        Pointer to the SpecFrame.
2052 *     status
2053 *        Pointer to the inherited status variable.
2054 
2055 *  Returned Value:
2056 *     Pointer to a null-terminated character string containing the requested
2057 *     information.
2058 
2059 *  Notes:
2060 *     -  A NULL pointer will be returned if this function is invoked with the
2061 *     global error status set, or if it should fail for any reason.
2062 */
2063 
2064 /* Local Variables: */
2065    astDECLARE_GLOBALS            /* Declare the thread specific global data */
2066    AstSpecFrame *this;           /* Pointer to SpecFrame structure */
2067    AstStdOfRestType sor;         /* Code identifying standard of rest */
2068    AstSystemType system;         /* Code identifying type of coordinates */
2069    const char *sor_string;       /* Pointer to SOR description */
2070    const char *result;           /* Pointer to result string */
2071    double rf;                    /* Rest frequency */
2072    int nc;                       /* No. of characters added */
2073    int pos;                      /* Buffer position to enter text */
2074 
2075 /* Check the global error status. */
2076    if ( !astOK ) return NULL;
2077 
2078 /* Get a pointer to the structure holding thread-specific global data. */
2079    astGET_GLOBALS(this_frame);
2080 
2081 /* Initialise. */
2082    result = NULL;
2083 
2084 /* Obtain a pointer to the SpecFrame structure. */
2085    this = (AstSpecFrame *) this_frame;
2086 
2087 /* See if a Title string has been set. If so, use the parent astGetTitle
2088    method to obtain a pointer to it. */
2089    if ( astTestTitle( this ) ) {
2090       result = (*parent_gettitle)( this_frame, status );
2091 
2092 /* Otherwise, we will generate a default Title string. Obtain the values of the
2093    SpecFrame's attributes that determine what this string will be. */
2094    } else {
2095       system = astGetSystem( this );
2096       sor = astGetStdOfRest( this );
2097       sor_string = StdOfRestString( sor, status );
2098       rf = astGetRestFreq( this );
2099 
2100 /* Classify the coordinate system type and create an appropriate Title
2101    string.  (Note that when invoking the astFmtDecimalYr function we must
2102    use a separate sprintf on each occasion so as not to over-write its
2103    internal buffer before the result string has been used.) */
2104       if ( astOK ) {
2105          result = gettitle_buff;
2106 
2107 /* Begin with the system's default label. */
2108          pos = sprintf( gettitle_buff, "%s", SystemLabel( system, status ) );
2109          gettitle_buff[ 0 ] = toupper( gettitle_buff[ 0 ] );
2110 
2111 /* Append the standard of rest in parentheses, if set. */
2112          if( astTestStdOfRest( this ) ) {
2113             nc = sprintf( gettitle_buff+pos, " (%s)", sor_string );
2114             pos += nc;
2115          }
2116 
2117 /* Append the rest frequency if relevant. */
2118          if( !ABS_SYSTEM(system) && ( astTestRestFreq( this ) ||
2119                                       astGetUseDefs( this ) ) ) {
2120             pos += sprintf( gettitle_buff+pos, ", rest frequency = %g GHz", rf*1.0E-9 );
2121          }
2122       }
2123    }
2124 
2125 /* If an error occurred, clear the returned pointer value. */
2126    if ( !astOK ) result = NULL;
2127 
2128 /* Return the result. */
2129    return result;
2130 }
2131 
GetSpecOriginCur(AstSpecFrame * this,int * status)2132 static double GetSpecOriginCur( AstSpecFrame *this, int *status ) {
2133 /*
2134 *  Name:
2135 *     GetSpecOriginCur
2136 
2137 *  Purpose:
2138 *     Obtain the SpecOrigin attribute for a SpecFrame in current units.
2139 
2140 *  Type:
2141 *     Private function.
2142 
2143 *  Synopsis:
2144 *     #include "timeframe.h"
2145 *     double GetSpecOriginCur( AstSpecFrame *this, int *status )
2146 
2147 *  Class Membership:
2148 *     SpecFrame virtual function
2149 
2150 *  Description:
2151 *     This function returns the SpecOrigin attribute for a SpecFrame, in
2152 *     the current units of the SpecFrame. The protected astGetSpecOrigin
2153 *     method can be used to obtain the time origin in the default units of
2154 *     the SpecFrame's System.
2155 
2156 *  Parameters:
2157 *     this
2158 *        Pointer to the SpecFrame.
2159 *     status
2160 *        Pointer to the inherited status variable.
2161 
2162 *  Returned Value:
2163 *     The SpecOrigin value, in the units, system and rest frame specified
2164 *     by the current values of the Unit, System and StdOfRest attributes
2165 *     within "this".
2166 
2167 *  Notes:
2168 *     - AST__BAD is returned if this function is invoked with
2169 *     the global error status set or if it should fail for any reason.
2170 */
2171 
2172 /* Local Variables: */
2173    AstMapping *map;
2174    const char *cur;
2175    const char *def;
2176    double result;
2177    double defval;
2178 
2179 /* Initialise. */
2180    result = AST__BAD;
2181 
2182 /* Check the global error status. */
2183    if ( !astOK ) return result;
2184 
2185 /* Get the value in the default units */
2186    result = astGetSpecOrigin( this );
2187 
2188 /* If SpecOrigin is non-zero and non-BAD we convert it to the current units.*/
2189    if( result != 0.0 && result != AST__BAD ) {
2190 
2191 /* Get the default units for the SpecFrame's System. */
2192       def = DefUnit( astGetSystem( this ), "astGetSpecOrigin", "SpecFrame", status );
2193 
2194 /* Get the current units from the SpecFrame. */
2195       cur = astGetUnit( this, 0 );
2196 
2197 /* If the units differ, get a Mapping from default to current units. */
2198       if( cur && def ){
2199          if( strcmp( cur, def ) ) {
2200             map = astUnitMapper( def, cur, NULL, NULL );
2201 
2202 /* Report an error if the units are incompatible. */
2203             if( !map ) {
2204                astError( AST__BADUN, "%s(%s): The current units (%s) are not suitable "
2205                          "for a SpecFrame.", status, "astGetSpecOrigin", astGetClass( this ),
2206                          cur );
2207 
2208 /* Otherwise, transform the stored origin value.*/
2209             } else {
2210                defval = result;
2211                astTran1( map, 1, &defval, 1, &result );
2212                map = astAnnul( map );
2213             }
2214          }
2215       }
2216    }
2217 
2218 /* Return the result. */
2219    return result;
2220 }
2221 
2222 
GetUnit(AstFrame * this_frame,int axis,int * status)2223 static const char *GetUnit( AstFrame *this_frame, int axis, int *status ) {
2224 /*
2225 *  Name:
2226 *     GetUnit
2227 
2228 *  Purpose:
2229 *     Obtain a pointer to the Unit string for a SpecFrame's axis.
2230 
2231 *  Type:
2232 *     Private function.
2233 
2234 *  Synopsis:
2235 *     #include "specframe.h"
2236 *     const char *GetUnit( AstFrame *this_frame, int axis )
2237 
2238 *  Class Membership:
2239 *     SpecFrame member function (over-rides the astGetUnit method inherited
2240 *     from the Frame class).
2241 
2242 *  Description:
2243 *     This function returns a pointer to the Unit string for a specified axis
2244 *     of a SpecFrame. If the Unit attribute has not been set for the axis, a
2245 *     pointer to a suitable default string is returned instead.
2246 
2247 *  Parameters:
2248 *     this
2249 *        Pointer to the SpecFrame.
2250 *     axis
2251 *        The number of the axis (zero-based) for which information is required.
2252 
2253 *  Returned Value:
2254 *     A pointer to a null-terminated string containing the Unit value.
2255 
2256 *  Notes:
2257 *     -  A NULL pointer will be returned if this function is invoked with the
2258 *     global error status set, or if it should fail for any reason.
2259 */
2260 
2261 /* Local Variables: */
2262    AstSpecFrame *this;           /* Pointer to the SpecFrame structure */
2263    AstSystemType system;         /* The SpecFrame's System value */
2264    const char *result;           /* Pointer value to return */
2265 
2266 /* Check the global error status. */
2267    if ( !astOK ) return NULL;
2268 
2269 /* Obtain a pointer to the SpecFrame structure. */
2270    this = (AstSpecFrame *) this_frame;
2271 
2272 /* Validate the axis index. */
2273    astValidateAxis( this, axis, 1, "astGetUnit" );
2274 
2275 /* If a value has been set for the Unit attribute, use the parent
2276    GetUnit method to return a pointer to the required Unit string. */
2277    if( astTestUnit( this, axis ) ){
2278       result = (*parent_getunit)( this_frame, axis, status );
2279 
2280 /* Otherwise, identify the spectral coordinate system described by the
2281    SpecFrame. */
2282    } else {
2283       system = astGetSystem( this );
2284 
2285 /* Return a string describing the default units. */
2286       result = DefUnit( system, "astGetUnit", astGetClass( this ), status );
2287    }
2288 
2289 /* If an error occurred, clear the returned value. */
2290    if ( !astOK ) result = NULL;
2291 
2292 /* Return the result. */
2293    return result;
2294 }
2295 
astInitSpecFrameVtab_(AstSpecFrameVtab * vtab,const char * name,int * status)2296 void astInitSpecFrameVtab_(  AstSpecFrameVtab *vtab, const char *name, int *status ) {
2297 /*
2298 *+
2299 *  Name:
2300 *     astInitSpecFrameVtab
2301 
2302 *  Purpose:
2303 *     Initialise a virtual function table for a SpecFrame.
2304 
2305 *  Type:
2306 *     Protected function.
2307 
2308 *  Synopsis:
2309 *     #include "specframe.h"
2310 *     void astInitSpecFrameVtab( AstSpecFrameVtab *vtab, const char *name )
2311 
2312 *  Class Membership:
2313 *     SpecFrame vtab initialiser.
2314 
2315 *  Description:
2316 *     This function initialises the component of a virtual function
2317 *     table which is used by the SpecFrame class.
2318 
2319 *  Parameters:
2320 *     vtab
2321 *        Pointer to the virtual function table. The components used by
2322 *        all ancestral classes will be initialised if they have not already
2323 *        been initialised.
2324 *     name
2325 *        Pointer to a constant null-terminated character string which contains
2326 *        the name of the class to which the virtual function table belongs (it
2327 *        is this pointer value that will subsequently be returned by the Object
2328 *        astClass function).
2329 *-
2330 */
2331 
2332 /* Local Variables: */
2333    astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
2334    AstFrameVtab *frame;          /* Pointer to Frame component of Vtab */
2335    AstObjectVtab *object;        /* Pointer to Object component of Vtab */
2336 
2337 /* Check the local error status. */
2338    if ( !astOK ) return;
2339 
2340 /* Get a pointer to the thread specific global data structure. */
2341    astGET_GLOBALS(NULL);
2342 
2343 /* Initialize the component of the virtual function table used by the
2344    parent class. */
2345    astInitFrameVtab( (AstFrameVtab *) vtab, name );
2346 
2347 /* Store a unique "magic" value in the virtual function table. This
2348    will be used (by astIsASpecFrame) to determine if an object belongs
2349    to this class.  We can conveniently use the address of the (static)
2350    class_check variable to generate this unique value. */
2351    vtab->id.check = &class_check;
2352    vtab->id.parent = &(((AstFrameVtab *) vtab)->id);
2353 
2354 /* Initialise member function pointers. */
2355 /* ------------------------------------ */
2356 /* Store pointers to the member functions (implemented here) that
2357    provide virtual methods for this class. */
2358    vtab->GetRefPos = GetRefPos;
2359    vtab->SetRefPos = SetRefPos;
2360 
2361    vtab->ClearAlignStdOfRest = ClearAlignStdOfRest;
2362    vtab->TestAlignStdOfRest = TestAlignStdOfRest;
2363    vtab->GetAlignStdOfRest = GetAlignStdOfRest;
2364    vtab->SetAlignStdOfRest = SetAlignStdOfRest;
2365 
2366    vtab->ClearSourceVRF = ClearSourceVRF;
2367    vtab->TestSourceVRF = TestSourceVRF;
2368    vtab->GetSourceVRF = GetSourceVRF;
2369    vtab->SetSourceVRF = SetSourceVRF;
2370 
2371    vtab->ClearSourceSys = ClearSourceSys;
2372    vtab->TestSourceSys = TestSourceSys;
2373    vtab->GetSourceSys = GetSourceSys;
2374    vtab->SetSourceSys = SetSourceSys;
2375 
2376    vtab->ClearRefDec = ClearRefDec;
2377    vtab->TestRefDec = TestRefDec;
2378    vtab->GetRefDec = GetRefDec;
2379    vtab->SetRefDec = SetRefDec;
2380 
2381    vtab->ClearRefRA = ClearRefRA;
2382    vtab->TestRefRA = TestRefRA;
2383    vtab->GetRefRA = GetRefRA;
2384    vtab->SetRefRA = SetRefRA;
2385 
2386    vtab->ClearRestFreq = ClearRestFreq;
2387    vtab->TestRestFreq = TestRestFreq;
2388    vtab->GetRestFreq = GetRestFreq;
2389    vtab->SetRestFreq = SetRestFreq;
2390 
2391    vtab->ClearStdOfRest = ClearStdOfRest;
2392    vtab->TestStdOfRest = TestStdOfRest;
2393    vtab->GetStdOfRest = GetStdOfRest;
2394    vtab->SetStdOfRest = SetStdOfRest;
2395 
2396    vtab->ClearSourceVel = ClearSourceVel;
2397    vtab->TestSourceVel = TestSourceVel;
2398    vtab->GetSourceVel = GetSourceVel;
2399    vtab->SetSourceVel = SetSourceVel;
2400 
2401    vtab->ClearSpecOrigin = ClearSpecOrigin;
2402    vtab->TestSpecOrigin = TestSpecOrigin;
2403    vtab->GetSpecOrigin = GetSpecOrigin;
2404    vtab->SetSpecOrigin = SetSpecOrigin;
2405 
2406    vtab->TestAlignSpecOffset = TestAlignSpecOffset;
2407    vtab->SetAlignSpecOffset = SetAlignSpecOffset;
2408    vtab->GetAlignSpecOffset = GetAlignSpecOffset;
2409    vtab->ClearAlignSpecOffset = ClearAlignSpecOffset;
2410 
2411 /* Save the inherited pointers to methods that will be extended, and
2412    replace them with pointers to the new member functions. */
2413    object = (AstObjectVtab *) vtab;
2414    frame = (AstFrameVtab *) vtab;
2415    parent_getobjsize = object->GetObjSize;
2416    object->GetObjSize = GetObjSize;
2417 
2418    parent_clearattrib = object->ClearAttrib;
2419    object->ClearAttrib = ClearAttrib;
2420    parent_getattrib = object->GetAttrib;
2421    object->GetAttrib = GetAttrib;
2422    parent_setattrib = object->SetAttrib;
2423    object->SetAttrib = SetAttrib;
2424    parent_testattrib = object->TestAttrib;
2425    object->TestAttrib = TestAttrib;
2426 
2427    parent_getdomain = frame->GetDomain;
2428    frame->GetDomain = GetDomain;
2429 
2430    parent_getsystem = frame->GetSystem;
2431    frame->GetSystem = GetSystem;
2432    parent_setsystem = frame->SetSystem;
2433    frame->SetSystem = SetSystem;
2434    parent_clearsystem = frame->ClearSystem;
2435    frame->ClearSystem = ClearSystem;
2436 
2437    parent_getalignsystem = frame->GetAlignSystem;
2438    frame->GetAlignSystem = GetAlignSystem;
2439 
2440    parent_getlabel = frame->GetLabel;
2441    frame->GetLabel = GetLabel;
2442 
2443    parent_getsymbol = frame->GetSymbol;
2444    frame->GetSymbol = GetSymbol;
2445 
2446    parent_gettitle = frame->GetTitle;
2447    frame->GetTitle = GetTitle;
2448 
2449    parent_clearunit = frame->ClearUnit;
2450    frame->ClearUnit = ClearUnit;
2451 
2452    parent_getunit = frame->GetUnit;
2453    frame->GetUnit = GetUnit;
2454 
2455    parent_setunit = frame->SetUnit;
2456    frame->SetUnit = SetUnit;
2457 
2458    parent_match = frame->Match;
2459    frame->Match = Match;
2460 
2461    parent_overlay = frame->Overlay;
2462    frame->Overlay = Overlay;
2463 
2464    parent_subframe = frame->SubFrame;
2465    frame->SubFrame = SubFrame;
2466 
2467 /* Store replacement pointers for methods which will be over-ridden by new
2468    member functions implemented here. */
2469    frame->GetActiveUnit = GetActiveUnit;
2470    frame->TestActiveUnit = TestActiveUnit;
2471    frame->ValidateSystem = ValidateSystem;
2472    frame->SystemString = SystemString;
2473    frame->SystemCode = SystemCode;
2474 
2475 /* Declare the copy constructor, destructor and class dump
2476    function. */
2477    astSetCopy( vtab, Copy );
2478    astSetDelete( vtab, Delete );
2479    astSetDump( vtab, Dump, "SpecFrame",
2480                "Description of spectral coordinate system" );
2481 
2482 /* If we have just initialised the vtab for the current class, indicate
2483    that the vtab is now initialised, and store a pointer to the class
2484    identifier in the base "object" level of the vtab. */
2485    if( vtab == &class_vtab ) {
2486       class_init = 1;
2487       astSetVtabClassIdentifier( vtab, &(vtab->id) );
2488    }
2489 }
2490 
MakeSpecMapping(AstSpecFrame * target,AstSpecFrame * result,AstSpecFrame * align_frm,int report,AstMapping ** map,int * status)2491 static int MakeSpecMapping( AstSpecFrame *target, AstSpecFrame *result,
2492                             AstSpecFrame *align_frm, int report,
2493                             AstMapping **map, int *status ) {
2494 /*
2495 *  Name:
2496 *     MakeSpecMapping
2497 
2498 *  Purpose:
2499 *     Generate a Mapping between two SpecFrames.
2500 
2501 *  Type:
2502 *     Private function.
2503 
2504 *  Synopsis:
2505 *     #include "specframe.h"
2506 *     int MakeSpecMapping( AstSpecFrame *target, AstSpecFrame *result,
2507 *                          AstSpecFrame *align_frm, int report,
2508 *                          AstMapping **map, int *status ) {
2509 
2510 *  Class Membership:
2511 *     SpecFrame member function.
2512 
2513 *  Description:
2514 *     This function takes two SpecFrames and generates a Mapping that
2515 *     converts between them, taking account of differences in their
2516 *     coordinate systems, rest frequency, standard of rest, etc.
2517 *
2518 *     In order to cut down the number of transformations to be considered,
2519 *     the scheme works by first converting from the target frame to an
2520 *     "alignment" Frame, using the attributes of the target to define the
2521 *     transformation. A transformation is then found from the alignment
2522 *     frame to the required result Frame,  using the attributes of the
2523 *     result to define the transformation. The alignment Frame is
2524 *     described by the AlignSystem and AlignStdOfRest attributes of the
2525 *     "align_frm" SpecFrame.
2526 *
2527 *     Thus, different forms of alignment can be obtained by suitable
2528 *     choice of the attributes of "align_frm". For instance, to compare the
2529 *     radio velocity dispersion of two lines at different rest frequencies,
2530 *     you would set "system=radio velocity" and (probably) "stdofrest=local
2531 *     group" in "align_frm". On the other hand if you wanted to re-calibrate
2532 *     an existing radio velocity Frame within a FrameSet to use a different
2533 *     rest frequency, you would make the SpecFrame the current Frame and then
2534 *     set the rest frequency attribute for the FrameSet. The "integrity
2535 *     checking" system in the FrameSet class would then get the Mapping
2536 *     between the original and the modified SpecFrames. In this case, the
2537 *     "alignment system" needs to be "frequency" since you want the original
2538 *     and modified SpecFrames to be aligned in frequency, not radio velocity.
2539 
2540 *  Parameters:
2541 *     target
2542 *        Pointer to the first SpecFrame.
2543 *     result
2544 *        Pointer to the second SpecFrame.
2545 *     align_frm
2546 *        A SpecFrame defining the system and standard of rest in which to
2547 *        align the target and result SpecFrames.
2548 *     report
2549 *        Should errors be reported if no match is possible? These reports
2550 *        will describe why no match was possible.
2551 *     map
2552 *        Pointer to a location which is to receive a pointer to the
2553 *        returned Mapping. The forward transformation of this Mapping
2554 *        will convert from "target" coordinates to "result"
2555 *        coordinates, and the inverse transformation will convert in
2556 *        the opposite direction (all coordinate values in radians).
2557 *     status
2558 *        Pointer to the inherited status variable.
2559 
2560 *  Returned Value:
2561 *     Non-zero if the Mapping could be generated, or zero if the two
2562 *     SpecFrames are sufficiently un-related that no meaningful Mapping
2563 *     can be produced (albeit an "unmeaningful" Mapping will be returned
2564 *     in this case, which will need to be annulled).
2565 
2566 *  Notes:
2567 *     A value of zero is returned if this function is invoked with the
2568 *     global error status set or if it should fail for any reason.
2569 */
2570 
2571 /* Local Constants: */
2572 #define MAX_ARGS 1               /* Max arguments for an SpecMap conversion */
2573 
2574 /* Local Variables: */
2575    AstMapping *map1;             /* Intermediate Mapping */
2576    AstMapping *map2;             /* Intermediate Mapping */
2577    AstMapping *umap1;            /* First Units Mapping */
2578    AstMapping *umap2;            /* Second Units Mapping */
2579    AstSpecMap *specmap;          /* Pointer to SpecMap */
2580    AstShiftMap *sm;              /* ShiftMap pointer */
2581    AstSpecFrame *align_target;   /* Alignment Frame with target properties */
2582    AstSpecFrame *align_result;   /* Alignment Frame with result properties */
2583    AstSystemType serr;           /* Erroneous system */
2584    AstSystemType align_system;   /* Code to identify alignment system */
2585    AstSystemType target_system;  /* Code to identify target system */
2586    AstSystemType result_system;  /* Code to identify result system */
2587    const char *uerr;             /* Erroneous units */
2588    const char *ures;             /* Results units */
2589    const char *utarg;            /* Target units */
2590    const char *vmess;            /* Text for use in error messages */
2591    double args[ MAX_ARGS ];      /* Conversion argument array */
2592    double target_rf;             /* Target rest frequency (Hz) */
2593    double result_rf;             /* Result rest frequency (Hz) */
2594    double target_origin;         /* Target origin */
2595    double result_origin;         /* Result origin */
2596    int match;                    /* Mapping can be generated? */
2597    int step2;                    /* Perform the 2nd step in the Mapping? */
2598    int step3;                    /* Perform the 3rd step in the Mapping? */
2599    int step4;                    /* Perform the 4th step in the Mapping? */
2600    int step5;                    /* Perform the 5th step in the Mapping? */
2601    int step6;                    /* Perform the 6th step in the Mapping? */
2602    int step7;                    /* Perform the 7th step in the Mapping? */
2603 
2604 /* Check the global error status. */
2605    if ( !astOK ) return 0;
2606 
2607 /* Initialise the returned values. */
2608    match = 1;
2609    *map = NULL;
2610 
2611 /* Create an initial (null) SpecMap. This is a 1D Mapping which converts
2612    spectral axis values between different systems and standard of rest.
2613    The axis units used by the SpecMap class match the default units used
2614    by this class. Any discrepancy between units is taken into account at
2615    the end of this function, once the total SpecMap has been created. */
2616    specmap = astSpecMap( 1, 0, "", status );
2617 
2618 /* Define local macros as shorthand for adding spectral coordinate
2619    conversions to this SpecMap.  Each macro simply stores details of
2620    the additional arguments in the "args" array and then calls
2621    astSpecAdd. The macros differ in the number of additional argument
2622    values. */
2623 #define TRANSFORM_0(cvt) \
2624         astSpecAdd( specmap, cvt, NULL );
2625 
2626 #define TRANSFORM_1(cvt,arg0) \
2627         args[ 0 ] = arg0; \
2628         astSpecAdd( specmap, cvt, args );
2629 
2630 /* Get all the necessary attributes from the result, target and alignment
2631    Frames. */
2632    target_rf = astGetRestFreq( target );
2633    result_rf = astGetRestFreq( result );
2634 
2635    target_system = astGetSystem( target );
2636    result_system = astGetSystem( result );
2637    align_system = astGetSystem( align_frm );
2638 
2639 /* Define text for error messages.*/
2640    vmess = "convert between spectral systems";
2641 
2642 /* Verify that values for the standard of rest have been set if required
2643    (i.e if the UseDefs attribute of either SpecFrame is false). */
2644    VerifyAttrs( result, vmess, "StdOfRest", "astMatch", status );
2645    VerifyAttrs( target, vmess, "StdOfRest", "astMatch", status );
2646 
2647 /* There are two different strategies for alignment. I'll use the Source
2648    rest frame as an example, although the same argument applies to other
2649    rest frames. In the first strategy, all "Source" rest frames are
2650    considered equal. That is, if two SpecFrames respresent (for example)
2651    frequencies in the source frame, then the SpecFrames are aligned using
2652    a UnitMap even if the details of the two source rest frames differ.
2653    This is usually what users want to see when (for instance) aligning
2654    plots of two spectra which both represent source frequencies but where
2655    the source frames details differ. In the second strategy, "Source"
2656    rest frames are aligned using a SpecMap that takes into account any
2657    differences in the properties of the source rest frames. This is what
2658    should happen when changes are made to the properties of a SpecFrame
2659    within a FrameSet. For instance, if the user changes the SourceVel
2660    attribute of the current Frame (assumed here to be a SpecFrame) in a
2661    FrameSet, then the process of restoring the integrity of the FrameSet
2662    (see frameset.c for details of integrity restoration) should cause the
2663    base->current Mapping in the FrameSet to be modified to reflect the
2664    new SourceVel value.
2665 
2666    So if the current call to this function is part of the process of
2667    restoring a FrameSet's integrity following changes to the FrameSet's
2668    current Frame, then we want to retain the properties of the supplied
2669    alignment Frame. So we use clones of the supplied alignment Frame. */
2670    if( astGetFrameFlags( target ) & AST__INTFLAG ) {
2671       align_target = astClone( align_frm );
2672       align_result = astClone( align_frm );
2673 
2674 /* Buf if we are not restoring the integrity of a FrameSet, we want
2675    to ignore any differences in the properties that define the available
2676    rest frames. So create copies of the alignment Frame in which the
2677    properies defining the available rest frames are the same as in the
2678    target and result Frames. */
2679    } else {
2680       align_target = astCopy( align_frm );
2681       astSetEpoch( align_target, astGetEpoch( target ) );
2682       astSetRefRA( align_target, astGetRefRA( target ) );
2683       astSetRefDec( align_target, astGetRefDec( target ) );
2684       astSetSourceVRF( align_target, astGetSourceVRF( target ) );
2685       astSetSourceVel( align_target, astGetSourceVel( target ) );
2686       astSetObsLat( align_target, astGetObsLat( target ) );
2687       astSetObsLon( align_target, astGetObsLon( target ) );
2688       astSetObsAlt( align_target, astGetObsAlt( target ) );
2689 
2690       align_result = astCopy( align_frm );
2691       astSetEpoch( align_result, astGetEpoch( result ) );
2692       astSetRefRA( align_result, astGetRefRA( result ) );
2693       astSetRefDec( align_result, astGetRefDec( result ) );
2694       astSetSourceVRF( align_result, astGetSourceVRF( result ) );
2695       astSetSourceVel( align_result, astGetSourceVel( result ) );
2696       astSetObsLat( align_result, astGetObsLat( result ) );
2697       astSetObsLon( align_result, astGetObsLon( result ) );
2698       astSetObsAlt( align_result, astGetObsAlt( result ) );
2699    }
2700 
2701 /* The supported spectral coordinate systems fall into two groups;
2702    "relative", and "absolute". The relative systems define each axis
2703    value with respect to the rest frequency, whereas the absolute systems
2704    have axis values which do not depend on the rest frequency. In order
2705    to convert an axis value from a system in one group to a system in the
2706    other group, the rest frequency must be known. However, the rest
2707    frequency is not necessary in order to convert axis values between two
2708    systems belonging to the same group.  Determine if the alignment system
2709    is absolute or relative. If absolute, we ignore the system of the supplied
2710    "align_frm" and align in frequency, since aligning in any absolute system
2711    will automatically ensure that all the other absolute systems are aligned.
2712    Similarly, aligning in any relative system will automatically ensure that
2713    all the other relative systems are aligned. Doing this cuts down the
2714    complexity of the conversion process since we do not need to check every
2715    possible alignment system. */
2716    align_system = ( ABS_SYSTEM( align_system ) ) ? AST__FREQ : AST__VREL;
2717 
2718 /* The total Mapping is made up of the following steps in series:
2719 
2720   0) Convert from an offset value to an absolute value (if SpecOrigin set)
2721   1) Convert target units to default units for the targets system
2722   2) Convert from target system in target SOR to frequency in target SOR
2723   3) Convert from freq in target SOR to freq in alignment SOR
2724   4) Convert from freq in alignment SOR to alignment system in alignment SOR
2725   5) Convert from alignment system in alignment SOR to freq in alignment SOR
2726   6) Convert from freq in alignment SOR to freq in result SOR
2727   7) Convert from freq in result SOR to result system in result SOR
2728   8) Convert default units for the result system to results unit
2729   9) Convert from an absolute value to an offset value (if SpecOrigin set)
2730 
2731    Steps 1,2,3,4 are performed using the attributes of the target (rest
2732    frequency, reference farem, etc), whilst steps 5,6,7,8 are performed
2733    using the attributes of the target (rest frequency, reference frame,
2734    etc). It is necessary to go from target system to alignment system
2735    via frequency because SOR conversion can only be performed in the
2736    frequency domain.
2737 
2738    Some of these steps may not be necessary. Initially assume all steps
2739    are necessary (we leave steps 0, 1, 8 and 9 out of this process and
2740    implement them once all other steps have been done). */
2741    step2 = 1;
2742    step3 = 1;
2743    step4 = 1;
2744    step5 = 1;
2745    step6 = 1;
2746    step7 = 1;
2747 
2748 /* Step 2 is not necessary if the target system is frequency. */
2749    if( target_system == AST__FREQ ) step2 = 0;
2750 
2751 /* Step 3 is not necessary if the alignment SOR is the same as the target
2752    SOR. */
2753    if( EqualSor( target, align_target, status ) ) step3 = 0;
2754 
2755 /* Step 6 is not necessary if the alignment SOR is the same as the result
2756    SOR. */
2757    if( EqualSor( result, align_result, status ) ) step6 = 0;
2758 
2759 /* Step 7 is not necessary if the result system is frequency. */
2760    if( result_system == AST__FREQ ) step7 = 0;
2761 
2762 /* Steps 4 and 5 are not necessary if the alignment system is frequency,
2763    or if the target and result rest frequencies are equal. */
2764    if( align_system == AST__FREQ || result_rf == target_rf ) step4 = step5 = 0;
2765 
2766 /* Steps 3 and 6 are not necessary if steps 4 and 5 are not necessary, and
2767    the target sor equals the result sor. */
2768    if( !step4 && !step5 && EqualSor( target, result, status ) ) step3 = step6 = 0;
2769 
2770 /* Steps 2 and 7 are not necessary if steps 3, 4, 5 and 6 are not necessary,
2771    and the target sor equals the result sor, and the target and results
2772    systems are equal (if the systems are relative they must also have equal
2773    rest frequencies). */
2774    if( !step3 && !step4 && !step5 && !step6 && EqualSor( target, result, status ) &&
2775        target_system == result_system ) {
2776       if( !ABS_SYSTEM( target_system ) && result_rf == target_rf ) step2 = step7 = 0;
2777    }
2778 
2779 
2780 /* Now we know which steps are needed, let's do them (we delay unit
2781    conversion to the end)... */
2782 
2783 /* Step 2: target system in target rest frame to frequency in target rest
2784    frame. */
2785    if( step2 ) {
2786       if( target_system != AST__FREQ ) {
2787 
2788 /* If the target system is absolute, we can convert directly to frequency. */
2789          if ( target_system == AST__ENERGY ) {
2790             TRANSFORM_0( "ENTOFR" )
2791 
2792          } else if ( target_system == AST__WAVENUM ) {
2793             TRANSFORM_0( "WNTOFR" )
2794 
2795          } else if ( target_system == AST__WAVELEN ) {
2796             TRANSFORM_0( "WVTOFR" )
2797 
2798          } else if ( target_system == AST__AIRWAVE ) {
2799             TRANSFORM_0( "AWTOFR" )
2800 
2801 /* If the target target_system is relative, we first need to convert to
2802    apparent radial velocity, and then to frequency using the rest frequency. */
2803          } else {
2804 
2805             if ( target_system == AST__VRADIO ) {
2806                TRANSFORM_0( "VRTOVL" )
2807 
2808             } else if ( target_system == AST__VOPTICAL ) {
2809                TRANSFORM_0( "VOTOVL" )
2810 
2811             } else if ( target_system == AST__REDSHIFT ) {
2812                TRANSFORM_0( "ZOTOVL" )
2813 
2814             } else if ( target_system == AST__BETA ) {
2815                TRANSFORM_0( "BTTOVL" )
2816             }
2817 
2818             VerifyAttrs( target, vmess, "RestFreq", "astMatch", status );
2819             TRANSFORM_1( "VLTOFR", target_rf )
2820          }
2821       }
2822    }
2823 
2824 /* Step 3: frequency in target rest frame to frequency in alignment rest
2825    frame. */
2826    if( step3 ) match = SorConvert( target, align_target, specmap, status );
2827 
2828 /* Step 4: frequency in alignment rest frame to alignment system in alignment
2829    rest frame. The alignment will be either relativistic velocity or
2830    frequency. */
2831    if( step4 ) {
2832       if( align_system == AST__VREL ) {
2833          VerifyAttrs( target, vmess, "RestFreq", "astMatch", status );
2834          TRANSFORM_1( "FRTOVL", target_rf )
2835       }
2836    }
2837 
2838 /* Step 5: Alignment system in alignment rest frame to frequency in alignment
2839    rest frame (from now on use the attributes of the result SpecFrame to
2840    define the conversion parameters). */
2841    if( step5 ) {
2842       if( align_system == AST__VREL ) {
2843          VerifyAttrs( result, vmess, "RestFreq", "astMatch", status );
2844          TRANSFORM_1( "VLTOFR", result_rf )
2845       }
2846    }
2847 
2848 /* Step 6: frequency in alignment rest frame to frequency in result rest
2849    frame. */
2850    if( step6 ) match = SorConvert( align_result, result, specmap, status );
2851 
2852 /* Step 7: frequency in result rest frame to result system in result rest
2853    frame. */
2854    if( step7 ) {
2855       if( result_system != AST__FREQ ) {
2856 
2857 /* If the results system is absolute, we can convert directly. */
2858          if ( result_system == AST__ENERGY ) {
2859             TRANSFORM_0( "FRTOEN" )
2860 
2861          } else if ( result_system == AST__WAVENUM ) {
2862             TRANSFORM_0( "FRTOWN" )
2863 
2864          } else if ( result_system == AST__WAVELEN ) {
2865             TRANSFORM_0( "FRTOWV" )
2866 
2867          } else if ( result_system == AST__AIRWAVE ) {
2868             TRANSFORM_0( "FRTOAW" )
2869 
2870 /* If the result system is relative, we first need to convert to apparent
2871    radial velocity from frequency using the rest frequency. Report an error
2872    if the rest frequency is undefined. */
2873          } else {
2874             VerifyAttrs( result, vmess, "RestFreq", "astMatch", status );
2875             TRANSFORM_1( "FRTOVL", result_rf )
2876 
2877 /* Now convert from apparent radial velocity to the required result system. */
2878             if ( result_system == AST__VRADIO ) {
2879                TRANSFORM_0( "VLTOVR" )
2880 
2881             } else if ( result_system == AST__VOPTICAL ) {
2882                TRANSFORM_0( "VLTOVO" )
2883 
2884             } else if ( result_system == AST__REDSHIFT ) {
2885                TRANSFORM_0( "VLTOZO" )
2886 
2887             } else if ( result_system == AST__BETA ) {
2888                TRANSFORM_0( "VLTOBT" )
2889             }
2890          }
2891       }
2892    }
2893 
2894 /* The SpecMap created above class assumes that the axis values supplied to
2895    its Transform method are in units which correspond to the default units
2896    for its class (the returned values also use these units). However,
2897    the Unit attributes of the supplied Frames may have been set to some
2898    non-default value, and so we need to add Mappings before and after the
2899    SpecMap which convert to and from the default units. Find the Mapping
2900    from the target Frame Units to the default Units for the target's system. */
2901    utarg = astGetUnit( target, 0 );
2902    umap1 = astUnitMapper( utarg, SpecMapUnit( target_system, "MakeSpecMap",
2903                                               "SpecFrame", status ), NULL, NULL );
2904 
2905 /* Find the Mapping from the default Units for the result's system to the
2906    Units of the result Frame. */
2907    ures = astGetUnit( result, 0 );
2908    umap2 = astUnitMapper( SpecMapUnit( result_system, "MakeSpecMap",
2909                                        "SpecFrame", status ), ures, NULL, NULL );
2910 
2911 /* If both units Mappings were created OK, sandwich the SpecMap between
2912    them. */
2913    if( umap1 && umap2 ) {
2914       map1 = (AstMapping *) astCmpMap( umap1, specmap, 1, "", status );
2915       map2 = (AstMapping *) astCmpMap( map1, umap2, 1, "", status );
2916       map1 = astAnnul( map1 );
2917 
2918 /* If the simplified SpecMap is a UnitMap, and the target and result
2919    units are the same, we do not need to know the mapping between units.
2920    Otherwise, report an error and indicate that we cannot convert between
2921    the Frames. */
2922    } else {
2923       map2 = astSimplify( specmap );
2924       if( !astIsAUnitMap( map2 ) || strcmp( ures, utarg ) ) {
2925          match = 0;
2926          if( astOK && report ) {
2927             if( !umap1 ) {
2928                uerr = utarg;
2929                serr = astGetSystem( target );
2930             } else {
2931                uerr = ures;
2932                serr = astGetSystem( result );
2933             }
2934             astError( AST__BADUN, "astMatch(SpecFrame): Inappropriate units (%s) "
2935                       "specified for a %s axis.", status, uerr, SystemLabel( serr, status ) );
2936          }
2937       }
2938    }
2939 
2940 /* Step 0: offset to absolute value in target system. Prepend the Maping created
2941    above with a ShiftMap that does the required shift of origin. */
2942    target_origin = GetSpecOriginCur( target, status );
2943    if( target_origin != 0.0 ) {
2944       sm = astShiftMap( 1, &target_origin, "", status );
2945       map1 = (AstMapping *) astCmpMap( sm, map2, 1, "", status );
2946       sm = astAnnul( sm );
2947    } else {
2948       map1 = astClone( map2 );
2949    }
2950    map2 = astAnnul( map2 );
2951 
2952 /* Step 9: absolute value to offset in result system. If we are aligning in the
2953    offset system, use the transformed target origin as the new zero point.
2954    Otherwise use the origin from the result frame. First get the origin for the
2955    result system. */
2956    if( astGetAlignSpecOffset( target ) && astGetAlignSpecOffset( result ) ) {
2957       result_origin = 0.0;
2958       astTran1( map1, 1, &result_origin, 1, &result_origin );
2959    } else {
2960       result_origin = GetSpecOriginCur( result, status );
2961    }
2962 
2963 /* Now create the ShiftMap and apend it to the end of the Maping. */
2964    if( result_origin != 0.0 ) {
2965       result_origin = -result_origin;
2966       sm = astShiftMap( 1, &result_origin, "", status );
2967       map2 = (AstMapping *) astCmpMap( map1, sm, 1, "", status );
2968       sm = astAnnul( sm );
2969    } else {
2970       map2 = astClone( map1 );
2971    }
2972    map1 = astAnnul( map1 );
2973 
2974 /* Return the simplified Mapping. */
2975    *map = astSimplify( map2 );
2976 
2977 /* Annul remaining resources. */
2978    map2 = astAnnul( map2 );
2979    specmap = astAnnul( specmap );
2980    if( umap1 ) umap1 = astAnnul( umap1 );
2981    if( umap2 ) umap2 = astAnnul( umap2 );
2982    align_result = astAnnul( align_result );
2983    align_target = astAnnul( align_target );
2984 
2985 /* If an error occurred, annul the returned Mapping and clear the returned
2986    values. */
2987    if ( !astOK ) {
2988       *map = astAnnul( *map );
2989       match = 0;
2990    }
2991 
2992 /* Return the result. */
2993    return match;
2994 
2995 /* Undefine macros local to this function. */
2996 #undef MAX_ARGS
2997 #undef TRANSFORM_0
2998 #undef TRANSFORM_1
2999 }
3000 
Match(AstFrame * template_frame,AstFrame * target,int matchsub,int ** template_axes,int ** target_axes,AstMapping ** map,AstFrame ** result,int * status)3001 static int Match( AstFrame *template_frame, AstFrame *target, int matchsub,
3002                   int **template_axes, int **target_axes, AstMapping **map,
3003                   AstFrame **result, int *status ) {
3004 /*
3005 *  Name:
3006 *     Match
3007 
3008 *  Purpose:
3009 *     Determine if conversion is possible between two coordinate systems.
3010 
3011 *  Type:
3012 *     Private function.
3013 
3014 *  Synopsis:
3015 *     #include "specframe.h"
3016 *     int Match( AstFrame *template, AstFrame *target, int matchsub,
3017 *                int **template_axes, int **target_axes,
3018 *                AstMapping **map, AstFrame **result, int *status )
3019 
3020 *  Class Membership:
3021 *     SpecFrame member function (over-rides the protected astMatch method
3022 *     inherited from the Frame class).
3023 
3024 *  Description:
3025 *     This function matches a "template" SpecFrame to a "target" Frame and
3026 *     determines whether it is possible to convert coordinates between them.
3027 *     If it is, a mapping that performs the transformation is returned along
3028 *     with a new Frame that describes the coordinate system that results when
3029 *     this mapping is applied to the "target" coordinate system. In addition,
3030 *     information is returned to allow the axes in this "result" Frame to be
3031 *     associated with the corresponding axes in the "target" and "template"
3032 *     Frames from which they are derived.
3033 
3034 *  Parameters:
3035 *     template
3036 *        Pointer to the template SpecFrame. This describes the coordinate
3037 *        system (or set of possible coordinate systems) into which we wish to
3038 *        convert our coordinates.
3039 *     target
3040 *        Pointer to the target Frame. This describes the coordinate system in
3041 *        which we already have coordinates.
3042 *     matchsub
3043 *        If zero then a match only occurs if the template is of the same
3044 *        class as the target, or of a more specialised class. If non-zero
3045 *        then a match can occur even if this is not the case.
3046 *     template_axes
3047 *        Address of a location where a pointer to int will be returned if the
3048 *        requested coordinate conversion is possible. This pointer will point
3049 *        at a dynamically allocated array of integers with one element for each
3050 *        axis of the "result" Frame (see below). It must be freed by the caller
3051 *        (using astFree) when no longer required.
3052 *
3053 *        For each axis in the result Frame, the corresponding element of this
3054 *        array will return the index of the template SpecFrame axis from
3055 *        which it is derived. If it is not derived from any template
3056 *        SpecFrame axis, a value of -1 will be returned instead.
3057 *     target_axes
3058 *        Address of a location where a pointer to int will be returned if the
3059 *        requested coordinate conversion is possible. This pointer will point
3060 *        at a dynamically allocated array of integers with one element for each
3061 *        axis of the "result" Frame (see below). It must be freed by the caller
3062 *        (using astFree) when no longer required.
3063 *
3064 *        For each axis in the result Frame, the corresponding element of this
3065 *        array will return the index of the target Frame axis from which it
3066 *        is derived. If it is not derived from any target Frame axis, a value
3067 *        of -1 will be returned instead.
3068 *     map
3069 *        Address of a location where a pointer to a new Mapping will be
3070 *        returned if the requested coordinate conversion is possible. If
3071 *        returned, the forward transformation of this Mapping may be used to
3072 *        convert coordinates between the "target" Frame and the "result"
3073 *        Frame (see below) and the inverse transformation will convert in the
3074 *        opposite direction.
3075 *     result
3076 *        Address of a location where a pointer to a new Frame will be returned
3077 *        if the requested coordinate conversion is possible. If returned, this
3078 *        Frame describes the coordinate system that results from applying the
3079 *        returned Mapping (above) to the "target" coordinate system. In
3080 *        general, this Frame will combine attributes from (and will therefore
3081 *        be more specific than) both the target and the template Frames. In
3082 *        particular, when the template allows the possibility of transformaing
3083 *        to any one of a set of alternative coordinate systems, the "result"
3084 *        Frame will indicate which of the alternatives was used.
3085 *     status
3086 *        Pointer to the inherited status variable.
3087 
3088 *  Returned Value:
3089 *     A non-zero value is returned if the requested coordinate conversion is
3090 *     possible. Otherwise zero is returned (this will not in itself result in
3091 *     an error condition).
3092 
3093 *  Notes:
3094 *     -  A value of zero will be returned if this function is invoked with the
3095 *     global error status set, or if it should fail for any reason.
3096 
3097 *  Implementation Notes:
3098 *     This implementation addresses the matching of a SpecFrame class
3099 *     object to any other class of Frame. A SpecFrame will match any class
3100 *     of SpecFrame (i.e. possibly from a derived class) but will not match
3101 *     a less specialised class of Frame.
3102 */
3103 
3104 /* Local Variables: */
3105    AstFrame *frame0;             /* Pointer to Frame underlying axis 0 */
3106    AstSpecFrame *template;       /* Pointer to template SpecFrame structure */
3107    int iaxis0;                   /* Axis index underlying axis 0 */
3108    int iaxis;                    /* Axis index */
3109    int match;                    /* Coordinate conversion possible? */
3110    int target_axis0;             /* Index of SpecFrame axis in the target */
3111    int target_naxes;             /* Number of target axes */
3112 
3113 /* Initialise the returned values. */
3114    *template_axes = NULL;
3115    *target_axes = NULL;
3116    *map = NULL;
3117    *result = NULL;
3118    match = 0;
3119 
3120 /* Check the global error status. */
3121    if ( !astOK ) return match;
3122 
3123 /* Obtain a pointer to the template SpecFrame structure. */
3124    template = (AstSpecFrame *) template_frame;
3125 
3126 /* Obtain the number of axes in the target Frame. */
3127    target_naxes = astGetNaxes( target );
3128 
3129 /* The first criterion for a match is that the template matches as a
3130    Frame class object. This ensures that the number of axes (1) and
3131    domain, etc. of the target Frame are suitable. Invoke the parent
3132    "astMatch" method to verify this. */
3133    match = (*parent_match)( template_frame, target, matchsub,
3134                             template_axes, target_axes, map, result, status );
3135 
3136 /* If a match was found, annul the returned objects, which are not
3137    needed, but keep the memory allocated for the axis association
3138    arrays, which we will re-use. */
3139    if ( astOK && match ) {
3140       *map = astAnnul( *map );
3141       *result = astAnnul( *result );
3142    }
3143 
3144 /* If OK so far, obtain pointers to the primary Frames which underlie
3145    all target axes. Stop when a SpecFrame axis is found. */
3146    if ( match && astOK ) {
3147       match = 0;
3148       for( iaxis = 0; iaxis < target_naxes; iaxis++ ) {
3149          astPrimaryFrame( target, iaxis, &frame0, &iaxis0 );
3150          if( astIsASpecFrame( frame0 ) ) {
3151             frame0 = astAnnul( frame0 );
3152             target_axis0 = iaxis;
3153             match = 1;
3154             break;
3155          } else {
3156             frame0 = astAnnul( frame0 );
3157          }
3158       }
3159 
3160    }
3161 
3162 /* Check at least one SpecFrame axis was found it the target. Store the
3163    axis associataions. */
3164    if( match && astOK ) {
3165       (*template_axes)[ 0 ] = 0;
3166       (*target_axes)[ 0 ] = target_axis0;
3167 
3168 /* Use the target's "astSubFrame" method to create a new Frame (the
3169    result Frame) with copies of the target axes in the required
3170    order. This process also overlays the template attributes on to the
3171    target Frame and returns a Mapping between the target and result
3172    Frames which effects the required coordinate conversion. */
3173       match = astSubFrame( target, template, 1, *target_axes, *template_axes,
3174                            map, result );
3175 
3176    }
3177 
3178 /* If an error occurred, or conversion to the result Frame's
3179    coordinate system was not possible, then free all memory, annul the
3180    returned objects, and reset the returned value. */
3181    if ( !astOK || !match ) {
3182       *template_axes = astFree( *template_axes );
3183       *target_axes = astFree( *target_axes );
3184       if( *map ) *map = astAnnul( *map );
3185       if( *result ) *result = astAnnul( *result );
3186       match = 0;
3187    }
3188 
3189 
3190 /* Return the result. */
3191    return match;
3192 }
3193 
OriginStdOfRest(AstSpecFrame * this,AstStdOfRestType newsor,const char * method,int * status)3194 static void OriginStdOfRest( AstSpecFrame *this, AstStdOfRestType newsor,
3195                              const char *method, int *status ){
3196 /*
3197 *  Name:
3198 *     OriginStdOfRest
3199 
3200 *  Purpose:
3201 *     Convert the SpecOrigin in a SpecFrame to a new rest frame.
3202 
3203 *  Type:
3204 *     Private function.
3205 
3206 *  Synopsis:
3207 *     #include "specframe.h"
3208 *     void OriginStdOfRest( AstSpecFrame *this, AstStdOfRestType newsor,
3209 *                           const char *method, int *status )
3210 
3211 *  Class Membership:
3212 *     SpecFrame member function
3213 
3214 *  Description:
3215 *     This function converts the value of the SpecOrigin attribute stored
3216 *     within a supplied SpecFrame from the rest frame currently associated
3217 *     with the SpecFrame, to the new rest frame indicated by "newsor".
3218 
3219 *  Parameters:
3220 *     this
3221 *        Point to the SpecFrame. On entry, the SpecOrigin value is
3222 *        assumed to refer to the re st frame given by the astGetStdOfRest
3223 *        method. On exit, the SpecOrigin value refers to the rest frame
3224 *        supplied in "newsor". The StdOfRest attribute of the SpecFrame
3225 *        should then be modified in order to keep things consistent.
3226 *     newsor
3227 *        The rest frame to which the SpecOrigin value stored within "this"
3228 *        should refer on exit.
3229 *     method
3230 *        Pointer to a string holding the name of the method to be
3231 *        included in any error messages.
3232 *     status
3233 *        Pointer to the inherited status variable.
3234 
3235 */
3236 
3237 
3238 /* Local Variables: */
3239    AstSpecFrame *sf;
3240    AstFrameSet *fs;
3241    double origin;
3242    double neworigin;
3243 
3244 /* Check the global error status. */
3245    if ( !astOK ) return;
3246 
3247 /* Do nothing if the SpecOrigin attribute has not been assigned a value. */
3248    if( astTestSpecOrigin( this ) ) {
3249 
3250 /* Do nothing if the rest frame will not change. */
3251       if( newsor != astGetStdOfRest( this ) ) {
3252 
3253 /* Save the original SpecOrigin value (in the current SpecFrame units) and then
3254    clear it. */
3255          origin = GetSpecOriginCur( this, status );
3256          astClearSpecOrigin( this );
3257 
3258 /* Take a copy of the SpecFrame and set the new StdOfRest. */
3259          sf = astCopy( this );
3260          astSetStdOfRest( sf, newsor );
3261 
3262 /* Create a Mapping to perform the rest frame change, then use it to convert
3263    the value to the new rest frame. */
3264          fs = astConvert( this, sf, "" );
3265          neworigin = AST__BAD;
3266          if( fs ) {
3267             astTran1( fs, 1, &origin, 1, &neworigin );
3268             fs = astAnnul( fs );
3269          }
3270 
3271 /* If succesful, convert from the current units to the default units, and store
3272    in "this". */
3273          if( neworigin != AST__BAD ) {
3274             astSetSpecOrigin( this, ToUnits( this, astGetUnit( this, 0 ), neworigin,
3275                               method, status ) );
3276 
3277          } else if( astOK ) {
3278             astError( AST__ATSER, "%s(%s): Cannot convert the SpecOrigin "
3279                       "value to a different rest frame.", status, method,
3280                       astGetClass( this ) );
3281          }
3282       }
3283    }
3284 }
3285 
OriginSystem(AstSpecFrame * this,AstSystemType oldsys,const char * method,int * status)3286 static void OriginSystem( AstSpecFrame *this, AstSystemType oldsys,
3287                           const char *method, int *status ){
3288 /*
3289 *  Name:
3290 *     OriginSystem
3291 
3292 *  Purpose:
3293 *     Convert the SpecOrigin in a SpecFrame to a new System.
3294 
3295 *  Type:
3296 *     Private function.
3297 
3298 *  Synopsis:
3299 *     #include "specframe.h"
3300 *     void OriginSystem( AstSpecFrame *this, AstSystemType oldsys,
3301 *                        const char *method, int *status )
3302 
3303 *  Class Membership:
3304 *     SpecFrame member function
3305 
3306 *  Description:
3307 *     This function converts the value of the SpecOrigin attribute stored
3308 *     within a supplied SpecFrame from its original System, etc, to the
3309 *     System, etc, currently associated with the SpecFrame.
3310 
3311 *  Parameters:
3312 *     this
3313 *        Point to the SpecFrame. On entry, the SpecOrigin value is
3314 *        assumed to refer to the System given by "oldsys", etc. On exit, the
3315 *        SpecOrigin value refers to the System returned by the astGetSystem
3316 *        method, etc.
3317 *     oldsys
3318 *        The System to which the SpecOrigin value stored within "this"
3319 *        refers on entry.
3320 *     method
3321 *        A string containing the method name for error messages.
3322 *     status
3323 *        Pointer to the inherited status variable.
3324 
3325 */
3326 
3327 /* Local Variables: */
3328    AstSpecFrame *sf1;
3329    AstSpecFrame *sf2;
3330    AstFrameSet *fs;
3331    double origin;
3332    double neworigin;
3333 
3334 /* Check the global error status. */
3335    if ( !astOK ) return;
3336 
3337 /* Do nothing if the SpecOrigin attribute has not been assigned a value. */
3338    if( astTestSpecOrigin( this ) ) {
3339 
3340 /* Do nothing if the System will not change. */
3341       if( oldsys != astGetSystem( this ) ) {
3342 
3343 /* Note the original SpecOrigin value, in the SpecFrame's default units. */
3344          origin = astGetSpecOrigin( this );
3345 
3346 /* Take a copy of the original SpecFrame and ensure the Units, SpecOrigin and
3347    AlignSpecOffset attributes are cleared. */
3348          sf1 = astCopy( this );
3349          astClearUnit( sf1, 0 );
3350          astClearSpecOrigin( sf1 );
3351          astClearAlignSpecOffset( sf1 );
3352 
3353 /* Take another copy of the SpecFrame and set the old system. */
3354          sf2 = astCopy( sf1 );
3355          astSetSystem( sf2, oldsys );
3356 
3357 /* Create a Mapping to perform the rest frame change, then use it to convert
3358    the value to the current system. */
3359          fs = astConvert( sf2, sf1, "" );
3360          neworigin = AST__BAD;
3361          if( fs ) {
3362             astTran1( fs, 1, &origin, 1, &neworigin );
3363             fs = astAnnul( fs );
3364          }
3365 
3366 /* Free resources */
3367          sf1 = astAnnul( sf1 );
3368          sf2 = astAnnul( sf2 );
3369 
3370 /* If succesful, store it in "this". */
3371          if( neworigin != AST__BAD ) {
3372             astSetSpecOrigin( this, neworigin );
3373 
3374          } else if( astOK ) {
3375             astError( AST__ATSER, "%s(%s): Cannot convert the SpecOrigin "
3376                       "value to a different spectral system.", status, method,
3377                       astGetClass( this ) );
3378          }
3379       }
3380    }
3381 }
3382 
3383 
Overlay(AstFrame * template,const int * template_axes,AstFrame * result,int * status)3384 static void Overlay( AstFrame *template, const int *template_axes,
3385                      AstFrame *result, int *status ) {
3386 /*
3387 *  Name:
3388 *     Overlay
3389 
3390 *  Purpose:
3391 *     Overlay the attributes of a template SpecFrame on to another Frame.
3392 
3393 *  Type:
3394 *     Private function.
3395 
3396 *  Synopsis:
3397 *     #include "specframe.h"
3398 *     void Overlay( AstFrame *template, const int *template_axes,
3399 *                   AstFrame *result, int *status )
3400 
3401 *  Class Membership:
3402 *     SpecFrame member function (over-rides the protected astOverlay method
3403 *     inherited from the Frame class).
3404 
3405 *  Description:
3406 *     This function overlays attributes of a SpecFrame (the "template") on to
3407 *     another Frame, so as to over-ride selected attributes of that second
3408 *     Frame. Normally only those attributes which have been specifically set
3409 *     in the template will be transferred. This implements a form of
3410 *     defaulting, in which a Frame acquires attributes from the template, but
3411 *     retains its original attributes (as the default) if new values have not
3412 *     previously been explicitly set in the template.
3413 *
3414 *     Note that if the result Frame is a SpecFrame and a change of spectral
3415 *     coordinate system occurs as a result of overlaying its System
3416 *     attribute, then some of its original attribute values may no
3417 *     longer be appropriate (e.g. the Title, or attributes describing
3418 *     its axes). In this case, these will be cleared before overlaying
3419 *     any new values.
3420 
3421 *  Parameters:
3422 *     template
3423 *        Pointer to the template SpecFrame, for which values should have been
3424 *        explicitly set for any attribute which is to be transferred.
3425 *     template_axes
3426 *        Pointer to an array of int, with one element for each axis of the
3427 *        "result" Frame (see below). For each axis in the result frame, the
3428 *        corresponding element of this array should contain the (zero-based)
3429 *        index of the template axis to which it corresponds. This array is used
3430 *        to establish from which template axis any axis-dependent attributes
3431 *        should be obtained.
3432 *
3433 *        If any axis in the result Frame is not associated with a template
3434 *        axis, the corresponding element of this array should be set to -1.
3435 *
3436 *        If a NULL pointer is supplied, the template and result axis
3437 *        indicies are assumed to be identical.
3438 *     result
3439 *        Pointer to the Frame which is to receive the new attribute values.
3440 *     status
3441 *        Pointer to the inherited status variable.
3442 
3443 *  Returned Value:
3444 *     void
3445 
3446 *  Notes:
3447 *     -  In general, if the result Frame is not from the same class as the
3448 *     template SpecFrame, or from a class derived from it, then attributes may
3449 *     exist in the template SpecFrame which do not exist in the result Frame.
3450 *     In this case, these attributes will not be transferred.
3451 */
3452 
3453 /* Local Variables: */
3454    AstFrame *templt;             /* Copy of supplied template Frame */
3455    AstSystemType new_system;     /* Code identifying new cordinates */
3456    AstSystemType old_system;     /* Code identifying old coordinates */
3457    const char *method;           /* Pointer to method string */
3458    const char *new_class;        /* Pointer to template class string */
3459    const char *old_class;        /* Pointer to result class string */
3460    int specframe;                /* Result Frame is a SpecFrame? */
3461 
3462 /* Check the global error status. */
3463    if ( !astOK ) return;
3464 
3465 /* Initialise strings used in error messages. */
3466    new_class = astGetClass( template );
3467    old_class = astGetClass( result );
3468    method = "astOverlay";
3469 
3470 /* Get the old and new systems. */
3471    old_system = astGetSystem( result );
3472    new_system = astGetSystem( template );
3473 
3474 /* It may be necessary to make temporary changes to the template Frame
3475    below. In order to ensure that we make no permanent changes to the
3476    supplied frame, we will, if necessary, take a deep copy of the
3477    supplied Frame, storing a pointer to the copy in "templt". If it is
3478    not necessary to make any changes to the template, we still want
3479    "templt" to hold a usable pointer, so we initialise it now to hold a
3480    clone of the supplied pointer. This pointer will be replaced by a
3481    pointer to a deep copy (if required) below. */
3482    templt = astClone( template );
3483 
3484 /* If the result Frame is a SpecFrame, we must test to see if overlaying its
3485    System attribute will change the type of coordinate system it describes.
3486    Determine the value of this attribute for the result and template
3487    SpecFrames. */
3488    specframe = astIsASpecFrame( result );
3489    if( specframe ) {
3490 
3491 /* If the coordinate system will change, any value already set for the result
3492    SpecFrame's Title will no longer be appropriate, so clear it. */
3493       if ( new_system != old_system ) {
3494          astClearTitle( result );
3495 
3496 /* If the systems have the same default units, we can retain the current
3497    Unit value. */
3498          if( strcmp( DefUnit( new_system, method, new_class, status ),
3499                      DefUnit( old_system, method, old_class, status ) ) ) {
3500             astClearUnit( result, 0 );
3501          }
3502 
3503 /* If necessary, clear inappropriate values for all those axis attributes
3504    whose access functions are over-ridden by this class (these access functions
3505    will then provide suitable defaults appropriate to the new coordinate system
3506    instead). */
3507          astClearLabel( result, 0 );
3508          astClearSymbol( result, 0 );
3509       }
3510 
3511 /* If the result Frame is not a SpecFrame, we must temporarily clear the
3512    System and AlignSystem values since the values used by this class
3513    are only appropriate to this class. Use a deep copy to avoid the danger
3514    of making any permanent changes to the suppied Frame. */
3515    } else {
3516       if( astTestSystem( template ) ) {
3517          templt = astAnnul( templt );
3518          templt = astCopy( template );
3519          astClearSystem( templt );
3520          astClearAlignSystem( templt );
3521       }
3522    }
3523 
3524 /* Invoke the parent class astOverlay method to transfer attributes inherited
3525    from the parent class. */
3526    (*parent_overlay)( templt, template_axes, result, status );
3527 
3528 /* Check if the result Frame is a SpecFrame or from a class derived from
3529    SpecFrame. If not, we cannot transfer SpecFrame attributes to it as it is
3530    insufficiently specialised. In this case simply omit these attributes. */
3531    if ( specframe && astOK ) {
3532 
3533 /* Define macros that test whether an attribute is set in the template and,
3534    if so, transfers its value to the result. */
3535 #define OVERLAY(attribute) \
3536    if ( astTest##attribute( template ) ) { \
3537       astSet##attribute( result, astGet##attribute( template ) ); \
3538    }
3539 
3540 /* Use the macro to transfer each SpecFrame attribute in turn. Note,
3541    SourceVRF must be overlayed before SourceVel. Otherwise the stored value
3542    for SourceVel would be changed from the default SourceVRF to the specified
3543    SourceVRF when SourceVRF was overlayed. */
3544       OVERLAY(AlignStdOfRest)
3545       OVERLAY(AlignSpecOffset);
3546       OVERLAY(RefDec)
3547       OVERLAY(RefRA)
3548       OVERLAY(RestFreq)
3549       OVERLAY(SourceSys)
3550       OVERLAY(SourceVRF)
3551       OVERLAY(SourceVel)
3552       OVERLAY(StdOfRest)
3553       OVERLAY(SpecOrigin)
3554    }
3555 
3556 /* Free resources */
3557    templt = astAnnul( templt );
3558 
3559 /* Undefine macros local to this function. */
3560 #undef OVERLAY
3561 }
3562 
SetAttrib(AstObject * this_object,const char * setting,int * status)3563 static void SetAttrib( AstObject *this_object, const char *setting, int *status ) {
3564 /*
3565 *  Name:
3566 *     SetAttrib
3567 
3568 *  Purpose:
3569 *     Set an attribute value for a SpecFrame.
3570 
3571 *  Type:
3572 *     Private function.
3573 
3574 *  Synopsis:
3575 *     #include "specframe.h"
3576 *     void SetAttrib( AstObject *this, const char *setting, int *status )
3577 
3578 *  Class Membership:
3579 *     SpecFrame member function (extends the astSetAttrib method inherited from
3580 *     the Mapping class).
3581 
3582 *  Description:
3583 *     This function assigns an attribute value for a SpecFrame, the attribute
3584 *     and its value being specified by means of a string of the form:
3585 *
3586 *        "attribute= value "
3587 *
3588 *     Here, "attribute" specifies the attribute name and should be in lower
3589 *     case with no white space present. The value to the right of the "="
3590 *     should be a suitable textual representation of the value to be assigned
3591 *     and this will be interpreted according to the attribute's data type.
3592 *     White space surrounding the value is only significant for string
3593 *     attributes.
3594 
3595 *  Parameters:
3596 *     this
3597 *        Pointer to the SpecFrame.
3598 *     setting
3599 *        Pointer to a null terminated string specifying the new attribute
3600 *        value.
3601 *     status
3602 *        Pointer to the inherited status variable.
3603 
3604 *  Returned Value:
3605 *     void
3606 
3607 *  Notes:
3608 *     This protected method is intended to be invoked by the Object astSet
3609 *     method and makes additional attributes accessible to it.
3610 */
3611 
3612 /* Local Vaiables: */
3613    AstMapping *umap;             /* Mapping between units */
3614    AstSpecFrame *this;           /* Pointer to the SpecFrame structure */
3615    AstStdOfRestType sor;         /* Standard of rest type code */
3616    AstSystemType sys;            /* Spectral system type code */
3617    char *a;                      /* Pointer to next character */
3618    char *new_setting;            /* Pointer value to new attribute setting */
3619    double dval;                  /* Double atribute value */
3620    double dtemp;                 /* Temporary double atribute value */
3621    int ival;                     /* Integer attribute value */
3622    int len;                      /* Length of setting string */
3623    int ulen;                     /* Used length of setting string */
3624    int namelen;                  /* Length of attribute name in setting */
3625    int nc;                       /* Number of characters read by astSscanf */
3626    int off;                      /* Offset of attribute value */
3627 
3628 /* Check the global error status. */
3629    if ( !astOK ) return;
3630 
3631 /* Obtain a pointer to the SpecFrame structure. */
3632    this = (AstSpecFrame *) this_object;
3633 
3634 /* Create an FK5 J2000 SkyFrame which will be used for formatting and
3635    unformatting sky positions, etc. */
3636    LOCK_MUTEX2
3637    if( !skyframe ) {
3638       astBeginPM;
3639       skyframe = astSkyFrame( "system=FK5,equinox=J2000", status );
3640       astEndPM;
3641    }
3642    UNLOCK_MUTEX2
3643 
3644 /* Obtain the length of the setting string. */
3645    len = strlen( setting );
3646 
3647 /* Obtain the used length of the setting string. */
3648    ulen = astChrLen( setting );
3649 
3650 /* Test for each recognised attribute in turn, using "astSscanf" to parse the
3651    setting string and extract the attribute value (or an offset to it in the
3652    case of string values). In each case, use the value set in "nc" to check
3653    that the entire string was matched. Once a value has been obtained, use the
3654    appropriate method to set it. */
3655 
3656 /* First look for axis attributes defined by the Frame class. Since a
3657    SpecFrame has only 1 axis, we allow these attributes to be specified
3658    without a trailing "(axis)" string. */
3659    if ( !strncmp( setting, "direction=", 10 ) ||
3660         !strncmp( setting, "bottom=", 7 ) ||
3661         !strncmp( setting, "top=", 4 ) ||
3662         !strncmp( setting, "format=", 7 ) ||
3663         !strncmp( setting, "label=", 6 ) ||
3664         !strncmp( setting, "symbol=", 7 ) ||
3665         !strncmp( setting, "unit=", 5 ) ) {
3666 
3667 /* Create a new setting string from the original by appending the string
3668    "(1)" to the end of the attribute name and then use the parent SetAttrib
3669    method. */
3670       new_setting = astMalloc( len + 4 );
3671       if( new_setting ) {
3672          memcpy( new_setting, setting, len + 1 );
3673          a = strchr( new_setting, '=' );
3674          namelen = a - new_setting;
3675          memcpy( a, "(1)", 4 );
3676          a += 3;
3677          strcpy( a, setting + namelen );
3678          (*parent_setattrib)( this_object, new_setting, status );
3679          new_setting = astFree( new_setting );
3680       }
3681 
3682 /* AlignStdOfRest. */
3683 /* --------------- */
3684    } else if ( nc = 0,
3685         ( 0 == astSscanf( setting, "alignstdofrest=%n%*s %n", &off, &nc ) )
3686         && ( nc >= len ) ) {
3687 
3688 /* Convert the string to a StdOfRest code before use. */
3689       sor = StdOfRestCode( setting + off, status );
3690       if ( sor != AST__BADSOR ) {
3691          astSetAlignStdOfRest( this, sor );
3692 
3693 /* Report an error if the string value wasn't recognised. */
3694       } else {
3695          astError( AST__ATTIN, "astSetAttrib(%s): Invalid standard of rest "
3696                    "description \"%s\".", status, astGetClass( this ), setting+off );
3697       }
3698 
3699 /* GeoLat. */
3700 /* ------- */
3701 /* Retained for backward compatibility with older versions of AST in which
3702    SpecFrame had GeoLon/Lat attributes (now ObsLon/Lat are used instead). */
3703    } else if ( nc = 0,
3704               ( 0 == astSscanf( setting, "geolat=%n%*s %n", &off, &nc ) )
3705               && ( nc >= 7 ) ) {
3706       new_setting = astStore( NULL, setting, len + 1 );
3707       new_setting[ 0 ] = 'o';
3708       new_setting[ 1 ] = 'b';
3709       new_setting[ 2 ] = 's';
3710       astSetAttrib( this, new_setting );
3711       new_setting = astFree( new_setting );
3712 
3713 /* GeoLon. */
3714 /* ------- */
3715    } else if ( nc = 0,
3716               ( 0 == astSscanf( setting, "geolon=%n%*s %n", &off, &nc ) )
3717               && ( nc >= 7 ) ) {
3718       new_setting = astStore( NULL, setting, len + 1 );
3719       new_setting[ 0 ] = 'o';
3720       new_setting[ 1 ] = 'b';
3721       new_setting[ 2 ] = 's';
3722       astSetAttrib( this, new_setting );
3723       new_setting = astFree( new_setting );
3724 
3725 /* RefDec. */
3726 /* ------- */
3727    } else if ( nc = 0,
3728               ( 0 == astSscanf( setting, "refdec=%n%*s %n", &off, &nc ) )
3729               && ( nc >= 7 ) ) {
3730 
3731 /* Convert the string to a radians value before use. */
3732       ival = astUnformat( skyframe, 1, setting + off, &dval );
3733       if ( ival == ulen - off  ) {
3734          astSetRefDec( this, dval );
3735 
3736 /* Report an error if the string value wasn't recognised. */
3737       } else {
3738          astError( AST__ATTIN, "astSetAttrib(%s): Invalid reference "
3739                    "declination \"%s\".", status, astGetClass( this ), setting + off );
3740       }
3741 
3742 /* RefRA. */
3743 /* ------ */
3744    } else if ( nc = 0,
3745               ( 0 == astSscanf( setting, "refra=%n%*s %n", &off, &nc ) )
3746               && ( nc >= 6 ) ) {
3747 
3748 /* Convert the string to a radians value before use. */
3749       ival = astUnformat( skyframe, 0, setting + off, &dval );
3750       if ( ival == ulen - off  ) {
3751          astSetRefRA( this, dval );
3752 
3753 /* Report an error if the string value wasn't recognised. */
3754       } else {
3755          astError( AST__ATTIN, "astSetAttrib(%s): Invalid reference right "
3756                    "ascension \"%s\".", status, astGetClass( this ), setting + off );
3757       }
3758 
3759 /* AlignSpecOffset. */
3760 /* ---------------- */
3761    } else if ( nc = 0,
3762              ( 1 == astSscanf( setting, "alignspecoffset= %d %n", &ival, &nc ) )
3763                && ( nc >= len ) ) {
3764       astSetAlignSpecOffset( this, ival );
3765 
3766 /* RestFreq. */
3767 /* --------- */
3768 /* Without any units indication - assume GHz. Convert to Hz for storage. */
3769    } else if ( nc = 0,
3770         ( 1 == astSscanf( setting, "restfreq= %lg %n", &dval, &nc ) )
3771         && ( nc >= len ) ) {
3772       astSetRestFreq( this, dval*1.0E9 );
3773 
3774 /* With units indication. */
3775    } else if ( nc = 0,
3776         ( 1 == astSscanf( setting, "restfreq= %lg %n%*s %n", &dval, &off, &nc ) )
3777         && ( nc >= len ) ) {
3778 
3779 /* Is there a Mapping from the supplied units to Hz? If so, use the
3780    Mapping to convert the supplied value to Hz. */
3781       if( ( umap = astUnitMapper( setting + off, "Hz", NULL, NULL ) ) ) {
3782          astTran1( umap, 1, &dval, 1, &dtemp );
3783          umap = astAnnul( umap );
3784 
3785 /* Otherwise, if there is a Mapping from the supplied units to metre,
3786    assume the supplied unit is a vacuum wavelength. */
3787       } else if( ( umap = astUnitMapper( setting + off, "m", NULL, NULL ) ) ) {
3788 
3789 /* Convert the supplied wavelength to metres. */
3790          astTran1( umap, 1, &dval, 1, &dtemp );
3791          umap = astAnnul( umap );
3792 
3793 /* Convert the wavelength (m) to frequency (Hz). */
3794          if( dtemp != AST__BAD && dtemp != 0.0 ) {
3795             dtemp = AST__C/dtemp;
3796          } else if( astOK ) {
3797             astError( AST__ATTIN, "astSetAttrib(%s): Invalid rest wavelength "
3798                    "\"%g %s\" supplied.", status, astGetClass( this ), dval, setting + off );
3799          }
3800 
3801 /* Otherwise, if there is a Mapping from the supplied units to Joule,
3802    assume the supplied unit is an energy. */
3803       } else if( ( umap = astUnitMapper( setting + off, "J", NULL, NULL ) ) ) {
3804 
3805 /* Convert the supplied energy to Joules. */
3806          astTran1( umap, 1, &dval, 1, &dtemp );
3807          umap = astAnnul( umap );
3808 
3809 /* Convert the energy (J) to frequency (Hz). */
3810          if( dtemp != AST__BAD ) {
3811             dtemp *= 1.0/AST__H;
3812          } else if( astOK ) {
3813             astError( AST__ATTIN, "astSetAttrib(%s): Invalid rest energy "
3814                    "\"%g %s\" supplied.", status, astGetClass( this ), dval, setting + off );
3815          }
3816 
3817 /* Otherwise report an error. */
3818       } else if( astOK ) {
3819          astError( AST__ATTIN, "astSetAttrib(%s): Rest frequency given in an "
3820                    "unsupported system of units \"%g %s\".", status,
3821                    astGetClass( this ), dval, setting + off );
3822       }
3823 
3824 /* Set the rest frequency. */
3825       astSetRestFreq( this, dtemp );
3826 
3827 /* SourceVel. */
3828 /* ---------- */
3829    } else if ( nc = 0,
3830         ( 1 == astSscanf( setting, "sourcevel= %lg %n", &dval, &nc ) )
3831         && ( nc >= len ) ) {
3832 
3833 /* Convert from km/s to m/s if the SourceVel value is a velocity. */
3834          if( astGetSourceSys( this ) == AST__VREL ||
3835              astGetSourceSys( this ) == AST__VRADIO ||
3836              astGetSourceSys( this ) == AST__VOPTICAL ) dval *= 1.0E3;
3837 
3838 /* Store the value */
3839       astSetSourceVel( this, dval );
3840 
3841 /* SourceVRF */
3842 /* --------- */
3843    } else if ( nc = 0,
3844         ( 0 == astSscanf( setting, "sourcevrf=%n%*s %n", &off, &nc ) )
3845         && ( nc >= len ) ) {
3846 
3847 /* Convert the string to a StdOfRest code before use. */
3848       sor = StdOfRestCode( setting + off, status );
3849       if ( sor != AST__BADSOR ) {
3850          astSetSourceVRF( this, sor );
3851 
3852 /* Report an error if the string value wasn't recognised. */
3853       } else {
3854          astError( AST__ATTIN, "astSetAttrib(%s): Invalid standard of rest "
3855                    "description \"%s\".", status, astGetClass( this ), setting+off );
3856       }
3857 
3858 /* SourceSys */
3859 /* --------- */
3860    } else if ( nc = 0,
3861         ( 0 == astSscanf( setting, "sourcesys=%n%*s %n", &off, &nc ) )
3862         && ( nc >= len ) ) {
3863 
3864 /* Convert the string to a System code before use. */
3865       sys = SystemCode( (AstFrame *) this, setting + off, status );
3866       astSetSourceSys( this, sys );
3867 
3868 /* StdOfRest. */
3869 /* ---------- */
3870    } else if ( nc = 0,
3871               ( 0 == astSscanf( setting, "stdofrest=%n%*s %n", &off, &nc ) )
3872               && ( nc >= len ) ) {
3873 
3874 /* Convert the string to a StdOfRest code before use. */
3875       sor = StdOfRestCode( setting + off, status );
3876       if ( sor != AST__BADSOR ) {
3877          astSetStdOfRest( this, sor );
3878 
3879 /* Report an error if the string value wasn't recognised. */
3880       } else {
3881          astError( AST__ATTIN, "astSetAttrib(%s): Invalid standard of rest "
3882                    "description \"%s\".", status, astGetClass( this ), setting + off );
3883       }
3884 
3885 /* SpecOrigin */
3886 /* ---------- */
3887 
3888 /* Floating-point without any units indication - assume the current Unit
3889    value. Convert from current units to default units for current system. */
3890    } else if ( nc = 0,
3891         ( 1 == astSscanf( setting, "specorigin= %lg %n", &dval, &nc ) )
3892         && ( nc >= len ) ) {
3893 
3894       astSetSpecOrigin( this, ToUnits( this, astGetUnit( this, 0 ), dval,
3895                                        "astSetSpecOrigin", status ) );
3896 
3897 /* Floating-point with units. Convert the supplied value to the default units
3898    for the SpecFrame's System. */
3899    } else if ( nc = 0,
3900         ( 1 == astSscanf( setting, "specorigin= %lg %n%*s %n", &dval, &off, &nc ) )
3901         && ( nc >= len ) ) {
3902       astSetSpecOrigin( this, ToUnits( this, setting + off, dval, "astSetSpecOrigin", status ) );
3903 
3904 /* Pass any unrecognised setting to the parent method for further
3905    interpretation. */
3906    } else {
3907       (*parent_setattrib)( this_object, setting, status );
3908    }
3909 }
3910 
SetRefPos(AstSpecFrame * this,AstSkyFrame * frm,double lon,double lat,int * status)3911 static void SetRefPos( AstSpecFrame *this, AstSkyFrame *frm, double lon,
3912                        double lat, int *status ){
3913 /*
3914 *++
3915 *  Name:
3916 c     astSetRefPos
3917 f     AST_SETREFPOS
3918 
3919 *  Purpose:
3920 *     Set the reference position in a specified celestial coordinate system.
3921 
3922 *  Type:
3923 *     Public virtual function.
3924 
3925 *  Synopsis:
3926 c     #include "specframe.h"
3927 c     void astSetRefPos( AstSpecFrame *this, AstSkyFrame *frm, double lon,
3928 c                        double lat )
3929 f     CALL AST_SETREFPOS( THIS, FRM, LON, LAT, STATUS )
3930 
3931 *  Class Membership:
3932 *     Frame method.
3933 
3934 *  Description:
3935 c     This function
3936 f     This routine
3937 *     sets the reference position (see attributes RefRA and RefDec) using
3938 *     axis values (in radians) supplied within the celestial coordinate
3939 *     system represented by a supplied SkyFrame.
3940 
3941 *  Parameters:
3942 c     this
3943 f     THIS = INTEGER (Given)
3944 *        Pointer to the SpecFrame.
3945 c     frm
3946 f     FRM = INTEGER (Given)
3947 *        Pointer to the SkyFrame which defines the celestial coordinate
3948 *        system in which the longitude and latitude values are supplied.
3949 c        If NULL
3950 f        If AST__NULL
3951 *        is supplied, then the supplied longitude and latitude values are
3952 *        assumed to be FK5 J2000 RA and Dec values.
3953 c     lon
3954 f     LON = DOUBLE PRECISION (Given)
3955 *        The longitude of the reference point, in the coordinate system
3956 *        represented by the supplied SkyFrame (radians).
3957 c     lat
3958 f     LAT = DOUBLE PRECISION (Given)
3959 *        The latitude of the reference point, in the coordinate system
3960 *        represented by the supplied SkyFrame (radians).
3961 f     STATUS = INTEGER (Given and Returned)
3962 f        The global status.
3963 
3964 *--
3965 */
3966 
3967 /* Local Variables: */
3968    AstFrameSet *fs;            /* Conversion FrameSet */
3969    AstFrame *fb;               /* Base Frame */
3970    AstFrame *fc;               /* Current Frame */
3971    double xin[ 1 ];            /* Axis 1 values */
3972    double yin[ 1 ];            /* Axis 2 values */
3973    double xout[ 1 ];           /* Axis 1 values */
3974    double yout[ 1 ];           /* Axis 2 values */
3975 
3976 /* Check the global error status. */
3977    if ( !astOK ) return;
3978 
3979 /* If no SkyFrame was supplied, just store the supplied RefRA and RefDec
3980    values. */
3981    if( !frm ) {
3982       astSetRefRA( this, lon );
3983       astSetRefDec( this, lat );
3984 
3985 /* Otherwise, convert the supplied values from the requested system. */
3986    } else {
3987 
3988 /* Create an FK5 J2000 SkyFrame which will be used for formatting and
3989    unformatting sky positions, etc. */
3990       LOCK_MUTEX2
3991       if( !skyframe ) {
3992          astBeginPM;
3993          skyframe = astSkyFrame( "system=FK5,equinox=J2000", status );
3994          astEndPM;
3995       }
3996       UNLOCK_MUTEX2
3997 
3998 /* Find the Mapping from the supplied SkyFrame, to the SkyFrame which
3999    describes the internal format in which the RefRA and RefDec attribute
4000    values are stored. */
4001       fs = astFindFrame( frm, skyframe, "" );
4002 
4003 /* If alignment was possible, use the Mapping to transform the supplied
4004    axis values, checking to see if the axes of the supplied SkyFrame have
4005    been permuted. */
4006       if( fs ) {
4007 
4008 /* Find the longitude axis in the Base Frame, and store the supplied
4009    longitude and latitude values. */
4010          fb = astGetFrame( fs, AST__BASE );
4011          if( astGetLonAxis( fb ) == 0 ) {
4012             xin[ 0 ] = lon;
4013             yin[ 0 ] = lat;
4014          } else {
4015             xin[ 0 ] = lat;
4016             yin[ 0 ] = lon;
4017          }
4018          astTran2( fs, 1, xin, yin, 1, xout, yout );
4019 
4020 /* Store the corresponding RefRA and RefDec values. */
4021          fc = astGetFrame( fs, AST__CURRENT );
4022          if( astGetLonAxis( fc ) == 0 ) {
4023             astSetRefRA( this, xout[ 0 ] );
4024             astSetRefDec( this, yout[ 0 ] );
4025          } else {
4026             astSetRefRA( this, yout[ 0 ] );
4027             astSetRefDec( this, xout[ 0 ] );
4028          }
4029 
4030 /* Annul object references. */
4031          fc = astAnnul( fc );
4032          fb = astAnnul( fb );
4033          fs = astAnnul( fs );
4034       }
4035    }
4036 }
4037 
SetStdOfRest(AstSpecFrame * this,AstStdOfRestType value,int * status)4038 static void SetStdOfRest( AstSpecFrame *this, AstStdOfRestType value, int *status ) {
4039 /*
4040 *+
4041 *  Name:
4042 *     astSetStdOfRest
4043 
4044 *  Purpose:
4045 *     Set the StdOfRest attribute for a SpecFrame.
4046 
4047 *  Type:
4048 *     Protected function.
4049 
4050 *  Synopsis:
4051 *     #include "specframe.h"
4052 *     void astSetStdOfRest( AstSpecFrame *this, AstStdOfRestType value )
4053 
4054 *  Class Membership:
4055 *     SpecFrame virtual function
4056 
4057 *  Description:
4058 *     This function set a new value for the StdOfRest attribute for a
4059 *     SpecFrame.
4060 
4061 *  Parameters:
4062 *     this
4063 *        Pointer to the SpecFrame.
4064 *     value
4065 *        The new value.
4066 
4067 *-
4068 */
4069 
4070 /* Check the global error status. */
4071    if ( !astOK ) return;
4072 
4073 
4074 /* Validate the StdOfRest value being set and report an error if necessary. */
4075    if( value < FIRST_SOR || value > LAST_SOR ) {
4076       astError( AST__ATTIN, "%s(%s): Bad value (%d) given for StdOfRest attribute.", status,
4077                             "astSetStdOfRest", astGetClass( this ), (int) value );
4078 
4079 /* Otherwise set the new StdOfRest */
4080    } else {
4081 
4082 /* Modify the SpecOrigin value stored in the SpecFrame structure to refer
4083    to the new rest frame. */
4084       OriginStdOfRest( this, value, "astSetStdOfRest", status );
4085 
4086 /* Store the new value for the rest frame in the SpecFrame structure. */
4087       this->stdofrest = value;
4088 
4089    }
4090 }
4091 
SetSystem(AstFrame * this_frame,AstSystemType newsys,int * status)4092 static void SetSystem( AstFrame *this_frame, AstSystemType newsys, int *status ) {
4093 /*
4094 *  Name:
4095 *     SetSystem
4096 
4097 *  Purpose:
4098 *     Set the System attribute for a SpecFrame.
4099 
4100 *  Type:
4101 *     Private function.
4102 
4103 *  Synopsis:
4104 *     #include "specframe.h"
4105 *     void SetSystem( AstFrame *this_frame, AstSystemType newsys, int *status )
4106 
4107 *  Class Membership:
4108 *     SpecFrame member function (over-rides the astSetSystem protected
4109 *     method inherited from the Frame class).
4110 
4111 *  Description:
4112 *     This function sets the System attribute for a SpecFrame.
4113 
4114 *  Parameters:
4115 *     this
4116 *        Pointer to the SpecFrame.
4117 *     newsys
4118 *        The new System value to be stored.
4119 *     status
4120 *        Pointer to the inherited status variable.
4121 
4122 */
4123 
4124 /* Local Variables: */
4125    AstSpecFrame *this;           /* Pointer to SpecFrame structure */
4126    AstSystemType oldsys;         /* Original System value */
4127 
4128 /* Check the global error status. */
4129    if ( !astOK ) return;
4130 
4131 /* Obtain a pointer to the SpecFrame structure. */
4132    this = (AstSpecFrame *) this_frame;
4133 
4134 /* Save the original System value */
4135    oldsys = astGetSystem( this_frame );
4136 
4137 /* Use the parent SetSystem method to store the new System value. */
4138    (*parent_setsystem)( this_frame, newsys, status );
4139 
4140 /* If the system has changed... */
4141    if( oldsys != newsys ) {
4142 
4143 /* Changing the System value will in general require the Units to change
4144    as well. If the user has previously specified the units to be used with
4145    the new system, then re-instate them (they are stored in the "usedunits"
4146    array in the SpecFrame structure). Otherwise, clear the units so that
4147    the default units will eb used with the new System. */
4148       if( (int) newsys < this->nuunits && this->usedunits &&
4149           this->usedunits[ (int) newsys ] ) {
4150          astSetUnit( this, 0, this->usedunits[ (int) newsys ] );
4151       } else {
4152          astClearUnit( this, 0 );
4153       }
4154 
4155 /* Modify the stored SpecOrigin. */
4156       OriginSystem( this, oldsys, "astSetSystem", status );
4157 
4158 /* Also, clear all attributes which have system-specific defaults. */
4159       astClearLabel( this_frame, 0 );
4160       astClearSymbol( this_frame, 0 );
4161       astClearTitle( this_frame );
4162    }
4163 }
4164 
SetUnit(AstFrame * this_frame,int axis,const char * value,int * status)4165 static void SetUnit( AstFrame *this_frame, int axis, const char *value, int *status ) {
4166 /*
4167 *  Name:
4168 *     SetUnit
4169 
4170 *  Purpose:
4171 *     Set a pointer to the Unit string for a SpecFrame's axis.
4172 
4173 *  Type:
4174 *     Private function.
4175 
4176 *  Synopsis:
4177 *     #include "specframe.h"
4178 *     void SetUnit( AstFrame *this_frame, int axis, const char *value )
4179 
4180 *  Class Membership:
4181 *     SpecFrame member function (over-rides the astSetUnit method inherited
4182 *     from the Frame class).
4183 
4184 *  Description:
4185 *     This function stores a pointer to the Unit string for a specified axis
4186 *     of a SpecFrame. It also stores the string in the "usedunits" array
4187 *     in the SpecFrame structure, in the element associated with the
4188 *     current System.
4189 
4190 *  Parameters:
4191 *     this
4192 *        Pointer to the SpecFrame.
4193 *     axis
4194 *        The number of the axis (zero-based) for which information is required.
4195 *     unit
4196 *        The new string to store.
4197 */
4198 
4199 /* Local Variables: */
4200    AstSpecFrame *this;           /* Pointer to the SpecFrame structure */
4201    int i;                        /* Loop counter */
4202    int system;                   /* The SpecFrame's System value */
4203 
4204 /* Check the global error status. */
4205    if ( !astOK ) return;
4206 
4207 /* Obtain a pointer to the SpecFrame structure. */
4208    this = (AstSpecFrame *) this_frame;
4209 
4210 /* Validate the axis index. */
4211    astValidateAxis( this, axis, 1, "astSetUnit" );
4212 
4213 /* Store the supplied value as the UsedUnit for the current System. First
4214    ensure the array is big enough. Free any previous value stored for the
4215    current system. */
4216    system = (int) astGetSystem( this );
4217    if( system >= this->nuunits ) {
4218       this->usedunits = astGrow( this->usedunits, system + 1,
4219                                  sizeof(char *) );
4220       if( astOK ) {
4221          for( i = this->nuunits; i < system + 1; i++ ) this->usedunits[ i ] = NULL;
4222          this->nuunits = system + 1;
4223       }
4224    }
4225 
4226 /* Now store a copy of the value, if it is different to the stored string. */
4227    if( astOK && ( !this->usedunits[ system ] ||
4228                   strcmp( this->usedunits[ system ], value ) ) ) {
4229       this->usedunits[ system ] = astStore( this->usedunits[ system ],
4230                                             value, strlen( value ) + 1 );
4231    }
4232 
4233 /* Now use the parent SetUnit method to store the value in the Axis
4234    structure */
4235    (*parent_setunit)( this_frame, axis, value, status );
4236 
4237 }
4238 
SorConvert(AstSpecFrame * this,AstSpecFrame * that,AstSpecMap * specmap,int * status)4239 static int SorConvert( AstSpecFrame *this, AstSpecFrame *that,
4240                        AstSpecMap *specmap, int *status ) {
4241 /*
4242 *  Name:
4243 *     SorConvert
4244 
4245 *  Purpose:
4246 *     Add a conversion to a SpecMap which transforms between two
4247 *     standards of rest.
4248 
4249 *  Type:
4250 *     Private function.
4251 
4252 *  Synopsis:
4253 *     #include "specframe.h"
4254 *     int SorConvert( AstSpecFrame *this, AstSpecFrame *that,
4255 *                     AstSpecMap *specmap, int *status )
4256 
4257 *  Class Membership:
4258 *     SpecFrame member function.
4259 
4260 *  Description:
4261 *     This function adds a conversion to a SpecMap which transforms
4262 *     frequencies from the standard of rest specified by "this" to
4263 *     the standard of rest specified by "that". Note the conversion is
4264 *     always between frequency in the two rest frames no matter what the
4265 *     System attributes of the two SpecFrames may be (which are ignored).
4266 
4267 *  Parameters:
4268 *     this
4269 *        The SpecFrame which defines the input rest frame.
4270 *     that
4271 *        The SpecFrame which defines the output rest frame.
4272 *     specmap
4273 *        The SpecMap to which the conversion is to be added.
4274 *     status
4275 *        Pointer to the inherited status variable.
4276 
4277 *  Returned Value:
4278 *     Zero is returned if the conversion could not be performed. One is
4279 *     returned otherwise.
4280 
4281 */
4282 
4283 /* Local Constants: */
4284 #define MAX_ARGS 7               /* Max arguments for an SpecMap conversion */
4285 
4286 /* Local Variables: */
4287    AstStdOfRestType from;        /* Input standard of rest */
4288    AstStdOfRestType to;          /* Output standard of rest */
4289    const char *vmess;            /* Text for use in error messages */
4290    double args[ MAX_ARGS ];      /* Conversion argument array */
4291    double dec;                   /* DEC of source (radians, FK5 J2000) */
4292    double epoch;                 /* Epoch of observation (MJD) */
4293    double alt;                   /* Observers geodetic altitude (radians) */
4294    double lat;                   /* Observers geodetic latitude (radians) */
4295    double lon;                   /* Observers geodetic longitude (radians) */
4296    double ra;                    /* RA of source (radians, FK5 J2000) */
4297    int result;                   /* Returned value */
4298 
4299 /* Initialise */
4300    result = 1;
4301 
4302 /* Check the global error status. */
4303    if ( !astOK ) return result;
4304 
4305 /* No conversion is required if the rest frames are equal. */
4306    if( !EqualSor( this, that, status ) ) {
4307 
4308 /* Define local macros as shorthand for adding spectral coordinate
4309    conversions to the SpecMap.  Each macro simply stores details of
4310    the additional arguments in the "args" array and then calls
4311    astSpecAdd. The macros differ in the number of additional argument
4312    values. */
4313 #define TRANSFORM_2(cvt,arg0,arg1) \
4314         args[ 0 ] = arg0; \
4315         args[ 1 ] = arg1; \
4316         astSpecAdd( specmap, cvt, args );
4317 
4318 #define TRANSFORM_3(cvt,arg0,arg1,arg2) \
4319         args[ 0 ] = arg0; \
4320         args[ 1 ] = arg1; \
4321         args[ 2 ] = arg2; \
4322         astSpecAdd( specmap, cvt, args );
4323 
4324 #define TRANSFORM_6(cvt,arg0,arg1,arg2,arg3,arg4,arg5) \
4325         args[ 0 ] = arg0; \
4326         args[ 1 ] = arg1; \
4327         args[ 2 ] = arg2; \
4328         args[ 3 ] = arg3; \
4329         args[ 4 ] = arg4; \
4330         args[ 5 ] = arg5; \
4331         astSpecAdd( specmap, cvt, args );
4332 
4333 /* A string for use in error messages. */
4334       vmess = "convert between different standards of rest";
4335 
4336 /* Get the required values from "this". */
4337       from = astGetStdOfRest( this );
4338       ra = astGetRefRA( this );
4339       dec = astGetRefDec( this );
4340       lon = astGetObsLon( this );
4341       lat = astGetObsLat( this );
4342       alt = astGetObsAlt( this );
4343       epoch = astGetEpoch( this );
4344 
4345 /* Verify that the reference RA and DEC can be used (they are needed by all
4346    the conversions used below). */
4347       VerifyAttrs( this, vmess, "RefRA RefDec", "astMatch", status );
4348 
4349 /* Convert from the "this" rest frame to heliographic. */
4350       if( from == AST__TPSOR ) {
4351          VerifyAttrs( this, vmess, "ObsLon ObsLat ObsAlt Epoch", "astMatch", status );
4352          TRANSFORM_6( "TPF2HL", lon, lat, alt, epoch, ra, dec )
4353 
4354       } else if( from == AST__GESOR ) {
4355          VerifyAttrs( this, vmess, "Epoch", "astMatch", status );
4356          TRANSFORM_3( "GEF2HL", epoch, ra, dec )
4357 
4358       } else if( from == AST__BYSOR ) {
4359          VerifyAttrs( this, vmess, "Epoch", "astMatch", status );
4360          TRANSFORM_3( "BYF2HL", epoch, ra, dec )
4361 
4362       } else if( from == AST__LKSOR ) {
4363          TRANSFORM_2( "LKF2HL", ra, dec )
4364 
4365       } else if( from == AST__LDSOR ) {
4366          TRANSFORM_2( "LDF2HL", ra, dec )
4367 
4368       } else if( from == AST__LGSOR ) {
4369          TRANSFORM_2( "LGF2HL", ra, dec )
4370 
4371       } else if( from == AST__GLSOR ) {
4372          TRANSFORM_2( "GLF2HL", ra, dec )
4373 
4374       } else if( from == AST__SCSOR ) {
4375          TRANSFORM_3( "USF2HL", ConvertSourceVel( this, AST__HLSOR, AST__VREL, status ),
4376                       ra, dec )
4377       }
4378 
4379 /* Now go from heliocentric to the "to" frame. */
4380       to = astGetStdOfRest( that );
4381       ra = astGetRefRA( that );
4382       dec = astGetRefDec( that );
4383       lon = astGetObsLon( that );
4384       lat = astGetObsLat( that );
4385       alt = astGetObsAlt( that );
4386       epoch = astGetEpoch( that );
4387       VerifyAttrs( that, vmess, "RefRA RefDec", "astMatch", status );
4388 
4389       if( to == AST__TPSOR ) {
4390          VerifyAttrs( that, vmess, "ObsLon ObsLat ObsAlt Epoch", "astMatch", status );
4391          TRANSFORM_6( "HLF2TP", lon, lat, alt, epoch, ra, dec )
4392 
4393       } else if( to == AST__GESOR ) {
4394          VerifyAttrs( that, vmess, "Epoch", "astMatch", status );
4395          TRANSFORM_3( "HLF2GE", epoch, ra, dec )
4396 
4397       } else if( to == AST__BYSOR ) {
4398          VerifyAttrs( that, vmess, "Epoch", "astMatch", status );
4399          TRANSFORM_3( "HLF2BY", epoch, ra, dec )
4400 
4401       } else if( to == AST__LKSOR ) {
4402          TRANSFORM_2( "HLF2LK", ra, dec )
4403 
4404       } else if( to == AST__LDSOR ) {
4405          TRANSFORM_2( "HLF2LD", ra, dec )
4406 
4407       } else if( to == AST__LGSOR ) {
4408          TRANSFORM_2( "HLF2LG", ra, dec )
4409 
4410       } else if( to == AST__GLSOR ) {
4411          TRANSFORM_2( "HLF2GL", ra, dec )
4412 
4413       } else if( to == AST__SCSOR ) {
4414          TRANSFORM_3( "HLF2US", ConvertSourceVel( that, AST__HLSOR, AST__VREL, status ),
4415                       ra, dec )
4416       }
4417    }
4418 
4419 /* Return the result. */
4420    return result;
4421 
4422 /* Undefine macros local to this function. */
4423 #undef MAX_ARGS
4424 #undef TRANSFORM_2
4425 #undef TRANSFORM_3
4426 #undef TRANSFORM_6
4427 }
4428 
SpecMapUnit(AstSystemType system,const char * method,const char * class,int * status)4429 static const char *SpecMapUnit( AstSystemType system, const char *method,
4430                                 const char *class, int *status ){
4431 /*
4432 *  Name:
4433 *     SpecMapUnit
4434 
4435 *  Purpose:
4436 *     Return the default units for a spectral coordinate system type used
4437 *     by the SpecMap class.
4438 
4439 *  Type:
4440 *     Private function.
4441 
4442 *  Synopsis:
4443 *     #include "specframe.h"
4444 *     const char *SpecMapUnit( AstSystemType system, const char *method,
4445 *                              const char *class, int *status )
4446 
4447 *  Class Membership:
4448 *     SpecFrame member function.
4449 
4450 *  Description:
4451 *     This function returns a textual representation of the
4452 *     units used by the SpecMap class for the specified spectral
4453 *     coordinate system. In general, the SpecMap class uses SI units
4454 *     (m/s, Hz, m, etc), but this class (SpecFrame) has default units
4455 *     more appropriate to astronomers (km/s, GHz, Angstroms, etc).
4456 
4457 *  Parameters:
4458 *     system
4459 *        The spectral coordinate system.
4460 *     method
4461 *        Pointer to a string holding the name of the calling method.
4462 *        This is only for use in constructing error messages.
4463 *     class
4464 *        Pointer to a string holding the name of the supplied object class.
4465 *        This is only for use in constructing error messages.
4466 *     status
4467 *        Pointer to the inherited status variable.
4468 
4469 *  Returned Value:
4470 *     A string describing the default units. This string follows the
4471 *     units syntax described in FITS WCS paper I "Representations of world
4472 *     coordinates in FITS" (Greisen & Calabretta).
4473 
4474 *  Notes:
4475 *     - A NULL pointer is returned if this function is invoked with
4476 *     the global error status set or if it should fail for any reason.
4477 */
4478 
4479 /* Local Variables: */
4480    const char *result;           /* Value to return */
4481 
4482 /* Initialize */
4483    result = NULL;
4484 
4485 /* Check the global error status. */
4486    if ( !astOK ) return result;
4487 
4488 /* Get an identifier for the default units. */
4489    if( system == AST__FREQ ) {
4490       result = "Hz";
4491    } else if( system == AST__ENERGY ) {
4492       result = "J";
4493    } else if( system == AST__WAVENUM ) {
4494       result = "1/m";
4495    } else if( system == AST__WAVELEN ) {
4496       result = "m";
4497    } else if( system == AST__AIRWAVE ) {
4498       result = "m";
4499    } else if( system == AST__VRADIO ) {
4500       result = "m/s";
4501    } else if( system == AST__VOPTICAL ) {
4502       result = "m/s";
4503    } else if( system == AST__REDSHIFT ) {
4504       result = "";
4505    } else if( system == AST__BETA ) {
4506       result = "";
4507    } else if( system == AST__VREL ) {
4508       result = "m/s";
4509 
4510 /* Report an error if the coordinate system was not recognised. */
4511    } else {
4512       astError( AST__SCSIN, "%s(%s): Corrupt %s contains illegal System "
4513                 "identification code (%d).", status, method, class, class,
4514                 (int) system );
4515    }
4516 
4517 /* Return the result. */
4518    return result;
4519 }
4520 
StdOfRestCode(const char * sor,int * status)4521 static AstStdOfRestType StdOfRestCode( const char *sor, int *status ) {
4522 /*
4523 *  Name:
4524 *     StdOfRestCode
4525 
4526 *  Purpose:
4527 *     Convert a string into a standard of rest type code.
4528 
4529 *  Type:
4530 *     Private function.
4531 
4532 *  Synopsis:
4533 *     #include "specframe.h"
4534 *     AstStdOfRestType StdOfRestCode( const char *sor )
4535 
4536 *  Class Membership:
4537 *     SpecFrame member function.
4538 
4539 *  Description:
4540 *     This function converts a string used for the external description of
4541 *     a standard of rest into a SpecFrame standard of rest type code
4542 *     (StdOfRest attribute value). It is the inverse of the
4543 *     StdOfRestString function.
4544 
4545 *  Parameters:
4546 *     sor
4547 *        Pointer to a constant null-terminated string containing the
4548 *        external description of the standard of rest.
4549 
4550 *  Returned Value:
4551 *     The StdOfRest type code.
4552 
4553 *  Notes:
4554 *     - A value of AST__BADSOR is returned if the standard of rest
4555 *     description was not recognised. This does not produce an error.
4556 *     - A value of AST__BADSOR is also returned if this function
4557 *     is invoked with the global error status set or if it should fail
4558 *     for any reason.
4559 */
4560 
4561 /* Local Variables: */
4562    AstStdOfRestType result;      /* Result value to return */
4563 
4564 /* Initialise. */
4565    result = AST__BADSOR;
4566 
4567 /* Check the global error status. */
4568    if ( !astOK ) return result;
4569 
4570 /* Match the "sor" string against each possibility and assign the
4571    result. */
4572    if ( astChrMatch( "TOPO", sor ) || astChrMatch( "TOPOCENT", sor ) || astChrMatch( "TOPOCENTRIC", sor ) ) {
4573       result = AST__TPSOR;
4574 
4575    } else if ( astChrMatch( "GEO", sor ) || astChrMatch( "GEOCENTR", sor ) || astChrMatch( "GEOCENTRIC", sor ) ) {
4576       result = AST__GESOR;
4577 
4578    } else if ( astChrMatch( "BARY", sor ) || astChrMatch( "BARYCENT", sor ) || astChrMatch( "BARYCENTRIC", sor ) ) {
4579       result = AST__BYSOR;
4580 
4581    } else if ( astChrMatch( "HELIO", sor ) || astChrMatch( "HELIOCEN", sor ) || astChrMatch( "HELIOCENTRIC", sor ) ) {
4582       result = AST__HLSOR;
4583 
4584    } else if ( astChrMatch( "LSRK", sor ) || astChrMatch( "LSR", sor ) ) {
4585       result = AST__LKSOR;
4586 
4587    } else if ( astChrMatch( "LSRD", sor ) ) {
4588       result = AST__LDSOR;
4589 
4590    } else if ( astChrMatch( "GAL", sor ) || astChrMatch( "GALACTOC", sor ) || astChrMatch( "GALACTIC", sor ) ) {
4591       result = AST__GLSOR;
4592 
4593    } else if ( astChrMatch( "LG", sor ) || astChrMatch( "LOCALGRP", sor ) ||
4594                astChrMatch( "LOCAL_GROUP", sor ) || astChrMatch( "LOCAL-GROUP", sor ) ) {
4595       result = AST__LGSOR;
4596 
4597    } else if ( astChrMatch( "SOURCE", sor ) || astChrMatch( "SRC", sor ) ) {
4598       result = AST__SCSOR;
4599 
4600    }
4601 
4602 /* Return the result. */
4603    return result;
4604 }
4605 
StdOfRestString(AstStdOfRestType sor,int * status)4606 static const char *StdOfRestString( AstStdOfRestType sor, int *status ) {
4607 /*
4608 *  Name:
4609 *     StdOfRestString
4610 
4611 *  Purpose:
4612 *     Convert a standard of rest type code into a string.
4613 
4614 *  Type:
4615 *     Private function.
4616 
4617 *  Synopsis:
4618 *     #include "specframe.h"
4619 *     const char *StdOfRestString( AstStdOfRestType sor, int *status )
4620 
4621 *  Class Membership:
4622 *     SpecFrame member function.
4623 
4624 *  Description:
4625 *     This function converts a SpecFrame standard of rest type code
4626 *     (StdOfRest attribute value) into a string suitable for use as an
4627 *     external representation of the standard of rest type.
4628 
4629 *  Parameters:
4630 *     sor
4631 *        The standard of rest type code.
4632 *     status
4633 *        Pointer to the inherited status variable.
4634 
4635 *  Returned Value:
4636 *     Pointer to a constant null-terminated string containing the
4637 *     textual equivalent of the type code supplied.
4638 
4639 *  Notes:
4640 *     - A NULL pointer value is returned if the standard of rest
4641 *     code was not recognised. This does not produce an error.
4642 *     - A NULL pointer value is also returned if this function is
4643 *     invoked with the global error status set or if it should fail
4644 *     for any reason.
4645 */
4646 
4647 /* Local Variables: */
4648    const char *result;           /* Pointer value to return */
4649 
4650 /* Initialise. */
4651    result = NULL;
4652 
4653 /* Check the global error status. */
4654    if ( !astOK ) return result;
4655 
4656 /* Match the "sor" value against each possibility and convert to a
4657    string pointer. (Where possible, return the same string as would be
4658    used in the FITS WCS representation of the standard of rest). */
4659    switch ( sor ) {
4660 
4661    case AST__TPSOR:
4662       result = "Topocentric";
4663       break;
4664 
4665    case AST__GESOR:
4666       result = "Geocentric";
4667       break;
4668 
4669    case AST__BYSOR:
4670       result = "Barycentric";
4671       break;
4672 
4673    case AST__HLSOR:
4674       result = "Heliocentric";
4675       break;
4676 
4677    case AST__LDSOR:
4678       result = "LSRD";
4679       break;
4680 
4681    case AST__LKSOR:
4682       result = "LSRK";
4683       break;
4684 
4685    case AST__LGSOR:
4686       result = "Local_group";
4687       break;
4688 
4689    case AST__GLSOR:
4690       result = "Galactic";
4691       break;
4692 
4693    case AST__SCSOR:
4694       result = "Source";
4695       break;
4696 
4697    }
4698 
4699 /* Return the result pointer. */
4700    return result;
4701 }
4702 
SubFrame(AstFrame * target_frame,AstFrame * template,int result_naxes,const int * target_axes,const int * template_axes,AstMapping ** map,AstFrame ** result,int * status)4703 static int SubFrame( AstFrame *target_frame, AstFrame *template,
4704                      int result_naxes, const int *target_axes,
4705                      const int *template_axes, AstMapping **map,
4706                      AstFrame **result, int *status ) {
4707 /*
4708 *  Name:
4709 *     SubFrame
4710 
4711 *  Purpose:
4712 *     Select axes from a SpecFrame and convert to the new coordinate
4713 *     system.
4714 
4715 *  Type:
4716 *     Private function.
4717 
4718 *  Synopsis:
4719 *     #include "specframe.h"
4720 *     int SubFrame( AstFrame *target, AstFrame *template,
4721 *                   int result_naxes, const int *target_axes,
4722 *                   const int *template_axes, AstMapping **map,
4723 *                   AstFrame **result, int *status )
4724 
4725 *  Class Membership:
4726 *     SpecFrame member function (over-rides the protected astSubFrame
4727 *     method inherited from the Frame class).
4728 
4729 *  Description:
4730 *     This function selects a requested sub-set (or super-set) of the axes
4731 *     from a "target" SpecFrame and creates a new Frame with copies of
4732 *     the selected axes assembled in the requested order. It then
4733 *     optionally overlays the attributes of a "template" Frame on to the
4734 *     result. It returns both the resulting Frame and a Mapping that
4735 *     describes how to convert between the coordinate systems described by
4736 *     the target and result Frames. If necessary, this Mapping takes
4737 *     account of any differences in the Frames' attributes due to the
4738 *     influence of the template.
4739 
4740 *  Parameters:
4741 *     target
4742 *        Pointer to the target SpecFrame, from which axes are to be
4743 *        selected.
4744 *     template
4745 *        Pointer to the template Frame, from which new attributes for the
4746 *        result Frame are to be obtained. Optionally, this may be NULL, in
4747 *        which case no overlaying of template attributes will be performed.
4748 *     result_naxes
4749 *        Number of axes to be selected from the target Frame. This number may
4750 *        be greater than or less than the number of axes in this Frame (or
4751 *        equal).
4752 *     target_axes
4753 *        Pointer to an array of int with result_naxes elements, giving a list
4754 *        of the (zero-based) axis indices of the axes to be selected from the
4755 *        target SpecFrame. The order in which these are given determines
4756 *        the order in which the axes appear in the result Frame. If any of the
4757 *        values in this array is set to -1, the corresponding result axis will
4758 *        not be derived from the target Frame, but will be assigned default
4759 *        attributes instead.
4760 *     template_axes
4761 *        Pointer to an array of int with result_naxes elements. This should
4762 *        contain a list of the template axes (given as zero-based axis indices)
4763 *        with which the axes of the result Frame are to be associated. This
4764 *        array determines which axes are used when overlaying axis-dependent
4765 *        attributes of the template on to the result. If any element of this
4766 *        array is set to -1, the corresponding result axis will not receive any
4767 *        template attributes.
4768 *
4769 *        If the template argument is given as NULL, this array is not used and
4770 *        a NULL pointer may also be supplied here.
4771 *     map
4772 *        Address of a location to receive a pointer to the returned Mapping.
4773 *        The forward transformation of this Mapping will describe how to
4774 *        convert coordinates from the coordinate system described by the target
4775 *        SpecFrame to that described by the result Frame. The inverse
4776 *        transformation will convert in the opposite direction.
4777 *     result
4778 *        Address of a location to receive a pointer to the result Frame.
4779 *     status
4780 *        Pointer to the inherited status variable.
4781 
4782 *  Returned Value:
4783 *     A non-zero value is returned if coordinate conversion is possible
4784 *     between the target and the result Frame. Otherwise zero is returned and
4785 *     *map and *result are returned as NULL (but this will not in itself
4786 *     result in an error condition). In general, coordinate conversion should
4787 *     always be possible if no template Frame is supplied but may not always
4788 *     be possible otherwise.
4789 
4790 *  Notes:
4791 *     -  A value of zero will be returned if this function is invoked with the
4792 *     global error status set, or if it should fail for any reason.
4793 
4794 *  Implementation Notes:
4795 *     -  This implementation addresses the selection of axes from a
4796 *     SpecFrame object. This results in another object of the same class
4797 *     only if the single SpecFrame axis is selected exactly once.
4798 *     Otherwise, the result is a Frame class object which inherits the
4799 *     SpecFrame's axis information (if appropriate) but none of the other
4800 *     properties of a SpecFrame.
4801 *     -  In the event that a SpecFrame results, the returned Mapping will
4802 *     take proper account of the relationship between the target and result
4803 *     coordinate systems.
4804 *     -  In the event that a Frame class object results, the returned Mapping
4805 *     will only represent a selection/permutation of axes.
4806 
4807 *  Implementation Deficiencies:
4808 *     -  Any axis selection is currently permitted. Probably this should be
4809 *     restricted so that each axis can only be selected once. The
4810 *     astValidateAxisSelection method will do this but currently there are bugs
4811 *     in the CmpFrame class that cause axis selections which will not pass this
4812 *     test. Install the validation when these are fixed.
4813 */
4814 
4815 /* Local Variables: */
4816    AstSpecFrame *target;      /* Pointer to the SpecFrame structure */
4817    AstSpecFrame *temp;        /* Pointer to copy of target SpecFrame */
4818    AstSpecFrame *align_frm;   /* Frame in which to align the SpecFrames */
4819    int match;                 /* Coordinate conversion is possible? */
4820    int report;                /* Report errors if SpecFrames cannot be aligned? */
4821 
4822 /* Initialise the returned values. */
4823    *map = NULL;
4824    *result = NULL;
4825    match = 0;
4826 
4827 /* Check the global error status. */
4828    if ( !astOK ) return match;
4829 
4830 /* Obtain a pointer to the target SpecFrame structure. */
4831    target = (AstSpecFrame *) target_frame;
4832 
4833 /* Result is a SpecFrame. */
4834 /* -------------------------- */
4835 /* Check if the result Frame is to have one axis obtained by selecting
4836    the single target SpecFrame axis. If so, the result will also be
4837    a SpecFrame. */
4838    if ( ( result_naxes == 1 ) && ( target_axes[ 0 ] == 0 ) ) {
4839 
4840 /* Form the result from a copy of the target. */
4841       *result = astCopy( target );
4842 
4843 /* Initialise a flag to indicate that MakeSpecMapping should not report
4844    errors if no Mapping can be created. */
4845       report = 0;
4846 
4847 /* If required, overlay the template attributes on to the result SpecFrame.
4848    Also get the system and standard of rest in which to align the two
4849    SpecFrames. These are the values from the template (if there is a
4850    template). */
4851       if ( template ) {
4852          astOverlay( template, template_axes, *result );
4853          if( astIsASpecFrame( template ) ) {
4854             align_frm = astCopy( template );
4855 
4856 /* Since we now know that both the template and target are SpecFrames, it
4857    should usually be possible to convert betwen them. If conversion is
4858    *not* possible (fpr instance if no rest frequency is availalbe, etc)
4859    then the user will probably be interested in knowing the reason why
4860    conversion is not possible. Therefore, indicate that MakeSpecMapping
4861    should report errors if no Mapping can be created. */
4862             report = 1;
4863 
4864          } else {
4865             align_frm = astCopy( target );
4866          }
4867 
4868 /* If no template was supplied, align in the System and StdOfRest of the
4869    target. */
4870       } else {
4871          VerifyAttrs( target, "convert between different spectral systems",
4872                    "StdOfRest", "astMatch", status );
4873          align_frm = astCopy( target );
4874       }
4875 
4876 /* The MakeSpecMapping function uses the System and StdOfRest attributes to
4877    define the alignment frame. But the AlignSystem and AlignStdOfRest
4878    attributes should be used for this purpose. Therefore, copy the values
4879    of the AlignSystem and AlignStdOfRest attributes to the System and
4880    StdOfRest attribute. */
4881       astSetSystem( align_frm, astGetAlignSystem( align_frm ) );
4882       astSetStdOfRest( align_frm, astGetAlignStdOfRest( align_frm ) );
4883 
4884 /* Generate a Mapping that takes account of changes in the sky coordinate
4885    system (equinox, epoch, etc.) between the target SpecFrame and the result
4886    SpecFrame. If this Mapping can be generated, set "match" to indicate that
4887    coordinate conversion is possible. If the template is a specframe,
4888    report errors if a match is not possible. */
4889       match = ( MakeSpecMapping( target, (AstSpecFrame *) *result,
4890                 align_frm, report, map, status ) != 0 );
4891 
4892 /* Free resources. */
4893       align_frm = astAnnul( align_frm );
4894 
4895 /* Result is not a SpecFrame. */
4896 /* ------------------------------ */
4897 /* In this case, we select axes as if the target were from the Frame
4898    class.  However, since the resulting data will then be separated
4899    from their enclosing SpecFrame, default attribute values may differ
4900    if the methods for obtaining them were over-ridden by the SpecFrame
4901    class. To overcome this, we ensure that these values are explicitly
4902    set for the result Frame (rather than relying on their defaults). */
4903    } else {
4904 
4905 /* Make a temporary copy of the target SpecFrame. We will explicitly
4906    set the attribute values in this copy so as not to modify the original. */
4907       temp = astCopy( target );
4908 
4909 /* Define a macro to test if an attribute is set. If not, set it
4910    explicitly to its default value. */
4911 #define SET(attribute) \
4912    if ( !astTest##attribute( temp ) ) { \
4913       astSet##attribute( temp, astGet##attribute( temp ) ); \
4914    }
4915 
4916 /* Set attribute values which apply to the Frame as a whole and which
4917    we want to retain, but whose defaults are over-ridden by the
4918    SpecFrame class. */
4919       SET(Domain)
4920       SET(Title)
4921 
4922 /* Define a macro to test if an attribute is set for axis zero (the only
4923    axis of a SpecFrame). If not, set it explicitly to its default value. */
4924 #define SET_AXIS(attribute) \
4925    if ( !astTest##attribute( temp, 0 ) ) { \
4926       astSet##attribute( temp, 0, \
4927                          astGet##attribute( temp, 0 ) ); \
4928    }
4929 
4930 /* Use this macro to set explicit values for all the axis attributes
4931    for which the SpecFrame class over-rides the default value. */
4932       SET_AXIS(Label)
4933       SET_AXIS(Symbol)
4934       SET_AXIS(Unit)
4935 
4936 /* Clear attributes which have an extended range of values allowed by
4937    this class. */
4938       astClearSystem( temp );
4939       astClearAlignSystem( temp );
4940 
4941 /* Invoke the astSubFrame method inherited from the Frame class to
4942    produce the result Frame by selecting the required set of axes and
4943    overlaying the template Frame's attributes. */
4944       match = (*parent_subframe)( (AstFrame *) temp, template,
4945                                   result_naxes, target_axes, template_axes,
4946                                   map, result, status );
4947 
4948 /* Delete the temporary copy of the target SpecFrame. */
4949       temp = astDelete( temp );
4950    }
4951 
4952 /* If an error occurred or no match was found, annul the returned
4953    objects and reset the returned result. */
4954    if ( !astOK || !match ) {
4955       if( *map ) *map = astAnnul( *map );
4956       if( *result ) *result = astAnnul( *result );
4957       match = 0;
4958    }
4959 
4960 /* Return the result. */
4961    return match;
4962 
4963 /* Undefine macros local to this function. */
4964 #undef SET
4965 #undef SET_AXIS
4966 }
4967 
SystemCode(AstFrame * this,const char * system,int * status)4968 static AstSystemType SystemCode( AstFrame *this, const char *system, int *status ) {
4969 /*
4970 *  Name:
4971 *     SystemCode
4972 
4973 *  Purpose:
4974 *     Convert a string into a coordinate system type code.
4975 
4976 *  Type:
4977 *     Private function.
4978 
4979 *  Synopsis:
4980 *     #include "specframe.h"
4981 *     AstSystemType SystemCode( AstFrame *this, const char *system, int *status )
4982 
4983 *  Class Membership:
4984 *     SpecFrame member function (over-rides the astSystemCode method
4985 *     inherited from the Frame class).
4986 
4987 *  Description:
4988 *     This function converts a string used for the external
4989 *     description of a coordinate system into a SpecFrame
4990 *     coordinate system type code (System attribute value). It is the
4991 *     inverse of the astSystemString function.
4992 
4993 *  Parameters:
4994 *     this
4995 *        The Frame.
4996 *     system
4997 *        Pointer to a constant null-terminated string containing the
4998 *        external description of the sky coordinate system.
4999 *     status
5000 *        Pointer to the inherited status variable.
5001 
5002 *  Returned Value:
5003 *     The System type code.
5004 
5005 *  Notes:
5006 *     - A value of AST__BADSYSTEM is returned if the sky coordinate
5007 *     system description was not recognised. This does not produce an
5008 *     error.
5009 *     - A value of AST__BADSYSTEM is also returned if this function
5010 *     is invoked with the global error status set or if it should fail
5011 *     for any reason.
5012 */
5013 
5014 /* Local Variables: */
5015    AstSystemType result;      /* Result value to return */
5016 
5017 /* Initialise. */
5018    result = AST__BADSYSTEM;
5019 
5020 /* Check the global error status. */
5021    if ( !astOK ) return result;
5022 
5023 /* Match the "system" string against each possibility and assign the
5024    result. */
5025    if ( astChrMatch( "FREQ", system ) ) {
5026       result = AST__FREQ;
5027 
5028    } else if ( astChrMatch( "ENER", system ) || astChrMatch( "ENERGY", system ) ) {
5029       result = AST__ENERGY;
5030 
5031    } else if ( astChrMatch( "WAVN", system ) || astChrMatch( "WAVENUM", system ) ) {
5032       result = AST__WAVENUM;
5033 
5034    } else if ( astChrMatch( "WAVE", system ) || astChrMatch( "WAVELEN", system ) ) {
5035       result = AST__WAVELEN;
5036 
5037    } else if ( astChrMatch( "AWAV", system ) || astChrMatch( "AIRWAVE", system ) ) {
5038       result = AST__AIRWAVE;
5039 
5040    } else if ( astChrMatch( "VRAD", system ) || astChrMatch( "VRADIO", system ) ) {
5041       result = AST__VRADIO;
5042 
5043    } else if ( astChrMatch( "VOPT", system ) || astChrMatch( "VOPTICAL", system ) ) {
5044       result = AST__VOPTICAL;
5045 
5046    } else if ( astChrMatch( "ZOPT", system ) || astChrMatch( "REDSHIFT", system ) ) {
5047       result = AST__REDSHIFT;
5048 
5049    } else if ( astChrMatch( "BETA", system ) ) {
5050       result = AST__BETA;
5051 
5052    } else if ( astChrMatch( "VELO", system ) || astChrMatch( "VREL", system ) ) {
5053       result = AST__VREL;
5054 
5055    }
5056 
5057 /* Return the result. */
5058    return result;
5059 }
5060 
SystemLabel(AstSystemType system,int * status)5061 static const char *SystemLabel( AstSystemType system, int *status ) {
5062 /*
5063 *  Name:
5064 *     SystemLabel
5065 
5066 *  Purpose:
5067 *     Return a label for a coordinate system type code.
5068 
5069 *  Type:
5070 *     Private function.
5071 
5072 *  Synopsis:
5073 *     #include "specframe.h"
5074 *     const char *SystemLabel( AstSystemType system, int *status )
5075 
5076 *  Class Membership:
5077 *     SpecFrame member function.
5078 
5079 *  Description:
5080 *     This function converts a SpecFrame coordinate system type code
5081 *     (System attribute value) into a descriptive string for human readers.
5082 
5083 *  Parameters:
5084 *     system
5085 *        The coordinate system type code.
5086 *     status
5087 *        Pointer to the inherited status variable.
5088 
5089 *  Returned Value:
5090 *     Pointer to a constant null-terminated string containing the
5091 *     textual equivalent of the type code supplied.
5092 
5093 *  Notes:
5094 *     - A NULL pointer value is returned if the sky coordinate system
5095 *     code was not recognised. This does not produce an error.
5096 *     - A NULL pointer value is also returned if this function is
5097 *     invoked with the global error status set or if it should fail
5098 *     for any reason.
5099 */
5100 
5101 /* Local Variables: */
5102    const char *result;           /* Pointer value to return */
5103 
5104 /* Initialise. */
5105    result = NULL;
5106 
5107 /* Check the global error status. */
5108    if ( !astOK ) return result;
5109 
5110 /* Match the "system" value against each possibility and convert to a
5111    string pointer. */
5112    switch ( system ) {
5113 
5114    case AST__FREQ:
5115       result = "frequency";
5116       break;
5117 
5118    case AST__ENERGY:
5119       result = "energy";
5120       break;
5121 
5122    case AST__WAVENUM:
5123       result = "wave-number";
5124       break;
5125 
5126    case AST__WAVELEN:
5127       result = "wavelength";
5128       break;
5129 
5130    case AST__AIRWAVE:
5131       result = "wavelength in air";
5132       break;
5133 
5134    case AST__VRADIO:
5135       result = "radio velocity";
5136       break;
5137 
5138    case AST__VOPTICAL:
5139       result = "optical velocity";
5140       break;
5141 
5142    case AST__REDSHIFT:
5143       result = "redshift";
5144       break;
5145 
5146    case AST__BETA:
5147       result = "beta factor";
5148       break;
5149 
5150    case AST__VREL:
5151       result = "apparent radial velocity";
5152       break;
5153    }
5154 
5155 /* Return the result pointer. */
5156    return result;
5157 }
5158 
SystemString(AstFrame * this,AstSystemType system,int * status)5159 static const char *SystemString( AstFrame *this, AstSystemType system, int *status ) {
5160 /*
5161 *  Name:
5162 *     SystemString
5163 
5164 *  Purpose:
5165 *     Convert a coordinate system type code into a string.
5166 
5167 *  Type:
5168 *     Private function.
5169 
5170 *  Synopsis:
5171 *     #include "specframe.h"
5172 *     const char *SystemString( AstFrame *this, AstSystemType system, int *status )
5173 
5174 *  Class Membership:
5175 *     SpecFrame member function (over-rides the astSystemString method
5176 *     inherited from the Frame class).
5177 
5178 *  Description:
5179 *     This function converts a SpecFrame coordinate system type code
5180 *     (System attribute value) into a string suitable for use as an
5181 *     external representation of the coordinate system type.
5182 
5183 *  Parameters:
5184 *     this
5185 *        The Frame.
5186 *     system
5187 *        The coordinate system type code.
5188 *     status
5189 *        Pointer to the inherited status variable.
5190 
5191 *  Returned Value:
5192 *     Pointer to a constant null-terminated string containing the
5193 *     textual equivalent of the type code supplied.
5194 
5195 *  Notes:
5196 *     - A NULL pointer value is returned if the sky coordinate system
5197 *     code was not recognised. This does not produce an error.
5198 *     - A NULL pointer value is also returned if this function is
5199 *     invoked with the global error status set or if it should fail
5200 *     for any reason.
5201 */
5202 
5203 /* Local Variables: */
5204    const char *result;           /* Pointer value to return */
5205 
5206 /* Initialise. */
5207    result = NULL;
5208 
5209 /* Check the global error status. */
5210    if ( !astOK ) return result;
5211 
5212 /* Match the "system" value against each possibility and convert to a
5213    string pointer. (Where possible, return the same string as would be
5214    used in the FITS WCS representation of the coordinate system). */
5215    switch ( system ) {
5216 
5217    case AST__FREQ:
5218       result = "FREQ";
5219       break;
5220 
5221    case AST__ENERGY:
5222       result = "ENER";
5223       break;
5224 
5225    case AST__WAVENUM:
5226       result = "WAVN";
5227       break;
5228 
5229    case AST__WAVELEN:
5230       result = "WAVE";
5231       break;
5232 
5233    case AST__AIRWAVE:
5234       result = "AWAV";
5235       break;
5236 
5237    case AST__VRADIO:
5238       result = "VRAD";
5239       break;
5240 
5241    case AST__VOPTICAL:
5242       result = "VOPT";
5243       break;
5244 
5245    case AST__REDSHIFT:
5246       result = "ZOPT";
5247       break;
5248 
5249    case AST__BETA:
5250       result = "BETA";
5251       break;
5252 
5253    case AST__VREL:
5254       result = "VELO";
5255       break;
5256    }
5257 
5258 /* Return the result pointer. */
5259    return result;
5260 }
5261 
TestActiveUnit(AstFrame * this_frame,int * status)5262 static int TestActiveUnit( AstFrame *this_frame, int *status ) {
5263 /*
5264 *  Name:
5265 *     TestActiveUnit
5266 
5267 *  Purpose:
5268 *     Test the ActiveUnit flag for a SpecFrame.
5269 
5270 *  Type:
5271 *     Private function.
5272 
5273 *  Synopsis:
5274 *     #include "specframe.h"
5275 *     int TestActiveUnit( AstFrame *this_frame, int *status )
5276 
5277 *  Class Membership:
5278 *     SpecFrame member function (over-rides the astTestActiveUnit protected
5279 *     method inherited from the Frame class).
5280 
5281 *  Description:
5282 *    This function test the value of the ActiveUnit flag for a SpecFrame,
5283 *    which is always "unset".
5284 
5285 *  Parameters:
5286 *     this
5287 *        Pointer to the SpecFrame.
5288 *     status
5289 *        Pointer to the inherited status variable.
5290 
5291 *  Returned Value:
5292 *     The result of the test (0).
5293 
5294 */
5295    return 0;
5296 }
5297 
TestAttrib(AstObject * this_object,const char * attrib,int * status)5298 static int TestAttrib( AstObject *this_object, const char *attrib, int *status ) {
5299 /*
5300 *  Name:
5301 *     TestAttrib
5302 
5303 *  Purpose:
5304 *     Test if a specified attribute value is set for a SpecFrame.
5305 
5306 *  Type:
5307 *     Private function.
5308 
5309 *  Synopsis:
5310 *     #include "specframe.h"
5311 *     int TestAttrib( AstObject *this, const char *attrib, int *status )
5312 
5313 *  Class Membership:
5314 *     SpecFrame member function (over-rides the astTestAttrib protected
5315 *     method inherited from the Frame class).
5316 
5317 *  Description:
5318 *     This function returns a boolean result (0 or 1) to indicate whether
5319 *     a value has been set for one of a SpecFrame's attributes.
5320 
5321 *  Parameters:
5322 *     this
5323 *        Pointer to the SpecFrame.
5324 *     attrib
5325 *        Pointer to a null terminated string specifying the attribute
5326 *        name.  This should be in lower case with no surrounding white
5327 *        space.
5328 *     status
5329 *        Pointer to the inherited status variable.
5330 
5331 *  Returned Value:
5332 *     One if a value has been set, otherwise zero.
5333 
5334 *  Notes:
5335 *     - This function uses one-based axis numbering so that it is
5336 *     suitable for external (public) use.
5337 *     - A value of zero will be returned if this function is invoked
5338 *     with the global status set, or if it should fail for any reason.
5339 */
5340 
5341 /* Local Variables: */
5342    AstSpecFrame *this;           /* Pointer to the SpecFrame structure */
5343    char *new_attrib;             /* Pointer value to new attribute name */
5344    int len;                      /* Length of attrib string */
5345    int result;                   /* Result value to return */
5346 
5347 /* Initialise. */
5348    result = 0;
5349 
5350 /* Check the global error status. */
5351    if ( !astOK ) return result;
5352 
5353 /* Obtain a pointer to the SpecFrame structure. */
5354    this = (AstSpecFrame *) this_object;
5355 
5356 /* Obtain the length of the attrib string. */
5357    len = strlen( attrib );
5358 
5359 /* Check the attribute name and test the appropriate attribute. */
5360 
5361 /* First look for axis attributes defined by the Frame class. Since a
5362    SpecFrame has only 1 axis, we allow these attributes to be specified
5363    without a trailing "(axis)" string. */
5364    if ( !strcmp( attrib, "direction" ) ||
5365         !strcmp( attrib, "bottom" ) ||
5366         !strcmp( attrib, "top" ) ||
5367         !strcmp( attrib, "format" ) ||
5368         !strcmp( attrib, "label" ) ||
5369         !strcmp( attrib, "symbol" ) ||
5370         !strcmp( attrib, "unit" ) ) {
5371 
5372 /* Create a new attribute name from the original by appending the string
5373    "(1)" and then use the parent TestAttrib method. */
5374       new_attrib = astMalloc( len + 4 );
5375       if( new_attrib ) {
5376          memcpy( new_attrib, attrib, len );
5377          memcpy( new_attrib + len, "(1)", 4 );
5378          result = (*parent_testattrib)( this_object, new_attrib, status );
5379          new_attrib = astFree( new_attrib );
5380       }
5381 
5382 /* AlignStdOfRest. */
5383 /* --------------- */
5384    } else if ( !strcmp( attrib, "alignstdofrest" ) ) {
5385       result = astTestAlignStdOfRest( this );
5386 
5387 /* GeoLat. */
5388 /* ------- */
5389 /* Retained for backward compatibility with older versions of AST in which
5390    SpecFrame had GeoLon/Lat attributes (now ObsLon/Lat are used instead). */
5391    } else if ( !strcmp( attrib, "geolat" ) ) {
5392       result = astTestAttrib( this, "obslat" );
5393 
5394 /* GeoLon. */
5395 /* ------- */
5396    } else if ( !strcmp( attrib, "geolon" ) ) {
5397       result = astTestAttrib( this, "obslon" );
5398 
5399 /* RefDec. */
5400 /* ------- */
5401    } else if ( !strcmp( attrib, "refdec" ) ) {
5402       result = astTestRefDec( this );
5403 
5404 /* RefRA. */
5405 /* ------ */
5406    } else if ( !strcmp( attrib, "refra" ) ) {
5407       result = astTestRefRA( this );
5408 
5409 /* RestFreq. */
5410 /* --------- */
5411    } else if ( !strcmp( attrib, "restfreq" ) ) {
5412       result = astTestRestFreq( this );
5413 
5414 /* SourceVel. */
5415 /* ---------- */
5416    } else if ( !strcmp( attrib, "sourcevel" ) ) {
5417       result = astTestSourceVel( this );
5418 
5419 /* SourceVRF */
5420 /* --------- */
5421    } else if ( !strcmp( attrib, "sourcevrf" ) ) {
5422       result = astTestSourceVRF( this );
5423 
5424 /* SourceSys */
5425 /* --------- */
5426    } else if ( !strcmp( attrib, "sourcesys" ) ) {
5427       result = astTestSourceSys( this );
5428 
5429 /* StdOfRest. */
5430 /* ---------- */
5431    } else if ( !strcmp( attrib, "stdofrest" ) ) {
5432       result = astTestStdOfRest( this );
5433 
5434 /* SpecOrigin. */
5435 /* --------- */
5436    } else if ( !strcmp( attrib, "specorigin" ) ) {
5437       result = astTestSpecOrigin( this );
5438 
5439 /* AlignSpecOffset */
5440 /* --------------- */
5441    } else if ( !strcmp( attrib, "alignspecoffset" ) ) {
5442       result = astTestAlignSpecOffset( this );
5443 
5444 /* If the attribute is not recognised, pass it on to the parent method
5445    for further interpretation. */
5446    } else {
5447       result = (*parent_testattrib)( this_object, attrib, status );
5448    }
5449 
5450 /* Return the result, */
5451    return result;
5452 }
5453 
ToUnits(AstSpecFrame * this,const char * oldunit,double oldval,const char * method,int * status)5454 static double ToUnits( AstSpecFrame *this, const char *oldunit, double oldval,
5455                        const char *method, int *status ){
5456 /*
5457 *
5458 *  Name:
5459 *     ToUnits
5460 
5461 *  Purpose:
5462 *     Convert a supplied spectral value to the default units of the supplied
5463 *     SpecFrame.
5464 
5465 *  Type:
5466 *     Private function.
5467 
5468 *  Synopsis:
5469 *     #include "timeframe.h"
5470 *     double ToUnits( AstSpecFrame *this, const char *oldunit, double oldval,
5471 *                     const char *method, int *status )
5472 
5473 *  Class Membership:
5474 *     SpecFrame member function
5475 
5476 *  Description:
5477 *     This function converts the supplied value from the supplied units to
5478 *     the default units associated with the supplied SpecFrame's System.
5479 
5480 *  Parameters:
5481 *     this
5482 *        Pointer to the SpecFrame.
5483 *     oldunit
5484 *        The units in which "oldval" is supplied.
5485 *     oldval
5486 *        The value to be converted.
5487 *     method
5488 *        Pointer to a string holding the name of the method to be
5489 *        included in any error messages.
5490 *     status
5491 *        Pointer to the inherited status variable.
5492 
5493 *  Returned Value:
5494 *     The converted value.
5495 
5496 */
5497 
5498 /* Local Variables: */
5499    AstMapping *map;
5500    const char *defunit;
5501    double result;
5502 
5503 /* Initialise. */
5504    result = AST__BAD;
5505 
5506 /* Check the global error status. */
5507    if ( !astOK ) return result;
5508 
5509 /* Get default units associated with the System attribute of the supplied
5510    SpecFrame, and find a Mapping from the old units to the default. */
5511    defunit = DefUnit( astGetSystem( this ), method, "SpecFrame", status );
5512    map = astUnitMapper( oldunit, defunit, NULL, NULL );
5513    if( map ) {
5514 
5515 /* Use the Mapping to convert the supplied value. */
5516       astTran1( map, 1, &oldval, 1, &result );
5517 
5518 /* Free resources. */
5519       map = astAnnul( map );
5520 
5521 /* Report an error if no conversion is possible. */
5522    } else if( astOK ){
5523       astError( AST__BADUN, "%s(%s): Cannot convert the supplied attribute "
5524                 "value from units of %s to %s.", status, method, astGetClass( this ),
5525                  oldunit, defunit );
5526    }
5527 
5528 /* Return the result */
5529    return result;
5530 }
5531 
5532 
ValidateSystem(AstFrame * this,AstSystemType system,const char * method,int * status)5533 static int ValidateSystem( AstFrame *this, AstSystemType system, const char *method, int *status ) {
5534 /*
5535 *
5536 *  Name:
5537 *     ValidateSystem
5538 
5539 *  Purpose:
5540 *     Validate a value for a Frame's System attribute.
5541 
5542 *  Type:
5543 *     Protected virtual function.
5544 
5545 *  Synopsis:
5546 *     #include "specframe.h"
5547 *     int ValidateSystem( AstFrame *this, AstSystemType system,
5548 *                         const char *method, int *status )
5549 
5550 *  Class Membership:
5551 *     SpecFrame member function (over-rides the astValidateSystem method
5552 *     inherited from the Frame class).
5553 
5554 *  Description:
5555 *     This function checks the validity of the supplied system value.
5556 *     If the value is valid, it is returned unchanged. Otherwise, an
5557 *     error is reported and a value of AST__BADSYSTEM is returned.
5558 
5559 *  Parameters:
5560 *     this
5561 *        Pointer to the Frame.
5562 *     system
5563 *        The system value to be checked.
5564 *     method
5565 *        Pointer to a constant null-terminated character string
5566 *        containing the name of the method that invoked this function
5567 *        to validate an axis index. This method name is used solely
5568 *        for constructing error messages.
5569 *     status
5570 *        Pointer to the inherited status variable.
5571 
5572 *  Returned Value:
5573 *     The validated system value.
5574 
5575 *  Notes:
5576 *     - A value of AST__BADSYSTEM will be returned if this function is invoked
5577 *     with the global error status set, or if it should fail for any
5578 *     reason.
5579 */
5580 
5581 /* Local Variables: */
5582    AstSystemType result;              /* Validated system value */
5583 
5584 /* Initialise. */
5585    result = AST__BADSYSTEM;
5586 
5587 /* Check the global error status. */
5588    if ( !astOK ) return result;
5589 
5590 /* If the value is out of bounds, report an error. */
5591    if ( system < FIRST_SYSTEM || system > LAST_SYSTEM ) {
5592          astError( AST__AXIIN, "%s(%s): Bad value (%d) given for the System "
5593                    "or AlignSystem attribute of a %s.", status, method,
5594                    astGetClass( this ), (int) system, astGetClass( this ) );
5595 
5596 /* Otherwise, return the supplied value. */
5597    } else {
5598       result = system;
5599    }
5600 
5601 /* Return the result. */
5602    return result;
5603 }
5604 
VerifyAttrs(AstSpecFrame * this,const char * purp,const char * attrs,const char * method,int * status)5605 static void VerifyAttrs( AstSpecFrame *this, const char *purp,
5606                          const char *attrs, const char *method, int *status ) {
5607 /*
5608 *  Name:
5609 *     VerifyAttrs
5610 
5611 *  Purpose:
5612 *     Verify that usable attribute values are available.
5613 
5614 *  Type:
5615 *     Private function.
5616 
5617 *  Synopsis:
5618 *     #include "specframe.h"
5619 *     void VerifyAttrs( AstSpecFrame *this, const char *purp,
5620 *                       const char *attrs, const char *method, int *status  )
5621 
5622 *  Class Membership:
5623 *     SpecFrame member function
5624 
5625 *  Description:
5626 *     This function tests each attribute listed in "attrs". It returns
5627 *     without action if 1) an explicit value has been set for each attribute
5628 *     or 2) the UseDefs attribute of the supplied SpecFrame is non-zero.
5629 *
5630 *     If UseDefs is zero (indicating that default values should not be
5631 *     used for attributes), and any of the named attributes does not have
5632 *     an explicitly set value, then an error is reported.
5633 
5634 *  Parameters:
5635 *     this
5636 *        Pointer to the SpecFrame.
5637 *     purp
5638 *        Pointer to a text string containing a message which will be
5639 *        included in any error report. This shouldindicate the purpose
5640 *        for which the attribute value is required.
5641 *     attrs
5642 *        A string holding a space separated list of attribute names.
5643 *     method
5644 *        A string holding the name of the calling method for use in error
5645 *        messages.
5646 *     status
5647 *        Pointer to the inherited status variable.
5648 
5649 */
5650 
5651 /* Local Variables: */
5652    const char *a;
5653    const char *desc;
5654    const char *p;
5655    int len;
5656    int set;
5657    int state;
5658 
5659 /* Check inherited status */
5660    if( !astOK ) return;
5661 
5662 /* If the SpecFrame has a non-zero value for its UseDefs attribute, then
5663    all attributes are assumed to have usable values, since the defaults
5664    will be used if no explicit value has been set. So we only need to do
5665    any checks if UseDefs is zero. */
5666    if( !astGetUseDefs( this ) ) {
5667 
5668 /* Stop compiler warnings about uninitialised variables */
5669       a = NULL;
5670       desc = NULL;
5671       len = 0;
5672       set = 0;
5673 
5674 /* Loop round the "attrs" string identifying the start and length of each
5675    non-blank word in the string. */
5676       state = 0;
5677       p = attrs;
5678       while( 1 ) {
5679          if( state == 0 ) {
5680             if( !isspace( *p ) ) {
5681                a = p;
5682                len = 1;
5683                state = 1;
5684             }
5685          } else {
5686             if( isspace( *p ) || !*p ) {
5687 
5688 /* The end of a word has just been reached. Compare it to each known
5689    attribute value. Get a flag indicating if the attribute has a set
5690    value, and a string describing the attribute.*/
5691                if( len > 0 ) {
5692 
5693                   if( !strncmp( "ObsLat", a, len ) ) {
5694                      set = astTestObsLat( this );
5695                      desc = "observer's latitude";
5696 
5697                   } else if( !strncmp( "ObsLon", a, len ) ) {
5698                      set = astTestObsLon( this );
5699                      desc = "observer's longitude";
5700 
5701                   } else if( !strncmp( "ObsAlt", a, len ) ) {
5702                      set = astTestObsAlt( this );
5703                      desc = "observer's altitude";
5704 
5705                   } else if( !strncmp( "RefRA", a, len ) ) {
5706                      set = astTestRefRA( this );
5707                      desc = "source RA";
5708 
5709                   } else if( !strncmp( "RefDec", a, len ) ) {
5710                      set = astTestRefDec( this );
5711                      desc = "source Dec";
5712 
5713                   } else if( !strncmp( "RestFreq", a, len ) ) {
5714                      set = astTestRestFreq( this );
5715                      desc = "rest frequency";
5716 
5717                   } else if( !strncmp( "SourceVel", a, len ) ) {
5718                      set = astTestSourceVel( this );
5719                      desc = "source velocity";
5720 
5721                   } else if( !strncmp( "StdOfRest", a, len ) ) {
5722                      set = astTestStdOfRest( this );
5723                      desc = "spectral standard of rest";
5724 
5725                   } else if( !strncmp( "Epoch", a, len ) ) {
5726                      set = astTestEpoch( this );
5727                      desc = "epoch of observation";
5728 
5729                   } else {
5730                      astError( AST__INTER, "VerifyAttrs(SpecFrame): "
5731                                "Unknown attribute name \"%.*s\" supplied (AST "
5732                                "internal programming error).", status, len, a );
5733                   }
5734 
5735 /* If the attribute does not have a set value, report an error. */
5736                   if( !set && astOK ) {
5737                      astError( AST__NOVAL, "%s(%s): Cannot %s.", status, method,
5738                                astGetClass( this ), purp );
5739                      astError( AST__NOVAL, "No value has been set for "
5740                                "the AST \"%.*s\" attribute (%s).", status, len, a,
5741                                desc );
5742                   }
5743 
5744 /* Continue the word search algorithm. */
5745                }
5746                len = 0;
5747                state = 0;
5748             } else {
5749                len++;
5750             }
5751          }
5752          if( !*(p++) ) break;
5753       }
5754    }
5755 }
5756 
5757 /* Functions which access class attributes. */
5758 /* ---------------------------------------- */
5759 /*
5760 *att++
5761 *  Name:
5762 *     AlignSpecOffset
5763 
5764 *  Purpose:
5765 *     Align SpecFrames using the offset coordinate system?
5766 
5767 *  Type:
5768 *     Public attribute.
5769 
5770 *  Synopsis:
5771 *     Integer (boolean).
5772 
5773 *  Description:
5774 *     This attribute is a boolean value which controls how a SpecFrame
5775 *     behaves when it is used (by
5776 c     astFindFrame or astConvert) as a template to match another (target)
5777 f     AST_FINDFRAME or AST_CONVERT) as a template to match another (target)
5778 *     SpecFrame. It determines whether alignment occurs between the offset
5779 *     values defined by the current value of the SpecOffset attribute, or
5780 *     between the corresponding absolute spectral values.
5781 *
5782 *     The default value of zero results in the two SpecFrames being aligned
5783 *     so that a given absolute spectral value in one is mapped to the same
5784 *     absolute value in the other. A non-zero value results in the SpecFrames
5785 *     being aligned so that a given offset value in one is mapped to the same
5786 *     offset value in the other.
5787 
5788 *  Applicability:
5789 *     SpecFrame
5790 *        All SpecFrames have this attribute.
5791 *att--
5792 */
5793 astMAKE_CLEAR(SpecFrame,AlignSpecOffset,alignspecoffset,-INT_MAX)
5794 astMAKE_GET(SpecFrame,AlignSpecOffset,int,0,( ( this->alignspecoffset != -INT_MAX ) ?
5795                                    this->alignspecoffset : 0 ))
5796 astMAKE_SET(SpecFrame,AlignSpecOffset,int,alignspecoffset,( value != 0 ))
5797 astMAKE_TEST(SpecFrame,AlignSpecOffset,( this->alignspecoffset != -INT_MAX ))
5798 
5799 
5800 
5801 /*
5802 *att++
5803 *  Name:
5804 *     AlignStdOfRest
5805 
5806 *  Purpose:
5807 *     Standard of rest to use when aligning SpecFrames.
5808 
5809 *  Type:
5810 *     Public attribute.
5811 
5812 *  Synopsis:
5813 *     String.
5814 
5815 *  Description:
5816 *     This attribute controls how a SpecFrame behaves when it is used (by
5817 c     astFindFrame or astConvert) as a template to match another (target)
5818 f     AST_FINDFRAME or AST_CONVERT) as a template to match another (target)
5819 *     SpecFrame. It identifies the standard of rest in which alignment is
5820 *     to occur. See the StdOfRest attribute for a desription of the values
5821 *     which may be assigned to this attribute. The default AlignStdOfRest
5822 *     value is "Helio" (heliographic).
5823 *
5824 c     When astFindFrame or astConvert is used on two SpecFrames (potentially
5825 f     When AST_FindFrame or AST_CONVERT is used on two SpecFrames (potentially
5826 *     describing different spectral coordinate systems), it returns a Mapping
5827 *     which can be used to transform a position in one SpecFrame into the
5828 *     corresponding position in the other. The Mapping is made up of the
5829 *     following steps in the indicated order:
5830 *
5831 *     - Map values from the system used by the target (wavelength,
5832 *     apparent radial velocity, etc) to the system specified by the
5833 *     AlignSystem attribute, using the target's rest frequency if necessary.
5834 *
5835 *     - Map these values from the target's standard of rest to the standard of
5836 *     rest specified by the AlignStdOfRest attribute, using the Epoch, ObsLat,
5837 *     ObsLon, ObsAlt, RefDec and RefRA attributes of the target to define the
5838 *     two standards of rest.
5839 *
5840 *     - Map these values from the standard of rest specified by the
5841 *     AlignStdOfRest attribute, to the template's standard of rest, using the
5842 *     Epoch, ObsLat, ObsLon, ObsAlt, RefDec and RefRA attributes of the
5843 *     template to define the two standards of rest.
5844 *
5845 *     - Map these values from the system specified by the AlignSystem
5846 *     attribute, to the system used by the template, using the template's
5847 *     rest frequency if necessary.
5848 
5849 *  Applicability:
5850 *     SpecFrame
5851 *        All SpecFrames have this attribute.
5852 
5853 *att--
5854 */
5855 /* The AlignStdOfRest value has a value of AST__BADSOR when not set yielding
5856    a default of AST__HLSOR. */
5857 astMAKE_TEST(SpecFrame,AlignStdOfRest,( this->alignstdofrest != AST__BADSOR ))
astMAKE_CLEAR(SpecFrame,AlignStdOfRest,alignstdofrest,AST__BADSOR)5858 astMAKE_CLEAR(SpecFrame,AlignStdOfRest,alignstdofrest,AST__BADSOR)
5859 astMAKE_GET(SpecFrame,AlignStdOfRest,AstStdOfRestType,AST__BADSOR,(
5860             ( this->alignstdofrest == AST__BADSOR ) ? AST__HLSOR : this->alignstdofrest ) )
5861 
5862 /* Validate the AlignStdOfRest value being set and report an error if necessary. */
5863 astMAKE_SET(SpecFrame,AlignStdOfRest,AstStdOfRestType,alignstdofrest,(
5864             ( ( value >= FIRST_SOR ) && ( value <= LAST_SOR ) ) ?
5865                  value :
5866                  ( astError( AST__ATTIN, "%s(%s): Bad value (%d) "
5867                              "given for AlignStdOfRest attribute.", status,
5868                              "astSetAlignStdOfRest", astGetClass( this ), (int) value ),
5869 
5870 /* Leave the value unchanged on error. */
5871                                             this->alignstdofrest ) ) )
5872 
5873 /*
5874 *att++
5875 *  Name:
5876 *     RefDec
5877 
5878 *  Purpose:
5879 *     The declination of the reference point
5880 
5881 *  Type:
5882 *     Public attribute.
5883 
5884 *  Synopsis:
5885 *     String.
5886 
5887 *  Description:
5888 *     This attribute specifies the FK5 J2000.0 declination of a reference
5889 *     point on the sky. See the description of attribute RefRA for details.
5890 *     The default RefDec is "0:0:0".
5891 
5892 *  Applicability:
5893 *     SpecFrame
5894 *        All SpecFrames have this attribute.
5895 
5896 *att--
5897 */
5898 /* The reference declination (FK5 J2000, radians). Clear the RefDec value by
5899    setting it to AST__BAD, which results in a default value of zero. Any
5900    value is acceptable. */
5901 astMAKE_CLEAR(SpecFrame,RefDec,refdec,AST__BAD)
5902 astMAKE_GET(SpecFrame,RefDec,double,0.0,((this->refdec!=AST__BAD)?this->refdec:0.0))
5903 astMAKE_SET(SpecFrame,RefDec,double,refdec,value)
5904 astMAKE_TEST(SpecFrame,RefDec,( this->refdec != AST__BAD ))
5905 
5906 /*
5907 *att++
5908 *  Name:
5909 *     RefRA
5910 
5911 *  Purpose:
5912 *     The right ascension of the reference point
5913 
5914 *  Type:
5915 *     Public attribute.
5916 
5917 *  Synopsis:
5918 *     String.
5919 
5920 *  Description:
5921 *     This attribute, together with the RefDec attribute, specifies the FK5
5922 *     J2000.0 coordinates of a reference point on the sky. For 1-dimensional
5923 *     spectra, this should normally be the position of the source. For
5924 *     spectral data with spatial coverage (spectral cubes, etc), this should
5925 *     be close to centre of the spatial coverage. It is used to define the
5926 *     correction for Doppler shift to be applied when using the
5927 c     astFindFrame or astConvert
5928 f     AST_FINDFRAME or AST_CONVERT
5929 *     method to convert between different standards of rest.
5930 *
5931 *     The SpecFrame class assumes this velocity correction is spatially
5932 *     invariant. If a single SpecFrame is used (for instance, as a
5933 *     component of a CmpFrame) to describe spectral values at different
5934 *     points on the sky, then it is assumes that the doppler shift at any
5935 *     spatial position is the same as at the reference position. The
5936 *     maximum velocity error introduced by this assumption is of the order
5937 *     of V*SIN(FOV), where FOV is the angular field of view, and V is the
5938 *     relative velocity of the two standards of rest. As an example, when
5939 *     correcting from the observers rest frame (i.e. the topocentric rest
5940 *     frame) to the kinematic local standard of rest the maximum value of V
5941 *     is about 20 km/s, so for 5 arc-minute field of view the maximum velocity
5942 *     error introduced by the correction will be about 0.03 km/s. As another
5943 *     example, the maximum error when correcting from the observers rest frame
5944 *     to the local group is about 5 km/s over a 1 degree field of view.
5945 *
5946 *     The RefRA and RefDec attributes are stored internally in radians, but
5947 *     are converted to and from a string for access. The format "hh:mm:ss.ss"
5948 *     is used for RefRA, and "dd:mm:ss.s" is used for RefDec. The methods
5949 c     astSetRefPos and astGetRefPos may be used to access the values of
5950 f     AST_SETREFPOS and AST_GETREFPOS may be used to access the value of
5951 *     these attributes directly as unformatted values in radians.
5952 *
5953 *     The default for RefRA is "0:0:0".
5954 
5955 *  Applicability:
5956 *     SpecFrame
5957 *        All SpecFrames have this attribute.
5958 
5959 *att--
5960 */
5961 /* The reference right ascension (FK5 J2000, radians). Clear the RefRA value
5962    by setting it to AST__BAD, which gives a default value of 0.0. Any
5963    value is acceptable. */
5964 astMAKE_CLEAR(SpecFrame,RefRA,refra,AST__BAD)
5965 astMAKE_GET(SpecFrame,RefRA,double,0.0,((this->refra!=AST__BAD)?this->refra:0.0))
5966 astMAKE_SET(SpecFrame,RefRA,double,refra,value)
5967 astMAKE_TEST(SpecFrame,RefRA,( this->refra != AST__BAD ))
5968 
5969 
5970 /*
5971 *att++
5972 *  Name:
5973 *     RestFreq
5974 
5975 *  Purpose:
5976 *     The rest frequency.
5977 
5978 *  Type:
5979 *     Public attribute.
5980 
5981 *  Synopsis:
5982 *     Floating point.
5983 
5984 *  Description:
5985 *     This attribute specifies the frequency corresponding to zero
5986 *     velocity. It is used when converting between between velocity-based
5987 *     coordinate systems and and other coordinate systems (such as frequency,
5988 *     wavelength, energy, etc). The default value is 1.0E5 GHz.
5989 *
5990 *     When setting a new value for this attribute, the new value can be
5991 *     supplied either directly as a frequency, or indirectly as a wavelength
5992 *     or energy, in which case the supplied value is converted to a frequency
5993 *     before being stored. The nature of the supplied value is indicated by
5994 *     appending text to the end of the numerical value indicating the units in
5995 *     which the value is supplied. If the units are not specified, then the
5996 *     supplied value is assumed to be a frequency in units of GHz. If the
5997 *     supplied unit is a unit of frequency, the supplied value is assumed to
5998 *     be a frequency in the given units. If the supplied unit is a unit of
5999 *     length, the supplied value is assumed to be a (vacuum) wavelength. If
6000 *     the supplied unit is a unit of energy, the supplied value is assumed to
6001 *     be an energy. For instance, the following strings all result in
6002 *     a rest frequency of around 1.4E14 Hz being used: "1.4E5", "1.4E14 Hz",
6003 *     "1.4E14 s**-1", "1.4E5 GHz", "2.14E-6 m", "21400 Angstrom", "9.28E-20 J",
6004 *     "9.28E-13 erg", "0.58 eV", etc.
6005 *
6006 *     When getting the value of this attribute, the returned value is
6007 *     always a frequency in units of GHz.
6008 
6009 *  Applicability:
6010 *     SpecFrame
6011 *        All SpecFrames have this attribute.
6012 
6013 *att--
6014 */
6015 /* The rest frequency (Hz). Clear the RestFreq value by setting it to AST__BAD,
6016    which gives 1.0E14 as the default value. Any value is acceptable. */
6017 astMAKE_CLEAR(SpecFrame,RestFreq,restfreq,AST__BAD)
6018 astMAKE_GET(SpecFrame,RestFreq,double,1.0E14,((this->restfreq!=AST__BAD)?this->restfreq:1.0E14))
6019 astMAKE_SET(SpecFrame,RestFreq,double,restfreq,value)
6020 astMAKE_TEST(SpecFrame,RestFreq,( this->restfreq != AST__BAD ))
6021 
6022 /*
6023 *att++
6024 *  Name:
6025 *     SourceVel
6026 
6027 *  Purpose:
6028 *     The source velocity.
6029 
6030 *  Type:
6031 *     Public attribute.
6032 
6033 *  Synopsis:
6034 *     Floating point.
6035 
6036 *  Description:
6037 *     This attribute (together with SourceSys, SourceVRF, RefRA and RefDec)
6038 *     defines the "Source" standard of rest (see attribute StdOfRest). This is
6039 *     a rest frame which is moving towards the position given by RefRA and
6040 *     RefDec at a  velocity given by SourceVel. A positive value means
6041 *     the source is moving away from the observer. When a new value is
6042 *     assigned to this attribute, the supplied value is assumed to refer
6043 *     to the spectral system specified by the SourceSys attribute. For
6044 *     instance, the SourceVel value may be supplied as a radio velocity, a
6045 *     redshift, a beta factor, etc. Similarly, when the current value of
6046 *     the SourceVel attribute is obtained, the returned value will refer
6047 *     to the spectral system specified by the SourceSys value. If the
6048 *     SourceSys value is changed, any value previously stored for the SourceVel
6049 *     attribute will be changed automatically from the old spectral system
6050 *     to the new spectral system.
6051 *
6052 *     When setting a value for SourceVel, the value should be supplied in the
6053 *     rest frame specified by the SourceVRF attribute. Likewise, when getting
6054 *     the value of SourceVel, it will be returned in the rest frame specified
6055 *     by the SourceVRF attribute.
6056 *
6057 *     The default SourceVel value is zero.
6058 
6059 *  Applicability:
6060 *     SpecFrame
6061 *        All SpecFrames have this attribute.
6062 
6063 *  Notes:
6064 *     - It is important to set an appropriate value for SourceVRF and
6065 *     SourceSys before setting a value for SourceVel. If a new value is later
6066 *     set for SourceVRF or SourceSys, the value stored for SourceVel will
6067 *     simultaneously be changed to the new standard of rest or spectral
6068 *     system.
6069 
6070 *att--
6071 */
6072 /* The source velocity (velocities are stored internally in m/s). Clear it
6073    by setting it to AST__BAD, which returns a default value of zero. Any
6074    value is acceptable. */
6075 astMAKE_CLEAR(SpecFrame,SourceVel,sourcevel,AST__BAD)
6076 astMAKE_SET(SpecFrame,SourceVel,double,sourcevel,value)
6077 astMAKE_TEST(SpecFrame,SourceVel,( this->sourcevel != AST__BAD ))
6078 astMAKE_GET(SpecFrame,SourceVel,double,0.0,((this->sourcevel!=AST__BAD)?this->sourcevel:0.0))
6079 
6080 /*
6081 *att++
6082 *  Name:
6083 *     SourceVRF
6084 
6085 *  Purpose:
6086 *     Rest frame in which the source velocity is stored.
6087 
6088 *  Type:
6089 *     Public attribute.
6090 
6091 *  Synopsis:
6092 *     String.
6093 
6094 *  Description:
6095 *     This attribute identifies the rest frame in which the source
6096 *     velocity or redshift is stored (the source velocity or redshift is
6097 *     accessed using attribute SourceVel). When setting a new value for the
6098 *     SourceVel attribute, the source velocity or redshift should be supplied
6099 *     in the rest frame indicated by this attribute. Likewise, when getting
6100 *     the value of the SourceVel attribute, the velocity or redshift will be
6101 *     returned in this rest frame.
6102 *
6103 *     If the value of SourceVRF is changed, the value stored for SourceVel
6104 *     will be converted from the old to the new rest frame.
6105 *
6106 *     The values which can be supplied are the same as for the StdOfRest
6107 *     attribute (except that SourceVRF cannot be set to "Source"). The
6108 *     default value is "Helio".
6109 
6110 *  Applicability:
6111 *     SpecFrame
6112 *        All SpecFrames have this attribute.
6113 
6114 *att--
6115 */
6116 /* The SourceVRF value has a value of AST__BADSOR when not set yielding
6117    a default of AST__HLSOR. */
6118 astMAKE_TEST(SpecFrame,SourceVRF,( this->sourcevrf != AST__BADSOR ))
6119 astMAKE_GET(SpecFrame,SourceVRF,AstStdOfRestType,AST__BADSOR,(
6120             ( this->sourcevrf == AST__BADSOR ) ? AST__HLSOR : this->sourcevrf ) )
6121 
6122 /* When clearing SourceVRF, convert the SourceVel value to heliocentric
6123   (but only if set)*/
6124 astMAKE_CLEAR(SpecFrame,SourceVRF,sourcevrf,((astTestSourceVel( this )?
6125 astSetSourceVel( this, ConvertSourceVel( this, AST__HLSOR, astGetSourceSys( this ), status ) ),NULL:NULL),AST__BADSOR))
6126 
6127 /* Validate the SourceVRF value being set and report an error if necessary.
6128    If OK, convert the stored SourceVel value into the new rest frame (but
6129 only if set)*/
6130 astMAKE_SET(SpecFrame,SourceVRF,AstStdOfRestType,sourcevrf,(
6131             ( ( value >= FIRST_SOR ) && ( value <= LAST_SOR ) && value != AST__SCSOR ) ?
6132                  (astTestSourceVel( this )?
6133                  astSetSourceVel( this, ConvertSourceVel( this, value, astGetSourceSys( this ), status )),NULL:NULL), value:
6134                  ( astError( AST__ATTIN, "%s(%s): Bad value (%d) "
6135                              "given for SourceVRF attribute.", status,
6136                              "astSetSourceVRF", astGetClass( this ), (int) value ),
6137 
6138 /* Leave the value unchanged on error. */
6139                                             this->sourcevrf ) ) )
6140 
6141 /*
6142 *att++
6143 *  Name:
6144 *     SourceSys
6145 
6146 *  Purpose:
6147 *     Spectral system in which the source velocity is stored.
6148 
6149 *  Type:
6150 *     Public attribute.
6151 
6152 *  Synopsis:
6153 *     String.
6154 
6155 *  Description:
6156 *     This attribute identifies the spectral system in which the
6157 *     SourceVel attribute value (the source velocity) is supplied and
6158 *     returned. It can be one of the following:
6159 *
6160 *        - "VRAD" or "VRADIO": Radio velocity (km/s)
6161 *        - "VOPT" or "VOPTICAL": Optical velocity (km/s)
6162 *        - "ZOPT" or "REDSHIFT": Redshift (dimensionless)
6163 *        - "BETA": Beta factor (dimensionless)
6164 *        - "VELO" or "VREL": Apparent radial ("relativistic") velocity (km/s)
6165 *
6166 *     When setting a new value for the SourceVel attribute, the source
6167 *     velocity should be supplied in the spectral system indicated
6168 *     by this attribute. Likewise, when getting the value of the SourceVel
6169 *     attribute, the velocity will be returned in this spectral system.
6170 *
6171 *     If the value of SourceSys is changed, the value stored for SourceVel
6172 *     will be converted from the old to the new spectral systems.
6173 *
6174 *     The default value is "VELO" (apparent radial velocity).
6175 
6176 *  Applicability:
6177 *     SpecFrame
6178 *        All SpecFrames have this attribute.
6179 
6180 *att--
6181 */
6182 /* The SourceSys value has a value of AST__BADSYS when not set yielding
6183    a default of AST__VREL. */
6184 astMAKE_TEST(SpecFrame,SourceSys,( this->sourcesys != AST__BADSYSTEM ))
6185 astMAKE_GET(SpecFrame,SourceSys,AstSystemType,AST__BADSYSTEM,(
6186             ( this->sourcesys == AST__BADSYSTEM ) ? AST__VREL : this->sourcesys ) )
6187 
6188 /* When clearing SourceSys, convert the SourceVel value to relativistic
6189    velocity (but only if set) */
6190 astMAKE_CLEAR(SpecFrame,SourceSys,sourcesys,((astTestSourceVel( this )?
6191 astSetSourceVel( this, ConvertSourceVel( this, astGetSourceVRF( this ),
6192                                          AST__VREL, status ) ),NULL:NULL),AST__BADSYSTEM))
6193 
6194 /* Validate the SourceSys value being set and report an error if necessary.
6195    If OK, convert the stored SourceVel value into the new rest frame (but
6196    only if set)*/
6197 astMAKE_SET(SpecFrame,SourceSys,AstSystemType,sourcesys,(
6198             ( ( value == AST__VREL ) || ( value == AST__BETA ) ||
6199               ( value == AST__VRADIO ) || ( value == AST__REDSHIFT ) ||
6200               ( value == AST__VOPTICAL ) ) ?
6201               (astTestSourceVel( this )?
6202                astSetSourceVel( this, ConvertSourceVel( this, astGetSourceVRF( this ),
6203                                                         value, status )),NULL:NULL),
6204                                                         value:
6205                  ( astError( AST__ATTIN, "%s(%s): Bad value (%d) "
6206                              "given for SourceSys attribute.", status,
6207                              "astSetSourceSys", astGetClass( this ), (int) value ),
6208 
6209 /* Leave the value unchanged on error. */
6210                                             this->sourcesys ) ) )
6211 
6212 /*
6213 *att++
6214 *  Name:
6215 *     StdOfRest
6216 
6217 *  Purpose:
6218 *     Standard of rest.
6219 
6220 *  Type:
6221 *     Public attribute.
6222 
6223 *  Synopsis:
6224 *     String.
6225 
6226 *  Description:
6227 *     This attribute identifies the standard of rest to which the spectral
6228 *     axis values of a SpecFrame refer, and may take any of the values
6229 *     listed in the "Standards of Rest" section (below).
6230 *
6231 *     The default StdOfRest value is "Helio".
6232 
6233 *  Applicability:
6234 *     SpecFrame
6235 *        All SpecFrames have this attribute.
6236 
6237 *  Standards of Rest:
6238 *     The SpecFrame class supports the following StdOfRest values (all are
6239 *     case-insensitive):
6240 *
6241 *     - "Topocentric", "Topocent" or "Topo": The observers rest-frame (assumed
6242 *     to be on the surface of the earth). Spectra recorded in this standard of
6243 *     rest suffer a Doppler shift which varies over the course of a day
6244 *     because of the rotation of the observer around the axis of the earth.
6245 *     This standard of rest must be qualified using the ObsLat, ObsLon,
6246 *     ObsAlt, Epoch, RefRA and RefDec attributes.
6247 *
6248 *     - "Geocentric", "Geocentr" or "Geo": The rest-frame of the earth centre.
6249 *     Spectra recorded in this standard of rest suffer a Doppler shift which
6250 *     varies over the course of a year because of the rotation of the earth
6251 *     around the Sun. This standard of rest must be qualified using the Epoch,
6252 *     RefRA and RefDec attributes.
6253 *
6254 *     - "Barycentric", "Barycent" or "Bary": The rest-frame of the solar-system
6255 *     barycentre. Spectra recorded in this standard of rest suffer a Doppler
6256 *     shift which depends both on the velocity of the Sun through the Local
6257 *     Standard of Rest, and on the movement of the planets through the solar
6258 *     system. This standard of rest must be qualified using the Epoch, RefRA
6259 *     and RefDec attributes.
6260 *
6261 *     - "Heliocentric", "Heliocen" or "Helio": The rest-frame of the Sun.
6262 *     Spectra recorded in this standard of rest suffer a Doppler shift which
6263 *     depends on the velocity of the Sun through the Local Standard of Rest.
6264 *     This standard of rest must be qualified using the RefRA and RefDec
6265 *     attributes.
6266 *
6267 *     - "LSRK", "LSR": The rest-frame of the kinematical Local Standard of
6268 *     Rest. Spectra recorded in this standard of rest suffer a Doppler shift
6269 *     which depends on the velocity of the kinematical Local Standard of Rest
6270 *     through the galaxy. This standard of rest must be qualified using the
6271 *     RefRA and RefDec attributes.
6272 *
6273 *     - "LSRD": The rest-frame of the dynamical Local Standard of Rest. Spectra
6274 *     recorded in this standard of rest suffer a Doppler shift which depends
6275 *     on the velocity of the dynamical Local Standard of Rest through the
6276 *     galaxy.  This standard of rest must be qualified using the RefRA and
6277 *     RefDec attributes.
6278 *
6279 *     - "Galactic", "Galactoc" or "Gal": The rest-frame of the galactic centre.
6280 *     Spectra recorded in this standard of rest suffer a Doppler shift which
6281 *     depends on the velocity of the galactic centre through the local group.
6282 *     This standard of rest must be qualified using the RefRA and RefDec
6283 *     attributes.
6284 *
6285 *     - "Local_group", "Localgrp" or "LG": The rest-frame of the local group.
6286 *     This standard of rest must be qualified using the RefRA and RefDec
6287 *     attributes.
6288 *
6289 *     - "Source", or "src": The rest-frame of the source. This standard of
6290 *     rest must be qualified using the RefRA, RefDec and SourceVel attributes.
6291 *
6292 *     Where more than one alternative System value is shown above, the
6293 *     first of these will be returned when an enquiry is made.
6294 *att--
6295 */
6296 /* The StdOfRest value has a value of AST__BADSOR when not set yielding
6297    a default of AST__HLSOR. */
6298 astMAKE_TEST(SpecFrame,StdOfRest,( this->stdofrest != AST__BADSOR ))
6299 astMAKE_GET(SpecFrame,StdOfRest,AstStdOfRestType,AST__BADSOR,(
6300             ( this->stdofrest == AST__BADSOR ) ? AST__HLSOR : this->stdofrest ) )
6301 
6302 /*
6303 *att++
6304 *  Name:
6305 *     SpecOrigin
6306 
6307 *  Purpose:
6308 *     The zero point for SpecFrame axis values
6309 
6310 *  Type:
6311 *     Public attribute.
6312 
6313 *  Synopsis:
6314 *     Floating point.
6315 
6316 *  Description:
6317 *     This specifies the origin from which all spectral values are measured.
6318 *     The default value (zero) results in the SpecFrame describing
6319 *     absolute spectral values in the system given by the System attribute
6320 *     (e.g. frequency, velocity, etc). If a SpecFrame is to be used to
6321 *     describe offset from some origin, the SpecOrigin attribute
6322 *     should be set to hold the required origin value. The SpecOrigin value
6323 *     stored inside the SpecFrame structure is modified whenever SpecFrame
6324 *     attribute values are changed so that it refers to the original spectral
6325 *     position.
6326 *
6327 *     When setting a new value for this attribute, the supplied value is assumed
6328 *     to be in the system, units and standard of rest described by the SpecFrame.
6329 *     Likewise, when getting the value of this attribute, the value is returned
6330 *     in the system, units and standard of rest described by the SpecFrame. If
6331 *     any of these attributes are changed, then any previously stored SpecOrigin
6332 *     value will also be changed so that refers to the new system, units or
6333 *     standard of rest.
6334 
6335 *  Applicability:
6336 *     SpecFrame
6337 *        All SpecFrames have this attribute.
6338 
6339 *att--
6340 */
6341 /* The spec origin, stored internally in the default units associated
6342    with the current System value. Clear the SpecOrigin value  by setting it
6343    to AST__BAD, which gives 0.0 as the default value. Any value is acceptable. */
6344 astMAKE_CLEAR(SpecFrame,SpecOrigin,specorigin,AST__BAD)
6345 astMAKE_GET(SpecFrame,SpecOrigin,double,0.0,((this->specorigin!=AST__BAD)?this->specorigin:0.0))
6346 astMAKE_SET(SpecFrame,SpecOrigin,double,specorigin,value)
6347 astMAKE_TEST(SpecFrame,SpecOrigin,( this->specorigin != AST__BAD ))
6348 
6349 /* Copy constructor. */
6350 /* ----------------- */
6351 static void Copy( const AstObject *objin, AstObject *objout, int *status ) {
6352 /*
6353 *  Name:
6354 *     Copy
6355 
6356 *  Purpose:
6357 *     Copy constructor for SpecFrame objects.
6358 
6359 *  Type:
6360 *     Private function.
6361 
6362 *  Synopsis:
6363 *     void Copy( const AstObject *objin, AstObject *objout, int *status )
6364 
6365 *  Description:
6366 *     This function implements the copy constructor for SpecFrame objects.
6367 
6368 *  Parameters:
6369 *     objin
6370 *        Pointer to the object to be copied.
6371 *     objout
6372 *        Pointer to the object being constructed.
6373 *     status
6374 *        Pointer to the inherited status variable.
6375 
6376 *  Notes:
6377 *     -  This constructor makes a deep copy.
6378 */
6379 
6380 /* Local Variables: */
6381    AstSpecFrame *in;             /* Pointer to input SpecFrame */
6382    AstSpecFrame *out;            /* Pointer to output SpecFrame */
6383    char *usedunit;               /* Pointer to an element of usedunits array */
6384    int i;                        /* Loop count */
6385    int nused;                    /* Size of "usedunits" array */
6386 
6387 /* Check the global error status. */
6388    if ( !astOK ) return;
6389 
6390 /* Obtain pointers to the input and output SpecFrames. */
6391    in = (AstSpecFrame *) objin;
6392    out = (AstSpecFrame *) objout;
6393 
6394 /* Nullify the pointers stored in the output object since these will
6395    currently be pointing at the input data (since the output is a simple
6396    byte-for-byte copy of the input). Otherwise, the input data could be
6397    freed by accidient if the output object is deleted due to an error
6398    occuring in this function. */
6399    out->usedunits = NULL;
6400 
6401 /* Store the last used units in the output SpecMap. */
6402    if( in && in->usedunits ) {
6403       nused = in->nuunits;
6404       out->usedunits = astMalloc( nused*sizeof( char * ) );
6405       if( out->usedunits ) {
6406          for( i = 0; i < nused; i++ ) {
6407             usedunit = in->usedunits[ i ];
6408             if( usedunit ) {
6409                out->usedunits[ i ] = astStore( NULL, usedunit,
6410                                                strlen( usedunit ) + 1 );
6411             } else {
6412                out->usedunits[ i ] = NULL;
6413             }
6414          }
6415       }
6416    }
6417 
6418 /* If an error has occurred, free the output resources. */
6419    if( !astOK ) Delete( (AstObject *) out, status );
6420 
6421 }
6422 
6423 /* Destructor. */
6424 /* ----------- */
Delete(AstObject * obj,int * status)6425 static void Delete( AstObject *obj, int *status ) {
6426 /*
6427 *  Name:
6428 *     Delete
6429 
6430 *  Purpose:
6431 *     Destructor for SpecFrame objects.
6432 
6433 *  Type:
6434 *     Private function.
6435 
6436 *  Synopsis:
6437 *     void Delete( AstObject *obj, int *status )
6438 
6439 *  Description:
6440 *     This function implements the destructor for SpecFrame objects.
6441 
6442 *  Parameters:
6443 *     obj
6444 *        Pointer to the object to be deleted.
6445 *     status
6446 *        Pointer to the inherited status variable.
6447 
6448 *  Notes:
6449 *     This function attempts to execute even if the global error status is
6450 *     set.
6451 */
6452 
6453 /* Local Variables: */
6454    AstSpecFrame *this;
6455    int i;
6456 
6457 /* Release the memory referred to in the SpecFrame structure. */
6458    this = (AstSpecFrame *) obj;
6459    if( this && this->usedunits ) {
6460       for( i = 0; i < this->nuunits; i++ ) {
6461          this->usedunits[ i ] = astFree( this->usedunits[ i ] );
6462       }
6463       this->usedunits = astFree( this->usedunits );
6464    }
6465 }
6466 
6467 /* Dump function. */
6468 /* -------------- */
Dump(AstObject * this_object,AstChannel * channel,int * status)6469 static void Dump( AstObject *this_object, AstChannel *channel, int *status ) {
6470 /*
6471 *  Name:
6472 *     Dump
6473 
6474 *  Purpose:
6475 *     Dump function for SpecFrame objects.
6476 
6477 *  Type:
6478 *     Private function.
6479 
6480 *  Synopsis:
6481 *     void Dump( AstObject *this, AstChannel *channel, int *status )
6482 
6483 *  Description:
6484 *     This function implements the Dump function which writes out data
6485 *     for the SpecFrame class to an output Channel.
6486 
6487 *  Parameters:
6488 *     this
6489 *        Pointer to the SpecFrame whose data are being written.
6490 *     channel
6491 *        Pointer to the Channel to which the data are being written.
6492 *     status
6493 *        Pointer to the inherited status variable.
6494 */
6495 
6496 /* Local Variables: */
6497    AstSpecFrame *this;           /* Pointer to the SpecFrame structure */
6498    AstStdOfRestType sor;         /* StdOfRest attribute value */
6499    AstSystemType sys;            /* Spectral system value */
6500    char buff[ 20 ];              /* Buffer for item name */
6501    char comm[ 50 ];              /* Buffer for comment */
6502    const char *sval;             /* Pointer to string value */
6503    double dval;                  /* Double value */
6504    int i;                        /* Loop count */
6505    int ival;                     /* int value */
6506    int j;                        /* Loop count */
6507    int set;                      /* Attribute value set? */
6508 
6509 /* Check the global error status. */
6510    if ( !astOK ) return;
6511 
6512 /* Obtain a pointer to the SpecFrame structure. */
6513    this = (AstSpecFrame *) this_object;
6514 
6515 /* Write out values representing the instance variables for the
6516    SpecFrame class.  Accompany these with appropriate comment strings,
6517    possibly depending on the values being written.*/
6518 
6519 /* In the case of attributes, we first use the appropriate (private)
6520    Test...  member function to see if they are set. If so, we then use
6521    the (private) Get... function to obtain the value to be written
6522    out.
6523 
6524    For attributes which are not set, we use the astGet... method to
6525    obtain the value instead. This will supply a default value
6526    (possibly provided by a derived class which over-rides this method)
6527    which is more useful to a human reader as it corresponds to the
6528    actual default attribute value.  Since "set" will be zero, these
6529    values are for information only and will not be read back. */
6530 
6531 /* StdOfRest. */
6532 /* ---------- */
6533    set = TestStdOfRest( this, status );
6534    sor = set ? GetStdOfRest( this, status ) : astGetStdOfRest( this );
6535 
6536 /* If set, convert explicitly to a string for the external
6537    representation. */
6538    sval = "";
6539    if ( set ) {
6540       if ( astOK ) {
6541          sval = StdOfRestString( sor, status );
6542 
6543 /* Report an error if the StdOfRest value was not recognised. */
6544          if ( !sval ) {
6545             astError( AST__SCSIN,
6546                      "%s(%s): Corrupt %s contains invalid standard of rest "
6547                      "identification code (%d).", status, "astWrite",
6548                      astGetClass( channel ), astGetClass( this ), (int) sor );
6549          }
6550       }
6551 
6552 /* If not set, use astGetAttrib which returns a string value using
6553    (possibly over-ridden) methods. */
6554    } else {
6555       sval = astGetAttrib( this_object, "stdofrest" );
6556    }
6557 
6558 /* Write out the value. */
6559    astWriteString( channel, "SoR", set, 1, sval, "Standard of rest" );
6560 
6561 /* AlignStdOfRest. */
6562 /* --------------- */
6563    set = TestAlignStdOfRest( this, status );
6564    sor = set ? GetAlignStdOfRest( this, status ) : astGetAlignStdOfRest( this );
6565 
6566 /* If set, convert explicitly to a string for the external representation. */
6567    if ( set ) {
6568       if ( astOK ) {
6569          sval = StdOfRestString( sor, status );
6570 
6571 /* Report an error if the StdOfRest value was not recognised. */
6572          if ( !sval ) {
6573             astError( AST__SCSIN,
6574                      "%s(%s): Corrupt %s contains invalid alignment standard "
6575                      "of rest identification code (%d).", status, "astWrite",
6576                      astGetClass( channel ), astGetClass( this ), (int) sor );
6577          }
6578       }
6579 
6580 /* If not set, use astGetAttrib which returns a string value using
6581    (possibly over-ridden) methods. */
6582    } else {
6583       sval = astGetAttrib( this_object, "alignstdofrest" );
6584    }
6585 
6586 /* Write out the value. */
6587    astWriteString( channel, "AlSoR", set, 0, sval, "Alignment standard of rest" );
6588 
6589 /* RefRA. */
6590 /* ------ */
6591    set = TestRefRA( this, status );
6592    dval = set ? GetRefRA( this, status ) : astGetRefRA( this );
6593    astWriteDouble( channel, "RefRA", set, 0, dval, "Reference RA (rads, FK5 J2000)" );
6594 
6595 /* RefDec. */
6596 /* ------- */
6597    set = TestRefDec( this, status );
6598    dval = set ? GetRefDec( this, status ) : astGetRefDec( this );
6599    astWriteDouble( channel, "RefDec", set, 0, dval, "Reference Dec (rads, FK5 J2000)" );
6600 
6601 /* RestFreq. */
6602 /* --------- */
6603    set = TestRestFreq( this, status );
6604    dval = set ? GetRestFreq( this, status ) : astGetRestFreq( this );
6605    astWriteDouble( channel, "RstFrq", set, 0, dval, "Rest frequency (Hz)" );
6606 
6607 /* SourceVel. */
6608 /* ---------- */
6609    set = TestSourceVel( this, status );
6610    dval = set ? GetSourceVel( this, status ) : astGetSourceVel( this );
6611    astWriteDouble( channel, "SrcVel", set, 0, dval, "Source velocity (m/s)" );
6612 
6613 /* SourceVRF. */
6614 /* ---------- */
6615    set = TestSourceVRF( this, status );
6616    sor = set ? GetSourceVRF( this, status ) : astGetSourceVRF( this );
6617 
6618 /* If set, convert explicitly to a string for the external representation. */
6619    if ( set ) {
6620       if ( astOK ) {
6621          sval = StdOfRestString( sor, status );
6622 
6623 /* Report an error if the value was not recognised. */
6624          if ( !sval ) {
6625             astError( AST__SCSIN,
6626                      "%s(%s): Corrupt %s contains invalid source velocity "
6627                      "rest frame identification code (%d).", status, "astWrite",
6628                      astGetClass( channel ), astGetClass( this ), (int) sor );
6629          }
6630       }
6631 
6632 /* If not set, use astGetAttrib which returns a string value using
6633    (possibly over-ridden) methods. */
6634    } else {
6635       sval = astGetAttrib( this_object, "sourcevrf" );
6636    }
6637 
6638 /* Write out the value. */
6639    astWriteString( channel, "SrcVRF", set, 0, sval, "Source velocity rest frame" );
6640 
6641 /* SourceSys. */
6642 /* ---------- */
6643    set = TestSourceSys( this, status );
6644    sys = set ? GetSourceSys( this, status ) : astGetSourceSys( this );
6645 
6646 /* If set, convert explicitly to a string for the external representation. */
6647    if ( set ) {
6648       if ( astOK ) {
6649          sval = SystemString( (AstFrame *) this, sys, status );
6650 
6651 /* Report an error if the value was not recognised. */
6652          if ( !sval ) {
6653             astError( AST__SCSIN,
6654                      "%s(%s): Corrupt %s contains invalid source velocity "
6655                      "spectral system identification code (%d).", status, "astWrite",
6656                      astGetClass( channel ), astGetClass( this ), (int) sys );
6657          }
6658       }
6659 
6660 /* If not set, use astGetAttrib which returns a string value using
6661    (possibly over-ridden) methods. */
6662    } else {
6663       sval = astGetAttrib( this_object, "sourcesys" );
6664    }
6665 
6666 /* Write out the value. */
6667    astWriteString( channel, "SrcSys", set, 0, sval, "Source velocity spectral system" );
6668 
6669 /* AlignSpecOffset. */
6670 /* ---------------- */
6671    set = TestAlignSpecOffset( this, status );
6672    ival = set ? GetAlignSpecOffset( this, status ) : astGetAlignSpecOffset( this );
6673    astWriteInt( channel, "AlSpOf", set, 0, ival,
6674                 ival ? "Align in offset coords" :
6675                        "Align in system coords" );
6676 
6677 /* UsedUnits */
6678 /* --------- */
6679    if( this->usedunits ) {
6680       for( i = 0; i < this->nuunits; i++ ) {
6681          if( this->usedunits[ i ] ) {
6682             sprintf( buff, "U%s", astSystemString( this, (AstSystemType) i ));
6683             for( j = 2; j < strlen( buff ); j++ ) buff[ j ] = tolower( buff[ j ] );
6684             sprintf( comm, "Preferred units for %s", SystemLabel( (AstSystemType) i, status ) );
6685             astWriteString( channel, buff, 1, 0, this->usedunits[ i ], comm );
6686          }
6687       }
6688    }
6689 
6690 /* SpecOrigin. */
6691 /* ----------- */
6692    set = TestSpecOrigin( this, status );
6693    dval = set ? GetSpecOrigin( this, status ) : astGetSpecOrigin( this );
6694    if( dval != AST__BAD ) {
6695       astWriteDouble( channel, "SpOrg", set, 0, dval, "Spec offset" );
6696    }
6697 
6698 }
6699 
6700 /* Standard class functions. */
6701 /* ========================= */
6702 /* Implement the astIsASpecFrame and astCheckSpecFrame functions using the
6703    macros defined for this purpose in the "object.h" header file. */
astMAKE_ISA(SpecFrame,Frame)6704 astMAKE_ISA(SpecFrame,Frame)
6705 astMAKE_CHECK(SpecFrame)
6706 
6707 AstSpecFrame *astSpecFrame_( const char *options, int *status, ...) {
6708 /*
6709 *+
6710 *  Name:
6711 *     astSpecFrame
6712 
6713 *  Purpose:
6714 *     Create a SpecFrame.
6715 
6716 *  Type:
6717 *     Protected function.
6718 
6719 *  Synopsis:
6720 *     #include "specframe.h"
6721 *     AstSpecFrame *astSpecFrame( const char *options, int *status, ... )
6722 
6723 *  Class Membership:
6724 *     SpecFrame constructor.
6725 
6726 *  Description:
6727 *     This function creates a new SpecFrame and optionally initialises its
6728 *     attributes.
6729 
6730 *  Parameters:
6731 *     options
6732 *        Pointer to a null terminated string containing an optional
6733 *        comma-separated list of attribute assignments to be used for
6734 *        initialising the new SpecFrame. The syntax used is the same as for the
6735 *        astSet method and may include "printf" format specifiers identified
6736 *        by "%" symbols in the normal way.
6737 *     status
6738 *        Pointer to the inherited status variable.
6739 *     ...
6740 *        If the "options" string contains "%" format specifiers, then an
6741 *        optional list of arguments may follow it in order to supply values to
6742 *        be substituted for these specifiers. The rules for supplying these
6743 *        are identical to those for the astSet method (and for the C "printf"
6744 *        function).
6745 
6746 *  Returned Value:
6747 *     A pointer to the new SpecFrame.
6748 
6749 *  Notes:
6750 *     -  A NULL pointer will be returned if this function is invoked with the
6751 *     global error status set, or if it should fail for any reason.
6752 *-
6753 
6754 *  Implementation Notes:
6755 *     - This function implements the basic SpecFrame constructor which
6756 *     is available via the protected interface to the SpecFrame class.
6757 *     A public interface is provided by the astSpecFrameId_ function.
6758 */
6759 
6760 /* Local Variables: */
6761    astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
6762    AstMapping *um;               /* Mapping from default to actual units */
6763    AstSpecFrame *new;            /* Pointer to new SpecFrame */
6764    AstSystemType s;              /* System */
6765    const char *u;                /* Units string */
6766    va_list args;                 /* Variable argument list */
6767 
6768 /* Get a pointer to the thread specific global data structure. */
6769    astGET_GLOBALS(NULL);
6770 
6771 /* Check the global status. */
6772    if ( !astOK ) return NULL;
6773 
6774 /* Initialise the SpecFrame, allocating memory and initialising the virtual
6775    function table as well if necessary. */
6776    new = astInitSpecFrame( NULL, sizeof( AstSpecFrame ), !class_init,
6777                            &class_vtab, "SpecFrame" );
6778 
6779 /* If successful, note that the virtual function table has been initialised. */
6780    if ( astOK ) {
6781       class_init = 1;
6782 
6783 /* Obtain the variable argument list and pass it along with the options string
6784    to the astVSet method to initialise the new SpecFrame's attributes. */
6785       va_start( args, status );
6786       astVSet( new, options, NULL, args );
6787       va_end( args );
6788 
6789 /* Check the Units are appropriate for the System. */
6790       u = astGetUnit( new, 0 );
6791       s = astGetSystem( new );
6792       um = astUnitMapper( DefUnit( s, "astSpecFrame", "SpecFrame", status ),
6793                           u, NULL, NULL );
6794       if( um ) {
6795          um = astAnnul( um );
6796       } else {
6797          astError( AST__BADUN, "astSpecFrame: Inappropriate units (%s) "
6798                    "specified for a %s axis.", status, u, SystemLabel( s, status ) );
6799       }
6800 
6801 /* If an error occurred, clean up by deleting the new object. */
6802       if ( !astOK ) new = astDelete( new );
6803    }
6804 
6805 /* Return a pointer to the new SpecFrame. */
6806    return new;
6807 }
6808 
astInitSpecFrame_(void * mem,size_t size,int init,AstSpecFrameVtab * vtab,const char * name,int * status)6809 AstSpecFrame *astInitSpecFrame_( void *mem, size_t size, int init,
6810                                  AstSpecFrameVtab *vtab, const char *name, int *status ) {
6811 /*
6812 *+
6813 *  Name:
6814 *     astInitSpecFrame
6815 
6816 *  Purpose:
6817 *     Initialise a SpecFrame.
6818 
6819 *  Type:
6820 *     Protected function.
6821 
6822 *  Synopsis:
6823 *     #include "specframe.h"
6824 *     AstSpecFrame *astInitSpecFrame( void *mem, size_t size, int init,
6825 *                                     AstFrameVtab *vtab, const char *name )
6826 
6827 *  Class Membership:
6828 *     SpecFrame initialiser.
6829 
6830 *  Description:
6831 *     This function is provided for use by class implementations to
6832 *     initialise a new SpecFrame object. It allocates memory (if
6833 *     necessary) to accommodate the SpecFrame plus any additional data
6834 *     associated with the derived class. It then initialises a
6835 *     SpecFrame structure at the start of this memory. If the "init"
6836 *     flag is set, it also initialises the contents of a virtual function
6837 *     table for a SpecFrame at the start of the memory passed via the
6838 *     "vtab" parameter.
6839 
6840 *  Parameters:
6841 *     mem
6842 *        A pointer to the memory in which the SpecFrame is to be
6843 *	created. This must be of sufficient size to accommodate the
6844 *	SpecFrame data (sizeof(SpecFrame)) plus any data used by
6845 *	the derived class. If a value of NULL is given, this function
6846 *	will allocate the memory itself using the "size" parameter to
6847 *	determine its size.
6848 *     size
6849 *        The amount of memory used by the SpecFrame (plus derived
6850 *	class data). This will be used to allocate memory if a value of
6851 *	NULL is given for the "mem" parameter. This value is also stored
6852 *	in the SpecFrame structure, so a valid value must be supplied
6853 *	even if not required for allocating memory.
6854 *     init
6855 *        A logical flag indicating if the SpecFrame's virtual function
6856 *	table is to be initialised. If this value is non-zero, the
6857 *	virtual function table will be initialised by this function.
6858 *     vtab
6859 *        Pointer to the start of the virtual function table to be
6860 *	associated with the new SpecFrame.
6861 *     name
6862 *        Pointer to a constant null-terminated character string which
6863 *	contains the name of the class to which the new object belongs
6864 *	(it is this pointer value that will subsequently be returned by
6865 *	the astGetClass method).
6866 
6867 *  Returned Value:
6868 *     A pointer to the new SpecFrame.
6869 
6870 *  Notes:
6871 *     -  A null pointer will be returned if this function is invoked with the
6872 *     global error status set, or if it should fail for any reason.
6873 *-
6874 */
6875 
6876 /* Local Variables: */
6877    AstSpecFrame *new;        /* Pointer to the new SpecFrame */
6878 
6879 /* Check the global status. */
6880    if ( !astOK ) return NULL;
6881 
6882 /* If necessary, initialise the virtual function table. */
6883    if ( init ) astInitSpecFrameVtab( vtab, name );
6884 
6885 /* Initialise a 1D Frame structure (the parent class) as the first component
6886    within the SpecFrame structure, allocating memory if necessary. */
6887    new = (AstSpecFrame *) astInitFrame( mem, size, 0,
6888                                         (AstFrameVtab *) vtab, name, 1 );
6889 
6890    if ( astOK ) {
6891 
6892 /* Initialise the SpecFrame data. */
6893 /* ----------------------------- */
6894 /* Initialise all attributes to their "undefined" values. */
6895       new->alignstdofrest = AST__BADSOR;
6896       new->refdec = AST__BAD;
6897       new->refra = AST__BAD;
6898       new->restfreq = AST__BAD;
6899       new->sourcevel = AST__BAD;
6900       new->sourcevrf = AST__BADSOR;
6901       new->sourcesys = AST__BADSYSTEM;
6902       new->stdofrest = AST__BADSOR;
6903       new->nuunits = 0;
6904       new->usedunits = NULL;
6905       new->specorigin = AST__BAD;
6906       new->alignspecoffset = -INT_MAX;
6907 
6908 /* If an error occurred, clean up by deleting the new object. */
6909       if ( !astOK ) new = astDelete( new );
6910 
6911    }
6912 
6913 /* Return a pointer to the new object. */
6914    return new;
6915 }
6916 
astLoadSpecFrame_(void * mem,size_t size,AstSpecFrameVtab * vtab,const char * name,AstChannel * channel,int * status)6917 AstSpecFrame *astLoadSpecFrame_( void *mem, size_t size,
6918                                  AstSpecFrameVtab *vtab,
6919                                  const char *name, AstChannel *channel, int *status ) {
6920 /*
6921 *+
6922 *  Name:
6923 *     astLoadSpecFrame
6924 
6925 *  Purpose:
6926 *     Load a SpecFrame.
6927 
6928 *  Type:
6929 *     Protected function.
6930 
6931 *  Synopsis:
6932 *     #include "specframe.h"
6933 *     AstSpecFrame *astLoadSpecFrame( void *mem, size_t size,
6934 *                                      AstSpecFrameVtab *vtab,
6935 *                                      const char *name, AstChannel *channel )
6936 
6937 *  Class Membership:
6938 *     SpecFrame loader.
6939 
6940 *  Description:
6941 *     This function is provided to load a new SpecFrame using data read
6942 *     from a Channel. It first loads the data used by the parent class
6943 *     (which allocates memory if necessary) and then initialises a
6944 *     SpecFrame structure in this memory, using data read from the
6945 *     input Channel.
6946 
6947 *  Parameters:
6948 *     mem
6949 *        A pointer to the memory into which the SpecFrame is to be
6950 *        loaded.  This must be of sufficient size to accommodate the
6951 *        SpecFrame data (sizeof(SpecFrame)) plus any data used by
6952 *        derived classes. If a value of NULL is given, this function
6953 *        will allocate the memory itself using the "size" parameter to
6954 *        determine its size.
6955 *     size
6956 *        The amount of memory used by the SpecFrame (plus derived class
6957 *        data).  This will be used to allocate memory if a value of
6958 *        NULL is given for the "mem" parameter. This value is also
6959 *        stored in the SpecFrame structure, so a valid value must be
6960 *        supplied even if not required for allocating memory.
6961 *
6962 *        If the "vtab" parameter is NULL, the "size" value is ignored
6963 *        and sizeof(AstSpecFrame) is used instead.
6964 *     vtab
6965 *        Pointer to the start of the virtual function table to be
6966 *        associated with the new SpecFrame. If this is NULL, a pointer
6967 *        to the (static) virtual function table for the SpecFrame class
6968 *        is used instead.
6969 *     name
6970 *        Pointer to a constant null-terminated character string which
6971 *        contains the name of the class to which the new object
6972 *        belongs (it is this pointer value that will subsequently be
6973 *        returned by the astGetClass method).
6974 *
6975 *        If the "vtab" parameter is NULL, the "name" value is ignored
6976 *        and a pointer to the string "SpecFrame" is used instead.
6977 
6978 *  Returned Value:
6979 *     A pointer to the new SpecFrame.
6980 
6981 *  Notes:
6982 *     - A null pointer will be returned if this function is invoked
6983 *     with the global error status set, or if it should fail for any
6984 *     reason.
6985 *-
6986 */
6987 
6988 /* Local Variables: */
6989    astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
6990    AstSpecFrame *new;            /* Pointer to the new SpecFrame */
6991    char buff[ 20 ];              /* Buffer for item name */
6992    char *sval;                   /* Pointer to string value */
6993    double obslat;                /* Value for ObsLat attribute */
6994    double obslon;                /* Get a pointer to the thread specific global data structure. */
6995 
6996 /* Value for ObsLon attribute */
6997    int i;                        /* Loop count */
6998    int j;                        /* Loop count */
6999    int nc;                       /* String length */
7000    int sys;                      /* System value */
7001 
7002    astGET_GLOBALS(channel);
7003 
7004 /* Initialise. */
7005    new = NULL;
7006 
7007 /* Check the global error status. */
7008    if ( !astOK ) return new;
7009 
7010 /* If a NULL virtual function table has been supplied, then this is
7011    the first loader to be invoked for this SpecFrame. In this case the
7012    SpecFrame belongs to this class, so supply appropriate values to be
7013    passed to the parent class loader (and its parent, etc.). */
7014    if ( !vtab ) {
7015       size = sizeof( AstSpecFrame );
7016       vtab = &class_vtab;
7017       name = "SpecFrame";
7018 
7019 /* If required, initialise the virtual function table for this class. */
7020       if ( !class_init ) {
7021          astInitSpecFrameVtab( vtab, name );
7022          class_init = 1;
7023       }
7024    }
7025 
7026 /* Invoke the parent class loader to load data for all the ancestral
7027    classes of the current one, returning a pointer to the resulting
7028    partly-built SpecFrame. */
7029    new = astLoadFrame( mem, size, (AstFrameVtab *) vtab, name,
7030                        channel );
7031    if ( astOK ) {
7032 
7033 /* Read input data. */
7034 /* ================ */
7035 /* Request the input Channel to read all the input data appropriate to
7036    this class into the internal "values list". */
7037        astReadClassData( channel, "SpecFrame" );
7038 
7039 /* Now read each individual data item from this list and use it to
7040    initialise the appropriate instance variable(s) for this class. */
7041 
7042 /* In the case of attributes, we first read the "raw" input value,
7043    supplying the "unset" value as the default. If a "set" value is
7044    obtained, we then use the appropriate (private) Set... member
7045    function to validate and set the value properly. */
7046 
7047 /* StdOfRest. */
7048 /* ---------- */
7049 /* Set the default and read the external representation as a string. */
7050        new->stdofrest = AST__BADSOR;
7051        sval = astReadString( channel, "sor", NULL );
7052 
7053 /* If a value was read, convert from a string to a StdOfRest code. */
7054        if ( sval ) {
7055           if ( astOK ) {
7056              new->stdofrest = StdOfRestCode( sval, status );
7057 
7058 /* Report an error if the value wasn't recognised. */
7059              if ( new->stdofrest == AST__BADSOR ) {
7060                 astError( AST__ATTIN,
7061                           "astRead(%s): Invalid standard of rest description "
7062                           "\"%s\".", status, astGetClass( channel ), sval );
7063              }
7064           }
7065 
7066 /* Free the string value. */
7067           sval = astFree( sval );
7068        }
7069 
7070 /* AlignStdOfRest. */
7071 /* --------------- */
7072 /* Set the default and read the external representation as a string. */
7073        new->alignstdofrest = AST__BADSOR;
7074        sval = astReadString( channel, "alsor", NULL );
7075 
7076 /* If a value was read, convert from a string to a StdOfRest code. */
7077        if ( sval ) {
7078           if ( astOK ) {
7079              new->alignstdofrest = StdOfRestCode( sval, status );
7080 
7081 /* Report an error if the value wasn't recognised. */
7082              if ( new->alignstdofrest == AST__BADSOR ) {
7083                 astError( AST__ATTIN,
7084                           "astRead(%s): Invalid alignment standard of rest "
7085                           "description \"%s\".", status, astGetClass( channel ), sval );
7086              }
7087           }
7088 
7089 /* Free the string value. */
7090           sval = astFree( sval );
7091        }
7092 
7093 /* GeoLat. */
7094 /* ------- */
7095 /* Retained for backward compatibility with older versions of AST in
7096    which SpecFrame had a GeoLat attribute (now ObsLat is used instead). */
7097       if( !astTestObsLat( new ) ) {
7098          obslat = astReadDouble( channel, "geolat", AST__BAD );
7099          if ( obslat != AST__BAD ) astSetObsLat( new, obslat );
7100       }
7101 
7102 /* GeoLon. */
7103 /* ------- */
7104 /* Retained for backward compatibility with older versions of AST in
7105    which SpecFrame had a GeoLon attribute (now ObsLon is used instead). */
7106       if( !astTestObsLon( new ) ) {
7107          obslon = astReadDouble( channel, "geolon", AST__BAD );
7108          if ( obslon != AST__BAD ) astSetObsLon( new, obslon );
7109       }
7110 
7111 /* RefRA. */
7112 /* ------ */
7113       new->refra = astReadDouble( channel, "refra", AST__BAD );
7114       if ( TestRefRA( new, status ) ) SetRefRA( new, new->refra, status );
7115 
7116 /* RefDec. */
7117 /* ------- */
7118       new->refdec = astReadDouble( channel, "refdec", AST__BAD );
7119       if ( TestRefDec( new, status ) ) SetRefDec( new, new->refdec, status );
7120 
7121 /* RestFreq. */
7122 /* --------- */
7123       new->restfreq = astReadDouble( channel, "rstfrq", AST__BAD );
7124       if ( TestRestFreq( new, status ) ) SetRestFreq( new, new->restfreq, status );
7125 
7126 /* AlignSpecOffset */
7127 /* --------------- */
7128       new->alignspecoffset = astReadInt( channel, "alspof", -INT_MAX );
7129       if ( TestAlignSpecOffset( new, status ) ) SetAlignSpecOffset( new, new->alignspecoffset, status );
7130 
7131 /* SourceVel. */
7132 /* ---------- */
7133       new->sourcevel = astReadDouble( channel, "srcvel", AST__BAD );
7134       if ( TestSourceVel( new, status ) ) SetSourceVel( new, new->sourcevel, status );
7135 
7136 /* SourceVRF */
7137 /* --------- */
7138 /* Set the default and read the external representation as a string. */
7139        new->sourcevrf = AST__BADSOR;
7140        sval = astReadString( channel, "srcvrf", NULL );
7141 
7142 /* If a value was read, convert from a string to a StdOfRest code. */
7143        if ( sval ) {
7144           if ( astOK ) {
7145              new->sourcevrf = StdOfRestCode( sval, status );
7146 
7147 /* Report an error if the value wasn't recognised. */
7148              if ( new->sourcevrf == AST__BADSOR ) {
7149                 astError( AST__ATTIN,
7150                           "astRead(%s): Invalid source velocity rest frame "
7151                           "description \"%s\".", status, astGetClass( channel ), sval );
7152              }
7153           }
7154 
7155 /* Free the string value. */
7156           sval = astFree( sval );
7157        }
7158 
7159 /* SourceSys */
7160 /* --------- */
7161 /* Set the default and read the external representation as a string. */
7162        new->sourcesys = AST__BADSYSTEM;
7163        sval = astReadString( channel, "srcsys", NULL );
7164 
7165 /* If a value was read, convert from a string to a System code. */
7166        if ( sval ) {
7167           if ( astOK ) {
7168              new->sourcesys = SystemCode( (AstFrame *) new, sval, status );
7169 
7170 /* Report an error if the value wasn't recognised. */
7171              if ( new->sourcesys == AST__BADSYSTEM ) {
7172                 astError( AST__ATTIN,
7173                           "astRead(%s): Invalid source velocity spectral system "
7174                           "description \"%s\".", status, astGetClass( channel ), sval );
7175              }
7176           }
7177 
7178 /* Free the string value. */
7179           sval = astFree( sval );
7180        }
7181 
7182 /* UsedUnits */
7183 /* --------- */
7184       new->nuunits = 0;
7185       new->usedunits = NULL;
7186       for( sys = FIRST_SYSTEM; sys <= LAST_SYSTEM; sys++ ) {
7187          nc = sprintf( buff, "u%s", astSystemString( new, (AstSystemType) sys ));
7188          for( j = 0; j < nc; j++ ) buff[ j ] = tolower( buff[ j ] );
7189          sval = astReadString( channel, buff, NULL );
7190          if( sval ) {
7191             if( (int) sys >= new->nuunits ) {
7192                new->usedunits = astGrow( new->usedunits, sys + 1,
7193                                           sizeof(char *) );
7194                if( astOK ) {
7195                   for( i = new->nuunits; i < sys + 1; i++ ) new->usedunits[ i ] = NULL;
7196                   new->nuunits = sys + 1;
7197                }
7198             } else {
7199                new->usedunits[ sys ] = astFree( new->usedunits[ sys ] );
7200             }
7201             if( astOK ) {
7202                new->usedunits[ sys ] = astStore( new->usedunits[ sys ],
7203                                                  sval, strlen( sval ) + 1 );
7204             }
7205             sval = astFree( sval);
7206          }
7207       }
7208 
7209 /* SpecOrigin. */
7210 /* --------- */
7211       new->specorigin = astReadDouble( channel, "sporg", AST__BAD );
7212       if ( TestSpecOrigin( new, status ) ) SetSpecOrigin( new, new->specorigin, status );
7213 
7214 
7215 /* If an error occurred, clean up by deleting the new SpecFrame. */
7216        if ( !astOK ) new = astDelete( new );
7217    }
7218 
7219 /* Return the new SpecFrame pointer. */
7220    return new;
7221 }
7222 
7223 /* Virtual function interfaces. */
7224 /* ============================ */
7225 /* These provide the external interface to the virtual functions defined by
7226    this class. Each simply checks the global error status and then locates and
7227    executes the appropriate member function, using the function pointer stored
7228    in the object's virtual function table (this pointer is located using the
7229    astMEMBER macro defined in "object.h").
7230 
7231    Note that the member function may not be the one defined here, as it may
7232    have been over-ridden by a derived class. However, it should still have the
7233    same interface. */
astGetRefPos_(AstSpecFrame * this,AstSkyFrame * frm,double * lon,double * lat,int * status)7234 void astGetRefPos_( AstSpecFrame *this, AstSkyFrame *frm, double *lon,
7235                     double *lat, int *status ){
7236    if ( !astOK ) return;
7237    (**astMEMBER(this,SpecFrame,GetRefPos))(this,frm,lon,lat, status );
7238 }
astSetRefPos_(AstSpecFrame * this,AstSkyFrame * frm,double lon,double lat,int * status)7239 void astSetRefPos_( AstSpecFrame *this, AstSkyFrame *frm, double lon,
7240                     double lat, int *status ){
7241    if ( !astOK ) return;
7242    (**astMEMBER(this,SpecFrame,SetRefPos))(this,frm,lon,lat, status );
7243 }
7244 
astSetStdOfRest_(AstSpecFrame * this,AstStdOfRestType value,int * status)7245 void astSetStdOfRest_( AstSpecFrame *this, AstStdOfRestType value, int *status ){
7246    if ( !astOK ) return;
7247    (**astMEMBER(this,SpecFrame,SetStdOfRest))(this,value, status );
7248 }
7249 
astClearStdOfRest_(AstSpecFrame * this,int * status)7250 void astClearStdOfRest_( AstSpecFrame *this, int *status ){
7251    if ( !astOK ) return;
7252    (**astMEMBER(this,SpecFrame,ClearStdOfRest))(this, status );
7253 }
7254 
7255 
7256 
7257 /* Special public interface functions. */
7258 /* =================================== */
7259 /* These provide the public interface to certain special functions
7260    whose public interface cannot be handled using macros (such as
7261    astINVOKE) alone. In general, they are named after the
7262    corresponding protected version of the function, but with "Id"
7263    appended to the name. */
7264 
7265 /* Public Interface Function Prototypes. */
7266 /* ------------------------------------- */
7267 /* The following functions have public prototypes only (i.e. no
7268    protected prototypes), so we must provide local prototypes for use
7269    within this module. */
7270 AstSpecFrame *astSpecFrameId_( const char *, ... );
7271 
7272 /* Special interface function implementations. */
7273 /* ------------------------------------------- */
astSpecFrameId_(const char * options,...)7274 AstSpecFrame *astSpecFrameId_( const char *options, ... ) {
7275 /*
7276 *++
7277 *  Name:
7278 c     astSpecFrame
7279 f     AST_SPECFRAME
7280 
7281 *  Purpose:
7282 *     Create a SpecFrame.
7283 
7284 *  Type:
7285 *     Public function.
7286 
7287 *  Synopsis:
7288 c     #include "specframe.h"
7289 c     AstSpecFrame *astSpecFrame( const char *options, ... )
7290 f     RESULT = AST_SPECFRAME( OPTIONS, STATUS )
7291 
7292 *  Class Membership:
7293 *     SpecFrame constructor.
7294 
7295 *  Description:
7296 *     This function creates a new SpecFrame and optionally initialises
7297 *     its attributes.
7298 *
7299 *     A SpecFrame is a specialised form of one-dimensional Frame which
7300 *     represents various coordinate systems used to describe positions within
7301 *     an electro-magnetic spectrum. The particular coordinate system to be
7302 *     used is specified by setting the SpecFrame's System attribute (the
7303 *     default is wavelength) qualified, as necessary, by other attributes
7304 *     such as the rest frequency, the standard of rest, the epoch of
7305 *     observation, etc (see the description of the System attribute for
7306 *     details).
7307 *
7308 *     By setting a value for thr SpecOrigin attribute, a SpecFrame can be made
7309 *     to represent offsets from a given spectral position, rather than absolute
7310 
7311 *  Parameters:
7312 c     options
7313 f     OPTIONS = CHARACTER * ( * ) (Given)
7314 c        Pointer to a null-terminated string containing an optional
7315 c        comma-separated list of attribute assignments to be used for
7316 c        initialising the new SpecFrame. The syntax used is identical to
7317 c        that for the astSet function and may include "printf" format
7318 c        specifiers identified by "%" symbols in the normal way.
7319 c        If no initialisation is required, a zero-length string may be
7320 c        supplied.
7321 f        A character string containing an optional comma-separated
7322 f        list of attribute assignments to be used for initialising the
7323 f        new SpecFrame. The syntax used is identical to that for the
7324 f        AST_SET routine. If no initialisation is required, a blank
7325 f        value may be supplied.
7326 c     ...
7327 c        If the "options" string contains "%" format specifiers, then
7328 c        an optional list of additional arguments may follow it in
7329 c        order to supply values to be substituted for these
7330 c        specifiers. The rules for supplying these are identical to
7331 c        those for the astSet function (and for the C "printf"
7332 c        function).
7333 f     STATUS = INTEGER (Given and Returned)
7334 f        The global status.
7335 
7336 *  Returned Value:
7337 c     astSpecFrame()
7338 f     AST_SPECFRAME = INTEGER
7339 *        A pointer to the new SpecFrame.
7340 
7341 *  Examples:
7342 c     frame = astSpecFrame( "" );
7343 f     FRAME = AST_SPECFRAME( ' ', STATUS )
7344 *        Creates a SpecFrame to describe the default wavelength spectral
7345 *        coordinate system. The RestFreq attribute (rest frequency) is
7346 *        unspecified, so it will not be possible to align this SpecFrame
7347 *        with another SpecFrame on the basis of a velocity-based system. The
7348 *        standard of rest is also unspecified. This means that alignment
7349 *        will be possible with other SpecFrames, but no correction will be
7350 *        made for Doppler shift caused by change of rest frame during the
7351 *        alignment.
7352 c     frame = astSpecFrame( "System=VELO, RestFreq=1.0E15, StdOfRest=LSRK" );
7353 f     FRAME = AST_SPECFRAME( 'System=VELO, RestFreq=1.0E15, StdOfRest=LSRK', STATUS )
7354 *        Creates a SpecFrame describing a apparent radial velocity ("VELO") axis
7355 *        with rest frequency 1.0E15 Hz (about 3000 Angstroms), measured
7356 *        in the kinematic Local Standard of Rest ("LSRK"). Since the
7357 *        source position has not been specified (using attributes RefRA and
7358 *        RefDec), it will only be possible to align this SpecFrame with
7359 *        other SpecFrames which are also measured in the LSRK standard of
7360 *        rest.
7361 
7362 *  Notes:
7363 *     - When conversion between two SpecFrames is requested (as when
7364 c     supplying SpecFrames to astConvert),
7365 f     supplying SpecFrames AST_CONVERT),
7366 *     account will be taken of the nature of the spectral coordinate systems
7367 *     they represent, together with any qualifying rest frequency, standard
7368 *     of rest, epoch values, etc. The AlignSystem and AlignStdOfRest
7369 *     attributes will also be taken into account. The results will therefore
7370 *     fully reflect the relationship between positions measured in the two
7371 *     systems. In addition, any difference in the Unit attributes of the two
7372 *     systems will also be taken into account.
7373 *     - A null Object pointer (AST__NULL) will be returned if this
7374 c     function is invoked with the AST error status set, or if it
7375 f     function is invoked with STATUS set to an error value, or if it
7376 *     should fail for any reason.
7377 *--
7378 
7379 *  Implementation Notes:
7380 *     - This function implements the external (public) interface to
7381 *     the astSpecFrame constructor function. It returns an ID value
7382 *     (instead of a true C pointer) to external users, and must be
7383 *     provided because astSpecFrame_ has a variable argument list which
7384 *     cannot be encapsulated in a macro (where this conversion would
7385 *     otherwise occur).
7386 *     - The variable argument list also prevents this function from
7387 *     invoking astSpecFrame_ directly, so it must be a
7388 *     re-implementation of it in all respects, except for the final
7389 *     conversion of the result to an ID value.
7390 */
7391 
7392 /* Local Variables: */
7393    astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
7394    AstMapping *um;               /* Mapping from default to actual units */
7395    AstSpecFrame *new;            /* Pointer to new SpecFrame */
7396    AstSystemType s;              /* System */
7397    const char *u;                /* Units string */
7398    va_list args;                 /* Variable argument list */
7399 
7400    int *status;                  /* Pointer to inherited status value */
7401    astGET_GLOBALS(NULL);         /* Get a pointer to the thread specific global data structure. */
7402 
7403 /* Get a pointer to the inherited status value. */
7404    status = astGetStatusPtr;
7405 
7406 /* Check the global status. */
7407    if ( !astOK ) return NULL;
7408 
7409 /* Initialise the SpecFrame, allocating memory and initialising the virtual
7410    function table as well if necessary. */
7411    new = astInitSpecFrame( NULL, sizeof( AstSpecFrame ), !class_init,
7412                            &class_vtab, "SpecFrame" );
7413 
7414 /* If successful, note that the virtual function table has been initialised. */
7415    if ( astOK ) {
7416       class_init = 1;
7417 
7418 /* Obtain the variable argument list and pass it along with the options string
7419    to the astVSet method to initialise the new SpecFrame's attributes. */
7420       va_start( args, options );
7421       astVSet( new, options, NULL, args );
7422       va_end( args );
7423 
7424 /* Check the Units are appropriate for the System. */
7425       u = astGetUnit( new, 0 );
7426       s = astGetSystem( new );
7427       um = astUnitMapper( DefUnit( s, "astSpecFrame", "SpecFrame", status ),
7428                           u, NULL, NULL );
7429       if( um ) {
7430          um = astAnnul( um );
7431       } else {
7432          astError( AST__BADUN, "astSpecFrame: Inappropriate units (%s) "
7433                    "specified for a %s axis.", status, u, SystemLabel( s, status ) );
7434       }
7435 
7436 /* If an error occurred, clean up by deleting the new object. */
7437       if ( !astOK ) new = astDelete( new );
7438    }
7439 
7440 /* Return an ID value for the new SpecFrame. */
7441    return astMakeId( new );
7442 }
7443 
7444