 
Geometric Event Finding Hands-On Lesson, using MPO (Python)
===========================================================================
 
   March 01, 2023
 
 
Overview
--------------------------------------------------------
 
   This lesson illustrates how the Geometry Finder (GF) subsystem of the
   SpiceyPy Toolkit can be used to find time intervals when specified
   geometric conditions are satisfied.
 
   In this lesson the student is asked to construct a program that finds
   the time intervals, within a specified time range, when the BepiColombo
   Mercury Planetary Orbiter (MPO) is visible from ESA's deep space station
   in New Norcia. Possible occultation of the spacecraft by Mercury is to
   be considered.
 
 
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.
 
 
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 FK and IK are text kernels, the information provided in
   them can be viewed using any text editor, while the meta information
   provided in binary kernels -- SPKs and CKs -- can be viewed using
   ``commnt'' or ``spacit'' utility programs located in ``cspice/exe'' of
   Toolkit installation tree.
 
 
Tutorials
 
   The following SPICE tutorials 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
      Lunar-Earth PCK   Lunar and Earth Orientation Data
      GF                The SPICE Geometry Finder (GF) subsystem
      DSK               Detailed Target Shape (Topography) Data
 
   These tutorials are available from the NAIF server at JPL:
 
      http://naif.jpl.nasa.gov/naif/tutorials.html
 
 
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
      ---------------  -----------------------------------------
      cells.req        Cell/window initialization
      frames.req       Using reference frames
      gf.req           The SPICE geometry finder (GF) subsystem
      kernel.req       Loading SPICE kernels
      naif_ids.req     Body and reference frame names
      pck.req          Obtaining planetary constants data
      spk.req          Computing positions and velocities
      time.req         UTC to ET time conversion
      windows.req      The SPICE window data type
 
 
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 SpiceyPy functions perform functions of interest, as well
   as the names of the source files that contain these functions.
 
 
 
SpiceyPy API Documentation
 
   A SpiceyPy function's parameters specification is available using the
   built-in Python help system. A more detailed specification of the API
   can be found in the CSPICE HTML API documentation page located under
   ``cspice/doc/html/cspice''.
 
   For example, the Python help function
 
      >>> import spiceypy
      >>> help(spiceypy.str2et)
 
   describes of the str2et function's parameters, while the document
 
      cspice/doc/html/cspice/str2et_c.html
 
   describes extensively the str2et functionality.
 
 
Kernels Used
--------------------------------------------------------
 
   The following kernels are used in examples provided in this lesson:
 
      1.  Solar System Ephemeris SPK, subsetted to cover only the time
          range of interest:
 
             de432s.bsp
 
      2.  ESA stations SPK:
 
             estrack_v04.bsp
 
      3.  ESA stations frame definitions:
 
             estrack_v04.tf
 
      4.  EARTH_FIXED/ITRF93 frame connection:
 
             earthfixeditrf93.tf
 
      5.  Binary PCK for Earth:
 
             earth_070425_370426_predict.bpc
 
      6.  BepiColombo MPO Spacecraft Trajectory SPK, subsetted to cover
          only the time range of interest:
 
             bc_mpo_mlt_50037_20260314_20280529_v05.bsp
 
      7.  Generic LSK:
 
             naif0012.tls
 
      8.  Generic PCK:
 
             pck00011.tpc
 
      9.  Low-resolution Mercury DSK:
 
             mercury_lowres.bds
 
 
   These SPICE kernels are included in the lesson package.
 
 
SpiceyPy 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    viewpr     spiceypy.furnsh  spiceypy.rpd     1-8
                         spiceypy.wninsd  spiceypy.str2et
                         spiceypy.gfposc  spiceypy.timout
                         spiceypy.unload  spiceypy.wncard
                                          spiceypy.wnfetd
 
         2    visibl     spiceypy.furnsh  spiceypy.rpd     1-9
                         spiceypy.wninsd  spiceypy.str2et
                         spiceypy.gfposc  spiceypy.timout
                         spiceypy.gfoclt  spiceypy.wndifd
                         spiceypy.unload  spiceypy.wncard
                                          spiceypy.wnfetd
 
              extra (*)  spiceypy.gfdist  spiceypy.repmc   1,6-8
                         spiceypy.kclear  spiceypy.repmf
 
 
         (*) Additional APIs and kernels used in Extra Credit tasks.
 
   Use the Python built-in help system on the various functions listed
   above for the API parameters' description, and refer to the headers of
   their corresponding CSPICE versions for detailed interface
   specifications.
 
 
Find View Periods
===========================================================================
 
 
Task Statement
--------------------------------------------------------
 
   Write a program that finds the set of time intervals, within the time
   range
 
      2027 JAN 03 TDB
      2027 JAN 06 TDB
 
   when the BepiColombo Mercury Planetary Orbiter (MPO) is visible from
   ESA's New Norcia station. These time intervals are frequently called
   ``view periods.''
 
   The spacecraft is considered visible if its apparent position (that is,
   its position corrected for light time and stellar aberration) has
   elevation of at least 6 degrees in the topocentric reference frame
   NEW_NORCIA_TOPO. In this exercise, we ignore the possibility of
   occultation of the spacecraft by Mercury.
 
   Use a search step size that ensures that no view periods of duration 5
   minutes or longer will be missed by the search.
 
   Display the start and stop times of these intervals using TDB calendar
   dates and millisecond precision.
 
 
