Main Page
Most Used CSPICE APIs

CSPICE APIs for accessing SPICE kernel data

CSPICE APIs for computing derived geometry

CSPICE APIs for checking geometric conditions

CSPICE APIs for finding times of geometric events

CSPICE APIs for coordinate conversions

CSPICE APIs for operations with 3D vectors and matrices

Miscellaneous CSPICE APIs


Top

Loading and Unloading SPICE Kernels

APIs:

  • furnsh_c - loads an individual kernel or a collection of kernels.
  • unload_c - unloads an individual kernel or a collection of kernels.

Brief Example:

   Generic LSK and PCK files listed in a meta-kernel named
   "mykernels.furnsh" and containing

      \begindata
         KERNELS_TO_LOAD = (
                            '/kernels/gen/lsk/naif0008.tls'
                            '/kernels/gen/pck/pck00008.tps'
                           )
      \begintext

   are loaded with a single call to furnsh_c:

      furnsh_c ( "mykernels.furnsh" );

Reference documents:

  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Converting between UTC and Ephemeris Time (LSK)

APIs:

  • furnsh_c - loads an individual kernel or a collection of kernels.
  • str2et_c - converts a time string to ET seconds past J2000.
  • timout_c - converts ET seconds past J2000 to a time string.

Brief Example:

   The following example loads an LSK file, converts a UTC string to ET
   seconds, adds 1 day and converts TDB seconds back to a UTC string in
   ISO DOY format:

      SpiceDouble   et;
      SpiceChar     utc[32];

      /*
         load LSK file
      */
      furnsh_c  ( "naif0008.tls" );

      /*
         convert UTC to ET
      */
      str2et_c  ( "2005 DEC 31 12:00", &et );

      /*
         add 1 day to ET and convert it back to UTC
      */
      timout_c ( et+spd_c(), "YYYY-DOYTHR:MN:SC.### ::RND", 32, utc );

Reference documents:

  • time.req - reference on time systems and conversions supported in SPICE
  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Converting between Ephemeris Time and Spacecraft Clock (SCLK)

APIs:

  • furnsh_c - loads individual kernels and collections of kernels.
  • scs2e_c - converts an SCLK string to ET seconds past J2000.
  • sce2s_c - converts ET seconds past J2000 to SCLK string.
  • sct2e_c - converts an encoded SCLK to ET seconds past J2000.
  • sce2c_c - converts ET seconds past J2000 to encoded SCLK.
  • scencd_c - converts an SCLK string to encoded SCLK.
  • scdecd_c - converts an encoded SCLK to SCLK string.

Brief Example:

   The following example loads an LSK file and an MGS SCLK file, then
   converts an MGS SCLK string to ET and encoded SCLK:

      SpiceChar   * sclk;
      SpiceDouble   et;
      SpiceDouble   sclkdp;
      SpiceInt      scid = -94;

      /*
         load LSK and SCLK files
      */
      furnsh_c ( "naif0008.tls" );
      furnsh_c ( "mgs_sclkscet_00061.tsc" );

      /*
         convert string SCLK to ET and to encoded SCLK
      */
      sclk = "820584056:121";
      scs2e_c  ( scid, sclk, &et );
      scencd_c ( scid, sclk, &sclkdp );

Reference documents:

  • sclk.req - reference on on-board clock (SCLK) implementation in SPICE
  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Constants and Orientation for Natural Bodies (PCK)

APIs:

  • furnsh_c - loads individual kernels and collections of kernels.
  • bodfnd_c - determines whether values exist for some item for a body in the kernel pool.
  • bodvrd_c - retrieves from the kernel pool the values of an item associated with a body.
  • pxform_c - returns the 3x3 matrix rotating a position vector one frame to another.
  • sxform_c - returns the 6x6 matrix rotating a state vector from one frame to another.
  • pckfrm_c - finds the set of reference frame class ID codes of all frames in a binary PCK file.
  • pckcov_c - finds the coverage window for a reference frame in a binary PCK file.

Brief Example:

   The following example retrieves radii for Mars and computes
   orientation of the Mars body-fixed frame:

      SpiceDouble   et;
      SpiceDouble   mat[3][3];
      SpiceDouble   radii[3];
      SpiceInt      n;

      /*
         load LSK and PCK files
      */
      furnsh_c ( "naif0008.tls" );
      furnsh_c ( "pck00008.tpc" );

      /*
         retrieve Mars radii
      */
      bodvrd_c ( "MARS", "RADII", 3, &n, radii );

      /*
         convert UTC to ET
      */
      str2et_c ( "2005-DEC-28 12:00", &et );

      /*
         compute Mars orientation relative to the J2000 frame
      */
      pxform_c ( "J2000", "IAU_MARS", et, mat );

Reference documents:

  • pck.req - reference on Planetary Constants Kernel (PCK) subsystem of SPICE
  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Computing Transformations Between Frames (FK)

APIs:

  • furnsh_c - loads individual kernels and collections of kernels.
  • pxform_c - returns the 3x3 matrix rotating a position vector one frame to another.
  • sxform_c - returns the 6x6 matrix rotating a state vector from one frame to another.
  • pxfrm2_c - returns the 3x3 matrix rotating a position vector from one frame at a specified epoch to another frame at a different epoch.

Brief Example:

   The following example computes orientation of the CASSINI ISS WAC
   camera relative to the Saturn body-fixed frame:

      SpiceDouble   et;
      SpiceDouble   mat[3][3];

      /*
         load kernels: LSK, PCK, CASSINI SCLK, FK and CK
      */
      furnsh_c ( "naif0008.tls" );
      furnsh_c ( "pck00008.tpc" );
      furnsh_c ( "cas00101.tsc" );
      furnsh_c ( "cas_v39.tf" );
      furnsh_c ( "05362_06002ra.bc" );

      /*
         convert UTC to ET
      */
      str2et_c ( "2005-DEC-28 12:00", &et );

      /*
         compute orientation of CASSINI relative to J2000 frame
      */
      pxform_c ( "IAU_SATURN", "CASSINI_ISS_WAC", et, mat );

Reference documents:

  • frames.req - reference on Frames Kernel (FK) and Frames subsystem of SPICE
  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Computing Positions of Spacecraft and Natural Bodies (SPK)

APIs:

  • furnsh_c - loads individual kernels and collections of kernels.
  • spkezr_c - returns the state of a target body relative to an observing body.
  • spkpos_c - returns the position of a target body relative to an observing body.
  • spkcpo_c - returns the state of a target body relative to a constant-position observer location.
  • spkcpt_c - returns the state of a constant-position target location relative to an observing body.
  • spkcvo_c - returns the state of a target body relative to a constant-velocity observer location.
  • spkcvt_c - returns the state of a constant-velocity target location relative to an observing body.
  • spkobj_c - finds the set of ID codes for all objects in a specified SPK file.
  • spkcov_c - finds the coverage window for a specified object in a specified SPK file.

