🛰️航天仿真算法库 SpaceAST 0.0.1
载入中...
搜索中...
未找到
轨道计算

详细描述

提供轨道根数、轨道参数计算及转换相关接口。

class  ast::IOrbitDesigner
 
class  ast::BaseOrbitDesigner
 
class  ast::CircularOrbitDesigner
 
class  ast::CriticallyInclinedOrbitDesigner
 
class  ast::RepeatingOrbitDesigner
 
class  ast::RepeatingSunSyncOrbitDesigner
 
class  ast::SimpleOrbitDesigner
 
class  ast::StationaryOrbitDesigner
 
class  ast::SunSynchronousOrbitDesigner
 
class  ast::CartState
 
class  ast::OrbElem
 
class  ast::ModOrbElem
 
class  ast::EquinElem
 
class  ast::ModEquinElem
 
class  ast::DelaunayElem
 

函数

double ast::aRAANRate (double gm, double j2, double rb, double a, double ecc, double inc)
 
double ast::aArgPeriRate (double gm, double j2, double rb, double a, double ecc, double inc)
 
double ast::aMeanAnomalyRate (double gm, double j2, double rb, double a, double ecc, double inc)
 
double ast::aMeanArgLatRate (double gm, double j2, double rb, double a, double ecc, double inc)
 
double ast::aJ2Period (double gm, double j2, double rb, double a, double ecc, double inc)
 
errc_t ast::aSunSynchronousInclination (double gm, double j2, double rb, double bodyMeanMotion, double a, double ecc, double &inc)
 
errc_t ast::aSunSynchronousSemiMajorAxis (double gm, double j2, double rb, double bodyMeanMotion, double inc, double ecc, double &semiMajorAxis)
 
A_ALWAYS_INLINE double ast::aSunSynchronousInclination (double gm, double j2, double rb, double bodyMeanMotion, double a, double ecc)
 
A_ALWAYS_INLINE double ast::aSunSynchronousSemiMajorAxis (double gm, double j2, double rb, double bodyMeanMotion, double inc, double ecc)
 
errc_t ast::coe2rv (const double *coe, double gm, double *pos, double *vel)
 
errc_t ast::coe2mee (const double *coe, double *mee)
 
void ast::ee2rv (const double *ee, double gm, double *pos, double *vel)
 
void ast::mee2rv (const double *mee, double gm, double *pos, double *vel)
 
errc_t ast::rv2mee (const double *pos, const double *vel, double gm, double *mee)
 
errc_t ast::mee2coe (const double *mee, double *coe)
 
void ast::rv2ee (const double *pos, const double *vel, double gm, double *ee)
 
errc_t ast::rv2moe (const double *pos, const double *vel, double gm, double *moe)
 
errc_t ast::rv2coe (const double *pos, const double *vel, double gm, double *coe)
 
errc_t ast::ee2moe (const double *ee, double *moe)
 
errc_t ast::moe2ee (const double *moe, double *ee)
 
errc_t ast::moe2coe (const double *moe, double *coe)
 
void ast::coe2moe (const double *coe, double *moe)
 
errc_t ast::moe2rv (const double *moe, double gm, double *pos, double *vel)
 
void ast::moe2mee (const double *moe, double *mee)
 
void ast::coe2ee (const double *coe, double *ee)
 
void ast::ee2coe (const double *ee, double *coe)
 
void ast::ee2mee (const double *ee, double *mee)
 
void ast::mee2ee (const double *mee, double *ee)
 
void ast::mee2moe (const double *mee, double *moe)
 
errc_t ast::coe2dela (const double *coe, double gm, double *dela)
 
errc_t ast::dela2coe (const double *dela, double gm, double *coe)
 
void ast::aModEquinElemToCart (const ModEquinElem &mee, double gm, Vector3d &pos, Vector3d &vel)
 
errc_t ast::aCartToModEquinElem (const Vector3d &pos, const Vector3d &vel, double gm, ModEquinElem &mee)
 
errc_t ast::aOrbElemToModEquinElem (const OrbElem &elem, ModEquinElem &mee)
 
errc_t ast::aModEquinElemToOrbElem (const ModEquinElem &mee, OrbElem &elem)
 
errc_t ast::aCartToModOrbElem (const Vector3d &pos, const Vector3d &vel, double gm, ModOrbElem &modOrb)
 
errc_t ast::aCartToOrbElem (const Vector3d &pos, const Vector3d &vel, double gm, OrbElem &elem)
 
errc_t ast::aEquinElemToModOrb (const EquinElem &equinElem, ModOrbElem &modOrb)
 
errc_t ast::aModOrbToEquinElem (const ModOrbElem &modOrb, EquinElem &equinElem)
 
errc_t ast::aModOrbElemToCart (const ModOrbElem &modOrb, double gm, Vector3d &pos, Vector3d &vel)
 
errc_t ast::aOrbElemToCart (const OrbElem &elem, double gm, Vector3d &pos, Vector3d &vel)
 
void ast::aCartToEquinElem (const Vector3d &pos, const Vector3d &vel, double gm, EquinElem &equinElem)
 
void ast::aEquinElemToCart (const EquinElem &equinElem, double gm, Vector3d &pos, Vector3d &vel)
 
errc_t ast::aOrbElemToDelaunay (const OrbElem &elem, double gm, DelaunayElem &delaunay)
 
errc_t ast::aDelaunayToOrbElem (const DelaunayElem &delaunay, double gm, OrbElem &elem)
 
A_ALWAYS_INLINE ModOrbElem ast::aCartToModOrbElem (const Vector3d &r, const Vector3d &v, double gm)
 
double ast::aApoAltToApoRad (double apogeeAlt, double bodyRadius)
 
double ast::aApoAltToMeanMotion (double apogeeAlt, double eccentricity, double bodyRadius, double gm)
 
double ast::aApoAltToPeriAlt (double apogeeAlt, double eccentricity, double bodyRadius)
 
double ast::aApoAltToPeriRad (double apogeeAlt, double eccentricity, double bodyRadius)
 
double ast::aApoAltToPeriod (double apogeeAlt, double eccentricity, double bodyRadius, double gm)
 
double ast::aApoAltToSMA (double apogeeAlt, double eccentricity, double bodyRadius)
 
double ast::aApoRadToApoAlt (double apogeeRad, double bodyRadius)
 
double ast::aApoRadToMeanMotion (double apogeeRad, double eccentricity, double gm)
 
double ast::aApoRadToPeriAlt (double apogeeRad, double eccentricity, double bodyRadius)
 
double ast::aApoRadToPeriod (double apogeeRad, double eccentricity, double gm)
 
double ast::aApoRadToPeriRad (double apogeeRad, double eccentricity)
 
double ast::aApoRadToEcc (double apogeeRad, double semiMajorAxis)
 
double ast::aApoRadToSMA (double apogeeRad, double eccentricity)
 
double ast::aEccToMean (double eccAnomaly, double eccentricity)
 
double ast::aEccToTrue (double eccAnomaly, double eccentricity)
 
double ast::aEccToTimePastAscNode (double eccAnomaly, double argPeri, double semiMajorAxis, double eccentricity, double gm)
 
double ast::aEccToTimePastPeri (double eccAnomaly, double semiMajorAxis, double eccentricity, double gm)
 
double ast::aMeanMotionToApoAlt (double meanMotion, double eccentricity, double bodyRadius, double gm)
 
