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 18package org.apache.commons.math.ode.nonstiff; 19 20import org.apache.commons.math.util.FastMath; 21 22 23/** 24 * This class implements the 5(4) Dormand-Prince integrator for Ordinary 25 * Differential Equations. 26 27 * <p>This integrator is an embedded Runge-Kutta integrator 28 * of order 5(4) used in local extrapolation mode (i.e. the solution 29 * is computed using the high order formula) with stepsize control 30 * (and automatic step initialization) and continuous output. This 31 * method uses 7 functions evaluations per step. However, since this 32 * is an <i>fsal</i>, the last evaluation of one step is the same as 33 * the first evaluation of the next step and hence can be avoided. So 34 * the cost is really 6 functions evaluations per step.</p> 35 * 36 * <p>This method has been published (whithout the continuous output 37 * that was added by Shampine in 1986) in the following article : 38 * <pre> 39 * A family of embedded Runge-Kutta formulae 40 * J. R. Dormand and P. J. Prince 41 * Journal of Computational and Applied Mathematics 42 * volume 6, no 1, 1980, pp. 19-26 43 * </pre></p> 44 * 45 * @version $Revision: 990655 $ $Date: 2010-08-29 23:49:40 +0200 (dim. 29 août 2010) $ 46 * @since 1.2 47 */ 48 49public class DormandPrince54Integrator extends EmbeddedRungeKuttaIntegrator { 50 51 /** Integrator method name. */ 52 private static final String METHOD_NAME = "Dormand-Prince 5(4)"; 53 54 /** Time steps Butcher array. */ 55 private static final double[] STATIC_C = { 56 1.0/5.0, 3.0/10.0, 4.0/5.0, 8.0/9.0, 1.0, 1.0 57 }; 58 59 /** Internal weights Butcher array. */ 60 private static final double[][] STATIC_A = { 61 {1.0/5.0}, 62 {3.0/40.0, 9.0/40.0}, 63 {44.0/45.0, -56.0/15.0, 32.0/9.0}, 64 {19372.0/6561.0, -25360.0/2187.0, 64448.0/6561.0, -212.0/729.0}, 65 {9017.0/3168.0, -355.0/33.0, 46732.0/5247.0, 49.0/176.0, -5103.0/18656.0}, 66 {35.0/384.0, 0.0, 500.0/1113.0, 125.0/192.0, -2187.0/6784.0, 11.0/84.0} 67 }; 68 69 /** Propagation weights Butcher array. */ 70 private static final double[] STATIC_B = { 71 35.0/384.0, 0.0, 500.0/1113.0, 125.0/192.0, -2187.0/6784.0, 11.0/84.0, 0.0 72 }; 73 74 /** Error array, element 1. */ 75 private static final double E1 = 71.0 / 57600.0; 76 77 // element 2 is zero, so it is neither stored nor used 78 79 /** Error array, element 3. */ 80 private static final double E3 = -71.0 / 16695.0; 81 82 /** Error array, element 4. */ 83 private static final double E4 = 71.0 / 1920.0; 84 85 /** Error array, element 5. */ 86 private static final double E5 = -17253.0 / 339200.0; 87 88 /** Error array, element 6. */ 89 private static final double E6 = 22.0 / 525.0; 90 91 /** Error array, element 7. */ 92 private static final double E7 = -1.0 / 40.0; 93 94 /** Simple constructor. 95 * Build a fifth order Dormand-Prince integrator with the given step bounds 96 * @param minStep minimal step (must be positive even for backward 97 * integration), the last step can be smaller than this 98 * @param maxStep maximal step (must be positive even for backward 99 * integration) 100 * @param scalAbsoluteTolerance allowed absolute error 101 * @param scalRelativeTolerance allowed relative error 102 */ 103 public DormandPrince54Integrator(final double minStep, final double maxStep, 104 final double scalAbsoluteTolerance, 105 final double scalRelativeTolerance) { 106 super(METHOD_NAME, true, STATIC_C, STATIC_A, STATIC_B, new DormandPrince54StepInterpolator(), 107 minStep, maxStep, scalAbsoluteTolerance, scalRelativeTolerance); 108 } 109 110 /** Simple constructor. 111 * Build a fifth order Dormand-Prince integrator with the given step bounds 112 * @param minStep minimal step (must be positive even for backward 113 * integration), the last step can be smaller than this 114 * @param maxStep maximal step (must be positive even for backward 115 * integration) 116 * @param vecAbsoluteTolerance allowed absolute error 117 * @param vecRelativeTolerance allowed relative error 118 */ 119 public DormandPrince54Integrator(final double minStep, final double maxStep, 120 final double[] vecAbsoluteTolerance, 121 final double[] vecRelativeTolerance) { 122 super(METHOD_NAME, true, STATIC_C, STATIC_A, STATIC_B, new DormandPrince54StepInterpolator(), 123 minStep, maxStep, vecAbsoluteTolerance, vecRelativeTolerance); 124 } 125 126 /** {@inheritDoc} */ 127 @Override 128 public int getOrder() { 129 return 5; 130 } 131 132 /** {@inheritDoc} */ 133 @Override 134 protected double estimateError(final double[][] yDotK, 135 final double[] y0, final double[] y1, 136 final double h) { 137 138 double error = 0; 139 140 for (int j = 0; j < mainSetDimension; ++j) { 141 final double errSum = E1 * yDotK[0][j] + E3 * yDotK[2][j] + 142 E4 * yDotK[3][j] + E5 * yDotK[4][j] + 143 E6 * yDotK[5][j] + E7 * yDotK[6][j]; 144 145 final double yScale = FastMath.max(FastMath.abs(y0[j]), FastMath.abs(y1[j])); 146 final double tol = (vecAbsoluteTolerance == null) ? 147 (scalAbsoluteTolerance + scalRelativeTolerance * yScale) : 148 (vecAbsoluteTolerance[j] + vecRelativeTolerance[j] * yScale); 149 final double ratio = h * errSum / tol; 150 error += ratio * ratio; 151 152 } 153 154 return FastMath.sqrt(error / mainSetDimension); 155 156 } 157 158} 159