Brief Example:

   The following example computes the geometric state (position and
   velocity) of MGS relative to Mars in the J2000 reference frame:

      SpiceDouble   et;
      SpiceDouble   state[6];
      SpiceDouble   lt;

      /*
         load kernels: LSK and MGS and DE SPKs
      */
      furnsh_c ( "naif0008.tls" );
      furnsh_c ( "de405.bsp" );
      furnsh_c ( "mgs_ext22.bsp" );

      /*
         convert UTC to ET
      */
      str2et_c ( "2006 JAN 31 01:00", &et );

      /*
         compute geometric state of MGS relative to Mars
      */
      spkezr_c ( "MGS", et, "J2000", "NONE", "MARS", state,
                 &lt );

Reference documents:

  • spk.req - reference on Spacecraft and Planet Ephemeris (SPK) subsystem of SPICE
  • frames.req - reference on Frames Kernel (FK) and Frames subsystem of SPICE
  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Computing Orientation for Spacecraft and Instruments (CK)

APIs:

  • furnsh_c - loads individual kernels and collections of kernels.
  • pxform_c - returns the 3x3 matrix rotating a position vector one frame to another.
  • sxform_c - returns the 6x6 matrix rotating a state vector from one frame to another.
  • ckobj_c - finds the set of ID codes for all objects in a specified CK file.
  • ckcov_c - finds the coverage window for a specified object in a specified CK file.
  • ckgp_c - gets pointing for a specified CK ID at a specified SCLK time.
  • ckgpav_c - gets pointing and angular velocity for a specified CK ID at a specified SCLK time.

Brief Example:

   The following example computes orientation of the CASSINI spacecraft
   provided in CK files using the Frames subsystem routine pxform_c:

      SpiceDouble   et;
      SpiceDouble   mat[3][3];

      /*
         load kernels: LSK, CASSINI SCLK, FK and CK
      */
      furnsh_c ( "naif0008.tls" );
      furnsh_c ( "cas00101.tsc" );
      furnsh_c ( "cas_v39.tf" );
      furnsh_c ( "05362_06002ra.bc" );

      /*
         convert UTC to ET
      */
      str2et_c ( "2005-DEC-28 12:00", &et );

      /*
         compute orientation of CASSINI relative to J2000 frame
      */
      pxform_c ( "J2000", "CASSINI_SC_COORD", et, mat );

Reference documents:

  • ck.req - reference on spacecraft orientation (CK) subsystem of SPICE
  • frames.req - reference on Frames Kernel (FK) and Frames subsystem of SPICE
  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Retrieving Instrument Parameters (IK)

APIs:

  • furnsh_c - loads individual kernels and collections of kernels.
  • getfvn_c - returns the field-of-view (FOV) configuration for an instrument specified by name.
  • getfov_c - returns the field-of-view (FOV) configuration for an instrument specified by ID.
  • gdpool_c - returns the double precision value of a kernel variable from the kernel pool.
  • gipool_c - returns the integer value of a kernel variable from the kernel pool.
  • gcpool_c - returns the character value of a kernel variable from the kernel pool.

Brief Example:

   The following example retrieves the M01 THEMIS IR camera focal
   length and FOV parameters:

      SpiceBoolean  found;
      SpiceChar     frame[32];
      SpiceChar     shape[32];
      SpiceDouble   fl;
      SpiceDouble   bsight[3];
      SpiceDouble   bounds[4][3];
      SpiceInt      n;

      /*
         load M01 FK and THEMIS IK
      */
      furnsh_c ( "m01_v28.tf" );
      furnsh_c ( "m01_themis_v32.ti" );

      /*
         retrieve M01 THEMIS IR pixel size
      */
      gdpool_c ( "INS-53031_FOCAL_LENGTH", 0, 1, &n, &fl,
                 &found );

      /*
         retrieve M01 THEMIS IR FOV parameters
      */
      getfvn_c ( "M01_THEMIS_IR", 4, 32, 32, shape, frame, bsight, &n,
                 bounds );

Reference documents:

  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Computing surface coordinates using digital shape (DSK)

APIs:

  • furnsh_c - loads individual kernels and collections of kernels.
  • latsrf_c - maps an array of planetocentric longitude/latitude coordinate pairs to surface points on a body, modeled as an ellipsoid or a digital shape (DSK).
  • srfnrm_c - maps an array of surface points on a body, modeled as an ellipsoid or a digital shape (DSK), to the corresponding outward surface normal vectors.
  • dskz02_c - returns plate model size parameters (plate count and vertex count) for a type 2 DSK segment.
  • dskp02_c - returns triangular plates from a type 2 DSK segment.
  • dskv02_c - returns vertices from a type 2 DSK segment.
  • dskobj_c - returns the set of body ID codes of all objects for which data are provided in a DSK file.
  • dsksrf_c - returns the set of surface ID codes for all surfaces associated with a body in a DSK file.

Brief Example:

   The following example computes surface points and normals for
   for a one-degree longitude/latitude grid on Phobos modeled by
   a DSK:

      SpiceInt                n;
      SpiceInt                i;
      SpiceInt                j;

      SpiceDouble             grid   [360*179][2];
      SpiceDouble             srfpts [360*179][3];
      SpiceDouble             normls [360*179][3];


      /*
         Load generic PCK and Phobos DSK.
      */
      furnsh_c ( "pck00010.tpc" );
      furnsh_c ( "phobos512.bds" );

      /*
         Generate longitude/latitude pairs of the grid (without
         poles).
      */
      n = 0;
      for ( i = 1;  i <= 360;  i++ ) {
         for ( j = 1;  j <= 179;  j++ ) {
            grid[n][0] = rpd_c() *   i;
            grid[n][1] = rpd_c() * ( j - 90 );
            ++n;
         }
      }

      /*
         Compute surface points based on DSK.
      */
      latsrf_c ( "DSK/UNPRIORITIZED", "PHOBOS", 0.0, "IAU_PHOBOS",
                 n, grid, srfpts );

      /*
         Compute surface normals based on DSK.
      */
      srfnrm_c ( "DSK/UNPRIORITIZED", "PHOBOS", 0.0, "IAU_PHOBOS",
                 n, srfpts, normls );

Reference documents:

  • dsk.req - reference on Digital Shape Kernels (DSK) subsystem of SPICE
  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Mapping Between Object Names and NAIF IDs

APIs:

  • furnsh_c - loads individual kernels and collections of kernels.
  • bodc2n_c - translates the NAIF integer code of a body into a common name for that body.
  • bodn2c_c - translates the name of a body or object to the corresponding NAIF integer ID code.

