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_spk14a

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


Abstract


   CSPICE_SPK14A adds data to a type 14 SPK segment associated with
   `handle'.

I/O


   Given:

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

               help, handle
                  LONG = Scalar

      ncsets   the number of Chebyshev coefficient sets and epochs to be stored
               in the segment.

               help, ncsets
                  LONG = Scalar

      coeffs   a time-ordered array of `ncsets' of Chebyshev polynomial
               coefficients to be placed in the segment of the SPK file.

               help, coeffs
                  DOUBLE = Array[ncsets*SETSZ]

               Each set has size SETSZ = 2 + 6*(CHBDEG+1), where CHBDEG is the
               degree of the Chebyshev polynomials used to represent the
               ephemeris information.

               These sets are used to compute the state vector, which
               consists of position components X, Y, Z and velocity
               components dx/dt, dy/dt, dz/dt, of a body relative to
               a center of motion.

               See the -Particulars section for details on how to store
               the coefficient sets in the array.

      epochs   contains the initial epochs (ephemeris seconds past J2000)
               corresponding to the Chebyshev coefficients in `coeffs'.

               help, epochs
                  DOUBLE = Array[ncsets]

               The I'th epoch is associated with the I'th Chebyshev
               coefficient set. The epochs must form a strictly increasing
               sequence.

   the call:

      cspice_spk14a, handle, ncsets, coeffs, epochs

   stores the ephemeris data in a segment in the SPK file associated with
   `handle'.

   See the -Particulars section for details about the structure of a
   type 14 SPK segment.

Parameters


   None.

