Main Page
In-situ Sensing Hands-On Lesson, using CASSINI (FORTRAN)

Table of Contents


   In-situ Sensing Hands-On Lesson, using CASSINI (FORTRAN)
      Overview
      Note About HTML Links
      References
         Tutorials
         Required Readings
         The Permuted Index
         API Documentation
      Kernels Used
      SPICE Modules Used

   Step-1: ``UTC to ET''
      ``UTC to ET'' Task Statement
      ``UTC to ET'' Hints
      ``UTC to ET'' Solution Steps
      ``UTC to ET'' Code

   Step-2: ``SCLK to ET''
      ``SCLK to ET'' Task Statement
      ``SCLK to ET'' Hints
      ``SCLK to ET'' Solution Steps
      ``SCLK to ET'' Code

   Step-3: ``Spacecraft State''
      ``Spacecraft State'' Task Statement
      ``Spacecraft State'' Hints
      ``Spacecraft State'' Solution Steps
      ``Spacecraft State'' Code

   Step-4: ``Sun Direction''
      ``Sun Direction'' Task Statement
      ``Sun Direction'' Hints
      ``Sun Direction'' Solution Steps
      ``Sun Direction'' Code

   Step-5: ``Sub-Spacecraft Point''
      ``Sub-Spacecraft Point'' Task Statement
      ``Sub-Spacecraft Point'' Hints
      ``Sub-Spacecraft Point'' Solution Steps
      ``Sub-Spacecraft Point'' Code

   Step-6: ``Spacecraft Velocity''
      ``Spacecraft Velocity'' Task Statement
      ``Spacecraft Velocity'' Hints
      ``Spacecraft Velocity'' Solution Steps
      ``Spacecraft Velocity'' Code Program ``scvel.f'':




Top

In-situ Sensing Hands-On Lesson, using CASSINI (FORTRAN)





March 01, 2023



Top

Overview




In this lesson you will develop a simple program illustrating how SPICE can be used to compute various kinds of geometry information applicable to the experiments carried out by an in-situ instrument flown on an interplanetary spacecraft.



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 ``toolkit/'' 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 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 ``toolkit/exe'' of Toolkit installation tree.



Top

Tutorials



The following SPICE tutorials serve as references for the discussions in this lesson:

 
   Name              Lesson steps/routines it describes
   ----------------  -----------------------------------------------
   Time              UTC to ET and SCLK to ET
   Loading Kernels   Loading SPICE kernels
   SCLK              SCLK to ET time conversion
   SPK               Computing positions and velocities
   Frames            Computing transformations between frames
These tutorials are available 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 ``toolkit/doc'' directory in the SPICE Toolkit installation tree.

   Name             Lesson steps/routines that it describes
   ---------------  -----------------------------------------
   kernel.req       Loading SPICE kernels
   naif_ids.req     Body and reference frame names
   spk.req          Computing positions and velocities
   sclk.req         SCLK to ET time conversion
   time.req         UTC to ET time conversion


Top

The Permuted Index



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

This text document provides a simple mechanism by which users can discover which SPICE routines perform functions of interest, as well as the names of the source files that contain these routines. It is particularly useful for FORTRAN programmers because some of the routines are entry points; the names of these routines do not translate directly into the name of the respective source files that contain them.



Top

API Documentation



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

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

   toolkit/src/spicelib/str2et.for
Since some of the FORTRAN routines are entry points they may be part of a source file that has different name. The ``Permuted Index'' document mentioned above can be used to locate the name of their source file.



Top

Kernels Used




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

   #  FILE NAME                 TYPE DESCRIPTION
   -- ------------------------- ---- -----------------------------------
   1  naif0008.tls              LSK  Generic LSK
   2  cas00084.tsc              SCLK Cassini SCLK
   3  030201AP_SK_SM546_T45.bsp SPK  Cassini Spacecraft SPK
   4  981005_PLTEPH-DE405S.bsp  SPK  Planetary Ephemeris SPK
   5  04135_04171pc_psiv2.bc    CK   Cassini Spacecraft CK
   6  cas_v37.tf                FK   Cassini FK
   7  020514_SE_SAT105.bsp      SPK  Saturnian Satellite Ephemeris SPK
   8  cpck05Mar2004.tpc         PCK  Cassini project PCK
These SPICE kernels are included in the lesson package.



Top

SPICE Modules Used




This section provides a complete list of the routines 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   ROUTINES   FUNCTIONS  KERNELS
   ------- ---------  ---------  ---------  ----------
      1    convrt     FURNSH                1
                      STR2ET
 
      2    sclket     FURNSH                1,2
                      STR2ET
                      SCS2E
 
      3    getsta     FURNSH                1-4
                      STR2ET
                      SCS2E
                      SPKEZR
 
      4    soldir     FURNSH                1-6
                      STR2ET
                      SCS2E
                      SPKEZR
                      SPKPOS
                      VHATIP
 
      5    sscpnt     FURNSH     DPR        1-8
                      STR2ET
                      SCS2E
                      SPKEZR
                      SPKPOS
                      VHATIP
                      SUBPNT
                      RECLAT
                      PXFORM
                      MXV
 
      6    scvel      FURNSH     DPR        1-8
                      STR2ET
                      SCS2E
                      SPKEZR
                      SPKPOS
                      VHATIP
                      SUBPNT
                      RECLAT
                      PXFORM
                      MXV
                      VPACK
                      VHAT
Refer to the headers of the various routines listed above, as detailed interface specifications are provided with the source code.



Top

Step-1: ``UTC to ET''







Top

``UTC to ET'' Task Statement




Write a program that computes and prints the Ephemeris Time (ET) corresponding to ``2004-06-11T19:32:00'' UTC, as the number of ephemeris seconds past J2000, .



Top

``UTC to ET'' Hints




Find out what SPICE kernel(s) is(are) needed to support this conversion. Reference the ``time.req'' and/or ``Time'' tutorial.

Find out what routine should be called to load necessary kernel(s). Reference the ``kernel.req'' and/or ``Loading Kernels'' tutorial.