Brief Example:

   The following example uses bodn2c_c to get the NAIF ID for Phobos
   (built into the SPICE system) and for the M01 THEMIS IR camera
   (defined in the M01 FK file):

      SpiceBoolean  found;
      Spiceint      phobos_id;
      Spiceint      themis_ir_id;

      /*
         load FK defining M01 name-ID mappings
      */
      furnsh_c ( "m01_v28.tf" );

      /*
         get NAIF ID for Phobos
      */
      bodn2c_c ( "PHOBOS", &phobos_id, &found );

      /*
         get NAIF ID for THEMIS IR camera
      */
      bodn2c_c ( "M01_THEMIS_IR", &themis_ir_id, &found );

Reference documents:

  • naif_ids.req - reference on object names and IDs supported in SPICE
  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Mapping objects to their preferred reference frames

APIs:

  • furnsh_c - loads individual kernels and collections of kernels.
  • cnmfrm_c - returns the preferred reference frame for an object specified by name.
  • cidfrm_c - returns the preferred reference frame for an object specified by ID.

Brief Example:

   The following example loads a text kernel that sets the preferred
   reference frame for the Earth to 'ITRF93' and calls cnmfrm_c to
   retrieve the name of this frame:

      SpiceInt                frcode;
      SpiceChar               frname[32];
      SpiceBoolean            found;

      /*
         load FK associating 'ITRF93' with the Earth
      */
      furnsh_c ( "earth_assoc_itrf93.tf" );

      /*
         retrieve the name of the preferred frame for the Earth
      */
      cnmfrm_c ( "EARTH", 32, &frcode, frname, &found );

Reference documents:

  • frames.req - reference on Frames Kernel (FK) and Frames subsystem of SPICE
  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Mapping between surface names and NAIF IDs

APIs:

  • furnsh_c - loads individual kernels and collections of kernels.
  • srfcss_c - translates a surface ID code, together with a body name, to the corresponding surface name.
  • srfs2c_c - translates a surface string, together with a body name, to the corresponding surface ID code.
  • srfc2s_c - translates a surface ID code, together with a body ID code, to the corresponding surface name.
  • srfscc_c - translates a surface string, together with a body ID code, to the corresponding surface ID code.

Brief Example:

   The following example uses srfs2c_c to get the surface ID for
   the 67P/Churyumov-Gerasimenko DSK surface named
   "ROS_CG_M001_OMSDLAM_N_V1":

      SpiceInt                surfid;
      SpiceBoolean            found;

      /*
         Load ROS FK defining surface name-ID mappings
      */
      furnsh_c ( "ROS_V25.TF" );

      /*
         Get NAIF ID for the surface "ROS_CG_M001_OMSDLAM_N_V1"
      */
      srfs2c_c ( "ROS_CG_M001_OMSDLAM_N_V1", "67P/C-G",
                 &surfid, &found );

Reference documents:

  • dsk.req - reference on Digital Shape Kernels (DSK) subsystem of SPICE
  • naif_ids.req - reference on object names and IDs supported in SPICE
  • kernel.req - reference on loading/unloading kernels and text kernel format

Top

Computing Planetocentric, Planetodetic, and Planetographic Coordinates

APIs:

  • reclat_c - converts from rectangular to planetocentric coordinates.
  • latrec_c - converts from planetocentric to rectangular coordinates.
  • srfrec_c - converts from planetocentric lat/lon of a surface point on a body to rectangular coordinates.
  • recgeo_c - converts from rectangular to geodetic coordinates.
  • georec_c - converts from geodetic to rectangular coordinates.
  • recpgr_c - converts from rectangular to planetographic coordinates.
  • pgrrec_c - converts from planetographic to rectangular coordinates.

Brief Example:

   The following example computes the planetocentric, planetodetic and
   planetographic coordinates for a point given as a Cartesian vector
   with respect to Mars in the Mars body-fixed frame (IAU_MARS):

      SpiceDouble   vec[3];
      SpiceDouble   radii[3];
      SpiceDouble   f;
      SpiceDouble   pcr,   pclon, pclat;
      SpiceDouble   pdalt, pdlon, pdlat;
      SpiceDouble   pgalt, pglon, pglat;
      SpiceInt      n;

      /*
         load PCK file
      */
      furnsh_c( "pck00008.tpc" );

      /*
         retrieve Mars radii
      */
      bodvrd_c( "MARS", "RADII", 3, &n, radii );

      /*
         package 3D vector
      */
      vpack_c ( 3373.850, -351.034, -117.267, vec );

      /*
         compute planetocentric coordinates
      */
      reclat_c( vec, &pcr, &pclon, &pclat );

      /*
         compute planetodetic coordinates
      */
      recgeo_c( vec, radii[0], f, &pdlon, &pdlat,
                &pdalt );

      /*
         compute planetographic coordinates
      */
      recpgr_c( "MARS", vec, radii[0], f, &pglon, &pglat,
                &pgalt );

Reference documents:

  • pck.req - reference on Planetary Constants Kernel (PCK) subsystem of SPICE

Top

Computing azimuth and elevation

APIs:

  • azlcpo_c - returns azimuth and elevation a target body relative to a constant-position observer location.
  • azlrec_c - converts from azimuth and elevation to rectangular coordinates.
  • recazl_c - converts from rectangular coordinates to azimuth and elevation.

Brief Example:

   The following example computes the azimuth/elevation coordinates
   of the Sun as seen from a location on Mars specified by its
   planetocentric coordinates:

      SpiceDouble   r;
      SpiceDouble   lon;
      SpiceDouble   lat;
      SpiceDouble   obspos[3];
      SpiceDouble   et;
      SpiceDouble   azlsta[6];
      SpiceDouble   ltim;

      /*
         load kernels: LSK, PCK, and satellite ephemeris SPK
      */
      furnsh_c( "naif0012.tls" );
      furnsh_c( "pck00010.tpc" );
      furnsh_c( "mar097s.bsp"  );

      /*
         convert planetocentric coordinates of the observer
         to rectangular coordinates.
      */
      r   = 3390.42;
      lon =  175.30;
      lat =  -14.59;

      latrec_c( r, lon * rpd_c(), lat * rpd_c(), obspos );

      /*
         convert UTC to ET
      */
      str2et_c ( "2022 JAN 31 01:00", &et );

      /*
         compute azimuth/elevation coordinates with the azimuth
         increasing clockwise about and the elevation positive
         towards +Z axis of the local topocentric reference frame.
      */
      azlcpo_c ( "ELLIPSOID", "SUN", et, "CN+S", SPICEFALSE, SPICETRUE,
                 obspos, "MARS", "IAU_MARS", azlsta, &ltim );

Top

Computing Surface Intercept Point

