astro/basic/orbit_test.go
starainrt 3ffdbe0034
feat: 扩展天文计算能力
- 新增日食、月食、本地可见性、中心线、半影区域、SVG 图示与沙罗周期信息
- 新增行星冲合、留、方照、物理星历、视直径、相位、亮肢角、轨道节点等计算
- 新增木星伽利略卫星位置、现象与接触事件计算
- 新增恒星星表、星座判定、自行修正与观测辅助能力
- 新增 coord、formula、orbit、sundial、lite/sun、lite/moon 等扩展包
- 完善农历年号、月相英文别名、视差角、大气质量、折射、日晷与双星计算
- 增加 NASA、JPL Horizons、IMCCE 等回归测试数据与基线测试
- 重构基础算法文件组织,补充大量公开 API 注释和语义回归测试
- 更新中文和英文 README,补充示例、精度说明、SVG 配图
2026-05-01 22:38:44 +08:00

127 lines
3.7 KiB
Go

package basic
import (
"math"
"testing"
)
func TestOrbitCircularStateAtEpoch(t *testing.T) {
elements := OrbitElements{
EpochJD: orbitReferenceJD,
A: 2,
E: 0,
I: 0,
Omega: 0,
W: 0,
M0: 0,
}
vector := OrbitHeliocentricXYZJ2000(orbitReferenceJD, elements)
assertSameFloat(t, "x", vector[0], 2)
assertSameFloat(t, "y", vector[1], 0)
assertSameFloat(t, "z", vector[2], 0)
lon, lat, distance := OrbitHeliocentricEclipticJ2000(orbitReferenceJD, elements)
assertSameFloat(t, "lon", lon, 0)
assertSameFloat(t, "lat", lat, 0)
assertSameFloat(t, "distance", distance, 2)
}
func TestOrbitQuarterPeriodOnCircularOrbit(t *testing.T) {
elements := OrbitElements{
EpochJD: orbitReferenceJD,
A: 1,
E: 0,
I: 0,
Omega: 0,
W: 0,
M0: 0,
}
quarterPeriodDays := 90 / OrbitMeanMotion(elements)
vector := OrbitHeliocentricXYZJ2000(orbitReferenceJD+quarterPeriodDays, elements)
if math.Abs(vector[0]) > 1e-10 || math.Abs(vector[1]-1) > 1e-10 || math.Abs(vector[2]) > 1e-10 {
t.Fatalf("quarter-period vector mismatch: got %+v", vector)
}
}
func TestOrbitMeanAndTrueAnomalyCircularMatch(t *testing.T) {
elements := OrbitElements{
EpochJD: orbitReferenceJD,
A: 1.523679,
E: 0,
I: 1.85,
Omega: 49.5,
W: 286.5,
M0: 123.4,
}
jd := orbitReferenceJD + 123.456
meanAnomaly := OrbitMeanAnomaly(jd, elements)
trueAnomaly := OrbitTrueAnomaly(jd, elements)
if math.Abs(meanAnomaly-trueAnomaly) > 1e-12 {
t.Fatalf("circular mean/true anomaly mismatch: mean=%.18f true=%.18f", meanAnomaly, trueAnomaly)
}
}
func TestPerihelionFormMatchesClassicalEllipticState(t *testing.T) {
classical := OrbitElements{
EpochJD: 2457305.5,
A: 3.462249489765068,
E: 0.6409081306555051,
I: 7.040294906760007,
Omega: 50.13557380441372,
W: 12.79824973415729,
M0: 8.859927418758764,
}
perihelion := OrbitElements{
Q: 1.243265641416762,
E: classical.E,
I: classical.I,
Omega: classical.Omega,
W: classical.W,
TpJD: 2457247.588657863465,
}
jd := 2457308.5
classicalVector := OrbitHeliocentricXYZJ2000(jd, classical)
perihelionVector := OrbitHeliocentricXYZJ2000(jd, perihelion)
for i := range classicalVector {
if math.Abs(classicalVector[i]-perihelionVector[i]) > 1e-11 {
t.Fatalf("component %d mismatch: classical=%.18f perihelion=%.18f", i, classicalVector[i], perihelionVector[i])
}
}
}
func TestParabolicPerihelionStateAtTp(t *testing.T) {
elements := OrbitElements{Q: 1, E: 1, I: 0, Omega: 0, W: 0, TpJD: orbitReferenceJD}
vector := OrbitHeliocentricXYZJ2000(orbitReferenceJD, elements)
assertSameFloat(t, "x", vector[0], 1)
assertSameFloat(t, "y", vector[1], 0)
assertSameFloat(t, "z", vector[2], 0)
if !math.IsNaN(OrbitMeanAnomaly(orbitReferenceJD, elements)) {
t.Fatalf("parabolic mean anomaly should be NaN")
}
}
func TestHyperbolicOrbitProducesFiniteState(t *testing.T) {
elements := OrbitElements{Q: 0.5, E: 1.2, I: 12, Omega: 30, W: 45, TpJD: orbitReferenceJD}
vector := OrbitHeliocentricXYZJ2000(orbitReferenceJD+20, elements)
for i, component := range vector {
if math.IsNaN(component) || math.IsInf(component, 0) {
t.Fatalf("component %d not finite: %.18f", i, component)
}
}
}
func TestOrbitInvalidEllipticElementsReturnNaN(t *testing.T) {
elements := OrbitElements{EpochJD: orbitReferenceJD, A: 1, E: 1.1}
if !math.IsNaN(OrbitMeanMotion(elements)) {
t.Fatalf("expected NaN mean motion for invalid elements")
}
vector := OrbitHeliocentricXYZJ2000(orbitReferenceJD, elements)
for i, component := range vector {
if !math.IsNaN(component) {
t.Fatalf("component %d expected NaN, got %.18f", i, component)
}
}
}