 
In-situ Sensing Hands-On Lesson (C)
===========================================================================
 
   May 11, 2009
 
 
Overview
--------------------------------------------------------
 
   In this lesson you will develop a simple program illustrating how SPICE
   can be used to compute various kinds of geometry information applicable
   to the experiments carried out by an in-situ instrument flown on an
   interplanetary spacecraft.
 
 
Note About HTML Links
--------------------------------------------------------
 
   The HTML version of this lesson contains links pointing to various HTML
   documents provided with the Toolkit. All of these links are relative
   and, in order to function, require this document to be in a certain
   location in the Toolkit HTML documentation directory tree.
 
   In order for the links to be resolved, create a subdirectory called
   ``lessons'' under the ``doc/html'' directory of the Toolkit tree and
   copy this document to that subdirectory before loading it into a Web
   browser.
 
 
References
--------------------------------------------------------
 
   This section lists SPICE documents referred to in this lesson.
 
   Of these documents, the ``Tutorials'' contains the highest level
   descriptions with the least number of details while the ``Required
   Reading'' documents contain much more detailed specifications. The most
   complete specifications are provided in the ``Headers'' -- the comments
   in the top section of the source file.
 
   In some cases the lesson explanations also refer to the information
   provided in the meta-data area of the kernels used in the lesson
   examples. It is especially true in case of the FK and IK files, which
   often contain comprehensive descriptions of the frames, instrument FOVs,
   etc. Since both FK and IK are text kernels, the information provided in
   them can be viewed using any text editor, while the meta information
   provided in binary kernels -- SPKs and CKs -- can be viewed using
   ``commnt'' or ``spacit'' utility programs located in ``cspice/exe'' of
   Toolkit installation tree.
 
 
Tutorials
 
   The following SPICE tutorials are referenced in this lesson:
 
      Name             Lesson steps/routines that it describes
      ---------------  -----------------------------------------
      Time             UTC to ET and SCLK to ET
      Loading Kernels  Loading SPICE kernels
      SCLK             SCLK to ET time conversion
      SPK              Computing positions and velocities
      Frames           Computing transformations between frames
 
   These tutorials are available in printed form and as MS Office or PDF
   files from NAIF server at JPL:
 
      http://naif.jpl.nasa.gov/naif/tutorials.html
 
 
Required Reading Documents
 
   The Toolkit includes a set of Required Reading documents. You find these
   documents in the documentation directory, ``cspice/doc'', directory of
   the Toolkit installation trees.
 
      Name             Lesson steps/routines that it describes
      ---------------  -----------------------------------------
      time.req         UTC to ET time conversion
      kernel.req       Loading SPICE kernels
      sclk.req         SCLK to ET time conversion
      naif_ids.req     Body and reference frame names
      spk.req          Computing positions and velocities
 
   Another very useful document, also distributed with the Toolkit, is
   ``Permuted Index'', called ``spicelib.idx'' for FORTRAN or
   ``cspice.idx'' for C, IDL, and MATLAB also found in the ``doc''
   directory.
 
   This text document provides an easy way to find what SPICE routine(s)
   performs a particular function of interest and the name of the source
   file that contains this function (this is especially useful for FORTRAN
   because some of the routines are entry points and, therefore, their name
   is different from the name of the source file in which they are
   located.)
 
 
Source Code Header Comments
 
   The most detailed specification of a given SPICE FORTRAN or C routine is
   contained in the header section of its source code. The source code is
   distributed with the Toolkit and is located under
   ``toolkit/src/spicelib'' in FORTRAN and under ``cspice/src/cspice'' in C
   Toolkits.
 
   For example the source code of the STR2ET/str2et_c routine is
 
      toolkit/src/spicelib/str2et.for
 
   in the FORTRAN Toolkit and in
 
      cspice/src/cspice/str2et_c.c
 
   in the C Toolkit.
 
   Since some of the FORTRAN routines are entry points they are usually
   part of a source file that has different name. The ``Permuted Index''
   document mentioned above can be used to locate the name of their source
   file.
 
 
Kernels Used
--------------------------------------------------------
 
   The kernels used in this lessons:
 
      File Name                 Type Description
      ------------------------- ---- --------------------------
      naif0008.tls              LSK  Generic LSK
      cpck05Mar2004.tpc         PCK  Cassini project PCK
      cas00084.tsc              SCLK Cassini SCLK
      020514_SE_SAT105.bsp      SPK  Saturnian Satellite Ephemeris SPK
      030201AP_SK_SM546_T45.bsp SPK  Cassini Spacecraft SPK
      981005_PLTEPH-DE405S.bsp  SPK  Planetary Ephemeris SPK
      sat128.bsp                SPK  Saturnian Satellite Ephemeris SPK
      04135_04171pc_psiv2.bc    CK   Cassini Spacecraft CK
      cas_v37.tf                FK   Cassini FK
 
   These SPICE kernels are included in the lesson package available from
   the NAIF server at JPL:
 
      ftp://naif.jpl.nasa.gov/pub/naif/toolkit_docs/Lessons/
 
 
CSPICE Routines Used
--------------------------------------------------------
 
   The CSPICE routines demonstrated in the lesson:
 
      Name        Function that it performs
      ----------  ---------------------------------------------------
      furnsh_c    Loads kernels, individually or listed in meta-kernel
      str2et_c    Converts UTC to ET
      scs2e_c     Converts SCLK to ET
      spkezr_c    Computes states (position & velocity)
      spkpos_c    Computes positions
      vhat_c      Find unit vector along a 3d vector
      subpnt_c    Computes body-fixed coordinates of sub-observer point
      reclat_c    Converts rectangular coordinated to latitudinal
      vsub_c      Subtracts 3d vectors
      pxform_c    Computes 3x3 matrix rotating vectors between frames
      mxv_c       Multiplies 3d vector by 3x3 matrix
      vpack_c     Packs 3 number into a 3d vector
 
   The most detailed documentation source for these routines are their
   headers.
 
 
