1 /*
2 *class++
3 *  Name:
4 *     DSBSpecFrame
5 
6 *  Purpose:
7 *     Dual sideband spectral coordinate system description.
8 
9 *  Constructor Function:
10 c     astDSBSpecFrame
11 f     AST_DSBSPECFRAME
12 
13 *  Description:
14 *     A DSBSpecFrame is a specialised form of SpecFrame which represents
15 *     positions in a spectrum obtained using a dual sideband instrument.
16 *     Such an instrument produces a spectrum in which each point contains
17 *     contributions from two distinctly different frequencies, one from
18 *     the "lower side band" (LSB) and one from the "upper side band" (USB).
19 *     Corresponding LSB and USB frequencies are connected by the fact
20 *     that they are an equal distance on either side of a fixed central
21 *     frequency known as the "Local Oscillator" (LO) frequency.
22 *
23 *     When quoting a position within such a spectrum, it is necessary to
24 *     indicate whether the quoted position is the USB position or the
25 *     corresponding LSB position. The SideBand attribute provides this
26 *     indication. Another option that the SideBand attribute provides is
27 *     to represent a spectral position by its topocentric offset from the
28 *     LO frequency.
29 *
30 *     In practice, the LO frequency is specified by giving the distance
31 *     from the LO frequency to some "central" spectral position. Typically
32 *     this central position is that of some interesting spectral feature.
33 *     The distance from this central position to the LO frequency is known
34 *     as the "intermediate frequency" (IF). The value supplied for IF can
35 *     be a signed value in order to indicate whether the LO frequency is
36 *     above or below the central position.
37 
38 *  Inheritance:
39 *     The DSBSpecFrame class inherits from the SpecFrame class.
40 
41 *  Attributes:
42 *     In addition to those attributes common to all SpecFrames, every
43 *     DSBSpecFrame also has the following attributes:
44 *
45 *     - AlignSideBand: Should alignment occur between sidebands?
46 *     - DSBCentre: The central position of interest.
47 *     - IF: The intermediate frequency used to define the LO frequency.
48 *     - ImagFreq: The image sideband equivalent of the rest frequency.
49 *     - SideBand: Indicates which sideband the DSBSpecFrame represents.
50 
51 *  Functions:
52 c     The DSBSpecFrame class does not define any new functions beyond those
53 f     The DSBSpecFrame class does not define any new routines beyond those
54 *     which are applicable to all SpecFrames.
55 
56 *  Copyright:
57 *     Copyright (C) 1997-2006 Council for the Central Laboratory of the
58 *     Research Councils
59 
60 *  Licence:
61 *     This program is free software: you can redistribute it and/or
62 *     modify it under the terms of the GNU Lesser General Public
63 *     License as published by the Free Software Foundation, either
64 *     version 3 of the License, or (at your option) any later
65 *     version.
66 *
67 *     This program is distributed in the hope that it will be useful,
68 *     but WITHOUT ANY WARRANTY; without even the implied warranty of
69 *     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
70 *     GNU Lesser General Public License for more details.
71 *
72 *     You should have received a copy of the GNU Lesser General
73 *     License along with this program.  If not, see
74 *     <http://www.gnu.org/licenses/>.
75 
76 *  Authors:
77 *     DSB: David Berry (Starlink)
78 
79 *  History:
80 *     5-AUG-2004 (DSB):
81 *        Original version.
82 *     7-OCT-2004 (DSB):
83 *        Fixed SetAttrib code which assigns values to SideBand. Previously
84 *        all supplied values were ignored, leaving SideBand unchanged.
85 *     2-SEP-2005 (DSB):
86 *        Allow conversion in any Domain within TopoMap (sometimes
87 *        SpecFrames have a new Domain set which is not equal to SPECTRUM").
88 *     12-SEP-2005 (DSB):
89 *        Set all attributes required to described the RestFreq value
90 *        before determining Mapping from RestFreq to ImagFreq in
91 *        GetImageFreq.
92 *     2-DEC-2005 (DSB):
93 *        Change default Domain from SPECTRUM to DSBSPECTRUM
94 *     3-APR-2006 (DSB):
95 *        Fix memory leak in astLoadDSBSpecFrame.
96 *     6-OCT-2006 (DSB):
97 *        Guard against annulling null pointers in subFrame.
98 *     27-OCT-2006 (DSB):
99 *        Added AlignSideBand attribute.
100 *     31-OCT-2006 (DSB):
101 *        Use AlignSideBand attribute in SubFrame only if we are not
102 *        currently restoring a FrameSet's integrity.
103 *     31-JAN-2007 (DSB):
104 *        Modified so that a DSBSpecFrame can be used as a template to find a
105 *        DSBSpecFrame (or SpecFrame) contained within a CmpFrame. This
106 *        involves changes in Match.
107 *     1-MAY-2007 (DSB):
108 *        The default for AlignSideband has been changed from 1 to 0.
109 *     8-MAY-2007 (DSB):
110 *        Correct initialisation of alignsideband in astInitDSBSpecFrame_.
111 *     19-OCT-2007 (DSB):
112 *        Ignore SideBand alignment if the AlignSideBand attribute is zero
113 *        in either the target or the template.
114 *     16-JAN-2007 (DSB):
115 *        Modify SubFrame so that DSBSpecFrames are aligned in the
116 *        observed sideband (LSB or USB) rather than always being aligned
117 *        in the USB.
118 *     12-FEB-2010 (DSB):
119 *        Report an error if the local oscillator frequency looks silly
120 *        (specifically, if it less than the absolute intermediate frequency).
121 *     29-APR-2011 (DSB):
122 *        Prevent astFindFrame from matching a subclass template against a
123 *        superclass target.
124 *class--
125 
126 *  Implementation Deficiencies:
127 *     - The default values for System and StdOfRest inherited from the
128 *     SpecFrame class are "Wave" and "Heliocentric". These are not
129 *     usually what is wanted for a DSB instrument. Defaults such as
130 *     "Freq" and "Topo" may be more appropriate. However, changing the
131 *     defaults inherited from SpecFrame may cause problems in the
132 *     astConvert algorithm. The astConvertX function in frame.c includes
133 *     the following implementation deficiency warning: "One likely
134 *     problem is with attributes which default in both the source and
135 *     destination Frames. This means they also default in the common
136 *     coordinate system. If these default values were to differ when
137 *     matching different target Frames, however, we would be in trouble,
138 *     because the common coordinate system would not then be remaining
139 *     constant. The longer-term solution to this is probably to provide
140 *     some mechanism to "fix" all attribute values for a Frame, by taking
141 *     any attributes that are un-set and explicitly setting a firm value
142 *     (equal to the default) so they cannot then change". So the defaults
143 *     should probably be left unchanged until this fix is made.
144 
145 */
146 
147 /* Module Macros. */
148 /* ============== */
149 /* Set the name of the class we are implementing. This indicates to
150    the header files that define class interfaces that they should make
151    "protected" symbols available. */
152 #define astCLASS DSBSpecFrame
153 
154 #define BADSB -9999
155 #define FIRST_SB -1
156 #define LSB -1
157 #define LO  0
158 #define USB 1
159 #define LAST_SB 1
160 
161 /* Include files. */
162 /* ============== */
163 /* Interface definitions. */
164 /* ---------------------- */
165 
166 #include "globals.h"             /* Thread-safe global data access */
167 #include "error.h"               /* Error reporting facilities */
168 #include "memory.h"              /* Memory management facilities */
169 #include "object.h"              /* Base Object class */
170 #include "channel.h"             /* I/O channels */
171 #include "specframe.h"           /* Spectral frames (parent class) */
172 #include "unit.h"                /* Unit handling */
173 #include "cmpmap.h"              /* Compound Mappings */
174 #include "unitmap.h"             /* Unit Mappings */
175 #include "winmap.h"              /* Window Mappings */
176 #include "dsbspecframe.h"        /* Interface definition for this class */
177 #include "globals.h"             /* Thread-safe global data access */
178 
179 /* Error code definitions. */
180 /* ----------------------- */
181 #include "ast_err.h"             /* AST error codes */
182 
183 /* C header files. */
184 /* --------------- */
185 #include <stdarg.h>
186 #include <stddef.h>
187 #include <stdio.h>
188 #include <string.h>
189 #include <ctype.h>
190 
191 /* Module Variables. */
192 /* ================= */
193 
194 /* Address of this static variable is used as a unique identifier for
195    member of this class. */
196 static int class_check;
197 
198 /* Pointers to parent class methods which are extended by this class. */
199 static const char *(* parent_getattrib)( AstObject *, const char *, int * );
200 static const char *(* parent_getlabel)( AstFrame *, int, int * );
201 static int (* parent_match)( AstFrame *, AstFrame *, int, int **, int **, AstMapping **, AstFrame **, int * );
202 static int (* parent_subframe)( AstFrame *, AstFrame *, int, const int *, const int *, AstMapping **, AstFrame **, int * );
203 static int (* parent_testattrib)( AstObject *, const char *, int * );
204 static void (* parent_clearattrib)( AstObject *, const char *, int * );
205 static void (* parent_setattrib)( AstObject *, const char *, int * );
206 static void (* parent_overlay)( AstFrame *, const int *, AstFrame *, int * );
207 static const char *(* parent_getdomain)( AstFrame *, int * );
208 
209 /* Define macros for accessing each item of thread specific global data. */
210 #ifdef THREAD_SAFE
211 
212 /* Define how to initialise thread-specific globals. */
213 #define GLOBAL_inits \
214    globals->Class_Init = 0; \
215    globals->GetAttrib_Buff[ 0 ] = 0; \
216    globals->GetLabel_Buff[ 0 ] = 0; \
217 
218 /* Create the function that initialises global data for this module. */
219 astMAKE_INITGLOBALS(DSBSpecFrame)
220 
221 /* Define macros for accessing each item of thread specific global data. */
222 #define class_init astGLOBAL(DSBSpecFrame,Class_Init)
223 #define class_vtab astGLOBAL(DSBSpecFrame,Class_Vtab)
224 #define getattrib_buff astGLOBAL(DSBSpecFrame,GetAttrib_Buff)
225 #define getlabel_buff astGLOBAL(DSBSpecFrame,GetLabel_Buff)
226 
227 
228 
229 /* If thread safety is not needed, declare and initialise globals at static
230    variables. */
231 #else
232 
233 /* Define the thread-specific globals for this class. */
234 
235 /* Buffer returned by GetAttrib. */
236 static char getattrib_buff[ 101 ];
237 
238 /* Default Label string buffer */
239 static char getlabel_buff[ 101 ];
240 
241 
242 /* Define the class virtual function table and its initialisation flag
243    as static variables. */
244 static AstDSBSpecFrameVtab class_vtab;   /* Virtual function table */
245 static int class_init = 0;       /* Virtual function table initialised? */
246 
247 #endif
248 
249 /* External Interface Function Prototypes. */
250 /* ======================================= */
251 /* The following functions have public prototypes only (i.e. no
252    protected prototypes), so we must provide local prototypes for use
253    within this module. */
254 AstDSBSpecFrame *astDSBSpecFrameId_( const char *, ... );
255 
256 /* Prototypes for Private Member Functions. */
257 /* ======================================== */
258 
259 static AstMapping *TopoMap( AstDSBSpecFrame *, int, const char *, int * );
260 static AstMapping *ToLOMapping( AstDSBSpecFrame *, const char *, int * )__attribute__((unused));
261 static AstMapping *ToLSBMapping( AstDSBSpecFrame *, const char *, int * );
262 static AstMapping *ToUSBMapping( AstDSBSpecFrame *, const char *, int * );
263 static const char *GetAttrib( AstObject *, const char *, int * );
264 static const char *GetLabel( AstFrame *, int, int * );
265 static double GetImagFreq( AstDSBSpecFrame *, int * );
266 static double GetLO( AstDSBSpecFrame *, const char *, const char *, int * );
267 static int Match( AstFrame *, AstFrame *, int, int **, int **, AstMapping **, AstFrame **, int * );
268 static int SubFrame( AstFrame *, AstFrame *, int, const int *, const int *, AstMapping **, AstFrame **, int * );
269 static int TestAttrib( AstObject *, const char *, int * );
270 static void ClearAttrib( AstObject *, const char *, int * );
271 static void Dump( AstObject *, AstChannel *, int * );
272 static void Overlay( AstFrame *, const int *, AstFrame *, int * );
273 static void SetAttrib( AstObject *, const char *, int * );
274 static void VerifyAttrs( AstDSBSpecFrame *, const char *, const char *, const char *, int * );
275 static const char *GetDomain( AstFrame *, int * );
276 
277 static double GetIF( AstDSBSpecFrame *, int * );
278 static int TestIF( AstDSBSpecFrame *, int * );
279 static void ClearIF( AstDSBSpecFrame *, int * );
280 static void SetIF( AstDSBSpecFrame *, double, int * );
281 
282 static double GetDSBCentre( AstDSBSpecFrame *, int * );
283 static int TestDSBCentre( AstDSBSpecFrame *, int * );
284 static void ClearDSBCentre( AstDSBSpecFrame *, int * );
285 static void SetDSBCentre( AstDSBSpecFrame *, double, int * );
286 
287 static int GetSideBand( AstDSBSpecFrame *, int * );
288 static int TestSideBand( AstDSBSpecFrame *, int * );
289 static void ClearSideBand( AstDSBSpecFrame *, int * );
290 static void SetSideBand( AstDSBSpecFrame *, int, int * );
291 
292 static int GetAlignSideBand( AstDSBSpecFrame *, int * );
293 static int TestAlignSideBand( AstDSBSpecFrame *, int * );
294 static void ClearAlignSideBand( AstDSBSpecFrame *, int * );
295 static void SetAlignSideBand( AstDSBSpecFrame *, int, int * );
296 
297 
298 /* Member functions. */
299 /* ================= */
ClearAttrib(AstObject * this_object,const char * attrib,int * status)300 static void ClearAttrib( AstObject *this_object, const char *attrib, int *status ) {
301 /*
302 *  Name:
303 *     ClearAttrib
304 
305 *  Purpose:
306 *     Clear an attribute value for a DSBSpecFrame.
307 
308 *  Type:
309 *     Private function.
310 
311 *  Synopsis:
312 *     #include "dsbspecframe.h"
313 *     void ClearAttrib( AstObject *this, const char *attrib, int *status )
314 
315 *  Class Membership:
316 *     DSBSpecFrame member function (over-rides the astClearAttrib protected
317 *     method inherited from the SpecFrame class).
318 
319 *  Description:
320 *     This function clears the value of a specified attribute for a
321 *     DSBSpecFrame, so that the default value will subsequently be used.
322 
323 *  Parameters:
324 *     this
325 *        Pointer to the DSBSpecFrame.
326 *     attrib
327 *        Pointer to a null-terminated string specifying the attribute
328 *        name.  This should be in lower case with no surrounding white
329 *        space.
330 *     status
331 *        Pointer to the inherited status variable.
332 */
333 
334 /* Local Variables: */
335    AstDSBSpecFrame *this;        /* Pointer to the DSBSpecFrame structure */
336 
337 /* Check the global error status. */
338    if ( !astOK ) return;
339 
340 /* Obtain a pointer to the DSBSpecFrame structure. */
341    this = (AstDSBSpecFrame *) this_object;
342 
343 /* Check the attribute name and clear the appropriate attribute. */
344 
345 /* DSBCentre. */
346 /* ---------- */
347    if ( !strcmp( attrib, "dsbcentre" ) ) {
348       astClearDSBCentre( this );
349 
350 /* IF */
351 /* -- */
352    } else if ( !strcmp( attrib, "if" ) ) {
353       astClearIF( this );
354 
355 /* SideBand */
356 /* -------- */
357    } else if ( !strcmp( attrib, "sideband" ) ) {
358       astClearSideBand( this );
359 
360 /* AlignSideBand */
361 /* ------------- */
362    } else if ( !strcmp( attrib, "alignsideband" ) ) {
363       astClearAlignSideBand( this );
364 
365 /* Read-only attributes. */
366 /* --------------------- */
367 /* Test if the attribute name matches any of the read-only attributes
368    of this class. If it does, then report an error. */
369    } else if ( !strcmp( attrib, "imagfreq" ) ) {
370       astError( AST__NOWRT, "astClear: Invalid attempt to clear the \"%s\" "
371                 "value for a %s.", status, attrib, astGetClass( this ) );
372       astError( AST__NOWRT, "This is a read-only attribute." , status);
373 
374 /* If the attribute is not recognised, pass it on to the parent method
375    for further interpretation. */
376    } else {
377       (*parent_clearattrib)( this_object, attrib, status );
378    }
379 }
380 
381 
GetAttrib(AstObject * this_object,const char * attrib,int * status)382 static const char *GetAttrib( AstObject *this_object, const char *attrib, int *status ) {
383 /*
384 *  Name:
385 *     GetAttrib
386 
387 *  Purpose:
388 *     Get the value of a specified attribute for a DSBSpecFrame.
389 
390 *  Type:
391 *     Private function.
392 
393 *  Synopsis:
394 *     #include "dsbspecframe.h"
395 *     const char *GetAttrib( AstObject *this, const char *attrib, int *status )
396 
397 *  Class Membership:
398 *     DSBSpecFrame member function (over-rides the protected astGetAttrib
399 *     method inherited from the SpecFrame class).
400 
401 *  Description:
402 *     This function returns a pointer to the value of a specified
403 *     attribute for a DSBSpecFrame, formatted as a character string.
404 
405 *  Parameters:
406 *     this
407 *        Pointer to the DSBSpecFrame.
408 *     attrib
409 *        Pointer to a null-terminated string containing the name of
410 *        the attribute whose value is required. This name should be in
411 *        lower case, with all white space removed.
412 *     status
413 *        Pointer to the inherited status variable.
414 
415 *  Returned Value:
416 *     - Pointer to a null-terminated string containing the attribute
417 *     value.
418 
419 *  Notes:
420 *     - The returned string pointer may point at memory allocated
421 *     within the DSBSpecFrame, or at static memory. The contents of the
422 *     string may be over-written or the pointer may become invalid
423 *     following a further invocation of the same function or any
424 *     modification of the DSBSpecFrame. A copy of the string should
425 *     therefore be made if necessary.
426 *     - A NULL pointer will be returned if this function is invoked
427 *     with the global error status set, or if it should fail for any
428 *     reason.
429 */
430 
431 /* Local Variables: */
432    astDECLARE_GLOBALS            /* Declare the thread specific global data */
433    AstDSBSpecFrame *this;        /* Pointer to the DSBSpecFrame structure */
434    AstMapping *tmap;             /* Ptr to Mapping from topofreq to this */
435    const char *result;           /* Pointer value to return */
436    double dval;                  /* Attribute value */
437    double dtemp;                 /* Attribute value */
438    int ival;                     /* Attribute value */
439 
440 /* Initialise. */
441    result = NULL;
442 
443 /* Check the global error status. */
444    if ( !astOK ) return result;
445 
446 /* Get a pointer to the structure holding thread-specific global data. */
447    astGET_GLOBALS(this_object);
448 
449 /* Obtain a pointer to the SpecFrame structure. */
450    this = (AstDSBSpecFrame *) this_object;
451 
452 /* Compare "attrib" with each recognised attribute name in turn,
453    obtaining the value of the required attribute. If necessary, write
454    the value into "getattrib_buff" as a null-terminated string in an appropriate
455    format.  Set "result" to point at the result string. */
456 
457 /* DSBCentre */
458 /* --------- */
459    if ( !strcmp( attrib, "dsbcentre" ) ) {
460 
461 /* Get the value as topocentric frequency in Hz. */
462       dval = astGetDSBCentre( this );
463 
464 /* Find the Mapping from topocentric frequency in Hz to the spectral system
465    described by this SpecFrame. */
466       tmap = TopoMap( this, 0, "astGetAttrib", status );
467       if ( astOK ) {
468 
469 /* Transform the internal value from topocentric frequency into the required
470    system. */
471          astTran1( tmap, 1, &dval, 1, &dtemp );
472          if( dtemp == AST__BAD ) {
473             astError( AST__INTER, "astGetAttrib(%s): Cannot convert DSBCentre "
474                       "value from topocentric frequency to the required "
475                       "system.", status, astGetClass( this ) );
476          } else {
477 
478 /* Format it. */
479             (void) sprintf( getattrib_buff, "%.*g", DBL_DIG, dtemp );
480             result = getattrib_buff;
481          }
482          tmap = astAnnul( tmap );
483       }
484 
485 /* IF */
486 /* -- */
487    } else if ( !strcmp( attrib, "if" ) ) {
488       dval = astGetIF( this );
489       if ( astOK ) {
490          (void) sprintf( getattrib_buff, "%.*g", DBL_DIG, dval*1.0E-9 );
491          result = getattrib_buff;
492       }
493 
494 /* ImagFreq */
495 /* -------- */
496    } else if ( !strcmp( attrib, "imagfreq" ) ) {
497       dval = astGetImagFreq( this );
498       if ( astOK ) {
499          (void) sprintf( getattrib_buff, "%.*g", DBL_DIG, dval*1.0E-9 );
500          result = getattrib_buff;
501       }
502 
503 /* SideBand */
504 /* -------- */
505    } else if ( !strcmp( attrib, "sideband" ) ) {
506       ival = astGetSideBand( this );
507       if ( astOK ) {
508          result = ( ival == USB ) ? "USB" : (( ival == LO ) ? "LO" : "LSB" );
509       }
510 
511 /* AlignSideBand */
512 /* ------------- */
513    } else if ( !strcmp( attrib, "alignsideband" ) ) {
514       ival = astGetAlignSideBand( this ) ? 1 : 0;
515       if ( astOK ) {
516          (void) sprintf( getattrib_buff, "%d", ival );
517          result = getattrib_buff;
518       }
519 
520 /* If the attribute name was not recognised, pass it on to the parent
521    method for further interpretation. */
522    } else {
523       result = (*parent_getattrib)( this_object, attrib, status );
524    }
525 
526 /* Return the result. */
527    return result;
528 
529 }
530 
GetDomain(AstFrame * this_frame,int * status)531 static const char *GetDomain( AstFrame *this_frame, int *status ) {
532 /*
533 *  Name:
534 *     GetDomain
535 
536 *  Purpose:
537 *     Obtain a pointer to the Domain attribute string for a DSBSpecFrame.
538 
539 *  Type:
540 *     Private function.
541 
542 *  Synopsis:
543 *     #include "dsbspecframe.h"
544 *     const char *GetDomain( AstFrame *this, int *status )
545 
546 *  Class Membership:
547 *     DSBSpecFrame member function (over-rides the astGetDomain protected
548 *     method inherited from the SpecFrame class).
549 
550 *  Description:
551 *    This function returns a pointer to the Domain attribute string
552 *    for a DSBSpecFrame.
553 
554 *  Parameters:
555 *     this
556 *        Pointer to the DSBSpecFrame.
557 *     status
558 *        Pointer to the inherited status variable.
559 
560 *  Returned Value:
561 *     A pointer to a constant null-terminated string containing the
562 *     Domain value.
563 
564 *  Notes:
565 *     - The returned pointer or the string it refers to may become
566 *     invalid following further invocation of this function or
567 *     modification of the DSBSpecFrame.
568 *     - A NULL pointer is returned if this function is invoked with
569 *     the global error status set or if it should fail for any reason.
570 */
571 
572 /* Local Variables: */
573    AstDSBSpecFrame *this;        /* Pointer to DSBSpecFrame structure */
574    const char *result;           /* Pointer value to return */
575 
576 /* Initialise. */
577    result = NULL;
578 
579 /* Check the global error status. */
580    if ( !astOK ) return result;
581 
582 /* Obtain a pointer to the DSBSpecFrame structure. */
583    this = (AstDSBSpecFrame *) this_frame;
584 
585 /* If a Domain attribute string has been set, invoke the parent method
586    to obtain a pointer to it. */
587    if ( astTestDomain( this ) ) {
588       result = (*parent_getdomain)( this_frame, status );
589 
590 /* Otherwise, provide a pointer to a suitable default string. */
591    } else {
592       result = "DSBSPECTRUM";
593    }
594 
595 /* Return the result. */
596    return result;
597 }
598 
GetImagFreq(AstDSBSpecFrame * this,int * status)599 static double GetImagFreq( AstDSBSpecFrame *this, int *status ) {
600 /*
601 *+
602 *  Name:
603 *     astGetImagFreq
604 
605 *  Purpose:
606 *     Get the value of the ImagFreq attribute.
607 
608 *  Type:
609 *     Protected virtual function.
610 
611 *  Synopsis:
612 *     #include "dsbspecframe.h"
613 *     double GetImagFreq( AstDSBSpecFrame *this )
614 
615 *  Class Membership:
616 *     DSBSpecFrame method.
617 
618 *  Description:
619 *     This function returns the image sideband frequency corresponding to
620 *     the rest frequency.
621 
622 *  Parameters:
623 *     this
624 *        Pointer to the Frame.
625 
626 *  Returned Value:
627 *     The required frequency, in Hz.
628 
629 *  Notes:
630 *     - A value of AST__BAD will be returned if this function is invoked
631 *     with the global error status set, or if it should fail for any
632 *     reason.
633 *-
634 */
635 
636 /* Local Variables: */
637    AstDSBSpecFrame *rf_frame;/* DSBSpecFrame describing the rest frequency */
638    AstMapping *map;          /* Pointer to "Observed to Image" mapping */
639    double result;            /* The returned frequency */
640    double rf;                /* Rest frequency in observed sideband */
641    int sb;                   /* SideBand value */
642 
643 /* Check the global error status. */
644    if ( !astOK ) return AST__BAD;
645 
646 /* The RestFreq attribute is an observed sideband frequency in the
647    source's standard of rest, measured in Hz. Temporaily set attributes
648    to these values. Create a copy of the supplied DSBSpecFrame and set
649    its attributes to these values.  */
650    rf_frame = astCopy( this );
651    astSetStdOfRest( rf_frame, AST__SCSOR );
652    astSetSystem( rf_frame, AST__FREQ );
653    astSetUnit( rf_frame, 0, "Hz" );
654    astSetC( rf_frame, "SideBand", "observed" );
655 
656 /* Create a Mapping which transforms positions from the observed to the
657    image sideband. */
658    sb = astGetSideBand( rf_frame );
659    if( sb == USB ) {
660       map = ToLSBMapping( rf_frame, "astGetImagFreq", status );
661 
662    } else if( sb == LSB ) {
663       map = ToUSBMapping( rf_frame, "astGetImagFreq", status );
664 
665    } else {
666       map = NULL;
667       astError( AST__INTER, "astGetImagFreq(%s): Illegal sideband value "
668                 "(%d) encountered (internal AST programming error).", status,
669                 astGetClass( this ), sb );
670    }
671 
672 /* Get the rest frequency in Hz, and transform it using the above Mapping. */
673    rf = astGetRestFreq( rf_frame );
674    astTran1( map, 1, &rf, 1, &result );
675 
676 /* Free resources */
677    map = astAnnul( map );
678    rf_frame = astAnnul( rf_frame );
679 
680 /* If an error has occurrred, return AST__BAD. */
681    if( !astOK ) result = AST__BAD;
682 
683 /* Return the result. */
684    return result;
685 
686 }
687 
GetLabel(AstFrame * this,int axis,int * status)688 static const char *GetLabel( AstFrame *this, int axis, int *status ) {
689 /*
690 *  Name:
691 *     GetLabel
692 
693 *  Purpose:
694 *     Access the Label string for a DSBSpecFrame axis.
695 
696 *  Type:
697 *     Private function.
698 
699 *  Synopsis:
700 *     #include "dsbspecframe.h"
701 *     const char *GetLabel( AstFrame *this, int axis, int *status )
702 
703 *  Class Membership:
704 *     DSBSpecFrame member function (over-rides the astGetLabel method inherited
705 *     from the SpecFrame class).
706 
707 *  Description:
708 *     This function returns a pointer to the Label string for a specified axis
709 *     of a DSBSpecFrame.
710 
711 *  Parameters:
712 *     this
713 *        Pointer to the SpecFrame.
714 *     axis
715 *        Axis index (zero-based) identifying the axis for which information is
716 *        required.
717 *     status
718 *        Pointer to the inherited status variable.
719 
720 *  Returned Value:
721 *     Pointer to a constant null-terminated character string containing the
722 *     requested information.
723 
724 *  Notes:
725 *     -  A NULL pointer will be returned if this function is invoked with the
726 *     global error status set, or if it should fail for any reason.
727 */
728 
729 /* Local Variables: */
730    astDECLARE_GLOBALS            /* Declare the thread specific global data */
731    const char *result;           /* Pointer to label string */
732 
733 /* Check the global error status. */
734    if ( !astOK ) return NULL;
735 
736 /* Get a pointer to the structure holding thread-specific global data. */
737    astGET_GLOBALS(this);
738 
739 /* Initialise. */
740    result = NULL;
741 
742 /* Validate the axis index. */
743    astValidateAxis( this, axis, 1, "astGetLabel" );
744 
745 /* Invoke the parent astGetLabel method to obtain a pointer to it. */
746    result = (*parent_getlabel)( this, axis, status );
747 
748 /* Check if this is a default value. If so, append a string indicating
749    the sideband. */
750    if ( !astTestLabel( this, axis ) ) {
751 
752 /* If OK, supply a pointer to a suitable default label string. */
753       sprintf( getlabel_buff, "%s (%s)", result, astGetAttrib( this, "sideband" ) );
754       result = getlabel_buff;
755    }
756 
757 /* Return the result. */
758    return result;
759 
760 }
761 
GetLO(AstDSBSpecFrame * this,const char * check_msg,const char * method,int * status)762 static double GetLO( AstDSBSpecFrame *this, const char *check_msg,
763                      const char *method, int *status ) {
764 /*
765 *  Name:
766 *     GetLO
767 
768 *  Purpose:
769 *     Get the Local Oscillator frequency.
770 
771 *  Type:
772 *     Private function.
773 
774 *  Synopsis:
775 *     #include "dsbspecframe.h"
776 *     double GetLO( AstDSBSpecFrame *this, const char *check_msg,
777 *                   const char *method, int *status )
778 
779 *  Class Membership:
780 *     DSBSpecFrame method.
781 
782 *  Description:
783 *     This function returns the local oscillator frequency in topocentric
784 *     frequency.
785 
786 *  Parameters:
787 *     this
788 *        Pointer to the Frame.
789 *     check_msg
790 *        If not NULL, an error will be reported if either the DSBCentre
791 *        or IF attribute has not been set to an explicit value. In this
792 *        case, the error message will include the supplied text.
793 *     method
794 *        The name of the calling method - used in error messages.
795 *     status
796 *        Pointer to the inherited status value.
797 
798 *  Returned Value:
799 *     The local oscillator frequency, in Hz.
800 
801 *  Notes:
802 *     - An error is reported if the local oscillator frequency looks
803 *     un-physical (specifically, if it is less than the absolute value of
804 *     the intermediate frequency).
805 *     - A value of AST__BAD will be returned if this function is invoked
806 *     with the global error status set, or if it should fail for any
807 *     reason.
808 *-
809 */
810 
811 /* Local Variables: */
812    double f_if;              /* Intermediate frequency (topo,HZ) */
813    double result;            /* The returned frequency */
814 
815 /* Check the global error status. */
816    if ( !astOK ) return AST__BAD;
817 
818 /* If required, check that explicit values have been assigned to the required
819    attributes (i.e. report an error if a default value would be used for
820    either attribute). */
821    if( check_msg ) VerifyAttrs( this, check_msg, "IF DSBCentre", method,
822                                 status );
823 
824 /* The local oscillator is the sum of the intermediate frequency and the
825    observation centre frequency. */
826    f_if = astGetIF( this );
827    result = astGetDSBCentre( this ) + f_if;
828 
829 /* Check the local  oscillator frequency is no smaller than the absolute
830    intermediate frequency. */
831    if( result < fabs( f_if ) && astOK ) {
832       astError( AST__ATTIN, "%s(%s): The local oscillator frequency (%g Hz) "
833                 "is too low (less than the intermediate frequency: %g Hz).",
834                 status, method, astGetClass( this ), result, fabs( f_if ) );
835       astError( AST__ATTIN, "   This could be caused by a bad value for"
836                 " either the IF attribute (currently %g Hz) or the DSBCentre "
837                 "attribute (currently %g Hz).", status, f_if,
838                 astGetDSBCentre( this ) );
839    }
840 
841 /* If an error has occurrred, return AST__BAD. */
842    if( !astOK ) result = AST__BAD;
843 
844 /* Return the result. */
845    return result;
846 }
847 
astInitDSBSpecFrameVtab_(AstDSBSpecFrameVtab * vtab,const char * name,int * status)848 void astInitDSBSpecFrameVtab_(  AstDSBSpecFrameVtab *vtab, const char *name, int *status ) {
849 /*
850 *+
851 *  Name:
852 *     astInitDSBSpecFrameVtab
853 
854 *  Purpose:
855 *     Initialise a virtual function table for a DSBSpecFrame.
856 
857 *  Type:
858 *     Protected function.
859 
860 *  Synopsis:
861 *     #include "dsbspecframe.h"
862 *     void astInitDSBSpecFrameVtab( AstDSBSpecFrameVtab *vtab, const char *name )
863 
864 *  Class Membership:
865 *     DSBSpecFrame vtab initialiser.
866 
867 *  Description:
868 *     This function initialises the component of a virtual function
869 *     table which is used by the DSBSpecFrame class.
870 
871 *  Parameters:
872 *     vtab
873 *        Pointer to the virtual function table. The components used by
874 *        all ancestral classes will be initialised if they have not already
875 *        been initialised.
876 *     name
877 *        Pointer to a constant null-terminated character string which contains
878 *        the name of the class to which the virtual function table belongs (it
879 *        is this pointer value that will subsequently be returned by the Object
880 *        astClass function).
881 *-
882 */
883 
884 /* Local Variables: */
885    astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
886    AstObjectVtab *object;        /* Pointer to Object component of Vtab */
887    AstFrameVtab *frame;          /* Pointer to Frame component of Vtab */
888 
889 /* Check the local error status. */
890    if ( !astOK ) return;
891 
892 /* Get a pointer to the thread specific global data structure. */
893    astGET_GLOBALS(NULL);
894 
895 /* Initialize the component of the virtual function table used by the
896    parent class. */
897    astInitSpecFrameVtab( (AstSpecFrameVtab *) vtab, name );
898 
899 /* Store a unique "magic" value in the virtual function table. This
900    will be used (by astIsADSBSpecFrame) to determine if an object belongs
901    to this class.  We can conveniently use the address of the (static)
902    class_check variable to generate this unique value. */
903    vtab->id.check = &class_check;
904    vtab->id.parent = &(((AstSpecFrameVtab *) vtab)->id);
905 
906 /* Initialise member function pointers. */
907 /* ------------------------------------ */
908 /* Store pointers to the member functions (implemented here) that provide
909    virtual methods for this class. */
910    vtab->ClearDSBCentre = ClearDSBCentre;
911    vtab->TestDSBCentre = TestDSBCentre;
912    vtab->GetDSBCentre = GetDSBCentre;
913    vtab->SetDSBCentre = SetDSBCentre;
914 
915    vtab->ClearIF = ClearIF;
916    vtab->TestIF = TestIF;
917    vtab->GetIF = GetIF;
918    vtab->SetIF = SetIF;
919 
920    vtab->ClearSideBand = ClearSideBand;
921    vtab->TestSideBand = TestSideBand;
922    vtab->GetSideBand = GetSideBand;
923    vtab->SetSideBand = SetSideBand;
924 
925    vtab->ClearAlignSideBand = ClearAlignSideBand;
926    vtab->TestAlignSideBand = TestAlignSideBand;
927    vtab->GetAlignSideBand = GetAlignSideBand;
928    vtab->SetAlignSideBand = SetAlignSideBand;
929 
930    vtab->GetImagFreq = GetImagFreq;
931 
932 /* Save the inherited pointers to methods that will be extended, and
933    replace them with pointers to the new member functions. */
934    object = (AstObjectVtab *) vtab;
935    frame = (AstFrameVtab *) vtab;
936 
937    parent_clearattrib = object->ClearAttrib;
938    object->ClearAttrib = ClearAttrib;
939 
940    parent_getattrib = object->GetAttrib;
941    object->GetAttrib = GetAttrib;
942 
943    parent_setattrib = object->SetAttrib;
944    object->SetAttrib = SetAttrib;
945 
946    parent_testattrib = object->TestAttrib;
947    object->TestAttrib = TestAttrib;
948 
949    parent_getdomain = frame->GetDomain;
950    frame->GetDomain = GetDomain;
951 
952    parent_overlay = frame->Overlay;
953    frame->Overlay = Overlay;
954 
955    parent_match = frame->Match;
956    frame->Match = Match;
957 
958    parent_subframe = frame->SubFrame;
959    frame->SubFrame = SubFrame;
960 
961    parent_getlabel = frame->GetLabel;
962    frame->GetLabel = GetLabel;
963 
964 /* Declare the class delete function.*/
965    astSetDump( vtab, Dump, "DSBSpecFrame", "Dual sideband spectral axis" );
966 
967 /* If we have just initialised the vtab for the current class, indicate
968    that the vtab is now initialised, and store a pointer to the class
969    identifier in the base "object" level of the vtab. */
970    if( vtab == &class_vtab ) {
971       class_init = 1;
972       astSetVtabClassIdentifier( vtab, &(vtab->id) );
973    }
974 }
975 
Match(AstFrame * template_frame,AstFrame * target,int matchsub,int ** template_axes,int ** target_axes,AstMapping ** map,AstFrame ** result,int * status)976 static int Match( AstFrame *template_frame, AstFrame *target, int matchsub,
977                   int **template_axes, int **target_axes, AstMapping **map,
978                   AstFrame **result, int *status ) {
979 /*
980 *  Name:
981 *     Match
982 
983 *  Purpose:
984 *     Determine if conversion is possible between two coordinate systems.
985 
986 *  Type:
987 *     Private function.
988 
989 *  Synopsis:
990 *     #include "dsbspecframe.h"
991 *     int Match( AstFrame *template, AstFrame *target, int matchsub,
992 *                int **template_axes, int **target_axes,
993 *                AstMapping **map, AstFrame **result, int *status )
994 
995 *  Class Membership:
996 *     DSBSpecFrame member function (over-rides the protected astMatch method
997 *     inherited from the SpecFrame class).
998 
999 *  Description:
1000 *     This function matches a "template" DSBSpecFrame to a "target" Frame and
1001 *     determines whether it is possible to convert coordinates between them.
1002 *     If it is, a mapping that performs the transformation is returned along
1003 *     with a new Frame that describes the coordinate system that results when
1004 *     this mapping is applied to the "target" coordinate system. In addition,
1005 *     information is returned to allow the axes in this "result" Frame to be
1006 *     associated with the corresponding axes in the "target" and "template"
1007 *     Frames from which they are derived.
1008 
1009 *  Parameters:
1010 *     template
1011 *        Pointer to the template DSBSpecFrame. This describes the coordinate
1012 *        system (or set of possible coordinate systems) into which we wish to
1013 *        convert our coordinates.
1014 *     target
1015 *        Pointer to the target Frame. This describes the coordinate system in
1016 *        which we already have coordinates.
1017 *     matchsub
1018 *        If zero then a match only occurs if the template is of the same
1019 *        class as the target, or of a more specialised class. If non-zero
1020 *        then a match can occur even if this is not the case.
1021 *     template_axes
1022 *        Address of a location where a pointer to int will be returned if the
1023 *        requested coordinate conversion is possible. This pointer will point
1024 *        at a dynamically allocated array of integers with one element for each
1025 *        axis of the "result" Frame (see below). It must be freed by the caller
1026 *        (using astFree) when no longer required.
1027 *
1028 *        For each axis in the result Frame, the corresponding element of this
1029 *        array will return the index of the template DSBSpecFrame axis from
1030 *        which it is derived. If it is not derived from any template
1031 *        DSBSpecFrame axis, a value of -1 will be returned instead.
1032 *     target_axes
1033 *        Address of a location where a pointer to int will be returned if the
1034 *        requested coordinate conversion is possible. This pointer will point
1035 *        at a dynamically allocated array of integers with one element for each
1036 *        axis of the "result" Frame (see below). It must be freed by the caller
1037 *        (using astFree) when no longer required.
1038 *
1039 *        For each axis in the result Frame, the corresponding element of this
1040 *        array will return the index of the target Frame axis from which it
1041 *        is derived. If it is not derived from any target Frame axis, a value
1042 *        of -1 will be returned instead.
1043 *     map
1044 *        Address of a location where a pointer to a new Mapping will be
1045 *        returned if the requested coordinate conversion is possible. If
1046 *        returned, the forward transformation of this Mapping may be used to
1047 *        convert coordinates between the "target" Frame and the "result"
1048 *        Frame (see below) and the inverse transformation will convert in the
1049 *        opposite direction.
1050 *     result
1051 *        Address of a location where a pointer to a new Frame will be returned
1052 *        if the requested coordinate conversion is possible. If returned, this
1053 *        Frame describes the coordinate system that results from applying the
1054 *        returned Mapping (above) to the "target" coordinate system. In
1055 *        general, this Frame will combine attributes from (and will therefore
1056 *        be more specific than) both the target and the template Frames. In
1057 *        particular, when the template allows the possibility of transformaing
1058 *        to any one of a set of alternative coordinate systems, the "result"
1059 *        Frame will indicate which of the alternatives was used.
1060 *     status
1061 *        Pointer to the inherited status variable.
1062 
1063 *  Returned Value:
1064 *     A non-zero value is returned if the requested coordinate conversion is
1065 *     possible. Otherwise zero is returned (this will not in itself result in
1066 *     an error condition).
1067 
1068 *  Notes:
1069 *     -  A value of zero will be returned if this function is invoked with the
1070 *     global error status set, or if it should fail for any reason.
1071 
1072 *  Implementation Notes:
1073 *     This implementation addresses the matching of a DSBSpecFrame class
1074 *     object to any other class of Frame. A DSBSpecFrame will match any class
1075 *     of DSBSpecFrame (i.e. possibly from a derived class) but will not match
1076 *     a less specialised class of Frame (except for a SpecFrame).
1077 */
1078 
1079 /* Local Variables: */
1080    AstDSBSpecFrame *template;    /* Pointer to template DSBSpecFrame structure */
1081    AstFrame *frame0;             /* Pointer to Frame underlying axis 0 */
1082    int iaxis0;                   /* Axis index underlying axis 0 */
1083    int match;                    /* Coordinate conversion possible? */
1084 
1085 /* Initialise the returned values. */
1086    *template_axes = NULL;
1087    *target_axes = NULL;
1088    *map = NULL;
1089    *result = NULL;
1090    match = 0;
1091 
1092 /* Check the global error status. */
1093    if ( !astOK ) return match;
1094 
1095 /* Obtain a pointer to the template DSBSpecFrame structure. */
1096    template = (AstDSBSpecFrame *) template_frame;
1097 
1098 /* The first criterion for a match is that the template matches as a
1099    SpecFrame class object. This ensures that the number of axes (1) and
1100    domain, class, etc. of the target Frame are suitable. Invoke the parent
1101    "astMatch" method to verify this. */
1102    match = (*parent_match)( template_frame, target, matchsub,
1103                             template_axes, target_axes, map, result, status );
1104 
1105 /* If a match was found, the target Frame must be (or contain) a SpecFrame,
1106    but this target SpecFrame may be a simple SpecFrame rather than a
1107    DSBSpecFrame. We use the returned objects directly if the target
1108    SpecFrame is not a DSBSpecFrame. So if a DSBSpecFrame and a base
1109    SpecFrame are aligned, this will result in the DSBSpecFrame behaving as
1110    a normal SpecFrame. */
1111    if ( astOK && match ) {
1112 
1113 /* Get the primary Frame associated with the matching target axis. */
1114       astPrimaryFrame( target, (*target_axes)[ 0 ], &frame0, &iaxis0 );
1115 
1116 /* Skip this next section, thus retaining the values returned by the
1117    parent Match method above, if the target axis is not a DSBSpecFrame. */
1118       if( astIsADSBSpecFrame( frame0 ) ) {
1119 
1120 /* Annul the returned objects, which are not needed, but keep the axis
1121    association arrays which already hold the correct values. */
1122          *map = astAnnul( *map );
1123          *result = astAnnul( *result );
1124 
1125 /* Use the target's "astSubFrame" method to create a new Frame (the
1126    result Frame) with a copy of of the target axis. This process also
1127    overlays the template attributes on to the target Frame and returns a
1128    Mapping between the target and result Frames which effects the required
1129    coordinate conversion. */
1130          match = astSubFrame( target, template, 1, *target_axes, *template_axes,
1131                               map, result );
1132       }
1133 
1134 /* Free resources. */
1135       frame0 = astAnnul( frame0 );
1136 
1137    }
1138 
1139 /* If an error occurred, or conversion to the result Frame's coordinate
1140    system was not possible, then free all memory, annul the returned
1141    objects, and reset the returned value. */
1142    if ( !astOK || !match ) {
1143       if( *template_axes ) *template_axes = astFree( *template_axes );
1144       if( *target_axes ) *target_axes = astFree( *target_axes );
1145       if( *map ) *map = astAnnul( *map );
1146       if( *result ) *result = astAnnul( *result );
1147       match = 0;
1148    }
1149 
1150 /* Return the result. */
1151    return match;
1152 }
1153 
Overlay(AstFrame * template,const int * template_axes,AstFrame * result,int * status)1154 static void Overlay( AstFrame *template, const int *template_axes,
1155                      AstFrame *result, int *status ) {
1156 /*
1157 *  Name:
1158 *     Overlay
1159 
1160 *  Purpose:
1161 *     Overlay the attributes of a template DSBSpecFrame on to another Frame.
1162 
1163 *  Type:
1164 *     Private function.
1165 
1166 *  Synopsis:
1167 *     #include "specframe.h"
1168 *     void Overlay( AstFrame *template, const int *template_axes,
1169 *                   AstFrame *result, int *status )
1170 
1171 *  Class Membership:
1172 *     DSBSpecFrame member function (over-rides the protected astOverlay method
1173 *     inherited from the SpecFrame class).
1174 
1175 *  Description:
1176 *     This function overlays attributes of a DSBSpecFrame (the "template") on to
1177 *     another Frame, so as to over-ride selected attributes of that second
1178 *     Frame. Normally only those attributes which have been specifically set
1179 *     in the template will be transferred. This implements a form of
1180 *     defaulting, in which a Frame acquires attributes from the template, but
1181 *     retains its original attributes (as the default) if new values have not
1182 *     previously been explicitly set in the template.
1183 *
1184 *     Note that if the result Frame is a DSBSpecFrame and a change of spectral
1185 *     coordinate system occurs as a result of overlaying its System
1186 *     attribute, then some of its original attribute values may no
1187 *     longer be appropriate (e.g. the Title, or attributes describing
1188 *     its axes). In this case, these will be cleared before overlaying
1189 *     any new values.
1190 
1191 *  Parameters:
1192 *     template
1193 *        Pointer to the template DSBSpecFrame, for which values should have been
1194 *        explicitly set for any attribute which is to be transferred.
1195 *     template_axes
1196 *        Pointer to an array of int, with one element for each axis of the
1197 *        "result" Frame (see below). For each axis in the result frame, the
1198 *        corresponding element of this array should contain the (zero-based)
1199 *        index of the template axis to which it corresponds. This array is used
1200 *        to establish from which template axis any axis-dependent attributes
1201 *        should be obtained.
1202 *
1203 *        If any axis in the result Frame is not associated with a template
1204 *        axis, the corresponding element of this array should be set to -1.
1205 *
1206 *        If a NULL pointer is supplied, the template and result axis
1207 *        indicies are assumed to be identical.
1208 *     result
1209 *        Pointer to the Frame which is to receive the new attribute values.
1210 *     status
1211 *        Pointer to the inherited status variable.
1212 
1213 *  Returned Value:
1214 *     void
1215 
1216 *  Notes:
1217 *     -  In general, if the result Frame is not from the same class as the
1218 *     template DSBSpecFrame, or from a class derived from it, then attributes may
1219 *     exist in the template DSBSpecFrame which do not exist in the result Frame.
1220 *     In this case, these attributes will not be transferred.
1221 */
1222 
1223 /* Check the global error status. */
1224    if ( !astOK ) return;
1225 
1226 /* Invoke the parent class astOverlay method to transfer attributes inherited
1227    from the parent class. */
1228    (*parent_overlay)( template, template_axes, result, status );
1229 
1230 /* Check if the result Frame is a DSBSpecFrame or from a class derived from
1231    DSBSpecFrame. If not, we cannot transfer DSBSpecFrame attributes to it as it is
1232    insufficiently specialised. In this case simply omit these attributes. */
1233    if( astIsADSBSpecFrame( result ) && astOK ) {
1234 
1235 /* Define macros that test whether an attribute is set in the template and,
1236    if so, transfers its value to the result. */
1237 #define OVERLAY(attribute) \
1238    if ( astTest##attribute( template ) ) { \
1239       astSet##attribute( result, astGet##attribute( template ) ); \
1240    }
1241 
1242 /* Use the macro to transfer each DSBSpecFrame attribute in turn. */
1243       OVERLAY(DSBCentre)
1244       OVERLAY(IF)
1245       OVERLAY(SideBand)
1246       OVERLAY(AlignSideBand)
1247    }
1248 
1249 /* Undefine macros local to this function. */
1250 #undef OVERLAY
1251 }
1252 
SetAttrib(AstObject * this_object,const char * setting,int * status)1253 static void SetAttrib( AstObject *this_object, const char *setting, int *status ) {
1254 /*
1255 *  Name:
1256 *     SetAttrib
1257 
1258 *  Purpose:
1259 *     Set an attribute value for a DSBSpecFrame.
1260 
1261 *  Type:
1262 *     Private function.
1263 
1264 *  Synopsis:
1265 *     #include "dsbspecframe.h"
1266 *     void SetAttrib( AstObject *this, const char *setting )
1267 
1268 *  Class Membership:
1269 *     DSBSpecFrame member function (over-rides the astSetAttrib protected
1270 *     method inherited from the SpecFrame class).
1271 
1272 *  Description:
1273 *     This function assigns an attribute value for a DSBSpecFrame, the
1274 *     attribute and its value being specified by means of a string of
1275 *     the form:
1276 *
1277 *        "attribute= value "
1278 *
1279 *     Here, "attribute" specifies the attribute name and should be in
1280 *     lower case with no white space present. The value to the right
1281 *     of the "=" should be a suitable textual representation of the
1282 *     value to be assigned and this will be interpreted according to
1283 *     the attribute's data type.  White space surrounding the value is
1284 *     only significant for string attributes.
1285 
1286 *  Parameters:
1287 *     this
1288 *        Pointer to the DSBSpecFrame.
1289 *     setting
1290 *        Pointer to a null-terminated string specifying the new attribute
1291 *        value.
1292 */
1293 
1294 /* Local Variables: */
1295    AstDSBSpecFrame *this;        /* Pointer to the DSBSpecFrame structure */
1296    AstMapping *tmap;             /* Ptr to Mapping from this to topofreq */
1297    AstMapping *umap;             /* Ptr to Mapping between units */
1298    double dtemp;                 /* Attribute value */
1299    double dval;                  /* Attribute value */
1300    int ival;                     /* Attribute value */
1301    int len;                      /* Length of setting string */
1302    int nc;                       /* Used length */
1303    int off;                      /* Offset to start of string */
1304 
1305 /* Check the global error status. */
1306    if ( !astOK ) return;
1307 
1308 /* Obtain a pointer to the DSBSpecFrame structure. */
1309    this = (AstDSBSpecFrame *) this_object;
1310 
1311 /* Obtain the length of the setting string. */
1312    len = (int) strlen( setting );
1313 
1314 /* Test for each recognised attribute in turn, using "astSscanf" to parse the
1315    setting string and extract the attribute value (or an offset to it in the
1316    case of string values). In each case, use the value set in "nc" to check
1317    that the entire string was matched. Once a value has been obtained, use the
1318    appropriate method to set it. */
1319 
1320 /* DSBCentre */
1321 /* --------- */
1322    if ( strstr( setting, "dsbcentre=" ) ) {
1323 
1324 /* Without any units indication - assume it is supplied in the system of
1325    the DSBSpecFrame. */
1326       int ok = 0;
1327       if( nc = 0,
1328            ( 1 == astSscanf( setting, "dsbcentre= %lg %n", &dval, &nc ) )
1329            && ( nc >= len ) ) {
1330          ok = 1;
1331 
1332 /* With units indication. Is there a Mapping from the supplied units to the
1333    units used by the DSBSpecFrame? If so, use the Mapping to convert the
1334    supplied value to the required units. */
1335       } else if ( nc = 0,
1336            ( 1 == astSscanf( setting, "dsbcentre= %lg %n%*s %n", &dval, &off, &nc ) )
1337            && ( nc >= len ) ) {
1338 
1339          if( ( umap = astUnitMapper( setting + off, astGetUnit( this, 0 ), NULL, NULL ) ) ) {
1340             astTran1( umap, 1, &dval, 1, &dtemp );
1341             dval = dtemp;
1342             umap = astAnnul( umap );
1343             if( astOK && dval != AST__BAD ) ok = 1;
1344 
1345 /* Otherwise report an error. */
1346          } else if( astOK ) {
1347             astError( AST__ATTIN, "astSetAttrib(%s): Value supplied for "
1348                       "attribute \"DSBCentre\" (%s) uses units which are "
1349                       "inappropriate for the current spectral system (%s).", status,
1350                        astGetClass( this ), setting + 10,
1351                        astGetTitle( this ) );
1352          }
1353       }
1354 
1355 /* Convert the value from the supplied system to topocentric frequency in
1356    Hx, and store. */
1357       if( ok ) {
1358 
1359 /* Find the Mapping from the spectral system described by this SpecFrame to
1360    topocentric frequency in Hz. */
1361          tmap = TopoMap( this, 1, "astSetAttrib", status );
1362          if ( astOK ) {
1363 
1364 /* Transform the supplied value to topocentric frequency. */
1365             astTran1( tmap, 1, &dval, 1, &dtemp );
1366             if( dtemp == AST__BAD ) {
1367                astError( AST__ATTIN, "astSetAttrib(%s): The setting \"%s\" is "
1368                          "invalid for a %s.", status, astGetClass( this ), setting,
1369                          astGetClass( this ) );
1370             } else {
1371 
1372 /* Store it. */
1373                astSetDSBCentre( this, dtemp );
1374 
1375             }
1376             tmap = astAnnul( tmap );
1377          }
1378 
1379       } else if( astOK ) {
1380          astError( AST__ATTIN, "astSetAttrib(%s): The setting \"%s\" is "
1381                    "invalid for a %s.", status, astGetClass( this ), setting,
1382                    astGetClass( this ) );
1383       }
1384 
1385 /* IF */
1386 /* -- */
1387 /* Without any units indication - assume GHz. Convert to Hz for storage. */
1388    } else if ( nc = 0,
1389         ( 1 == astSscanf( setting, "if= %lg %n", &dval, &nc ) )
1390         && ( nc >= len ) ) {
1391       astSetIF( this, dval*1.0E9 );
1392 
1393 /* With units indication. */
1394    } else if ( nc = 0,
1395         ( 1 == astSscanf( setting, "if= %lg %n%*s %n", &dval, &off, &nc ) )
1396         && ( nc >= len ) ) {
1397 
1398 /* Is there a Mapping from the supplied units to Hz? If so, use the
1399    Mapping to convert the supplied value to Hz. */
1400       if( ( umap = astUnitMapper( setting + off, "Hz", NULL, NULL ) ) ) {
1401          astTran1( umap, 1, &dval, 1, &dtemp );
1402          umap = astAnnul( umap );
1403 
1404 /* Otherwise report an error. */
1405       } else if( astOK ) {
1406          astError( AST__ATTIN, "astSetAttrib(%s): Intermediate frequency given "
1407                    "in an inappropriate system of units \"%g %s\".", status,
1408                    astGetClass( this ), dval, setting + off );
1409       }
1410 
1411 /* Set the intermediate frequency. */
1412       astSetIF( this, dtemp );
1413 
1414 /* SideBand */
1415 /* -------- */
1416    } else if ( nc = 0,
1417                ( 0 == astSscanf( setting, "sideband= %n%*s %n", &ival, &nc ) )
1418                && ( nc >= len ) ) {
1419 
1420       if( astChrMatch( "usb", setting+ival ) ) {
1421          astSetSideBand( this, USB );
1422 
1423       } else if( astChrMatch( "lsb", setting+ival ) ) {
1424          astSetSideBand( this, LSB );
1425 
1426       } else if( astChrMatch( "lo", setting+ival ) ) {
1427          astSetSideBand( this, LO );
1428 
1429       } else if( astChrMatch( "observed", setting+ival ) ) {
1430          astSetSideBand( this, ( astGetIF( this ) > 0 ) ? LSB : USB );
1431 
1432       } else if( astChrMatch( "image", setting+ival ) ) {
1433          astSetSideBand( this, ( astGetIF( this ) <= 0 ) ? LSB : USB );
1434 
1435       } else {
1436          astError( AST__ATTIN, "astSetAttrib(%s): The setting \"%s\" is "
1437                    "invalid for a %s.", status, astGetClass( this ), setting,
1438                    astGetClass( this ) );
1439       }
1440 
1441 /* AlignSideBand */
1442 /* ------------- */
1443    } else if ( nc = 0,
1444         ( 1 == astSscanf( setting, "alignsideband= %d %n", &ival, &nc ) )
1445         && ( nc >= len ) ) {
1446       astSetAlignSideBand( this, ival );
1447 
1448 /* Read-only attributes. */
1449 /* --------------------- */
1450 /* Define a macro to see if the setting string matches any of the
1451    read-only attributes of this class. */
1452 #define MATCH(attrib) \
1453         ( nc = 0, ( 0 == astSscanf( setting, attrib "=%*[^\n]%n", &nc ) ) && \
1454                   ( nc >= len ) )
1455 
1456 /* Use this macro to report an error if a read-only attribute has been
1457    specified. */
1458    } else if ( MATCH( "imagfreq" ) ) {
1459       astError( AST__NOWRT, "astSet: The setting \"%s\" is invalid for a %s.", status,
1460                 setting, astGetClass( this ) );
1461       astError( AST__NOWRT, "This is a read-only attribute." , status);
1462 
1463 /* Pass any unrecognised setting to the parent method for further
1464    interpretation. */
1465    } else {
1466       (*parent_setattrib)( this_object, setting, status );
1467    }
1468 }
1469 
SubFrame(AstFrame * target_frame,AstFrame * template,int result_naxes,const int * target_axes,const int * template_axes,AstMapping ** map,AstFrame ** result,int * status)1470 static int SubFrame( AstFrame *target_frame, AstFrame *template,
1471                      int result_naxes, const int *target_axes,
1472                      const int *template_axes, AstMapping **map,
1473                      AstFrame **result, int *status ) {
1474 /*
1475 *  Name:
1476 *     SubFrame
1477 
1478 *  Purpose:
1479 *     Select axes from a DSBSpecFrame and convert to the new coordinate
1480 *     system.
1481 
1482 *  Type:
1483 *     Private function.
1484 
1485 *  Synopsis:
1486 *     #include "dsbspecframe.h"
1487 *     int SubFrame( AstFrame *target, AstFrame *template,
1488 *                   int result_naxes, const int *target_axes,
1489 *                   const int *template_axes, AstMapping **map,
1490 *                   AstFrame **result, int *status )
1491 
1492 *  Class Membership:
1493 *     DSBSpecFrame member function (over-rides the protected astSubFrame
1494 *     method inherited from the SpecFrame class).
1495 
1496 *  Description:
1497 *     This function selects a requested sub-set (or super-set) of the axes
1498 *     from a "target" DSBSpecFrame and creates a new Frame with copies of
1499 *     the selected axes assembled in the requested order. It then
1500 *     optionally overlays the attributes of a "template" Frame on to the
1501 *     result. It returns both the resulting Frame and a Mapping that
1502 *     describes how to convert between the coordinate systems described by
1503 *     the target and result Frames. If necessary, this Mapping takes
1504 *     account of any differences in the Frames' attributes due to the
1505 *     influence of the template.
1506 
1507 *  Parameters:
1508 *     target
1509 *        Pointer to the target DSBSpecFrame, from which axes are to be
1510 *        selected.
1511 *     template
1512 *        Pointer to the template Frame, from which new attributes for the
1513 *        result Frame are to be obtained. Optionally, this may be NULL, in
1514 *        which case no overlaying of template attributes will be performed.
1515 *     result_naxes
1516 *        Number of axes to be selected from the target Frame. This number may
1517 *        be greater than or less than the number of axes in this Frame (or
1518 *        equal).
1519 *     target_axes
1520 *        Pointer to an array of int with result_naxes elements, giving a list
1521 *        of the (zero-based) axis indices of the axes to be selected from the
1522 *        target DSBSpecFrame. The order in which these are given determines
1523 *        the order in which the axes appear in the result Frame. If any of the
1524 *        values in this array is set to -1, the corresponding result axis will
1525 *        not be derived from the target Frame, but will be assigned default
1526 *        attributes instead.
1527 *     template_axes
1528 *        Pointer to an array of int with result_naxes elements. This should
1529 *        contain a list of the template axes (given as zero-based axis indices)
1530 *        with which the axes of the result Frame are to be associated. This
1531 *        array determines which axes are used when overlaying axis-dependent
1532 *        attributes of the template on to the result. If any element of this
1533 *        array is set to -1, the corresponding result axis will not receive any
1534 *        template attributes.
1535 *
1536 *        If the template argument is given as NULL, this array is not used and
1537 *        a NULL pointer may also be supplied here.
1538 *     map
1539 *        Address of a location to receive a pointer to the returned Mapping.
1540 *        The forward transformation of this Mapping will describe how to
1541 *        convert coordinates from the coordinate system described by the target
1542 *        DSBSpecFrame to that described by the result Frame. The inverse
1543 *        transformation will convert in the opposite direction.
1544 *     result
1545 *        Address of a location to receive a pointer to the result Frame.
1546 *     status
1547 *        Pointer to the inherited status variable.
1548 
1549 *  Returned Value:
1550 *     A non-zero value is returned if coordinate conversion is possible
1551 *     between the target and the result Frame. Otherwise zero is returned and
1552 *     *map and *result are returned as NULL (but this will not in itself
1553 *     result in an error condition). In general, coordinate conversion should
1554 *     always be possible if no template Frame is supplied but may not always
1555 *     be possible otherwise.
1556 
1557 *  Notes:
1558 *     -  A value of zero will be returned if this function is invoked with the
1559 *     global error status set, or if it should fail for any reason.
1560 
1561 *  Implementation Notes:
1562 *     -  This implementation addresses the selection of axes from a
1563 *     DSBSpecFrame object. This results in another object of the same class
1564 *     only if the single DSBSpecFrame axis is selected exactly once.
1565 *     Otherwise, the result is a Frame class object which inherits the
1566 *     DSBSpecFrame's axis information (if appropriate) but none of the other
1567 *     properties of a DSBSpecFrame.
1568 *     -  In the event that a DSBSpecFrame results, the returned Mapping will
1569 *     take proper account of the relationship between the target and result
1570 *     coordinate systems.
1571 *     -  In the event that a Frame class object results, the returned Mapping
1572 *     will only represent a selection/permutation of axes.
1573 
1574 *  Implementation Deficiencies:
1575 *     -  Any axis selection is currently permitted. Probably this should be
1576 *     restricted so that each axis can only be selected once. The
1577 *     astValidateAxisSelection method will do this but currently there are bugs
1578 *     in the CmpFrame class that cause axis selections which will not pass this
1579 *     test. Install the validation when these are fixed.
1580 */
1581 
1582 /* Local Variables: */
1583    AstDSBSpecFrame *dsbresult;/* Pointer to the DSBSpecFrame result Frame */
1584    AstDSBSpecFrame *dsbtarget;/* Pointer to the DSBSpecFrame target Frame */
1585    AstMapping *map1;          /* Intermediate Mapping */
1586    AstMapping *map2;          /* Intermediate Mapping */
1587    AstMapping *map3;          /* Intermediate Mapping */
1588    int alignsb;               /* Use sidebands to align the Frames? */
1589    int match;                 /* Coordinate conversion is possible? */
1590    int obs_sb;                /* The observed sideband value */
1591    int old_sb;                /* The original Sideband value */
1592 
1593 /* Initialise the returned values. */
1594    *map = NULL;
1595    *result = NULL;
1596    match = 0;
1597 
1598 /* Check the global error status. */
1599    if ( !astOK ) return match;
1600 
1601 /* Invoke the astSubFrame method inherited from the parent SpecFrame
1602    class. This will (if possible) create a result Frame which is a
1603    DSBSpecFrame (since the supplied target Frame is a DSBSpecFrame).
1604    However, the Mapping from target to result Frame will take no account
1605    of any differences in the values of the attributes specific to the
1606    DSBSpecFrame class. */
1607    match = (*parent_subframe)( target_frame, template, result_naxes,
1608                                target_axes, template_axes, map, result, status );
1609 
1610 /* If a match occurred, and the result and template Frames are both
1611    DSBSpecFrames, we now modify the Mapping to take account of
1612    DSBSpecFrame-specific attributes. */
1613    if( match && template && astIsADSBSpecFrame( template ) &&
1614                             astIsADSBSpecFrame( *result ) ) {
1615 
1616 /* Get pointers to the two DSBSpecFrames */
1617       dsbtarget = (AstDSBSpecFrame *) target_frame;
1618 
1619 /* See whether alignment occurs between sidebands. If the current call to
1620    this function is part of the process of restoring a FrameSet's integrity
1621    following changes to the FrameSet's current Frame, then we ignore the
1622    setting of the AlignSideBand attributes and use 1. This ensures that
1623    when the SideBand attribute (for instance) is changed via a FrameSet
1624    pointer, the Mappings within the FrameSet are modified to produce
1625    frequencies in the new SideBand. In most other cases, astronomers
1626    usually want to align the DSBSpecFrames as if they were basic SpecFrames
1627    (that is, ignoring the setting of the SideBand attribute). */
1628       if( astGetFrameFlags( target_frame ) & AST__INTFLAG ) {
1629          alignsb = 1;
1630       } else {
1631          alignsb = astGetAlignSideBand( dsbtarget ) &&
1632                    astGetAlignSideBand( (AstDSBSpecFrame *) template );
1633       }
1634 
1635 /* If we are aligning the sidebands we need to modify the Mapping
1636    returned above by the parent SubFrame method. The existing Mapping
1637    will convert between the spectral systems represented by the two
1638    DSBSpecFrames but will not take account of any difference in
1639    sidebands. */
1640       if( alignsb ) {
1641 
1642 /* We assume that alignment occurs in the observed sideband. Determine
1643    which side band is the observed sideband in the target. */
1644          old_sb = astGetSideBand( dsbtarget );
1645          astSetC( dsbtarget, "SideBand", "observed" );
1646          obs_sb = astGetSideBand( dsbtarget );
1647          astSetSideBand( dsbtarget, old_sb );
1648 
1649 /* Create a Mapping which transforms positions from the target to an exact
1650    copy of the target in which the SideBand attribute is set to the
1651    observed (USB or LSB) sideband. This will be a UnitMap if the target
1652    already represents the observed sideband. */
1653          if( obs_sb == USB ) {
1654             map1 = ToUSBMapping( dsbtarget, "astSubFrame", status );
1655 
1656          } else if( obs_sb == LSB ) {
1657             map1 = ToLSBMapping( dsbtarget, "astSubFrame", status );
1658 
1659          } else {
1660             map1 = NULL;
1661             astError( AST__INTER, "astGetImagFreq(%s): Illegal sideband value "
1662                       "(%d) encountered (internal AST programming error).", status,
1663                       astGetClass( target_frame ), obs_sb );
1664          }
1665 
1666 /* Determine which side band is the observed sideband in the result. */
1667          dsbresult = (AstDSBSpecFrame *) *result;
1668          old_sb = astGetSideBand( dsbresult );
1669          astSetC( dsbresult, "SideBand", "observed" );
1670          obs_sb = astGetSideBand( dsbresult );
1671          astSetSideBand( dsbresult, old_sb );
1672 
1673 /* Create a Mapping which transforms positions from the result to an exact
1674    copy of the result in which the SideBand attribute is set to the
1675    obserfed sideband. This will be a UnitMap if the target already represents
1676    the observed sideband. */
1677          if( obs_sb == USB ) {
1678             map2 = ToUSBMapping( dsbresult, "astSubFrame", status );
1679 
1680          } else if( obs_sb == LSB ) {
1681             map2 = ToLSBMapping( dsbresult, "astSubFrame", status );
1682 
1683          } else {
1684             map2 = NULL;
1685             astError( AST__INTER, "astGetImagFreq(%s): Illegal sideband value "
1686                       "(%d) encountered (internal AST programming error).", status,
1687                       astGetClass( target_frame ), obs_sb );
1688          }
1689 
1690 /* Invert it to get the mapping from the observed sideband to the result. */
1691          astInvert( map2 );
1692 
1693 /* Form a Mapping which first maps target values to the observed sideband,
1694    then applies the Mapping returned by the parent SubFrame method in
1695    order to convert between spectral systems, and then converts from the
1696    observed sideband to the SideBand of the result. */
1697          map3 = (AstMapping *) astCmpMap( map1, *map, 1, "", status );
1698          map1 = astAnnul( map1 );
1699          *map = astAnnul( *map );
1700          map1 = (AstMapping *) astCmpMap( map3, map2, 1, "", status );
1701          map3 = astAnnul( map3 );
1702          map2 = astAnnul( map2 );
1703 
1704 /* Returned the simplified Mapping. */
1705          *map = astSimplify( map1 );
1706          map1 = astAnnul( map1 );
1707       }
1708    }
1709 
1710 /* If an error occurred or no match was found, annul the returned
1711    objects and reset the returned result. */
1712    if ( !astOK || !match ) {
1713       if( *map ) *map = astAnnul( *map );
1714       if( *result ) *result = astAnnul( *result );
1715       match = 0;
1716    }
1717 
1718 /* Return the result. */
1719    return match;
1720 
1721 }
1722 
TestAttrib(AstObject * this_object,const char * attrib,int * status)1723 static int TestAttrib( AstObject *this_object, const char *attrib, int *status ) {
1724 /*
1725 *  Name:
1726 *     TestAttrib
1727 
1728 *  Purpose:
1729 *     Test if a specified attribute value is set for a DSBSpecFrame.
1730 
1731 *  Type:
1732 *     Private function.
1733 
1734 *  Synopsis:
1735 *     #include "dsbspecframe.h"
1736 *     int TestAttrib( AstObject *this, const char *attrib, int *status )
1737 
1738 *  Class Membership:
1739 *     DSBSpecFrame member function (over-rides the astTestAttrib protected
1740 *     method inherited from the SpecFrame class).
1741 
1742 *  Description:
1743 *     This function returns a boolean result (0 or 1) to indicate whether
1744 *     a value has been set for one of a DSBSpecFrame's attributes.
1745 
1746 *  Parameters:
1747 *     this
1748 *        Pointer to the DSBSpecFrame.
1749 *     attrib
1750 *        Pointer to a null-terminated string specifying the attribute
1751 *        name.  This should be in lower case with no surrounding white
1752 *        space.
1753 *     status
1754 *        Pointer to the inherited status variable.
1755 
1756 *  Returned Value:
1757 *     One if a value has been set, otherwise zero.
1758 
1759 *  Notes:
1760 *     - A value of zero will be returned if this function is invoked
1761 *     with the global status set, or if it should fail for any reason.
1762 */
1763 
1764 /* Local Variables: */
1765    AstDSBSpecFrame *this;        /* Pointer to the DSBSpecFrame structure */
1766    int result;                   /* Result value to return */
1767 
1768 /* Initialise. */
1769    result = 0;
1770 
1771 /* Check the global error status. */
1772    if ( !astOK ) return result;
1773 
1774 /* Obtain a pointer to the DSBSpecFrame structure. */
1775    this = (AstDSBSpecFrame *) this_object;
1776 
1777 /* Check the attribute name and test the appropriate attribute. */
1778 
1779 /* DSBCentre */
1780 /* --------- */
1781    if ( !strcmp( attrib, "dsbcentre" ) ) {
1782       result = astTestDSBCentre( this );
1783 
1784 /* IF */
1785 /* -- */
1786    } else if ( !strcmp( attrib, "if" ) ) {
1787       result = astTestIF( this );
1788 
1789 /* SideBand */
1790 /* -------- */
1791    } else if ( !strcmp( attrib, "sideband" ) ) {
1792       result = astTestSideBand( this );
1793 
1794 /* AlignSideBand */
1795 /* ------------- */
1796    } else if ( !strcmp( attrib, "alignsideband" ) ) {
1797       result = astTestAlignSideBand( this );
1798 
1799 /* Read-only attributes. */
1800 /* --------------------- */
1801 /* Test if the attribute name matches any of the read-only attributes
1802    of this class. If it does, then return zero. */
1803    } else if ( !strcmp( attrib, "imagfreq" ) ) {
1804       result = 0;
1805 
1806 /* If the attribute is not recognised, pass it on to the parent method
1807    for further interpretation. */
1808    } else {
1809       result = (*parent_testattrib)( this_object, attrib, status );
1810    }
1811 
1812 /* Return the result, */
1813    return result;
1814 }
1815 
ToLOMapping(AstDSBSpecFrame * this,const char * method,int * status)1816 static AstMapping *ToLOMapping( AstDSBSpecFrame *this, const char *method, int *status ){
1817 /*
1818 *  Name:
1819 *     ToLOMapping
1820 
1821 *  Purpose:
1822 *     Create a Mapping which transforms a DSBSpecFrame to offset from the
1823 *     local oscillator frequency.
1824 
1825 *  Type:
1826 *     Private function.
1827 
1828 *  Synopsis:
1829 *     #include "dsbspecframe.h"
1830 *     AstMapping *ToLOMapping( AstDSBSpecFrame *this, const char *method, int *status )
1831 
1832 *  Class Membership:
1833 *     DSBSpecFrame member function
1834 
1835 *  Description:
1836 *     This function returns a pointer to a new Mapping which transforms
1837 *     positions in the supplied DSBSpecFrame into an offset from the local
1838 *     oscillator frequency. This will be a UnitMap if the DSBSpecFrame
1839 *     already represents offset from the local oscillator frequency.
1840 
1841 *  Parameters:
1842 *     this
1843 *        Pointer to the DSBSpecFrame.
1844 *     method
1845 *        Pointer to a null-terminated string containing the name of the
1846 *        public invoking method. This is only used in the construction of
1847 *        error messages.
1848 *     status
1849 *        Pointer to the inherited status variable.
1850 
1851 *  Returned Value:
1852 *     Pointer to a new Mapping.
1853 
1854 *  Notes:
1855 *     - A NULL pointer will be returned if this function is invoked
1856 *     with the global error status set, or if it should fail for any
1857 *     reason.
1858 */
1859 
1860 /* Local Variables: */
1861    AstMapping *fmap;         /* LSB to USB (topo freq) */
1862    AstMapping *map1;         /* This to USB (topo freq) */
1863    AstMapping *map2;         /* This (LSB) to This (USB) */
1864    AstMapping *result;       /* Pointer to the returned Mapping */
1865    AstMapping *tmap;         /* This to topocentric freq */
1866    double f_lo;              /* Local oscillator freq (topo Hz) */
1867    double f_in_a;            /* First LSB or USB freq */
1868    double f_in_b;            /* Second LSB or USB freq */
1869    double f_out_a;           /* First LO freq */
1870    double f_out_b;           /* Second LO freq */
1871    int sb;                   /* SideBand value */
1872 
1873 /* Initialise. */
1874    result = NULL;
1875 
1876 /* Check the global error status. */
1877    if ( !astOK ) return result;
1878 
1879 /* If the DSBSpecFrame already represents LO offset, return a UnitMap.*/
1880    sb = astGetSideBand( this );
1881    if( sb == LO ) {
1882       result = (AstMapping *) astUnitMap( 1, "", status );
1883 
1884 /* If the DSBSpecFrame represents the USB or LSB, create a suitable WinMap. */
1885    } else {
1886 
1887 /* Find the Mapping from the spectral system described by this SpecFrame to
1888    topocentric frequency in Hz. */
1889       tmap = TopoMap( this, 1, method, status );
1890 
1891 /* Calculate the local oscillator frequency (topocentric in Hertz). */
1892       f_lo = GetLO( this, "create a Mapping to upper sideband",
1893                     "astGetImagFreq", status );
1894 
1895 /* Create a 1D WinMap which converts f_in to f_out. */
1896       if( sb == LSB ) {
1897          f_in_a = 0.0;
1898          f_in_b = f_lo;
1899          f_out_a = f_lo;
1900          f_out_b = 0.0;
1901       } else {
1902          f_in_a = 0.0;
1903          f_in_b = -f_lo;
1904          f_out_a = f_lo;
1905          f_out_b = 0.0;
1906       }
1907 
1908       fmap = (AstMapping *) astWinMap( 1, &f_in_a, &f_in_b, &f_out_a, &f_out_b, "", status );
1909 
1910 /* Construct the Mapping: input to f_in, f_in to f_out, f_out to input */
1911       map1 = (AstMapping *) astCmpMap( tmap, fmap, 1, "", status );
1912       astInvert( tmap );
1913       map2 = (AstMapping *) astCmpMap( map1, tmap, 1, "", status );
1914 
1915 /* Simplify */
1916       result = astSimplify( map2 );
1917 
1918 /* Free resources */
1919       tmap = astAnnul( tmap );
1920       fmap = astAnnul( fmap );
1921       map1 = astAnnul( map1 );
1922       map2 = astAnnul( map2 );
1923    }
1924 
1925 /* Return NULL if an error has occurred. */
1926    if( !astOK ) result = astAnnul( result );
1927 
1928 /* Return the result. */
1929    return result;
1930 
1931 }
1932 
ToLSBMapping(AstDSBSpecFrame * this,const char * method,int * status)1933 static AstMapping *ToLSBMapping( AstDSBSpecFrame *this, const char *method, int *status ){
1934 /*
1935 *  Name:
1936 *     ToLSBMapping
1937 
1938 *  Purpose:
1939 *     Create a Mapping which transforms a DSBSpecFrame to the lower
1940 *     sideband.
1941 
1942 *  Type:
1943 *     Private function.
1944 
1945 *  Synopsis:
1946 *     #include "dsbspecframe.h"
1947 *     AstMapping *ToLSBMapping( AstDSBSpecFrame *this, const char *method, int *status )
1948 
1949 *  Class Membership:
1950 *     DSBSpecFrame member function
1951 
1952 *  Description:
1953 *     This function returns a pointer to a new Mapping which transforms
1954 *     positions in the supplied DSBSpecFrame to the lower sideband. This
1955 *     will be a UnitMap if the DSBSpecFrame already represents the lower
1956 *     sideband.
1957 
1958 *  Parameters:
1959 *     this
1960 *        Pointer to the DSBSpecFrame.
1961 *     method
1962 *        Pointer to a null-terminated string containing the name of the
1963 *        public invoking method. This is only used in the construction of
1964 *        error messages.
1965 *     status
1966 *        Pointer to the inherited status variable.
1967 
1968 *  Returned Value:
1969 *     Pointer to a new Mapping.
1970 
1971 *  Notes:
1972 *     - A NULL pointer will be returned if this function is invoked
1973 *     with the global error status set, or if it should fail for any
1974 *     reason.
1975 */
1976 
1977 /* Local Variables: */
1978    AstMapping *fmap;         /* LSB to USB (topo freq) */
1979    AstMapping *map1;         /* This to USB (topo freq) */
1980    AstMapping *map2;         /* This (LSB) to This (USB) */
1981    AstMapping *result;       /* Pointer to the returned Mapping */
1982    AstMapping *tmap;         /* This to topocentric freq */
1983    double f_lo;              /* Local oscillator freq (topo Hz) */
1984    double f_out_a;           /* First LSB freq */
1985    double f_out_b;           /* Second LSB freq */
1986    double f_in_a;            /* First USB or LO freq */
1987    double f_in_b;            /* Second USB or LO freq */
1988    int sb;                   /* SideBand value */
1989 
1990 /* Initialise. */
1991    result = NULL;
1992 
1993 /* Check the global error status. */
1994    if ( !astOK ) return result;
1995 
1996 /* If the DSBSpecFrame already represents the LSB, return a UnitMap.*/
1997    sb = astGetSideBand( this );
1998    if( sb == LSB ) {
1999       result = (AstMapping *) astUnitMap( 1, "", status );
2000 
2001 /* If the DSBSpecFrame represents the USB or LO offset, create a suitable
2002    WinMap. */
2003    } else {
2004 
2005 /* Find the Mapping from the spectral system described by this SpecFrame to
2006    topocentric frequency in Hz. */
2007       tmap = TopoMap( this, 1, method, status );
2008 
2009 /* Calculate the local oscillator frequency (topocentric in Hertz). */
2010       f_lo = GetLO( this, "create a Mapping to lower sideband",
2011                     "astGetImagFreq", status );
2012 
2013 /* Create a 1D WinMap which converts USB or LO to LSB. */
2014       if( sb == USB ) {
2015          f_in_a = 0.0;
2016          f_in_b = 2*f_lo;
2017          f_out_a = 2*f_lo;
2018          f_out_b = 0.0;
2019       } else {
2020          f_in_a = 0.0;
2021          f_in_b = f_lo;
2022          f_out_a = f_lo;
2023          f_out_b = 0.0;
2024       }
2025 
2026       fmap = (AstMapping *) astWinMap( 1, &f_in_a, &f_in_b, &f_out_a, &f_out_b, "", status );
2027 
2028 /* Construct the Mapping: input to f_in, f_in to f_out, f_out to input */
2029       map1 = (AstMapping *) astCmpMap( tmap, fmap, 1, "", status );
2030       astInvert( tmap );
2031       map2 = (AstMapping *) astCmpMap( map1, tmap, 1, "", status );
2032 
2033 /* Simplify */
2034       result = astSimplify( map2 );
2035 
2036 /* Free resources */
2037       tmap = astAnnul( tmap );
2038       fmap = astAnnul( fmap );
2039       map1 = astAnnul( map1 );
2040       map2 = astAnnul( map2 );
2041    }
2042 
2043 /* Return NULL if an error has occurred. */
2044    if( !astOK ) result = astAnnul( result );
2045 
2046 /* Return the result. */
2047    return result;
2048 
2049 }
2050 
TopoMap(AstDSBSpecFrame * this,int forward,const char * method,int * status)2051 static AstMapping *TopoMap( AstDSBSpecFrame *this, int forward,
2052                             const char *method, int *status ){
2053 /*
2054 *  Name:
2055 *     TopoMap
2056 
2057 *  Purpose:
2058 *     Create a Mapping which transforms a DSBSpecFrame to topocentric
2059 *     frequency.
2060 
2061 *  Type:
2062 *     Private function.
2063 
2064 *  Synopsis:
2065 *     #include "dsbspecframe.h"
2066 *     AstMapping *TopoMap( AstDSBSpecFrame *this, int forward,
2067 *                          const char *method, int *status )
2068 
2069 *  Class Membership:
2070 *     DSBSpecFrame member function
2071 
2072 *  Description:
2073 *     This function returns a pointer to a new Mapping which transforms
2074 *     positions in the supplied DSBSpecFrame to the corresponding
2075 *     topocentric frequency values in Hz (or the inverse of this).
2076 
2077 *  Parameters:
2078 *     this
2079 *        Pointer to the DSBSpecFrame.
2080 *     forward
2081 *        If zero, the calcuated Mapping is inverted before being returned.
2082 *     method
2083 *        Pointer to a null-terminated string containing the name of the
2084 *        public invoking method. This is only used in the construction of
2085 *        error messages.
2086 *     status
2087 *        Pointer to the inherited status variable.
2088 
2089 *  Returned Value:
2090 *     Pointer to a new Mapping.
2091 
2092 *  Notes:
2093 *     - A NULL pointer will be returned if this function is invoked
2094 *     with the global error status set, or if it should fail for any
2095 *     reason.
2096 */
2097 
2098 /* Local Variables: */
2099    AstMapping *result;       /* The returned Mapping */
2100    AstFrameSet *fs;          /* FrameSet connecting tf1 and tf2 */
2101    AstSpecFrame *tf1;        /* SpecFrame corresponding to this DSBSpecFrame */
2102    AstSpecFrame *tf2;        /* Topocentric frequency SpecFrame */
2103    int template_axis;        /* The axis to overlay */
2104 
2105 /* Initialise. */
2106    result = NULL;
2107 
2108 /* Check the global error status. */
2109    if ( !astOK ) return result;
2110 
2111 /* Make a SpecFrame and then overlay the SpecFrame attributes of this
2112    DSBSpecFrame onto the new SpecFrame. This means it inherits the current
2113    values of things like ObsLon and ObsLat. */
2114    tf1 = astSpecFrame( "", status );
2115    template_axis = 0;
2116    (*parent_overlay)( (AstFrame *) this, &template_axis, (AstFrame *) tf1, status );
2117 
2118 /* Copy this new SpecFrame and set its attributes to describe topocentric
2119    frequency in Hz. Ensure that alignment occurs in the topocentric Frame. */
2120    astSetAlignStdOfRest( tf1, AST__TPSOR);
2121    tf2 = astCopy( tf1 );
2122    astSetSystem( tf2, AST__FREQ );
2123    astSetStdOfRest( tf2, AST__TPSOR );
2124    astSetUnit( tf2, 0, "Hz" );
2125 
2126 /* Find the Mapping from the spectral system described by this SpecFrame to
2127    topocentric frequency in Hz. */
2128    fs = astConvert( tf1, tf2, "" );
2129    if ( astOK ) {
2130       if( !fs ) {
2131          astError( AST__INTER, "%s(%s): Cannot convert DSBCentre "
2132                    "value from the supplied system to topocentric frequency "
2133                    "(internal AST programming error).", status, method,
2134                    astGetClass( this ) );
2135       } else {
2136          result = astGetMapping( fs, AST__BASE, AST__CURRENT );
2137          if( !forward ) astInvert( result );
2138       }
2139       fs = astAnnul( fs );
2140    }
2141 
2142 /* Free resources */
2143    tf1 = astAnnul( tf1 );
2144    tf2 = astAnnul( tf2 );
2145 
2146 /* Annul the result if an error has occurred. */
2147    if( !astOK ) result = astAnnul( result );
2148 
2149 /* Return the result. */
2150    return result;
2151 
2152 }
2153 
ToUSBMapping(AstDSBSpecFrame * this,const char * method,int * status)2154 static AstMapping *ToUSBMapping( AstDSBSpecFrame *this, const char *method, int *status ){
2155 /*
2156 *  Name:
2157 *     ToUSBMapping
2158 
2159 *  Purpose:
2160 *     Create a Mapping which transforms a DSBSpecFrame to the upper
2161 *     sideband.
2162 
2163 *  Type:
2164 *     Private function.
2165 
2166 *  Synopsis:
2167 *     #include "dsbspecframe.h"
2168 *     AstMapping *ToUSBMapping( AstDSBSpecFrame *this, const char *method, int *status )
2169 
2170 *  Class Membership:
2171 *     DSBSpecFrame member function
2172 
2173 *  Description:
2174 *     This function returns a pointer to a new Mapping which transforms
2175 *     positions in the supplied DSBSpecFrame to the upper sideband. This
2176 *     will be a UnitMap if the DSBSpecFrame already represents the upper
2177 *     sideband.
2178 
2179 *  Parameters:
2180 *     this
2181 *        Pointer to the DSBSpecFrame.
2182 *     method
2183 *        Pointer to a null-terminated string containing the name of the
2184 *        public invoking method. This is only used in the construction of
2185 *        error messages.
2186 *     status
2187 *        Pointer to the inherited status variable.
2188 
2189 *  Returned Value:
2190 *     Pointer to a new Mapping.
2191 
2192 *  Notes:
2193 *     - A NULL pointer will be returned if this function is invoked
2194 *     with the global error status set, or if it should fail for any
2195 *     reason.
2196 */
2197 
2198 /* Local Variables: */
2199    AstMapping *fmap;         /* LSB to USB (topo freq) */
2200    AstMapping *map1;         /* This to USB (topo freq) */
2201    AstMapping *map2;         /* This (LSB) to This (USB) */
2202    AstMapping *result;       /* Pointer to the returned Mapping */
2203    AstMapping *tmap;         /* This to topocentric freq */
2204    double f_lo;              /* Local oscillator freq (topo Hz) */
2205    double f_in_a;            /* First LSB or LO freq */
2206    double f_in_b;            /* Second LSB or LO freq */
2207    double f_out_a;           /* First USB freq */
2208    double f_out_b;           /* Second USB freq */
2209    int sb;                   /* SideBand value */
2210 
2211 /* Initialise. */
2212    result = NULL;
2213 
2214 /* Check the global error status. */
2215    if ( !astOK ) return result;
2216 
2217 /* If the DSBSpecFrame already represents the USB, return a UnitMap.*/
2218    sb = astGetSideBand( this );
2219    if( sb == USB ) {
2220       result = (AstMapping *) astUnitMap( 1, "", status );
2221 
2222 /* If the DSBSpecFrame represents the LSB, or LO offset, create a suitable
2223    WinMap. */
2224    } else {
2225 
2226 /* Find the Mapping from the spectral system described by this SpecFrame to
2227    topocentric frequency in Hz. */
2228       tmap = TopoMap( this, 1, method, status );
2229 
2230 /* Calculate the local oscillator frequency (topocentric in Hertz). */
2231       f_lo = GetLO( this, "create a Mapping to upper sideband",
2232                     "astGetImagFreq", status );
2233 
2234 /* Create a 1D WinMap which converts f_in to f_out. */
2235       if( sb == LSB ) {
2236          f_in_a = 0.0;
2237          f_in_b = 2*f_lo;
2238          f_out_a = 2*f_lo;
2239          f_out_b = 0.0;
2240       } else {
2241          f_in_a = 0.0;
2242          f_in_b = -f_lo;
2243          f_out_a = f_lo;
2244          f_out_b = 0.0;
2245       }
2246 
2247       fmap = (AstMapping *) astWinMap( 1, &f_in_a, &f_in_b, &f_out_a, &f_out_b, "", status );
2248 
2249 /* Construct the Mapping: input to f_in, f_in to f_out, f_out to input */
2250       map1 = (AstMapping *) astCmpMap( tmap, fmap, 1, "", status );
2251       astInvert( tmap );
2252       map2 = (AstMapping *) astCmpMap( map1, tmap, 1, "", status );
2253 
2254 /* Simplify */
2255       result = astSimplify( map2 );
2256 
2257 /* Free resources */
2258       tmap = astAnnul( tmap );
2259       fmap = astAnnul( fmap );
2260       map1 = astAnnul( map1 );
2261       map2 = astAnnul( map2 );
2262    }
2263 
2264 /* Return NULL if an error has occurred. */
2265    if( !astOK ) result = astAnnul( result );
2266 
2267 /* Return the result. */
2268    return result;
2269 
2270 }
2271 
VerifyAttrs(AstDSBSpecFrame * this,const char * purp,const char * attrs,const char * method,int * status)2272 static void VerifyAttrs( AstDSBSpecFrame *this, const char *purp,
2273                          const char *attrs, const char *method, int *status ) {
2274 /*
2275 *  Name:
2276 *     VerifyAttrs
2277 
2278 *  Purpose:
2279 *     Verify that usable attribute values are available.
2280 
2281 *  Type:
2282 *     Private function.
2283 
2284 *  Synopsis:
2285 *     #include "dsbspecframe.h"
2286 *     void VerifyAttrs( AstDSBSpecFrame *this, const char *purp,
2287 *                       const char *attrs, const char *method, int *status  )
2288 
2289 *  Class Membership:
2290 *     DSBSpecFrame member function
2291 
2292 *  Description:
2293 *     This function tests each attribute listed in "attrs". It returns
2294 *     without action if 1) an explicit value has been set for each attribute
2295 *     or 2) the UseDefs attribute of the supplied DSBSpecFrame is non-zero.
2296 *
2297 *     If UseDefs is zero (indicating that default values should not be
2298 *     used for attributes), and any of the named attributes does not have
2299 *     an explicitly set value, then an error is reported.
2300 
2301 *  Parameters:
2302 *     this
2303 *        Pointer to the DSBSpecFrame.
2304 *     purp
2305 *        Pointer to a text string containing a message which will be
2306 *        included in any error report. This shouldindicate the purpose
2307 *        for which the attribute value is required.
2308 *     attrs
2309 *        A string holding a space separated list of attribute names.
2310 *     method
2311 *        A string holding the name of the calling method for use in error
2312 *        messages.
2313 *     status
2314 *        Pointer to the inherited status variable.
2315 
2316 */
2317 
2318 /* Local Variables: */
2319    const char *a;
2320    const char *desc;
2321    const char *p;
2322    int len;
2323    int set;
2324    int state;
2325 
2326 /* Check inherited status */
2327    if( !astOK ) return;
2328 
2329 /* If the DSBSpecFrame has a non-zero value for its UseDefs attribute, then
2330    all attributes are assumed to have usable values, since the defaults
2331    will be used if no explicit value has been set. So we only need to do
2332    any checks if UseDefs is zero. */
2333    if( !astGetUseDefs( this ) ) {
2334 
2335 /* Initialise variables to avoid compiler warnings. */
2336       a = NULL;
2337       desc = NULL;
2338       len = 0;
2339       set = 0;
2340 
2341 /* Loop round the "attrs" string identifying the start and length of each
2342    non-blank word in the string. */
2343       state = 0;
2344       p = attrs;
2345       while( 1 ) {
2346          if( state == 0 ) {
2347             if( !isspace( *p ) ) {
2348                a = p;
2349                len = 1;
2350                state = 1;
2351             }
2352          } else {
2353             if( isspace( *p ) || !*p ) {
2354 
2355 /* The end of a word has just been reached. Compare it to each known
2356    attribute value. Get a flag indicating if the attribute has a set
2357    value, and a string describing the attribute.*/
2358                if( len > 0 ) {
2359 
2360                   if( !strncmp( "DSBCentre", a, len ) ) {
2361                      set = astTestDSBCentre( this );
2362                      desc = "central position of interest";
2363 
2364                   } else if( !strncmp( "IF", a, len ) ) {
2365                      set = astTestIF( this );
2366                      desc = "intermediate frequency";
2367 
2368                   } else {
2369                      astError( AST__INTER, "VerifyAttrs(DSBSpecFrame): "
2370                                "Unknown attribute name \"%.*s\" supplied (AST "
2371                                "internal programming error).", status, len, a );
2372                   }
2373 
2374 /* If the attribute does not have a set value, report an error. */
2375                   if( !set && astOK ) {
2376                      astError( AST__NOVAL, "%s(%s): Cannot %s.", status, method,
2377                                astGetClass( this ), purp );
2378                      astError( AST__NOVAL, "No value has been set for "
2379                                "the AST \"%.*s\" attribute (%s).", status, len, a,
2380                                desc );
2381                   }
2382 
2383 /* Continue the word search algorithm. */
2384                }
2385                len = 0;
2386                state = 0;
2387             } else {
2388                len++;
2389             }
2390          }
2391          if( !*(p++) ) break;
2392       }
2393    }
2394 }
2395 
2396 
2397 /* Functions which access class attributes. */
2398 /* ---------------------------------------- */
2399 /* Implement member functions to access the attributes associated with
2400    this class using the macros defined for this purpose in the
2401    "object.h" file. For a description of each attribute, see the class
2402    interface (in the associated .h file). */
2403 
2404 /*
2405 *att++
2406 *  Name:
2407 *     ImagFreq
2408 
2409 *  Purpose:
2410 *     The image sideband equivalent of the rest frequency.
2411 
2412 *  Type:
2413 *     Public attribute.
2414 
2415 *  Synopsis:
2416 *     Floating point, read-only.
2417 
2418 *  Description:
2419 *     This is a read-only attribute giving the frequency which
2420 *     corresponds to the rest frequency but is in the opposite sideband.
2421 
2422 *     The value is calculated by first transforming the rest frequency
2423 *     (given by the RestFreq attribute) from the standard of rest of the
2424 *     source (given by the SourceVel and SourceVRF attributes) to the
2425 *     standard of rest of the observer (i.e. the topocentric standard of
2426 *     rest). The resulting topocentric frequency is assumed to be in the
2427 *     same sideband as the value given for the DSBCentre attribute (the
2428 *     "observed" sideband), and is transformed to the other sideband (the
2429 *     "image" sideband). The new frequency is converted back to the standard
2430 *     of rest of the source, and the resulting value is returned as the
2431 *     attribute value, in units of GHz.
2432 
2433 *  Applicability:
2434 *     DSBSpecFrame
2435 *        All DSBSpecFrames have this attribute.
2436 
2437 *att--
2438 */
2439 
2440 
2441 /*
2442 *att++
2443 *  Name:
2444 *     DSBCentre
2445 
2446 *  Purpose:
2447 *     The central position of interest in a dual sideband spectrum.
2448 
2449 *  Type:
2450 *     Public attribute.
2451 
2452 *  Synopsis:
2453 *     Floating point.
2454 
2455 *  Description:
2456 *     This attribute specifies the central position of interest in a dual
2457 *     sideband spectrum. Its sole use is to determine the local oscillator
2458 *     frequency (the frequency which marks the boundary between the lower
2459 *     and upper sidebands). See the description of the IF (intermediate
2460 *     frequency) attribute for details of how the local oscillator frequency
2461 *     is calculated. The sideband containing this central position is
2462 *     referred to as the "observed" sideband, and the other sideband as
2463 *     the "image" sideband.
2464 *
2465 *     The value is accessed as a position in the spectral system
2466 *     represented by the SpecFrame attributes inherited by this class, but
2467 *     is stored internally as topocentric frequency. Thus, if the System
2468 *     attribute of the DSBSpecFrame is set to "VRAD", the Unit attribute
2469 *     set to "m/s" and the StdOfRest attribute set to "LSRK", then values
2470 *     for the DSBCentre attribute should be supplied as radio velocity in
2471 *     units of "m/s" relative to the kinematic LSR (alternative units may
2472 *     be used by appending a suitable units string to the end of the value).
2473 *     This value is then converted to topocentric frequency and stored. If
2474 *     (say) the Unit attribute is subsequently changed to "km/s" before
2475 *     retrieving the current value of the DSBCentre attribute, the stored
2476 *     topocentric frequency will be converted back to LSRK radio velocity,
2477 *     this time in units of "km/s", before being returned.
2478 *
2479 *     The default value for this attribute is 30 GHz.
2480 
2481 *  Applicability:
2482 *     DSBSpecFrame
2483 *        All DSBSpecFrames have this attribute.
2484 
2485 *  Note:
2486 *     - The attributes which define the transformation to or from topocentric
2487 *     frequency should be assigned their correct values before accessing
2488 *     this attribute. These potentially include System, Unit, StdOfRest,
2489 *     ObsLon, ObsLat, ObsAlt, Epoch, RefRA, RefDec and RestFreq.
2490 
2491 *att--
2492 */
2493 /* The central frequency (topocentric frequency in Hz). */
astMAKE_CLEAR(DSBSpecFrame,DSBCentre,dsbcentre,AST__BAD)2494 astMAKE_CLEAR(DSBSpecFrame,DSBCentre,dsbcentre,AST__BAD)
2495 astMAKE_GET(DSBSpecFrame,DSBCentre,double,3.0E10,((this->dsbcentre!=AST__BAD)?this->dsbcentre:3.0E10))
2496 astMAKE_SET(DSBSpecFrame,DSBCentre,double,dsbcentre,value)
2497 astMAKE_TEST(DSBSpecFrame,DSBCentre,( this->dsbcentre != AST__BAD ))
2498 
2499 
2500 /*
2501 *att++
2502 *  Name:
2503 *     IF
2504 
2505 *  Purpose:
2506 *     The intermediate frequency in a dual sideband spectrum.
2507 
2508 *  Type:
2509 *     Public attribute.
2510 
2511 *  Synopsis:
2512 *     Floating point.
2513 
2514 *  Description:
2515 *     This attribute specifies the (topocentric) intermediate frequency in
2516 *     a dual sideband spectrum. Its sole use is to determine the local
2517 *     oscillator (LO) frequency (the frequency which marks the boundary
2518 *     between the lower and upper sidebands). The LO frequency is
2519 *     equal to the sum of the centre frequency and the intermediate
2520 *     frequency. Here, the "centre frequency" is the topocentric
2521 *     frequency in Hz corresponding to the current value of the DSBCentre
2522 *     attribute. The value of the IF attribute may be positive or
2523 *     negative: a positive value results in the LO frequency being above
2524 *     the central frequency, whilst a negative IF value results in the LO
2525 *     frequency being below the central frequency. The sign of the IF
2526 *     attribute value determines the default value for the SideBand
2527 *     attribute.
2528 *
2529 *     When setting a new value for this attribute, the units in which the
2530 *     frequency value is supplied may be indicated by appending a suitable
2531 *     string to the end of the formatted value. If the units are not
2532 *     specified, then the supplied value is assumed to be in units of GHz.
2533 *     For instance, the following strings all result in an IF of 4 GHz being
2534 *     used: "4.0", "4.0 GHz", "4.0E9 Hz", etc.
2535 *
2536 *     When getting the value of this attribute, the returned value is
2537 *     always in units of GHz. The default value for this attribute is 4 GHz.
2538 
2539 *  Applicability:
2540 *     DSBSpecFrame
2541 *        All DSBSpecFrames have this attribute.
2542 
2543 *att--
2544 */
2545 /* The intermediate frequency (topocentric in Hz). */
2546 astMAKE_CLEAR(DSBSpecFrame,IF,ifr,AST__BAD)
2547 astMAKE_GET(DSBSpecFrame,IF,double,4.0E9,((this->ifr!=AST__BAD)?this->ifr:4.0E9))
2548 astMAKE_SET(DSBSpecFrame,IF,double,ifr,value)
2549 astMAKE_TEST(DSBSpecFrame,IF,( this->ifr != AST__BAD ))
2550 
2551 /*
2552 *att++
2553 *  Name:
2554 *     SideBand
2555 
2556 *  Purpose:
2557 *     Indicates which sideband a dual sideband spectrum represents.
2558 
2559 *  Type:
2560 *     Public attribute.
2561 
2562 *  Synopsis:
2563 *     String.
2564 
2565 *  Description:
2566 *     This attribute indicates whether the DSBSpecFrame currently
2567 *     represents its lower or upper sideband, or an offset from the local
2568 *     oscillator frequency. When querying the current value, the returned
2569 *     string is always one of "usb" (for upper sideband), "lsb" (for lower
2570 *     sideband), or "lo" (for offset from the local oscillator frequency).
2571 *     When setting a new value, any of the strings "lsb", "usb", "observed",
2572 *     "image" or "lo" may be supplied (case insensitive). The "observed"
2573 *     sideband is which ever sideband (upper or lower) contains the central
2574 *     spectral position given by attribute DSBCentre, and the "image"
2575 *     sideband is the other sideband. It is the sign of the IF attribute
2576 *     which determines if the observed sideband is the upper or lower
2577 *     sideband. The default value for SideBand is the observed sideband.
2578 
2579 *  Applicability:
2580 *     DSBSpecFrame
2581 *        All DSBSpecFrames have this attribute.
2582 
2583 *att--
2584 */
2585 /* Protected access to the SideBand attribute uses BADSB to indicate
2586    "unset". Other negative values mean "LSB", zero means "LO" and
2587    positive values mean "USB". */
2588 astMAKE_CLEAR(DSBSpecFrame,SideBand,sideband,BADSB)
2589 astMAKE_SET(DSBSpecFrame,SideBand,int,sideband,((value<0)?LSB:((value==0)?LO:USB)))
2590 astMAKE_TEST(DSBSpecFrame,SideBand,( this->sideband != BADSB ))
2591 astMAKE_GET(DSBSpecFrame,SideBand,int,USB,(this->sideband == BADSB ? ((astGetIF( this )>0)?LSB:USB):this->sideband))
2592 
2593 /*
2594 *att++
2595 *  Name:
2596 *     AlignSideBand
2597 
2598 *  Purpose:
2599 *     Should the SideBand attribute be taken into account when aligning
2600 *     this DSBSpecFrame with another DSBSpecFrame?
2601 
2602 *  Type:
2603 *     Public attribute.
2604 
2605 *  Synopsis:
2606 *     Integer (boolean).
2607 
2608 *  Description:
2609 *     This attribute controls how a DSBSpecFrame behaves when an attempt
2610 *     is made to align it with another DSBSpecFrame using
2611 c     astFindFrame or astConvert.
2612 f     AST_FINDFRAME or AST_CONVERT.
2613 *     If both DSBSpecFrames have a non-zero value for AlignSideBand, the
2614 *     value of the SideBand attribute in each DSBSpecFrame is used so that
2615 *     alignment occurs between sidebands. That is, if one DSBSpecFrame
2616 *     represents USB and the other represents LSB then
2617 c     astFindFrame and astConvert
2618 f     AST_FINDFRAME and AST_CONVERT
2619 *     will recognise that the DSBSpecFrames represent different sidebands
2620 *     and will take this into account when constructing the Mapping that
2621 *     maps positions in one DSBSpecFrame into the other. If AlignSideBand
2622 *     in either DSBSpecFrame is set to zero, then the values of the SideBand
2623 *     attributes are ignored. In the above example, this would result in a
2624 *     frequency in the first DSBSpecFrame being mapped onto the same
2625 *     frequency in the second DSBSpecFrame, even though those frequencies
2626 *     refer to different sidebands. In other words, if either AlignSideBand
2627 *     attribute is zero, then the two DSBSpecFrames aligns like basic
2628 *     SpecFrames. The default value for AlignSideBand is zero.
2629 *
2630 c     When astFindFrame or astConvert
2631 f     When AST_FINDFRAME or AST_CONVERT
2632 *     is used on two DSBSpecFrames (potentially describing different spectral
2633 *     coordinate systems and/or sidebands), it returns a Mapping which can be
2634 *     used to transform a position in one DSBSpecFrame into the corresponding
2635 *     position in the other. The Mapping is made up of the following steps in
2636 *     the indicated order:
2637 *
2638 *     - If both DSBSpecFrames have a value of 1 for the AlignSideBand
2639 *     attribute, map values from the target's current sideband (given by its
2640 *     SideBand attribute) to the observed sideband (whether USB or LSB). If
2641 *     the target already represents the observed sideband, this step will
2642 *     leave the values unchanged. If either of the two DSBSpecFrames have a
2643 *     value of zero for its AlignSideBand attribute, then this step is omitted.
2644 *
2645 *     - Map the values from the spectral system of the target to the spectral
2646 *     system of the template. This Mapping takes into account all the
2647 *     inherited SpecFrame attributes such as System, StdOfRest, Unit, etc.
2648 *
2649 *     - If both DSBSpecFrames have a value of 1 for the AlignSideBand
2650 *     attribute, map values from the result's observed sideband to the
2651 *     result's current sideband (given by its SideBand attribute). If the
2652 *     result already represents the observed sideband, this step will leave
2653 *     the values unchanged. If either of the two DSBSpecFrames have a value
2654 *     of zero for its AlignSideBand attribute, then this step is omitted.
2655 
2656 *  Applicability:
2657 *     DSBSpecFrame
2658 *        All DSBSpecFrames have this attribute.
2659 
2660 *att--
2661 */
2662 /* The AlignSideBand value has a value of -1 when not set yielding a
2663    default of 0. */
2664 astMAKE_TEST(DSBSpecFrame,AlignSideBand,( this->alignsideband != -1 ))
2665 astMAKE_CLEAR(DSBSpecFrame,AlignSideBand,alignsideband,-1)
2666 astMAKE_GET(DSBSpecFrame,AlignSideBand,int,-1,((this->alignsideband==-1)?0:this->alignsideband) )
2667 astMAKE_SET(DSBSpecFrame,AlignSideBand,int,alignsideband,(value?1:0))
2668 
2669 /* Copy constructor. */
2670 /* ----------------- */
2671 /* None needed */
2672 
2673 /* Destructor. */
2674 /* ----------- */
2675 /* None needed */
2676 
2677 /* Dump function. */
2678 /* -------------- */
2679 static void Dump( AstObject *this_object, AstChannel *channel, int *status ) {
2680 /*
2681 *  Name:
2682 *     Dump
2683 
2684 *  Purpose:
2685 *     Dump function for DSBSpecFrame objects.
2686 
2687 *  Type:
2688 *     Private function.
2689 
2690 *  Synopsis:
2691 *     void Dump( AstObject *this, AstChannel *channel, int *status )
2692 
2693 *  Description:
2694 *     This function implements the Dump function which writes out data
2695 *     for the DSBSpecFrame class to an output Channel.
2696 
2697 *  Parameters:
2698 *     this
2699 *        Pointer to the DSBSpecFrame whose data are being written.
2700 *     channel
2701 *        Pointer to the Channel to which the data are being written.
2702 *     status
2703 *        Pointer to the inherited status variable.
2704 */
2705 
2706 /* Local Variables: */
2707    AstDSBSpecFrame *this;        /* Pointer to the DSBSpecFrame structure */
2708    const char *cval;             /* Attribute value */
2709    const char *comm;             /* Attribute comment */
2710    double dval;                  /* Attribute value */
2711    int ival;                     /* Attribute value */
2712    int set;                      /* Is attribute set? */
2713 
2714 /* Check the global error status. */
2715    if ( !astOK ) return;
2716 
2717 /* Obtain a pointer to the DSBSpecFrame structure. */
2718    this = (AstDSBSpecFrame *) this_object;
2719 
2720 /* In the case of attributes, we first use the appropriate (private)
2721    Test...  member function to see if they are set. If so, we then use
2722    the (private) Get... function to obtain the value to be written
2723    out.
2724 
2725    For attributes which are not set, we use the astGet... method to
2726    obtain the value instead. This will supply a default value
2727    (possibly provided by a derived class which over-rides this method)
2728    which is more useful to a human reader as it corresponds to the
2729    actual default attribute value.  Since "set" will be zero, these
2730    values are for information only and will not be read back. */
2731 
2732 /* DSBCentre */
2733 /* --------- */
2734    set = TestDSBCentre( this, status );
2735    dval = set ? GetDSBCentre( this, status ) : astGetDSBCentre( this );
2736    astWriteDouble( channel, "DSBCen", set, 1, dval, "Central frequency (Hz topo)" );
2737 
2738 /* IF */
2739 /* -- */
2740    set = TestIF( this, status );
2741    dval = set ? GetIF( this, status ) : astGetIF( this );
2742    astWriteDouble( channel, "IF", set, 1, dval, "Intermediate frequency (Hz)" );
2743 
2744 /* SideBand */
2745 /* -------- */
2746    set = TestSideBand( this, status );
2747    ival = set ? GetSideBand( this, status ) : astGetSideBand( this );
2748    if( ival == LSB ) {
2749       cval = "LSB";
2750       comm = "Represents lower sideband";
2751 
2752    } else if( ival == LO ) {
2753       cval = "LO";
2754       comm = "Represents offset from LO frequency";
2755 
2756    } else {
2757       cval = "USB";
2758       comm = "Represents upper sideband";
2759    }
2760    astWriteString( channel, "SideBn", set, 1, cval, comm );
2761 
2762 /* AlignSideBand */
2763 /* ------------- */
2764    set = TestAlignSideBand( this, status );
2765    ival = set ? GetAlignSideBand( this, status ) : astGetAlignSideBand( this );
2766    astWriteInt( channel, "AlSdBn", set, 1, ival, "Align sidebands?" );
2767 }
2768 
2769 /* Standard class functions. */
2770 /* ========================= */
2771 /* Implement the astIsADSBSpecFrame and astCheckDSBSpecFrame functions using the macros
2772    defined for this purpose in the "object.h" header file. */
astMAKE_ISA(DSBSpecFrame,SpecFrame)2773 astMAKE_ISA(DSBSpecFrame,SpecFrame)
2774 astMAKE_CHECK(DSBSpecFrame)
2775 
2776 AstDSBSpecFrame *astDSBSpecFrame_( const char *options, int *status, ...) {
2777 /*
2778 *++
2779 *  Name:
2780 c     astDSBSpecFrame
2781 f     AST_DSBSPECFRAME
2782 
2783 *  Purpose:
2784 *     Create a DSBSpecFrame.
2785 
2786 *  Type:
2787 *     Public function.
2788 
2789 *  Synopsis:
2790 c     #include "dsbspecframe.h"
2791 c     AstDSBSpecFrame *astDSBSpecFrame( const char *options, ... )
2792 f     RESULT = AST_DSBSPECFRAME( OPTIONS, STATUS )
2793 
2794 *  Class Membership:
2795 *     DSBSpecFrame constructor.
2796 
2797 *  Description:
2798 *     This function creates a new DSBSpecFrame and optionally initialises its
2799 *     attributes.
2800 *
2801 *     A DSBSpecFrame is a specialised form of SpecFrame which represents
2802 *     positions in a spectrum obtained using a dual sideband instrument.
2803 *     Such an instrument produces a spectrum in which each point contains
2804 *     contributions from two distinctly different frequencies, one from
2805 *     the "lower side band" (LSB) and one from the "upper side band" (USB).
2806 *     Corresponding LSB and USB frequencies are connected by the fact
2807 *     that they are an equal distance on either side of a fixed central
2808 *     frequency known as the "Local Oscillator" (LO) frequency.
2809 *
2810 *     When quoting a position within such a spectrum, it is necessary to
2811 *     indicate whether the quoted position is the USB position or the
2812 *     corresponding LSB position. The SideBand attribute provides this
2813 *     indication. Another option that the SideBand attribute provides is
2814 *     to represent a spectral position by its topocentric offset from the
2815 *     LO frequency.
2816 *
2817 *     In practice, the LO frequency is specified by giving the distance
2818 *     from the LO frequency to some "central" spectral position. Typically
2819 *     this central position is that of some interesting spectral feature.
2820 *     The distance from this central position to the LO frequency is known
2821 *     as the "intermediate frequency" (IF). The value supplied for IF can
2822 *     be a signed value in order to indicate whether the LO frequency is
2823 *     above or below the central position.
2824 
2825 *  Parameters:
2826 c     options
2827 f     OPTIONS = CHARACTER * ( * ) (Given)
2828 c        Pointer to a null-terminated string containing an optional
2829 c        comma-separated list of attribute assignments to be used for
2830 c        initialising the new DSBSpecFrame. The syntax used is identical to
2831 c        that for the astSet function and may include "printf" format
2832 c        specifiers identified by "%" symbols in the normal way.
2833 f        A character string containing an optional comma-separated
2834 f        list of attribute assignments to be used for initialising the
2835 f        new DSBSpecFrame. The syntax used is identical to that for the
2836 f        AST_SET routine.
2837 c     ...
2838 c        If the "options" string contains "%" format specifiers, then
2839 c        an optional list of additional arguments may follow it in
2840 c        order to supply values to be substituted for these
2841 c        specifiers. The rules for supplying these are identical to
2842 c        those for the astSet function (and for the C "printf"
2843 c        function).
2844 f     STATUS = INTEGER (Given and Returned)
2845 f        The global status.
2846 
2847 *  Returned Value:
2848 c     astDSBSpecFrame()
2849 f     AST_DSBSPECFRAME = INTEGER
2850 *        A pointer to the new DSBSpecFrame.
2851 
2852 *  Notes:
2853 *     - A null Object pointer (AST__NULL) will be returned if this
2854 c     function is invoked with the AST error status set, or if it
2855 f     function is invoked with STATUS set to an error value, or if it
2856 *     should fail for any reason.
2857 *--
2858 */
2859 
2860 /* Local Variables: */
2861    astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
2862    AstDSBSpecFrame *new;        /* Pointer to new DSBSpecFrame */
2863    va_list args;                /* Variable argument list */
2864 
2865 /* Get a pointer to the thread specific global data structure. */
2866    astGET_GLOBALS(NULL);
2867 
2868 /* Check the global status. */
2869    if ( !astOK ) return NULL;
2870 
2871 /* Initialise the DSBSpecFrame, allocating memory and initialising the
2872    virtual function table as well if necessary. */
2873    new = astInitDSBSpecFrame( NULL, sizeof( AstDSBSpecFrame ), !class_init, &class_vtab,
2874                               "DSBSpecFrame" );
2875 
2876 /* If successful, note that the virtual function table has been
2877    initialised. */
2878    if ( astOK ) {
2879       class_init = 1;
2880 
2881 /* Obtain the variable argument list and pass it along with the options string
2882    to the astVSet method to initialise the new DSBSpecFrame's attributes. */
2883       va_start( args, status );
2884       astVSet( new, options, NULL, args );
2885       va_end( args );
2886 
2887 /* If an error occurred, clean up by deleting the new object. */
2888       if ( !astOK ) new = astDelete( new );
2889    }
2890 
2891 /* Return a pointer to the new DSBSpecFrame. */
2892    return new;
2893 }
2894 
astDSBSpecFrameId_(const char * options,...)2895 AstDSBSpecFrame *astDSBSpecFrameId_( const char *options, ... ) {
2896 /*
2897 *  Name:
2898 *     astDSBSpecFrameId_
2899 
2900 *  Purpose:
2901 *     Create a DSBSpecFrame.
2902 
2903 *  Type:
2904 *     Private function.
2905 
2906 *  Synopsis:
2907 *     #include "dsbspecframe.h"
2908 *     AstDSBSpecFrame *astDSBSpecFrameId_( const char *options, ... )
2909 
2910 *  Class Membership:
2911 *     DSBSpecFrame constructor.
2912 
2913 *  Description:
2914 *     This function implements the external (public) interface to the
2915 *     astDSBSpecFrame constructor function. It returns an ID value (instead
2916 *     of a true C pointer) to external users, and must be provided
2917 *     because astDSBSpecFrame_ has a variable argument list which cannot be
2918 *     encapsulated in a macro (where this conversion would otherwise
2919 *     occur).
2920 *
2921 *     The variable argument list also prevents this function from
2922 *     invoking astDSBSpecFrame_ directly, so it must be a re-implementation
2923 *     of it in all respects, except for the final conversion of the
2924 *     result to an ID value.
2925 
2926 *  Parameters:
2927 *     As for astDSBSpecFrame_.
2928 
2929 *  Returned Value:
2930 *     The ID value associated with the new DSBSpecFrame.
2931 */
2932 
2933 /* Local Variables: */
2934    astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
2935    AstDSBSpecFrame *new;        /* Pointer to new DSBSpecFrame */
2936    va_list args;                /* Variable argument list */
2937 
2938    int *status;                  /* Pointer to inherited status value */
2939 
2940 /* Get a pointer to the inherited status value. */
2941    status = astGetStatusPtr;
2942 
2943 /* Get a pointer to the thread specific global data structure. */
2944    astGET_GLOBALS(NULL);
2945 
2946 /* Check the global status. */
2947    if ( !astOK ) return NULL;
2948 
2949 /* Initialise the DSBSpecFrame, allocating memory and initialising the
2950    virtual function table as well if necessary. */
2951    new = astInitDSBSpecFrame( NULL, sizeof( AstDSBSpecFrame ), !class_init, &class_vtab,
2952                               "DSBSpecFrame" );
2953 
2954 /* If successful, note that the virtual function table has been
2955    initialised. */
2956    if ( astOK ) {
2957       class_init = 1;
2958 
2959 /* Obtain the variable argument list and pass it along with the options string
2960    to the astVSet method to initialise the new DSBSpecFrame's attributes. */
2961       va_start( args, options );
2962       astVSet( new, options, NULL, args );
2963       va_end( args );
2964 
2965 /* If an error occurred, clean up by deleting the new object. */
2966       if ( !astOK ) new = astDelete( new );
2967    }
2968 
2969 /* Return an ID value for the new DSBSpecFrame. */
2970    return astMakeId( new );
2971 }
2972 
astInitDSBSpecFrame_(void * mem,size_t size,int init,AstDSBSpecFrameVtab * vtab,const char * name,int * status)2973 AstDSBSpecFrame *astInitDSBSpecFrame_( void *mem, size_t size, int init,
2974                                        AstDSBSpecFrameVtab *vtab, const char *name, int *status ) {
2975 /*
2976 *+
2977 *  Name:
2978 *     astInitDSBSpecFrame
2979 
2980 *  Purpose:
2981 *     Initialise a DSBSpecFrame.
2982 
2983 *  Type:
2984 *     Protected function.
2985 
2986 *  Synopsis:
2987 *     #include "dsbspecframe.h"
2988 *     AstDSBSpecFrame *astInitDSBSpecFrame( void *mem, size_t size, int init,
2989 *                               AstDSBSpecFrameVtab *vtab, const char *name )
2990 
2991 *  Class Membership:
2992 *     DSBSpecFrame initialiser.
2993 
2994 *  Description:
2995 *     This function is provided for use by class implementations to initialise
2996 *     a new DSBSpecFrame object. It allocates memory (if necessary) to accommodate
2997 *     the DSBSpecFrame plus any additional data associated with the derived class.
2998 *     It then initialises a DSBSpecFrame structure at the start of this memory. If
2999 *     the "init" flag is set, it also initialises the contents of a virtual
3000 *     function table for a DSBSpecFrame at the start of the memory passed via the
3001 *     "vtab" parameter.
3002 
3003 *  Parameters:
3004 *     mem
3005 *        A pointer to the memory in which the DSBSpecFrame is to be initialised.
3006 *        This must be of sufficient size to accommodate the DSBSpecFrame data
3007 *        (sizeof(DSBSpecFrame)) plus any data used by the derived class. If a value
3008 *        of NULL is given, this function will allocate the memory itself using
3009 *        the "size" parameter to determine its size.
3010 *     size
3011 *        The amount of memory used by the DSBSpecFrame (plus derived class data).
3012 *        This will be used to allocate memory if a value of NULL is given for
3013 *        the "mem" parameter. This value is also stored in the DSBSpecFrame
3014 *        structure, so a valid value must be supplied even if not required for
3015 *        allocating memory.
3016 *     init
3017 *        A logical flag indicating if the DSBSpecFrame's virtual function table is
3018 *        to be initialised. If this value is non-zero, the virtual function
3019 *        table will be initialised by this function.
3020 *     vtab
3021 *        Pointer to the start of the virtual function table to be associated
3022 *        with the new DSBSpecFrame.
3023 *     name
3024 *        Pointer to a constant null-terminated character string which contains
3025 *        the name of the class to which the new object belongs (it is this
3026 *        pointer value that will subsequently be returned by the astGetClass
3027 *        method).
3028 
3029 *  Returned Value:
3030 *     A pointer to the new DSBSpecFrame.
3031 
3032 *  Notes:
3033 *     -  A null pointer will be returned if this function is invoked with the
3034 *     global error status set, or if it should fail for any reason.
3035 *-
3036 */
3037 
3038 /* Local Variables: */
3039    AstDSBSpecFrame *new;        /* Pointer to new DSBSpecFrame */
3040 
3041 /* Check the global status. */
3042    if ( !astOK ) return NULL;
3043 
3044 /* If necessary, initialise the virtual function table. */
3045    if ( init ) astInitDSBSpecFrameVtab( vtab, name );
3046 
3047 /* Initialise. */
3048    new = NULL;
3049 
3050 /* Initialise a SpecFrame structure (the parent class) as the first component
3051    within the DSBSpecFrame structure, allocating memory if necessary. Specify that
3052    the SpecFrame should be defined in both the forward and inverse directions. */
3053    new = (AstDSBSpecFrame *) astInitSpecFrame( mem, size, 0,
3054                                        (AstSpecFrameVtab *) vtab, name );
3055    if ( astOK ) {
3056 
3057 /* Initialise the DSBSpecFrame data. */
3058 /* --------------------------------- */
3059       new->dsbcentre = AST__BAD;
3060       new->ifr = AST__BAD;
3061       new->sideband = BADSB;
3062       new->alignsideband = -1;
3063 
3064 /* If an error occurred, clean up by deleting the new DSBSpecFrame. */
3065       if ( !astOK ) new = astDelete( new );
3066    }
3067 
3068 /* Return a pointer to the new DSBSpecFrame. */
3069    return new;
3070 }
3071 
astLoadDSBSpecFrame_(void * mem,size_t size,AstDSBSpecFrameVtab * vtab,const char * name,AstChannel * channel,int * status)3072 AstDSBSpecFrame *astLoadDSBSpecFrame_( void *mem, size_t size,
3073                            AstDSBSpecFrameVtab *vtab, const char *name,
3074                            AstChannel *channel, int *status ) {
3075 /*
3076 *+
3077 *  Name:
3078 *     astLoadDSBSpecFrame
3079 
3080 *  Purpose:
3081 *     Load a DSBSpecFrame.
3082 
3083 *  Type:
3084 *     Protected function.
3085 
3086 *  Synopsis:
3087 *     #include "dsbspecframe.h"
3088 *     AstDSBSpecFrame *astLoadDSBSpecFrame( void *mem, size_t size,
3089 *                               AstDSBSpecFrameVtab *vtab, const char *name,
3090 *                               AstChannel *channel )
3091 
3092 *  Class Membership:
3093 *     DSBSpecFrame loader.
3094 
3095 *  Description:
3096 *     This function is provided to load a new DSBSpecFrame using data read
3097 *     from a Channel. It first loads the data used by the parent class
3098 *     (which allocates memory if necessary) and then initialises a
3099 *     DSBSpecFrame structure in this memory, using data read from the input
3100 *     Channel.
3101 *
3102 *     If the "init" flag is set, it also initialises the contents of a
3103 *     virtual function table for a DSBSpecFrame at the start of the memory
3104 *     passed via the "vtab" parameter.
3105 
3106 
3107 *  Parameters:
3108 *     mem
3109 *        A pointer to the memory into which the DSBSpecFrame is to be
3110 *        loaded.  This must be of sufficient size to accommodate the
3111 *        DSBSpecFrame data (sizeof(DSBSpecFrame)) plus any data used by derived
3112 *        classes. If a value of NULL is given, this function will
3113 *        allocate the memory itself using the "size" parameter to
3114 *        determine its size.
3115 *     size
3116 *        The amount of memory used by the DSBSpecFrame (plus derived class
3117 *        data).  This will be used to allocate memory if a value of
3118 *        NULL is given for the "mem" parameter. This value is also
3119 *        stored in the DSBSpecFrame structure, so a valid value must be
3120 *        supplied even if not required for allocating memory.
3121 *
3122 *        If the "vtab" parameter is NULL, the "size" value is ignored
3123 *        and sizeof(AstDSBSpecFrame) is used instead.
3124 *     vtab
3125 *        Pointer to the start of the virtual function table to be
3126 *        associated with the new DSBSpecFrame. If this is NULL, a pointer
3127 *        to the (static) virtual function table for the DSBSpecFrame class
3128 *        is used instead.
3129 *     name
3130 *        Pointer to a constant null-terminated character string which
3131 *        contains the name of the class to which the new object
3132 *        belongs (it is this pointer value that will subsequently be
3133 *        returned by the astGetClass method).
3134 *
3135 *        If the "vtab" parameter is NULL, the "name" value is ignored
3136 *        and a pointer to the string "DSBSpecFrame" is used instead.
3137 
3138 *  Returned Value:
3139 *     A pointer to the new DSBSpecFrame.
3140 
3141 *  Notes:
3142 *     - A null pointer will be returned if this function is invoked
3143 *     with the global error status set, or if it should fail for any
3144 *     reason.
3145 *-
3146 */
3147 
3148 /* Local Constants. */
3149    astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
3150 #define KEY_LEN 50               /* Maximum length of a keyword */
3151 
3152 /* Local Variables: */
3153    AstDSBSpecFrame *new;         /* Pointer to the new DSBSpecFrame */
3154    char *text;                   /* Pointer to string value */
3155 
3156 /* Get a pointer to the thread specific global data structure. */
3157    astGET_GLOBALS(channel);
3158 
3159 /* Initialise. */
3160    new = NULL;
3161 
3162 /* Check the global error status. */
3163    if ( !astOK ) return new;
3164 
3165 /* If a NULL virtual function table has been supplied, then this is
3166    the first loader to be invoked for this DSBSpecFrame. In this case the
3167    DSBSpecFrame belongs to this class, so supply appropriate values to be
3168    passed to the parent class loader (and its parent, etc.). */
3169    if ( !vtab ) {
3170       size = sizeof( AstDSBSpecFrame );
3171       vtab = &class_vtab;
3172       name = "DSBSpecFrame";
3173 
3174 /* If required, initialise the virtual function table for this class. */
3175       if ( !class_init ) {
3176          astInitDSBSpecFrameVtab( vtab, name );
3177          class_init = 1;
3178       }
3179    }
3180 
3181 /* Invoke the parent class loader to load data for all the ancestral
3182    classes of the current one, returning a pointer to the resulting
3183    partly-built DSBSpecFrame. */
3184    new = astLoadSpecFrame( mem, size, (AstSpecFrameVtab *) vtab, name,
3185                            channel );
3186    if ( astOK ) {
3187 
3188 /* Read input data. */
3189 /* ================ */
3190 /* Request the input Channel to read all the input data appropriate to
3191    this class into the internal "values list". */
3192        astReadClassData( channel, "DSBSpecFrame" );
3193 
3194 /* Now read each individual data item from this list and use it to
3195    initialise the appropriate instance variable(s) for this class. */
3196 
3197 /* In the case of attributes, we first read the "raw" input value,
3198    supplying the "unset" value as the default. If a "set" value is
3199    obtained, we then use the appropriate (private) Set... member
3200    function to validate and set the value properly. */
3201 
3202 /* DSBCentre */
3203 /* --------- */
3204       new->dsbcentre = astReadDouble( channel, "dsbcen", AST__BAD );
3205       if ( TestDSBCentre( new, status ) ) SetDSBCentre( new, new->dsbcentre, status );
3206 
3207 /* IF */
3208 /* -- */
3209       new->ifr = astReadDouble( channel, "if", AST__BAD );
3210       if ( TestIF( new, status ) ) SetIF( new, new->ifr, status );
3211 
3212 /* SideBand */
3213 /* -------- */
3214       text = astReadString( channel, "sidebn", " " );
3215       if( astOK ) {
3216          if( !strcmp( text, " " ) ) {
3217             new->sideband = BADSB;
3218          } else if( !strcmp( text, "USB" ) ) {
3219             new->sideband = USB;
3220          } else if( !strcmp( text, "LSB" ) ) {
3221             new->sideband = LSB;
3222          } else if( !strcmp( text, "LO" ) ) {
3223             new->sideband = LO;
3224          } else {
3225             astError( AST__ATTIN, "astRead(%s): Invalid SideBand description "
3226                       "\"%s\".", status, astGetClass( channel ), text );
3227          }
3228          if ( TestSideBand( new, status ) ) SetSideBand( new, new->sideband, status );
3229          text = astFree( text );
3230       }
3231 
3232 /* AlignSideBand */
3233 /* ------------- */
3234       new->alignsideband = astReadInt( channel, "alsdbn", -1 );
3235       if( TestAlignSideBand( new, status ) ) SetAlignSideBand( new, new->alignsideband, status );
3236 
3237 /* If an error occurred, clean up by deleting the new DSBSpecFrame. */
3238       if ( !astOK ) new = astDelete( new );
3239    }
3240 
3241 /* Return the new DSBSpecFrame pointer. */
3242    return new;
3243 }
3244 
3245 /* Virtual function interfaces. */
3246 /* ============================ */
3247 /* These provide the external interface to the virtual functions defined by
3248    this class. Each simply checks the global error status and then locates and
3249    executes the appropriate member function, using the function pointer stored
3250    in the object's virtual function table (this pointer is located using the
3251    astMEMBER macro defined in "object.h").
3252 
3253    Note that the member function may not be the one defined here, as it may
3254    have been over-ridden by a derived class. However, it should still have the
3255    same interface. */
3256 
astGetImagFreq_(AstDSBSpecFrame * this,int * status)3257 double astGetImagFreq_( AstDSBSpecFrame *this, int *status ) {
3258    if ( !astOK ) return AST__BAD;
3259    return (**astMEMBER(this,DSBSpecFrame,GetImagFreq))( this, status );
3260 }
3261 
3262 
3263 
3264 
3265 
3266 
3267