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 * https://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 /*
19 * This is not the original file distributed by the Apache Software Foundation
20 * It has been modified by the Hipparchus project
21 */
22
23 package org.hipparchus.ode.nonstiff;
24
25 import org.hipparchus.Field;
26 import org.hipparchus.CalculusFieldElement;
27 import org.hipparchus.exception.MathIllegalArgumentException;
28 import org.hipparchus.linear.Array2DRowFieldMatrix;
29 import org.hipparchus.linear.FieldMatrix;
30 import org.hipparchus.ode.FieldEquationsMapper;
31 import org.hipparchus.ode.FieldODEStateAndDerivative;
32 import org.hipparchus.util.FastMath;
33
34
35 /**
36 * This class implements explicit Adams-Bashforth integrators for Ordinary
37 * Differential Equations.
38 *
39 * <p>Adams-Bashforth methods (in fact due to Adams alone) are explicit
40 * multistep ODE solvers. This implementation is a variation of the classical
41 * one: it uses adaptive stepsize to implement error control, whereas
42 * classical implementations are fixed step size. The value of state vector
43 * at step n+1 is a simple combination of the value at step n and of the
44 * derivatives at steps n, n-1, n-2 ... Depending on the number k of previous
45 * steps one wants to use for computing the next value, different formulas
46 * are available:</p>
47 * <ul>
48 * <li>k = 1: y<sub>n+1</sub> = y<sub>n</sub> + h y'<sub>n</sub></li>
49 * <li>k = 2: y<sub>n+1</sub> = y<sub>n</sub> + h (3y'<sub>n</sub>-y'<sub>n-1</sub>)/2</li>
50 * <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>
51 * <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>
52 * <li>...</li>
53 * </ul>
54 *
55 * <p>A k-steps Adams-Bashforth method is of order k.</p>
56 *
57 * <p> There must be sufficient time for the {@link #setStarterIntegrator(org.hipparchus.ode.FieldODEIntegrator)
58 * starter integrator} to take several steps between the the last reset event, and the end
59 * of integration, otherwise an exception may be thrown during integration. The user can
60 * adjust the end date of integration, or the step size of the starter integrator to
61 * ensure a sufficient number of steps can be completed before the end of integration.
62 * </p>
63 *
64 * <p><strong>Implementation details</strong></p>
65 *
66 * <p>We define scaled derivatives s<sub>i</sub>(n) at step n as:
67 * \[
68 * \left\{\begin{align}
69 * s_1(n) &= h y'_n \text{ for first derivative}\\
70 * s_2(n) &= \frac{h^2}{2} y_n'' \text{ for second derivative}\\
71 * s_3(n) &= \frac{h^3}{6} y_n''' \text{ for third derivative}\\
72 * &\cdots\\
73 * s_k(n) &= \frac{h^k}{k!} y_n^{(k)} \text{ for } k^\mathrm{th} \text{ derivative}
74 * \end{align}\right.
75 * \]</p>
76 *
77 * <p>The definitions above use the classical representation with several previous first
78 * derivatives. Lets define
79 * \[
80 * q_n = [ s_1(n-1) s_1(n-2) \ldots s_1(n-(k-1)) ]^T
81 * \]
82 * (we omit the k index in the notation for clarity). With these definitions,
83 * Adams-Bashforth methods can be written:
84 * \[
85 * \left\{\begin{align}
86 * k = 1: & y_{n+1} = y_n + s_1(n) \\
87 * k = 2: & y_{n+1} = y_n + \frac{3}{2} s_1(n) + [ \frac{-1}{2} ] q_n \\
88 * k = 3: & y_{n+1} = y_n + \frac{23}{12} s_1(n) + [ \frac{-16}{12} \frac{5}{12} ] q_n \\
89 * k = 4: & y_{n+1} = y_n + \frac{55}{24} s_1(n) + [ \frac{-59}{24} \frac{37}{24} \frac{-9}{24} ] q_n \\
90 * & \cdots
91 * \end{align}\right.
92 * \]
93 *
94 * <p>Instead of using the classical representation with first derivatives only (y<sub>n</sub>,
95 * s<sub>1</sub>(n) and q<sub>n</sub>), our implementation uses the Nordsieck vector with
96 * higher degrees scaled derivatives all taken at the same step (y<sub>n</sub>, s<sub>1</sub>(n)
97 * and r<sub>n</sub>) where r<sub>n</sub> is defined as:
98 * \[
99 * r_n = [ s_2(n), s_3(n) \ldots s_k(n) ]^T
100 * \]
101 * (here again we omit the k index in the notation for clarity)
102 * </p>
103 *
104 * <p>Taylor series formulas show that for any index offset i, s<sub>1</sub>(n-i) can be
105 * computed from s<sub>1</sub>(n), s<sub>2</sub>(n) ... s<sub>k</sub>(n), the formula being exact
106 * for degree k polynomials.
107 * \[
108 * s_1(n-i) = s_1(n) + \sum_{j\gt 0} (j+1) (-i)^j s_{j+1}(n)
109 * \]
110 * The previous formula can be used with several values for i to compute the transform between
111 * classical representation and Nordsieck vector. The transform between r<sub>n</sub>
112 * and q<sub>n</sub> resulting from the Taylor series formulas above is:
113 * \[
114 * q_n = s_1(n) u + P r_n
115 * \]
116 * where u is the [ 1 1 ... 1 ]<sup>T</sup> vector and P is the (k-1)×(k-1) matrix built
117 * with the \((j+1) (-i)^j\) terms with i being the row number starting from 1 and j being
118 * the column number starting from 1:
119 * \[
120 * P=\begin{bmatrix}
121 * -2 & 3 & -4 & 5 & \ldots \\
122 * -4 & 12 & -32 & 80 & \ldots \\
123 * -6 & 27 & -108 & 405 & \ldots \\
124 * -8 & 48 & -256 & 1280 & \ldots \\
125 * & & \ldots\\
126 * \end{bmatrix}
127 * \]
128 *
129 * <p>Using the Nordsieck vector has several advantages:</p>
130 * <ul>
131 * <li>it greatly simplifies step interpolation as the interpolator mainly applies
132 * Taylor series formulas,</li>
133 * <li>it simplifies step changes that occur when discrete events that truncate
134 * the step are triggered,</li>
135 * <li>it allows to extend the methods in order to support adaptive stepsize.</li>
136 * </ul>
137 *
138 * <p>The Nordsieck vector at step n+1 is computed from the Nordsieck vector at step n as follows:
139 * <ul>
140 * <li>y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li>
141 * <li>s<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, y<sub>n+1</sub>)</li>
142 * <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>
143 * </ul>
144 * <p>where A is a rows shifting matrix (the lower left part is an identity matrix):</p>
145 * <pre>
146 * [ 0 0 ... 0 0 | 0 ]
147 * [ ---------------+---]
148 * [ 1 0 ... 0 0 | 0 ]
149 * A = [ 0 1 ... 0 0 | 0 ]
150 * [ ... | 0 ]
151 * [ 0 0 ... 1 0 | 0 ]
152 * [ 0 0 ... 0 1 | 0 ]
153 * </pre>
154 *
155 * <p>The P<sup>-1</sup>u vector and the P<sup>-1</sup> A P matrix do not depend on the state,
156 * they only depend on k and therefore are precomputed once for all.</p>
157 *
158 * @param <T> the type of the field elements
159 */
160 public class AdamsBashforthFieldIntegrator<T extends CalculusFieldElement<T>> extends AdamsFieldIntegrator<T> {
161
162 /** Integrator method name. */
163 private static final String METHOD_NAME = "Adams-Bashforth";
164
165 /**
166 * Build an Adams-Bashforth integrator with the given order and step control parameters.
167 * @param field field to which the time and state vector elements belong
168 * @param nSteps number of steps of the method excluding the one being computed
169 * @param minStep minimal step (sign is irrelevant, regardless of
170 * integration direction, forward or backward), the last step can
171 * be smaller than this
172 * @param maxStep maximal step (sign is irrelevant, regardless of
173 * integration direction, forward or backward), the last step can
174 * be smaller than this
175 * @param scalAbsoluteTolerance allowed absolute error
176 * @param scalRelativeTolerance allowed relative error
177 * @exception MathIllegalArgumentException if order is 1 or less
178 */
179 public AdamsBashforthFieldIntegrator(final Field<T> field, final int nSteps,
180 final double minStep, final double maxStep,
181 final double scalAbsoluteTolerance,
182 final double scalRelativeTolerance)
183 throws MathIllegalArgumentException {
184 super(field, METHOD_NAME, nSteps, nSteps, minStep, maxStep,
185 scalAbsoluteTolerance, scalRelativeTolerance);
186 }
187
188 /**
189 * Build an Adams-Bashforth integrator with the given order and step control parameters.
190 * @param field field to which the time and state vector elements belong
191 * @param nSteps number of steps of the method excluding the one being computed
192 * @param minStep minimal step (sign is irrelevant, regardless of
193 * integration direction, forward or backward), the last step can
194 * be smaller than this
195 * @param maxStep maximal step (sign is irrelevant, regardless of
196 * integration direction, forward or backward), the last step can
197 * be smaller than this
198 * @param vecAbsoluteTolerance allowed absolute error
199 * @param vecRelativeTolerance allowed relative error
200 * @exception IllegalArgumentException if order is 1 or less
201 */
202 public AdamsBashforthFieldIntegrator(final Field<T> field, final int nSteps,
203 final double minStep, final double maxStep,
204 final double[] vecAbsoluteTolerance,
205 final double[] vecRelativeTolerance)
206 throws IllegalArgumentException {
207 super(field, METHOD_NAME, nSteps, nSteps, minStep, maxStep,
208 vecAbsoluteTolerance, vecRelativeTolerance);
209 }
210
211 /** {@inheritDoc} */
212 @Override
213 protected double errorEstimation(final T[] previousState, final T predictedTime,
214 final T[] predictedState, final T[] predictedScaled,
215 final FieldMatrix<T> predictedNordsieck) {
216
217 final StepsizeHelper helper = getStepSizeHelper();
218 double error = 0;
219 for (int i = 0; i < helper.getMainSetDimension(); ++i) {
220 final double tol = helper.getTolerance(i, FastMath.abs(predictedState[i].getReal()));
221
222 // apply Taylor formula from high order to low order,
223 // for the sake of numerical accuracy
224 double variation = 0;
225 int sign = predictedNordsieck.getRowDimension() % 2 == 0 ? -1 : 1;
226 for (int k = predictedNordsieck.getRowDimension() - 1; k >= 0; --k) {
227 variation += sign * predictedNordsieck.getEntry(k, i).getReal();
228 sign = -sign;
229 }
230 variation -= predictedScaled[i].getReal();
231
232 final double ratio = (predictedState[i].getReal() - previousState[i].getReal() + variation) / tol;
233 error += ratio * ratio;
234
235 }
236
237 return FastMath.sqrt(error / helper.getMainSetDimension());
238
239 }
240
241 /** {@inheritDoc} */
242 @Override
243 protected AdamsFieldStateInterpolator<T> finalizeStep(final T stepSize, final T[] predictedY,
244 final T[] predictedScaled, final Array2DRowFieldMatrix<T> predictedNordsieck,
245 final boolean isForward,
246 final FieldODEStateAndDerivative<T> globalPreviousState,
247 final FieldODEStateAndDerivative<T> globalCurrentState,
248 final FieldEquationsMapper<T> equationsMapper) {
249 return new AdamsFieldStateInterpolator<>(getStepSize(), globalCurrentState,
250 predictedScaled, predictedNordsieck, isForward,
251 getStepStart(), globalCurrentState, equationsMapper);
252 }
253
254 }