Step-1: ``UTC to ET''
===========================================================================
 
 
``UTC to ET'' Task Statement
--------------------------------------------------------
 
   Write a program that computes and prints the Ephemeris Time (ET)
   corresponding to ``2004-06-11T19:32:00'' UTC, as the number of ephemeris
   seconds past J2000, .
 
 
``UTC to ET'' Hints
--------------------------------------------------------
 
   Find out what SPICE kernel(s) is(are) needed to support this conversion.
   Reference the ``time.req'' and/or ``Time'' tutorial.
 
   Find necessary kernel(s) on the NAIF's FTP site.
 
   Find out what routine should be called to load necessary kernel(s).
   Reference the ``kernel.req'' and/or ``Loading Kernels'' tutorial.
 
   Find the ``loader'' routine calling sequence specification. Look at the
   ``time.req'' and that routine's source code header. This routine may be
   an entry point, in which case there will be no source file with the same
   name. To find out in which source file this entry point is, search for
   its name in the ``Permuted Index''.
 
   Find the routine(s) used to convert time between UTC and ET. Look at the
   ``time.req'' and/or ``Time'' tutorial.
 
   Find the ``converter'' routine(s) calling sequence specification. Look
   in the ``time.req'' and the routine's source code header.
 
   Put all calls together in a program, add variable declarations (the
   routine header's ``Declarations'' and ``Examples'' sections are a good
   place to look for declaration specification and examples) and output
   print statements. Compile it and link it against CSPICE.
 
 
``UTC to ET'' Solution Steps
--------------------------------------------------------
 
   Only one kernel file is needed to support this conversion -- an LSK file
   ``naif0008.tls''.
 
   As any other SPICE kernel this file can be loaded by the furnsh_c
   routine. For that, the name of the file can put be provided as a sole
   argument of this routine:
 
      #include "SpiceUsr.h"
      ...
      SpiceChar         * lskfle = "naif0008.tls";
 
      furnsh_c ( lskfle );
 
   or it can be listed in a meta-kernel:
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           )
      \begintext
 
   the name of which, let's call it ``spice_example.tm'', can be then
   provided as a sole argument of the furnsh_c routine:
 
      #include "SpiceUsr.h"
      ...
      SpiceChar         * mkfile = "spice_example.tm";
 
      furnsh_c ( mkfile );
 
   While the second option seems to involve a bit more work -- it requires
   making an extra file -- it is a much better way to go if you plan to
   load more kernels as you extend the program. With the meta-kernel
   approach simply adding more kernels to the list in KERNEL_TO_LOAD
   without changing the program code will accomplish that.
 
   The highest level CSPICE time routine converting UTC to ET is str2et_c
   (``cspice/src/cspice/str2et_c.c'').
 
   It has two arguments -- input time string representing UTC in a variety
   of formats (see str2et_c header's section ``Particulars'' for the
   complete description of input time formats) and output DP number of ET
   seconds past J2000. A call to str2et_c converting a given UTC to ET
   could look like this:
 
      SpiceChar         * utc = "2004-06-11T19:32:00";
      SpiceDouble         et;
      ...
      str2et_c ( utc, &et );
 
   By combining furnsh_c and str2et_c calls and required declarations and
   by adding a simple print statement, one would get a complete program
   that prints ET for the given UTC epoch.
 
   The program's source code then needs to be compiled and linked against
   CSPICE. Assuming that the program was saved in a file called
 
   "spice_example.c", this can be done with the following command on a Sun
   workstation:
 
   Sun C:
 
      cc -c -Xc -o spice_example spice_example.c cspice.a -lm
 
   gcc:
 
      gcc -c -ansi -Wall -o spice_example spice_example.c cspice.a -lm
 
   The command assumes that the library file(s) ``cspice.a'' and the CSPICE
   include files *.h are located in current directory, which may not be the
   case.
 
   When you run the program's executable, ``spice_example'', it produces
   the following output (the output below was generated by this program
   compiled with gcc on a PC running Linux; your output may differ slightly
   in its format and numeric precision):
 
      > ./spice_example
      utc  = 2004-06-11T19:32:00
      et   =     140254384.184625
      >
 
 
``UTC to ET'' Code
--------------------------------------------------------
 
   Program ``spice_example.c'':
 
         #include <stdio.h>
         #include "SpiceUsr.h"
 
         int main()
      {
         SpiceChar             * mkfile;
         SpiceChar             * utc;
         SpiceDouble             et;
 
 
         mkfile = "spice_example.tm";
         furnsh_c ( mkfile   );
 
         utc    = "2004-06-11T19:32:00";
         str2et_c ( utc, &et );
 
         printf ( "utc  = %s     \n", utc );
         printf ( "et   = %20.6f \n", et  );
 
         return ( 0 );
      }
 
   Meta-kernel file ``spice_example.tm'':
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           )
      \begintext
 
 
Step-2: ``SCLK to ET''
===========================================================================
 
 
``SCLK to ET'' Task Statement
--------------------------------------------------------
 
   Extend the program from Step-1 to compute and print ET for the following
   CASSINI on-board clock epoch ``1465674964.105''.
 
 
