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    *      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  
18  package org.apache.commons.math3.ode.nonstiff;
19  
20  import java.util.Arrays;
21  
22  import org.apache.commons.math3.exception.DimensionMismatchException;
23  import org.apache.commons.math3.exception.MaxCountExceededException;
24  import org.apache.commons.math3.exception.NoBracketingException;
25  import org.apache.commons.math3.exception.NumberIsTooSmallException;
26  import org.apache.commons.math3.linear.Array2DRowRealMatrix;
27  import org.apache.commons.math3.linear.RealMatrixPreservingVisitor;
28  import org.apache.commons.math3.ode.EquationsMapper;
29  import org.apache.commons.math3.ode.ExpandableStatefulODE;
30  import org.apache.commons.math3.ode.sampling.NordsieckStepInterpolator;
31  import org.apache.commons.math3.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   * <h3>Implementation details</h3>
60   *
61   * <p>We define scaled derivatives s<sub>i</sub>(n) at step n as:
62   * <pre>
63   * s<sub>1</sub>(n) = h y'<sub>n</sub> for first derivative
64   * s<sub>2</sub>(n) = h<sup>2</sup>/2 y''<sub>n</sub> for second derivative
65   * s<sub>3</sub>(n) = h<sup>3</sup>/6 y'''<sub>n</sub> for third derivative
66   * ...
67   * s<sub>k</sub>(n) = h<sup>k</sup>/k! y<sup>(k)</sup><sub>n</sub> for k<sup>th</sup> derivative
68   * </pre></p>
69   *
70   * <p>The definitions above use the classical representation with several previous first
71   * derivatives. Lets define
72   * <pre>
73   *   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>
74   * </pre>
75   * (we omit the k index in the notation for clarity). With these definitions,
76   * Adams-Moulton methods can be written:
77   * <ul>
78   *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n+1)</li>
79   *   <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>
80   *   <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>
81   *   <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>
82   *   <li>...</li>
83   * </ul></p>
84   *
85   * <p>Instead of using the classical representation with first derivatives only (y<sub>n</sub>,
86   * s<sub>1</sub>(n+1) and q<sub>n+1</sub>), our implementation uses the Nordsieck vector with
87   * higher degrees scaled derivatives all taken at the same step (y<sub>n</sub>, s<sub>1</sub>(n)
88   * and r<sub>n</sub>) where r<sub>n</sub> is defined as:
89   * <pre>
90   * r<sub>n</sub> = [ s<sub>2</sub>(n), s<sub>3</sub>(n) ... s<sub>k</sub>(n) ]<sup>T</sup>
91   * </pre>
92   * (here again we omit the k index in the notation for clarity)
93   * </p>
94   *
95   * <p>Taylor series formulas show that for any index offset i, s<sub>1</sub>(n-i) can be
96   * computed from s<sub>1</sub>(n), s<sub>2</sub>(n) ... s<sub>k</sub>(n), the formula being exact
97   * for degree k polynomials.
98   * <pre>
99   * s<sub>1</sub>(n-i) = s<sub>1</sub>(n) + &sum;<sub>j</sub> j (-i)<sup>j-1</sup> s<sub>j</sub>(n)
100  * </pre>
101  * The previous formula can be used with several values for i to compute the transform between
102  * classical representation and Nordsieck vector. The transform between r<sub>n</sub>
103  * and q<sub>n</sub> resulting from the Taylor series formulas above is:
104  * <pre>
105  * q<sub>n</sub> = s<sub>1</sub>(n) u + P r<sub>n</sub>
106  * </pre>
107  * where u is the [ 1 1 ... 1 ]<sup>T</sup> vector and P is the (k-1)&times;(k-1) matrix built
108  * with the j (-i)<sup>j-1</sup> terms:
109  * <pre>
110  *        [  -2   3   -4    5  ... ]
111  *        [  -4  12  -32   80  ... ]
112  *   P =  [  -6  27 -108  405  ... ]
113  *        [  -8  48 -256 1280  ... ]
114  *        [          ...           ]
115  * </pre></p>
116  *
117  * <p>Using the Nordsieck vector has several advantages:
118  * <ul>
119  *   <li>it greatly simplifies step interpolation as the interpolator mainly applies
120  *   Taylor series formulas,</li>
121  *   <li>it simplifies step changes that occur when discrete events that truncate
122  *   the step are triggered,</li>
123  *   <li>it allows to extend the methods in order to support adaptive stepsize.</li>
124  * </ul></p>
125  *
126  * <p>The predicted Nordsieck vector at step n+1 is computed from the Nordsieck vector at step
127  * n as follows:
128  * <ul>
129  *   <li>Y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li>
130  *   <li>S<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, Y<sub>n+1</sub>)</li>
131  *   <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>
132  * </ul>
133  * where A is a rows shifting matrix (the lower left part is an identity matrix):
134  * <pre>
135  *        [ 0 0   ...  0 0 | 0 ]
136  *        [ ---------------+---]
137  *        [ 1 0   ...  0 0 | 0 ]
138  *    A = [ 0 1   ...  0 0 | 0 ]
139  *        [       ...      | 0 ]
140  *        [ 0 0   ...  1 0 | 0 ]
141  *        [ 0 0   ...  0 1 | 0 ]
142  * </pre>
143  * From this predicted vector, the corrected vector is computed as follows:
144  * <ul>
145  *   <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>
146  *   <li>s<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, y<sub>n+1</sub>)</li>
147  *   <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>
148  * </ul>
149  * where the upper case Y<sub>n+1</sub>, S<sub>1</sub>(n+1) and R<sub>n+1</sub> represent the
150  * predicted states whereas the lower case y<sub>n+1</sub>, s<sub>n+1</sub> and r<sub>n+1</sub>
151  * represent the corrected states.</p>
152  *
153  * <p>The P<sup>-1</sup>u vector and the P<sup>-1</sup> A P matrix do not depend on the state,
154  * they only depend on k and therefore are precomputed once for all.</p>
155  *
156  * @version $Id: AdamsMoultonIntegrator.java 1463684 2013-04-02 19:04:13Z luc $
157  * @since 2.0
158  */
159 public class AdamsMoultonIntegrator extends AdamsIntegrator {
160 
161     /** Integrator method name. */
162     private static final String METHOD_NAME = "Adams-Moulton";
163 
164     /**
165      * Build an Adams-Moulton integrator with the given order and error control parameters.
166      * @param nSteps number of steps of the method excluding the one being computed
167      * @param minStep minimal step (sign is irrelevant, regardless of
168      * integration direction, forward or backward), the last step can
169      * be smaller than this
170      * @param maxStep maximal step (sign is irrelevant, regardless of
171      * integration direction, forward or backward), the last step can
172      * be smaller than this
173      * @param scalAbsoluteTolerance allowed absolute error
174      * @param scalRelativeTolerance allowed relative error
175      * @exception NumberIsTooSmallException if order is 1 or less
176      */
177     public AdamsMoultonIntegrator(final int nSteps,
178                                   final double minStep, final double maxStep,
179                                   final double scalAbsoluteTolerance,
180                                   final double scalRelativeTolerance)
181         throws NumberIsTooSmallException {
182         super(METHOD_NAME, nSteps, nSteps + 1, minStep, maxStep,
183               scalAbsoluteTolerance, scalRelativeTolerance);
184     }
185 
186     /**
187      * Build an Adams-Moulton integrator with the given order and error control parameters.
188      * @param nSteps number of steps of the method excluding the one being computed
189      * @param minStep minimal step (sign is irrelevant, regardless of
190      * integration direction, forward or backward), the last step can
191      * be smaller than this
192      * @param maxStep maximal step (sign is irrelevant, regardless of
193      * integration direction, forward or backward), the last step can
194      * be smaller than this
195      * @param vecAbsoluteTolerance allowed absolute error
196      * @param vecRelativeTolerance allowed relative error
197      * @exception IllegalArgumentException if order is 1 or less
198      */
199     public AdamsMoultonIntegrator(final int nSteps,
200                                   final double minStep, final double maxStep,
201                                   final double[] vecAbsoluteTolerance,
202                                   final double[] vecRelativeTolerance)
203         throws IllegalArgumentException {
204         super(METHOD_NAME, nSteps, nSteps + 1, minStep, maxStep,
205               vecAbsoluteTolerance, vecRelativeTolerance);
206     }
207 
208 
209     /** {@inheritDoc} */
210     @Override
211     public void integrate(final ExpandableStatefulODE equations,final double t)
212         throws NumberIsTooSmallException, DimensionMismatchException,
213                MaxCountExceededException, NoBracketingException {
214 
215         sanityChecks(equations, t);
216         setEquations(equations);
217         final boolean forward = t > equations.getTime();
218 
219         // initialize working arrays
220         final double[] y0   = equations.getCompleteState();
221         final double[] y    = y0.clone();
222         final double[] yDot = new double[y.length];
223         final double[] yTmp = new double[y.length];
224         final double[] predictedScaled = new double[y.length];
225         Array2DRowRealMatrix nordsieckTmp = null;
226 
227         // set up two interpolators sharing the integrator arrays
228         final NordsieckStepInterpolator interpolator = new NordsieckStepInterpolator();
229         interpolator.reinitialize(y, forward,
230                                   equations.getPrimaryMapper(), equations.getSecondaryMappers());
231 
232         // set up integration control objects
233         initIntegration(equations.getTime(), y0, t);
234 
235         // compute the initial Nordsieck vector using the configured starter integrator
236         start(equations.getTime(), y, t);
237         interpolator.reinitialize(stepStart, stepSize, scaled, nordsieck);
238         interpolator.storeTime(stepStart);
239 
240         double hNew = stepSize;
241         interpolator.rescale(hNew);
242 
243         isLastStep = false;
244         do {
245 
246             double error = 10;
247             while (error >= 1.0) {
248 
249                 stepSize = hNew;
250 
251                 // predict a first estimate of the state at step end (P in the PECE sequence)
252                 final double stepEnd = stepStart + stepSize;
253                 interpolator.setInterpolatedTime(stepEnd);
254                 final ExpandableStatefulODE expandable = getExpandable();
255                 final EquationsMapper primary = expandable.getPrimaryMapper();
256                 primary.insertEquationData(interpolator.getInterpolatedState(), yTmp);
257                 int index = 0;
258                 for (final EquationsMapper secondary : expandable.getSecondaryMappers()) {
259                     secondary.insertEquationData(interpolator.getInterpolatedSecondaryState(index), yTmp);
260                     ++index;
261                 }
262 
263                 // evaluate a first estimate of the derivative (first E in the PECE sequence)
264                 computeDerivatives(stepEnd, yTmp, yDot);
265 
266                 // update Nordsieck vector
267                 for (int j = 0; j < y0.length; ++j) {
268                     predictedScaled[j] = stepSize * yDot[j];
269                 }
270                 nordsieckTmp = updateHighOrderDerivativesPhase1(nordsieck);
271                 updateHighOrderDerivativesPhase2(scaled, predictedScaled, nordsieckTmp);
272 
273                 // apply correction (C in the PECE sequence)
274                 error = nordsieckTmp.walkInOptimizedOrder(new Corrector(y, predictedScaled, yTmp));
275 
276                 if (error >= 1.0) {
277                     // reject the step and attempt to reduce error by stepsize control
278                     final double factor = computeStepGrowShrinkFactor(error);
279                     hNew = filterStep(stepSize * factor, forward, false);
280                     interpolator.rescale(hNew);
281                 }
282             }
283 
284             // evaluate a final estimate of the derivative (second E in the PECE sequence)
285             final double stepEnd = stepStart + stepSize;
286             computeDerivatives(stepEnd, yTmp, yDot);
287 
288             // update Nordsieck vector
289             final double[] correctedScaled = new double[y0.length];
290             for (int j = 0; j < y0.length; ++j) {
291                 correctedScaled[j] = stepSize * yDot[j];
292             }
293             updateHighOrderDerivativesPhase2(predictedScaled, correctedScaled, nordsieckTmp);
294 
295             // discrete events handling
296             System.arraycopy(yTmp, 0, y, 0, y.length);
297             interpolator.reinitialize(stepEnd, stepSize, correctedScaled, nordsieckTmp);
298             interpolator.storeTime(stepStart);
299             interpolator.shift();
300             interpolator.storeTime(stepEnd);
301             stepStart = acceptStep(interpolator, y, yDot, t);
302             scaled    = correctedScaled;
303             nordsieck = nordsieckTmp;
304 
305             if (!isLastStep) {
306 
307                 // prepare next step
308                 interpolator.storeTime(stepStart);
309 
310                 if (resetOccurred) {
311                     // some events handler has triggered changes that
312                     // invalidate the derivatives, we need to restart from scratch
313                     start(stepStart, y, t);
314                     interpolator.reinitialize(stepStart, stepSize, scaled, nordsieck);
315 
316                 }
317 
318                 // stepsize control for next step
319                 final double  factor     = computeStepGrowShrinkFactor(error);
320                 final double  scaledH    = stepSize * factor;
321                 final double  nextT      = stepStart + scaledH;
322                 final boolean nextIsLast = forward ? (nextT >= t) : (nextT <= t);
323                 hNew = filterStep(scaledH, forward, nextIsLast);
324 
325                 final double  filteredNextT      = stepStart + hNew;
326                 final boolean filteredNextIsLast = forward ? (filteredNextT >= t) : (filteredNextT <= t);
327                 if (filteredNextIsLast) {
328                     hNew = t - stepStart;
329                 }
330 
331                 interpolator.rescale(hNew);
332             }
333 
334         } while (!isLastStep);
335 
336         // dispatch results
337         equations.setTime(stepStart);
338         equations.setCompleteState(y);
339 
340         resetInternalState();
341 
342     }
343 
344     /** Corrector for current state in Adams-Moulton method.
345      * <p>
346      * This visitor implements the Taylor series formula:
347      * <pre>
348      * 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>
349      * </pre>
350      * </p>
351      */
352     private class Corrector implements RealMatrixPreservingVisitor {
353 
354         /** Previous state. */
355         private final double[] previous;
356 
357         /** Current scaled first derivative. */
358         private final double[] scaled;
359 
360         /** Current state before correction. */
361         private final double[] before;
362 
363         /** Current state after correction. */
364         private final double[] after;
365 
366         /** Simple constructor.
367          * @param previous previous state
368          * @param scaled current scaled first derivative
369          * @param state state to correct (will be overwritten after visit)
370          */
371         public Corrector(final double[] previous, final double[] scaled, final double[] state) {
372             this.previous = previous;
373             this.scaled   = scaled;
374             this.after    = state;
375             this.before   = state.clone();
376         }
377 
378         /** {@inheritDoc} */
379         public void start(int rows, int columns,
380                           int startRow, int endRow, int startColumn, int endColumn) {
381             Arrays.fill(after, 0.0);
382         }
383 
384         /** {@inheritDoc} */
385         public void visit(int row, int column, double value) {
386             if ((row & 0x1) == 0) {
387                 after[column] -= value;
388             } else {
389                 after[column] += value;
390             }
391         }
392 
393         /**
394          * End visiting the Nordsieck vector.
395          * <p>The correction is used to control stepsize. So its amplitude is
396          * considered to be an error, which must be normalized according to
397          * error control settings. If the normalized value is greater than 1,
398          * the correction was too large and the step must be rejected.</p>
399          * @return the normalized correction, if greater than 1, the step
400          * must be rejected
401          */
402         public double end() {
403 
404             double error = 0;
405             for (int i = 0; i < after.length; ++i) {
406                 after[i] += previous[i] + scaled[i];
407                 if (i < mainSetDimension) {
408                     final double yScale = FastMath.max(FastMath.abs(previous[i]), FastMath.abs(after[i]));
409                     final double tol = (vecAbsoluteTolerance == null) ?
410                                        (scalAbsoluteTolerance + scalRelativeTolerance * yScale) :
411                                        (vecAbsoluteTolerance[i] + vecRelativeTolerance[i] * yScale);
412                     final double ratio  = (after[i] - before[i]) / tol;
413                     error += ratio * ratio;
414                 }
415             }
416 
417             return FastMath.sqrt(error / mainSetDimension);
418 
419         }
420     }
421 
422 }