1 /* Copyright 2002-2022 CS GROUP 2 * Licensed to CS GROUP (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 p = pva.getPosition(); 167 final double r2 = p.getNormSq(); 168 final double r = FastMath.sqrt(r2); 169 final Vector3D keplerianAcceleration = new Vector3D(-mu / (r * r2), p); 170 171 // Check if acceleration is null or relatively close to 0 compared to the keplerain acceleration 172 final Vector3D a = pva.getAcceleration(); 173 if (a == null || a.getNorm() < 1e-9 * keplerianAcceleration.getNorm()) { 174 return false; 175 } 176 177 final Vector3D nonKeplerianAcceleration = a.subtract(keplerianAcceleration); 178 179 if ( nonKeplerianAcceleration.getNorm() > 1e-9 * keplerianAcceleration.getNorm()) { 180 // we have a relevant acceleration, we can compute derivatives 181 return true; 182 } else { 183 // the provided acceleration is either too small to be reliable (probably even 0), or NaN 184 return false; 185 } 186 } 187 188 /** Get the orbit type. 189 * @return orbit type 190 */ 191 public abstract OrbitType getType(); 192 193 /** Ensure the defining frame is a pseudo-inertial frame. 194 * @param frame frame to check 195 * @exception IllegalArgumentException if frame is not a {@link 196 * Frame#isPseudoInertial pseudo-inertial frame} 197 */ 198 private static void ensurePseudoInertialFrame(final Frame frame) 199 throws IllegalArgumentException { 200 if (!frame.isPseudoInertial()) { 201 throw new OrekitIllegalArgumentException(OrekitMessages.NON_PSEUDO_INERTIAL_FRAME, 202 frame.getName()); 203 } 204 } 205 206 /** Get the frame in which the orbital parameters are defined. 207 * @return frame in which the orbital parameters are defined 208 */ 209 public Frame getFrame() { 210 return frame; 211 } 212 213 /** Get the semi-major axis. 214 * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p> 215 * @return semi-major axis (m) 216 */ 217 public abstract double getA(); 218 219 /** Get the semi-major axis derivative. 220 * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p> 221 * <p> 222 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 223 * </p> 224 * @return semi-major axis derivative (m/s) 225 * @see #hasDerivatives() 226 * @since 9.0 227 */ 228 public abstract double getADot(); 229 230 /** Get the first component of the equinoctial eccentricity vector derivative. 231 * @return first component of the equinoctial eccentricity vector derivative 232 */ 233 public abstract double getEquinoctialEx(); 234 235 /** Get the first component of the equinoctial eccentricity vector. 236 * <p> 237 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 238 * </p> 239 * @return first component of the equinoctial eccentricity vector 240 * @see #hasDerivatives() 241 * @since 9.0 242 */ 243 public abstract double getEquinoctialExDot(); 244 245 /** Get the second component of the equinoctial eccentricity vector derivative. 246 * @return second component of the equinoctial eccentricity vector derivative 247 */ 248 public abstract double getEquinoctialEy(); 249 250 /** Get the second component of the equinoctial eccentricity vector. 251 * <p> 252 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 253 * </p> 254 * @return second component of the equinoctial eccentricity vector 255 * @see #hasDerivatives() 256 * @since 9.0 257 */ 258 public abstract double getEquinoctialEyDot(); 259 260 /** Get the first component of the inclination vector. 261 * @return first component of the inclination vector 262 */ 263 public abstract double getHx(); 264 265 /** Get the first component of the inclination vector derivative. 266 * <p> 267 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 268 * </p> 269 * @return first component of the inclination vector derivative 270 * @see #hasDerivatives() 271 * @since 9.0 272 */ 273 public abstract double getHxDot(); 274 275 /** Get the second component of the inclination vector. 276 * @return second component of the inclination vector 277 */ 278 public abstract double getHy(); 279 280 /** Get the second component of the inclination vector derivative. 281 * <p> 282 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 283 * </p> 284 * @return second component of the inclination vector derivative 285 * @see #hasDerivatives() 286 * @since 9.0 287 */ 288 public abstract double getHyDot(); 289 290 /** Get the eccentric longitude argument. 291 * @return E + ω + Ω eccentric longitude argument (rad) 292 */ 293 public abstract double getLE(); 294 295 /** Get the eccentric longitude argument derivative. 296 * <p> 297 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 298 * </p> 299 * @return d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s) 300 * @see #hasDerivatives() 301 * @since 9.0 302 */ 303 public abstract double getLEDot(); 304 305 /** Get the true longitude argument. 306 * @return v + ω + Ω true longitude argument (rad) 307 */ 308 public abstract double getLv(); 309 310 /** Get the true longitude argument derivative. 311 * <p> 312 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 313 * </p> 314 * @return d(v + ω + Ω)/dt true longitude argument derivative (rad/s) 315 * @see #hasDerivatives() 316 * @since 9.0 317 */ 318 public abstract double getLvDot(); 319 320 /** Get the mean longitude argument. 321 * @return M + ω + Ω mean longitude argument (rad) 322 */ 323 public abstract double getLM(); 324 325 /** Get the mean longitude argument derivative. 326 * <p> 327 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 328 * </p> 329 * @return d(M + ω + Ω)/dt mean longitude argument derivative (rad/s) 330 * @see #hasDerivatives() 331 * @since 9.0 332 */ 333 public abstract double getLMDot(); 334 335 // Additional orbital elements 336 337 /** Get the eccentricity. 338 * @return eccentricity 339 */ 340 public abstract double getE(); 341 342 /** Get the eccentricity derivative. 343 * <p> 344 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 345 * </p> 346 * @return eccentricity derivative 347 * @see #hasDerivatives() 348 * @since 9.0 349 */ 350 public abstract double getEDot(); 351 352 /** Get the inclination. 353 * @return inclination (rad) 354 */ 355 public abstract double getI(); 356 357 /** Get the inclination derivative. 358 * <p> 359 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}. 360 * </p> 361 * @return inclination derivative (rad/s) 362 * @see #hasDerivatives() 363 * @since 9.0 364 */ 365 public abstract double getIDot(); 366 367 /** Check if orbit includes derivatives. 368 * @return true if orbit includes derivatives 369 * @see #getADot() 370 * @see #getEquinoctialExDot() 371 * @see #getEquinoctialEyDot() 372 * @see #getHxDot() 373 * @see #getHyDot() 374 * @see #getLEDot() 375 * @see #getLvDot() 376 * @see #getLMDot() 377 * @see #getEDot() 378 * @see #getIDot() 379 * @since 9.0 380 */ 381 public boolean hasDerivatives() { 382 return !Double.isNaN(getADot()); 383 } 384 385 /** Get the central acceleration constant. 386 * @return central acceleration constant 387 */ 388 public double getMu() { 389 return mu; 390 } 391 392 /** Get the Keplerian period. 393 * <p>The Keplerian period is computed directly from semi major axis 394 * and central acceleration constant.</p> 395 * @return Keplerian period in seconds, or positive infinity for hyperbolic orbits 396 */ 397 public double getKeplerianPeriod() { 398 final double a = getA(); 399 return (a < 0) ? Double.POSITIVE_INFINITY : 2.0 * FastMath.PI * a * FastMath.sqrt(a / mu); 400 } 401 402 /** Get the Keplerian mean motion. 403 * <p>The Keplerian mean motion is computed directly from semi major axis 404 * and central acceleration constant.</p> 405 * @return Keplerian mean motion in radians per second 406 */ 407 public double getKeplerianMeanMotion() { 408 final double absA = FastMath.abs(getA()); 409 return FastMath.sqrt(mu / absA) / absA; 410 } 411 412 /** Get the date of orbital parameters. 413 * @return date of the orbital parameters 414 */ 415 public AbsoluteDate getDate() { 416 return date; 417 } 418 419 /** Get the {@link TimeStampedPVCoordinates} in a specified frame. 420 * @param outputFrame frame in which the position/velocity coordinates shall be computed 421 * @return pvCoordinates in the specified output frame 422 * @see #getPVCoordinates() 423 */ 424 public TimeStampedPVCoordinates getPVCoordinates(final Frame outputFrame) { 425 if (pvCoordinates == null) { 426 pvCoordinates = initPVCoordinates(); 427 } 428 429 // If output frame requested is the same as definition frame, 430 // PV coordinates are returned directly 431 if (outputFrame == frame) { 432 return pvCoordinates; 433 } 434 435 // Else, PV coordinates are transformed to output frame 436 final Transform t = frame.getTransformTo(outputFrame, date); 437 return t.transformPVCoordinates(pvCoordinates); 438 } 439 440 /** {@inheritDoc} */ 441 public TimeStampedPVCoordinates getPVCoordinates(final AbsoluteDate otherDate, final Frame otherFrame) { 442 return shiftedBy(otherDate.durationFrom(getDate())).getPVCoordinates(otherFrame); 443 } 444 445 446 /** Get the {@link TimeStampedPVCoordinates} in definition frame. 447 * @return pvCoordinates in the definition frame 448 * @see #getPVCoordinates(Frame) 449 */ 450 public TimeStampedPVCoordinates getPVCoordinates() { 451 if (pvCoordinates == null) { 452 pvCoordinates = initPVCoordinates(); 453 } 454 return pvCoordinates; 455 } 456 457 /** Compute the position/velocity coordinates from the canonical parameters. 458 * @return computed position/velocity coordinates 459 */ 460 protected abstract TimeStampedPVCoordinates initPVCoordinates(); 461 462 /** Get a time-shifted orbit. 463 * <p> 464 * The orbit can be slightly shifted to close dates. The shifting model is a 465 * Keplerian one if no derivatives are available in the orbit, or Keplerian 466 * plus quadratic effect of the non-Keplerian acceleration if derivatives are 467 * available. Shifting is <em>not</em> intended as a replacement for proper 468 * orbit propagation but should be sufficient for small time shifts or coarse 469 * accuracy. 470 * </p> 471 * @param dt time shift in seconds 472 * @return a new orbit, shifted with respect to the instance (which is immutable) 473 */ 474 public abstract Orbit shiftedBy(double dt); 475 476 /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters. 477 * <p> 478 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with 479 * respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter 480 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot. 481 * </p> 482 * @param type type of the position angle to use 483 * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix 484 * is larger than 6x6, only the 6x6 upper left corner will be modified 485 */ 486 public void getJacobianWrtCartesian(final PositionAngle type, final double[][] jacobian) { 487 488 final double[][] cachedJacobian; 489 synchronized (this) { 490 switch (type) { 491 case MEAN : 492 if (jacobianMeanWrtCartesian == null) { 493 // first call, we need to compute the Jacobian and cache it 494 jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian(); 495 } 496 cachedJacobian = jacobianMeanWrtCartesian; 497 break; 498 case ECCENTRIC : 499 if (jacobianEccentricWrtCartesian == null) { 500 // first call, we need to compute the Jacobian and cache it 501 jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian(); 502 } 503 cachedJacobian = jacobianEccentricWrtCartesian; 504 break; 505 case TRUE : 506 if (jacobianTrueWrtCartesian == null) { 507 // first call, we need to compute the Jacobian and cache it 508 jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian(); 509 } 510 cachedJacobian = jacobianTrueWrtCartesian; 511 break; 512 default : 513 throw new OrekitInternalError(null); 514 } 515 } 516 517 // fill the user provided array 518 for (int i = 0; i < cachedJacobian.length; ++i) { 519 System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length); 520 } 521 522 } 523 524 /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters. 525 * <p> 526 * Element {@code jacobian[i][j]} is the derivative of Cartesian coordinate i of the orbit with 527 * respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate 528 * x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters. 529 * </p> 530 * @param type type of the position angle to use 531 * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix 532 * is larger than 6x6, only the 6x6 upper left corner will be modified 533 */ 534 public void getJacobianWrtParameters(final PositionAngle type, final double[][] jacobian) { 535 536 final double[][] cachedJacobian; 537 synchronized (this) { 538 switch (type) { 539 case MEAN : 540 if (jacobianWrtParametersMean == null) { 541 // first call, we need to compute the Jacobian and cache it 542 jacobianWrtParametersMean = createInverseJacobian(type); 543 } 544 cachedJacobian = jacobianWrtParametersMean; 545 break; 546 case ECCENTRIC : 547 if (jacobianWrtParametersEccentric == null) { 548 // first call, we need to compute the Jacobian and cache it 549 jacobianWrtParametersEccentric = createInverseJacobian(type); 550 } 551 cachedJacobian = jacobianWrtParametersEccentric; 552 break; 553 case TRUE : 554 if (jacobianWrtParametersTrue == null) { 555 // first call, we need to compute the Jacobian and cache it 556 jacobianWrtParametersTrue = createInverseJacobian(type); 557 } 558 cachedJacobian = jacobianWrtParametersTrue; 559 break; 560 default : 561 throw new OrekitInternalError(null); 562 } 563 } 564 565 // fill the user-provided array 566 for (int i = 0; i < cachedJacobian.length; ++i) { 567 System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length); 568 } 569 570 } 571 572 /** Create an inverse Jacobian. 573 * @param type type of the position angle to use 574 * @return inverse Jacobian 575 */ 576 private double[][] createInverseJacobian(final PositionAngle type) { 577 578 // get the direct Jacobian 579 final double[][] directJacobian = new double[6][6]; 580 getJacobianWrtCartesian(type, directJacobian); 581 582 // invert the direct Jacobian 583 final RealMatrix matrix = MatrixUtils.createRealMatrix(directJacobian); 584 final DecompositionSolver solver = new QRDecomposition(matrix).getSolver(); 585 return solver.getInverse().getData(); 586 587 } 588 589 /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters. 590 * <p> 591 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with 592 * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter 593 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot. 594 * </p> 595 * @return 6x6 Jacobian matrix 596 * @see #computeJacobianEccentricWrtCartesian() 597 * @see #computeJacobianTrueWrtCartesian() 598 */ 599 protected abstract double[][] computeJacobianMeanWrtCartesian(); 600 601 /** Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters. 602 * <p> 603 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with 604 * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter 605 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot. 606 * </p> 607 * @return 6x6 Jacobian matrix 608 * @see #computeJacobianMeanWrtCartesian() 609 * @see #computeJacobianTrueWrtCartesian() 610 */ 611 protected abstract double[][] computeJacobianEccentricWrtCartesian(); 612 613 /** Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters. 614 * <p> 615 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with 616 * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter 617 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot. 618 * </p> 619 * @return 6x6 Jacobian matrix 620 * @see #computeJacobianMeanWrtCartesian() 621 * @see #computeJacobianEccentricWrtCartesian() 622 */ 623 protected abstract double[][] computeJacobianTrueWrtCartesian(); 624 625 /** Add the contribution of the Keplerian motion to parameters derivatives 626 * <p> 627 * This method is used by integration-based propagators to evaluate the part of Keplerian 628 * motion to evolution of the orbital state. 629 * </p> 630 * @param type type of the position angle in the state 631 * @param gm attraction coefficient to use 632 * @param pDot array containing orbital state derivatives to update (the Keplerian 633 * part must be <em>added</em> to the array components, as the array may already 634 * contain some non-zero elements corresponding to non-Keplerian parts) 635 */ 636 public abstract void addKeplerContribution(PositionAngle type, double gm, double[] pDot); 637 638 /** Fill a Jacobian half row with a single vector. 639 * @param a coefficient of the vector 640 * @param v vector 641 * @param row Jacobian matrix row 642 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 643 */ 644 protected static void fillHalfRow(final double a, final Vector3D v, final double[] row, final int j) { 645 row[j] = a * v.getX(); 646 row[j + 1] = a * v.getY(); 647 row[j + 2] = a * v.getZ(); 648 } 649 650 /** Fill a Jacobian half row with a linear combination of vectors. 651 * @param a1 coefficient of the first vector 652 * @param v1 first vector 653 * @param a2 coefficient of the second vector 654 * @param v2 second vector 655 * @param row Jacobian matrix row 656 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 657 */ 658 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2, 659 final double[] row, final int j) { 660 row[j] = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX()); 661 row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY()); 662 row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ()); 663 } 664 665 /** Fill a Jacobian half row with a linear combination of vectors. 666 * @param a1 coefficient of the first vector 667 * @param v1 first vector 668 * @param a2 coefficient of the second vector 669 * @param v2 second vector 670 * @param a3 coefficient of the third vector 671 * @param v3 third vector 672 * @param row Jacobian matrix row 673 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 674 */ 675 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2, 676 final double a3, final Vector3D v3, 677 final double[] row, final int j) { 678 row[j] = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX()); 679 row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY()); 680 row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ()); 681 } 682 683 /** Fill a Jacobian half row with a linear combination of vectors. 684 * @param a1 coefficient of the first vector 685 * @param v1 first vector 686 * @param a2 coefficient of the second vector 687 * @param v2 second vector 688 * @param a3 coefficient of the third vector 689 * @param v3 third vector 690 * @param a4 coefficient of the fourth vector 691 * @param v4 fourth vector 692 * @param row Jacobian matrix row 693 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 694 */ 695 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2, 696 final double a3, final Vector3D v3, final double a4, final Vector3D v4, 697 final double[] row, final int j) { 698 row[j] = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX(), a4, v4.getX()); 699 row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY(), a4, v4.getY()); 700 row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ(), a4, v4.getZ()); 701 } 702 703 /** Fill a Jacobian half row with a linear combination of vectors. 704 * @param a1 coefficient of the first vector 705 * @param v1 first vector 706 * @param a2 coefficient of the second vector 707 * @param v2 second vector 708 * @param a3 coefficient of the third vector 709 * @param v3 third vector 710 * @param a4 coefficient of the fourth vector 711 * @param v4 fourth vector 712 * @param a5 coefficient of the fifth vector 713 * @param v5 fifth vector 714 * @param row Jacobian matrix row 715 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 716 */ 717 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2, 718 final double a3, final Vector3D v3, final double a4, final Vector3D v4, 719 final double a5, final Vector3D v5, 720 final double[] row, final int j) { 721 final double[] a = new double[] { 722 a1, a2, a3, a4, a5 723 }; 724 row[j] = MathArrays.linearCombination(a, new double[] { 725 v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX() 726 }); 727 row[j + 1] = MathArrays.linearCombination(a, new double[] { 728 v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY() 729 }); 730 row[j + 2] = MathArrays.linearCombination(a, new double[] { 731 v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ() 732 }); 733 } 734 735 /** Fill a Jacobian half row with a linear combination of vectors. 736 * @param a1 coefficient of the first vector 737 * @param v1 first vector 738 * @param a2 coefficient of the second vector 739 * @param v2 second vector 740 * @param a3 coefficient of the third vector 741 * @param v3 third vector 742 * @param a4 coefficient of the fourth vector 743 * @param v4 fourth vector 744 * @param a5 coefficient of the fifth vector 745 * @param v5 fifth vector 746 * @param a6 coefficient of the sixth vector 747 * @param v6 sixth vector 748 * @param row Jacobian matrix row 749 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set) 750 */ 751 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2, 752 final double a3, final Vector3D v3, final double a4, final Vector3D v4, 753 final double a5, final Vector3D v5, final double a6, final Vector3D v6, 754 final double[] row, final int j) { 755 final double[] a = new double[] { 756 a1, a2, a3, a4, a5, a6 757 }; 758 row[j] = MathArrays.linearCombination(a, new double[] { 759 v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX(), v6.getX() 760 }); 761 row[j + 1] = MathArrays.linearCombination(a, new double[] { 762 v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY(), v6.getY() 763 }); 764 row[j + 2] = MathArrays.linearCombination(a, new double[] { 765 v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ(), v6.getZ() 766 }); 767 } 768 769 }