``SCLK to ET'' Hints
--------------------------------------------------------
 
   Find out what additional (to those already loaded in Step-1) SPICE
   kernel(s) is(are) needed to support SCLK to ET conversion. Look at the
   ``sclk.req'' and/or ``SCLK'' tutorial.
 
   Find necessary kernel(s) on the NAIF's FTP site.
 
   Modify the program or meta-kernel to load this(these) kernels.
 
   Find the routine(s) needed to convert time between SCLK and ET. Look at
   the ``sclk.req'' and/or ``Time'' and ``SCLK'' tutorials.
 
   Find the ``converter'' routine's calling sequence specification. Look in
   the ``sclk.req'' and the routine's source code header.
 
   Look at ``naif_ids.req'' and the comments in the additional kernel(s)
   that you have loaded for information on proper values of input arguments
   of this routine.
 
   Add calls to the ``converter'' routine(s), necessary variable
   declarations (the routine header's ``Declarations'' and ``Examples''
   sections are a good place to look for declaration specification and
   examples), and output print statements to the program. Re-compile and
   re-link it against CSPICE.
 
 
``SCLK to ET'' Solution Steps
--------------------------------------------------------
 
   A CASSINI SCLK file is needed additionally to the LSK file loaded in the
   Step-1 to support this conversion.
 
   No code change is needed in the loading portion of the program if a
   meta-kernel approach was used in the Step-1. The program will load the
   file if it will be added to the list of kernels in the KERNELS_TO_LOAD
   variable:
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           'kernels/sclk/cas00084.tsc'
                           )
      \begintext
 
   The highest level CSPICE routine converting SCLK to ET is scs2e_c
   (``cspice/src/cspice/scs2e_c.c'').
 
   It has three arguments -- NAIF ID for CASSINI s/c (-82 as described by
   ``naif_ids.req'' document), input time string representing CASSINI SCLK,
   and output DP number of ET seconds past J2000. A call to str2et_c
   converting given SCLK to ET could look like this:
 
      SpiceChar          * sclk = "1465674964.105";
      SpiceInt             scid = -82;
      ...
      scs2e_c ( scid, sclk, &et );
 
   By adding the scs2e_c call, required declarations and a simple print
   statement, one would get a complete program that prints ET for the given
   SCLK epoch.
 
   The program's source code then needs to be re-compiled and re-linked
   against CSPICE. It can be done using the same compile command as in
   Step-1:
 
      Sun C:
 
      cc -c -Xc -o spice_example spice_example.c cspice.a -lm
 
      gcc:
 
      gcc -c -ansi -Wall -o spice_example spice_example.c cspice.a -lm
 
   When you run the program's executable, ``spice_example'', it produces
   the following output (the output below was generated by this program
   compiled with gcc on a PC running Linux; your output may differ slightly
   in its format and numeric precision):
 
      > ./spice_example
      utc  = 2004-06-11T19:32:00
      et   =     140254384.184625
      sclk = 1465674964.105
      et   =     140254384.183426
      >
 
 
``SCLK to ET'' Code
--------------------------------------------------------
 
   Program ``spice_example.c'':
 
         #include <stdio.h>
         #include "SpiceUsr.h"
 
         int main()
      {
         SpiceChar             * mkfile;
         SpiceChar             * utc;
         SpiceChar             * sclk;
 
         SpiceDouble             et;
 
         SpiceInt                scid;
 
 
         mkfile  =  "spice_example.tm";
         furnsh_c ( mkfile   );
 
         utc     =  "2004-06-11T19:32:00";
         str2et_c ( utc, &et );
 
         printf ( "utc  = %s     \n", utc );
         printf ( "et   = %20.6f \n", et  );
 
         scid   = -82;
         sclk   = "1465674964.105";
         scs2e_c ( scid, sclk, &et );
 
         printf ( "sclk = %s     \n", sclk );
         printf ( "et   = %20.6f \n", et  );
 
 
         return ( 0 );
      }
 
   Meta-kernel file ``spice_example.tm'':
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           'kernels/sclk/cas00084.tsc'
                           )
      \begintext
 
 
Step-3: ``Spacecraft State''
===========================================================================
 
 
``Spacecraft State'' Task Statement
--------------------------------------------------------
 
   Extend the program from Step-2 to compute geometric state -- position
   and velocity -- of the CASSINI spacecraft with respect to the Sun in the
   Ecliptic frame at the epoch specified by SCLK time from Step-2.
 
 
``Spacecraft State'' Hints
--------------------------------------------------------
 
   Find out what additional (to those already loaded in Steps-1&2) SPICE
   kernel(s) is(are) needed to support state computation. Look at the
   ``spk.req'' and/or ``SPK'' tutorial.
 
   Find necessary kernel(s) on the NAIF's FTP site.
 
   Verify that the kernels contain enough data to compute the state of
   interest. Use ``brief'' utility program located under ``toolkit/exe''
   directory for that.
 
   Modify the meta-kernel to load this(these) kernels.
 
   Determine the routine(s) needed to compute states. Look at the
   ``spk.req'' and/or ``SPK'' tutorial presentation.
 
   Find the the routine(s) calling sequence specification. Look in the
   ``spk.req'' and the routine's source code header.
 
   Reference the ``naif_ids.req'' and ``frames.req'' and the routine(s)
   header ``Inputs'' and ``Particulars'' sections to determine proper
   values of the input arguments of this routine.
 
   Add calls to the routine(s), necessary variable declarations and output
   print statements to the program. Re-compile and re-link it against
   CSPICE.
 
 
