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_spkcpo

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


Abstract


   CSPICE_SPKCPO returns the state of a specified target relative to
   an "observer," where the observer has constant position in a
   specified reference frame. The observer's position is provided
   by the calling program rather than by loaded SPK files.

I/O


   Given:

      target   scalar string name of a target body.

               help, target
                  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 target is earth.

               Case and leading and trailing blanks are not significant
               in the string `target'.


      et       scalar double precision ephemeris time at which the state of the
               target relative to the observer is to be computed.

               help, et
                  DOUBLE = Scalar

               `et' is expressed as seconds past J2000 TDB. `et' refers to
               time at the observer's location.


      outref   scalar string name of the reference frame with respect to which
               the output state is expressed.

               help, outref
                  STRING = Scalar

               When `outref' is time-dependent (non-inertial), its
               orientation relative to the J2000 frame is evaluated in
               the manner commanded by the input argument `refloc' (see
               description below).

               Case and leading and trailing blanks are not significant
               in the string `outref'.


      refloc   scalar string name indicating the output reference frame
               evaluation locus: this is the location associated with the epoch
               at which this routine is to evaluate the orientation, relative
               to the J2000 frame, of the output frame `outref'.

               help, refloc
                  STRING = Scalar

               The values and meanings of `refloc' are:

                  'OBSERVER'  Evaluate `outref' at the observer's
                              epoch `et'.

                              Normally the locus 'OBSERVER' should
                              be selected when `outref' is centered
                              at the observer.


                  'TARGET'    Evaluate `outref' at the target epoch;
                              letting `ltime' be the one-way light time
                              between the target and observer, the
                              target epoch is

                                 et-ltime  if reception aberration
                                           corrections are used

                                 et+ltime  if transmission aberration
                                           corrections are used

                                 et        if no aberration corrections
                                           are used

                              Normally the locus 'TARGET' should
                              be selected when `outref' is centered
                              at the target object.


                  'CENTER'    Evaluate the frame `outref' at the epoch
                              associated its center. This epoch,
                              which we'll call `etctr', is determined
                              as follows:

                                 Let `ltctr' be the one-way light time
                                 between the observer and the center
                                 of `outref'. Then `etctr' is

                                    et-ltctr  if reception
                                              aberration corrections
                                              are used

                                    et+ltctr  if transmission
                                              aberration corrections
                                              are used

                                    et        if no aberration
                                              corrections are used


                              The locus 'CENTER' should be selected
                              when the user intends to obtain
                              results compatible with those produced
                              by cspice_spkezr.

               When `outref' is inertial, all choices of `refloc'
               yield the same results.

               Case and leading and trailing blanks are not
               significant in the string `refloc'.

      abcorr   scalar string name indicating the aberration corrections to be
               applied to the observer-target state to account for one-way
               light time and stellar aberration.

               help, abcorr
                  STRING = Scalar

               `abcorr' may be any of the following:

                  'NONE'     Apply no correction. Return the
                             geometric state of the target
                             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. 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.

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


               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 and stellar
                             aberration corrections.


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

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


      obspos   double precision 3-vector defining the fixed (constant)
               geometric position of an observer relative to its center of
               motion `obsctr', expressed in the reference frame `obsref'.

               help, obspos
                  DOUBLE = Array[3]

               Units are always km.

      obsctr   scalar string name of the center of motion of `obspos'.

               help, obsctr
                  STRING = Scalar

               The ephemeris of `obsctr' is provided by loaded SPK files.

               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 center of motion.

               Case and leading and trailing blanks are not
               significant in the string `obsctr'.


      obsref   scalar string name of the reference frame relative to which the
               input position `obspos' is expressed.

               help, obsref
                  STRING = Scalar

               The observer has constant position relative to its center of
               motion in this reference frame.

               Case and leading and trailing blanks are not
               significant in the string `obsref'.

   the call:

      cspice_spkcpo, target, et,     outref, refloc, $
                     abcorr, obspos, obsctr,         $
                     obsref, state,  ltime

   returns:

      state    a double precision Cartesian 6-vector representing the position
               and velocity of the target relative to the specified observer.

               help, state
                  DOUBLE = Array[6]

               `state' is corrected for the specified aberrations and is
               expressed with respect to the reference frame specified by
               `outref'. The first three components of `state' 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 `state' 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 `state' is the derivative
               with respect to time of the position component of
               `state'.

               Units are always km and km/sec.

               When `state' is expressed in a time-dependent
               (non-inertial) output frame, the orientation of that
               frame relative to the J2000 frame is evaluated in the
               manner indicated by the input argument `refloc' (see
               description above).

      ltime    scalar double precision one-way light time between the observer
               and target in seconds.

               help, ltime
                  DOUBLE = Scalar

               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.

               Please note, CSPICE documentation and source code
               uniformly uses the variable name `lt' to designate
               the light-time between an observer and target. IDL
               uses "lt" as the less-than numeric comparison
               operator and so does not allow "lt" as a variable name.
               Therefore, Icy documentation uses the name `ltime'
               for the light-time value.

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) Compute apparent solar azimuth and elevation as seen from a
      specified surface point on the earth.

      Task Description
      ================

      In this example we'll use the location of the DSN station
      DSS-14 as our surface point.

      We'll perform the solar azimuth and elevation computation two
      ways:

         - Using a station frame kernel to provide the
           specification of a topocentric reference frame
           centered at DSS-14.

         - Computing inline the transformation from the earth-fixed,
           earth-centered frame ITRF93 to a topocentric frame
           centered at DSS-14.

      Note that results of the two computations will differ
      slightly. There are three sources of the differences:

         1) The station position is time-dependent due to tectonic
            plate motion, and epochs of the station positions used
            to specify the axes of the topocentric frame are
            different in the two cases. This gives rise to different
            orientations of the frame's axes relative to the frame
            ITRF93.

         2) The two computations use different earth radii; this
            results in computation of different geodetic latitudes
            of the station. This difference also affects the
            topocentric frame orientation relative to ITRF93.

         3) The station movement between `et' and the epoch at which
            the DSS-14_TOPO frame is specified contributes a very
            small offset---on the order of 10 cm---to the station-sun
            position vector, expressed in the ITRF93 frame.


      Kernels
      =======

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


         KPL/MK

         File name: spkcpo_ex1.tm

         This is the meta-kernel file for the header code example for
         the subroutine cspice_spkcpo. These kernel files can be found on
         the NAIF website.

         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
            ---------                        --------
            de421.bsp                        Planetary ephemeris
            pck00010.tpc                     Planet orientation and
                                             radii
            naif0010.tls                     Leapseconds
            earth_720101_070426.bpc          Earth historical
                                             binary PCK
            earthstns_itrf93_050714.bsp      DSN station SPK
            earth_topo_050714.tf             DSN station FK
            mgs_moc_v20.ti                   MGS MOC instrument
                                             parameters
            mgs_sclkscet_00061.tsc           MGS SCLK coefficients
            mgs_sc_ext12.bc                  MGS s/c bus attitude
            mgs_ext12_ipng_mgs95j.bsp        MGS ephemeris

         \begindata

         KERNELS_TO_LOAD = ( 'de421.bsp',
                             'pck00010.tpc',
                             'naif0010.tls',
                             'earth_720101_070426.bpc',
                             'earthstns_itrf93_050714.bsp',
                             'earth_topo_050714.tf',
                             'mgs_moc_v20.ti',
                             'mgs_sclkscet_00061.tsc',
                             'mgs_sc_ext12.bc',
                             'mgs_ext12_ipng_mgs95j.bsp'  )

         \begintext

         End of meta-kernel.


      Example code begins here.


      ;;
      ;; This program uses cspice_spkcpo to compute solar azimuth
      ;; and elevation at a given surface point on the earth.
      ;;
      PRO spkcpo_ex1

         ;;
         ;; Local parameters.
         ;;
         META   =  'spkcpo_ex1.tm'
         TIMFMT =  'YYYY MON DD HR:MN:SC.###### UTC'
         TIMFM2 =  'YYYY MON DD HR:MN:SC.###### TDB ::TDB'
         TIMLEN =  41

         ;;
         ;; Initial values
         ;;
         z = [ 0.0D, 0.0, 1.0 ]


         ;;
         ;; Load SPICE kernels.
         ;;
         cspice_furnsh, META

         ;;
         ;; Convert the observation time to seconds past J2000 TDB.
         ;;
         obstim = '2003 OCT 13 06:00:00.000000 UTC';

         cspice_str2et, obstim, et

         ;;
         ;; Set the target, observer center, observer frame, and
         ;; observer position relative to its center.
         ;;
         target = 'SUN'
         obsctr = 'EARTH'
         obsref = 'ITRF93'

         ;;
         ;; Set the position of DSS-14 relative to the earth's
         ;; center at the J2000 epoch, expressed in the
         ;; ITRF93 reference frame. Values come from the
         ;; earth station SPK specified in the meta-kernel.
         ;;
         ;; The actual station velocity is non-zero due
         ;; to tectonic plate motion; we ignore the motion
         ;; in this example. See the routine cspice_spkcvo for an
         ;; example in which the plate motion is accounted for.
         ;;
         obspos =  [ -2353.6213656676991D0,                                   $
                     -4641.3414911499403D0,                                   $
                      3677.0523293197439D0 ]

         ;;
         ;; Find the apparent state of the sun relative
         ;; to the station in the DSS-14_TOPO reference frame.
         ;; Evaluate the output frame's orientation, that is the
         ;; orientation of the DSS-14_TOPO frame relative to the
         ;; J2000 frame, at the observation epoch. This
         ;; correction is obtained by setting `refloc' to
         ;; 'OBSERVER'.
         ;;
         outref = 'DSS-14_TOPO'
         abcorr = 'CN+S'

         refloc = 'OBSERVER'

         ;;
         ;; Compute the observer-target state.
         ;;
         cspice_spkcpo, target, et,     outref, refloc, abcorr,              $
                        obspos, obsctr, obsref, state0, lt0

         ;;
         ;; Compute planetocentric coordinates of the
         ;; observer-target position in the local
         ;; topocentric reference frame DSS-14_TOPO.
         ;;
         cspice_reclat, state0[0:2], r, lon, lat

         ;;
         ;; Compute solar azimuth. The latitude we've
         ;; already computed is the elevation. Express
         ;; both angles in degrees.
         ;;
         el =   lat * cspice_dpr()
         az = - lon * cspice_dpr()

         if ( az lt 0.d0 ) then begin
            az +=  360.0d0
         endif

         ;;
         ;; Display the computed state, light time, and angles.
         ;;
         cspice_timout, et-lt0, TIMFMT, TIMLEN, emitim

         print, ' Frame evaluation locus:     ', refloc
         print, ' '
         print, ' Target:                     ', target
         print, ' Observation time:           ', obstim
         print, ' Observer center:            ', obsctr
         print, ' Observer frame:             ', obsref
         print, ' Emission time:              ', emitim
         print, ' Output reference frame:     ', outref
         print, ' Aberration correction:      ', abcorr
         print, ' '
         print, ' Observer-target position (km):'
         print, format = '(3F20.8)', state0[0:2]
         print, ' Observer-target velocity (km/s):'
         print, format = '(3F20.8)', state0[3:5]
         print, format='(" Light time (s):        ",F20.8)', lt0

         print, ' '
         print, format = '(" Solar azimuth (deg):     ", F20.8)', az
         print, format = '(" Solar elevation (deg):   ", F20.8)', el

         ;;
         ;; For an arbitrary surface point, we might not
         ;; have a frame kernel available. In this case
         ;; we can look up the state in the observer frame
         ;; using cspice_spkcpo and then convert the state to
         ;; the local topocentric frame. We'll first
         ;; create the transformation matrix for converting
         ;; vectors in the observer frame to the topocentric
         ;; frame.
         ;;
         ;; First step: find the geodetic (planetodetic)
         ;; coordinates of the observer. We need the
         ;; equatorial radius and flattening coefficient
         ;; of the reference ellipsoid.
         ;;
         cspice_bodvrd, 'EARTH', 'RADII', 3, radii

         re = radii[0]
         rp = radii[2]

         f  = ( re - rp ) / re

         cspice_recgeo, obspos, re, f, obslon, obslat, obsalt

         ;;
         ;; Find the outward surface normal on the reference
         ;; ellipsoid at the observer's longitude and latitude.
         ;;
         cspice_latrec, 1.D0, obslon, obslat, normal

         ;;
         ;; The topocentric frame has its +Z axis aligned
         ;; with `normal' and its +X axis pointed north.
         ;; The north direction is aligned with the component
         ;; of the ITRF93 +Z axis orthogonal to the topocentric
         ;; +Z axis.
         ;;
         cspice_twovec, normal, 3, z, 1, xform

         outref = 'ITRF93'
         abcorr = 'CN+S'

         refloc = 'OBSERVER'

         ;;
         ;; Compute the observer-target state.
         ;;
         cspice_spkcpo, target, et,     outref, refloc, abcorr,              $
                        obspos, obsctr, obsref, state1, lt1

         ;;
         ;; Convert the position to the topocentric frame.
         ;;
         cspice_mxv, xform, state1[0:2], topvec

         ;;
         ;; Compute azimuth and elevation.
         ;;
         cspice_reclat, topvec, r, lon, lat

         el =   lat * cspice_dpr()
         az = - lon * cspice_dpr()

         if ( az lt 0.d0 ) then begin
            az +=  360.d0
         endif

         print, ' '
         print, ' AZ/EL computed without frame kernel:'
         print, ' Distance between last two '
         print, format='(" positions (km):   ", F20.8)',                     $
                                  cspice_vdist( state0[0:2], topvec )

         print, ' '
         print, format = '(" Solar azimuth (deg):     ", F20.8)', az
         print, format = '(" Solar elevation (deg):   ", F20.8)', el

         ;;
         ;; 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:


       Frame evaluation locus:     OBSERVER

       Target:                     SUN
       Observation time:           2003 OCT 13 06:00:00.000000 UTC
       Observer center:            EARTH
       Observer frame:             ITRF93
       Emission time:              2003 OCT 13 05:51:42.068322 UTC
       Output reference frame:     DSS-14_TOPO
       Aberration correction:      CN+S

       Observer-target position (km):
         62512272.82074845   58967494.42513601 -122059095.46751881
       Observer-target velocity (km/s):
             2475.97326517      -9870.26706232      -3499.90809969
       Light time (s):                497.93167797

       Solar azimuth (deg):             316.67141599
       Solar elevation (deg):           -54.85253168

       AZ/EL computed without frame kernel:
       Distance between last two
       positions (km):             3.07056970

       Solar azimuth (deg):             316.67141786
       Solar elevation (deg):           -54.85253216


Particulars


   This routine computes observer-target states for observers whose
   trajectories are not provided by SPK files.

   Observers supported by this routine must have constant position
   with respect to a specified center of motion, expressed in a
   caller-specified reference frame. The state of the center of
   motion relative to the target must be computable using
   loaded SPK data.

   For applications in which the observer has constant, non-zero velocity
   relative to its center of motion, the Icy routine

      cspice_spkcvo     { SPK, constant velocity observer state }

   can be used.

   This routine is suitable for computing states of target ephemeris
   objects, as seen from landmarks on the surface of an extended
   object, in cases where no SPK data are available for those
   landmarks.

   This routine's treatment of the output reference frame differs
   from that of the principal SPK API routines

      cspice_spkezr
      cspice_spkez
      cspice_spkpos
      cspice_spkezp

   which require both observer and target ephemerides to be provided
   by loaded SPK files:

      The SPK API routines listed above evaluate the orientation of the
      output reference frame (with respect to the J2000 frame) at an
      epoch corrected for one-way light time between the observer and
      the center of the output frame. When the center of the output
      frame is not the target (for example, when the target is on the
      surface of Mars and the output frame is centered at Mars'
      center), the epoch of evaluation may not closely match the
      light-time corrected epoch associated with the target itself. A
      similar problem may occur when the observer is a surface point on
      an extended body and the output frame is centered at the body
      center: the listed routines will correct the orientation of the
      output frame for one-way light time between the frame center and
      the observer.

      This routine allows the caller to dictate how the orientation
      of the output reference frame is to be evaluated. The caller
      passes to this routine an input string called the output
      frame's evaluation "locus." This string specifies the location
      associated with the output frame's evaluation epoch. The three
      possible values of the locus are

         'TARGET'
         'OBSERVER'
         'CENTER'

      The choice of locus has an effect when aberration corrections
      are used and the output frame is non-inertial.

      When the locus is 'TARGET' and light time corrections are
      used, the orientation of the output frame is evaluated at the
      epoch obtained by correcting the observation epoch `et' for
      one-way light time `ltime'. The evaluation epoch will be either
      et-ltime or et+ltime for reception or transmission corrections
      respectively.

      For remote sensing applications where the target is a surface
      point on an extended object, and the orientation of that
      object should be evaluated at the emission time, the locus
      'TARGET' should be used.

      When the output frame's orientation should be evaluated at
      the observation epoch `et', which is the case when the
      output frame is centered at the observer, the locus
      'OBSERVER' should be used.

      The locus option 'CENTER' is provided for compatibility
      with existing SPK state computation APIs such as cspice_spkezr.

      Note that the output frame evaluation locus does not affect
      the computation of light time between the target and
      observer.


   The SPK routines that compute observer-target states for
   combinations of objects having ephemerides provided by the SPK
   system and objects having constant position or constant velocity
   are

      cspice_spkcpo {SPK, Constant position observer}
      cspice_spkcpt {SPK, Constant position target}
      cspice_spkcvo {SPK, Constant velocity observer}
      cspice_spkcvt {SPK, Constant velocity target}

