AdamsNordsieckTransformer.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 java.util.HashMap;
  20. import java.util.Map;

  21. import org.hipparchus.fraction.BigFraction;
  22. import org.hipparchus.linear.Array2DRowFieldMatrix;
  23. import org.hipparchus.linear.Array2DRowRealMatrix;
  24. import org.hipparchus.linear.ArrayFieldVector;
  25. import org.hipparchus.linear.FieldDecompositionSolver;
  26. import org.hipparchus.linear.FieldLUDecomposition;
  27. import org.hipparchus.linear.FieldMatrix;
  28. import org.hipparchus.linear.MatrixUtils;
  29. import org.hipparchus.linear.QRDecomposition;
  30. import org.hipparchus.linear.RealMatrix;

  31. /** Transformer to Nordsieck vectors for Adams integrators.
  32.  * <p>This class is used by {@link AdamsBashforthIntegrator Adams-Bashforth} and
  33.  * {@link AdamsMoultonIntegrator Adams-Moulton} integrators to convert between
  34.  * classical representation with several previous first derivatives and Nordsieck
  35.  * representation with higher order scaled derivatives.</p>
  36.  *
  37.  * <p>We define scaled derivatives s<sub>i</sub>(n) at step n as:
  38.  * \[
  39.  *   \left\{\begin{align}
  40.  *   s_1(n) &amp;= h y'_n \text{ for first derivative}\\
  41.  *   s_2(n) &amp;= \frac{h^2}{2} y_n'' \text{ for second derivative}\\
  42.  *   s_3(n) &amp;= \frac{h^3}{6} y_n''' \text{ for third derivative}\\
  43.  *   &amp;\cdots\\
  44.  *   s_k(n) &amp;= \frac{h^k}{k!} y_n^{(k)} \text{ for } k^\mathrm{th} \text{ derivative}
  45.  *   \end{align}\right.
  46.  * \]</p>
  47.  *
  48.  * <p>With the previous definition, the classical representation of multistep methods
  49.  * uses first derivatives only, i.e. it handles y<sub>n</sub>, s<sub>1</sub>(n) and
  50.  * q<sub>n</sub> where q<sub>n</sub> is defined as:
  51.  * \[
  52.  *   q_n = [ s_1(n-1) s_1(n-2) \ldots s_1(n-(k-1)) ]^T
  53.  * \]
  54.  * (we omit the k index in the notation for clarity).</p>
  55.  *
  56.  * <p>Another possible representation uses the Nordsieck vector with
  57.  * higher degrees scaled derivatives all taken at the same step, i.e it handles y<sub>n</sub>,
  58.  * s<sub>1</sub>(n) and r<sub>n</sub>) where r<sub>n</sub> is defined as:
  59.  * \[
  60.  * r_n = [ s_2(n), s_3(n) \ldots s_k(n) ]^T
  61.  * \]
  62.  * (here again we omit the k index in the notation for clarity)
  63.  * </p>
  64.  *
  65.  * <p>Taylor series formulas show that for any index offset i, s<sub>1</sub>(n-i) can be
  66.  * computed from s<sub>1</sub>(n), s<sub>2</sub>(n) ... s<sub>k</sub>(n), the formula being exact
  67.  * for degree k polynomials.
  68.  * \[
  69.  * s_1(n-i) = s_1(n) + \sum_{j\gt 0} (j+1) (-i)^j s_{j+1}(n)
  70.  * \]
  71.  * The previous formula can be used with several values for i to compute the transform between
  72.  * classical representation and Nordsieck vector at step end. The transform between r<sub>n</sub>
  73.  * and q<sub>n</sub> resulting from the Taylor series formulas above is:
  74.  * \[
  75.  * q_n = s_1(n) u + P r_n
  76.  * \]
  77.  * where u is the [ 1 1 ... 1 ]<sup>T</sup> vector and P is the (k-1)&times;(k-1) matrix built
  78.  * with the \((j+1) (-i)^j\) terms with i being the row number starting from 1 and j being
  79.  * the column number starting from 1:
  80.  * \[
  81.  *   P=\begin{bmatrix}
  82.  *   -2  &amp;  3 &amp;   -4 &amp;    5 &amp; \ldots \\
  83.  *   -4  &amp; 12 &amp;  -32 &amp;   80 &amp; \ldots \\
  84.  *   -6  &amp; 27 &amp; -108 &amp;  405 &amp; \ldots \\
  85.  *   -8  &amp; 48 &amp; -256 &amp; 1280 &amp; \ldots \\
  86.  *       &amp;    &amp;  \ldots\\
  87.  *    \end{bmatrix}
  88.  * \]
  89.  *
  90.  * <p>Changing -i into +i in the formula above can be used to compute a similar transform between
  91.  * classical representation and Nordsieck vector at step start. The resulting matrix is simply
  92.  * the absolute value of matrix P.</p>
  93.  *
  94.  * <p>For {@link AdamsBashforthIntegrator Adams-Bashforth} method, the Nordsieck vector
  95.  * at step n+1 is computed from the Nordsieck vector at step n as follows:
  96.  * <ul>
  97.  *   <li>y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li>
  98.  *   <li>s<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, y<sub>n+1</sub>)</li>
  99.  *   <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>
  100.  * </ul>
  101.  * <p>where A is a rows shifting matrix (the lower left part is an identity matrix):</p>
  102.  * <pre>
  103.  *        [ 0 0   ...  0 0 | 0 ]
  104.  *        [ ---------------+---]
  105.  *        [ 1 0   ...  0 0 | 0 ]
  106.  *    A = [ 0 1   ...  0 0 | 0 ]
  107.  *        [       ...      | 0 ]
  108.  *        [ 0 0   ...  1 0 | 0 ]
  109.  *        [ 0 0   ...  0 1 | 0 ]
  110.  * </pre>
  111.  *
  112.  * <p>For {@link AdamsMoultonIntegrator Adams-Moulton} method, the predicted Nordsieck vector
  113.  * at step n+1 is computed from the Nordsieck vector at step n as follows:
  114.  * <ul>
  115.  *   <li>Y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li>
  116.  *   <li>S<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, Y<sub>n+1</sub>)</li>
  117.  *   <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>
  118.  * </ul>
  119.  * From this predicted vector, the corrected vector is computed as follows:
  120.  * <ul>
  121.  *   <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>
  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> = R<sub>n+1</sub> + (s<sub>1</sub>(n+1) - S<sub>1</sub>(n+1)) P<sup>-1</sup> u</li>
  124.  * </ul>
  125.  * <p>where the upper case Y<sub>n+1</sub>, S<sub>1</sub>(n+1) and R<sub>n+1</sub> represent the
  126.  * predicted states whereas the lower case y<sub>n+1</sub>, s<sub>n+1</sub> and r<sub>n+1</sub>
  127.  * represent the corrected states.</p>
  128.  *
  129.  * <p>We observe that both methods use similar update formulas. In both cases a P<sup>-1</sup>u
  130.  * vector and a P<sup>-1</sup> A P matrix are used that do not depend on the state,
  131.  * they only depend on k. This class handles these transformations.</p>
  132.  *
  133.  */
  134. public class AdamsNordsieckTransformer {

  135.     /** Cache for already computed coefficients. */
  136.     private static final Map<Integer, AdamsNordsieckTransformer> CACHE = new HashMap<>();

  137.     /** Update matrix for the higher order derivatives h<sup>2</sup>/2 y'', h<sup>3</sup>/6 y''' ... */
  138.     private final Array2DRowRealMatrix update;

  139.     /** Update coefficients of the higher order derivatives wrt y'. */
  140.     private final double[] c1;

  141.     /** Simple constructor.
  142.      * @param n number of steps of the multistep method
  143.      * (excluding the one being computed)
  144.      */
  145.     private AdamsNordsieckTransformer(final int n) {

  146.         final int rows = n - 1;

  147.         // compute exact coefficients
  148.         FieldMatrix<BigFraction> bigP = buildP(rows);
  149.         FieldDecompositionSolver<BigFraction> pSolver =
  150.                 new FieldLUDecomposition<>(bigP).getSolver();

  151.         BigFraction[] u = new BigFraction[rows];
  152.         Arrays.fill(u, BigFraction.ONE);
  153.         BigFraction[] bigC1 = pSolver.solve(new ArrayFieldVector<>(u, false)).toArray();

  154.         // update coefficients are computed by combining transform from
  155.         // Nordsieck to multistep, then shifting rows to represent step advance
  156.         // then applying inverse transform
  157.         BigFraction[][] shiftedP = bigP.getData();
  158.         for (int i = shiftedP.length - 1; i > 0; --i) {
  159.             // shift rows
  160.             shiftedP[i] = shiftedP[i - 1];
  161.         }
  162.         shiftedP[0] = new BigFraction[rows];
  163.         Arrays.fill(shiftedP[0], BigFraction.ZERO);
  164.         FieldMatrix<BigFraction> bigMSupdate =
  165.             pSolver.solve(new Array2DRowFieldMatrix<>(shiftedP, false));

  166.         // convert coefficients to double
  167.         update         = MatrixUtils.bigFractionMatrixToRealMatrix(bigMSupdate);
  168.         c1             = new double[rows];
  169.         for (int i = 0; i < rows; ++i) {
  170.             c1[i] = bigC1[i].doubleValue();
  171.         }

  172.     }

  173.     /** Get the Nordsieck transformer for a given number of steps.
  174.      * @param nSteps number of steps of the multistep method
  175.      * (excluding the one being computed)
  176.      * @return Nordsieck transformer for the specified number of steps
  177.      */
  178.     public static AdamsNordsieckTransformer getInstance(final int nSteps) { // NOPMD - PMD false positive
  179.         synchronized(CACHE) {
  180.             AdamsNordsieckTransformer t = CACHE.get(nSteps);
  181.             if (t == null) {
  182.                 t = new AdamsNordsieckTransformer(nSteps);
  183.                 CACHE.put(nSteps, t);
  184.             }
  185.             return t;
  186.         }
  187.     }

  188.     /** Build the P matrix.
  189.      * <p>The P matrix general terms are shifted \((j+1) (-i)^j\) terms
  190.      * with i being the row number starting from 1 and j being the column
  191.      * number starting from 1:
  192.      * <pre>
  193.      *        [  -2   3   -4    5  ... ]
  194.      *        [  -4  12  -32   80  ... ]
  195.      *   P =  [  -6  27 -108  405  ... ]
  196.      *        [  -8  48 -256 1280  ... ]
  197.      *        [          ...           ]
  198.      * </pre></p>
  199.      * @param rows number of rows of the matrix
  200.      * @return P matrix
  201.      */
  202.     private FieldMatrix<BigFraction> buildP(final int rows) {

  203.         final BigFraction[][] pData = new BigFraction[rows][rows];

  204.         for (int i = 1; i <= pData.length; ++i) {
  205.             // build the P matrix elements from Taylor series formulas
  206.             final BigFraction[] pI = pData[i - 1];
  207.             final int factor = -i;
  208.             int aj = factor;
  209.             for (int j = 1; j <= pI.length; ++j) {
  210.                 pI[j - 1] = new BigFraction(aj * (j + 1));
  211.                 aj *= factor;
  212.             }
  213.         }

  214.         return new Array2DRowFieldMatrix<>(pData, false);

  215.     }

  216.     /** Initialize the high order scaled derivatives at step start.
  217.      * @param h step size to use for scaling
  218.      * @param t first steps times
  219.      * @param y first steps states
  220.      * @param yDot first steps derivatives
  221.      * @return Nordieck vector at start of first step (h<sup>2</sup>/2 y''<sub>n</sub>,
  222.      * h<sup>3</sup>/6 y'''<sub>n</sub> ... h<sup>k</sup>/k! y<sup>(k)</sup><sub>n</sub>)
  223.      */

  224.     public Array2DRowRealMatrix initializeHighOrderDerivatives(final double h, final double[] t,
  225.                                                                final double[][] y,
  226.                                                                final double[][] yDot) {

  227.         // using Taylor series with di = ti - t0, we get:
  228.         //  y(ti)  - y(t0)  - di y'(t0) =   di^2 / h^2 s2 + ... +   di^k     / h^k sk + O(h^k)
  229.         //  y'(ti) - y'(t0)             = 2 di   / h^2 s2 + ... + k di^(k-1) / h^k sk + O(h^(k-1))
  230.         // we write these relations for i = 1 to i= 1+n/2 as a set of n + 2 linear
  231.         // equations depending on the Nordsieck vector [s2 ... sk rk], so s2 to sk correspond
  232.         // to the appropriately truncated Taylor expansion, and rk is the Taylor remainder.
  233.         // The goal is to have s2 to sk as accurate as possible considering the fact the sum is
  234.         // truncated and we don't want the error terms to be included in s2 ... sk, so we need
  235.         // to solve also for the remainder
  236.         final double[][] a     = new double[c1.length + 1][c1.length + 1];
  237.         final double[][] b     = new double[c1.length + 1][y[0].length];
  238.         final double[]   y0    = y[0];
  239.         final double[]   yDot0 = yDot[0];
  240.         for (int i = 1; i < y.length; ++i) {

  241.             final double di    = t[i] - t[0];
  242.             final double ratio = di / h;
  243.             double dikM1Ohk    =  1 / h;

  244.             // linear coefficients of equations
  245.             // y(ti) - y(t0) - di y'(t0) and y'(ti) - y'(t0)
  246.             final double[] aI    = a[2 * i - 2];
  247.             final double[] aDotI = (2 * i - 1) < a.length ? a[2 * i - 1] : null;
  248.             for (int j = 0; j < aI.length; ++j) {
  249.                 dikM1Ohk *= ratio;
  250.                 aI[j]     = di      * dikM1Ohk;
  251.                 if (aDotI != null) {
  252.                     aDotI[j]  = (j + 2) * dikM1Ohk;
  253.                 }
  254.             }

  255.             // expected value of the previous equations
  256.             final double[] yI    = y[i];
  257.             final double[] yDotI = yDot[i];
  258.             final double[] bI    = b[2 * i - 2];
  259.             final double[] bDotI = (2 * i - 1) < b.length ? b[2 * i - 1] : null;
  260.             for (int j = 0; j < yI.length; ++j) {
  261.                 bI[j]    = yI[j] - y0[j] - di * yDot0[j];
  262.                 if (bDotI != null) {
  263.                     bDotI[j] = yDotI[j] - yDot0[j];
  264.                 }
  265.             }

  266.         }

  267.         // solve the linear system to get the best estimate of the Nordsieck vector [s2 ... sk],
  268.         // with the additional terms s(k+1) and c grabbing the parts after the truncated Taylor expansion
  269.         final QRDecomposition decomposition = new QRDecomposition(new Array2DRowRealMatrix(a, false));
  270.         final RealMatrix x = decomposition.getSolver().solve(new Array2DRowRealMatrix(b, false));

  271.         // extract just the Nordsieck vector [s2 ... sk]
  272.         final Array2DRowRealMatrix truncatedX = new Array2DRowRealMatrix(x.getRowDimension() - 1, x.getColumnDimension());
  273.         for (int i = 0; i < truncatedX.getRowDimension(); ++i) {
  274.             for (int j = 0; j < truncatedX.getColumnDimension(); ++j) {
  275.                 truncatedX.setEntry(i, j, x.getEntry(i, j));
  276.             }
  277.         }
  278.         return truncatedX;

  279.     }

  280.     /** Update the high order scaled derivatives for Adams integrators (phase 1).
  281.      * <p>The complete update of high order derivatives has a form similar to:
  282.      * \[
  283.      * r_{n+1} = (s_1(n) - s_1(n+1)) P^{-1} u + P^{-1} A P r_n
  284.      * \]
  285.      * this method computes the P<sup>-1</sup> A P r<sub>n</sub> part.</p>
  286.      * @param highOrder high order scaled derivatives
  287.      * (h<sup>2</sup>/2 y'', ... h<sup>k</sup>/k! y(k))
  288.      * @return updated high order derivatives
  289.      * @see #updateHighOrderDerivativesPhase2(double[], double[], Array2DRowRealMatrix)
  290.      */
  291.     public Array2DRowRealMatrix updateHighOrderDerivativesPhase1(final Array2DRowRealMatrix highOrder) {
  292.         return update.multiply(highOrder);
  293.     }

  294.     /** Update the high order scaled derivatives Adams integrators (phase 2).
  295.      * <p>The complete update of high order derivatives has a form similar to:
  296.      * \[
  297.      * r_{n+1} = (s_1(n) - s_1(n+1)) P^{-1} u + P^{-1} A P r_n
  298.      * \]
  299.      * this method computes the (s<sub>1</sub>(n) - s<sub>1</sub>(n+1)) P<sup>-1</sup> u part.</p>
  300.      * <p>Phase 1 of the update must already have been performed.</p>
  301.      * @param start first order scaled derivatives at step start
  302.      * @param end first order scaled derivatives at step end
  303.      * @param highOrder high order scaled derivatives, will be modified
  304.      * (h<sup>2</sup>/2 y'', ... h<sup>k</sup>/k! y(k))
  305.      * @see #updateHighOrderDerivativesPhase1(Array2DRowRealMatrix)
  306.      */
  307.     public void updateHighOrderDerivativesPhase2(final double[] start,
  308.                                                  final double[] end,
  309.                                                  final Array2DRowRealMatrix highOrder) {
  310.         final double[][] data = highOrder.getDataRef();
  311.         for (int i = 0; i < data.length; ++i) {
  312.             final double[] dataI = data[i];
  313.             final double c1I = c1[i];
  314.             for (int j = 0; j < dataI.length; ++j) {
  315.                 dataI[j] += c1I * (start[j] - end[j]);
  316.             }
  317.         }
  318.     }

  319. }