Learning Goals
--------------------------------------------------------
 
   Exposure to SPICE GF event finding routines. Familiarity with SPICE
   windows and routines that manipulate them. Exposure to SPICE time
   parsing and output formatting routines.
 
 
Approach
--------------------------------------------------------
 
 
Solution steps
 
   A possible solution could consist of the following steps:
 
   Preparation:
 
       1.   Decide what SPICE kernels are necessary. Use the SPICE summary
            tool BRIEF to examine the coverage of the binary kernels and
            verify the availability of required data.
 
       2.   Create a meta-kernel listing the SPICE kernels to be loaded.
            (Hint: consult a programming example tutorial, or the
            Introduction to Kernels tutorial, for a reminder of how to do
            this.)
 
            Name the meta-kernel 'viewpr.tm'.
 
   Next, write a program that performs the following steps:
 
       1.   Use spiceypy.furnsh to load the meta-kernel.
 
       2.   Create confinement and output SpiceyPy windows using
            stypes.SPICEDOUBLE_CELL.
 
       3.   Insert the given time bounds into the confinement window using
            spiceypy.wninsd.
 
       4.   Select a step size for searching for visibility state
            transitions: in this case, each target rise or set event is a
            state transition.
 
            The step size must be large enough so the search proceeds with
            reasonable speed, but small enough so that no visibility
            transition events---that is, target rise or set events---are
            missed.
 
       5.   Use the GF routine spiceypy.gfposc to find the window of times,
            within the confinement window, during which the BepiColombo MPO
            spacecraft is above the elevation limit as seen from ESA's New
            Norcia station, in the reference frame NEW_NORCIA_TOPO.
 
            Use light time and stellar aberration corrections for the
            apparent position of the spacecraft as seen from the station.
 
       6.   Fetch and display the contents of the result window. Use
            spiceypy.wnfetd to extract from the result window the start and
            stop times of each time interval. Display each of the intervals
            in the result window as a pair of start and stop times. Express
            each time as a TDB calendar date using the routine
            spiceypy.timout.
 
   You may find it useful to consult the references listed above. In
   particular, the header of the SPICE GF function spiceypy.gfposc contains
   pertinent documentation.
 
 
Solution
--------------------------------------------------------
 
 
Solution Meta-Kernel
 
   The meta-kernel we created for the solution to this exercise is named
   'viewpr.tm'. Its contents follow:
 
      KPL/MK
 
         This is the meta-kernel used in the solution of the tasks in the
         Geometric Event Finding Hands On Lesson.
 
         The names and contents of the kernels referenced by this
         meta-kernel are as follows:
 
            1.  Solar System Ephemeris SPK, subsetted to cover only the
                time range of interest:
 
                   de432s.bsp
 
            2.  ESA stations SPK:
 
                   estrack_v04.bsp
 
            3.  ESA stations frame definitions:
 
                   estrack_v04.tf
 
            4.  EARTH_FIXED/ITRF93 frame connection:
 
                   earthfixeditrf93.tf
 
            5.  Binary PCK for Earth:
 
                   earth_070425_370426_predict.bpc
 
            6.  BepiColombo MPO Spacecraft Trajectory SPK, subsetted to
                cover only the time range of interest:
 
                   bc_mpo_mlt_50037_20260314_20280529_v05.bsp
 
            7.  Generic LSK:
 
                   naif0012.tls
 
            8.  Generic PCK:
 
                  pck00011.tpc
 
 
      \begindata
 
         KERNELS_TO_LOAD = (
 
            'kernels/spk/de432s.bsp'
            'kernels/spk/estrack_v04.bsp'
            'kernels/fk/estrack_v04.tf'
            'kernels/fk/earthfixeditrf93.tf'
            'kernels/pck/earth_070425_370426_predict.bpc'
            'kernels/lsk/naif0012.tls'
            'kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp'
            'kernels/pck/pck00011.tpc'
 
                           )
 
      \begintext
 
 
 