Find the ``loader'' routine calling sequence specification. Look at the ``time.req'' and that routine's source code header. This routine may be an entry point, in which case there will be no source file with the same name. To find out in which source file this entry point is, search for its name in the ``Permuted Index''.

Find the routine(s) used to convert time between UTC and ET. Look at the ``time.req'' and/or ``Time'' tutorial.

Find the ``converter'' routine(s) calling sequence specification. Look in the ``time.req'' and the routine's source code header.

Put all calls together in a program, add variable declarations (the routine header's ``Declarations'' and ``Examples'' sections are a good place to look for declaration specification and examples) and output print statements. Compile it and link it against SPICELIB.



Top

``UTC to ET'' Solution Steps




Only one kernel file is needed to support this conversion -- an LSK file ``naif0008.tls''.

As any other SPICE kernel this file can be loaded by the FURNSH routine (entry in ``toolkit/src/spicelib/keeper.f''). For that, the name of the file can be provided as a sole argument of this routine:

   CHARACTER*(256)      LSKFLE
   ...
   LSKFLE = 'naif0008.tls'
   ...
   CALL FURNSH( LSKFLE )
or it can be listed in a meta-kernel:

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        )
   \begintext
the name of which, let's call it ``convrt.tm'', can be then provided as a sole argument of the FURNSH routine:

         CHARACTER*(256)      MKFILE
 
         ...
 
         MKFILE = 'convrt.tm'
         CALL FURNSH( MKFILE )
While the second option seems to involve a bit more work -- it requires making an extra file -- it is a much better way to go if you plan to load more kernels as you extend the program. With the meta-kernel approach simply adding more kernels to the list in KERNEL_TO_LOAD without changing the program code will accomplish that.

The highest level SPICE time routine converting UTC to ET is STR2ET (``toolkit/src/spicelib/str2et.f'').

It has two arguments -- input time string representing UTC in a variety of formats (see STR2ET header's section ``Particulars'' for the complete description of input time formats) and output DP number of ET seconds past J2000. A call to STR2ET converting a given UTC to ET could look like this:

         CHARACTER*(32)       UTC
         DOUBLE PRECISION     ET
 
         ...
 
         UTC = '2004-06-11T19:32:00'
         CALL STR2ET( UTC, ET )
By combining FURNSH and STR2ET calls and required declarations and by adding a simple print statement, one would get a complete program that prints ET for the given UTC epoch.

The program's source code then needs to be compiled and linked against SPICELIB.

When you run the program's executable, ``convrt'', it produces the following output:

   > ./convrt
   UTC       = 2004-06-11T19:32:00
   ET        =     140254384.184625


Top

``UTC to ET'' Code




Program ``convrt.f'':

         PROGRAM CONVRT
         IMPLICIT NONE
 
         CHARACTER*(256)      MKFILE
         CHARACTER*(32)       UTC
         DOUBLE PRECISION     ET
 
         MKFILE = 'convrt.tm'
         CALL FURNSH( MKFILE )
 
         UTC = '2004-06-11T19:32:00'
         CALL STR2ET( UTC, ET )
 
         WRITE(*,'(2A)')      'UTC       = ', UTC
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         END
Meta-kernel file ``convrt.tm'':

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        )
   \begintext


Top

Step-2: ``SCLK to ET''







Top

``SCLK to ET'' Task Statement




Extend the program from Step-1 to compute and print ET for the following CASSINI on-board clock epoch ``1465674964.105''.



Top

``SCLK to ET'' Hints




Find out what additional (to those already loaded in Step-1) SPICE kernel(s) is(are) needed to support SCLK to ET conversion. Look at the ``sclk.req'' and/or ``SCLK'' tutorial.

Modify the program or meta-kernel to load this(these) kernels.

Find the routine(s) needed to convert time between SCLK and ET. Look at the ``sclk.req'' and/or ``Time'' and ``SCLK'' tutorials.

Find the ``converter'' routine's calling sequence specification. Look in the ``sclk.req'' and the routine's source code header.

Look at ``naif_ids.req'' and the comments in the additional kernel(s) that you have loaded for information on proper values of input arguments of this routine.

Add calls to the ``converter'' routine(s), necessary variable declarations (the routine header's ``Declarations'' and ``Examples'' sections are a good place to look for declaration specification and examples), and output print statements to the program. Re-compile and re-link it against SPICELIB.



Top

``SCLK to ET'' Solution Steps




A CASSINI SCLK file is needed additionally to the LSK file loaded in the Step-1 to support this conversion.

No code change is needed in the loading portion of the program if a meta-kernel approach was used in the Step-1. The program will load the file if it will be added to the list of kernels in the KERNELS_TO_LOAD variable:

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
      cas00084.tsc                Cassini SCLK.
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        'kernels/sclk/cas00084.tsc'
                        )
   \begintext
The highest level SPICE routine converting SCLK to ET is SCS2E (``toolkit/src/spicelib/scs2e.f'').

It has three arguments -- NAIF ID for CASSINI s/c (-82 as described by ``naif_ids.req'' document), input time string representing CASSINI SCLK, and output DP number of ET seconds past J2000. A call to STR2ET converting given SCLK to ET could look like this:

         CHARACTER*(32)       SCLK
 
         ...
 
         INTEGER              SCID
 
         ...
 
         SCID = -82
         SCLK = '1465674964.105'
         CALL SCS2E( SCID, SCLK, ET )
By adding the SCS2E call, required declarations and a simple print statement, one would get a complete program that prints ET for the given SCLK epoch.

When you run the program's executable, ``sclket'', it produces the following output:

   > ./sclket
   UTC       = 2004-06-11T19:32:00
   ET        =     140254384.184625
   SCLK      = 1465674964.105
   ET        =     140254384.183426


Top

``SCLK to ET'' Code




