void spkw09_c ( SpiceInt handle,
SpiceInt body,
SpiceInt center,
ConstSpiceChar * frame,
SpiceDouble first,
SpiceDouble last,
ConstSpiceChar * segid,
SpiceInt degree,
SpiceInt n,
ConstSpiceDouble states[][6],
ConstSpiceDouble epochs[] )
Write a type 9 segment to an SPK file.
NAIF_IDS
SPC
SPK
TIME
EPHEMERIS
FILES
Variable I/O Description
-------- --- --------------------------------------------------
handle I Handle of an SPK file open for writing.
body I NAIF code for an ephemeris object.
center I NAIF code for center of motion of body.
frame I Reference frame name.
first I Start time of interval covered by segment.
last I End time of interval covered by segment.
segid I Segment identifier.
degree I Degree of interpolating polynomials.
n I Number of states.
states I Array of states.
epochs I Array of epochs corresponding to states.
maxdeg P Maximum allowed degree of interpolating polynomial.
handle is the file handle of an SPK file that has been
opened for writing.
body is the NAIF integer code for an ephemeris object
whose state relative to another body is described
by the segment to be created.
center is the NAIF integer code for the center of motion
of the object identified by body.
frame is the NAIF name for a reference frame
relative to which the state information for body
is specified.
first,
last are, respectively, the start and stop times of
the time interval over which the segment defines
the state of body.
segid is the segment identifier. An SPK segment
identifier may contain up to 40 characters.
degree is the degree of the Lagrange polynomials used to
interpolate the states. All components of the
state vectors are interpolated by polynomials of
fixed degree.
n is the number of states in the input state vector
array.
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.
epochs is an array of epochs corresponding to the members
of the state array. The epochs are specified as
seconds past J2000, TDB.
None. See $Particulars for a description of the effect of this
routine.
MAXDEG is the maximum allowed degree of the interpolating
polynomial. If the value of MAXDEG is increased,
the CSPICE routine spkpvn_ must be changed
accordingly. In particular, the size of the
record passed to spkrNN_ and spkeNN_ must be
increased, and comments describing the record size
must be changed.
The current value of MAXDEG is 15.
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.
2) If the last non-blank character of segid occurs past index 40,
the error SPICE(SEGIDTOOLONG) is signaled.
3) If segid contains any nonprintable characters, the error
SPICE(NONPRINTABLECHARS) is signaled.
4) If degree is not at least 1 or is greater than MAXDEG, the
error SPICE(INVALIDDEGREE) is signaled.
5) If the number of states n is not at least degree+1, the error
SPICE(TOOFEWSTATES) will be signaled.
6) If first is greater than or equal to last then the error
SPICE(BADDESCRTIMES) will be signaled.
7) If the elements of the array epochs are not in strictly
increasing order, the error SPICE(TIMESOUTOFORDER) will be
signaled.
8) If the first epoch epochs[0] is greater than first, the error
SPICE(BADDESCRTIMES) will be signaled.
9) If the last epoch epochs[n] is less than last, the error
SPICE(BADDESCRTIMES) will be signaled.
10) If either the input frame or segment ID string pointer is null,
the error SPICE(NULLPOINTER) is signaled.
11) If either the input frame or segment ID string is empty,
the error SPICE(EMPTYSTRING) is signaled.
A new type 9 SPK segment is written to the SPK file attached
to handle.
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.
Suppose that you have states and are prepared to produce
a segment of type 09 in an SPK file.
The following code fragment could be used to add the new segment
to a previously opened SPK file attached to HANDLE. The file must
have been opened with write access.
#include "SpiceUsr.h"
.
.
.
/.
Create a segment identifier.
./
#define SEGID "MY_SAMPLE_SPK_TYPE_9_SEGMENT"
/.
Write the segment.
./
spkw09_c ( handle, body, center, frame,
first, last, segid, degree,
n, states, epochs );
None.
None.
K.R. Gehringer (JPL)
N.J. Bachman (JPL)
J.M. Lynch (JPL)
W.L. Taber (JPL)
-CSPICE Version 1.0.0, 21-JUN-1999 (KRG) (NJB) (JML) (WLT)
write spk type_9 ephemeris data segment
Link to routine spkw09_c source file spkw09_c.c
|