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

Procedure
Abstract
Required_Reading
Keywords
Brief_I/O
Detailed_Input
Detailed_Output
Parameters
Exceptions
Files
Particulars
Examples
Restrictions
Literature_References
Author_and_Institution
Version
Index_Entries

Procedure

   void spk14e_c ( SpiceInt   handle ) 

Abstract

 
   End the type 14 SPK segment currently being written to the SPK 
   file associated with handle. See also spk14a_c and spk14b_c. 
 

Required_Reading

 
   None. 
 

Keywords

 
   SPK 
 

Brief_I/O

 
   VARIABLE  I/O  DESCRIPTION 
   --------  ---  -------------------------------------------------- 
   handle     I   The handle of an SPK file open for writing. 
 

Detailed_Input

 
   handle   is the file handle of an SPK file that has been 
            opened for writing, and to which a type 14 segment is 
            being written. 
 

Detailed_Output

 
   None.    The type 14 segment in the SPK file associated with 
            handle will be ended, making the addition of the data 
            to the file permanent. 
 
            See the Particulars section for details about the 
            structure of a type 14 SPK segment. 
 

Parameters

 
   None. 
 

Exceptions

 
   None. 
 

Files

 
   See the argument handle. 
 

Particulars

 
   This routine ends a type 14 SPK segment which is being written to 
   the SPK file associated with handle. Ending the SPK segment is a 
   necessary step in the process of making the data a permanent part 
   of the SPK file. 
 
   This routine is one of a set of three routines for creating and 
   adding data to type 14 SPK segments. These routines are: 
 
      spk14b_c: Begin a type 14 SPK segment. This routine must be 
                called before any data may be added to a type 14 
                segment. 
 
      spk14a_c: Add data to a type 14 SPK segment. This routine may be 
                called any number of times after a call to SPK14B to 
                add type 14 records to the SPK segment that was 
                started. 
 
      spk14e_c: 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 a body relative to a
   center of motion. 
 
   The ephemeris data supplied to the type 14 SPK writer routines 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   | 
      ------------------------------------------------ 
 

