HermiteRuleFactory.java

  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.  * This is not the original file distributed by the Apache Software Foundation
  19.  * It has been modified by the Hipparchus project
  20.  */
  21. package org.hipparchus.analysis.integration.gauss;

  22. import org.hipparchus.exception.MathIllegalArgumentException;
  23. import org.hipparchus.util.FastMath;
  24. import org.hipparchus.util.Pair;

  25. /**
  26.  * Factory that creates a
  27.  * <a href="http://en.wikipedia.org/wiki/Gauss-Hermite_quadrature">
  28.  * Gauss-type quadrature rule using Hermite polynomials</a>
  29.  * of the first kind.
  30.  * Such a quadrature rule allows the calculation of improper integrals
  31.  * of a function
  32.  * <p>
  33.  *  \(f(x) e^{-x^2}\)
  34.  * </p>
  35.  * <p>
  36.  * Recurrence relation and weights computation follow
  37.  * <a href="http://en.wikipedia.org/wiki/Abramowitz_and_Stegun">
  38.  * Abramowitz and Stegun, 1964</a>.
  39.  * </p>
  40.  *
  41.  */
  42. public class HermiteRuleFactory extends AbstractRuleFactory {

  43.     /** √π. */
  44.     private static final double SQRT_PI = 1.77245385090551602729;

  45.     /** Empty constructor.
  46.      * <p>
  47.      * This constructor is not strictly necessary, but it prevents spurious
  48.      * javadoc warnings with JDK 18 and later.
  49.      * </p>
  50.      * @since 3.0
  51.      */
  52.     public HermiteRuleFactory() { // NOPMD - unnecessary constructor added intentionally to make javadoc happy
  53.         // nothing to do
  54.     }

  55.     /** {@inheritDoc} */
  56.     @Override
  57.     protected Pair<double[], double[]> computeRule(int numberOfPoints)
  58.         throws MathIllegalArgumentException {

  59.         if (numberOfPoints == 1) {
  60.             // Break recursion.
  61.             return new Pair<>(new double[] { 0 } , new double[] { SQRT_PI });
  62.         }

  63.         // find nodes as roots of Hermite polynomial
  64.         final double[] points = findRoots(numberOfPoints, new Hermite(numberOfPoints)::ratio);
  65.         enforceSymmetry(points);

  66.         // compute weights
  67.         final double[] weights = new double[numberOfPoints];
  68.         final Hermite hm1 = new Hermite(numberOfPoints - 1);
  69.         for (int i = 0; i < numberOfPoints; i++) {
  70.             final double y = hm1.hNhNm1(points[i])[0];
  71.             weights[i] = SQRT_PI / (numberOfPoints * y * y);
  72.         }

  73.         return new Pair<>(points, weights);

  74.     }

  75.     /** Hermite polynomial, normalized to avoid overflow.
  76.      * <p>
  77.      * The regular Hermite polynomials and associated weights are given by:
  78.      *   <pre>
  79.      *     H₀(x)   = 1
  80.      *     H₁(x)   = 2 x
  81.      *     Hₙ₊₁(x) = 2x Hₙ(x) - 2n Hₙ₋₁(x), and H'ₙ(x) = 2n Hₙ₋₁(x)
  82.      *     wₙ(xᵢ) = [2ⁿ⁻¹ n! √π]/[n Hₙ₋₁(xᵢ)]²
  83.      *   </pre>
  84.      * </p>
  85.      * <p>
  86.      * In order to avoid overflow with normalize the polynomials hₙ(x) = Hₙ(x) / √[2ⁿ n!]
  87.      * so the recurrence relations and weights become:
  88.      *   <pre>
  89.      *     h₀(x)   = 1
  90.      *     h₁(x)   = √2 x
  91.      *     hₙ₊₁(x) = [√2 x hₙ(x) - √n hₙ₋₁(x)]/√(n+1), and h'ₙ(x) = 2n hₙ₋₁(x)
  92.      *     uₙ(xᵢ) = √π/[n Nₙ₋₁(xᵢ)²]
  93.      *   </pre>
  94.      * </p>
  95.      */
  96.     private static class Hermite {

  97.         /** √2. */
  98.         private static final double SQRT2 = FastMath.sqrt(2);

  99.         /** Degree. */
  100.         private final int degree;

  101.         /** Simple constructor.
  102.          * @param degree polynomial degree
  103.          */
  104.         Hermite(int degree) {
  105.             this.degree = degree;
  106.         }

  107.         /** Compute ratio H(x)/H'(x).
  108.          * @param x point at which ratio must be computed
  109.          * @return ratio H(x)/H'(x)
  110.          */
  111.         public double ratio(double x) {
  112.             double[] h = hNhNm1(x);
  113.             return h[0] / (h[1] * 2 * degree);
  114.         }

  115.         /** Compute Nₙ(x) and Nₙ₋₁(x).
  116.          * @param x point at which polynomials are evaluated
  117.          * @return array containing Nₙ(x) at index 0 and Nₙ₋₁(x) at index 1
  118.          */
  119.         private double[] hNhNm1(final double x) {
  120.             double[] h = { SQRT2 * x, 1 };
  121.             double sqrtN = 1;
  122.             for (int n = 1; n < degree; n++) {
  123.                 // apply recurrence relation hₙ₊₁(x) = [√2 x hₙ(x) - √n hₙ₋₁(x)]/√(n+1)
  124.                 final double sqrtNp = FastMath.sqrt(n + 1);
  125.                 final double hp = (h[0] * x * SQRT2 - h[1] * sqrtN) / sqrtNp;
  126.                 h[1]  = h[0];
  127.                 h[0]  = hp;
  128.                 sqrtN = sqrtNp;
  129.             }
  130.             return h;
  131.         }

  132.     }

  133. }