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