// 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 [pso2] = CL_kp_anomConvertCir(type_anom1,type_anom2,ex,ey,pso1) // Conversion of anomaly (circular or equinoctical elements) // // Calling Sequence // [pso2] = CL_kp_anomConvertCir(type_anom1, type_anom2, ex, ey, pso1) // // Description // //

Converts anomaly from one type to another, for circular or equinoctial orbital elements

//

//

The components of the eccentricity vector ex, ey are of the form:

//

ex = e*cos(A)

//

ey = e*sin(A)

//

Where:

//

A = w (argument of perigee) for circular orbital elements, or

//

A = W + w (RAAN + argument of perigee) for equinoctial orbital elements.

//

//

The "pso" quantity that can be converted is:

//

- pso = A + M if type_anom = "M"

//

- pso = A + v if type_anom = "v"

//

- pso = A + E if type_anom = "E"

//

//

Note:

//

- Only elliptical orbits are handled (eccentricity < 1)

//

//
// // Parameters // type_anom1: (string) Type of input anomaly (1x1) // type_anom2: (string) Type of output anomaly (1x1) // ex: x component of eccentricity vector (1xN or 1x1) // ey: y component of eccentricity vector (1xN or 1x1) // pso1: Input anomaly (1xN or 1x1) // pso2: Output anomaly (1xN) // // Authors // CNES - DCT/SB // // See Also // CL_kp_anomConvert // // Examples // // Mean argument of latitude to true argument of latitude: // ex = 0.1; // ey = 0.2; // psoM = [-0.1, 10.2]; // rad // psov = CL_kp_anomConvertCir("M", "v", ex, ey, psoM) // // // Consistency test: // psoM - CL_kp_anomConvertCir("v", "M", ex, ey, psov) // => 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 value for type_anom1"); end if (type_anom2 <> "M" & type_anom2 <> "v" & type_anom2 <> "E") CL__error("Invalid value for type_anom2"); end // check size / resize [ex, ey, pso1] = CL__checkInputs(ex, 1, ey, 1, pso1, 1); // Check validity of eccentricity // only ecc in [0,1[ allowed (otherwise: error) ecc = sqrt(ex.^2 + ey.^2); if (find(ecc >= 1) <> []) CL__error("Invalid eccentricity") end if (type_anom2 == type_anom1) pso2 = pso1; else // atan(ey,ex) = pom (circular elements) or pom + gom (equinoctical elements) // NB: A is 0 if both ex and ey are 0. A = atan(ey, ex); // Using CL_kp_anomConvert (usual Kepler's equation) is sufficient. No need for specific function. pso2 = CL_kp_anomConvert(type_anom1, type_anom2, ecc, pso1 - A) + A; end endfunction