Examples

 
   Assume we have the following for each of the examples that 
   follow. 
 
      handle   is the handle of an SPK file opened with write 
               access. 
 
      segid    is a character string of no more than 40 characters 
               which provides a pedigree for the data in the SPK 
               segment we will create. 
 
      body     is the NAIF ID code for the body whose ephemeris 
               is to be placed into the file. 
 
      center   is the center of motion for the ephemeris of body. 
 
      reffrm   is the name of the SPICE reference frame for the 
               ephemeris. 
 
      first    is the starting epoch, in seconds past J2000, for 
               the ephemeris data to be placed into the segment. 
 
      last     is the ending epoch, in seconds past J2000, for 
               the ephemeris data to be placed into the segment. 
 
   Example 1: 
 
      For this example, we also assume that: 
 
         n        is the number of type 14 records that we want to 
                  put into a segment in an SPK file. 
 
         recrds   contains n type 14 records packaged for the SPK 
                  file. 
 
         etstrt   contains the initial epochs for each of the 
                  records contained in RECRDS, where 
 
                     etstrt[i] < etstrt[i+1], i = 0, n-2 
 
                     etstrt[1] <= first, etstrt[n-1] < last 
 
                     etstrt[i+1], i = 0, n-2, is the ending epoch for 
                     record i as well as the initial epoch for record 
                     i+1. 
 
      Then the following code fragment demonstrates how to create a 
      type 14 SPK segment if all of the data for the segment is 
      available at one time. 
 
         #include "SpiceUsr.h"
            .
            .
            .
         
         #define SPK  "example.bsp"
            
         /.
         If the segment is to be appended to an existing file, open
         that file for "append" access.  Otherwise, create a new file.
         ./
         
         if ( exists_c(SPK) )
         {
            spkopa_c ( SPK, &handle );
         } 
         else
         {
            /.
            New files are supplied with an internal file name.  
            Comment area space may be reserved at this time; the
            units are characters.
            ./
            ifname = "Sample type 14 SPK file.";
            ncomch = 1024;
            
            spkopn_c ( SPK, ifname, ncomch, &handle );
         }
          
          
         /.
         Begin the segment. 
         ./
         spk14b_c ( handle, segid, body, center, reffrm, 
                    first,  last,  chbdeg               );
       
         /.
         Add the data to the segment all at once. 
         ./
         spk14a_c ( handle, n, recrds, etstrt ); 
      
         /.
         End the segment, making the segment a permanent addition 
         to the SPK file. 
         ./
         spk14e_c ( handle ); 
         
             .
             .
             .
         /.
         After all segments have been loaded, close the SPK file.
         ./
         spkcls_c ( handle );
         
 
   Example 2: 
 
      In this example we want to add type 14 SPK records, as described
      above in the Particulars section, to the segments being written
      as they are generated.  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
      one at a time.
 
      For this example, assume that we want to generate n type 14 SPK 
      records, one for each of n time intervals, writing them all to 
      the same segment in the SPK file. Let 
 
         n        be the number of type 14 records that we want to 
                  generate and put into a segment in an SPK file. 
 
         record   be an array with enough room to hold a single type 
                  14 record, i.e. record should have dimension at 
                  least 6 * (chbdeg + 1 ) + 2. 
 
         start    be an array of n times that are the beginning 
                  epochs for each of the intervals of interest. The 
                  times should be in increasing order and the start 
                  time for the first interval should equal the 
                  starting time for the segment. 
 
                     start[i] < start[i+1], i = 0, n-2 
 
                     start[0] = first 
 
         stop     be an array of n times that are the ending epochs 
                  for each of the intervals of interest. The times 
                  should be in increasing order and the stop time for 
                  interval i should equal the start time for interval 
                  i+1, i.e., we want to have continuous coverage in 
                  time across all of the records. Also, the stop time 
                  for the last interval should equal the ending time 
                  for the segment. 
 
                     stop[i]   < stop [i+1], i = 0, n-2 
 
                     stop[i]   = start[i+1], i = 0, n-2 
 
                     stop[n-1] = last 
 

         genrec( time1, time2, record ) 
 
                  be a subroutine that generates a type 14 SPK record 
                  for a time interval specified by time1 and time2. 
 

      Then the following code fragment demonstrates how to create a 
      type 14 SPK segment if all of the data for the segment is not 
      available at one time. 
 
         #include "SpiceUsr.h"
            .
            .
            .
        
         /.
         Begin the segment. 
         ./
         spk14b_c ( handle, segid, body, center, reffrm, 
                    first,  last,  chbdeg                ); 
 
        
         /.
         Generate the records and write them to the segment in the 
         SPK file one at at time. 
         ./   
         
         for ( i = 0;  i < n;  i++ )
         {
            genrec   ( start[i],    stop[i], record  ); 
            spk14a_c ( handle,   1, record,  start+i );
         }
 
         /.
         End the segment, making the segment a permanent addition 
         to the SPK file. 
         ./   
         spk14e_c ( handle );
 

Restrictions

 
   1) The type 14 SPK segment being closed must have been started by 
      the routine spk14b_c, the routine which begins a type 14 SPK 
      segment. 
 

Literature_References

 
   None. 
 

Author_and_Institution

 
   N.J. Bachman        (JPL)
   K.R. Gehringer      (JPL) 
 

Version

 
   -CSPICE Version 1.0.1, 16-JAN-2003 (EDW)

       Trivial correction to the header.

   -CSPICE Version 1.0.0, 29-JUL-1999 (NJB) (KRG)

Index_Entries

 
   end a type_14 spk segment 
 

Link to routine spk14e_c source file spk14e_c.c

Wed Apr  5 17:54:43 2017