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.math4.legacy.ode.nonstiff; 019 020import java.util.Arrays; 021 022import org.apache.commons.math4.legacy.exception.DimensionMismatchException; 023import org.apache.commons.math4.legacy.exception.MaxCountExceededException; 024import org.apache.commons.math4.legacy.exception.NoBracketingException; 025import org.apache.commons.math4.legacy.exception.NumberIsTooSmallException; 026import org.apache.commons.math4.legacy.linear.Array2DRowRealMatrix; 027import org.apache.commons.math4.legacy.linear.RealMatrixPreservingVisitor; 028import org.apache.commons.math4.legacy.ode.EquationsMapper; 029import org.apache.commons.math4.legacy.ode.ExpandableStatefulODE; 030import org.apache.commons.math4.legacy.ode.sampling.NordsieckStepInterpolator; 031import org.apache.commons.math4.core.jdkmath.JdkMath; 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 * <p><b>Implementation details</b></p> 060 * 061 * <p>We define scaled derivatives s<sub>i</sub>(n) at step n as: 062 * <div style="white-space: pre"><code> 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 * </code></div> 069 * 070 * <p>The definitions above use the classical representation with several previous first 071 * derivatives. Lets define 072 * <div style="white-space: pre"><code> 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 * </code></div> 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> 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 * <div style="white-space: pre"><code> 090 * r<sub>n</sub> = [ s<sub>2</sub>(n), s<sub>3</sub>(n) ... s<sub>k</sub>(n) ]<sup>T</sup> 091 * </code></div> 092 * (here again we omit the k index in the notation for clarity) 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 * <div style="white-space: pre"><code> 098 * s<sub>1</sub>(n-i) = s<sub>1</sub>(n) + ∑<sub>j>0</sub> (j+1) (-i)<sup>j</sup> s<sub>j+1</sub>(n) 099 * </code></div> 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 * <div style="white-space: pre"><code> 104 * q<sub>n</sub> = s<sub>1</sub>(n) u + P r<sub>n</sub> 105 * </code></div> 106 * where u is the [ 1 1 ... 1 ]<sup>T</sup> vector and P is the (k-1)×(k-1) matrix built 107 * with the (j+1) (-i)<sup>j</sup> terms with i being the row number starting from 1 and j being 108 * the column number starting from 1: 109 * <pre> 110 * [ -2 3 -4 5 ... ] 111 * [ -4 12 -32 80 ... ] 112 * P = [ -6 27 -108 405 ... ] 113 * [ -8 48 -256 1280 ... ] 114 * [ ... ] 115 * </pre> 116 * 117 * <p>Using the Nordsieck vector has several advantages: 118 * <ul> 119 * <li>it greatly simplifies step interpolation as the interpolator mainly applies 120 * Taylor series formulas,</li> 121 * <li>it simplifies step changes that occur when discrete events that truncate 122 * the step are triggered,</li> 123 * <li>it allows to extend the methods in order to support adaptive stepsize.</li> 124 * </ul> 125 * 126 * <p>The predicted Nordsieck vector at step n+1 is computed from the Nordsieck vector at step 127 * n as follows: 128 * <ul> 129 * <li>Y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n) + u<sup>T</sup> r<sub>n</sub></li> 130 * <li>S<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, Y<sub>n+1</sub>)</li> 131 * <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> 132 * </ul> 133 * where A is a rows shifting matrix (the lower left part is an identity matrix): 134 * <pre> 135 * [ 0 0 ... 0 0 | 0 ] 136 * [ ---------------+---] 137 * [ 1 0 ... 0 0 | 0 ] 138 * A = [ 0 1 ... 0 0 | 0 ] 139 * [ ... | 0 ] 140 * [ 0 0 ... 1 0 | 0 ] 141 * [ 0 0 ... 0 1 | 0 ] 142 * </pre> 143 * From this predicted vector, the corrected vector is computed as follows: 144 * <ul> 145 * <li>y<sub>n+1</sub> = y<sub>n</sub> + S<sub>1</sub>(n+1) + [ -1 +1 -1 +1 ... ±1 ] r<sub>n+1</sub></li> 146 * <li>s<sub>1</sub>(n+1) = h f(t<sub>n+1</sub>, y<sub>n+1</sub>)</li> 147 * <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> 148 * </ul> 149 * where the upper case Y<sub>n+1</sub>, S<sub>1</sub>(n+1) and R<sub>n+1</sub> represent the 150 * predicted states whereas the lower case y<sub>n+1</sub>, s<sub>n+1</sub> and r<sub>n+1</sub> 151 * represent the corrected states. 152 * 153 * <p>The P<sup>-1</sup>u vector and the P<sup>-1</sup> A P matrix do not depend on the state, 154 * they only depend on k and therefore are precomputed once for all.</p> 155 * 156 * @since 2.0 157 */ 158public 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 /** {@inheritDoc} */ 208 @Override 209 public void integrate(final ExpandableStatefulODE equations,final double t) 210 throws NumberIsTooSmallException, DimensionMismatchException, 211 MaxCountExceededException, NoBracketingException { 212 213 sanityChecks(equations, t); 214 setEquations(equations); 215 final boolean forward = t > equations.getTime(); 216 217 // initialize working arrays 218 final double[] y0 = equations.getCompleteState(); 219 final double[] y = y0.clone(); 220 final double[] yDot = new double[y.length]; 221 final double[] yTmp = new double[y.length]; 222 final double[] predictedScaled = new double[y.length]; 223 Array2DRowRealMatrix nordsieckTmp = null; 224 225 // set up two interpolators sharing the integrator arrays 226 final NordsieckStepInterpolator interpolator = new NordsieckStepInterpolator(); 227 interpolator.reinitialize(y, forward, 228 equations.getPrimaryMapper(), equations.getSecondaryMappers()); 229 230 // set up integration control objects 231 initIntegration(equations.getTime(), y0, t); 232 233 // compute the initial Nordsieck vector using the configured starter integrator 234 start(equations.getTime(), y, t); 235 interpolator.reinitialize(stepStart, stepSize, scaled, nordsieck); 236 interpolator.storeTime(stepStart); 237 238 double hNew = stepSize; 239 interpolator.rescale(hNew); 240 241 isLastStep = false; 242 do { 243 244 double error = 10; 245 while (error >= 1.0) { 246 247 stepSize = hNew; 248 249 // predict a first estimate of the state at step end (P in the PECE sequence) 250 final double stepEnd = stepStart + stepSize; 251 interpolator.setInterpolatedTime(stepEnd); 252 final ExpandableStatefulODE expandable = getExpandable(); 253 final EquationsMapper primary = expandable.getPrimaryMapper(); 254 primary.insertEquationData(interpolator.getInterpolatedState(), yTmp); 255 int index = 0; 256 for (final EquationsMapper secondary : expandable.getSecondaryMappers()) { 257 secondary.insertEquationData(interpolator.getInterpolatedSecondaryState(index), yTmp); 258 ++index; 259 } 260 261 // evaluate a first estimate of the derivative (first E in the PECE sequence) 262 computeDerivatives(stepEnd, yTmp, yDot); 263 264 // update Nordsieck vector 265 for (int j = 0; j < y0.length; ++j) { 266 predictedScaled[j] = stepSize * yDot[j]; 267 } 268 nordsieckTmp = updateHighOrderDerivativesPhase1(nordsieck); 269 updateHighOrderDerivativesPhase2(scaled, predictedScaled, nordsieckTmp); 270 271 // apply correction (C in the PECE sequence) 272 error = nordsieckTmp.walkInOptimizedOrder(new Corrector(y, predictedScaled, yTmp)); 273 274 if (error >= 1.0) { 275 // reject the step and attempt to reduce error by stepsize control 276 final double factor = computeStepGrowShrinkFactor(error); 277 hNew = filterStep(stepSize * factor, forward, false); 278 interpolator.rescale(hNew); 279 } 280 } 281 282 // evaluate a final estimate of the derivative (second E in the PECE sequence) 283 final double stepEnd = stepStart + stepSize; 284 computeDerivatives(stepEnd, yTmp, yDot); 285 286 // update Nordsieck vector 287 final double[] correctedScaled = new double[y0.length]; 288 for (int j = 0; j < y0.length; ++j) { 289 correctedScaled[j] = stepSize * yDot[j]; 290 } 291 updateHighOrderDerivativesPhase2(predictedScaled, correctedScaled, nordsieckTmp); 292 293 // discrete events handling 294 System.arraycopy(yTmp, 0, y, 0, y.length); 295 interpolator.reinitialize(stepEnd, stepSize, correctedScaled, nordsieckTmp); 296 interpolator.storeTime(stepStart); 297 interpolator.shift(); 298 interpolator.storeTime(stepEnd); 299 stepStart = acceptStep(interpolator, y, yDot, t); 300 scaled = correctedScaled; 301 nordsieck = nordsieckTmp; 302 303 if (!isLastStep) { 304 305 // prepare next step 306 interpolator.storeTime(stepStart); 307 308 if (resetOccurred) { 309 // some events handler has triggered changes that 310 // invalidate the derivatives, we need to restart from scratch 311 start(stepStart, y, t); 312 interpolator.reinitialize(stepStart, stepSize, scaled, nordsieck); 313 } 314 315 // stepsize control for next step 316 final double factor = computeStepGrowShrinkFactor(error); 317 final double scaledH = stepSize * factor; 318 final double nextT = stepStart + scaledH; 319 final boolean nextIsLast = forward ? (nextT >= t) : (nextT <= t); 320 hNew = filterStep(scaledH, forward, nextIsLast); 321 322 final double filteredNextT = stepStart + hNew; 323 final boolean filteredNextIsLast = forward ? (filteredNextT >= t) : (filteredNextT <= t); 324 if (filteredNextIsLast) { 325 hNew = t - stepStart; 326 } 327 328 interpolator.rescale(hNew); 329 } 330 } while (!isLastStep); 331 332 // dispatch results 333 equations.setTime(stepStart); 334 equations.setCompleteState(y); 335 336 resetInternalState(); 337 } 338 339 /** Corrector for current state in Adams-Moulton method. 340 * <p> 341 * This visitor implements the Taylor series formula: 342 * <pre> 343 * Y<sub>n+1</sub> = y<sub>n</sub> + s<sub>1</sub>(n+1) + [ -1 +1 -1 +1 ... ±1 ] r<sub>n+1</sub> 344 * </pre> 345 * </p> 346 */ 347 private final class Corrector implements RealMatrixPreservingVisitor { 348 349 /** Previous state. */ 350 private final double[] previous; 351 352 /** Current scaled first derivative. */ 353 private final double[] scaled; 354 355 /** Current state before correction. */ 356 private final double[] before; 357 358 /** Current state after correction. */ 359 private final double[] after; 360 361 /** Simple constructor. 362 * @param previous previous state 363 * @param scaled current scaled first derivative 364 * @param state state to correct (will be overwritten after visit) 365 */ 366 Corrector(final double[] previous, final double[] scaled, final double[] state) { 367 this.previous = previous; 368 this.scaled = scaled; 369 this.after = state; 370 this.before = state.clone(); 371 } 372 373 /** {@inheritDoc} */ 374 @Override 375 public void start(int rows, int columns, 376 int startRow, int endRow, int startColumn, int endColumn) { 377 Arrays.fill(after, 0.0); 378 } 379 380 /** {@inheritDoc} */ 381 @Override 382 public void visit(int row, int column, double value) { 383 if ((row & 0x1) == 0) { 384 after[column] -= value; 385 } else { 386 after[column] += value; 387 } 388 } 389 390 /** 391 * End visiting the Nordsieck vector. 392 * <p>The correction is used to control stepsize. So its amplitude is 393 * considered to be an error, which must be normalized according to 394 * error control settings. If the normalized value is greater than 1, 395 * the correction was too large and the step must be rejected.</p> 396 * @return the normalized correction, if greater than 1, the step 397 * must be rejected 398 */ 399 @Override 400 public double end() { 401 402 double error = 0; 403 for (int i = 0; i < after.length; ++i) { 404 after[i] += previous[i] + scaled[i]; 405 if (i < mainSetDimension) { 406 final double yScale = JdkMath.max(JdkMath.abs(previous[i]), JdkMath.abs(after[i])); 407 final double tol = (vecAbsoluteTolerance == null) ? 408 (scalAbsoluteTolerance + scalRelativeTolerance * yScale) : 409 (vecAbsoluteTolerance[i] + vecRelativeTolerance[i] * yScale); 410 final double ratio = (after[i] - before[i]) / tol; // (corrected-predicted)/tol 411 error += ratio * ratio; 412 } 413 } 414 415 return JdkMath.sqrt(error / mainSetDimension); 416 } 417 } 418}