1 /*
2 
3 -Procedure spk14a_c ( SPK, add data to a type 14 segment )
4 
5 -Abstract
6 
7    Add data to a type 14 SPK segment associated with handle. See
8    also spk14b_c and spk14e_c.
9 
10 -Disclaimer
11 
12    THIS SOFTWARE AND ANY RELATED MATERIALS WERE CREATED BY THE
13    CALIFORNIA INSTITUTE OF TECHNOLOGY (CALTECH) UNDER A U.S.
14    GOVERNMENT CONTRACT WITH THE NATIONAL AERONAUTICS AND SPACE
15    ADMINISTRATION (NASA). THE SOFTWARE IS TECHNOLOGY AND SOFTWARE
16    PUBLICLY AVAILABLE UNDER U.S. EXPORT LAWS AND IS PROVIDED "AS-IS"
17    TO THE RECIPIENT WITHOUT WARRANTY OF ANY KIND, INCLUDING ANY
18    WARRANTIES OF PERFORMANCE OR MERCHANTABILITY OR FITNESS FOR A
19    PARTICULAR USE OR PURPOSE (AS SET FORTH IN UNITED STATES UCC
20    SECTIONS 2312-2313) OR FOR ANY PURPOSE WHATSOEVER, FOR THE
21    SOFTWARE AND RELATED MATERIALS, HOWEVER USED.
22 
23    IN NO EVENT SHALL CALTECH, ITS JET PROPULSION LABORATORY, OR NASA
24    BE LIABLE FOR ANY DAMAGES AND/OR COSTS, INCLUDING, BUT NOT
25    LIMITED TO, INCIDENTAL OR CONSEQUENTIAL DAMAGES OF ANY KIND,
26    INCLUDING ECONOMIC DAMAGE OR INJURY TO PROPERTY AND LOST PROFITS,
27    REGARDLESS OF WHETHER CALTECH, JPL, OR NASA BE ADVISED, HAVE
28    REASON TO KNOW, OR, IN FACT, SHALL KNOW OF THE POSSIBILITY.
29 
30    RECIPIENT BEARS ALL RISK RELATING TO QUALITY AND PERFORMANCE OF
31    THE SOFTWARE AND ANY RELATED MATERIALS, AND AGREES TO INDEMNIFY
32    CALTECH AND NASA FOR ALL THIRD-PARTY CLAIMS RESULTING FROM THE
33    ACTIONS OF RECIPIENT IN THE USE OF THE SOFTWARE.
34 
35 -Required_Reading
36 
37    SPK
38 
39 -Keywords
40 
41    SPK
42 
43 */
44 
45    #include "SpiceUsr.h"
46    #include "SpiceZfc.h"
47    #include "SpiceZmc.h"
48    #undef    spk14a_c
49 
50 
spk14a_c(SpiceInt handle,SpiceInt ncsets,ConstSpiceDouble coeffs[],ConstSpiceDouble epochs[])51    void spk14a_c ( SpiceInt           handle,
52                    SpiceInt           ncsets,
53                    ConstSpiceDouble   coeffs [],
54                    ConstSpiceDouble   epochs []  )
55 
56 /*
57 
58 -Brief_I/O
59 
60    VARIABLE  I/O  DESCRIPTION
61    --------  ---  --------------------------------------------------
62    handle     I   The handle of an SPK file open for writing.
63    ncsets     I   The number of coefficient sets and epochs.
64    coeffs     I   The collection of coefficient sets.
65    epochs     I   The epochs associated with the coefficient sets.
66 
67 -Detailed_Input
68 
69    handle   is the file handle of an SPK file that has been
70             opened for writing.
71 
72    ncsets   is the number of Chebyshev coefficient sets and epochs
73             to be stored in the segment.
74 
75    coeffs   contains a time-ordered array of Chebyshev coefficient
76             sets for computing the state vector of a body, packed one
77             after the other into an array.  A state vector contains
78             the position, X, Y, Z coordinates, and the velocities,
79             dX/dt, dY/dt, dZ/dt, for the position of a body relative
80             to a center of motion.
81 
82             See the Particulars section for details on how to store
83             the coefficient sets in the array.
84 
85    epochs   contains the initial epochs (ephemeris seconds past
86             J2000) corresponding to the Chebyshev coefficients in
87             COEFFS. The I'th epoch is associated with the I'th
88             Chebyshev coefficient set. The epochs must form a
89             strictly increasing sequence.
90 
91 -Detailed_Output
92 
93    None.    The ephemeris data is stored in a segment in the SPK file
94             associated with handle.
95 
96             See the Particulars section for details about the
97             structure of a type 14 SPK segment.
98 
99 -Parameters
100 
101    None.
102 
103 -Particulars
104 
105    This routine adds data to a type 14 SPK segment that is associated
106    with the input argument handle. The segment must have been started
107    by a call to the routine spk14b_c, the routine which begins a type
108    14 SPK segment.
109 
110    This routine is one of a set of three routines for creating and
111    adding data to type 14 SPK segments. These routines are:
112 
113       spk14b_c:  Begin a type 14 SPK segment. This routine must be
114                  called before any data may be added to a type 14
115                  segment.
116 
117       spk14a_c:  Add data to a type 14 SPK segment. This routine may be
118                  called any number of times after a call to spk14b_c to
119                  add type 14 records to the SPK segment that was
120                  started.
121 
122       spk14e_c:  End a type 14 SPK segment. This routine is called to
123                  make the type 14 segment a permanent addition to the
124                  SPK file. Once this routine is called, no further type
125                  14 records may be added to the segment. A new segment
126                  must be started.
127 
128    A type 14 SPK segment consists of coefficient sets for fixed order
129    Chebyshev polynomials over consecutive time intervals, where the
130    time intervals need not all be of the same length. The Chebyshev
131    polynomials represent the position, X, Y, and Z coordinates, and
132    the velocities, dX/dt, dY/dt, and dZ/dt, of body relative to
133    center.
134 
135    The ephemeris data supplied to the type 14 SPK writer is packed
136    into an array as a sequence of logical records,
137 
138       -----------------------------------------------------
139       | Record 1 | Record 2 | ... | Record N-1 | Record N |
140       -----------------------------------------------------
141 
142    with each record has the following format.
143 
144          ------------------------------------------------
145          |  The midpoint of the approximation interval  |
146          ------------------------------------------------
147          |  The radius of the approximation interval    |
148          ------------------------------------------------
149          |  CHBDEG+1 coefficients for the X coordinate  |
150          ------------------------------------------------
151          |  CHBDEG+1 coefficients for the Y coordinate  |
152          ------------------------------------------------
153          |  CHBDEG+1 coefficients for the Z coordinate  |
154          ------------------------------------------------
155          |  CHBDEG+1 coefficients for the X velocity    |
156          ------------------------------------------------
157          |  CHBDEG+1 coefficients for the Y velocity    |
158          ------------------------------------------------
159          |  CHBDEG+1 coefficients for the Z velocity    |
160          ------------------------------------------------
161 
162 
163 -Examples
164 
165    Assume we have the following for each of the examples that
166    follow.
167 
168       handle   is the handle of an SPK file opened with write
169                access.
170 
171       segid    is a character string of no more than 40 characters
172                which provides a pedigree for the data in the SPK
173                segment we will create.
174 
175       body     is the NAIF ID code for the body whose ephemeris
176                is to be placed into the file.
177 
178       center   is the center of motion for the ephemeris of body.
179 
180       reffrm   is the name of the SPICE reference frame for the
181                ephemeris.
182 
183       first    is the starting epoch, in seconds past J2000, for
184                the ephemeris data to be placed into the segment.
185 
186       last     is the ending epoch, in seconds past J2000, for
187                the ephemeris data to be placed into the segment.
188 
189    Example 1:
190 
191       For this example, we also assume that:
192 
193          n        is the number of type 14 records that we want to
194                   put into a segment in an SPK file.
195 
196          recrds   contains n type 14 records packaged for the SPK
197                   file.
198 
199          etstrt   contains the initial epochs for each of the
200                   records contained in RECRDS, where
201 
202                      etstrt[i] < etstrt[i+1], i = 0, n-2
203 
204                      etstrt[1] <= first, etstrt[n-1] < last
205 
206                      etstrt[i+1], i = 0, n-2, is the ending epoch for
207                      record i as well as the initial epoch for record
208                      i+1.
209 
210       Then the following code fragment demonstrates how to create a
211       type 14 SPK segment if all of the data for the segment is
212       available at one time.
213 
214          #include "SpiceUsr.h"
215             .
216             .
217             .
218 
219          #define SPK  "example.bsp"
220 
221          /.
222          If the segment is to be appended to an existing file, open
223          that file for "append" access.  Otherwise, create a new file.
224          ./
225 
226          if ( exists_c(SPK) )
227          {
228             spkopa_c ( SPK, &handle );
229          }
230          else
231          {
232             /.
233             New files are supplied with an internal file name.
234             Comment area space may be reserved at this time; the
235             units are characters.
236             ./
237             ifname = "Sample type 14 SPK file.";
238             ncomch = 1024;
239 
240             spkopn_c ( SPK, ifname, ncomch, &handle );
241          }
242 
243 
244          /.
245          Begin the segment.
246          ./
247          spk14b_c ( handle, segid, body, center, reffrm,
248                     first,  last,  chbdeg               );
249 
250          /.
251          Add the data to the segment all at once.
252          ./
253          spk14a_c ( handle, n, recrds, etstrt );
254 
255          /.
256          End the segment, making the segment a permanent addition
257          to the SPK file.
258          ./
259          spk14e_c ( handle );
260 
261              .
262              .
263              .
264          /.
265          After all segments have been loaded, close the SPK file.
266          ./
267          spkcls_c ( handle );
268 
269 
270    Example 2:
271 
272       In this example we want to add type 14 SPK records, as described
273       above in the Particulars section, to the segments being written
274       as they are generated.  The ability to write the records in this
275       way is useful if computer memory is limited. It may also be
276       convenient from a programming perspective to write the records
277       one at a time.
278 
279       For this example, assume that we want to generate n type 14 SPK
280       records, one for each of n time intervals, writing them all to
281       the same segment in the SPK file. Let
282 
283          n        be the number of type 14 records that we want to
284                   generate and put into a segment in an SPK file.
285 
286          record   be an array with enough room to hold a single type
287                   14 record, i.e. record should have dimension at
288                   least 6 * (chbdeg + 1 ) + 2.
289 
290          start    be an array of n times that are the beginning
291                   epochs for each of the intervals of interest. The
292                   times should be in increasing order and the start
293                   time for the first interval should equal the
294                   starting time for the segment.
295 
296                      start[i] < start[i+1], i = 0, n-2
297 
298                      start[0] = first
299 
300          stop     be an array of n times that are the ending epochs
301                   for each of the intervals of interest. The times
302                   should be in increasing order and the stop time for
303                   interval i should equal the start time for interval
304                   i+1, i.e., we want to have continuous coverage in
305                   time across all of the records. Also, the stop time
306                   for the last interval should equal the ending time
307                   for the segment.
308 
309                      stop[i]   < stop [i+1], i = 0, n-2
310 
311                      stop[i]   = start[i+1], i = 0, n-2
312 
313                      stop[n-1] = last
314 
315 
316          genrec( time1, time2, record )
317 
318                   be a subroutine that generates a type 14 SPK record
319                   for a time interval specified by time1 and time2.
320 
321 
322       Then the following code fragment demonstrates how to create a
323       type 14 SPK segment if all of the data for the segment is not
324       available at one time.
325 
326          #include "SpiceUsr.h"
327             .
328             .
329             .
330 
331          /.
332          Begin the segment.
333          ./
334          spk14b_c ( handle, segid, body, center, reffrm,
335                     first,  last,  chbdeg                );
336 
337 
338          /.
339          Generate the records and write them to the segment in the
340          SPK file one at at time.
341          ./
342 
343          for ( i = 0;  i < n;  i++ )
344          {
345             genrec   ( start[i],     stop[i],  record   );
346             spk14a_c ( handle,   1,  record,   start+i  );
347          }
348 
349          /.
350          End the segment, making the segment a permanent addition
351          to the SPK file.
352          ./
353          spk14e_c ( handle );
354 
355 
356 
357 -Restrictions
358 
359    1) The type 14 SPK segment to which we are adding data must have
360       been started by the routine SPK14B, the routine which begins a
361       type 14 SPK segment.
362 
363 -Exceptions
364 
365    1) If the number of coefficient sets and epochs is not positive,
366       the error SPICE(INVALIDARGUMENT) will be signaled.
367 
368 -Files
369 
370    None.
371 
372 -Author_and_Institution
373 
374    N.J. Bachman        (JPL)
375    K.R. Gehringer      (JPL)
376 
377 -Literature_References
378 
379    None.
380 
381 -Version
382 
383    -CSPICE Version 1.0.0, 29-JUL-1999 (NJB) (KRG)
384 
385 -Index_Entries
386 
387    add data to a type_14 spk segment
388 
389 -&
390 */
391 
392 { /* Begin spk14a_c */
393 
394 
395 
396    /*
397    Participate in error tracing.
398    */
399    chkin_c ( "spk14a_c" );
400 
401 
402    spk14a_ ( ( integer    * ) &handle,
403              ( integer    * ) &ncsets,
404              ( doublereal * ) coeffs,
405              ( doublereal * ) epochs  );
406 
407 
408    chkout_c ( "spk14a_c" );
409 
410 } /* End spk14a_c */
411