// 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 kepculated by CEA, CNRS and INRIA at the following URL
// 'http://www.cecill.info'.
function [pos,vel,jacob] = CL_oe_kep2car(kep, mu)
// Classical Keplerian to cartesian orbital elements
//
// Calling Sequence
// [pos, vel, jacob] = CL_oe_kep2car(kep [, mu])
//
// Description
//
//
Converts classical Keplerian orbital elements to cartesian orbital elements.
//
The transformation jacobian is optionally computed.
//
See Orbital elements for more details.
//
//
//
// Parameters
// kep: Classical Keplerian orbital elements [sma;ecc;inc;pom;gom;anm] [m,rad] (6xN)
// mu : (optional) Gravitational constant. [m^3/s^2] (default value is %CL_mu)
// pos: Position vector [x;y;z] [m] (3xN)
// vel: Velocity vector [vx;vy;vz] [m/s] (3xN)
// jacob: (optional) Transformation jacobian (See Orbital elements for more details) (6x6xN)
//
// Authors
// CNES - DCT/SB
//
// See also
// CL_oe_car2kep
//
// Examples
// // // Example 1
// kep = [7000.e3; 0.1; 0.5; 1; 2; 3];
// [pos, vel] = CL_oe_kep2car(kep);
//
// // Example 2
// kep = [7000.e3; 0.1; 0.5; 1; 2; 3];
// [pos, vel, jacob1] = CL_oe_kep2car(kep);
// [kep2, jacob2] = CL_oe_car2kep(pos, vel);
// kep2 - kep // => zero
// jacob2 * jacob1 // => identity
// Declarations:
EPS_ORB = CL__dataGetEnv("epsOrb", internal=%t);
// Code:
if (~exists("mu", "local")); mu = CL__dataGetEnv("mu"); end
// Handle [] cases
if (kep == [])
pos = [];
vel = [];
jacob = [];
return;
end
// Check validity of input
[isvalid,orbit_type] = CL__oe_isValid("kep",kep);
if (~isvalid); CL__error("Invalid input"); end;
if (find(orbit_type == 3) <> []);
CL__error("Invalid orbital element (parabolic orbit)");
end
// Flag to compute jacobian
cjac = %f;
if (argn(1) == 3)
cjac = %t;
end
N = size(kep,2);
pos = zeros(3,N);
vel = zeros(3,N);
jacob = [];
if (cjac)
jacob = zeros(6,6,N);
end
// Elliptic orbits
I = find(orbit_type == 1);
if (I <> [])
[pos(:,I), vel(:,I), jacob_tmp] = CL__oe_kep2car_ell(kep(:,I),mu,cjac);
if (cjac)
jacob(:,:,I) = jacob_tmp;
end
end
// Hyperbolic orbits
I = find(orbit_type == 2);
if (I <> [])
[pos(:,I), vel(:,I), jacob_tmp] = CL__oe_kep2car_hyp(kep(:,I),mu,cjac);
if (cjac)
jacob(:,:,I) = jacob_tmp;
end
end
// Parabolic orbits
// I = find(orbit_type == 3);
// if (I <> [])
// [pos(:,I), vel(:,I), jacob_tmp] = CL__oe_kep2car_parab(kep(:,I),mu,cjac);
// if (cjac)
// jacob(:,:,I) = jacob_tmp;
// end
// end
endfunction