🛰️航天仿真算法库 SpaceAST 0.0.1
载入中...
搜索中...
未找到
OrbitElement.hpp
浏览该文件的文档.
1
20
21#pragma once
22
23#include "AstGlobal.h"
24#include "AstCore/Vector.hpp"
25#include "AstCore/OrbitParam.hpp"
26#include "AstCore/TimePoint.hpp"
27#include "AstUtil/Constants.h"
28#include <string>
29
30AST_NAMESPACE_BEGIN
31
46{
47public:
48 Vector3d pos_;
49 Vector3d vel_;
50public:
53 static CartState Zero() {return CartState{Vector3d::Zero(), Vector3d::Zero()};}
54
55 A_DEF_POD_ITERABLE(double)
56
57
59 const Vector3d& pos() const {return pos_;}
60 Vector3d& pos() {return pos_;}
61
64 const Vector3d& vel() const {return vel_;}
65 Vector3d& vel() {return vel_;}
66
67 double& x() {return pos_.x();}
68 double& y() {return pos_.y();}
69 double& z() {return pos_.z();}
70 double& vx() {return vel_.x();}
71 double& vy() {return vel_.y();}
72 double& vz() {return vel_.z();}
73
74 double x() const {return pos_.x();}
75 double y() const {return pos_.y();}
76 double z() const {return pos_.z();}
77 double vx() const {return vel_.x();}
78 double vy() const {return vel_.y();}
79 double vz() const {return vel_.z();}
80
82 AST_CORE_API
83 std::string toString() const;
84};
85
86
89{
90public:
91 double a_;
92 double e_;
93 double i_;
94 double raan_;
95 double argper_;
96 double trueA_;
97public:
99 double getMeanMotion(double gm) const {return aSMAToMeanMotion(getSMA(), gm);}
100
101 double getSMA() const {return a_;}
102
103 double getA() const {return a_;}
104
105 double getE() const {return e_;}
106
107 double getI() const {return i_;}
108
109 double getRAAN() const {return raan_;}
110
111 double getArgPer() const {return argper_;}
112
113 double getTrueA() const {return trueA_;}
114
116 AST_CORE_API
117 std::string toString() const;
118public:
119 A_DEF_POD_ITERABLE(double)
120 AST_DEF_ACCESS_METHOD(double, a)
121 AST_DEF_ACCESS_METHOD(double, e)
122 AST_DEF_ACCESS_METHOD(double, i)
123 AST_DEF_ACCESS_METHOD(double, raan)
124 AST_DEF_ACCESS_METHOD(double, argper)
125 AST_DEF_ACCESS_METHOD(double, trueA)
126};
127
128
131{
132public:
133 double rp_;
134 double e_;
135 double i_;
136 double raan_;
137 double argper_;
138 double trueA_;
139public:
141 double getSMA() const {return rp_ / (1 - e_);}
142
144 double getMeanMotion(double gm) const {return aPeriRadToMeanMotion(getPeriRad(), getEcc(), gm);}
145
147 double getPeriod(double gm) const {return kTwoPI / getMeanMotion(gm);}
148
150 double getApoRad() const{return aPeriRadToApoRad(getPeriRad(), getEcc());}
151
153 double getApoAlt(double bodyRadius) const{ return aPeriRadToApoAlt(getPeriRad(), getEcc(), bodyRadius);}
154
156 double getPeriRad() const {return rp_;}
157
159 double getPeriAlt(double bodyRadius) const{ return getPeriRad() - bodyRadius;}
160
162 double getEcc() const {return e_;}
163
165 double getInc() const {return i_;}
166
168 double getRAAN() const {return raan_;}
169
171 double getArgPeri() const {return argper_;}
172
174 double getTrueAnomaly() const {return trueA_;}
175
177 double getMeanAnomaly() const {return aTrueToMean(getTrueAnomaly(), getEcc());}
178
180 double getEccAnomaly() const {return aTrueToEcc(getTrueAnomaly(), getEcc());}
181
183 double getArgLat() const {return aTrueToArgLat(getTrueAnomaly(), getArgPeri());}
184
186 double getTimePastPeri(double gm) const{return aTrueToTimePastPeri(getTrueAnomaly(), getSMA(), getEcc(), gm);}
187
189 double getTimePastAscNode(double gm) const{return aTrueToTimePastAscNode(getTrueAnomaly(), getArgPeri(), getSMA(), getEcc(), gm);}
190
192 TimePoint getTimeOfAscNodePassage(const TimePoint& stateEpoch, double gm) const
193 {
194 return stateEpoch - getTimePastAscNode(gm);
195 }
196
198 double getLAN(Axes* inertialAxes, const TimePoint& timeOfAscNodePassage, Axes* bodyFixedAxes) const
199 {
200 return aRAANToLAN(getRAAN(), inertialAxes, timeOfAscNodePassage, bodyFixedAxes);
201 }
202
203public:
205 double getA() const {return getSMA();}
206
208 double getP() const {return rp_ * (1 + e_);}
209
210 double getE() const {return getEcc();}
211
212 double getI() const {return getInc();}
213
214 double getArgPer() const {return getArgPeri();}
215
216 double getTrueA() const {return getTrueAnomaly();}
217
219 AST_CORE_API
220 std::string toString() const;
221public:
222 A_DEF_POD_ITERABLE(double)
223 AST_DEF_ACCESS_METHOD(double, rp)
224 AST_DEF_ACCESS_METHOD(double, e)
225 AST_DEF_ACCESS_METHOD(double, i)
226 AST_DEF_ACCESS_METHOD(double, raan)
227 AST_DEF_ACCESS_METHOD(double, argper)
228 AST_DEF_ACCESS_METHOD(double, trueA)
229};
230
231
234{
235public:
236 double a_;
237 double h_;
238 double k_;
239 double p_;
240 double q_;
241 double lambda_;
242public:
243 A_DEF_POD_ITERABLE(double)
244 AST_DEF_ACCESS_METHOD(double, a)
245 AST_DEF_ACCESS_METHOD(double, h)
246 AST_DEF_ACCESS_METHOD(double, k)
247 AST_DEF_ACCESS_METHOD(double, p)
248 AST_DEF_ACCESS_METHOD(double, q)
249 AST_DEF_ACCESS_METHOD(double, lambda)
250};
251
254{
255public:
256 double p_;
257 double f_;
258 double g_;
259 double h_;
260 double k_;
261 double L_;
262public:
263 A_DEF_POD_ITERABLE(double)
264 AST_DEF_ACCESS_METHOD(double, p)
265 AST_DEF_ACCESS_METHOD(double, f)
266 AST_DEF_ACCESS_METHOD(double, g)
267 AST_DEF_ACCESS_METHOD(double, h)
268 AST_DEF_ACCESS_METHOD(double, k)
269 AST_DEF_ACCESS_METHOD(double, L)
270};
271
272
273
283{
284
285public:
286 double L_;
287 double G_;
288 double H_;
289 double l_;
290 double g_;
291 double h_;
292public:
293 A_DEF_POD_ITERABLE(double)
294 AST_DEF_ACCESS_METHOD(double, L)
295 AST_DEF_ACCESS_METHOD(double, G)
296 AST_DEF_ACCESS_METHOD(double, H)
297 AST_DEF_ACCESS_METHOD(double, l)
298 AST_DEF_ACCESS_METHOD(double, g)
299 AST_DEF_ACCESS_METHOD(double, h)
300};
301
302
309AST_CORE_CAPI errc_t coe2rv(const double* coe, double gm, double* pos, double* vel);
310
315AST_CORE_CAPI errc_t coe2mee(const double* coe, double* mee);
316
322AST_CORE_CAPI void ee2rv(const double* ee, double gm, double* pos, double* vel);
323
329AST_CORE_CAPI void mee2rv(const double* mee, double gm, double* pos, double* vel);
330
337AST_CORE_CAPI errc_t rv2mee(const double* pos, const double* vel, double gm, double* mee);
338
343AST_CORE_CAPI errc_t mee2coe(const double* mee, double* coe);
344
350AST_CORE_CAPI void rv2ee(const double* pos, const double* vel, double gm, double* ee);
351
358AST_CORE_CAPI errc_t rv2moe(const double* pos, const double* vel, double gm, double* moe);
359
366AST_CORE_CAPI errc_t rv2coe(const double* pos, const double* vel, double gm, double* coe);
367
372AST_CORE_CAPI errc_t ee2moe(const double* ee, double* moe);
373
378AST_CORE_CAPI errc_t moe2ee(const double* moe, double* ee);
379
383AST_CORE_CAPI errc_t moe2coe(const double* moe, double* coe);
384
388AST_CORE_CAPI void coe2moe(const double* coe, double* moe);
389
395AST_CORE_CAPI errc_t moe2rv(const double* moe, double gm, double* pos, double* vel);
396
400AST_CORE_CAPI void moe2mee(const double* moe, double* mee);
401
405AST_CORE_CAPI void coe2ee(const double* coe, double* ee);
406
410AST_CORE_CAPI void ee2coe(const double* ee, double* coe);
411
415AST_CORE_CAPI void ee2mee(const double* ee, double* mee);
416
420AST_CORE_CAPI void mee2ee(const double* mee, double* ee);
421
425AST_CORE_CAPI void mee2moe(const double* mee, double* moe);
426
427
432AST_CORE_CAPI errc_t coe2dela(const double* coe, double gm, double* dela);
433
434
438AST_CORE_CAPI errc_t dela2coe(const double* dela, double gm, double* coe);
439
440
446AST_CORE_CAPI void aModEquinElemToCart(
447 const ModEquinElem& mee,
448 double gm,
449 Vector3d& pos,
450 Vector3d& vel
451);
452
459AST_CORE_CAPI errc_t aCartToModEquinElem(
460 const Vector3d& pos,
461 const Vector3d& vel,
462 double gm,
463 ModEquinElem& mee
464);
465
470AST_CORE_CAPI errc_t aOrbElemToModEquinElem(
471 const OrbElem& elem,
472 ModEquinElem& mee
473);
474
479AST_CORE_CAPI errc_t aModEquinElemToOrbElem(
480 const ModEquinElem& mee,
481 OrbElem& elem
482);
483
490AST_CORE_CAPI errc_t aCartToModOrbElem(
491 const Vector3d& pos,
492 const Vector3d& vel,
493 double gm,
494 ModOrbElem& modOrb);
495
502AST_CORE_CAPI
503errc_t aCartToOrbElem (
504 const Vector3d& pos,
505 const Vector3d& vel,
506 double gm,
507 OrbElem& elem);
508
513AST_CORE_CAPI
514errc_t aEquinElemToModOrb (
515 const EquinElem& equinElem,
516 ModOrbElem& modOrb);
517
522AST_CORE_CAPI
523errc_t aModOrbToEquinElem (
524 const ModOrbElem& modOrb,
525 EquinElem& equinElem);
526
532AST_CORE_CAPI
533errc_t aModOrbElemToCart (
534 const ModOrbElem& modOrb,
535 double gm,
536 Vector3d& pos,
537 Vector3d& vel);
538
545AST_CORE_CAPI
546errc_t aOrbElemToCart (
547 const OrbElem& elem,
548 double gm,
549 Vector3d& pos,
550 Vector3d& vel);
551
557AST_CORE_CAPI
559 const Vector3d& pos,
560 const Vector3d& vel,
561 double gm,
562 EquinElem& equinElem);
563
569AST_CORE_CAPI
571 const EquinElem& equinElem,
572 double gm,
573 Vector3d& pos,
574 Vector3d& vel);
575
576
581AST_CORE_CAPI
582errc_t aOrbElemToDelaunay(
583 const OrbElem& elem,
584 double gm,
585 DelaunayElem& delaunay);
586
587
592AST_CORE_CAPI
593errc_t aDelaunayToOrbElem(
594 const DelaunayElem& delaunay,
595 double gm,
596 OrbElem& elem);
597
598
604A_ALWAYS_INLINE
605ModOrbElem aCartToModOrbElem(const Vector3d& r, const Vector3d& v, double gm)
606{
607 ModOrbElem modOrbElem;
608 aCartToModOrbElem(r, v, gm, modOrbElem);
609 return modOrbElem;
610}
611
612
616AST_NAMESPACE_END
617
618AST_DECL_TYPE_ALIAS(CartState)
Unit h
小时
定义 Unit.cpp:439
Unit L
定义 Unit.cpp:458
Unit g
定义 Unit.cpp:443
轴系类
定义 Axes.hpp:70
直角坐标
定义 OrbitElement.hpp:46
Vector3d vel_
速度
定义 OrbitElement.hpp:49
Vector3d pos_
位置
定义 OrbitElement.hpp:48
const Vector3d & vel() const
获取速度
定义 OrbitElement.hpp:64
static CartState Zero()
获取零状态
定义 OrbitElement.hpp:53
德洛奈根数
定义 OrbitElement.hpp:283
春分点根数
定义 OrbitElement.hpp:234
double lambda_
mean longitude = M + raan + argper
定义 OrbitElement.hpp:241
double k_
e*cos(argper + raan)
定义 OrbitElement.hpp:238
double q_
tan(i/2)*cos(raan)
定义 OrbitElement.hpp:240
double a_
semimajor axis length
定义 OrbitElement.hpp:236
double p_
tan(i/2)*sin(raan)
定义 OrbitElement.hpp:239
double h_
e*sin(argper + raan) omegabar=argper + raan
定义 OrbitElement.hpp:237
改进春分点轨道根数, 180度奇异
定义 OrbitElement.hpp:254
double h_
h = tan(i/2)cos(RAAN)
定义 OrbitElement.hpp:259
double k_
k = tan(i/2)sin(RAAN)
定义 OrbitElement.hpp:260
double L_
L = RAAN + argper + trueA
定义 OrbitElement.hpp:261
double g_
g = e*sin(argper+RAAN)
定义 OrbitElement.hpp:258
double f_
f = e*cos(argper+RAAN)
定义 OrbitElement.hpp:257
double p_
p = a(1-e^2) 半通径
定义 OrbitElement.hpp:256
修正轨道根数
定义 OrbitElement.hpp:131
double rp_
近拱点半径
定义 OrbitElement.hpp:133
double getMeanMotion(double gm) const
计算平均角速度
定义 OrbitElement.hpp:144
double getA() const
计算半长轴
定义 OrbitElement.hpp:205
double getLAN(Axes *inertialAxes, const TimePoint &timeOfAscNodePassage, Axes *bodyFixedAxes) const
计算升交点经度
定义 OrbitElement.hpp:198
double trueA_
真近点角
定义 OrbitElement.hpp:138
double getSMA() const
计算半长轴
定义 OrbitElement.hpp:141
double getArgLat() const
计算纬度幅角
定义 OrbitElement.hpp:183
double getEccAnomaly() const
计算偏近点角
定义 OrbitElement.hpp:180
double getRAAN() const
计算升交点赤经
定义 OrbitElement.hpp:168
double getTimePastAscNode(double gm) const
计算过升交点后经过的时间
定义 OrbitElement.hpp:189
double getPeriAlt(double bodyRadius) const
计算近拱点高度
定义 OrbitElement.hpp:159
double getP() const
计算半通径
定义 OrbitElement.hpp:208
double getPeriod(double gm) const
计算周期
定义 OrbitElement.hpp:147
double e_
偏心率
定义 OrbitElement.hpp:134
double getApoRad() const
计算远拱点半径
定义 OrbitElement.hpp:150
double argper_
近拱点角
定义 OrbitElement.hpp:137
double getTimePastPeri(double gm) const
计算过近地点后经过的时间
定义 OrbitElement.hpp:186
double getMeanAnomaly() const
计算平近点角
定义 OrbitElement.hpp:177
double getPeriRad() const
计算近拱点半径
定义 OrbitElement.hpp:156
double getApoAlt(double bodyRadius) const
计算远拱点高度
定义 OrbitElement.hpp:153
double raan_
升交点赤经
定义 OrbitElement.hpp:136
double getEcc() const
计算偏心率
定义 OrbitElement.hpp:162
double getArgPeri() const
计算近拱点幅角
定义 OrbitElement.hpp:171
double getInc() const
计算轨道倾角
定义 OrbitElement.hpp:165
double i_
轨道倾角
定义 OrbitElement.hpp:135
TimePoint getTimeOfAscNodePassage(const TimePoint &stateEpoch, double gm) const
计算过升交点时刻
定义 OrbitElement.hpp:192
double getTrueAnomaly() const
计算真近点角
定义 OrbitElement.hpp:174
经典轨道根数
定义 OrbitElement.hpp:89
double a_
长半轴
定义 OrbitElement.hpp:91
double argper_
近拱点角
定义 OrbitElement.hpp:95
double trueA_
真近点角
定义 OrbitElement.hpp:96
double raan_
升交点赤经
定义 OrbitElement.hpp:94
double i_
轨道倾角
定义 OrbitElement.hpp:93
double getMeanMotion(double gm) const
计算平均角速度变化率
定义 OrbitElement.hpp:99
double e_
偏心率
定义 OrbitElement.hpp:92
绝对时间点
定义 TimePoint.hpp:106
errc_t aCartToModEquinElem(const Vector3d &pos, const Vector3d &vel, double gm, ModEquinElem &mee)
直角坐标转换为改进春分点轨道根数(类引用版本)
定义 OrbitElement.cpp:1036
errc_t moe2rv(const double *moe, double gm, double *pos, double *vel)
修正轨道根数转换为直角坐标
定义 OrbitElement.cpp:680
double aPeriRadToApoRad(double perigeeRad, double eccentricity)
近地点半径转换为远地点半径
定义 OrbitParam.cpp:316
double aPeriRadToMeanMotion(double perigeeRad, double eccentricity, double gm)
近地点半径转换为平均角速度
定义 OrbitParam.cpp:321
void mee2ee(const double *mee, double *ee)
改进春分点轨道根数转换为春分点根数
定义 OrbitElement.cpp:857
double aRAANToLAN(double raan, Axes *inertialAxes, const TimePoint &timeOfAscNodePassage, Axes *bodyFixedAxes)
升交点赤经转换为升交点经度
定义 OrbitParam.cpp:511
errc_t aModOrbToEquinElem(const ModOrbElem &modOrb, EquinElem &equinElem)
修正轨道根数转换为春分点根数(类引用版本)
定义 OrbitElement.cpp:1098
errc_t mee2coe(const double *mee, double *coe)
改进春分点轨道根数转换为经典轨道根数
定义 OrbitElement.cpp:267
void rv2ee(const double *pos, const double *vel, double gm, double *ee)
直角坐标转换为春分点根数
定义 OrbitElement.cpp:329
void coe2moe(const double *coe_, double *moe_)
经典轨道根数转换为修正轨道根数
定义 OrbitElement.cpp:669
double aSMAToMeanMotion(double semiMajorAxis, double gm)
长半轴转换为平均角速度
定义 OrbitParam.cpp:384
errc_t coe2rv(const double *coe, double gm, double *pos, double *vel)
经典轨道根数转换为直角坐标
定义 OrbitElement.cpp:69
errc_t ee2moe(const double *ee, double *moe)
春分点根数转换为修正轨道根数
定义 OrbitElement.cpp:545
errc_t moe2ee(const double *moe, double *ee)
修正轨道根数转换为春分点根数
定义 OrbitElement.cpp:622
double aTrueToMean(double trueAnomaly, double eccentricity)
真近点角转换为平近点角
定义 OrbitParam.cpp:484
errc_t rv2moe(const double *pos, const double *vel, double gm, double *moe)
直角坐标转换为修正轨道根数
定义 OrbitElement.cpp:397
void aModEquinElemToCart(const ModEquinElem &mee, double gm, Vector3d &pos, Vector3d &vel)
改进春分点轨道根数转换为直角坐标(类引用版本)
定义 OrbitElement.cpp:1024
errc_t coe2mee(const double *coe, double *mee)
经典轨道根数转换为改进春分点轨道根数
定义 OrbitElement.cpp:104
errc_t aDelaunayToOrbElem(const DelaunayElem &delaunay, double gm, OrbElem &elem)
修正轨道根数转换为经典轨道根数(类引用版本)
定义 OrbitElement.cpp:1153
errc_t moe2coe(const double *moe_, double *coe_)
修正轨道根数转换为经典轨道根数
定义 OrbitElement.cpp:656
double aTrueToTimePastAscNode(double trueAnomaly, double argPeri, double semiMajorAxis, double eccentricity, double gm)
真近点角转换为过升交点后时间
定义 OrbitParam.cpp:489
errc_t aOrbElemToDelaunay(const OrbElem &elem, double gm, DelaunayElem &delaunay)
经典轨道根数转换为修正轨道根数(类引用版本)
定义 OrbitElement.cpp:1148
void ee2rv(const double *ee, double gm, double *pos, double *vel)
春分点根数转换为直角坐标
定义 OrbitElement.cpp:127
void mee2rv(const double *mee, double gm, double *pos, double *vel)
改进春分点轨道根数转换为直角坐标
定义 OrbitElement.cpp:181
errc_t rv2coe(const double *pos, const double *vel, double gm, double *coe)
直角坐标转换为经典轨道根数
定义 OrbitElement.cpp:461
void aCartToEquinElem(const Vector3d &pos, const Vector3d &vel, double gm, EquinElem &equinElem)
直角坐标转换为春分点根数(类引用版本)
定义 OrbitElement.cpp:1129
void aEquinElemToCart(const EquinElem &equinElem, double gm, Vector3d &pos, Vector3d &vel)
春分点根数转换为直角坐标(类引用版本)
定义 OrbitElement.cpp:1139
void ee2mee(const double *ee, double *mee)
春分点根数转换为改进春分点轨道根数
定义 OrbitElement.cpp:821
errc_t aCartToModOrbElem(const Vector3d &pos, const Vector3d &vel, double gm, ModOrbElem &modOrb)
直角坐标转换为修正轨道根数(类引用版本)
定义 OrbitElement.cpp:1066
errc_t aModEquinElemToOrbElem(const ModEquinElem &mee, OrbElem &elem)
改进春分点轨道根数转换为经典轨道根数(类引用版本)
定义 OrbitElement.cpp:1057
void mee2moe(const double *mee, double *moe)
改进春分点轨道根数转换为修正轨道根数
定义 OrbitElement.cpp:876
errc_t aOrbElemToModEquinElem(const OrbElem &elem, ModEquinElem &mee)
经典轨道根数转换为改进春分点轨道根数(类引用版本)
定义 OrbitElement.cpp:1048
double aTrueToEcc(double f, double e)
真近点角转换为偏近点角
定义 OrbitParam.cpp:457
double aTrueToArgLat(double trueAnomaly, double argPeri)
真近点角转换为纬度幅角
定义 OrbitParam.hpp:653
errc_t aOrbElemToCart(const OrbElem &elem, double gm, Vector3d &pos, Vector3d &vel)
经典轨道根数转换为直角坐标(类引用版本)
定义 OrbitElement.cpp:1118
errc_t aEquinElemToModOrb(const EquinElem &equinElem, ModOrbElem &modOrb)
春分点根数转换为修正轨道根数(类引用版本)
定义 OrbitElement.cpp:1089
errc_t rv2mee(const double *pos_, const double *vel_, double gm, double *mee)
直角坐标转换为改进春分点轨道根数
定义 OrbitElement.cpp:204
errc_t aModOrbElemToCart(const ModOrbElem &modOrb, double gm, Vector3d &pos, Vector3d &vel)
修正轨道根数转换为直角坐标(类引用版本)
定义 OrbitElement.cpp:1107
void coe2ee(const double *coe, double *ee)
经典轨道根数转换为春分点根数
定义 OrbitElement.cpp:730
errc_t aCartToOrbElem(const Vector3d &pos, const Vector3d &vel, double gm, OrbElem &elem)
直角坐标转换为经典轨道根数(类引用版本)
定义 OrbitElement.cpp:1077
void ee2coe(const double *ee, double *coe)
春分点根数转换为经典轨道根数
定义 OrbitElement.cpp:748
double aTrueToTimePastPeri(double trueAnomaly, double semiMajorAxis, double eccentricity, double gm)
真近点角转换为过近心点后时间
定义 OrbitParam.cpp:494
errc_t coe2dela(const double *coeIn, double gm, double *delaOut)
经典轨道根数转换为德洛奈根数
定义 OrbitElement.cpp:931
void moe2mee(const double *moe, double *mee)
修正轨道根数转换为改进春分点轨道根数
定义 OrbitElement.cpp:713
double aPeriRadToApoAlt(double perigeeRad, double eccentricity, double bodyRadius)
近地点半径转换为远地点高度
定义 OrbitParam.cpp:312
errc_t dela2coe(const double *delaIn, double gm, double *coeOut)
德洛奈根数转换为经典轨道根数
定义 OrbitElement.cpp:985
constexpr double kTwoPI
2*PI
定义 Constants.h:54