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 org.apache.commons.math3.Field;
21 import org.apache.commons.math3.RealFieldElement;
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.Array2DRowFieldMatrix;
27 import org.apache.commons.math3.linear.FieldMatrix;
28 import org.apache.commons.math3.ode.FieldExpandableODE;
29 import org.apache.commons.math3.ode.FieldODEState;
30 import org.apache.commons.math3.ode.FieldODEStateAndDerivative;
31 import org.apache.commons.math3.util.MathArrays;
32 
33 
34 /**
35  * This class implements explicit Adams-Bashforth integrators for Ordinary
36  * Differential Equations.
37  *
38  * <p>Adams-Bashforth methods (in fact due to Adams alone) are explicit
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, n-1, n-2 ... Depending on the number k of previous
44  * steps one wants to use for computing the next value, different formulas
45  * are available:</p>
46  * <ul>
47  *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + h y'<sub>n</sub></li>
48  *   <li>k = 2: y<sub>n+1</sub> = y<sub>n</sub> + h (3y'<sub>n</sub>-y'<sub>n-1</sub>)/2</li>
49  *   <li>k = 3: y<sub>n+1</sub> = y<sub>n</sub> + h (23y'<sub>n</sub>-16y'<sub>n-1</sub>+5y'<sub>n-2</sub>)/12</li>
50  *   <li>k = 4: y<sub>n+1</sub> = y<sub>n</sub> + h (55y'<sub>n</sub>-59y'<sub>n-1</sub>+37y'<sub>n-2</sub>-9y'<sub>n-3</sub>)/24</li>
51  *   <li>...</li>
52  * </ul>
53  *
54  * <p>A k-steps Adams-Bashforth method is of order k.</p>
55  *
56  * <h3>Implementation details</h3>
57  *
58  * <p>We define scaled derivatives s<sub>i</sub>(n) at step n as:
59  * <pre>
60  * s<sub>1</sub>(n) = h y'<sub>n</sub> for first derivative
61  * s<sub>2</sub>(n) = h<sup>2</sup>/2 y''<sub>n</sub> for second derivative
62  * s<sub>3</sub>(n) = h<sup>3</sup>/6 y'''<sub>n</sub> for third derivative
63  * ...
64  * s<sub>k</sub>(n) = h<sup>k</sup>/k! y<sup>(k)</sup><sub>n</sub> for k<sup>th</sup> derivative
65  * </pre></p>
66  *
67  * <p>The definitions above use the classical representation with several previous first
68  * derivatives. Lets define
69  * <pre>
70  *   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>
71  * </pre>
72  * (we omit the k index in the notation for clarity). With these definitions,
73  * Adams-Bashforth methods can be written:
74  * <ul>
75  *   <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n)</li>
76  *   <li>k = 2: y<sub>n+1</sub> = y<sub>n</sub> + 3/2 s<sub>1</sub>(n) + [ -1/2 ] q<sub>n</sub></li>
77  *   <li>k = 3: y<sub>n+1</sub> = y<sub>n</sub> + 23/12 s<sub>1</sub>(n) + [ -16/12 5/12 ] q<sub>n</sub></li>
78  *   <li>k = 4: y<sub>n+1</sub> = y<sub>n</sub> + 55/24 s<sub>1</sub>(n) + [ -59/24 37/24 -9/24 ] q<sub>n</sub></li>
79  *   <li>...</li>
80  * </ul></p>
81  *
82  * <p>Instead of using the classical representation with first derivatives only (y<sub>n</sub>,
83  * s<sub>1</sub>(n) and q<sub>n</sub>), our implementation uses the Nordsieck vector with
84  * higher degrees scaled derivatives all taken at the same step (y<sub>n</sub>, s<sub>1</sub>(n)
85  * and r<sub>n</sub>) where r<sub>n</sub> is defined as:
86  * <pre>
87  * r<sub>n</sub> = [ s<sub>2</sub>(n), s<sub>3</sub>(n) ... s<sub>k</sub>(n) ]<sup>T</sup>
88  * </pre>
89  * (here again we omit the k index in the notation for clarity)
90  * </p>
91  *
92  * <p>Taylor series formulas show that for any index offset i, s<sub>1</sub>(n-i) can be
93  * computed from s<sub>1</sub>(n), s<sub>2</sub>(n) ... s<sub>k</sub>(n), the formula being exact
94  * for degree k polynomials.
95  * <pre>
96  * 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)
97  * </pre>
98  * The previous formula can be used with several values for i to compute the transform between
99  * classical representation and Nordsieck vector. The transform between r<sub>n</sub>
100  * and q<sub>n</sub> resulting from the Taylor series formulas above is:
101  * <pre>
102  * q<sub>n</sub> = s<sub>1</sub>(n) u + P r<sub>n</sub>
103  * </pre>
104  * where u is the [ 1 1 ... 1 ]<sup>T</sup> vector and P is the (k-1)&times;(k-1) matrix built
105  * with the (j+1) (-i)<sup>j</sup> terms with i being the row number starting from 1 and j being
106  * the column number starting from 1:
107  * <pre>
108  *        [  -2   3   -4    5  ... ]
109  *        [  -4  12  -32   80  ... ]
110  *   P =  [  -6  27 -108  405  ... ]
111  *        [  -8  48 -256 1280  ... ]
112  *        [          ...           ]
113  * </pre></p>
114  *
115  * <p>Using the Nordsieck vector has several advantages:
116  * <ul>
117  *   <li>it greatly simplifies step interpolation as the interpolator mainly applies
118  *   Taylor series formulas,</li>
119  *   <li>it simplifies step changes that occur when discrete events that truncate
120  *   the step are triggered,</li>
121  *   <li>it allows to extend the methods in order to support adaptive stepsize.</li>
122  * </ul></p>
123  *
124  * <p>The Nordsieck vector at step n+1 is computed from the Nordsieck vector at step n as follows:
125  * <ul>
126  *   <li>y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li>
127  *   <li>s<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, y<sub>n+1</sub>)</li>
128  *   <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>
129  * </ul>
130  * where A is a rows shifting matrix (the lower left part is an identity matrix):
131  * <pre>
132  *        [ 0 0   ...  0 0 | 0 ]
133  *        [ ---------------+---]
134  *        [ 1 0   ...  0 0 | 0 ]
135  *    A = [ 0 1   ...  0 0 | 0 ]
136  *        [       ...      | 0 ]
137  *        [ 0 0   ...  1 0 | 0 ]
138  *        [ 0 0   ...  0 1 | 0 ]
139  * </pre></p>
140  *
141  * <p>The P<sup>-1</sup>u vector and the P<sup>-1</sup> A P matrix do not depend on the state,
142  * they only depend on k and therefore are precomputed once for all.</p>
143  *
144  * @param <T> the type of the field elements
145  * @since 3.6
146  */
147 public class AdamsBashforthFieldIntegrator<T extends RealFieldElement<T>> extends AdamsFieldIntegrator<T> {
148 
149     /** Integrator method name. */
150     private static final String METHOD_NAME = "Adams-Bashforth";
151 
152     /**
153      * Build an Adams-Bashforth integrator with the given order and step control parameters.
154      * @param field field to which the time and state vector elements belong
155      * @param nSteps number of steps of the method excluding the one being computed
156      * @param minStep minimal step (sign is irrelevant, regardless of
157      * integration direction, forward or backward), the last step can
158      * be smaller than this
159      * @param maxStep maximal step (sign is irrelevant, regardless of
160      * integration direction, forward or backward), the last step can
161      * be smaller than this
162      * @param scalAbsoluteTolerance allowed absolute error
163      * @param scalRelativeTolerance allowed relative error
164      * @exception NumberIsTooSmallException if order is 1 or less
165      */
AdamsBashforthFieldIntegrator(final Field<T> field, final int nSteps, final double minStep, final double maxStep, final double scalAbsoluteTolerance, final double scalRelativeTolerance)166     public AdamsBashforthFieldIntegrator(final Field<T> field, final int nSteps,
167                                          final double minStep, final double maxStep,
168                                          final double scalAbsoluteTolerance,
169                                          final double scalRelativeTolerance)
170         throws NumberIsTooSmallException {
171         super(field, METHOD_NAME, nSteps, nSteps, minStep, maxStep,
172               scalAbsoluteTolerance, scalRelativeTolerance);
173     }
174 
175     /**
176      * Build an Adams-Bashforth integrator with the given order and step control parameters.
177      * @param field field to which the time and state vector elements belong
178      * @param nSteps number of steps of the method excluding the one being computed
179      * @param minStep minimal step (sign is irrelevant, regardless of
180      * integration direction, forward or backward), the last step can
181      * be smaller than this
182      * @param maxStep maximal step (sign is irrelevant, regardless of
183      * integration direction, forward or backward), the last step can
184      * be smaller than this
185      * @param vecAbsoluteTolerance allowed absolute error
186      * @param vecRelativeTolerance allowed relative error
187      * @exception IllegalArgumentException if order is 1 or less
188      */
AdamsBashforthFieldIntegrator(final Field<T> field, final int nSteps, final double minStep, final double maxStep, final double[] vecAbsoluteTolerance, final double[] vecRelativeTolerance)189     public AdamsBashforthFieldIntegrator(final Field<T> field, final int nSteps,
190                                          final double minStep, final double maxStep,
191                                          final double[] vecAbsoluteTolerance,
192                                          final double[] vecRelativeTolerance)
193         throws IllegalArgumentException {
194         super(field, METHOD_NAME, nSteps, nSteps, minStep, maxStep,
195               vecAbsoluteTolerance, vecRelativeTolerance);
196     }
197 
198     /** Estimate error.
199      * <p>
200      * Error is estimated by interpolating back to previous state using
201      * the state Taylor expansion and comparing to real previous state.
202      * </p>
203      * @param previousState state vector at step start
204      * @param predictedState predicted state vector at step end
205      * @param predictedScaled predicted value of the scaled derivatives at step end
206      * @param predictedNordsieck predicted value of the Nordsieck vector at step end
207      * @return estimated normalized local discretization error
208      */
errorEstimation(final T[] previousState, final T[] predictedState, final T[] predictedScaled, final FieldMatrix<T> predictedNordsieck)209     private T errorEstimation(final T[] previousState,
210                               final T[] predictedState,
211                               final T[] predictedScaled,
212                               final FieldMatrix<T> predictedNordsieck) {
213 
214         T error = getField().getZero();
215         for (int i = 0; i < mainSetDimension; ++i) {
216             final T yScale = predictedState[i].abs();
217             final T tol = (vecAbsoluteTolerance == null) ?
218                           yScale.multiply(scalRelativeTolerance).add(scalAbsoluteTolerance) :
219                           yScale.multiply(vecRelativeTolerance[i]).add(vecAbsoluteTolerance[i]);
220 
221             // apply Taylor formula from high order to low order,
222             // for the sake of numerical accuracy
223             T variation = getField().getZero();
224             int sign = predictedNordsieck.getRowDimension() % 2 == 0 ? -1 : 1;
225             for (int k = predictedNordsieck.getRowDimension() - 1; k >= 0; --k) {
226                 variation = variation.add(predictedNordsieck.getEntry(k, i).multiply(sign));
227                 sign      = -sign;
228             }
229             variation = variation.subtract(predictedScaled[i]);
230 
231             final T ratio  = predictedState[i].subtract(previousState[i]).add(variation).divide(tol);
232             error = error.add(ratio.multiply(ratio));
233 
234         }
235 
236         return error.divide(mainSetDimension).sqrt();
237 
238     }
239 
240     /** {@inheritDoc} */
241     @Override
integrate(final FieldExpandableODE<T> equations, final FieldODEState<T> initialState, final T finalTime)242     public FieldODEStateAndDerivative<T> integrate(final FieldExpandableODE<T> equations,
243                                                    final FieldODEState<T> initialState,
244                                                    final T finalTime)
245         throws NumberIsTooSmallException, DimensionMismatchException,
246                MaxCountExceededException, NoBracketingException {
247 
248         sanityChecks(initialState, finalTime);
249         final T   t0 = initialState.getTime();
250         final T[] y  = equations.getMapper().mapState(initialState);
251         setStepStart(initIntegration(equations, t0, y, finalTime));
252         final boolean forward = finalTime.subtract(initialState.getTime()).getReal() > 0;
253 
254         // compute the initial Nordsieck vector using the configured starter integrator
255         start(equations, getStepStart(), finalTime);
256 
257         // reuse the step that was chosen by the starter integrator
258         FieldODEStateAndDerivative<T> stepStart = getStepStart();
259         FieldODEStateAndDerivative<T> stepEnd   =
260                         AdamsFieldStepInterpolator.taylor(stepStart,
261                                                           stepStart.getTime().add(getStepSize()),
262                                                           getStepSize(), scaled, nordsieck);
263 
264         // main integration loop
265         setIsLastStep(false);
266         do {
267 
268             T[] predictedY = null;
269             final T[] predictedScaled = MathArrays.buildArray(getField(), y.length);
270             Array2DRowFieldMatrix<T> predictedNordsieck = null;
271             T error = getField().getZero().add(10);
272             while (error.subtract(1.0).getReal() >= 0.0) {
273 
274                 // predict a first estimate of the state at step end
275                 predictedY = stepEnd.getState();
276 
277                 // evaluate the derivative
278                 final T[] yDot = computeDerivatives(stepEnd.getTime(), predictedY);
279 
280                 // predict Nordsieck vector at step end
281                 for (int j = 0; j < predictedScaled.length; ++j) {
282                     predictedScaled[j] = getStepSize().multiply(yDot[j]);
283                 }
284                 predictedNordsieck = updateHighOrderDerivativesPhase1(nordsieck);
285                 updateHighOrderDerivativesPhase2(scaled, predictedScaled, predictedNordsieck);
286 
287                 // evaluate error
288                 error = errorEstimation(y, predictedY, predictedScaled, predictedNordsieck);
289 
290                 if (error.subtract(1.0).getReal() >= 0.0) {
291                     // reject the step and attempt to reduce error by stepsize control
292                     final T factor = computeStepGrowShrinkFactor(error);
293                     rescale(filterStep(getStepSize().multiply(factor), forward, false));
294                     stepEnd = AdamsFieldStepInterpolator.taylor(getStepStart(),
295                                                                 getStepStart().getTime().add(getStepSize()),
296                                                                 getStepSize(),
297                                                                 scaled,
298                                                                 nordsieck);
299 
300                 }
301             }
302 
303             // discrete events handling
304             setStepStart(acceptStep(new AdamsFieldStepInterpolator<T>(getStepSize(), stepEnd,
305                                                                       predictedScaled, predictedNordsieck, forward,
306                                                                       getStepStart(), stepEnd,
307                                                                       equations.getMapper()),
308                                     finalTime));
309             scaled    = predictedScaled;
310             nordsieck = predictedNordsieck;
311 
312             if (!isLastStep()) {
313 
314                 System.arraycopy(predictedY, 0, y, 0, y.length);
315 
316                 if (resetOccurred()) {
317                     // some events handler has triggered changes that
318                     // invalidate the derivatives, we need to restart from scratch
319                     start(equations, getStepStart(), finalTime);
320                 }
321 
322                 // stepsize control for next step
323                 final T       factor     = computeStepGrowShrinkFactor(error);
324                 final T       scaledH    = getStepSize().multiply(factor);
325                 final T       nextT      = getStepStart().getTime().add(scaledH);
326                 final boolean nextIsLast = forward ?
327                                            nextT.subtract(finalTime).getReal() >= 0 :
328                                            nextT.subtract(finalTime).getReal() <= 0;
329                 T hNew = filterStep(scaledH, forward, nextIsLast);
330 
331                 final T       filteredNextT      = getStepStart().getTime().add(hNew);
332                 final boolean filteredNextIsLast = forward ?
333                                                    filteredNextT.subtract(finalTime).getReal() >= 0 :
334                                                    filteredNextT.subtract(finalTime).getReal() <= 0;
335                 if (filteredNextIsLast) {
336                     hNew = finalTime.subtract(getStepStart().getTime());
337                 }
338 
339                 rescale(hNew);
340                 stepEnd = AdamsFieldStepInterpolator.taylor(getStepStart(), getStepStart().getTime().add(getStepSize()),
341                                                             getStepSize(), scaled, nordsieck);
342 
343             }
344 
345         } while (!isLastStep());
346 
347         final FieldODEStateAndDerivative<T> finalState = getStepStart();
348         setStepStart(null);
349         setStepSize(null);
350         return finalState;
351 
352     }
353 
354 }
355