| 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 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-spacecraft and Sub-solar Points (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 a Triaxial Ellipsoid (fovint) Task Statement Learning Goals Approach Solution Solution Meta-Kernel Solution Source Code Solution Sample Output Extra Credit 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
1. Generic LSK:
naif0012.tls
2. ExoMars-16 TGO SCLK:
em16_tgo_step_20160414.tsc
3. Solar System Ephemeris SPK, subsetted to cover only the time
range of interest:
de430.bsp
4. Martian Satellite Ephemeris SPK, subsetted to cover only the
time range of interest:
mar085.bsp
5. ExoMars-16 TGO Spacecraft Trajectory SPK, subsetted to cover
only the time range of interest:
em16_tgo_mlt_20171205_20230115_v01.bsp
6. ExoMars-16 TGO FK:
em16_tgo_v07.tf
7. ExoMars-16 TGO Spacecraft CK, subsetted to cover only the time
range of interest::
em16_tgo_sc_slt_npo_20171205_20230115_s20160414_v01.bc
8. Generic PCK:
pck00010.tpc
9. NOMAD IK:
em16_tgo_nomad_v02.ti
These 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
sce2s_c
2 getsta furnsh_c vnorm_c 1,3-6
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-6,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_c
Refer 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/naif0012.tls',
'kernels/sclk/em16_tgo_step_20160414.tsc'
)
\begintext
Solution Source Code
#include <stdio.h>
/*
Standard CSPICE User Include File
*/
#include "SpiceUsr.h"
/*
Local Parameters
*/
#define METAKR "convtm.tm"
#define SCLKID -143
#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: 2018 JUN 11 19:32:00
ET Seconds Past J2000: 582017589.185
Calendar ET (etcal_c): 2018 JUN 11 19:33:09.184
Calendar ET (timout_c): 2018-JUN-11T19:33:09
Spacecraft Clock Time: 1/0070841719.26698
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: 2458281.3146896
Earliest UTC convertible to SCLK: 2016-03-13T21:34:13.194
UTC time from spacecraft clock: 2018-06-11T19:32:00.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.
\begindata
KERNELS_TO_LOAD = (
'kernels/lsk/naif0012.tls',
'kernels/spk/de430.bsp',
'kernels/spk/mar085.bsp',
'kernels/spk/em16_tgo_mlt_20171205_20230115_v01.bsp',
'kernels/fk/em16_tgo_v07.tf'
)
\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.
Since the SPICE body/ID mapping for TGO is not
yet included in the standard library, we will
need the frame kernel where the mapping is
defined.
*/
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 Mars as seen from
ExoMars-16 TGO in the J2000 frame. All of the ephemeris
readers return states in units of kilometers and
kilometers per second.
*/
spkezr_c ( "MARS", et, "J2000", "LT+S",
"TGO", state, <ime );
printf ( " Apparent state of Mars as seen "
"from ExoMars-16 TGO 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
ExoMars-16 TGO 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",
"TGO", pos, <ime );
printf ( " Apparent position of Earth as seen "
"from ExoMars-16 TGO 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] );
/*
We need only display LTIME, as it is precisely the
light time in which we are interested.
*/
printf ( " One way light time between ExoMars-16 "
"TGO and the apparent\n" );
printf ( " position of Earth (seconds): "
"%16.3f\n", ltime );
/*
Compute the apparent position of the Sun as seen
from Mars in the J2000 frame.
*/
spkpos_c ( "SUN", et, "J2000", "LT+S",
"MARS", pos, <ime );
printf ( " Apparent position of Sun as seen "
"from Mars 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 Mars. The above SPKPOS call gives us
the apparent distance, so we need to adjust our
aberration correction appropriately.
*/
spkpos_c ( "SUN", et, "J2000", "NONE",
"MARS", 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 "
"Mars body centers:\n" );
printf ( " (AU): %16.3f\n", dist );
return(0);
}
Solution Sample Output
Converting UTC Time: 2018 JUN 11 19:32:00
ET seconds past J2000: 582017589.185
Apparent state of Mars as seen from ExoMars-16 TGO in the J2000
frame (km, km/s):
X = 2911.822
Y = -2033.802
Z = -1291.701
VX = 1.310
VY = -0.056
VZ = 3.104
Apparent position of Earth as seen from ExoMars-16 TGO in the
J2000 frame (km):
X = -49609884.080
Y = 57070665.862
Z = 30304236.930
One way light time between ExoMars-16 TGO and the apparent
position of Earth (seconds): 271.738
Apparent position of Sun as seen from Mars in the
J2000 frame (km):
X = -24712734.289
Y = 194560532.943
Z = 89906636.789
Actual distance between Sun and Mars body centers:
(AU): 1.442
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
spkezr_c ( "MARS", et, "J2000", "LT+S",
"-143", state, <ime );
http://naif.jpl.nasa.gov/pub/naif/generic_kernels/spk/satellites/
Actual position of Jupiter as seen from Mars in the
J2000 frame (km):
X = -536521483.296
Y = -384722940.462
Z = -145930841.439
Additional kernels used in this task:
a. EDM lander FK:
em16_emd_v00.tf
b. EDM landing site SPK:
em16_edm_sot_landing_site_20161020_21000101_v01.bsp
c. Generic PCK, where the Mars orientation constants are
provided:
pck00010.tpc
Apparent position of EDM Landing Site (EDM_LANDING_SITE, NAIF
ID -117900) as seen from ExoMars-16 TGO in the J2000 frame
(km):
X = -131.716
Y = -2168.989
Z = 208.792
Actual (geometric) position of Sun as seen from Mars in the
J2000 frame (km):
X = -24730875.201
Y = 194558449.560
Z = 89906170.855
Light-time corrected position of Sun as seen from Mars in the
J2000 frame (km):
X = -24730866.489
Y = 194558445.246
Z = 89906168.754
Apparent position of Sun as seen from Mars in the
J2000 frame (km):
X = -24712734.289
Y = 194560532.943
Z = 89906636.789
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/naif0012.tls',
'kernels/sclk/em16_tgo_step_20160414.tsc',
'kernels/spk/de430.bsp',
'kernels/spk/mar085.bsp',
'kernels/spk/em16_tgo_mlt_20171205_20230115_v01.bsp',
'kernels/fk/em16_tgo_v07.tf',
'kernels/ck/em16_tgo_sc_slt_npo_20171205_20230115_s20160414_v01.bc',
'kernels/pck/pck00010.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 ExoMars-16 TGO
The necessary ephemerides
A planetary constants file (PCK)
A spacecraft orientation kernel for ExoMars-16 TGO (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 Mars as seen from
ExoMars-16 TGO in the J2000 reference frame.
*/
spkezr_c ( "MARS", et, "J2000", "LT+S",
"TGO", state, <ime );
/*
Now obtain the transformation from the inertial
J2000 frame to the non-inertial body-fixed IAU_MARS
frame. Since we want the apparent position, we
need to subtract ltime from et.
*/
sxform_c ( "J2000", "IAU_MARS", et-ltime, sform );
/*
Now rotate the apparent J2000 state into IAU_MARS
with the following matrix multiplication:
*/
mxvg_c ( sform, state, 6, 6, bfixst );
/*
Display the results.
*/
printf ( " Apparent state of Mars as seen "
"from ExoMars-16 TGO in the IAU_MARS\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 ( "MARS", et, "IAU_MARS", "LT+S",
"TGO", state, <ime );
/*
Display the results.
*/
printf ( " Apparent state of Mars as seen "
"from ExoMars-16 TGO in the IAU_MARS\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_MARS 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 Mars as seen from the orbiter
and the nominal instrument view direction. First,
compute the apparent position of Mars as seen from
ExoMars-16 TGO in the J2000 frame.
*/
spkpos_c ( "MARS", et, "J2000", "LT+S",
"TGO", pos, <ime );
/*
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 -Y axis
of the TGO_SPACECRAFT frame defined there.
*/
bsight[0] = 0.0;
bsight[1] = -1.0;
bsight[2] = 0.0;
/*
Now compute the rotation matrix from TGO_SPACECRAFT
into J2000.
*/
pxform_c ( "TGO_SPACECRAFT", "J2000", et, pform );
/*
And multiply the result to obtain the nominal instrument
view direction 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 Mars and the\n" );
printf ( " ExoMars-16 TGO nominal "
"instrument view direction (degrees): \n" );
printf ( " %16.3f\n", sep );
/*
Or alternatively we can work in the spacecraft
frame directly.
*/
spkpos_c ( "MARS", et, "TGO_SPACECRAFT", "LT+S",
"TGO", pos, <ime );
/*
The nominal instrument view direction is the -Y-axis
in the TGO_SPACECRAFT frame.
*/
bsight[0] = 0.0;
bsight[1] = -1.0;
bsight[2] = 0.0;
/*
Lastly compute the angular separation.
*/
convrt_c ( vsep_c(bsight, pos), "RADIANS",
"DEGREES", &sep );
printf ( " Angular separation between the "
"apparent position of Mars and the\n" );
printf ( " ExoMars-16 TGO nominal "
"instrument view direction computed\n" );
printf ( " using vectors in the "
"TGO_SPACECRAFT frame (degrees):\n" );
printf ( " %16.3f\n", sep );
return(0);
}
Solution Sample Output
Converting UTC Time: 2018 JUN 11 19:32:00
ET seconds past J2000: 582017589.185
Apparent state of Mars as seen from ExoMars-16 TGO in the IAU_MARS
body-fixed frame (km, km/s):
X = -2843.464125
Y = 2235.459544
Z = 1095.894969
VX = 0.311443
VY = -1.151929
VZ = 3.082123
Apparent state of Mars as seen from ExoMars-16 TGO in the IAU_MARS
body-fixed frame (km, km/s) obtained using spkezr_c directly:
X = -2843.464125
Y = 2235.459544
Z = 1095.894969
VX = 0.311443
VY = -1.151929
VZ = 3.082123
Angular separation between the apparent position of Mars and the
ExoMars-16 TGO nominal instrument view direction (degrees):
5.438
Angular separation between the apparent position of Mars and the
ExoMars-16 TGO nominal instrument view direction computed
using vectors in the TGO_SPACECRAFT frame (degrees):
5.438
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
Segment No.: 1
Object: -143000
Interval Begin ET Interval End ET AV
------------------------ ------------------------ ---
2018-JUN-11 00:01:09.184 2018-JUN-12 06:28:03.102 Y
2018-JUN-12 06:58:03.102 2018-JUN-12 18:15:43.102 Y
2018-JUN-12 18:45:43.102 2018-JUN-13 04:03:23.102 Y
2018-JUN-13 04:33:23.102 2018-JUN-13 07:59:43.102 Y
2018-JUN-13 08:29:43.102 2018-JUN-13 12:01:09.184 Y
Object: -143000
Interval Begin ET Interval End ET AV
------------------------ ------------------------ ---
2018-JUN-11 00:01:09.184 2018-JUN-13 12:01:09.184 Y
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/naif0012.tls',
'kernels/spk/de430.bsp',
'kernels/spk/mar085.bsp',
'kernels/spk/em16_tgo_mlt_20171205_20230115_v01.bsp',
'kernels/fk/em16_tgo_v07.tf',
'kernels/pck/pck00010.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)
A frames kernel (TF) with the TGO ID/name mapping
*/
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 ExoMars-16 TGO
on Mars.
*/
subpnt_c ( "NEAR POINT: ELLIPSOID",
"MARS", et, "IAU_MARS", "LT+S",
"TGO", spoint, &trgepc, srfvec );
printf ( " Apparent sub-observer point of "
"ExoMars-16 TGO on Mars in the\n" );
printf ( " IAU_MARS 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 Mars
as seen from ExoMars-16 TGO.
*/
subslr_c ( "NEAR POINT: ELLIPSOID",
"MARS", et, "IAU_MARS", "LT+S",
"TGO", spoint, &trgepc, srfvec );
printf ( " Apparent sub-solar point on Mars "
"as seen from ExoMars-16 TGO in\n" );
printf ( " the IAU_MARS 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: 2018 JUN 11 19:32:00
ET seconds past J2000: 582017589.185
Apparent sub-observer point of ExoMars-16 TGO on Mars in the
IAU_MARS frame (km):
X = 2554.165
Y = -2008.010
Z = -983.240
ALT = 385.045
Apparent sub-solar point on Mars as seen from ExoMars-16 TGO in
the IAU_MARS frame (km):
X = 487.589
Y = -3348.610
Z = -286.697
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 Mars as seen from ExoMars-16 TGO in
the IAU_MARS frame using the 'Near Point: ellipsoid' method
(km):
X = 487.589
Y = -3348.610
Z = -286.697
Apparent sub-solar point on Mars as seen from ExoMars-16 TGO in
the IAU_MARS frame using the 'Intercept: ellipsoid' method
(km):
X = 487.547
Y = -3348.322
Z = -290.077
Apparent sub-spacecraft point of ExoMars-16 TGO on Phobos in
the IAU_PHOBOS frame using the 'Near Point: ellipsoid' method
(km):
X = 12.059
Y = 4.173
Z = -0.675
Planetocentric coordinates of the sub-spacecraft point on
Phobos (degrees, km):
LAT = -3.030
LON = 19.088
R = 12.779
Planetographic coordinates of the sub-spacecraft point on
Phobos (degrees, km):
LAT = -6.267
LON = 340.912
ALT = -0.202
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/naif0012.tls',
'kernels/sclk/em16_tgo_step_20160414.tsc',
'kernels/spk/de430.bsp',
'kernels/spk/mar085.bsp',
'kernels/spk/em16_tgo_mlt_20171205_20230115_v01.bsp',
'kernels/fk/em16_tgo_v07.tf',
'kernels/ck/em16_tgo_sc_slt_npo_20171205_20230115_s20160414_v01.bc',
'kernels/pck/pck00010.tpc',
'kernels/ik/em16_tgo_nomad_v02.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",
"ExoMars-16 TGO NOMAD LNO Nadir 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 lnonid;
SpiceInt marsid;
SpiceInt sc;
SpiceBoolean found;
/*
Load the kernels that this program requires. We
will need:
A leapseconds kernel.
A SCLK kernel for ExoMars-16 TGO.
Any necessary ephemerides.
The ExoMars-16 TGO frame kernel.
An ExoMars-16 TGO C-kernel.
A PCK file with Mars constants.
The ExoMars-16 TGO NOMAD 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 NOMAD LNO Nadir aperture. To do this we will need the
ID code for TGO_NOMAD_LNO_NAD.
*/
bodn2c_c ( "TGO_NOMAD_LNO_NAD", &lnonid, &found );
/*
Stop the program if the code was not found.
*/
if ( !found )
{
printf ( "Unable to locate the ID code for "
"TGO_NOMAD_LNO_NAD\n" );
exit ( EXIT_FAILURE );
}
/*
Now retrieve the field of view parameters.
*/
getfov_c ( lnonid, 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 Mars.
*/
sincpt_c ( "Ellipsoid", "MARS", et, "IAU_MARS",
"LT+S", "TGO", 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_MARS frame of the intersection.
*/
printf ( " Position vector of surface intercept "
"in the IAU_MARS 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", "MARS", et, "IAU_MARS",
"LT+S", "TGO", 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 Mars.
*/
bodn2c_c ( "MARS", &marsid, &found );
/*
The ID code for MARS 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 Mars.\n" );
exit ( EXIT_FAILURE );
}
/*
Compute local solar time corresponding to the TDB light
time corrected epoch at the intercept.
*/
et2lst_c ( trgepc,
marsid,
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: 2018 JUN 11 19:32:00
ET seconds past J2000: 582017589.185
Vector: Boundary Corner 1
Position vector of surface intercept in the IAU_MARS frame (km):
X = 2535.004
Y = -2028.528
Z = -990.594
Planetocentric coordinates of the intercept (degrees):
LAT = -16.967
LON = -38.667
Phase angle (degrees): 48.207
Solar incidence angle (degrees): 43.872
Emission angle (degrees): 4.798
Vector: Boundary Corner 2
Position vector of surface intercept in the IAU_MARS frame (km):
X = 2525.056
Y = -2042.075
Z = -988.196
Planetocentric coordinates of the intercept (degrees):
LAT = -16.925
LON = -38.963
Phase angle (degrees): 50.707
Solar incidence angle (degrees): 43.586
Emission angle (degrees): 7.432
Vector: Boundary Corner 3
Position vector of surface intercept in the IAU_MARS frame (km):
X = 2525.201
Y = -2042.104
Z = -987.770
Planetocentric coordinates of the intercept (degrees):
LAT = -16.917
LON = -38.962
Phase angle (degrees): 50.708
Solar incidence angle (degrees): 43.585
Emission angle (degrees): 7.413
Vector: Boundary Corner 4
Position vector of surface intercept in the IAU_MARS frame (km):
X = 2535.149
Y = -2028.558
Z = -990.170
Planetocentric coordinates of the intercept (degrees):
LAT = -16.960
LON = -38.666
Phase angle (degrees): 48.208
Solar incidence angle (degrees): 43.871
Emission angle (degrees): 4.769
Vector: ExoMars-16 TGO NOMAD LNO Nadir Boresight
Position vector of surface intercept in the IAU_MARS frame (km):
X = 2530.122
Y = -2035.307
Z = -989.188
Planetocentric coordinates of the intercept (degrees):
LAT = -16.942
LON = -38.814
Phase angle (degrees): 49.457
Solar incidence angle (degrees): 43.729
Emission angle (degrees): 6.086
Local Solar Time at boresight intercept (24 Hour Clock):
14:51:36
Extra Credit
|