void spkw18_c ( SpiceInt handle,
SpiceSPK18Subtype subtyp,
SpiceInt body,
SpiceInt center,
ConstSpiceChar * frame,
SpiceDouble first,
SpiceDouble last,
ConstSpiceChar * segid,
SpiceInt degree,
SpiceInt n,
const void * packts,
ConstSpiceDouble epochs[] )
Write a type 18 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.
subtyp I SPK type 18 subtype code.
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.
subtyp is an integer code indicating the subtype of the
the segment to be created.
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 nominal degree of the polynomials used to
interpolate the states contained in the input
packets. All components of the state vectors are
interpolated by polynomials of the specified
degree, except near the segment boundaries,
or if the total number of states in the segment
is too few to allow interpolation using the
specified degree.
n is the number of packets in the input packet
array.
packts contains a time-ordered array of data packets
representing geometric states of body relative
to center, specified relative to frame. The
packet structure depends on the segment subtype
as follows:
Type 0 (indicated by code S18TP0):
x, y, z, dx/dt, dy/dt, dz/dt,
vx, vy, vz, dvx/dt, dvy/dt, dvz/dt
where x, y, z represent Cartesian position
components and vx, vy, vz represent Cartesian
velocity components. Note well: vx, vy, and
vz *are not necessarily equal* to the time
derivatives of x, y, and z. This packet
structure mimics that of the Rosetta/MEX orbit
file from which the data are taken.
Type 1 (indicated by code S18TP1):
x, y, z, dx/dt, dy/dt, dz/dt
where x, y, z represent Cartesian position
components and vx, vy, vz represent Cartesian
velocity components.
Position units are kilometers, velocity units
are kilometers per second, and acceleration units
are kilometers per second per second.
epochs is an array of epochs corresponding to the members
of the packets 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.
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 window size implied by DEGREE is odd, the error
SPICE(INVALIDDEGREE) is signaled.
6) If the number of packets n is not at least 1,
the error SPICE(TOOFEWSTATES) will be signaled.
7) If first is greater than or equal to last then the error
SPICE(BADDESCRTIMES) will be signaled.
8) If the elements of the array epochs are not in strictly
increasing order, the error SPICE(TIMESOUTOFORDER) will be
signaled.
9) If the first epoch epochs[0] is greater than first, the error
SPICE(BADDESCRTIMES) will be signaled.
10) If the last epoch epochs[n-1] is less than last, the error
SPICE(BADDESCRTIMES) will be signaled.
11) If either the input frame or segment ID string pointer is null,
the error SPICE(NULLPOINTER) is signaled.
12) If either the input frame or segment ID string is empty,
the error SPICE(EMPTYSTRING) is signaled.
A new type 18 SPK segment is written to the SPK file attached
to HANDLE.
This routine writes an SPK type 18 data segment to the open SPK
file according to the format described in the type 18 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 18 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_18_SEGMENT"
/.
Write the segment.
./
spkw18_c ( handle, subtyp, body, center,
frame, first, last, segid,
degree, n, states, epochs );
None.
None.
N.J. Bachman (JPL)
-CSPICE Version 1.0.1, 29-APR-2003 (NJB)
Description of error condition arising from invalid window
size was corrected.
-CSPICE Version 1.0.0, 16-AUG-2002 (NJB)
write spk type_18 ephemeris data segment
Link to routine spkw18_c source file spkw18_c.c
|