Skip to content

Commit 23fe424

Browse files
author
tteil
committed
add python unit test
1 parent 01794f5 commit 23fe424

File tree

1 file changed

+111
-0
lines changed

1 file changed

+111
-0
lines changed
Lines changed: 111 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,111 @@
1+
# ISC License
2+
#
3+
# Copyright (c) 2025, Laboratory for Atmospheric and Space Physics, University of Colorado at Boulder
4+
#
5+
# Permission to use, copy, modify, and/or distribute this software for any
6+
# purpose with or without fee is hereby granted, provided that the above
7+
# copyright notice and this permission notice appear in all copies.
8+
#
9+
# THE SOFTWARE IS PROVIDED "AS IS" AND THE AUTHOR DISCLAIMS ALL WARRANTIES
10+
# WITH REGARD TO THIS SOFTWARE INCLUDING ALL IMPLIED WARRANTIES OF
11+
# MERCHANTABILITY AND FITNESS. IN NO EVENT SHALL THE AUTHOR BE LIABLE FOR
12+
# ANY SPECIAL, DIRECT, INDIRECT, OR CONSEQUENTIAL DAMAGES OR ANY DAMAGES
13+
# WHATSOEVER RESULTING FROM LOSS OF USE, DATA OR PROFITS, WHETHER IN AN
14+
# ACTION OF CONTRACT, NEGLIGENCE OR OTHER TORTIOUS ACTION, ARISING OUT OF
15+
# OR IN CONNECTION WITH THE USE OR PERFORMANCE OF THIS SOFTWARE.
16+
17+
from Basilisk.utilities import orbitalMotion as oe_py
18+
from Basilisk.architecture.orbitalMotion import OrbitalMotion as oe_cpp
19+
from Basilisk.architecture.orbitalMotion import ClassicalElements as elements_cpp
20+
import numpy as np
21+
22+
23+
def test_rv_oe():
24+
gravitational_parameter = 1e10
25+
position = np.array([1e5, 1e6, 1e2])
26+
velocity = np.array([10, 20, 30])
27+
28+
elements = oe_py.rv2elem(gravitational_parameter, position, velocity)
29+
r, v = oe_py.elem2rv(gravitational_parameter, elements)
30+
np.testing.assert_almost_equal(position, r, 8)
31+
np.testing.assert_almost_equal(velocity, v, 8)
32+
33+
gravitational_parameter = 1e10
34+
position = np.array([1e5, 1e6, 0])
35+
h = np.array([0.0, 0.0, 1.0])
36+
velocity = np.sqrt(gravitational_parameter / np.linalg.norm(position)) * np.cross(
37+
h, position / np.linalg.norm(position)
38+
)
39+
elements = oe_py.rv2elem(gravitational_parameter, position, velocity)
40+
np.testing.assert_almost_equal(elements.a, np.linalg.norm(position), 8)
41+
np.testing.assert_almost_equal(elements.e, 0, 8)
42+
np.testing.assert_almost_equal(elements.i, 0, 8)
43+
44+
45+
def test_compare_conversions():
46+
gravitational_parameter = 1e10
47+
position = np.array([1e5, 1e6, 1e2])
48+
velocity = np.array([10, 20, 30])
49+
50+
# state vector to orbital elements
51+
elements_py = oe_py.rv2elem(gravitational_parameter, position, velocity)
52+
elements_cpp = oe_cpp.cartesianStateToElements(gravitational_parameter, position, velocity)
53+
np.testing.assert_almost_equal(elements_py.a, elements_cpp.semiMajorAxis, 8)
54+
np.testing.assert_almost_equal(elements_py.e, elements_cpp.eccentricity, 8)
55+
np.testing.assert_almost_equal(elements_py.i, elements_cpp.inclination, 8)
56+
np.testing.assert_almost_equal(elements_py.Omega, elements_cpp.rightAscensionAscendingNode, 8)
57+
np.testing.assert_almost_equal(elements_py.omega, elements_cpp.argPeriapsis, 8)
58+
np.testing.assert_almost_equal(elements_py.f, elements_cpp.trueAnomaly, 8)
59+
60+
# orbital elements to state vectors
61+
r_py, v_py = oe_py.elem2rv(gravitational_parameter, elements_py)
62+
state_cpp = oe_cpp.elementsToCartesianState(gravitational_parameter, elements_cpp)
63+
np.testing.assert_almost_equal(r_py, np.array(state_cpp.position).flatten(), 8)
64+
np.testing.assert_almost_equal(v_py, np.array(state_cpp.velocity).flatten(), 8)
65+
66+
eccentric_anomaly = 0.5
67+
# eccentric to true anomaly and back
68+
true_anomaly_py = oe_py.E2f(eccentric_anomaly, elements_py.e)
69+
true_anomaly_cpp = oe_cpp.eccentricToTrueAnomaly(eccentric_anomaly, elements_py.e)
70+
np.testing.assert_almost_equal(true_anomaly_py, true_anomaly_cpp, 8)
71+
72+
eccentric_anomaly_py = oe_py.f2E(true_anomaly_py, elements_py.e)
73+
eccentric_anomaly_cpp = oe_cpp.trueToEccentricAnomaly(true_anomaly_cpp, elements_py.e)
74+
np.testing.assert_almost_equal(eccentric_anomaly, eccentric_anomaly_py, 8)
75+
np.testing.assert_almost_equal(eccentric_anomaly, eccentric_anomaly_cpp, 8)
76+
77+
# eccentric to mean anomaly
78+
mean_anomaly_py = oe_py.E2M(eccentric_anomaly, elements_py.e)
79+
mean_anomaly_cpp = oe_cpp.eccentricToMeanAnomaly(eccentric_anomaly, elements_py.e)
80+
np.testing.assert_almost_equal(mean_anomaly_py, mean_anomaly_cpp, 8)
81+
82+
eccentric_anomaly_py = oe_py.M2E(mean_anomaly_py, elements_py.e)
83+
eccentric_anomaly_cpp = oe_cpp.meanToEccentricAnomaly(mean_anomaly_cpp, elements_py.e)
84+
np.testing.assert_almost_equal(eccentric_anomaly_py, eccentric_anomaly_cpp, 8)
85+
np.testing.assert_almost_equal(eccentric_anomaly_py, eccentric_anomaly, 8)
86+
np.testing.assert_almost_equal(eccentric_anomaly_cpp, eccentric_anomaly, 8)
87+
88+
# true to hyperbolic anomaly
89+
e_hyperbolic = 1.2
90+
true_anomaly = 0.234
91+
hyper_anomaly_py = oe_py.f2H(true_anomaly, e_hyperbolic)
92+
hyper_anomaly_cpp = oe_cpp.trueToHyperbolicAnomaly(true_anomaly, e_hyperbolic)
93+
np.testing.assert_almost_equal(hyper_anomaly_py, hyper_anomaly_cpp, 8)
94+
95+
true_anomaly_py = oe_py.H2f(hyper_anomaly_py, e_hyperbolic)
96+
true_anomaly_cpp = oe_cpp.hyperbolicToTrueAnomaly(hyper_anomaly_cpp, e_hyperbolic)
97+
np.testing.assert_almost_equal(true_anomaly_py, true_anomaly_cpp, 8)
98+
np.testing.assert_almost_equal(true_anomaly_py, true_anomaly, 8)
99+
np.testing.assert_almost_equal(true_anomaly_cpp, true_anomaly, 8)
100+
101+
# true to mean anomaly
102+
true_anomaly = 0.234
103+
mean_anomaly_py = oe_py.E2M(oe_py.f2E(true_anomaly, elements_py.e), elements_py.e)
104+
mean_anomaly_cpp = oe_cpp.trueToMeanAnomaly(true_anomaly, elements_py.e)
105+
np.testing.assert_almost_equal(mean_anomaly_py, mean_anomaly_cpp, 8)
106+
107+
true_anomaly_py = oe_py.E2f(oe_py.M2E(mean_anomaly_py, elements_py.e), elements_py.e)
108+
true_anomaly_cpp = oe_cpp.meanToTrueAnomaly(mean_anomaly_cpp, elements_py.e)
109+
np.testing.assert_almost_equal(true_anomaly_py, true_anomaly_cpp, 8)
110+
np.testing.assert_almost_equal(true_anomaly_py, true_anomaly, 8)
111+
np.testing.assert_almost_equal(true_anomaly_cpp, true_anomaly, 8)

0 commit comments

Comments
 (0)