Index of Functions: A  B  C  D  E  F  G  H  I  J  K  L  M  N  O  P  Q  R  S  T  U  V  W  X 
Index Page
cspice_pxfrm2

Table of contents
Abstract
I/O
Parameters
Examples
Particulars
Exceptions
Files
Restrictions
Required_Reading
Literature_References
Author_and_Institution
Version
Index_Entries

Abstract


   CSPICE_PXFRM2 returns the 3x3 matrix that transforms position
   vectors from one specified frame at a specified epoch to another
   specified frame at another specified epoch.

I/O


   Given:

      from     name of a reference frame recognized by Mice that
               corresponds to the input `etfrom'.

               [1,c1] = size(from); char = class(from)

                  or

               [1,1] = size(from); cell = class(from)

      to       name of a reference frame recognized by Mice that
               corresponds to the desired output at `etto''.

               [1,c2] = size(to); char = class(to)

                  or

               [1,1] = size(to); cell = class(to)

      etfrom   epoch(s) in ephemeris seconds past the epoch of J2000 (TDB)
               corresponding to the `from' reference frame.

               [1,n] = size(etfrom); double = class(etfrom)

      etto     epoch(s) in ephemeris seconds past the epoch of J2000 (TDB)
               that corresponds to the `to' reference frame.

               [1,n] = size(etto); double = class(etto)

   the call:

      [rotate] = cspice_pxfrm2( from, to, etfrom, etto )

   returns:

      rotate   operator(s) that transform position vector(s) from the
               reference frame `from' at epoch `etfrom' to the frame `to' at
               epoch `etto'.

               If [1,1] = size(etfrom) then [3,3]   = size(rotate)
               If [1,n] = size(etfrom) then [3,3,n] = size(rotate)
                                             double = class(rotate)

               If (x, y, z) is a position relative to the reference
               frame `from' at time `etfrom' then the vector ( x', y',
               z') is the same position relative to the frame `to' at
               epoch `etto'. Here the vector ( x', y', z' ) is defined
               by the equation:

                    -   -       -        -     -  -
                   | x'  |     |          |   | x  |
                   | y'  |  =  |  rotate  |   | y  |
                   | z'  |     |          |   | z  |
                    -   -       -        -     -  -

              `rotate' returns with the same vectorization measure (N)
               as `etfrom' and `etto'.

Parameters


   None.

Examples


   Any numerical results shown for this example may differ between
   platforms as the results depend on the SPICE kernels used as input
   and the machine specific arithmetic implementation.

   1) Suppose that MGS has taken a picture of Mars at time `etrec' with
      the MOC narrow angle camera. We want to know the latitude and
      longitude associated with two pixels projected to Mars'
      surface: the boresight and one along the boundary of the
      field of view (FOV). Due to light time, the photons taken in
      the picture left Mars at time `etemit', when Mars was at a
      different state than at time `etrec'.

      In order to solve this problem, we could use the cspice_sincpt
      routine for both pixels, but this would be slow. Instead, we
      will assume that the light time for each pixel is the same. We
      will call cspice_sincpt once to get the light time and surface point
      associated with the boresight. Then, we will rotate the first
      FOV boundary vector from the camera frame at `etrec' to the
      body-fixed Mars frame at `etemit', and call the faster routine
      cspice_surfpt to retrieve the surface point for the FOV boundary
      vector.

      This example problem could be extended to find the latitude
      and longitude associated with every pixel in an instrument's
      field of view, but this example is simplified to only solve
      for two pixels: the boresight and one along the boundary of
      the field of view.

      Assumptions:

         1)  The light times from the surface points in the camera's
             field of view to the camera are equal.

         2)  The camera offset from the center of gravity of the
             spacecraft is zero. If the data are more accurate
             and precise, this assumption can be easily discarded.

         3)  An ellipsoid shape model for the target body is
             sufficient.

         4)  The boundary field of view vector returned from cspice_getfov
             is associated with a boundary field of view pixel. If
             this example were extended to include a geometric camera
             model, this assumption would not be needed since the
             direction vectors associated with each pixel would be
             calculated from the geometric camera model.

      Use the meta-kernel shown below to load the required SPICE
      kernels.


        KPL/MK

         File name: pxfrm2_ex1.tm

         This is the meta-kernel file for the example problem for
         the subroutine cspice_pxfrm2. These kernel files can be found in
         the NAIF archives.

         In order for an application to use this meta-kernel, the
         kernels referenced here must be present in the user's
         current working directory.

         The names and contents of the kernels referenced
         by this meta-kernel are as follows:

            File name                     Contents
            ---------                     --------
            de421.bsp                     Planetary ephemeris
            pck00009.tpc                  Planet orientation and
                                          radii
            naif0009.tls                  Leapseconds
            mgs_ext12_ipng_mgs95j.bsp     MGS ephemeris
            mgs_moc_v20.ti                MGS MOC instrument
                                          parameters
            mgs_sclkscet_00061.tsc        MGS SCLK coefficients
            mgs_sc_ext12.bc               MGS s/c bus attitude

          \begindata

          KERNELS_TO_LOAD = ( 'de421.bsp',
                              'pck00009.tpc',
                              'naif0009.tls',
                              'mgs_ext12_ipng_mgs95j.bsp',
                              'mgs_moc_v20.ti',
                              'mgs_sclkscet_00061.tsc',
                              'mgs_sc_ext12.bc' )

          \begintext

          End of meta-kernel.


      Example code begins here.


      function pxfrm2_ex1()

         % --------------------- Program Setup ---------------------

         metakr = 'pxfrm2_ex1.tm';
         camera = 'MGS_MOC_NA';
         NCORNR = 4;
         ABCORR = 'CN+S';

         %
         % Load kernels
         %
         cspice_furnsh( metakr );

         %
         % Convert the time the picture was taken from a
         % UTC time string to seconds past J2000, TDB.
         %
         etrec = cspice_str2et( '2003 OCT 13 06:00:00 UTC' );

         %
         % Assume the one-way light times from different
         % surface points on Mars to MGS within the camera's
         % FOV are equal. This means the photons that make
         % up different pixels were all emitted from Mars at
         % `etemit' and received by MGS at `etrec'.  It would be
         % slow to process images using cspice_sincpt for every
         % pixel.  Instead, we will use cspice_sincpt on the
         % boresight pixel and use cspice_surfpt for the first FOV
         % boundary pixel.  If this example program were extended
         % to include all of the camera's pixels, cspice_surfpt would
         % be used for the remaining pixels.
         %
         % Get the MGS MOC Narrow angle camera (MGS_MOC_NA)
         % ID code. Then look up the field of view (FOV)
         % parameters by calling cspice_getfov.
         %
         [camid, found] = cspice_bodn2c( camera );

         if ( ~found )
             txt = sprintf( ['SPICE(NOTRANSLATION) Could not find ', ...
                             'ID code for instrument %s.'], camera );
             error( txt )
         end

         %
         % cspice_getfov will return the name of the camera-fixed frame
         % in the string OBSREF, the camera boresight vector in
         % the array BSIGHT, and the FOV corner vectors in the
         % array BOUNDS.
         %
         [shape, obsref, bsight, bounds] = cspice_getfov( camid, NCORNR);

         fprintf( '\nObservation Reference Frame:  %s\n', obsref );

         %
         % ----------- Boresight Surface Intercept -----------
         %
         % Retrieve the time, surface intercept point, and vector
         % from MGS to the boresight surface intercept point
         % in IAU_MARS coordinates.
         %
         [ spoint, etemit, srfvec, found ] = ...
                 cspice_sincpt( 'Ellipsoid', 'Mars', etrec,  'IAU_MARS', ...
                                 ABCORR,     'MGS' , obsref,  bsight );

         if ( ~found )
             txt = sprintf( ['SPICE(NOINTERCEPT)' ...
                             'Intercept not found for boresight vector.' ]);
             error( txt )
         end

         %
         % Convert the intersection point of the boresight
         % vector and Mars from rectangular into latitudinal
         % coordinates. Convert radians to degrees.
         %
         [ radius, lon, lat ] = cspice_reclat( spoint );

         lon = lon * cspice_dpr;
         lat = lat * cspice_dpr;

         fprintf( ['\n'                                         ...
                   'Boresight surface intercept coordinates:\n' ...
                   '    Radius    (km) :  %f\n'                 ...
                   '    Latitude  (deg):  %f\n'                 ...
                   '    Longitude (deg):  %f\n' ],              ...
                    radius, lat, lon );

         %
         % ------ 1st Boundary FOV Surface Intercept (cspice_surfpt) -----
         %
         % Now we will transform the first FOV corner vector into the
         % IAU_MARS frame so the surface intercept point can be
         % calculated using cspice_surfpt, which is faster than
         % cspice_subpnt.
         %
         % If this example program were extended to include all
         % of the pixels in the camera's FOV, a few steps, such as
         % finding the rotation matrix from the camera frame to the
         % IAU_MARS frame, looking up the semi-axis values for Mars,
         % and finding the position of MGS with respect to Mars could
         % be done once and used for every pixel.
         %
         % Find the rotation matrix from the ray's reference
         % frame at the time the photons were received (etrec)
         % to IAU_MARS at the time the photons were emitted
         % (etemit).
         %
         [rotate] = cspice_pxfrm2( obsref, 'IAU_MARS', etrec, etemit );

         %
         % Look up the semi-axis values for Mars.
         %
         radii = cspice_bodvrd( 'MARS', 'RADII', 3 );

         %
         % Find the position of the center of Mars with respect
         % to MGS.  The position of the observer with respect
         % to Mars is required for the call to cspice_surfpt.  Note:
         % the apparent position of MGS with respect to Mars is
         % not the same as the negative of Mars with respect to MGS.
         %
         pos_mgs_wrt_mars = spoint - srfvec;

         %
         % The first boundary FOV pixel must be rotated into the
         % IAU_MARS reference frame.
         %
         bndvec = rotate * bounds(:,1);

         %
         % Calculate the surface point of the boundary FOV
         % vector.
         %
         [surface_point, found] = cspice_surfpt ( pos_mgs_wrt_mars, ...
                                                  bndvec, radii(1), ...
                                                  radii(2), radii(3) );

         if ( ~found )
             txt = 'SPICE(NOTFOUND)Could not calculate surface point.';
             error( txt )
         end

         surf_point_using_surfpt = surface_point;

         %
         % Convert the intersection point of the boundary
         % FOV vector and Mars from rectangular into
         % latitudinal coordinates. Convert radians
         % to degrees.
         %
         [ radius, lon, lat ] = cspice_reclat( surface_point );

         lon = lon * cspice_dpr;
         lat = lat * cspice_dpr;

         fprintf( ['\n'                                        ...
                   'Boundary vector surface intercept\n'       ...
                   'coordinates using cspice_surfpt:\n'        ...
                   '    Radius    (km) :  %f\n'                ...
                   '    Latitude  (deg):  %f\n'                ...
                   '    Longitude (deg):  %f\n'                ...
                   '    Emit time using boresight LT (s):  %10.9f\n'], ...
                    radius, lat, lon, etemit );

         %
         % ------ 1st Boundary FOV Surface Intercept Verification ----
         %
         % For verification only, we will calculate the surface
         % intercept coordinates for the first boundary vector
         % using cspice_sincpt and compare to the faster
         % cspice_surfpt method.
         %
         [ surface_point, etemit, srfvec, found ] = ...
                  cspice_sincpt( 'Ellipsoid', 'Mars', etrec,  'IAU_MARS', ...
                                  ABCORR,     'MGS' , obsref,  bounds(:,1) );

         if ( ~found )
             txt = sprintf( ['SPICE(NOTFOUND)' ...
                             'Intercept not found for ' ...
                             'the boundary FOV vector.' ]);
             error( txt )
         end

         %
         % Convert the intersection point of the first boundary
         % vector and Mars from rectangular into latitudinal
         % coordinates. Convert radians to degrees.
         %
         [ radius, lon, lat ] = cspice_reclat( surface_point );

         lon = lon * cspice_dpr;
         lat = lat * cspice_dpr;

         fprintf( ['\n'                                              ...
                   'Boundary vector surface intercept\n'             ...
                   'coordinates using cspice_sincpt:\n'              ...
                   '    Radius    (km) :  %f\n'                      ...
                   '    Latitude  (deg):  %f\n'                      ...
                   '    Longitude (deg):  %f\n'                      ...
                   '    Emit time using boresight LT (s):  %10.9f\n\n'], ...
                    radius, lat, lon, etemit );

         distance = cspice_vdist ( surf_point_using_surfpt, surface_point );

         fprintf( [ 'Distance between surface points of the first\n' ...
                    'boundary vector using cspice_surfpt and\n'      ...
                    'cspice_sincpt:\n'                               ...
                    '    Distance (mm):   %f\n' ],                   ...
                    distance*10^6 );

         %
         % It's always good form to unload kernels after use,
         % particularly in MATLAB due to data persistence.
         %
         cspice_kclear


      When this program was executed on a Mac/Intel/Octave5.x/64-bit
      platform, the output was:


      Observation Reference Frame:  MGS_MOC_NA

      Boresight surface intercept coordinates:
          Radius    (km) :  3384.940410
          Latitude  (deg):  -48.479580
          Longitude (deg):  -123.436454

      Boundary vector surface intercept
      coordinates using cspice_surfpt:
          Radius    (km) :  3384.941136
          Latitude  (deg):  -48.477482
          Longitude (deg):  -123.474080
          Emit time using boresight LT (s):  119296864.181059480

      Boundary vector surface intercept
      coordinates using cspice_sincpt:
          Radius    (km) :  3384.941136
          Latitude  (deg):  -48.477482
          Longitude (deg):  -123.474079
          Emit time using boresight LT (s):  119296864.181059465

      Distance between surface points of the first
      boundary vector using cspice_surfpt and
      cspice_sincpt:
          Distance (mm):   32.139880