``Spacecraft State'' Solution Steps
--------------------------------------------------------
 
   A CASSINI spacecraft trajectory SPK and generic planetary ephemeris SPK
   files are needed to support computation of the state of interest.
 
   The file names can be added to the meta-kernel to get them loaded into
   the program:
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           'kernels/sclk/cas00084.tsc'
                           'kernels/spk/020514_SE_SAT105.bsp'
                           'kernels/spk/030201AP_SK_SM546_T45.bsp'
                           'kernels/spk/981005_PLTEPH-DE405S.bsp'
                           'kernels/spk/sat128.bsp'
                           )
      \begintext
 
   The highest level CSPICE routine computing states is spkezr_c
   (``cspice/src/cspice/spkezr_c.c'').
 
   We are interested in computing CASSINI position and velocity with
   respect to the Sun, therefore the target and observer names should be
   set to 'CASSINI' and 'Sun' (both names can be found in
   ``naif_ids.req'').
 
   The state should be in ecliptic frame, therefore the name of the frame
   in which the state should be computed is 'ECLIPJ2000' (see
   ``frames.req'' document.)
 
   Since we need only the geometric position, the `abcorr' argument of the
   routine should be set to "NONE" (see aberration correction discussion in
   the (``cspice/src/cspice/spkezr_c.c''). header).
 
   Putting it all together, we get:
 
      SpiceChar           * target;
      SpiceChar           * frame;
      SpiceChar           * corrtn;
      SpiceChar           * observ;
      SpiceDouble           state [6];
      SpiceDouble           ltime;
      ...
      target = "CASSINI";
      frame  = "ECLIPJ2000";
      corrtn = "NONE";
      observ = "SUN";
 
      spkezr_c ( target, et, frame, corrtn, observ, state, &ltime );
 
   The updated program with added calls, required declarations and simple
   print statements produces the following output (the output below was
   generated by this program compiled with gcc on a PC running Linux; your
   output may differ slightly in its format and numeric precision):
 
      > ./spice_example
      utc  = 2004-06-11T19:32:00
      et   =     140254384.184625
      sclk = 1465674964.105
      et   =     140254384.183426
      state = -3.765991e+08 1.294488e+09 -7.064853e+06 -5.164226e+00
      8.017189e-01 4.060306e-02
      >
 
 
``Spacecraft State'' Code
--------------------------------------------------------
 
   Program ``spice_example.c'':
 
 
         #include <stdio.h>
         #include "SpiceUsr.h"
 
         int main()
      {
         SpiceChar             * mkfile;
         SpiceChar             * utc;
         SpiceChar             * sclk;
         SpiceChar             * target;
         SpiceChar             * frame;
         SpiceChar             * corrtn;
         SpiceChar             * observ;
 
         SpiceDouble             et;
         SpiceDouble             state [6];
         SpiceDouble             ltime;
 
         SpiceInt                scid;
 
 
         mkfile  =  "spice_example.tm";
         furnsh_c ( mkfile );
 
         utc     =  "2004-06-11T19:32:00";
         str2et_c ( utc, &et );
 
         printf ( "utc  = %s     \n", utc );
         printf ( "et   = %20.6f \n", et  );
 
         scid   = -82;
         sclk   = "1465674964.105";
         scs2e_c ( scid, sclk, &et );
 
         printf ( "sclk = %s     \n", sclk );
         printf ( "et   = %20.6f \n", et  );
 
         target = "CASSINI";
         frame  = "ECLIPJ2000";
         corrtn = "NONE";
         observ = "SUN";
 
         spkezr_c ( target, et, frame, corrtn, observ, state, &ltime );
 
         printf ( "state = %e %e %e %e %e %e\n",
                  state[0], state[1], state[2],
                  state[3], state[4], state[5] );
 
         return ( 0 );
      }
 
   Meta-kernel file ``spice_example.tm'':
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           'kernels/sclk/cas00084.tsc'
                           'kernels/spk/020514_SE_SAT105.bsp'
                           'kernels/spk/030201AP_SK_SM546_T45.bsp'
                           'kernels/spk/981005_PLTEPH-DE405S.bsp'
                           'kernels/spk/sat128.bsp'
                           )
      \begintext
 
 
Step-4: ``Sun Direction''
===========================================================================
 
 
``Sun Direction'' Task Statement
--------------------------------------------------------
 
   Extend the program from Step-3 to compute apparent direction of the Sun
   in the INMS frame at the epoch specified by SCLK time from Step-2.
 
 
``Sun Direction'' Hints
--------------------------------------------------------
 
   Determine the additional SPICE kernels needed to support the direction
   computation, knowing that they should provide the s/c and instrument
   frame orientation. Retrieve these kernels from the NAIF's FTP site.
 
   Verify that the orientation data in the kernels have adequate coverage
   to support computation of the direction of interest. Use ``ckbrief''
   utility program located under ``toolkit/exe'' directory for that.
 
   Modify the meta-kernel to load this(these) kernels.
 
   Determine the proper input arguments for the spkpos_c call to calculate
   the direction (which is the position portion of the output state). Look
   through the Frames Kernel find the name of the frame to used.
 
   Add calls to the routine(s), necessary variable declarations and output
   print statements to the program. Re-compile and re-link it against
   CSPICE.
 
 
``Sun Direction'' Solution Steps
--------------------------------------------------------
 
   A CASSINI spacecraft orientation CK file, providing s/c orientation with
   respect to an inertial frame, and CASSINI FK file, providing orientation
   of the INMS frame with respect to the s/ frame, are needed additionally
   to already loaded kernels to support computation of this direction.
 
   The file names can be added to the meta-kernel to get them loaded into
   the program:
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           'kernels/sclk/cas00084.tsc'
                           'kernels/spk/020514_SE_SAT105.bsp'
                           'kernels/spk/030201AP_SK_SM546_T45.bsp'
                           'kernels/spk/981005_PLTEPH-DE405S.bsp'
                           'kernels/spk/sat128.bsp'
                           'kernels/ck/04135_04171pc_psiv2.bc'
                           'kernels/fk/cas_v37.tf'
                           )
      \begintext
 
   The same highest level CSPICE routine computing positions, spkpos_c, can
   be used to compute this direction.
 
   Since this is the direction of the Sun as seen from the s/c, the target
   argument should be set to 'Sun' and the observer argument should be set
   to "CASSINI" The name of the INMS frame is "CASSINI_INMS", the
   definition and description of this frame are provided in the CASSINI FK
   file, ``cassini_v02.tf''.
 
   Since the apparent, or ``as seen'', position is sought for, the `abcorr'
   argument of the routine should be set to "LT+S" (see aberration
   correction discussion in the (``cspice/src/cspice/spkpos_c.c'') header).
 
   If desired, the position can then be turned into a unit vector using
   vhat_c function (``cspice/src/cspice/vhat_c.c''). Putting it all
   together, we get:
 
      SpiceDouble          sundir [3];
      ...
      target = "SUN";
      frame  = "CASSINI_INMS";
      corrtn = "LT+S";
      observ = "CASSINI";
 
      spkpos_c ( target, et, frame, corrtn, observ, sundir, &ltime );
      vhat_c   ( sundir, sundir );
 
   The updated program with added calls, required declarations and simple
   print statements produces the following output (the output below was
   generated by this program compiled with gcc on a PC running Linux; your
   output may differ slightly in its format and numeric precision):
 
      > ./spice_example
      utc  = 2004-06-11T19:32:00
      et   =     140254384.184625
      sclk = 1465674964.105
      et   =     140254384.183426
      state = -3.765991e+08 1.294488e+09 -7.064853e+06 -5.164226e+00
      8.017189e-01 4.060306e-02
      sundir = -2.902040e-01 8.816312e-01 3.721667e-01
      >
 
 
