🛰️航天仿真算法库 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};
81
82
85{
86public:
87 double a_;
88 double e_;
89 double i_;
90 double raan_;
91 double argper_;
92 double trueA_;
93public:
95 double getMeanMotion(double gm) const {return aSMAToMeanMotion(getSMA(), gm);}
96
97 double getSMA() const {return a_;}
98
99 double getA() const {return a_;}
100
101 double getE() const {return e_;}
102
103 double getI() const {return i_;}
104
105 double getRAAN() const {return raan_;}
106
107 double getArgPer() const {return argper_;}
108
109 double getTrueA() const {return trueA_;}
110
112 AST_CORE_API
113 std::string toString() const;
114public:
115 A_DEF_POD_ITERABLE(double)
116 AST_DEF_ACCESS_METHOD(double, a)
117 AST_DEF_ACCESS_METHOD(double, e)
118 AST_DEF_ACCESS_METHOD(double, i)
119 AST_DEF_ACCESS_METHOD(double, raan)
120 AST_DEF_ACCESS_METHOD(double, argper)
121 AST_DEF_ACCESS_METHOD(double, trueA)
122};
123
124
127{
128public:
129 double rp_;
130 double e_;
131 double i_;
132 double raan_;
133 double argper_;
134 double trueA_;
135public:
137 double getSMA() const {return rp_ / (1 - e_);}
138
140 double getMeanMotion(double gm) const {return aPeriRadToMeanMotion(getPeriRad(), getEcc(), gm);}
141
143 double getPeriod(double gm) const {return kTwoPI / getMeanMotion(gm);}
144
146 double getApoRad() const{return aPeriRadToApoRad(getPeriRad(), getEcc());}
147
149 double getApoAlt(double bodyRadius) const{ return aPeriRadToApoAlt(getPeriRad(), getEcc(), bodyRadius);}
150
152 double getPeriRad() const {return rp_;}
153
155 double getPeriAlt(double bodyRadius) const{ return getPeriRad() - bodyRadius;}
156
158 double getEcc() const {return e_;}
159
161 double getInc() const {return i_;}
162
164 double getRAAN() const {return raan_;}
165
167 double getArgPeri() const {return argper_;}
168
170 double getTrueAnomaly() const {return trueA_;}
171
173 double getMeanAnomaly() const {return aTrueToMean(getTrueAnomaly(), getEcc());}
174
176 double getEccAnomaly() const {return aTrueToEcc(getTrueAnomaly(), getEcc());}
177
179 double getArgLat() const {return aTrueToArgLat(getTrueAnomaly(), getArgPeri());}
180
182 double getTimePastPeri(double gm) const{return aTrueToTimePastPeri(getTrueAnomaly(), getSMA(), getEcc(), gm);}
183
185 double getTimePastAscNode(double gm) const{return aTrueToTimePastAscNode(getTrueAnomaly(), getArgPeri(), getSMA(), getEcc(), gm);}
186
188 TimePoint getTimeOfAscNodePassage(const TimePoint& stateEpoch, double gm) const
189 {
190 return stateEpoch - getTimePastAscNode(gm);
191 }
192
194 double getLAN(Axes* inertialAxes, const TimePoint& timeOfAscNodePassage, Axes* bodyFixedAxes) const
195 {
196 return aRAANToLAN(getRAAN(), inertialAxes, timeOfAscNodePassage, bodyFixedAxes);
197 }
198
199public:
201 double getA() const {return getSMA();}
202
204 double getP() const {return rp_ * (1 + e_);}
205
206 double getE() const {return getEcc();}
207
208 double getI() const {return getInc();}
209
210 double getArgPer() const {return getArgPeri();}
211
212 double getTrueA() const {return getTrueAnomaly();}
213
215 AST_CORE_API
216 std::string toString() const;
217public:
218 A_DEF_POD_ITERABLE(double)
219 AST_DEF_ACCESS_METHOD(double, rp)
220 AST_DEF_ACCESS_METHOD(double, e)
221 AST_DEF_ACCESS_METHOD(double, i)
222 AST_DEF_ACCESS_METHOD(double, raan)
223 AST_DEF_ACCESS_METHOD(double, argper)
224 AST_DEF_ACCESS_METHOD(double, trueA)
225};
226
227
230{
231public:
232 double a_;
233 double h_;
234 double k_;
235 double p_;
236 double q_;
237 double lambda_;
238public:
239 A_DEF_POD_ITERABLE(double)
240 AST_DEF_ACCESS_METHOD(double, a)
241 AST_DEF_ACCESS_METHOD(double, h)
242 AST_DEF_ACCESS_METHOD(double, k)
243 AST_DEF_ACCESS_METHOD(double, p)
244 AST_DEF_ACCESS_METHOD(double, q)
245 AST_DEF_ACCESS_METHOD(double, lambda)
246};
247
250{
251public:
252 double p_;
253 double f_;
254 double g_;
255 double h_;
256 double k_;
257 double L_;
258public:
259 A_DEF_POD_ITERABLE(double)
260 AST_DEF_ACCESS_METHOD(double, p)
261 AST_DEF_ACCESS_METHOD(double, f)
262 AST_DEF_ACCESS_METHOD(double, g)
263 AST_DEF_ACCESS_METHOD(double, h)
264 AST_DEF_ACCESS_METHOD(double, k)
265 AST_DEF_ACCESS_METHOD(double, L)
266};
267
268
269
279{
280
281public:
282 double L_;
283 double G_;
284 double H_;
285 double l_;
286 double g_;
287 double h_;
288public:
289 A_DEF_POD_ITERABLE(double)
290 AST_DEF_ACCESS_METHOD(double, L)
291 AST_DEF_ACCESS_METHOD(double, G)
292 AST_DEF_ACCESS_METHOD(double, H)
293 AST_DEF_ACCESS_METHOD(double, l)
294 AST_DEF_ACCESS_METHOD(double, g)
295 AST_DEF_ACCESS_METHOD(double, h)
296};
297
298
305AST_CORE_CAPI errc_t coe2rv(const double* coe, double gm, double* pos, double* vel);
306
311AST_CORE_CAPI errc_t coe2mee(const double* coe, double* mee);
312
318AST_CORE_CAPI void ee2rv(const double* ee, double gm, double* pos, double* vel);
319
325AST_CORE_CAPI void mee2rv(const double* mee, double gm, double* pos, double* vel);
326
333AST_CORE_CAPI errc_t rv2mee(const double* pos, const double* vel, double gm, double* mee);
334
339AST_CORE_CAPI errc_t mee2coe(const double* mee, double* coe);
340
346AST_CORE_CAPI void rv2ee(const double* pos, const double* vel, double gm, double* ee);
347
354AST_CORE_CAPI errc_t rv2moe(const double* pos, const double* vel, double gm, double* moe);
355
362AST_CORE_CAPI errc_t rv2coe(const double* pos, const double* vel, double gm, double* coe);
363
368AST_CORE_CAPI errc_t ee2moe(const double* ee, double* moe);
369
374AST_CORE_CAPI errc_t moe2ee(const double* moe, double* ee);
375
379AST_CORE_CAPI errc_t moe2coe(const double* moe, double* coe);
380
384AST_CORE_CAPI void coe2moe(const double* coe, double* moe);
385
391AST_CORE_CAPI errc_t moe2rv(const double* moe, double gm, double* pos, double* vel);
392
396AST_CORE_CAPI void moe2mee(const double* moe, double* mee);
397
401AST_CORE_CAPI void coe2ee(const double* coe, double* ee);
402
406AST_CORE_CAPI void ee2coe(const double* ee, double* coe);
407
411AST_CORE_CAPI void ee2mee(const double* ee, double* mee);
412
416AST_CORE_CAPI void mee2ee(const double* mee, double* ee);
417
421AST_CORE_CAPI void mee2moe(const double* mee, double* moe);
422
423
428AST_CORE_CAPI errc_t coe2dela(const double* coe, double gm, double* dela);
429
430
434AST_CORE_CAPI errc_t dela2coe(const double* dela, double gm, double* coe);
435
436
442AST_CORE_CAPI void aModEquinElemToCart(
443 const ModEquinElem& mee,
444 double gm,
445 Vector3d& pos,
446 Vector3d& vel
447);
448
455AST_CORE_CAPI errc_t aCartToModEquinElem(
456 const Vector3d& pos,
457 const Vector3d& vel,
458 double gm,
459 ModEquinElem& mee
460);
461
466AST_CORE_CAPI errc_t aOrbElemToModEquinElem(
467 const OrbElem& elem,
468 ModEquinElem& mee
469);
470
475AST_CORE_CAPI errc_t aModEquinElemToOrbElem(
476 const ModEquinElem& mee,
477 OrbElem& elem
478);
479
486AST_CORE_CAPI errc_t aCartToModOrbElem(
487 const Vector3d& pos,
488 const Vector3d& vel,
489 double gm,
490 ModOrbElem& modOrb);
491
498AST_CORE_CAPI
499errc_t aCartToOrbElem (
500 const Vector3d& pos,
501 const Vector3d& vel,
502 double gm,
503 OrbElem& elem);
504
509AST_CORE_CAPI
510errc_t aEquinElemToModOrb (
511 const EquinElem& equinElem,
512 ModOrbElem& modOrb);
513
518AST_CORE_CAPI
519errc_t aModOrbToEquinElem (
520 const ModOrbElem& modOrb,
521 EquinElem& equinElem);
522
528AST_CORE_CAPI
529errc_t aModOrbElemToCart (
530 const ModOrbElem& modOrb,
531 double gm,
532 Vector3d& pos,
533 Vector3d& vel);
534
541AST_CORE_CAPI
542errc_t aOrbElemToCart (
543 const OrbElem& elem,
544 double gm,
545 Vector3d& pos,
546 Vector3d& vel);
547
553AST_CORE_CAPI
555 const Vector3d& pos,
556 const Vector3d& vel,
557 double gm,
558 EquinElem& equinElem);
559
565AST_CORE_CAPI
567 const EquinElem& equinElem,
568 double gm,
569 Vector3d& pos,
570 Vector3d& vel);
571
572
577AST_CORE_CAPI
578errc_t aOrbElemToDelaunay(
579 const OrbElem& elem,
580 double gm,
581 DelaunayElem& delaunay);
582
583
588AST_CORE_CAPI
589errc_t aDelaunayToOrbElem(
590 const DelaunayElem& delaunay,
591 double gm,
592 OrbElem& elem);
593
594
600A_ALWAYS_INLINE
601ModOrbElem aCartToModOrbElem(const Vector3d& r, const Vector3d& v, double gm)
602{
603 ModOrbElem modOrbElem;
604 aCartToModOrbElem(r, v, gm, modOrbElem);
605 return modOrbElem;
606}
607
608
612AST_NAMESPACE_END
613
614AST_DECL_TYPE_ALIAS(CartState)
Unit h
小时
定义 Unit.cpp:429
Unit L
定义 Unit.cpp:448
Unit g
定义 Unit.cpp:433
轴系类
定义 Axes.hpp:69
直角坐标
定义 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:279
春分点根数
定义 OrbitElement.hpp:230
double lambda_
mean longitude = M + raan + argper
定义 OrbitElement.hpp:237
double k_
e*cos(argper + raan)
定义 OrbitElement.hpp:234
double q_
tan(i/2)*cos(raan)
定义 OrbitElement.hpp:236
double a_
semimajor axis length
定义 OrbitElement.hpp:232
double p_
tan(i/2)*sin(raan)
定义 OrbitElement.hpp:235
double h_
e*sin(argper + raan) omegabar=argper + raan
定义 OrbitElement.hpp:233
改进春分点轨道根数, 180度奇异
定义 OrbitElement.hpp:250
double h_
h = tan(i/2)cos(RAAN)
定义 OrbitElement.hpp:255
double k_
k = tan(i/2)sin(RAAN)
定义 OrbitElement.hpp:256
double L_
L = RAAN + argper + trueA
定义 OrbitElement.hpp:257
double g_
g = e*sin(argper+RAAN)
定义 OrbitElement.hpp:254
double f_
f = e*cos(argper+RAAN)
定义 OrbitElement.hpp:253
double p_
p = a(1-e^2) 半通径
定义 OrbitElement.hpp:252
修正轨道根数
定义 OrbitElement.hpp:127
double rp_
近拱点半径
定义 OrbitElement.hpp:129
double getMeanMotion(double gm) const
计算平均角速度
定义 OrbitElement.hpp:140
double getA() const
计算半长轴
定义 OrbitElement.hpp:201
double getLAN(Axes *inertialAxes, const TimePoint &timeOfAscNodePassage, Axes *bodyFixedAxes) const
计算升交点经度
定义 OrbitElement.hpp:194
double trueA_
真近点角
定义 OrbitElement.hpp:134
double getSMA() const
计算半长轴
定义 OrbitElement.hpp:137
double getArgLat() const
计算纬度幅角
定义 OrbitElement.hpp:179
double getEccAnomaly() const
计算偏近点角
定义 OrbitElement.hpp:176
double getRAAN() const
计算升交点赤经
定义 OrbitElement.hpp:164
double getTimePastAscNode(double gm) const
计算过升交点后经过的时间
定义 OrbitElement.hpp:185
double getPeriAlt(double bodyRadius) const
计算近拱点高度
定义 OrbitElement.hpp:155
double getP() const
计算半通径
定义 OrbitElement.hpp:204
double getPeriod(double gm) const
计算周期
定义 OrbitElement.hpp:143
double e_
偏心率
定义 OrbitElement.hpp:130
double getApoRad() const
计算远拱点半径
定义 OrbitElement.hpp:146
double argper_
近拱点角
定义 OrbitElement.hpp:133
double getTimePastPeri(double gm) const
计算过近地点后经过的时间
定义 OrbitElement.hpp:182
double getMeanAnomaly() const
计算平近点角
定义 OrbitElement.hpp:173
double getPeriRad() const
计算近拱点半径
定义 OrbitElement.hpp:152
double getApoAlt(double bodyRadius) const
计算远拱点高度
定义 OrbitElement.hpp:149
double raan_
升交点赤经
定义 OrbitElement.hpp:132
double getEcc() const
计算偏心率
定义 OrbitElement.hpp:158
double getArgPeri() const
计算近拱点幅角
定义 OrbitElement.hpp:167
double getInc() const
计算轨道倾角
定义 OrbitElement.hpp:161
double i_
轨道倾角
定义 OrbitElement.hpp:131
TimePoint getTimeOfAscNodePassage(const TimePoint &stateEpoch, double gm) const
计算过升交点时刻
定义 OrbitElement.hpp:188
double getTrueAnomaly() const
计算真近点角
定义 OrbitElement.hpp:170
经典轨道根数
定义 OrbitElement.hpp:85
double a_
长半轴
定义 OrbitElement.hpp:87
double argper_
近拱点角
定义 OrbitElement.hpp:91
double trueA_
真近点角
定义 OrbitElement.hpp:92
double raan_
升交点赤经
定义 OrbitElement.hpp:90
double i_
轨道倾角
定义 OrbitElement.hpp:89
double getMeanMotion(double gm) const
计算平均角速度变化率
定义 OrbitElement.hpp:95
double e_
偏心率
定义 OrbitElement.hpp:88
绝对时间点
定义 TimePoint.hpp:106
errc_t aCartToModEquinElem(const Vector3d &pos, const Vector3d &vel, double gm, ModEquinElem &mee)
直角坐标转换为改进春分点轨道根数(类引用版本)
定义 OrbitElement.cpp:1028
errc_t moe2rv(const double *moe, double gm, double *pos, double *vel)
修正轨道根数转换为直角坐标
定义 OrbitElement.cpp:672
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:849
double aRAANToLAN(double raan, Axes *inertialAxes, const TimePoint &timeOfAscNodePassage, Axes *bodyFixedAxes)
升交点赤经转换为升交点经度
定义 OrbitParam.cpp:511
errc_t aModOrbToEquinElem(const ModOrbElem &modOrb, EquinElem &equinElem)
修正轨道根数转换为春分点根数(类引用版本)
定义 OrbitElement.cpp:1090
errc_t mee2coe(const double *mee, double *coe)
改进春分点轨道根数转换为经典轨道根数
定义 OrbitElement.cpp:259
void rv2ee(const double *pos, const double *vel, double gm, double *ee)
直角坐标转换为春分点根数
定义 OrbitElement.cpp:321
void coe2moe(const double *coe_, double *moe_)
经典轨道根数转换为修正轨道根数
定义 OrbitElement.cpp:661
double aSMAToMeanMotion(double semiMajorAxis, double gm)
长半轴转换为平均角速度
定义 OrbitParam.cpp:384
errc_t coe2rv(const double *coe, double gm, double *pos, double *vel)
经典轨道根数转换为直角坐标
定义 OrbitElement.cpp:61
errc_t ee2moe(const double *ee, double *moe)
春分点根数转换为修正轨道根数
定义 OrbitElement.cpp:537
errc_t moe2ee(const double *moe, double *ee)
修正轨道根数转换为春分点根数
定义 OrbitElement.cpp:614
double aTrueToMean(double trueAnomaly, double eccentricity)
真近点角转换为平近点角
定义 OrbitParam.cpp:484
errc_t rv2moe(const double *pos, const double *vel, double gm, double *moe)
直角坐标转换为修正轨道根数
定义 OrbitElement.cpp:389
void aModEquinElemToCart(const ModEquinElem &mee, double gm, Vector3d &pos, Vector3d &vel)
改进春分点轨道根数转换为直角坐标(类引用版本)
定义 OrbitElement.cpp:1016
errc_t coe2mee(const double *coe, double *mee)
经典轨道根数转换为改进春分点轨道根数
定义 OrbitElement.cpp:96
errc_t aDelaunayToOrbElem(const DelaunayElem &delaunay, double gm, OrbElem &elem)
修正轨道根数转换为经典轨道根数(类引用版本)
定义 OrbitElement.cpp:1145
errc_t moe2coe(const double *moe_, double *coe_)
修正轨道根数转换为经典轨道根数
定义 OrbitElement.cpp:648
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:1140
void ee2rv(const double *ee, double gm, double *pos, double *vel)
春分点根数转换为直角坐标
定义 OrbitElement.cpp:119
void mee2rv(const double *mee, double gm, double *pos, double *vel)
改进春分点轨道根数转换为直角坐标
定义 OrbitElement.cpp:173
errc_t rv2coe(const double *pos, const double *vel, double gm, double *coe)
直角坐标转换为经典轨道根数
定义 OrbitElement.cpp:453
void aCartToEquinElem(const Vector3d &pos, const Vector3d &vel, double gm, EquinElem &equinElem)
直角坐标转换为春分点根数(类引用版本)
定义 OrbitElement.cpp:1121
void aEquinElemToCart(const EquinElem &equinElem, double gm, Vector3d &pos, Vector3d &vel)
春分点根数转换为直角坐标(类引用版本)
定义 OrbitElement.cpp:1131
void ee2mee(const double *ee, double *mee)
春分点根数转换为改进春分点轨道根数
定义 OrbitElement.cpp:813
errc_t aCartToModOrbElem(const Vector3d &pos, const Vector3d &vel, double gm, ModOrbElem &modOrb)
直角坐标转换为修正轨道根数(类引用版本)
定义 OrbitElement.cpp:1058
errc_t aModEquinElemToOrbElem(const ModEquinElem &mee, OrbElem &elem)
改进春分点轨道根数转换为经典轨道根数(类引用版本)
定义 OrbitElement.cpp:1049
void mee2moe(const double *mee, double *moe)
改进春分点轨道根数转换为修正轨道根数
定义 OrbitElement.cpp:868
errc_t aOrbElemToModEquinElem(const OrbElem &elem, ModEquinElem &mee)
经典轨道根数转换为改进春分点轨道根数(类引用版本)
定义 OrbitElement.cpp:1040
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:1110
errc_t aEquinElemToModOrb(const EquinElem &equinElem, ModOrbElem &modOrb)
春分点根数转换为修正轨道根数(类引用版本)
定义 OrbitElement.cpp:1081
errc_t rv2mee(const double *pos_, const double *vel_, double gm, double *mee)
直角坐标转换为改进春分点轨道根数
定义 OrbitElement.cpp:196
errc_t aModOrbElemToCart(const ModOrbElem &modOrb, double gm, Vector3d &pos, Vector3d &vel)
修正轨道根数转换为直角坐标(类引用版本)
定义 OrbitElement.cpp:1099
void coe2ee(const double *coe, double *ee)
经典轨道根数转换为春分点根数
定义 OrbitElement.cpp:722
errc_t aCartToOrbElem(const Vector3d &pos, const Vector3d &vel, double gm, OrbElem &elem)
直角坐标转换为经典轨道根数(类引用版本)
定义 OrbitElement.cpp:1069
void ee2coe(const double *ee, double *coe)
春分点根数转换为经典轨道根数
定义 OrbitElement.cpp:740
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:923
void moe2mee(const double *moe, double *mee)
修正轨道根数转换为改进春分点轨道根数
定义 OrbitElement.cpp:705
double aPeriRadToApoAlt(double perigeeRad, double eccentricity, double bodyRadius)
近地点半径转换为远地点高度
定义 OrbitParam.cpp:312
errc_t dela2coe(const double *delaIn, double gm, double *coeOut)
德洛奈根数转换为经典轨道根数
定义 OrbitElement.cpp:977
constexpr double kTwoPI
2*PI
定义 Constants.h:54