Program ``sclket.f'':

         PROGRAM SCLKET
         IMPLICIT NONE
 
         CHARACTER*(256)      MKFILE
         CHARACTER*(32)       UTC
         CHARACTER*(32)       SCLK
         DOUBLE PRECISION     ET
         INTEGER              SCID
 
 
         MKFILE = 'sclket.tm'
         CALL FURNSH( MKFILE )
 
         UTC = '2004-06-11T19:32:00'
         CALL STR2ET( UTC, ET )
 
         WRITE(*,'(2A)')      'UTC       = ', UTC
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         SCID = -82
         SCLK = '1465674964.105'
         CALL SCS2E( SCID, SCLK, ET )
 
         WRITE(*,'(2A)')      'SCLK      = ', SCLK
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         END
Meta-kernel file ``sclket.tm'':

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
      cas00084.tsc                Cassini SCLK.
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        'kernels/sclk/cas00084.tsc'
                        )
   \begintext


Top

Step-3: ``Spacecraft State''







Top

``Spacecraft State'' Task Statement




Extend the program from Step-2 to compute geometric state -- position and velocity -- of the CASSINI spacecraft with respect to the Sun in the Ecliptic frame at the epoch specified by SCLK time from Step-2.



Top

``Spacecraft State'' Hints




Find out what additional (to those already loaded in Steps-1&2) SPICE kernel(s) is(are) needed to support state computation. Look at the ``spk.req'' and/or ``SPK'' tutorial.

Verify that the kernels contain enough data to compute the state of interest. Use ``brief'' utility program located under ``toolkit/exe'' directory for that.

Modify the meta-kernel to load this(these) kernels.

Determine the routine(s) needed to compute states. Look at the ``spk.req'' and/or ``SPK'' tutorial presentation.

Find the the routine(s) calling sequence specification. Look in the ``spk.req'' and the routine's source code header.

Reference the ``naif_ids.req'' and ``frames.req'' and the routine(s) header ``Inputs'' and ``Particulars'' sections to determine proper values of the input arguments of this routine.

Add calls to the routine(s), necessary variable declarations and output print statements to the program. Re-compile and re-link it against SPICELIB.



Top

``Spacecraft State'' Solution Steps




A CASSINI spacecraft trajectory SPK and generic planetary ephemeris SPK files are needed to support computation of the state of interest.

The file names can be added to the meta-kernel to get them loaded into the program:

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
      cas00084.tsc                Cassini SCLK.
      030201AP_SK_SM546_T45.bsp   Cassini Spacecraft SPK.
      981005_PLTEPH-DE405S.bsp    Planetary Ephemeris SPK.
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        'kernels/sclk/cas00084.tsc'
                        'kernels/spk/030201AP_SK_SM546_T45.bsp'
                        'kernels/spk/981005_PLTEPH-DE405S.bsp'
                        )
   \begintext
The highest level SPICE routine computing states is SPKEZR (``toolkit/src/spicelib/spkezr.f'').

We are interested in computing CASSINI position and velocity with respect to the Sun, therefore the target and observer names should be set to 'CASSINI' and 'Sun' (both names can be found in ``naif_ids.req'').

The state should be in ecliptic frame, therefore the name of the frame in which the state should be computed is 'ECLIPJ2000' (see ``frames.req'' document.)

Since we need only the geometric position, the ABCORR argument of the routine should be set to 'NONE' (see aberration correction discussion in the (``toolkit/src/spicelib/spkezr.f''). header).

Putting it all together, we get:

         CHARACTER*(32)       TARGET
         CHARACTER*(32)       FRAME
         CHARACTER*(32)       CORRTN
         CHARACTER*(32)       OBSERV
 
         ...
 
         DOUBLE PRECISION     STATE  ( 6 )
         DOUBLE PRECISION     LT
 
         ...
 
         TARGET = 'CASSINI'
         FRAME  = 'ECLIPJ2000'
         CORRTN = 'NONE'
         OBSERV = 'SUN'
         CALL SPKEZR( TARGET, ET, FRAME, CORRTN, OBSERV, STATE, LT )
When you run the program's executable, ``getsta'', it produces the following output:

   > ./getsta
   UTC       = 2004-06-11T19:32:00
   ET        =     140254384.184625
   SCLK      = 1465674964.105
   ET        =     140254384.183426
    X        =    -376599061.916539
    Y        =    1294487780.929154
    Z        =      -7064853.054698
   VX        =            -5.164226
   VY        =             0.801719
   VZ        =             0.040603


Top

``Spacecraft State'' Code




Program ``getsta.f'':

         PROGRAM GETSTA
         IMPLICIT NONE
 
         CHARACTER*(256)      MKFILE
         CHARACTER*(32)       UTC
         CHARACTER*(32)       SCLK
         CHARACTER*(32)       TARGET
         CHARACTER*(32)       FRAME
         CHARACTER*(32)       CORRTN
         CHARACTER*(32)       OBSERV
         DOUBLE PRECISION     ET
         DOUBLE PRECISION     STATE  ( 6 )
         DOUBLE PRECISION     LT
         INTEGER              SCID
 
 
         MKFILE = 'getsta.tm'
         CALL FURNSH( MKFILE )
 
         UTC = '2004-06-11T19:32:00'
         CALL STR2ET( UTC, ET )
 
         WRITE(*,'(2A)')      'UTC       = ', UTC
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         SCID = -82
         SCLK = '1465674964.105'
         CALL SCS2E( SCID, SCLK, ET )
 
         WRITE(*,'(2A)')      'SCLK      = ', SCLK
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         TARGET = 'CASSINI'
         FRAME  = 'ECLIPJ2000'
         CORRTN = 'NONE'
         OBSERV = 'SUN'
         CALL SPKEZR( TARGET, ET, FRAME, CORRTN, OBSERV, STATE, LT )
 
         WRITE(*,'(A,F20.6)') ' X        = ', STATE(1)
         WRITE(*,'(A,F20.6)') ' Y        = ', STATE(2)
         WRITE(*,'(A,F20.6)') ' Z        = ', STATE(3)
         WRITE(*,'(A,F20.6)') 'VX        = ', STATE(4)
         WRITE(*,'(A,F20.6)') 'VY        = ', STATE(5)
         WRITE(*,'(A,F20.6)') 'VZ        = ', STATE(6)
 
         END
