gonum.org/v1/gonum@v0.15.1-0.20240517103525-f853624cb1bb/integrate/quad/quad.go (about) 1 // Copyright ©2015 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 "sync" 10 ) 11 12 // FixedLocationer computes a set of quadrature locations and weights and stores 13 // them in-place into x and weight respectively. The number of points generated is equal to 14 // the len(x). The weights and locations should be chosen such that 15 // 16 // int_min^max f(x) dx ≈ \sum_i w_i f(x_i) 17 type FixedLocationer interface { 18 FixedLocations(x, weight []float64, min, max float64) 19 } 20 21 // FixedLocationSingler wraps the FixedLocationSingle method. 22 type FixedLocationSingler interface { 23 // FixedLocationSingle returns the location and weight for 24 // element k in a fixed quadrature rule with n total samples 25 // and integral bounds from min to max. 26 FixedLocationSingle(n, k int, min, max float64) (x, weight float64) 27 } 28 29 // Fixed approximates the integral of the function f from min to max using a fixed 30 // n-point quadrature rule. During evaluation, f will be evaluated n times using 31 // the weights and locations specified by rule. That is, Fixed estimates 32 // 33 // int_min^max f(x) dx ≈ \sum_i w_i f(x_i) 34 // 35 // If rule is nil, an acceptable default is chosen, otherwise it is 36 // assumed that the properties of the integral match the assumptions of rule. 37 // For example, Legendre assumes that the integration bounds are finite. If 38 // rule is also a FixedLocationSingler, the quadrature points are computed 39 // individually rather than as a unit. 40 // 41 // If concurrent <= 0, f is evaluated serially, while if concurrent > 0, f 42 // may be evaluated with at most concurrent simultaneous evaluations. 43 // 44 // min must be less than or equal to max, and n must be positive, otherwise 45 // Fixed will panic. 46 func Fixed(f func(float64) float64, min, max float64, n int, rule FixedLocationer, concurrent int) float64 { 47 // TODO(btracey): When there are Hermite polynomial quadrature, add an additional 48 // example to the documentation comment that talks about weight functions. 49 if n <= 0 { 50 panic("quad: non-positive number of locations") 51 } 52 if min > max { 53 panic("quad: min > max") 54 } 55 if min == max { 56 return 0 57 } 58 intfunc := f 59 // If rule is non-nil it is assumed that the function and the constraints 60 // of rule are aligned. If it is nil, wrap the function and do something 61 // reasonable. 62 // TODO(btracey): Replace wrapping with other quadrature rules when 63 // we have rules that support infinite-bound integrals. 64 if rule == nil { 65 // int_a^b f(x)dx = int_u^-1(a)^u^-1(b) f(u(t))u'(t)dt 66 switch { 67 case math.IsInf(max, 1) && math.IsInf(min, -1): 68 // u(t) = (t/(1-t^2)) 69 min = -1 70 max = 1 71 intfunc = func(x float64) float64 { 72 v := 1 - x*x 73 return f(x/v) * (1 + x*x) / (v * v) 74 } 75 case math.IsInf(max, 1): 76 // u(t) = a + t / (1-t) 77 a := min 78 min = 0 79 max = 1 80 intfunc = func(x float64) float64 { 81 v := 1 - x 82 return f(a+x/v) / (v * v) 83 } 84 case math.IsInf(min, -1): 85 // u(t) = a - (1-t)/t 86 a := max 87 min = 0 88 max = 1 89 intfunc = func(x float64) float64 { 90 return f(a-(1-x)/x) / (x * x) 91 } 92 } 93 rule = Legendre{} 94 } 95 singler, isSingler := rule.(FixedLocationSingler) 96 97 var xs, weights []float64 98 if !isSingler { 99 xs = make([]float64, n) 100 weights = make([]float64, n) 101 rule.FixedLocations(xs, weights, min, max) 102 } 103 104 if concurrent > n { 105 concurrent = n 106 } 107 108 if concurrent <= 0 { 109 var integral float64 110 // Evaluate in serial. 111 if isSingler { 112 for k := 0; k < n; k++ { 113 x, weight := singler.FixedLocationSingle(n, k, min, max) 114 integral += weight * intfunc(x) 115 } 116 return integral 117 } 118 for i, x := range xs { 119 integral += weights[i] * intfunc(x) 120 } 121 return integral 122 } 123 124 // Evaluate concurrently 125 tasks := make(chan int) 126 127 // Launch distributor 128 go func() { 129 for i := 0; i < n; i++ { 130 tasks <- i 131 } 132 close(tasks) 133 }() 134 135 var mux sync.Mutex 136 var integral float64 137 var wg sync.WaitGroup 138 wg.Add(concurrent) 139 for i := 0; i < concurrent; i++ { 140 // Launch workers 141 go func() { 142 defer wg.Done() 143 var subIntegral float64 144 for k := range tasks { 145 var x, weight float64 146 if isSingler { 147 x, weight = singler.FixedLocationSingle(n, k, min, max) 148 } else { 149 x = xs[k] 150 weight = weights[k] 151 } 152 f := intfunc(x) 153 subIntegral += f * weight 154 } 155 mux.Lock() 156 integral += subIntegral 157 mux.Unlock() 158 }() 159 } 160 wg.Wait() 161 return integral 162 }