Vulcan
Aerospace Engineering Utilities Built on Janus
Loading...
Searching...
No Matches
TransferMechanics.hpp
Go to the documentation of this file.
1// Vulcan Transfer Mechanics
2// Orbital maneuver delta-v computations
3#pragma once
4
5#include <cmath>
6#include <janus/janus.hpp>
7#include <tuple>
8#include <utility>
12
14
27template <typename Scalar>
28std::pair<Scalar, Scalar> hohmann_delta_v(const Scalar &r1, const Scalar &r2,
29 double mu = constants::earth::mu) {
30 // Transfer orbit semi-major axis
31 const Scalar a_t = (r1 + r2) / 2.0;
32
33 // Velocities using vis-viva
34 const Scalar v1 = quantities::circular_velocity(r1, mu);
35 const Scalar v2 = quantities::circular_velocity(r2, mu);
36 const Scalar v_t1 = quantities::velocity(r1, a_t, mu);
37 const Scalar v_t2 = quantities::velocity(r2, a_t, mu);
38
39 // Delta-v magnitudes (absolute value for outward transfer)
40 const Scalar dv1 = janus::abs(v_t1 - v1);
41 const Scalar dv2 = janus::abs(v2 - v_t2);
42
43 return {dv1, dv2};
44}
45
55template <typename Scalar>
56Scalar hohmann_total_delta_v(const Scalar &r1, const Scalar &r2,
57 double mu = constants::earth::mu) {
58 auto [dv1, dv2] = hohmann_delta_v(r1, r2, mu);
59 return dv1 + dv2;
60}
61
71template <typename Scalar>
72Scalar hohmann_transfer_time(const Scalar &r1, const Scalar &r2,
73 double mu = constants::earth::mu) {
74 const Scalar a_t = (r1 + r2) / 2.0;
75 return M_PI * janus::sqrt(a_t * a_t * a_t / mu);
76}
77
91template <typename Scalar>
92std::tuple<Scalar, Scalar, Scalar>
93bielliptic_delta_v(const Scalar &r1, const Scalar &r2, const Scalar &r_b,
94 double mu = constants::earth::mu) {
95 // First transfer orbit: r1 to r_b
96 const Scalar a1 = (r1 + r_b) / 2.0;
97 // Second transfer orbit: r_b to r2
98 const Scalar a2 = (r_b + r2) / 2.0;
99
100 const Scalar v1 = quantities::circular_velocity(r1, mu);
101 const Scalar v2 = quantities::circular_velocity(r2, mu);
102
103 const Scalar v_t1_peri = quantities::velocity(r1, a1, mu);
104 const Scalar v_t1_apo = quantities::velocity(r_b, a1, mu);
105 const Scalar v_t2_apo = quantities::velocity(r_b, a2, mu);
106 const Scalar v_t2_peri = quantities::velocity(r2, a2, mu);
107
108 const Scalar dv1 = janus::abs(v_t1_peri - v1);
109 const Scalar dv2 = janus::abs(v_t2_apo - v_t1_apo);
110 const Scalar dv3 = janus::abs(v2 - v_t2_peri);
111
112 return {dv1, dv2, dv3};
113}
114
126template <typename Scalar>
127Scalar plane_change_delta_v(const Scalar &v, const Scalar &delta_i) {
128 return 2.0 * v * janus::sin(janus::abs(delta_i) / 2.0);
129}
130
142template <typename Scalar>
143Scalar combined_maneuver_delta_v(const Scalar &v1, const Scalar &v2,
144 const Scalar &delta_i) {
145 const Scalar cos_di = janus::cos(delta_i);
146 return janus::sqrt(v1 * v1 + v2 * v2 - 2.0 * v1 * v2 * cos_di);
147}
148
149} // namespace vulcan::orbital::transfer
constexpr double mu
Gravitational parameter (GM) [m^3/s^2].
Definition Constants.hpp:13
Scalar velocity(const Scalar &r, const Scalar &a, double mu=constants::earth::mu)
Orbital velocity (vis-viva equation).
Definition OrbitalQuantities.hpp:35
Scalar circular_velocity(const Scalar &r, double mu=constants::earth::mu)
Circular orbit velocity at given radius.
Definition OrbitalQuantities.hpp:75
Definition TransferMechanics.hpp:13
Scalar combined_maneuver_delta_v(const Scalar &v1, const Scalar &v2, const Scalar &delta_i)
Combined plane change and altitude change.
Definition TransferMechanics.hpp:143
Scalar plane_change_delta_v(const Scalar &v, const Scalar &delta_i)
Simple plane change delta-v.
Definition TransferMechanics.hpp:127
std::pair< Scalar, Scalar > hohmann_delta_v(const Scalar &r1, const Scalar &r2, double mu=constants::earth::mu)
Hohmann transfer delta-v.
Definition TransferMechanics.hpp:28
Scalar hohmann_total_delta_v(const Scalar &r1, const Scalar &r2, double mu=constants::earth::mu)
Total Hohmann transfer delta-v.
Definition TransferMechanics.hpp:56
std::tuple< Scalar, Scalar, Scalar > bielliptic_delta_v(const Scalar &r1, const Scalar &r2, const Scalar &r_b, double mu=constants::earth::mu)
Bielliptic transfer delta-v.
Definition TransferMechanics.hpp:93
Scalar hohmann_transfer_time(const Scalar &r1, const Scalar &r2, double mu=constants::earth::mu)
Hohmann transfer time.
Definition TransferMechanics.hpp:72