Meta-kernel file ``getsta.tm'':

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
      cas00084.tsc                Cassini SCLK.
      030201AP_SK_SM546_T45.bsp   Cassini Spacecraft SPK.
      981005_PLTEPH-DE405S.bsp    Planetary Ephemeris SPK.
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        'kernels/sclk/cas00084.tsc'
                        'kernels/spk/030201AP_SK_SM546_T45.bsp'
                        'kernels/spk/981005_PLTEPH-DE405S.bsp'
                        )
   \begintext


Top

Step-4: ``Sun Direction''







Top

``Sun Direction'' Task Statement




Extend the program from Step-3 to compute apparent direction of the Sun in the INMS frame at the epoch specified by SCLK time from Step-2.



Top

``Sun Direction'' Hints




Determine the additional SPICE kernels needed to support the direction computation, knowing that they should provide the s/c and instrument frame orientation.

Verify that the orientation data in the kernels have adequate coverage to support computation of the direction of interest. Use ``ckbrief'' utility program located under ``toolkit/exe'' directory for that.

Modify the meta-kernel to load this(these) kernels.

Determine the proper input arguments for the SPKPOS call to calculate the direction (which is the position portion of the output state). Look through the Frames Kernel find the name of the frame to used.

Add calls to the routine(s), necessary variable declarations and output print statements to the program. Re-compile and re-link it against SPICELIB.



Top

``Sun Direction'' Solution Steps




A CASSINI spacecraft orientation CK file, providing s/c orientation with respect to an inertial frame, and CASSINI FK file, providing orientation of the INMS frame with respect to the s/c frame, are needed additionally to already loaded kernels to support computation of this direction.

The file names can be added to the meta-kernel to get them loaded into the program:

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
      cas00084.tsc                Cassini SCLK.
      030201AP_SK_SM546_T45.bsp   Cassini Spacecraft SPK.
      981005_PLTEPH-DE405S.bsp    Planetary Ephemeris SPK.
      04135_04171pc_psiv2.bc      Cassini Spacecraft CK.
      cas_v37.tf                  Cassini FK.
 
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        'kernels/sclk/cas00084.tsc'
                        'kernels/spk/030201AP_SK_SM546_T45.bsp'
                        'kernels/spk/981005_PLTEPH-DE405S.bsp'
                        'kernels/ck/04135_04171pc_psiv2.bc'
                        'kernels/fk/cas_v37.tf'
                        )
   \begintext
The same highest level SPICE routine computing positions, SPKPOS, can be used to compute this direction.

Since this is the direction of the Sun as seen from the s/c, the target argument should be set to 'Sun' and the observer argument should be set to 'CASSINI'. The name of the INMS frame is 'CASSINI_INMS', the definition and description of this frame are provided in the CASSINI FK file, ``cas_v37.tf''.

Since the apparent, or ``as seen'', position is sought for, the ABCORR argument of the routine should be set to 'LT+S' (see aberration correction discussion in the (``toolkit/src/spicelib/spkpos.f'') header).

If desired, the position can then be turned into a unit vector using VHATIP routine (``toolkit/src/spicelib/vhatip.f''). Putting it all together, we get:

         DOUBLE PRECISION     SUNDIR ( 3 )
 
         ...
 
         TARGET = 'SUN'
         FRAME  = 'CASSINI_INMS'
         CORRTN = 'LT+S'
         OBSERV = 'CASSINI'
         CALL SPKPOS( TARGET, ET, FRAME, CORRTN, OBSERV, SUNDIR, LT )
         CALL VHATIP( SUNDIR )
When you run the program's executable, ``soldir'', it produces the following output:

   > ./soldir
   UTC       = 2004-06-11T19:32:00
   ET        =     140254384.184625
   SCLK      = 1465674964.105
   ET        =     140254384.183426
    X        =    -376599061.916539
    Y        =    1294487780.929154
    Z        =      -7064853.054698
   VX        =            -5.164226
   VY        =             0.801719
   VZ        =             0.040603
   SUNDIR(X) =            -0.290204
   SUNDIR(Y) =             0.881631
   SUNDIR(Z) =             0.372167


Top

``Sun Direction'' Code




Program ``soldir.f'':

         PROGRAM SOLDIR
         IMPLICIT NONE
 
         CHARACTER*(256)      MKFILE
         CHARACTER*(32)       UTC
         CHARACTER*(32)       SCLK
         CHARACTER*(32)       TARGET
         CHARACTER*(32)       FRAME
         CHARACTER*(32)       CORRTN
         CHARACTER*(32)       OBSERV
         DOUBLE PRECISION     ET
         DOUBLE PRECISION     STATE  ( 6 )
         DOUBLE PRECISION     SUNDIR ( 3 )
         DOUBLE PRECISION     LT
         INTEGER              SCID
 
 
         MKFILE = 'soldir.tm'
         CALL FURNSH( MKFILE )
 
         UTC = '2004-06-11T19:32:00'
         CALL STR2ET( UTC, ET )
 
         WRITE(*,'(2A)')      'UTC       = ', UTC
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         SCID = -82
         SCLK = '1465674964.105'
         CALL SCS2E( SCID, SCLK, ET )
 
         WRITE(*,'(2A)')      'SCLK      = ', SCLK
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         TARGET = 'CASSINI'
         FRAME  = 'ECLIPJ2000'
         CORRTN = 'NONE'
         OBSERV = 'SUN'
         CALL SPKEZR( TARGET, ET, FRAME, CORRTN, OBSERV, STATE, LT )
 
         WRITE(*,'(A,F20.6)') ' X        = ', STATE(1)
         WRITE(*,'(A,F20.6)') ' Y        = ', STATE(2)
         WRITE(*,'(A,F20.6)') ' Z        = ', STATE(3)
         WRITE(*,'(A,F20.6)') 'VX        = ', STATE(4)
         WRITE(*,'(A,F20.6)') 'VY        = ', STATE(5)
         WRITE(*,'(A,F20.6)') 'VZ        = ', STATE(6)
 
         TARGET = 'SUN'
         FRAME  = 'CASSINI_INMS'
         CORRTN = 'LT+S'
         OBSERV = 'CASSINI'
         CALL SPKPOS( TARGET, ET, FRAME, CORRTN, OBSERV, SUNDIR, LT )
         CALL VHATIP( SUNDIR )
 
         WRITE(*,'(A,F20.6)') 'SUNDIR(X) = ', SUNDIR(1)
         WRITE(*,'(A,F20.6)') 'SUNDIR(Y) = ', SUNDIR(2)
         WRITE(*,'(A,F20.6)') 'SUNDIR(Z) = ', SUNDIR(3)
 
         END
