Remote Sensing Hands-On Lesson (C)





October 14, 2004



Overview




In this lesson you will develop a series of simple programs that demonstrate the usage of CSPICE to compute a variety of different geometric quantities applicable to experiments carried out by a remote sensing instrument flown on an interplanetary spacecraft. This particular lesson focuses on a framing camera flying on the Cassini spacecraft, but many of the concepts are easily extended and generalized to other scenarios.



References






Tutorials



The following SPICE tutorials are referred to by the discussions in this lesson:

   Name             Lesson steps/functions it describes
   ---------------  -----------------------------------------
   Time             Time Conversion
   SCLK and LSK     Time Conversion
   SPK              Obtaining Ephemeris Data
   Frames           Reference Frames
   Using Frames     Reference Frames
   PCK              Planetary Constants Data
   CK               Spacecraft Orientation Data
These tutorials are available from the NAIF ftp server at JPL:

   ftp://naif.jpl.nasa.gov/pub/naif/toolkit_docs/Tutorials


Required Readings



The Required Reading documents are provided with the Toolkit and are located under the ``cspice/doc'' directory in the C installation tree.

   Name             Lesson steps/functions that it describes
   ---------------  -----------------------------------------
   time.req         Time Conversion
   sclk.req         SCLK Time Conversion
   spk.req          Obtaining Ephemeris Data
   frames.req       Using Reference Frames
   pck.req          Obtaining Planetary Constants Data
   ck.req           Obtaining Spacecraft Orientation Data
   naif_ids.req     Determining Body ID Codes


The Permuted Index



Another useful document distributed with the Toolkit is the permuted index. This is located under the ``cspice/doc'' directory in the C installation tree. This text document provides a simple mechanism to discover what CSPICE functions perform a particular function of interest as well as the name of the source module that contains the function.



Source Code Headers



The most detailed specification of a given CSPICE function is contained in the header section of its source code. The source code is distributed with the Toolkit and is located under ``cspice/src/cspice'' in the C versions. For example the header of str2et_c is contained in the file:

   cspice/src/cspice/str2et_c.c


Kernels Used




The programs that are produced in the course of this lesson will compute geometry for the Cassini orbiter. The following CASSINI SPICE kernels will be used:

   #  FILE NAME                 TYPE  DESCRIPTION
   -- ------------------------- ----  ------------------------
   1  naif0007.tls              LSK   Generic LSK
   2  cas00084.tsc              SCLK  Cassini SCLK
   3  sat128.bsp                SPK   Saturnian Satellite Ephemeris
   4  981005_PLTEPH-DE405S.bsp  SPK   Solar System Ephemeris
   5  020514_SE_SAT105.bsp      SPK   Saturnian Satellite Ephemeris
   6  030201AP_SK_SM546_T45.bsp SPK   Cassini Spacecraft SPK
   7  cas_v37.tf                FK    Cassini FK
   8  04135_04171pc_psiv2.bc    CK    Cassini Spacecraft CK
   9  cpck05Mar2004.tpc         PCK   Cassini Project PCK
   10 cas_iss_v09.ti            IK    ISS Instrument Kernel


CSPICE Modules Used




This section provides a complete summary of the functions, and the kernels that are suggested for usage in each of the exercises in this tutorial. (You may wish to not look at this list unless/until you ``get stuck'' while working on your own.)

   CHAPTER EXERCISE   FUNCTIONS  NON-VOID   KERNELS
   ------- ---------  ---------  ---------  -------
     1     convtm     furnsh_c              1,2
                      prompt_c
                      str2et_c
                      etcal_c
                      timout_c
                      sce2c_c
                      sce2s_c
 
     2     getsta     furnsh_c   vnorm_c    1,3-6
                      prompt_c
                      str2et_c
                      spkezr_c
                      spkpos_c
                      convrt_c
 
     3     xform      furnsh_c   vsep_c     1-9
                      prompt_c
                      str2et_c
                      spkezr_c
                      sxform_c
                      mxvg_c
                      spkpos_c
                      pxform_c
                      mxv_c
                      convrt_c
 
     4     subpts     furnsh_c              1,3-6,9
                      prompt_c
                      str2et_c
                      subpt_c
                      subsol_c
 
     5     fovint     furnsh_c   dpr_c      1-10
                      prompt_c
                      str2et_c
                      bodn2c_c
                      byebye_c
                      getfov_c
                      srfxpt_c
                      reclat_c
 
     6     angles     furnsh_c   dpr_c      1-10
                      prompt_c
                      str2et_c
                      bodn2c_c
                      byebye_c
                      getfov_c
                      srfxpt_c
                      reclat_c
                      illum_c
                      et2lst_c
Refer to the headers of the various functions listed above, as detailed interface specifications are provided with the source code.



Time Conversion (convtm)







Task Statement




Write a program that prompts the user for an input UTC time string, converts it to the following time systems and output formats:

and displays the results. Use the program to convert "2004 jun 11 19:32:00" UTC into these alternate systems.



Learning Goals




Familiarity with the various time conversion and parsing functions available in the Toolkit. Exposure to source code headers and their usage in learning to call functions.



Approach




The solution to the problem can be broken down into a series of simple steps:

You may find it useful to consult the permuted index, the headers of various source modules, and the ``Time Required Reading'' and ``SCLK Required Reading'' documents.

When completing the ``calendar format'' step above, consider using one of two possible methods: etcal_c or timout_c.



Solution






Solution Meta-Kernel