Solution Code
 
   The example program below shows one possible solution.
 
 
      #
      # Solution viewpr
      #
      from __future__ import print_function
      import spiceypy.utils.support_types as stypes
      import spiceypy
 
      def viewpr():
          #
          # Local Parameters
          #
          METAKR = 'viewpr.tm'
          TDBFMT = 'YYYY MON DD HR:MN:SC.### (TDB) ::TDB'
          MAXIVL = 1000
          MAXWIN = 2 * MAXIVL
 
          #
          # Load the meta-kernel.
          #
          spiceypy.furnsh( METAKR )
 
          #
          # Assign the inputs for our search.
          #
          # Since we're interested in the apparent location of the
          # target, we use light time and stellar aberration
          # corrections. We use the "converged Newtonian" form
          # of the light time correction because this choice may
          # increase the accuracy of the occultation times we'll
          # compute using gfoclt.
          #
          srfpt  = 'NEW_NORCIA'
          obsfrm = 'NEW_NORCIA_TOPO'
          target = 'MPO'
          abcorr = 'CN+S'
          start  = '2027 JAN 03 TDB'
          stop   = '2027 JAN 06 TDB'
          elvlim =  6.0
 
          #
          # The elevation limit above has units of degrees; we convert
          # this value to radians for computation using SPICE routines.
          # We'll store the equivalent value in radians in revlim.
          #
          revlim = spiceypy.rpd() * elvlim
 
          #
          # Since SPICE doesn't directly support the AZ/EL coordinate
          # system, we use the equivalent constraint
          #
          #    latitude > revlim
          #
          # in the latitudinal coordinate system, where the reference
          # frame is topocentric and is centered at the viewing location.
          #
          crdsys = 'LATITUDINAL'
          coord  = 'LATITUDE'
          relate = '>'
 
          #
          # The adjustment value only applies to absolute extrema
          # searches; simply give it an initial value of zero
          # for this inequality search.
          #
          adjust = 0.0
 
          #
          # stepsz is the step size, measured in seconds, used to search
          # for times bracketing a state transition. Since we don't expect
          # any events of interest to be shorter than five minutes, and
          # since the separation between events is well over 5 minutes,
          # we'll use this value as our step size. Units are seconds.
          #
          stepsz = 300.0
 
          #
          # Display a banner for the output report:
          #
          print( '\n{:s}\n'.format(
                 'Inputs for target visibility search:' )  )
 
          print( '   Target                       = '
                 '{:s}'.format( target )  )
          print( '   Observation surface location = '
                 '{:s}'.format( srfpt  )  )
          print( '   Observer\'s reference frame   = '
                 '{:s}'.format( obsfrm )  )
          print( '   Elevation limit (degrees)    = '
                 '{:f}'.format( elvlim )  )
          print( '   Aberration correction        = '
                 '{:s}'.format( abcorr )  )
          print( '   Step size (seconds)          = '
                 '{:f}'.format( stepsz )  )
 
          #
          # Convert the start and stop times to ET.
          #
          etbeg = spiceypy.str2et( start )
          etend = spiceypy.str2et( stop  )
 
          #
          # Display the search interval start and stop times
          # using the format shown below.
          #
          #    2004 MAY 06 20:15:00.000 (TDB)
          #
          timstr = spiceypy.timout( etbeg, TDBFMT )
          print( '   Start time                   = '
                 '{:s}'.format(timstr) )
 
          timstr = spiceypy.timout( etend, TDBFMT )
          print( '   Stop time                    = '
                 '{:s}'.format(timstr) )
 
          print( ' ' )
 
          #
          # Initialize the "confinement" window with the interval
          # over which we'll conduct the search.
          #
          cnfine = stypes.SPICEDOUBLE_CELL(2)
          spiceypy.wninsd( etbeg, etend, cnfine )
 
          #
          # In the call below, the maximum number of window
          # intervals gfposc can store internally is set to MAXIVL.
          # We set the cell size to MAXWIN to achieve this.
          #
          riswin = stypes.SPICEDOUBLE_CELL( MAXWIN )
 
          #
          # Now search for the time period, within our confinement
          # window, during which the apparent target has elevation
          # at least equal to the elevation limit.
          #
          spiceypy.gfposc( target, obsfrm, abcorr, srfpt,
                           crdsys, coord,  relate, revlim,
                           adjust, stepsz, MAXIVL, cnfine, riswin )
 
          #
          # The function wncard returns the number of intervals
          # in a SPICE window.
          #
          winsiz = spiceypy.wncard( riswin )
 
          if winsiz == 0:
 
              print( 'No events were found.' )
 
          else:
 
              #
              # Display the visibility time periods.
              #
              print( 'Visibility times of {0:s} '
                     'as seen from {1:s}:\n'.format(
                      target, srfpt )                )
 
              for  i  in  range(winsiz):
                  #
                  # Fetch the start and stop times of
                  # the ith interval from the search result
                  # window riswin.
                  #
                  [intbeg, intend] = spiceypy.wnfetd( riswin, i )
 
                  #
                  # Convert the rise time to a TDB calendar string.
                  #
                  timstr = spiceypy.timout( intbeg, TDBFMT )
 
                  #
                  # Write the string to standard output.
                  #
                  if  i  ==  0:
 
                      print( 'Visibility or window start time:'
                             '  {:s}'.format( timstr )          )
                  else:
 
                      print( 'Visibility start time:          '
                             '  {:s}'.format( timstr )          )
 
                  #
                  # Convert the set time to a TDB calendar string.
                  #
                  timstr = spiceypy.timout( intend, TDBFMT )
 
                  #
                  # Write the string to standard output.
                  #
                  if  i  ==  (winsiz-1):
 
                      print( 'Visibility or window stop time: '
                             '  {:s}'.format( timstr )          )
                  else:
 
                      print( 'Visibility stop time:           '
                             '  {:s}'.format( timstr )          )
 
                  print( ' ' )
 
          spiceypy.unload( METAKR )
 
      if __name__ == '__main__':
          viewpr()
 
 
