1 /* 2 * Licensed to the Hipparchus project 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 Hipparchus project 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 package org.hipparchus.analysis.integration.gauss; 18 19 import org.hipparchus.CalculusFieldElement; 20 import org.hipparchus.Field; 21 import org.hipparchus.exception.MathIllegalArgumentException; 22 import org.hipparchus.util.MathArrays; 23 import org.hipparchus.util.Pair; 24 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><p> 35 * Recurrence relation and weights computation follow 36 * <a href="http://en.wikipedia.org/wiki/Abramowitz_and_Stegun"> 37 * Abramowitz and Stegun, 1964</a>. 38 * </p><p> 39 * The coefficients of the standard Hermite polynomials grow very rapidly. 40 * In order to avoid overflows, each Hermite polynomial is normalized with 41 * respect to the underlying scalar product. 42 * @param <T> Type of the number used to represent the points and weights of 43 * the quadrature rules. 44 * @since 2.0 45 */ 46 public class FieldHermiteRuleFactory<T extends CalculusFieldElement<T>> extends FieldAbstractRuleFactory<T> { 47 48 /** Simple constructor 49 * @param field field to which rule coefficients belong 50 */ 51 public FieldHermiteRuleFactory(final Field<T> field) { 52 super(field); 53 } 54 55 /** {@inheritDoc} */ 56 @Override 57 protected Pair<T[], T[]> computeRule(int numberOfPoints) 58 throws MathIllegalArgumentException { 59 60 final Field<T> field = getField(); 61 final T sqrtPi = field.getZero().getPi().sqrt(); 62 63 if (numberOfPoints == 1) { 64 // Break recursion. 65 final T[] points = MathArrays.buildArray(field, numberOfPoints); 66 final T[] weights = MathArrays.buildArray(field, numberOfPoints); 67 points[0] = field.getZero(); 68 weights[0] = sqrtPi; 69 return new Pair<>(points, weights); 70 } 71 72 // find nodes as roots of Hermite polynomial 73 final T[] points = findRoots(numberOfPoints, new Hermite<>(field, numberOfPoints)::ratio); 74 enforceSymmetry(points); 75 76 // compute weights 77 final T[] weights = MathArrays.buildArray(field, numberOfPoints); 78 final Hermite<T> hm1 = new Hermite<>(field, numberOfPoints - 1); 79 for (int i = 0; i < numberOfPoints; i++) { 80 final T y = hm1.hNhNm1(points[i])[0]; 81 weights[i] = sqrtPi.divide(y.square().multiply(numberOfPoints)); 82 } 83 84 return new Pair<>(points, weights); 85 86 } 87 88 /** Hermite polynomial, normalized to avoid overflow. 89 * <p> 90 * The regular Hermite polynomials and associated weights are given by: 91 * <pre> 92 * H₀(x) = 1 93 * H₁(x) = 2 x 94 * Hₙ₊₁(x) = 2x Hₙ(x) - 2n Hₙ₋₁(x), and H'ₙ(x) = 2n Hₙ₋₁(x) 95 * wₙ(xᵢ) = [2ⁿ⁻¹ n! √π]/[n Hₙ₋₁(xᵢ)]² 96 * </pre> 97 * </p> 98 * <p> 99 * In order to avoid overflow with normalize the polynomials hₙ(x) = Hₙ(x) / √[2ⁿ n!] 100 * so the recurrence relations and weights become: 101 * <pre> 102 * h₀(x) = 1 103 * h₁(x) = √2 x 104 * hₙ₊₁(x) = [√2 x hₙ(x) - √n hₙ₋₁(x)]/√(n+1), and h'ₙ(x) = 2n hₙ₋₁(x) 105 * uₙ(xᵢ) = √π/[n Nₙ₋₁(xᵢ)²] 106 * </pre> 107 * </p> 108 * @param <T> Type of the field elements. 109 */ 110 private static class Hermite<T extends CalculusFieldElement<T>> { 111 112 /** √2. */ 113 private final T sqrt2; 114 115 /** Degree. */ 116 private final int degree; 117 118 /** Simple constructor. 119 * @param field field to which rule coefficients belong 120 * @param degree polynomial degree 121 */ 122 Hermite(Field<T> field, int degree) { 123 this.sqrt2 = field.getZero().newInstance(2).sqrt(); 124 this.degree = degree; 125 } 126 127 /** Compute ratio H(x)/H'(x). 128 * @param x point at which ratio must be computed 129 * @return ratio H(x)/H'(x) 130 */ 131 public T ratio(T x) { 132 T[] h = hNhNm1(x); 133 return h[0].divide(h[1].multiply(2 * degree)); 134 } 135 136 /** Compute Nₙ(x) and Nₙ₋₁(x). 137 * @param x point at which polynomials are evaluated 138 * @return array containing Nₙ(x) at index 0 and Nₙ₋₁(x) at index 1 139 */ 140 private T[] hNhNm1(final T x) { 141 T[] h = MathArrays.buildArray(x.getField(), 2); 142 h[0] = sqrt2.multiply(x); 143 h[1] = x.getField().getOne(); 144 T sqrtN = x.getField().getOne(); 145 for (int n = 1; n < degree; n++) { 146 // apply recurrence relation hₙ₊₁(x) = [√2 x hₙ(x) - √n hₙ₋₁(x)]/√(n+1) 147 final T sqrtNp = x.getField().getZero().newInstance(n + 1).sqrt(); 148 final T hp = (h[0].multiply(x).multiply(sqrt2).subtract(h[1].multiply(sqrtN))).divide(sqrtNp); 149 h[1] = h[0]; 150 h[0] = hp; 151 sqrtN = sqrtNp; 152 } 153 return h; 154 } 155 156 } 157 158 }