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