J2解析预报器
2026/1/21大约 2 分钟
J2解析预报器
概述
J2解析预报器考虑地球扁率J2项的解析预报器,提供比二体预报器更高精度的轨道预报。该方法在保持较高计算效率的同时,能够考虑地球非球形引力场的主要摄动效应。
核心概念
J2摄动项
J2项是地球引力场球谐展开中的主要项,表示地球的扁率效应。J2摄动主要影响:
- 轨道面进动:轨道平面绕地球自转轴进动
- 近地点进动:近地点在轨道面内旋转
预报原理
J2解析预报器基于平均轨道根数的长期变化率:
- 计算J2摄动引起的轨道根数长期变化率
- 对轨道根数进行线性预报
- 计算对应的位置和速度向量
主要接口
类接口
class J2Analytical : public J2J4Analytical构造函数:
J2Analytical(const ModOrbElem& modOrbElem, const TimePoint& epoch, double gm, double j2, double re);
J2Analytical(const Vector3d& pos, const Vector3d& vel, const TimePoint& epoch, double gm, double j2, double re);函数接口
AST_CORE_API errc_t aJ2AnalyticalProp(double duration, double gm, double j2, double re, Vector3d& r, Vector3d& v);
AST_CORE_API errc_t aJ2AnalyticalProp(double duration, double gm, double j2, double re, ModOrbElem& modOrbElem);参数说明:
duration:传播时间间隔(秒)gm:引力常数(m³/s²)j2:J2项系数(无量纲)re:天体半径(米)r:输入/输出位置向量(米)v:输入/输出速度向量(米/秒)modOrbElem:输入/输出改进轨道根数
用法示例
1. 位置速度传播
#include "AstCore/J2Analytical.hpp"
#include "AstCore/TimePoint.hpp"
#include "AstMath/Vector.hpp"
#include "AstUtil/Literals.hpp"
#include "AstUtil/Constants.h"
#include <iostream>
#include <iomanip>
AST_USING_NAMESPACE
using namespace _AST literals;
int main()
{
// 地球参数
const double gm_earth = kEarthGrav;
const double j2_earth = 1.08262668e-3; // 地球J2系数
const double re_earth = kEarthRadius; // 地球半径
// 初始轨道状态 (近地轨道,倾角45度)
Vector3d position{7000_km, 0.0, 0.0};
Vector3d velocity{0.0, 7.5_km_s, 1.0_km_s}; // 包含Z分量产生倾角
std::cout << std::fixed << std::setprecision(3);
std::cout << "J2解析传播器 - 位置速度传播示例" << std::endl;
std::cout << "==================================" << std::endl;
std::cout << "初始轨道状态:" << std::endl;
std::cout << "位置: (" << position[0]/1000.0 << ", "
<< position[1]/1000.0 << ", "
<< position[2]/1000.0 << ") km" << std::endl;
std::cout << "速度: (" << velocity[0]/1000.0 << ", "
<< velocity[1]/1000.0 << ", "
<< velocity[2]/1000.0 << ") km/s" << std::endl;
std::cout << std::endl;
// 传播时间间隔 (1天)
double duration = 24 * 3600.0; // 秒
// 执行J2解析传播
errc_t result = aJ2AnalyticalProp(duration, gm_earth, j2_earth, re_earth, position, velocity);
if (result == eNoError) {
std::cout << "传播后轨道状态 (" << duration/3600.0 << " 小时后):" << std::endl;
std::cout << "位置: (" << position[0]/1000.0 << ", "
<< position[1]/1000.0 << ", "
<< position[2]/1000.0 << ") km" << std::endl;
std::cout << "速度: (" << velocity[0]/1000.0 << ", "
<< velocity[1]/1000.0 << ", "
<< velocity[2]/1000.0 << ") km/s" << std::endl;
// 计算轨道参数变化
double initial_radius = 7000.0; // km
double final_radius = std::sqrt(position[0]*position[0] + position[1]*position[1] + position[2]*position[2]) / 1000.0;
std::cout << std::endl;
std::cout << "轨道参数变化:" << std::endl;
std::cout << "初始轨道半径: " << initial_radius << " km" << std::endl;
std::cout << "最终轨道半径: " << final_radius << " km" << std::endl;
std::cout << "半径变化: " << (final_radius - initial_radius) << " km" << std::endl;
} else {
std::cout << "传播失败,错误码: " << result << std::endl;
}
return 0;
}2. 轨道根数传播
#include "AstCore/J2Analytical.hpp"
#include "AstCore/TimePoint.hpp"
#include "AstCore/OrbitElement.hpp"
#include "AstMath/Vector.hpp"
#include "AstUtil/Literals.hpp"
#include "AstUtil/Constants.h"
#include <iostream>
#include <iomanip>
#include <cmath>
AST_USING_NAMESPACE
using namespace _AST literals;
int main()
{
// 地球参数
const double gm_earth = kEarthGrav;
const double j2_earth = 1.08262668e-3;
const double re_earth = kEarthRadius;
// 创建改进轨道根数 (太阳同步轨道示例)
ModOrbElem orb_elem;
orb_elem.rp_ = 7000_km; // 近拱点半径
orb_elem.e_ = 0.001; // 偏心率
orb_elem.i_ = 98.0 * kDegToRad; // 倾角 (太阳同步轨道)
orb_elem.raan_ = 0.0; // 升交点赤经
orb_elem.argper_ = 0.0; // 近拱点角
orb_elem.trueA_ = 0.0; // 真近点角
std::cout << std::fixed << std::setprecision(6);
std::cout << "J2解析传播器 - 轨道根数传播示例" << std::endl;
std::cout << "==================================" << std::endl;
std::cout << "初始轨道根数:" << std::endl;
std::cout << "半长轴: " << orb_elem.getA() / 1000.0 << " km" << std::endl;
std::cout << "偏心率: " << orb_elem.e_ << std::endl;
std::cout << "倾角: " << orb_elem.i_ * kRadToDeg << " 度" << std::endl;
std::cout << "升交点赤经: " << orb_elem.raan_ * kRadToDeg << " 度" << std::endl;
std::cout << "近拱点角: " << orb_elem.argper_ * kRadToDeg << " 度" << std::endl;
std::cout << "真近点角: " << orb_elem.trueA_ * kRadToDeg << " 度" << std::endl;
std::cout << std::endl;
// 传播多个时间点,观察J2摄动效应
double duration_day = 1.0; // 天
double duration_sec = duration_day * 24 * 3600.0;
std::cout << "轨道根数随时间变化 (" << duration_day << " 天):" << std::endl;
std::cout << "时间(天)\t升交点赤经(度)\t近地点幅角(度)" << std::endl;
std::cout << "------------------------------------------------" << std::endl;
for (int i = 0; i <= 10; i++) {
double time = i * duration_sec / 10.0;
// 复制初始轨道根数
ModOrbElem current_orb = orb_elem;
// 执行J2解析传播
errc_t result = aJ2AnalyticalProp(time, gm_earth, j2_earth, re_earth, current_orb);
if (result == eNoError) {
std::cout << std::setw(6) << time/(24*3600.0) << "\t"
<< std::setw(12) << current_orb.raan_ * kRadToDeg << "\t"
<< std::setw(12) << current_orb.argper_ * kRadToDeg << std::endl;
} else {
std::cout << "传播失败,错误码: " << result << std::endl;
break;
}
}
// 计算J2引起的升交点退行率
double semi_major_axis = orb_elem.getA();
double eccentricity = orb_elem.e_;
double inclination = orb_elem.i_;
// 升交点退行率公式: dΩ/dt = -3/2 * J2 * n * (R/p)^2 * cos(i)
double n = std::sqrt(gm_earth / (semi_major_axis * semi_major_axis * semi_major_axis)); // 平均运动
double p = semi_major_axis * (1 - eccentricity * eccentricity); // 半通径
double raan_rate = -1.5 * j2_earth * n * std::pow(re_earth / p, 2) * std::cos(inclination);
std::cout << std::endl;
std::cout << "J2摄动理论值:" << std::endl;
std::cout << "升交点退行率: " << raan_rate * kRadToDeg * 86400.0 << " 度/天" << std::endl;
return 0;
}3. 类接口使用
#include "AstCore/J2Analytical.hpp"
#include "AstCore/J2Analytical.hpp"
#include "AstCore/TimePoint.hpp"
#include "AstMath/Vector.hpp"
#include "AstUtil/Literals.hpp"
#include "AstUtil/Constants.h"
#include <iostream>
#include <iomanip>
AST_USING_NAMESPACE
using namespace _AST literals;
int main()
{
// 地球参数
const double gm_earth = kEarthGrav;
const double j2_earth = 1.08262668e-3;
const double re_earth = kEarthRadius;
// 创建时间点
TimePoint epoch = TimePoint::FromUTC(2026, 1, 1, 0, 0, 0);
std::cout << std::fixed << std::setprecision(6);
std::cout << "J2解析传播器 - 类接口使用示例" << std::endl;
std::cout << "==================================" << std::endl;
// 示例1: 使用轨道根数构造
std::cout << "\n示例1: 使用轨道根数构造传播器" << std::endl;
{
ModOrbElem orb_elem;
orb_elem.rp_ = 7000_km;
orb_elem.e_ = 0.01;
orb_elem.i_ = 45.0 * kDegToRad;
orb_elem.raan_ = 0.0;
orb_elem.argper_ = 0.0;
orb_elem.trueA_ = 0.0;
J2Analytical propagator1(orb_elem, epoch, gm_earth, j2_earth, re_earth);
std::cout << "使用轨道根数成功创建J2传播器" << std::endl;
}
// 示例2: 使用位置速度构造
std::cout << "\n示例2: 使用位置速度构造传播器" << std::endl;
{
Vector3d position{7000_km, 0.0, 0.0};
Vector3d velocity{0.0, 7.5_km_s, 1.0_km_s};
J2Analytical propagator2(position, velocity, epoch, gm_earth, j2_earth, re_earth);
std::cout << "使用位置速度成功创建J2传播器" << std::endl;
// 演示继承的基类功能
std::cout << "传播器类型: J2Analytical" << std::endl;
std::cout << "继承自: J2J4Analytical" << std::endl;
}
// 示例3: 完整的传播流程
std::cout << "\n示例3: 完整的传播流程演示" << std::endl;
{
Vector3d position{7000_km, 0.0, 0.0};
Vector3d velocity{0.0, 7.5_km_s, 0.0};
// 创建传播器
J2Analytical propagator(position, velocity, epoch, gm_earth, j2_earth, re_earth);
// 传播多个时间点
double duration = 24 * 3600.0; // 1天
std::cout << "时间传播结果:" << std::endl;
std::cout << "时间(小时)\tX(km)\t\tY(km)\t\tZ(km)" << std::endl;
std::cout << "----------------------------------------------" << std::endl;
for (int i = 0; i <= 4; i++) {
double time = i * duration / 4.0;
// 注意:这里需要调用实际的传播方法
// 由于J2Analytical类的具体传播方法未在头文件中显示,
// 这里仅演示类的基本使用
std::cout << std::setw(8) << time/3600.0 << "\t"
<< std::setw(10) << "待实现" << "\t"
<< std::setw(10) << "待实现" << "\t"
<< std::setw(10) << "待实现" << std::endl;
}
}
std::cout << "\n类接口使用演示完成" << std::endl;
std::cout << "注意:实际传播功能需要调用类的具体方法" << std::endl;
return 0;
}精度与限制
精度特点
- 计算精度:考虑J2项长期摄动,精度优于二体传播器
- 适用范围:适用于需要考虑地球扁率效应的近地轨道
- 计算效率:较高,适合需要较高效率的应用
使用限制
- 仅考虑J2项摄动,忽略高阶项和其他摄动力
- 对于高精度轨道预报需求,建议使用HPOP等数值方法
- 长时间传播时,由于模型简化,误差会累积
地球J2参数
对于地球轨道计算,推荐使用以下参数:
- J2系数:1.08262668e-3
- 地球半径:6378137.0米
- 引力常数:3.986004418e14 m³/s²
注意事项
- 输入的位置和速度向量应在同一惯性坐标系下,例如J2000、ICRF、TOD等
- J2系数和天体半径应与中心天体匹配
- 传播时间间隔不宜过长,避免数值误差累积
API参考
///
/// @file J2Analytical.hpp
/// @brief ~
/// @details ~
/// @author axel
/// @date 2025-12-31
/// @copyright 版权所有 (C) 2025-present, ast项目.
///
/// ast项目(https://github.com/space-ast/ast)
/// 本项目基于 Apache 2.0 开源许可证分发。
/// 您可在遵守许可证条款的前提下使用、修改和分发本软件。
/// 许可证全文请见:
///
/// http://www.apache.org/licenses/LICENSE-2.0
///
/// 重要须知:
/// 软件按"现有状态"提供,无任何明示或暗示的担保条件。
/// 除非法律要求或书面同意,作者与贡献者不承担任何责任。
/// 使用本软件所产生的风险,需由您自行承担。
#pragma once
#include "AstGlobal.h"
#include "AstCore/OrbitElement.hpp"
#include "AstCore/TimePoint.hpp"
#include "J2J4Analytical.hpp"
AST_NAMESPACE_BEGIN
/*!
@ingroup Core
@defgroup Propagator 轨道预报
@brief 提供二体、J2、J4、HPOP等轨道预报器
*/
/*!
@addtogroup Propagator
@{
*/
class ModOrbElem; ///< 改进轨道根数
class TimePoint; ///< 时间点
class AST_CORE_API J2Analytical : public J2J4Analytical
{
public:
/// @brief 构造函数
/// @details ~
/// @param modOrbElem 改进轨道根数
/// @param epoch 参考时间点
/// @param gm 引力常数
/// @param j2 J2项
/// @param re 天体半径
J2Analytical(const ModOrbElem& modOrbElem, const TimePoint& epoch, double gm, double j2, double re);
/// @brief 构造函数
/// @details ~
/// @param pos 位置向量
/// @param vel 速度向量
/// @param epoch 参考时间点
/// @param gm 引力常数
/// @param j2 J2项
/// @param re 天体半径
J2Analytical(const Vector3d& pos, const Vector3d& vel, const TimePoint& epoch, double gm, double j2, double re);
protected:
double j2_; ///< J2项
double re_; ///< 天体半径
};
/// @brief J2轨道预报
/// @details
/// @param duration 时间间隔
/// @param gm 引力常数
/// @param j2 J2项
/// @param re 天体半径
/// @param r 位置向量
/// @param v 速度向量
/// @return 错误码
AST_CORE_API errc_t aJ2AnalyticalProp(double duration, double gm, double j2, double re, Vector3d& r, Vector3d& v);
/// @brief J2轨道预报
/// @details
/// @param duration 时间间隔
/// @param gm 引力常数
/// @param j2 J2项
/// @param re 天体半径
/// @param modOrbElem 改进轨道根数
/// @return 错误码
AST_CORE_API errc_t aJ2AnalyticalProp(double duration, double gm, double j2, double re, ModOrbElem& modOrbElem);
/*! @} */
AST_NAMESPACE_END