double ast::aMeanMotionToApoRad (double meanMotion, double eccentricity, double gm)
 
double ast::aMeanMotionToPeriAlt (double meanMotion, double eccentricity, double bodyRadius, double gm)
 
double ast::aMeanMotionToPeriRad (double meanMotion, double eccentricity, double gm)
 
double ast::aMeanMotionToPeriod (double meanMotn)
 
double ast::aMeanMotionToSMA (double meanMotn, double gm)
 
double ast::aMeanToEcc (double meanAnomaly, double eccentricity, double eps=1e-14, int maxIter=100)
 
double ast::aMeanToTimePastAscNode (double meanAnomaly, double argPeri, double semiMajorAxis, double eccentricity, double gm)
 
double ast::aMeanToTimePastPeri (double meanAnomaly, double semiMajorAxis, double gm)
 
double ast::aMeanToTrue (double meanAnomaly, double eccentricity, double eps=1e-14, int maxIter=100)
 
double ast::aPeriAltToApoAlt (double perigeeAlt, double eccentricity, double bodyRadius)
 
double ast::aPeriAltToApoRad (double perigeeAlt, double eccentricity, double bodyRadius)
 
double ast::aPeriAltToMeanMotion (double perigeeAlt, double eccentricity, double bodyRadius, double gm)
 
double ast::aPeriAltToPeriRad (double perigeeAlt, double bodyRadius)
 
double ast::aPeriAltToPeriod (double perigeeAlt, double eccentricity, double bodyRadius, double gm)
 
double ast::aPeriAltToSMA (double perigeeAlt, double eccentricity, double bodyRadius)
 
double ast::aPeriRadToApoAlt (double perigeeRad, double eccentricity, double bodyRadius)
 
double ast::aPeriRadToApoRad (double perigeeRad, double eccentricity)
 
double ast::aPeriRadToMeanMotion (double perigeeRad, double eccentricity, double gm)
 
double ast::aPeriRadToPeriAlt (double perigeeRad, double bodyRadius)
 
double ast::aPeriRadToPeriod (double perigeeRad, double eccentricity, double gm)
 
double ast::aPeriRadToSMA (double perigeeRad, double eccentricity)
 
double ast::aPeriRadToEcc (double perigeeRad, double semiMajorAxis)
 
double ast::aPeriodToApoAlt (double period, double eccentricity, double bodyRadius, double gm)
 
double ast::aPeriodToApoRad (double period, double eccentricity, double gm)
 
double ast::aPeriodToMeanMotion (double period)
 
double ast::aPeriodToPeriAlt (double period, double eccentricity, double bodyRadius, double gm)
 
double ast::aPeriodToPeriRad (double period, double eccentricity, double gm)
 
double ast::aPeriodToSMA (double period, double gm)
 
double ast::aPeriRadApoRadToEcc (double perigeeRad, double apogeeRad)
 
double ast::aPeriAltApoAltToEcc (double perigeeAlt, double apogeeAlt, double bodyRadius)
 
double ast::aSMAToApoAlt (double semiMajorAxis, double eccentricity, double bodyRadius)
 
double ast::aSMAToApoRad (double semiMajorAxis, double eccentricity)
 
double ast::aSMAToMeanMotion (double semiMajorAxis, double gm)
 
double ast::aSMAToPeriAlt (double semiMajorAxis, double eccentricity, double bodyRadius)
 
double ast::aSMAToPeriRad (double semiMajorAxis, double eccentricity)
 
double ast::aSMAToPeriod (double semiMajorAxis, double gm)
 
double ast::aSMAToSMinAx (double semiMajorAxis, double eccentricity)
 
double ast::aSMAToSParam (double semiMajorAxis, double eccentricity)
 
double ast::aSMinAxToSMA (double semiminorAxis, double eccentricity)
 
double ast::aTimePastAscNodeToEcc (double TimePastAscNode, double argPeri, double semiMajorAxis, double eccentricity, double gm, double eps=1e-14, int maxIter=100)
 
double ast::aTimePastAscNodeToMean (double TimePastAscNode, double argPeri, double semiMajorAxis, double eccentricity, double gm)
 
double ast::aTimePastAscNodeToTimePastPeri (double TimePastAscNode, double argPeri, double semiMajorAxis, double eccentricity, double gm)
 
double ast::aTimePastAscNodeToTrue (double TimePastAscNode, double argPeri, double semiMajorAxis, double eccentricity, double gm, double eps=1e-14, int maxIter=100)
 
double ast::aTimePastPeriToEcc (double TimePastPeri, double semiMajorAxis, double eccentricity, double gm, double eps=1e-14, int maxIter=100)
 
double ast::aTimePastPeriToMean (double TimePastPeri, double semiMajorAxis, double gm)
 
double ast::aTimePastPeriToTrue (double TimePastPeri, double semiMajorAxis, double eccentricity, double gm, double eps=1e-14, int maxIter=100)
 
double ast::aTimePastPeriToTimePastAscNode (double TimePastPeri, double argPeri, double semiMajorAxis, double eccentricity, double gm)
 
double ast::aTrueToEcc (double trueAnomaly, double eccentricity)
 
double ast::aTrueToMean (double trueAnomaly, double eccentricity)
 
double ast::aTrueToTimePastAscNode (double trueAnomaly, double argPeri, double semiMajorAxis, double eccentricity, double gm)
 
double ast::aTrueToTimePastPeri (double trueAnomaly, double semiMajorAxis, double eccentricity, double gm)
 
double ast::aTrueToArgLat (double trueAnomaly, double argPeri)
 
double ast::aArgLatToTrue (double argLat, double argPeri)
 
double ast::aTrueToTrueLong (double trueAnomaly, double argPeri, double raan)
 
double ast::aArgPeriToLongPeri (double argPeri, double raan)
 
double ast::aRAANToLAN (double raan, Axes *inertialAxes, const TimePoint &timeOfAscNodePassage, Axes *bodyFixedAxes)
 
double ast::aLANToRAAN (double lan, Axes *bodyFixedAxes, const TimePoint &timeOfAscNodePassage, Axes *inertialAxes)
 
double ast::aRepeatGroundTrackSMA (int daysToRepeat, int revsToRepeat, double gm, double bodyRotRate)
 
double ast::aEccToFlat (double eccentricity)
 
double ast::aFlatToEcc (double flatFactor)
 

函数说明

◆ aApoAltToApoRad()

double ast::aApoAltToApoRad ( double apogeeAlt,
double bodyRadius )
inline

远地点高度转换为远地点半径

轨道参数与其常用缩写

轨道参数 全称 常用缩写
远地点高度 apoapsis altitude apoAlt
近地点高度 periapsis altitude periAlt
远地点半径 apoapsis radius apoRad
近地点半径 periapsis radius periRad
轨道周期 period period
平均角速度 mean motion meanMotion
偏心率 eccentricity ecc
长半轴 semi-major axis sma
近地点幅角 argument of perigee argPeri
升交点赤经 right ascension of ascending node raan
真近点角 true anomaly trueAnomaly
平近点角 mean anomaly meanAnomaly
偏近点角 eccentric anomaly eccAnomaly
纬度幅角 argument of latitude argLatAOP
真近点经度 true longitude trueLong
近地点经度 longitude of perigee longPeri
过近心点后时间 time past periapsis TPP,timePastPeri
过升交点后时间 time past ascending node TPAN,timePastAscNode
扁率 flattening factor flatFactor
引力参数 gravitational parameter gm
天体半径 body radius bodyRadius
天体自转角速度 body rotation rate bodyRotRate

