Index of Functions: A  B  C  D  E  F  G  H  I  J  K  L  M  N  O  P  Q  R  S  T  U  V  W  X 
Index Page
cspice_spkezr

Table of contents
Abstract
I/O
Parameters
Examples
Particulars
Exceptions
Files
Restrictions
Required_Reading
Literature_References
Author_and_Institution
Version
Index_Entries


Abstract


   CSPICE_SPKEZR returns the state (position and velocity) of
   a target body relative to an observing body, optionally
   corrected for light  time (planetary aberration) and stellar
   aberration.

I/O


   Given:

      targ     the name of a target body.

               help, targ
                  STRING = Scalar

               Optionally, you may supply the integer ID code for the object
               as an integer string. For example both 'MOON' and '301' are
               legitimate strings that indicate the moon is the target body.

               The target and observer define a state vector whose
               position component points from the observer to the
               target.

      et       the ephemeris time, expressed as seconds past J2000 TDB, at
               which the state of the target body relative to the observer is
               to be computed, or an N-vector of times.

               help, et
                  DOUBLE = Scalar   or   DOUBLE = Array[N]

               `et' refers to time at the observer's location.

      ref      the name of the reference frame relative to which the output
               state vector should be expressed.

               help, ref
                  STRING = Scalar

               This may be any frame supported by the SPICE system, including
               built-in frames (documented in the Frames Required Reading) and
               frames defined by a loaded frame kernel (FK).

               When `ref' designates a non-inertial frame, the
               orientation of the frame is evaluated at an epoch
               dependent on the selected aberration correction.
               See the description of the output state vector `starg'
               for details.

      abcorr   indicates the aberration corrections to be applied to the state
               of the target body to account for one-way light time and stellar
               aberration.

               help, abcorr
                  STRING = Scalar

               See the discussion in the -Particulars section for
               recommendations on how to choose aberration corrections.

               `abcorr' may be any of the following:

                  'NONE'     Apply no correction. Return the
                             geometric state of the target body
                             relative to the observer.

               The following values of `abcorr' apply to the
               "reception" case in which photons depart from the
               target's location at the light-time corrected epoch
               et-ltime and *arrive* at the observer's location at `et':

                  'LT'       Correct for one-way light time (also
                             called "planetary aberration") using a
                             Newtonian formulation. This correction
                             yields the state of the target at the
                             moment it emitted photons arriving at
                             the observer at `et'.

                             The light time correction uses an
                             iterative solution of the light time
                             equation (see -Particulars for details).
                             The solution invoked by the 'LT' option
                             uses one iteration.

                  'LT+S'     Correct for one-way light time and
                             stellar aberration using a Newtonian
                             formulation. This option modifies the
                             state obtained with the 'LT' option to
                             account for the observer's velocity
                             relative to the solar system
                             barycenter. The result is the apparent
                             state of the target---the position and
                             velocity of the target as seen by the
                             observer.

                  'CN'       Converged Newtonian light time
                             correction. In solving the light time
                             equation, the 'CN' correction iterates
                             until the solution converges (three
                             iterations on all supported platforms).
                             Whether the 'CN+S' solution is
                             substantially more accurate than the
                             'LT' solution depends on the geometry
                             of the participating objects and on the
                             accuracy of the input data. In all
                             cases this routine will execute more
                             slowly when a converged solution is
                             computed. See the -Particulars section
                             below for a discussion of precision of
                             light time corrections.

                  'CN+S'     Converged Newtonian light time
                             correction and stellar aberration
                             correction.

               The following values of `abcorr' apply to the
               "transmission" case in which photons *depart* from
               the observer's location at `et' and arrive at the
               target's location at the light-time corrected epoch
               et+ltime:

                  'XLT'      "Transmission" case: correct for
                             one-way light time using a Newtonian
                             formulation. This correction yields the
                             state of the target at the moment it
                             receives photons emitted from the
                             observer's location at `et'.

                  'XLT+S'    "Transmission" case: correct for
                             one-way light time and stellar
                             aberration using a Newtonian
                             formulation  This option modifies the
                             state obtained with the 'XLT' option to
                             account for the observer's velocity
                             relative to the solar system
                             barycenter. The position component of
                             the computed target state indicates the
                             direction that photons emitted from the
                             observer's location must be "aimed" to
                             hit the target.

                  'XCN'      "Transmission" case: converged
                             Newtonian light time correction.

                  'XCN+S'    "Transmission" case: converged
                             Newtonian light time correction and
                             stellar aberration correction.

               Neither special nor general relativistic effects are
               accounted for in the aberration corrections applied
               by this routine.

               Case and blanks are not significant in the string
               `abcorr'.

      obs      the name of an observing body.

               help, obs
                  STRING = Scalar

               Optionally, you may supply the ID code of the object as an
               integer string. For example, both 'EARTH' and '399' are
               legitimate strings to supply to indicate the observer is Earth.

   the call:

      cspice_spkezr, targ, et, ref, abcorr, obs, starg, ltime

   returns:

      starg    a Cartesian state vector representing the position and velocity
               of the target body relative to the specified observer, or an
               N-vector of vectors.

               help, starg
                  DOUBLE = Array[6]   or   DOUBLE = Array[6,N]

               `starg' is corrected for the specified aberrations, and is
               expressed with respect to the reference frame specified by
               `ref'. The first three components of `starg' represent the x-,
               y- and z-components of the target's position; the last three
               components form the corresponding velocity vector.

               The position component of `starg' points from the
               observer's location at `et' to the aberration-corrected
               location of the target. Note that the sense of the
               position vector is independent of the direction of
               radiation travel implied by the aberration
               correction.

               The velocity component of `starg' is the derivative
               with respect to time of the position component of
               `starg'.

               Units are always km and km/sec.

               Non-inertial frames are treated as follows: letting
               `ltcent' be the one-way light time between the observer
               and the central body associated with the frame, the
               orientation of the frame is evaluated at et-ltcent,
               et+ltcent, or `et' depending on whether the requested
               aberration correction is, respectively, for received
               radiation, transmitted radiation, or is omitted.
               `ltcent' is computed using the method indicated by
               `abcorr'.

               `starg' returns with the same measure of vectorization (N)
               as `et'.

      ltime    the one-way light time between the observer and target in
               seconds, or an N-vector of times.

               help, ltime
                  DOUBLE = Scalar   or   DOUBLE = Array[N]

               If the target state is corrected for aberrations, then `ltime'
               is the one-way light time between the observer and the light
               time corrected target location.

               `ltime' returns with the same measure of vectorization (N)
               as `et'.

Parameters


   None.

Examples


   Any numerical results shown for this example may differ between
   platforms as the results depend on the SPICE kernels used as input
   and the machine specific arithmetic implementation.

   1) Load a planetary SPK, and look up the state vector of Mars
      as seen from the Earth in the J2000 frame with aberration
      corrections 'LT+S' (ligth time plus stellar aberration) at
      different epochs.

      Use the meta-kernel shown below to load the required SPICE
      kernels.


         KPL/MK

         File: spkezr_ex1.tm

         This meta-kernel is intended to support operation of SPICE
         example programs. The kernels shown here should not be
         assumed to contain adequate or correct versions of data
         required by SPICE-based user applications.

         In order for an application to use this meta-kernel, the
         kernels referenced here must be present in the user's
         current working directory.

         The names and contents of the kernels referenced
         by this meta-kernel are as follows:

            File name                        Contents
            ---------                        --------
            de430.bsp                        Planetary ephemeris
            mar097.bsp                       Mars satellite ephemeris
            naif0011.tls                     Leapseconds


         \begindata

            KERNELS_TO_LOAD = ( 'de430.bsp',
                                'mar097.bsp',
                                'naif0011.tls' )

         \begintext

         End of meta-kernel


      Example code begins here.


      PRO spkezr_ex1

         ;;
         ;; Load a set of kernels: an SPK file, a PCK
         ;; file and a leapseconds file. Use a meta
         ;; kernel for convenience.
         ;;
         cspice_furnsh, 'spkezr_ex1.tm'

         ;;
         ;; Define parameters for a state lookup:
         ;;
         ;; Return the state vector of Mars (499) as seen from
         ;; Earth (399) in the J2000 frame
         ;; using aberration correction LT+S (light time plus
         ;; stellar aberration) at the epoch
         ;; July 4, 2003 11:00 AM PST.
         ;;
         target   = 'Mars'
         epoch    = 'July 4, 2003 11:00 AM PST'
         frame    = 'J2000'
         abcorr   = 'LT+S'
         observer = 'Earth'

         ;;
         ;; Convert the epoch to ephemeris time.
         ;;
         cspice_str2et, epoch, et

         ;;
         ;; Look-up the state for the defined parameters.
         ;;
         cspice_spkezr, target, et, frame, abcorr, observer, $
                     state, ltime
         ;;
         ;; Output...
         ;;
         print, 'The position of    : ', target
         print, 'As observed from   : ', observer
         print, 'In reference frame : ', frame
         print, 'At epoch           : ', epoch
         print

         ;;
         ;; The first three entries of state contain the
         ;; X, Y, Z position components. The final three contain
         ;; the Vx, Vy, Vz velocity components.
         ;;
         print, 'Scalar'
         print, 'R (kilometers)     : ', state[0:2]
         print, 'V (kilometers/sec) : ', state[3:5]
         print, 'Light time (secs)  : ', ltime
         print, ' between observer'
         print, ' and target'
         print

         ;;
         ;; Create a vector of et's, starting at 'epoch'
         ;; in steps of 100000 ephemeris seconds.
         ;;
         SIZE   = 5
         vec_et = dindgen( SIZE )*100000. + et

         ;;
         ;; Look up the 'state' vectors and light time values
         ;; 'ltime'  corresponding to the vector of input
         ;; ephemeris time 'vec_et'.
         ;;
         ;;
         cspice_spkezr, target  , vec_et, frame, abcorr, $
                        observer, state , ltime

         print, 'Vector'
         cspice_et2utc, vec_et, 'C', 3, vec_epoch

         ;;
         ;; When called with a vector 'et', cspice_spkezr returns
         ;; 'state' as an 6xN matrix. Extract the ith state from the
         ;; matrix as:
         ;;
         ;;   state_i = state[*,i]
         ;;
         for i=0, 4 do begin

            print, 'At epoch (UTC)     : ', vec_epoch[i]
            print, 'R (kilometers)     : ', (state[*,i])[0:2]
            print, 'V (kilometers/sec) : ', (state[*,i])[3:5]
            print, 'Light time (secs)  : ', ltime[i]
            print, ' between observer'
            print, ' and target'
            print

         endfor

         ;;
         ;; It's always good form to unload kernels after use,
         ;; particularly in IDL due to data persistence.
         ;;
         cspice_kclear

      END


      When this program was executed on a Mac/Intel/IDL8.x/64-bit
      platform, the output was:


      The position of    : Mars
      As observed from   : Earth
      In reference frame : J2000
      At epoch           : July 4, 2003 11:00 AM PST

      Scalar
      R (kilometers)     :        73822235.      -27127919.      -18741306.
      V (kilometers/sec) :       -6.8085133       7.5139962       3.0012985
      Light time (secs)  :        269.68988
       between observer
       and target

      Vector
      At epoch (UTC)     : 2003 JUL 04 19:00:00.000
      R (kilometers)     :        73822235.      -27127919.      -18741306.
      V (kilometers/sec) :       -6.8085133       7.5139962       3.0012985
      Light time (secs)  :        269.68988
       between observer
       and target

      At epoch (UTC)     : 2003 JUL 05 22:46:40.000
      R (kilometers)     :        73140185.      -26390525.      -18446763.
      V (kilometers/sec) :       -6.8312194       7.2341558       2.8896967
      Light time (secs)  :        266.56404
       between observer
       and target

      At epoch (UTC)     : 2003 JUL 07 02:33:20.000
      R (kilometers)     :        72456240.      -25681031.      -18163339.
      V (kilometers/sec) :       -6.8464804       6.9560179       2.7789284
      Light time (secs)  :        263.48035
       between observer
       and target

      At epoch (UTC)     : 2003 JUL 08 06:20:00.000
      R (kilometers)     :        71771127.      -24999260.      -17890947.
      V (kilometers/sec) :       -6.8546121       6.6797297       2.6690806
      Light time (secs)  :        260.43952
       between observer
       and target

      At epoch (UTC)     : 2003 JUL 09 10:06:40.000
      R (kilometers)     :        71085544.      -24345021.      -17629491.
      V (kilometers/sec) :       -6.8559458       6.4053554       2.5602008
      Light time (secs)  :        257.44220
       between observer
       and target


Particulars


   This routine is part of the user interface to the SPICE ephemeris
   system. It allows you to retrieve state information for any
   ephemeris object relative to any other in a reference frame that
   is convenient for further computations.

   This routine is identical in function to the routine cspice_spkez except
   that it allows you to refer to ephemeris objects by name (via a
   character string).

   Please refer to the Aberration Corrections Required Reading
   abcorr.req for detailed information describing the nature and
   calculation of the applied corrections.

Exceptions


   1)  If name of target or observer cannot be translated to its NAIF
       ID code, the error SPICE(IDCODENOTFOUND) is signaled by a
       routine in the call tree of this routine.

   2)  If the reference frame `ref' is not a recognized reference
       frame, the error SPICE(UNKNOWNFRAME) is signaled by a routine
       in the call tree of this routine.

   3)  If the loaded kernels provide insufficient data to compute the
       requested state vector, an error is signaled by a routine in
       the call tree of this routine.

   4)  If an error occurs while reading an SPK or other kernel file,
       the error  is signaled by a routine in the call tree
       of this routine.

   5)  If any of the input arguments, `targ', `et', `ref', `abcorr'
       or `obs', is undefined, an error is signaled by the IDL error
       handling system.

   6)  If any of the input arguments, `targ', `et', `ref', `abcorr'
       or `obs', is not of the expected type, or it does not have the
       expected dimensions and size, an error is signaled by the Icy
       interface.

   7)  If any of the output arguments, `starg' or `ltime', is not a
       named variable, an error is signaled by the Icy interface.

