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
018package org.apache.commons.math3.ode.nonstiff;
019
020import java.util.Arrays;
021
022import org.apache.commons.math3.exception.DimensionMismatchException;
023import org.apache.commons.math3.exception.MaxCountExceededException;
024import org.apache.commons.math3.exception.NoBracketingException;
025import org.apache.commons.math3.exception.NumberIsTooSmallException;
026import org.apache.commons.math3.linear.Array2DRowRealMatrix;
027import org.apache.commons.math3.linear.RealMatrixPreservingVisitor;
028import org.apache.commons.math3.ode.EquationsMapper;
029import org.apache.commons.math3.ode.ExpandableStatefulODE;
030import org.apache.commons.math3.ode.sampling.NordsieckStepInterpolator;
031import org.apache.commons.math3.util.FastMath;
032
033
034/**
035 * This class implements implicit Adams-Moulton integrators for Ordinary
036 * Differential Equations.
037 *
038 * <p>Adams-Moulton methods (in fact due to Adams alone) are implicit
039 * multistep ODE solvers. This implementation is a variation of the classical
040 * one: it uses adaptive stepsize to implement error control, whereas
041 * classical implementations are fixed step size. The value of state vector
042 * at step n+1 is a simple combination of the value at step n and of the
043 * derivatives at steps n+1, n, n-1 ... Since y'<sub>n+1</sub> is needed to
044 * compute y<sub>n+1</sub>, another method must be used to compute a first
045 * estimate of y<sub>n+1</sub>, then compute y'<sub>n+1</sub>, then compute
046 * a final estimate of y<sub>n+1</sub> using the following formulas. Depending
047 * on the number k of previous steps one wants to use for computing the next
048 * value, different formulas are available for the final estimate:</p>
049 * <ul>
050 *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + h y'<sub>n+1</sub></li>
051 *   <li>k = 2: y<sub>n+1</sub> = y<sub>n</sub> + h (y'<sub>n+1</sub>+y'<sub>n</sub>)/2</li>
052 *   <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>
053 *   <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>
054 *   <li>...</li>
055 * </ul>
056 *
057 * <p>A k-steps Adams-Moulton method is of order k+1.</p>
058 *
059 * <h3>Implementation details</h3>
060 *
061 * <p>We define scaled derivatives s<sub>i</sub>(n) at step n as:
062 * <pre>
063 * s<sub>1</sub>(n) = h y'<sub>n</sub> for first derivative
064 * s<sub>2</sub>(n) = h<sup>2</sup>/2 y''<sub>n</sub> for second derivative
065 * s<sub>3</sub>(n) = h<sup>3</sup>/6 y'''<sub>n</sub> for third derivative
066 * ...
067 * s<sub>k</sub>(n) = h<sup>k</sup>/k! y<sup>(k)</sup><sub>n</sub> for k<sup>th</sup> derivative
068 * </pre></p>
069 *
070 * <p>The definitions above use the classical representation with several previous first
071 * derivatives. Lets define
072 * <pre>
073 *   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>
074 * </pre>
075 * (we omit the k index in the notation for clarity). With these definitions,
076 * Adams-Moulton methods can be written:
077 * <ul>
078 *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n+1)</li>
079 *   <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>
080 *   <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>
081 *   <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>
082 *   <li>...</li>
083 * </ul></p>
084 *
085 * <p>Instead of using the classical representation with first derivatives only (y<sub>n</sub>,
086 * s<sub>1</sub>(n+1) and q<sub>n+1</sub>), our implementation uses the Nordsieck vector with
087 * higher degrees scaled derivatives all taken at the same step (y<sub>n</sub>, s<sub>1</sub>(n)
088 * and r<sub>n</sub>) where r<sub>n</sub> is defined as:
089 * <pre>
090 * r<sub>n</sub> = [ s<sub>2</sub>(n), s<sub>3</sub>(n) ... s<sub>k</sub>(n) ]<sup>T</sup>
091 * </pre>
092 * (here again we omit the k index in the notation for clarity)
093 * </p>
094 *
095 * <p>Taylor series formulas show that for any index offset i, s<sub>1</sub>(n-i) can be
096 * computed from s<sub>1</sub>(n), s<sub>2</sub>(n) ... s<sub>k</sub>(n), the formula being exact
097 * for degree k polynomials.
098 * <pre>
099 * s<sub>1</sub>(n-i) = s<sub>1</sub>(n) + &sum;<sub>j&gt;0</sub> (j+1) (-i)<sup>j</sup> s<sub>j+1</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+1) (-i)<sup>j</sup> terms with i being the row number starting from 1 and j being
109 * the column number starting from 1:
110 * <pre>
111 *        [  -2   3   -4    5  ... ]
112 *        [  -4  12  -32   80  ... ]
113 *   P =  [  -6  27 -108  405  ... ]
114 *        [  -8  48 -256 1280  ... ]
115 *        [          ...           ]
116 * </pre></p>
117 *
118 * <p>Using the Nordsieck vector has several advantages:
119 * <ul>
120 *   <li>it greatly simplifies step interpolation as the interpolator mainly applies
121 *   Taylor series formulas,</li>
122 *   <li>it simplifies step changes that occur when discrete events that truncate
123 *   the step are triggered,</li>
124 *   <li>it allows to extend the methods in order to support adaptive stepsize.</li>
125 * </ul></p>
126 *
127 * <p>The predicted Nordsieck vector at step n+1 is computed from the Nordsieck vector at step
128 * n as follows:
129 * <ul>
130 *   <li>Y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li>
131 *   <li>S<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, Y<sub>n+1</sub>)</li>
132 *   <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>
133 * </ul>
134 * where A is a rows shifting matrix (the lower left part is an identity matrix):
135 * <pre>
136 *        [ 0 0   ...  0 0 | 0 ]
137 *        [ ---------------+---]
138 *        [ 1 0   ...  0 0 | 0 ]
139 *    A = [ 0 1   ...  0 0 | 0 ]
140 *        [       ...      | 0 ]
141 *        [ 0 0   ...  1 0 | 0 ]
142 *        [ 0 0   ...  0 1 | 0 ]
143 * </pre>
144 * From this predicted vector, the corrected vector is computed as follows:
145 * <ul>
146 *   <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>
147 *   <li>s<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, y<sub>n+1</sub>)</li>
148 *   <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>
149 * </ul>
150 * where the upper case Y<sub>n+1</sub>, S<sub>1</sub>(n+1) and R<sub>n+1</sub> represent the
151 * predicted states whereas the lower case y<sub>n+1</sub>, s<sub>n+1</sub> and r<sub>n+1</sub>
152 * represent the corrected states.</p>
153 *
154 * <p>The P<sup>-1</sup>u vector and the P<sup>-1</sup> A P matrix do not depend on the state,
155 * they only depend on k and therefore are precomputed once for all.</p>
156 *
157 * @since 2.0
158 */
159public 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    /** {@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                final ExpandableStatefulODE expandable = getExpandable();
254                final EquationsMapper primary = expandable.getPrimaryMapper();
255                primary.insertEquationData(interpolator.getInterpolatedState(), yTmp);
256                int index = 0;
257                for (final EquationsMapper secondary : expandable.getSecondaryMappers()) {
258                    secondary.insertEquationData(interpolator.getInterpolatedSecondaryState(index), yTmp);
259                    ++index;
260                }
261
262                // evaluate a first estimate of the derivative (first E in the PECE sequence)
263                computeDerivatives(stepEnd, yTmp, yDot);
264
265                // update Nordsieck vector
266                for (int j = 0; j < y0.length; ++j) {
267                    predictedScaled[j] = stepSize * yDot[j];
268                }
269                nordsieckTmp = updateHighOrderDerivativesPhase1(nordsieck);
270                updateHighOrderDerivativesPhase2(scaled, predictedScaled, nordsieckTmp);
271
272                // apply correction (C in the PECE sequence)
273                error = nordsieckTmp.walkInOptimizedOrder(new Corrector(y, predictedScaled, yTmp));
274
275                if (error >= 1.0) {
276                    // reject the step and attempt to reduce error by stepsize control
277                    final double factor = computeStepGrowShrinkFactor(error);
278                    hNew = filterStep(stepSize * factor, forward, false);
279                    interpolator.rescale(hNew);
280                }
281            }
282
283            // evaluate a final estimate of the derivative (second E in the PECE sequence)
284            final double stepEnd = stepStart + stepSize;
285            computeDerivatives(stepEnd, yTmp, yDot);
286
287            // update Nordsieck vector
288            final double[] correctedScaled = new double[y0.length];
289            for (int j = 0; j < y0.length; ++j) {
290                correctedScaled[j] = stepSize * yDot[j];
291            }
292            updateHighOrderDerivativesPhase2(predictedScaled, correctedScaled, nordsieckTmp);
293
294            // discrete events handling
295            System.arraycopy(yTmp, 0, y, 0, y.length);
296            interpolator.reinitialize(stepEnd, stepSize, correctedScaled, nordsieckTmp);
297            interpolator.storeTime(stepStart);
298            interpolator.shift();
299            interpolator.storeTime(stepEnd);
300            stepStart = acceptStep(interpolator, y, yDot, t);
301            scaled    = correctedScaled;
302            nordsieck = nordsieckTmp;
303
304            if (!isLastStep) {
305
306                // prepare next step
307                interpolator.storeTime(stepStart);
308
309                if (resetOccurred) {
310                    // some events handler has triggered changes that
311                    // invalidate the derivatives, we need to restart from scratch
312                    start(stepStart, y, t);
313                    interpolator.reinitialize(stepStart, stepSize, scaled, nordsieck);
314
315                }
316
317                // stepsize control for next step
318                final double  factor     = computeStepGrowShrinkFactor(error);
319                final double  scaledH    = stepSize * factor;
320                final double  nextT      = stepStart + scaledH;
321                final boolean nextIsLast = forward ? (nextT >= t) : (nextT <= t);
322                hNew = filterStep(scaledH, forward, nextIsLast);
323
324                final double  filteredNextT      = stepStart + hNew;
325                final boolean filteredNextIsLast = forward ? (filteredNextT >= t) : (filteredNextT <= t);
326                if (filteredNextIsLast) {
327                    hNew = t - stepStart;
328                }
329
330                interpolator.rescale(hNew);
331            }
332
333        } while (!isLastStep);
334
335        // dispatch results
336        equations.setTime(stepStart);
337        equations.setCompleteState(y);
338
339        resetInternalState();
340
341    }
342
343    /** Corrector for current state in Adams-Moulton method.
344     * <p>
345     * This visitor implements the Taylor series formula:
346     * <pre>
347     * 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>
348     * </pre>
349     * </p>
350     */
351    private class Corrector implements RealMatrixPreservingVisitor {
352
353        /** Previous state. */
354        private final double[] previous;
355
356        /** Current scaled first derivative. */
357        private final double[] scaled;
358
359        /** Current state before correction. */
360        private final double[] before;
361
362        /** Current state after correction. */
363        private final double[] after;
364
365        /** Simple constructor.
366         * @param previous previous state
367         * @param scaled current scaled first derivative
368         * @param state state to correct (will be overwritten after visit)
369         */
370        Corrector(final double[] previous, final double[] scaled, final double[] state) {
371            this.previous = previous;
372            this.scaled   = scaled;
373            this.after    = state;
374            this.before   = state.clone();
375        }
376
377        /** {@inheritDoc} */
378        public void start(int rows, int columns,
379                          int startRow, int endRow, int startColumn, int endColumn) {
380            Arrays.fill(after, 0.0);
381        }
382
383        /** {@inheritDoc} */
384        public void visit(int row, int column, double value) {
385            if ((row & 0x1) == 0) {
386                after[column] -= value;
387            } else {
388                after[column] += value;
389            }
390        }
391
392        /**
393         * End visiting the Nordsieck vector.
394         * <p>The correction is used to control stepsize. So its amplitude is
395         * considered to be an error, which must be normalized according to
396         * error control settings. If the normalized value is greater than 1,
397         * the correction was too large and the step must be rejected.</p>
398         * @return the normalized correction, if greater than 1, the step
399         * must be rejected
400         */
401        public double end() {
402
403            double error = 0;
404            for (int i = 0; i < after.length; ++i) {
405                after[i] += previous[i] + scaled[i];
406                if (i < mainSetDimension) {
407                    final double yScale = FastMath.max(FastMath.abs(previous[i]), FastMath.abs(after[i]));
408                    final double tol    = (vecAbsoluteTolerance == null) ?
409                                          (scalAbsoluteTolerance + scalRelativeTolerance * yScale) :
410                                          (vecAbsoluteTolerance[i] + vecRelativeTolerance[i] * yScale);
411                    final double ratio  = (after[i] - before[i]) / tol; // (corrected-predicted)/tol
412                    error += ratio * ratio;
413                }
414            }
415
416            return FastMath.sqrt(error / mainSetDimension);
417
418        }
419    }
420
421}