peri- 近 apo- 远 ge- 地球 apsis 拱点

periapsis: 近拱点 perigee: 近地点,绕地球(Geo) 转 perihelion: 近日点,绕太阳(Helios) 转 perilune或Periselene: 近月点,绕月球(Luna/Selene) 转 perijove: 近木点,绕木星(Jove/Jupiter) 转 periastron: 近星点,绕恒星(Aster)转(通常在双星系统中使用)

参数
apogeeAlt远地点高度 [m]
bodyRadius中心天体半径 [m]
返回
远地点半径 [m]

◆ aApoAltToMeanMotion()

AST_CORE_CAPI double ast::aApoAltToMeanMotion ( double apogeeAlt,
double eccentricity,
double bodyRadius,
double gm )

远地点高度转换为平均角速度

参数
apogeeAlt远地点高度 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
gm引力参数 [m^3/s^2]
返回
平均角速度 [rad/s]

◆ aApoAltToPeriAlt()

AST_CORE_CAPI double ast::aApoAltToPeriAlt ( double apogeeAlt,
double eccentricity,
double bodyRadius )

远地点高度转换为近地点高度

参数
apogeeAlt远地点高度 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
返回
近地点高度 [m]

◆ aApoAltToPeriod()

AST_CORE_CAPI double ast::aApoAltToPeriod ( double apogeeAlt,
double eccentricity,
double bodyRadius,
double gm )

远地点高度转换为轨道周期

参数
apogeeAlt远地点高度 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
gm引力参数 [m^3/s^2]
返回
轨道周期 [s]

◆ aApoAltToPeriRad()

AST_CORE_CAPI double ast::aApoAltToPeriRad ( double apogeeAlt,
double eccentricity,
double bodyRadius )

远地点高度转换为近地点半径

参数
apogeeAlt远地点高度 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
返回
近地点半径 [m]

◆ aApoAltToSMA()

AST_CORE_CAPI double ast::aApoAltToSMA ( double apogeeAlt,
double eccentricity,
double bodyRadius )

远地点高度转换为长半轴

参数
apogeeAlt远地点高度 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
返回
长半轴 [m]

◆ aApoRadToApoAlt()

AST_CORE_CAPI double ast::aApoRadToApoAlt ( double apogeeRad,
double bodyRadius )

远地点半径转换为远地点高度

参数
apogeeRad远地点半径 [m]
bodyRadius中心天体半径 [m]
返回
远地点高度 [m]

◆ aApoRadToEcc()

double ast::aApoRadToEcc ( double apogeeRad,
double semiMajorAxis )
inline

远地点半径转换为偏心率

参数
apogeeRad远地点半径 [m]
semiMajorAxis长半轴 [m]
返回
偏心率

◆ aApoRadToMeanMotion()

AST_CORE_CAPI double ast::aApoRadToMeanMotion ( double apogeeRad,
double eccentricity,
double gm )

远地点半径转换为平均角速度

参数
apogeeRad远地点半径 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
平均角速度 [rad/s]

◆ aApoRadToPeriAlt()

AST_CORE_CAPI double ast::aApoRadToPeriAlt ( double apogeeRad,
double eccentricity,
double bodyRadius )

远地点半径转换为近地点高度

参数
apogeeRad远地点半径 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
返回
近地点高度 [m]

◆ aApoRadToPeriod()

AST_CORE_CAPI double ast::aApoRadToPeriod ( double apogeeRad,
double eccentricity,
double gm )

远地点半径转换为轨道周期

参数
apogeeRad远地点半径 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
轨道周期 [s]

◆ aApoRadToPeriRad()

AST_CORE_CAPI double ast::aApoRadToPeriRad ( double apogeeRad,
double eccentricity )

远地点半径转换为近地点半径

参数
apogeeRad远地点半径 [m]
eccentricity偏心率
返回
近地点半径 [m]

◆ aApoRadToSMA()

AST_CORE_CAPI double ast::aApoRadToSMA ( double apogeeRad,
double eccentricity )

远地点半径转换为长半轴

参数
apogeeRad远地点半径 [m]
eccentricity偏心率
返回
长半轴 [m]

◆ aArgLatToTrue()

double ast::aArgLatToTrue ( double argLat,
double argPeri )
inline

纬度幅角转换为真近点角

参数
argLat纬度幅角 [rad]
argPeri近地点幅角 [rad]
返回
真近点角 [rad]

◆ aArgPeriRate()

AST_CORE_CAPI double ast::aArgPeriRate ( double gm,
double j2,
double rb,
double a,
double ecc,
double inc )

计算轨道近点幅角变化率

参数
gm天体的引力参数
j2天体的J2项
rb天体的半径
a轨道的半长轴
ecc轨道的偏心率
inc轨道的倾角
返回
轨道近点幅角变化率

◆ aArgPeriToLongPeri()

AST_CORE_CAPI double ast::aArgPeriToLongPeri ( double argPeri,
double raan )

近地点幅角转换为近地点经度

参数
argPeri近地点幅角 [rad]
raan升交点赤经 [rad]
返回
近地点经度 [rad]

◆ aCartToEquinElem()

AST_CORE_CAPI void ast::aCartToEquinElem ( const Vector3d & pos,
const Vector3d & vel,
double gm,
EquinElem & equinElem )

直角坐标转换为春分点根数(类引用版本)

参数
pos位置矢量 [m]
vel速度矢量 [m/s]
gm引力参数 [m^3/s^2]
equinElem输出春分点根数

◆ aCartToModEquinElem()

AST_CORE_CAPI errc_t ast::aCartToModEquinElem ( const Vector3d & pos,
const Vector3d & vel,
double gm,
ModEquinElem & mee )

直角坐标转换为改进春分点轨道根数(类引用版本)

参数
pos位置矢量 [m]
vel速度矢量 [m/s]
gm引力参数 [m^3/s^2]
mee输出改进春分点轨道根数
返回
错误码,成功返回eNoError

◆ aCartToModOrbElem() [1/2]

AST_CORE_CAPI errc_t ast::aCartToModOrbElem ( const Vector3d & pos,
const Vector3d & vel,
double gm,
ModOrbElem & modOrb )

直角坐标转换为修正轨道根数(类引用版本)

参数
pos位置矢量 [m]
vel速度矢量 [m/s]
gm引力参数 [m^3/s^2]
modOrb输出修正轨道根数
返回
错误码,成功返回eNoError

◆ aCartToModOrbElem() [2/2]

A_ALWAYS_INLINE ModOrbElem ast::aCartToModOrbElem ( const Vector3d & r,
const Vector3d & v,
double gm )

直角坐标转换为修正轨道根数(类引用版本)

参数
pos位置矢量 [m]
vel速度矢量 [m/s]
gm引力参数 [m^3/s^2]
返回
修正轨道根数

◆ aCartToOrbElem()

AST_CORE_CAPI errc_t ast::aCartToOrbElem ( const Vector3d & pos,
const Vector3d & vel,
double gm,
OrbElem & elem )

直角坐标转换为经典轨道根数(类引用版本)