``Sun Direction'' Code
--------------------------------------------------------
 
   Program ``spice_example.c'':
 
 
         #include <stdio.h>
         #include "SpiceUsr.h"
 
         int main()
      {
         SpiceChar             * mkfile;
         SpiceChar             * utc;
         SpiceChar             * sclk;
         SpiceChar             * target;
         SpiceChar             * frame;
         SpiceChar             * corrtn;
         SpiceChar             * observ;
 
         SpiceDouble             et;
         SpiceDouble             state  [6];
         SpiceDouble             sundir [3];
         SpiceDouble             ltime;
 
         SpiceInt                scid;
 
 
         mkfile  =  "spice_example.tm";
         furnsh_c ( mkfile );
 
         utc     =  "2004-06-11T19:32:00";
         str2et_c ( utc, &et );
 
         printf ( "utc  = %s     \n", utc );
         printf ( "et   = %20.6f \n", et  );
 
         scid   = -82;
         sclk   = "1465674964.105";
         scs2e_c ( scid, sclk, &et );
 
         printf ( "sclk = %s     \n", sclk );
         printf ( "et   = %20.6f \n", et  );
 
         target = "CASSINI";
         frame  = "ECLIPJ2000";
         corrtn = "NONE";
         observ = "SUN";
 
         spkezr_c ( target, et, frame, corrtn, observ, state, &ltime );
 
         printf ( "state = %e %e %e %e %e %e\n",
                  state[0], state[1], state[2],
                  state[3], state[4], state[5] );
 
         target = "SUN";
         frame  = "CASSINI_INMS";
         corrtn = "LT+S";
         observ = "CASSINI";
 
         spkpos_c ( target, et, frame, corrtn, observ, sundir, &ltime );
         vhat_c   ( sundir, sundir );
 
         printf ( "sundir = %e %e %e\n",
                  sundir[0], sundir[1], sundir[2] );
 
         return ( 0 );
      }
 
   Meta-kernel file ``spice_example.tm'':
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           'kernels/sclk/cas00084.tsc'
                           'kernels/spk/020514_SE_SAT105.bsp'
                           'kernels/spk/030201AP_SK_SM546_T45.bsp'
                           'kernels/spk/981005_PLTEPH-DE405S.bsp'
                           'kernels/spk/sat128.bsp'
                           'kernels/ck/04135_04171pc_psiv2.bc'
                           'kernels/fk/cas_v37.tf'
                           )
      \begintext
 
 
Step-5: ``Sub-Spacecraft Point''
===========================================================================
 
 
``Sub-Spacecraft Point'' Task Statement
--------------------------------------------------------
 
   Extend the program from Step-4 to compute planetocentric longitude and
   and latitude of the sub-spacecraft point on Phoebe, and the direction
   from the spacecraft to that point in the INMS frame.
 
 
``Sub-Spacecraft Point'' Hints
--------------------------------------------------------
 
   Find the CSPICE routine that computes sub-observer point coordinates.
   Use ``Most Used CSPICE APIs'' or ``subpt'' cookbook program for that.
 
   Refer to the routine's header to determine the additional kernels needed
   for this direction computation. Get these kernels from the NAIF's FTP
   site. Modify the meta-kernel to load this(these) kernels.
 
   Determine the proper input arguments for the routine. Refer to the
   routine's header for that information.
 
   Convert the surface point Cartesian vector returned by this routine to
   latitudinal coordinates. Use ``Permuted Index'' to find the routine that
   does this conversion. Refer to the routine's header for input/output
   argument specifications.
 
   Since the Cartesian vector from the spacecraft to the sub-spacecraft
   point is computed in the Phoebe body-fixed frame, it should be
   transformed into the instrument frame get the direction we are looking
   for. Refer to ``frames.req'' and/or ``Frames'' tutorial to determine the
   name of the routine computing transformations and use it to compute
   transformation from Phoebe body-fixed to the INMS frame.
 
   Using ``Permuted Index'' find the routine that multiplies 3x3 matrix by
   3d vector and use it to rotate the vector to the instrument frame.
 
   Add calls to the routine(s), necessary variable declarations and output
   print statements to the program. Re-compile and re-link it against
   CSPICE.
 
 