The meta-kernel we created for the solution to this exercise is named 'convtm.mk'. Its contents follow:

   KPL/MK
 
   This is the meta-kernel used in the solution of the ``Time
   Conversion'' task in the Remote Sensing Hands On Lesson.
 
      \begindata
      KERNELS_TO_LOAD = ( 'kernels/lsk/naif0007.tls',
                          'kernels/sclk/cas00084.tsc' )
      \begintext
 
 


Solution Source Code



A sample solution to the problem follows:

      #include <stdio.h>
 
      /*
      Standard CSPICE User Include File
      */
      #include "SpiceUsr.h"
 
      /*
      Local Parameters
      */
 
      #define METAKR "convtm.mk"
      #define SCLKID -82
      #define STRLEN 50
 
      int main (void)
      {
 
         /*
         Local Variables
         */
         SpiceChar               calet  [STRLEN];
         SpiceChar               sclkst [STRLEN];
         SpiceChar               utctim [STRLEN];
 
         SpiceDouble             et;
 
         /*
         Load the kernels this program requires.
         Both the spacecraft clock kernel and a
         leapseconds kernel should be listed in
         the meta-kernel.
         */
         furnsh_c ( METAKR );
 
         /*
         Prompt the user for the input time string.
         */
         prompt_c ( "Input UTC Time: ", STRLEN, utctim );
 
         printf ( "Converting UTC Time: %s\n", utctim );
 
         /*
         Convert utctim to ET.
         */
         str2et_c ( utctim, &et );
 
         printf ( "   ET Seconds Past 2000: %16.3f\n", et );
 
         /*
         Now convert ET to a calendar time
         string.  This can be accomplished in two
         ways.
         */
         etcal_c ( et, STRLEN, calet );
 
         printf ( "   Calendar ET (etcal_c): %s\n", calet );
 
         /*
         Or use timout_c for finer control over the
         output format.  The picture below was built
         by examining the header of timout_c.
         */
         timout_c ( et,     "YYYY-MON-DDTHR:MN:SC ::TDB",
                    STRLEN, calet                         );
 
         printf ( "   Calendar ET (timout_c): %s\n", calet );
 
         /*
         Convert ET to spacecraft clock time.
         */
         sce2s_c ( SCLKID, et, STRLEN, sclkst );
 
         printf ( "   Spacecraft Clock Time: %s\n", sclkst );
 
         return(0);
      }


Solution Sample Output



After compiling the program, execute it:

   Converting UTC Time: 2004 jun 11 19:32:00
      ET Seconds Past 2000:    140254384.185
      Calendar ET (etcal_c): 2004 JUN 11 19:33:04.184
      Calendar ET (timout_c): 2004-JUN-11T19:33:04
      Spacecraft Clock Time: 1/1465674964.105


Obtaining Target States and Positions (getsta)







Task Statement




Write a program that prompts the user for an input UTC time string, computes the following quantities at that epoch:

and displays the results. Use the program to compute these quantities at "2004 jun 11 19:32:00" UTC.



Learning Goals




Understand the anatomy of an spkezr_c call. Discover the difference between spkezr_c and spkpos_c. Familiarity with the Toolkit utility ``brief''. Exposure to unit conversion with CSPICE.



Approach




The solution to the problem can be broken down into a series of simple steps:

You may find it useful to consult the permuted index, the headers of various source modules, and the ``SPK Required Reading'' document.

When deciding which SPK files to load, the Toolkit utility ``brief'' may be of some use.

``brief'' is located in the ``cspice/exe'' directory for C toolkits. Consult its user's guide available in ``cspice/doc/brief.ug'' for details.



Solution






Solution Meta-Kernel



The meta-kernel we created for the solution to this exercise is named 'getsta.mk'. Its contents follow:

   KPL/MK
 
   This is the meta-kernel used in the solution of the
   ``Obtaining Target States and Positions'' task in the
   Remote Sensing Hands On Lesson.
 
      \begindata
      KERNELS_TO_LOAD = ( 'kernels/lsk/naif0007.tls',
                          'kernels/spk/sat128.bsp'
                          'kernels/spk/981005_PLTEPH-DE405S.bsp',
                          'kernels/spk/020514_SE_SAT105.bsp',
                          'kernels/spk/030201AP_SK_SM546_T45.bsp' )
      \begintext


Solution Source Code



