AdamsMoultonIntegrator.java

  1. /*
  2.  * Licensed to the Hipparchus project under one or more
  3.  * contributor license agreements.  See the NOTICE file distributed with
  4.  * this work for additional information regarding copyright ownership.
  5.  * The Hipparchus project 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.  *      https://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.hipparchus.ode.nonstiff;

  18. import java.util.Arrays;

  19. import org.hipparchus.exception.MathIllegalArgumentException;
  20. import org.hipparchus.exception.MathIllegalStateException;
  21. import org.hipparchus.linear.Array2DRowRealMatrix;
  22. import org.hipparchus.linear.RealMatrix;
  23. import org.hipparchus.linear.RealMatrixPreservingVisitor;
  24. import org.hipparchus.ode.EquationsMapper;
  25. import org.hipparchus.ode.LocalizedODEFormats;
  26. import org.hipparchus.ode.ODEStateAndDerivative;
  27. import org.hipparchus.ode.nonstiff.interpolators.AdamsStateInterpolator;
  28. import org.hipparchus.util.FastMath;


  29. /**
  30.  * This class implements implicit Adams-Moulton integrators for Ordinary
  31.  * Differential Equations.
  32.  *
  33.  * <p>Adams-Moulton methods (in fact due to Adams alone) are implicit
  34.  * multistep ODE solvers. This implementation is a variation of the classical
  35.  * one: it uses adaptive stepsize to implement error control, whereas
  36.  * classical implementations are fixed step size. The value of state vector
  37.  * at step n+1 is a simple combination of the value at step n and of the
  38.  * derivatives at steps n+1, n, n-1 ... Since y'<sub>n+1</sub> is needed to
  39.  * compute y<sub>n+1</sub>, another method must be used to compute a first
  40.  * estimate of y<sub>n+1</sub>, then compute y'<sub>n+1</sub>, then compute
  41.  * a final estimate of y<sub>n+1</sub> using the following formulas. Depending
  42.  * on the number k of previous steps one wants to use for computing the next
  43.  * value, different formulas are available for the final estimate:</p>
  44.  * <ul>
  45.  *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + h y'<sub>n+1</sub></li>
  46.  *   <li>k = 2: y<sub>n+1</sub> = y<sub>n</sub> + h (y'<sub>n+1</sub>+y'<sub>n</sub>)/2</li>
  47.  *   <li>k = 3: y<sub>n+1</sub> = y<sub>n</sub> + h (5y'<sub>n+1</sub>+8y'<sub>n</sub>-y'<sub>n-1</sub>)/12</li>
  48.  *   <li>k = 4: y<sub>n+1</sub> = y<sub>n</sub> + h (9y'<sub>n+1</sub>+19y'<sub>n</sub>-5y'<sub>n-1</sub>+y'<sub>n-2</sub>)/24</li>
  49.  *   <li>...</li>
  50.  * </ul>
  51.  *
  52.  * <p>A k-steps Adams-Moulton method is of order k+1.</p>
  53.  *
  54.  * <p> There must be sufficient time for the {@link #setStarterIntegrator(org.hipparchus.ode.ODEIntegrator)
  55.  * starter integrator} to take several steps between the the last reset event, and the end
  56.  * of integration, otherwise an exception may be thrown during integration. The user can
  57.  * adjust the end date of integration, or the step size of the starter integrator to
  58.  * ensure a sufficient number of steps can be completed before the end of integration.
  59.  * </p>
  60.  *
  61.  * <p><strong>Implementation details</strong></p>
  62.  *
  63.  * <p>We define scaled derivatives s<sub>i</sub>(n) at step n as:
  64.  * \[
  65.  *   \left\{\begin{align}
  66.  *   s_1(n) &amp;= h y'_n \text{ for first derivative}\\
  67.  *   s_2(n) &amp;= \frac{h^2}{2} y_n'' \text{ for second derivative}\\
  68.  *   s_3(n) &amp;= \frac{h^3}{6} y_n''' \text{ for third derivative}\\
  69.  *   &amp;\cdots\\
  70.  *   s_k(n) &amp;= \frac{h^k}{k!} y_n^{(k)} \text{ for } k^\mathrm{th} \text{ derivative}
  71.  *   \end{align}\right.
  72.  * \]</p>
  73.  *
  74.  * <p>The definitions above use the classical representation with several previous first
  75.  * derivatives. Lets define
  76.  * \[
  77.  *   q_n = [ s_1(n-1) s_1(n-2) \ldots s_1(n-(k-1)) ]^T
  78.  * \]
  79.  * (we omit the k index in the notation for clarity). With these definitions,
  80.  * Adams-Moulton methods can be written:</p>
  81.  * <ul>
  82.  *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n+1)</li>
  83.  *   <li>k = 2: y<sub>n+1</sub> = y<sub>n</sub> + 1/2 s<sub>1</sub>(n+1) + [ 1/2 ] q<sub>n+1</sub></li>
  84.  *   <li>k = 3: y<sub>n+1</sub> = y<sub>n</sub> + 5/12 s<sub>1</sub>(n+1) + [ 8/12 -1/12 ] q<sub>n+1</sub></li>
  85.  *   <li>k = 4: y<sub>n+1</sub> = y<sub>n</sub> + 9/24 s<sub>1</sub>(n+1) + [ 19/24 -5/24 1/24 ] q<sub>n+1</sub></li>
  86.  *   <li>...</li>
  87.  * </ul>
  88.  *
  89.  * <p>Instead of using the classical representation with first derivatives only (y<sub>n</sub>,
  90.  * s<sub>1</sub>(n+1) and q<sub>n+1</sub>), our implementation uses the Nordsieck vector with
  91.  * higher degrees scaled derivatives all taken at the same step (y<sub>n</sub>, s<sub>1</sub>(n)
  92.  * and r<sub>n</sub>) where r<sub>n</sub> is defined as:
  93.  * \[
  94.  * r_n = [ s_2(n), s_3(n) \ldots s_k(n) ]^T
  95.  * \]
  96.  * (here again we omit the k index in the notation for clarity)
  97.  * </p>
  98.  *
  99.  * <p>Taylor series formulas show that for any index offset i, s<sub>1</sub>(n-i) can be
  100.  * computed from s<sub>1</sub>(n), s<sub>2</sub>(n) ... s<sub>k</sub>(n), the formula being exact
  101.  * for degree k polynomials.
  102.  * \[
  103.  * s_1(n-i) = s_1(n) + \sum_{j\gt 0} (j+1) (-i)^j s_{j+1}(n)
  104.  * \]
  105.  * The previous formula can be used with several values for i to compute the transform between
  106.  * classical representation and Nordsieck vector. The transform between r<sub>n</sub>
  107.  * and q<sub>n</sub> resulting from the Taylor series formulas above is:
  108.  * \[
  109.  * q_n = s_1(n) u + P r_n
  110.  * \]
  111.  * where u is the [ 1 1 ... 1 ]<sup>T</sup> vector and P is the (k-1)&times;(k-1) matrix built
  112.  * with the \((j+1) (-i)^j\) terms with i being the row number starting from 1 and j being
  113.  * the column number starting from 1:
  114.  * \[
  115.  *   P=\begin{bmatrix}
  116.  *   -2  &amp;  3 &amp;   -4 &amp;    5 &amp; \ldots \\
  117.  *   -4  &amp; 12 &amp;  -32 &amp;   80 &amp; \ldots \\
  118.  *   -6  &amp; 27 &amp; -108 &amp;  405 &amp; \ldots \\
  119.  *   -8  &amp; 48 &amp; -256 &amp; 1280 &amp; \ldots \\
  120.  *       &amp;    &amp;  \ldots\\
  121.  *    \end{bmatrix}
  122.  * \]
  123.  *
  124.  * <p>Using the Nordsieck vector has several advantages:</p>
  125.  * <ul>
  126.  *   <li>it greatly simplifies step interpolation as the interpolator mainly applies
  127.  *   Taylor series formulas,</li>
  128.  *   <li>it simplifies step changes that occur when discrete events that truncate
  129.  *   the step are triggered,</li>
  130.  *   <li>it allows to extend the methods in order to support adaptive stepsize.</li>
  131.  * </ul>
  132.  *
  133.  * <p>The predicted Nordsieck vector at step n+1 is computed from the Nordsieck vector at step
  134.  * n as follows:
  135.  * <ul>
  136.  *   <li>Y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li>
  137.  *   <li>S<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, Y<sub>n+1</sub>)</li>
  138.  *   <li>R<sub>n+1</sub> = (s<sub>1</sub>(n) - S<sub>1</sub>(n+1)) P<sup>-1</sup> u + P<sup>-1</sup> A P r<sub>n</sub></li>
  139.  * </ul>
  140.  * where A is a rows shifting matrix (the lower left part is an identity matrix):
  141.  * <pre>
  142.  *        [ 0 0   ...  0 0 | 0 ]
  143.  *        [ ---------------+---]
  144.  *        [ 1 0   ...  0 0 | 0 ]
  145.  *    A = [ 0 1   ...  0 0 | 0 ]
  146.  *        [       ...      | 0 ]
  147.  *        [ 0 0   ...  1 0 | 0 ]
  148.  *        [ 0 0   ...  0 1 | 0 ]
  149.  * </pre>
  150.  * From this predicted vector, the corrected vector is computed as follows:
  151.  * <ul>
  152.  *   <li>y<sub>n+1</sub> = y<sub>n</sub> + S<sub>1</sub>(n+1) + [ -1 +1 -1 +1 ... &plusmn;1 ] r<sub>n+1</sub></li>
  153.  *   <li>s<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, y<sub>n+1</sub>)</li>
  154.  *   <li>r<sub>n+1</sub> = R<sub>n+1</sub> + (s<sub>1</sub>(n+1) - S<sub>1</sub>(n+1)) P<sup>-1</sup> u</li>
  155.  * </ul>
  156.  * <p>where the upper case Y<sub>n+1</sub>, S<sub>1</sub>(n+1) and R<sub>n+1</sub> represent the
  157.  * predicted states whereas the lower case y<sub>n+1</sub>, s<sub>n+1</sub> and r<sub>n+1</sub>
  158.  * represent the corrected states.</p>
  159.  *
  160.  * <p>The P<sup>-1</sup>u vector and the P<sup>-1</sup> A P matrix do not depend on the state,
  161.  * they only depend on k and therefore are precomputed once for all.</p>
  162.  *
  163.  */
  164. public class AdamsMoultonIntegrator extends AdamsIntegrator {

  165.     /** Name of integration scheme. */
  166.     public static final String METHOD_NAME = "Adams-Moulton";

  167.     /**
  168.      * Build an Adams-Moulton integrator with the given order and error control parameters.
  169.      * @param nSteps number of steps of the method excluding the one being computed
  170.      * @param minStep minimal step (sign is irrelevant, regardless of
  171.      * integration direction, forward or backward), the last step can
  172.      * be smaller than this
  173.      * @param maxStep maximal step (sign is irrelevant, regardless of
  174.      * integration direction, forward or backward), the last step can
  175.      * be smaller than this
  176.      * @param scalAbsoluteTolerance allowed absolute error
  177.      * @param scalRelativeTolerance allowed relative error
  178.      * @exception MathIllegalArgumentException if order is 1 or less
  179.      */
  180.     public AdamsMoultonIntegrator(final int nSteps,
  181.                                   final double minStep, final double maxStep,
  182.                                   final double scalAbsoluteTolerance,
  183.                                   final double scalRelativeTolerance)
  184.         throws MathIllegalArgumentException {
  185.         super(METHOD_NAME, nSteps, nSteps + 1, minStep, maxStep,
  186.               scalAbsoluteTolerance, scalRelativeTolerance);
  187.     }

  188.     /**
  189.      * Build an Adams-Moulton integrator with the given order and error control parameters.
  190.      * @param nSteps number of steps of the method excluding the one being computed
  191.      * @param minStep minimal step (sign is irrelevant, regardless of
  192.      * integration direction, forward or backward), the last step can
  193.      * be smaller than this
  194.      * @param maxStep maximal step (sign is irrelevant, regardless of
  195.      * integration direction, forward or backward), the last step can
  196.      * be smaller than this
  197.      * @param vecAbsoluteTolerance allowed absolute error
  198.      * @param vecRelativeTolerance allowed relative error
  199.      * @exception IllegalArgumentException if order is 1 or less
  200.      */
  201.     public AdamsMoultonIntegrator(final int nSteps,
  202.                                   final double minStep, final double maxStep,
  203.                                   final double[] vecAbsoluteTolerance,
  204.                                   final double[] vecRelativeTolerance)
  205.         throws IllegalArgumentException {
  206.         super(METHOD_NAME, nSteps, nSteps + 1, minStep, maxStep,
  207.               vecAbsoluteTolerance, vecRelativeTolerance);
  208.     }

  209.     /** {@inheritDoc} */
  210.     @Override
  211.     protected double errorEstimation(final double[] previousState, final double predictedTime,
  212.                                      final double[] predictedState,
  213.                                      final double[] predictedScaled,
  214.                                      final RealMatrix predictedNordsieck) {
  215.         final double error = predictedNordsieck.walkInOptimizedOrder(new Corrector(previousState, predictedScaled, predictedState));
  216.         if (Double.isNaN(error)) {
  217.             throw new MathIllegalStateException(LocalizedODEFormats.NAN_APPEARING_DURING_INTEGRATION,
  218.                                                 predictedTime);
  219.         }
  220.         return error;
  221.     }

  222.     /** {@inheritDoc} */
  223.     @Override
  224.     protected AdamsStateInterpolator finalizeStep(final double stepSize, final double[] predictedState,
  225.                                                   final double[] predictedScaled, final Array2DRowRealMatrix predictedNordsieck,
  226.                                                   final boolean isForward,
  227.                                                   final ODEStateAndDerivative globalPreviousState,
  228.                                                   final ODEStateAndDerivative globalCurrentState,
  229.                                                   final EquationsMapper equationsMapper) {

  230.         final double[] correctedYDot = computeDerivatives(globalCurrentState.getTime(), predictedState);

  231.         // update Nordsieck vector
  232.         final double[] correctedScaled = new double[predictedState.length];
  233.         for (int j = 0; j < correctedScaled.length; ++j) {
  234.             correctedScaled[j] = getStepSize() * correctedYDot[j];
  235.         }
  236.         updateHighOrderDerivativesPhase2(predictedScaled, correctedScaled, predictedNordsieck);

  237.         final ODEStateAndDerivative updatedStepEnd =
  238.                         equationsMapper.mapStateAndDerivative(globalCurrentState.getTime(),
  239.                                                               predictedState, correctedYDot);
  240.         return new AdamsStateInterpolator(getStepSize(), updatedStepEnd,
  241.                                           correctedScaled, predictedNordsieck, isForward,
  242.                                           getStepStart(), updatedStepEnd,
  243.                                           equationsMapper);

  244.     }

  245.     /** Corrector for current state in Adams-Moulton method.
  246.      * <p>
  247.      * This visitor implements the Taylor series formula:
  248.      * <pre>
  249.      * Y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n+1) + [ -1 +1 -1 +1 ... &plusmn;1 ] r<sub>n+1</sub>
  250.      * </pre>
  251.      * </p>
  252.      */
  253.     private class Corrector implements RealMatrixPreservingVisitor {

  254.         /** Previous state. */
  255.         private final double[] previous;

  256.         /** Current scaled first derivative. */
  257.         private final double[] scaled;

  258.         /** Current state before correction. */
  259.         private final double[] before;

  260.         /** Current state after correction. */
  261.         private final double[] after;

  262.         /** Simple constructor.
  263.          * <p>
  264.          * All arrays will be stored by reference to caller arrays.
  265.          * </p>
  266.          * @param previous previous state
  267.          * @param scaled current scaled first derivative
  268.          * @param state state to correct (will be overwritten after visit)
  269.          */
  270.         Corrector(final double[] previous, final double[] scaled, final double[] state) {
  271.             this.previous = previous; // NOPMD - array reference storage is intentional and documented here
  272.             this.scaled   = scaled;   // NOPMD - array reference storage is intentional and documented here
  273.             this.after    = state;    // NOPMD - array reference storage is intentional and documented here
  274.             this.before   = state.clone();
  275.         }

  276.         /** {@inheritDoc} */
  277.         @Override
  278.         public void start(int rows, int columns,
  279.                           int startRow, int endRow, int startColumn, int endColumn) {
  280.             Arrays.fill(after, 0.0);
  281.         }

  282.         /** {@inheritDoc} */
  283.         @Override
  284.         public void visit(int row, int column, double value) {
  285.             if ((row & 0x1) == 0) {
  286.                 after[column] -= value;
  287.             } else {
  288.                 after[column] += value;
  289.             }
  290.         }

  291.         /**
  292.          * End visiting the Nordsieck vector.
  293.          * <p>The correction is used to control stepsize. So its amplitude is
  294.          * considered to be an error, which must be normalized according to
  295.          * error control settings. If the normalized value is greater than 1,
  296.          * the correction was too large and the step must be rejected.</p>
  297.          * @return the normalized correction, if greater than 1, the step
  298.          * must be rejected
  299.          */
  300.         @Override
  301.         public double end() {

  302.             final StepsizeHelper helper = getStepSizeHelper();
  303.             double error = 0;
  304.             for (int i = 0; i < after.length; ++i) {
  305.                 after[i] += previous[i] + scaled[i];
  306.                 if (i < helper.getMainSetDimension()) {
  307.                     final double tol   = helper.getTolerance(i, FastMath.max(FastMath.abs(previous[i]), FastMath.abs(after[i])));
  308.                     final double ratio = (after[i] - before[i]) / tol; // (corrected-predicted)/tol
  309.                     error += ratio * ratio;
  310.                 }
  311.             }

  312.             return FastMath.sqrt(error / helper.getMainSetDimension());

  313.         }
  314.     }

  315. }