Solution Sample Output
 
   Numerical results shown for this example may differ across platforms
   since the results depend on the SPICE kernels used as input and on the
   host platform's arithmetic implementation.
 
   Execute the program. The output is:
 
 
      Inputs for target visibility search:
 
         Target                       = MPO
         Observation surface location = NEW_NORCIA
         Observer's reference frame   = NEW_NORCIA_TOPO
         Elevation limit (degrees)    = 6.000000
         Aberration correction        = CN+S
         Step size (seconds)          = 300.000000
         Start time                   = 2027 JAN 03 00:00:00.000 (TDB)
         Stop time                    = 2027 JAN 06 00:00:00.000 (TDB)
 
      Visibility times of MPO as seen from NEW_NORCIA:
 
      Visibility or window start time:  2027 JAN 03 00:00:00.000 (TDB)
      Visibility stop time:             2027 JAN 03 10:58:25.063 (TDB)
 
      Visibility start time:            2027 JAN 03 21:55:08.488 (TDB)
      Visibility stop time:             2027 JAN 04 11:01:14.279 (TDB)
 
      Visibility start time:            2027 JAN 04 21:58:41.333 (TDB)
      Visibility stop time:             2027 JAN 05 11:04:00.020 (TDB)
 
      Visibility start time:            2027 JAN 05 22:02:18.477 (TDB)
      Visibility or window stop time:   2027 JAN 06 00:00:00.000 (TDB)
 
 
 
Find Times when Target is Visible
===========================================================================
 
 
Task Statement
--------------------------------------------------------
 
   Extend the program of the previous chapter to find times when the
   BepiColombo MPO orbiter is:
 
       --   Above the elevation limit in the NEW_NORCIA_TOPO topocentric
            reference frame.
 
       --   and is not occulted by Mercury
 
   Finding time intervals that satisfy the second condition requires a
   search for occultations of the spacecraft by Mercury. Perform this
   search twice: once using an ellipsoidal shape model for Mercury, and
   once using a DSK shape model.
 
   Compute the final results twice as well, using the results of both
   occultation searches.
 
   For each of the two shape model cases, store the set of time intervals
   when the spacecraft is visible in a SpiceyPy window. We'll call this the
   ``result window.''
 
   Display each of the intervals in the result window as a pair of start
   and stop times. Express each time as a TDB calendar date using the same
   format as in the previous program.
 
 
Learning Goals
--------------------------------------------------------
 
   Familiarity with the GF occultation finding routine spiceypy.gfoclt.
   Experience with Digital Shape Kernel (DSK) shape models. Further
   experience with the SpiceyPy window functions.
 
 
Approach
--------------------------------------------------------
 
 
Solution steps
 
   A possible solution would consist of the following steps:
 
       1.   Use the meta-kernel from the previous chapter as the starting
            point. Add more kernels to it as needed.
 
            Name the meta-kernel 'visibl.tm'.
 
       2.   Include the code from the program of the previous chapter in a
            new source file; modify this code to create the new program.
 
       3.   Your program will need additional windows to capture the
            results of occultation searches performed using both
            ellipsoidal and DSK shape models. Additional windows will be
            needed to compute the set differences of the elevation search
            (``view period'') window and each of the occultation search
            windows. Further details are provided below.
 
            Create additional output SpiceyPy windows using
            stypes.SPICEDOUBLE_CELL.
 
       4.   The remaining steps can be performed twice: once using an
            ellipsoidal shape model for Mercury, and once using a DSK
            Mercury shape model. Alternatively, two copies of the entire
            solution program can be created: one for each shape model.
 
       5.   Search for occultations of the BepiColombo MPO orbiter as seen
            from New Norcia station using spiceypy.gfoclt. Use as the
            confinement window for this search the result window from the
            elevation search performed by spiceypy.gfposc.
 
            Since occultations occur when the apparent BepiColombo MPO
            spacecraft position is behind the apparent figure of Mercury,
            light time correction must be performed for the occultation
            search. To improve accuracy of the occultation state
            determination, use ``converged Newtonian'' light time
            correction.
 
       6.   Use the SpiceyPy window subtraction routine spiceypy.wndifd to
            subtract the window of times when the spacecraft is occulted
            from the window of times when the spacecraft is above the
            elevation limit. The difference window is the final result.
 
       7.   Modify the code to display the contents of the difference
            window.
 
   This completes the assignment.
 
 
Solution
--------------------------------------------------------
 
 
Solution Meta-Kernel
 
   The meta-kernel we created for the solution to this exercise is named
   'visibl.tm'. Its contents follow:
 
      KPL/MK
 
         This is the meta-kernel used in the solution of the tasks in the
         Geometric Event Finding Hands On Lesson.
 
         The names and contents of the kernels referenced by this
         meta-kernel are as follows:
 
            1.  Solar System Ephemeris SPK, subsetted to cover only the
                time range of interest:
 
                   de432s.bsp
 
            2.  ESA stations SPK:
 
                   estrack_v04.bsp
 
            3.  ESA stations frame definitions:
 
                   estrack_v04.tf
 
            4.  EARTH_FIXED/ITRF93 frame connection:
 
                   earthfixeditrf93.tf
 
            5.  Binary PCK for Earth:
 
                   earth_070425_370426_predict.bpc
 
            6.  BepiColombo MPO Spacecraft Trajectory SPK, subsetted to
                cover only the time range of interest:
 
                   bc_mpo_mlt_50037_20260314_20280529_v05.bsp
 
            7.  Generic LSK:
 
                   naif0012.tls
 
            8.  Generic PCK:
 
                  pck00011.tpc
 
            9. Low-resolution Mercury DSK:
 
                  mercury_lowres.bds
 
      \begindata
 
         KERNELS_TO_LOAD = (
 
            'kernels/spk/de432s.bsp'
            'kernels/spk/estrack_v04.bsp'
            'kernels/fk/estrack_v04.tf'
            'kernels/fk/earthfixeditrf93.tf'
            'kernels/pck/earth_070425_370426_predict.bpc'
            'kernels/lsk/naif0012.tls'
            'kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp'
            'kernels/pck/pck00011.tpc'
            'kernels/dsk/mercury_lowres.bds'
 
                           )
 
      \begintext
 
 
 
