1 /* Copyright 2002-2019 CS Systèmes d'Information 2 * Licensed to CS Systèmes d'Information (CS) under one or more 3 * contributor license agreements. See the NOTICE file distributed with 4 * this work for additional information regarding copyright ownership. 5 * CS licenses this file to You under the Apache License, Version 2.0 6 * (the "License"); you may not use this file except in compliance with 7 * the License. You may obtain a copy of the License at 8 * 9 * http://www.apache.org/licenses/LICENSE-2.0 10 * 11 * Unless required by applicable law or agreed to in writing, software 12 * distributed under the License is distributed on an "AS IS" BASIS, 13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. 14 * See the License for the specific language governing permissions and 15 * limitations under the License. 16 */ 17 package org.orekit.orbits; 18 19 import java.io.Serializable; 20 21 import org.hipparchus.geometry.euclidean.threed.Vector3D; 22 import org.hipparchus.linear.DecompositionSolver; 23 import org.hipparchus.linear.MatrixUtils; 24 import org.hipparchus.linear.QRDecomposition; 25 import org.hipparchus.linear.RealMatrix; 26 import org.hipparchus.util.FastMath; 27 import org.hipparchus.util.MathArrays; 28 import org.orekit.errors.OrekitIllegalArgumentException; 29 import org.orekit.errors.OrekitInternalError; 30 import org.orekit.errors.OrekitMessages; 31 import org.orekit.frames.Frame; 32 import org.orekit.frames.Transform; 33 import org.orekit.time.AbsoluteDate; 34 import org.orekit.time.TimeInterpolable; 35 import org.orekit.time.TimeShiftable; 36 import org.orekit.time.TimeStamped; 37 import org.orekit.utils.PVCoordinates; 38 import org.orekit.utils.PVCoordinatesProvider; 39 import org.orekit.utils.TimeStampedPVCoordinates; 40 41 /** 42 * This class handles orbital parameters. 43 44 * <p> 45 * For user convenience, both the Cartesian and the equinoctial elements 46 * are provided by this class, regardless of the canonical representation 47 * implemented in the derived class (which may be classical Keplerian 48 * elements for example). 49 * </p> 50 * <p> 51 * The parameters are defined in a frame specified by the user. It is important 52 * to make sure this frame is consistent: it probably is inertial and centered 53 * on the central body. This information is used for example by some 54 * force models. 55 * </p> 56 * <p> 57 * Instance of this class are guaranteed to be immutable. 58 * </p> 59 * @author Luc Maisonobe 60 * @author Guylaine Prat 61 * @author Fabien Maussion 62 * @author Véronique Pommier-Maurussane 63 */ 64 public abstract class Orbit 65 implements TimeStamped, TimeShiftable<Orbit>, TimeInterpolable<Orbit>, 66 Serializable, PVCoordinatesProvider { 67 68 /** Serializable UID. */ 69 private static final long serialVersionUID = 438733454597999578L; 70 71 /** Frame in which are defined the orbital parameters. */ 72 private final Frame frame; 73 74 /** Date of the orbital parameters. */ 75 private final AbsoluteDate date; 76 77 /** Value of mu used to compute position and velocity (m³/s²). */ 78 private final double mu; 79 80 /** Computed PVCoordinates. */ 81 private transient TimeStampedPVCoordinates pvCoordinates; 82 83 /** Jacobian of the orbital parameters with mean angle with respect to the Cartesian coordinates. */ 84 private transient double[][] jacobianMeanWrtCartesian; 85 86 /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with mean angle. */ 87 private transient double[][] jacobianWrtParametersMean; 88 89 /** Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian coordinates. */ 90 private transient double[][] jacobianEccentricWrtCartesian; 91 92 /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with eccentric angle. */ 93 private transient double[][] jacobianWrtParametersEccentric; 94 95 /** Jacobian of the orbital parameters with true angle with respect to the Cartesian coordinates. */ 96 private transient double[][] jacobianTrueWrtCartesian; 97 98 /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with true angle. */ 99 private transient double[][] jacobianWrtParametersTrue; 100 101 /** Default constructor. 102 * Build a new instance with arbitrary default elements. 103 * @param frame the frame in which the parameters are defined 104 * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame}) 105 * @param date date of the orbital parameters 106 * @param mu central attraction coefficient (m^3/s^2) 107 * @exception IllegalArgumentException if frame is not a {@link 108 * Frame#isPseudoInertial pseudo-inertial frame} 109 */ 110 protected Orbit(final Frame frame, final AbsoluteDate date, final double mu) 111 throws IllegalArgumentException { 112 ensurePseudoInertialFrame(frame); 113 this.date = date; 114 this.mu = mu; 115 this.pvCoordinates = null; 116 this.frame = frame; 117 jacobianMeanWrtCartesian = null; 118 jacobianWrtParametersMean = null; 119 jacobianEccentricWrtCartesian = null; 120 jacobianWrtParametersEccentric = null; 121 jacobianTrueWrtCartesian = null; 122 jacobianWrtParametersTrue = null; 123 } 124 125 /** Set the orbit from Cartesian parameters. 126 * 127 * <p> The acceleration provided in {@code pvCoordinates} is accessible using 128 * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods 129 * use {@code mu} and the position to compute the acceleration, including 130 * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}. 131 * 132 * @param pvCoordinates the position and velocity in the inertial frame 133 * @param frame the frame in which the {@link TimeStampedPVCoordinates} are defined 134 * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame}) 135 * @param mu central attraction coefficient (m^3/s^2) 136 * @exception IllegalArgumentException if frame is not a {@link 137 * Frame#isPseudoInertial pseudo-inertial frame} 138 */ 139 protected Orbit(final TimeStampedPVCoordinates pvCoordinates, final Frame frame, final double mu) 140 throws IllegalArgumentException { 141 ensurePseudoInertialFrame(frame); 142 this.date = pvCoordinates.getDate(); 143 this.mu = mu; 144 if (pvCoordinates.getAcceleration().getNormSq() == 0) { 145 // the acceleration was not provided, 146 // compute it from Newtonian attraction 147 final double r2 = pvCoordinates.getPosition().getNormSq(); 148 final double r3 = r2 * FastMath.sqrt(r2); 149 this.pvCoordinates = new TimeStampedPVCoordinates(pvCoordinates.getDate(), 150 pvCoordinates.getPosition(), 151 pvCoordinates.getVelocity(), 152 new Vector3D(-mu / r3, pvCoordinates.getPosition())); 153 } else { 154 this.pvCoordinates = pvCoordinates; 155 } 156 this.frame = frame; 157 } 158 159 /** Check if Cartesian coordinates include non-Keplerian acceleration. 160 * @param pva Cartesian coordinates 161 * @param mu central attraction coefficient 162 * @return true if Cartesian coordinates include non-Keplerian acceleration 163 */ 164 protected static boolean hasNonKeplerianAcceleration(final PVCoordinates pva, final double mu) { 165 166 final Vector3D a = pva.getAcceleration(); 167 if (a == null) { 168 return false; 169 } 170 171 final Vector3D p = pva.getPosition(); 172 final double r2 = p.getNormSq(); 173 final double r = FastMath.sqrt(r2); 174 final Vector3D keplerianAcceleration = new Vector3D(-mu / (r * r2), p); 175 if (a.getNorm() > 1.0e-9 * keplerianAcceleration.getNorm()) { 176 // we have a relevant acceleration, we can compute derivatives 177 return true; 178 } else { 179 // the provided acceleration is either too small to be reliable (probably even 0), or NaN 180 return false; 181 } 182 183 } 184 185 /** Get the orbit type. 186 * @return orbit type 187 */ 188 public abstract OrbitType getType(); 189 190 /** Ensure the defining frame is a pseudo-inertial frame. 191 * @param frame frame to check 192 * @exception IllegalArgumentException if frame is not a {@link 193 * Frame#isPseudoInertial pseudo-inertial frame} 194 */ 195 private static void ensurePseudoInertialFrame(final Frame frame) 196 throws IllegalArgumentException { 197 if (!frame.isPseudoInertial()) { 198 throw new OrekitIllegalArgumentException(OrekitMessages.NON_PSEUDO_INERTIAL_FRAME, 199 frame.getName()); 200 } 201 } 202 203 /** Get the frame in which the orbital parameters are defined. 204 * @return frame in which the orbital parameters are defined 205 */ 206 public Frame getFrame() { 207 return frame; 208 } 209 210 /** Get the semi-major axis. 211 * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p> 212 * @return semi-major axis (m) 213 */ 214 public abstract double getA(); 215 216 /** Get the semi-major axis derivative. 217 * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p> 218 * <p> 219 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 220 * </p> 221 * @return semi-major axis derivative (m/s) 222 * @see #hasDerivatives() 223 * @since 9.0 224 */ 225 public abstract double getADot(); 226 227 /** Get the first component of the equinoctial eccentricity vector derivative. 228 * @return first component of the equinoctial eccentricity vector derivative 229 */ 230 public abstract double getEquinoctialEx(); 231 232 /** Get the first component of the equinoctial eccentricity vector. 233 * <p> 234 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 235 * </p> 236 * @return first component of the equinoctial eccentricity vector 237 * @see #hasDerivatives() 238 * @since 9.0 239 */ 240 public abstract double getEquinoctialExDot(); 241 242 /** Get the second component of the equinoctial eccentricity vector derivative. 243 * @return second component of the equinoctial eccentricity vector derivative 244 */ 245 public abstract double getEquinoctialEy(); 246 247 /** Get the second component of the equinoctial eccentricity vector. 248 * <p> 249 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 250 * </p> 251 * @return second component of the equinoctial eccentricity vector 252 * @see #hasDerivatives() 253 * @since 9.0 254 */ 255 public abstract double getEquinoctialEyDot(); 256 257 /** Get the first component of the inclination vector. 258 * @return first component of the inclination vector 259 */ 260 public abstract double getHx(); 261 262 /** Get the first component of the inclination vector derivative. 263 * <p> 264 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 265 * </p> 266 * @return first component of the inclination vector derivative 267 * @see #hasDerivatives() 268 * @since 9.0 269 */ 270 public abstract double getHxDot(); 271 272 /** Get the second component of the inclination vector. 273 * @return second component of the inclination vector 274 */ 275 public abstract double getHy(); 276 277 /** Get the second component of the inclination vector derivative. 278 * <p> 279 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 280 * </p> 281 * @return second component of the inclination vector derivative 282 * @see #hasDerivatives() 283 * @since 9.0 284 */ 285 public abstract double getHyDot(); 286 287 /** Get the eccentric longitude argument. 288 * @return E + ω + Ω eccentric longitude argument (rad) 289 */ 290 public abstract double getLE(); 291 292 /** Get the eccentric longitude argument derivative. 293 * <p> 294 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 295 * </p> 296 * @return d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s) 297 * @see #hasDerivatives() 298 * @since 9.0 299 */ 300 public abstract double getLEDot(); 301 302 /** Get the true longitude argument. 303 * @return v + ω + Ω true longitude argument (rad) 304 */ 305 public abstract double getLv(); 306 307 /** Get the true longitude argument derivative. 308 * <p> 309 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 310 * </p> 311 * @return d(v + ω + Ω)/dt true longitude argument derivative (rad/s) 312 * @see #hasDerivatives() 313 * @since 9.0 314 */ 315 public abstract double getLvDot(); 316 317 /** Get the mean longitude argument. 318 * @return M + ω + Ω mean longitude argument (rad) 319 */ 320 public abstract double getLM(); 321 322 /** Get the mean longitude argument derivative. 323 * <p> 324 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 325 * </p> 326 * @return d(M + ω + Ω)/dt mean longitude argument derivative (rad/s) 327 * @see #hasDerivatives() 328 * @since 9.0 329 */ 330 public abstract double getLMDot(); 331 332 // Additional orbital elements 333 334 /** Get the eccentricity. 335 * @return eccentricity 336 */ 337 public abstract double getE(); 338 339 /** Get the eccentricity derivative. 340 * <p> 341 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 342 * </p> 343 * @return eccentricity derivative 344 * @see #hasDerivatives() 345 * @since 9.0 346 */ 347 public abstract double getEDot(); 348 349 /** Get the inclination. 350 * @return inclination (rad) 351 */ 352 public abstract double getI(); 353 354 /** Get the inclination derivative. 355 * <p> 356 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 357 * </p> 358 * @return inclination derivative (rad/s) 359 * @see #hasDerivatives() 360 * @since 9.0 361 */ 362 public abstract double getIDot(); 363 364 /** Check if orbit includes derivatives. 365 * @return true if orbit includes derivatives 366 * @see #getADot() 367 * @see #getEquinoctialExDot() 368 * @see #getEquinoctialEyDot() 369 * @see #getHxDot() 370 * @see #getHyDot() 371 * @see #getLEDot() 372 * @see #getLvDot() 373 * @see #getLMDot() 374 * @see #getEDot() 375 * @see #getIDot() 376 * @since 9.0 377 */ 378 public boolean hasDerivatives() { 379 return !Double.isNaN(getADot()); 380 } 381 382 /** Get the central acceleration constant. 383 * @return central acceleration constant 384 */ 385 public double getMu() { 386 return mu; 387 } 388 389 /** Get the Keplerian period. 390 * <p>The Keplerian period is computed directly from semi major axis 391 * and central acceleration constant.</p> 392 * @return Keplerian period in seconds, or positive infinity for hyperbolic orbits 393 */ 394 public double getKeplerianPeriod() { 395 final double a = getA(); 396 return (a < 0) ? Double.POSITIVE_INFINITY : 2.0 * FastMath.PI * a * FastMath.sqrt(a / mu); 397 } 398 399 /** Get the Keplerian mean motion. 400 * <p>The Keplerian mean motion is computed directly from semi major axis 401 * and central acceleration constant.</p> 402 * @return Keplerian mean motion in radians per second 403 */ 404 public double getKeplerianMeanMotion() { 405 final double absA = FastMath.abs(getA()); 406 return FastMath.sqrt(mu / absA) / absA; 407 } 408 409 /** Get the date of orbital parameters. 410 * @return date of the orbital parameters 411 */ 412 public AbsoluteDate getDate() { 413 return date; 414 } 415 416 /** Get the {@link TimeStampedPVCoordinates} in a specified frame. 417 * @param outputFrame frame in which the position/velocity coordinates shall be computed 418 * @return pvCoordinates in the specified output frame 419 * @see #getPVCoordinates() 420 */ 421 public TimeStampedPVCoordinates getPVCoordinates(final Frame outputFrame) { 422 if (pvCoordinates == null) { 423 pvCoordinates = initPVCoordinates(); 424 } 425 426 // If output frame requested is the same as definition frame, 427 // PV coordinates are returned directly 428 if (outputFrame == frame) { 429 return pvCoordinates; 430 } 431 432 // Else, PV coordinates are transformed to output frame 433 final Transform t = frame.getTransformTo(outputFrame, date); 434 return t.transformPVCoordinates(pvCoordinates); 435 } 436 437 /** {@inheritDoc} */ 438 public TimeStampedPVCoordinates getPVCoordinates(final AbsoluteDate otherDate, final Frame otherFrame) { 439 return shiftedBy(otherDate.durationFrom(getDate())).getPVCoordinates(otherFrame); 440 } 441 442 443 /** Get the {@link TimeStampedPVCoordinates} in definition frame. 444 * @return pvCoordinates in the definition frame 445 * @see #getPVCoordinates(Frame) 446 */ 447 public TimeStampedPVCoordinates getPVCoordinates() { 448 if (pvCoordinates == null) { 449 pvCoordinates = initPVCoordinates(); 450 } 451 return pvCoordinates; 452 } 453 454 /** Compute the position/velocity coordinates from the canonical parameters. 455 * @return computed position/velocity coordinates 456 */ 457 protected abstract TimeStampedPVCoordinates initPVCoordinates(); 458 459 /** Get a time-shifted orbit. 460 * <p> 461 * The orbit can be slightly shifted to close dates. The shifting model is a 462 * Keplerian one if no derivatives are available in the orbit, or Keplerian 463 * plus quadratic effect of the non-Keplerian acceleration if derivatives are 464 * available. Shifting is <em>not</em> intended as a replacement for proper 465 * orbit propagation but should be sufficient for small time shifts or coarse 466 * accuracy. 467 * </p> 468 * @param dt time shift in seconds 469 * @return a new orbit, shifted with respect to the instance (which is immutable) 470 */ 471 public abstract Orbit shiftedBy(double dt); 472 473 /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters. 474 * <p> 475 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with 476 * respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter 477 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot. 478 * </p> 479 * @param type type of the position angle to use 480 * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix 481 * is larger than 6x6, only the 6x6 upper left corner will be modified 482 */ 483 public void getJacobianWrtCartesian(final PositionAngle type, final double[][] jacobian) { 484 485 final double[][] cachedJacobian; 486 synchronized (this) { 487 switch (type) { 488 case MEAN : 489 if (jacobianMeanWrtCartesian == null) { 490 // first call, we need to compute the Jacobian and cache it 491 jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian(); 492 } 493 cachedJacobian = jacobianMeanWrtCartesian; 494 break; 495 case ECCENTRIC : 496 if (jacobianEccentricWrtCartesian == null) { 497 // first call, we need to compute the Jacobian and cache it 498 jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian(); 499 } 500 cachedJacobian = jacobianEccentricWrtCartesian; 501 break; 502 case TRUE : 503 if (jacobianTrueWrtCartesian == null) { 504 // first call, we need to compute the Jacobian and cache it 505 jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian(); 506 } 507 cachedJacobian = jacobianTrueWrtCartesian; 508 break; 509 default : 510 throw new OrekitInternalError(null); 511 } 512 } 513 514 // fill the user provided array 515 for (int i = 0; i < cachedJacobian.length; ++i) { 516 System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length); 517 } 518 519 } 520 521 /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters. 522 * <p> 523 * Element {@code jacobian[i][j]} is the derivative of Cartesian coordinate i of the orbit with 524 * respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate 525 * x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters. 526 * </p> 527 * @param type type of the position angle to use 528 * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix 529 * is larger than 6x6, only the 6x6 upper left corner will be modified 530 */ 531 public void getJacobianWrtParameters(final PositionAngle type, final double[][] jacobian) { 532 533 final double[][] cachedJacobian; 534 synchronized (this) { 535 switch (type) { 536 case MEAN : 537 if (jacobianWrtParametersMean == null) { 538 // first call, we need to compute the Jacobian and cache it 539 jacobianWrtParametersMean = createInverseJacobian(type); 540 } 541 cachedJacobian = jacobianWrtParametersMean; 542 break; 543 case ECCENTRIC : 544 if (jacobianWrtParametersEccentric == null) { 545 // first call, we need to compute the Jacobian and cache it 546 jacobianWrtParametersEccentric = createInverseJacobian(type); 547 } 548 cachedJacobian = jacobianWrtParametersEccentric; 549 break; 550 case TRUE : 551 if (jacobianWrtParametersTrue == null) { 552 // first call, we need to compute the Jacobian and cache it 553 jacobianWrtParametersTrue = createInverseJacobian(type); 554 } 555 cachedJacobian = jacobianWrtParametersTrue; 556 break; 557 default : 558 throw new OrekitInternalError(null); 559 } 560 } 561 562 // fill the user-provided array 563 for (int i = 0; i < cachedJacobian.length; ++i) { 564 System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length); 565 } 566 567 } 568 569 /** Create an inverse Jacobian. 570 * @param type type of the position angle to use 571 * @return inverse Jacobian 572 */ 573 private double[][] createInverseJacobian(final PositionAngle type) { 574 575 // get the direct Jacobian 576 final double[][] directJacobian = new double[6][6]; 577 getJacobianWrtCartesian(type, directJacobian); 578 579 // invert the direct Jacobian 580 final RealMatrix matrix = MatrixUtils.createRealMatrix(directJacobian); 581 final DecompositionSolver solver = new QRDecomposition(matrix).getSolver(); 582 return solver.getInverse().getData(); 583 584 } 585 586 /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters. 587 * <p> 588 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with 589 * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter 590 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot. 591 * </p> 592 * @return 6x6 Jacobian matrix 593 * @see #computeJacobianEccentricWrtCartesian() 594 * @see #computeJacobianTrueWrtCartesian() 595 */ 596 protected abstract double[][] computeJacobianMeanWrtCartesian(); 597 598 /** Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters. 599 * <p> 600 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with 601 * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter 602 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot. 603 * </p> 604 * @return 6x6 Jacobian matrix 605 * @see #computeJacobianMeanWrtCartesian() 606 * @see #computeJacobianTrueWrtCartesian() 607 */ 608 protected abstract double[][] computeJacobianEccentricWrtCartesian(); 609 610 /** Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters. 611 * <p> 612 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with 613 * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter 614 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot. 615 * </p> 616 * @return 6x6 Jacobian matrix 617 * @see #computeJacobianMeanWrtCartesian() 618 * @see #computeJacobianEccentricWrtCartesian() 619 */ 620 protected abstract double[][] computeJacobianTrueWrtCartesian(); 621 622 /** Add the contribution of the Keplerian motion to parameters derivatives 623 * <p> 624 * This method is used by integration-based propagators to evaluate the part of Keplerian 625 * motion to evolution of the orbital state. 626 * </p> 627 * @param type type of the position angle in the state 628 * @param gm attraction coefficient to use 629 * @param pDot array containing orbital state derivatives to update (the Keplerian 630 * part must be <em>added</em> to the array components, as the array may already 631 * contain some non-zero elements corresponding to non-Keplerian parts) 632 */ 633 public abstract void addKeplerContribution(PositionAngle type, double gm, double[] pDot); 634 635 /** Fill a Jacobian half row with a single vector. 636 * @param a coefficient of the vector 637 * @param v vector 638 * @param row Jacobian matrix row 639 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 640 */ 641 protected static void fillHalfRow(final double a, final Vector3D v, final double[] row, final int j) { 642 row[j] = a * v.getX(); 643 row[j + 1] = a * v.getY(); 644 row[j + 2] = a * v.getZ(); 645 } 646 647 /** Fill a Jacobian half row with a linear combination of vectors. 648 * @param a1 coefficient of the first vector 649 * @param v1 first vector 650 * @param a2 coefficient of the second vector 651 * @param v2 second vector 652 * @param row Jacobian matrix row 653 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 654 */ 655 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2, 656 final double[] row, final int j) { 657 row[j] = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX()); 658 row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY()); 659 row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ()); 660 } 661 662 /** Fill a Jacobian half row with a linear combination of vectors. 663 * @param a1 coefficient of the first vector 664 * @param v1 first vector 665 * @param a2 coefficient of the second vector 666 * @param v2 second vector 667 * @param a3 coefficient of the third vector 668 * @param v3 third vector 669 * @param row Jacobian matrix row 670 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 671 */ 672 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2, 673 final double a3, final Vector3D v3, 674 final double[] row, final int j) { 675 row[j] = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX()); 676 row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY()); 677 row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ()); 678 } 679 680 /** Fill a Jacobian half row with a linear combination of vectors. 681 * @param a1 coefficient of the first vector 682 * @param v1 first vector 683 * @param a2 coefficient of the second vector 684 * @param v2 second vector 685 * @param a3 coefficient of the third vector 686 * @param v3 third vector 687 * @param a4 coefficient of the fourth vector 688 * @param v4 fourth vector 689 * @param row Jacobian matrix row 690 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 691 */ 692 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2, 693 final double a3, final Vector3D v3, final double a4, final Vector3D v4, 694 final double[] row, final int j) { 695 row[j] = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX(), a4, v4.getX()); 696 row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY(), a4, v4.getY()); 697 row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ(), a4, v4.getZ()); 698 } 699 700 /** Fill a Jacobian half row with a linear combination of vectors. 701 * @param a1 coefficient of the first vector 702 * @param v1 first vector 703 * @param a2 coefficient of the second vector 704 * @param v2 second vector 705 * @param a3 coefficient of the third vector 706 * @param v3 third vector 707 * @param a4 coefficient of the fourth vector 708 * @param v4 fourth vector 709 * @param a5 coefficient of the fifth vector 710 * @param v5 fifth vector 711 * @param row Jacobian matrix row 712 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 713 */ 714 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2, 715 final double a3, final Vector3D v3, final double a4, final Vector3D v4, 716 final double a5, final Vector3D v5, 717 final double[] row, final int j) { 718 final double[] a = new double[] { 719 a1, a2, a3, a4, a5 720 }; 721 row[j] = MathArrays.linearCombination(a, new double[] { 722 v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX() 723 }); 724 row[j + 1] = MathArrays.linearCombination(a, new double[] { 725 v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY() 726 }); 727 row[j + 2] = MathArrays.linearCombination(a, new double[] { 728 v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ() 729 }); 730 } 731 732 /** Fill a Jacobian half row with a linear combination of vectors. 733 * @param a1 coefficient of the first vector 734 * @param v1 first vector 735 * @param a2 coefficient of the second vector 736 * @param v2 second vector 737 * @param a3 coefficient of the third vector 738 * @param v3 third vector 739 * @param a4 coefficient of the fourth vector 740 * @param v4 fourth vector 741 * @param a5 coefficient of the fifth vector 742 * @param v5 fifth vector 743 * @param a6 coefficient of the sixth vector 744 * @param v6 sixth vector 745 * @param row Jacobian matrix row 746 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 747 */ 748 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2, 749 final double a3, final Vector3D v3, final double a4, final Vector3D v4, 750 final double a5, final Vector3D v5, final double a6, final Vector3D v6, 751 final double[] row, final int j) { 752 final double[] a = new double[] { 753 a1, a2, a3, a4, a5, a6 754 }; 755 row[j] = MathArrays.linearCombination(a, new double[] { 756 v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX(), v6.getX() 757 }); 758 row[j + 1] = MathArrays.linearCombination(a, new double[] { 759 v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY(), v6.getY() 760 }); 761 row[j + 2] = MathArrays.linearCombination(a, new double[] { 762 v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ(), v6.getZ() 763 }); 764 } 765 766 }