参数
pos位置矢量 [m]
vel速度矢量 [m/s]
gm引力参数 [m^3/s^2]
elem输出经典轨道根数
返回
错误码,成功返回eNoError

◆ aDelaunayToOrbElem()

AST_CORE_CAPI errc_t ast::aDelaunayToOrbElem ( const DelaunayElem & delaunay,
double gm,
OrbElem & elem )

修正轨道根数转换为经典轨道根数(类引用版本)

参数
delaunay修正轨道根数
gm引力参数 [m^3/s^2]
elem输出经典轨道根数

◆ aEccToFlat()

AST_CORE_CAPI double ast::aEccToFlat ( double eccentricity)

偏心率转换为扁率

参数
eccentricity偏心率
返回
扁率

◆ aEccToMean()

AST_CORE_CAPI double ast::aEccToMean ( double eccAnomaly,
double eccentricity )

偏近点角转换为平近点角

参数
eccAnomaly偏近点角 [rad]
eccentricity偏心率
返回
平近点角 [rad]

◆ aEccToTimePastAscNode()

AST_CORE_CAPI double ast::aEccToTimePastAscNode ( double eccAnomaly,
double argPeri,
double semiMajorAxis,
double eccentricity,
double gm )

偏近点角转换为过升交点后时间

参数
eccAnomaly偏近点角 [rad]
argPeri近地点幅角 [rad]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
过升交点后时间 [s]

◆ aEccToTimePastPeri()

AST_CORE_CAPI double ast::aEccToTimePastPeri ( double eccAnomaly,
double semiMajorAxis,
double eccentricity,
double gm )

偏近点角转换为过近心点后时间

参数
eccAnomaly偏近点角 [rad]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
过近心点后时间 [s]

◆ aEccToTrue()

AST_CORE_CAPI double ast::aEccToTrue ( double eccAnomaly,
double eccentricity )

偏近点角转换为真近点角

参数
eccAnomaly偏近点角 [rad]
eccentricity偏心率
返回
真近点角 [rad]

◆ aEquinElemToCart()

AST_CORE_CAPI void ast::aEquinElemToCart ( const EquinElem & equinElem,
double gm,
Vector3d & pos,
Vector3d & vel )

春分点根数转换为直角坐标(类引用版本)

参数
equinElem春分点根数
gm引力参数 [m^3/s^2]
pos输出位置矢量 [m]
vel输出速度矢量 [m/s]

◆ aEquinElemToModOrb()

AST_CORE_CAPI errc_t ast::aEquinElemToModOrb ( const EquinElem & equinElem,
ModOrbElem & modOrb )

春分点根数转换为修正轨道根数(类引用版本)

参数
equinElem春分点根数
modOrb输出修正轨道根数
返回
错误码,成功返回eNoError

◆ aFlatToEcc()

AST_CORE_CAPI double ast::aFlatToEcc ( double flatFactor)

扁率转换为偏心率

参数
flatFactor扁率
返回
偏心率

◆ aJ2Period()

AST_CORE_CAPI double ast::aJ2Period ( double gm,
double j2,
double rb,
double a,
double ecc,
double inc )

计算轨道在J2长期项影响下的轨道周期/交点周期

参数
gm天体的引力参数
j2天体的J2项
rb天体的半径
a轨道的半长轴
ecc轨道的偏心率
inc轨道的倾角
返回
J2长期项影响下的轨道周期/交点周期

◆ aLANToRAAN()

AST_CORE_CAPI double ast::aLANToRAAN ( double lan,
Axes * bodyFixedAxes,
const TimePoint & timeOfAscNodePassage,
Axes * inertialAxes )

升交点经度转换为升交点赤经

参数
lan升交点经度 [rad]
bodyFixedAxes天体固连系
timeOfAscNodePassage过升交点时刻
inertialAxes定义升交点的惯性系
返回
升交点赤经 [rad]

◆ aMeanAnomalyRate()

AST_CORE_CAPI double ast::aMeanAnomalyRate ( double gm,
double j2,
double rb,
double a,
double ecc,
double inc )

计算轨道平近点角变化率

参数
gm天体的引力参数
j2天体的J2项
rb天体的半径
a轨道的半长轴
ecc轨道的偏心率
inc轨道的倾角
返回
轨道平近点角变化率

◆ aMeanArgLatRate()

AST_CORE_CAPI double ast::aMeanArgLatRate ( double gm,
double j2,
double rb,
double a,
double ecc,
double inc )

计算维度幅角平均变化率(平近点角变化率 + 近点幅角变化率)

参数
gm天体的引力参数
j2天体的J2项
rb天体的半径
a轨道的半长轴
ecc轨道的偏心率
inc轨道的倾角
返回
轨道平近点幅角变化率

◆ aMeanMotionToApoAlt()

AST_CORE_CAPI double ast::aMeanMotionToApoAlt ( double meanMotion,
double eccentricity,
double bodyRadius,
double gm )

平均角速度转换为远地点高度

参数
meanMotion平均角速度 [rad/s]
eccentricity偏心率
bodyRadius中心天体半径 [m]
gm引力参数 [m^3/s^2]
返回
远地点高度 [m]

◆ aMeanMotionToApoRad()

AST_CORE_CAPI double ast::aMeanMotionToApoRad ( double meanMotion,
double eccentricity,
double gm )

平均角速度转换为远地点半径

参数
meanMotion平均角速度 [rad/s]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
远地点半径 [m]

◆ aMeanMotionToPeriAlt()

AST_CORE_CAPI double ast::aMeanMotionToPeriAlt ( double meanMotion,
double eccentricity,
double bodyRadius,
double gm )

平均角速度转换为近地点高度

参数
meanMotion平均角速度 [rad/s]
eccentricity偏心率
bodyRadius中心天体半径 [m]
gm引力参数 [m^3/s^2]
返回
近地点高度 [m]

◆ aMeanMotionToPeriod()

AST_CORE_CAPI double ast::aMeanMotionToPeriod ( double meanMotn)

平均角速度转换为轨道周期

参数
meanMotn平均角速度 [rad/s]
返回
轨道周期 [s]

◆ aMeanMotionToPeriRad()

AST_CORE_CAPI double ast::aMeanMotionToPeriRad ( double meanMotion,
double eccentricity,
double gm )

平均角速度转换为近地点半径

参数
meanMotion平均角速度 [rad/s]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
近地点半径 [m]

◆ aMeanMotionToSMA()

AST_CORE_CAPI double ast::aMeanMotionToSMA ( double meanMotn,
double gm )

平均角速度转换为长半轴

参数
meanMotn平均角速度 [rad/s]
gm引力参数 [m^3/s^2]
返回
长半轴 [m]

◆ aMeanToEcc()

AST_CORE_CAPI double ast::aMeanToEcc ( double meanAnomaly,
double eccentricity,
double eps = 1e-14,
int maxIter = 100 )

平近点角转换为偏近点角(使用牛顿迭代法)

参数
meanAnomaly平近点角 [rad]
eccentricity偏心率
eps迭代精度,默认1e-14
maxIter最大迭代次数,默认100
返回
偏近点角 [rad]

◆ aMeanToTimePastAscNode()

AST_CORE_CAPI double ast::aMeanToTimePastAscNode ( double meanAnomaly,
double argPeri,
double semiMajorAxis,
double eccentricity,
double gm )

平近点角转换为过升交点后时间