APIs:

  • sincpt_c - computes the surface intercept point of the ray on a body, modeled as an ellipsoid or a digital shape (DSK), at a specified epoch.
  • dskxv_c - computes ray-surface intercepts for a set of rays, using data provided by multiple loaded DSK segments.
  • dskxsi_c - computes a ray-surface intercept using data provided by multiple loaded DSK segments and returns information about the source of the data defining the surface on which the intercept was found.

Brief Example:

   The following example computes the surface intercept point of the
   MGS MOC NA camera boresight with the surface of Mars modeled as
   an ellipsoid:

      SpiceBoolean  found;
      SpiceChar     frame[32];
      SpiceChar     shape[32];
      SpiceDouble   et;
      SpiceDouble   bsight[3];
      SpiceDouble   bounds [4][3];
      SpiceDouble   spoint[3];
      SpiceDouble   trgepc;
      SpiceDouble   srfvec[3];
      SpiceInt      n;

      /*
         load kernels: LSK, PCK, MGS SCLK, planet/satellite
         ephemeris SPK, MGS spacecraft SPK, MGS spacecraft
         CK, MGS FK, and MOC IK
      */
      furnsh_c ( "naif0008.tls" );
      furnsh_c ( "mars_iau2000_v0.tpc" );
      furnsh_c ( "mgs_sclkscet_00061.tsc" );
      furnsh_c ( "mar063.bsp" );
      furnsh_c ( "mgs_ext22.bsp" );
      furnsh_c ( "mgs_sc_ext22.bc" );
      furnsh_c ( "mgs_v10.tf" );
      furnsh_c ( "mgs_moc_v20.ti" );

      /*
         retrieve MOC NA camera boresight direction
      */
      getfvn_c ( "MGS_MOC_NA", 4, 32, 32, shape, frame, bsight, &n,
                 bounds );

      /*
         convert UTC to ET
      */
      str2et_c ( "2006 JAN 31 01:00", &et );

      /*
         compute surface intercept point
      */
      sincpt_c ( "Ellipsoid", "MARS", et, "IAU_MARS",
                 "CN+S", "MGS", frame, bsight,
                 spoint, &trgepc, srfvec, &found );

.

Top

Computing point on a ray nearest to the body's surface

APIs:

  • tangpt_c - computes the "tangent point", the point on a ray nearest to the target's surface, modeled as an ellipsoid, at a specified epoch.

Brief Example:

   The following example computes the tangent point on the MAVEN
   IUVS instrument boresight nearest to the surface of Mars modeled
   as an ellipsoid:

      SpiceChar            frame  [32];
      SpiceChar            shape  [32];
      SpiceDouble          bsight [3];
      SpiceDouble          bounds [12][3];
      SpiceDouble          tanpt  [3];
      SpiceDouble          alt;
      SpiceDouble          range;
      SpiceDouble          srfpt  [3];
      SpiceDouble          trgepc;
      SpiceDouble          srfvec [3];
      SpiceInt             n;

      /*
         load kernels: LSK, PCK, MAVEN SCLK, planet/satellite
         ephemeris SPK, MAVEN spacecraft SPK, MAVEN FK, IUVS
         IK, and MAVEN spacecraft, IPP, and IUVS CKs.
      */
      furnsh_c ( "naif0012.tls"                       );
      furnsh_c ( "pck00010.tpc"                       );
      furnsh_c ( "mvn_sclkscet_00086.tsc"             );
      furnsh_c ( "mar097s.bsp"                        );
      furnsh_c ( "maven_orb_rec_201001_210101_v1.bsp" );
      furnsh_c ( "maven_v09.tf"                       );
      furnsh_c ( "maven_iuvs_v11.ti"                  );
      furnsh_c ( "mvn_sc_rel_201005_201011_v01.bc"    );
      furnsh_c ( "mvn_app_rel_201005_201011_v01.bc"   );
      furnsh_c ( "mvn_iuvs_rem_201001_201231_v01.bc"  );

      /*
         retrieve IUVS boresight direction
      */
      getfvn_c ( 'MAVEN_IUVS', 12, 32, 32, shape, frame, bsight, &n,
                 bounds );

      /*
         convert UTC to ET
      */
      str2et_c ( "2020-10-11 16:01:43", &et );

      /*
         compute tangent point
      */
      tangpt_c ( "ELLIPSOID", "MARS", et, "IAU_MARS", "CN+S",
                 "TANGENT POINT", "MAVEN", frame, bsight,
                 tanpt, &alt, &range, srfpt, &trgepc, srfvec );

Top

Computing Sub-observer and Sub-solar Points

APIs:

  • subpnt_c - computes the sub-observer point on a body, modeled as an ellipsoid or a digital shape (DSK), at a particular epoch.
  • subslr_c - computes the sub-solar point on a body, modeled as an ellipsoid or a digital shape (DSK), as seen by an observer at a particular epoch.

Brief Example:

   The following example computes the sub-spacecraft and sub-solar
   points on Mars, modeled as an ellipsoid, for MGS:

      SpiceDouble   et;
      SpiceDouble   subsc[3];
      SpiceDouble   subsolar[3];
      SpiceDouble   srfvec[3];
      SpiceDouble   trgepc;

      /*
         load kernels: LSK, PCK, planet/satellite SPK
         and MGS spacecraft SPK
      */
      furnsh_c ( "naif0008.tls" );
      furnsh_c ( "mars_iau2000_v0.tpc" );
      furnsh_c ( "mar063.bsp" );
      furnsh_c ( "mgs_ext22.bsp" );

      /*
         convert UTC to ET
      */
      str2et_c ( "2006 JAN 31 01:00", &et );

      /*
         compute sub-spacecraft point
      */
      subpnt_c ( "Near point: ellipsoid",
                 "MARS", et, "IAU_MARS", "LT+S", "MGS",
                 subsc, &trgepc, srfvec );

      /*
         compute sub-solar point
      */
      subslr_c ( "Near point: ellipsoid",
                 "MARS", et, "IAU_MARS", "LT+S", "MGS",
                 subsolar, &trgepc, srfvec );

Top

Computing Illumination Angles

APIs:

  • ilumin_c - computes the illumination angles at a specified surface point of a target body, modeled as an ellipsoid or a digital shape (DSK), as seen from an observer body, illuminated by the Sun.
  • illumg_c - computes the illumination angles at a specified surface point of a target body, modeled as an ellipsoid or a digital shape (DSK), as seen from an observer body, illuminated by a user specified body.
  • illumf_c - computes the illumination angles at a specified surface point of a target body, modeled as an ellipsoid or a digital shape (DSK), as seen from an observer body, illuminated by a user specified body, with flags indicating whether the point is visible from the observer and whether it is illuminated.
  • phaseq_c - computes the apparent phase angle between the centers of target, observer, and illuminator ephemeris objects.

