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