Main Page
Remote Sensing Hands-On Lesson, using MPO (C)

Table of Contents


   Remote Sensing Hands-On Lesson, using MPO (C)
      Overview
      Note About HTML Links
      References
         Tutorials
         Required Readings
         The Permuted Index
         API Documentation
      Kernels Used
      CSPICE Modules Used

   Time Conversion (convtm)
      Task Statement
      Learning Goals
      Approach
      Solution
         Solution Meta-Kernel
         Solution Source Code
         Solution Sample Output
      Extra Credit
         Task statements and questions
         Solutions and answers

   Obtaining Target States and Positions (getsta)
      Task Statement
      Learning Goals
      Approach
      Solution
         Solution Meta-Kernel
         Solution Source Code
         Solution Sample Output
      Extra Credit
         Task statements and questions
         Solutions and answers

   Spacecraft Orientation and Reference Frames (xform)
      Task Statement
      Learning Goals
      Approach
      Solution
         Solution Meta-Kernel
         Solution Source Code
         Solution Sample Output
      Extra Credit
         Task statements and questions
         Solutions and answers

   Computing Sub-s/c and Sub-solar Points on an Ellipsoid and a DSK (subpts)
      Task Statement
      Learning Goals
      Approach
      Solution
         Solution Meta-Kernel
         Solution Source Code
         Solution Sample Output
      Extra Credit
         Task statements and questions
         Solutions and answers

   Intersecting Vectors with an Ellipsoid and a DSK (fovint)
      Task Statement
      Learning Goals
      Approach
      Solution
         Solution Meta-Kernel
         Solution Source Code
         Solution Sample Output
      Extra Credit




Top

Remote Sensing Hands-On Lesson, using MPO (C)





March 01, 2023



Top

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 spectrometer flying on the BepiColombo MPO spacecraft, but many of the concepts are easily extended and generalized to other scenarios.



Top

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, if not done already by installing the lessons package under the Toolkit's ``doc/html'' directory, create a subdirectory called ``lessons'' under the ``doc/html'' directory of the ``cspice/'' tree and copy this document to that subdirectory before loading it into a Web browser.



Top

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 ``API Documentation''.

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 the 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.



Top

Tutorials



The following SPICE tutorials serve as references for 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
   DSK               Detailed Target Shape (Topography) Data
These tutorials are available from the NAIF server at JPL:

   https://naif.jpl.nasa.gov/naif/tutorials.html


Top

Required Readings



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

   Name             Lesson steps/functions that it describes
   ---------------  -----------------------------------------
   ck.req           Obtaining spacecraft orientation data
   dsk.req          Obtaining detailed body shape data
   frames.req       Using reference frames
   naif_ids.req     Determining body ID codes
   pck.req          Obtaining planetary constants data
   sclk.req         SCLK time conversion
   spk.req          Obtaining ephemeris data
   time.req         Time conversion


Top

The Permuted Index



Another useful document distributed with the Toolkit is the permuted index. It is located under the ``cspice/doc'' directory in the C installation tree.

This text document provides a simple mechanism by which users can discover which CSPICE functions perform functions of interest, as well as the names of the source files that contain these functions.



Top

API Documentation



The most detailed specification of a given SPICE C routine is contained in the header section of its source code. The source code is distributed with the Toolkit and is located under the ``cspice/src/cspice'' path.

For example the path of the source code of the str2et_c routine is

   cspice/src/cspice/str2et_c.c


Top

Kernels Used




The following kernels are used in examples provided in this lesson:

   1.  Generic LSK:
 
          naif0012.tls
 
   2.  BepiColombo MPO SCLK:
 
          bc_mpo_step_20230117.tsc
 
   3.  Solar System Ephemeris SPK, subsetted to cover only the time
       range of interest:
 
          de432s.bsp
 
   4.  BepiColombo MPO Spacecraft Trajectory SPK, subsetted to cover
       only the time range of interest:
 
          bc_mpo_mlt_50037_20260314_20280529_v05.bsp
 
   5.  BepiColombo MPO FK:
 
          bc_mpo_v32.tf
 
   6.  BepiColombo MPO Spacecraft CK, subsetted to cover only the time
       range of interest:
 
          bc_mpo_sc_slt_50028_20260314_20280529_f20181127_v03.bc
 
   7.  Generic PCK:
 
          pck00011.tpc
 
   8.  Low-resolution Mercury DSK:
 
          mercury_lowres.bds
 
   9.  SIMBIO-SYS IK:
 
          bc_mpo_simbio-sys_v08.ti
 
These SPICE kernels are included in the lesson package.

In addition to these kernels, the extra credit exercises require the following kernels:

   #  FILE NAME       TYPE DESCRIPTION
   -- --------------- ---- ---------------------------------------------
   10 jup365_2027.bsp SPK  Generic Jovian Satellite Ephemeris SPK
These SPICE kernels are available from the NAIF server at JPL, in the ``satellites/a_old_versions'' subdurectory:

   https://naif.jpl.nasa.gov/pub/naif/generic_kernels/spk/


Top

CSPICE Modules Used




This section provides a complete list of the functions and kernels that are suggested for usage in each of the exercises in this lesson. (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
                      sce2s_c
 
           extra (*)  unload_c   unitim_c   1,2
                      sct2e_c
                      et2utc_c
                      scs2e_c
 
      2    getsta     furnsh_c   vnorm_c    1,3,4
                      prompt_c
                      str2et_c
                      spkezr_c
                      spkpos_c
                      convrt_c
 
           extra (*)  kclear_c              1,4,10
                      unload_c
 
      3    xform      furnsh_c   vsep_c     1-7
                      prompt_c
                      str2et_c
                      spkezr_c
                      sxform_c
                      mxvg_c
                      spkpos_c
                      pxform_c
                      mxv_c
                      convrt_c
 
           extra (*)  kclear_c              1-7
                      unload_c
 
      4    subpts     furnsh_c   vnorm_c    1,3-4,7,8
                      prompt_c
                      str2et_c
                      subpnt_c
                      subslr_c
 
           extra (*)  kclear_c   dpr_c      1,3-4,7,10
                      reclat_c
                      bodvrd_c
                      recpgr_c
 
      5    fovint     furnsh_c   dpr_c      1-9
                      prompt_c
                      str2et_c
                      getfvn_c
                      bodn2c_c
                      sincpt_c
                      reclat_c
                      illumf_c
                      et2lst_c
 
 
      (*) Additional APIs and kernels used in Extra Credit tasks.
Refer to the headers of the various functions listed above, as detailed interface specifications are provided with the source code.



Top

Time Conversion (convtm)







Top

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:

    1. Ephemeris Time (ET) in seconds past J2000

    2. Calendar Ephemeris Time

    3. Spacecraft Clock Time

and displays the results. Use the program to convert "2027 JAN 05 02:04:36" UTC into these alternate systems.



Top

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.



Top

Approach




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

    -- Decide which SPICE kernels are necessary. Prepare a meta-kernel listing the kernels and load it into the program.

    -- Prompt the user for an input UTC time string.

    -- Convert the input time string into ephemeris time expressed as seconds past J2000 TDB. Display the result.

    -- Convert ephemeris time into a calendar format. Display the result.

    -- Convert ephemeris time into a spacecraft clock string. Display the result.

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

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



Top

Solution






Top

Solution Meta-Kernel



The meta-kernel we created for the solution to this exercise is named 'convtm.tm'. 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.
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
         1. Generic LSK:
 
               naif0012.tls
 
         2. BepiColombo MPO SCLK:
 
               bc_mpo_step_20230117.tsc
 
 
   \begindata
 
    KERNELS_TO_LOAD = (
 
    'kernels/lsk/naif0012.tls',
    'kernels/sclk/bc_mpo_step_20230117.tsc'
 
                      )
 
   \begintext


Top

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.tm"
      #define SCLKID -121
      #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 J2000: %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);
      }


Top

Solution Sample Output



After compiling the program, execute it:

   Input UTC Time: 2027 JAN 05 02:04:36
   Converting UTC Time: 2027 JAN 05 02:04:36
      ET Seconds Past J2000:    852386745.184
      Calendar ET (etcal_c): 2027 JAN 05 02:05:45.184
      Calendar ET (timout_c): 2027-JAN-05T02:05:45
      Spacecraft Clock Time: 1/0863834674:28127