Solution Code
 
 
      #
      # Solution visibl
      #
      from __future__ import print_function
 
      #
      # SpiceyPy package:
      #
      import spiceypy.utils.support_types as stypes
      import spiceypy
 
      def visibl():
          #
          # Local Parameters
          #
          METAKR = 'visibl.tm'
          TDBFMT = 'YYYY MON DD HR:MN:SC.### TDB ::TDB'
          MAXIVL = 1000
          MAXWIN = 2 * MAXIVL
 
          #
          # Load the meta-kernel.
          #
          spiceypy.furnsh( METAKR )
 
          #
          # Assign the inputs for our search.
          #
          # Since we're interested in the apparent location of the
          # target, we use light time and stellar aberration
          # corrections. We use the "converged Newtonian" form
          # of the light time correction because this choice may
          # increase the accuracy of the occultation times we'll
          # compute using gfoclt.
          #
          srfpt  = 'NEW_NORCIA'
          obsfrm = 'NEW_NORCIA_TOPO'
          target = 'MPO'
          abcorr = 'CN+S'
          start  = '2027 JAN 03 TDB'
          stop   = '2027 JAN 06 TDB'
          elvlim =  6.0
 
          #
          # The elevation limit above has units of degrees; we convert
          # this value to radians for computation using SPICE routines.
          # We'll store the equivalent value in radians in revlim.
          #
          revlim = spiceypy.rpd() * elvlim
 
          #
          # We model the target shape as a point. We either model the
          # blocking body's shape as an ellipsoid, or we represent
          # its shape using actual topographic data. No body-fixed
          # reference frame is required for the target since its
          # orientation is not used.
          #
          back   = target
          bshape = 'POINT'
          bframe = ' '
          front  = 'MERCURY'
          fshape = 'ELLIPSOID'
          fframe = 'IAU_MERCURY'
 
          #
          # The occultation type should be set to 'ANY' for a point
          # target.
          #
          occtyp = 'any'
 
          #
          # Since SPICE doesn't directly support the AZ/EL coordinate
          # system, we use the equivalent constraint
          #
          #    latitude > revlim
          #
          # in the latitudinal coordinate system, where the reference
          # frame is topocentric and is centered at the viewing location.
          #
          crdsys = 'LATITUDINAL'
          coord  = 'LATITUDE'
          relate = '>'
 
          #
          # The adjustment value only applies to absolute extrema
          # searches; simply give it an initial value of zero
          # for this inequality search.
          #
          adjust = 0.0
 
          #
          # stepsz is the step size, measured in seconds, used to search
          # for times bracketing a state transition. Since we don't expect
          # any events of interest to be shorter than five minutes, and
          # since the separation between events is well over 5 minutes,
          # we'll use this value as our step size. Units are seconds.
          #
          stepsz = 300.0
 
          #
          # Display a banner for the output report:
          #
          print( '\n{:s}\n'.format(
                 'Inputs for target visibility search:' )  )
 
          print( '   Target                       = '
                 '{:s}'.format( target )  )
          print( '   Observation surface location = '
                 '{:s}'.format( srfpt  )  )
          print( '   Observer\'s reference frame   = '
                 '{:s}'.format( obsfrm )  )
          print( '   Blocking body                = '
                 '{:s}'.format( front  )  )
          print( '   Blocker\'s reference frame    = '
                 '{:s}'.format( fframe )  )
          print( '   Elevation limit (degrees)    = '
                 '{:f}'.format( elvlim )  )
          print( '   Aberration correction        = '
                 '{:s}'.format( abcorr )  )
          print( '   Step size (seconds)          = '
                 '{:f}'.format( stepsz )  )
 
          #
          # Convert the start and stop times to ET.
          #
          etbeg = spiceypy.str2et( start )
          etend = spiceypy.str2et( stop  )
 
          #
          # Display the search interval start and stop times
          # using the format shown below.
          #
          #    2004 MAY 06 20:15:00.000 (TDB)
          #
          btmstr = spiceypy.timout( etbeg, TDBFMT )
          print( '   Start time                   = '
                 '{:s}'.format(btmstr) )
 
          etmstr = spiceypy.timout( etend, TDBFMT )
          print( '   Stop time                    = '
                 '{:s}'.format(etmstr) )
 
          print( ' ' )
 
          #
          # Initialize the "confinement" window with the interval
          # over which we'll conduct the search.
          #
          cnfine = stypes.SPICEDOUBLE_CELL(2)
          spiceypy.wninsd( etbeg, etend, cnfine )
 
          #
          # In the call below, the maximum number of window
          # intervals gfposc can store internally is set to MAXIVL.
          # We set the cell size to MAXWIN to achieve this.
          #
          riswin = stypes.SPICEDOUBLE_CELL( MAXWIN )
 
          #
          # Now search for the time period, within our confinement
          # window, during which the apparent target has elevation
          # at least equal to the elevation limit.
          #
          spiceypy.gfposc( target, obsfrm, abcorr, srfpt,
                           crdsys, coord,  relate, revlim,
                           adjust, stepsz, MAXIVL, cnfine, riswin )
 
          #
          # Now find the times when the apparent target is above
          # the elevation limit and is not occulted by the
          # blocking body (Mercury). We'll find the window of times when
          # the target is above the elevation limit and *is* occulted,
          # then subtract that window from the view period window
          # riswin found above.
          #
          # For this occultation search, we can use riswin as
          # the confinement window because we're not interested in
          # occultations that occur when the target is below the
          # elevation limit.
          #
          # Find occultations within the view period window.
          #
          print( ' Searching using ellipsoid target shape model...' )
 
          eocwin = stypes.SPICEDOUBLE_CELL( MAXWIN )
 
          fshape = 'ELLIPSOID'
 
          spiceypy.gfoclt( occtyp, front,  fshape,  fframe,
                           back,   bshape, bframe,  abcorr,
                           srfpt,  stepsz, riswin,  eocwin )
          print( ' Done.' )
 
          #
          # Subtract the occultation window from the view period
          # window: this yields the time periods when the target
          # is visible.
          #
          evswin = spiceypy.wndifd( riswin, eocwin )
 
          #
          #  Repeat the search using low-resolution DSK data
          # for the front body.
          #
          print( ' Searching using DSK target shape model...' )
 
          docwin = stypes.SPICEDOUBLE_CELL( MAXWIN )
 
          fshape = 'DSK/UNPRIORITIZED'
 
          spiceypy.gfoclt( occtyp, front,  fshape,  fframe,
                           back,   bshape, bframe,  abcorr,
                           srfpt,  stepsz, riswin,  docwin )
          print( ' Done.\n' )
 
          dvswin = spiceypy.wndifd( riswin, docwin )
 
          #
          # The function wncard returns the number of intervals
          # in a SPICE window.
          #
          winsiz = spiceypy.wncard( evswin )
 
          if winsiz == 0:
 
              print( 'No events were found.' )
 
          else:
              #
              # Display the visibility time periods.
              #
              print( 'Visibility start and stop times of '
                     '{0:s} as seen from {1:s}\n'
                     'using both ellipsoidal and DSK '
                     'target shape models:\n'.format(
                         target, srfpt )                 )
 
              for  i  in  range(winsiz):
                  #
                  # Fetch the start and stop times of
                  # the ith interval from the ellipsoid
                  # search result window evswin.
                  #
                  [intbeg, intend] = spiceypy.wnfetd( evswin, i )
 
                  #
                  # Convert the rise time to TDB calendar strings.
                  # Write the results.
                  #
                  btmstr = spiceypy.timout( intbeg, TDBFMT )
                  etmstr = spiceypy.timout( intend, TDBFMT )
 
                  print( ' Ell: {:s} : {:s}'.format( btmstr, etmstr ) )
 
                  #
                  # Fetch the start and stop times of
                  # the ith interval from the DSK
                  # search result window dvswin.
                  #
                  [dintbg, dinten] = spiceypy.wnfetd( dvswin, i )
 
                  #
                  # Convert the rise time to TDB calendar strings.
                  # Write the results.
                  #
                  btmstr = spiceypy.timout( dintbg, TDBFMT )
                  etmstr = spiceypy.timout( dinten, TDBFMT )
 
                  print( ' DSK: {:s} : {:s}\n'.format( btmstr, etmstr ) )
              #
              # End of result display loop.
              #
 
          spiceypy.unload( METAKR )
 
      if __name__ == '__main__':
          visibl()
 
 