Examples


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

   1) This example demonstrates how to create an SPK type 14 kernel
      containing only one segment, given a set of Chebyshev
      coefficients and their associated epochs, if all of the data
      for the segment is available at one time.


      Example code begins here.


      PRO spk14a_ex1

         ;;
         ;; Create an SPK type 14. First create a segment
         ;; identifier.
         ;;
         body      = 3
         center    = 10
         ref       = 'J2000'
         segid     = 'SPK type 14 test segment'
         SPK14     = 'spk14a_ex1.bsp'
         N_RECORDS = 4
         CHBDEG    = 2

         ;;
         ;; Define the coefficients.
         ;;
         discreteEpochs = [ 100.d, 200.d, 300.d, 400.d, 500.d ]

         ChebyRecords14 = [                            $
                           150.0d,                     $
                           50.0d,                      $
                           1.0101d, 1.0102d, 1.0103d,  $
                           1.0201d, 1.0202d, 1.0203d,  $
                           1.0301d, 1.0302d, 1.0303d,  $
                           1.0401d, 1.0402d, 1.0403d,  $
                           1.0501d, 1.0502d, 1.0503d,  $
                           1.0601d, 1.0602d, 1.0603d   $
                          ,                            $
                           250.0d,                     $
                           50.0d,                      $
                           2.0101d, 2.0102d, 2.0103d,  $
                           2.0201d, 2.0202d, 2.0203d,  $
                           2.0301d, 2.0302d, 2.0303d,  $
                           2.0401d, 2.0402d, 2.0403d,  $
                           2.0501d, 2.0502d, 2.0503d,  $
                           2.0601d, 2.0602d, 2.0603d   $
                          ,                            $
                           350.0d,                     $
                           50.0d,                      $
                           3.0101d, 3.0102d, 3.0103d,  $
                           3.0201d, 3.0202d, 3.0203d,  $
                           3.0301d, 3.0302d, 3.0303d,  $
                           3.0401d, 3.0402d, 3.0403d,  $
                           3.0501d, 3.0502d, 3.0503d,  $
                           3.0601d, 3.0602d, 3.0603d   $
                          ,                            $
                           450.0d,                     $
                           50.0d,                      $
                           4.0101d, 4.0102d, 4.0103d,  $
                           4.0201d, 4.0202d, 4.0203d,  $
                           4.0301d, 4.0302d, 4.0303d,  $
                           4.0401d, 4.0402d, 4.0403d,  $
                           4.0501d, 4.0502d, 4.0503d,  $
                           4.0601d, 4.0602d, 4.0603d   $
                        ]

         ;;
         ;; Open a new SPK file.
         ;;
         cspice_spkopn, SPK14, 'Type 14 SPK internal file name.', 4, handle

         ;;
         ;; Create a type 14 segment.
         ;;
         ;; Begin the segment
         ;;
         cspice_spk14b, handle,                      $
                        segid,                       $
                        body,                        $
                        center,                      $
                        ref,                         $
                        discreteEpochs[0],           $
                        discreteEpochs[N_RECORDS],   $
                        CHBDEG

         ;;
         ;; Add the data.
         ;;
         cspice_spk14a, handle,           $
                        N_RECORDS,        $
                        ChebyRecords14,   $
                        discreteEpochs

         ;;
         ;; End the segment.
         ;;
         cspice_spk14e, handle

         ;;
         ;; SAFELY close the SPK file.
         ;;
         cspice_spkcls, handle

      END


      When this program is executed, no output is presented on
      screen. After run completion, a new SPK type 14 exists in
      the output directory.

   2) This example demonstrates how to add type 14 SPK records to the
      segment being written, one at a time. The ability to write the
      records in this way is useful if computer memory is limited. It
      may also be convenient from a programming perspective to write
      the records this way.


      Example code begins here.


      PRO spk14a_ex2

         ;;
         ;; Create an SPK type 14. First create a segment
         ;; identifier.
         ;;
         body      = 3
         center    = 10
         ref       = 'J2000'
         segid     = 'SPK type 14 test segment'
         SPK14     = 'spk14a_ex2.bsp'
         N_RECORDS = 4
         CHBDEG    = 2

         ;;
         ;; Define the coefficients.
         ;;
         discreteEpochs = [ 100.d, 200.d, 300.d, 400.d, 500.d ]

         ChebyRecords14 = [                            $
                          [150.0d,                     $
                           50.0d,                      $
                           1.0101d, 1.0102d, 1.0103d,  $
                           1.0201d, 1.0202d, 1.0203d,  $
                           1.0301d, 1.0302d, 1.0303d,  $
                           1.0401d, 1.0402d, 1.0403d,  $
                           1.0501d, 1.0502d, 1.0503d,  $
                           1.0601d, 1.0602d, 1.0603d]  $
                          ,                            $
                          [250.0d,                     $
                           50.0d,                      $
                           2.0101d, 2.0102d, 2.0103d,  $
                           2.0201d, 2.0202d, 2.0203d,  $
                           2.0301d, 2.0302d, 2.0303d,  $
                           2.0401d, 2.0402d, 2.0403d,  $
                           2.0501d, 2.0502d, 2.0503d,  $
                           2.0601d, 2.0602d, 2.0603d]  $
                          ,                            $
                          [350.0d,                     $
                           50.0d,                      $
                           3.0101d, 3.0102d, 3.0103d,  $
                           3.0201d, 3.0202d, 3.0203d,  $
                           3.0301d, 3.0302d, 3.0303d,  $
                           3.0401d, 3.0402d, 3.0403d,  $
                           3.0501d, 3.0502d, 3.0503d,  $
                           3.0601d, 3.0602d, 3.0603d]  $
                          ,                            $
                          [450.0d,                     $
                           50.0d,                      $
                           4.0101d, 4.0102d, 4.0103d,  $
                           4.0201d, 4.0202d, 4.0203d,  $
                           4.0301d, 4.0302d, 4.0303d,  $
                           4.0401d, 4.0402d, 4.0403d,  $
                           4.0501d, 4.0502d, 4.0503d,  $
                           4.0601d, 4.0602d, 4.0603d]  $
                        ]

         ;;
         ;; Open a new SPK file.
         ;;
         cspice_spkopn, SPK14, 'Type 14 SPK internal file name.', 4, handle

         ;;
         ;; Create a type 14 segment.
         ;;
         ;; Begin the segment
         ;;
         cspice_spk14b, handle,                      $
                        segid,                       $
                        body,                        $
                        center,                      $
                        ref,                         $
                        discreteEpochs[0],           $
                        discreteEpochs[N_RECORDS],   $
                        CHBDEG

         ;;
         ;;  Write the records to the segment in the
         ;;  SPK file one at at time.
         ;;
         for i = 0, N_RECORDS-1 do begin

            cspice_spk14a, handle, 1, ChebyRecords14[*,i], [discreteEpochs[i]]

         endfor

         ;;
         ;; End the segment.
         ;;
         cspice_spk14e, handle

         ;;
         ;; SAFELY close the SPK file.
         ;;
         cspice_spkcls, handle

      END


      When this program is executed, no output is presented on
      screen. After run completion, a new SPK type 14 exists in
      the output directory.