Meta-kernel file ``soldir.tm'':

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
      cas00084.tsc                Cassini SCLK.
      030201AP_SK_SM546_T45.bsp   Cassini Spacecraft SPK.
      981005_PLTEPH-DE405S.bsp    Planetary Ephemeris SPK.
      04135_04171pc_psiv2.bc      Cassini Spacecraft CK.
      cas_v37.tf                  Cassini FK.
 
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        'kernels/sclk/cas00084.tsc'
                        'kernels/spk/030201AP_SK_SM546_T45.bsp'
                        'kernels/spk/981005_PLTEPH-DE405S.bsp'
                        'kernels/ck/04135_04171pc_psiv2.bc'
                        'kernels/fk/cas_v37.tf'
                        )
   \begintext


Top

Step-5: ``Sub-Spacecraft Point''







Top

``Sub-Spacecraft Point'' Task Statement




Extend the program from Step-4 to compute planetocentric longitude and and latitude of the sub-spacecraft point on Phoebe, and the direction from the spacecraft to that point in the INMS frame.



Top

``Sub-Spacecraft Point'' Hints




Find the SPICE routine that computes sub-observer point coordinates. Use ``Most Used SPICE APIs'' or ``subpt'' cookbook program for that.

Refer to the routine's header to determine the additional kernels needed for this direction computation. Modify the meta-kernel to load this(these) kernels.

Determine the proper input arguments for the routine. Refer to the routine's header for that information.

Convert the surface point Cartesian vector returned by this routine to latitudinal coordinates. Use ``Permuted Index'' to find the routine that does this conversion. Refer to the routine's header for input/output argument specifications.

Since the Cartesian vector from the spacecraft to the sub-spacecraft point is computed in the Phoebe body-fixed frame, it should be transformed into the instrument frame get the direction we are looking for. Refer to ``frames.req'' and/or ``Frames'' tutorial to determine the name of the routine computing transformations and use it to compute transformation from Phoebe body-fixed to the INMS frame.

Using ``Permuted Index'' find the routine that multiplies 3x3 matrix by 3d vector and use it to rotate the vector to the instrument frame.

Add calls to the routine(s), necessary variable declarations and output print statements to the program. Re-compile and re-link it against SPICELIB.



Top

``Sub-Spacecraft Point'' Solution Steps




The SUBPNT routine (``toolkit/src/spicelib/subpnt.f'') can be used to compute the sub-observer point and the vector from the observer to that point with a single call. To determine this point as the closest point on the Phoebe ellipsoid, the METHOD argument has to be set to 'NEAR POINT: ELLIPSOID'. For our case the TARGET is 'PHOEBE', the target body-fixed frame is 'IAU_PHOEBE', and the observer is 'CASSINI'.

Since the s/c is close to Phoebe, light time does not need to be taken into account and, therefore, the ABCORR argument can be set to 'NONE'.

In order for SUBPNT to compute the nearest point location, a PCK file containing Phoebe radii has to be loaded into the program (see ``Files'' section of the routine's header.) All other files required for this computation are already being loaded by the program. With PCK file name added to it, the updated meta-kernel will look like this:

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
      cas00084.tsc                Cassini SCLK.
      020514_SE_SAT105.bsp        Saturnian Satellite Ephemeris SPK.
      030201AP_SK_SM546_T45.bsp   Cassini Spacecraft SPK.
      981005_PLTEPH-DE405S.bsp    Planetary Ephemeris SPK.
      04135_04171pc_psiv2.bc      Cassini Spacecraft CK.
      cas_v37.tf                  Cassini FK.
      cpck05Mar2004.tpc           Cassini project PCK.
 
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        'kernels/sclk/cas00084.tsc'
                        'kernels/spk/020514_SE_SAT105.bsp'
                        'kernels/spk/030201AP_SK_SM546_T45.bsp'
                        'kernels/spk/981005_PLTEPH-DE405S.bsp'
                        'kernels/ck/04135_04171pc_psiv2.bc'
                        'kernels/fk/cas_v37.tf'
                        'kernels/pck/cpck05Mar2004.tpc'
                        )
   \begintext
The sub-spacecraft point Cartesian vector can be converted to planetocentric radius, longitude and latitude using the RECLAT routine (``toolkit/src/spicelib/reclat.f'').

The vector from the spacecraft to the sub-spacecraft point returned by SUBPNT has to be rotated from the body-fixed frame to the instrument frame. The name of the routine that computes 3x3 matrices rotating vectors from one frame to another is PXFORM (``toolkit/src/spicelib/pxform.f'').

In our case the FROM argument should be set to 'IAU_PHOEBE' and the TO argument should be set to 'CASSINI_INMS'

The vector should be then multiplied by this matrix to rotate it to the instrument frame. The MXV routine performs that function (``toolkit/src/spicelib/mxv.f'')

After applying the rotation, normalize the resultant vector using the VHATIP routine.

For output the longitude and latitude angles returned by RECLAT in radians can be converted to degrees by multiplying by DPR function (``toolkit/src/spicelib/dpr.f'').