Exceptions


   1)  If either the name of the center of motion or the target
       cannot be translated to its NAIF ID code, an error is signaled
       by a routine in the call tree of this routine.

   2)  If the reference frame `outref' is unrecognized, an error is
       signaled by a routine in the call tree of this routine.

   3)  If the reference frame `obsref' is unrecognized, an error is
       signaled by a routine in the call tree of this routine.

   4)  If the frame evaluation locus `refloc' is not recognized, an
       error is signaled by a routine in the call tree of this
       routine.

   5)  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.

   6)  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.

   7)  If the aberration correction `abcorr' is not recognized, an
       error is signaled by a routine in the call tree of this
       routine.

   8)  If any of the input arguments, `target', `et', `outref',
       `refloc', `abcorr', `obspos', `obsctr' or `obsref', is
       undefined, an error is signaled by the IDL error handling
       system.

   9)  If any of the input arguments, `target', `et', `outref',
       `refloc', `abcorr', `obspos', `obsctr' or `obsref', is not of
       the expected type, or it does not have the expected dimensions
       and size, an error is signaled by the Icy interface.

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

Files


   Appropriate kernels must be loaded by the calling program before
   this routine is called.

   The following data are required:

   -  SPK data: ephemeris data for the observer center and target
      must be loaded. If aberration corrections are used, the
      states of the observer center and target relative to the
      solar system barycenter must be calculable from the
      available ephemeris data. Typically ephemeris data are made
      available by loading one or more SPK files using cspice_furnsh.

   The following data may be required:

   -  PCK data: if the target frame is a PCK frame, rotation data
      for the target frame must be loaded. These may be provided
      in a text or binary PCK file.

   -  Frame data: if a frame definition not built into SPICE is
      required, for example to convert the observer-target state
      to the output frame, that definition must be available in
      the kernel pool. Typically frame definitions are supplied
      by loading a frame kernel using cspice_furnsh.

   -  Additional kernels: if any frame used in this routine's
      state computation is a CK frame, then at least one CK and
      corresponding SCLK kernel is required. If dynamic frames
      are used, additional SPK, PCK, CK, or SCLK kernels may be
      required.

   In all cases, kernel data are normally loaded once per program
   run, NOT every time this routine is called.

Restrictions


   1)  This routine may not be suitable for work with stars or other
       objects having large distances from the observer, due to loss
       of precision in position vectors.

Required_Reading


   FRAMES.REQ
   ICY.REQ
   PCK.REQ
   SPK.REQ
   TIME.REQ

Literature_References


   None.

Author_and_Institution


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

Version


   -Icy Version 1.0.1, 01-JUN-2021 (JDR)

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

       Edited to the header to comply with NAIF standard. Added
       example's task description.

       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.0.0, 09-APR-2012 (EDW)

Index_Entries


   state relative to constant_position_observer
   state relative to constant_position surface_point
   state relative to surface_point on extended_object
   state relative to landmark on extended_object



Fri Dec 31 18:43:07 2021