/* $Procedure DVSEP ( Derivative of separation angle ) */ doublereal dvsep_(doublereal *s1, doublereal *s2) { /* System generated locals */ doublereal ret_val; /* Local variables */ logical safe; extern doublereal vdot_(doublereal *, doublereal *); doublereal numr; extern /* Subroutine */ int chkin_(char *, ftnlen); doublereal denom; extern /* Subroutine */ int dvhat_(doublereal *, doublereal *); extern doublereal dpmax_(void); extern /* Subroutine */ int vcrss_(doublereal *, doublereal *, doublereal *); extern doublereal vnorm_(doublereal *); extern logical vzero_(doublereal *); doublereal u1[6], u2[6]; extern /* Subroutine */ int sigerr_(char *, ftnlen), chkout_(char *, ftnlen), setmsg_(char *, ftnlen); doublereal pcross[3]; extern logical return_(void); /* $ Abstract */ /* Calculate the time derivative of the separation angle between */ /* two input states, S1 and S2. */ /* $ Disclaimer */ /* THIS SOFTWARE AND ANY RELATED MATERIALS WERE CREATED BY THE */ /* CALIFORNIA INSTITUTE OF TECHNOLOGY (CALTECH) UNDER A U.S. */ /* GOVERNMENT CONTRACT WITH THE NATIONAL AERONAUTICS AND SPACE */ /* ADMINISTRATION (NASA). THE SOFTWARE IS TECHNOLOGY AND SOFTWARE */ /* PUBLICLY AVAILABLE UNDER U.S. EXPORT LAWS AND IS PROVIDED "AS-IS" */ /* TO THE RECIPIENT WITHOUT WARRANTY OF ANY KIND, INCLUDING ANY */ /* WARRANTIES OF PERFORMANCE OR MERCHANTABILITY OR FITNESS FOR A */ /* PARTICULAR USE OR PURPOSE (AS SET FORTH IN UNITED STATES UCC */ /* SECTIONS 2312-2313) OR FOR ANY PURPOSE WHATSOEVER, FOR THE */ /* SOFTWARE AND RELATED MATERIALS, HOWEVER USED. */ /* IN NO EVENT SHALL CALTECH, ITS JET PROPULSION LABORATORY, OR NASA */ /* BE LIABLE FOR ANY DAMAGES AND/OR COSTS, INCLUDING, BUT NOT */ /* LIMITED TO, INCIDENTAL OR CONSEQUENTIAL DAMAGES OF ANY KIND, */ /* INCLUDING ECONOMIC DAMAGE OR INJURY TO PROPERTY AND LOST PROFITS, */ /* REGARDLESS OF WHETHER CALTECH, JPL, OR NASA BE ADVISED, HAVE */ /* REASON TO KNOW, OR, IN FACT, SHALL KNOW OF THE POSSIBILITY. */ /* RECIPIENT BEARS ALL RISK RELATING TO QUALITY AND PERFORMANCE OF */ /* THE SOFTWARE AND ANY RELATED MATERIALS, AND AGREES TO INDEMNIFY */ /* CALTECH AND NASA FOR ALL THIRD-PARTY CLAIMS RESULTING FROM THE */ /* ACTIONS OF RECIPIENT IN THE USE OF THE SOFTWARE. */ /* $ Required_Reading */ /* None. */ /* $ Keywords */ /* GEOMETRY */ /* DERIVATIVES */ /* $ Declarations */ /* $ Brief_I/O */ /* VARIABLE I/O DESCRIPTION */ /* -------- --- -------------------------------------------------- */ /* S1 I State vector of the first body. */ /* S2 I State vector of the second body. */ /* $ Detailed_Input */ /* S1 the state vector of the first target body as seen from */ /* the observer. */ /* S2 the state vector of the second target body as seen from */ /* the observer. */ /* An implicit assumption exists that both states lie in the same */ /* reference frame with the same observer for the same epoch. If this */ /* is not the case, the numerical result has no meaning. */ /* $ Detailed_Output */ /* The function returns the double precision value of the time */ /* derivative of the angular separation between S1 and S2. */ /* $ Parameters */ /* None. */ /* $ Exceptions */ /* 1) The error SPICE(NUMERICOVERFLOW) signals if the inputs S1, S2 */ /* define states with an angular separation rate ~ DPMAX(). */ /* 2) If called in RETURN mode, the return has value 0. */ /* 3) Linear dependent position components of S1 and S1 constitutes */ /* a non-error exception. The function returns 0 for this case. */ /* $ Files */ /* None. */ /* $ Particulars */ /* In this discussion, the notation */ /* < V1, V2 > */ /* indicates the dot product of vectors V1 and V2. The notation */ /* V1 x V2 */ /* indicates the cross product of vectors V1 and V2. */ /* To start out, note that we need consider only unit vectors, */ /* since the angular separation of any two non-zero vectors */ /* equals the angular separation of the corresponding unit vectors. */ /* Call these vectors U1 and U2; let their velocities be V1 and V2. */ /* For unit vectors having angular separation */ /* THETA */ /* the identity */ /* || U1 x U1 || = ||U1|| * ||U2|| * sin(THETA) (1) */ /* reduces to */ /* || U1 x U2 || = sin(THETA) (2) */ /* and the identity */ /* | < U1, U2 > | = || U1 || * || U2 || * cos(THETA) (3) */ /* reduces to */ /* | < U1, U2 > | = cos(THETA) (4) */ /* Since THETA is an angular separation, THETA is in the range */ /* 0 : Pi */ /* Then letting s be +1 if cos(THETA) > 0 and -1 if cos(THETA) < 0, */ /* we have for any value of THETA other than 0 or Pi */ /* 2 1/2 */ /* cos(THETA) = s * ( 1 - sin (THETA) ) (5) */ /* or */ /* 2 1/2 */ /* < U1, U2 > = s * ( 1 - sin (THETA) ) (6) */ /* At this point, for any value of THETA other than 0 or Pi, */ /* we can differentiate both sides with respect to time (T) */ /* to obtain */ /* 2 -1/2 */ /* < U1, V2 > + < V1, U2 > = s * (1/2)(1 - sin (THETA)) */ /* * (-2) sin(THETA)*cos(THETA) */ /* * d(THETA)/dT (7a) */ /* Using equation (5), and noting that s = 1/s, we can cancel */ /* the cosine terms on the right hand side */ /* -1 */ /* < U1, V2 > + < V1, U2 > = (1/2)(cos(THETA)) */ /* * (-2) sin(THETA)*cos(THETA) */ /* * d(THETA)/dT (7b) */ /* With (7b) reducing to */ /* < U1, V2 > + < V1, U2 > = - sin(THETA) * d(THETA)/dT (8) */ /* Using equation (2) and switching sides, we obtain */ /* || U1 x U2 || * d(THETA)/dT = - < U1, V2 > - < V1, U2 > (9) */ /* or, provided U1 and U2 are linearly independent, */ /* d(THETA)/dT = ( - < U1, V2 > - < V1, U2 > ) / ||U1 x U2|| (10) */ /* Note for times when U1 and U2 have angular separation 0 or Pi */ /* radians, the derivative of angular separation with respect to */ /* time doesn't exist. (Consider the graph of angular separation */ /* with respect to time; typically the graph is roughly v-shaped at */ /* the singular points.) */ /* $ Examples */ /* PROGRAM DVSEP_T */ /* IMPLICIT NONE */ /* DOUBLE PRECISION ET */ /* DOUBLE PRECISION LT */ /* DOUBLE PRECISION DSEPT */ /* DOUBLE PRECISION STATEE (6) */ /* DOUBLE PRECISION STATEM (6) */ /* INTEGER STRLEN */ /* PARAMETER ( STRLEN = 64 ) */ /* CHARACTER*(STRLEN) BEGSTR */ /* DOUBLE PRECISION DVSEP */ /* C */ /* C Load kernels. */ /* C */ /* CALL FURNSH ('standard.tm') */ /* C */ /* C An arbitrary time. */ /* C */ /* BEGSTR = 'JAN 1 2009' */ /* CALL STR2ET( BEGSTR, ET ) */ /* C */ /* C Calculate the state vectors sun to Moon, sun to earth at ET. */ /* C */ /* C */ /* CALL SPKEZR ( 'EARTH', ET, 'J2000', 'NONE', 'SUN', */ /* . STATEE, LT) */ /* CALL SPKEZR ( 'MOON', ET, 'J2000', 'NONE', 'SUN', */ /* . STATEM, LT) */ /* C */ /* C Calculate the time derivative of the angular separation of */ /* C the earth and Moon as seen from the sun at ET. */ /* C */ /* DSEPT = DVSEP( STATEE, STATEM ) */ /* WRITE(*,*) 'Time derivative of angular separation: ', DSEPT */ /* END */ /* The program compiled on OS X with g77 outputs (radians/sec): */ /* Time derivative of angular separation: 3.81211936E-09 */ /* $ Restrictions */ /* None. */ /* $ Literature_References */ /* None. */ /* $ Author_and_Institution */ /* E.D. Wright (JPL) */ /* N.J. Bachman (JPL) */ /* $ Version */ /* - SPICELIB Version 1.0.1, 15-MAR-2010 (EDW) */ /* Trivial header format clean-up. */ /* - SPICELIB Version 1.0.1, 31-MAR-2009 (EDW) */ /* -& */ /* $ Index_Entries */ /* time derivative of angular separation */ /* -& */ /* SPICELIB functions */ /* Local variables */ if (return_()) { ret_val = 0.; return ret_val; } chkin_("DVSEP", (ftnlen)5); /* Compute the unit vectors and corresponding time derivatives */ /* for the input state vectors. */ dvhat_(s1, u1); dvhat_(s2, u2); /* Calculate the cross product vector of U1 and U2. As both vectors */ /* have magnitude one, the magnitude of the cross product equals */ /* sin(THETA), with THETA the angle between S1 and S2. */ vcrss_(u1, u2, pcross); /* Now calculate the time derivate of the angular separation between */ /* S1 and S2. */ /* The routine needs to guard against both division by zero */ /* and numeric overflow. Before carrying out the division */ /* indicated by equation (10), the routine should verify that */ /* || U1 x U2 || > fudge factor * | numerator | / DPMAX() */ /* A fudge factor of 10.D0 should suffice. */ /* Note that the inequality is strict. */ /* Handle the parallel and anti-parallel cases. */ if (vzero_(pcross)) { ret_val = 0.; chkout_("DVSEP", (ftnlen)5); return ret_val; } /* Now check for possible overflow. */ numr = vdot_(u1, &u2[3]) + vdot_(&u1[3], u2); denom = vnorm_(pcross); safe = denom > abs(numr) * 10. / dpmax_(); if (! safe) { ret_val = 0.; setmsg_("Numerical overflow event.", (ftnlen)25); sigerr_("SPICE(NUMERICOVERFLOW)", (ftnlen)22); chkout_("DVSEP", (ftnlen)5); return ret_val; } ret_val = -numr / denom; chkout_("DVSEP", (ftnlen)5); return ret_val; } /* dvsep_ */
/* $Procedure ZZSTELAB ( Private --- stellar aberration correction ) */ /* Subroutine */ int zzstelab_(logical *xmit, doublereal *accobs, doublereal * vobs, doublereal *starg, doublereal *scorr, doublereal *dscorr) { /* System generated locals */ integer i__1; doublereal d__1, d__2; /* Builtin functions */ double sqrt(doublereal); integer s_rnge(char *, integer, char *, integer); /* Local variables */ extern /* Subroutine */ int vadd_(doublereal *, doublereal *, doublereal * ); doublereal dphi, rhat[3]; extern /* Subroutine */ int vhat_(doublereal *, doublereal *); extern doublereal vdot_(doublereal *, doublereal *); extern /* Subroutine */ int vequ_(doublereal *, doublereal *); doublereal term1[3], term2[3], term3[3], c__, lcacc[3]; integer i__; doublereal s; extern /* Subroutine */ int chkin_(char *, ftnlen); doublereal saoff[6] /* was [3][2] */, drhat[3]; extern /* Subroutine */ int dvhat_(doublereal *, doublereal *); doublereal ptarg[3], evobs[3], srhat[6], vphat[3], vtarg[3]; extern /* Subroutine */ int vlcom_(doublereal *, doublereal *, doublereal *, doublereal *, doublereal *), vperp_(doublereal *, doublereal *, doublereal *); extern doublereal vnorm_(doublereal *); extern logical vzero_(doublereal *); extern /* Subroutine */ int vlcom3_(doublereal *, doublereal *, doublereal *, doublereal *, doublereal *, doublereal *, doublereal *), cleard_(integer *, doublereal *); doublereal vp[3]; extern doublereal clight_(void); doublereal dptmag, ptgmag, eptarg[3], dvphat[3], lcvobs[3]; extern /* Subroutine */ int qderiv_(integer *, doublereal *, doublereal *, doublereal *, doublereal *), sigerr_(char *, ftnlen), chkout_( char *, ftnlen), setmsg_(char *, ftnlen); doublereal svphat[6]; extern logical return_(void); extern /* Subroutine */ int vminus_(doublereal *, doublereal *); doublereal sgn, dvp[3], svp[6]; /* $ Abstract */ /* SPICE Private routine intended solely for the support of SPICE */ /* routines. Users should not call this routine directly due */ /* to the volatile nature of this routine. */ /* Return the state (position and velocity) of a target body */ /* relative to an observing body, optionally corrected for light */ /* time (planetary aberration) and stellar aberration. */ /* $ Disclaimer */ /* THIS SOFTWARE AND ANY RELATED MATERIALS WERE CREATED BY THE */ /* CALIFORNIA INSTITUTE OF TECHNOLOGY (CALTECH) UNDER A U.S. */ /* GOVERNMENT CONTRACT WITH THE NATIONAL AERONAUTICS AND SPACE */ /* ADMINISTRATION (NASA). THE SOFTWARE IS TECHNOLOGY AND SOFTWARE */ /* PUBLICLY AVAILABLE UNDER U.S. EXPORT LAWS AND IS PROVIDED "AS-IS" */ /* TO THE RECIPIENT WITHOUT WARRANTY OF ANY KIND, INCLUDING ANY */ /* WARRANTIES OF PERFORMANCE OR MERCHANTABILITY OR FITNESS FOR A */ /* PARTICULAR USE OR PURPOSE (AS SET FORTH IN UNITED STATES UCC */ /* SECTIONS 2312-2313) OR FOR ANY PURPOSE WHATSOEVER, FOR THE */ /* SOFTWARE AND RELATED MATERIALS, HOWEVER USED. */ /* IN NO EVENT SHALL CALTECH, ITS JET PROPULSION LABORATORY, OR NASA */ /* BE LIABLE FOR ANY DAMAGES AND/OR COSTS, INCLUDING, BUT NOT */ /* LIMITED TO, INCIDENTAL OR CONSEQUENTIAL DAMAGES OF ANY KIND, */ /* INCLUDING ECONOMIC DAMAGE OR INJURY TO PROPERTY AND LOST PROFITS, */ /* REGARDLESS OF WHETHER CALTECH, JPL, OR NASA BE ADVISED, HAVE */ /* REASON TO KNOW, OR, IN FACT, SHALL KNOW OF THE POSSIBILITY. */ /* RECIPIENT BEARS ALL RISK RELATING TO QUALITY AND PERFORMANCE OF */ /* THE SOFTWARE AND ANY RELATED MATERIALS, AND AGREES TO INDEMNIFY */ /* CALTECH AND NASA FOR ALL THIRD-PARTY CLAIMS RESULTING FROM THE */ /* ACTIONS OF RECIPIENT IN THE USE OF THE SOFTWARE. */ /* $ Required_Reading */ /* SPK */ /* $ Keywords */ /* EPHEMERIS */ /* $ Declarations */ /* $ Brief_I/O */ /* Variable I/O Description */ /* -------- --- -------------------------------------------------- */ /* XMIT I Reception/transmission flag. */ /* ACCOBS I Observer acceleration relative to SSB. */ /* VOBS I Observer velocity relative to to SSB. */ /* STARG I State of target relative to observer. */ /* SCORR O Stellar aberration correction for position. */ /* DSCORR O Stellar aberration correction for velocity. */ /* $ Detailed_Input */ /* XMIT is a logical flag which is set to .TRUE. for the */ /* "transmission" case in which photons *depart* from */ /* the observer's location at an observation epoch ET */ /* and arrive at the target's location at the light-time */ /* corrected epoch ET+LT, where LT is the one-way light */ /* time between observer and target; XMIT is set to */ /* .FALSE. for "reception" case in which photons depart */ /* from the target's location at the light-time */ /* corrected epoch ET-LT and *arrive* at the observer's */ /* location at ET. */ /* Note that the observation epoch is not used in this */ /* routine. */ /* XMIT must be consistent with any light time */ /* corrections used for the input state STARG: if that */ /* state has been corrected for "reception" light time; */ /* XMIT must be .FALSE.; otherwise XMIT must be .TRUE. */ /* ACCOBS is the geometric acceleration of the observer */ /* relative to the solar system barycenter. Units are */ /* km/sec**2. ACCOBS must be expressed relative to */ /* an inertial reference frame. */ /* VOBS is the geometric velocity of the observer relative to */ /* the solar system barycenter. VOBS must be expressed */ /* relative to the same inertial reference frame as */ /* ACCOBS. Units are km/sec. */ /* STARG is the Cartesian state of the target relative to the */ /* observer. Normally STARG has been corrected for */ /* one-way light time, but this is not required. STARG */ /* must be expressed relative to the same inertial */ /* reference frame as ACCOBS. Components are */ /* (x, y, z, dx, dy, dz). Units are km and km/sec. */ /* $ Detailed_Output */ /* SCORR is the stellar aberration correction for the position */ /* component of STARG. Adding SCORR to this position */ /* vector produces the input observer-target position, */ /* corrected for stellar aberration. */ /* The reference frame of SCORR is the common frame */ /* relative to which the inputs ACCOBS, VOBS, and STARG */ /* are expressed. Units are km. */ /* DSCORR is the stellar aberration correction for the velocity */ /* component of STARG. Adding DSCORR to this velocity */ /* vector produces the input observer-target velocity, */ /* corrected for stellar aberration. */ /* The reference frame of DSCORR is the common frame */ /* relative to which the inputs ACCOBS, VOBS, and STARG */ /* are expressed. Units are km/s. */ /* $ Parameters */ /* None. */ /* $ Exceptions */ /* 1) If attempt to divide by zero occurs, the error */ /* SPICE(DIVIDEBYZERO) will be signaled. This case may occur */ /* due to uninitialized inputs. */ /* 2) Loss of precision will occur for geometric cases in which */ /* VOBS is nearly parallel to the position component of STARG. */ /* $ Files */ /* None. */ /* $ Particulars */ /* This routine computes a Newtonian estimate of the stellar */ /* aberration correction of an input state. Normally the input state */ /* has already been corrected for one-way light time. */ /* Since stellar aberration corrections are typically "small" */ /* relative to the magnitude of the input observer-target position */ /* and velocity, this routine avoids loss of precision by returning */ /* the corrections themselves rather than the corrected state */ /* vector. This allows the caller to manipulate (for example, */ /* interpolate) the corrections with greater accuracy. */ /* $ Examples */ /* See SPICELIB routine SPKACS. */ /* $ Restrictions */ /* None. */ /* $ Literature_References */ /* SPK Required Reading. */ /* $ Author_and_Institution */ /* N.J. Bachman (JPL) */ /* $ Version */ /* - SPICELIB Version 2.0.0, 15-APR-2014 (NJB) */ /* Added RETURN test and discovery check-in. */ /* Check for division by zero was added. This */ /* case might occur due to uninitialized inputs. */ /* - SPICELIB Version 1.0.1, 12-FEB-2009 (NJB) */ /* Minor updates were made to the inline documentation. */ /* - SPICELIB Version 1.0.0, 17-JAN-2008 (NJB) */ /* -& */ /* Note for the maintenance programmer */ /* =================================== */ /* The source code of the test utility T_ZZSTLABN must be */ /* kept in sync with the source code of this routine. That */ /* routine uses a value of SEPLIM that forces the numeric */ /* branch of the velocity computation to be taken in all */ /* cases. See the documentation of that routine for details. */ /* SPICELIB functions */ /* Local parameters */ /* Let PHI be the (non-negative) rotation angle of the stellar */ /* aberration correction; then SEPLIM is a limit on how close PHI */ /* may be to zero radians while stellar aberration velocity is */ /* computed analytically. When sin(PHI) is less than SEPLIM, the */ /* velocity must be computed numerically. */ /* Let TDELTA be the time interval, measured in seconds, */ /* used for numerical differentiation of the stellar */ /* aberration correction, when this is necessary. */ /* Local variables */ /* Use discovery check-in. */ if (return_()) { return 0; } /* In the discussion below, the dot product of vectors X and Y */ /* is denoted by */ /* <X,Y> */ /* The speed of light is denoted by the lower case letter "c." BTW, */ /* variable names used here are case-sensitive: upper case "C" */ /* represents a different quantity which is unrelated to the speed */ /* of light. */ /* Variable names ending in "HAT" denote unit vectors. Variable */ /* names starting with "D" denote derivatives with respect to time. */ /* We'll compute the correction SCORR and its derivative with */ /* respect to time DSCORR for the reception case. In the */ /* transmission case, we perform the same computation with the */ /* negatives of the observer velocity and acceleration. */ /* In the code below, we'll store the position and velocity portions */ /* of the input observer-target state STARG in the variables PTARG */ /* and VTARG, respectively. */ /* Let VP be the component of VOBS orthogonal to PTARG. VP */ /* is defined as */ /* VOBS - < VOBS, RHAT > RHAT (1) */ /* where RHAT is the unit vector */ /* PTARG/||PTARG|| */ /* Then */ /* ||VP||/c (2) */ /* is the magnitude of */ /* s = sin( phi ) (3) */ /* where phi is the stellar aberration correction angle. We'll */ /* need the derivative with respect to time of (2). */ /* Differentiating (1) with respect to time yields the */ /* velocity DVP, where, letting */ /* DRHAT = d(RHAT) / dt */ /* VPHAT = VP / ||VP|| */ /* DVPMAG = d( ||VP|| ) / dt */ /* we have */ /* DVP = d(VP)/dt */ /* = ACCOBS - ( ( <VOBS,DRHAT> + <ACCOBS, RHAT> )*RHAT */ /* + <VOBS,RHAT> * DRHAT ) (4) */ /* and */ /* DVPMAG = < DVP, VPHAT > (5) */ /* Now we can find the derivative with respect to time of */ /* the stellar aberration angle phi: */ /* ds/dt = d(sin(phi))/dt = d(phi)/dt * cos(phi) (6) */ /* Using (2) and (5), we have for positive phi, */ /* ds/dt = (1/c)*DVPMAG = (1/c)*<DVP, VPHAT> (7) */ /* Then for positive phi */ /* d(phi)/dt = (1/cos(phi)) * (1/c) * <DVP, VPHAT> (8) */ /* Equation (8) is well-defined as along as VP is non-zero: */ /* if VP is the zero vector, VPHAT is undefined. We'll treat */ /* the singular and near-singular cases separately. */ /* The aberration correction itself is a rotation by angle phi */ /* from RHAT towards VP, so the corrected vector is */ /* ( sin(phi)*VPHAT + cos(phi)*RHAT ) * ||PTARG|| */ /* and we can express the offset of the corrected vector from */ /* PTARG, which is the output SCORR, as */ /* SCORR = */ /* ( sin(phi)*VPHAT + (cos(phi)-1)*RHAT ) * ||PTARG|| (9) */ /* Let DPTMAG be defined as */ /* DPTMAG = d ( ||PTARG|| ) / dt (10) */ /* Then the derivative with respect to time of SCORR is */ /* DSCORR = */ /* ( sin(phi)*DVPHAT */ /* + cos(phi)*d(phi)/dt * VPHAT */ /* + (cos(phi) - 1) * DRHAT */ /* + ( -sin(phi)*d(phi)/dt ) * RHAT ) * ||PTARG|| */ /* + ( sin(phi)*VPHAT + (cos(phi)-1)*RHAT ) * DPTMAG (11) */ /* Computations begin here: */ /* Split STARG into position and velocity components. Compute */ /* RHAT */ /* DRHAT */ /* VP */ /* DPTMAG */ if (*xmit) { vminus_(vobs, lcvobs); vminus_(accobs, lcacc); } else { vequ_(vobs, lcvobs); vequ_(accobs, lcacc); } vequ_(starg, ptarg); vequ_(&starg[3], vtarg); dvhat_(starg, srhat); vequ_(srhat, rhat); vequ_(&srhat[3], drhat); vperp_(lcvobs, rhat, vp); dptmag = vdot_(vtarg, rhat); /* Compute sin(phi) and cos(phi), which we'll call S and C */ /* respectively. Note that phi is always close to zero for */ /* realistic inputs (for which ||VOBS|| << CLIGHT), so the */ /* cosine term is positive. */ s = vnorm_(vp) / clight_(); /* Computing MAX */ d__1 = 0., d__2 = 1 - s * s; c__ = sqrt((max(d__1,d__2))); if (c__ == 0.) { /* C will be used as a divisor later (in the computation */ /* of DPHI), so we'll put a stop to the problem here. */ chkin_("ZZSTELAB", (ftnlen)8); setmsg_("Cosine of the aberration angle is 0; this cannot occur for " "realistic observer velocities. This case can arise due to un" "initialized inputs. This cosine value is used as a divisor i" "n a later computation, so it must not be equal to zero.", ( ftnlen)234); sigerr_("SPICE(DIVIDEBYZERO)", (ftnlen)19); chkout_("ZZSTELAB", (ftnlen)8); return 0; } /* Compute the unit vector VPHAT and the stellar */ /* aberration correction. We avoid relying on */ /* VHAT's exception handling for the zero vector. */ if (vzero_(vp)) { cleard_(&c__3, vphat); } else { vhat_(vp, vphat); } /* Now we can use equation (9) to obtain the stellar */ /* aberration correction SCORR: */ /* SCORR = */ /* ( sin(phi)*VPHAT + (cos(phi)-1)*RHAT ) * ||PTARG|| */ ptgmag = vnorm_(ptarg); d__1 = ptgmag * s; d__2 = ptgmag * (c__ - 1.); vlcom_(&d__1, vphat, &d__2, rhat, scorr); /* Now we use S as an estimate of PHI to decide if we're */ /* going to differentiate the stellar aberration correction */ /* analytically or numerically. */ /* Note that S is non-negative by construction, so we don't */ /* need to use the absolute value of S here. */ if (s >= 1e-6) { /* This is the analytic case. */ /* Compute DVP---the derivative of VP with respect to time. */ /* Recall equation (4): */ /* DVP = d(VP)/dt */ /* = ACCOBS - ( ( <VOBS,DRHAT> + <ACCOBS, RHAT> )*RHAT */ /* + <VOBS,RHAT> * DRHAT ) */ d__1 = -vdot_(lcvobs, drhat) - vdot_(lcacc, rhat); d__2 = -vdot_(lcvobs, rhat); vlcom3_(&c_b7, lcacc, &d__1, rhat, &d__2, drhat, dvp); vhat_(vp, vphat); /* Now we can compute DVPHAT, the derivative of VPHAT: */ vequ_(vp, svp); vequ_(dvp, &svp[3]); dvhat_(svp, svphat); vequ_(&svphat[3], dvphat); /* Compute the DPHI, the time derivative of PHI, using equation 8: */ /* d(phi)/dt = (1/cos(phi)) * (1/c) * <DVP, VPHAT> */ dphi = 1. / (c__ * clight_()) * vdot_(dvp, vphat); /* At long last we've assembled all of the "ingredients" required */ /* to compute DSCORR: */ /* DSCORR = */ /* ( sin(phi)*DVPHAT */ /* + cos(phi)*d(phi)/dt * VPHAT */ /* + (cos(phi) - 1) * DRHAT */ /* + ( -sin(phi)*d(phi)/dt ) * RHAT ) * ||PTARG|| */ /* + ( sin(phi)*VPHAT + (cos(phi)-1)*RHAT ) * DPTMAG */ d__1 = c__ * dphi; vlcom_(&s, dvphat, &d__1, vphat, term1); d__1 = c__ - 1.; d__2 = -s * dphi; vlcom_(&d__1, drhat, &d__2, rhat, term2); vadd_(term1, term2, term3); d__1 = dptmag * s; d__2 = dptmag * (c__ - 1.); vlcom3_(&ptgmag, term3, &d__1, vphat, &d__2, rhat, dscorr); } else { /* This is the numeric case. We're going to differentiate */ /* the stellar aberration correction offset vector using */ /* a quadratic estimate. */ for (i__ = 1; i__ <= 2; ++i__) { /* Set the sign of the time offset. */ if (i__ == 1) { sgn = -1.; } else { sgn = 1.; } /* Estimate the observer's velocity relative to the */ /* solar system barycenter at the current epoch. We use */ /* the local copies of the input velocity and acceleration */ /* to make a linear estimate. */ d__1 = sgn * 1.; vlcom_(&c_b7, lcvobs, &d__1, lcacc, evobs); /* Estimate the observer-target vector. We use the */ /* observer-target state velocity to make a linear estimate. */ d__1 = sgn * 1.; vlcom_(&c_b7, starg, &d__1, &starg[3], eptarg); /* Let RHAT be the unit observer-target position. */ /* Compute the component of the observer's velocity */ /* that is perpendicular to the target position; call */ /* this vector VP. Also compute the unit vector in */ /* the direction of VP. */ vhat_(eptarg, rhat); vperp_(evobs, rhat, vp); if (vzero_(vp)) { cleard_(&c__3, vphat); } else { vhat_(vp, vphat); } /* Compute the sine and cosine of the correction */ /* angle. */ s = vnorm_(vp) / clight_(); /* Computing MAX */ d__1 = 0., d__2 = 1 - s * s; c__ = sqrt((max(d__1,d__2))); /* Compute the vector offset of the correction. */ ptgmag = vnorm_(eptarg); d__1 = ptgmag * s; d__2 = ptgmag * (c__ - 1.); vlcom_(&d__1, vphat, &d__2, rhat, &saoff[(i__1 = i__ * 3 - 3) < 6 && 0 <= i__1 ? i__1 : s_rnge("saoff", i__1, "zzstelab_", ( ftnlen)597)]); } /* Now compute the derivative. */ qderiv_(&c__3, saoff, &saoff[3], &c_b7, dscorr); } /* At this point the correction offset SCORR and its derivative */ /* with respect to time DSCORR are both set. */ return 0; } /* zzstelab_ */