A sample solution to the problem follows:

      #include <stdio.h>
 
      /*
      Standard CSPICE User Include File
      */
      #include "SpiceUsr.h"
 
      /*
      Local Parameters
      */
 
      #define METAKR "getsta.mk"
      #define STRLEN 50
 
      int main (void)
      {
         /*
         Local Variables
         */
         SpiceChar               utctim [STRLEN];
 
         SpiceDouble             dist;
         SpiceDouble             et;
         SpiceDouble             ltime;
         SpiceDouble             pos   [3];
         SpiceDouble             state [6];
 
         /*
         Load the kernels that this program requires.  We
         will need a leapseconds kernel to convert input
         UTC time strings into ET.  We also will need the
         necessary SPK files with coverage for the bodies
         in which we are interested.
         */
         furnsh_c ( METAKR );
 
         /*
         Prompt the user for the input time string.
         */
         prompt_c (  "Input UTC Time: ", STRLEN, utctim );
 
         printf ( "Converting UTC Time: %s\n", utctim  );
 
         /*
         Convert utctim to ET.
         */
         str2et_c ( utctim, &et );
 
         printf ( "   ET Seconds Past 2000: %16.3f\n", et );
 
         /*
         Compute the apparent state of Phoebe as seen from
         CASSINI in the J2000 frame.  All of the ephemeris
         readers return states in units of kilometers and
         kilometers per second.
         */
         spkezr_c ( "PHOEBE", et,    "J2000", "LT+S",
                    "CASSINI",  state, <ime          );
 
         printf ( "   Apparent State of Phoebe as seen "
                  "from CASSINI in the J2000\n"        );
         printf ( "      frame (km, km/s):\n"          );
         printf ( "      X = %16.3f\n", state[0]       );
         printf ( "      Y = %16.3f\n", state[1]       );
         printf ( "      Z = %16.3f\n", state[2]       );
         printf ( "     VX = %16.3f\n", state[3]       );
         printf ( "     VY = %16.3f\n", state[4]       );
         printf ( "     VZ = %16.3f\n", state[5]       );
 
         /*
         Compute the apparent position of Earth as seen from
         CASSINI in the J2000 frame.  Note: We could have
         continued using spkezr_c and simply ignored the
         velocity components.
         */
         spkpos_c ( "EARTH", et,  "J2000", "LT+S",
                    "CASSINI",   pos, <ime             );
 
         printf ( "   Apparent position of Earth as "
                  "seen from CASSINI in the J2000\n"     );
         printf ( "      frame (km): \n"                 );
         printf ( "      X = %16.3f\n", pos[0]           );
         printf ( "      Y = %16.3f\n", pos[1]           );
         printf ( "      Z = %16.3f\n", pos[2]           );
 
 
         /*
         We need only display LTIME, as it is precisely the
         light time in which we are interested.
         */
         printf ( "   One way light time between CASSINI and "
                  "the apparent position\n"                );
         printf ( "      of Earth (seconds): %16.3f\n", ltime  );
 
         /*
         Compute the apparent position of the Sun as seen
         from Phoebe in the J2000 frame.
         */
         spkpos_c ( "SUN",  et,  "J2000", "LT+S",
                    "PHOEBE", pos, <ime                );
 
         printf ( "   Apparent position of Sun as seen "
                  "from Phoebe in the\n"                 );
         printf ( "      J2000 frame (km): \n"           );
         printf ( "      X = %16.3f\n", pos[0]           );
         printf ( "      Y = %16.3f\n", pos[1]           );
         printf ( "      Z = %16.3f\n", pos[2]           );
 
         /*
         Now we need to compute the actual distance between
         the Sun and Phoebe.  The above SPKPOS call gives us
         the apparent distance, so we need to adjust our
         aberration correction appropriately.
         */
         spkpos_c ( "SUN",  et,  "J2000", "NONE",
                    "PHOEBE", pos, <ime                );
 
         /*
         Compute the distance between the body centers in
         kilometers.
         */
         dist = vnorm_c ( pos );
 
         /*
         Convert this value to AU using convrt_c.
         */
         convrt_c ( dist, "KM", "AU", &dist );
 
         printf ( "   Actual distance between Sun and "
                  "Phoebe\n"                             );
         printf ( "      (AU): %16.3f\n", dist           );
 
         return(0);
      }


Solution Sample Output



After compiling the program, execute it:

   Converting UTC Time: 2004 jun 11 19:32:00
      ET Seconds Past 2000:    140254384.185
      Apparent State of Phoebe as seen from CASSINI in the J2000
         frame (km, km/s):
         X =         -119.921
         Y =         2194.139
         Z =          -57.639
        VX =           -5.980
        VY =           -2.119
        VZ =           -0.295
      Apparent position of Earth as seen from CASSINI in the J2000
         frame (km):
         X =    353019393.123
         Y =  -1328180352.140
         Z =   -568134171.697
      One way light time between CASSINI and the apparent position
         of Earth (seconds):         4960.427
      Apparent position of Sun as seen from Phoebe in the
         J2000 frame (km):
         X =    376551465.272
         Y =  -1190495630.303
         Z =   -508438699.110
      Actual distance between Sun and Phoebe
         (AU):            9.012


Spacecraft Orientation and Reference Frames (xform)







Task Statement




Write a program that prompts the user for an input time string, computes the following at the epoch of interest:

and displays the results. Use the program to compute these quantities at the epoch "2004 jun 11 19:32:00" UTC.



Learning Goals




Familiarity with the different types of kernels involved in chaining reference frames together, both inertial and non-inertial. Discover some of the matrix and vector math functions. Understand the difference between pxform_c and sxform_c.



Approach




The solution to the problem can be broken down into a series of simple steps:

HINT: Several of the steps above may be compressed into a single using CSPICE functions with which you are already familiar. The ``long-way'' presented above is intended to facilitate the introduction of the functions pxform_c and sxform_c.

You may find it useful to consult the permuted index, the headers of various source modules, and the following toolkit documentation:

This particular example makes use of many of the different types of SPICE kernels. You should spend a few moments thinking about which kernels you will need and what data they provide.



Solution






Solution Meta-Kernel



The meta-kernel we created for the solution to this exercise is named 'xform.mk'. Its contents follow:

   KPL/MK
 
   This is the meta-kernel used in the solution of the ``Spacecraft
   Orientation and Reference Frames'' task in the Remote Sensing
   Hands On Lesson.
 
      \begindata
      KERNELS_TO_LOAD = ( 'kernels/lsk/naif0007.tls',
                          'kernels/sclk/cas00084.tsc',
                          'kernels/spk/sat128.bsp'
                          'kernels/spk/981005_PLTEPH-DE405S.bsp',
                          'kernels/spk/020514_SE_SAT105.bsp',
                          'kernels/spk/030201AP_SK_SM546_T45.bsp',
                          'kernels/fk/cas_v37.tf',
                          'kernels/ck/04135_04171pc_psiv2.bc',
                          'kernels/pck/cpck05Mar2004.tpc' )
      \begintext


