PHSRM Remote Sensing Hands-On Lesson (FORTRAN) |
Table of ContentsPHSRM Remote Sensing Hands-On Lesson (FORTRAN) Overview Note About HTML Links References Tutorials Required Readings The Permuted Index Source Code Header Comments Kernels Used SPICE Modules Used Time Conversion (convtm) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output Obtaining Target States and Positions (getsta) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output Spacecraft Orientation and Reference Frames (xform) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output Computing Sub-spacecraft and Sub-solar Points (subpts) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output Intersecting Vectors with a Triaxial Ellipsoid (fovint) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output PHSRM Remote Sensing Hands-On Lesson (FORTRAN)
Overview
Note About HTML Links
In order for the links to be resolved, 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. ReferencesTutorials
Name Lesson steps/routines 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 DataThese tutorials are available from the NAIF ftp server at JPL:
http://naif.jpl.nasa.gov/naif/tutorials.html Required Readings
Name Lesson steps/routines that it describes --------------- ----------------------------------------- time.req Time Conversion sclk.req SCLK Time Conversion spk.req Obtaining Ephemeris Data frames.req Using Reference Frames pck.req Obtaining Planetary Constants Data ck.req Obtaining Spacecraft Orientation Data naif_ids.req Determining Body ID Codes The Permuted Index
This text document provides a simple mechanism to discover what SPICE routines perform a particular function of interest as well as the name of the source module that contains the routine. This is particularly useful for FORTRAN programmers because some of the routines are entry points and, therefore, the name does not translate directly into the name of the source module that contains them. Source Code Header Comments
For example the source code of the STR2ET/str2et_c routine is
toolkit/src/spicelib/str2et.forin the FORTRAN Toolkit and in
cspice/src/cspice/str2et_c.cin the C Toolkit. Since some of the FORTRAN routines are entry points they are usually 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. Kernels Used
# FILE NAME TYPE DESCRIPTION -- ------------------------------- ---- ------------------------ 1 naif0009.tls LSK Generic LSK 2 phsrm_201008031645.tsc SCLK PHSRM SCLK 3 de421xs.bsp SPK Solar System Ephemeris 4 phobos_kiam_101231_v00.bsp SPK Phobos Ephemeris 5 phsrm_130114_130114_130214_nom2.bsp SPK PHSRM Spacecraft SPK 6 phsrm_v00.tf FK PHSRM FK 7 phsrm_sc_test2_111108_130214_v00.bc CK PHSRM Spacecraft CK 8 pck00009.tpc PCK Generic PCK 9 phsrm_tsns_v00.ti IK PHSRM TSNS IKThese SPICE kernels are included in the lesson package available from the PHSRM server at IKI:
http://spice.ikiweb.ru/PHSRM/kernels SPICE Modules Used
CHAPTER EXERCISE ROUTINES FUNCTIONS KERNELS ------- --------- --------- --------- --------- 1 convtm FURNSH 1,2 PROMPT STR2ET ETCAL TIMOUT SCE2C SCE2S 2 getsta FURNSH VNORM 1,3-6 PROMPT STR2ET SPKEZR SPKPOS CONVRT 3 xform FURNSH VSEP 1-8 PROMPT STR2ET SPKEZR SXFORM MXVG SPKPOS PXFORM MXV CONVRT 4 subpts FURNSH 1,3-6,8 PROMPT STR2ET SUBPT SUBSOL 5 fovint FURNSH DPR 1-9 PROMPT STR2ET BODN2C BYEBYE GETFOV SINCPT RECLAT ILUMIN ET2LSTRefer to the headers of the various routines listed above, as detailed interface specifications are provided with the source code. Time Conversion (convtm)Task Statement
Learning Goals
Approach
When completing the ``calendar format'' step above, consider using one of two possible methods: ETCAL or 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. \begindata KERNELS_TO_LOAD = ( 'kernels/lsk/naif0009.tls', 'kernels/sclk/phsrm_201008031645.tsc' ) \begintext Solution Source Code
PROGRAM CONVTM IMPLICIT NONE C C Local Parameters C C The name of the meta-kernel that lists the kernels C to load into the program. C CHARACTER*(*) METAKR PARAMETER ( METAKR = 'convtm.tm' ) C C The spacecraft clock ID code for PHSRM. C INTEGER SCLKID PARAMETER ( SCLKID = -555 ) C C The length of various string variables. C INTEGER STRLEN PARAMETER ( STRLEN = 50 ) C C Local Variables C CHARACTER*(STRLEN) CALET CHARACTER*(STRLEN) SCLKST CHARACTER*(STRLEN) UTCTIM DOUBLE PRECISION ET C C Load the kernels this program requires. C Both the spacecraft clock kernel and a C leapseconds kernel should be listed C in the meta-kernel. C CALL FURNSH ( METAKR ) C C Prompt the user for the input time string. C CALL PROMPT ( 'Input UTC Time: ', UTCTIM ) WRITE (*,*) 'Converting UTC Time: ', UTCTIM C C Convert UTCTIM to ET. C CALL STR2ET ( UTCTIM, ET ) WRITE (*,'(A,F16.3)') ' ET Seconds Past J2000: ', ET C C Now convert ET to a formal calendar time C string. This can be accomplished in two C ways. C CALL ETCAL ( ET, CALET ) WRITE (*,*) ' Calendar ET (ETCAL): ', CALET C C Or use TIMOUT for finer control over the C output format. The picture below was built C by examining the header of TIMOUT. C CALL TIMOUT ( ET, 'YYYY-MON-DDTHR:MN:SC ::TDB', CALET ) WRITE (*,*) ' Calendar ET (TIMOUT): ', CALET C C Convert ET to spacecraft clock time. C CALL SCE2S ( SCLKID, ET, SCLKST ) WRITE (*,*) ' Spacecraft Clock Time: ', SCLKST END Solution Sample Output
Converting UTC Time: 2013-02-10 20:40:00 ET Seconds Past J2000: 413800866.185 Calendar ET (ETCAL): 2013 FEB 10 20:41:06.185 Calendar ET (TIMOUT): 2013-FEB-10T20:41:06 Spacecraft Clock Time: 1/0461:09600000 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 ``toolkit/exe'' directory for FORTRAN toolkits. Consult its user's guide available in ``toolkit/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. \begindata KERNELS_TO_LOAD = ( 'kernels/lsk/naif0009.tls', 'kernels/spk/de421xs.bsp', 'kernels/spk/phobos_kiam_101231_v00.bsp', 'kernels/spk/phsrm_130114_130114_130214_nom2.bsp' 'kernels/fk/phsrm_v00.tf' ) \begintext Solution Source Code
PROGRAM GETSTA IMPLICIT NONE C C SPICELIB Functions C DOUBLE PRECISION VNORM C C Local Parameters C C C The name of the meta-kernel that lists the kernels C to load into the program. C CHARACTER*(*) METAKR PARAMETER ( METAKR = 'getsta.tm' ) C C The length of various string variables. C INTEGER STRLEN PARAMETER ( STRLEN = 50 ) C C Local Variables C CHARACTER*(STRLEN) UTCTIM DOUBLE PRECISION DIST DOUBLE PRECISION ET DOUBLE PRECISION LTIME DOUBLE PRECISION POS ( 3 ) DOUBLE PRECISION STATE ( 6 ) C C Load the kernels that this program requires. We C will need a leapseconds kernel to convert input C UTC time strings into ET. We also will need the C necessary SPK files with coverage for the bodies C in which we are interested. C CALL FURNSH ( METAKR ) C C Prompt the user for the input time string. C CALL PROMPT ( 'Input UTC Time: ', UTCTIM ) WRITE (*,*) 'Converting UTC Time: ', UTCTIM C C Convert UTCTIM to ET. C CALL STR2ET ( UTCTIM, ET ) WRITE (*,'(A,F16.3)') ' ET seconds past J2000: ', ET C C Compute the apparent state of Phobos as seen from C PHSRM in the J2000 frame. All of the ephemeris C readers return states in units of kilometers and C kilometers per second. C CALL SPKEZR ( 'PHOBOS', ET, 'J2000', 'LT+S', . 'PHSRM', STATE, LTIME ) WRITE (*,*) ' Apparent state of Phobos as seen from ' .// 'PHSRM in the J2000' WRITE (*,*) ' frame (km, km/s):' WRITE (*,'(A,F16.3)') ' X = ', STATE(1) WRITE (*,'(A,F16.3)') ' Y = ', STATE(2) WRITE (*,'(A,F16.3)') ' Z = ', STATE(3) WRITE (*,'(A,F16.3)') ' VX = ', STATE(4) WRITE (*,'(A,F16.3)') ' VY = ', STATE(5) WRITE (*,'(A,F16.3)') ' VZ = ', STATE(6) C C Compute the apparent position of Earth as seen from C PHSRM in the J2000 frame. Note: We could have continued C using SPKEZR and simply ignored the velocity components. C CALL SPKPOS ( 'EARTH', ET, 'J2000', 'LT+S', . 'PHSRM', POS, LTIME ) WRITE (*,*) ' Apparent position of Earth as seen from ' .// 'PHSRM in the J2000' WRITE (*,*) ' frame (km):' WRITE (*,'(A,F16.3)') ' X = ', POS(1) WRITE (*,'(A,F16.3)') ' Y = ', POS(2) WRITE (*,'(A,F16.3)') ' Z = ', POS(3) C C We need only display LTIME, as it is precisely the light C time in which we are interested. C WRITE (*,*) ' One way light time between PHSRM and the ' .// 'apparent position' WRITE (*,'(A,F16.3)') ' of Earth (seconds): ', LTIME C C Compute the apparent position of the Sun as seen from C Phobos in the J2000 frame. C CALL SPKPOS ( 'SUN', ET, 'J2000', 'LT+S', . 'PHOBOS', POS, LTIME ) WRITE (*,*) ' Apparent position of Sun as seen from ' .// 'Phobos in the' WRITE (*,*) ' J2000 frame (km):' WRITE (*,'(A,F16.3)') ' X = ', POS(1) WRITE (*,'(A,F16.3)') ' Y = ', POS(2) WRITE (*,'(A,F16.3)') ' Z = ', POS(3) C C Now we need to compute the actual distance between the Sun C and Phobos. The above SPKPOS call gives us the apparent C distance, so we need to adjust our aberration correction C appropriately. C CALL SPKPOS ( 'SUN', ET, 'J2000', 'NONE', . 'PHOBOS', POS, LTIME ) C C Compute the distance between the body centers in C kilometers. C DIST = VNORM(POS) C C Convert this value to AU using CONVRT. C CALL CONVRT ( DIST, 'KM', 'AU', DIST ) WRITE (*,*) ' Actual distance between Sun and Phobos body ' .// 'centers: ' WRITE (*,'(A,F16.3)') ' (AU):', DIST END Solution Sample Output
Converting UTC Time: 2013-02-10 20:40:00 ET seconds past J2000: 413800866.185 Apparent state of Phobos as seen from PHSRM in the J2000 frame (km, km/s): X = 43.646 Y = 5.698 Z = -20.731 VX = 0.008 VY = -0.002 VZ = -0.005 Apparent position of Earth as seen from PHSRM in the J2000 frame (km): X = -318212438.123 Y = 123091882.078 Z = 59808626.721 One way light time between PHSRM and the apparent position of Earth (seconds): 1155.441 Apparent position of Sun as seen from Phobos in the J2000 frame (km): X = -201733820.746 Y = 39838264.775 Z = 23717367.433 Actual distance between Sun and Phobos body centers: (AU): 1.384 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. \begindata KERNELS_TO_LOAD = ( 'kernels/lsk/naif0009.tls', 'kernels/sclk/phsrm_201008031645.tsc', 'kernels/spk/de421xs.bsp', 'kernels/spk/phobos_kiam_101231_v00.bsp', 'kernels/spk/phsrm_130114_130114_130214_nom2.bsp', 'kernels/fk/phsrm_v00.tf', 'kernels/ck/phsrm_sc_test2_111108_130214_v00.bc', 'kernels/pck/pck00009.tpc' ) \begintext Solution Source Code
PROGRAM XFORM IMPLICIT NONE C C SPICELIB Functions C DOUBLE PRECISION VSEP C C Local Parameters C C C The name of the meta-kernel that lists the kernels C to load into the program. C CHARACTER*(*) METAKR PARAMETER ( METAKR = 'xform.tm' ) C C The length of various string variables. C INTEGER STRLEN PARAMETER ( STRLEN = 50 ) C C Local Variables C CHARACTER*(STRLEN) UTCTIM DOUBLE PRECISION ET DOUBLE PRECISION LTIME DOUBLE PRECISION STATE ( 6 ) DOUBLE PRECISION BFIXST ( 6 ) DOUBLE PRECISION POS ( 3 ) DOUBLE PRECISION SXFMAT ( 6, 6 ) DOUBLE PRECISION PFORM ( 3, 3 ) DOUBLE PRECISION BSIGHT ( 3 ) DOUBLE PRECISION SEP C C Load the kernels that this program requires. We C will need: C C A leapseconds kernel C A spacecraft clock kernel for PHSRM C The necessary ephemerides C A planetary constants file (PCK) C A spacecraft orientation kernel for PHSRM (CK) C A frame kernel (TF) C CALL FURNSH ( METAKR ) C C Prompt the user for the input time string. C CALL PROMPT ( 'Input UTC Time: ', UTCTIM ) WRITE (*,*) 'Converting UTC Time: ', UTCTIM C C Convert UTCTIM to ET. C CALL STR2ET ( UTCTIM, ET ) WRITE (*,'(A,F16.3)') ' ET seconds past J2000: ', ET C C Compute the apparent state of Phobos as seen from PHSRM C in the J2000 reference frame. C CALL SPKEZR ( 'PHOBOS', ET, 'J2000', 'LT+S', . 'PHSRM', STATE, LTIME ) C C Now obtain the transformation from the inertial C J2000 frame to the non-inertial, body-fixed IAU_PHOBOS C frame. Since we'll use this transformation to produce C the apparent state in the IAU_PHOBOS reference frame, C we need to correct the orientation of this frame for C one-way light time; hence we subtract LTIME from ET C in the call below. C CALL SXFORM ( 'J2000', 'IAU_PHOBOS', ET-LTIME, SXFMAT ) C C Now transform the apparent J2000 state into IAU_PHOBOS C with the following matrix multiplication: C CALL MXVG ( SXFMAT, STATE, 6, 6, BFIXST ) C C Display the results. C WRITE (*,*) ' Apparent state of Phobos as seen from ' .// 'PHSRM in the IAU_PHOBOS' WRITE (*,*) ' body-fixed frame (km, km/s):' WRITE (*,'(A,F19.6)') ' X = ', BFIXST(1) WRITE (*,'(A,F19.6)') ' Y = ', BFIXST(2) WRITE (*,'(A,F19.6)') ' Z = ', BFIXST(3) WRITE (*,'(A,F19.6)') ' VX = ', BFIXST(4) WRITE (*,'(A,F19.6)') ' VY = ', BFIXST(5) WRITE (*,'(A,F19.6)') ' VZ = ', BFIXST(6) C C It is worth pointing out, all of the above could have C been done with a single call to SPKEZR: C CALL SPKEZR ( 'PHOBOS', ET, 'IAU_PHOBOS', 'LT+S', . 'PHSRM', STATE, LTIME ) C C Display the results. C WRITE (*,*) ' Apparent state of Phobos as seen from PHSRM ' .// 'in the IAU_PHOBOS' WRITE (*,*) ' body-fixed frame (km, km/s) ' .// 'obtained using SPKEZR directly:' WRITE (*,'(A,F19.6)') ' X = ', STATE(1) WRITE (*,'(A,F19.6)') ' Y = ', STATE(2) WRITE (*,'(A,F19.6)') ' Z = ', STATE(3) WRITE (*,'(A,F19.6)') ' VX = ', STATE(4) WRITE (*,'(A,F19.6)') ' VY = ', STATE(5) WRITE (*,'(A,F19.6)') ' VZ = ', STATE(6) C C Note that the velocity found by using SPKEZR C to compute the state in the IAU_PHOBOS frame differs C at the few mm/second level from that found previously C by calling SPKEZR and then SXFORM. Computing velocity C via a single call to SPKEZR as we've done immediately C above is slightly more accurate because it accounts for C the effect of the rate of change of light time on the C apparent angular velocity of the target's body-fixed C reference frame. C C Now we are to compute the angular separation between C the apparent position of the Sun as seen from the C spacecraft and the normal vector of the solar C arrays. First, compute the apparent position of C the Sun as seen from PHSRM in the J2000 frame. C CALL SPKPOS ( 'SUN', ET, 'J2000', 'LT+S', . 'PHSRM', POS, LTIME ) C C Now set the direction of the solar array normal C at this same epoch. From reading the frame kernel C we know that the solar array normal is nominally the C +X axis of the PHSRM_SPACECRAFT frame defined there. C BSIGHT(1) = 1.0D0 BSIGHT(2) = 0.0D0 BSIGHT(3) = 0.0D0 C C Now compute the rotation matrix from PHSRM_SPACECRAFT into C J2000. C CALL PXFORM ( 'PHSRM_SPACECRAFT', 'J2000', ET, PFORM ) C C And multiply the result to obtain the nominal solar array C normal in the J2000 reference frame. C CALL MXV ( PFORM, BSIGHT, BSIGHT ) C C Lastly compute the angular separation. C CALL CONVRT ( VSEP(BSIGHT, POS), 'RADIANS', . 'DEGREES', SEP ) WRITE (*,*) ' Angular separation between the ' .// 'apparent position of' WRITE (*,*) ' Sun and the PHSRM ' .// 'solar array normal (degrees): ' WRITE (*,'(A,F19.3)') ' ', SEP C C Or, alternately we can work in the spacecraft C frame directly. C CALL SPKPOS ( 'SUN', ET, 'PHSRM_SPACECRAFT', 'LT+S', . 'PHSRM', POS, LTIME ) C C The solar array normal is the X-axis in the C PHSRM_SPACECRAFT frame. C BSIGHT(1) = 1.0D0 BSIGHT(2) = 0.0D0 BSIGHT(3) = 0.0D0 C C Lastly compute the angular separation. C CALL CONVRT ( VSEP(BSIGHT, POS), 'RADIANS', . 'DEGREES', SEP ) WRITE (*,*) ' Angular separation between the ' .// 'apparent position of' WRITE (*,*) ' Sun and the PHSRM ' .// 'solar array normal computed ' WRITE (*,*) ' using vectors in the PHSRM_SPACECRAFT ' .// 'frame (degrees): ' WRITE (*,'(A,F19.3)') ' ', SEP END Solution Sample Output
Converting UTC Time: 2013-02-10 20:40:00 ET seconds past J2000: 413800866.185 Apparent state of Phobos as seen from PHSRM in the IAU_PHOBOS body-fixed frame (km, km/s): X = 32.950479 Y = -35.796727 Z = -0.250178 VX = -0.004507 VY = -0.016033 VZ = -0.000054 Apparent state of Phobos as seen from PHSRM in the IAU_PHOBOS body-fixed frame (km, km/s) obtained using SPKEZR directly: X = 32.950479 Y = -35.796727 Z = -0.250178 VX = -0.004507 VY = -0.016033 VZ = -0.000054 Angular separation between the apparent position of Sun and the PHSRM solar array normal (degrees): 25.807 Angular separation between the apparent position of Sun and the PHSRM solar array normal computed using vectors in the PHSRM_SPACECRAFT frame (degrees): 25.807 Computing Sub-spacecraft and Sub-solar Points (subpts)Task Statement
Learning Goals
Approach
One point worth considering: Which method do you want to use to compute the sub-solar (or sub-observer) point? SolutionSolution Meta-Kernel
KPL/MK This is the meta-kernel used in the solution of the ``Computing Sub-spacecraft and Sub-solar Points'' task in the Remote Sensing Hands On Lesson. \begindata KERNELS_TO_LOAD = ( 'kernels/lsk/naif0009.tls', 'kernels/spk/de421xs.bsp', 'kernels/spk/phobos_kiam_101231_v00.bsp', 'kernels/spk/phsrm_130114_130114_130214_nom2.bsp', 'kernels/pck/pck00009.tpc' 'kernels/fk/phsrm_v00.tf' ) \begintext Solution Source Code
PROGRAM SUBPTS IMPLICIT NONE C C SPICELIB functions C DOUBLE PRECISION VNORM C C Local Parameters C C C The name of the meta-kernel that lists the kernels C to load into the program. C CHARACTER*(*) METAKR PARAMETER ( METAKR = 'subpts.tm' ) C C The length of various string variables. C INTEGER STRLEN PARAMETER ( STRLEN = 50 ) C C Local Variables C CHARACTER*(STRLEN) UTCTIM DOUBLE PRECISION ET DOUBLE PRECISION SPOINT ( 3 ) DOUBLE PRECISION SRFVEC ( 3 ) DOUBLE PRECISION TRGEPC C C Load the kernels that this program requires. We C will need: C C A leapseconds kernel C The necessary ephemerides C A planetary constants file (PCK) C CALL FURNSH ( METAKR ) C C Prompt the user for the input time string. C CALL PROMPT ( 'Input UTC Time: ', UTCTIM ) WRITE (*,*) 'Converting UTC Time: ', UTCTIM C C Convert UTCTIM to ET. C CALL STR2ET ( UTCTIM, ET ) WRITE (*,'(A,F16.3)') ' ET seconds past J2000: ', ET C C Compute the apparent sub-observer point of PHSRM on Phobos. C CALL SUBPNT ( 'NEAR POINT: ELLIPSOID', . 'PHOBOS', ET, 'IAU_PHOBOS', 'LT+S', . 'PHSRM', SPOINT, TRGEPC, SRFVEC ) WRITE (*,*) ' Apparent sub-observer point of PHSRM ' .// 'on Phobos in the' WRITE (*,*) ' IAU_PHOBOS frame (km):' WRITE (*,'(A,F16.3)') ' X = ', SPOINT(1) WRITE (*,'(A,F16.3)') ' Y = ', SPOINT(2) WRITE (*,'(A,F16.3)') ' Z = ', SPOINT(3) WRITE (*,'(A,F16.3)') ' ALT = ', VNORM(SRFVEC) C C Compute the apparent sub-solar point on Phobos as seen C from PHSRM. C CALL SUBSLR ( 'NEAR POINT: ELLIPSOID', . 'PHOBOS', ET, 'IAU_PHOBOS', 'LT+S', . 'PHSRM', SPOINT, TRGEPC, SRFVEC ) WRITE (*,*) ' Apparent sub-solar point on Phobos as ' .// 'seen from PHSRM in' WRITE (*,*) ' the IAU_PHOBOS frame (km):' WRITE (*,'(A,F16.3)') ' X = ', SPOINT(1) WRITE (*,'(A,F16.3)') ' Y = ', SPOINT(2) WRITE (*,'(A,F16.3)') ' Z = ', SPOINT(3) END Solution Sample Output
Converting UTC Time: 2013-02-10 20:40:00 ET seconds past J2000: 413800866.185 Apparent sub-observer point of PHSRM on Phobos in the IAU_PHOBOS frame (km): X = -9.501 Y = 7.897 Z = 0.040 ALT = 36.446 Apparent sub-solar point on Phobos as seen from PHSRM in the IAU_PHOBOS frame (km): X = -7.904 Y = 8.282 Z = -2.986 Intersecting Vectors with a Triaxial Ellipsoid (fovint)Task Statement
At each point of intersection compute the following:
Use this program to compute values at the epoch:
Learning Goals
Approach
SolutionSolution Meta-Kernel
KPL/MK This is the meta-kernel used in the solution of the ``Intersecting Vectors with a Triaxial Ellipsoid'' task in the Remote Sensing Hands On Lesson. \begindata KERNELS_TO_LOAD = ( 'kernels/lsk/naif0009.tls', 'kernels/sclk/phsrm_201008031645.tsc', 'kernels/spk/de421xs.bsp', 'kernels/spk/phobos_kiam_101231_v00.bsp', 'kernels/spk/phsrm_130114_130114_130214_nom2.bsp', 'kernels/fk/phsrm_v00.tf', 'kernels/ck/phsrm_sc_test2_111108_130214_v00.bc', 'kernels/pck/pck00009.tpc', 'kernels/ik/phsrm_tsns_v00.ti' ) \begintext Solution Source Code
PROGRAM FOVINT IMPLICIT NONE C C SPICELIB functions C DOUBLE PRECISION DPR C C Local Parameters C C C The name of the meta-kernel that lists the kernels C to load into the program. C CHARACTER*(*) METAKR PARAMETER ( METAKR = 'fovint.tm' ) C C The length of various string variables. C INTEGER STRLEN PARAMETER ( STRLEN = 50 ) C C The maximum number of boundary corner vectors C we can retrieve. C INTEGER BCVLEN PARAMETER ( BCVLEN = 5 ) C C Local Variables C CHARACTER*(STRLEN) AMPM CHARACTER*(STRLEN) INSFRM CHARACTER*(STRLEN) SHAPE CHARACTER*(STRLEN) TIME CHARACTER*(STRLEN) UTCTIM CHARACTER*(STRLEN) VECNAM ( BCVLEN ) DOUBLE PRECISION BOUNDS ( 3, BCVLEN ) DOUBLE PRECISION BSIGHT ( 3 ) DOUBLE PRECISION EMISSN DOUBLE PRECISION ET DOUBLE PRECISION LAT DOUBLE PRECISION LON DOUBLE PRECISION PHASE DOUBLE PRECISION POINT ( 3 ) DOUBLE PRECISION RADIUS DOUBLE PRECISION SOLAR DOUBLE PRECISION SRFVEC ( 3 ) DOUBLE PRECISION TRGEPC INTEGER HR INTEGER I INTEGER MN INTEGER N INTEGER NACID INTEGER PHOEID INTEGER SC LOGICAL FOUND C C Load the kernels that this program requires. We C will need: C C A leapseconds kernel. C A SCLK kernel for PHSRM. C Any necessary ephemerides. C The PHSRM frame kernel. C A PHSRM C-kernel. C A PCK file with Phobos constants. C The PHSRM I-kernel. C CALL FURNSH ( METAKR ) C C Prompt the user for the input time string. C CALL PROMPT ( 'Input UTC Time: ', UTCTIM ) WRITE (*,*) 'Converting UTC Time: ', UTCTIM C C Convert UTCTIM to ET. C CALL STR2ET ( UTCTIM, ET ) WRITE (*,'(A,F16.3)') ' ET seconds past J2000: ', ET C C Now we need to obtain the FOV configuration of the NAC 1 C camera. To do this we will need the ID code for C PHSRM_TSNS_NAC_1. C CALL BODN2C ( 'PHSRM_TSNS_NAC_1', NACID, FOUND ) C C Stop the program if the code was not found. C IF ( .NOT. FOUND ) THEN WRITE (*,*) 'Unable to locate the ID code for ' . // 'PHSRM_TSNS_NAC_1' CALL BYEBYE ( 'FAILURE' ) END IF C C Now retrieve the field of view parameters. C CALL GETFOV ( NACID, BCVLEN, SHAPE, INSFRM, . BSIGHT, N, BOUNDS ) C C Rather than treat BSIGHT as a separate vector, C copy it into the last slot of BOUNDS. C CALL MOVED ( BSIGHT, 3, BOUNDS(1,5) ) C C Define names for each of the vectors for display C purposes. C VECNAM (1) = 'Boundary Corner 1' VECNAM (2) = 'Boundary Corner 2' VECNAM (3) = 'Boundary Corner 3' VECNAM (4) = 'Boundary Corner 4' VECNAM (5) = 'PHSRM TSNS NAC 1 Boresight' C C Now perform the same set of calculations for each C vector listed in the BOUNDS array. C DO I = 1, 5 C C Call SINCPT to determine coordinates of the C intersection of this vector with the surface C of Phobos. C CALL SINCPT ( 'Ellipsoid', 'PHOBOS', ET, . 'IAU_PHOBOS', 'LT+S', 'PHSRM', . INSFRM, BOUNDS(1,I), POINT, . TRGEPC, SRFVEC, FOUND ) C C Check the found flag. Display a message if the point C of intersection was not found, otherwise continue with C the calculations. C WRITE (*,*) 'Vector: ', VECNAM(I) IF ( .NOT. FOUND ) THEN WRITE (*,*) 'No intersection point found at ' . // 'this epoch for this vector.' ELSE C C Now, we have discovered a point of intersection. C Start by displaying the position vector in the C IAU_PHOBOS frame of the intersection. C WRITE (*,*) ' Position vector of ' . // 'surface intercept in ' . // 'the IAU_PHOBOS frame (km):' WRITE (*,'(A,F16.3)') ' X = ', POINT(1) WRITE (*,'(A,F16.3)') ' Y = ', POINT(2) WRITE (*,'(A,F16.3)') ' Z = ', POINT(3) C C Display the planetocentric latitude and longitude C of the intercept. C CALL RECLAT ( POINT, RADIUS, LON, LAT ) WRITE (*,*) ' Planetocentric coordinates of the ' . // 'intercept (degrees):' WRITE (*,'(A,F16.3)') ' LAT = ', LAT * DPR() WRITE (*,'(A,F16.3)') ' LON = ', LON * DPR() C C Compute the illumination angles at this C point. C CALL ILUMIN ( 'Ellipsoid', 'PHOBOS', ET, . 'IAU_PHOBOS', 'LT+S', 'PHSRM', . POINT, TRGEPC, SRFVEC, . PHASE, SOLAR, EMISSN ) WRITE (*,'(A,F16.3)') ' Phase angle (degrees):' . // ' ', PHASE * DPR() WRITE (*,'(A,F16.3)') ' Solar incidence angle ' . // '(degrees): ', SOLAR * DPR() WRITE (*,'(A,F16.3)') ' Emission angle (degree' . // 's): ', EMISSN* DPR() END IF WRITE (*,*) ' ' END DO C C Lastly compute the local solar time at the boresight C intersection. C IF ( FOUND ) THEN C C Get Phobos ID. C CALL BODN2C ( 'PHOBOS', PHOEID, FOUND ) C C Stop the program if the code was not found. C IF ( .NOT. FOUND ) THEN WRITE (*,*) 'Unable to locate the ID code for ' . // 'PHOBOS' CALL BYEBYE ( 'FAILURE' ) END IF C C Compute local time corresponding to the TDB light time C corrected epoch at the intercept. C CALL ET2LST ( TRGEPC, . PHOEID, . LON, . 'PLANETOCENTRIC', . HR, . MN, . SC, . TIME, . AMPM ) WRITE (*,*) ' Local Solar Time at boresight ' . // 'intercept (24 Hour Clock): ' WRITE (*,*) ' ', TIME ELSE WRITE (*,*) ' No boresight intercept to compute ' . // 'local solar time.' END IF END Solution Sample Output
Converting UTC Time: 2013-02-10 20:40:00 ET seconds past J2000: 413800866.185 Vector: Boundary Corner 1 Position vector of surface intercept in the IAU_PHOBOS frame (km): X = -7.873 Y = 9.060 Z = -0.203 Planetocentric coordinates of the intercept (degrees): LAT = -0.971 LON = 130.989 Phase angle (degrees): 26.338 Solar incidence angle (degrees): 22.465 Emission angle (degrees): 12.132 Vector: Boundary Corner 2 Position vector of surface intercept in the IAU_PHOBOS frame (km): X = -7.889 Y = 9.044 Z = 0.340 Planetocentric coordinates of the intercept (degrees): LAT = 1.622 LON = 131.101 Phase angle (degrees): 25.527 Solar incidence angle (degrees): 26.790 Emission angle (degrees): 12.108 Vector: Boundary Corner 3 Position vector of surface intercept in the IAU_PHOBOS frame (km): X = -8.355 Y = 8.748 Z = 0.316 Planetocentric coordinates of the intercept (degrees): LAT = 1.498 LON = 133.684 Phase angle (degrees): 25.279 Solar incidence angle (degrees): 26.523 Emission angle (degrees): 8.968 Vector: Boundary Corner 4 Position vector of surface intercept in the IAU_PHOBOS frame (km): X = -8.338 Y = 8.764 Z = -0.226 Planetocentric coordinates of the intercept (degrees): LAT = -1.068 LON = 133.574 Phase angle (degrees): 26.096 Solar incidence angle (degrees): 22.151 Emission angle (degrees): 9.076 Vector: PHSRM TSNS NAC 1 Boresight Position vector of surface intercept in the IAU_PHOBOS frame (km): X = -8.120 Y = 8.909 Z = 0.057 Planetocentric coordinates of the intercept (degrees): LAT = 0.270 LON = 132.344 Phase angle (degrees): 25.807 Solar incidence angle (degrees): 24.456 Emission angle (degrees): 10.241 Local Solar Time at boresight intercept (24 Hour Clock): 12:34:36 |