1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17 package org.orekit.propagation.analytical.gnss;
18
19 import org.hipparchus.Field;
20 import org.hipparchus.analysis.differentiation.Gradient;
21 import org.hipparchus.analysis.differentiation.GradientField;
22 import org.hipparchus.analysis.differentiation.UnivariateDerivative2;
23 import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
24 import org.hipparchus.geometry.euclidean.threed.Vector3D;
25 import org.hipparchus.linear.MatrixUtils;
26 import org.hipparchus.linear.QRDecomposition;
27 import org.hipparchus.linear.RealMatrix;
28 import org.hipparchus.linear.RealVector;
29 import org.hipparchus.util.FastMath;
30 import org.hipparchus.util.FieldSinCos;
31 import org.hipparchus.util.SinCos;
32 import org.orekit.attitudes.Attitude;
33 import org.orekit.attitudes.AttitudeProvider;
34 import org.orekit.frames.Frame;
35 import org.orekit.orbits.CartesianOrbit;
36 import org.orekit.orbits.FieldKeplerianAnomalyUtility;
37 import org.orekit.orbits.KeplerianAnomalyUtility;
38 import org.orekit.orbits.KeplerianOrbit;
39 import org.orekit.orbits.Orbit;
40 import org.orekit.orbits.PositionAngleType;
41 import org.orekit.propagation.AbstractMatricesHarvester;
42 import org.orekit.propagation.SpacecraftState;
43 import org.orekit.propagation.analytical.AbstractAnalyticalPropagator;
44 import org.orekit.propagation.analytical.gnss.data.FieldGnssOrbitalElements;
45 import org.orekit.propagation.analytical.gnss.data.GNSSOrbitalElements;
46 import org.orekit.time.AbsoluteDate;
47 import org.orekit.utils.DoubleArrayDictionary;
48 import org.orekit.utils.FieldPVCoordinates;
49 import org.orekit.utils.PVCoordinates;
50 import org.orekit.utils.ParameterDriver;
51 import org.orekit.utils.TimeSpanMap.Span;
52
53 import java.util.ArrayList;
54 import java.util.Collections;
55 import java.util.List;
56
57
58
59
60
61
62
63
64 public class GNSSPropagator extends AbstractAnalyticalPropagator {
65
66
67
68
69 private static final int MAX_ITER = 100;
70
71
72
73
74 private static final double TOL_P = 1.0e-6;
75
76
77
78
79 private static final double TOL_V = 1.0e-9;
80
81
82
83
84 private static final int FREE_PARAMETERS = 6;
85
86
87
88
89 private static final double EPS = 1.0e-12;
90
91
92 private GNSSOrbitalElements<?> orbitalElements;
93
94
95 private final Frame eci;
96
97
98 private final Frame ecef;
99
100
101
102
103
104
105
106
107
108 GNSSPropagator(final GNSSOrbitalElements<?> orbitalElements, final Frame eci,
109 final Frame ecef, final AttitudeProvider provider, final double mass) {
110 super(provider);
111
112 this.orbitalElements = orbitalElements;
113
114 this.eci = eci;
115
116 this.ecef = ecef;
117
118 final Orbit orbit = propagateOrbit(orbitalElements.getDate());
119 final Attitude attitude = provider.getAttitude(orbit, orbit.getDate(), orbit.getFrame());
120
121
122 super.resetInitialState(new SpacecraftState(orbit, attitude).withMass(mass));
123
124 }
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140 GNSSPropagator(final SpacecraftState initialState, final GNSSOrbitalElements<?> nonKeplerianElements,
141 final Frame ecef, final AttitudeProvider provider, final double mass) {
142 this(buildOrbitalElements(initialState, nonKeplerianElements, ecef, provider, mass),
143 initialState.getFrame(), ecef, provider, initialState.getMass());
144 }
145
146
147
148
149
150
151 public Frame getECI() {
152 return eci;
153 }
154
155
156
157
158
159
160
161 public Frame getECEF() {
162 return ecef;
163 }
164
165
166
167
168
169
170 public double getMU() {
171 return orbitalElements.getMu();
172 }
173
174
175
176
177
178 public GNSSOrbitalElements<?> getOrbitalElements() {
179 return orbitalElements;
180 }
181
182
183
184
185 @Override
186 protected AbstractMatricesHarvester createHarvester(final String stmName, final RealMatrix initialStm,
187 final DoubleArrayDictionary initialJacobianColumns) {
188
189
190 final GnssHarvester harvester = new GnssHarvester(this, stmName, initialStm, initialJacobianColumns);
191
192
193 addAdditionalDataProvider(harvester);
194
195
196 return harvester;
197
198 }
199
200
201
202
203 @Override
204 protected List<String> getJacobiansColumnsNames() {
205 final List<String> columnsNames = new ArrayList<>();
206 for (final ParameterDriver driver : orbitalElements.getParametersDrivers()) {
207 if (driver.isSelected() && !columnsNames.contains(driver.getNamesSpanMap().getFirstSpan().getData())) {
208
209
210 for (Span<String> span = driver.getNamesSpanMap().getFirstSpan(); span != null; span = span.next()) {
211 columnsNames.add(span.getData());
212 }
213 }
214 }
215 Collections.sort(columnsNames);
216 return columnsNames;
217 }
218
219
220 @Override
221 public Orbit propagateOrbit(final AbsoluteDate date) {
222
223 final PVCoordinates pvaInECEF = propagateInEcef(date);
224
225 final PVCoordinates pvaInECI = ecef.getTransformTo(eci, date).transformPVCoordinates(pvaInECEF);
226
227 return new CartesianOrbit(pvaInECI, eci, date, getMU());
228 }
229
230
231
232
233
234
235
236
237
238
239 public PVCoordinates propagateInEcef(final AbsoluteDate date) {
240
241 final UnivariateDerivative2 tk = new UnivariateDerivative2(getTk(date), 1.0, 0.0);
242
243 final UnivariateDerivative2 ak = tk.multiply(orbitalElements.getADot()).add(orbitalElements.getSma());
244
245 final UnivariateDerivative2 nA = tk.multiply(orbitalElements.getDeltaN0Dot() * 0.5).
246 add(orbitalElements.getDeltaN0()).
247 add(orbitalElements.getMeanMotion0());
248
249 final UnivariateDerivative2 mk = tk.multiply(nA).add(orbitalElements.getM0());
250
251 final UnivariateDerivative2 e = tk.newInstance(orbitalElements.getE());
252 final UnivariateDerivative2 ek = FieldKeplerianAnomalyUtility.ellipticMeanToEccentric(e, mk);
253
254 final UnivariateDerivative2 vk = FieldKeplerianAnomalyUtility.ellipticEccentricToTrue(e, ek);
255
256 final UnivariateDerivative2 phik = vk.add(orbitalElements.getPa());
257 final FieldSinCos<UnivariateDerivative2> cs2phi = FastMath.sinCos(phik.multiply(2));
258
259 final UnivariateDerivative2 dphik = cs2phi.cos().multiply(orbitalElements.getCuc()).add(cs2phi.sin().multiply(orbitalElements.getCus()));
260
261 final UnivariateDerivative2 drk = cs2phi.cos().multiply(orbitalElements.getCrc()).add(cs2phi.sin().multiply(orbitalElements.getCrs()));
262
263 final UnivariateDerivative2 dik = cs2phi.cos().multiply(orbitalElements.getCic()).add(cs2phi.sin().multiply(orbitalElements.getCis()));
264
265 final FieldSinCos<UnivariateDerivative2> csuk = FastMath.sinCos(phik.add(dphik));
266
267 final UnivariateDerivative2 rk = ek.cos().multiply(-orbitalElements.getE()).add(1).multiply(ak).add(drk);
268
269 final UnivariateDerivative2 ik = tk.multiply(orbitalElements.getIDot()).add(orbitalElements.getI0()).add(dik);
270 final FieldSinCos<UnivariateDerivative2> csik = FastMath.sinCos(ik);
271
272 final UnivariateDerivative2 xk = csuk.cos().multiply(rk);
273 final UnivariateDerivative2 yk = csuk.sin().multiply(rk);
274
275 final double thetaDot = orbitalElements.getAngularVelocity();
276 final FieldSinCos<UnivariateDerivative2> csomk =
277 FastMath.sinCos(tk.multiply(orbitalElements.getOmegaDot() - thetaDot).
278 add(orbitalElements.getOmega0() - thetaDot * orbitalElements.getTime()));
279
280 final FieldVector3D<UnivariateDerivative2> positionWithDerivatives =
281 new FieldVector3D<>(xk.multiply(csomk.cos()).subtract(yk.multiply(csomk.sin()).multiply(csik.cos())),
282 xk.multiply(csomk.sin()).add(yk.multiply(csomk.cos()).multiply(csik.cos())),
283 yk.multiply(csik.sin()));
284 return new PVCoordinates(positionWithDerivatives);
285 }
286
287
288
289
290
291
292
293 private double getTk(final AbsoluteDate date) {
294 final double cycleDuration = orbitalElements.getCycleDuration();
295
296 double tk = date.durationFrom(orbitalElements.getDate());
297
298 while (tk > 0.5 * cycleDuration) {
299 tk -= cycleDuration;
300 }
301 while (tk < -0.5 * cycleDuration) {
302 tk += cycleDuration;
303 }
304
305 return tk;
306 }
307
308
309 @Override
310 public Frame getFrame() {
311 return eci;
312 }
313
314
315 @Override
316 protected double getMass(final AbsoluteDate date) {
317 return getInitialState().getMass();
318 }
319
320
321 @Override
322 public void resetInitialState(final SpacecraftState state) {
323 orbitalElements = buildOrbitalElements(state, orbitalElements, ecef, getAttitudeProvider(), state.getMass());
324 final Orbit orbit = propagateOrbit(orbitalElements.getDate());
325 final Attitude attitude = getAttitudeProvider().getAttitude(orbit, orbit.getDate(), orbit.getFrame());
326 super.resetInitialState(new SpacecraftState(orbit, attitude).withMass(state.getMass()));
327 }
328
329
330 @Override
331 protected void resetIntermediateState(final SpacecraftState state, final boolean forward) {
332 resetInitialState(state);
333 }
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350 private static GNSSOrbitalElements<?> buildOrbitalElements(final SpacecraftState initialState,
351 final GNSSOrbitalElements<?> nonKeplerianElements,
352 final Frame ecef, final AttitudeProvider provider,
353 final double mass) {
354
355
356 final Frame frozenEcef = ecef.getFrozenFrame(initialState.getFrame(), initialState.getDate(), "frozen");
357 final KeplerianOrbit orbit = approximateInitialOrbit(initialState, nonKeplerianElements, frozenEcef);
358
359
360 final PVCoordinates targetPV = initialState.getPVCoordinates();
361 final FieldGnssOrbitalElements<Gradient, ?> gElements = convert(nonKeplerianElements, orbit);
362 for (int i = 0; i < MAX_ITER; ++i) {
363
364
365 final FieldGnssPropagator<Gradient> gPropagator =
366 new FieldGnssPropagator<>(gElements, initialState.getFrame(), ecef, provider,
367 gElements.getMu().newInstance(mass));
368 final FieldPVCoordinates<Gradient> gPV = gPropagator.getInitialState().getPVCoordinates();
369
370
371 final RealMatrix jacobian = MatrixUtils.createRealMatrix(FREE_PARAMETERS, FREE_PARAMETERS);
372 jacobian.setRow(0, gPV.getPosition().getX().getGradient());
373 jacobian.setRow(1, gPV.getPosition().getY().getGradient());
374 jacobian.setRow(2, gPV.getPosition().getZ().getGradient());
375 jacobian.setRow(3, gPV.getVelocity().getX().getGradient());
376 jacobian.setRow(4, gPV.getVelocity().getY().getGradient());
377 jacobian.setRow(5, gPV.getVelocity().getZ().getGradient());
378
379
380 final RealVector residuals = MatrixUtils.createRealVector(FREE_PARAMETERS);
381 residuals.setEntry(0, targetPV.getPosition().getX() - gPV.getPosition().getX().getValue());
382 residuals.setEntry(1, targetPV.getPosition().getY() - gPV.getPosition().getY().getValue());
383 residuals.setEntry(2, targetPV.getPosition().getZ() - gPV.getPosition().getZ().getValue());
384 residuals.setEntry(3, targetPV.getVelocity().getX() - gPV.getVelocity().getX().getValue());
385 residuals.setEntry(4, targetPV.getVelocity().getY() - gPV.getVelocity().getY().getValue());
386 residuals.setEntry(5, targetPV.getVelocity().getZ() - gPV.getVelocity().getZ().getValue());
387 final RealVector correction = new QRDecomposition(jacobian, EPS).getSolver().solve(residuals);
388
389
390 gElements.setSma(gElements.getSma().add(correction.getEntry(0)));
391 gElements.setE(gElements.getE().add(correction.getEntry(1)));
392 gElements.setI0(gElements.getI0().add(correction.getEntry(2)));
393 gElements.setPa(gElements.getPa().add(correction.getEntry(3)));
394 gElements.setOmega0(gElements.getOmega0().add(correction.getEntry(4)));
395 gElements.setM0(gElements.getM0().add(correction.getEntry(5)));
396
397 final double deltaP = FastMath.sqrt(residuals.getEntry(0) * residuals.getEntry(0) +
398 residuals.getEntry(1) * residuals.getEntry(1) +
399 residuals.getEntry(2) * residuals.getEntry(2));
400 final double deltaV = FastMath.sqrt(residuals.getEntry(3) * residuals.getEntry(3) +
401 residuals.getEntry(4) * residuals.getEntry(4) +
402 residuals.getEntry(5) * residuals.getEntry(5));
403
404 if (deltaP < TOL_P && deltaV < TOL_V) {
405 break;
406 }
407
408 }
409
410 return gElements.toNonField();
411
412 }
413
414
415
416
417
418
419
420
421 private static KeplerianOrbit approximateInitialOrbit(final SpacecraftState initialState,
422 final GNSSOrbitalElements<?> nonKeplerianElements,
423 final Frame frozenEcef) {
424
425
426
427 final PVCoordinates pv = initialState.getPVCoordinates(frozenEcef);
428 final Vector3D p = pv.getPosition();
429 final Vector3D v = pv.getVelocity();
430
431
432 final double rk = p.getNorm();
433
434
435 final Vector3D normal = pv.getMomentum().normalize();
436 final double cosIk = normal.getZ();
437 final double ik = Vector3D.angle(normal, Vector3D.PLUS_K);
438
439
440 final double q = FastMath.hypot(normal.getX(), normal.getY());
441 final double cos = -normal.getY() / q;
442 final double sin = normal.getX() / q;
443 final double xk = p.getX() * cos + p.getY() * sin;
444 final double yk = (p.getY() * cos - p.getX() * sin) / cosIk;
445
446
447 final double uk = FastMath.atan2(yk, xk);
448
449
450 double phi = uk;
451 for (int i = 0; i < MAX_ITER; ++i) {
452 final double previous = phi;
453 final SinCos cs2Phi = FastMath.sinCos(2 * phi);
454 phi = uk - (cs2Phi.cos() * nonKeplerianElements.getCuc() + cs2Phi.sin() * nonKeplerianElements.getCus());
455 if (FastMath.abs(phi - previous) <= EPS) {
456 break;
457 }
458 }
459 final SinCos cs2phi = FastMath.sinCos(2 * phi);
460
461
462
463 final double i0 = ik - (cs2phi.cos() * nonKeplerianElements.getCic() + cs2phi.sin() * nonKeplerianElements.getCis());
464 final double om0 = FastMath.atan2(sin, cos) +
465 nonKeplerianElements.getAngularVelocity() * nonKeplerianElements.getTime();
466
467
468 final double mu = initialState.getOrbit().getMu();
469 final double rV2OMu = rk * v.getNormSq() / mu;
470 final double sma = rk / (2 - rV2OMu);
471 final double eCosE = rV2OMu - 1;
472 final double eSinE = Vector3D.dotProduct(p, v) / FastMath.sqrt(mu * sma);
473 final double e = FastMath.hypot(eCosE, eSinE);
474 final double eccentricAnomaly = FastMath.atan2(eSinE, eCosE);
475 final double aop = phi - eccentricAnomaly;
476 final double meanAnomaly = KeplerianAnomalyUtility.ellipticEccentricToMean(e, eccentricAnomaly);
477
478 return new KeplerianOrbit(sma, e, i0, aop, om0, meanAnomaly, PositionAngleType.MEAN,
479 PositionAngleType.MEAN, frozenEcef,
480 initialState.getDate(), mu);
481
482 }
483
484
485
486
487
488
489
490 private static FieldGnssOrbitalElements<Gradient, ?> convert(final GNSSOrbitalElements<?> elements,
491 final KeplerianOrbit orbit) {
492
493 final Field<Gradient> field = GradientField.getField(FREE_PARAMETERS);
494 final FieldGnssOrbitalElements<Gradient, ?> gElements = elements.toField(field);
495
496
497 gElements.setSma(Gradient.variable(FREE_PARAMETERS, 0, orbit.getA()));
498 gElements.setE(Gradient.variable(FREE_PARAMETERS, 1, orbit.getE()));
499 gElements.setI0(Gradient.variable(FREE_PARAMETERS, 2, orbit.getI()));
500 gElements.setPa(Gradient.variable(FREE_PARAMETERS, 3, orbit.getPerigeeArgument()));
501 gElements.setOmega0(Gradient.variable(FREE_PARAMETERS, 4, orbit.getRightAscensionOfAscendingNode()));
502 gElements.setM0(Gradient.variable(FREE_PARAMETERS, 5, orbit.getMeanAnomaly()));
503
504 return gElements;
505
506 }
507
508 }