Solution Source Code



A sample solution to the problem follows:

      #include <stdio.h>
 
      /*
      Standard CSPICE User Include File
      */
      #include "SpiceUsr.h"
 
      /*
      Local Parameters
      */
 
      #define METAKR "xform.mk"
      #define STRLEN 50
 
      int main (void)
      {
 
         /*
         Local Variables
         */
         SpiceChar               utctim [STRLEN];
 
         SpiceDouble             et;
         SpiceDouble             ltime;
         SpiceDouble             state  [6];
         SpiceDouble             bfixst [6];
         SpiceDouble             pos    [3];
         SpiceDouble             sform  [6][6];
         SpiceDouble             pform  [3][3];
         SpiceDouble             bsight [3];
         SpiceDouble             sep;
 
 
         /*
         Load the kernels that this program requires.  We
         will need:
 
            A leapseconds kernel
            A spacecraft clock kernel for CASSINI
            The necessary ephemerides
            A planetary constants file (PCK)
            A spacecraft orientation kernel for CASSINI (CK)
            A frame kernel (TF)
         */
         furnsh_c ( METAKR );
 
         /*
         Prompt the user for the input time string.
         */
         prompt_c ( "Input UTC Time: ", STRLEN, utctim );
 
         printf ( "Converting UTC Time: %s\n", utctim );
 
         /*
         Convert utctim to ET.
         */
         str2et_c ( utctim, &et );
 
         printf ( "   ET Seconds Past 2000: %16.3f\n", et );
 
         /*
         Compute the apparent state of Phoebe as seen from
         CASSINI in the J2000 reference frame.
         */
         spkezr_c ( "PHOEBE", et,    "J2000", "LT+S",
                    "CASSINI",  state, <ime              );
 
         /*
         Now obtain the transformation from the inertial
         J2000 frame to the non-inertial body-fixed IAU_PHOEBE
         frame.  Since we want the apparent position, we
         need to subtract ltime from et.
         */
         sxform_c ( "J2000", "IAU_PHOEBE", et-ltime, sform );
 
         /*
         Now rotate the apparent J2000 state into IAU_PHOEBE
         with the following matrix multiplication:
         */
         mxvg_c ( sform, state, 6, 6, bfixst );
 
         /*
         Display the results.
         */
         printf ( "   Apparent state of Phoebe as seen "
                  "from CASSINI in the IAU_PHOEBE\n"    );
         printf ( "      body-fixed frame (km, km/s):\n");
         printf ( "      X = %19.6f\n", bfixst[0]       );
         printf ( "      Y = %19.6f\n", bfixst[1]       );
         printf ( "      Z = %19.6f\n", bfixst[2]       );
         printf ( "     VX = %19.6f\n", bfixst[3]       );
         printf ( "     VY = %19.6f\n", bfixst[4]       );
         printf ( "     VZ = %19.6f\n", bfixst[5]       );
 
         /*
         It is worth pointing out, all of the above could
         have been done with a single use of spkezr_c:
         */
         spkezr_c ( "PHOEBE", et,    "IAU_PHOEBE", "LT+S",
                    "CASSINI",  state, <ime           );
 
         /*
         Display the results.
         */
         printf ( "   Apparent state of Phoebe as seen "
                  "from CASSINI in the IAU_PHOEBE\n"    );
         printf ( "      body-fixed frame (km, km/s) "
                  "obtained\n"                          );
         printf ( "      using spkezr_c directly:\n"    );
         printf ( "      X = %19.6f\n", state[0]        );
         printf ( "      Y = %19.6f\n", state[1]        );
         printf ( "      Z = %19.6f\n", state[2]        );
         printf ( "     VX = %19.6f\n", state[3]        );
         printf ( "     VY = %19.6f\n", state[4]        );
         printf ( "     VZ = %19.6f\n", state[5]        );
 
         /*
         Now we are to compute the angular separation between
         the apparent position of the Earth as seen from the
         orbiter and the nominal boresight of the high gain
         antenna.  First, compute the apparent position of
         the Earth as seen from CASSINI in the J2000 frame.
         */
         spkpos_c ( "EARTH", et,  "J2000", "LT+S",
                    "CASSINI",   pos, <ime            );
 
         /*
         Now compute the location of the antenna boresight
         at this same epoch.  From reading the frame kernel
         we know that the antenna boresight is nominally the
         +Z axis of the CASSINI_HGA frame defined there.
         */
         bsight[0] = 0.0;
         bsight[1] = 0.0;
         bsight[2] = 1.0;
 
         /*
         Now compute the rotation matrix from CASSINI_HGA into
         J2000.
         */
         pxform_c ( "CASSINI_HGA", "J2000", et, pform  );
 
         /*
         And multiply the result to obtain the nominal
         antenna boresight in the J2000 reference frame.
         */
         mxv_c ( pform, bsight, bsight );
 
         /*
         Lastly compute the angular separation.
         */
         convrt_c ( vsep_c(bsight, pos), "RADIANS",
                    "DEGREES",           &sep          );
 
         printf ( "   Angular separation between the "
                  "apparent position of\n"             );
         printf ( "      Earth and the CASSINI high "
                  "gain antenna boresight (degrees):\n");
         printf ( "      %16.3f\n", sep                );
 
         /*
         Or alternatively we can work in the antenna
         frame directly.
         */
         spkpos_c ( "EARTH", et,  "CASSINI_HGA", "LT+S",
                    "CASSINI",   pos, <ime              );
 
         /*
         The antenna boresight is the Z-axis in the
         CASSINI_HGA frame.
         */
         bsight[0] = 0.0;
         bsight[1] = 0.0;
         bsight[2] = 1.0;
 
         /*
         Lastly compute the angular separation.
         */
         convrt_c ( vsep_c(bsight, pos), "RADIANS",
                    "DEGREES",           &sep          );
 
         printf ( "   Angular separation between the "
                  "apparent position of\n"             );
         printf ( "      Earth and the CASSINI high "
                  "gain antenna boresight computed\n"  );
         printf ( "      using vectors in the CASSINI_HGA "
                  "frame (degrees):\n"                 );
         printf ( "      %16.3f\n", sep                );
 
         return(0);
      }