``Sub-Spacecraft Point'' Solution Steps
--------------------------------------------------------
 
   The subpnt_c routine (``cspice/src/cspice/subpnt_c.c'') can be used to
   compute the sub-observer point and the vector from the observer to that
   point with a single call. To determine this point as the closest point
   on the Phoebe ellipsoid, the `method' argument has to be set to "NEAR
   POINT: ELLIPSOID". For our case the `target' is "PHOEBE", the target
   body-fixed frame is "IAU_PHOEBE", and the observer is "CASSINI".
 
   Since the s/c is close to Phoebe, light time does not need to be taken
   into account and, therefore, the `abcorr' argument can be set to "NONE".
 
   In order for subpnt_c to compute the nearest point location, a PCK file
   containing Phoebe radii has to be loaded into the program (see ``Files''
   section of the routine's header.) All other files required for this
   computation are already being loaded by the program. With PCK file name
   added to it, the updated meta-kernel will look like this:
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           'kernels/sclk/cas00084.tsc'
                           'kernels/spk/020514_SE_SAT105.bsp'
                           'kernels/spk/030201AP_SK_SM546_T45.bsp'
                           'kernels/spk/981005_PLTEPH-DE405S.bsp'
                           'kernels/spk/sat128.bsp'
                           'kernels/ck/04135_04171pc_psiv2.bc'
                           'kernels/fk/cas_v37.tf'
                           'kernels/pck/cpck05Mar2004.tpc'
                           )
      \begintext
 
   The sub-spacecraft point Cartesian vector can be converted to
   planetocentric radius, longitude and latitude using the reclat_c routine
   (``cspice/src/cspice/reclat_c.c'').
 
   The vector from the spacecraft to the sub-spacecraft point returned by
   subpnt_c has to be rotated from the body-fixed frame to the instrument
   frame. The name of the routine that computes 3x3 matrices rotating
   vectors from one frame to another is pxform_c
   (``cspice/src/cspice/pxform_c.c''). In our case the `from' argument
   should be set to "IAU_PHOEBE" and the `to' argument should be set to
   "CASSINI_INMS"
 
   The vector should be then multiplied by this matrix to rotate it to the
   instrument frame. The mxv_c routine performs that function
   (``cspice/src/cspice/mxv_c.c'')
 
   After applying the rotation, normalize the resultant vector using the
   vhat_c function.
 
   For output the longitude and latitude angles returned by reclat_c in
   radians can be converted to degrees by multiplying by dpr_c function
   (``cspice/src/cspice/dpr_c.c''). Putting it all together, we get:
 
         SpiceChar             * method;
         SpiceChar             * fromfr;
         SpiceChar             * tofr;
 
         SpiceDouble             spoint [3];
         SpiceDouble             trgepc;
         SpiceDouble             srfvec [3];
         SpiceDouble             srad;
         SpiceDouble             slon;
         SpiceDouble             slat;
         SpiceDouble             sbpdir [3];
         SpiceDouble             m2imat [3][3];
 
         ...
 
         method = "NEAR POINT: ELLIPSOID";
         target = "PHOEBE";
         frame  = "IAU_PHOEBE";
         corrtn = "NONE";
         observ = "CASSINI";
 
         subpnt_c ( method, target, et, frame, corrtn, observ,
                                            spoint, &trgepc, srfvec );
 
         reclat_c ( spoint, &srad, &slon, &slat );
 
         fromfr = "IAU_PHOEBE";
         tofr   = "CASSINI_INMS";
 
         pxform_c ( fromfr, tofr, et, m2imat );
 
         mxv_c  ( m2imat, srfvec, sbpdir );
         vhat_c ( sbpdir, sbpdir );
 
         printf ( "lon    =  %e \n", slon * dpr_c()  );
         printf ( "lat    =  %e \n", slat * dpr_c()  );
         printf ( "sbpdir =  %e %e %e \n",
                   sbpdir[0], sbpdir[1], sbpdir[2]       );
 
   The updated program with added calls, required declarations and simple
   print statements produces the following output (the output below was
   generated by this program compiled with gcc on a PC running Linux; your
   output may differ slightly in its format and numeric precision):
 
      > ./spice_example
      utc  = 2004-06-11T19:32:00
      et   =     140254384.184625
      sclk = 1465674964.105
      et   =     140254384.183426
      state = -3.765991e+08 1.294488e+09 -7.064853e+06 -5.164226e+00
      8.017189e-01 4.060306e-02
      sundir = -2.902040e-01 8.816312e-01 3.721667e-01
      lon    =  2.342316e+01
      lat    =  3.709797e+00
      sbpdir =  -7.762071e-04 -9.998732e-01 -1.590546e-02
      >
 
 
``Sub-Spacecraft Point'' Code
--------------------------------------------------------
 
   Program ``spice_example.c'':
 
 
         #include <stdio.h>
         #include "SpiceUsr.h"
 
         int main()
      {
         SpiceChar             * mkfile;
         SpiceChar             * utc;
         SpiceChar             * sclk;
         SpiceChar             * target;
         SpiceChar             * frame;
         SpiceChar             * corrtn;
         SpiceChar             * observ;
         SpiceChar             * method;
         SpiceChar             * fromfr;
         SpiceChar             * tofr;
 
         SpiceDouble             et;
         SpiceDouble             state  [6];
         SpiceDouble             sundir [3];
         SpiceDouble             ltime;
         SpiceDouble             spoint [3];
         SpiceDouble             trgepc;
         SpiceDouble             srfvec [3];
         SpiceDouble             srad;
         SpiceDouble             slon;
         SpiceDouble             slat;
         SpiceDouble             sbpdir [3];
         SpiceDouble             m2imat [3][3];
 
         SpiceInt                scid;
 
 
         mkfile  =  "spice_example.tm";
         furnsh_c ( mkfile );
 
         utc     =  "2004-06-11T19:32:00";
         str2et_c ( utc, &et );
 
         printf ( "utc  = %s     \n", utc );
         printf ( "et   = %20.6f \n", et  );
 
         scid   = -82;
         sclk   = "1465674964.105";
         scs2e_c ( scid, sclk, &et );
 
         printf ( "sclk = %s     \n", sclk );
         printf ( "et   = %20.6f \n", et  );
 
         target = "CASSINI";
         frame  = "ECLIPJ2000";
         corrtn = "NONE";
         observ = "SUN";
 
         spkezr_c ( target, et, frame, corrtn, observ, state, &ltime );
 
         printf ( "state = %e %e %e %e %e %e\n",
                  state[0], state[1], state[2],
                  state[3], state[4], state[5] );
 
         target = "SUN";
         frame  = "CASSINI_INMS";
         corrtn = "LT+S";
         observ = "CASSINI";
 
         spkpos_c ( target, et, frame, corrtn, observ, sundir, &ltime );
         vhat_c   ( sundir, sundir );
 
         printf ( "sundir = %e %e %e\n",
                  sundir[0], sundir[1], sundir[2] );
 
         method = "NEAR POINT: ELLIPSOID";
         target = "PHOEBE";
         frame  = "IAU_PHOEBE";
         corrtn = "NONE";
         observ = "CASSINI";
 
         subpnt_c ( method, target, et, frame, corrtn, observ,
                                            spoint, &trgepc, srfvec );
 
         reclat_c ( spoint, &srad, &slon, &slat );
 
         fromfr = "IAU_PHOEBE";
         tofr   = "CASSINI_INMS";
 
         pxform_c ( fromfr, tofr, et, m2imat );
 
         mxv_c  ( m2imat, srfvec, sbpdir );
         vhat_c ( sbpdir, sbpdir );
 
         printf ( "lon    =  %e \n", slon * dpr_c()  );
         printf ( "lat    =  %e \n", slat * dpr_c()  );
         printf ( "sbpdir =  %e %e %e \n",
                   sbpdir[0], sbpdir[1], sbpdir[2]       );
 
         return ( 0 );
      }
 
   Meta-kernel file ``spice_example.tm'':
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           'kernels/sclk/cas00084.tsc'
                           'kernels/spk/020514_SE_SAT105.bsp'
                           'kernels/spk/030201AP_SK_SM546_T45.bsp'
                           'kernels/spk/981005_PLTEPH-DE405S.bsp'
                           'kernels/spk/sat128.bsp'
                           'kernels/ck/04135_04171pc_psiv2.bc'
                           'kernels/fk/cas_v37.tf'
                           'kernels/pck/cpck05Mar2004.tpc'
                           )
      \begintext
 
 
Step-6: ``Spacecraft Velocity''
===========================================================================
 
 
``Spacecraft Velocity'' Task Statement
--------------------------------------------------------
 
   Extend the program from Step-5 to compute the spacecraft velocity with
   respect to Phoebe in the INMS frame.
 
 
``Spacecraft Velocity'' Hints
--------------------------------------------------------
 
   Compute velocity of the spacecraft with respect to Phoebe in some
   inertial frame, for example J2000. Recall that velocity is the last
   three components of the state vector returned by spkezr_c.
 
   Since the velocity vector is computed in the inertial frame, it should
   be rotated to the instrument frame. Look at the previous step the
   routine that compute necessary rotation and rotate vectors.
 
   Add calls to the routine(s), necessary variable declarations and output
   print statements to the program. Re-compile and re-link it against
   CSPICE.
 
 
``Spacecraft Velocity'' Solution Steps
--------------------------------------------------------
 
   All kernels required for computations in this step are already being
   loaded by the program, therefore, the meta-kernel does not need to be
   changed.
 
   The spacecraft velocity vector is the last three components of the state
   returned by spkezr_c. To compute velocity of CASSINI with respect to
   Phoebe in the J2000 inertial frame the spkezr_c arguments should be set
   to "CASSINI" (TARG), "PHOEBE" (OBS), "J2000" (REF) and "NONE" (ABCORR).
 
   For convenience the velocity can be copied from the output state in to a
   3d vector using the vpack_c routine (``cspice/src/cspice/vpack_c.c'').
 
   The computed velocity vector has to be rotated from the J2000 frame to
   the instrument frame. The pxform_c routine used in the previous step can
   be used to compute the rotation matrix needed for that. In this case the
   frame name arguments should be set to "J2000" (FROM) and "CASSINI_INMS"
   (TO).
 
   As in the previous step the difference vector should be then multiplied
   by this rotation matrix using the mxv_c routine. After applying the
   rotation, normalize the resultant vector using the vhat_c routine.
 
   Putting it all together, we get:
 
         SpiceDouble             scvdir [3];
         SpiceDouble             j2imat [3][3];
 
         ...
 
         target = "CASSINI";
         frame  = "J2000";
         corrtn = "NONE";
         observ = "PHOEBE";
 
         spkezr_c ( target, et, frame, corrtn, observ, state, &ltime );
         vpack_c  ( state[3], state[4], state[5], scvdir );
 
         fromfr = "J2000";
         tofr   = "CASSINI_INMS";
         pxform_c ( fromfr, tofr, et, j2imat );
 
         mxv_c  ( j2imat, scvdir, scvdir );
         vhat_c ( scvdir, scvdir );
 
         printf ( "scvdir = %e %e %e \n",
                   scvdir[0], scvdir[1], scvdir[2] );
 
   The updated program with added calls, required declarations and simple
   print statements produces the following output (the output below was
   generated by this program compiled with gcc on a PC running Linux; your
   output may differ slightly in its format and numeric precision):
 
      > ./spice_example
      utc  = 2004-06-11T19:32:00
      et   =     140254384.184625
      sclk = 1465674964.105
      et   =     140254384.183426
      state = -3.765991e+08 1.294488e+09 -7.064853e+06 -5.164226e+00
      8.017189e-01 4.060306e-02
      sundir = -2.902040e-01 8.816312e-01 3.721667e-01
      lon    =  2.342316e+01
      lat    =  3.709797e+00
      sbpdir =  -7.762071e-04 -9.998732e-01 -1.590546e-02
      scvdir = 3.957849e-01 -2.928077e-01 8.704125e-01
      >
 
   Note that computing the spacecraft velocity in the instrument frame by a
   single call to spkezr_c by specifying "CASSINI_INMS" in the `ref'
   argument returns an incorrect result. Such computation will take into
   account the spacecraft angular velocity from the CK files, which should
   not be considered in this case.
 
 