Files


   This routine computes states using SPK files that have been loaded into
   the SPICE system, normally via the kernel loading interface routine
   cspice_furnsh. See the routine cspice_furnsh and the SPK and KERNEL
   Required Reading for further information on loading (and unloading)
   kernels.

   If the output state `starg' is to be expressed relative to a
   non-inertial frame, or if any of the ephemeris data used to
   compute `starg' are expressed relative to a non-inertial frame in
   the SPK files providing those data, additional kernels may be
   needed to enable the reference frame transformations required to
   compute the state. These additional kernels may be C-kernels, PCK
   files or frame kernels. Any such kernels must already be loaded
   at the time this routine is called.

Restrictions


   None.

Required_Reading


   ABCORR.REQ
   FRAMES.REQ
   ICY.REQ
   NAIF_IDS.REQ
   SPK.REQ
   TIME.REQ

Literature_References


   None.

Author_and_Institution


   N.J. Bachman        (JPL)
   J. Diaz del Rio     (ODC Space)
   E.D. Wright         (JPL)

Version


   -Icy Version 1.1.6, 02-NOV-2021 (JDR)

       Edited the header to comply with NAIF standard. Added example's
       problem statement and meta-kernel.

       Updated -Particulars section.

       Added -Parameters, -Exceptions, -Files, -Restrictions,
       -Literature_References and -Author_and_Institution sections.

       Removed reference to the routine's corresponding CSPICE header from
       -Abstract section.

       Added arguments' type and size information in the -I/O section.

   -Icy Version 1.1.5, 22-DEC-2015 (EDW)

       Corrected typo in -Version section. The 07-NOV-2013
       version entry lacked a version ID.

       Particulars updated to refer to Aberration Corrections
       Required Reading document.

   -Icy Version 1.1.4, 07-JUL-2014 (NJB) (EDW)

       Discussion of light time corrections was updated. Assertions
       that converged light time corrections are unlikely to be
       useful were removed.

   -Icy Version 1.1.3, 07-NOV-2013 (EDW)

       Added aberration algorithm explanation to -Particulars section.

   -Icy Version 1.1.2, 22-DEC-2008 (EDW)

       Header edits performed to improve argument descriptions.
       These descriptions should now closely match the descriptions
       in the corresponding CSPICE routine.

   -Icy Version 1.1.1, 29-SEP-2007 (EDW)

       Replaced the comment fragment in the -I/O section

          "return with the same order"

       with

          "return with the same measure of vectorization"

   -Icy Version 1.1.0, 01-AUG-2004 (EDW)

        Added capability to process vector 'et' input
        returning a matrix 'starg' and vector 'ltime'
        on output.

   -Icy Version 1.0.0, 16-JUN-2003 (EDW)

Index_Entries


   using body names get target state relative to an observer
   get state relative to observer corrected for aberrations
   read ephemeris data
   read trajectory data



Fri Dec 31 18:43:07 2021