Solution Sample Output



After compiling the program, execute it:

   Converting UTC Time: 2004 jun 11 19:32:00
      ET Seconds Past 2000:    140254384.185
      Apparent state of Phoebe as seen from CASSINI in the IAU_PHOEBE
         body-fixed frame (km, km/s):
         X =        -1982.639762
         Y =         -934.530471
         Z =         -166.562595
        VX =            3.970729
        VY =           -3.812531
        VZ =           -2.371665
      Apparent state of Phoebe as seen from CASSINI in the IAU_PHOEBE
         body-fixed frame (km, km/s) obtained
         using spkezr_c directly:
         X =        -1982.639762
         Y =         -934.530471
         Z =         -166.562595
        VX =            3.970729
        VY =           -3.812531
        VZ =           -2.371665
      Angular separation between the apparent position of
         Earth and the CASSINI high gain antenna boresight (degrees):
                   71.924
      Angular separation between the apparent position of
         Earth and the CASSINI high gain antenna boresight computed
         using vectors in the CASSINI_HGA frame (degrees):
                   71.924


Computing Sub-spacecraft and Sub-solar Points (subpts)







Task Statement




Write a program that prompts the user for an input UTC time string, computes the following quantities at that epoch:

and displays the results. Use the program to compute these quantities at "2004 jun 11 19:32:00" UTC.



Learning Goals




Discover higher level geometry calculation functions in CSPICE and their usage as it relates to CASSINI.



Approach




This particular problem is more of an exercise in searching the permuted index to find the appropriate functions and then reading their headers to understand how to call them.

One point worth considering: Which method do you want to use to compute the sub-solar (or sub-observer) point?



Solution






Solution Meta-Kernel



The meta-kernel we created for the solution to this exercise is named 'subpts.mk'. Its contents follow:

   KPL/MK
 
   This is the meta-kernel used in the solution of the
   ``Computing Sub-spacecraft and Sub-solar Points'' task
   in the Remote Sensing Hands On Lesson.
 
      \begindata
      KERNELS_TO_LOAD = ( 'kernels/lsk/naif0007.tls',
                          'kernels/spk/sat128.bsp'
                          'kernels/spk/981005_PLTEPH-DE405S.bsp',
                          'kernels/spk/020514_SE_SAT105.bsp',
                          'kernels/spk/030201AP_SK_SM546_T45.bsp',
                          'kernels/pck/cpck05Mar2004.tpc' )
      \begintext


Solution Source Code



A sample solution to the problem follows:

      #include <stdio.h>
 
      /*
      Standard CSPICE User Include File
      */
      #include "SpiceUsr.h"
 
      /*
      Local Parameters
      */
 
      #define METAKR "subpts.mk"
      #define STRLEN 50
 
      int main (void)
      {
 
         /*
         Local Variables
         */
         SpiceChar               utctim [STRLEN];
 
         SpiceDouble             alt;
         SpiceDouble             et;
         SpiceDouble             spoint [3];
 
         /*
         Load the kernels that this program requires.  We
         will need:
 
            A leapseconds kernel
            The necessary ephemerides
            A planetary constants file (PCK)
         */
         furnsh_c ( METAKR );
 
         /*
         Prompt the user for the input time string.
         */
         prompt_c ( "Input UTC Time: ", STRLEN, utctim );
 
         printf ( "Converting UTC Time: %s\n", utctim );
 
         /*
         Convert utctim to ET.
         */
         str2et_c ( utctim, &et );
 
         printf ( "   ET Seconds Past 2000: %16.3f\n", et );
 
         /*
         Compute the apparent sub-observer point of CASSINI
         on Phoebe.
         */
         subpt_c ( "NEAR POINT", "PHOEBE", et,  "LT+S",
                   "CASSINI", spoint, &alt                 );
 
         printf ( "   Apparent Sub-Observer point of CASSINI "
                  "on Phoebe in IAU_PHOEBE\n"              );
         printf ( "      (km):\n"                          );
         printf ( "      X = %16.3f\n", spoint[0]          );
         printf ( "      Y = %16.3f\n", spoint[1]          );
         printf ( "      Z = %16.3f\n", spoint[2]          );
         printf ( "    ALT = %16.3f\n", alt                );
 
         /*
         Compute the apparent sub-solar point on Phoebe
         as seen from CASSINI.
         */
         subsol_c ( "NEAR POINT", "PHOEBE", et, "LT+S",
                    "CASSINI",        spoint               );
 
         printf ( "   Apparent Sub-Solar point on Phoebe "
                  "as seen from CASSINI in IAU_PHOEBE\n"   );
         printf ( "      (km):\n"                          );
         printf ( "      X = %16.3f\n", spoint[0]          );
         printf ( "      Y = %16.3f\n", spoint[1]          );
         printf ( "      Z = %16.3f\n", spoint[2]          );
 
         return(0);
      }


Solution Sample Output



