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