// 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 [anom2] = CL_kp_anomConvert(type_anom1,type_anom2,ecc,anom1) // Conversion of anomaly // // Calling Sequence // [anom2] = CL_kp_anomConvert(type_anom1, type_anom2, ecc, anom1) // // Description // //

Converts anomaly from one type to another.

//

Available types are:

//

- "M" : mean anomaly

//

- "v" : true anomaly

//

- "E" : eccentric anomaly

//

//

Elliptic and hyperbolic orbits are supported.

//

See Orbital elements for more details // on orbital elements.

//

//
// // Parameters // type_anom1: (string) Type of input anomaly ("M", "E", "v") (1x1) // type_anom2: (string) Type of output anomaly ("M", "E", "v") (1x1) // ecc : Eccentricity (1xN or 1x1) // anom1: Input anomaly of type type_anom1 (1xN or 1x1) // anom2: Output anomaly of type type_anom2 (1xN) // // Authors // CNES - DCT/SB // // See Also // CL_kp_M2E // CL_kp_E2M // CL_kp_M2v // CL_kp_v2M // // Examples // // Mean anomaly to true anomaly : // ecc = 0.3; // M = [-0.1, 10.2]; // rad // v = CL_kp_anomConvert("M", "v", ecc, M) // // // Consistency test: // M - CL_kp_anomConvert("v", "M", ecc, v) // => 0 // Argument checking if (typeof(type_anom1) <> "string" | typeof(type_anom2) <> "string") CL__error("Invalid type for type_anom1 or type_anom2"); end if (type_anom1 <> "M" & type_anom1 <> "v" & type_anom1 <> "E" ) CL__error("Invalid anomaly type (type_anom1)"); end if (type_anom2 <> "M" & type_anom2 <> "v" & type_anom2 <> "E" ) CL__error("Invalid anomaly type (type_anom2)"); end // Note: // Checks would not be necessary if type_anom1 <> type_anom2 // (as already done in lower level functions) if (find (ecc < 0) <> []) CL__error("Invalid eccentricity"); end // check size / resize [ecc, anom1] = CL__checkInputs(ecc, 1, anom1, 1); // parabolic orbits (ecc == 1) => result is %nan I = find(ecc == 1); ecc(I) = %nan; anom1(I) = %nan; // necessary if type_anom2 == type_anom1 // Conversion if (type_anom2 == type_anom1) // no change of type anom2 = anom1; elseif (type_anom1 == "M") // Mean anomaly if (type_anom2 == "v") // From mean to true anomaly anom2 = CL_kp_M2v(ecc, anom1); elseif (type_anom2 == "E") // From mean to eccentric anomaly anom2 = CL_kp_M2E(ecc, anom1); end elseif (type_anom1 == "v") // True anomaly if (type_anom2 == "M") // From true to mean anomaly anom2 = CL_kp_v2M(ecc, anom1); elseif (type_anom2 == "E") // From true to eccentric anomaly anom2 = CL_kp_v2E(ecc, anom1); end elseif (type_anom1 == "E") // Eccentric anomaly if (type_anom2 == "M") // From eccentric to mean anomaly anom2 = CL_kp_E2M(ecc, anom1); elseif (type_anom2 == "v") // From eccentric to true anomaly anom2 = CL_kp_E2v(ecc, anom1); end end endfunction