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 */
017package org.apache.commons.math3.analysis.integration;
018
019import org.apache.commons.math3.analysis.UnivariateFunction;
020import org.apache.commons.math3.analysis.integration.gauss.GaussIntegratorFactory;
021import org.apache.commons.math3.analysis.integration.gauss.GaussIntegrator;
022import org.apache.commons.math3.exception.MathIllegalArgumentException;
023import org.apache.commons.math3.exception.MaxCountExceededException;
024import org.apache.commons.math3.exception.NotStrictlyPositiveException;
025import org.apache.commons.math3.exception.NumberIsTooSmallException;
026import org.apache.commons.math3.exception.TooManyEvaluationsException;
027import org.apache.commons.math3.exception.util.LocalizedFormats;
028import org.apache.commons.math3.util.FastMath;
029
030/**
031 * This algorithm divides the integration interval into equally-sized
032 * sub-interval and on each of them performs a
033 * <a href="http://mathworld.wolfram.com/Legendre-GaussQuadrature.html">
034 * Legendre-Gauss</a> quadrature.
035 * Because of its <em>non-adaptive</em> nature, this algorithm can
036 * converge to a wrong value for the integral (for example, if the
037 * function is significantly different from zero toward the ends of the
038 * integration interval).
039 * In particular, a change of variables aimed at estimating integrals
040 * over infinite intervals as proposed
041 * <a href="http://en.wikipedia.org/w/index.php?title=Numerical_integration#Integrals_over_infinite_intervals">
042 *  here</a> should be avoided when using this class.
043 *
044 * @since 3.1
045 */
046
047public class IterativeLegendreGaussIntegrator
048    extends BaseAbstractUnivariateIntegrator {
049    /** Factory that computes the points and weights. */
050    private static final GaussIntegratorFactory FACTORY
051        = new GaussIntegratorFactory();
052    /** Number of integration points (per interval). */
053    private final int numberOfPoints;
054
055    /**
056     * Builds an integrator with given accuracies and iterations counts.
057     *
058     * @param n Number of integration points.
059     * @param relativeAccuracy Relative accuracy of the result.
060     * @param absoluteAccuracy Absolute accuracy of the result.
061     * @param minimalIterationCount Minimum number of iterations.
062     * @param maximalIterationCount Maximum number of iterations.
063     * @throws NotStrictlyPositiveException if minimal number of iterations
064     * or number of points are not strictly positive.
065     * @throws NumberIsTooSmallException if maximal number of iterations
066     * is smaller than or equal to the minimal number of iterations.
067     */
068    public IterativeLegendreGaussIntegrator(final int n,
069                                            final double relativeAccuracy,
070                                            final double absoluteAccuracy,
071                                            final int minimalIterationCount,
072                                            final int maximalIterationCount)
073        throws NotStrictlyPositiveException, NumberIsTooSmallException {
074        super(relativeAccuracy, absoluteAccuracy, minimalIterationCount, maximalIterationCount);
075        if (n <= 0) {
076            throw new NotStrictlyPositiveException(LocalizedFormats.NUMBER_OF_POINTS, n);
077        }
078       numberOfPoints = n;
079    }
080
081    /**
082     * Builds an integrator with given accuracies.
083     *
084     * @param n Number of integration points.
085     * @param relativeAccuracy Relative accuracy of the result.
086     * @param absoluteAccuracy Absolute accuracy of the result.
087     * @throws NotStrictlyPositiveException if {@code n < 1}.
088     */
089    public IterativeLegendreGaussIntegrator(final int n,
090                                            final double relativeAccuracy,
091                                            final double absoluteAccuracy)
092        throws NotStrictlyPositiveException {
093        this(n, relativeAccuracy, absoluteAccuracy,
094             DEFAULT_MIN_ITERATIONS_COUNT, DEFAULT_MAX_ITERATIONS_COUNT);
095    }
096
097    /**
098     * Builds an integrator with given iteration counts.
099     *
100     * @param n Number of integration points.
101     * @param minimalIterationCount Minimum number of iterations.
102     * @param maximalIterationCount Maximum number of iterations.
103     * @throws NotStrictlyPositiveException if minimal number of iterations
104     * is not strictly positive.
105     * @throws NumberIsTooSmallException if maximal number of iterations
106     * is smaller than or equal to the minimal number of iterations.
107     * @throws NotStrictlyPositiveException if {@code n < 1}.
108     */
109    public IterativeLegendreGaussIntegrator(final int n,
110                                            final int minimalIterationCount,
111                                            final int maximalIterationCount)
112        throws NotStrictlyPositiveException, NumberIsTooSmallException {
113        this(n, DEFAULT_RELATIVE_ACCURACY, DEFAULT_ABSOLUTE_ACCURACY,
114             minimalIterationCount, maximalIterationCount);
115    }
116
117    /** {@inheritDoc} */
118    @Override
119    protected double doIntegrate()
120        throws MathIllegalArgumentException, TooManyEvaluationsException, MaxCountExceededException {
121        // Compute first estimate with a single step.
122        double oldt = stage(1);
123
124        int n = 2;
125        while (true) {
126            // Improve integral with a larger number of steps.
127            final double t = stage(n);
128
129            // Estimate the error.
130            final double delta = FastMath.abs(t - oldt);
131            final double limit =
132                FastMath.max(getAbsoluteAccuracy(),
133                             getRelativeAccuracy() * (FastMath.abs(oldt) + FastMath.abs(t)) * 0.5);
134
135            // check convergence
136            if (getIterations() + 1 >= getMinimalIterationCount() &&
137                delta <= limit) {
138                return t;
139            }
140
141            // Prepare next iteration.
142            final double ratio = FastMath.min(4, FastMath.pow(delta / limit, 0.5 / numberOfPoints));
143            n = FastMath.max((int) (ratio * n), n + 1);
144            oldt = t;
145            incrementCount();
146        }
147    }
148
149    /**
150     * Compute the n-th stage integral.
151     *
152     * @param n Number of steps.
153     * @return the value of n-th stage integral.
154     * @throws TooManyEvaluationsException if the maximum number of evaluations
155     * is exceeded.
156     */
157    private double stage(final int n)
158        throws TooManyEvaluationsException {
159        // Function to be integrated is stored in the base class.
160        final UnivariateFunction f = new UnivariateFunction() {
161                /** {@inheritDoc} */
162                public double value(double x)
163                    throws MathIllegalArgumentException, TooManyEvaluationsException {
164                    return computeObjectiveValue(x);
165                }
166            };
167
168        final double min = getMin();
169        final double max = getMax();
170        final double step = (max - min) / n;
171
172        double sum = 0;
173        for (int i = 0; i < n; i++) {
174            // Integrate over each sub-interval [a, b].
175            final double a = min + i * step;
176            final double b = a + step;
177            final GaussIntegrator g = FACTORY.legendreHighPrecision(numberOfPoints, a, b);
178            sum += g.integrate(f);
179        }
180
181        return sum;
182    }
183}