Putting it all together, we get:

   C
   C     SPICELIB functions
   C
         DOUBLE PRECISION     DPR
 
         ...
 
         CHARACTER*(32)       METHOD
         CHARACTER*(32)       FROMFR
         CHARACTER*(32)       TOFR
 
         ...
 
         DOUBLE PRECISION     SPOINT ( 3 )
         DOUBLE PRECISION     TRGEPC
         DOUBLE PRECISION     SRFVEC ( 3 )
         DOUBLE PRECISION     SRAD
         DOUBLE PRECISION     SLON
         DOUBLE PRECISION     SLAT
         DOUBLE PRECISION     SBPDIR ( 3 )
         DOUBLE PRECISION     M2IMAT ( 3, 3 )
 
         ...
 
         METHOD = 'NEAR POINT: ELLIPSOID'
         TARGET = 'PHOEBE'
         FRAME  = 'IAU_PHOEBE'
         CORRTN = 'NONE'
         OBSERV = 'CASSINI'
         CALL SUBPNT ( METHOD, TARGET, ET, FRAME, CORRTN, OBSERV,
        .              SPOINT, TRGEPC, SRFVEC )
 
         CALL RECLAT( SPOINT, SRAD, SLON, SLAT )
 
         FROMFR = 'IAU_PHOEBE'
         TOFR   = 'CASSINI_INMS'
         CALL PXFORM( FROMFR, TOFR, ET, M2IMAT )
 
         CALL MXV   ( M2IMAT, SRFVEC, SBPDIR )
         CALL VHATIP( SBPDIR )
 
         WRITE(*,'(A,F20.6)') 'LON       = ', SLON*DPR()
         WRITE(*,'(A,F20.6)') 'LAT       = ', SLAT*DPR()
When you run the program's executable, ``sscpnt'', it produces the following output:

   > ./sscpnt
   UTC       = 2004-06-11T19:32:00
   ET        =     140254384.184625
   SCLK      = 1465674964.105
   ET        =     140254384.183426
    X        =    -376599061.916539
    Y        =    1294487780.929154
    Z        =      -7064853.054698
   VX        =            -5.164226
   VY        =             0.801719
   VZ        =             0.040603
   SUNDIR(X) =            -0.290204
   SUNDIR(Y) =             0.881631
   SUNDIR(Z) =             0.372167
   LON       =            23.423158
   LAT       =             3.709797
   SBPDIR(X) =            -0.000776
   SBPDIR(Y) =            -0.999873
   SBPDIR(Z) =            -0.015905


Top

``Sub-Spacecraft Point'' Code




Program ``sccpnt.f'':

         PROGRAM SSCPNT
         IMPLICIT NONE
 
   C
   C     SPICELIB functions
   C
         DOUBLE PRECISION     DPR
 
         CHARACTER*(256)      MKFILE
         CHARACTER*(32)       UTC
         CHARACTER*(32)       SCLK
         CHARACTER*(32)       TARGET
         CHARACTER*(32)       FRAME
         CHARACTER*(32)       CORRTN
         CHARACTER*(32)       OBSERV
         CHARACTER*(32)       METHOD
         CHARACTER*(32)       FROMFR
         CHARACTER*(32)       TOFR
         DOUBLE PRECISION     ET
         DOUBLE PRECISION     STATE  ( 6 )
         DOUBLE PRECISION     SUNDIR ( 3 )
         DOUBLE PRECISION     LT
         DOUBLE PRECISION     SPOINT ( 3 )
         DOUBLE PRECISION     TRGEPC
         DOUBLE PRECISION     SRFVEC ( 3 )
         DOUBLE PRECISION     SRAD
         DOUBLE PRECISION     SLON
         DOUBLE PRECISION     SLAT
         DOUBLE PRECISION     SBPDIR ( 3 )
         DOUBLE PRECISION     M2IMAT ( 3, 3 )
         INTEGER              SCID
 
 
         MKFILE = 'sscpnt.tm'
         CALL FURNSH( MKFILE )
 
         UTC = '2004-06-11T19:32:00'
         CALL STR2ET( UTC, ET )
 
         WRITE(*,'(2A)')      'UTC       = ', UTC
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         SCID = -82
         SCLK = '1465674964.105'
         CALL SCS2E( SCID, SCLK, ET )
 
         WRITE(*,'(2A)')      'SCLK      = ', SCLK
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         TARGET = 'CASSINI'
         FRAME  = 'ECLIPJ2000'
         CORRTN = 'NONE'
         OBSERV = 'SUN'
         CALL SPKEZR( TARGET, ET, FRAME, CORRTN, OBSERV, STATE, LT )
 
         WRITE(*,'(A,F20.6)') ' X        = ', STATE(1)
         WRITE(*,'(A,F20.6)') ' Y        = ', STATE(2)
         WRITE(*,'(A,F20.6)') ' Z        = ', STATE(3)
         WRITE(*,'(A,F20.6)') 'VX        = ', STATE(4)
         WRITE(*,'(A,F20.6)') 'VY        = ', STATE(5)
         WRITE(*,'(A,F20.6)') 'VZ        = ', STATE(6)
 
         TARGET = 'SUN'
         FRAME  = 'CASSINI_INMS'
         CORRTN = 'LT+S'
         OBSERV = 'CASSINI'
         CALL SPKPOS( TARGET, ET, FRAME, CORRTN, OBSERV, SUNDIR, LT )
         CALL VHATIP( SUNDIR )
 
         WRITE(*,'(A,F20.6)') 'SUNDIR(X) = ', SUNDIR(1)
         WRITE(*,'(A,F20.6)') 'SUNDIR(Y) = ', SUNDIR(2)
         WRITE(*,'(A,F20.6)') 'SUNDIR(Z) = ', SUNDIR(3)
 
         METHOD = 'NEAR POINT: ELLIPSOID'
         TARGET = 'PHOEBE'
         FRAME  = 'IAU_PHOEBE'
         CORRTN = 'NONE'
         OBSERV = 'CASSINI'
         CALL SUBPNT ( METHOD, TARGET, ET, FRAME, CORRTN, OBSERV,
        .              SPOINT, TRGEPC, SRFVEC )
 
         CALL RECLAT( SPOINT, SRAD, SLON, SLAT )
 
         FROMFR = 'IAU_PHOEBE'
         TOFR   = 'CASSINI_INMS'
         CALL PXFORM( FROMFR, TOFR, ET, M2IMAT )
 
         CALL MXV   ( M2IMAT, SRFVEC, SBPDIR )
         CALL VHATIP( SBPDIR )
 
         WRITE(*,'(A,F20.6)') 'LON       = ', SLON*DPR()
         WRITE(*,'(A,F20.6)') 'LAT       = ', SLAT*DPR()
         WRITE(*,'(A,F20.6)') 'SBPDIR(X) = ', SBPDIR(1)
         WRITE(*,'(A,F20.6)') 'SBPDIR(Y) = ', SBPDIR(2)
         WRITE(*,'(A,F20.6)') 'SBPDIR(Z) = ', SBPDIR(3)
 
         END