After compiling the program, execute it:

   Converting UTC Time: 2004 jun 11 19:32:00
      ET Seconds Past 2000:    140254384.185
      Apparent Sub-Observer point of CASSINI on Phoebe in IAU_PHOEBE
         (km):
         X =          104.498
         Y =           45.269
         Z =            7.383
       ALT =         2084.116
      Apparent Sub-Solar point on Phoebe as seen from CASSINI in IAU_PHOEBE
         (km):
         X =           78.681
         Y =           76.879
         Z =          -21.885


Intersecting Vectors with a Triaxial Ellipsoid (fovint)







Task Statement




Write a program that prompts the user for an input UTC time string and computes the intersection of the CASSINI ISS NAC camera boresight with the surface of Phoebe and presents it in the following coordinates:

If this intersection is found, the program displays the results of the above computations, otherwise it indicates no intersection has occurred. Use this program to compute values at the following epochs:



Learning Goals




Understand how field of view parameters are retrieved from instrument kernels. Learn how various standard planetary constants are retrieved from text PCKs. Discover how to compute the intersection of field of view vectors with triaxial ellipsoidal target bodies.



Approach




This problem can be broken down into several simple, small steps:

It may be useful to consult the CASSINI ISS instrument kernel to determine the name of the NAC camera as well as its configuration. This exercise may make use of some of the concepts and (loosely) code from the ``Spacecraft Orientation and Reference Frames'' task.



Solution






Solution Meta-Kernel



The meta-kernel we created for the solution to this exercise is named 'fovint.mk'. Its contents follow:

   KPL/MK
 
   This is the meta-kernel used in the solution of the
   ``Intersecting Vectors with a Triaxial Ellipsoid'' task
   in the Remote Sensing Hands On Lesson.
 
      \begindata
      KERNELS_TO_LOAD = ( 'kernels/lsk/naif0007.tls',
                          'kernels/sclk/cas00084.tsc',
                          'kernels/spk/sat128.bsp'
                          'kernels/spk/981005_PLTEPH-DE405S.bsp',
                          'kernels/spk/020514_SE_SAT105.bsp',
                          'kernels/spk/030201AP_SK_SM546_T45.bsp',
                          'kernels/fk/cas_v37.tf',
                          'kernels/ck/04135_04171pc_psiv2.bc',
                          'kernels/pck/cpck05Mar2004.tpc',
                          'kernels/ik/cas_iss_v09.ti' )
      \begintext


Solution Source Code



A sample solution to the problem follows:

      #include <stdio.h>
 
      /*
      Standard CSPICE User Include File
      */
      #include "SpiceUsr.h"
      #include <stdlib.h>
 
      /*
      Local Parameters
      */
 
      #define METAKR "fovint.mk"
      #define STRLEN 50
      #define BCVLEN 4
 
      int main (void)
      {
 
         /*
         Local Variables
         */
         SpiceChar               frame  [STRLEN];
         SpiceChar               shape  [STRLEN];
         SpiceChar               utctim [STRLEN];
 
         SpiceDouble             lat;
         SpiceDouble             lon;
         SpiceDouble             bounds [BCVLEN][3];
         SpiceDouble             bsight [3];
         SpiceDouble             dist;
         SpiceDouble             et;
         SpiceDouble             obspos [3];
         SpiceDouble             point  [3];
         SpiceDouble             radius;
         SpiceDouble             trgepc;
 
         SpiceInt                n;
         SpiceInt                nacid;
 
         SpiceBoolean            found;
 
         /*
         Load the kernels that this program requires.  We
         will need:
 
            A leapseconds kernel.
            A SCLK kernel for CASSINI.
            Any necessary ephemerides.
            The CASSINI frame kernel.
            A CASSINI C-kernel.
            A PCK file with Phoebe constants.
            The CASSINI ISS I-kernel.
         */
         furnsh_c ( METAKR );
 
         /*
         Prompt the user for the input time string.
         */
         prompt_c ( "Input UTC Time: ", STRLEN, utctim );
 
         printf ( "Converting UTC Time: %s\n", utctim );
 
         /*
         Convert utctim to ET.
         */
         str2et_c ( utctim, &et );
 
         printf ( "   ET Seconds Past 2000: %16.3f\n", et );
 
         /*
         Now we need to obtain the FOV configuration of
         the ISS NAC camera. To do this we will need the
         ID code for CASSINI_ISS_NAC.
         */
         bodn2c_c ( "CASSINI_ISS_NAC", &nacid, &found );
 
         /*
         Stop the program if the code was not found.
         */
         if ( !found )
         {
            printf ( "Unable to locate the ID code for "
                     "CASSINI_ISS_NAC.\n"                  );
            exit   ( EXIT_FAILURE );
         }
 
         /*
         Now retrieve the field of view parameters.
         */
         getfov_c ( nacid, BCVLEN, STRLEN, STRLEN,
                    shape, frame,  bsight, &n,      bounds );
 
         /*
         Call srfxpt_c to determine coordinates of the
         intersection of this vector with the surface
         of Phoebe.
         */
         srfxpt_c ( "Ellipsoid",
                    "PHOEBE",  et, "LT+S",
                    "CASSINI", frame, bsight,
                    point,  &dist, &trgepc, obspos, &found );
 
         /*
         Check the found flag.  Display a message if the
         point of intersection was not found and stop.
         */
         if ( !found )
         {
            printf ( "No intersection point found at "
                     "this epoch.\n"                   );
            exit   ( EXIT_SUCCESS );
         }
 
         /*
         Now, we have discovered a point of intersection.
         Start by displaying the position vector in the
         IAU_PHOEBE frame of the intersection.
         */
         printf ( "   Position vector of CASSINI NAC camera "
                  "boresight surface intercept\n"         );
         printf ( "      in the IAU_PHOEBE frame "
                  "(km):\n"                               );
         printf ( "      X = %16.3f\n", point[0]          );
         printf ( "      Y = %16.3f\n", point[1]          );
         printf ( "      Z = %16.3f\n", point[2]          );
 
         /*
         Now express the coordinates of this point in
         planetocentric latitude and longitude.
         */
         reclat_c ( point, &radius, &lon, &lat );
 
         /*
         Convert the angles to degrees for displaying.
         */
         printf ( "   Planetocentric coordinates of the "
                  "intercept (degrees):\n"              );
         printf ( "      LAT = %16.3f\n", lat * dpr_c() );
         printf ( "      LON = %16.3f\n", lon * dpr_c() );
 
         return(0);
      }