Particulars


   cspice_pxfrm2 is most commonly used to transform a position between
   time-dependent reference frames.

   For more examples of where to use cspice_pxfrm2, please see:

         cspice_sincpt
         cspice_surfpt
         cspice_subslr
         cspice_ilumin

Exceptions


   1)  If sufficient information has not been supplied via loaded
       SPICE kernels to compute the transformation between the
       two frames, an error is signaled by a routine
       in the call tree of this routine.

   2)  If either frame `from' or `to' is not recognized, the error
       SPICE(UNKNOWNFRAME) is signaled by a routine in the call tree
       of this routine.

   3)  If any of the input arguments, `from', `to', `etfrom' or
       `etto', is undefined, an error is signaled by the Matlab error
       handling system.

   4)  If any of the input arguments, `from', `to', `etfrom' or
       `etto', is not of the expected type, or it does not have the
       expected dimensions and size, an error is signaled by the Mice
       interface.

   5)  If the input vectorizable arguments `etfrom' and `etto' do not
       have the same measure of vectorization (N), an error is
       signaled by the Mice interface.

Files


   Appropriate kernels must be loaded by the calling program before
   this routine is called. Kernels that may be required include
   SPK files, PCK files, frame kernels, C-kernels, and SCLK kernels.

   Such kernel data are normally loaded once per program
   run, NOT every time this routine is called.

Restrictions


   None.

Required_Reading


   FRAMES
   MICE.REQ

Literature_References


   None.

Author_and_Institution


   J. Diaz del Rio     (ODC Space)
   S.C. Krening        (JPL)
   E.D. Wright         (JPL)

Version


   -Mice Version 1.1.0, 26-NOV-2021 (EDW) (JDR)

       Added -Parameters, -Exceptions, -Files, -Restrictions,
       -Literature_References and -Author_and_Institution sections, and
       completed -Particulars section. Edited the header to comply with NAIF
       standard.

       Added code example's meta-kernel.

       Eliminated use of "lasterror" in rethrow.

       Removed reference to the function's corresponding CSPICE header from
       -Required_Reading section.

   -Mice Version 1.0.0, 12-MAR-2012 (SCK)

Index_Entries


   Position transformation matrix for different epochs


Fri Dec 31 18:44:26 2021