Meta-kernel file ``sscpnt.tm'':

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
      cas00084.tsc                Cassini SCLK.
      020514_SE_SAT105.bsp        Saturnian Satellite Ephemeris SPK.
      030201AP_SK_SM546_T45.bsp   Cassini Spacecraft SPK.
      981005_PLTEPH-DE405S.bsp    Planetary Ephemeris SPK.
      04135_04171pc_psiv2.bc      Cassini Spacecraft CK.
      cas_v37.tf                  Cassini FK.
      cpck05Mar2004.tpc           Cassini project PCK.
 
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        'kernels/sclk/cas00084.tsc'
                        'kernels/spk/020514_SE_SAT105.bsp'
                        'kernels/spk/030201AP_SK_SM546_T45.bsp'
                        'kernels/spk/981005_PLTEPH-DE405S.bsp'
                        'kernels/ck/04135_04171pc_psiv2.bc'
                        'kernels/fk/cas_v37.tf'
                        'kernels/pck/cpck05Mar2004.tpc'
                        )
   \begintext


Top

Step-6: ``Spacecraft Velocity''







Top

``Spacecraft Velocity'' Task Statement




Extend the program from Step-5 to compute the spacecraft velocity with respect to Phoebe in the INMS frame.



Top

``Spacecraft Velocity'' Hints




Compute velocity of the spacecraft with respect to Phoebe in some inertial frame, for example J2000. Recall that velocity is the last three components of the state vector returned by SPKEZR.

Since the velocity vector is computed in the inertial frame, it should be rotated to the instrument frame. Look at the previous step the routine that compute necessary rotation and rotate vectors.

Add calls to the routine(s), necessary variable declarations and output print statements to the program. Re-compile and re-link it against SPICELIB.



Top

``Spacecraft Velocity'' Solution Steps




All kernels required for computations in this step are already being loaded by the program, therefore, the meta-kernel does not need to be changed.

The spacecraft velocity vector is the last three components of the state returned by SPKEZR. To compute velocity of CASSINI with respect to Phoebe in the J2000 inertial frame the SPKEZR arguments should be set to 'CASSINI' (TARG), 'PHOEBE' (OBS), 'J2000' (REF) and 'NONE' (ABCORR).

For convenience the velocity can be copied from the output state in to a 3d vector using the VPACK routine (``toolkit/src/spicelib/vpack.f'').

The computed velocity vector has to be rotated from the J2000 frame to the instrument frame. The PXFORM routine used in the previous step can be used to compute the rotation matrix needed for that. In this case the frame name arguments should be set to 'J2000' (FROM) and 'CASSINI_INMS' (TO).

As in the previous step the difference vector should be then multiplied by this rotation matrix using the MXV routine. After applying the rotation, normalize the resultant vector using the VHAT routine.

Putting it all together, we get:

         DOUBLE PRECISION     SCVDIR ( 3 )
         DOUBLE PRECISION     TMPDIR ( 3 )
         DOUBLE PRECISION     J2IMAT ( 3, 3 )
 
         ...
 
         TARGET = 'CASSINI'
         FRAME  = 'J2000'
         CORRTN = 'NONE'
         OBSERV = 'PHOEBE'
         CALL SPKEZR( TARGET, ET, FRAME, CORRTN, OBSERV, STATE, LT )
         CALL VPACK ( STATE(4), STATE(5), STATE(6), SCVDIR )
 
         FROMFR = 'J2000'
         TOFR   = 'CASSINI_INMS'
         CALL PXFORM( FROMFR, TOFR, ET, J2IMAT )
 
         CALL MXV   ( J2IMAT, SCVDIR, TMPDIR )
         CALL VHAT  ( TMPDIR, SCVDIR )
When you run the program's executable, ``scvel'', it produces the following output:

   > ./scvel
   UTC       = 2004-06-11T19:32:00
   ET        =     140254384.184625
   SCLK      = 1465674964.105
   ET        =     140254384.183426
    X        =    -376599061.916539
    Y        =    1294487780.929154
    Z        =      -7064853.054698
   VX        =            -5.164226
   VY        =             0.801719
   VZ        =             0.040603
   SUNDIR(X) =            -0.290204
   SUNDIR(Y) =             0.881631
   SUNDIR(Z) =             0.372167
   LON       =            23.423158
   LAT       =             3.709797
   SBPDIR(X) =            -0.000776
   SBPDIR(Y) =            -0.999873
   SBPDIR(Z) =            -0.015905
   SCVDIR(X) =             0.395785
   SCVDIR(Y) =            -0.292808
   SCVDIR(Z) =             0.870413
Note that computing the spacecraft velocity in the instrument frame by a single call to SPKEZR by specifying 'CASSINI_INMS' in the REF argument returns an incorrect result. Such computation will take into account the spacecraft angular velocity from the CK files, which should not be considered in this case.



Top