Solution Sample Output
 
   Numerical results shown for this example may differ across platforms
   since the results depend on the SPICE kernels used as input and on the
   host platform's arithmetic implementation.
 
   Execute the program. The output is:
 
 
      Inputs for target visibility search:
 
         Target                       = MPO
         Observation surface location = NEW_NORCIA
         Observer's reference frame   = NEW_NORCIA_TOPO
         Blocking body                = MERCURY
         Blocker's reference frame    = IAU_MERCURY
         Elevation limit (degrees)    = 6.000000
         Aberration correction        = CN+S
         Step size (seconds)          = 300.000000
         Start time                   = 2027 JAN 03 00:00:00.000 TDB
         Stop time                    = 2027 JAN 06 00:00:00.000 TDB
 
       Searching using ellipsoid target shape model...
       Done.
       Searching using DSK target shape model...
       Done.
 
      Visibility start and stop times of MPO as seen from NEW_NORCIA
      using both ellipsoidal and DSK target shape models:
 
       Ell: 2027 JAN 03 00:00:00.000 TDB : 2027 JAN 03 01:28:03.419 TDB
       DSK: 2027 JAN 03 00:00:00.000 TDB : 2027 JAN 03 01:28:03.202 TDB
 
       Ell: 2027 JAN 03 02:00:42.993 TDB : 2027 JAN 03 03:49:50.750 TDB
       DSK: 2027 JAN 03 02:00:43.226 TDB : 2027 JAN 03 03:49:50.589 TDB
 
       Ell: 2027 JAN 03 04:22:20.803 TDB : 2027 JAN 03 06:11:38.050 TDB
       DSK: 2027 JAN 03 04:22:21.167 TDB : 2027 JAN 03 06:11:37.927 TDB
 
       Ell: 2027 JAN 03 06:43:58.528 TDB : 2027 JAN 03 08:33:25.499 TDB
       DSK: 2027 JAN 03 06:43:58.803 TDB : 2027 JAN 03 08:33:25.452 TDB
 
       Ell: 2027 JAN 03 09:05:36.070 TDB : 2027 JAN 03 10:55:12.991 TDB
       DSK: 2027 JAN 03 09:05:36.483 TDB : 2027 JAN 03 10:55:13.005 TDB
 
       Ell: 2027 JAN 03 21:55:08.488 TDB : 2027 JAN 03 22:44:11.490 TDB
       DSK: 2027 JAN 03 21:55:08.488 TDB : 2027 JAN 03 22:44:11.836 TDB
 
       Ell: 2027 JAN 03 23:15:19.552 TDB : 2027 JAN 04 01:05:59.339 TDB
       DSK: 2027 JAN 03 23:15:20.564 TDB : 2027 JAN 04 01:05:59.788 TDB
 
       Ell: 2027 JAN 04 01:36:56.572 TDB : 2027 JAN 04 03:27:47.253 TDB
       DSK: 2027 JAN 04 01:36:56.903 TDB : 2027 JAN 04 03:27:47.794 TDB
 
       Ell: 2027 JAN 04 03:58:33.411 TDB : 2027 JAN 04 05:49:35.238 TDB
       DSK: 2027 JAN 04 03:58:33.685 TDB : 2027 JAN 04 05:49:35.857 TDB
 
       Ell: 2027 JAN 04 06:20:10.165 TDB : 2027 JAN 04 08:11:23.310 TDB
       DSK: 2027 JAN 04 06:20:10.819 TDB : 2027 JAN 04 08:11:23.843 TDB
 
       Ell: 2027 JAN 04 08:41:46.813 TDB : 2027 JAN 04 10:33:11.480 TDB
       DSK: 2027 JAN 04 08:41:47.399 TDB : 2027 JAN 04 10:33:12.291 TDB
 
       Ell: 2027 JAN 04 21:58:41.333 TDB : 2027 JAN 04 22:22:13.911 TDB
       DSK: 2027 JAN 04 21:58:41.333 TDB : 2027 JAN 04 22:22:13.969 TDB
 
       Ell: 2027 JAN 04 22:51:24.368 TDB : 2027 JAN 05 00:44:02.576 TDB
       DSK: 2027 JAN 04 22:51:24.088 TDB : 2027 JAN 05 00:44:02.498 TDB
 
       Ell: 2027 JAN 05 01:13:00.256 TDB : 2027 JAN 05 03:05:51.406 TDB
       DSK: 2027 JAN 05 01:13:00.056 TDB : 2027 JAN 05 03:05:51.377 TDB
 
       Ell: 2027 JAN 05 03:34:36.025 TDB : 2027 JAN 05 05:27:40.260 TDB
       DSK: 2027 JAN 05 03:34:36.194 TDB : 2027 JAN 05 05:27:40.400 TDB
 
       Ell: 2027 JAN 05 05:56:11.727 TDB : 2027 JAN 05 07:49:29.298 TDB
       DSK: 2027 JAN 05 05:56:11.995 TDB : 2027 JAN 05 07:49:29.743 TDB
 
       Ell: 2027 JAN 05 08:17:47.213 TDB : 2027 JAN 05 10:11:18.377 TDB
       DSK: 2027 JAN 05 08:17:47.173 TDB : 2027 JAN 05 10:11:18.893 TDB
 
       Ell: 2027 JAN 05 10:39:22.574 TDB : 2027 JAN 05 11:04:00.020 TDB
       DSK: 2027 JAN 05 10:39:22.690 TDB : 2027 JAN 05 11:04:00.020 TDB
 
       Ell: 2027 JAN 05 22:27:17.183 TDB : 2027 JAN 06 00:00:00.000 TDB
       DSK: 2027 JAN 05 22:27:17.436 TDB : 2027 JAN 06 00:00:00.000 TDB
 
 
 
