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 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 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.I_DOT_INDEX] = getMU().newInstance(orbitalElements.getIDot());
140 parameters[GNSSOrbitalElements.OMEGA_DOT_INDEX] = getMU().newInstance(orbitalElements.getOmegaDot());
141 parameters[GNSSOrbitalElements.CUC_INDEX] = getMU().newInstance(orbitalElements.getCuc());
142 parameters[GNSSOrbitalElements.CUS_INDEX] = getMU().newInstance(orbitalElements.getCus());
143 parameters[GNSSOrbitalElements.CRC_INDEX] = getMU().newInstance(orbitalElements.getCrc());
144 parameters[GNSSOrbitalElements.CRS_INDEX] = getMU().newInstance(orbitalElements.getCrs());
145 parameters[GNSSOrbitalElements.CIC_INDEX] = getMU().newInstance(orbitalElements.getCic());
146 parameters[GNSSOrbitalElements.CIS_INDEX] = getMU().newInstance(orbitalElements.getCis());
147 return parameters;
148 }
149
150
151 @Override
152 public List<ParameterDriver> getParametersDrivers() {
153 return orbitalElements.getParametersDrivers();
154 }
155
156
157
158
159
160
161 public Frame getECI() {
162 return eci;
163 }
164
165
166
167
168
169
170
171 public Frame getECEF() {
172 return ecef;
173 }
174
175
176
177
178
179
180 public T getMU() {
181 return orbitalElements.getMu();
182 }
183
184
185 @Override
186 public FieldOrbit<T> propagateOrbit(final FieldAbsoluteDate<T> date,
187 final T[] parameters) {
188
189 final FieldPVCoordinates<T> pvaInECEF = propagateInEcef(date, parameters);
190
191 final FieldPVCoordinates<T> pvaInECI = ecef.getTransformTo(eci, date).transformPVCoordinates(pvaInECEF);
192
193 return new FieldCartesianOrbit<>(pvaInECI, eci, date, getMU());
194 }
195
196
197
198
199
200
201
202
203
204
205
206 public FieldPVCoordinates<T> propagateInEcef(final FieldAbsoluteDate<T> date, final T[] parameters) {
207
208 final FieldUnivariateDerivative2<T> tk = new FieldUnivariateDerivative2<>(getTk(date),
209 date.getField().getOne(),
210 date.getField().getZero());
211
212
213 final FieldUnivariateDerivative2<T> ak = tk.multiply(orbitalElements.getADot()).add(orbitalElements.getSma());
214
215 final FieldUnivariateDerivative2<T> nA = tk.multiply(orbitalElements.getDeltaN0Dot().multiply(0.5)).
216 add(orbitalElements.getDeltaN0()).
217 add(orbitalElements.getMeanMotion0());
218
219 final FieldUnivariateDerivative2<T> mk = tk.multiply(nA).add(orbitalElements.getM0());
220
221 final FieldUnivariateDerivative2<T> e = tk.newInstance(orbitalElements.getE());
222 final FieldUnivariateDerivative2<T> ek = FieldKeplerianAnomalyUtility.ellipticMeanToEccentric(e, mk);
223
224 final FieldUnivariateDerivative2<T> vk = FieldKeplerianAnomalyUtility.ellipticEccentricToTrue(e, ek);
225
226 final FieldUnivariateDerivative2<T> phik = vk.add(orbitalElements.getPa());
227 final FieldSinCos<FieldUnivariateDerivative2<T>> cs2phi = FastMath.sinCos(phik.multiply(2));
228
229 final FieldUnivariateDerivative2<T> dphik = cs2phi.cos().multiply(parameters[GNSSOrbitalElements.CUC_INDEX]).
230 add(cs2phi.sin().multiply(parameters[GNSSOrbitalElements.CUS_INDEX]));
231
232 final FieldUnivariateDerivative2<T> drk = cs2phi.cos().multiply(parameters[GNSSOrbitalElements.CRC_INDEX]).
233 add(cs2phi.sin().multiply(parameters[GNSSOrbitalElements.CRS_INDEX]));
234
235 final FieldUnivariateDerivative2<T> dik = cs2phi.cos().multiply(parameters[GNSSOrbitalElements.CIC_INDEX]).
236 add(cs2phi.sin().multiply(parameters[GNSSOrbitalElements.CIS_INDEX]));
237
238 final FieldSinCos<FieldUnivariateDerivative2<T>> csuk = FastMath.sinCos(phik.add(dphik));
239
240 final FieldUnivariateDerivative2<T> rk = ek.cos().multiply(e.negate()).add(1).multiply(ak).add(drk);
241
242 final FieldUnivariateDerivative2<T> ik = tk.multiply(parameters[GNSSOrbitalElements.I_DOT_INDEX]).
243 add(orbitalElements.getI0()).add(dik);
244 final FieldSinCos<FieldUnivariateDerivative2<T>> csik = FastMath.sinCos(ik);
245
246 final FieldUnivariateDerivative2<T> xk = csuk.cos().multiply(rk);
247 final FieldUnivariateDerivative2<T> yk = csuk.sin().multiply(rk);
248
249 final FieldSinCos<FieldUnivariateDerivative2<T>> csomk =
250 FastMath.sinCos(tk.multiply(parameters[GNSSOrbitalElements.OMEGA_DOT_INDEX].
251 subtract(orbitalElements.getAngularVelocity())).
252 add(orbitalElements.getOmega0()).
253 subtract(parameters[GNSSOrbitalElements.TIME_INDEX].multiply(orbitalElements.getAngularVelocity())));
254
255 final FieldVector3D<FieldUnivariateDerivative2<T>> positionWithDerivatives =
256 new FieldVector3D<>(xk.multiply(csomk.cos()).subtract(yk.multiply(csomk.sin()).multiply(csik.cos())),
257 xk.multiply(csomk.sin()).add(yk.multiply(csomk.cos()).multiply(csik.cos())),
258 yk.multiply(csik.sin()));
259 return new FieldPVCoordinates<>(positionWithDerivatives);
260
261 }
262
263
264
265
266
267
268
269 private T getTk(final FieldAbsoluteDate<T> date) {
270
271 T tk = date.durationFrom(orbitalElements.getDate());
272
273 while (tk.getReal() > 0.5 * orbitalElements.getCycleDuration()) {
274 tk = tk.subtract(orbitalElements.getCycleDuration());
275 }
276 while (tk.getReal() < -0.5 * orbitalElements.getCycleDuration()) {
277 tk = tk.add(orbitalElements.getCycleDuration());
278 }
279
280 return tk;
281 }
282
283
284 @Override
285 public Frame getFrame() {
286 return eci;
287 }
288
289
290 @Override
291 protected T getMass(final FieldAbsoluteDate<T> date) {
292 return getInitialState().getMass();
293 }
294
295
296 @Override
297 public void resetInitialState(final FieldSpacecraftState<T> state) {
298 orbitalElements = buildOrbitalElements(state, orbitalElements, ecef, getAttitudeProvider(), state.getMass());
299 final FieldOrbit<T> orbit = propagateOrbit(orbitalElements.getDate(), defaultParameters());
300 final FieldAttitude<T> attitude = getAttitudeProvider().getAttitude(orbit, orbit.getDate(), orbit.getFrame());
301 super.resetInitialState(new FieldSpacecraftState<>(orbit, attitude).withMass(state.getMass()));
302 }
303
304
305 @Override
306 protected void resetIntermediateState(final FieldSpacecraftState<T> state, final boolean forward) {
307 resetInitialState(state);
308 }
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325 private static <T extends CalculusFieldElement<T>> FieldGnssOrbitalElements<T, ?>
326 buildOrbitalElements(final FieldSpacecraftState<T> initialState,
327 final FieldGnssOrbitalElements<T, ?> nonKeplerianElements,
328 final Frame ecef, final AttitudeProvider provider,
329 final T mass) {
330
331 final Field<T> field = initialState.getDate().getField();
332
333
334 final Frame frozenEcef = ecef.getFrozenFrame(initialState.getFrame(),
335 initialState.getDate().toAbsoluteDate(),
336 "frozen");
337 final FieldKeplerianOrbit<T> orbit = approximateInitialOrbit(initialState, nonKeplerianElements, frozenEcef);
338
339
340 final FieldPVCoordinates<T> targetPV = initialState.getPVCoordinates();
341 final FieldGnssOrbitalElements<FieldGradient<T>, ?> gElements = convert(nonKeplerianElements, orbit);
342 for (int i = 0; i < MAX_ITER; ++i) {
343
344
345 final FieldGnssPropagator<FieldGradient<T>> gPropagator =
346 new FieldGnssPropagator<>(gElements, initialState.getFrame(), ecef, provider,
347 gElements.getMu().newInstance(mass));
348 final FieldPVCoordinates<FieldGradient<T>> gPV = gPropagator.getInitialState().getPVCoordinates();
349
350
351 final FieldMatrix<T> jacobian = MatrixUtils.createFieldMatrix(field, FREE_PARAMETERS, FREE_PARAMETERS);
352 jacobian.setRow(0, gPV.getPosition().getX().getGradient());
353 jacobian.setRow(1, gPV.getPosition().getY().getGradient());
354 jacobian.setRow(2, gPV.getPosition().getZ().getGradient());
355 jacobian.setRow(3, gPV.getVelocity().getX().getGradient());
356 jacobian.setRow(4, gPV.getVelocity().getY().getGradient());
357 jacobian.setRow(5, gPV.getVelocity().getZ().getGradient());
358
359
360 final FieldVector<T> residuals = MatrixUtils.createFieldVector(field, FREE_PARAMETERS);
361 residuals.setEntry(0, targetPV.getPosition().getX().subtract(gPV.getPosition().getX().getValue()));
362 residuals.setEntry(1, targetPV.getPosition().getY().subtract(gPV.getPosition().getY().getValue()));
363 residuals.setEntry(2, targetPV.getPosition().getZ().subtract(gPV.getPosition().getZ().getValue()));
364 residuals.setEntry(3, targetPV.getVelocity().getX().subtract(gPV.getVelocity().getX().getValue()));
365 residuals.setEntry(4, targetPV.getVelocity().getY().subtract(gPV.getVelocity().getY().getValue()));
366 residuals.setEntry(5, targetPV.getVelocity().getZ().subtract(gPV.getVelocity().getZ().getValue()));
367 final FieldVector<T> correction = new FieldQRDecomposition<>(jacobian, field.getZero().newInstance(EPS)).
368 getSolver().
369 solve(residuals);
370
371
372 gElements.setSma(gElements.getSma().add(correction.getEntry(0)));
373 gElements.setE(gElements.getE().add(correction.getEntry(1)));
374 gElements.setI0(gElements.getI0().add(correction.getEntry(2)));
375 gElements.setPa(gElements.getPa().add(correction.getEntry(3)));
376 gElements.setOmega0(gElements.getOmega0().add(correction.getEntry(4)));
377 gElements.setM0(gElements.getM0().add(correction.getEntry(5)));
378
379 final double deltaP = FastMath.sqrt(residuals.getEntry(0).getReal() * residuals.getEntry(0).getReal() +
380 residuals.getEntry(1).getReal() * residuals.getEntry(1).getReal() +
381 residuals.getEntry(2).getReal() * residuals.getEntry(2).getReal());
382 final double deltaV = FastMath.sqrt(residuals.getEntry(3).getReal() * residuals.getEntry(3).getReal() +
383 residuals.getEntry(4).getReal() * residuals.getEntry(4).getReal() +
384 residuals.getEntry(5).getReal() * residuals.getEntry(5).getReal());
385
386 if (deltaP < TOL_P && deltaV < TOL_V) {
387 break;
388 }
389
390 }
391
392 return gElements.changeField(FieldGradient::getValue);
393
394 }
395
396
397
398
399
400
401
402
403 private static <T extends CalculusFieldElement<T>> FieldKeplerianOrbit<T>
404 approximateInitialOrbit(final FieldSpacecraftState<T> initialState,
405 final FieldGnssOrbitalElements<T, ?> nonKeplerianElements,
406 final Frame frozenEcef) {
407
408
409
410 final FieldPVCoordinates<T> pv = initialState.getPVCoordinates(frozenEcef);
411 final FieldVector3D<T> p = pv.getPosition();
412 final FieldVector3D<T> v = pv.getVelocity();
413
414
415 final T rk = p.getNorm();
416
417
418 final FieldVector3D<T> normal = pv.getMomentum().normalize();
419 final T cosIk = normal.getZ();
420 final T ik = FieldVector3D.angle(normal, Vector3D.PLUS_K);
421
422
423 final T q = FastMath.hypot(normal.getX(), normal.getY());
424 final T cos = normal.getY().negate().divide(q);
425 final T sin = normal.getX().divide(q);
426 final T xk = p.getX().multiply(cos).add(p.getY().multiply(sin));
427 final T yk = p.getY().multiply(cos).subtract(p.getX().multiply(sin)).divide(cosIk);
428
429
430 final T uk = FastMath.atan2(yk, xk);
431
432
433 T phi = uk;
434 for (int i = 0; i < MAX_ITER; ++i) {
435 final T previous = phi;
436 final FieldSinCos<T> cs2Phi = FastMath.sinCos(phi.multiply(2));
437 phi = uk.subtract(cs2Phi.cos().multiply(nonKeplerianElements.getCuc()).add(cs2Phi.sin().multiply(nonKeplerianElements.getCus())));
438 if (FastMath.abs(phi.subtract(previous).getReal()) <= EPS) {
439 break;
440 }
441 }
442 final FieldSinCos<T> cs2phi = FastMath.sinCos(phi.multiply(2));
443
444
445
446 final T i0 = ik.subtract(cs2phi.cos().multiply(nonKeplerianElements.getCic()).add(cs2phi.sin().multiply(nonKeplerianElements.getCis())));
447 final T om0 = FastMath.atan2(sin, cos).
448 add(nonKeplerianElements.getAngularVelocity() * nonKeplerianElements.getTime());
449
450
451 final T mu = initialState.getOrbit().getMu();
452 final T rV2OMu = rk.multiply(v.getNormSq()).divide(mu);
453 final T sma = rk.divide(rV2OMu.negate().add(2));
454 final T eCosE = rV2OMu.subtract(1);
455 final T eSinE = FieldVector3D.dotProduct(p, v).divide(FastMath.sqrt(mu.multiply(sma)));
456 final T e = FastMath.hypot(eCosE, eSinE);
457 final T eccentricAnomaly = FastMath.atan2(eSinE, eCosE);
458 final T aop = phi.subtract(eccentricAnomaly);
459 final T meanAnomaly = FieldKeplerianAnomalyUtility.ellipticEccentricToMean(e, eccentricAnomaly);
460
461 return new FieldKeplerianOrbit<>(sma, e, i0, aop, om0, meanAnomaly, PositionAngleType.MEAN,
462 PositionAngleType.MEAN, frozenEcef,
463 initialState.getDate(), mu);
464
465 }
466
467
468
469
470
471
472
473 private static <T extends CalculusFieldElement<T>> FieldGnssOrbitalElements<FieldGradient<T>, ?>
474 convert(final FieldGnssOrbitalElements<T, ?> elements, final FieldKeplerianOrbit<T> orbit) {
475
476 final FieldGnssOrbitalElements<FieldGradient<T>, ?> gElements =
477 elements.changeField(t -> FieldGradient.constant(FREE_PARAMETERS, t));
478
479
480 gElements.setSma(FieldGradient.variable(FREE_PARAMETERS, 0, orbit.getA()));
481 gElements.setE(FieldGradient.variable(FREE_PARAMETERS, 1, orbit.getE()));
482 gElements.setI0(FieldGradient.variable(FREE_PARAMETERS, 2, orbit.getI()));
483 gElements.setPa(FieldGradient.variable(FREE_PARAMETERS, 3, orbit.getPerigeeArgument()));
484 gElements.setOmega0(FieldGradient.variable(FREE_PARAMETERS, 4, orbit.getRightAscensionOfAscendingNode()));
485 gElements.setM0(FieldGradient.variable(FREE_PARAMETERS, 5, orbit.getMeanAnomaly()));
486
487 return gElements;
488
489 }
490
491 }