``Spacecraft Velocity'' Code Program ``scvel.f'':




         PROGRAM SCVEL
         IMPLICIT NONE
 
   C
   C     SPICELIB functions
   C
         DOUBLE PRECISION     DPR
 
         CHARACTER*(256)      MKFILE
         CHARACTER*(32)       UTC
         CHARACTER*(32)       SCLK
         CHARACTER*(32)       TARGET
         CHARACTER*(32)       FRAME
         CHARACTER*(32)       CORRTN
         CHARACTER*(32)       OBSERV
         CHARACTER*(32)       METHOD
         CHARACTER*(32)       FROMFR
         CHARACTER*(32)       TOFR
         DOUBLE PRECISION     ET
         DOUBLE PRECISION     STATE  ( 6 )
         DOUBLE PRECISION     SUNDIR ( 3 )
         DOUBLE PRECISION     LT
         DOUBLE PRECISION     SPOINT ( 3 )
         DOUBLE PRECISION     TRGEPC
         DOUBLE PRECISION     SRFVEC ( 3 )
         DOUBLE PRECISION     SRAD
         DOUBLE PRECISION     SLON
         DOUBLE PRECISION     SLAT
         DOUBLE PRECISION     SBPDIR ( 3 )
         DOUBLE PRECISION     M2IMAT ( 3, 3 )
         DOUBLE PRECISION     SCVDIR ( 3 )
         DOUBLE PRECISION     TMPDIR ( 3 )
         DOUBLE PRECISION     J2IMAT ( 3, 3 )
         INTEGER              SCID
 
 
         MKFILE = 'scvel.tm'
         CALL FURNSH( MKFILE )
 
         UTC = '2004-06-11T19:32:00'
         CALL STR2ET( UTC, ET )
 
         WRITE(*,'(2A)')      'UTC       = ', UTC
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         SCID = -82
         SCLK = '1465674964.105'
         CALL SCS2E( SCID, SCLK, ET )
 
         WRITE(*,'(2A)')      'SCLK      = ', SCLK
         WRITE(*,'(A,F20.6)') 'ET        = ', ET
 
         TARGET = 'CASSINI'
         FRAME  = 'ECLIPJ2000'
         CORRTN = 'NONE'
         OBSERV = 'SUN'
         CALL SPKEZR( TARGET, ET, FRAME, CORRTN, OBSERV, STATE, LT )
 
         WRITE(*,'(A,F20.6)') ' X        = ', STATE(1)
         WRITE(*,'(A,F20.6)') ' Y        = ', STATE(2)
         WRITE(*,'(A,F20.6)') ' Z        = ', STATE(3)
         WRITE(*,'(A,F20.6)') 'VX        = ', STATE(4)
         WRITE(*,'(A,F20.6)') 'VY        = ', STATE(5)
         WRITE(*,'(A,F20.6)') 'VZ        = ', STATE(6)
 
         TARGET = 'SUN'
         FRAME  = 'CASSINI_INMS'
         CORRTN = 'LT+S'
         OBSERV = 'CASSINI'
         CALL SPKPOS( TARGET, ET, FRAME, CORRTN, OBSERV, SUNDIR, LT )
         CALL VHATIP( SUNDIR )
 
         WRITE(*,'(A,F20.6)') 'SUNDIR(X) = ', SUNDIR(1)
         WRITE(*,'(A,F20.6)') 'SUNDIR(Y) = ', SUNDIR(2)
         WRITE(*,'(A,F20.6)') 'SUNDIR(Z) = ', SUNDIR(3)
 
         METHOD = 'NEAR POINT: ELLIPSOID'
         TARGET = 'PHOEBE'
         FRAME  = 'IAU_PHOEBE'
         CORRTN = 'NONE'
         OBSERV = 'CASSINI'
         CALL SUBPNT ( METHOD, TARGET, ET, FRAME, CORRTN, OBSERV,
        .              SPOINT, TRGEPC, SRFVEC )
 
         CALL RECLAT( SPOINT, SRAD, SLON, SLAT )
 
         FROMFR = 'IAU_PHOEBE'
         TOFR   = 'CASSINI_INMS'
         CALL PXFORM( FROMFR, TOFR, ET, M2IMAT )
 
         CALL MXV   ( M2IMAT, SRFVEC, SBPDIR )
         CALL VHATIP( SBPDIR )
 
         WRITE(*,'(A,F20.6)') 'LON       = ', SLON*DPR()
         WRITE(*,'(A,F20.6)') 'LAT       = ', SLAT*DPR()
         WRITE(*,'(A,F20.6)') 'SBPDIR(X) = ', SBPDIR(1)
         WRITE(*,'(A,F20.6)') 'SBPDIR(Y) = ', SBPDIR(2)
         WRITE(*,'(A,F20.6)') 'SBPDIR(Z) = ', SBPDIR(3)
 
         TARGET = 'CASSINI'
         FRAME  = 'J2000'
         CORRTN = 'NONE'
         OBSERV = 'PHOEBE'
         CALL SPKEZR( TARGET, ET, FRAME, CORRTN, OBSERV, STATE, LT )
         CALL VPACK ( STATE(4), STATE(5), STATE(6), SCVDIR )
 
         FROMFR = 'J2000'
         TOFR   = 'CASSINI_INMS'
         CALL PXFORM( FROMFR, TOFR, ET, J2IMAT )
 
         CALL MXV   ( J2IMAT, SCVDIR, TMPDIR )
         CALL VHAT  ( TMPDIR, SCVDIR )
 
         WRITE(*,'(A,F20.6)') 'SCVDIR(X) = ', SCVDIR(1)
         WRITE(*,'(A,F20.6)') 'SCVDIR(Y) = ', SCVDIR(2)
         WRITE(*,'(A,F20.6)') 'SCVDIR(Z) = ', SCVDIR(3)
 
         END
Meta-kernel file ``scvel.tm'':

   KPL/MK
 
      The names and contents of the kernels referenced by this
      meta-kernel are as follows:
 
 
      File Name                   Description
      --------------------------  ----------------------------------
      naif0008.tls                Generic LSK.
      cas00084.tsc                Cassini SCLK.
      020514_SE_SAT105.bsp        Saturnian Satellite Ephemeris SPK.
      030201AP_SK_SM546_T45.bsp   Cassini Spacecraft SPK.
      981005_PLTEPH-DE405S.bsp    Planetary Ephemeris SPK.
      04135_04171pc_psiv2.bc      Cassini Spacecraft CK.
      cas_v37.tf                  Cassini FK.
      cpck05Mar2004.tpc           Cassini project PCK.
 
 
   \begindata
      KERNELS_TO_LOAD = (
                        'kernels/lsk/naif0008.tls'
                        'kernels/sclk/cas00084.tsc'
                        'kernels/spk/020514_SE_SAT105.bsp'
                        'kernels/spk/030201AP_SK_SM546_T45.bsp'
                        'kernels/spk/981005_PLTEPH-DE405S.bsp'
                        'kernels/ck/04135_04171pc_psiv2.bc'
                        'kernels/fk/cas_v37.tf'
                        'kernels/pck/cpck05Mar2004.tpc'
                        )
   \begintext