gonum.org/v1/gonum@v0.14.0/integrate/quad/hermite.go (about)

     1  // Copyright ©2016 The Gonum Authors. All rights reserved.
     2  // Use of this source code is governed by a BSD-style
     3  // license that can be found in the LICENSE file.
     4  
     5  package quad
     6  
     7  import (
     8  	"math"
     9  
    10  	"gonum.org/v1/gonum/floats"
    11  	"gonum.org/v1/gonum/mathext"
    12  )
    13  
    14  // Hermite generates sample locations and weights for performing quadrature with
    15  // a squared-exponential weight
    16  //
    17  //	int_-inf^inf e^(-x^2) f(x) dx .
    18  type Hermite struct{}
    19  
    20  func (h Hermite) FixedLocations(x, weight []float64, min, max float64) {
    21  	// TODO(btracey): Implement the case where x > 20, x < 200 so that we don't
    22  	// need to store all of that data.
    23  
    24  	// Algorithm adapted from Chebfun http://www.chebfun.org/.
    25  	//
    26  	// References:
    27  	// Algorithm:
    28  	// G. H. Golub and J. A. Welsch, "Calculation of Gauss quadrature rules",
    29  	// Math. Comp. 23:221-230, 1969.
    30  	// A. Glaser, X. Liu and V. Rokhlin, "A fast algorithm for the
    31  	// calculation of the roots of special functions", SIAM Journal
    32  	// on Scientific Computing", 29(4):1420-1438:, 2007.
    33  	// A. Townsend, T. Trogdon, and S.Olver, Fast computation of Gauss quadrature
    34  	// nodes and weights on the whole real line, IMA J. Numer. Anal., 36: 337–358,
    35  	// 2016. http://arxiv.org/abs/1410.5286
    36  
    37  	if len(x) != len(weight) {
    38  		panic("hermite: slice length mismatch")
    39  	}
    40  	if min >= max {
    41  		panic("hermite: min >= max")
    42  	}
    43  	if !math.IsInf(min, -1) || !math.IsInf(max, 1) {
    44  		panic("hermite: non-infinite bound")
    45  	}
    46  	h.locations(x, weight)
    47  }
    48  
    49  func (h Hermite) locations(x, weights []float64) {
    50  	n := len(x)
    51  	switch {
    52  	case 0 < n && n <= 200:
    53  		copy(x, xCacheHermite[n-1])
    54  		copy(weights, wCacheHermite[n-1])
    55  	case n > 200:
    56  		h.locationsAsy(x, weights)
    57  	}
    58  }
    59  
    60  // Algorithm adapted from Chebfun http://www.chebfun.org/. Specific code
    61  // https://github.com/chebfun/chebfun/blob/development/hermpts.m.
    62  
    63  // Original Copyright Notice:
    64  
    65  /*
    66  Copyright (c) 2015, The Chancellor, Masters and Scholars of the University
    67  of Oxford, and the Chebfun Developers. All rights reserved.
    68  
    69  Redistribution and use in source and binary forms, with or without
    70  modification, are permitted provided that the following conditions are met:
    71      * Redistributions of source code must retain the above copyright
    72        notice, this list of conditions and the following disclaimer.
    73      * Redistributions in binary form must reproduce the above copyright
    74        notice, this list of conditions and the following disclaimer in the
    75        documentation and/or other materials provided with the distribution.
    76      * Neither the name of the University of Oxford nor the names of its
    77        contributors may be used to endorse or promote products derived from
    78        this software without specific prior written permission.
    79  
    80  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
    81  ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
    82  WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
    83  DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR
    84  ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
    85  (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
    86  LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND
    87  ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
    88  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
    89  SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
    90  */
    91  
    92  // locationAsy returns the node locations and weights of a Hermite quadrature rule
    93  // with len(x) points.
    94  func (h Hermite) locationsAsy(x, w []float64) {
    95  	// A. Townsend, T. Trogdon, and S.Olver, Fast computation of Gauss quadrature
    96  	// nodes and weights the whole real line, IMA J. Numer. Anal.,
    97  	// 36: 337–358, 2016. http://arxiv.org/abs/1410.5286
    98  
    99  	// Find the positive locations and weights.
   100  	n := len(x)
   101  	l := n / 2
   102  	xa := x[l:]
   103  	wa := w[l:]
   104  	for i := range xa {
   105  		xa[i], wa[i] = h.locationsAsy0(i, n)
   106  	}
   107  	// Flip around zero -- copy the negative x locations with the corresponding
   108  	// weights.
   109  	if n%2 == 0 {
   110  		l--
   111  	}
   112  	for i, v := range xa {
   113  		x[l-i] = -v
   114  	}
   115  	for i, v := range wa {
   116  		w[l-i] = v
   117  	}
   118  	sumW := floats.Sum(w)
   119  	c := math.SqrtPi / sumW
   120  	floats.Scale(c, w)
   121  }
   122  
   123  // locationsAsy0 returns the location and weight for location i in an n-point
   124  // quadrature rule. The rule is symmetric, so i should be <= n/2 + n%2.
   125  func (h Hermite) locationsAsy0(i, n int) (x, w float64) {
   126  	const convTol = 1e-16
   127  	const convIter = 20
   128  	theta0 := h.hermiteInitialGuess(i, n)
   129  	t0 := theta0 / math.Sqrt(2*float64(n)+1)
   130  	theta0 = math.Acos(t0)
   131  	sqrt2np1 := math.Sqrt(2*float64(n) + 1)
   132  	var vali, dvali float64
   133  	for k := 0; k < convIter; k++ {
   134  		vali, dvali = h.hermpolyAsyAiry(i, n, theta0)
   135  		dt := -vali / (math.Sqrt2 * sqrt2np1 * dvali * math.Sin(theta0))
   136  		theta0 -= dt
   137  		if math.Abs(theta0) < convTol {
   138  			break
   139  		}
   140  	}
   141  	x = sqrt2np1 * math.Cos(theta0)
   142  	ders := x*vali + math.Sqrt2*dvali
   143  	w = math.Exp(-x*x) / (ders * ders)
   144  	return x, w
   145  }
   146  
   147  // hermpolyAsyAiry evaluates the Hermite polynomials using the Airy asymptotic
   148  // formula in theta-space.
   149  func (h Hermite) hermpolyAsyAiry(i, n int, t float64) (valVec, dvalVec float64) {
   150  	musq := 2*float64(n) + 1
   151  	cosT := math.Cos(t)
   152  	sinT := math.Sin(t)
   153  	sin2T := 2 * cosT * sinT
   154  	eta := 0.5*t - 0.25*sin2T
   155  	chi := -math.Pow(3*eta/2, 2.0/3)
   156  	phi := math.Pow(-chi/(sinT*sinT), 1.0/4)
   157  	cnst := 2 * math.SqrtPi * math.Pow(musq, 1.0/6) * phi
   158  	airy0 := real(mathext.AiryAi(complex(math.Pow(musq, 2.0/3)*chi, 0)))
   159  	airy1 := real(mathext.AiryAiDeriv(complex(math.Pow(musq, 2.0/3)*chi, 0)))
   160  
   161  	// Terms in 12.10.43:
   162  	const (
   163  		a1 = 15.0 / 144
   164  		b1 = -7.0 / 5 * a1
   165  		a2 = 5.0 * 7 * 9 * 11.0 / 2.0 / 144.0 / 144.0
   166  		b2 = -13.0 / 11 * a2
   167  		a3 = 7.0 * 9 * 11 * 13 * 15 * 17 / 6.0 / 144.0 / 144.0 / 144.0
   168  		b3 = -19.0 / 17 * a3
   169  	)
   170  
   171  	// Pre-compute terms.
   172  	cos2T := cosT * cosT
   173  	cos3T := cos2T * cosT
   174  	cos4T := cos3T * cosT
   175  	cos5T := cos4T * cosT
   176  	cos7T := cos5T * cos2T
   177  	cos9T := cos7T * cos2T
   178  
   179  	chi2 := chi * chi
   180  	chi3 := chi2 * chi
   181  	chi4 := chi3 * chi
   182  	chi5 := chi4 * chi
   183  
   184  	phi6 := math.Pow(phi, 6)
   185  	phi12 := phi6 * phi6
   186  	phi18 := phi12 * phi6
   187  
   188  	// u polynomials in 12.10.9.
   189  	u1 := (cos3T - 6*cosT) / 24.0
   190  	u2 := (-9*cos4T + 249*cos2T + 145) / 1152.0
   191  	u3 := (-4042*cos9T + 18189*cos7T - 28287*cos5T - 151995*cos3T - 259290*cosT) / 414720.0
   192  
   193  	val := airy0
   194  	B0 := -(phi6*u1 + a1) / chi2
   195  	val += B0 * airy1 / math.Pow(musq, 4.0/3)
   196  	A1 := (phi12*u2 + b1*phi6*u1 + b2) / chi3
   197  	val += A1 * airy0 / (musq * musq)
   198  	B1 := -(phi18*u3 + a1*phi12*u2 + a2*phi6*u1 + a3) / chi5
   199  	val += B1 * airy1 / math.Pow(musq, 4.0/3+2)
   200  	val *= cnst
   201  
   202  	// Derivative.
   203  	eta = 0.5*t - 0.25*sin2T
   204  	chi = -math.Pow(3*eta/2, 2.0/3)
   205  	phi = math.Pow(-chi/(sinT*sinT), 1.0/4)
   206  	cnst = math.Sqrt2 * math.SqrtPi * math.Pow(musq, 1.0/3) / phi
   207  
   208  	// v polynomials in 12.10.10.
   209  	v1 := (cos3T + 6*cosT) / 24
   210  	v2 := (15*cos4T - 327*cos2T - 143) / 1152
   211  	v3 := (259290*cosT + 238425*cos3T - 36387*cos5T + 18189*cos7T - 4042*cos9T) / 414720
   212  
   213  	C0 := -(phi6*v1 + b1) / chi
   214  	dval := C0 * airy0 / math.Pow(musq, 2.0/3)
   215  	dval += airy1
   216  	C1 := -(phi18*v3 + b1*phi12*v2 + b2*phi6*v1 + b3) / chi4
   217  	dval += C1 * airy0 / math.Pow(musq, 2.0/3+2)
   218  	D1 := (phi12*v2 + a1*phi6*v1 + a2) / chi3
   219  	dval += D1 * airy1 / (musq * musq)
   220  	dval *= cnst
   221  	return val, dval
   222  }
   223  
   224  // hermiteInitialGuess returns the initial guess for node i in an n-point Hermite
   225  // quadrature rule. The rule is symmetric, so i should be <= n/2 + n%2.
   226  func (h Hermite) hermiteInitialGuess(i, n int) float64 {
   227  	// There are two different formulas for the initial guesses of the hermite
   228  	// quadrature locations. The first uses the Gatteschi formula and is good
   229  	// near x = sqrt(n+0.5)
   230  	//  [1] L. Gatteschi, Asymptotics and bounds for the zeros of Laguerre
   231  	//  polynomials: a survey, J. Comput. Appl. Math., 144 (2002), pp. 7-27.
   232  	// The second is the Tricomi initial guesses, good near x = 0. This is
   233  	// equation 2.1 in [1] and is originally from
   234  	//  [2] F. G. Tricomi, Sugli zeri delle funzioni di cui si conosce una
   235  	//  rappresentazione asintotica, Ann. Mat. Pura Appl. 26 (1947), pp. 283-300.
   236  
   237  	// If the number of points is odd, there is a quadrature point at 1, which
   238  	// has an initial guess of 0.
   239  	if n%2 == 1 {
   240  		if i == 0 {
   241  			return 0
   242  		}
   243  		i--
   244  	}
   245  
   246  	m := n / 2
   247  	a := -0.5
   248  	if n%2 == 1 {
   249  		a = 0.5
   250  	}
   251  	nu := 4*float64(m) + 2*a + 2
   252  
   253  	// Find the split between Gatteschi guesses and Tricomi guesses.
   254  	p := 0.4985 + math.SmallestNonzeroFloat64
   255  	pidx := int(math.Floor(p * float64(n)))
   256  
   257  	// Use the Tricomi initial guesses in the first half where x is nearer to zero.
   258  	// Note: zeros of besselj(+/-.5,x) are integer and half-integer multiples of pi.
   259  	if i < pidx {
   260  		rhs := math.Pi * (4*float64(m) - 4*(float64(i)+1) + 3) / nu
   261  		tnk := math.Pi / 2
   262  		for k := 0; k < 7; k++ {
   263  			val := tnk - math.Sin(tnk) - rhs
   264  			dval := 1 - math.Cos(tnk)
   265  			dTnk := val / dval
   266  			tnk -= dTnk
   267  			if math.Abs(dTnk) < 1e-14 {
   268  				break
   269  			}
   270  		}
   271  		vc := math.Cos(tnk / 2)
   272  		t := vc * vc
   273  		return math.Sqrt(nu*t - (5.0/(4.0*(1-t)*(1-t))-1.0/(1-t)-1+3*a*a)/3/nu)
   274  	}
   275  
   276  	// Use Gatteschi guesses in the second half where x is nearer to sqrt(n+0.5)
   277  	i = m - (i + 1)
   278  	var ar float64
   279  	if i < len(airyRtsExact) {
   280  		ar = airyRtsExact[i]
   281  	} else {
   282  		t := 3.0 / 8 * math.Pi * (4*(float64(i)+1) - 1)
   283  		ar = math.Pow(t, 2.0/3) * (1 +
   284  			5.0/48*math.Pow(t, -2) -
   285  			5.0/36*math.Pow(t, -4) +
   286  			77125.0/82944*math.Pow(t, -6) -
   287  			108056875.0/6967296*math.Pow(t, -8) +
   288  			162375596875.0/334430208*math.Pow(t, -10))
   289  	}
   290  	r := nu + math.Pow(2, 2.0/3)*ar*math.Pow(nu, 1.0/3) +
   291  		0.2*math.Pow(2, 4.0/3)*ar*ar*math.Pow(nu, -1.0/3) +
   292  		(11.0/35-a*a-12.0/175*ar*ar*ar)/nu +
   293  		(16.0/1575*ar+92.0/7875*math.Pow(ar, 4))*math.Pow(2, 2.0/3)*math.Pow(nu, -5.0/3) -
   294  		(15152.0/3031875*math.Pow(ar, 5)+1088.0/121275*ar*ar)*math.Pow(2, 1.0/3)*math.Pow(nu, -7.0/3)
   295  	if r < 0 {
   296  		ar = 0
   297  	} else {
   298  		ar = math.Sqrt(r)
   299  	}
   300  	return ar
   301  }
   302  
   303  // airyRtsExact are the first airy roots.
   304  var airyRtsExact = []float64{
   305  	-2.338107410459762,
   306  	-4.087949444130970,
   307  	-5.520559828095555,
   308  	-6.786708090071765,
   309  	-7.944133587120863,
   310  	-9.022650853340979,
   311  	-10.040174341558084,
   312  	-11.008524303733260,
   313  	-11.936015563236262,
   314  	-12.828776752865757,
   315  }