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