Solution Sample Output



After compiling the program, execute it:

   Converting UTC Time: 2004 jun 11 19:32:00
      ET Seconds Past 2000:    140254384.185
      Position vector of CASSINI NAC camera boresight surface intercept
         in the IAU_PHOEBE frame (km):
         X =           86.390
         Y =           72.089
         Z =            8.255
      Planetocentric coordinates of the intercept (degrees):
         LAT =            4.196
         LON =           39.844


Computing Illumination Angles and Local Time (angles)







Task Statement




Write a program that prompts the user for an input time string and computes the intersection of the CASSINI NAC camera boresight and field of view boundary vectors with the surface of Phoebe. At these points of intersection, if they exist, compute the following:

Additionally compute the local solar time at the intercept of the camera boresight with the surface of Phoebe.

Display the results of the above computations if an intersection occurs, otherwise indicate the absence of an intersection. Use this program to compute values at the epoch "2004-01-12T4:15.24.000" UTC.



Learning Goals




Discover another high level geometry function and another time conversion function in CSPICE. Reinforce the concepts introduced in the previous task.



Approach




Making use of the code you wrote for the previous task is probably the fastest means to an end. A significant percentage of the task is devoted to similar computations.

This problem can be broken down into several steps:

For each vector in the set of boundary corner vectors, and for the boresight vector, perform the following operations:

At this point, if a boresight intercept was located, then proceed.



Solution






Solution Meta-Kernel



The meta-kernel we created for the solution to this exercise is named 'angles.mk'. Its contents follow:

   KPL/MK
   This is the meta-kernel used in the solution of the
   ``Computing Illumination Angles and Local Time'' task
   in the Remote Sensing Hands On Lesson.
 
      \begindata
      KERNELS_TO_LOAD = ( 'kernels/lsk/naif0007.tls',
                          'kernels/sclk/cas00084.tsc',
                          'kernels/spk/sat128.bsp'
                          'kernels/spk/981005_PLTEPH-DE405S.bsp',
                          'kernels/spk/020514_SE_SAT105.bsp',
                          'kernels/spk/030201AP_SK_SM546_T45.bsp',
                          'kernels/fk/cas_v37.tf',
                          'kernels/ck/04135_04171pc_psiv2.bc',
                          'kernels/pck/cpck05Mar2004.tpc',
                          'kernels/ik/cas_iss_v09.ti' )
      \begintext


Solution Source Code