Brief Example:

   The following example computes the illumination angles for a point
   specified using planetocentric coordinates, observed by MGS:

      SpiceDouble   r   = 3390.42;
      SpiceDouble   lon =  175.30;
      SpiceDouble   lat =  -14.59;
      SpiceDouble   point[3];
      SpiceDouble   et;
      SpiceDouble   srfvec[3];
      SpiceDouble   trgepc;
      SpiceDouble   phase, solar, emissn;

      /*
         load kernels: LSK, PCK, planet/satellite SPK
         and MGS spacecraft SPK
      */
      furnsh_c ( "naif0008.tls" );
      furnsh_c ( "mars_iau2000_v0.tpc" );
      furnsh_c ( "mar063.bsp" );
      furnsh_c ( "mgs_ext22.bsp" );

      /*
         convert planetocentric r/lon/lat to Cartesian vector
      */
      latrec_c( r, lon * rpd_c(), lat * rpd_c(), point );

      /*
         convert UTC to ET
      */
      str2et_c ( "2006 JAN 31 01:00", &et );

      /*
         compute illumination angles, modeling Mars as
         an ellipsoid
      */
      ilumin_c ( "Ellipsoid", "MARS",  et, "IAU_MARS",
                 "LT+S", "MGS", point,
                 &trgepc, srfvec, &phase, &solar, &emissn );

Top

Computing limb and terminator

APIs:

  • limbpt_c - computes limb points on a body, modeled as an ellipsoid or a digital shape (DSK).
  • termpt_c - computes umbral or penumbral terminator points on a body, modeled as an ellipsoid or a digital shape (DSK).

Brief Example:

   The following example computes the apparent limb points on
   Phobos, modeled by a DSK, as seen from Mars:

      #define MAXN         10000

      SpiceDouble          et;
      SpiceDouble          z       [3] = { 0.0, 0.0, 1.0 };
      SpiceInt             npts    [MAXN];
      SpiceDouble          points  [MAXN][3];
      SpiceDouble          trgeps  [MAXN];
      SpiceDouble          tangts  [MAXN][3];

      /*
         Load kernels: LSK, PCK, planet/satellite SPK
         and Phobos DSK
      */
      furnsh_c ( "naif0011.tls"  );
      furnsh_c ( "pck00010.tpc"  );
      furnsh_c ( "de430.bsp"     );
      furnsh_c ( "mar097.bsp"    );
      furnsh_c ( "phobos512.bds" );

      /*
         Convert UTC to ET
      */
      str2et_c ( "2008 aug 11 00:00:00", &et );

      /*
         Compute limb points in 360 cutting half-planes
         spaced 1 degree apart about the Mars-Phobos
         direction.
      */
      limbpt_c ( "TANGENT/DSK/UNPRIORITIZED", "PHOBOS", et,
               "IAU_PHOBOS", "CN+S", "CENTER", "MARS", z,
               1.0/dpr_c(), 360, 1.0e-4, 1.0e-7, MAXN,
               npts, points, trgeps, tangts );

Top

Computing angular separation between two bodies

APIs:

  • trgsep_c - Computes the angular separation between two spherical or point bodies at a specified epoch.

Brief Example:

   The following example computes the apparent angular separation
   of the Sun and Moon modeled as spheres as observed from the Earth:

      SpiceDouble          et;
      SpiceDouble          angle;

      /*
         Load kernels: LSK, PCK, and planetary ephemeris SPK
      */
      furnsh_c ( "naif0012.tls"  );
      furnsh_c ( "pck00010.tpc"  );
      furnsh_c ( "de421.bsp"     );

      /*
         Convert UTC to ET
      */
      str2et_c ( "2021-12-04T07:33:00", &et );

      /*
         Compute angular separation
      */
      angle = trgsep_c ( et,
                        "SUN",   "SPHERE", "IAU_SUN",
                        "MOON",  "SPHERE", "IAU_MOON",
                        "EARTH", "CN+S" );

Top

Computing and Propagating Orbital Elements

APIs:

  • conics_c - determines the state of an orbiting body from a set orbital elements.
  • oscelt_c - determines the set of orbital elements corresponding to the state of a body.
  • getelm_c - parses the "lines" of a two-line element set, returning the elements in units suitable for use in SPICE APIs.
  • evsgp4_c - evaluates NORAD two-line elements for earth orbiting spacecraft using algorithms by D. Vallado.

Brief Example:

   The following example computes the set of orbital elements for the
   state of the Mars Express spacecraft at a given time:

      SpiceInt      n;
      SpiceDouble   gm;
      SpiceDouble   et;
      SpiceDouble   state[6];
      SpiceDouble   lt;
      SpiceDouble   elts[8];

      /*
         load kernels: LSK, MEX trajectory SPK, and gravity PCK
      */
      furnsh_c ( "naif0008.tls" );
      furnsh_c ( "ORMM__050901000000_00165.BSP" );
      furnsh_c ( "DE403-MASSES.TPC" );

      /*
         retrieve GM for Mars
      */
      bodvrd_c ( "MARS", "GM", 1, &n, gm );

      /*
         convert UTC to ET
      */
      str2et_c ( "2005 SEP 02 04:50:45", &et );

      /*
         compute state of MEX at given UTC
      */
      spkezr_c ( "MEX", et, "MARSIAU", "NONE", "MARS",
                 state, &lt );

      /*
         compute orbital elements
      */
      oscelt_c ( state, et, gm, elts );

Top

Checking for in Field-Of-View (FOV) conditions

APIs:

  • fovray_c - determines if a specified ray is within the FOV of a specified instrument at a given time.
  • fovtrg_c - determines if a specified ephemeris object is within the FOV of a specified instrument at a given time.

Brief Example:

   The following example determines if Phoebe is in the CASSINI
   ISS NAC FOV on 2004-06-11 18:00 UTC:

      SpiceDouble   et;
      SpiceBoolean  visibl;

      /*
         Load kernels.
      */
      furnsh_c( "cas_2004_v17.tm" );

      /*
         Convert UTC to ET.
      */
      str2et_c( "2004-06-11 18:00", &et );

      /*
         Check for in-FOV condition.
      */
      fovtrg_c ( "CASSINI_ISS_NAC", "PHOEBE",
                 "ELLIPSOID", "IAU_PHOEBE", "CN+S",
                 "CASSINI", et, &visibl );

Top

Checking for occultation conditions

APIs:

  • occult_c - determines the occultation condition (not occulted, partially, etc.) of one target relative to another target as seen by an observer at a given time, with targets modeled as points, ellipsoids, or digital shapes (DSK).

