001    /*
002     * Licensed to the Apache Software Foundation (ASF) under one or more
003     * contributor license agreements.  See the NOTICE file distributed with
004     * this work for additional information regarding copyright ownership.
005     * The ASF licenses this file to You under the Apache License, Version 2.0
006     * (the "License"); you may not use this file except in compliance with
007     * the License.  You may obtain a copy of the License at
008     *
009     *      http://www.apache.org/licenses/LICENSE-2.0
010     *
011     * Unless required by applicable law or agreed to in writing, software
012     * distributed under the License is distributed on an "AS IS" BASIS,
013     * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
014     * See the License for the specific language governing permissions and
015     * limitations under the License.
016     */
017    
018    package org.apache.commons.math3.ode.nonstiff;
019    
020    import java.util.Arrays;
021    
022    import org.apache.commons.math3.exception.DimensionMismatchException;
023    import org.apache.commons.math3.exception.MaxCountExceededException;
024    import org.apache.commons.math3.exception.NoBracketingException;
025    import org.apache.commons.math3.exception.NumberIsTooSmallException;
026    import org.apache.commons.math3.linear.Array2DRowRealMatrix;
027    import org.apache.commons.math3.linear.RealMatrixPreservingVisitor;
028    import org.apache.commons.math3.ode.ExpandableStatefulODE;
029    import org.apache.commons.math3.ode.sampling.NordsieckStepInterpolator;
030    import org.apache.commons.math3.util.FastMath;
031    
032    
033    /**
034     * This class implements implicit Adams-Moulton integrators for Ordinary
035     * Differential Equations.
036     *
037     * <p>Adams-Moulton methods (in fact due to Adams alone) are implicit
038     * multistep ODE solvers. This implementation is a variation of the classical
039     * one: it uses adaptive stepsize to implement error control, whereas
040     * classical implementations are fixed step size. The value of state vector
041     * at step n+1 is a simple combination of the value at step n and of the
042     * derivatives at steps n+1, n, n-1 ... Since y'<sub>n+1</sub> is needed to
043     * compute y<sub>n+1</sub>, another method must be used to compute a first
044     * estimate of y<sub>n+1</sub>, then compute y'<sub>n+1</sub>, then compute
045     * a final estimate of y<sub>n+1</sub> using the following formulas. Depending
046     * on the number k of previous steps one wants to use for computing the next
047     * value, different formulas are available for the final estimate:</p>
048     * <ul>
049     *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + h y'<sub>n+1</sub></li>
050     *   <li>k = 2: y<sub>n+1</sub> = y<sub>n</sub> + h (y'<sub>n+1</sub>+y'<sub>n</sub>)/2</li>
051     *   <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>
052     *   <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>
053     *   <li>...</li>
054     * </ul>
055     *
056     * <p>A k-steps Adams-Moulton method is of order k+1.</p>
057     *
058     * <h3>Implementation details</h3>
059     *
060     * <p>We define scaled derivatives s<sub>i</sub>(n) at step n as:
061     * <pre>
062     * s<sub>1</sub>(n) = h y'<sub>n</sub> for first derivative
063     * s<sub>2</sub>(n) = h<sup>2</sup>/2 y''<sub>n</sub> for second derivative
064     * s<sub>3</sub>(n) = h<sup>3</sup>/6 y'''<sub>n</sub> for third derivative
065     * ...
066     * s<sub>k</sub>(n) = h<sup>k</sup>/k! y<sup>(k)</sup><sub>n</sub> for k<sup>th</sup> derivative
067     * </pre></p>
068     *
069     * <p>The definitions above use the classical representation with several previous first
070     * derivatives. Lets define
071     * <pre>
072     *   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>
073     * </pre>
074     * (we omit the k index in the notation for clarity). With these definitions,
075     * Adams-Moulton methods can be written:
076     * <ul>
077     *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n+1)</li>
078     *   <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>
079     *   <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>
080     *   <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>
081     *   <li>...</li>
082     * </ul></p>
083     *
084     * <p>Instead of using the classical representation with first derivatives only (y<sub>n</sub>,
085     * s<sub>1</sub>(n+1) and q<sub>n+1</sub>), our implementation uses the Nordsieck vector with
086     * higher degrees scaled derivatives all taken at the same step (y<sub>n</sub>, s<sub>1</sub>(n)
087     * and r<sub>n</sub>) where r<sub>n</sub> is defined as:
088     * <pre>
089     * r<sub>n</sub> = [ s<sub>2</sub>(n), s<sub>3</sub>(n) ... s<sub>k</sub>(n) ]<sup>T</sup>
090     * </pre>
091     * (here again we omit the k index in the notation for clarity)
092     * </p>
093     *
094     * <p>Taylor series formulas show that for any index offset i, s<sub>1</sub>(n-i) can be
095     * computed from s<sub>1</sub>(n), s<sub>2</sub>(n) ... s<sub>k</sub>(n), the formula being exact
096     * for degree k polynomials.
097     * <pre>
098     * 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)
099     * </pre>
100     * The previous formula can be used with several values for i to compute the transform between
101     * classical representation and Nordsieck vector. The transform between r<sub>n</sub>
102     * and q<sub>n</sub> resulting from the Taylor series formulas above is:
103     * <pre>
104     * q<sub>n</sub> = s<sub>1</sub>(n) u + P r<sub>n</sub>
105     * </pre>
106     * where u is the [ 1 1 ... 1 ]<sup>T</sup> vector and P is the (k-1)&times;(k-1) matrix built
107     * with the j (-i)<sup>j-1</sup> terms:
108     * <pre>
109     *        [  -2   3   -4    5  ... ]
110     *        [  -4  12  -32   80  ... ]
111     *   P =  [  -6  27 -108  405  ... ]
112     *        [  -8  48 -256 1280  ... ]
113     *        [          ...           ]
114     * </pre></p>
115     *
116     * <p>Using the Nordsieck vector has several advantages:
117     * <ul>
118     *   <li>it greatly simplifies step interpolation as the interpolator mainly applies
119     *   Taylor series formulas,</li>
120     *   <li>it simplifies step changes that occur when discrete events that truncate
121     *   the step are triggered,</li>
122     *   <li>it allows to extend the methods in order to support adaptive stepsize.</li>
123     * </ul></p>
124     *
125     * <p>The predicted Nordsieck vector at step n+1 is computed from the Nordsieck vector at step
126     * n as follows:
127     * <ul>
128     *   <li>Y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li>
129     *   <li>S<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, Y<sub>n+1</sub>)</li>
130     *   <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>
131     * </ul>
132     * where A is a rows shifting matrix (the lower left part is an identity matrix):
133     * <pre>
134     *        [ 0 0   ...  0 0 | 0 ]
135     *        [ ---------------+---]
136     *        [ 1 0   ...  0 0 | 0 ]
137     *    A = [ 0 1   ...  0 0 | 0 ]
138     *        [       ...      | 0 ]
139     *        [ 0 0   ...  1 0 | 0 ]
140     *        [ 0 0   ...  0 1 | 0 ]
141     * </pre>
142     * From this predicted vector, the corrected vector is computed as follows:
143     * <ul>
144     *   <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>
145     *   <li>s<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, y<sub>n+1</sub>)</li>
146     *   <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>
147     * </ul>
148     * where the upper case Y<sub>n+1</sub>, S<sub>1</sub>(n+1) and R<sub>n+1</sub> represent the
149     * predicted states whereas the lower case y<sub>n+1</sub>, s<sub>n+1</sub> and r<sub>n+1</sub>
150     * represent the corrected states.</p>
151     *
152     * <p>The P<sup>-1</sup>u vector and the P<sup>-1</sup> A P matrix do not depend on the state,
153     * they only depend on k and therefore are precomputed once for all.</p>
154     *
155     * @version $Id: AdamsMoultonIntegrator.java 1416643 2012-12-03 19:37:14Z tn $
156     * @since 2.0
157     */
158    public class AdamsMoultonIntegrator extends AdamsIntegrator {
159    
160        /** Integrator method name. */
161        private static final String METHOD_NAME = "Adams-Moulton";
162    
163        /**
164         * Build an Adams-Moulton integrator with the given order and error control parameters.
165         * @param nSteps number of steps of the method excluding the one being computed
166         * @param minStep minimal step (sign is irrelevant, regardless of
167         * integration direction, forward or backward), the last step can
168         * be smaller than this
169         * @param maxStep maximal step (sign is irrelevant, regardless of
170         * integration direction, forward or backward), the last step can
171         * be smaller than this
172         * @param scalAbsoluteTolerance allowed absolute error
173         * @param scalRelativeTolerance allowed relative error
174         * @exception NumberIsTooSmallException if order is 1 or less
175         */
176        public AdamsMoultonIntegrator(final int nSteps,
177                                      final double minStep, final double maxStep,
178                                      final double scalAbsoluteTolerance,
179                                      final double scalRelativeTolerance)
180            throws NumberIsTooSmallException {
181            super(METHOD_NAME, nSteps, nSteps + 1, minStep, maxStep,
182                  scalAbsoluteTolerance, scalRelativeTolerance);
183        }
184    
185        /**
186         * Build an Adams-Moulton integrator with the given order and error control parameters.
187         * @param nSteps number of steps of the method excluding the one being computed
188         * @param minStep minimal step (sign is irrelevant, regardless of
189         * integration direction, forward or backward), the last step can
190         * be smaller than this
191         * @param maxStep maximal step (sign is irrelevant, regardless of
192         * integration direction, forward or backward), the last step can
193         * be smaller than this
194         * @param vecAbsoluteTolerance allowed absolute error
195         * @param vecRelativeTolerance allowed relative error
196         * @exception IllegalArgumentException if order is 1 or less
197         */
198        public AdamsMoultonIntegrator(final int nSteps,
199                                      final double minStep, final double maxStep,
200                                      final double[] vecAbsoluteTolerance,
201                                      final double[] vecRelativeTolerance)
202            throws IllegalArgumentException {
203            super(METHOD_NAME, nSteps, nSteps + 1, minStep, maxStep,
204                  vecAbsoluteTolerance, vecRelativeTolerance);
205        }
206    
207    
208        /** {@inheritDoc} */
209        @Override
210        public void integrate(final ExpandableStatefulODE equations,final double t)
211            throws NumberIsTooSmallException, DimensionMismatchException,
212                   MaxCountExceededException, NoBracketingException {
213    
214            sanityChecks(equations, t);
215            setEquations(equations);
216            final boolean forward = t > equations.getTime();
217    
218            // initialize working arrays
219            final double[] y0   = equations.getCompleteState();
220            final double[] y    = y0.clone();
221            final double[] yDot = new double[y.length];
222            final double[] yTmp = new double[y.length];
223            final double[] predictedScaled = new double[y.length];
224            Array2DRowRealMatrix nordsieckTmp = null;
225    
226            // set up two interpolators sharing the integrator arrays
227            final NordsieckStepInterpolator interpolator = new NordsieckStepInterpolator();
228            interpolator.reinitialize(y, forward,
229                                      equations.getPrimaryMapper(), equations.getSecondaryMappers());
230    
231            // set up integration control objects
232            initIntegration(equations.getTime(), y0, t);
233    
234            // compute the initial Nordsieck vector using the configured starter integrator
235            start(equations.getTime(), y, t);
236            interpolator.reinitialize(stepStart, stepSize, scaled, nordsieck);
237            interpolator.storeTime(stepStart);
238    
239            double hNew = stepSize;
240            interpolator.rescale(hNew);
241    
242            isLastStep = false;
243            do {
244    
245                double error = 10;
246                while (error >= 1.0) {
247    
248                    stepSize = hNew;
249    
250                    // predict a first estimate of the state at step end (P in the PECE sequence)
251                    final double stepEnd = stepStart + stepSize;
252                    interpolator.setInterpolatedTime(stepEnd);
253                    System.arraycopy(interpolator.getInterpolatedState(), 0, yTmp, 0, y0.length);
254    
255                    // evaluate a first estimate of the derivative (first E in the PECE sequence)
256                    computeDerivatives(stepEnd, yTmp, yDot);
257    
258                    // update Nordsieck vector
259                    for (int j = 0; j < y0.length; ++j) {
260                        predictedScaled[j] = stepSize * yDot[j];
261                    }
262                    nordsieckTmp = updateHighOrderDerivativesPhase1(nordsieck);
263                    updateHighOrderDerivativesPhase2(scaled, predictedScaled, nordsieckTmp);
264    
265                    // apply correction (C in the PECE sequence)
266                    error = nordsieckTmp.walkInOptimizedOrder(new Corrector(y, predictedScaled, yTmp));
267    
268                    if (error >= 1.0) {
269                        // reject the step and attempt to reduce error by stepsize control
270                        final double factor = computeStepGrowShrinkFactor(error);
271                        hNew = filterStep(stepSize * factor, forward, false);
272                        interpolator.rescale(hNew);
273                    }
274                }
275    
276                // evaluate a final estimate of the derivative (second E in the PECE sequence)
277                final double stepEnd = stepStart + stepSize;
278                computeDerivatives(stepEnd, yTmp, yDot);
279    
280                // update Nordsieck vector
281                final double[] correctedScaled = new double[y0.length];
282                for (int j = 0; j < y0.length; ++j) {
283                    correctedScaled[j] = stepSize * yDot[j];
284                }
285                updateHighOrderDerivativesPhase2(predictedScaled, correctedScaled, nordsieckTmp);
286    
287                // discrete events handling
288                System.arraycopy(yTmp, 0, y, 0, y.length);
289                interpolator.reinitialize(stepEnd, stepSize, correctedScaled, nordsieckTmp);
290                interpolator.storeTime(stepStart);
291                interpolator.shift();
292                interpolator.storeTime(stepEnd);
293                stepStart = acceptStep(interpolator, y, yDot, t);
294                scaled    = correctedScaled;
295                nordsieck = nordsieckTmp;
296    
297                if (!isLastStep) {
298    
299                    // prepare next step
300                    interpolator.storeTime(stepStart);
301    
302                    if (resetOccurred) {
303                        // some events handler has triggered changes that
304                        // invalidate the derivatives, we need to restart from scratch
305                        start(stepStart, y, t);
306                        interpolator.reinitialize(stepStart, stepSize, scaled, nordsieck);
307    
308                    }
309    
310                    // stepsize control for next step
311                    final double  factor     = computeStepGrowShrinkFactor(error);
312                    final double  scaledH    = stepSize * factor;
313                    final double  nextT      = stepStart + scaledH;
314                    final boolean nextIsLast = forward ? (nextT >= t) : (nextT <= t);
315                    hNew = filterStep(scaledH, forward, nextIsLast);
316    
317                    final double  filteredNextT      = stepStart + hNew;
318                    final boolean filteredNextIsLast = forward ? (filteredNextT >= t) : (filteredNextT <= t);
319                    if (filteredNextIsLast) {
320                        hNew = t - stepStart;
321                    }
322    
323                    interpolator.rescale(hNew);
324                }
325    
326            } while (!isLastStep);
327    
328            // dispatch results
329            equations.setTime(stepStart);
330            equations.setCompleteState(y);
331    
332            resetInternalState();
333    
334        }
335    
336        /** Corrector for current state in Adams-Moulton method.
337         * <p>
338         * This visitor implements the Taylor series formula:
339         * <pre>
340         * 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>
341         * </pre>
342         * </p>
343         */
344        private class Corrector implements RealMatrixPreservingVisitor {
345    
346            /** Previous state. */
347            private final double[] previous;
348    
349            /** Current scaled first derivative. */
350            private final double[] scaled;
351    
352            /** Current state before correction. */
353            private final double[] before;
354    
355            /** Current state after correction. */
356            private final double[] after;
357    
358            /** Simple constructor.
359             * @param previous previous state
360             * @param scaled current scaled first derivative
361             * @param state state to correct (will be overwritten after visit)
362             */
363            public Corrector(final double[] previous, final double[] scaled, final double[] state) {
364                this.previous = previous;
365                this.scaled   = scaled;
366                this.after    = state;
367                this.before   = state.clone();
368            }
369    
370            /** {@inheritDoc} */
371            public void start(int rows, int columns,
372                              int startRow, int endRow, int startColumn, int endColumn) {
373                Arrays.fill(after, 0.0);
374            }
375    
376            /** {@inheritDoc} */
377            public void visit(int row, int column, double value) {
378                if ((row & 0x1) == 0) {
379                    after[column] -= value;
380                } else {
381                    after[column] += value;
382                }
383            }
384    
385            /**
386             * End visiting the Nordsieck vector.
387             * <p>The correction is used to control stepsize. So its amplitude is
388             * considered to be an error, which must be normalized according to
389             * error control settings. If the normalized value is greater than 1,
390             * the correction was too large and the step must be rejected.</p>
391             * @return the normalized correction, if greater than 1, the step
392             * must be rejected
393             */
394            public double end() {
395    
396                double error = 0;
397                for (int i = 0; i < after.length; ++i) {
398                    after[i] += previous[i] + scaled[i];
399                    if (i < mainSetDimension) {
400                        final double yScale = FastMath.max(FastMath.abs(previous[i]), FastMath.abs(after[i]));
401                        final double tol = (vecAbsoluteTolerance == null) ?
402                                           (scalAbsoluteTolerance + scalRelativeTolerance * yScale) :
403                                           (vecAbsoluteTolerance[i] + vecRelativeTolerance[i] * yScale);
404                        final double ratio  = (after[i] - before[i]) / tol;
405                        error += ratio * ratio;
406                    }
407                }
408    
409                return FastMath.sqrt(error / mainSetDimension);
410    
411            }
412        }
413    
414    }