Particulars


   This routine adds data to a type 14 SPK segment that is associated
   with the input argument handle. The segment must have been started
   by a call to the routine cspice_spk14b, the routine which begins a type
   14 SPK segment.

   This routine is one of a set of three routines for creating and
   adding data to type 14 SPK segments. These routines are:

      cspice_spk14b: Begin a type 14 SPK segment. This routine must be
                     called before any data may be added to a type 14
                     segment.

      cspice_spk14a: Add data to a type 14 SPK segment. This routine may be
                     called any number of times after a call to cspice_spk14b
                     to add type 14 records to the SPK segment that was
                     started.

      cspice_spk14e: End a type 14 SPK segment. This routine is called to
                     make the type 14 segment a permanent addition to the
                     SPK file. Once this routine is called, no further type
                     14 records may be added to the segment. A new segment
                     must be started.

   A type 14 SPK segment consists of coefficient sets for fixed order
   Chebyshev polynomials over consecutive time intervals, where the
   time intervals need not all be of the same length. The Chebyshev
   polynomials represent the position, X, Y, and Z coordinates, and
   the velocities, dX/dt, dY/dt, and dZ/dt, of body relative to
   center.

   The ephemeris data supplied to the type 14 SPK writer is packed
   into an array as a sequence of logical records,

      -----------------------------------------------------
      | Record 1 | Record 2 | ... | Record N-1 | Record N |
      -----------------------------------------------------

   with each record has the following format.

         ------------------------------------------------
         |  The midpoint of the approximation interval  |
         ------------------------------------------------
         |  The radius of the approximation interval    |
         ------------------------------------------------
         |  chbdeg+1 coefficients for the X coordinate  |
         ------------------------------------------------
         |  chbdeg+1 coefficients for the Y coordinate  |
         ------------------------------------------------
         |  chbdeg+1 coefficients for the Z coordinate  |
         ------------------------------------------------
         |  chbdeg+1 coefficients for the X velocity    |
         ------------------------------------------------
         |  chbdeg+1 coefficients for the Y velocity    |
         ------------------------------------------------
         |  chbdeg+1 coefficients for the Z velocity    |
         ------------------------------------------------

Exceptions


   1)  If the number of coefficient sets and epochs is not positive,
       the error SPICE(INVALIDARGUMENT) is signaled by a routine in
       the call tree of this routine.

   2)  If any of the input arguments, `handle', `ncsets', `coeffs' or
       `epochs', is undefined, an error is signaled by the IDL error
       handling system.

   3)  If any of the input arguments, `handle', `ncsets', `coeffs' 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


   None.

Restrictions


   1)  The type 14 SPK segment to which we are adding data must have
       been started by the routine cspice_spk14b, the routine which begins a
       type 14 SPK segment.

Required_Reading


   ICY.REQ
   SPK.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.0.1, 27-AUG-2021 (JDR) (NJB)

       Edited the header to comply with NAIF standard. Added example's problem
       statement.

       Extended "coeffs" argument description to provide the size of the
       Chebyshev polynomials sets.

       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


   add data to a type_14 SPK segment



Fri Dec 31 18:43:07 2021