Brief Example:

   The following example determines if the Sun is occulted by Saturn as
   seen from CASSINI on 2005-05-03 06:00 UTC, with both Sun and Saturn
   modeled as ellipsoids:

      SpiceDouble   et;
      SpiceInt      ocltid;

      /*
         Load kernels.
      */
      furnsh_c( "cas_2005_v17.tm" );

      /*
         Convert UTC to ET.
      */
      str2et_c ( "2005-05-03 06:00", &et );

      /*
         Check for in-FOV condition.
      */
      occult_c ( "SATURN", "ELLIPSOID", "IAU_SATURN",
                 "SUN",    "ELLIPSOID", "IAU_SUN",
                 "CN",     "CASSINI",   et,  &ocltid );

      if ( ocltid > 0 ) {
         printf ( "CASSINI IS IN SATURN SHADOW.\n" );
      }

Top

Finding times of events satisfying numerical constraints

APIs:

  • gfdist_c - determines time intervals when a specified constraint on observer-target distance is met.
  • gfilum_c - determines time intervals when a specified constraint on the observed phase, solar incidence, or emission angle at a surface point is met.
  • gfpa_c - determines time intervals when a specified constraint on the phase angle between the illuminator, target, and observer body centers is met.
  • gfposc_c - determines time intervals when a coordinate of an observer-target position vector satisfies a numerical constraint.
  • gfrr_c - determines time intervals when a specified constraint on the observer-target range rate is met.
  • gfsep_c - determines time intervals when the angular separation between the position vectors of two target bodies relative to an observer satisfies a numerical relationship.
  • gfsntc_c - determines time intervals when a coordinate of a ray-surface intercept position vector satisfies a numerical constraint.
  • gfsubc_c - determines time intervals when a coordinate of a sub-observer point position vector satisfies a numerical constraint.

Brief Example:

   The following example determines time intervals between Jan 1 and
   and April 1, 2007 when the distance between the Moon and the Earth
   was greater than 400,000 km.

      #define  MAXWIN  200

      SPICEDOUBLE_CELL      ( cnfine, MAXWIN );
      SPICEDOUBLE_CELL      ( result, MAXWIN );

      SpiceDouble             adjust;
      SpiceDouble             et0;
      SpiceDouble             et1;
      SpiceDouble             refval;
      SpiceDouble             step;

      /*
         Load kernels.
      */
      furnsh_c( "naif0008.tls" );
      furnsh_c( "de421.bsp"    );

      /*
         Store the time bounds in the confinement window.
      */
      str2et_c ( "2007 JAN 1", &et0 );
      str2et_c ( "2007 APR 1", &et1 );

      wninsd_c ( et0, et1, &cnfine );

      /*
         Set search parameters. Use a step size of 1 day (in units of
         seconds).
      */
      step   = spd_c();
      refval = 4.e5;
      adjust = 0.0;

      /*
         Perform search.
      */
      gfdist_c ( "MOON", "NONE", "EARTH", ">",     refval,
                 adjust, step,   MAXWIN,  &cnfine, &result );

Reference documents:

  • gf.req - reference on the Geometry Finder (GF) subsystem of SPICE

Top

Finding times of events satisfying FOV constraints

APIs:

  • gfrfov_c - determines time intervals when a specified ray intersects the space bounded by the field-of-view (FOV) of a specified instrument.
  • gftfov_c - determines time intervals when a specified ephemeris object intersects the space bounded by the field-of-view (FOV) of a specified instrument.

Brief Example:

   The following example determines time intervals between 09:00 and
   11:15 UTC on Jun 11, 2004 when Saturn's satellite Phoebe was within
   the FOV of the Cassini narrow angle camera (CASSINI_ISS_NAC).

      #define  MAXWIN         10000

      SPICEDOUBLE_CELL ( cnfine, MAXWIN );
      SPICEDOUBLE_CELL ( result, MAXWIN );

      SpiceDouble             et0;
      SpiceDouble             et1;
      SpiceDouble             stepsz;

      /*
         Load kernels.
      */
      furnsh_c( "naif0009.tls"               );
      furnsh_c( "cpck05Mar2004.tpc"          );
      furnsh_c( "981005_PLTEPH-DE405S.bsp"   );
      furnsh_c( "020514_SE_SAT105.bsp"       );
      furnsh_c( "030201AP_SK_SM546_T45.bsp"  );
      furnsh_c( "cas_v37.tf"                 );
      furnsh_c( "04135_04171pc_psiv2.bc"     );
      furnsh_c( "cas00084.tsc"               );
      furnsh_c( "cas_iss_v09.ti"             );

      /*
         Store the time bounds in the confinement window.
      */
      str2et_c ( "2004 JUN 11 09:00", &et0 );
      str2et_c ( "2004 JUN 11 11:15", &et1 );

      wninsd_c ( et0, et1, &cnfine );

      /*
         Set search parameters. Use a step size of 10 seconds.
      */
      stepsz = 10.0;

      /*
         Perform search.
      */
      gftfov_c ( "CASSINI_ISS_NAC", "PHOEBE", "ELLIPSOID",
                 "IAU_PHOEBE", "LT+S", "CASSINI",
                 stepsz, &cnfine, &result );

Reference documents:

  • gf.req - reference on the Geometry Finder (GF) subsystem of SPICE

Top

Finding times of occultations

APIs:

  • gfoclt_c - determines time intervals when an observer sees one target occulted by another, with targets modeled as points, ellipsoids, or digital shapes (DSK).

Brief Example:

   The following example determines time intervals in December, 2001
   when the Sun is occulted by the Moon (solar eclipses) as seen from
   the center of the Earth, with both Sun and Moon modeled as
   ellipsoids:

      #define MAXWIN          200

      SPICEDOUBLE_CELL      ( cnfine, MAXWIN );
      SPICEDOUBLE_CELL      ( result, MAXWIN );

      SpiceDouble             et0;
      SpiceDouble             et1;
      SpiceDouble             step;

      /*
         Load kernels.
      */
      furnsh_c( "naif0008.tls"  );
      furnsh_c( "de421.bsp"     );
      furnsh_c( "pck00008.tpc"  );

      /*
         Store the time bounds in the confinement window.
      */
      str2et_c ( "2001 DEC 01", &et0 );
      str2et_c ( "2002 JAN 01", &et1 );

      wninsd_c ( et0, et1, &cnfine );

      /*
         Set search parameters. Select a 3-minute step.
      */
      step = 180.0;

      /*
         Perform search.
      */
      gfoclt_c ( "ANY",
                 "MOON",    "ELLIPSOID",  "IAU_MOON",
                 "SUN",     "ELLIPSOID",  "IAU_SUN",
                 "LT",      "EARTH",      step,
                 &cnfine,   &result                 );

Reference documents:

  • gf.req - reference on the Geometry Finder (GF) subsystem of SPICE

Top

Converting from and to rectangular coordinates

