AdamsMoultonFieldIntegrator.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 java.util.Arrays;

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


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

  157.     /** Integrator method name. */
  158.     private static final String METHOD_NAME = "Adams-Moulton";

  159.     /**
  160.      * Build an Adams-Moulton integrator with the given order and error 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 NumberIsTooSmallException if order is 1 or less
  172.      */
  173.     public AdamsMoultonFieldIntegrator(final Field<T> field, final int nSteps,
  174.                                        final double minStep, final double maxStep,
  175.                                        final double scalAbsoluteTolerance,
  176.                                        final double scalRelativeTolerance)
  177.         throws NumberIsTooSmallException {
  178.         super(field, METHOD_NAME, nSteps, nSteps + 1, minStep, maxStep,
  179.               scalAbsoluteTolerance, scalRelativeTolerance);
  180.     }

  181.     /**
  182.      * Build an Adams-Moulton integrator with the given order and error 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 AdamsMoultonFieldIntegrator(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 + 1, minStep, maxStep,
  201.               vecAbsoluteTolerance, vecRelativeTolerance);
  202.     }

  203.     /** {@inheritDoc} */
  204.     @Override
  205.     public FieldODEStateAndDerivative<T> integrate(final FieldExpandableODE<T> equations,
  206.                                                    final FieldODEState<T> initialState,
  207.                                                    final T finalTime)
  208.         throws NumberIsTooSmallException, DimensionMismatchException,
  209.                MaxCountExceededException, NoBracketingException {

  210.         sanityChecks(initialState, finalTime);
  211.         final T   t0 = initialState.getTime();
  212.         final T[] y  = equations.getMapper().mapState(initialState);
  213.         setStepStart(initIntegration(equations, t0, y, finalTime));
  214.         final boolean forward = finalTime.subtract(initialState.getTime()).getReal() > 0;

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

  217.         // reuse the step that was chosen by the starter integrator
  218.         FieldODEStateAndDerivative<T> stepStart = getStepStart();
  219.         FieldODEStateAndDerivative<T> stepEnd   =
  220.                         AdamsFieldStepInterpolator.taylor(stepStart,
  221.                                                           stepStart.getTime().add(getStepSize()),
  222.                                                           getStepSize(), scaled, nordsieck);

  223.         // main integration loop
  224.         setIsLastStep(false);
  225.         do {

  226.             T[] predictedY = null;
  227.             final T[] predictedScaled = MathArrays.buildArray(getField(), y.length);
  228.             Array2DRowFieldMatrix<T> predictedNordsieck = null;
  229.             T error = getField().getZero().add(10);
  230.             while (error.subtract(1.0).getReal() >= 0.0) {

  231.                 // predict a first estimate of the state at step end (P in the PECE sequence)
  232.                 predictedY = stepEnd.getState();

  233.                 // evaluate a first estimate of the derivative (first E in the PECE sequence)
  234.                 final T[] yDot = computeDerivatives(stepEnd.getTime(), predictedY);

  235.                 // update Nordsieck vector
  236.                 for (int j = 0; j < predictedScaled.length; ++j) {
  237.                     predictedScaled[j] = getStepSize().multiply(yDot[j]);
  238.                 }
  239.                 predictedNordsieck = updateHighOrderDerivativesPhase1(nordsieck);
  240.                 updateHighOrderDerivativesPhase2(scaled, predictedScaled, predictedNordsieck);

  241.                 // apply correction (C in the PECE sequence)
  242.                 error = predictedNordsieck.walkInOptimizedOrder(new Corrector(y, predictedScaled, predictedY));

  243.                 if (error.subtract(1.0).getReal() >= 0.0) {
  244.                     // reject the step and attempt to reduce error by stepsize control
  245.                     final T factor = computeStepGrowShrinkFactor(error);
  246.                     rescale(filterStep(getStepSize().multiply(factor), forward, false));
  247.                     stepEnd = AdamsFieldStepInterpolator.taylor(getStepStart(),
  248.                                                                 getStepStart().getTime().add(getStepSize()),
  249.                                                                 getStepSize(),
  250.                                                                 scaled,
  251.                                                                 nordsieck);
  252.                 }
  253.             }

  254.             // evaluate a final estimate of the derivative (second E in the PECE sequence)
  255.             final T[] correctedYDot = computeDerivatives(stepEnd.getTime(), predictedY);

  256.             // update Nordsieck vector
  257.             final T[] correctedScaled = MathArrays.buildArray(getField(), y.length);
  258.             for (int j = 0; j < correctedScaled.length; ++j) {
  259.                 correctedScaled[j] = getStepSize().multiply(correctedYDot[j]);
  260.             }
  261.             updateHighOrderDerivativesPhase2(predictedScaled, correctedScaled, predictedNordsieck);

  262.             // discrete events handling
  263.             stepEnd = new FieldODEStateAndDerivative<>(stepEnd.getTime(), predictedY, correctedYDot);
  264.             setStepStart(acceptStep(new AdamsFieldStepInterpolator<>(getStepSize(), stepEnd,
  265.                                                                       correctedScaled, predictedNordsieck, forward,
  266.                                                                       getStepStart(), stepEnd,
  267.                                                                       equations.getMapper()),
  268.                                     finalTime));
  269.             scaled    = correctedScaled;
  270.             nordsieck = predictedNordsieck;

  271.             if (!isLastStep()) {

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

  273.                 if (resetOccurred()) {
  274.                     // some events handler has triggered changes that
  275.                     // invalidate the derivatives, we need to restart from scratch
  276.                     start(equations, getStepStart(), finalTime);
  277.                 }

  278.                 // stepsize control for next step
  279.                 final T  factor     = computeStepGrowShrinkFactor(error);
  280.                 final T  scaledH    = getStepSize().multiply(factor);
  281.                 final T  nextT      = getStepStart().getTime().add(scaledH);
  282.                 final boolean nextIsLast = forward ?
  283.                                            nextT.subtract(finalTime).getReal() >= 0 :
  284.                                            nextT.subtract(finalTime).getReal() <= 0;
  285.                 T hNew = filterStep(scaledH, forward, nextIsLast);

  286.                 final T  filteredNextT      = getStepStart().getTime().add(hNew);
  287.                 final boolean filteredNextIsLast = forward ?
  288.                                                    filteredNextT.subtract(finalTime).getReal() >= 0 :
  289.                                                    filteredNextT.subtract(finalTime).getReal() <= 0;
  290.                 if (filteredNextIsLast) {
  291.                     hNew = finalTime.subtract(getStepStart().getTime());
  292.                 }

  293.                 rescale(hNew);
  294.                 stepEnd = AdamsFieldStepInterpolator.taylor(getStepStart(), getStepStart().getTime().add(getStepSize()),
  295.                                                             getStepSize(), scaled, nordsieck);
  296.             }
  297.         } while (!isLastStep());

  298.         final FieldODEStateAndDerivative<T> finalState = getStepStart();
  299.         setStepStart(null);
  300.         setStepSize(null);
  301.         return finalState;
  302.     }

  303.     /** Corrector for current state in Adams-Moulton method.
  304.      * <p>
  305.      * This visitor implements the Taylor series formula:
  306.      * <pre>
  307.      * 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>
  308.      * </pre>
  309.      * </p>
  310.      */
  311.     private final class Corrector implements FieldMatrixPreservingVisitor<T> {

  312.         /** Previous state. */
  313.         private final T[] previous;

  314.         /** Current scaled first derivative. */
  315.         private final T[] scaled;

  316.         /** Current state before correction. */
  317.         private final T[] before;

  318.         /** Current state after correction. */
  319.         private final T[] after;

  320.         /** Simple constructor.
  321.          * @param previous previous state
  322.          * @param scaled current scaled first derivative
  323.          * @param state state to correct (will be overwritten after visit)
  324.          */
  325.         Corrector(final T[] previous, final T[] scaled, final T[] state) {
  326.             this.previous = previous;
  327.             this.scaled   = scaled;
  328.             this.after    = state;
  329.             this.before   = state.clone();
  330.         }

  331.         /** {@inheritDoc} */
  332.         @Override
  333.         public void start(int rows, int columns,
  334.                           int startRow, int endRow, int startColumn, int endColumn) {
  335.             Arrays.fill(after, getField().getZero());
  336.         }

  337.         /** {@inheritDoc} */
  338.         @Override
  339.         public void visit(int row, int column, T value) {
  340.             if ((row & 0x1) == 0) {
  341.                 after[column] = after[column].subtract(value);
  342.             } else {
  343.                 after[column] = after[column].add(value);
  344.             }
  345.         }

  346.         /**
  347.          * End visiting the Nordsieck vector.
  348.          * <p>The correction is used to control stepsize. So its amplitude is
  349.          * considered to be an error, which must be normalized according to
  350.          * error control settings. If the normalized value is greater than 1,
  351.          * the correction was too large and the step must be rejected.</p>
  352.          * @return the normalized correction, if greater than 1, the step
  353.          * must be rejected
  354.          */
  355.         @Override
  356.         public T end() {

  357.             T error = getField().getZero();
  358.             for (int i = 0; i < after.length; ++i) {
  359.                 after[i] = after[i].add(previous[i].add(scaled[i]));
  360.                 if (i < mainSetDimension) {
  361.                     final T yScale = RealFieldElement.max(previous[i].abs(), after[i].abs());
  362.                     final T tol = (vecAbsoluteTolerance == null) ?
  363.                                   yScale.multiply(scalRelativeTolerance).add(scalAbsoluteTolerance) :
  364.                                   yScale.multiply(vecRelativeTolerance[i]).add(vecAbsoluteTolerance[i]);
  365.                     final T ratio  = after[i].subtract(before[i]).divide(tol); // (corrected-predicted)/tol
  366.                     error = error.add(ratio.multiply(ratio));
  367.                 }
  368.             }

  369.             return error.divide(mainSetDimension).sqrt();
  370.         }
  371.     }
  372. }