Remote Sensing Hands-On Lesson (C) |
Table of ContentsRemote Sensing Hands-On Lesson (C) Overview Note About HTML Links References Tutorials Required Readings The Permuted Index Source Code Header Comments Kernels Used CSPICE 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 Remote Sensing Hands-On Lesson (C)
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/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 DataThese tutorials are available from the NAIF ftp server at JPL:
http://naif.jpl.nasa.gov/naif/tutorials.html Required Readings
Name Lesson steps/functions 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 CSPICE functions perform a particular function of interest as well as the name of the source module that contains the function. 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 naif0008.tls LSK Generic LSK 2 cas00084.tsc SCLK Cassini SCLK 3 981005_PLTEPH-DE405S.bsp SPK Solar System Ephemeris 4 020514_SE_SAT105.bsp SPK Saturnian Satellite Ephemeris 5 030201AP_SK_SM546_T45.bsp SPK Cassini Spacecraft SPK 6 cas_v37.tf FK Cassini FK 7 04135_04171pc_psiv2.bc CK Cassini Spacecraft CK 8 cpck05Mar2004.tpc PCK Cassini Project PCK 9 cas_iss_v09.ti IK ISS Instrument KernelThese SPICE kernels are included in the lesson package available from the NAIF server at JPL:
ftp://naif.jpl.nasa.gov/pub/naif/toolkit_docs/Lessons/ CSPICE Modules Used
CHAPTER EXERCISE FUNCTIONS NON-VOID KERNELS ------- --------- --------- --------- ------- 1 convtm furnsh_c 1,2 prompt_c str2et_c etcal_c timout_c sce2c_c sce2s_c 2 getsta furnsh_c vnorm_c 1,3-5 prompt_c str2et_c spkezr_c spkpos_c convrt_c 3 xform furnsh_c vsep_c 1-8 prompt_c str2et_c spkezr_c sxform_c mxvg_c spkpos_c pxform_c mxv_c convrt_c 4 subpts furnsh_c 1,3-5,8 prompt_c str2et_c subpt_c subsol_c 5 fovint furnsh_c dpr_c 1-9 prompt_c str2et_c bodn2c_c getfov_c sincpt_c reclat_c ilumin_c et2lst_cRefer to the headers of the various functions 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_c or timout_c. 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/naif0008.tls', 'kernels/sclk/cas00084.tsc' ) \begintext Solution Source Code
#include <stdio.h> /* Standard CSPICE User Include File */ #include "SpiceUsr.h" /* Local Parameters */ #define METAKR "convtm.tm" #define SCLKID -82 #define STRLEN 50 int main (void) { /* Local Variables */ SpiceChar calet [STRLEN]; SpiceChar sclkst [STRLEN]; SpiceChar utctim [STRLEN]; SpiceDouble et; /* Load the kernels this program requires. Both the spacecraft clock kernel and a leapseconds kernel should be listed in the meta-kernel. */ furnsh_c ( METAKR ); /* Prompt the user for the input time string. */ prompt_c ( "Input UTC Time: ", STRLEN, utctim ); printf ( "Converting UTC Time: %s\n", utctim ); /* Convert utctim to ET. */ str2et_c ( utctim, &et ); printf ( " ET Seconds Past J2000: %16.3f\n", et ); /* Now convert ET to a calendar time string. This can be accomplished in two ways. */ etcal_c ( et, STRLEN, calet ); printf ( " Calendar ET (etcal_c): %s\n", calet ); /* Or use timout_c for finer control over the output format. The picture below was built by examining the header of timout_c. */ timout_c ( et, "YYYY-MON-DDTHR:MN:SC ::TDB", STRLEN, calet ); printf ( " Calendar ET (timout_c): %s\n", calet ); /* Convert ET to spacecraft clock time. */ sce2s_c ( SCLKID, et, STRLEN, sclkst ); printf ( " Spacecraft Clock Time: %s\n", sclkst ); return(0); } Solution Sample Output
Converting UTC Time: 2004 jun 11 19:32:00 ET Seconds Past J2000: 140254384.185 Calendar ET (etcal_c): 2004 JUN 11 19:33:04.184 Calendar ET (timout_c): 2004-JUN-11T19:33:04 Spacecraft Clock Time: 1/1465674964.105 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. \begindata KERNELS_TO_LOAD = ( 'kernels/lsk/naif0008.tls', 'kernels/spk/981005_PLTEPH-DE405S.bsp', 'kernels/spk/020514_SE_SAT105.bsp', 'kernels/spk/030201AP_SK_SM546_T45.bsp' ) \begintext Solution Source Code
#include <stdio.h> /* Standard CSPICE User Include File */ #include "SpiceUsr.h" /* Local Parameters */ #define METAKR "getsta.tm" #define STRLEN 50 int main (void) { /* Local Variables */ SpiceChar utctim [STRLEN]; SpiceDouble dist; SpiceDouble et; SpiceDouble ltime; SpiceDouble pos [3]; SpiceDouble state [6]; /* Load the kernels that this program requires. We will need a leapseconds kernel to convert input UTC time strings into ET. We also will need the necessary SPK files with coverage for the bodies in which we are interested. */ furnsh_c ( METAKR ); /* Prompt the user for the input time string. */ prompt_c ( "Input UTC Time: ", STRLEN, utctim ); printf ( "Converting UTC Time: %s\n", utctim ); /* Convert utctim to ET. */ str2et_c ( utctim, &et ); printf ( " ET seconds past J2000: %16.3f\n", et ); /* Compute the apparent state of Phoebe as seen from CASSINI in the J2000 frame. All of the ephemeris readers return states in units of kilometers and kilometers per second. */ spkezr_c ( "PHOEBE", et, "J2000", "LT+S", "CASSINI", state, <ime ); printf ( " Apparent state of Phoebe as seen " "from CASSINI in the J2000\n" ); printf ( " frame (km, km/s):\n" ); printf ( " X = %16.3f\n", state[0] ); printf ( " Y = %16.3f\n", state[1] ); printf ( " Z = %16.3f\n", state[2] ); printf ( " VX = %16.3f\n", state[3] ); printf ( " VY = %16.3f\n", state[4] ); printf ( " VZ = %16.3f\n", state[5] ); /* Compute the apparent position of Earth as seen from CASSINI in the J2000 frame. Note: We could have continued using spkezr_c and simply ignored the velocity components. */ spkpos_c ( "EARTH", et, "J2000", "LT+S", "CASSINI", pos, <ime ); printf ( " Apparent position of Earth as " "seen from CASSINI in the J2000\n" ); printf ( " frame (km): \n" ); printf ( " X = %16.3f\n", pos[0] ); printf ( " Y = %16.3f\n", pos[1] ); printf ( " Z = %16.3f\n", pos[2] ); /* We need only display LTIME, as it is precisely the light time in which we are interested. */ printf ( " One way light time between CASSINI and " "the apparent position\n" ); printf ( " of Earth (seconds): %16.3f\n", ltime ); /* Compute the apparent position of the Sun as seen from Phoebe in the J2000 frame. */ spkpos_c ( "SUN", et, "J2000", "LT+S", "PHOEBE", pos, <ime ); printf ( " Apparent position of Sun as seen " "from Phoebe in the\n" ); printf ( " J2000 frame (km): \n" ); printf ( " X = %16.3f\n", pos[0] ); printf ( " Y = %16.3f\n", pos[1] ); printf ( " Z = %16.3f\n", pos[2] ); /* Now we need to compute the actual distance between the Sun and Phoebe. The above SPKPOS call gives us the apparent distance, so we need to adjust our aberration correction appropriately. */ spkpos_c ( "SUN", et, "J2000", "NONE", "PHOEBE", pos, <ime ); /* Compute the distance between the body centers in kilometers. */ dist = vnorm_c ( pos ); /* Convert this value to AU using convrt_c. */ convrt_c ( dist, "KM", "AU", &dist ); printf ( " Actual distance between Sun and " "Phoebe body centers:\n" ); printf ( " (AU): %16.3f\n", dist ); return(0); } Solution Sample Output
Converting UTC Time: 2004 JUN 11 19:32:00 ET seconds past J2000: 140254384.185 Apparent state of Phoebe as seen from CASSINI in the J2000 frame (km, km/s): X = -119.921 Y = 2194.139 Z = -57.639 VX = -5.980 VY = -2.119 VZ = -0.295 Apparent position of Earth as seen from CASSINI in the J2000 frame (km): X = 353019393.123 Y = -1328180352.140 Z = -568134171.697 One way light time between CASSINI and the apparent position of Earth (seconds): 4960.427 Apparent position of Sun as seen from Phoebe in the J2000 frame (km): X = 376551465.272 Y = -1190495630.303 Z = -508438699.110 Actual distance between Sun and Phoebe body centers: (AU): 9.012 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/naif0008.tls', 'kernels/sclk/cas00084.tsc', 'kernels/spk/981005_PLTEPH-DE405S.bsp', 'kernels/spk/020514_SE_SAT105.bsp', 'kernels/spk/030201AP_SK_SM546_T45.bsp', 'kernels/fk/cas_v37.tf', 'kernels/ck/04135_04171pc_psiv2.bc', 'kernels/pck/cpck05Mar2004.tpc' ) \begintext Solution Source Code
#include <stdio.h> /* Standard CSPICE User Include File */ #include "SpiceUsr.h" /* Local Parameters */ #define METAKR "xform.tm" #define STRLEN 50 int main (void) { /* Local Variables */ SpiceChar utctim [STRLEN]; SpiceDouble et; SpiceDouble ltime; SpiceDouble state [6]; SpiceDouble bfixst [6]; SpiceDouble pos [3]; SpiceDouble sform [6][6]; SpiceDouble pform [3][3]; SpiceDouble bsight [3]; SpiceDouble sep; /* Load the kernels that this program requires. We will need: A leapseconds kernel A spacecraft clock kernel for CASSINI The necessary ephemerides A planetary constants file (PCK) A spacecraft orientation kernel for CASSINI (CK) A frame kernel (TF) */ furnsh_c ( METAKR ); /* Prompt the user for the input time string. */ prompt_c ( "Input UTC Time: ", STRLEN, utctim ); printf ( "Converting UTC Time: %s\n", utctim ); /* Convert utctim to ET. */ str2et_c ( utctim, &et ); printf ( " ET seconds past J2000: %16.3f\n", et ); /* Compute the apparent state of Phoebe as seen from CASSINI in the J2000 reference frame. */ spkezr_c ( "PHOEBE", et, "J2000", "LT+S", "CASSINI", state, <ime ); /* Now obtain the transformation from the inertial J2000 frame to the non-inertial body-fixed IAU_PHOEBE frame. Since we want the apparent position, we need to subtract ltime from et. */ sxform_c ( "J2000", "IAU_PHOEBE", et-ltime, sform ); /* Now rotate the apparent J2000 state into IAU_PHOEBE with the following matrix multiplication: */ mxvg_c ( sform, state, 6, 6, bfixst ); /* Display the results. */ printf ( " Apparent state of Phoebe as seen " "from CASSINI in the IAU_PHOEBE\n" ); printf ( " body-fixed frame (km, km/s):\n"); printf ( " X = %19.6f\n", bfixst[0] ); printf ( " Y = %19.6f\n", bfixst[1] ); printf ( " Z = %19.6f\n", bfixst[2] ); printf ( " VX = %19.6f\n", bfixst[3] ); printf ( " VY = %19.6f\n", bfixst[4] ); printf ( " VZ = %19.6f\n", bfixst[5] ); /* It is worth pointing out, all of the above could have been done with a single use of spkezr_c: */ spkezr_c ( "PHOEBE", et, "IAU_PHOEBE", "LT+S", "CASSINI", state, <ime ); /* Display the results. */ printf ( " Apparent state of Phoebe as seen " "from CASSINI in the IAU_PHOEBE\n" ); printf ( " body-fixed frame (km, km/s) " "obtained using spkezr_c directly:\n" ); printf ( " X = %19.6f\n", state[0] ); printf ( " Y = %19.6f\n", state[1] ); printf ( " Z = %19.6f\n", state[2] ); printf ( " VX = %19.6f\n", state[3] ); printf ( " VY = %19.6f\n", state[4] ); printf ( " VZ = %19.6f\n", state[5] ); /* Note that the velocity found by using spkezr_c to compute the state in the IAU_PHOEBE frame differs at the few mm/second level from that found previously by calling spkezr_c and then sxform_c. Computing velocity via a single call to spkezr_c as we've done immediately above is slightly more accurate because it accounts for the effect of the rate of change of light time on the apparent angular velocity of the target's body-fixed reference frame. Now we are to compute the angular separation between the apparent position of the Earth as seen from the orbiter and the nominal boresight of the high gain antenna. First, compute the apparent position of the Earth as seen from CASSINI in the J2000 frame. */ spkpos_c ( "EARTH", et, "J2000", "LT+S", "CASSINI", pos, <ime ); /* Now compute the location of the antenna boresight at this same epoch. From reading the frame kernel we know that the antenna boresight is nominally the +Z axis of the CASSINI_HGA frame defined there. */ bsight[0] = 0.0; bsight[1] = 0.0; bsight[2] = 1.0; /* Now compute the rotation matrix from CASSINI_HGA into J2000. */ pxform_c ( "CASSINI_HGA", "J2000", et, pform ); /* And multiply the result to obtain the nominal antenna boresight in the J2000 reference frame. */ mxv_c ( pform, bsight, bsight ); /* Lastly compute the angular separation. */ convrt_c ( vsep_c(bsight, pos), "RADIANS", "DEGREES", &sep ); printf ( " Angular separation between the " "apparent position of\n" ); printf ( " Earth and the CASSINI high " "gain antenna boresight (degrees):\n"); printf ( " %16.3f\n", sep ); /* Or alternatively we can work in the antenna frame directly. */ spkpos_c ( "EARTH", et, "CASSINI_HGA", "LT+S", "CASSINI", pos, <ime ); /* The antenna boresight is the Z-axis in the CASSINI_HGA frame. */ bsight[0] = 0.0; bsight[1] = 0.0; bsight[2] = 1.0; /* Lastly compute the angular separation. */ convrt_c ( vsep_c(bsight, pos), "RADIANS", "DEGREES", &sep ); printf ( " Angular separation between the " "apparent position of\n" ); printf ( " Earth and the CASSINI high " "gain antenna boresight computed\n" ); printf ( " using vectors in the CASSINI_HGA " "frame (degrees):\n" ); printf ( " %16.3f\n", sep ); return(0); } Solution Sample Output
Converting UTC Time: 2004 JUN 11 19:32:00 ET seconds past J2000: 140254384.185 Apparent state of Phoebe as seen from CASSINI in the IAU_PHOEBE body-fixed frame (km, km/s): X = -1982.639762 Y = -934.530471 Z = -166.562595 VX = 3.970833 VY = -3.812498 VZ = -2.371663 Apparent state of Phoebe as seen from CASSINI in the IAU_PHOEBE body-fixed frame (km, km/s) obtained using spkezr_c directly: X = -1982.639762 Y = -934.530471 Z = -166.562595 VX = 3.970832 VY = -3.812496 VZ = -2.371663 Angular separation between the apparent position of Earth and the CASSINI high gain antenna boresight (degrees): 71.924 Angular separation between the apparent position of Earth and the CASSINI high gain antenna boresight computed using vectors in the CASSINI_HGA frame (degrees): 71.924 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/naif0008.tls', 'kernels/spk/981005_PLTEPH-DE405S.bsp', 'kernels/spk/020514_SE_SAT105.bsp', 'kernels/spk/030201AP_SK_SM546_T45.bsp', 'kernels/pck/cpck05Mar2004.tpc' ) \begintext Solution Source Code
#include <stdio.h> /* Standard CSPICE User Include File */ #include "SpiceUsr.h" /* Local Parameters */ #define METAKR "subpts.tm" #define STRLEN 50 int main (void) { /* Local Variables */ SpiceChar utctim [STRLEN]; SpiceDouble et; SpiceDouble spoint [3]; SpiceDouble srfvec [3]; SpiceDouble trgepc; /* Load the kernels that this program requires. We will need: A leapseconds kernel The necessary ephemerides A planetary constants file (PCK) */ furnsh_c ( METAKR ); /* Prompt the user for the input time string. */ prompt_c ( "Input UTC Time: ", STRLEN, utctim ); printf ( "Converting UTC Time: %s\n", utctim ); /* Convert utctim to ET. */ str2et_c ( utctim, &et ); printf ( " ET seconds past J2000: %16.3f\n", et ); /* Compute the apparent sub-observer point of CASSINI on Phoebe. */ subpnt_c ( "NEAR POINT: ELLIPSOID", "PHOEBE", et, "IAU_PHOEBE", "LT+S", "CASSINI", spoint, &trgepc, srfvec ); printf ( " Apparent sub-observer point of CASSINI " "on Phoebe in the\n" ); printf ( " IAU_PHOEBE frame (km):\n" ); printf ( " X = %16.3f\n", spoint[0] ); printf ( " Y = %16.3f\n", spoint[1] ); printf ( " Z = %16.3f\n", spoint[2] ); printf ( " ALT = %16.3f\n", vnorm_c(srfvec) ); /* Compute the apparent sub-solar point on Phoebe as seen from CASSINI. */ subslr_c ( "NEAR POINT: ELLIPSOID", "PHOEBE", et, "IAU_PHOEBE", "LT+S", "CASSINI", spoint, &trgepc, srfvec ); printf ( " Apparent sub-solar point on Phoebe " "as seen from CASSINI in\n" ); printf ( " the IAU_PHOEBE frame (km):\n" ); printf ( " X = %16.3f\n", spoint[0] ); printf ( " Y = %16.3f\n", spoint[1] ); printf ( " Z = %16.3f\n", spoint[2] ); return(0); } Solution Sample Output
Converting UTC Time: 2004 JUN 11 19:32:00 ET seconds past J2000: 140254384.185 Apparent sub-observer point of CASSINI on Phoebe in the IAU_PHOEBE frame (km): X = 104.498 Y = 45.269 Z = 7.383 ALT = 2084.116 Apparent sub-solar point on Phoebe as seen from CASSINI in the IAU_PHOEBE frame (km): X = 78.681 Y = 76.879 Z = -21.885 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/naif0008.tls', 'kernels/sclk/cas00084.tsc', 'kernels/spk/981005_PLTEPH-DE405S.bsp', 'kernels/spk/020514_SE_SAT105.bsp', 'kernels/spk/030201AP_SK_SM546_T45.bsp', 'kernels/fk/cas_v37.tf', 'kernels/ck/04135_04171pc_psiv2.bc', 'kernels/pck/cpck05Mar2004.tpc', 'kernels/ik/cas_iss_v09.ti' ) \begintext Solution Source Code
#include <stdio.h> /* Standard CSPICE User Include File */ #include "SpiceUsr.h" #include <stdlib.h> /* Local Parameters */ #define METAKR "fovint.tm" #define STRLEN 50 #define BCVLEN 5 int main (void) { /* Local Variables */ SpiceChar ampm [STRLEN]; SpiceChar insfrm [STRLEN]; SpiceChar shape [STRLEN]; SpiceChar time [STRLEN]; SpiceChar utctim [STRLEN]; SpiceChar *vecnam[] = { "Boundary Corner 1", "Boundary Corner 2", "Boundary Corner 3", "Boundary Corner 4", "Cassini NAC Boresight" }; SpiceDouble bounds [BCVLEN][3]; SpiceDouble bsight [3]; SpiceDouble emissn; SpiceDouble et; SpiceDouble lat; SpiceDouble lon; SpiceDouble phase; SpiceDouble point [3]; SpiceDouble radius; SpiceDouble solar; SpiceDouble srfvec [3]; SpiceDouble trgepc; SpiceInt hr; SpiceInt i; SpiceInt mn; SpiceInt n; SpiceInt nacid; SpiceInt phoeid; SpiceInt sc; SpiceBoolean found; /* Load the kernels that this program requires. We will need: A leapseconds kernel. A SCLK kernel for CASSINI. Any necessary ephemerides. The CASSINI frame kernel. A CASSINI C-kernel. A PCK file with Phoebe constants. The CASSINI ISS I-kernel. */ furnsh_c ( METAKR ); /* Prompt the user for the input time string. */ prompt_c ( "Input UTC Time: ", STRLEN, utctim ); printf ( "Converting UTC Time: %s\n", utctim ); /* Convert utctim to ET. */ str2et_c ( utctim, &et ); printf ( " ET seconds past J2000: %16.3f\n", et ); /* Now we need to obtain the FOV configuration of the ISS NAC camera. To do this we will need the ID code for CASSINI_ISS_NAC. */ bodn2c_c ( "CASSINI_ISS_NAC", &nacid, &found ); /* Stop the program if the code was not found. */ if ( !found ) { printf ( "Unable to locate the ID code for " "CASSINI_ISS_NAC\n" ); exit ( EXIT_FAILURE ); } /* Now retrieve the field of view parameters. */ getfov_c ( nacid, BCVLEN, STRLEN, STRLEN, shape, insfrm, bsight, &n, bounds ); /* Rather than treat BSIGHT as a separate vector, copy it into the last slot of BOUNDS. */ for ( i=0; i<3; i++ ) { bounds[4][i] = bsight[i]; } /* Now perform the same set of calculations for each vector listed in the BOUNDS array. */ for ( i=0; i<5; i++ ) { /* Call sincpt_c to determine coordinates of the intersection of this vector with the surface of Phoebe. */ sincpt_c ( "Ellipsoid", "PHOEBE", et, "IAU_PHOEBE", "LT+S", "CASSINI", insfrm, bounds[i], point, &trgepc, srfvec, &found ); /* Check the found flag. Display a message if the point of intersection was not found, otherwise continue with the calculations. */ printf ( "Vector: %s\n", vecnam[i] ); if ( !found ) { printf ( "No intersection point found at " "this epoch for this vector.\n" ); } else { /* Now, we have discovered a point of intersection. Start by displaying the position vector in the IAU_PHOEBE frame of the intersection. */ printf ( " Position vector of surface intercept " "in the IAU_PHOEBE frame (km):\n" ); printf ( " X = %16.3f\n", point[0] ); printf ( " Y = %16.3f\n", point[1] ); printf ( " Z = %16.3f\n", point[2] ); /* Display the planetocentric latitude and longitude of the intercept. */ reclat_c ( point, &radius, &lon, &lat ); printf ( " Planetocentric coordinates of " "the intercept (degrees):\n" ); printf ( " LAT = %16.3f\n", lat * dpr_c() ); printf ( " LON = %16.3f\n", lon * dpr_c() ); /* Compute the illumination angles at this point. */ ilumin_c ( "Ellipsoid", "PHOEBE", et, "IAU_PHOEBE", "LT+S", "CASSINI", point, &trgepc, srfvec, &phase, &solar, &emissn ); printf ( " Phase angle (degrees): " "%16.3f\n", phase * dpr_c() ); printf ( " Solar incidence angle (degrees): " "%16.3f\n", solar * dpr_c() ); printf ( " Emission angle (degrees): " "%16.3f\n", emissn * dpr_c() ); } printf ( "\n" ); } /* Lastly compute the local solar time at the boresight intersection. */ if ( found ) { /* Get ID code of Phoebe. */ bodn2c_c ( "PHOEBE", &phoeid, &found ); /* The ID code for PHOEBE is built-in to the library. However, it is good programming practice to get in the habit of checking your found-flags. */ if ( !found ) { printf ( "Unable to locate the body ID code " "for Phoebe.\n" ); exit ( EXIT_FAILURE ); } /* Compute local solar time corresponding to the TDB light time corrected epoch at the intercept. */ et2lst_c ( trgepc, phoeid, lon, "PLANETOCENTRIC", STRLEN, STRLEN, &hr, &mn, &sc, time, ampm ); printf ( " Local Solar Time at boresight " "intercept (24 Hour Clock):\n" ); printf ( " %s\n", time ); } else { printf ( " No boresight intercept to " "compute local solar time.\n" ); } return(0); } Solution Sample Output
Converting UTC Time: 2004 JUN 11 19:32:00 ET seconds past J2000: 140254384.185 Vector: Boundary Corner 1 Position vector of surface intercept in the IAU_PHOEBE frame (km): X = 91.026 Y = 67.190 Z = 2.030 Planetocentric coordinates of the intercept (degrees): LAT = 1.028 LON = 36.433 Phase angle (degrees): 28.110 Solar incidence angle (degrees): 16.121 Emission angle (degrees): 14.627 Vector: Boundary Corner 2 Position vector of surface intercept in the IAU_PHOEBE frame (km): X = 89.991 Y = 66.726 Z = 14.733 Planetocentric coordinates of the intercept (degrees): LAT = 7.492 LON = 36.556 Phase angle (degrees): 27.894 Solar incidence angle (degrees): 22.894 Emission angle (degrees): 14.988 Vector: Boundary Corner 3 Position vector of surface intercept in the IAU_PHOEBE frame (km): X = 80.963 Y = 76.643 Z = 14.427 Planetocentric coordinates of the intercept (degrees): LAT = 7.373 LON = 43.430 Phase angle (degrees): 28.171 Solar incidence angle (degrees): 21.315 Emission angle (degrees): 21.977 Vector: Boundary Corner 4 Position vector of surface intercept in the IAU_PHOEBE frame (km): X = 81.997 Y = 77.106 Z = 1.699 Planetocentric coordinates of the intercept (degrees): LAT = 0.865 LON = 43.239 Phase angle (degrees): 28.385 Solar incidence angle (degrees): 13.882 Emission angle (degrees): 21.763 Vector: Cassini NAC Boresight Position vector of surface intercept in the IAU_PHOEBE frame (km): X = 86.390 Y = 72.089 Z = 8.255 Planetocentric coordinates of the intercept (degrees): LAT = 4.196 LON = 39.844 Phase angle (degrees): 28.139 Solar incidence angle (degrees): 18.247 Emission angle (degrees): 17.858 Local Solar Time at boresight intercept (24 Hour Clock): 11:31:50 |