APIs:

  • reclat_c - converts from rectangular to latitudinal coordinates.
  • latrec_c - converts from latitudinal to rectangular coordinates.
  • srfrec_c - converts from planetocentric lat/lon of a surface point on a body to rectangular coordinates.
  • recgeo_c - converts from rectangular to geodetic coordinates.
  • georec_c - converts from geodetic to rectangular coordinates.
  • recpgr_c - converts from rectangular to planetographic coordinates.
  • pgrrec_c - converts from planetographic to rectangular coordinates.
  • recrad_c - converts from rectangular coordinates to range, right ascension, and declination.
  • radrec_c - converts from range, right ascension, and declination to rectangular coordinates.
  • recsph_c - converts from rectangular to spherical coordinates.
  • sphrec_c - converts from spherical to rectangular coordinates.
  • reccyl_c - converts from rectangular to cylindrical coordinates.
  • cylrec_c - converts from cylindrical to rectangular coordinates.

Brief Example:

   The following example computes the planetocentric, planetodetic and
   planetographic coordinates for a point given as a Cartesian vector
   with respect to Mars in the Mars body-fixed frame (IAU_MARS):

      SpiceDouble   vec[3];
      SpiceDouble   radii[3];
      SpiceDouble   f;
      SpiceDouble   pcr,   pclon, pclat;
      SpiceDouble   pdalt, pdlon, pdlat;
      SpiceDouble   pgalt, pglon, pglat;
      SpiceInt      n;

      /*
         load PCK file
      */
      furnsh_c( "pck00008.tpc" );

      /*
         retrieve Mars radii
      */
      bodvrd_c( "MARS", "RADII", 3, &n, radii );

      /*
         pack 3D vector
      */
      vpack_c ( 3373.850, -351.034, -117.267, vec );

      /*
         compute planetocentric coordinates
      */
      reclat_c( vec, &pcr, &pclon, &pclat );

      /*
         compute planetodetic coordinates
      */
      recgeo_c( vec, radii[0], f, &pdlon, &pdlat,
                &pdalt );

      /*
         compute planetographic coordinates
      */
      recpgr_c( "MARS", vec, radii[0], f, &pglon, &pglat,
                &pgalt );

Top

Converting from and to Spherical Coordinates

APIs:

  • sphrec_c - converts from spherical to rectangular coordinates.
  • recsph_c - converts from rectangular to spherical coordinates.
  • sphcyl_c - converts from spherical to cylindrical coordinates.
  • cylsph_c - converts from cylindrical to spherical coordinates.
  • sphlat_c - converts from spherical to latitudinal coordinates.
  • latsph_c - converts from latitudinal to spherical coordinates.

Brief Example:

   The following example computes the spherical coordinates of a point
   specified as a Cartesian vector:

      SpiceDouble   vec[3];
      SpiceDouble   r;
      SpiceDouble   colat;
      SpiceDouble   lon;

      vpack_c  ( 1.0, 1.0, 1.0, vec );

      recsph_c ( vec, &r, &colat, &lon );

Top

Converting from and to cylindrical coordinates

APIs:

  • cylrec_c - converts from cylindrical to rectangular coordinates.
  • reccyl_c - converts from rectangular to cylindrical coordinates.
  • cylsph_c - converts from cylindrical to spherical coordinates.
  • sphcyl_c - converts from spherical to cylindrical coordinates.
  • cyllat_c - converts from cylindrical to latitudinal coordinates.
  • latcyl_c - converts from latitudinal to cylindrical coordinates.

Brief Example:

   The following example computes the cylindrical coordinates of a
   point specified as a Cartesian vector:

      SpiceDouble   vec[3];
      SpiceDouble   r;
      SpiceDouble   lon;
      SpiceDouble   z;

      vpack_c  ( 1.0, 1.0, 1.0, vec );

      reccyl_c ( vec, &r, &lon, &z );

Top

Converting from and to Latitudinal Coordinates

APIs:

  • latrec_c - converts from latitudinal to rectangular coordinates.
  • reclat_c - converts from rectangular to latitudinal coordinates.
  • latsph_c - converts from latitudinal to spherical coordinates.
  • sphlat_c - converts from spherical to latitudinal coordinates.
  • latcyl_c - converts from latitudinal to cylindrical coordinates.
  • cyllat_c - converts from cylindrical to latitudinal coordinates.

Brief Example:

   The following example converts the planetocentric coordinates of a
   point to a Cartesian vector:

      SpiceDouble   vec[3];
      SpiceDouble   r   = 3390.42;
      SpiceDouble   lon =  175.30;
      SpiceDouble   lat =  -14.59;

      latrec_c( r, lon * rpd_c(), lat * rpd_c(), vec );

Top

Converting from and to R, RA, and DEC

APIs:

  • radrec_c - converts from range, right ascension, and declination to rectangular coordinates.
  • recrad_c - converts from rectangular coordinates to range, right ascension, and declination.

Brief Example:

   The following example computes the unit vector corresponding to an
   RA and DEC given in degrees:

      SpiceDouble   ra  = 120.0;
      SpiceDouble   dec = -30.0;
      SpiceDouble   vec[3];

      radrec_c( 1.0, ra * rpd_c(), dec * rpd_c(), vec );

Top

Converting from and to Geodetic Coordinates

APIs:

  • georec_c - converts from geodetic to rectangular coordinates.
  • recgeo_c - converts from rectangular to geodetic coordinates.

Brief Example:

   The following example converts the areodetic coordinates of a
   landing site (Lon = 175.3 deg, Lat = -14.59 deg, Alt = -1.91 km),
   given relative to the IAU 2000 Mars ellipsoid, to a Cartesian
   vector:

      SpiceDouble   lon;
      SpiceDouble   lat;
      SpiceDouble   alt;
      SpiceDouble   radii[3];
      SpiceDouble   f;
      SpiceDouble   vec[3];

      SpiceInt      n;

      furnsh_c( "pck00008.tpc" );

      bodvrd_c( "MARS", "RADII", 3, &n, radii );

      lon =  175.30 * rpd_c();
      lat =  -14.59 * rpd_c();
      alt =   -1.91;

      f = (radii[0]-radii[2])/radii[0];

      georec_c( lon, lat, alt, radii[0], f, vec );

Top

Transforming states from one coordinate system to another

APIs:

  • xfmsta_c - transforms states between coordinate systems -- rectangular, cylindrical, latitudinal, spherical, geodetic, and planetographic.

Brief Example:

   The following example converts a state specified in the
   rectangular coordinate system (X, Y, Z, DX, DY, DZ) to
   the state specified in the latitudinal coordinate system
   (R, LONG, LAT, DR, DLONG, DLAT):

      SpiceDouble   istate[6];
      SpiceDouble   ostate[6];

      istate[0] = 1000.0;
      istate[1] = 2000.0;
      istate[2] = 3000.0;
      istate[3] =   20.0;
      istate[4] =  -40.0;
      istate[6] =  -60.0;

      xfmsta_c ( istate, "RECTANGULAR",
                 "LATITUDINAL", " ",
                 ostate );

      printf( "R     = %20.10f KM\n",    ostate[0] );
      printf( "LONG  = %20.10f RAD\n",   ostate[1] );
      printf( "LAT   = %20.10f RAD\n",   ostate[2] );
      printf( "DR    = %20.10f KM/S\n",  ostate[3] );
      printf( "DLONG = %20.10f RAD/S\n", ostate[4] );
      printf( "DLAT  = %20.10f RAD/S\n", ostate[5] );