参数
meanAnomaly平近点角 [rad]
argPeri近地点幅角 [rad]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
过升交点后时间 [s]

◆ aMeanToTimePastPeri()

AST_CORE_CAPI double ast::aMeanToTimePastPeri ( double meanAnomaly,
double semiMajorAxis,
double gm )

平近点角转换为过近心点后时间

参数
meanAnomaly平近点角 [rad]
semiMajorAxis长半轴 [m]
gm引力参数 [m^3/s^2]
返回
过近心点后时间 [s]

◆ aMeanToTrue()

AST_CORE_CAPI double ast::aMeanToTrue ( double meanAnomaly,
double eccentricity,
double eps = 1e-14,
int maxIter = 100 )

平近点角转换为真近点角

参数
meanAnomaly平近点角 [rad]
eccentricity偏心率
eps迭代精度,默认1e-14
maxIter最大迭代次数,默认100
返回
真近点角 [rad]

◆ aModEquinElemToCart()

AST_CORE_CAPI void ast::aModEquinElemToCart ( const ModEquinElem & mee,
double gm,
Vector3d & pos,
Vector3d & vel )

改进春分点轨道根数转换为直角坐标(类引用版本)

参数
mee改进春分点轨道根数
gm引力参数 [m^3/s^2]
pos输出位置矢量 [m]
vel输出速度矢量 [m/s]

◆ aModEquinElemToOrbElem()

AST_CORE_CAPI errc_t ast::aModEquinElemToOrbElem ( const ModEquinElem & mee,
OrbElem & elem )

改进春分点轨道根数转换为经典轨道根数(类引用版本)

参数
mee改进春分点轨道根数
elem输出经典轨道根数
返回
错误码,成功返回eNoError

◆ aModOrbElemToCart()

AST_CORE_CAPI errc_t ast::aModOrbElemToCart ( const ModOrbElem & modOrb,
double gm,
Vector3d & pos,
Vector3d & vel )

修正轨道根数转换为直角坐标(类引用版本)

参数
modOrb修正轨道根数
gm引力参数 [m^3/s^2]
pos输出位置矢量 [m]
vel输出速度矢量 [m/s]

◆ aModOrbToEquinElem()

AST_CORE_CAPI errc_t ast::aModOrbToEquinElem ( const ModOrbElem & modOrb,
EquinElem & equinElem )

修正轨道根数转换为春分点根数(类引用版本)

参数
modOrb修正轨道根数
equinElem输出春分点根数
返回
错误码,成功返回eNoError

◆ aOrbElemToCart()

AST_CORE_CAPI errc_t ast::aOrbElemToCart ( const OrbElem & elem,
double gm,
Vector3d & pos,
Vector3d & vel )

经典轨道根数转换为直角坐标(类引用版本)

参数
elem经典轨道根数
gm引力参数 [m^3/s^2]
pos输出位置矢量 [m]
vel输出速度矢量 [m/s]
返回
错误码,成功返回eNoError

◆ aOrbElemToDelaunay()

AST_CORE_CAPI errc_t ast::aOrbElemToDelaunay ( const OrbElem & elem,
double gm,
DelaunayElem & delaunay )

经典轨道根数转换为修正轨道根数(类引用版本)

参数
elem经典轨道根数
gm引力参数 [m^3/s^2]
delaunay输出修正轨道根数

◆ aOrbElemToModEquinElem()

AST_CORE_CAPI errc_t ast::aOrbElemToModEquinElem ( const OrbElem & elem,
ModEquinElem & mee )

经典轨道根数转换为改进春分点轨道根数(类引用版本)

参数
elem经典轨道根数
mee输出改进春分点轨道根数
返回
错误码,成功返回eNoError

◆ aPeriAltApoAltToEcc()

AST_CORE_CAPI double ast::aPeriAltApoAltToEcc ( double perigeeAlt,
double apogeeAlt,
double bodyRadius )

根据近地点和远地点高度计算偏心率

参数
perigeeAlt近地点高度 [m]
apogeeAlt远地点高度 [m]
bodyRadius中心天体半径 [m]
返回
偏心率

◆ aPeriAltToApoAlt()

AST_CORE_CAPI double ast::aPeriAltToApoAlt ( double perigeeAlt,
double eccentricity,
double bodyRadius )

近地点高度转换为远地点高度

参数
perigeeAlt近地点高度 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
返回
远地点高度 [m]

◆ aPeriAltToApoRad()

AST_CORE_CAPI double ast::aPeriAltToApoRad ( double perigeeAlt,
double eccentricity,
double bodyRadius )

近地点高度转换为远地点半径

参数
perigeeAlt近地点高度 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
返回
远地点半径 [m]

◆ aPeriAltToMeanMotion()

AST_CORE_CAPI double ast::aPeriAltToMeanMotion ( double perigeeAlt,
double eccentricity,
double bodyRadius,
double gm )

近地点高度转换为平均角速度

参数
perigeeAlt近地点高度 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
gm引力参数 [m^3/s^2]
返回
平均角速度 [rad/s]

◆ aPeriAltToPeriod()

AST_CORE_CAPI double ast::aPeriAltToPeriod ( double perigeeAlt,
double eccentricity,
double bodyRadius,
double gm )

近地点高度转换为轨道周期

参数
perigeeAlt近地点高度 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
gm引力参数 [m^3/s^2]
返回
轨道周期 [s]

◆ aPeriAltToPeriRad()

AST_CORE_CAPI double ast::aPeriAltToPeriRad ( double perigeeAlt,
double bodyRadius )

近地点高度转换为近地点半径

参数
perigeeAlt近地点高度 [m]
bodyRadius中心天体半径 [m]
返回
近地点半径 [m]

◆ aPeriAltToSMA()

AST_CORE_CAPI double ast::aPeriAltToSMA ( double perigeeAlt,
double eccentricity,
double bodyRadius )

近地点高度转换为长半轴

参数
perigeeAlt近地点高度 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
返回
长半轴 [m]

◆ aPeriodToApoAlt()

AST_CORE_CAPI double ast::aPeriodToApoAlt ( double period,
double eccentricity,
double bodyRadius,
double gm )

轨道周期转换为远地点高度

参数
period轨道周期 [s]
eccentricity偏心率
bodyRadius中心天体半径 [m]
gm引力参数 [m^3/s^2]
返回
远地点高度 [m]

◆ aPeriodToApoRad()

AST_CORE_CAPI double ast::aPeriodToApoRad ( double period,
double eccentricity,
double gm )

轨道周期转换为远地点半径

参数
period轨道周期 [s]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
远地点半径 [m]

◆ aPeriodToMeanMotion()

AST_CORE_CAPI double ast::aPeriodToMeanMotion ( double period)

轨道周期转换为平均角速度

参数
period轨道周期 [s]
返回
平均角速度 [rad/s]

◆ aPeriodToPeriAlt()

AST_CORE_CAPI double ast::aPeriodToPeriAlt ( double period,
double eccentricity,
double bodyRadius,
double gm )

轨道周期转换为近地点高度

参数
period轨道周期 [s]
eccentricity偏心率
bodyRadius中心天体半径 [m]
gm引力参数 [m^3/s^2]
返回
近地点高度 [m]

◆ aPeriodToPeriRad()

AST_CORE_CAPI double ast::aPeriodToPeriRad ( double period,
double eccentricity,
double gm )

轨道周期转换为近地点半径