Extra Credit
===========================================================================
 
   In this ``extra credit'' section you will be presented with more complex
   tasks, aimed at improving your understanding of the geometry event
   finding subsystem and particularly the spiceypy.gfposc and
   spiceypy.gfdist 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 to the questions
   asked in these tasks.
 
 
Task statements
--------------------------------------------------------
 
       1.   Write a program that finds the times, within the time range
 
            2027 JAN 03 TDB
            2027 JAN 04 TDB
 
            when the BepiColombo Mercury Planetary Orbiter (MPO) crosses
            Mercury' equator. Display the results using TDB calendar dates
            and millisecond precision.
 
       2.   Write a program that finds the times, within the time range
 
            2027 JAN 03 TDB
            2027 JAN 04 TDB
 
            when the BepiColombo Mercury Planetary Orbiter (MPO) is at
            periapsis. Display the results using TDB calendar dates and
            millisecond precision.
 
       3.   Write a program that finds the times, within the time range
 
            2027 JAN 03 TDB
            2027 JAN 04 TDB
 
            when the BepiColombo Mercury Planetary Orbiter (MPO) is at
            apoapsis. Display the results using TDB calendar dates and
            millisecond precision.
 
 
Solutions
--------------------------------------------------------
 
       1.   Solution for the equator crossing search, using spiceypy.gfposc
            for the BepiColombo Mercury Planetary Orbiter (MPO) spacecraft
            latitude in the Mercury body-fixed frame equal to 0 degrees:
 
 
      Inputs for equator crossing search:
 
         Target                       = MPO
         Observer                     = MERCURY
         Observer's reference frame   = IAU_MERCURY
         Latitude limit (degrees)     = 0.000000
         Aberration correction        = NONE
         Step size (seconds)          = 300.000000
         Start time                   = 2027 JAN 03 00:00:00.000 (TDB)
         Stop time                    = 2027 JAN 04 00:00:00.000 (TDB)
 
      MPO MERCURY equator crossing times:
 
       Equator crossing or start time:  2027 JAN 03 00:21:02.744 (TDB)
       Equator crossing time:           2027 JAN 03 01:34:39.885 (TDB)
       Equator crossing time:           2027 JAN 03 02:42:45.780 (TDB)
       Equator crossing time:           2027 JAN 03 03:56:22.917 (TDB)
       Equator crossing time:           2027 JAN 03 05:04:28.849 (TDB)
       Equator crossing time:           2027 JAN 03 06:18:06.007 (TDB)
       Equator crossing time:           2027 JAN 03 07:26:11.888 (TDB)
       Equator crossing time:           2027 JAN 03 08:39:49.081 (TDB)
       Equator crossing time:           2027 JAN 03 09:47:54.901 (TDB)
       Equator crossing time:           2027 JAN 03 11:01:32.084 (TDB)
       Equator crossing time:           2027 JAN 03 12:09:37.930 (TDB)
       Equator crossing time:           2027 JAN 03 13:23:15.135 (TDB)
       Equator crossing time:           2027 JAN 03 14:31:20.939 (TDB)
       Equator crossing time:           2027 JAN 03 15:44:58.128 (TDB)
       Equator crossing time:           2027 JAN 03 16:53:03.959 (TDB)
       Equator crossing time:           2027 JAN 03 18:06:41.171 (TDB)
       Equator crossing time:           2027 JAN 03 19:14:46.962 (TDB)
       Equator crossing time:           2027 JAN 03 20:28:24.150 (TDB)
       Equator crossing time:           2027 JAN 03 21:36:29.935 (TDB)
       Equator crossing time:           2027 JAN 03 22:50:07.114 (TDB)
       Equator crossing or stop time:   2027 JAN 03 23:58:12.944 (TDB)
 
       2.   Solution for the periapsis search, using spiceypy.gfdist for
            the BepiColombo Mercury Planetary Orbiter (MPO) spacecraft
            distance from Mercury at a local minimum:
 
 
      Inputs for periapsis search:
 
         Target                       = MPO
         Observer                     = MERCURY
         Aberration correction        = NONE
         Step size (seconds)          = 300.000000
         Start time                   = 2027 JAN 03 00:00:00.000 (TDB)
         Stop time                    = 2027 JAN 04 00:00:00.000 (TDB)
 
      MPO periapsis times:
 
       Periapsis or start time:         2027 JAN 03 00:18:33.937 (TDB)
       Periapsis time:                  2027 JAN 03 02:40:16.998 (TDB)
       Periapsis time:                  2027 JAN 03 05:01:59.964 (TDB)
       Periapsis time:                  2027 JAN 03 07:23:43.026 (TDB)
       Periapsis time:                  2027 JAN 03 09:45:25.991 (TDB)
       Periapsis time:                  2027 JAN 03 12:07:09.042 (TDB)
       Periapsis time:                  2027 JAN 03 14:28:52.095 (TDB)
       Periapsis time:                  2027 JAN 03 16:50:35.082 (TDB)
       Periapsis time:                  2027 JAN 03 19:12:18.042 (TDB)
       Periapsis time:                  2027 JAN 03 21:34:01.097 (TDB)
       Periapsis or end time:           2027 JAN 03 23:55:44.079 (TDB)
 
       3.   Solution for the apoapsis search, using spiceypy.gfdist for the
            BepiColombo Mercury Planetary Orbiter (MPO) spacecraft distance
            from Mercury at a local maximum:
 
 
      Inputs for apoapsis search:
 
         Target                       = MPO
         Observer                     = MERCURY
         Aberration correction        = NONE
         Step size (seconds)          = 300.000000
         Start time                   = 2027 JAN 03 00:00:00.000 (TDB)
         Stop time                    = 2027 JAN 04 00:00:00.000 (TDB)
 
      MPO apoapsis times:
 
       Apoapsis or start time:          2027 JAN 03 01:29:25.529 (TDB)
       Apoapsis time:                   2027 JAN 03 03:51:08.495 (TDB)
       Apoapsis time:                   2027 JAN 03 06:12:51.561 (TDB)
       Apoapsis time:                   2027 JAN 03 08:34:34.611 (TDB)
       Apoapsis time:                   2027 JAN 03 10:56:17.595 (TDB)
       Apoapsis time:                   2027 JAN 03 13:18:00.653 (TDB)
       Apoapsis time:                   2027 JAN 03 15:39:43.611 (TDB)
       Apoapsis time:                   2027 JAN 03 18:01:26.677 (TDB)
       Apoapsis time:                   2027 JAN 03 20:23:09.638 (TDB)
       Apoapsis or stop time:           2027 JAN 03 22:44:52.618 (TDB)
 