Top

Performing simple operations on 3D vectors

APIs:

  • vpack_c - packs three scalar components into a vector.
  • vupack_c - unpacks three scalar components from a vector.
  • vadd_c - adds two 3D vectors.
  • vsub_c - computes the difference between two 3D vectors.
  • vcrss_c - computes the cross product of two 3D vectors.
  • vdot_c - computes the dot product of two 3D vectors.
  • vrel_c - returns the relative difference between two 3D vectors
  • vscl_c - multiplies a scalar and a 3D vector.
  • vminus_c - negates a 3D vector.
  • vequ_c - makes one 3D vector equal to another.
  • vzero_c - indicates whether a 3D vector is the zero vector.
  • vsep_c - finds the separation angle between two 3D vectors.
  • vdist_c - returns the distance between two 3D vectors.
  • vnorm_c - computes the magnitude of a 3D vector.
  • vhat_c - finds the unit vector along a 3D vector.
  • ucrss_c - computes the normalized cross product of two 3D vectors.
  • unorm_c - normalizes a 3D vector and return its magnitude.

Brief Example:

   The following example prints the angular separation between two 3D
   vectors in degrees:

      SpiceDouble   a[3];
      SpiceDouble   b[3];

      vpack_c( 1.0, 0.0, 0.0, a );
      vpack_c( 1.0, 1.0, 0.0, b );

      printf( "angular separation, deg = %12.3f\n",
              vsep_c(a,b) * dpr_c() );

Top

Projecting, Combining and Rotating 3D vectors

APIs:

  • rotvec_c - transform a 3D vector to a new coordinate system rotated by an angle about X, Y, or Z.
  • vperp_c - finds the component of a 3D vector that is perpendicular to a second 3D vector.
  • vproj_c - finds the projection of one 3D vector onto another 3D vector.
  • vrotv_c - rotates a 3D vector about a specified axis 3D vector by a specified angle.
  • vlcom_c - computes the vector linear combination a*v1 + b*v2 of two 3D vectors.
  • vlcom3_c - computes the vector linear combination a*v1 + b*v2 + c*v3 of three 3D vectors.

Brief Example:

   The following example finds the projection of one 3D vector onto
   another 3D vector:

      SpiceDouble   a[3];
      SpiceDouble   b[3];
      SpiceDouble   c[3];

      vproj_c( a, b, c );

Top

Performing Simple Operations on 3x3 Matrices

APIs:

  • mxm_c - multiplies two 3x3 matrices.
  • mxmt_c - multiplies a 3x3 matrix and the transpose of another 3x3 matrix.
  • mxv_c - multiplies a 3x3 matrix with a 3D vector.
  • mtxm_c - multiplies the transpose of a 3x3 matrix and a 3x3 matrix.
  • mtxv_c - multiplies the transpose of a 3x3 matrix on the left with a 3D vector on the right.
  • vtmv_c - multiplies the transpose of a 3D vector, a 3x3 matrix, and a 3D vector.
  • xpose_c - transposes a 3x3 matrix.
  • mequ_c - sets one 3x3 matrix equal to another.
  • det_c - computes the determinant of a 3x3 matrix.
  • trace_c - returns the trace of a 3x3 matrix.

Brief Example:

   The following example multiples two matrices and converts the
   resulting matrix to a SPICE-style quaternion:

      SpiceDouble   mat1[3][3];
      SpiceDouble   mat2[3][3];
      SpiceDouble   mat3[3][3];
      SpiceDouble   q[4];

      mxm_c( mat1, mat2, mat3 );
      m2q_c( mat3, q );

Top

Creating and Converting Transformation Matrices

APIs:

  • rotate_c - calculates the 3x3 matrix for a rotation of an angle about the X, Y or Z axis.
  • rotmat_c - applies a rotation of an angle about the X, Y, or Z axis to a matrix.
  • twovec_c - builds the transformation to a frame based on two vectors.
  • eul2m_c - constructs a rotation matrix from a set of Euler angles.
  • m2eul_c - factors a matrix as a product of three rotations about specified axes.
  • raxisa_c - computes the axis of the rotation given by a matrix and the angle about that axis.
  • axisar_c - construct a rotation matrix that rotates vectors by an angle about an axis.
  • m2q_c - finds a unit quaternion corresponding to a specified rotation matrix.
  • q2m_c - find the rotation matrix corresponding to a specified unit quaternion.

Brief Example:

   The following example creates a transformation matrix from right
   ascension, declination, and twist given in degrees:

      SpiceDouble   ra    =   10.0;
      SpiceDouble   dec   =   20.0;
      SpiceDouble   twist = -135.0;
      SpiceDouble   mat[3][3];

      ra    = ra    * rpd_c();
      dec   = dec   * rpd_c();
      twist = twist * rpd_c();

      eul2m_c( twist, halfpi_c() - dec, ra, 3, 2, 3, mat );

Reference documents:

  • rotation.req - reference on different representations of rotation transformations

Top

Accessing Physical and Mathematical Constants

APIs:

  • halfpi_c - returns half the value of pi.
  • pi_c - returns the value of pi.
  • twopi_c - returns twice the value of pi.
  • dpr_c - returns the number of degrees per radian.
  • rpd_c - returns the number of radians per degree.
  • spd_c - returns the number of seconds in a day.
  • jyear_c - returns the number of seconds per Julian year.
  • tyear_c - returns the number of seconds per tropical year.
  • clight_c - returns the speed of light in vacuo (km/sec)
  • b1900_c - returns the Julian Date corresponding to Besselian date 1900.0.
  • b1950_c - returns the Julian Date corresponding to Besselian date 1950.0.
  • j1900_c - returns the Julian Date of 1899 DEC 31 12:00:00 (1900 JAN 0.5).
  • j1950_c - returns the Julian Date of 1950 JAN 01 00:00:00 (1950 JAN 1.0).
  • j2000_c - returns the Julian Date of 2000 JAN 01 12:00:00 (2000 JAN 1.5).
  • j2100_c - returns the Julian Date of 2100 JAN 01 12:00:00 (2100 JAN 1.5).

Brief Example:

   The following example sets angle to pi/2 and prints it in degrees:

      SpiceDouble   angle;

      angle = halfpi_c();

      printf( "angle, deg = %12.3f\n", angle * dpr_c() );