参数
period轨道周期 [s]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
近地点半径 [m]

◆ aPeriodToSMA()

AST_CORE_CAPI double ast::aPeriodToSMA ( double period,
double gm )

轨道周期转换为长半轴

参数
period轨道周期 [s]
gm引力参数 [m^3/s^2]
返回
长半轴 [m]

◆ aPeriRadApoRadToEcc()

AST_CORE_CAPI double ast::aPeriRadApoRadToEcc ( double perigeeRad,
double apogeeRad )

根据近地点和远地点半径计算偏心率

参数
perigeeRad近地点半径 [m]
apogeeRad远地点半径 [m]
返回
偏心率

◆ aPeriRadToApoAlt()

AST_CORE_CAPI double ast::aPeriRadToApoAlt ( double perigeeRad,
double eccentricity,
double bodyRadius )

近地点半径转换为远地点高度

参数
perigeeRad近地点半径 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
返回
远地点高度 [m]

◆ aPeriRadToApoRad()

AST_CORE_CAPI double ast::aPeriRadToApoRad ( double perigeeRad,
double eccentricity )

近地点半径转换为远地点半径

参数
perigeeRad近地点半径 [m]
eccentricity偏心率
返回
远地点半径 [m]

◆ aPeriRadToEcc()

double ast::aPeriRadToEcc ( double perigeeRad,
double semiMajorAxis )
inline

近地点半径转换为偏心率

参数
perigeeRad近地点半径 [m]
semiMajorAxis长半轴 [m]
返回
偏心率

◆ aPeriRadToMeanMotion()

AST_CORE_CAPI double ast::aPeriRadToMeanMotion ( double perigeeRad,
double eccentricity,
double gm )

近地点半径转换为平均角速度

参数
perigeeRad近地点半径 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
平均角速度 [rad/s]

◆ aPeriRadToPeriAlt()

AST_CORE_CAPI double ast::aPeriRadToPeriAlt ( double perigeeRad,
double bodyRadius )

近地点半径转换为近地点高度

参数
perigeeRad近地点半径 [m]
bodyRadius中心天体半径 [m]
返回
近地点高度 [m]

◆ aPeriRadToPeriod()

AST_CORE_CAPI double ast::aPeriRadToPeriod ( double perigeeRad,
double eccentricity,
double gm )

近地点半径转换为轨道周期

参数
perigeeRad近地点半径 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
轨道周期 [s]

◆ aPeriRadToSMA()

AST_CORE_CAPI double ast::aPeriRadToSMA ( double perigeeRad,
double eccentricity )

近地点半径转换为长半轴

参数
perigeeRad近地点半径 [m]
eccentricity偏心率
返回
长半轴 [m]

◆ aRAANRate()

AST_CORE_CAPI double ast::aRAANRate ( double gm,
double j2,
double rb,
double a,
double ecc,
double inc )

计算轨道升交点赤经漂移速率

参数
gm天体的引力参数
j2天体的J2项
rb天体的半径
a轨道的半长轴
ecc轨道的偏心率
inc轨道的倾角
返回
轨道升交点赤经漂移速率

◆ aRAANToLAN()

AST_CORE_CAPI double ast::aRAANToLAN ( double raan,
Axes * inertialAxes,
const TimePoint & timeOfAscNodePassage,
Axes * bodyFixedAxes )

升交点赤经转换为升交点经度

参数
raan升交点赤经 [rad]
inertialAxes定义升交点的惯性系
timeOfAscNodePassage过升交点时刻
bodyFixedAxes天体固连系
返回
升交点经度 [rad]

◆ aRepeatGroundTrackSMA()

AST_CORE_CAPI double ast::aRepeatGroundTrackSMA ( int daysToRepeat,
int revsToRepeat,
double gm,
double bodyRotRate )

计算地面轨迹重复所需的轨道长半轴

参数
daysToRepeat重复天数
revsToRepeat重复圈数
gm引力参数 [m^3/s^2]
bodyRotRate中心天体自转角速度 [rad/s]
返回
满足地面轨迹重复条件的轨道长半轴 [m]

◆ aSMAToApoAlt()

AST_CORE_CAPI double ast::aSMAToApoAlt ( double semiMajorAxis,
double eccentricity,
double bodyRadius )

长半轴转换为远地点高度

参数
semiMajorAxis长半轴 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
返回
远地点高度 [m]

◆ aSMAToApoRad()

AST_CORE_CAPI double ast::aSMAToApoRad ( double semiMajorAxis,
double eccentricity )

长半轴转换为远地点半径

参数
semiMajorAxis长半轴 [m]
eccentricity偏心率
返回
远地点半径 [m]

◆ aSMAToMeanMotion()

AST_CORE_CAPI double ast::aSMAToMeanMotion ( double semiMajorAxis,
double gm )

长半轴转换为平均角速度

参数
semiMajorAxis长半轴 [m]
gm引力参数 [m^3/s^2]
返回
平均角速度 [rad/s]

◆ aSMAToPeriAlt()

AST_CORE_CAPI double ast::aSMAToPeriAlt ( double semiMajorAxis,
double eccentricity,
double bodyRadius )

长半轴转换为近地点高度

参数
semiMajorAxis长半轴 [m]
eccentricity偏心率
bodyRadius中心天体半径 [m]
返回
近地点高度 [m]

◆ aSMAToPeriod()

AST_CORE_CAPI double ast::aSMAToPeriod ( double semiMajorAxis,
double gm )

长半轴转换为轨道周期

参数
semiMajorAxis长半轴 [m]
gm引力参数 [m^3/s^2]
返回
轨道周期 [s]

◆ aSMAToPeriRad()

double ast::aSMAToPeriRad ( double semiMajorAxis,
double eccentricity )
inline

长半轴转换为近地点半径

参数
semiMajorAxis长半轴 [m]
eccentricity偏心率
返回
近地点半径 [m]

◆ aSMAToSMinAx()

AST_CORE_CAPI double ast::aSMAToSMinAx ( double semiMajorAxis,
double eccentricity )

长半轴转换为短半轴

参数
semiMajorAxis长半轴 [m]
eccentricity偏心率
返回
短半轴 [m]

◆ aSMAToSParam()

AST_CORE_CAPI double ast::aSMAToSParam ( double semiMajorAxis,
double eccentricity )

长半轴转换为半通径参数

参数
semiMajorAxis长半轴 [m]
eccentricity偏心率
返回
半通径 [m]

◆ aSMinAxToSMA()

AST_CORE_CAPI double ast::aSMinAxToSMA ( double semiminorAxis,
double eccentricity )

短半轴转换为长半轴

参数
semiminorAxis短半轴 [m]
eccentricity偏心率
返回
长半轴 [m]

◆ aSunSynchronousInclination() [1/2]

A_ALWAYS_INLINE double ast::aSunSynchronousInclination ( double gm,
double j2,
double rb,
double bodyMeanMotion,
double a,
double ecc )

计算太阳同步轨道的倾角

参数
[in]gm天体引力系数
[in]j2天体J2项
[in]rb天体半径
[in]bodyMeanMotion天体平均公转角速度
[in]a轨道半长轴
[in]ecc轨道偏心率
返回
double 轨道倾角,如果计算失败,则为NaN

◆ aSunSynchronousInclination() [2/2]