A sample solution to the problem follows:

      #include <stdio.h>
 
      /*
      Standard CSPICE User Include File
      */
      #include "SpiceUsr.h"
      #include <stdlib.h>
 
      /*
      Local Parameters
      */
 
      #define METAKR "angles.mk"
      #define STRLEN 50
      #define BCVLEN 5
 
      int main (void)
      {
 
         /*
         Local Variables
         */
         SpiceChar               ampm   [STRLEN];
         SpiceChar               frame  [STRLEN];
         SpiceChar               shape  [STRLEN];
         SpiceChar               time   [STRLEN];
         SpiceChar               utctim [STRLEN];
         SpiceChar               *vecnam[] = {
                                   "Boundary Corner 1",
                                   "Boundary Corner 2",
                                   "Boundary Corner 3",
                                   "Boundary Corner 4",
                                   "Boresight"
         };
 
         SpiceDouble             lat;
         SpiceDouble             lon;
         SpiceDouble             bounds [BCVLEN][3];
         SpiceDouble             bsight [3];
         SpiceDouble             dist;
         SpiceDouble             emissn;
         SpiceDouble             et;
         SpiceDouble             obspos [3];
         SpiceDouble             phase;
         SpiceDouble             point  [3];
         SpiceDouble             radius;
         SpiceDouble             solar;
         SpiceDouble             trgepc;
 
 
         SpiceInt                hr;
         SpiceInt                i;
         SpiceInt                phoeid;
         SpiceInt                mn;
         SpiceInt                n;
         SpiceInt                sc;
         SpiceInt                nacid;
 
         SpiceBoolean            found;
 
         /*
         Load the kernels that this program requires.  We
         will need:
 
            A leapseconds kernel.
            A SCLK kernel for CASSINI.
            Any necessary ephemerides.
            The CASSINI frame kernel.
            A CASSINI C-kernel.
            A PCK file with Phoebe constants.
            The CASSINI ISS I-kernel.
         */
         furnsh_c ( METAKR );
 
         /*
         Prompt the user for the input time string.
         */
         prompt_c ( "Input UTC Time: ", STRLEN, utctim );
 
         printf ( "Converting UTC Time: %s\n", utctim );
 
         /*
         Convert utctim to ET.
         */
         str2et_c ( utctim, &et );
 
         printf ( "   ET Seconds Past 2000: %16.3f\n", et );
 
         /*
         Now we need to obtain the FOV configuration of
         the ISS NAC camera.  To do this we will need the
         ID code for CASSINI_ISS_NAC.
         */
         bodn2c_c ( "CASSINI_ISS_NAC", &nacid, &found );
 
         /*
         Stop the program if the code was not found.
         */
         if ( !found )
         {
            printf ( "Unable to locate the ID code for "
                     "CASSINI_ISS_NAC\n"                   );
            exit   ( EXIT_FAILURE );
         }
 
         /*
         Now retrieve the field of view parameters.
         */
         getfov_c ( nacid, BCVLEN, STRLEN, STRLEN,
                    shape, frame,  bsight, &n,      bounds );
 
         /*
         Rather than treat BSIGHT as a separate vector,
         copy it into the last slot of BOUNDS.
         */
         for ( i=0; i<3; i++ )
         {
            bounds[4][i] = bsight[i];
         }
 
         /*
         Now perform the same set of calculations for each
         vector listed in the BOUNDS array.
         */
         for ( i=0; i<5; i++ )
         {
 
            /*
            Call srfxpt_c to determine coordinates of the
            intersection of this vector with the surface
            of Phoebe.
            */
            srfxpt_c ( "Ellipsoid",
                       "PHOEBE",  et, "LT+S",
                       "CASSINI", frame, bounds[i],
                       point,  &dist, &trgepc, obspos, &found );
 
            /*
            Check the found flag.  Display a message if
            the point of intersection was not found,
            otherwise continue with the calculations.
            */
            printf ( "Vector: %s\n", vecnam[i] );
 
            if ( !found )
            {
               printf ( "No intersection point found at "
                        "this epoch for this vector.\n"   );
            }
            else
            {
               /*
               Display the planetocentric latitude and longitude
               of the intercept.
               */
               reclat_c ( point, &radius, &lon, &lat );
 
               printf ( "   Planetocentric coordinates of "
                        "the intercept (degrees):\n"  );
               printf ( "    LAT = %16.3f\n", lat * dpr_c() );
               printf ( "    LON = %16.3f\n", lon * dpr_c() );
 
               /*
               Compute the illumination angles at this
               point.
               */
               illum_c ( "PHOEBE", et,    "LT+S", "CASSINI",
                         point, &phase, &solar, &emissn  );
 
               printf ( "   Phase angle (degrees):           "
                        "%16.3f\n", phase * dpr_c()          );
               printf ( "   Solar incidence angle (degrees): "
                        "%16.3f\n", solar * dpr_c()          );
               printf ( "   Emission angle (degrees):        "
                        "%16.3f\n", emissn * dpr_c()         );
 
            }
 
            printf ( "\n" );
 
         }
 
         /*
         Lastly compute the local solar time at the
         boresight intersection.
         */
         if ( found )
         {
 
            /*
            Get ID code of Phoebe.
            */
            bodn2c_c ( "PHOEBE", &phoeid, &found );
 
            /*
            The ID code for PHOEBE is built-in to the library.
            However, it is good programming practice to get
            in the habit of checking your found-flags.
            */
            if ( !found )
            {
               printf ( "Unable to locate the body ID code "
                        "for Phoebe.\n"                      );
               exit   ( EXIT_FAILURE );
            }
 
            /*
            Call et2lst_c to compute local time.
            */
            et2lst_c ( et,
                       phoeid,
                       lon,
                       "PLANETOCENTRIC",
                       STRLEN,
                       STRLEN,
                       &hr,
                       &mn,
                       &sc,
                       time,
                       ampm             );
 
            printf ( "   Local Solar Time at boresight "
                     "intercept (24 Hour Clock):\n"      );
            printf ( "      %s\n", time                  );
         }
         else
         {
            printf ( "   No boresight intercept to "
                     "compute local solar time.\n"   );
         }
 
         return(0);
      }


Solution Sample Output



After compiling the program, execute it:

   Converting UTC Time: 2004 jun 11 19:32:00
      ET Seconds Past 2000:    140254384.185
   Vector: Boundary Corner 1
      Planetocentric coordinates of the intercept (degrees):
       LAT =            1.028
       LON =           36.433
      Phase angle (degrees):                     28.110
      Solar incidence angle (degrees):           16.121
      Emission angle (degrees):                  14.627
 
   Vector: Boundary Corner 2
      Planetocentric coordinates of the intercept (degrees):
       LAT =            7.492
       LON =           36.556
      Phase angle (degrees):                     27.894
      Solar incidence angle (degrees):           22.894
      Emission angle (degrees):                  14.988
 
   Vector: Boundary Corner 3
      Planetocentric coordinates of the intercept (degrees):
       LAT =            7.373
       LON =           43.430
      Phase angle (degrees):                     28.171
      Solar incidence angle (degrees):           21.315
      Emission angle (degrees):                  21.977
 
   Vector: Boundary Corner 4
      Planetocentric coordinates of the intercept (degrees):
       LAT =            0.865
       LON =           43.239
      Phase angle (degrees):                     28.385
      Solar incidence angle (degrees):           13.882
      Emission angle (degrees):                  21.763
 
   Vector: Boresight
      Planetocentric coordinates of the intercept (degrees):
       LAT =            4.196
       LON =           39.844
      Phase angle (degrees):                     28.140
      Solar incidence angle (degrees):           18.247
      Emission angle (degrees):                  17.858
 
      Local Solar Time at boresight intercept (24 Hour Clock):
         11:31:50