// Copyright (c) CNES 2008 // // This software is part of CelestLab, a CNES toolbox for Scilab // // This software is governed by the CeCILL license under French law and // abiding by the rules of distribution of free software. You can use, // modify and/ or redistribute the software under the terms of the CeCILL // license as circulated by CEA, CNRS and INRIA at the following URL // 'http://www.cecill.info'. function [vect_inertial] = CL_fr_tnw2inertial(pos_car,vel_car,vect_tnw) // "tnw" Local orbital frame to inertial frame vector transformation // // Calling Sequence // vect_inertial = CL_fr_tnw2inertial(pos_car,vel_car,vect_tnw) // // Description // //

Converts the coordinates of a vector from the "tnw" local orbital frame to the inertial frame.

//

The inertial frame is implicitly the frame relative to which the satellite position and velocity are defined. The "qsw" local frame is built using these position and velocity vectors.

//

// //

See Local frames for more details on the definition of local frames.

//
//
// // Parameters // pos_car: satellite's position vector relative to the inertial frame [m] (3xN) // vel_car: satellite's velocity vector relative to the inertial frame [m/s] (3xN) // vect_tnw: cartesian coordinates in the "tnw" local orbital frame (3xN) // vect_inertial: cartesian coordinates in the inertial frame (3xN) // // Bibliography // 1) Mecanique spatiale, CNES - Cepadues 1995, Tome I, section 10.2.2.3 (Definition du Repere orbital local) // 2) CNES - MSLIB FORTRAN 90, Volume O (mo_tnw_geo) // // Authors // CNES - DCT/SB // // See also // CL_fr_tnwMat // CL_fr_inertial2tnw // CL_fr_qswMat // CL_fr_inertial2qsw // CL_fr_qsw2inertial // // Examples // // Example 1 // pos_car = [[3500.e3; 2500.e3; 5800.e3], [4500.e3; 2100.e3; 6800.e3]]; // vel_car = [[1.e3; 3.e3; 7.e3], [2.e3; 3.e3; 6.e3]]; // vect_tnw = [[1;0;0], [1;1;1]]; // vect_inertial = CL_fr_tnw2inertial(pos_car, vel_car, vect_tnw) // Declarations: // Code: [lhs,rhs] = argn(); if (rhs <> 3) CL__error("Invalid number of input arguments"); end //vector t, tangent to velocity vector t = CL_unitVector(vel_car) //CL_cross product pos_car^vel w0 = CL_cross(pos_car,vel_car) //vector w, giving direction of angular momentum vector (perpendicular to the osculating orbit plane) w = CL_unitVector(w0) //CL_cross product: n = w^t n = CL_cross(w,t) vect_inertial = [ vect_tnw(1,:).*t(1,:) + vect_tnw(2,:).*n(1,:) + vect_tnw(3,:).*w(1,:) ; vect_tnw(1,:).*t(2,:) + vect_tnw(2,:).*n(2,:) + vect_tnw(3,:).*w(2,:) ; vect_tnw(1,:).*t(3,:) + vect_tnw(2,:).*n(3,:) + vect_tnw(3,:).*w(3,:)]; endfunction