errc_t ast::aSunSynchronousInclination ( double gm,
double j2,
double rb,
double bodyMeanMotion,
double a,
double ecc,
double & inc )

计算太阳同步轨道的倾角

参数
[in]gm天体引力系数
[in]j2天体J2项
[in]rb天体半径
[in]bodyMeanMotion天体平均公转角速度
[in]a轨道半长轴
[in]ecc轨道偏心率
[out]inc轨道倾角,如果计算失败,则为NaN
返回
errc_t 错误码

◆ aSunSynchronousSemiMajorAxis() [1/2]

A_ALWAYS_INLINE double ast::aSunSynchronousSemiMajorAxis ( double gm,
double j2,
double rb,
double bodyMeanMotion,
double inc,
double ecc )

计算太阳同步轨道的半长轴

参数
[in]gm天体引力系数
[in]j2天体J2项
[in]rb天体半径
[in]bodyMeanMotion天体平均公转角速度
[in]inc轨道倾角
[in]ecc轨道偏心率
返回
double 轨道半长轴,如果计算失败,则为NaN

◆ aSunSynchronousSemiMajorAxis() [2/2]

errc_t ast::aSunSynchronousSemiMajorAxis ( double gm,
double j2,
double rb,
double bodyMeanMotion,
double inc,
double ecc,
double & semiMajorAxis )

计算太阳同步轨道的半长轴

参数
[in]gm天体引力系数
[in]j2天体J2项
[in]rb天体半径
[in]bodyMeanMotion天体平均公转角速度
[in]inc轨道倾角
[in]ecc轨道偏心率
[out]semiMajorAxis轨道半长轴,如果计算失败,则为NaN
返回
errc_t 错误码

◆ aTimePastAscNodeToEcc()

AST_CORE_CAPI double ast::aTimePastAscNodeToEcc ( double TimePastAscNode,
double argPeri,
double semiMajorAxis,
double eccentricity,
double gm,
double eps = 1e-14,
int maxIter = 100 )

过升交点后时间转换为偏近点角

参数
TimePastAscNode过升交点后时间 [s]
argPeri近地点幅角 [rad]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
eps迭代精度,默认1e-14
maxIter最大迭代次数,默认100
返回
偏近点角 [rad]

◆ aTimePastAscNodeToMean()

AST_CORE_CAPI double ast::aTimePastAscNodeToMean ( double TimePastAscNode,
double argPeri,
double semiMajorAxis,
double eccentricity,
double gm )

过升交点后时间转换为平近点角

参数
TimePastAscNode过升交点后时间 [s]
argPeri近地点幅角 [rad]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
平近点角 [rad]

◆ aTimePastAscNodeToTimePastPeri()

AST_CORE_CAPI double ast::aTimePastAscNodeToTimePastPeri ( double TimePastAscNode,
double argPeri,
double semiMajorAxis,
double eccentricity,
double gm )

过升交点后时间转换为过近心点后时间

参数
TimePastAscNode过升交点后时间 [s]
argPeri近地点幅角 [rad]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
过近心点后时间 [s]

◆ aTimePastAscNodeToTrue()

AST_CORE_CAPI double ast::aTimePastAscNodeToTrue ( double TimePastAscNode,
double argPeri,
double semiMajorAxis,
double eccentricity,
double gm,
double eps = 1e-14,
int maxIter = 100 )

过升交点后时间转换为真近点角

参数
TimePastAscNode过升交点后时间 [s]
argPeri近地点幅角 [rad]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
eps迭代精度,默认1e-14
maxIter最大迭代次数,默认100
返回
真近点角 [rad]

◆ aTimePastPeriToEcc()

AST_CORE_CAPI double ast::aTimePastPeriToEcc ( double TimePastPeri,
double semiMajorAxis,
double eccentricity,
double gm,
double eps = 1e-14,
int maxIter = 100 )

过近心点后时间转换为偏近点角

参数
TimePastPeri过近心点后时间 [s]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
eps迭代精度,默认1e-14
maxIter最大迭代次数,默认100
返回
偏近点角 [rad]

◆ aTimePastPeriToMean()

AST_CORE_CAPI double ast::aTimePastPeriToMean ( double TimePastPeri,
double semiMajorAxis,
double gm )

过近心点后时间转换为平近点角

参数
TimePastPeri过近心点后时间 [s]
semiMajorAxis长半轴 [m]
gm引力参数 [m^3/s^2]
返回
平近点角 [rad]

◆ aTimePastPeriToTimePastAscNode()

AST_CORE_CAPI double ast::aTimePastPeriToTimePastAscNode ( double TimePastPeri,
double argPeri,
double semiMajorAxis,
double eccentricity,
double gm )

过近心点后时间转换为过升交点后时间

参数
TimePastPeri过近心点后时间 [s]
argPeri近地点幅角 [rad]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
过升交点后时间 [s]

◆ aTimePastPeriToTrue()

AST_CORE_CAPI double ast::aTimePastPeriToTrue ( double TimePastPeri,
double semiMajorAxis,
double eccentricity,
double gm,
double eps = 1e-14,
int maxIter = 100 )

过近心点后时间转换为真近点角

参数
TimePastPeri过近心点后时间 [s]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
eps迭代精度,默认1e-14
maxIter最大迭代次数,默认100
返回
真近点角 [rad]

◆ aTrueToArgLat()

double ast::aTrueToArgLat ( double trueAnomaly,
double argPeri )
inline

真近点角转换为纬度幅角

参数
trueAnomaly真近点角 [rad]
argPeri近地点幅角 [rad]
返回
纬度幅角 [rad]

◆ aTrueToEcc()

AST_CORE_CAPI double ast::aTrueToEcc ( double trueAnomaly,
double eccentricity )

真近点角转换为偏近点角

参数
trueAnomaly真近点角 [rad]
eccentricity偏心率
返回
偏近点角 [rad]

◆ aTrueToMean()

AST_CORE_CAPI double ast::aTrueToMean ( double trueAnomaly,
double eccentricity )

真近点角转换为平近点角

参数
trueAnomaly真近点角 [rad]
eccentricity偏心率
返回
平近点角 [rad]

◆ aTrueToTimePastAscNode()

AST_CORE_CAPI double ast::aTrueToTimePastAscNode ( double trueAnomaly,
double argPeri,
double semiMajorAxis,
double eccentricity,
double gm )

真近点角转换为过升交点后时间

参数
trueAnomaly真近点角 [rad]
argPeri近地点幅角 [rad]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
过升交点后时间 [s]

◆ aTrueToTimePastPeri()

AST_CORE_CAPI double ast::aTrueToTimePastPeri ( double trueAnomaly,
double semiMajorAxis,
double eccentricity,
double gm )

真近点角转换为过近心点后时间

参数
trueAnomaly真近点角 [rad]
semiMajorAxis长半轴 [m]
eccentricity偏心率
gm引力参数 [m^3/s^2]
返回
过近心点后时间 [s]

◆ aTrueToTrueLong()

AST_CORE_CAPI double ast::aTrueToTrueLong ( double trueAnomaly,
double argPeri,
double raan )

真近点角转换为真近点经度

参数
trueAnomaly真近点角 [rad]
argPeri近地点幅角 [rad]
raan升交点赤经 [rad]
返回
真近点经度 [rad]

◆ coe2dela()

AST_CORE_CAPI errc_t ast::coe2dela ( const double * coe,
double gm,
double * dela )

