_base.py 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518
  1. from scipy._lib._array_api import array_namespace, xp_size
  2. from functools import cached_property
  3. class Rule:
  4. """
  5. Base class for numerical integration algorithms (cubatures).
  6. Finds an estimate for the integral of ``f`` over the region described by two arrays
  7. ``a`` and ``b`` via `estimate`, and find an estimate for the error of this
  8. approximation via `estimate_error`.
  9. If a subclass does not implement its own `estimate_error`, then it will use a
  10. default error estimate based on the difference between the estimate over the whole
  11. region and the sum of estimates over that region divided into ``2^ndim`` subregions.
  12. See Also
  13. --------
  14. FixedRule
  15. Examples
  16. --------
  17. In the following, a custom rule is created which uses 3D Genz-Malik cubature for
  18. the estimate of the integral, and the difference between this estimate and a less
  19. accurate estimate using 5-node Gauss-Legendre quadrature as an estimate for the
  20. error.
  21. >>> import numpy as np
  22. >>> from scipy.integrate import cubature
  23. >>> from scipy.integrate._rules import (
  24. ... Rule, ProductNestedFixed, GenzMalikCubature, GaussLegendreQuadrature
  25. ... )
  26. >>> def f(x, r, alphas):
  27. ... # f(x) = cos(2*pi*r + alpha @ x)
  28. ... # Need to allow r and alphas to be arbitrary shape
  29. ... npoints, ndim = x.shape[0], x.shape[-1]
  30. ... alphas_reshaped = alphas[np.newaxis, :]
  31. ... x_reshaped = x.reshape(npoints, *([1]*(len(alphas.shape) - 1)), ndim)
  32. ... return np.cos(2*np.pi*r + np.sum(alphas_reshaped * x_reshaped, axis=-1))
  33. >>> genz = GenzMalikCubature(ndim=3)
  34. >>> gauss = GaussKronrodQuadrature(npoints=21)
  35. >>> # Gauss-Kronrod is 1D, so we find the 3D product rule:
  36. >>> gauss_3d = ProductNestedFixed([gauss, gauss, gauss])
  37. >>> class CustomRule(Rule):
  38. ... def estimate(self, f, a, b, args=()):
  39. ... return genz.estimate(f, a, b, args)
  40. ... def estimate_error(self, f, a, b, args=()):
  41. ... return np.abs(
  42. ... genz.estimate(f, a, b, args)
  43. ... - gauss_3d.estimate(f, a, b, args)
  44. ... )
  45. >>> rng = np.random.default_rng()
  46. >>> res = cubature(
  47. ... f=f,
  48. ... a=np.array([0, 0, 0]),
  49. ... b=np.array([1, 1, 1]),
  50. ... rule=CustomRule(),
  51. ... args=(rng.random((2,)), rng.random((3, 2, 3)))
  52. ... )
  53. >>> res.estimate
  54. array([[-0.95179502, 0.12444608],
  55. [-0.96247411, 0.60866385],
  56. [-0.97360014, 0.25515587]])
  57. """
  58. def estimate(self, f, a, b, args=()):
  59. r"""
  60. Calculate estimate of integral of `f` in rectangular region described by
  61. corners `a` and ``b``.
  62. Parameters
  63. ----------
  64. f : callable
  65. Function to integrate. `f` must have the signature::
  66. f(x : ndarray, \*args) -> ndarray
  67. `f` should accept arrays ``x`` of shape::
  68. (npoints, ndim)
  69. and output arrays of shape::
  70. (npoints, output_dim_1, ..., output_dim_n)
  71. In this case, `estimate` will return arrays of shape::
  72. (output_dim_1, ..., output_dim_n)
  73. a, b : ndarray
  74. Lower and upper limits of integration as rank-1 arrays specifying the left
  75. and right endpoints of the intervals being integrated over. Infinite limits
  76. are currently not supported.
  77. args : tuple, optional
  78. Additional positional args passed to ``f``, if any.
  79. Returns
  80. -------
  81. est : ndarray
  82. Result of estimation. If `f` returns arrays of shape ``(npoints,
  83. output_dim_1, ..., output_dim_n)``, then `est` will be of shape
  84. ``(output_dim_1, ..., output_dim_n)``.
  85. """
  86. raise NotImplementedError
  87. def estimate_error(self, f, a, b, args=()):
  88. r"""
  89. Estimate the error of the approximation for the integral of `f` in rectangular
  90. region described by corners `a` and `b`.
  91. If a subclass does not override this method, then a default error estimator is
  92. used. This estimates the error as ``|est - refined_est|`` where ``est`` is
  93. ``estimate(f, a, b)`` and ``refined_est`` is the sum of
  94. ``estimate(f, a_k, b_k)`` where ``a_k, b_k`` are the coordinates of each
  95. subregion of the region described by ``a`` and ``b``. In the 1D case, this
  96. is equivalent to comparing the integral over an entire interval ``[a, b]`` to
  97. the sum of the integrals over the left and right subintervals, ``[a, (a+b)/2]``
  98. and ``[(a+b)/2, b]``.
  99. Parameters
  100. ----------
  101. f : callable
  102. Function to estimate error for. `f` must have the signature::
  103. f(x : ndarray, \*args) -> ndarray
  104. `f` should accept arrays `x` of shape::
  105. (npoints, ndim)
  106. and output arrays of shape::
  107. (npoints, output_dim_1, ..., output_dim_n)
  108. In this case, `estimate` will return arrays of shape::
  109. (output_dim_1, ..., output_dim_n)
  110. a, b : ndarray
  111. Lower and upper limits of integration as rank-1 arrays specifying the left
  112. and right endpoints of the intervals being integrated over. Infinite limits
  113. are currently not supported.
  114. args : tuple, optional
  115. Additional positional args passed to `f`, if any.
  116. Returns
  117. -------
  118. err_est : ndarray
  119. Result of error estimation. If `f` returns arrays of shape
  120. ``(npoints, output_dim_1, ..., output_dim_n)``, then `est` will be
  121. of shape ``(output_dim_1, ..., output_dim_n)``.
  122. """
  123. est = self.estimate(f, a, b, args)
  124. refined_est = 0
  125. for a_k, b_k in _split_subregion(a, b):
  126. refined_est += self.estimate(f, a_k, b_k, args)
  127. return self.xp.abs(est - refined_est)
  128. class FixedRule(Rule):
  129. """
  130. A rule implemented as the weighted sum of function evaluations at fixed nodes.
  131. Attributes
  132. ----------
  133. nodes_and_weights : (ndarray, ndarray)
  134. A tuple ``(nodes, weights)`` of nodes at which to evaluate ``f`` and the
  135. corresponding weights. ``nodes`` should be of shape ``(num_nodes,)`` for 1D
  136. cubature rules (quadratures) and more generally for N-D cubature rules, it
  137. should be of shape ``(num_nodes, ndim)``. ``weights`` should be of shape
  138. ``(num_nodes,)``. The nodes and weights should be for integrals over
  139. :math:`[-1, 1]^n`.
  140. See Also
  141. --------
  142. GaussLegendreQuadrature, GaussKronrodQuadrature, GenzMalikCubature
  143. Examples
  144. --------
  145. Implementing Simpson's 1/3 rule:
  146. >>> import numpy as np
  147. >>> from scipy.integrate._rules import FixedRule
  148. >>> class SimpsonsQuad(FixedRule):
  149. ... @property
  150. ... def nodes_and_weights(self):
  151. ... nodes = np.array([-1, 0, 1])
  152. ... weights = np.array([1/3, 4/3, 1/3])
  153. ... return (nodes, weights)
  154. >>> rule = SimpsonsQuad()
  155. >>> rule.estimate(
  156. ... f=lambda x: x**2,
  157. ... a=np.array([0]),
  158. ... b=np.array([1]),
  159. ... )
  160. [0.3333333]
  161. """
  162. def __init__(self):
  163. self.xp = None
  164. @property
  165. def nodes_and_weights(self):
  166. raise NotImplementedError
  167. def estimate(self, f, a, b, args=()):
  168. r"""
  169. Calculate estimate of integral of `f` in rectangular region described by
  170. corners `a` and `b` as ``sum(weights * f(nodes))``.
  171. Nodes and weights will automatically be adjusted from calculating integrals over
  172. :math:`[-1, 1]^n` to :math:`[a, b]^n`.
  173. Parameters
  174. ----------
  175. f : callable
  176. Function to integrate. `f` must have the signature::
  177. f(x : ndarray, \*args) -> ndarray
  178. `f` should accept arrays `x` of shape::
  179. (npoints, ndim)
  180. and output arrays of shape::
  181. (npoints, output_dim_1, ..., output_dim_n)
  182. In this case, `estimate` will return arrays of shape::
  183. (output_dim_1, ..., output_dim_n)
  184. a, b : ndarray
  185. Lower and upper limits of integration as rank-1 arrays specifying the left
  186. and right endpoints of the intervals being integrated over. Infinite limits
  187. are currently not supported.
  188. args : tuple, optional
  189. Additional positional args passed to `f`, if any.
  190. Returns
  191. -------
  192. est : ndarray
  193. Result of estimation. If `f` returns arrays of shape ``(npoints,
  194. output_dim_1, ..., output_dim_n)``, then `est` will be of shape
  195. ``(output_dim_1, ..., output_dim_n)``.
  196. """
  197. nodes, weights = self.nodes_and_weights
  198. if self.xp is None:
  199. self.xp = array_namespace(nodes)
  200. return _apply_fixed_rule(f, a, b, nodes, weights, args, self.xp)
  201. class NestedFixedRule(FixedRule):
  202. r"""
  203. A cubature rule with error estimate given by the difference between two underlying
  204. fixed rules.
  205. If constructed as ``NestedFixedRule(higher, lower)``, this will use::
  206. estimate(f, a, b) := higher.estimate(f, a, b)
  207. estimate_error(f, a, b) := \|higher.estimate(f, a, b) - lower.estimate(f, a, b)|
  208. (where the absolute value is taken elementwise).
  209. Attributes
  210. ----------
  211. higher : Rule
  212. Higher accuracy rule.
  213. lower : Rule
  214. Lower accuracy rule.
  215. See Also
  216. --------
  217. GaussKronrodQuadrature
  218. Examples
  219. --------
  220. >>> from scipy.integrate import cubature
  221. >>> from scipy.integrate._rules import (
  222. ... GaussLegendreQuadrature, NestedFixedRule, ProductNestedFixed
  223. ... )
  224. >>> higher = GaussLegendreQuadrature(10)
  225. >>> lower = GaussLegendreQuadrature(5)
  226. >>> rule = NestedFixedRule(
  227. ... higher,
  228. ... lower
  229. ... )
  230. >>> rule_2d = ProductNestedFixed([rule, rule])
  231. """
  232. def __init__(self, higher, lower):
  233. self.higher = higher
  234. self.lower = lower
  235. self.xp = None
  236. @property
  237. def nodes_and_weights(self):
  238. if self.higher is not None:
  239. return self.higher.nodes_and_weights
  240. else:
  241. raise NotImplementedError
  242. @property
  243. def lower_nodes_and_weights(self):
  244. if self.lower is not None:
  245. return self.lower.nodes_and_weights
  246. else:
  247. raise NotImplementedError
  248. def estimate_error(self, f, a, b, args=()):
  249. r"""
  250. Estimate the error of the approximation for the integral of `f` in rectangular
  251. region described by corners `a` and `b`.
  252. Parameters
  253. ----------
  254. f : callable
  255. Function to estimate error for. `f` must have the signature::
  256. f(x : ndarray, \*args) -> ndarray
  257. `f` should accept arrays `x` of shape::
  258. (npoints, ndim)
  259. and output arrays of shape::
  260. (npoints, output_dim_1, ..., output_dim_n)
  261. In this case, `estimate` will return arrays of shape::
  262. (output_dim_1, ..., output_dim_n)
  263. a, b : ndarray
  264. Lower and upper limits of integration as rank-1 arrays specifying the left
  265. and right endpoints of the intervals being integrated over. Infinite limits
  266. are currently not supported.
  267. args : tuple, optional
  268. Additional positional args passed to `f`, if any.
  269. Returns
  270. -------
  271. err_est : ndarray
  272. Result of error estimation. If `f` returns arrays of shape
  273. ``(npoints, output_dim_1, ..., output_dim_n)``, then `est` will be
  274. of shape ``(output_dim_1, ..., output_dim_n)``.
  275. """
  276. nodes, weights = self.nodes_and_weights
  277. lower_nodes, lower_weights = self.lower_nodes_and_weights
  278. if self.xp is None:
  279. self.xp = array_namespace(nodes)
  280. error_nodes = self.xp.concat([nodes, lower_nodes], axis=0)
  281. error_weights = self.xp.concat([weights, -lower_weights], axis=0)
  282. return self.xp.abs(
  283. _apply_fixed_rule(f, a, b, error_nodes, error_weights, args, self.xp)
  284. )
  285. class ProductNestedFixed(NestedFixedRule):
  286. """
  287. Find the n-dimensional cubature rule constructed from the Cartesian product of 1-D
  288. `NestedFixedRule` quadrature rules.
  289. Given a list of N 1-dimensional quadrature rules which support error estimation
  290. using NestedFixedRule, this will find the N-dimensional cubature rule obtained by
  291. taking the Cartesian product of their nodes, and estimating the error by taking the
  292. difference with a lower-accuracy N-dimensional cubature rule obtained using the
  293. ``.lower_nodes_and_weights`` rule in each of the base 1-dimensional rules.
  294. Parameters
  295. ----------
  296. base_rules : list of NestedFixedRule
  297. List of base 1-dimensional `NestedFixedRule` quadrature rules.
  298. Attributes
  299. ----------
  300. base_rules : list of NestedFixedRule
  301. List of base 1-dimensional `NestedFixedRule` qudarature rules.
  302. Examples
  303. --------
  304. Evaluate a 2D integral by taking the product of two 1D rules:
  305. >>> import numpy as np
  306. >>> from scipy.integrate import cubature
  307. >>> from scipy.integrate._rules import (
  308. ... ProductNestedFixed, GaussKronrodQuadrature
  309. ... )
  310. >>> def f(x):
  311. ... # f(x) = cos(x_1) + cos(x_2)
  312. ... return np.sum(np.cos(x), axis=-1)
  313. >>> rule = ProductNestedFixed(
  314. ... [GaussKronrodQuadrature(15), GaussKronrodQuadrature(15)]
  315. ... ) # Use 15-point Gauss-Kronrod, which implements NestedFixedRule
  316. >>> a, b = np.array([0, 0]), np.array([1, 1])
  317. >>> rule.estimate(f, a, b) # True value 2*sin(1), approximately 1.6829
  318. np.float64(1.682941969615793)
  319. >>> rule.estimate_error(f, a, b)
  320. np.float64(2.220446049250313e-16)
  321. """
  322. def __init__(self, base_rules):
  323. for rule in base_rules:
  324. if not isinstance(rule, NestedFixedRule):
  325. raise ValueError("base rules for product need to be instance of"
  326. "NestedFixedRule")
  327. self.base_rules = base_rules
  328. self.xp = None
  329. @cached_property
  330. def nodes_and_weights(self):
  331. nodes = _cartesian_product(
  332. [rule.nodes_and_weights[0] for rule in self.base_rules]
  333. )
  334. if self.xp is None:
  335. self.xp = array_namespace(nodes)
  336. weights = self.xp.prod(
  337. _cartesian_product(
  338. [rule.nodes_and_weights[1] for rule in self.base_rules]
  339. ),
  340. axis=-1,
  341. )
  342. return nodes, weights
  343. @cached_property
  344. def lower_nodes_and_weights(self):
  345. nodes = _cartesian_product(
  346. [cubature.lower_nodes_and_weights[0] for cubature in self.base_rules]
  347. )
  348. if self.xp is None:
  349. self.xp = array_namespace(nodes)
  350. weights = self.xp.prod(
  351. _cartesian_product(
  352. [cubature.lower_nodes_and_weights[1] for cubature in self.base_rules]
  353. ),
  354. axis=-1,
  355. )
  356. return nodes, weights
  357. def _cartesian_product(arrays):
  358. xp = array_namespace(*arrays)
  359. arrays_ix = xp.meshgrid(*arrays, indexing='ij')
  360. result = xp.reshape(xp.stack(arrays_ix, axis=-1), (-1, len(arrays)))
  361. return result
  362. def _split_subregion(a, b, xp, split_at=None):
  363. """
  364. Given the coordinates of a region like a=[0, 0] and b=[1, 1], yield the coordinates
  365. of all subregions, which in this case would be::
  366. ([0, 0], [1/2, 1/2]),
  367. ([0, 1/2], [1/2, 1]),
  368. ([1/2, 0], [1, 1/2]),
  369. ([1/2, 1/2], [1, 1])
  370. """
  371. xp = array_namespace(a, b)
  372. if split_at is None:
  373. split_at = (a + b) / 2
  374. left = [xp.stack((a[i], split_at[i])) for i in range(a.shape[0])]
  375. right = [xp.stack((split_at[i], b[i])) for i in range(b.shape[0])]
  376. a_sub = _cartesian_product(left)
  377. b_sub = _cartesian_product(right)
  378. for i in range(a_sub.shape[0]):
  379. yield a_sub[i, ...], b_sub[i, ...]
  380. def _apply_fixed_rule(f, a, b, orig_nodes, orig_weights, args, xp):
  381. # Downcast nodes and weights to common dtype of a and b
  382. result_dtype = a.dtype
  383. orig_nodes = xp.astype(orig_nodes, result_dtype)
  384. orig_weights = xp.astype(orig_weights, result_dtype)
  385. # Ensure orig_nodes are at least 2D, since 1D cubature methods can return arrays of
  386. # shape (npoints,) rather than (npoints, 1)
  387. if orig_nodes.ndim == 1:
  388. orig_nodes = orig_nodes[:, None]
  389. rule_ndim = orig_nodes.shape[-1]
  390. a_ndim = xp_size(a)
  391. b_ndim = xp_size(b)
  392. if rule_ndim != a_ndim or rule_ndim != b_ndim:
  393. raise ValueError(f"rule and function are of incompatible dimension, nodes have"
  394. f"ndim {rule_ndim}, while limit of integration has ndim"
  395. f"a_ndim={a_ndim}, b_ndim={b_ndim}")
  396. lengths = b - a
  397. # The underlying rule is for the hypercube [-1, 1]^n.
  398. #
  399. # To handle arbitrary regions of integration, it's necessary to apply a linear
  400. # change of coordinates to map each interval [a[i], b[i]] to [-1, 1].
  401. nodes = (orig_nodes + 1) * (lengths * 0.5) + a
  402. # Also need to multiply the weights by a scale factor equal to the determinant
  403. # of the Jacobian for this coordinate change.
  404. weight_scale_factor = xp.prod(lengths, dtype=result_dtype) / 2**rule_ndim
  405. weights = orig_weights * weight_scale_factor
  406. f_nodes = f(nodes, *args)
  407. weights_reshaped = xp.reshape(weights, (-1, *([1] * (f_nodes.ndim - 1))))
  408. # f(nodes) will have shape (num_nodes, output_dim_1, ..., output_dim_n)
  409. # Summing along the first axis means estimate will shape (output_dim_1, ...,
  410. # output_dim_n)
  411. est = xp.sum(weights_reshaped * f_nodes, axis=0, dtype=result_dtype)
  412. return est