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