1 /* zzspkpa0.f -- translated by f2c (version 19980913).
2    You must link the resulting object file with the libraries:
3 	-lf2c -lm   (in that order)
4 */
5 
6 #include "f2c.h"
7 
8 /* Table of constant values */
9 
10 static integer c__0 = 0;
11 static integer c__9 = 9;
12 
13 /* $Procedure ZZSPKPA0 ( S/P Kernel, apparent position only ) */
zzspkpa0_(integer * targ,doublereal * et,char * ref,doublereal * sobs,char * abcorr,doublereal * ptarg,doublereal * lt,ftnlen ref_len,ftnlen abcorr_len)14 /* Subroutine */ int zzspkpa0_(integer *targ, doublereal *et, char *ref,
15 	doublereal *sobs, char *abcorr, doublereal *ptarg, doublereal *lt,
16 	ftnlen ref_len, ftnlen abcorr_len)
17 {
18     /* Initialized data */
19 
20     static logical first = TRUE_;
21     static char flags[5*9] = "NONE " "LT   " "LT+S " "CN   " "CN+S " "XLT  "
22 	    "XLT+S" "XCN  " "XCN+S";
23     static char prvcor[5] = "     ";
24 
25     /* System generated locals */
26     integer i__1;
27     doublereal d__1;
28 
29     /* Builtin functions */
30     integer s_cmp(char *, char *, ftnlen, ftnlen);
31     /* Subroutine */ int s_copy(char *, char *, ftnlen, ftnlen);
32 
33     /* Local variables */
34     char corr[5];
35     extern /* Subroutine */ int zzspkgp0_(integer *, doublereal *, char *,
36 	    integer *, doublereal *, doublereal *, ftnlen), vsub_(doublereal *
37 	    , doublereal *, doublereal *);
38     static logical xmit;
39     extern /* Subroutine */ int vequ_(doublereal *, doublereal *);
40     doublereal tpos[3];
41     integer i__, refid;
42     extern /* Subroutine */ int chkin_(char *, ftnlen), errch_(char *, char *,
43 	     ftnlen, ftnlen);
44     static logical usecn, uselt;
45     extern doublereal vnorm_(doublereal *);
46     extern logical failed_(void);
47     extern doublereal clight_(void);
48     extern integer isrchc_(char *, integer *, char *, ftnlen, ftnlen);
49     extern /* Subroutine */ int stelab_(doublereal *, doublereal *,
50 	    doublereal *), sigerr_(char *, ftnlen), chkout_(char *, ftnlen),
51 	    stlabx_(doublereal *, doublereal *, doublereal *);
52     integer ltsign;
53     extern /* Subroutine */ int ljucrs_(integer *, char *, char *, ftnlen,
54 	    ftnlen), setmsg_(char *, ftnlen);
55     integer maxitr;
56     extern /* Subroutine */ int irfnum_(char *, integer *, ftnlen);
57     extern logical return_(void);
58     static logical usestl;
59     extern logical odd_(integer *);
60 
61 /* $ Abstract */
62 
63 /*     Return the position of a target body relative to an observer, */
64 /*     optionally corrected for light time and stellar aberration. */
65 
66 /* $ Disclaimer */
67 
68 /*     THIS SOFTWARE AND ANY RELATED MATERIALS WERE CREATED BY THE */
69 /*     CALIFORNIA INSTITUTE OF TECHNOLOGY (CALTECH) UNDER A U.S. */
70 /*     GOVERNMENT CONTRACT WITH THE NATIONAL AERONAUTICS AND SPACE */
71 /*     ADMINISTRATION (NASA). THE SOFTWARE IS TECHNOLOGY AND SOFTWARE */
72 /*     PUBLICLY AVAILABLE UNDER U.S. EXPORT LAWS AND IS PROVIDED "AS-IS" */
73 /*     TO THE RECIPIENT WITHOUT WARRANTY OF ANY KIND, INCLUDING ANY */
74 /*     WARRANTIES OF PERFORMANCE OR MERCHANTABILITY OR FITNESS FOR A */
75 /*     PARTICULAR USE OR PURPOSE (AS SET FORTH IN UNITED STATES UCC */
76 /*     SECTIONS 2312-2313) OR FOR ANY PURPOSE WHATSOEVER, FOR THE */
77 /*     SOFTWARE AND RELATED MATERIALS, HOWEVER USED. */
78 
79 /*     IN NO EVENT SHALL CALTECH, ITS JET PROPULSION LABORATORY, OR NASA */
80 /*     BE LIABLE FOR ANY DAMAGES AND/OR COSTS, INCLUDING, BUT NOT */
81 /*     LIMITED TO, INCIDENTAL OR CONSEQUENTIAL DAMAGES OF ANY KIND, */
82 /*     INCLUDING ECONOMIC DAMAGE OR INJURY TO PROPERTY AND LOST PROFITS, */
83 /*     REGARDLESS OF WHETHER CALTECH, JPL, OR NASA BE ADVISED, HAVE */
84 /*     REASON TO KNOW, OR, IN FACT, SHALL KNOW OF THE POSSIBILITY. */
85 
86 /*     RECIPIENT BEARS ALL RISK RELATING TO QUALITY AND PERFORMANCE OF */
87 /*     THE SOFTWARE AND ANY RELATED MATERIALS, AND AGREES TO INDEMNIFY */
88 /*     CALTECH AND NASA FOR ALL THIRD-PARTY CLAIMS RESULTING FROM THE */
89 /*     ACTIONS OF RECIPIENT IN THE USE OF THE SOFTWARE. */
90 
91 /* $ Required_Reading */
92 
93 /*     SPK */
94 
95 /* $ Keywords */
96 
97 /*     EPHEMERIS */
98 
99 /* $ Declarations */
100 /* $ Brief_I/O */
101 
102 /*     Variable  I/O  Description */
103 /*     --------  ---  -------------------------------------------------- */
104 /*     TARG       I   Target body. */
105 /*     ET         I   Observer epoch. */
106 /*     REF        I   Inertial reference frame of observer's state. */
107 /*     SOBS       I   State of observer wrt. solar system barycenter. */
108 /*     ABCORR     I   Aberration correction flag. */
109 /*     PTARG      O   Position of target. */
110 /*     LT         O   One way light time between observer and target. */
111 
112 /* $ Detailed_Input */
113 
114 /*     TARG        is the NAIF ID code for a target body.  The target */
115 /*                 and observer define a position vector which points */
116 /*                 from the observer to the target. */
117 
118 /*     ET          is the ephemeris time, expressed as seconds past */
119 /*                 J2000 TDB, at which the position of the target body */
120 /*                 relative to the observer is to be computed.  ET */
121 /*                 refers to time at the observer's location. */
122 
123 /*     REF         is the inertial reference frame with respect to which */
124 /*                 the observer's state SOBS is expressed. REF must be */
125 /*                 recognized by the SPICE Toolkit.  The acceptable */
126 /*                 frames are listed in the Frames Required Reading, as */
127 /*                 well as in the SPICELIB routine CHGIRF. */
128 
129 /*                 Case and blanks are not significant in the string REF. */
130 
131 /*     SOBS        is the geometric (uncorrected) state of the observer */
132 /*                 relative to the solar system barycenter at epoch ET. */
133 /*                 SOBS is a 6-vector:  the first three components of */
134 /*                 SOBS represent a Cartesian position vector; the last */
135 /*                 three components represent the corresponding velocity */
136 /*                 vector.  SOBS is expressed relative to the inertial */
137 /*                 reference frame designated by REF. */
138 
139 /*                 Units are always km and km/sec. */
140 
141 
142 /*     ABCORR      indicates the aberration corrections to be applied to */
143 /*                 the position of the target body to account for */
144 /*                 one-way light time and stellar aberration.  See the */
145 /*                 discussion in the Particulars section for */
146 /*                 recommendations on how to choose aberration */
147 /*                 corrections. */
148 
149 /*                 ABCORR may be any of the following: */
150 
151 /*                    'NONE'     Apply no correction. Return the */
152 /*                               geometric position of the target body */
153 /*                               relative to the observer. */
154 
155 /*                 The following values of ABCORR apply to the */
156 /*                 "reception" case in which photons depart from the */
157 /*                 target's location at the light-time corrected epoch */
158 /*                 ET-LT and *arrive* at the observer's location at ET: */
159 
160 /*                    'LT'       Correct for one-way light time (also */
161 /*                               called "planetary aberration") using a */
162 /*                               Newtonian formulation. This correction */
163 /*                               yields the position of the target at the */
164 /*                               moment it emitted photons arriving at */
165 /*                               the observer at ET. */
166 
167 /*                               The light time correction involves */
168 /*                               iterative solution of the light time */
169 /*                               equation (see Particulars for details). */
170 /*                               The solution invoked by the 'LT' option */
171 /*                               uses one iteration. */
172 
173 /*                    'LT+S'     Correct for one-way light time and */
174 /*                               stellar aberration using a Newtonian */
175 /*                               formulation. This option modifies the */
176 /*                               position obtained with the 'LT' option */
177 /*                               to account for the observer's velocity */
178 /*                               relative to the solar system */
179 /*                               barycenter. The result is the apparent */
180 /*                               position of the target---the position */
181 /*                               of the target as seen by the observer. */
182 
183 /*                    'CN'       Converged Newtonian light time */
184 /*                               correction. In solving the light time */
185 /*                               equation, the 'CN' correction iterates */
186 /*                               until the solution converges (three */
187 /*                               iterations on all supported platforms). */
188 /*                               Whether the 'CN+S' solution is */
189 /*                               substantially more accurate than the */
190 /*                               'LT' solution depends on the geometry */
191 /*                               of the participating objects and on the */
192 /*                               accuracy of the input data. In all */
193 /*                               cases this routine will execute more */
194 /*                               slowly when a converged solution is */
195 /*                               computed. See the Particulars section */
196 /*                               of SPKEZR for a discussion of precision */
197 /*                               of light time corrections. */
198 
199 /*                    'CN+S'     Converged Newtonian light time */
200 /*                               correction and stellar aberration */
201 /*                               correction. */
202 
203 
204 /*                 The following values of ABCORR apply to the */
205 /*                 "transmission" case in which photons *depart* from */
206 /*                 the observer's location at ET and arrive at the */
207 /*                 target's location at the light-time corrected epoch */
208 /*                 ET+LT: */
209 
210 /*                    'XLT'      "Transmission" case:  correct for */
211 /*                               one-way light time using a Newtonian */
212 /*                               formulation. This correction yields the */
213 /*                               position of the target at the moment it */
214 /*                               receives photons emitted from the */
215 /*                               observer's location at ET. */
216 
217 /*                    'XLT+S'    "Transmission" case:  correct for */
218 /*                               one-way light time and stellar */
219 /*                               aberration using a Newtonian */
220 /*                               formulation  This option modifies the */
221 /*                               position obtained with the 'XLT' option */
222 /*                               to account for the observer's velocity */
223 /*                               relative to the solar system */
224 /*                               barycenter. The target position */
225 /*                               indicates the direction that photons */
226 /*                               emitted from the observer's location */
227 /*                               must be "aimed" to hit the target. */
228 
229 /*                    'XCN'      "Transmission" case:  converged */
230 /*                               Newtonian light time correction. */
231 
232 /*                    'XCN+S'    "Transmission" case:  converged */
233 /*                               Newtonian light time correction and */
234 /*                               stellar aberration correction. */
235 
236 /*                 Neither special nor general relativistic effects are */
237 /*                 accounted for in the aberration corrections applied */
238 /*                 by this routine. */
239 
240 /*                 Case and blanks are not significant in the string */
241 /*                 ABCORR. */
242 
243 /* $ Detailed_Output */
244 
245 /*     PTARG       is a Cartesian 3-vector representing the position of */
246 /*                 the target body relative to the specified observer. */
247 /*                 PTARG is corrected for the specified aberrations, and */
248 /*                 is expressed with respect to the specified inertial */
249 /*                 reference frame.  The components of PTARG represent */
250 /*                 the x-, y- and z-components of the target's position. */
251 
252 /*                 The vector PTARG points from the observer's position */
253 /*                 at ET to the aberration-corrected location of the */
254 /*                 target. Note that the sense of the position vector is */
255 /*                 independent of the direction of radiation travel */
256 /*                 implied by the aberration correction. */
257 
258 /*                 Units are always km. */
259 
260 /*     LT          is the one-way light time between the observer and */
261 /*                 target in seconds.  If the target position is */
262 /*                 corrected for aberrations, then LT is the one-way */
263 /*                 light time between the observer and the light time */
264 /*                 corrected target location. */
265 
266 /* $ Parameters */
267 
268 /*     None. */
269 
270 /* $ Exceptions */
271 
272 /*     1) If the value of ABCORR is not recognized, the error */
273 /*        'SPICE(SPKINVALIDOPTION)' is signaled. */
274 
275 /*     2) If the reference frame requested is not a recognized */
276 /*        inertial reference frame the error 'SPICE(BADFRAME)' is */
277 /*        signaled. */
278 
279 /*     3) If the position of the target relative to the solar system */
280 /*        barycenter cannot be computed, the error will be diagnosed */
281 /*        by routines in the call tree of this routine. */
282 
283 /* $ Files */
284 
285 
286 /*     This routine computes positions using SPK files that have been */
287 /*     loaded into the SPICE system, normally via the kernel loading */
288 /*     interface routine FURNSH.  Application programs typically load */
289 /*     kernels once before this routine is called, for example during */
290 /*     program initialization; kernels need not be loaded repeatedly. */
291 /*     See the routine FURNSH and the SPK and KERNEL Required Reading */
292 /*     for further information on loading (and unloading) kernels. */
293 
294 /*     If any of the ephemeris data used to compute PTARG are expressed */
295 /*     relative to a non-inertial frame in the SPK files providing those */
296 /*     data, additional kernels may be needed to enable the reference */
297 /*     frame transformations required to compute PTARG.  Normally */
298 /*     these additional kernels are PCK files or frame kernels.  Any */
299 /*     such kernels must already be loaded at the time this routine is */
300 /*     called. */
301 
302 /* $ Particulars */
303 
304 /*     In space science or engineering applications one frequently */
305 /*     wishes to know where to point a remote sensing instrument, such */
306 /*     as an optical camera or radio antenna, in order to observe or */
307 /*     otherwise receive radiation from a target.  This pointing problem */
308 /*     is complicated by the finite speed of light:  one needs to point */
309 /*     to where the target appears to be as opposed to where it actually */
310 /*     is at the epoch of observation.  We use the adjectives */
311 /*     "geometric," "uncorrected," or "true" to refer to an actual */
312 /*     position or state of a target at a specified epoch.  When a */
313 /*     geometric position or state vector is modified to reflect how it */
314 /*     appears to an observer, we describe that vector by any of the */
315 /*     terms "apparent," "corrected," "aberration corrected," or "light */
316 /*     time and stellar aberration corrected." */
317 
318 /*     The SPICE Toolkit can correct for two phenomena affecting the */
319 /*     apparent location of an object:  one-way light time (also called */
320 /*     "planetary aberration") and stellar aberration.  Correcting for */
321 /*     one-way light time is done by computing, given an observer and */
322 /*     observation epoch, where a target was when the observed photons */
323 /*     departed the target's location.  The vector from the observer to */
324 /*     this computed target location is called a "light time corrected" */
325 /*     vector.  The light time correction depends on the motion of the */
326 /*     target, but it is independent of the velocity of the observer */
327 /*     relative to the solar system barycenter. Relativistic effects */
328 /*     such as light bending and gravitational delay are not accounted */
329 /*     for in the light time correction performed by this routine. */
330 
331 /*     The velocity of the observer also affects the apparent location */
332 /*     of a target:  photons arriving at the observer are subject to a */
333 /*     "raindrop effect" whereby their velocity relative to the observer */
334 /*     is, using a Newtonian approximation, the photons' velocity */
335 /*     relative to the solar system barycenter minus the velocity of the */
336 /*     observer relative to the solar system barycenter.  This effect is */
337 /*     called "stellar aberration."  Stellar aberration is independent */
338 /*     of the motion of the target.  The stellar aberration formula used */
339 /*     by this routine is non- relativistic. */
340 
341 /*     Stellar aberration corrections are applied after light time */
342 /*     corrections:  the light time corrected target position vector is */
343 /*     used as an input to the stellar aberration correction. */
344 
345 /*     When light time and stellar aberration corrections are both */
346 /*     applied to a geometric position vector, the resulting position */
347 /*     vector indicates where the target "appears to be" from the */
348 /*     observer's location. */
349 
350 /*     As opposed to computing the apparent position of a target, one */
351 /*     may wish to compute the pointing direction required for */
352 /*     transmission of photons to the target.  This requires correction */
353 /*     of the geometric target position for the effects of light time and */
354 /*     stellar aberration, but in this case the corrections are computed */
355 /*     for radiation traveling from the observer to the target. */
356 
357 /*     The "transmission" light time correction yields the target's */
358 /*     location as it will be when photons emitted from the observer's */
359 /*     location at ET arrive at the target.  The transmission stellar */
360 /*     aberration correction is the inverse of the traditional stellar */
361 /*     aberration correction:  it indicates the direction in which */
362 /*     radiation should be emitted so that, using a Newtonian */
363 /*     approximation, the sum of the velocity of the radiation relative */
364 /*     to the observer and of the observer's velocity, relative to the */
365 /*     solar system barycenter, yields a velocity vector that points in */
366 /*     the direction of the light time corrected position of the target. */
367 
368 /*     The traditional aberration corrections applicable to observation */
369 /*     and those applicable to transmission are related in a simple way: */
370 /*     one may picture the geometry of the "transmission" case by */
371 /*     imagining the "observation" case running in reverse time order, */
372 /*     and vice versa. */
373 
374 /*     One may reasonably object to using the term "observer" in the */
375 /*     transmission case, in which radiation is emitted from the */
376 /*     observer's location.  The terminology was retained for */
377 /*     consistency with earlier documentation. */
378 
379 /*     Below, we indicate the aberration corrections to use for some */
380 /*     common applications: */
381 
382 /*        1) Find the apparent direction of a target for a remote-sensing */
383 /*           observation. */
384 
385 /*              Use 'LT+S' or 'CN+S: apply both light time and stellar */
386 /*              aberration corrections. */
387 
388 /*           Note that using light time corrections alone ('LT' or 'CN') */
389 /*           is generally not a good way to obtain an approximation to */
390 /*           an apparent target vector: since light time and stellar */
391 /*           aberration corrections often partially cancel each other, */
392 /*           it may be more accurate to use no correction at all than to */
393 /*           use light time alone. */
394 
395 
396 /*        2) Find the corrected pointing direction to radiate a signal */
397 /*           to a target. This computation is often applicable for */
398 /*           implementing communications sessions. */
399 
400 /*              Use 'XLT+S' or 'XCN+S: apply both light time and stellar */
401 /*              aberration corrections for transmission. */
402 
403 
404 /*        3) Compute the apparent position of a target body relative */
405 /*           to a star or other distant object. */
406 
407 /*              Use 'LT', 'CN', 'LT+S', or 'CN+S' as needed to match the */
408 /*              correction applied to the position of the distant */
409 /*              object. For example, if a star position is obtained from */
410 /*              a catalog, the position vector may not be corrected for */
411 /*              stellar aberration. In this case, to find the angular */
412 /*              separation of the star and the limb of a planet, the */
413 /*              vector from the observer to the planet should be */
414 /*              corrected for light time but not stellar aberration. */
415 
416 
417 /*        4) Obtain an uncorrected state vector derived directly from */
418 /*           data in an SPK file. */
419 
420 /*              Use 'NONE'. */
421 
422 
423 /*        5) Use a geometric position vector as a low-accuracy estimate */
424 /*           of the apparent position for an application where execution */
425 /*           speed is critical: */
426 
427 /*              Use 'NONE'. */
428 
429 
430 /*        6) While this routine cannot perform the relativistic */
431 /*           aberration corrections required to compute positions */
432 /*           with the highest possible accuracy, it can supply the */
433 /*           geometric positions required as inputs to these */
434 /*           computations: */
435 
436 /*              Use 'NONE', then apply high-accuracy aberration */
437 /*              corrections (not available in the SPICE Toolkit). */
438 
439 
440 /*     Below, we discuss in more detail how the aberration corrections */
441 /*     applied by this routine are computed. */
442 
443 
444 /*     Geometric case */
445 /*     ============== */
446 
447 /*        SPKAPO begins by computing the geometric position T(ET) of the */
448 /*        target body relative to the solar system barycenter (SSB). */
449 /*        Subtracting the geometric position of the observer O(ET) gives */
450 /*        the geometric position of the target body relative to the */
451 /*        observer. The one-way light time, LT, is given by */
452 
453 /*                  | T(ET) - O(ET) | */
454 /*           LT = ------------------- */
455 /*                          c */
456 
457 /*        The geometric relationship between the observer, target, and */
458 /*        solar system barycenter is as shown: */
459 
460 
461 /*           SSB ---> O(ET) */
462 /*            |      / */
463 /*            |     / */
464 /*            |    / */
465 /*            |   /  T(ET) - O(ET) */
466 /*            V  V */
467 /*           T(ET) */
468 
469 
470 /*        The returned position vector is */
471 
472 /*           T(ET) - O(ET) */
473 
474 
475 /*     Reception case */
476 /*     ============== */
477 
478 /*        When any of the options 'LT', 'CN', 'LT+S', 'CN+S' are */
479 /*        selected, SPKAPO computes the position of the target body at */
480 /*        epoch ET-LT, where LT is the one-way light time.  Let T(t) */
481 /*        and O(t) represent the positions of the target and observer */
482 /*        relative to the solar system barycenter at time t; then LT */
483 /*        is the solution of the */
484 /*        light-time equation */
485 
486 /*                  | T(ET-LT) - O(ET) | */
487 /*           LT = ------------------------                            (1) */
488 /*                           c */
489 
490 /*        The ratio */
491 
492 /*            | T(ET) - O(ET) | */
493 /*          ---------------------                                     (2) */
494 /*                    c */
495 
496 /*        is used as a first approximation to LT; inserting (2) into the */
497 /*        RHS of the light-time equation (1) yields the "one-iteration" */
498 /*        estimate of the one-way light time. Repeating the process */
499 /*        until the estimates of LT converge yields the "converged */
500 /*        Newtonian" light time estimate. */
501 
502 /*        Subtracting the geometric position of the observer O(ET) gives */
503 /*        the position of the target body relative to the observer: */
504 /*        T(ET-LT) - O(ET). */
505 
506 /*           SSB ---> O(ET) */
507 /*            | \     | */
508 /*            |  \    | */
509 /*            |   \   | T(ET-LT) - O(ET) */
510 /*            |    \  | */
511 /*            V     V V */
512 /*           T(ET)  T(ET-LT) */
513 
514 
515 /*        The light-time corrected position is the vector */
516 
517 /*           T(ET-LT) - O(ET) */
518 
519 /*        If correction for stellar aberration is requested, the target */
520 /*        position is rotated toward the solar system barycenter-relative */
521 /*        velocity vector of the observer. The magnitude of the rotation */
522 /*        depends on the magnitude of the observer's velocity relative */
523 /*        to the solar system barycenter and the angle between */
524 /*        this velocity and the observer-target vector.  The rotation */
525 /*        is computed as follows: */
526 
527 /*           Let r be the light time corrected vector from the observer */
528 /*           to the object, and v be the velocity of the observer with */
529 /*           respect to the solar system barycenter. Let w be the angle */
530 /*           between them. The aberration angle phi is given by */
531 
532 /*              sin(phi) = v sin(w) / c */
533 
534 /*           Let h be the vector given by the cross product */
535 
536 /*              h = r X v */
537 
538 /*           Rotate r by phi radians about h to obtain the apparent */
539 /*           position of the object. */
540 
541 
542 
543 /*     Transmission case */
544 /*     ================== */
545 
546 /*        When any of the options 'XLT', 'XCN', 'XLT+S', 'XCN+S' are */
547 /*        selected, SPKAPO computes the position of the target body T at */
548 /*        epoch ET+LT, where LT is the one-way light time.  LT is the */
549 /*        solution of the light-time equation */
550 
551 /*                  | T(ET+LT) - O(ET) | */
552 /*           LT = ------------------------                            (3) */
553 /*                            c */
554 
555 /*        Subtracting the geometric position of the observer, O(ET), */
556 /*        gives the position of the target body relative to the */
557 /*        observer: T(ET-LT) - O(ET). */
558 
559 /*                   SSB --> O(ET) */
560 /*                  / |    * */
561 /*                 /  |  *  T(ET+LT) - O(ET) */
562 /*                /   |* */
563 /*               /   *| */
564 /*              V  V  V */
565 /*          T(ET+LT)  T(ET) */
566 
567 
568 /*        The light-time corrected position is */
569 
570 /*           T(ET+LT) - O(ET) */
571 
572 /*        If correction for stellar aberration is requested, the target */
573 /*        position is rotated away from the solar system barycenter- */
574 /*        relative velocity vector of the observer.  The magnitude of the */
575 /*        rotation depends on the magnitude of the velocity and the */
576 /*        angle between the velocity and the observer-target vector. */
577 /*        The rotation is computed as in the reception case, but the */
578 /*        sign of the rotation angle is negated. */
579 
580 /*     Neither special nor general relativistic effects are accounted */
581 /*     for in the aberration corrections performed by this routine. */
582 
583 /* $ Examples */
584 
585 /*     In the following code fragment, SPKSSB and SPKAPO are used */
586 /*     to display the position of Io (body 501) as seen from the */
587 /*     Voyager 2 spacecraft (Body -32) at a series of epochs. */
588 
589 /*     Normally, one would call the high-level reader SPKPOS to obtain */
590 /*     position vectors.  The example below illustrates the interface */
591 /*     of this routine, but is not intended as a recommendation on */
592 /*     how to use the SPICE SPK subsystem. */
593 
594 /*     The use of integer ID codes is necessitated by the low-level */
595 /*     interface of this routine. */
596 
597 /*        IO    = 501 */
598 /*        VGR2  = -32 */
599 
600 /*        DO WHILE ( EPOCH .LE. END ) */
601 
602 /*           CALL SPKSSB (  VGR2,  EPOCH, 'J2000', STVGR2  ) */
603 /*           CALL SPKAPO (  IO,    EPOCH, 'J2000', STVGR2, */
604 /*       .                 'LT+S', STIO,   LT              ) */
605 
606 /*           CALL RECRAD (  STIO,  RANGE,  RA,     DEC     ) */
607 /*           WRITE (*,*)  RA * DPR(),  DEC * DPR() */
608 
609 /*           EPOCH = EPOCH + DELTA */
610 
611 /*        END DO */
612 
613 /* $ Restrictions */
614 
615 /*     1) The ephemeris files to be used by SPKAPO must be loaded */
616 /*        (normally by the SPICELIB kernel loader FURNSH) before */
617 /*        this routine is called. */
618 
619 /*     2) Unlike most other SPK position computation routines, this */
620 /*        routine requires that the input state be relative to an */
621 /*        inertial reference frame.  Non-inertial frames are not */
622 /*        supported by this routine. */
623 
624 /*     3) In a future version of this routine, the implementation */
625 /*        of the aberration corrections may be enhanced to improve */
626 /*        accuracy. */
627 
628 /* $ Literature_References */
629 
630 /*     SPK Required Reading. */
631 
632 /* $ Author_and_Institution */
633 
634 /*     N.J. Bachman    (JPL) */
635 /*     H.A. Neilan     (JPL) */
636 /*     B.V. Semenov    (JPL) */
637 /*     W.L. Taber      (JPL) */
638 /*     I.M. Underwood  (JPL) */
639 
640 /* $ Version */
641 
642 /* -    SPICELIB Version 2.3.0, 03-JUL-2014 (NJB) (BVS) */
643 
644 /*        Discussion of light time corrections was updated. Assertions */
645 /*        that converged light time corrections are unlikely to be */
646 /*        useful were removed. */
647 
648 /*     Last update was 21-SEP-2013 (BVS) */
649 
650 /*        Updated to call LJUCRS instead of CMPRSS/UCASE. */
651 
652 /* -    SPICELIB Version 2.2.0, 17-MAY-2010 (NJB) */
653 
654 /*        Bug fix: routine now returns immediately after */
655 /*        state lookup failure. */
656 
657 /* -    SPICELIB Version 2.1.0, 31-AUG-2005 (NJB) */
658 
659 /*        Updated to remove non-standard use of duplicate arguments */
660 /*        in VSUB call. */
661 
662 /* -    SPICELIB Version 2.0.1, 20-OCT-2003 (EDW) */
663 
664 /*        Added mention that LT returns in seconds. */
665 /*        Corrected spelling errors. */
666 
667 /* -    SPICELIB Version 2.0.0, 18-DEC-2001 (NJB) */
668 
669 /*        Updated to handle aberration corrections for transmission */
670 /*        of radiation.  Formerly, only the reception case was */
671 /*        supported.  The header was revised and expanded to explain */
672 /*        the functionality of this routine in more detail. */
673 
674 /* -    SPICELIB Version 1.0.0, 03-MAR-1999 (WLT) */
675 
676 /* -& */
677 /* $ Index_Entries */
678 
679 /*     apparent position from spk file */
680 /*     get apparent position */
681 
682 /* -& */
683 /* $ Revisions */
684 
685 /* -    SPICELIB Version 1.1.0, 31-AUG-2005 (NJB) */
686 
687 /*        Updated to remove non-standard use of duplicate arguments */
688 /*        in VSUB call. */
689 
690 /* -& */
691 
692 /*     SPICELIB functions */
693 
694 
695 /*     Local parameters */
696 
697 
698 /*     Indices of flags in the FLAGS array: */
699 
700 
701 /*     NAIF ID code for the solar system barycenter: */
702 
703 
704 /*     Local variables */
705 
706 
707 /*     Saved variables */
708 
709 
710 /*     Initial values */
711 
712 
713 /*     Standard SPICE error handling. */
714 
715     if (return_()) {
716 	return 0;
717     } else {
718 	chkin_("ZZSPKPA0", (ftnlen)8);
719     }
720     if (first || s_cmp(abcorr, prvcor, abcorr_len, (ftnlen)5) != 0) {
721 
722 /*        The aberration correction flag differs from the value it */
723 /*        had on the previous call, if any.  Analyze the new flag. */
724 
725 /*        Remove leading and embedded white space from the aberration */
726 /*        correction flag, then convert to upper case. */
727 
728 	ljucrs_(&c__0, abcorr, corr, abcorr_len, (ftnlen)5);
729 
730 /*        Locate the flag in our list of flags. */
731 
732 	i__ = isrchc_(corr, &c__9, flags, (ftnlen)5, (ftnlen)5);
733 	if (i__ == 0) {
734 	    setmsg_("Requested aberration correction was #.", (ftnlen)38);
735 	    errch_("#", abcorr, (ftnlen)1, abcorr_len);
736 	    sigerr_("SPICE(SPKINVALIDOPTION)", (ftnlen)23);
737 	    chkout_("ZZSPKPA0", (ftnlen)8);
738 	    return 0;
739 	}
740 
741 /*        The aberration correction flag is recognized; save it. */
742 
743 	s_copy(prvcor, abcorr, (ftnlen)5, abcorr_len);
744 
745 /*        Set logical flags indicating the attributes of the requested */
746 /*        correction. */
747 
748 	xmit = i__ > 5;
749 	uselt = i__ == 2 || i__ == 3 || i__ == 6 || i__ == 7;
750 	usestl = i__ > 1 && odd_(&i__);
751 	usecn = i__ == 4 || i__ == 5 || i__ == 8 || i__ == 9;
752 	first = FALSE_;
753     }
754 
755 /*     See if the reference frame is a recognized inertial frame. */
756 
757     irfnum_(ref, &refid, ref_len);
758     if (refid == 0) {
759 	setmsg_("The requested frame '#' is not a recognized inertial frame. "
760 		, (ftnlen)60);
761 	errch_("#", ref, (ftnlen)1, ref_len);
762 	sigerr_("SPICE(BADFRAME)", (ftnlen)15);
763 	chkout_("ZZSPKPA0", (ftnlen)8);
764 	return 0;
765     }
766 
767 /*     Determine the sign of the light time offset. */
768 
769     if (xmit) {
770 	ltsign = 1;
771     } else {
772 	ltsign = -1;
773     }
774 
775 /*     Find the geometric position of the target body with respect to the */
776 /*     solar system barycenter. Subtract the position of the observer */
777 /*     to get the relative position. Use this to compute the one-way */
778 /*     light time. */
779 
780     zzspkgp0_(targ, et, ref, &c__0, ptarg, lt, ref_len);
781     if (failed_()) {
782 	chkout_("ZZSPKPA0", (ftnlen)8);
783 	return 0;
784     }
785     vsub_(ptarg, sobs, tpos);
786     vequ_(tpos, ptarg);
787     *lt = vnorm_(ptarg) / clight_();
788 
789 /*     To correct for light time, find the position of the target body */
790 /*     at the current epoch minus the one-way light time. Note that */
791 /*     the observer remains where he is. */
792 
793     if (uselt) {
794 	maxitr = 1;
795     } else if (usecn) {
796 	maxitr = 3;
797     } else {
798 	maxitr = 0;
799     }
800     i__1 = maxitr;
801     for (i__ = 1; i__ <= i__1; ++i__) {
802 	d__1 = *et + ltsign * *lt;
803 	zzspkgp0_(targ, &d__1, ref, &c__0, ptarg, lt, ref_len);
804 	if (failed_()) {
805 	    chkout_("ZZSPKPA0", (ftnlen)8);
806 	    return 0;
807 	}
808 	vsub_(ptarg, sobs, tpos);
809 	vequ_(tpos, ptarg);
810 	*lt = vnorm_(ptarg) / clight_();
811     }
812 
813 /*     At this point, PTARG contains the geometric or light-time */
814 /*     corrected position of the target relative to the observer, */
815 /*     depending on the specified correction. */
816 
817 /*     If stellar aberration correction is requested, perform it now. */
818 
819     if (usestl) {
820 	if (xmit) {
821 
822 /*           This is the transmission case. */
823 
824 /*           Compute the position vector obtained by applying */
825 /*           "reception" stellar aberration to PTARG. */
826 
827 	    stlabx_(ptarg, &sobs[3], tpos);
828 	    vequ_(tpos, ptarg);
829 	} else {
830 
831 /*           This is the reception case. */
832 
833 /*           Compute the position vector obtained by applying */
834 /*           "reception" stellar aberration to PTARG. */
835 
836 	    stelab_(ptarg, &sobs[3], tpos);
837 	    vequ_(tpos, ptarg);
838 	}
839     }
840     chkout_("ZZSPKPA0", (ftnlen)8);
841     return 0;
842 } /* zzspkpa0_ */
843 
844