AdamsBashforthFieldIntegrator.java

  1. /*
  2.  * Licensed to the Apache Software Foundation (ASF) 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 ASF 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. /*
  18.  * This is not the original file distributed by the Apache Software Foundation
  19.  * It has been modified by the Hipparchus project
  20.  */

  21. package org.hipparchus.ode.nonstiff;

  22. import org.hipparchus.Field;
  23. import org.hipparchus.CalculusFieldElement;
  24. import org.hipparchus.exception.MathIllegalArgumentException;
  25. import org.hipparchus.linear.Array2DRowFieldMatrix;
  26. import org.hipparchus.linear.FieldMatrix;
  27. import org.hipparchus.ode.FieldEquationsMapper;
  28. import org.hipparchus.ode.FieldODEStateAndDerivative;
  29. import org.hipparchus.ode.nonstiff.interpolators.AdamsFieldStateInterpolator;
  30. import org.hipparchus.util.FastMath;


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

  157.     /** Name of integration scheme. */
  158.     public static final String METHOD_NAME = AdamsBashforthIntegrator.METHOD_NAME;

  159.     /**
  160.      * Build an Adams-Bashforth integrator with the given order and step control parameters.
  161.      * @param field field to which the time and state vector elements belong
  162.      * @param nSteps number of steps of the method excluding the one being computed
  163.      * @param minStep minimal step (sign is irrelevant, regardless of
  164.      * integration direction, forward or backward), the last step can
  165.      * be smaller than this
  166.      * @param maxStep maximal step (sign is irrelevant, regardless of
  167.      * integration direction, forward or backward), the last step can
  168.      * be smaller than this
  169.      * @param scalAbsoluteTolerance allowed absolute error
  170.      * @param scalRelativeTolerance allowed relative error
  171.      * @exception MathIllegalArgumentException if order is 1 or less
  172.      */
  173.     public AdamsBashforthFieldIntegrator(final Field<T> field, final int nSteps,
  174.                                          final double minStep, final double maxStep,
  175.                                          final double scalAbsoluteTolerance,
  176.                                          final double scalRelativeTolerance)
  177.         throws MathIllegalArgumentException {
  178.         super(field, METHOD_NAME, nSteps, nSteps, minStep, maxStep,
  179.               scalAbsoluteTolerance, scalRelativeTolerance);
  180.     }

  181.     /**
  182.      * Build an Adams-Bashforth integrator with the given order and step control parameters.
  183.      * @param field field to which the time and state vector elements belong
  184.      * @param nSteps number of steps of the method excluding the one being computed
  185.      * @param minStep minimal step (sign is irrelevant, regardless of
  186.      * integration direction, forward or backward), the last step can
  187.      * be smaller than this
  188.      * @param maxStep maximal step (sign is irrelevant, regardless of
  189.      * integration direction, forward or backward), the last step can
  190.      * be smaller than this
  191.      * @param vecAbsoluteTolerance allowed absolute error
  192.      * @param vecRelativeTolerance allowed relative error
  193.      * @exception IllegalArgumentException if order is 1 or less
  194.      */
  195.     public AdamsBashforthFieldIntegrator(final Field<T> field, final int nSteps,
  196.                                          final double minStep, final double maxStep,
  197.                                          final double[] vecAbsoluteTolerance,
  198.                                          final double[] vecRelativeTolerance)
  199.         throws IllegalArgumentException {
  200.         super(field, METHOD_NAME, nSteps, nSteps, minStep, maxStep,
  201.               vecAbsoluteTolerance, vecRelativeTolerance);
  202.     }

  203.     /** {@inheritDoc} */
  204.     @Override
  205.     protected double errorEstimation(final T[] previousState, final T predictedTime,
  206.                                      final T[] predictedState, final T[] predictedScaled,
  207.                                      final FieldMatrix<T> predictedNordsieck) {

  208.         final StepsizeHelper helper = getStepSizeHelper();
  209.         double error = 0;
  210.         for (int i = 0; i < helper.getMainSetDimension(); ++i) {
  211.             final double tol = helper.getTolerance(i, FastMath.abs(predictedState[i].getReal()));

  212.             // apply Taylor formula from high order to low order,
  213.             // for the sake of numerical accuracy
  214.             double variation = 0;
  215.             int sign = predictedNordsieck.getRowDimension() % 2 == 0 ? -1 : 1;
  216.             for (int k = predictedNordsieck.getRowDimension() - 1; k >= 0; --k) {
  217.                 variation += sign * predictedNordsieck.getEntry(k, i).getReal();
  218.                 sign       = -sign;
  219.             }
  220.             variation -= predictedScaled[i].getReal();

  221.             final double ratio  = (predictedState[i].getReal() - previousState[i].getReal() + variation) / tol;
  222.             error              += ratio * ratio;

  223.         }

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

  225.     }

  226.     /** {@inheritDoc} */
  227.     @Override
  228.     protected AdamsFieldStateInterpolator<T> finalizeStep(final T stepSize, final T[] predictedY,
  229.                                                           final T[] predictedScaled, final Array2DRowFieldMatrix<T> predictedNordsieck,
  230.                                                           final boolean isForward,
  231.                                                           final FieldODEStateAndDerivative<T> globalPreviousState,
  232.                                                           final FieldODEStateAndDerivative<T> globalCurrentState,
  233.                                                           final FieldEquationsMapper<T> equationsMapper) {
  234.         return new AdamsFieldStateInterpolator<>(getStepSize(), globalCurrentState,
  235.                                                  predictedScaled, predictedNordsieck, isForward,
  236.                                                  getStepStart(), globalCurrentState, equationsMapper);
  237.     }

  238. }