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_spkw09

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


Abstract


   CSPICE_SPKW09 writes a type 9 SPK segment to an SPK file

I/O


   Given:

      handle   the file handle of an SPK file that has been opened for writing.

               help, handle
                  LONG = Scalar

      body     the NAIF integer code for an ephemeris object whose state
               relative to another body is described by the segment to be
               created.

               help, body
                  LONG = Scalar

      center   the NAIF integer code for the center of motion of the object
               identified by `body'.

               help, center
                  LONG = Scalar

      frame    the NAIF name for a reference frame relative to which the state
               information for `body' is specified.

               help, frame
                  STRING = Scalar

      first,
      last     respectively, the start and stop times of the time interval over
               which the segment defines the state of `body'.

               help, first
                  DOUBLE = Scalar
               help, last
                  DOUBLE = Scalar

      segid    the segment identifier.

               help, segid
                  STRING = Scalar

               An SPK segment identifier may contain up to 40 characters.

      degree   the degree of the Lagrange polynomials used to interpolate the
               states.

               help, degree
                  LONG = Scalar

               All components of the state vectors are interpolated by
               polynomials of fixed degree.

      n        the number of states in the input state vector array.

               help, n
                  LONG = Scalar

      states   contains a time-ordered array of geometric states ( x, y, z,
               dx/dt, dy/dt, dz/dt, in kilometers and kilometers per second )
               of `body' relative to `center', specified relative to `frame'.

               help, states
                  DOUBLE = Array[6,n]

      epochs   an array of epochs corresponding to the members of the state
               array.

               help, epochs
                  DOUBLE = Array[n]

               The epochs are specified as seconds past J2000, TDB.

   the call:

      cspice_spkw09, handle, body,   center, frame,  first, last,            $
                     segid,  degree, n,      states, epochs

   writes to the SPK file referred to by `handle' a type 9 SPK segment
   containing the data listed in `states'.

Parameters


   MAXDEG      is the maximum allowed degree of the interpolating
               polynomial.

               The current value of MAXDEG is 15.

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) Suppose that you have a time-ordered array of geometric states
      of a new object that follows Phobos, with a delay of 1 hour,
      in its orbit around Mars and are prepared to produce a segment
      of type 09 in an SPK file. Create a new SPK file with this
      segment. Use an existing SPK to create the input data for the
      SPK segment.

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


         KPL/MK

         File: spkw09_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
            ---------                        --------
            mar097.bsp                       Mars satellite ephemeris
            naif0012.tls                     Leapseconds

         \begindata

            KERNELS_TO_LOAD = ( 'mar097.bsp',
                                'naif0012.tls' )

         \begintext

         End of meta-kernel


      Example code begins here.


      PRO spkw09_ex1

         ;;
         ;; Local parameters.
         ;;
         SPKNAM = 'spkw09_ex1.bsp'
         DEGREE = 3
         MARS = 499
         NEPOCS = 800
         NOBJ = 403

         ;;
         ;; Load the input SPK file.
         ;;
         cspice_furnsh, 'spkw09_ex1.tm'

         ;;
         ;; Convert the input UTC to ephemeris time
         ;;
         cspice_str2et, '2018 Apr 03 08:35', et

         ;;
         ;; Create the time-ordered array of geometric states,
         ;; at unequal time steps.
         ;;
         time   = et
         step   = 60.0
         delta  = 10.0
         epochs = DBLARR( NEPOCS )
         states = DBLARR( 6, NEPOCS )

         for i=0, NEPOCS - 1L do begin

            cspice_spkezr, 'PHOBOS', time, 'J2000', 'NONE', 'MARS', $
                            state,   ltime

            states[*,i] = state
            epochs[i] = time + 3600.0

            time = time + step +                            $
                   sin( cspice_halfpi( ) * i / 2.0 ) * delta

         endfor

         ;;
         ;; Open a new SPK file, with 5000 characters reserved
         ;; for comments.
         ;;
         ifname = 'Test SPK type 9 internal filename.'
         cspice_spkopn, SPKNAM, ifname, 5000, handle

         ;;
         ;; Create a segment identifier.
         ;;
         segid = 'MY_SAMPLE_SPK_TYPE_9_SEGMENT'

         ;;
         ;; Write the segment.
         ;;
         cspice_spkw09, handle,    NOBJ,              MARS, 'J2000', $
                        epochs[0], epochs[NEPOCS-1L], segid, DEGREE, $
                        NEPOCS,    states,            epochs

         ;;
         ;; Close the new SPK file.
         ;;
         cspice_spkcls, handle

         ;;
         ;; Compute the state of Phobos as seen from Mars,
         ;; 12 hours after the input UTC time.
         ;;
         et = et + 43200.0
         cspice_spkezr, 'PHOBOS', et, 'J2000', 'NONE', 'MARS', $
                         state,   ltime

         print, format='(A)', 'Phobos as seen from Mars'
         print, format='(A,F20.6)', '   Epoch       (s):', et
         print, format='(A,3F14.6)', '   Position   (km):', state[0:2]
         print, format='(A,3F14.6)', '   Velocity (km/s):', state[3:5]
         print

         ;;
         ;; Load the newly created kernel, and compute the state
         ;; of the new object as seen from Mars, 13 hours after
         ;; the input UTC time.
         ;;
         cspice_furnsh, SPKNAM
         et = et + 3600.0

         cspice_spkezr, '403', et, 'J2000', 'NONE', 'MARS',  $
                        state, ltime

         print, format='(A)', 'Object 403 as seen from Mars'
         print, format='(A,F20.6)', '   Epoch       (s):', et
         print, format='(A,3F14.6)', '   Position   (km):', state[0:2]
         print, format='(A,3F14.6)', '   Velocity (km/s):', state[3:5]

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


      Phobos as seen from Mars
         Epoch       (s):    576059769.185657
         Position   (km):  -7327.262770   2414.326550   5207.106376
         Velocity (km/s):     -0.942893     -1.894731     -0.396715

      Object 403 as seen from Mars
         Epoch       (s):    576063369.185657
         Position   (km):  -7327.262770   2414.326550   5207.106376
         Velocity (km/s):     -0.942893     -1.894731     -0.396715


      Note that after run completion, a new SPK file exists in
      the output directory.

Particulars


   This routine writes an SPK type 09 data segment to the open SPK
   file according to the format described in the type 09 section of
   the SPK Required Reading. The SPK file must have been opened with
   write access.

Exceptions


   If any of the following exceptions occur, this routine will return
   without creating a new segment.

   1)  If `frame' is not a recognized name, the error
       SPICE(INVALIDREFFRAME) is signaled by a routine in the call
       tree of this routine.

   2)  If the last non-blank character of `segid' occurs past index 40,
       the error SPICE(SEGIDTOOLONG) is signaled by a routine in the
       call tree of this routine.

   3)  If `segid' contains any nonprintable characters, the error
       SPICE(NONPRINTABLECHARS) is signaled by a routine in the call
       tree of this routine.

   4)  If `degree' is not at least 1 or is greater than MAXDEG, the
       error SPICE(INVALIDDEGREE) is signaled by a routine in the
       call tree of this routine.

   5)  If the number of states `n' is not at least degree+1, the error
       SPICE(TOOFEWSTATES) is signaled by a routine in the call tree
       of this routine.

   6)  If `first' is greater than or equal to `last', the error
       SPICE(BADDESCRTIMES) is signaled by a routine in the call tree
       of this routine.

   7)  If the elements of the array `epochs' are not in strictly
       increasing order, the error SPICE(TIMESOUTOFORDER) is
       signaled by a routine in the call tree of this routine.

   8)  If the first epoch, epochs[0], is greater than `first', the
       error SPICE(BADDESCRTIMES) is signaled by a routine in the
       call tree of this routine.

   9)  If the last epoch, epochs[n-1], is less than `last', the error
       SPICE(BADDESCRTIMES) is signaled by a routine in the call tree
       of this routine.

   10) If any of the input arguments, `handle', `body', `center',
       `frame', `first', `last', `segid', `degree', `n', `states' or
       `epochs', is undefined, an error is signaled by the IDL error
       handling system.

   11) If any of the input arguments, `handle', `body', `center',
       `frame', `first', `last', `segid', `degree', `n', `states' or
       `epochs', is not of the expected type, or it does not have the
       expected dimensions and size, an error is signaled by the Icy
       interface.

Files


   A new type 9 SPK segment is written to the SPK file attached
   to `handle'.

Restrictions


   None.

Required_Reading


   ICY.REQ
   SPK.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)

       Edited the header to comply with NAIF standard. Modified code example
       to use real input data and produce solution.

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

       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, 16-JUN-2003 (EDW)

Index_Entries


   write SPK type_9 ephemeris data segment



Fri Dec 31 18:43:07 2021