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.  *      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.apache.commons.math4.legacy.ode.nonstiff;

  18. import org.apache.commons.math4.legacy.core.Field;
  19. import org.apache.commons.math4.legacy.core.RealFieldElement;
  20. import org.apache.commons.math4.legacy.exception.DimensionMismatchException;
  21. import org.apache.commons.math4.legacy.exception.MaxCountExceededException;
  22. import org.apache.commons.math4.legacy.exception.NoBracketingException;
  23. import org.apache.commons.math4.legacy.exception.NumberIsTooSmallException;
  24. import org.apache.commons.math4.legacy.linear.Array2DRowFieldMatrix;
  25. import org.apache.commons.math4.legacy.linear.FieldMatrix;
  26. import org.apache.commons.math4.legacy.ode.FieldExpandableODE;
  27. import org.apache.commons.math4.legacy.ode.FieldODEState;
  28. import org.apache.commons.math4.legacy.ode.FieldODEStateAndDerivative;
  29. import org.apache.commons.math4.legacy.core.MathArrays;


  30. /**
  31.  * This class implements explicit Adams-Bashforth integrators for Ordinary
  32.  * Differential Equations.
  33.  *
  34.  * <p>Adams-Bashforth methods (in fact due to Adams alone) are explicit
  35.  * multistep ODE solvers. This implementation is a variation of the classical
  36.  * one: it uses adaptive stepsize to implement error control, whereas
  37.  * classical implementations are fixed step size. The value of state vector
  38.  * at step n+1 is a simple combination of the value at step n and of the
  39.  * derivatives at steps n, n-1, n-2 ... Depending on the number k of previous
  40.  * steps one wants to use for computing the next value, different formulas
  41.  * are available:</p>
  42.  * <ul>
  43.  *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + h y'<sub>n</sub></li>
  44.  *   <li>k = 2: y<sub>n+1</sub> = y<sub>n</sub> + h (3y'<sub>n</sub>-y'<sub>n-1</sub>)/2</li>
  45.  *   <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>
  46.  *   <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>
  47.  *   <li>...</li>
  48.  * </ul>
  49.  *
  50.  * <p>A k-steps Adams-Bashforth method is of order k.</p>
  51.  *
  52.  * <p><b>Implementation details</b></p>
  53.  *
  54.  * <p>We define scaled derivatives s<sub>i</sub>(n) at step n as:
  55.  * <div style="white-space: pre"><code>
  56.  * s<sub>1</sub>(n) = h y'<sub>n</sub> for first derivative
  57.  * s<sub>2</sub>(n) = h<sup>2</sup>/2 y''<sub>n</sub> for second derivative
  58.  * s<sub>3</sub>(n) = h<sup>3</sup>/6 y'''<sub>n</sub> for third derivative
  59.  * ...
  60.  * s<sub>k</sub>(n) = h<sup>k</sup>/k! y<sup>(k)</sup><sub>n</sub> for k<sup>th</sup> derivative
  61.  * </code></div>
  62.  *
  63.  * <p>The definitions above use the classical representation with several previous first
  64.  * derivatives. Lets define
  65.  * <div style="white-space: pre"><code>
  66.  *   q<sub>n</sub> = [ s<sub>1</sub>(n-1) s<sub>1</sub>(n-2) ... s<sub>1</sub>(n-(k-1)) ]<sup>T</sup>
  67.  * </code></div>
  68.  * (we omit the k index in the notation for clarity). With these definitions,
  69.  * Adams-Bashforth methods can be written:
  70.  * <ul>
  71.  *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n)</li>
  72.  *   <li>k = 2: y<sub>n+1</sub> = y<sub>n</sub> + 3/2 s<sub>1</sub>(n) + [ -1/2 ] q<sub>n</sub></li>
  73.  *   <li>k = 3: y<sub>n+1</sub> = y<sub>n</sub> + 23/12 s<sub>1</sub>(n) + [ -16/12 5/12 ] q<sub>n</sub></li>
  74.  *   <li>k = 4: y<sub>n+1</sub> = y<sub>n</sub> + 55/24 s<sub>1</sub>(n) + [ -59/24 37/24 -9/24 ] q<sub>n</sub></li>
  75.  *   <li>...</li>
  76.  * </ul>
  77.  *
  78.  * <p>Instead of using the classical representation with first derivatives only (y<sub>n</sub>,
  79.  * s<sub>1</sub>(n) and q<sub>n</sub>), our implementation uses the Nordsieck vector with
  80.  * higher degrees scaled derivatives all taken at the same step (y<sub>n</sub>, s<sub>1</sub>(n)
  81.  * and r<sub>n</sub>) where r<sub>n</sub> is defined as:
  82.  * <div style="white-space: pre"><code>
  83.  * r<sub>n</sub> = [ s<sub>2</sub>(n), s<sub>3</sub>(n) ... s<sub>k</sub>(n) ]<sup>T</sup>
  84.  * </code></div>
  85.  * (here again we omit the k index in the notation for clarity)
  86.  *
  87.  * <p>Taylor series formulas show that for any index offset i, s<sub>1</sub>(n-i) can be
  88.  * computed from s<sub>1</sub>(n), s<sub>2</sub>(n) ... s<sub>k</sub>(n), the formula being exact
  89.  * for degree k polynomials.
  90.  * <div style="white-space: pre"><code>
  91.  * s<sub>1</sub>(n-i) = s<sub>1</sub>(n) + &sum;<sub>j&gt;0</sub> (j+1) (-i)<sup>j</sup> s<sub>j+1</sub>(n)
  92.  * </code></div>
  93.  * The previous formula can be used with several values for i to compute the transform between
  94.  * classical representation and Nordsieck vector. The transform between r<sub>n</sub>
  95.  * and q<sub>n</sub> resulting from the Taylor series formulas above is:
  96.  * <div style="white-space: pre"><code>
  97.  * q<sub>n</sub> = s<sub>1</sub>(n) u + P r<sub>n</sub>
  98.  * </code></div>
  99.  * where u is the [ 1 1 ... 1 ]<sup>T</sup> vector and P is the (k-1)&times;(k-1) matrix built
  100.  * with the (j+1) (-i)<sup>j</sup> terms with i being the row number starting from 1 and j being
  101.  * the column number starting from 1:
  102.  * <pre>
  103.  *        [  -2   3   -4    5  ... ]
  104.  *        [  -4  12  -32   80  ... ]
  105.  *   P =  [  -6  27 -108  405  ... ]
  106.  *        [  -8  48 -256 1280  ... ]
  107.  *        [          ...           ]
  108.  * </pre>
  109.  *
  110.  * <p>Using the Nordsieck vector has several advantages:
  111.  * <ul>
  112.  *   <li>it greatly simplifies step interpolation as the interpolator mainly applies
  113.  *   Taylor series formulas,</li>
  114.  *   <li>it simplifies step changes that occur when discrete events that truncate
  115.  *   the step are triggered,</li>
  116.  *   <li>it allows to extend the methods in order to support adaptive stepsize.</li>
  117.  * </ul>
  118.  *
  119.  * <p>The Nordsieck vector at step n+1 is computed from the Nordsieck vector at step n as follows:
  120.  * <ul>
  121.  *   <li>y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li>
  122.  *   <li>s<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, y<sub>n+1</sub>)</li>
  123.  *   <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>
  124.  * </ul>
  125.  * where A is a rows shifting matrix (the lower left part is an identity matrix):
  126.  * <pre>
  127.  *        [ 0 0   ...  0 0 | 0 ]
  128.  *        [ ---------------+---]
  129.  *        [ 1 0   ...  0 0 | 0 ]
  130.  *    A = [ 0 1   ...  0 0 | 0 ]
  131.  *        [       ...      | 0 ]
  132.  *        [ 0 0   ...  1 0 | 0 ]
  133.  *        [ 0 0   ...  0 1 | 0 ]
  134.  * </pre>
  135.  *
  136.  * <p>The P<sup>-1</sup>u vector and the P<sup>-1</sup> A P matrix do not depend on the state,
  137.  * they only depend on k and therefore are precomputed once for all.</p>
  138.  *
  139.  * @param <T> the type of the field elements
  140.  * @since 3.6
  141.  */
  142. public class AdamsBashforthFieldIntegrator<T extends RealFieldElement<T>> extends AdamsFieldIntegrator<T> {

  143.     /** Integrator method name. */
  144.     private static final String METHOD_NAME = "Adams-Bashforth";

  145.     /**
  146.      * Build an Adams-Bashforth integrator with the given order and step control parameters.
  147.      * @param field field to which the time and state vector elements belong
  148.      * @param nSteps number of steps of the method excluding the one being computed
  149.      * @param minStep minimal step (sign is irrelevant, regardless of
  150.      * integration direction, forward or backward), the last step can
  151.      * be smaller than this
  152.      * @param maxStep maximal step (sign is irrelevant, regardless of
  153.      * integration direction, forward or backward), the last step can
  154.      * be smaller than this
  155.      * @param scalAbsoluteTolerance allowed absolute error
  156.      * @param scalRelativeTolerance allowed relative error
  157.      * @exception NumberIsTooSmallException if order is 1 or less
  158.      */
  159.     public AdamsBashforthFieldIntegrator(final Field<T> field, final int nSteps,
  160.                                          final double minStep, final double maxStep,
  161.                                          final double scalAbsoluteTolerance,
  162.                                          final double scalRelativeTolerance)
  163.         throws NumberIsTooSmallException {
  164.         super(field, METHOD_NAME, nSteps, nSteps, minStep, maxStep,
  165.               scalAbsoluteTolerance, scalRelativeTolerance);
  166.     }

  167.     /**
  168.      * Build an Adams-Bashforth integrator with the given order and step control parameters.
  169.      * @param field field to which the time and state vector elements belong
  170.      * @param nSteps number of steps of the method excluding the one being computed
  171.      * @param minStep minimal step (sign is irrelevant, regardless of
  172.      * integration direction, forward or backward), the last step can
  173.      * be smaller than this
  174.      * @param maxStep maximal step (sign is irrelevant, regardless of
  175.      * integration direction, forward or backward), the last step can
  176.      * be smaller than this
  177.      * @param vecAbsoluteTolerance allowed absolute error
  178.      * @param vecRelativeTolerance allowed relative error
  179.      * @exception IllegalArgumentException if order is 1 or less
  180.      */
  181.     public AdamsBashforthFieldIntegrator(final Field<T> field, final int nSteps,
  182.                                          final double minStep, final double maxStep,
  183.                                          final double[] vecAbsoluteTolerance,
  184.                                          final double[] vecRelativeTolerance)
  185.         throws IllegalArgumentException {
  186.         super(field, METHOD_NAME, nSteps, nSteps, minStep, maxStep,
  187.               vecAbsoluteTolerance, vecRelativeTolerance);
  188.     }

  189.     /** Estimate error.
  190.      * <p>
  191.      * Error is estimated by interpolating back to previous state using
  192.      * the state Taylor expansion and comparing to real previous state.
  193.      * </p>
  194.      * @param previousState state vector at step start
  195.      * @param predictedState predicted state vector at step end
  196.      * @param predictedScaled predicted value of the scaled derivatives at step end
  197.      * @param predictedNordsieck predicted value of the Nordsieck vector at step end
  198.      * @return estimated normalized local discretization error
  199.      */
  200.     private T errorEstimation(final T[] previousState,
  201.                               final T[] predictedState,
  202.                               final T[] predictedScaled,
  203.                               final FieldMatrix<T> predictedNordsieck) {

  204.         T error = getField().getZero();
  205.         for (int i = 0; i < mainSetDimension; ++i) {
  206.             final T yScale = predictedState[i].abs();
  207.             final T tol = (vecAbsoluteTolerance == null) ?
  208.                           yScale.multiply(scalRelativeTolerance).add(scalAbsoluteTolerance) :
  209.                           yScale.multiply(vecRelativeTolerance[i]).add(vecAbsoluteTolerance[i]);

  210.             // apply Taylor formula from high order to low order,
  211.             // for the sake of numerical accuracy
  212.             T variation = getField().getZero();
  213.             int sign = (predictedNordsieck.getRowDimension() & 1) == 0 ? -1 : 1;
  214.             for (int k = predictedNordsieck.getRowDimension() - 1; k >= 0; --k) {
  215.                 variation = variation.add(predictedNordsieck.getEntry(k, i).multiply(sign));
  216.                 sign      = -sign;
  217.             }
  218.             variation = variation.subtract(predictedScaled[i]);

  219.             final T ratio  = predictedState[i].subtract(previousState[i]).add(variation).divide(tol);
  220.             error = error.add(ratio.multiply(ratio));
  221.         }

  222.         return error.divide(mainSetDimension).sqrt();
  223.     }

  224.     /** {@inheritDoc} */
  225.     @Override
  226.     public FieldODEStateAndDerivative<T> integrate(final FieldExpandableODE<T> equations,
  227.                                                    final FieldODEState<T> initialState,
  228.                                                    final T finalTime)
  229.         throws NumberIsTooSmallException, DimensionMismatchException,
  230.                MaxCountExceededException, NoBracketingException {

  231.         sanityChecks(initialState, finalTime);
  232.         final T   t0 = initialState.getTime();
  233.         final T[] y  = equations.getMapper().mapState(initialState);
  234.         setStepStart(initIntegration(equations, t0, y, finalTime));
  235.         final boolean forward = finalTime.subtract(initialState.getTime()).getReal() > 0;

  236.         // compute the initial Nordsieck vector using the configured starter integrator
  237.         start(equations, getStepStart(), finalTime);

  238.         // reuse the step that was chosen by the starter integrator
  239.         FieldODEStateAndDerivative<T> stepStart = getStepStart();
  240.         FieldODEStateAndDerivative<T> stepEnd   =
  241.                         AdamsFieldStepInterpolator.taylor(stepStart,
  242.                                                           stepStart.getTime().add(getStepSize()),
  243.                                                           getStepSize(), scaled, nordsieck);

  244.         // main integration loop
  245.         setIsLastStep(false);
  246.         do {

  247.             T[] predictedY = null;
  248.             final T[] predictedScaled = MathArrays.buildArray(getField(), y.length);
  249.             Array2DRowFieldMatrix<T> predictedNordsieck = null;
  250.             T error = getField().getZero().add(10);
  251.             while (error.subtract(1.0).getReal() >= 0.0) {

  252.                 // predict a first estimate of the state at step end
  253.                 predictedY = stepEnd.getState();

  254.                 // evaluate the derivative
  255.                 final T[] yDot = computeDerivatives(stepEnd.getTime(), predictedY);

  256.                 // predict Nordsieck vector at step end
  257.                 for (int j = 0; j < predictedScaled.length; ++j) {
  258.                     predictedScaled[j] = getStepSize().multiply(yDot[j]);
  259.                 }
  260.                 predictedNordsieck = updateHighOrderDerivativesPhase1(nordsieck);
  261.                 updateHighOrderDerivativesPhase2(scaled, predictedScaled, predictedNordsieck);

  262.                 // evaluate error
  263.                 error = errorEstimation(y, predictedY, predictedScaled, predictedNordsieck);

  264.                 if (error.subtract(1.0).getReal() >= 0.0) {
  265.                     // reject the step and attempt to reduce error by stepsize control
  266.                     final T factor = computeStepGrowShrinkFactor(error);
  267.                     rescale(filterStep(getStepSize().multiply(factor), forward, false));
  268.                     stepEnd = AdamsFieldStepInterpolator.taylor(getStepStart(),
  269.                                                                 getStepStart().getTime().add(getStepSize()),
  270.                                                                 getStepSize(),
  271.                                                                 scaled,
  272.                                                                 nordsieck);
  273.                 }
  274.             }

  275.             // discrete events handling
  276.             setStepStart(acceptStep(new AdamsFieldStepInterpolator<>(getStepSize(), stepEnd,
  277.                                                                       predictedScaled, predictedNordsieck, forward,
  278.                                                                       getStepStart(), stepEnd,
  279.                                                                       equations.getMapper()),
  280.                                     finalTime));
  281.             scaled    = predictedScaled;
  282.             nordsieck = predictedNordsieck;

  283.             if (!isLastStep()) {

  284.                 System.arraycopy(predictedY, 0, y, 0, y.length);

  285.                 if (resetOccurred()) {
  286.                     // some events handler has triggered changes that
  287.                     // invalidate the derivatives, we need to restart from scratch
  288.                     start(equations, getStepStart(), finalTime);
  289.                 }

  290.                 // stepsize control for next step
  291.                 final T       factor     = computeStepGrowShrinkFactor(error);
  292.                 final T       scaledH    = getStepSize().multiply(factor);
  293.                 final T       nextT      = getStepStart().getTime().add(scaledH);
  294.                 final boolean nextIsLast = forward ?
  295.                                            nextT.subtract(finalTime).getReal() >= 0 :
  296.                                            nextT.subtract(finalTime).getReal() <= 0;
  297.                 T hNew = filterStep(scaledH, forward, nextIsLast);

  298.                 final T       filteredNextT      = getStepStart().getTime().add(hNew);
  299.                 final boolean filteredNextIsLast = forward ?
  300.                                                    filteredNextT.subtract(finalTime).getReal() >= 0 :
  301.                                                    filteredNextT.subtract(finalTime).getReal() <= 0;
  302.                 if (filteredNextIsLast) {
  303.                     hNew = finalTime.subtract(getStepStart().getTime());
  304.                 }

  305.                 rescale(hNew);
  306.                 stepEnd = AdamsFieldStepInterpolator.taylor(getStepStart(), getStepStart().getTime().add(getStepSize()),
  307.                                                             getStepSize(), scaled, nordsieck);
  308.             }
  309.         } while (!isLastStep());

  310.         final FieldODEStateAndDerivative<T> finalState = getStepStart();
  311.         setStepStart(null);
  312.         setStepSize(null);
  313.         return finalState;
  314.     }
  315. }