// 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 [dv,anv, dv_car] = CL_man_inclination(ai,ei,inci,pomi,incf, posman,icheck,mu) // Inclination maneuver (elliptical orbits) - DEPRECATED // // Calling Sequence // [dv,anv] = CL_man_inclination(ai,ei,inci,pomi,incf [,posman,icheck,mu]) // // Description // //

This function is deprecated.

//

Replacement function: CL_man_dvInc

//

// //

Computes the velocity increment needed to change the inclination.

//

The maneuver consists in making the orbit rotate around the line of nodes, thus changing the velocity direction at one of the nodes.

//

The output argument dv is the velocity increment vector in spherical coordinates in the "qsw" local orbital frame.

//

anv is the true anomaly at the maneuver position.

//

The optional flag posman can be used to define the position of the maneuver (maneuver at ascending node, descending node, or "best choice").

//

The optional flag icheck is used to enforce no checking // on the final inclination value. If the targetted inclination is less // than 0 or more than pi, the ascending node will rotate 180 degrees. // If icheck is true, the right ascension of the ascending node is not changed by the maneuver.

//
//
// // Parameters // ai : Semi major-axis [m] (1xN or 1x1) // ei: Eccentricity (1xN or 1x1) // inci : Initial inclination [rad] (1xN or 1x1) // pomi: Argument of periapsis [rad] (1xN or 1x1) // incf : Final inclination [rad] (1xN or 1x1) // posman: (optional) Flag specifying the position of the maneuver: 1 or "an" -> ascending node, -1 or "dn" -> descending node, 0 or "opt" -> minimal |dv|. Default is at the ascending node. (1xN or 1x1) // icheck: (optional, boolean) Flag specifying if incf must be checked in the standard range for inclination values ([0, pi]). Default is %t. (1x1) // mu : (optional) Gravitational constant. [m^3/s^2] (default value is %CL_mu) // dv : Velocity increment vector in spherical coordinates in the "qsw" frame [lambda;phi;|dv|] [rad,rad,m/s] (3xN) // anv: True anomaly at maneuver position [rad] (1xN) // // Authors // CNES - DCT/SB // // See also // CL_man_biElliptic // CL_man_hohmann // CL_man_sma // CL_man_hohmannG // CL_man_apsidesLine // // Examples // ai = 7200.e3; // ei = 0.1; // inci = 1; // pomi = 1.4; // incf = 1.1; // [dv,anv] = CL_man_inclination(ai,ei,inci,pomi,incf,posman=0) // // Check results : // anm = CL_kp_v2M(ei,anv); // kep = [ai ; ei ; inci ; pomi ; 0 ; anm]; // kep1 = CL_man_applyDv(kep,dv) // Declarations: // Code: CL__warnDeprecated(); // deprecated function if (~exists("posman", "local")); posman = 1; end if (~exists("icheck", "local")); icheck = %t; end if (~exists("mu", "local")); mu = CL__dataGetEnv("mu"); end // convert posman type to "real" if (typeof(posman) == "string") str = posman; posman = %nan * ones(str); posman(find(str == "an")) = 1; posman(find(str == "dn")) = -1; posman(find(str == "opt")) = 0; end // checks arguments sizes are OK / resizes [ai,ei,inci,pomi,incf,posman] = CL__checkInputs(ai,1,ei,1,inci,1,pomi,1,incf,1,posman,1); if (find(ai <= 0 | ei < 0 | ei >= 1 | inci < 0 | inci > %pi) <> []) CL__error("Invalid input arguments (orbital elements)"); end if (icheck & find(incf < 0 | incf > %pi ) <> []) CL__error("Invalid input arguments (final inclination)"); end if (find(posman <> 0 & posman <> 1 & posman <> -1) <> []) CL__error("Invalid value for ''posman''"); end // 1) maneuver anomaly // maneuver at ascending node (default) anv = -pomi; // posman = -1 => descending node I = find(posman == -1); anv(I) = -pomi(I) + %pi * ones(I); // posman = 0 => ascending node, except if ascending node // closer to periapsis than descending node. // Note: "posman == 0" tested twice to simplify code if (find(posman == 0) <> []) // I: such that perigee "close" to ascending node I = find(posman == 0 & abs(CL_rMod(pomi, -%pi, %pi)) < %pi/2); anv(I) = -pomi(I) + %pi * ones(I); // descending node end // final result in [0, 2*pi] anv = CL_rMod(anv, 0, 2*%pi); // 2) velocity increment // vcosslope = |v| * cos(slope) // slope = %pi/2 - angle(velocity vector, radius vector) vcosslope = (1+ei.*cos(anv)) .* sqrt(mu./(ai.*(1-ei.^2))); // change of inclination // NB: "modulo" not necessary (results would be the same) dinc = incf - inci; // sign of rotation around radius vector: // = 1 at ascending node (dinc>0 <=> rotation angle>0) // = -1 at descending node (dinc>0 <=> rotation angle<0) // NB: assumes anv is in [0, 2*pi] sgnrot = -1 * ones(anv); I = find(abs(CL_rMod(anv+pomi,-%pi,%pi)) < %pi/2); // asc node sgnrot(I) = 1; // cartesian coordinates of dv in "qsw": // obtained by: dv = v2 - v1, with: // v1 = v * [sin(slope); cos(slope); 0] // v2 = image of v1 by rotation around radius vector, // of angle: dinc * sgnrot dv_car = [ zeros(dinc); -2 * vcosslope .* (sin(dinc/2)).^2; vcosslope .* sin(dinc) .* sgnrot ]; dv = CL_co_car2sph(dv_car); endfunction