gauss.qbk 4.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111
  1. [/
  2. Copyright (c) 2017 John Maddock
  3. Use, modification and distribution are subject to the
  4. Boost Software License, Version 1.0. (See accompanying file
  5. LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  6. ]
  7. [section:gauss Gauss-Legendre quadrature]
  8. [heading Synopsis]
  9. `#include <boost/math/quadrature/gauss.hpp>`
  10. namespace boost{ namespace math{ namespace quadrature{
  11. template <class Real, unsigned Points, class ``__Policy`` = boost::math::policies::policy<> >
  12. struct gauss
  13. {
  14. static const RandomAccessContainer& abscissa();
  15. static const RandomAccessContainer& weights();
  16. template <class F>
  17. static auto integrate(F f, Real* pL1 = nullptr)->decltype(std::declval<F>()(std::declval<Real>()))
  18. template <class F>
  19. static auto integrate(F f, Real a, Real b, Real* pL1 = nullptr)->decltype(std::declval<F>()(std::declval<Real>()))
  20. };
  21. }}} // namespaces
  22. [heading description]
  23. The `gauss` class template performs "one shot" non-adaptive Gauss-Legendre integration on some arbitrary function /f/ using the
  24. number of evaluation points as specified by /Points/.
  25. This is intentionally a very simple quadrature routine, it obtains no estimate of the error, and is not adaptive, but is very efficient
  26. in simple cases that involve integrating smooth "bell like" functions and functions with rapidly convergent power series.
  27. static const RandomAccessContainer& abscissa();
  28. static const RandomAccessContainer& weights();
  29. These functions provide direct access to the abscissa and weights used to perform the quadrature: the return type depends on the
  30. /Points/ template parameter, but is always a RandomAccessContainer type. Note that only positive (or zero) abscissa and weights
  31. are stored.
  32. template <class F>
  33. static auto integrate(F f, Real* pL1 = nullptr)->decltype(std::declval<F>()(std::declval<Real>()))
  34. Integrates /f/ over (-1,1), and optionally sets `*pL1` to the L1 norm of the returned value: if this is substantially larger
  35. than the return value, then the sum was ill-conditioned. Note however, that no error estimate is available.
  36. template <class F>
  37. static auto integrate(F f, Real a, Real b, Real* pL1 = nullptr)->decltype(std::declval<F>()(std::declval<Real>()))
  38. Integrates /f/ over (a,b), and optionally sets `*pL1` to the L1 norm of the returned value: if this is substantially larger
  39. than the return value, then the sum was ill-conditioned. Note however, that no error estimate is available. This function
  40. supports both finite and infinite /a/ and /b/, as long as `a < b`.
  41. The Gaussian quadrature routine support both real and complex-valued quadrature.
  42. For example, the Lambert-W function admits the integral representation
  43. [expression ['W(z) = 1/2\u03A0 \u222B[sub -\u03A0][super \u03A0] ((1- v cot(v) )^2 + v^2)/(z + v csc(v) exp(-v cot(v))) dv]]
  44. so it can be effectively computed via Gaussian quadrature using the following code:
  45. Complex z{2, 3};
  46. auto lw = [&z](Real v)->Complex {
  47. using std::cos;
  48. using std::sin;
  49. using std::exp;
  50. Real sinv = sin(v);
  51. Real cosv = cos(v);
  52. Real cotv = cosv/sinv;
  53. Real cscv = 1/sinv;
  54. Real t = (1-v*cotv)*(1-v*cotv) + v*v;
  55. Real x = v*cscv*exp(-v*cotv);
  56. Complex den = z + x;
  57. Complex num = t*(z/pi<Real>());
  58. Complex res = num/den;
  59. return res;
  60. };
  61. boost::math::quadrature::gauss<Real, 30> integrator;
  62. Complex W = integrator.integrate(lw, (Real) 0, pi<Real>());
  63. [heading Choosing the number of points]
  64. Internally class `gauss` has pre-computed tables of abscissa and weights for 7, 15, 20, 25 and 30 points at up to 100-decimal
  65. digit precision. That means that using for example, `gauss<double, 30>::integrate` incurs absolutely zero set-up overhead from
  66. computing the abscissa/weight pairs. When using multiprecision types with less than 100 digits of precision, then there is a small
  67. initial one time cost, while the abscissa/weight pairs are constructed from strings.
  68. However, for types with higher precision, or numbers of points other than those given above, the abscissa/weight pairs are computed
  69. when first needed and then cached for future use, which does incur a noticeable overhead. If this is likely to be an issue, then
  70. * Defining BOOST_MATH_GAUSS_NO_COMPUTE_ON_DEMAND will result in a compile-time error, whenever a combination of number type
  71. and number of points is used which does not have pre-computed values.
  72. * There is a program [@../../tools/gauss_kronrod_constants.cpp gauss_kronrod_constants.cpp] which was used to provide the
  73. pre-computed values already in gauss.hpp. The program can be trivially modified to generate code and constants for other precisions
  74. and numbers of points.
  75. [heading Examples]
  76. [import ../../example/gauss_example.cpp]
  77. [gauss_example]
  78. [endsect] [/section:gauss Gauss-Legendre quadrature]