Top

Extra Credit




In this ``extra credit'' section you will be presented with more complex tasks, aimed at improving your understanding of time conversions, the Toolkit routines that deal with them, and some common errors that may happen during the execution of these conversions.

These ``extra credit'' tasks are provided as task statements, and unlike the regular tasks, no approach or solution source code is provided. In the next section, you will find the numeric solutions (when applicable) and answers to the questions asked in these tasks.



Top

Task statements and questions



    1. Extend your program to convert the input UTC time string to TDB Julian Date. Convert "2027 JAN 05 02:04:36" UTC.

    2. Remove the LSK from the original meta-kernel and run your program again, using the same inputs as before. Has anything changed? Why?

    3. Remove the SCLK from the original meta-kernel and run your program again, using the same inputs as before. Has anything changed? Why?

    4. Modify your program to perform conversion of UTC or ephemeris time, to a spacecraft clock string using the NAIF ID for the BepiColombo MPO SIMBIO-SYS HRIC channel. Convert "2027 JAN 05 02:04:36" UTC.

    5. Find the earliest UTC time that can be converted to BepiColombo MPO spacecraft clock.

    6. Extend your program to convert the spacecraft clock time obtained in the regular task back to UTC Time and present it in ISO calendar date format, with a resolution of milliseconds.

    7. Examine the contents of the generic LSK and the BepiColombo MPO SCLK kernels. Can you understand and explain what you see?



Top

Solutions and answers



    1. Two methods exist in order to convert ephemeris time to Julian Date: unitim_c and timout_c. The difference between them is the type of output produced by each method. unitim_c returns the double precision value of an input epoch, while timout_c returns the string representation of the ephemeris time in Julian Date format (when picture input is set to 'JULIAND.######### ::TDB'). Refer to the function header for further details. The solution for the requested input UTC string is:

      Julian Date TDB:   2461410.5873285
    2. When running the original program without the LSK kernel, an error is produced:

 
   =====================================================================
   ===========
 
   Toolkit version: N0067
 
   SPICE(NOLEAPSECONDS) --
 
   The variable that points to the leapseconds (DELTET/DELTA_AT) could n
   ot be
   located in the kernel pool. It is likely that the leapseconds kernel
   has not
   been loaded.
 
   A traceback follows.  The name of the highest level module is first.
   str2et_c --> STR2ET --> TTRANS
 
   Oh, by the way:  The SPICELIB error handling actions are USER-TAILORA
   BLE.  You
   can choose whether the Toolkit aborts or continues when errors occur,
    which
   error messages to output, and where to send the output.  Please read
   the ERROR
   "Required Reading" file, or see the routines ERRACT, ERRDEV, and ERRP
   RT.
 
   =====================================================================
   ===========
    This error is triggered by str2et_c because the variable that points to the leapseconds is not present in the kernel pool and therefore the program lacks data required to perform the requested UTC to ephemeris time conversion.

    By default, SPICE will report, as a minimum, a short descriptive message and a expanded form of this short message where more details about the error are provided. If this error message is not sufficient for you to understand what has happened, you could go to the ``Exceptions'' section in the SPICELIB or CSPICE headers of the function that has triggered the error and find out more information about the possible causes.

    3. When running the original program without the SCLK kernel, an error is produced:

 
   =====================================================================
   ===========
 
   Toolkit version: N0067
 
   SPICE(KERNELVARNOTFOUND) -- The Variable Was not Found in the Kernel
   Pool.
 
   Kernel variable SCLK_DATA_TYPE_121 was not found in the kernel pool.
 
   A traceback follows.  The name of the highest level module is first.
   sce2s_c --> SCE2S --> SCE2T --> SCTYPE --> SCTY01
 
   Oh, by the way:  The SPICELIB error handling actions are USER-TAILORA
   BLE.  You
   can choose whether the Toolkit aborts or continues when errors occur,
    which
   error messages to output, and where to send the output.  Please read
   the ERROR
   "Required Reading" file, or see the routines ERRACT, ERRDEV, and ERRP
   RT.
 
   =====================================================================
   ===========
    This error is triggered by sce2s_c. In this case the error message may not give you enough information to understand what has actually happened. Nevertheless, the expanded form of this short message clearly indicates that the SCLK kernel for the spacecraft ID -121 has not been loaded.

    The UTC string to ephemeris time conversion and the conversion of ephemeris time into a calendar format worked normally as these conversions only require the LSK kernel to be loaded.

    4. The first thing you need to do is to find out what the NAIF ID is for the SIMBIO-SYS HRIC channel. In order to do so, examine the BepiColombo MPO frames definitions kernel listed above and look for the ``BepiColombo MPO Mission NAIF ID Codes'' or for the ``BepiColombo MPO NAIF ID Codes to Name Mapping'' and there, for the NAIF ID given to MPO_SIMBIO-SYS_HRIC_FPA (which is -121610). Then replace in your code the SCLK ID -121 with -121610. After compiling and executing the program using the original meta-kernel, you will be getting the same error as in the previous task. Despite the error being exactly the same, this case is different. Generally, spacecraft clocks are associated with the spacecraft ID and not with its payload, sensors or structures IDs. Therefore, in order to do conversions from/to spacecraft clock for payload, sensors or spacecraft structures, the spacecraft ID must be used.

    Note that this does not need to be true for all missions or payloads, as SPICE does not restrict the SCLKs to spacecraft IDs only. Please refer to your mission's SCLK kernels for particulars.

    5. Use sct2e_c with the encoding of the MPO spacecraft clock time set to 0.0 ticks and convert the resulting ephemeris time to UTC using either timout_c or et2utc_c. The solution for the requested SCLK string is:

      Earliest UTC convertible to SCLK: 1999-08-22T00:00:05.204
    6. Use scs2e_c with the SCLK string obtained in the computations performed in the regular tasks and convert the resulting ephemeris time to UTC using either et2utc_c, with 'ISOC' format and 3 digits precision, or using timout_c using the time picture 'YYYY-MM-DDTHR:MN:SC.### ::RND'. The solution of the requested conversion is:

      Spacecraft Clock Time:          1/0863834674:28127
      UTC time from spacecraft clock: 2027-01-05T02:04:36.000


Top

Obtaining Target States and Positions (getsta)







Top

Task Statement




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

    1. The apparent state of Mercury as seen from BepiColombo MPO in the J2000 frame, in kilometers and kilometers/second. This vector itself is not of any particular interest, but it is a useful intermediate quantity in some geometry calculations.

    2. The apparent position of the Earth as seen from BepiColombo MPO in the J2000 frame, in kilometers.

    3. The one-way light time between BepiColombo MPO and the apparent position of Earth, in seconds.

    4. The apparent position of the Sun as seen from Mercury in the J2000 frame (J2000), in kilometers.

    5. The actual (geometric) distance between the Sun and Mercury, in astronomical units.

and displays the results. Use the program to compute these quantities at "2027 JAN 05 02:04:36" UTC.



Top

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.



Top

Approach




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

    -- Decide which SPICE kernels are necessary. Prepare a meta-kernel listing the kernels and load it into the program.

    -- Prompt the user for an input time string.

    -- Convert the input time string into ephemeris time expressed as seconds past J2000 TDB.

    -- Compute the state of Mercury relative to BepiColombo MPO in the J2000 reference frame, corrected for aberrations.

    -- Compute the position of Earth relative to BepiColombo MPO in the J2000 reference frame, corrected for aberrations. (The function in the library that computes this also returns the one-way light time between BepiColombo MPO and Earth.)

    -- Compute the position of the Sun relative to Mercury in the J2000 reference frame, corrected for aberrations.

    -- Compute the position of the Sun relative to Mercury without correcting for aberration.

    Compute the length of this vector. This provides the desired distance in kilometers.

    -- Convert the distance in kilometers into AU.

You may find it useful to consult the permuted index, the headers of various source modules, and the ``SPK Required Reading'' (spk.req) 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.



Top

Solution






Top

Solution Meta-Kernel