``Spacecraft Velocity'' Code Program ``spice_example.c'':
--------------------------------------------------------
 
         #include <stdio.h>
         #include "SpiceUsr.h"
 
         int main()
      {
         SpiceChar             * mkfile;
         SpiceChar             * utc;
         SpiceChar             * sclk;
         SpiceChar             * target;
         SpiceChar             * frame;
         SpiceChar             * corrtn;
         SpiceChar             * observ;
         SpiceChar             * method;
         SpiceChar             * fromfr;
         SpiceChar             * tofr;
 
         SpiceDouble             et;
         SpiceDouble             state  [6];
         SpiceDouble             sundir [3];
         SpiceDouble             ltime;
         SpiceDouble             spoint [3];
         SpiceDouble             trgepc;
         SpiceDouble             srfvec [3];
         SpiceDouble             srad;
         SpiceDouble             slon;
         SpiceDouble             slat;
         SpiceDouble             sbpdir [3];
         SpiceDouble             m2imat [3][3];
         SpiceDouble             scvdir [3];
         SpiceDouble             j2imat [3][3];
 
         SpiceInt                scid;
 
 
         mkfile  =  "spice_example.tm";
         furnsh_c ( mkfile );
 
         utc     =  "2004-06-11T19:32:00";
         str2et_c ( utc, &et );
 
         printf ( "utc  = %s     \n", utc );
         printf ( "et   = %20.6f \n", et  );
 
         scid   = -82;
         sclk   = "1465674964.105";
         scs2e_c ( scid, sclk, &et );
 
         printf ( "sclk = %s     \n", sclk );
         printf ( "et   = %20.6f \n", et  );
 
         target = "CASSINI";
         frame  = "ECLIPJ2000";
         corrtn = "NONE";
         observ = "SUN";
 
         spkezr_c ( target, et, frame, corrtn, observ, state, &ltime );
 
         printf ( "state = %e %e %e %e %e %e\n",
                  state[0], state[1], state[2],
                  state[3], state[4], state[5] );
 
         target = "SUN";
         frame  = "CASSINI_INMS";
         corrtn = "LT+S";
         observ = "CASSINI";
 
         spkpos_c ( target, et, frame, corrtn, observ, sundir, &ltime );
         vhat_c   ( sundir, sundir );
 
         printf ( "sundir = %e %e %e\n",
                  sundir[0], sundir[1], sundir[2] );
 
         method = "NEAR POINT: ELLIPSOID";
         target = "PHOEBE";
         frame  = "IAU_PHOEBE";
         corrtn = "NONE";
         observ = "CASSINI";
 
         subpnt_c ( method, target, et, frame, corrtn, observ,
                                            spoint, &trgepc, srfvec );
 
         reclat_c ( spoint, &srad, &slon, &slat );
 
         fromfr = "IAU_PHOEBE";
         tofr   = "CASSINI_INMS";
 
         pxform_c ( fromfr, tofr, et, m2imat );
 
         mxv_c  ( m2imat, srfvec, sbpdir );
         vhat_c ( sbpdir, sbpdir );
 
         printf ( "lon    =  %e \n", slon * dpr_c()  );
         printf ( "lat    =  %e \n", slat * dpr_c()  );
         printf ( "sbpdir =  %e %e %e \n",
                   sbpdir[0], sbpdir[1], sbpdir[2]       );
 
         target = "CASSINI";
         frame  = "J2000";
         corrtn = "NONE";
         observ = "PHOEBE";
 
         spkezr_c ( target, et, frame, corrtn, observ, state, &ltime );
         vpack_c  ( state[3], state[4], state[5], scvdir );
 
         fromfr = "J2000";
         tofr   = "CASSINI_INMS";
         pxform_c ( fromfr, tofr, et, j2imat );
 
         mxv_c  ( j2imat, scvdir, scvdir );
         vhat_c ( scvdir, scvdir );
 
         printf ( "scvdir = %e %e %e \n",
                  scvdir[0], scvdir[1], scvdir[2] );
 
 
         return ( 0 );
      }
 
   Meta-kernel file ``spice_example.tm'':
 
      \begindata
         KERNELS_TO_LOAD = (
                           'kernels/lsk/naif0008.tls'
                           'kernels/sclk/cas00084.tsc'
                           'kernels/spk/020514_SE_SAT105.bsp'
                           'kernels/spk/030201AP_SK_SM546_T45.bsp'
                           'kernels/spk/981005_PLTEPH-DE405S.bsp'
                           'kernels/spk/sat128.bsp'
                           'kernels/ck/04135_04171pc_psiv2.bc'
                           'kernels/fk/cas_v37.tf'
                           'kernels/pck/cpck05Mar2004.tpc'
                           )
      \begintext
 