经典轨道根数转换为德洛奈根数

参数
coe经典轨道根数 [长半轴, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
gm引力参数 [m^3/s^2]
del输出德洛奈根数, see DelaunayElem

◆ coe2ee()

AST_CORE_CAPI void ast::coe2ee ( const double * coe,
double * ee )

经典轨道根数转换为春分点根数

参数
coe经典轨道根数 [长半轴, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
ee输出春分点根数 [长半轴, h, k, p, q, 平经度]

◆ coe2mee()

AST_CORE_CAPI errc_t ast::coe2mee ( const double * coe,
double * mee )

经典轨道根数转换为改进春分点轨道根数

参数
coe经典轨道根数 [长半轴, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
mee输出改进春分点轨道根数 [半通径, f, g, h, k, L]
返回
错误码,成功返回eNoError

◆ coe2moe()

AST_CORE_CAPI void ast::coe2moe ( const double * coe,
double * moe )

经典轨道根数转换为修正轨道根数

参数
coe经典轨道根数 [长半轴, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
moe输出修正轨道根数 [近拱点半径, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]

◆ coe2rv()

AST_CORE_CAPI errc_t ast::coe2rv ( const double * coe,
double gm,
double * pos,
double * vel )

经典轨道根数转换为直角坐标

参数
coe经典轨道根数 [长半轴, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
gm引力参数 [m^3/s^2]
pos输出位置矢量 [m]
vel输出速度矢量 [m/s]
返回
错误码,成功返回eNoError

◆ dela2coe()

AST_CORE_CAPI errc_t ast::dela2coe ( const double * dela,
double gm,
double * coe )

德洛奈根数转换为经典轨道根数

参数
del德洛奈根数, see DelaunayElem
coe输出经典轨道根数 [长半轴, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]

◆ ee2coe()

AST_CORE_CAPI void ast::ee2coe ( const double * ee,
double * coe )

春分点根数转换为经典轨道根数

参数
ee春分点根数 [长半轴, h, k, p, q, 平经度]
coe输出经典轨道根数 [长半轴, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]

◆ ee2mee()

AST_CORE_CAPI void ast::ee2mee ( const double * ee,
double * mee )

春分点根数转换为改进春分点轨道根数

参数
ee春分点根数 [长半轴, h, k, p, q, 平经度]
mee输出改进春分点轨道根数 [半通径, f, g, h, k, L]

◆ ee2moe()

AST_CORE_CAPI errc_t ast::ee2moe ( const double * ee,
double * moe )

春分点根数转换为修正轨道根数

参数
ee春分点根数 [长半轴, h, k, p, q, 平经度]
moe输出修正轨道根数 [近拱点半径, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
返回
错误码,成功返回eNoError

◆ ee2rv()

AST_CORE_CAPI void ast::ee2rv ( const double * ee,
double gm,
double * pos,
double * vel )

春分点根数转换为直角坐标

参数
ee春分点根数 [长半轴, h, k, p, q, 平经度]
gm引力参数 [m^3/s^2]
pos输出位置矢量 [m]
vel输出速度矢量 [m/s]

◆ mee2coe()

AST_CORE_CAPI errc_t ast::mee2coe ( const double * mee,
double * coe )

改进春分点轨道根数转换为经典轨道根数

参数
mee改进春分点轨道根数 [半通径, f, g, h, k, L]
coe输出经典轨道根数 [长半轴, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
返回
错误码,成功返回eNoError

◆ mee2ee()

AST_CORE_CAPI void ast::mee2ee ( const double * mee,
double * ee )

改进春分点轨道根数转换为春分点根数

参数
mee改进春分点轨道根数 [半通径, f, g, h, k, L]
ee输出春分点根数 [长半轴, h, k, p, q, 平经度]

◆ mee2moe()

AST_CORE_CAPI void ast::mee2moe ( const double * mee,
double * moe )

改进春分点轨道根数转换为修正轨道根数

参数
mee改进春分点轨道根数 [半通径, f, g, h, k, L]
moe输出修正轨道根数 [近拱点半径, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]

◆ mee2rv()

AST_CORE_CAPI void ast::mee2rv ( const double * mee,
double gm,
double * pos,
double * vel )

改进春分点轨道根数转换为直角坐标

参数
mee改进春分点轨道根数 [半通径, f, g, h, k, L]
gm引力参数 [m^3/s^2]
pos输出位置矢量 [m]
vel输出速度矢量 [m/s]

◆ moe2coe()

AST_CORE_CAPI errc_t ast::moe2coe ( const double * moe,
double * coe )

修正轨道根数转换为经典轨道根数

参数
moe修正轨道根数 [近拱点半径, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
coe输出经典轨道根数 [长半轴, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]

◆ moe2ee()

AST_CORE_CAPI errc_t ast::moe2ee ( const double * moe,
double * ee )

修正轨道根数转换为春分点根数

参数
moe修正轨道根数 [近拱点半径, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
ee输出春分点根数 [长半轴, h, k, p, q, 平经度]
返回
错误码,成功返回eNoError

◆ moe2mee()

AST_CORE_CAPI void ast::moe2mee ( const double * moe,
double * mee )

修正轨道根数转换为改进春分点轨道根数

参数
moe修正轨道根数 [近拱点半径, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
mee输出改进春分点轨道根数 [半通径, f, g, h, k, L]

◆ moe2rv()

AST_CORE_CAPI errc_t ast::moe2rv ( const double * moe,
double gm,
double * pos,
double * vel )

修正轨道根数转换为直角坐标

参数
moe修正轨道根数 [近拱点半径, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
gm引力参数 [m^3/s^2]
pos输出位置矢量 [m]
vel输出速度矢量 [m/s]

◆ rv2coe()

AST_CORE_CAPI errc_t ast::rv2coe ( const double * pos,
const double * vel,
double gm,
double * coe )

直角坐标转换为经典轨道根数

参数
pos位置矢量 [m]
vel速度矢量 [m/s]
gm引力参数 [m^3/s^2]
coe输出经典轨道根数 [长半轴, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
返回
错误码,成功返回eNoError

◆ rv2ee()

AST_CORE_CAPI void ast::rv2ee ( const double * pos,
const double * vel,
double gm,
double * ee )

直角坐标转换为春分点根数

参数
pos位置矢量 [m]
vel速度矢量 [m/s]
gm引力参数 [m^3/s^2]
ee输出春分点根数 [长半轴, h, k, p, q, 平经度]

◆ rv2mee()

AST_CORE_CAPI errc_t ast::rv2mee ( const double * pos,
const double * vel,
double gm,
double * mee )

直角坐标转换为改进春分点轨道根数

参数
pos位置矢量 [m]
vel速度矢量 [m/s]
gm引力参数 [m^3/s^2]
mee输出改进春分点轨道根数 [半通径, f, g, h, k, L]
返回
错误码,成功返回eNoError

◆ rv2moe()

AST_CORE_CAPI errc_t ast::rv2moe ( const double * pos,
const double * vel,
double gm,
double * moe )

直角坐标转换为修正轨道根数

参数
pos位置矢量 [m]
vel速度矢量 [m/s]
gm引力参数 [m^3/s^2]
moe输出修正轨道根数 [近拱点半径, 偏心率, 轨道倾角, 升交点赤经, 近拱点角, 真近点角]
返回
错误码,成功返回eNoError