The meta-kernel we created for the solution to this exercise is named 'getsta.tm'. 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.
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
         1. Generic LSK:
 
               naif0012.tls
 
         2. Solar System Ephemeris SPK, subsetted to cover only
            the time range of interest:
 
               de432s.bsp
 
         3. BepiColombo MPO Spacecraft Trajectory SPK, subsetted
            to cover only the time range of interest:
 
               bc_mpo_mlt_50037_20260314_20280529_v05.bsp
 
 
   \begindata
 
    KERNELS_TO_LOAD = (
 
    'kernels/lsk/naif0012.tls',
    'kernels/spk/de432s.bsp',
    'kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp',
 
                        )
 
   \begintext


Top

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.tm"
   #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 J2000: %16.3f\n", et );
 
      /*
      Compute the apparent state of Mercury as seen from
      BepiColombo MPO in the J2000 frame.  All of the ephemeris
      readers return states in units of kilometers and
      kilometers per second.
      */
      spkezr_c ( "MERCURY", et,    "J2000", "LT+S",
                 "MPO",     state, &ltime              );
 
      printf ( "   Apparent state of Mercury as seen "
               "from BepiColombo MPO in the\n"         );
      printf ( "      J2000 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
      BepiColombo MPO 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",
                 "MPO",   pos, &ltime                  );
 
      printf ( "   Apparent position of Earth as seen "
               "from BepiColombo MPO 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]            );
 
 
      /*
      We need only display LTIME, as it is precisely the
      light time in which we are interested.
      */
      printf ( "   One way light time between BepiColombo "
               "MPO and the apparent\n"                );
      printf ( "      position of Earth (seconds): "
               "%16.3f\n", ltime                       );
 
      /*
      Compute the apparent position of the Sun as seen
      from Mercury in the J2000 frame.
      */
      spkpos_c ( "SUN",     et,  "J2000", "LT+S",
                 "MERCURY", pos, &ltime               );
 
      printf ( "   Apparent position of Sun as seen "
               "from Mercury 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 Mercury.  The above SPKPOS call gives us
      the apparent distance, so we need to adjust our
      aberration correction appropriately.
      */
      spkpos_c ( "SUN",     et,  "J2000", "NONE",
                 "MERCURY", pos, &ltime               );
 
      /*
      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 "
               "Mercury body centers:\n"             );
      printf ( "      (AU): %16.3f\n", dist          );
 
      return(0);
   }


Top

Solution Sample Output



After compiling the program, execute it:

   Input UTC Time: 2027 JAN 05 02:04:36
   Converting UTC Time: 2027 JAN 05 02:04:36
      ET seconds past J2000:    852386745.184
      Apparent state of Mercury as seen from BepiColombo MPO in the
         J2000 frame (km, km/s):
         X =         -683.207
         Y =        -1438.946
         Z =        -2427.819
        VX =            0.036
        VY =            2.360
        VZ =           -1.783
      Apparent position of Earth as seen from BepiColombo MPO in the
         J2000 frame (km):
         X =    -59257854.691
         Y =    185201786.218
         Z =     88178321.179
      One way light time between BepiColombo MPO and the apparent
         position of Earth (seconds):          712.193
      Apparent position of Sun as seen from Mercury in the
         J2000 frame (km):
         X =    -23429947.239
         Y =     54297427.572
         Z =     31434173.468
      Actual distance between Sun and Mercury body centers:
         (AU):            0.448


Top

Extra Credit




In this ``extra credit'' section you will be presented with more complex tasks, aimed at improving your understanding of state computations, particularly the application of the different light time and stellar aberration corrections available in the spkezr_c function, and some common errors that may happen when computing these states.

These ``extra credit'' tasks are provided as task statements, and unlike the regular tasks, no approach or solution source code is provided. In the next section, you will find the numeric solutions (when applicable) and answers to the questions asked in these tasks.



Top

Task statements and questions



    1. Remove the planetary ephemerides SPK from the original meta-kernel and run your program again, using the same inputs as before. Has anything changed? Why?

    2. Extend your program to compute the geometric position of Jupiter as seen from Mercury in the J2000 frame (J2000), in kilometers.

    3. Extend, or modify, your program to compute the position of the Sun as seen from Mercury in the J2000 frame (J2000), in kilometers, using the following light time and aberration corrections: NONE, LT and LT+S. Explain the differences.

    4. Examine the BepiColombo MPO frames definition kernel to find the SPICE ID/name definitions.



Top

Solutions and answers



    1. When running the original program without the planetary ephemerides SPK, an error is produced by spkezr_c:

 
   =====================================================================
   ===========
 
   Toolkit version: N0067
 
   SPICE(SPKINSUFFDATA) --
 
   Insufficient ephemeris data has been loaded to compute the state of -
   121
   (BEPICOLOMBO MPO) relative to 0 (SOLAR SYSTEM BARYCENTER) at the ephe
   meris
   epoch 2027 JAN 05 02:05:45.184.
 
   A traceback follows.  The name of the highest level module is first.
   spkezr_c --> SPKEZR --> SPKEZ --> SPKACS --> SPKGEO
 
   Oh, by the way:  The SPICELIB error handling actions are USER-TAILORA
   BLE.  You
   can choose whether the Toolkit aborts or continues when errors occur,
    which
   error messages to output, and where to send the output.  Please read
   the ERROR
   "Required Reading" file, or see the routines ERRACT, ERRDEV, and ERRP
   RT.
 
   =====================================================================
   ===========
    This error is generated when trying to compute the apparent state of Mercury as seen from BepiColombo MPO in the J2000 frame because despite the BepiColombo MPO ephemeris data being relative to Mercury, the state of the spacecraft with respect to the solar system barycenter is required to compute the light time and stellar aberrations. The loaded SPK data are enough to compute geometric states of BepiColombo MPO with respect to Mercury center, and geometric states of Mercury barycenter with respect to the Solar System Barycenter, but insufficient to compute the state of the spacecraft relative to the Solar System Barycenter because the SPK data needed to compute geometric states of Mercury center relative to its barycenter are no longer loaded. Run ``brief'' on the SPKs used in the original task to find out which ephemeris objects are available from those kernels. If you want to find out what is the 'center of motion' for the ephemeris object(s) included in an SPK, use the -c option when running ``brief'':

 
   BRIEF -- Version 4.1.0, September 17, 2021 -- Toolkit Version N0067
 
 
   Summary for: kernels/spk/de432s.bsp
 
   Bodies: MERCURY BARYCENTER (1) w.r.t. SOLAR SYSTEM BARYCENTER (0)
           VENUS BARYCENTER (2) w.r.t. SOLAR SYSTEM BARYCENTER (0)
           EARTH BARYCENTER (3) w.r.t. SOLAR SYSTEM BARYCENTER (0)
           MARS BARYCENTER (4) w.r.t. SOLAR SYSTEM BARYCENTER (0)
           JUPITER BARYCENTER (5) w.r.t. SOLAR SYSTEM BARYCENTER (0)
           SATURN BARYCENTER (6) w.r.t. SOLAR SYSTEM BARYCENTER (0)
           URANUS BARYCENTER (7) w.r.t. SOLAR SYSTEM BARYCENTER (0)
           NEPTUNE BARYCENTER (8) w.r.t. SOLAR SYSTEM BARYCENTER (0)
           PLUTO BARYCENTER (9) w.r.t. SOLAR SYSTEM BARYCENTER (0)
           SUN (10) w.r.t. SOLAR SYSTEM BARYCENTER (0)
           MERCURY (199) w.r.t. MERCURY BARYCENTER (1)
           VENUS (299) w.r.t. VENUS BARYCENTER (2)
           MOON (301) w.r.t. EARTH BARYCENTER (3)
           EARTH (399) w.r.t. EARTH BARYCENTER (3)
           Start of Interval (UTC)             End of Interval (UTC)
           -----------------------------       -------------------------
   ----
           2027-JAN-02 23:01:53.350            2027-JAN-08 00:59:37.932
 
 
   Summary for: kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp
 
   Body: BEPICOLOMBO MPO (-121) w.r.t. MERCURY (199)
         Start of Interval (UTC)             End of Interval (UTC)
         -----------------------------       ---------------------------
   --
         2027-JAN-02 23:01:53.350            2027-JAN-08 00:59:37.932
 
   Bodies: -121000 w.r.t. BEPICOLOMBO MPO (-121)
           -121540 w.r.t. BEPICOLOMBO MPO (-121)
           -121600 w.r.t. BEPICOLOMBO MPO (-121)
           Start of Interval (UTC)             End of Interval (UTC)
           -----------------------------       -------------------------
   ----
           2027-JAN-02 23:01:53.350            2027-JAN-08 00:59:37.932
 
 
    2. If you run your extended program with the original meta-kernel, the SPICE(SPKINSUFFDATA) error should be produced by the spkpos_c function because you have not loaded enough ephemeris data to compute the position of Jupiter with respect to Mercury. The loaded SPKs contain data for Mercury relative to the Solar System Barycenter, and for the Jupiter System Barycenter relative to the Solar System Barycenter, but the data for Jupiter relative to the Jupiter System Barycenter are missing:

 
      Additional kernels required for this task:
 
         1. Generic Jovian Satellite Ephemeris SPK:
 
               jup365_2027.bsp
 
      available in the NAIF server at:
 
   https://naif.jpl.nasa.gov/pub/naif/generic_kernels/spk/
   satellites/a_old_versions
 
    Download the relevant SPK, add it to the meta-kernel and run again your extended program. The solution for the input UTC time "2027 JAN 05 02:04:36" when using the downloaded Jovian Satellite Ephemeris SPK:

      Actual position of Jupiter as seen from Mercury in the
         J2000 frame (km):
         X =   -623644094.418
         Y =    532767093.112
         Z =    251130102.035
    3. When using 'NONE' aberration corrections, spkpos_c returns the geometric position of the target body relative to the observer. If 'LT' is used, the returned vector corresponds to the position of the target at the moment it emitted photons arriving at the observer at `et'. If 'LT+S' is used instead, the returned vector takes into account the observer's velocity relative to the solar system barycenter. The solution for the input UTC time "2027 JAN 05 02:04:36" is:

      Actual (geometric) position of Sun as seen from Mercury in the
         J2000 frame (km):
         X =    -23438490.402
         Y =     54294213.485
         Z =     31433347.025
      Light-time corrected position of Sun as seen from Mercury in the
         J2000 frame (km):
         X =    -23438492.550
         Y =     54294212.272
         Z =     31433346.550
      Apparent position of Sun as seen from Mercury in the
         J2000 frame (km):
         X =    -23430052.903
         Y =     54297381.156
         Z =     31434164.775


Top

Spacecraft Orientation and Reference Frames (xform)







Top

Task Statement




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

    1. The apparent state of Mercury as seen from BepiColombo MPO in the IAU_MERCURY body-fixed frame. This vector itself is not of any particular interest, but it is a useful intermediate quantity in some geometry calculations.

    2. The angular separation between the apparent position of Mercury as seen from BepiColombo MPO and the nominal instrument view direction.

    The nominal instrument view direction is not provided by any kernel variable, but it is indicated in the BepiColombo MPO frame kernel cited above in the section ``Kernels Used'' to be the +Z axis of the MPO_SPACECRAFT frame.

Use the program to compute these quantities at the epoch "2027 JAN 05 02:04:36" UTC.



Top

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.



Top

Approach




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

    -- Decide which SPICE kernels are necessary. Prepare a meta-kernel listing the kernels and load it into the program.

    -- Prompt the user for an input time string.

    -- Convert the input time string into ephemeris time expressed as seconds past J2000 TDB.

    -- Compute the state of Mercury relative to BepiColombo MPO in the J2000 reference frame, corrected for aberrations.

    -- Compute the state transformation matrix from J2000 to IAU_MERCURY at the epoch, adjusted for light time.

    -- Multiply the state of Mercury relative to BepiColombo MPO in the J2000 reference frame by the state transformation matrix computed in the previous step.

    -- Compute the position of Mercury relative to BepiColombo MPO in the J2000 reference frame, corrected for aberrations.

    -- Determine what the nominal instrument view direction of the BepiColombo MPO spacecraft is by examining the frame kernel's content.

    -- Compute the rotation matrix from the BepiColombo MPO spacecraft frame to J2000.

    -- Multiply the nominal instrument view direction expressed in the BepiColombo MPO spacecraft frame by the rotation matrix from the previous step.

    -- Compute the separation between the result of the previous step and the apparent position of Mercury relative to BepiColombo MPO in the J2000 frame.

HINT: Several of the steps above may be compressed into a single step 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:

    2. PCK Required Reading (pck.req)

    3. SPK Required Reading (spk.req)

    4. CK Required Reading (ck.req)

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.



Top

Solution






Top

Solution Meta-Kernel



The meta-kernel we created for the solution to this exercise is named 'xform.tm'. 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.
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
         1. Generic LSK:
 
               naif0012.tls
 
         2. BepiColombo MPO SCLK:
 
               bc_mpo_step_20230117.tsc
 
         3. Solar System Ephemeris SPK, subsetted to cover only
            the time range of interest:
 
               de432s.bsp
 
         4. BepiColombo MPO Spacecraft Trajectory SPK, subsetted
            to cover only the time range of interest:
 
               bc_mpo_mlt_50037_20260314_20280529_v05.bsp
 
         5. BepiColombo MPO FK:
 
               bc_mpo_v32.tf
 
         6. BepiColombo MPO Spacecraft CK, subsetted to cover only
            the time range of interest:
 
               bc_mpo_sc_slt_50028_20260314_20280529_f20181127_v03.bc
 
         7. Generic PCK:
 
               pck00011.tpc
 
 
   \begindata
 
    KERNELS_TO_LOAD = (
 
    'kernels/lsk/naif0012.tls',
    'kernels/sclk/bc_mpo_step_20230117.tsc',
    'kernels/spk/de432s.bsp',
    'kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp',
    'kernels/fk/bc_mpo_v32.tf',
    'kernels/ck/bc_mpo_sc_slt_50028_20260314_20280529_f20181127_v03.bc',
    'kernels/pck/pck00011.tpc'
 
                       )
 
   \begintext


Top

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.tm"
   #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 BepiColombo MPO
         The necessary ephemerides
         A planetary constants file (PCK)
         A spacecraft orientation kernel for BepiColombo MPO (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 J2000: %16.3f\n", et );
 
      /*
      Compute the apparent state of Mercury as seen from
      BepiColombo MPO in the J2000 reference frame.
      */
      spkezr_c ( "MERCURY", et,    "J2000", "LT+S",
                 "MPO",     state, &ltime           );
 
      /*
      Now obtain the transformation from the inertial
      J2000 frame to the non-inertial body-fixed IAU_MERCURY
      frame.  Since we want the apparent position, we
      need to subtract ltime from et.
      */
      sxform_c ( "J2000", "IAU_MERCURY", et-ltime, sform );
 
      /*
      Now rotate the apparent J2000 state into IAU_MERCURY
      with the following matrix multiplication:
      */
      mxvg_c ( sform, state, 6, 6, bfixst );
 
      /*
      Display the results.
      */
      printf ( "   Apparent state of Mercury as seen "
               "from BepiColombo MPO in the\n" );
      printf ( "      IAU_MERCURY 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 ( "MERCURY", et,    "IAU_MERCURY", "LT+S",
                 "MPO",     state, &ltime                 );
 
      /*
      Display the results.
      */
      printf ( "   Apparent state of Mercury as seen "
               "from BepiColombo MPO in the\n"         );
      printf ( "      IAU_MERCURY body-fixed frame "
               "(km, km/s) obtained using\n"           );
      printf ( "      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]          );
 
      /*
      Note that the velocity found by using spkezr_c
      to compute the state in the IAU_MERCURY frame differs
      at the few mm/second level from that found previously
      by calling spkezr_c and then sxform_c. Computing
      velocity via a single call to spkezr_c as we've
      done immediately above is slightly more accurate because
      it accounts for the effect of the rate of change of
      light time on the apparent angular velocity of the
      target's body-fixed reference frame.
 
      Now we are to compute the angular separation between
      the apparent position of Mercury as seen from the orbiter
      and the nominal instrument view direction.  First,
      compute the apparent position of Mercury as seen from
      BepiColombo MPO in the J2000 frame.
      */
      spkpos_c ( "MERCURY", et,  "J2000", "LT+S",
                 "MPO",     pos, &ltime            );
 
      /*
      Now compute the location of the nominal instrument view
      direction.  From reading the frame kernel we know that
      the instrument view direction is nominally the +Z axis
      of the MPO_SPACECRAFT frame defined there.
      */
      bsight[0] = 0.0;
      bsight[1] = 0.0;
      bsight[2] = 1.0;
 
      /*
      Now compute the rotation matrix from MPO_SPACECRAFT
      into J2000.
      */
      pxform_c ( "MPO_SPACECRAFT", "J2000", et, pform  );
 
      /*
      And multiply the result to obtain the nominal instrument
      view direction 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 Mercury and\n"      );
      printf ( "      the BepiColombo MPO nominal "
               "instrument view direction\n"             );
      printf ( "      (degrees):\n"                      );
      printf ( "      %16.3f\n", sep                     );
 
      /*
      Or alternatively we can work in the spacecraft
      frame directly.
      */
      spkpos_c ( "MERCURY", et,  "MPO_SPACECRAFT", "LT+S",
                 "MPO",  pos, &ltime                     );
 
      /*
      The nominal instrument view direction is the +Z-axis
      in the MPO_SPACECRAFT 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 Mercury and\n"    );
      printf ( "      the BepiColombo MPO nominal "
               "instrument view direction computed\n"   );
      printf ( "      using vectors in the "
               "MPO_SPACECRAFT frame (degrees):\n"      );
      printf ( "      %16.3f\n", sep                    );
 
      return(0);
   }


Top

Solution Sample Output



After compiling the program, execute it:

   Input UTC Time: 2027 JAN 05 02:04:36
   Converting UTC Time: 2027 JAN 05 02:04:36
      ET seconds past J2000:    852386745.184
      Apparent state of Mercury as seen from BepiColombo MPO in the
         IAU_MERCURY body-fixed frame (km, km/s):
         X =        -2354.697620
         Y =         -762.547549
         Z =        -1518.408470
        VX =            1.208589
        VY =            0.394259
        VZ =           -2.671125
      Apparent state of Mercury as seen from BepiColombo MPO in the
         IAU_MERCURY body-fixed frame (km, km/s) obtained using
         spkezr_c directly:
         X =        -2354.697620
         Y =         -762.547549
         Z =        -1518.408470
        VX =            1.208589
        VY =            0.394259
        VZ =           -2.671125
      Angular separation between the apparent position of Mercury and
         the BepiColombo MPO nominal instrument view direction
         (degrees):
                    0.009
      Angular separation between the apparent position of Mercury and
         the BepiColombo MPO nominal instrument view direction computed
         using vectors in the MPO_SPACECRAFT frame (degrees):
                    0.009


Top

Extra Credit




In this ``extra credit'' section you will be presented with more complex tasks, aimed at improving your understanding of frame transformations, and some common errors that may happen when computing them.

These ``extra credit'' tasks are provided as task statements, and unlike the regular tasks, no approach or solution source code is provided. In the next section, you will find the numeric solutions (when applicable) and answers to the questions asked in these tasks.



Top

Task statements and questions



    1. Run the original program using the input UTC time "2027 JAN 06 15:32:05". Explain what happens.

    2. Compute the angular separation between the apparent position of the Sun as seen from BepiColombo MPO and the nominal instrument view direction. Is the science deck illuminated?



Top

Solutions and answers



    1. When running the original software using as input the UTC time string "2027 JAN 06 15:32:05":

 
   =====================================================================
   ===========
 
   Toolkit version: N0067
 
   SPICE(NOFRAMECONNECT) --
 
   At epoch 8.5252159418408E+08 TDB (2027 JAN 06 15:33:14.184 TDB), ther
   e is
   insufficient information available to transform from reference frame
   -121000
   (MPO_SPACECRAFT) to reference frame 1 (J2000). MPO_SPACECRAFT is a CK
    frame; a
   CK file containing data for instrument or structure -121000 at the ep
   och shown
   above, as well as a corresponding SCLK kernel, must be loaded in orde
   r to use
   this frame. Failure to find required CK data could be due to one or m
   ore CK
   files not having been loaded, or to the epoch shown above lying withi
   n a
   coverage gap or beyond the coverage bounds of the loaded CK files. It
    is also
   possible that no loaded CK file has required angular velocity data fo
   r the
   input epoch, even if a loaded CK does have attitude data for that epo
   ch. You
   can use CKBRIEF with the -dump option to display coverage intervals o
   f a CK
   file.
 
   A traceback follows.  The name of the highest level module is first.
   pxform_c --> PXFORM --> REFCHG
 
   Oh, by the way:  The SPICELIB error handling actions are USER-TAILORA
   BLE.  You
   can choose whether the Toolkit aborts or continues when errors occur,
    which
   error messages to output, and where to send the output.  Please read
   the ERROR
   "Required Reading" file, or see the routines ERRACT, ERRDEV, and ERRP
   RT.
 
   =====================================================================
   ===========
    pxform_c returns the SPICE(NOFRAMECONNECT) error, which indicates that there are not sufficient data to perform the transformation from the MPO_SPACECRAFT frame to J2000 at the requested epoch. If you summarize the BepiColombo MPO spacecraft CK using the ``ckbrief'' utility program with the -dump option (display interpolation intervals boundaries) you will find that the CK contains gaps within its segment:

 
   CKBRIEF -- Version 6.1.0, June 27, 2014 -- Toolkit Version N0067
 
 
   Summary for: kernels/ck/bc_mpo_sc_slt_50028_20260314_20280529_f201811
   27_v03.bc
 
   Segment No.: 1
 
   Object:  -121000
     Interval Begin UTC       Interval End UTC         AV
     ------------------------ ------------------------ ---
     2027-JAN-02 23:01:53.350 2027-JAN-06 11:04:56.368 Y
     2027-JAN-06 11:08:00.779 2027-JAN-06 15:30:56.685 Y
     2027-JAN-06 15:33:04.016 2027-JAN-06 22:05:57.865 Y
     2027-JAN-06 22:10:03.746 2027-JAN-08 00:59:37.932 Y
 
 
    whereas if you had used ckbrief without -dump you would have gotten the following information (only CK segment begin/end times):

 
   CKBRIEF -- Version 6.1.0, June 27, 2014 -- Toolkit Version N0067
 
 
   Summary for: kernels/ck/bc_mpo_sc_slt_50028_20260314_20280529_f201811
   27_v03.bc
 
   Object:  -121000
     Interval Begin UTC       Interval End UTC         AV
     ------------------------ ------------------------ ---
     2027-JAN-02 23:01:53.350 2027-JAN-08 00:59:37.932 Y
 
 
    which has insufficient detail to reveal the problem.

    2. By computing the apparent position of the Sun as seen from BepiColombo MPO in the MPO_SPACECRAFT frame, and the angular separation between this vector and the nominal instrument view direction (+Z-axis of the MPO_SPACECRAFT frame), you will find whether the science deck is illuminated. The solution for the input UTC time "2027 JAN 05 02:04:36" is:

   Angular separation between the apparent position of the Sun and the
   BepiColombo MPO nominal instrument view direction (degrees):
       135.393
 
   Science Deck illumination:
      BepiColombo MPO Science Deck IS NOT illuminated.
    since the angular separation is greater than 90 degrees.



Top

Computing Sub-s/c and Sub-solar Points on an Ellipsoid and a DSK (subpts)







Top

Task Statement




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

    1. The apparent sub-observer point of BepiColombo MPO on Mercury, in the body fixed frame IAU_MERCURY, in kilometers, and the spacecraft altitude as the distance between the spacecraft and this point, in kilometers.

    2. The apparent sub-solar point on Mercury, as seen from BepiColombo MPO in the body fixed frame IAU_MERCURY, in kilometers.

The program computes each point twice: once using an ellipsoidal shape model and the

        near point/ellipsoid
definition, and once using a DSK shape model and the

        nadir/dsk/unprioritized
definition.

The program displays the results. Use the program to compute these quantities at "2027 JAN 05 02:04:36" UTC.



Top

Learning Goals




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



Top

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: how would the results change if the sub-solar and sub-observer points were computed using the

        intercept/ellipsoid
and

        intercept/dsk/unprioritized
definitions? Which definition is appropriate?



Top

Solution






Top

Solution Meta-Kernel



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

   KPL/MK
 
      This is the meta-kernel used in the solution of the
      ``Computing Sub-s/c and Sub-solar Points on an Ellipsoid
      and a DSK'' task in the Remote Sensing Hands On Lesson.
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
         1. Generic LSK:
 
               naif0012.tls
 
         2. Solar System Ephemeris SPK, subsetted to cover only
            the time range of interest:
 
               de432s.bsp
 
         3. BepiColombo MPO Spacecraft Trajectory SPK, subsetted
            to cover only the time range of interest:
 
               bc_mpo_mlt_50037_20260314_20280529_v05.bsp
 
         4. Generic PCK:
 
               pck00011.tpc
 
         5. Low-resolution Mercury DSK:
 
               mercury_lowres.bds
 
   \begindata
 
    KERNELS_TO_LOAD = (
 
    'kernels/lsk/naif0012.tls',
    'kernels/spk/de432s.bsp',
    'kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp',
    'kernels/pck/pck00011.tpc'
    'kernels/dsk/mercury_lowres.bds'
 
                      )
 
   \begintext


Top

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.tm"
   #define STRLEN 50
 
   int main (void)
   {
      /*
      Local Variables
      */
      SpiceChar             * method;
      SpiceChar               utctim [STRLEN];
 
      SpiceDouble             et;
      SpiceDouble             spoint [3];
      SpiceDouble             srfvec [3];
      SpiceDouble             trgepc;
 
      SpiceInt                i;
 
      /*
      Load the kernels that this program requires.  We
      will need:
 
         A leapseconds kernel
         The necessary ephemerides
         A planetary constants file (PCK)
         A DSK file containing Mercury shape data
      */
      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 J2000: %16.3f\n", et );
 
 
      for ( i = 0;  i < 2;  i++ )
      {
         if ( i == 0 )
         {
            /*
            Use the "near point" sub-point definition
            and an ellipsoidal model.
            */
            method = "NEAR POINT/Ellipsoid";
         }
         else
         {
            /*
            Use the "nadir" sub-point definition
            and a DSK model.
            */
            method = "NADIR/DSK/Unprioritized";
         }
 
         printf ( "\n Sub-point/target shape model: %s\n\n",
                  method                                      );
 
         /*
         Compute the apparent sub-observer point of BepiColombo MPO
         on Mercury.
         */
         subpnt_c ( method,
                    "MERCURY", et,     "IAU_MERCURY", "LT+S",
                    "MPO",     spoint, &trgepc,       srfvec  );
 
         printf ( "   Apparent sub-observer point of BepiColombo "
                  "MPO on Mercury\n"                          );
         printf ( "   in the IAU_MERCURY frame (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", vnorm_c(srfvec)       );
 
         /*
         Compute the apparent sub-solar point on Mercury
         as seen from BepiColombo MPO.
         */
         subslr_c ( method,
                    "MERCURY", et,     "IAU_MERCURY", "LT+S",
                    "MPO",     spoint, &trgepc,       srfvec  );
 
         printf ( "   Apparent sub-solar point on Mercury "
                  "as seen from BepiColombo\n"                );
         printf ( "   MPO in the IAU_MERCURY frame (km):\n"   );
         printf ( "      X = %16.3f\n", spoint[0]             );
         printf ( "      Y = %16.3f\n", spoint[1]             );
         printf ( "      Z = %16.3f\n", spoint[2]             );
      }
 
      printf ( "\n" );
 
      return(0);
   }


Top

Solution Sample Output



After compiling the program, execute it:

   Input UTC Time: 2027 JAN 05 02:04:36
    Converting UTC Time: 2027 JAN 05 02:04:36
      ET seconds past J2000:    852386745.184
 
    Sub-point/target shape model: NEAR POINT/Ellipsoid
 
      Apparent sub-observer point of BepiColombo MPO on Mercury
      in the IAU_MERCURY frame (km):
         X =         1978.726
         Y =          640.793
         Z =         1275.611
       ALT =          463.634
      Apparent sub-solar point on Mercury as seen from BepiColombo
      MPO in the IAU_MERCURY frame (km):
         X =         1526.831
         Y =         1903.936
         Z =           -1.436
 
    Sub-point/target shape model: NADIR/DSK/Unprioritized
 
      Apparent sub-observer point of BepiColombo MPO on Mercury
      in the IAU_MERCURY frame (km):
         X =         1979.558
         Y =          641.062
         Z =         1276.148
       ALT =          462.608
      Apparent sub-solar point on Mercury as seen from BepiColombo
      MPO in the IAU_MERCURY frame (km):
         X =         1525.673
         Y =         1902.492
         Z =           -1.434
 


Top

Extra Credit




In this ``extra credit'' section you will be presented with more complex tasks, aimed at improving your understanding of subpnt_c and subslr_c functions.

These ``extra credit'' tasks are provided as task statements, and unlike the regular tasks, no approach or solution source code is provided. In the next section, you will find the numeric solutions (when applicable) and answers to the questions asked in these tasks.



Top

Task statements and questions



    1. Recompute the apparent sub-solar point on Mercury as seen from BepiColombo MPO in the body fixed frame IAU_MERCURY in kilometers using the 'Intercept/ellipsoid' method at "2027 JAN 05 02:04:36". Explain the differences.

    2. Compute the geometric sub-spacecraft point of BepiColombo MPO on Europa in the body fixed frame IAU_EUROPA in kilometers using the 'Near point/ellipsoid' method at "2027 JAN 05 02:04:36". This point itself is not of any particular interest, but it is useful as input to the next ``extra credit'' task.

    3. Transform the sub-spacecraft Cartesian coordinates obtained in the previous task to planetocentric and planetographic coordinates. When computing planetographic coordinates, retrieve Europa' radii by calling bodvrd_c and use the first element of the returned radii values as Europa' equatorial radius. Explain why planetocentric and planetographic latitudes and longitudes are different. Explain why the planetographic altitude for a point on the surface of Europa is not zero and whether this is correct or not.



Top

Solutions and answers



    1. The differences observed are due to the computation method. The ``Intercept/ellipsoid'' method defines the sub-solar point as the target surface intercept of the line containing the Sun and the target's center, while the ``Near point/ellipsoid'' method defines the sub-solar point as the the nearest point on the target relative to the Sun. Since Mercury is not spherical, these two points are not the same:

      Apparent sub-solar point on Mercury as seen from BepiColombo MPO
      in the IAU_MERCURY frame using the 'Near Point: ellipsoid' method
      (km):
         X =         1526.828
         Y =         1903.939
         Z =           -1.435
 
      Apparent sub-solar point on Mercury as seen from BepiColombo MPO
      in the IAU_MERCURY frame using the 'Intercept: ellipsoid' method
      (km):
         X =         1526.828
         Y =         1903.939
         Z =           -1.438
    2. Download the relevant SPK prodiving ephemeris data for Europa, add it to the meta-kernel and run again your extended program:

 
      Additional kernels required for this task:
 
         1. Generic Jovian Satellite Ephemeris SPK:
 
               jup365_2027.bsp
 
      available in the NAIF server at:
 
   https://naif.jpl.nasa.gov/pub/naif/generic_kernels/spk/
   satellites/a_old_versions
 
    The geometric sub-spacecraft point of BepiColombo MPO on Europa in the body fixed frame IAU_EUROPA in kilometers at "2027 JAN 05 02:04:36" UTC epoch is:

      Geometric sub-spacecraft point of BepiColombo MPO on Europa in
      the IAU_EUROPA frame using the 'Near Point: ellipsoid' method
      (km):
         X =         -753.484
         Y =        -1366.703
         Z =          -24.296
    3. The sub-spacecraft point of BepiColombo MPO on Europa in planetocentric and planetographic coordinates at "2027 JAN 05 02:04:36" UTC epoch is:

      Planetocentric coordinates of the BepiColombo MPO
      sub-spacecraft point on Europa (degrees, km):
      LAT =           -0.892
      LON =         -118.869
      R   =         1560.835
 
      Planetographic coordinates of the BepiColombo MPO
      sub-spacecraft point on Europa (degrees, km):
      LAT =           -0.895
      LON =          118.869
      ALT =           -1.764
    The planetocentric and planetographic longitudes are different (``graphic'' = 360 - ``centric'') because planetographic longitudes on Europa are measured positive west as defined by the Europa' rotation direction.

    The planetocentric and planetographic latitudes are different because the planetocentric latitude was computed as the angle between the direction from the center of the body to the point and the equatorial plane, while the planetographic latitude was computed as the angle between the surface normal at the point and the equatorial plane.

    The planetographic altitude is non zero because it was computed using a different and incorrect Europa surface model: a spheroid with equal equatorial radii. The surface point computed by subpnt_c was computed by treating Europa as a triaxial ellipsoid with different equatorial radii. The planetographic latitude is also incorrect because it is based on the normal to the surface of the spheroid rather than the ellipsoid, In general planetographic coordinates cannot be used for bodies with shapes modeled as triaxial ellipsoids.



Top

Intersecting Vectors with an Ellipsoid and a DSK (fovint)







Top

Task Statement




Write a program that prompts the user for an input UTC time string and, for that time, computes the intersection of the BepiColombo MPO SIMBIO-SYS HRIC channel boresight and field of view (FOV) boundary vectors with the surface of Mercury. Compute each intercept twice: once with Mercury' shape modeled as an ellipsoid, and once with Mercury' shape modeled by DSK data. The program presents each point of intersection as

    1. A Cartesian vector in the IAU_MERCURY frame

    2. Planetocentric (latitudinal) coordinates in the IAU_MERCURY frame.

For each of the camera FOV boundary and boresight vectors, if an intersection is found, the program displays the results of the above computations, otherwise it indicates no intersection exists.

At each point of intersection compute the following:

    3. Phase angle

    4. Solar incidence angle

    5. Emission angle

These angles should be computed using both ellipsoidal and DSK shape models.

Additionally compute the local solar time at the intercept of the spectrometer aperture boresight with the surface of Mercury, using both ellipsoidal and DSK shape models.

Use this program to compute values at the UTC epoch:

    "2027 JAN 05 02:04:36"



Top

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 target bodies whose shapes are modeled as ellipsoids or provided by DSKs. Discover another high level geometry function and another time conversion function in CSPICE.



Top

Approach




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

    -- Decide which SPICE kernels are necessary. Prepare a meta-kernel listing the kernels and load it into the program. Remember, you will need to find a kernel with information about the BepiColombo MPO SIMBIO-SYS spectrometer.

    -- Prompt the user for an input time string.

    -- Convert the input time string into ephemeris time expressed as seconds past J2000 TDB.

    -- Retrieve the FOV (field of view) configuration for the BepiColombo MPO SIMBIO-SYS HRIC channel.

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

    -- Compute the intercept of the vector with Mercury modeled as an ellipsoid or using DSK data.

    -- If this intercept is found, convert the position vector of the intercept into planetocentric coordinates.

    Then compute the phase, solar incidence, and emission angles at the intercept. Otherwise indicate to the user no intercept was found for this vector.

    -- Compute the planetocentric longitude of the boresight intercept.

Finally

    -- Compute the local solar time at the boresight intercept longitude on a 24-hour clock. The input time for this computation should be the TDB observation epoch minus one-way light time from the boresight intercept to the spacecraft.

It may be useful to consult the BepiColombo MPO SIMBIO-SYS instrument kernel to determine the name of the SIMBIO-SYS HRIC channel 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.



Top

Solution






Top

Solution Meta-Kernel



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

   KPL/MK
 
      This is the meta-kernel used in the solution of the
      ``Intersecting Vectors with an Ellipsoid and a DSK'' task
      in the Remote Sensing Hands On Lesson.
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
         1. Generic LSK:
 
               naif0012.tls
 
         2. BepiColombo MPO SCLK:
 
               bc_mpo_step_20230117.tsc
 
         3. Solar System Ephemeris SPK, subsetted to cover only
            the time range of interest:
 
               de432s.bsp
 
         4. BepiColombo MPO Spacecraft Trajectory SPK, subsetted
            to cover only the time range of interest:
 
               bc_mpo_mlt_50037_20260314_20280529_v05.bsp
 
         5. BepiColombo MPO FK:
 
               bc_mpo_v32.tf
 
         6. BepiColombo MPO Spacecraft CK, subsetted to cover only
            the time range of interest:
 
               bc_mpo_sc_slt_50028_20260314_20280529_f20181127_v03.bc
 
         7. Generic PCK:
 
               pck00011.tpc
 
         8. SIMBIO-SYS IK:
 
               bc_mpo_simbio-sys_v08.ti
 
         9. Low-resolution Mercury DSK:
 
               mercury_lowres.bds
 
   \begindata
 
    KERNELS_TO_LOAD = (
 
    'kernels/lsk/naif0012.tls',
    'kernels/sclk/bc_mpo_step_20230117.tsc',
    'kernels/spk/de432s.bsp',
    'kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp',
    'kernels/fk/bc_mpo_v32.tf',
    'kernels/ck/bc_mpo_sc_slt_50028_20260314_20280529_f20181127_v03.bc',
    'kernels/pck/pck00011.tpc',
    'kernels/ik/bc_mpo_simbio-sys_v08.ti'
    'kernels/dsk/mercury_lowres.bds'
 
                      )
 
   \begintext


Top

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.tm"
   #define STRLEN 50
 
   /*
   BCVLEN is the maximum number of boundary corner
   vectors we can retrieve. We've extended this array by 1
   element to make room for the boresight vector.
   */
   #define BCVLEN 5
 
   int main (void)
   {
      /*
      Local Variables
      */
      SpiceBoolean            lit;
      SpiceBoolean            visibl;
 
      SpiceChar               ampm   [STRLEN];
      SpiceChar              *boolst [2] = { "false", "true" };
      SpiceChar               insfrm [STRLEN];
      SpiceChar              *method [2];
      SpiceChar               shape  [STRLEN];
      SpiceChar               time   [STRLEN];
      SpiceChar               utctim [STRLEN];
      SpiceChar              *vecnam[] = {
                                "Boundary Corner 1",
                                "Boundary Corner 2",
                                "Boundary Corner 3",
                                "Boundary Corner 4",
                                "MPO SIMBIO-SYS HRIC Boresight" };
 
      SpiceDouble             bounds [BCVLEN][3];
      SpiceDouble             bsight [3];
      SpiceDouble             emissn;
      SpiceDouble             et;
      SpiceDouble             lat;
      SpiceDouble             lon;
      SpiceDouble             phase;
      SpiceDouble             point  [3];
      SpiceDouble             radius;
      SpiceDouble             solar;
      SpiceDouble             srfvec [3];
      SpiceDouble             trgepc;
 
      SpiceInt                hr;
      SpiceInt                i;
      SpiceInt                j;
      SpiceInt                mn;
      SpiceInt                n;
      SpiceInt                marsid;
      SpiceInt                sc;
 
      SpiceBoolean            found;
 
      /*
      Load the kernels that this program requires.  We
      will need:
 
         A leapseconds kernel.
         A SCLK kernel for BepiColombo MPO.
         Any necessary ephemerides.
         The BepiColombo MPO frame kernel.
         An BepiColombo MPO C-kernel.
         A PCK file with Mercury constants.
         The BepiColombo MPO SIMBIO-SYS I-kernel.
         A DSK file containing Mercury shape data.
      */
      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 J2000: %16.3f\n\n", et );
 
      /*
      Now we need to obtain the FOV configuration of
      the SIMBIO-SYS HRIC channel.
      */
      getfvn_c ( "MPO_SIMBIO-SYS_HRIC_FPA", BCVLEN, STRLEN, STRLEN,
                 shape,  insfrm, 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];
      }
 
      /*
      Set values of "method" string that specify use of
      ellipsoidal and DSK (topographic) shape models.
 
      In this case, we can use the same methods for calls to both
      sincpt_c and ilumin_c. Note that some SPICE routines require
      different "method" inputs from those shown here. See the
      API documentation of each routine for details.
      */
      method[0] = "Ellipsoid";
      method[1] = "DSK/Unprioritized";
 
      /*
      Get Mercury ID. We'll use this ID code later, when we
      compute local solar time.
      */
      bodn2c_c ( "MERCURY", &marsid, &found );
 
      /*
      The ID code for MERCURY 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 Mercury.\n"                      );
         exit   ( EXIT_FAILURE );
      }
 
      /*
      Now perform the same set of calculations for each
      vector listed in the BOUNDS array. Use both
      ellipsoidal and detailed (DSK) shape models.
      */
      for ( i=0; i<5; i++ )
      {
         printf ( "Vector: %s\n", vecnam[i] );
 
         for ( j = 0;  j < 2;  j++ )
         {
            printf ( "\n Target shape model: %s\n\n", method[j] );
 
            /*
            Call sincpt_c to determine coordinates of the
            intersection of this vector with the surface
            of Mercury.
            */
            sincpt_c ( method[j], "MERCURY", et,     "IAU_MERCURY",
                       "LT+S",    "MPO",     insfrm, bounds[i],
                       point,     &trgepc,   srfvec, &found       );
 
            /*
            Check the found flag.  Display a message if
            the point of intersection was not found,
            otherwise continue with the calculations.
            */
            if ( !found )
            {
               printf ( "No intersection point found at "
                        "this epoch for this vector.\n"   );
            }
            else
            {
               /*
               Now, we have discovered a point of intersection.
               Start by displaying the position vector in the
               IAU_MERCURY frame of the intersection.
               */
               printf ( "  Position vector of surface intercept "
                        "in the IAU_MERCURY\n  frame (km):\n" );
               printf ( "     X   = %16.3f\n", point[0] );
               printf ( "     Y   = %16.3f\n", point[1] );
               printf ( "     Z   = %16.3f\n", point[2] );
 
               /*
               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.
               */
               illumf_c( method[j],     "MERCURY", "SUN",  et,
                         "IAU_MERCURY", "LT+S",    "MPO",  point,
                         &trgepc,       srfvec,    &phase, &solar,
                         &emissn,       &visibl,   &lit            );
 
               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 ( "  Observer visible: %s\n", boolst[visibl] );
               printf ( "  Sun visible:      %s\n", boolst[lit]    );
 
               if ( i == 4 )
               {
                  /*
                  Compute local solar time corresponding to the TDB
                  light time corrected epoch at the boresignt
                  intercept.
                  */
                  et2lst_c ( trgepc,
                             marsid,
                             lon,
                             "PLANETOCENTRIC",
                             STRLEN,
                             STRLEN,
                             &hr,
                             &mn,
                             &sc,
                             time,
                             ampm             );
 
                  printf ( "\n  Local Solar Time at boresight "
                           "intercept (24 Hour Clock):\n"      );
                  printf ( "     %s\n", time                   );
               }
               /*
               End of LST computation block.
               */
            }
            /*
            End of shape model loop.
            */
         }
         /*
         End of vector loop.
         */
         printf ( "\n" );
      }
 
      return(0);
   }


Top

Solution Sample Output



After compiling the program, execute it:

   Input UTC Time: 2027 JAN 05 02:04:36
   Converting UTC Time: 2027 JAN 05 02:04:36
     ET seconds past J2000:    852386745.184
 
   Vector: Boundary Corner 1
 
    Target shape model: Ellipsoid
 
     Position vector of surface intercept in the IAU_MERCURY
     frame (km):
        X   =         1973.717
        Y   =          645.436
        Z   =         1281.009
     Planetocentric coordinates of the intercept (degrees):
        LAT =           31.670
        LON =           18.109
     Phase angle (degrees):                     44.735
     Solar incidence angle (degrees):           44.622
     Emission angle (degrees):                   1.280
     Observer visible: true
     Sun visible:      true
 
    Target shape model: DSK/Unprioritized
 
     Position vector of surface intercept in the IAU_MERCURY
     frame (km):
        X   =         1974.257
        Y   =          645.602
        Z   =         1281.346
     Planetocentric coordinates of the intercept (degrees):
        LAT =           31.670
        LON =           18.108
     Phase angle (degrees):                     44.735
     Solar incidence angle (degrees):           46.703
     Emission angle (degrees):                   4.145
     Observer visible: true
     Sun visible:      true
 
   Vector: Boundary Corner 2
 
    Target shape model: Ellipsoid
 
     Position vector of surface intercept in the IAU_MERCURY
     frame (km):
        X   =         1979.643
        Y   =          647.354
        Z   =         1270.875
     Planetocentric coordinates of the intercept (degrees):
        LAT =           31.391
        LON =           18.108
     Phase angle (degrees):                     45.641
     Solar incidence angle (degrees):           44.447
     Emission angle (degrees):                   1.198
     Observer visible: true
     Sun visible:      true
 
    Target shape model: DSK/Unprioritized
 
     Position vector of surface intercept in the IAU_MERCURY
     frame (km):
        X   =         1980.449
        Y   =          647.601
        Z   =         1271.407
     Planetocentric coordinates of the intercept (degrees):
        LAT =           31.391
        LON =           18.108
     Phase angle (degrees):                     45.641
     Solar incidence angle (degrees):           43.796
     Emission angle (degrees):                   1.894
     Observer visible: true
     Sun visible:      true
 
   Vector: Boundary Corner 3
 
    Target shape model: Ellipsoid
 
     Position vector of surface intercept in the IAU_MERCURY
     frame (km):
        X   =         1983.307
        Y   =          636.037
        Z   =         1270.876
     Planetocentric coordinates of the intercept (degrees):
        LAT =           31.391
        LON =           17.781
     Phase angle (degrees):                     44.501
     Solar incidence angle (degrees):           44.666
     Emission angle (degrees):                   1.195
     Observer visible: true
     Sun visible:      true
 
    Target shape model: DSK/Unprioritized
 
     Position vector of surface intercept in the IAU_MERCURY
     frame (km):
        X   =         1984.034
        Y   =          636.285
        Z   =         1271.361
     Planetocentric coordinates of the intercept (degrees):
        LAT =           31.391
        LON =           17.781
     Phase angle (degrees):                     44.501
     Solar incidence angle (degrees):           45.429
     Emission angle (degrees):                   2.027
     Observer visible: true
     Sun visible:      true
 
   Vector: Boundary Corner 4
 
    Target shape model: Ellipsoid
 
     Position vector of surface intercept in the IAU_MERCURY
     frame (km):
        X   =         1977.381
        Y   =          634.119
        Z   =         1281.010
     Planetocentric coordinates of the intercept (degrees):
        LAT =           31.670
        LON =           17.780
     Phase angle (degrees):                     43.576
     Solar incidence angle (degrees):           44.840
     Emission angle (degrees):                   1.278
     Observer visible: true
     Sun visible:      true
 
    Target shape model: DSK/Unprioritized
 
     Position vector of surface intercept in the IAU_MERCURY
     frame (km):
        X   =         1978.158
        Y   =          634.384
        Z   =         1281.499
     Planetocentric coordinates of the intercept (degrees):
        LAT =           31.670
        LON =           17.781
     Phase angle (degrees):                     43.576
     Solar incidence angle (degrees):           45.349
     Emission angle (degrees):                   1.920
     Observer visible: true
     Sun visible:      true
 
   Vector: MPO SIMBIO-SYS HRIC Boresight
 
    Target shape model: Ellipsoid
 
     Position vector of surface intercept in the IAU_MERCURY
     frame (km):
        X   =         1978.524
        Y   =          640.740
        Z   =         1275.950
     Planetocentric coordinates of the intercept (degrees):
        LAT =           31.530
        LON =           17.944
     Phase angle (degrees):                     44.609
     Solar incidence angle (degrees):           44.644
     Emission angle (degrees):                   0.059
     Observer visible: true
     Sun visible:      true
 
     Local Solar Time at boresight intercept (24 Hour Clock):
        09:46:41
 
    Target shape model: DSK/Unprioritized
 
     Position vector of surface intercept in the IAU_MERCURY
     frame (km):
        X   =         1979.357
        Y   =          641.010
        Z   =         1276.487
     Planetocentric coordinates of the intercept (degrees):
        LAT =           31.530
        LON =           17.944
     Phase angle (degrees):                     44.609
     Solar incidence angle (degrees):           45.349
     Emission angle (degrees):                   1.138
     Observer visible: true
     Sun visible:      true
 
     Local Solar Time at boresight intercept (24 Hour Clock):
        09:46:41
 


Top

Extra Credit




There are no ``extra credit'' tasks for this step of the lesson.