View Javadoc
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  
18  package org.hipparchus.ode.nonstiff;
19  
20  import java.util.Arrays;
21  
22  import org.hipparchus.exception.MathIllegalArgumentException;
23  import org.hipparchus.exception.MathIllegalStateException;
24  import org.hipparchus.linear.Array2DRowRealMatrix;
25  import org.hipparchus.linear.RealMatrix;
26  import org.hipparchus.linear.RealMatrixPreservingVisitor;
27  import org.hipparchus.ode.EquationsMapper;
28  import org.hipparchus.ode.LocalizedODEFormats;
29  import org.hipparchus.ode.ODEStateAndDerivative;
30  import org.hipparchus.util.FastMath;
31  
32  
33  /**
34   * This class implements implicit Adams-Moulton integrators for Ordinary
35   * Differential Equations.
36   *
37   * <p>Adams-Moulton methods (in fact due to Adams alone) are implicit
38   * multistep ODE solvers. This implementation is a variation of the classical
39   * one: it uses adaptive stepsize to implement error control, whereas
40   * classical implementations are fixed step size. The value of state vector
41   * at step n+1 is a simple combination of the value at step n and of the
42   * derivatives at steps n+1, n, n-1 ... Since y'<sub>n+1</sub> is needed to
43   * compute y<sub>n+1</sub>, another method must be used to compute a first
44   * estimate of y<sub>n+1</sub>, then compute y'<sub>n+1</sub>, then compute
45   * a final estimate of y<sub>n+1</sub> using the following formulas. Depending
46   * on the number k of previous steps one wants to use for computing the next
47   * value, different formulas are available for the final estimate:</p>
48   * <ul>
49   *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + h y'<sub>n+1</sub></li>
50   *   <li>k = 2: y<sub>n+1</sub> = y<sub>n</sub> + h (y'<sub>n+1</sub>+y'<sub>n</sub>)/2</li>
51   *   <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>
52   *   <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>
53   *   <li>...</li>
54   * </ul>
55   *
56   * <p>A k-steps Adams-Moulton method is of order k+1.</p>
57   *
58   * <p> There must be sufficient time for the {@link #setStarterIntegrator(org.hipparchus.ode.ODEIntegrator)
59   * starter integrator} to take several steps between the the last reset event, and the end
60   * of integration, otherwise an exception may be thrown during integration. The user can
61   * adjust the end date of integration, or the step size of the starter integrator to
62   * ensure a sufficient number of steps can be completed before the end of integration.
63   * </p>
64   *
65   * <p><strong>Implementation details</strong></p>
66   *
67   * <p>We define scaled derivatives s<sub>i</sub>(n) at step n as:
68   * \[
69   *   \left\{\begin{align}
70   *   s_1(n) &amp;= h y'_n \text{ for first derivative}\\
71   *   s_2(n) &amp;= \frac{h^2}{2} y_n'' \text{ for second derivative}\\
72   *   s_3(n) &amp;= \frac{h^3}{6} y_n''' \text{ for third derivative}\\
73   *   &amp;\cdots\\
74   *   s_k(n) &amp;= \frac{h^k}{k!} y_n^{(k)} \text{ for } k^\mathrm{th} \text{ derivative}
75   *   \end{align}\right.
76   * \]</p>
77   *
78   * <p>The definitions above use the classical representation with several previous first
79   * derivatives. Lets define
80   * \[
81   *   q_n = [ s_1(n-1) s_1(n-2) \ldots s_1(n-(k-1)) ]^T
82   * \]
83   * (we omit the k index in the notation for clarity). With these definitions,
84   * Adams-Moulton methods can be written:</p>
85   * <ul>
86   *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n+1)</li>
87   *   <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>
88   *   <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>
89   *   <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>
90   *   <li>...</li>
91   * </ul>
92   *
93   * <p>Instead of using the classical representation with first derivatives only (y<sub>n</sub>,
94   * s<sub>1</sub>(n+1) and q<sub>n+1</sub>), our implementation uses the Nordsieck vector with
95   * higher degrees scaled derivatives all taken at the same step (y<sub>n</sub>, s<sub>1</sub>(n)
96   * and r<sub>n</sub>) where r<sub>n</sub> is defined as:
97   * \[
98   * r_n = [ s_2(n), s_3(n) \ldots s_k(n) ]^T
99   * \]
100  * (here again we omit the k index in the notation for clarity)
101  * </p>
102  *
103  * <p>Taylor series formulas show that for any index offset i, s<sub>1</sub>(n-i) can be
104  * computed from s<sub>1</sub>(n), s<sub>2</sub>(n) ... s<sub>k</sub>(n), the formula being exact
105  * for degree k polynomials.
106  * \[
107  * s_1(n-i) = s_1(n) + \sum_{j\gt 0} (j+1) (-i)^j s_{j+1}(n)
108  * \]
109  * The previous formula can be used with several values for i to compute the transform between
110  * classical representation and Nordsieck vector. The transform between r<sub>n</sub>
111  * and q<sub>n</sub> resulting from the Taylor series formulas above is:
112  * \[
113  * q_n = s_1(n) u + P r_n
114  * \]
115  * where u is the [ 1 1 ... 1 ]<sup>T</sup> vector and P is the (k-1)&times;(k-1) matrix built
116  * with the \((j+1) (-i)^j\) terms with i being the row number starting from 1 and j being
117  * the column number starting from 1:
118  * \[
119  *   P=\begin{bmatrix}
120  *   -2  &amp;  3 &amp;   -4 &amp;    5 &amp; \ldots \\
121  *   -4  &amp; 12 &amp;  -32 &amp;   80 &amp; \ldots \\
122  *   -6  &amp; 27 &amp; -108 &amp;  405 &amp; \ldots \\
123  *   -8  &amp; 48 &amp; -256 &amp; 1280 &amp; \ldots \\
124  *       &amp;    &amp;  \ldots\\
125  *    \end{bmatrix}
126  * \]
127  *
128  * <p>Using the Nordsieck vector has several advantages:</p>
129  * <ul>
130  *   <li>it greatly simplifies step interpolation as the interpolator mainly applies
131  *   Taylor series formulas,</li>
132  *   <li>it simplifies step changes that occur when discrete events that truncate
133  *   the step are triggered,</li>
134  *   <li>it allows to extend the methods in order to support adaptive stepsize.</li>
135  * </ul>
136  *
137  * <p>The predicted Nordsieck vector at step n+1 is computed from the Nordsieck vector at step
138  * n as follows:
139  * <ul>
140  *   <li>Y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li>
141  *   <li>S<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, Y<sub>n+1</sub>)</li>
142  *   <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>
143  * </ul>
144  * where A is a rows shifting matrix (the lower left part is an identity matrix):
145  * <pre>
146  *        [ 0 0   ...  0 0 | 0 ]
147  *        [ ---------------+---]
148  *        [ 1 0   ...  0 0 | 0 ]
149  *    A = [ 0 1   ...  0 0 | 0 ]
150  *        [       ...      | 0 ]
151  *        [ 0 0   ...  1 0 | 0 ]
152  *        [ 0 0   ...  0 1 | 0 ]
153  * </pre>
154  * From this predicted vector, the corrected vector is computed as follows:
155  * <ul>
156  *   <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>
157  *   <li>s<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, y<sub>n+1</sub>)</li>
158  *   <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>
159  * </ul>
160  * <p>where the upper case Y<sub>n+1</sub>, S<sub>1</sub>(n+1) and R<sub>n+1</sub> represent the
161  * predicted states whereas the lower case y<sub>n+1</sub>, s<sub>n+1</sub> and r<sub>n+1</sub>
162  * represent the corrected states.</p>
163  *
164  * <p>The P<sup>-1</sup>u vector and the P<sup>-1</sup> A P matrix do not depend on the state,
165  * they only depend on k and therefore are precomputed once for all.</p>
166  *
167  */
168 public class AdamsMoultonIntegrator extends AdamsIntegrator {
169 
170     /** Integrator method name. */
171     private static final String METHOD_NAME = "Adams-Moulton";
172 
173     /**
174      * Build an Adams-Moulton integrator with the given order and error control parameters.
175      * @param nSteps number of steps of the method excluding the one being computed
176      * @param minStep minimal step (sign is irrelevant, regardless of
177      * integration direction, forward or backward), the last step can
178      * be smaller than this
179      * @param maxStep maximal step (sign is irrelevant, regardless of
180      * integration direction, forward or backward), the last step can
181      * be smaller than this
182      * @param scalAbsoluteTolerance allowed absolute error
183      * @param scalRelativeTolerance allowed relative error
184      * @exception MathIllegalArgumentException if order is 1 or less
185      */
186     public AdamsMoultonIntegrator(final int nSteps,
187                                   final double minStep, final double maxStep,
188                                   final double scalAbsoluteTolerance,
189                                   final double scalRelativeTolerance)
190         throws MathIllegalArgumentException {
191         super(METHOD_NAME, nSteps, nSteps + 1, minStep, maxStep,
192               scalAbsoluteTolerance, scalRelativeTolerance);
193     }
194 
195     /**
196      * Build an Adams-Moulton integrator with the given order and error control parameters.
197      * @param nSteps number of steps of the method excluding the one being computed
198      * @param minStep minimal step (sign is irrelevant, regardless of
199      * integration direction, forward or backward), the last step can
200      * be smaller than this
201      * @param maxStep maximal step (sign is irrelevant, regardless of
202      * integration direction, forward or backward), the last step can
203      * be smaller than this
204      * @param vecAbsoluteTolerance allowed absolute error
205      * @param vecRelativeTolerance allowed relative error
206      * @exception IllegalArgumentException if order is 1 or less
207      */
208     public AdamsMoultonIntegrator(final int nSteps,
209                                   final double minStep, final double maxStep,
210                                   final double[] vecAbsoluteTolerance,
211                                   final double[] vecRelativeTolerance)
212         throws IllegalArgumentException {
213         super(METHOD_NAME, nSteps, nSteps + 1, minStep, maxStep,
214               vecAbsoluteTolerance, vecRelativeTolerance);
215     }
216 
217     /** {@inheritDoc} */
218     @Override
219     protected double errorEstimation(final double[] previousState, final double predictedTime,
220                                      final double[] predictedState,
221                                      final double[] predictedScaled,
222                                      final RealMatrix predictedNordsieck) {
223         final double error = predictedNordsieck.walkInOptimizedOrder(new Corrector(previousState, predictedScaled, predictedState));
224         if (Double.isNaN(error)) {
225             throw new MathIllegalStateException(LocalizedODEFormats.NAN_APPEARING_DURING_INTEGRATION,
226                                                 predictedTime);
227         }
228         return error;
229     }
230 
231     /** {@inheritDoc} */
232     @Override
233     protected AdamsStateInterpolator finalizeStep(final double stepSize, final double[] predictedState,
234                                                   final double[] predictedScaled, final Array2DRowRealMatrix predictedNordsieck,
235                                                   final boolean isForward,
236                                                   final ODEStateAndDerivative globalPreviousState,
237                                                   final ODEStateAndDerivative globalCurrentState,
238                                                   final EquationsMapper equationsMapper) {
239 
240         final double[] correctedYDot = computeDerivatives(globalCurrentState.getTime(), predictedState);
241 
242         // update Nordsieck vector
243         final double[] correctedScaled = new double[predictedState.length];
244         for (int j = 0; j < correctedScaled.length; ++j) {
245             correctedScaled[j] = getStepSize() * correctedYDot[j];
246         }
247         updateHighOrderDerivativesPhase2(predictedScaled, correctedScaled, predictedNordsieck);
248 
249         final ODEStateAndDerivative updatedStepEnd =
250                         equationsMapper.mapStateAndDerivative(globalCurrentState.getTime(),
251                                                               predictedState, correctedYDot);
252         return new AdamsStateInterpolator(getStepSize(), updatedStepEnd,
253                                           correctedScaled, predictedNordsieck, isForward,
254                                           getStepStart(), updatedStepEnd,
255                                           equationsMapper);
256 
257     }
258 
259     /** Corrector for current state in Adams-Moulton method.
260      * <p>
261      * This visitor implements the Taylor series formula:
262      * <pre>
263      * 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>
264      * </pre>
265      * </p>
266      */
267     private class Corrector implements RealMatrixPreservingVisitor {
268 
269         /** Previous state. */
270         private final double[] previous;
271 
272         /** Current scaled first derivative. */
273         private final double[] scaled;
274 
275         /** Current state before correction. */
276         private final double[] before;
277 
278         /** Current state after correction. */
279         private final double[] after;
280 
281         /** Simple constructor.
282          * <p>
283          * All arrays will be stored by reference to caller arrays.
284          * </p>
285          * @param previous previous state
286          * @param scaled current scaled first derivative
287          * @param state state to correct (will be overwritten after visit)
288          */
289         Corrector(final double[] previous, final double[] scaled, final double[] state) {
290             this.previous = previous; // NOPMD - array reference storage is intentional and documented here
291             this.scaled   = scaled;   // NOPMD - array reference storage is intentional and documented here
292             this.after    = state;    // NOPMD - array reference storage is intentional and documented here
293             this.before   = state.clone();
294         }
295 
296         /** {@inheritDoc} */
297         @Override
298         public void start(int rows, int columns,
299                           int startRow, int endRow, int startColumn, int endColumn) {
300             Arrays.fill(after, 0.0);
301         }
302 
303         /** {@inheritDoc} */
304         @Override
305         public void visit(int row, int column, double value) {
306             if ((row & 0x1) == 0) {
307                 after[column] -= value;
308             } else {
309                 after[column] += value;
310             }
311         }
312 
313         /**
314          * End visiting the Nordsieck vector.
315          * <p>The correction is used to control stepsize. So its amplitude is
316          * considered to be an error, which must be normalized according to
317          * error control settings. If the normalized value is greater than 1,
318          * the correction was too large and the step must be rejected.</p>
319          * @return the normalized correction, if greater than 1, the step
320          * must be rejected
321          */
322         @Override
323         public double end() {
324 
325             final StepsizeHelper helper = getStepSizeHelper();
326             double error = 0;
327             for (int i = 0; i < after.length; ++i) {
328                 after[i] += previous[i] + scaled[i];
329                 if (i < helper.getMainSetDimension()) {
330                     final double tol   = helper.getTolerance(i, FastMath.max(FastMath.abs(previous[i]), FastMath.abs(after[i])));
331                     final double ratio = (after[i] - before[i]) / tol; // (corrected-predicted)/tol
332                     error += ratio * ratio;
333                 }
334             }
335 
336             return FastMath.sqrt(error / helper.getMainSetDimension());
337 
338         }
339     }
340 
341 }