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
018 package org.apache.commons.math.ode.nonstiff;
019
020
021 /**
022 * This class implements the 5(4) Dormand-Prince integrator for Ordinary
023 * Differential Equations.
024
025 * <p>This integrator is an embedded Runge-Kutta integrator
026 * of order 5(4) used in local extrapolation mode (i.e. the solution
027 * is computed using the high order formula) with stepsize control
028 * (and automatic step initialization) and continuous output. This
029 * method uses 7 functions evaluations per step. However, since this
030 * is an <i>fsal</i>, the last evaluation of one step is the same as
031 * the first evaluation of the next step and hence can be avoided. So
032 * the cost is really 6 functions evaluations per step.</p>
033 *
034 * <p>This method has been published (whithout the continuous output
035 * that was added by Shampine in 1986) in the following article :
036 * <pre>
037 * A family of embedded Runge-Kutta formulae
038 * J. R. Dormand and P. J. Prince
039 * Journal of Computational and Applied Mathematics
040 * volume 6, no 1, 1980, pp. 19-26
041 * </pre></p>
042 *
043 * @version $Revision: 810196 $ $Date: 2009-09-01 15:47:46 -0400 (Tue, 01 Sep 2009) $
044 * @since 1.2
045 */
046
047 public class DormandPrince54Integrator extends EmbeddedRungeKuttaIntegrator {
048
049 /** Integrator method name. */
050 private static final String METHOD_NAME = "Dormand-Prince 5(4)";
051
052 /** Time steps Butcher array. */
053 private static final double[] STATIC_C = {
054 1.0/5.0, 3.0/10.0, 4.0/5.0, 8.0/9.0, 1.0, 1.0
055 };
056
057 /** Internal weights Butcher array. */
058 private static final double[][] STATIC_A = {
059 {1.0/5.0},
060 {3.0/40.0, 9.0/40.0},
061 {44.0/45.0, -56.0/15.0, 32.0/9.0},
062 {19372.0/6561.0, -25360.0/2187.0, 64448.0/6561.0, -212.0/729.0},
063 {9017.0/3168.0, -355.0/33.0, 46732.0/5247.0, 49.0/176.0, -5103.0/18656.0},
064 {35.0/384.0, 0.0, 500.0/1113.0, 125.0/192.0, -2187.0/6784.0, 11.0/84.0}
065 };
066
067 /** Propagation weights Butcher array. */
068 private static final double[] STATIC_B = {
069 35.0/384.0, 0.0, 500.0/1113.0, 125.0/192.0, -2187.0/6784.0, 11.0/84.0, 0.0
070 };
071
072 /** Error array, element 1. */
073 private static final double E1 = 71.0 / 57600.0;
074
075 // element 2 is zero, so it is neither stored nor used
076
077 /** Error array, element 3. */
078 private static final double E3 = -71.0 / 16695.0;
079
080 /** Error array, element 4. */
081 private static final double E4 = 71.0 / 1920.0;
082
083 /** Error array, element 5. */
084 private static final double E5 = -17253.0 / 339200.0;
085
086 /** Error array, element 6. */
087 private static final double E6 = 22.0 / 525.0;
088
089 /** Error array, element 7. */
090 private static final double E7 = -1.0 / 40.0;
091
092 /** Simple constructor.
093 * Build a fifth order Dormand-Prince integrator with the given step bounds
094 * @param minStep minimal step (must be positive even for backward
095 * integration), the last step can be smaller than this
096 * @param maxStep maximal step (must be positive even for backward
097 * integration)
098 * @param scalAbsoluteTolerance allowed absolute error
099 * @param scalRelativeTolerance allowed relative error
100 */
101 public DormandPrince54Integrator(final double minStep, final double maxStep,
102 final double scalAbsoluteTolerance,
103 final double scalRelativeTolerance) {
104 super(METHOD_NAME, true, STATIC_C, STATIC_A, STATIC_B, new DormandPrince54StepInterpolator(),
105 minStep, maxStep, scalAbsoluteTolerance, scalRelativeTolerance);
106 }
107
108 /** Simple constructor.
109 * Build a fifth order Dormand-Prince integrator with the given step bounds
110 * @param minStep minimal step (must be positive even for backward
111 * integration), the last step can be smaller than this
112 * @param maxStep maximal step (must be positive even for backward
113 * integration)
114 * @param vecAbsoluteTolerance allowed absolute error
115 * @param vecRelativeTolerance allowed relative error
116 */
117 public DormandPrince54Integrator(final double minStep, final double maxStep,
118 final double[] vecAbsoluteTolerance,
119 final double[] vecRelativeTolerance) {
120 super(METHOD_NAME, true, STATIC_C, STATIC_A, STATIC_B, new DormandPrince54StepInterpolator(),
121 minStep, maxStep, vecAbsoluteTolerance, vecRelativeTolerance);
122 }
123
124 /** {@inheritDoc} */
125 @Override
126 public int getOrder() {
127 return 5;
128 }
129
130 /** {@inheritDoc} */
131 @Override
132 protected double estimateError(final double[][] yDotK,
133 final double[] y0, final double[] y1,
134 final double h) {
135
136 double error = 0;
137
138 for (int j = 0; j < y0.length; ++j) {
139 final double errSum = E1 * yDotK[0][j] + E3 * yDotK[2][j] +
140 E4 * yDotK[3][j] + E5 * yDotK[4][j] +
141 E6 * yDotK[5][j] + E7 * yDotK[6][j];
142
143 final double yScale = Math.max(Math.abs(y0[j]), Math.abs(y1[j]));
144 final double tol = (vecAbsoluteTolerance == null) ?
145 (scalAbsoluteTolerance + scalRelativeTolerance * yScale) :
146 (vecAbsoluteTolerance[j] + vecRelativeTolerance[j] * yScale);
147 final double ratio = h * errSum / tol;
148 error += ratio * ratio;
149
150 }
151
152 return Math.sqrt(error / y0.length);
153
154 }
155
156 }