| Remote Sensing Hands-On Lesson, using MPO (Python) |
Table of ContentsRemote Sensing Hands-On Lesson, using MPO (Python) Overview Note About HTML Links References Tutorials Required Readings The Permuted Index SpiceyPy API Documentation Kernels Used SpiceyPy Modules Used Time Conversion (convtm) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output Extra Credit Task statements and questions Solutions and answers Obtaining Target States and Positions (getsta) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output Extra Credit Task statements and questions Solutions and answers Spacecraft Orientation and Reference Frames (xform) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output Extra Credit Task statements and questions Solutions and answers Computing Sub-s/c and Sub-solar Points on an Ellipsoid and a DSK (subpts) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output Extra Credit Task statements and questions Solutions and answers Intersecting Vectors with an Ellipsoid and a DSK (fovint) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output Extra Credit Remote Sensing Hands-On Lesson, using MPO (Python)
Overview
Note About HTML Links
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
Of these documents, the ``Tutorials'' contains the highest level descriptions with the least number of details while the ``Required Reading'' documents contain much more detailed specifications. The most complete specifications are provided in the ``API Documentation''. In some cases the lesson explanations also refer to the information provided in the meta-data area of the kernels used in the lesson examples. It is especially true in case of the FK and IK files, which often contain comprehensive descriptions of the frames, instrument FOVs, etc. Since both the FK and IK are text kernels, the information provided in them can be viewed using any text editor, while the meta information provided in binary kernels---SPKs and CKs---can be viewed using ``commnt'' or ``spacit'' utility programs located in ``cspice/exe'' of Toolkit installation tree. Tutorials
Name Lesson steps/functions it describes ---------------- ----------------------------------------------- Time Time Conversion SCLK and LSK Time Conversion SPK Obtaining Ephemeris Data Frames Reference Frames Using Frames Reference Frames PCK Planetary Constants Data CK Spacecraft Orientation Data DSK Detailed Target Shape (Topography) DataThese tutorials are available from the NAIF server at JPL:
https://naif.jpl.nasa.gov/naif/tutorials.html Required Readings
Name Lesson steps/functions that it describes --------------- ----------------------------------------- ck.req Obtaining spacecraft orientation data dsk.req Obtaining detailed body shape data frames.req Using reference frames naif_ids.req Determining body ID codes pck.req Obtaining planetary constants data sclk.req SCLK time conversion spk.req Obtaining ephemeris data time.req Time conversion The Permuted Index
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
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.htmldescribes extensively the str2et functionality. Kernels Used
1. Generic LSK:
naif0012.tls
2. BepiColombo MPO SCLK:
bc_mpo_step_20230117.tsc
3. Solar System Ephemeris SPK, subsetted to cover only the time
range of interest:
de432s.bsp
4. BepiColombo MPO Spacecraft Trajectory SPK, subsetted to cover
only the time range of interest:
bc_mpo_mlt_50037_20260314_20280529_v05.bsp
5. BepiColombo MPO FK:
bc_mpo_v32.tf
6. BepiColombo MPO Spacecraft CK, subsetted to cover only the time
range of interest:
bc_mpo_sc_slt_50028_20260314_20280529_f20181127_v03.bc
7. Generic PCK:
pck00011.tpc
8. Low-resolution Mercury DSK:
mercury_lowres.bds
9. SIMBIO-SYS IK:
bc_mpo_simbio-sys_v08.ti
These SPICE kernels are included in the lesson package.
In addition to these kernels, the extra credit exercises require the following kernels:
# FILE NAME TYPE DESCRIPTION -- --------------- ---- --------------------------------------------- 10 jup365_2027.bsp SPK Generic Jovian Satellite Ephemeris SPKThese SPICE kernels are available from the NAIF server at JPL, in the ``satellites/a_old_versions'' subdurectory:
https://naif.jpl.nasa.gov/pub/naif/generic_kernels/spk/ SpiceyPy Modules Used
CHAPTER EXERCISE FUNCTIONS NON-VOID KERNELS
------- --------- --------------- --------------- ----------
1 convtm spiceypy.furnsh spiceypy.str2et 1,2
spiceypy.unload spiceypy.etcal
spiceypy.timout
spiceypy.sce2s
extra (*) spiceypy.unitim 1,2
spiceypy.sct2e
spiceypy.et2utc
spiceypy.scs2e
2 getsta spiceypy.furnsh spiceypy.str2et 1,3,4
spiceypy.unload spiceypy.spkezr
spiceypy.spkpos
spiceypy.vnorm
spiceypy.convrt
extra (*) spiceypy.kclear 1,4,10
3 xform spiceypy.furnsh spiceypy.str2et 1-7
spiceypy.unload spiceypy.spkezr
spiceypy.sxform
spiceypy.mxvg
spiceypy.spkpos
spiceypy.pxform
spiceypy.mxv
spiceypy.convrt
spiceypy.vsep
extra (*) spiceypy.kclear 1-7
4 subpts spiceypy.furnsh spiceypy.str2et 1,3-4,7,8
spiceypy.unload spiceypy.subpnt
spiceypy.vnorm
spiceypy.subslr
extra (*) spiceypy.kclear spiceypy.reclat 1,3-4,7,10
spiceypy.dpr
spiceypy.bodvrd
spiceypy.recpgr
5 fovint spiceypy.furnsh spiceypy.str2et 1-9
spiceypy.no_found_check spiceypy.getfvn
spiceypy.unload spiceypy.bodn2c
spiceypy.sincpt
spiceypy.reclat
spiceypy.dpr
spiceypy.illumf
spiceypy.et2lst
(*) 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.
Time Conversion (convtm)Task Statement
Learning Goals
Approach
When completing the ``calendar format'' step above, consider using one of two possible methods: spiceypy.etcal or spiceypy.timout. SolutionSolution Meta-Kernel
KPL/MK
This is the meta-kernel used in the solution of the ``Time
Conversion'' task in the Remote Sensing Hands On Lesson.
The names and contents of the kernels referenced by this
meta-kernel are as follows:
1. Generic LSK:
naif0012.tls
2. BepiColombo MPO SCLK:
bc_mpo_step_20230117.tsc
\begindata
KERNELS_TO_LOAD = (
'kernels/lsk/naif0012.tls',
'kernels/sclk/bc_mpo_step_20230117.tsc'
)
\begintext
Solution Source Code
#
# Solution convtm
#
from __future__ import print_function
from builtins import input
import spiceypy
def convtm():
#
# Local Parameters
#
METAKR = 'convtm.tm'
SCLKID = -121
#
# Load the kernels this program requires.
# Both the spacecraft clock kernel and a
# leapseconds kernel should be listed in
# the meta-kernel.
#
spiceypy.furnsh( METAKR )
#
# Prompt the user for the input time string.
#
utctim = input( 'Input UTC Time: ' )
print( 'Converting UTC Time: {:s}'.format( utctim ) )
#
# Convert utctim to ET.
#
et = spiceypy.str2et( utctim )
print( ' ET Seconds Past J2000: {:16.3f}'.format( et ) )
#
# Now convert ET to a calendar time string.
# This can be accomplished in two ways.
#
calet = spiceypy.etcal( et )
print( ' Calendar ET (etcal): {:s}'.format( calet ) )
#
# Or use timout for finer control over the
# output format. The picture below was built
# by examining the header of timout.
#
calet = spiceypy.timout( et, 'YYYY-MON-DDTHR:MN:SC ::TDB' )
print( ' Calendar ET (timout): {:s}'.format( calet ) )
#
# Convert ET to spacecraft clock time.
#
sclkst = spiceypy.sce2s( SCLKID, et )
print( ' Spacecraft Clock Time: {:s}'.format( sclkst ) )
spiceypy.unload( METAKR )
if __name__ == '__main__':
convtm()
Solution Sample Output
Input UTC Time: 2027 JAN 05 02:04:36
Converting UTC Time: 2027 JAN 05 02:04:36
ET Seconds Past J2000: 852386745.184
Calendar ET (etcal): 2027 JAN 05 02:05:45.184
Calendar ET (timout): 2027-JAN-05T02:05:45
Spacecraft Clock Time: 1/0863834674:28127
Extra Credit
These ``extra credit'' tasks are provided as task statements, and unlike the regular tasks, no approach or solution source code is provided. In the next section, you will find the numeric solutions (when applicable) and answers to the questions asked in these tasks. Task statements and questions
Solutions and answers
Julian Date TDB: 2461410.5873285
Traceback (most recent call last):
File "convtm.py", line 73, in <module>
convtm()
File "convtm.py", line 36, in convtm
et = spiceypy.str2et( utctim )
File "/home/bsemenov/virtenvs/ndt_3.8/lib/python3.8/site-packages/s
piceypy/spiceypy.py", line 139, in with_errcheck
check_for_spice_error(f)
File "/home/bsemenov/virtenvs/ndt_3.8/lib/python3.8/site-packages/s
piceypy/spiceypy.py", line 122, in check_for_spice_error
raise dynamically_instantiate_spiceyerror(
spiceypy.utils.exceptions.SpiceNOLEAPSECONDS:
=====================================================================
===========
Toolkit version: CSPICE_N0067
SPICE(NOLEAPSECONDS) --
The variable that points to the leapseconds (DELTET/DELTA_AT) could n
ot be located in the kernel pool. It is likely that the leapseconds
kernel has not been loaded.
str2et_c --> STR2ET --> TTRANS
=====================================================================
===========
Traceback (most recent call last):
File "convtm.py", line 73, in <module>
convtm()
File "convtm.py", line 64, in convtm
sclkst = spiceypy.sce2s( SCLKID, et )
File "/home/bsemenov/virtenvs/ndt_3.8/lib/python3.8/site-packages/s
piceypy/spiceypy.py", line 139, in with_errcheck
check_for_spice_error(f)
File "/home/bsemenov/virtenvs/ndt_3.8/lib/python3.8/site-packages/s
piceypy/spiceypy.py", line 122, in check_for_spice_error
raise dynamically_instantiate_spiceyerror(
spiceypy.utils.exceptions.SpiceKERNELVARNOTFOUND:
=====================================================================
===========
Toolkit version: CSPICE_N0067
SPICE(KERNELVARNOTFOUND) --
The Variable Was not Found in the Kernel Pool.
Kernel variable SCLK_DATA_TYPE_121 was not found in the kernel pool.
sce2s_c --> SCE2S --> SCE2T --> SCTYPE --> SCTY01
=====================================================================
===========
Earliest UTC convertible to SCLK: 1999-08-22T00:00:05.204
Spacecraft Clock Time: 1/0863834674:28127
UTC time from spacecraft clock: 2027-01-05T02:04:36.000
Obtaining Target States and Positions (getsta)Task Statement
Learning Goals
Approach
When deciding which SPK files to load, the Toolkit utility ``brief'' may be of some use. ``brief'' is located in the ``cspice/exe'' directory for C toolkits. Consult its user's guide available in ``cspice/doc/brief.ug'' for details. SolutionSolution Meta-Kernel
KPL/MK
This is the meta-kernel used in the solution of the
``Obtaining Target States and Positions'' task in the
Remote Sensing Hands On Lesson.
The names and contents of the kernels referenced by this
meta-kernel are as follows:
1. Generic LSK:
naif0012.tls
2. Solar System Ephemeris SPK, subsetted to cover only
the time range of interest:
de432s.bsp
3. BepiColombo MPO Spacecraft Trajectory SPK, subsetted
to cover only the time range of interest:
bc_mpo_mlt_50037_20260314_20280529_v05.bsp
\begindata
KERNELS_TO_LOAD = (
'kernels/lsk/naif0012.tls',
'kernels/spk/de432s.bsp',
'kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp',
)
\begintext
Solution Source Code
#
# Solution getsta.py
#
from __future__ import print_function
from builtins import input
import spiceypy
def getsta():
#
# Local parameters
#
METAKR = 'getsta.tm'
#
# Load the kernels that this program requires. We
# will need a leapseconds kernel to convert input
# UTC time strings into ET. We also will need the
# necessary SPK files with coverage for the bodies
# in which we are interested.
#
spiceypy.furnsh( METAKR )
#
#Prompt the user for the input time string.
#
utctim = input( 'Input UTC Time: ' )
print( 'Converting UTC Time: {:s}'.format(utctim) )
#
#Convert utctim to ET.
#
et = spiceypy.str2et( utctim )
print( ' ET seconds past J2000: {:16.3f}'.format(et) )
#
# Compute the apparent state of Mercury as seen from
# BepiColombo MPO in the J2000 frame. All of the ephemeris
# readers return states in units of kilometers and
# kilometers per second.
#
[state, ltime] = spiceypy.spkezr( 'MERCURY', et, 'J2000',
'LT+S', 'MPO' )
print( ' Apparent state of Mercury as seen '
'from BepiColombo MPO in the\n'
' J2000 frame (km, km/s):' )
print( ' X = {:16.3f}'.format(state[0]) )
print( ' Y = {:16.3f}'.format(state[1]) )
print( ' Z = {:16.3f}'.format(state[2]) )
print( ' VX = {:16.3f}'.format(state[3]) )
print( ' VY = {:16.3f}'.format(state[4]) )
print( ' VZ = {:16.3f}'.format(state[5]) )
#
# Compute the apparent position of Earth as seen from
# BepiColombo MPO in the J2000 frame. Note: We could have
# continued using spkezr and simply ignored the
# velocity components.
#
[pos, ltime] = spiceypy.spkpos( 'EARTH', et, 'J2000',
'LT+S', 'MPO', )
print( ' Apparent position of Earth as '
'seen from BepiColombo MPO in the\n'
' J2000 frame (km):' )
print( ' X = {:16.3f}'.format(pos[0]) )
print( ' Y = {:16.3f}'.format(pos[1]) )
print( ' Z = {:16.3f}'.format(pos[2]) )
#
# We need only display LTIME, as it is precisely the
# light time in which we are interested.
#
print( ' One way light time between BepiColombo '
'MPO and the apparent\n'
' position of Earth (seconds):'
' {:16.3f}'.format(ltime) )
#
# Compute the apparent position of the Sun as seen from
# Mercury in the J2000 frame.
#
[pos, ltime] = spiceypy.spkpos( 'SUN', et, 'J2000',
'LT+S', 'MERCURY', )
print( ' Apparent position of Sun as '
'seen from Mercury in the\n'
' J2000 frame (km):' )
print( ' X = {:16.3f}'.format(pos[0]) )
print( ' Y = {:16.3f}'.format(pos[1]) )
print( ' Z = {:16.3f}'.format(pos[2]) )
#
# Now we need to compute the actual distance between
# the Sun and Mercury. The above spkpos call gives us
# the apparent distance, so we need to adjust our
# aberration correction appropriately.
#
[pos, ltime] = spiceypy.spkpos( 'SUN', et, 'J2000',
'NONE', 'MERCURY' )
#
# Compute the distance between the body centers in
# kilometers.
#
dist = spiceypy.vnorm( pos )
#
# Convert this value to AU using convrt.
#
dist = spiceypy.convrt( dist, 'KM', 'AU' )
print( ' Actual distance between Sun and '
'Mercury body centers:\n'
' (AU): {:16.3f}'.format(dist) )
spiceypy.unload( METAKR )
if __name__ == '__main__':
getsta()
Solution Sample Output
Input UTC Time: 2027 JAN 05 02:04:36
Converting UTC Time: 2027 JAN 05 02:04:36
ET seconds past J2000: 852386745.184
Apparent state of Mercury as seen from BepiColombo MPO in the
J2000 frame (km, km/s):
X = -683.207
Y = -1438.946
Z = -2427.819
VX = 0.036
VY = 2.360
VZ = -1.783
Apparent position of Earth as seen from BepiColombo MPO in the
J2000 frame (km):
X = -59257854.691
Y = 185201786.218
Z = 88178321.179
One way light time between BepiColombo MPO and the apparent
position of Earth (seconds): 712.193
Apparent position of Sun as seen from Mercury in the
J2000 frame (km):
X = -23429947.239
Y = 54297427.572
Z = 31434173.468
Actual distance between Sun and Mercury body centers:
(AU): 0.448
Extra Credit
These ``extra credit'' tasks are provided as task statements, and unlike the regular tasks, no approach or solution source code is provided. In the next section, you will find the numeric solutions (when applicable) and answers to the questions asked in these tasks. Task statements and questions
Solutions and answers
Traceback (most recent call last):
File "getsta.py", line 128, in <module>
getsta()
File "getsta.py", line 46, in getsta
[state, ltime] = spiceypy.spkezr( 'MERCURY', et, 'J2000',
File "/home/bsemenov/virtenvs/ndt_3.8/lib/python3.8/site-packages/s
piceypy/spiceypy.py", line 139, in with_errcheck
check_for_spice_error(f)
File "/home/bsemenov/virtenvs/ndt_3.8/lib/python3.8/site-packages/s
piceypy/spiceypy.py", line 122, in check_for_spice_error
raise dynamically_instantiate_spiceyerror(
spiceypy.utils.exceptions.SpiceSPKINSUFFDATA:
=====================================================================
===========
Toolkit version: CSPICE_N0067
SPICE(SPKINSUFFDATA) --
Insufficient ephemeris data has been loaded to compute the state of -
121 (BEPICOLOMBO MPO) relative to 0 (SOLAR SYSTEM BARYCENTER) at the
ephemeris epoch 2027 JAN 05 02:05:45.184.
spkezr_c --> SPKEZR --> SPKEZ --> SPKACS --> SPKGEO
=====================================================================
===========
BRIEF -- Version 4.1.0, September 17, 2021 -- Toolkit Version N0067
Summary for: kernels/spk/de432s.bsp
Bodies: MERCURY BARYCENTER (1) w.r.t. SOLAR SYSTEM BARYCENTER (0)
VENUS BARYCENTER (2) w.r.t. SOLAR SYSTEM BARYCENTER (0)
EARTH BARYCENTER (3) w.r.t. SOLAR SYSTEM BARYCENTER (0)
MARS BARYCENTER (4) w.r.t. SOLAR SYSTEM BARYCENTER (0)
JUPITER BARYCENTER (5) w.r.t. SOLAR SYSTEM BARYCENTER (0)
SATURN BARYCENTER (6) w.r.t. SOLAR SYSTEM BARYCENTER (0)
URANUS BARYCENTER (7) w.r.t. SOLAR SYSTEM BARYCENTER (0)
NEPTUNE BARYCENTER (8) w.r.t. SOLAR SYSTEM BARYCENTER (0)
PLUTO BARYCENTER (9) w.r.t. SOLAR SYSTEM BARYCENTER (0)
SUN (10) w.r.t. SOLAR SYSTEM BARYCENTER (0)
MERCURY (199) w.r.t. MERCURY BARYCENTER (1)
VENUS (299) w.r.t. VENUS BARYCENTER (2)
MOON (301) w.r.t. EARTH BARYCENTER (3)
EARTH (399) w.r.t. EARTH BARYCENTER (3)
Start of Interval (UTC) End of Interval (UTC)
----------------------------- -------------------------
----
2027-JAN-02 23:01:53.350 2027-JAN-08 00:59:37.932
Summary for: kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp
Body: BEPICOLOMBO MPO (-121) w.r.t. MERCURY (199)
Start of Interval (UTC) End of Interval (UTC)
----------------------------- ---------------------------
--
2027-JAN-02 23:01:53.350 2027-JAN-08 00:59:37.932
Bodies: -121000 w.r.t. BEPICOLOMBO MPO (-121)
-121540 w.r.t. BEPICOLOMBO MPO (-121)
-121600 w.r.t. BEPICOLOMBO MPO (-121)
Start of Interval (UTC) End of Interval (UTC)
----------------------------- -------------------------
----
2027-JAN-02 23:01:53.350 2027-JAN-08 00:59:37.932
Additional kernels required for this task:
1. Generic Jovian Satellite Ephemeris SPK:
jup365_2027.bsp
available in the NAIF server at:
https://naif.jpl.nasa.gov/pub/naif/generic_kernels/spk/
satellites/a_old_versions
Actual position of Jupiter as seen from Mercury in the
J2000 frame (km):
X = -623644094.418
Y = 532767093.112
Z = 251130102.035
Actual (geometric) position of Sun as seen from Mercury in the
J2000 frame (km):
X = -23438490.402
Y = 54294213.485
Z = 31433347.025
Light-time corrected position of Sun as seen from Mercury in the
J2000 frame (km):
X = -23438492.550
Y = 54294212.272
Z = 31433346.550
Apparent position of Sun as seen from Mercury in the
J2000 frame (km):
X = -23430052.903
Y = 54297381.156
Z = 31434164.775
Spacecraft Orientation and Reference Frames (xform)Task Statement
Learning Goals
Approach
You may find it useful to consult the permuted index, the headers of various source modules, and the following toolkit documentation:
SolutionSolution Meta-Kernel
KPL/MK
This is the meta-kernel used in the solution of the ``Spacecraft
Orientation and Reference Frames'' task in the Remote Sensing
Hands On Lesson.
The names and contents of the kernels referenced by this
meta-kernel are as follows:
1. Generic LSK:
naif0012.tls
2. BepiColombo MPO SCLK:
bc_mpo_step_20230117.tsc
3. Solar System Ephemeris SPK, subsetted to cover only
the time range of interest:
de432s.bsp
4. BepiColombo MPO Spacecraft Trajectory SPK, subsetted
to cover only the time range of interest:
bc_mpo_mlt_50037_20260314_20280529_v05.bsp
5. BepiColombo MPO FK:
bc_mpo_v32.tf
6. BepiColombo MPO Spacecraft CK, subsetted to cover only
the time range of interest:
bc_mpo_sc_slt_50028_20260314_20280529_f20181127_v03.bc
7. Generic PCK:
pck00011.tpc
\begindata
KERNELS_TO_LOAD = (
'kernels/lsk/naif0012.tls',
'kernels/sclk/bc_mpo_step_20230117.tsc',
'kernels/spk/de432s.bsp',
'kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp',
'kernels/fk/bc_mpo_v32.tf',
'kernels/ck/bc_mpo_sc_slt_50028_20260314_20280529_f20181127_v03.bc',
'kernels/pck/pck00011.tpc'
)
\begintext
Solution Source Code
#
# Solution xform.py
#
from __future__ import print_function
from builtins import input
import spiceypy
def xform():
#
# Local parameters
#
METAKR = 'xform.tm'
#
# Load the kernels that this program requires. We
# will need:
#
# A leapseconds kernel
# A spacecraft clock kernel for BepiColombo MPO
# The necessary ephemerides
# A planetary constants file (PCK)
# A spacecraft orientation kernel for BepiColombo MPO (CK)
# A frame kernel (TF)
#
spiceypy.furnsh( METAKR )
#
# Prompt the user for the input time string.
#
utctim = input( 'Input UTC Time: ' )
print( 'Converting UTC Time: {:s}'.format(utctim) )
#
#Convert utctim to ET.
#
et = spiceypy.str2et( utctim )
print( ' ET seconds past J2000: {:16.3f}'.format(et) )
#
# Compute the apparent state of Mercury as seen from
# BepiColombo MPO in the J2000 frame.
#
[state, ltime] = spiceypy.spkezr( 'MERCURY', et, 'J2000',
'LT+S', 'MPO' )
#
# Now obtain the transformation from the inertial
# J2000 frame to the non-inertial body-fixed IAU_MERCURY
# frame. Since we want the apparent position, we
# need to subtract ltime from et.
#
sform = spiceypy.sxform( 'J2000', 'IAU_MERCURY', et-ltime )
#
# Now rotate the apparent J2000 state into IAU_MERCURY
# with the following matrix multiplication:
#
bfixst = spiceypy.mxvg ( sform, state )
#
# Display the results.
#
print( ' Apparent state of Mercury as seen '
'from BepiColombo MPO in the\n'
' IAU_MERCURY body-fixed frame (km, km/s):' )
print( ' X = {:19.6f}'.format(bfixst[0]) )
print( ' Y = {:19.6f}'.format(bfixst[1]) )
print( ' Z = {:19.6f}'.format(bfixst[2]) )
print( ' VX = {:19.6f}'.format(bfixst[3]) )
print( ' VY = {:19.6f}'.format(bfixst[4]) )
print( ' VZ = {:19.6f}'.format(bfixst[5]) )
#
# It is worth pointing out, all of the above could
# have been done with a single use of spkezr:
#
[state, ltime] = spiceypy.spkezr(
'MERCURY', et, 'IAU_MERCURY',
'LT+S', 'MPO' )
#
# Display the results.
#
print( ' Apparent state of Mercury as seen '
'from BepiColombo MPO in the\n'
' IAU_MERCURY body-fixed frame '
'(km, km/s) obtained using\n'
' spkezr directly:' )
print( ' X = {:19.6f}'.format(state[0]) )
print( ' Y = {:19.6f}'.format(state[1]) )
print( ' Z = {:19.6f}'.format(state[2]) )
print( ' VX = {:19.6f}'.format(state[3]) )
print( ' VY = {:19.6f}'.format(state[4]) )
print( ' VZ = {:19.6f}'.format(state[5]) )
#
# Note that the velocity found by using spkezr
# to compute the state in the IAU_MERCURY frame differs
# at the few mm/second level from that found previously
# by calling spkezr and then sxform. Computing
# velocity via a single call to spkezr as we've
# done immediately above is slightly more accurate because
# it accounts for the effect of the rate of change of
# light time on the apparent angular velocity of the
# target's body-fixed reference frame.
#
# Now we are to compute the angular separation between
# the apparent position of Mercury as seen from the orbiter
# and the nominal instrument view direction. First,
# compute the apparent position of Mercury as seen from
# BepiColombo MPO in the J2000 frame.
#
[pos, ltime] = spiceypy.spkpos( 'MERCURY', et, 'J2000',
'LT+S', 'MPO' )
#
# Now compute the location of the nominal instrument view
# direction. From reading the frame kernel we know that
# the instrument view direction is nominally the +Z axis
# of the MPO_SPACECRAFT frame defined there.
#
bsight = [ 0.0, 0.0, 1.0 ]
#
# Now compute the rotation matrix from MPO_SPACECRAFT into
# J2000.
#
pform = spiceypy.pxform( 'MPO_SPACECRAFT', 'J2000', et )
#
# And multiply the result to obtain the nominal instrument
# view direction in the J2000 reference frame.
#
bsight = spiceypy.mxv( pform, bsight )
#
# Lastly compute the angular separation.
#
sep = spiceypy.convrt( spiceypy.vsep(bsight, pos),
'RADIANS', 'DEGREES' )
print( ' Angular separation between the '
'apparent position of Mercury and\n'
' the BepiColombo MPO nominal instrument '
'view direction\n'
' (degrees):\n'
' {:16.3f}'.format(sep) )
#
# Or alternatively we can work in the spacecraft
# frame directly.
#
[pos, ltime] = spiceypy.spkpos(
'MERCURY', et, 'MPO_SPACECRAFT',
'LT+S', 'MPO' )
#
# The nominal instrument view direction is the +Z-axis
# in the MPO_SPACECRAFT frame.
#
bsight = [ 0.0, 0.0, 1.0 ]
#
# Lastly compute the angular separation.
#
sep = spiceypy.convrt( spiceypy.vsep(bsight, pos),
'RADIANS', 'DEGREES' )
print( ' Angular separation between the '
'apparent position of Mercury and\n'
' the BepiColombo MPO nominal instrument '
'view direction computed\n'
' using vectors in the '
'MPO_SPACECRAFT frame (degrees):\n'
' {:16.3f}'.format(sep) )
spiceypy.unload( METAKR )
if __name__ == '__main__':
xform()
Solution Sample Output
Input UTC Time: 2027 JAN 05 02:04:36
Converting UTC Time: 2027 JAN 05 02:04:36
ET seconds past J2000: 852386745.184
Apparent state of Mercury as seen from BepiColombo MPO in the
IAU_MERCURY body-fixed frame (km, km/s):
X = -2354.697620
Y = -762.547549
Z = -1518.408470
VX = 1.208589
VY = 0.394259
VZ = -2.671125
Apparent state of Mercury as seen from BepiColombo MPO in the
IAU_MERCURY body-fixed frame (km, km/s) obtained using
spkezr directly:
X = -2354.697620
Y = -762.547549
Z = -1518.408470
VX = 1.208589
VY = 0.394259
VZ = -2.671125
Angular separation between the apparent position of Mercury and
the BepiColombo MPO nominal instrument view direction
(degrees):
0.009
Angular separation between the apparent position of Mercury and
the BepiColombo MPO nominal instrument view direction computed
using vectors in the MPO_SPACECRAFT frame (degrees):
0.009
Extra Credit
These ``extra credit'' tasks are provided as task statements, and unlike the regular tasks, no approach or solution source code is provided. In the next section, you will find the numeric solutions (when applicable) and answers to the questions asked in these tasks. Task statements and questions
Solutions and answers
Traceback (most recent call last):
File "xform.py", line 185, in <module>
xform()
File "xform.py", line 131, in xform
pform = spiceypy.pxform( 'MPO_SPACECRAFT', 'J2000', et )
File "/home/bsemenov/virtenvs/ndt_3.8/lib/python3.8/site-packages/s
piceypy/spiceypy.py", line 139, in with_errcheck
check_for_spice_error(f)
File "/home/bsemenov/virtenvs/ndt_3.8/lib/python3.8/site-packages/s
piceypy/spiceypy.py", line 122, in check_for_spice_error
raise dynamically_instantiate_spiceyerror(
spiceypy.utils.exceptions.SpiceNOFRAMECONNECT:
=====================================================================
===========
Toolkit version: CSPICE_N0067
SPICE(NOFRAMECONNECT) --
At epoch 8.5252159418408E+08 TDB (2027 JAN 06 15:33:14.184 TDB), ther
e is insufficient information available to transform from reference f
rame -121000 (MPO_SPACECRAFT) to reference frame 1 (J2000). MPO_SPACE
CRAFT is a CK frame; a CK file containing data for instrument or stru
cture -121000 at the epoch shown above, as well as a corresponding SC
LK kernel, must be loaded in order to use this frame. Failure to find
required CK data could be due to one or more CK files not having bee
n loaded, or to the epoch shown above lying within a coverage gap or
beyond the coverage bounds of the loaded CK files. It is also possibl
e that no loaded CK file has required angular velocity data for the i
nput epoch, even if a loaded CK does have attitude data for that epoc
h. You can use CKBRIEF with the -dump option to display coverage inte
rvals of a CK file.
pxform_c --> PXFORM --> REFCHG
=====================================================================
===========
CKBRIEF -- Version 6.1.0, June 27, 2014 -- Toolkit Version N0067
Summary for: kernels/ck/bc_mpo_sc_slt_50028_20260314_20280529_f201811
27_v03.bc
Segment No.: 1
Object: -121000
Interval Begin UTC Interval End UTC AV
------------------------ ------------------------ ---
2027-JAN-02 23:01:53.350 2027-JAN-06 11:04:56.368 Y
2027-JAN-06 11:08:00.779 2027-JAN-06 15:30:56.685 Y
2027-JAN-06 15:33:04.016 2027-JAN-06 22:05:57.865 Y
2027-JAN-06 22:10:03.746 2027-JAN-08 00:59:37.932 Y
CKBRIEF -- Version 6.1.0, June 27, 2014 -- Toolkit Version N0067
Summary for: kernels/ck/bc_mpo_sc_slt_50028_20260314_20280529_f201811
27_v03.bc
Object: -121000
Interval Begin UTC Interval End UTC AV
------------------------ ------------------------ ---
2027-JAN-02 23:01:53.350 2027-JAN-08 00:59:37.932 Y
Angular separation between the apparent position of the Sun and the
BepiColombo MPO nominal instrument view direction (degrees):
135.393
Science Deck illumination:
BepiColombo MPO Science Deck IS NOT illuminated.
Computing Sub-s/c and Sub-solar Points on an Ellipsoid and a DSK (subpts)Task Statement
near point/ellipsoid
definition, and once using a DSK shape model and the
nadir/dsk/unprioritized
definition.
The program displays the results. Use the program to compute these quantities at "2027 JAN 05 02:04:36" UTC. Learning Goals
Approach
One point worth considering: how would the results change if the sub-solar and sub-observer points were computed using the
intercept/ellipsoid
and
intercept/dsk/unprioritized
definitions? Which definition is appropriate?
SolutionSolution Meta-Kernel
KPL/MK
This is the meta-kernel used in the solution of the
``Computing Sub-s/c and Sub-solar Points on an Ellipsoid
and a DSK'' task in the Remote Sensing Hands On Lesson.
The names and contents of the kernels referenced by this
meta-kernel are as follows:
1. Generic LSK:
naif0012.tls
2. Solar System Ephemeris SPK, subsetted to cover only
the time range of interest:
de432s.bsp
3. BepiColombo MPO Spacecraft Trajectory SPK, subsetted
to cover only the time range of interest:
bc_mpo_mlt_50037_20260314_20280529_v05.bsp
4. Generic PCK:
pck00011.tpc
5. Low-resolution Mercury DSK:
mercury_lowres.bds
\begindata
KERNELS_TO_LOAD = (
'kernels/lsk/naif0012.tls',
'kernels/spk/de432s.bsp',
'kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp',
'kernels/pck/pck00011.tpc'
'kernels/dsk/mercury_lowres.bds'
)
\begintext
Solution Source Code
#
# Solution subpts.py
#
from __future__ import print_function
from builtins import input
#
# SpiceyPy package:
#
import spiceypy
def subpts():
#
# Local parameters
#
METAKR = 'subpts.tm'
#
# Load the kernels that this program requires. We
# will need:
#
# A leapseconds kernel
# The necessary ephemerides
# A planetary constants file (PCK)
# A DSK file containing Mercury shape data
#
spiceypy.furnsh( METAKR )
#
#Prompt the user for the input time string.
#
utctim = input( 'Input UTC Time: ' )
print( ' Converting UTC Time: {:s}'.format(utctim) )
#
#Convert utctim to ET.
#
et = spiceypy.str2et( utctim )
print( ' ET seconds past J2000: {:16.3f}'.format(et) )
for i in range(2):
if i == 0:
#
# Use the "near point" sub-point definition
# and an ellipsoidal model.
#
method = 'NEAR POINT/Ellipsoid'
else:
#
# Use the "nadir" sub-point definition
# and a DSK model.
#
method = 'NADIR/DSK/Unprioritized'
print( '\n Sub-point/target shape model: {:s}\n'.format(
method ) )
#
# Compute the apparent sub-observer point of BepiColombo MPO
# on Mercury.
#
[spoint, trgepc, srfvec] = spiceypy.subpnt(
method, 'MERCURY', et,
'IAU_MERCURY', 'LT+S', 'MPO' )
print( ' Apparent sub-observer point of BepiColombo '
'MPO on Mercury\n'
' in the IAU_MERCURY frame (km):' )
print( ' X = {:16.3f}'.format(spoint[0]) )
print( ' Y = {:16.3f}'.format(spoint[1]) )
print( ' Z = {:16.3f}'.format(spoint[2]) )
print( ' ALT = {:16.3f}'.format(spiceypy.vnorm(srfvec)) )
#
# Compute the apparent sub-solar point on Mercury
# as seen from BepiColombo MPO.
#
[spoint, trgepc, srfvec] = spiceypy.subslr(
method, 'MERCURY', et,
'IAU_MERCURY', 'LT+S', 'MPO' )
print( ' Apparent sub-solar point on Mercury '
'as seen from BepiColombo \n'
' MPO in the IAU_MERCURY frame (km):' )
print( ' X = {:16.3f}'.format(spoint[0]) )
print( ' Y = {:16.3f}'.format(spoint[1]) )
print( ' Z = {:16.3f}'.format(spoint[2]) )
#
# End of computation block for "method"
#
print( '' )
spiceypy.unload( METAKR )
if __name__ == '__main__':
subpts()
Solution Sample Output
Input UTC Time: 2027 JAN 05 02:04:36
Converting UTC Time: 2027 JAN 05 02:04:36
ET seconds past J2000: 852386745.184
Sub-point/target shape model: NEAR POINT/Ellipsoid
Apparent sub-observer point of BepiColombo MPO on Mercury
in the IAU_MERCURY frame (km):
X = 1978.726
Y = 640.793
Z = 1275.611
ALT = 463.634
Apparent sub-solar point on Mercury as seen from BepiColombo
MPO in the IAU_MERCURY frame (km):
X = 1526.831
Y = 1903.936
Z = -1.436
Sub-point/target shape model: NADIR/DSK/Unprioritized
Apparent sub-observer point of BepiColombo MPO on Mercury
in the IAU_MERCURY frame (km):
X = 1979.558
Y = 641.062
Z = 1276.148
ALT = 462.608
Apparent sub-solar point on Mercury as seen from BepiColombo
MPO in the IAU_MERCURY frame (km):
X = 1525.673
Y = 1902.492
Z = -1.434
Extra Credit
These ``extra credit'' tasks are provided as task statements, and unlike the regular tasks, no approach or solution source code is provided. In the next section, you will find the numeric solutions (when applicable) and answers to the questions asked in these tasks. Task statements and questions
Solutions and answers
Apparent sub-solar point on Mercury as seen from BepiColombo MPO
in the IAU_MERCURY frame using the 'Near Point: ellipsoid' method
(km):
X = 1526.828
Y = 1903.939
Z = -1.435
Apparent sub-solar point on Mercury as seen from BepiColombo MPO
in the IAU_MERCURY frame using the 'Intercept: ellipsoid' method
(km):
X = 1526.828
Y = 1903.939
Z = -1.438
Additional kernels required for this task:
1. Generic Jovian Satellite Ephemeris SPK:
jup365_2027.bsp
available in the NAIF server at:
https://naif.jpl.nasa.gov/pub/naif/generic_kernels/spk/
satellites/a_old_versions
Geometric sub-spacecraft point of BepiColombo MPO on Europa in
the IAU_EUROPA frame using the 'Near Point: ellipsoid' method
(km):
X = -753.484
Y = -1366.703
Z = -24.296
Planetocentric coordinates of the BepiColombo MPO
sub-spacecraft point on Europa (degrees, km):
LAT = -0.892
LON = -118.869
R = 1560.835
Planetographic coordinates of the BepiColombo MPO
sub-spacecraft point on Europa (degrees, km):
LAT = -0.895
LON = 118.869
ALT = -1.764
Intersecting Vectors with an Ellipsoid and a DSK (fovint)Task Statement
At each point of intersection compute the following:
Additionally compute the local solar time at the intercept of the spectrometer aperture boresight with the surface of Mercury, using both ellipsoidal and DSK shape models. Use this program to compute values at the UTC epoch:
Learning Goals
Approach
SolutionSolution Meta-Kernel
KPL/MK
This is the meta-kernel used in the solution of the
``Intersecting Vectors with an Ellipsoid and a DSK'' task
in the Remote Sensing Hands On Lesson.
The names and contents of the kernels referenced by this
meta-kernel are as follows:
1. Generic LSK:
naif0012.tls
2. BepiColombo MPO SCLK:
bc_mpo_step_20230117.tsc
3. Solar System Ephemeris SPK, subsetted to cover only
the time range of interest:
de432s.bsp
4. BepiColombo MPO Spacecraft Trajectory SPK, subsetted
to cover only the time range of interest:
bc_mpo_mlt_50037_20260314_20280529_v05.bsp
5. BepiColombo MPO FK:
bc_mpo_v32.tf
6. BepiColombo MPO Spacecraft CK, subsetted to cover only
the time range of interest:
bc_mpo_sc_slt_50028_20260314_20280529_f20181127_v03.bc
7. Generic PCK:
pck00011.tpc
8. SIMBIO-SYS IK:
bc_mpo_simbio-sys_v08.ti
9. Low-resolution Mercury DSK:
mercury_lowres.bds
\begindata
KERNELS_TO_LOAD = (
'kernels/lsk/naif0012.tls',
'kernels/sclk/bc_mpo_step_20230117.tsc',
'kernels/spk/de432s.bsp',
'kernels/spk/bc_mpo_mlt_50037_20260314_20280529_v05.bsp',
'kernels/fk/bc_mpo_v32.tf',
'kernels/ck/bc_mpo_sc_slt_50028_20260314_20280529_f20181127_v03.bc',
'kernels/pck/pck00011.tpc',
'kernels/ik/bc_mpo_simbio-sys_v08.ti'
'kernels/dsk/mercury_lowres.bds'
)
\begintext
Solution Source Code
#
# Solution fovint.py
#
from __future__ import print_function
from builtins import input
#
# SpiceyPy package:
#
import spiceypy
from spiceypy.utils.support_types import SpiceyError
def fovint():
#
# Local parameters
#
METAKR = 'fovint.tm'
ROOM = 4
#
# Load the kernels that this program requires. We
# will need:
#
# A leapseconds kernel.
# A SCLK kernel for BepiColombo MPO.
# Any necessary ephemerides.
# The BepiColombo MPO frame kernel.
# An BepiColombo MPO C-kernel.
# A PCK file with Mercury constants.
# The BepiColombo MPO SIMBIO-SYS I-kernel.
# A DSK file containing Mercury shape data.
#
spiceypy.furnsh( METAKR )
#
# Prompt the user for the input time string.
#
utctim = input( 'Input UTC Time: ' )
print( 'Converting UTC Time: {:s}'.format(utctim) )
#
# Convert utctim to ET.
#
et = spiceypy.str2et( utctim )
print( ' ET seconds past J2000: {:16.3f}\n'.format(et) )
#
# Now we need to obtain the FOV configuration of
# the SIMBIO-SYS HRIC channel.
#
[ shape, insfrm, bsight,
n, bounds ] = spiceypy.getfvn( 'MPO_SIMBIO-SYS_HRIC_FPA', ROOM
)
#
# `bounds' is a numpy array. We'll convert it to a list.
#
# Rather than treat BSIGHT as a separate vector,
# copy it into the last slot of BOUNDS.
#
bounds = bounds.tolist()
bounds.append( bsight )
#
# Set vector names to be used for output.
#
vecnam = [ 'Boundary Corner 1',
'Boundary Corner 2',
'Boundary Corner 3',
'Boundary Corner 4',
'MPO SIMBIO-SYS HRIC Boresight' ]
#
# Set values of "method" string that specify use of
# ellipsoidal and DSK (topographic) shape models.
#
# In this case, we can use the same methods for calls to both
# spiceypy.sincpt and spiceypy.ilumin. Note that some SPICE
# routines require different "method" inputs from those
# shown here. See the API documentation of each routine
# for details.
#
method = [ 'Ellipsoid', 'DSK/Unprioritized' ]
#
# Get ID code of Mercury. We'll use this ID code later, when we
# compute local solar time.
#
try:
marsid = spiceypy.bodn2c( 'MERCURY' )
except SpiceyError as exc:
#
# The ID code for MERCURY is built-in to the library.
# However, it is good programming practice to get
# in the habit of handling exceptions that may
# be thrown when a quantity is not found.
#
#
# Also see the call to spiceypy.sincpt below for an
# example of handling a "not found" condition that is
# not an error.
#
print( 'Exception message is: {:s}'.format( exc.message ) )
raise
#
# Now perform the same set of calculations for each
# vector listed in the BOUNDS array. Use both
# ellipsoidal and detailed (DSK) shape models.
#
for i in range(5):
#
# Call sincpt to determine coordinates of the
# intersection of this vector with the surface
# of Mercury.
#
print( 'Vector: {:s}\n'.format( vecnam[i] ) )
for j in range(2):
print ( ' Target shape model: {:s}\n'.format(
method[j] ) )
#
# Treat an "intercept not found" condition as a normal
# case rather than as an exception. This way we
# distinguish between this condition and other
# exceptions, in particular, any other "not found"
# exceptions.
#
with spiceypy.no_found_check():
[point, trgepc, srfvec, found ] = spiceypy.sincpt(
method[j], 'MERCURY', et,
'IAU_MERCURY', 'LT+S', 'MPO',
insfrm, bounds[i] )
if found:
#
# Now, we have discovered a point of intersection.
# Start by displaying the position vector in the
# IAU_MERCURY frame of the intersection.
#
print( ' Position vector of surface intercept '
'in the IAU_MERCURY frame (km):' )
print( ' X = {:16.3f}'.format( point[0] ) )
print( ' Y = {:16.3f}'.format( point[1] ) )
print( ' Z = {:16.3f}'.format( point[2] ) )
#
# Display the planetocentric latitude and longitude
# of the intercept.
#
[radius, lon, lat] = spiceypy.reclat( point )
print( ' Planetocentric coordinates of '
'the intercept (degrees):' )
print( ' LAT = {:16.3f}'.format(
lat * spiceypy.dpr() ) )
print( ' LON = {:16.3f}'.format(
lon * spiceypy.dpr() ) )
#
# Compute the illumination angles at this
# point.
#
[ trgepc, srfvec, phase, solar,
emissn, visibl, lit ] = spiceypy.illumf(
method[j], 'MERCURY', 'SUN', et,
'IAU_MERCURY', 'LT+S', 'MPO', point )
print( ' Phase angle (degrees): '
'{:16.3f}'.format( phase*spiceypy.dpr() ) )
print( ' Solar incidence angle (degrees): '
'{:16.3f}'.format( solar*spiceypy.dpr() ) )
print( ' Emission angle (degrees): '
'{:16.3f}'.format( emissn*spiceypy.dpr()) )
print( ' Observer visible: {:s}'.format(
str(visibl) ) )
print( ' Sun visible: {:s}'.format(
str(lit) ) )
if i == 4:
#
# Compute local solar time corresponding
# to the light time corrected TDB epoch
# at the boresight intercept.
#
[hr, mn, sc, time, ampm] = spiceypy.et2lst(
trgepc,
marsid,
lon,
'PLANETOCENTRIC' )
print( '\n Local Solar Time at boresight '
'intercept (24 Hour Clock):\n'
' {:s}'.format( time ) )
#
# End of LST computation block.
#
else:
#
# Display a message if no intercept was found.
# Otherwise, continue with the calculations.
#
print( 'No intersection point found at '
'this epoch for this vector.\n' )
#
# End of spiceypy.sincpt intercept handling block.
#
print( '' )
#
# End of target shape model loop.
#
#
# End of vector loop.
#
spiceypy.unload( METAKR )
if __name__ == '__main__':
fovint()
Solution Sample Output
Input UTC Time: 2027 JAN 05 02:04:36
Converting UTC Time: 2027 JAN 05 02:04:36
ET seconds past J2000: 852386745.184
Vector: Boundary Corner 1
Target shape model: Ellipsoid
Position vector of surface intercept in the IAU_MERCURY
frame (km):
X = 1973.717
Y = 645.436
Z = 1281.009
Planetocentric coordinates of the intercept (degrees):
LAT = 31.670
LON = 18.109
Phase angle (degrees): 44.735
Solar incidence angle (degrees): 44.622
Emission angle (degrees): 1.280
Observer visible: true
Sun visible: true
Target shape model: DSK/Unprioritized
Position vector of surface intercept in the IAU_MERCURY
frame (km):
X = 1974.257
Y = 645.602
Z = 1281.346
Planetocentric coordinates of the intercept (degrees):
LAT = 31.670
LON = 18.108
Phase angle (degrees): 44.735
Solar incidence angle (degrees): 46.703
Emission angle (degrees): 4.145
Observer visible: true
Sun visible: true
Vector: Boundary Corner 2
Target shape model: Ellipsoid
Position vector of surface intercept in the IAU_MERCURY
frame (km):
X = 1979.643
Y = 647.354
Z = 1270.875
Planetocentric coordinates of the intercept (degrees):
LAT = 31.391
LON = 18.108
Phase angle (degrees): 45.641
Solar incidence angle (degrees): 44.447
Emission angle (degrees): 1.198
Observer visible: true
Sun visible: true
Target shape model: DSK/Unprioritized
Position vector of surface intercept in the IAU_MERCURY
frame (km):
X = 1980.449
Y = 647.601
Z = 1271.407
Planetocentric coordinates of the intercept (degrees):
LAT = 31.391
LON = 18.108
Phase angle (degrees): 45.641
Solar incidence angle (degrees): 43.796
Emission angle (degrees): 1.894
Observer visible: true
Sun visible: true
Vector: Boundary Corner 3
Target shape model: Ellipsoid
Position vector of surface intercept in the IAU_MERCURY
frame (km):
X = 1983.307
Y = 636.037
Z = 1270.876
Planetocentric coordinates of the intercept (degrees):
LAT = 31.391
LON = 17.781
Phase angle (degrees): 44.501
Solar incidence angle (degrees): 44.666
Emission angle (degrees): 1.195
Observer visible: true
Sun visible: true
Target shape model: DSK/Unprioritized
Position vector of surface intercept in the IAU_MERCURY
frame (km):
X = 1984.034
Y = 636.285
Z = 1271.361
Planetocentric coordinates of the intercept (degrees):
LAT = 31.391
LON = 17.781
Phase angle (degrees): 44.501
Solar incidence angle (degrees): 45.429
Emission angle (degrees): 2.027
Observer visible: true
Sun visible: true
Vector: Boundary Corner 4
Target shape model: Ellipsoid
Position vector of surface intercept in the IAU_MERCURY
frame (km):
X = 1977.381
Y = 634.119
Z = 1281.010
Planetocentric coordinates of the intercept (degrees):
LAT = 31.670
LON = 17.780
Phase angle (degrees): 43.576
Solar incidence angle (degrees): 44.840
Emission angle (degrees): 1.278
Observer visible: true
Sun visible: true
Target shape model: DSK/Unprioritized
Position vector of surface intercept in the IAU_MERCURY
frame (km):
X = 1978.158
Y = 634.384
Z = 1281.499
Planetocentric coordinates of the intercept (degrees):
LAT = 31.670
LON = 17.781
Phase angle (degrees): 43.576
Solar incidence angle (degrees): 45.349
Emission angle (degrees): 1.920
Observer visible: true
Sun visible: true
Vector: MPO SIMBIO-SYS HRIC Boresight
Target shape model: Ellipsoid
Position vector of surface intercept in the IAU_MERCURY
frame (km):
X = 1978.524
Y = 640.740
Z = 1275.950
Planetocentric coordinates of the intercept (degrees):
LAT = 31.530
LON = 17.944
Phase angle (degrees): 44.609
Solar incidence angle (degrees): 44.644
Emission angle (degrees): 0.059
Observer visible: true
Sun visible: true
Local Solar Time at boresight intercept (24 Hour Clock):
09:46:41
Target shape model: DSK/Unprioritized
Position vector of surface intercept in the IAU_MERCURY
frame (km):
X = 1979.357
Y = 641.010
Z = 1276.487
Planetocentric coordinates of the intercept (degrees):
LAT = 31.530
LON = 17.944
Phase angle (degrees): 44.609
Solar incidence angle (degrees): 45.349
Emission angle (degrees): 1.138
Observer visible: true
Sun visible: true
Local Solar Time at boresight intercept (24 Hour Clock):
09:46:41
Extra Credit
|