// 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 [varargout] = CL_man_dvHohmann(ai,af,mu,res) // Hohmann transfer // // Calling Sequence // [deltav,dv1,dv2,anv1,anv2] = CL_man_dvHohmann(ai,af [,mu]) // man = CL_man_dvHohmann(ai,af [,mu], res="s") // // Description // //

Computes the maneuvers of a Hohmann transfer from an initial // circular orbit with semi major axis ai to a final circular // orbit with semi major axis af.

//

The output argument deltav is the sum of the norms // of the velocity increments (|dv1| + |dv2|).

//

Velocity increments are expressed in the "qsw" frame.

//

If the argument res is present and is equal to "s", all the output data are returned in a structure.

//
//
// // Parameters // ai : Semi-major axis of initial circular orbit. [m] (1xN or 1x1) // af : Semi-major axis of final circular orbit. [m] (1xN or 1x1) // mu : (optional) Gravitational constant. [m^3/s^2] (default value is %CL_mu) // res : (string, optional) Type of output: "d" or "s" for . Default is "d". // deltav : Som of norms of velocity increments [m/s] (1xN) // dv1: First velocity increment in cartesian coordnates in the "qsw" local orbital frame [m/s]. (3xN) // dv2: Second velocity increment in cartesian coordnates in the "qsw" local orbital frame. [m/s] (3xN) // anv1: True anomaly at the location of the first velocity increment (in the initial orbit): as the initial orbit is circular, anv1 is arbitrarily set to 0. (1xN) // anv2: True anomaly at the location of the second velocity increment (in the intermediate orbit). [rad] (1xN) // man: Structure containing all the output data. // // Authors // CNES - DCT/SB // // See also // CL_man_dvHohmannG // CL_man_dvBiElliptic // CL_man_dvSma // // Examples // // 7200km to 7000km : // ai = 7200.e3; // af = 7000.e3; // [deltav,dv1,dv2,anv1,anv2] = CL_man_dvHohmann(ai,af); // // // Check results: // kep = [ai ; 0 ; %pi/2 ; 0 ; 0 ; anv1]; // kep1 = CL_man_applyDvKep(kep,dv1); // kep1(6) = anv2; // kep2 = CL_man_applyDvKep(kep1,dv2) // Declarations: // Code: if (~exists("mu", "local")); mu = CL__dataGetEnv("mu"); end if (~exists("res", "local")); res = "d"; end // check "res" argument if (res <> "d" & res <> "s"); CL__error("Invalid value for argument ''res''"); end if (argn(1) > 1 & res == "s"); CL__error("Invalid number of output arguments"); end // checks arguments sizes are OK / resizes [ai,af] = CL__checkInputs(ai,1,af,1); if (find(ai <= 0 | af <= 0) <> []) CL__error("Invalid input arguments"); end // 1st man : // initial orbit: ai / ai // maneuver: ai -> af (opposite man position) // // 2nd man : // initial orbit: af / ai (maneuver at af) // maneuver: ai -> af (opposite man position) [dv1, anv1] = CL__man_dvRaps(ai, ai, af, mu); [dv2, anv2] = CL__man_dvRaps(af, ai, af, mu); // sum of dv norms deltav = CL_norm(dv1) + CL_norm(dv2); // output if (res == "d") varargout = list(deltav, dv1, dv2, anv1, anv2); else varargout(1) = struct("deltav",deltav, "dv1", dv1, "dv2", dv2, "anv1", anv1, "anv2", anv2); end endfunction