_numdiff.py 28 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761
  1. """Routines for numerical differentiation."""
  2. import functools
  3. import numpy as np
  4. from numpy.linalg import norm
  5. from scipy.sparse.linalg import LinearOperator
  6. from ..sparse import issparse, csc_matrix, csr_matrix, coo_matrix, find
  7. from ._group_columns import group_dense, group_sparse
  8. def _adjust_scheme_to_bounds(x0, h, num_steps, scheme, lb, ub):
  9. """Adjust final difference scheme to the presence of bounds.
  10. Parameters
  11. ----------
  12. x0 : ndarray, shape (n,)
  13. Point at which we wish to estimate derivative.
  14. h : ndarray, shape (n,)
  15. Desired absolute finite difference steps.
  16. num_steps : int
  17. Number of `h` steps in one direction required to implement finite
  18. difference scheme. For example, 2 means that we need to evaluate
  19. f(x0 + 2 * h) or f(x0 - 2 * h)
  20. scheme : {'1-sided', '2-sided'}
  21. Whether steps in one or both directions are required. In other
  22. words '1-sided' applies to forward and backward schemes, '2-sided'
  23. applies to center schemes.
  24. lb : ndarray, shape (n,)
  25. Lower bounds on independent variables.
  26. ub : ndarray, shape (n,)
  27. Upper bounds on independent variables.
  28. Returns
  29. -------
  30. h_adjusted : ndarray, shape (n,)
  31. Adjusted absolute step sizes. Step size decreases only if a sign flip
  32. or switching to one-sided scheme doesn't allow to take a full step.
  33. use_one_sided : ndarray of bool, shape (n,)
  34. Whether to switch to one-sided scheme. Informative only for
  35. ``scheme='2-sided'``.
  36. """
  37. if scheme == '1-sided':
  38. use_one_sided = np.ones_like(h, dtype=bool)
  39. elif scheme == '2-sided':
  40. h = np.abs(h)
  41. use_one_sided = np.zeros_like(h, dtype=bool)
  42. else:
  43. raise ValueError("`scheme` must be '1-sided' or '2-sided'.")
  44. if np.all((lb == -np.inf) & (ub == np.inf)):
  45. return h, use_one_sided
  46. h_total = h * num_steps
  47. h_adjusted = h.copy()
  48. lower_dist = x0 - lb
  49. upper_dist = ub - x0
  50. if scheme == '1-sided':
  51. x = x0 + h_total
  52. violated = (x < lb) | (x > ub)
  53. fitting = np.abs(h_total) <= np.maximum(lower_dist, upper_dist)
  54. h_adjusted[violated & fitting] *= -1
  55. forward = (upper_dist >= lower_dist) & ~fitting
  56. h_adjusted[forward] = upper_dist[forward] / num_steps
  57. backward = (upper_dist < lower_dist) & ~fitting
  58. h_adjusted[backward] = -lower_dist[backward] / num_steps
  59. elif scheme == '2-sided':
  60. central = (lower_dist >= h_total) & (upper_dist >= h_total)
  61. forward = (upper_dist >= lower_dist) & ~central
  62. h_adjusted[forward] = np.minimum(
  63. h[forward], 0.5 * upper_dist[forward] / num_steps)
  64. use_one_sided[forward] = True
  65. backward = (upper_dist < lower_dist) & ~central
  66. h_adjusted[backward] = -np.minimum(
  67. h[backward], 0.5 * lower_dist[backward] / num_steps)
  68. use_one_sided[backward] = True
  69. min_dist = np.minimum(upper_dist, lower_dist) / num_steps
  70. adjusted_central = (~central & (np.abs(h_adjusted) <= min_dist))
  71. h_adjusted[adjusted_central] = min_dist[adjusted_central]
  72. use_one_sided[adjusted_central] = False
  73. return h_adjusted, use_one_sided
  74. @functools.lru_cache()
  75. def _eps_for_method(x0_dtype, f0_dtype, method):
  76. """
  77. Calculates relative EPS step to use for a given data type
  78. and numdiff step method.
  79. Progressively smaller steps are used for larger floating point types.
  80. Parameters
  81. ----------
  82. f0_dtype: np.dtype
  83. dtype of function evaluation
  84. x0_dtype: np.dtype
  85. dtype of parameter vector
  86. method: {'2-point', '3-point', 'cs'}
  87. Returns
  88. -------
  89. EPS: float
  90. relative step size. May be np.float16, np.float32, np.float64
  91. Notes
  92. -----
  93. The default relative step will be np.float64. However, if x0 or f0 are
  94. smaller floating point types (np.float16, np.float32), then the smallest
  95. floating point type is chosen.
  96. """
  97. # the default EPS value
  98. EPS = np.finfo(np.float64).eps
  99. x0_is_fp = False
  100. if np.issubdtype(x0_dtype, np.inexact):
  101. # if you're a floating point type then over-ride the default EPS
  102. EPS = np.finfo(x0_dtype).eps
  103. x0_itemsize = np.dtype(x0_dtype).itemsize
  104. x0_is_fp = True
  105. if np.issubdtype(f0_dtype, np.inexact):
  106. f0_itemsize = np.dtype(f0_dtype).itemsize
  107. # choose the smallest itemsize between x0 and f0
  108. if x0_is_fp and f0_itemsize < x0_itemsize:
  109. EPS = np.finfo(f0_dtype).eps
  110. if method in ["2-point", "cs"]:
  111. return EPS**0.5
  112. elif method in ["3-point"]:
  113. return EPS**(1/3)
  114. else:
  115. raise RuntimeError("Unknown step method, should be one of "
  116. "{'2-point', '3-point', 'cs'}")
  117. def _compute_absolute_step(rel_step, x0, f0, method):
  118. """
  119. Computes an absolute step from a relative step for finite difference
  120. calculation.
  121. Parameters
  122. ----------
  123. rel_step: None or array-like
  124. Relative step for the finite difference calculation
  125. x0 : np.ndarray
  126. Parameter vector
  127. f0 : np.ndarray or scalar
  128. method : {'2-point', '3-point', 'cs'}
  129. Returns
  130. -------
  131. h : float
  132. The absolute step size
  133. Notes
  134. -----
  135. `h` will always be np.float64. However, if `x0` or `f0` are
  136. smaller floating point dtypes (e.g. np.float32), then the absolute
  137. step size will be calculated from the smallest floating point size.
  138. """
  139. # this is used instead of np.sign(x0) because we need
  140. # sign_x0 to be 1 when x0 == 0.
  141. sign_x0 = (x0 >= 0).astype(float) * 2 - 1
  142. rstep = _eps_for_method(x0.dtype, f0.dtype, method)
  143. if rel_step is None:
  144. abs_step = rstep * sign_x0 * np.maximum(1.0, np.abs(x0))
  145. else:
  146. # User has requested specific relative steps.
  147. # Don't multiply by max(1, abs(x0) because if x0 < 1 then their
  148. # requested step is not used.
  149. abs_step = rel_step * sign_x0 * np.abs(x0)
  150. # however we don't want an abs_step of 0, which can happen if
  151. # rel_step is 0, or x0 is 0. Instead, substitute a realistic step
  152. dx = ((x0 + abs_step) - x0)
  153. abs_step = np.where(dx == 0,
  154. rstep * sign_x0 * np.maximum(1.0, np.abs(x0)),
  155. abs_step)
  156. return abs_step
  157. def _prepare_bounds(bounds, x0):
  158. """
  159. Prepares new-style bounds from a two-tuple specifying the lower and upper
  160. limits for values in x0. If a value is not bound then the lower/upper bound
  161. will be expected to be -np.inf/np.inf.
  162. Examples
  163. --------
  164. >>> _prepare_bounds([(0, 1, 2), (1, 2, np.inf)], [0.5, 1.5, 2.5])
  165. (array([0., 1., 2.]), array([ 1., 2., inf]))
  166. """
  167. lb, ub = [np.asarray(b, dtype=float) for b in bounds]
  168. if lb.ndim == 0:
  169. lb = np.resize(lb, x0.shape)
  170. if ub.ndim == 0:
  171. ub = np.resize(ub, x0.shape)
  172. return lb, ub
  173. def group_columns(A, order=0):
  174. """Group columns of a 2-D matrix for sparse finite differencing [1]_.
  175. Two columns are in the same group if in each row at least one of them
  176. has zero. A greedy sequential algorithm is used to construct groups.
  177. Parameters
  178. ----------
  179. A : array_like or sparse matrix, shape (m, n)
  180. Matrix of which to group columns.
  181. order : int, iterable of int with shape (n,) or None
  182. Permutation array which defines the order of columns enumeration.
  183. If int or None, a random permutation is used with `order` used as
  184. a random seed. Default is 0, that is use a random permutation but
  185. guarantee repeatability.
  186. Returns
  187. -------
  188. groups : ndarray of int, shape (n,)
  189. Contains values from 0 to n_groups-1, where n_groups is the number
  190. of found groups. Each value ``groups[i]`` is an index of a group to
  191. which ith column assigned. The procedure was helpful only if
  192. n_groups is significantly less than n.
  193. References
  194. ----------
  195. .. [1] A. Curtis, M. J. D. Powell, and J. Reid, "On the estimation of
  196. sparse Jacobian matrices", Journal of the Institute of Mathematics
  197. and its Applications, 13 (1974), pp. 117-120.
  198. """
  199. if issparse(A):
  200. A = csc_matrix(A)
  201. else:
  202. A = np.atleast_2d(A)
  203. A = (A != 0).astype(np.int32)
  204. if A.ndim != 2:
  205. raise ValueError("`A` must be 2-dimensional.")
  206. m, n = A.shape
  207. if order is None or np.isscalar(order):
  208. rng = np.random.RandomState(order)
  209. order = rng.permutation(n)
  210. else:
  211. order = np.asarray(order)
  212. if order.shape != (n,):
  213. raise ValueError("`order` has incorrect shape.")
  214. A = A[:, order]
  215. if issparse(A):
  216. groups = group_sparse(m, n, A.indices, A.indptr)
  217. else:
  218. groups = group_dense(m, n, A)
  219. groups[order] = groups.copy()
  220. return groups
  221. def approx_derivative(fun, x0, method='3-point', rel_step=None, abs_step=None,
  222. f0=None, bounds=(-np.inf, np.inf), sparsity=None,
  223. as_linear_operator=False, args=(), kwargs={}):
  224. """Compute finite difference approximation of the derivatives of a
  225. vector-valued function.
  226. If a function maps from R^n to R^m, its derivatives form m-by-n matrix
  227. called the Jacobian, where an element (i, j) is a partial derivative of
  228. f[i] with respect to x[j].
  229. Parameters
  230. ----------
  231. fun : callable
  232. Function of which to estimate the derivatives. The argument x
  233. passed to this function is ndarray of shape (n,) (never a scalar
  234. even if n=1). It must return 1-D array_like of shape (m,) or a scalar.
  235. x0 : array_like of shape (n,) or float
  236. Point at which to estimate the derivatives. Float will be converted
  237. to a 1-D array.
  238. method : {'3-point', '2-point', 'cs'}, optional
  239. Finite difference method to use:
  240. - '2-point' - use the first order accuracy forward or backward
  241. difference.
  242. - '3-point' - use central difference in interior points and the
  243. second order accuracy forward or backward difference
  244. near the boundary.
  245. - 'cs' - use a complex-step finite difference scheme. This assumes
  246. that the user function is real-valued and can be
  247. analytically continued to the complex plane. Otherwise,
  248. produces bogus results.
  249. rel_step : None or array_like, optional
  250. Relative step size to use. If None (default) the absolute step size is
  251. computed as ``h = rel_step * sign(x0) * max(1, abs(x0))``, with
  252. `rel_step` being selected automatically, see Notes. Otherwise
  253. ``h = rel_step * sign(x0) * abs(x0)``. For ``method='3-point'`` the
  254. sign of `h` is ignored. The calculated step size is possibly adjusted
  255. to fit into the bounds.
  256. abs_step : array_like, optional
  257. Absolute step size to use, possibly adjusted to fit into the bounds.
  258. For ``method='3-point'`` the sign of `abs_step` is ignored. By default
  259. relative steps are used, only if ``abs_step is not None`` are absolute
  260. steps used.
  261. f0 : None or array_like, optional
  262. If not None it is assumed to be equal to ``fun(x0)``, in this case
  263. the ``fun(x0)`` is not called. Default is None.
  264. bounds : tuple of array_like, optional
  265. Lower and upper bounds on independent variables. Defaults to no bounds.
  266. Each bound must match the size of `x0` or be a scalar, in the latter
  267. case the bound will be the same for all variables. Use it to limit the
  268. range of function evaluation. Bounds checking is not implemented
  269. when `as_linear_operator` is True.
  270. sparsity : {None, array_like, sparse matrix, 2-tuple}, optional
  271. Defines a sparsity structure of the Jacobian matrix. If the Jacobian
  272. matrix is known to have only few non-zero elements in each row, then
  273. it's possible to estimate its several columns by a single function
  274. evaluation [3]_. To perform such economic computations two ingredients
  275. are required:
  276. * structure : array_like or sparse matrix of shape (m, n). A zero
  277. element means that a corresponding element of the Jacobian
  278. identically equals to zero.
  279. * groups : array_like of shape (n,). A column grouping for a given
  280. sparsity structure, use `group_columns` to obtain it.
  281. A single array or a sparse matrix is interpreted as a sparsity
  282. structure, and groups are computed inside the function. A tuple is
  283. interpreted as (structure, groups). If None (default), a standard
  284. dense differencing will be used.
  285. Note, that sparse differencing makes sense only for large Jacobian
  286. matrices where each row contains few non-zero elements.
  287. as_linear_operator : bool, optional
  288. When True the function returns an `scipy.sparse.linalg.LinearOperator`.
  289. Otherwise it returns a dense array or a sparse matrix depending on
  290. `sparsity`. The linear operator provides an efficient way of computing
  291. ``J.dot(p)`` for any vector ``p`` of shape (n,), but does not allow
  292. direct access to individual elements of the matrix. By default
  293. `as_linear_operator` is False.
  294. args, kwargs : tuple and dict, optional
  295. Additional arguments passed to `fun`. Both empty by default.
  296. The calling signature is ``fun(x, *args, **kwargs)``.
  297. Returns
  298. -------
  299. J : {ndarray, sparse matrix, LinearOperator}
  300. Finite difference approximation of the Jacobian matrix.
  301. If `as_linear_operator` is True returns a LinearOperator
  302. with shape (m, n). Otherwise it returns a dense array or sparse
  303. matrix depending on how `sparsity` is defined. If `sparsity`
  304. is None then a ndarray with shape (m, n) is returned. If
  305. `sparsity` is not None returns a csr_matrix with shape (m, n).
  306. For sparse matrices and linear operators it is always returned as
  307. a 2-D structure, for ndarrays, if m=1 it is returned
  308. as a 1-D gradient array with shape (n,).
  309. See Also
  310. --------
  311. check_derivative : Check correctness of a function computing derivatives.
  312. Notes
  313. -----
  314. If `rel_step` is not provided, it assigned as ``EPS**(1/s)``, where EPS is
  315. determined from the smallest floating point dtype of `x0` or `fun(x0)`,
  316. ``np.finfo(x0.dtype).eps``, s=2 for '2-point' method and
  317. s=3 for '3-point' method. Such relative step approximately minimizes a sum
  318. of truncation and round-off errors, see [1]_. Relative steps are used by
  319. default. However, absolute steps are used when ``abs_step is not None``.
  320. If any of the absolute or relative steps produces an indistinguishable
  321. difference from the original `x0`, ``(x0 + dx) - x0 == 0``, then a
  322. automatic step size is substituted for that particular entry.
  323. A finite difference scheme for '3-point' method is selected automatically.
  324. The well-known central difference scheme is used for points sufficiently
  325. far from the boundary, and 3-point forward or backward scheme is used for
  326. points near the boundary. Both schemes have the second-order accuracy in
  327. terms of Taylor expansion. Refer to [2]_ for the formulas of 3-point
  328. forward and backward difference schemes.
  329. For dense differencing when m=1 Jacobian is returned with a shape (n,),
  330. on the other hand when n=1 Jacobian is returned with a shape (m, 1).
  331. Our motivation is the following: a) It handles a case of gradient
  332. computation (m=1) in a conventional way. b) It clearly separates these two
  333. different cases. b) In all cases np.atleast_2d can be called to get 2-D
  334. Jacobian with correct dimensions.
  335. References
  336. ----------
  337. .. [1] W. H. Press et. al. "Numerical Recipes. The Art of Scientific
  338. Computing. 3rd edition", sec. 5.7.
  339. .. [2] A. Curtis, M. J. D. Powell, and J. Reid, "On the estimation of
  340. sparse Jacobian matrices", Journal of the Institute of Mathematics
  341. and its Applications, 13 (1974), pp. 117-120.
  342. .. [3] B. Fornberg, "Generation of Finite Difference Formulas on
  343. Arbitrarily Spaced Grids", Mathematics of Computation 51, 1988.
  344. Examples
  345. --------
  346. >>> import numpy as np
  347. >>> from scipy.optimize._numdiff import approx_derivative
  348. >>>
  349. >>> def f(x, c1, c2):
  350. ... return np.array([x[0] * np.sin(c1 * x[1]),
  351. ... x[0] * np.cos(c2 * x[1])])
  352. ...
  353. >>> x0 = np.array([1.0, 0.5 * np.pi])
  354. >>> approx_derivative(f, x0, args=(1, 2))
  355. array([[ 1., 0.],
  356. [-1., 0.]])
  357. Bounds can be used to limit the region of function evaluation.
  358. In the example below we compute left and right derivative at point 1.0.
  359. >>> def g(x):
  360. ... return x**2 if x >= 1 else x
  361. ...
  362. >>> x0 = 1.0
  363. >>> approx_derivative(g, x0, bounds=(-np.inf, 1.0))
  364. array([ 1.])
  365. >>> approx_derivative(g, x0, bounds=(1.0, np.inf))
  366. array([ 2.])
  367. """
  368. if method not in ['2-point', '3-point', 'cs']:
  369. raise ValueError("Unknown method '%s'. " % method)
  370. x0 = np.atleast_1d(x0)
  371. if x0.ndim > 1:
  372. raise ValueError("`x0` must have at most 1 dimension.")
  373. lb, ub = _prepare_bounds(bounds, x0)
  374. if lb.shape != x0.shape or ub.shape != x0.shape:
  375. raise ValueError("Inconsistent shapes between bounds and `x0`.")
  376. if as_linear_operator and not (np.all(np.isinf(lb))
  377. and np.all(np.isinf(ub))):
  378. raise ValueError("Bounds not supported when "
  379. "`as_linear_operator` is True.")
  380. def fun_wrapped(x):
  381. f = np.atleast_1d(fun(x, *args, **kwargs))
  382. if f.ndim > 1:
  383. raise RuntimeError("`fun` return value has "
  384. "more than 1 dimension.")
  385. return f
  386. if f0 is None:
  387. f0 = fun_wrapped(x0)
  388. else:
  389. f0 = np.atleast_1d(f0)
  390. if f0.ndim > 1:
  391. raise ValueError("`f0` passed has more than 1 dimension.")
  392. if np.any((x0 < lb) | (x0 > ub)):
  393. raise ValueError("`x0` violates bound constraints.")
  394. if as_linear_operator:
  395. if rel_step is None:
  396. rel_step = _eps_for_method(x0.dtype, f0.dtype, method)
  397. return _linear_operator_difference(fun_wrapped, x0,
  398. f0, rel_step, method)
  399. else:
  400. # by default we use rel_step
  401. if abs_step is None:
  402. h = _compute_absolute_step(rel_step, x0, f0, method)
  403. else:
  404. # user specifies an absolute step
  405. sign_x0 = (x0 >= 0).astype(float) * 2 - 1
  406. h = abs_step
  407. # cannot have a zero step. This might happen if x0 is very large
  408. # or small. In which case fall back to relative step.
  409. dx = ((x0 + h) - x0)
  410. h = np.where(dx == 0,
  411. _eps_for_method(x0.dtype, f0.dtype, method) *
  412. sign_x0 * np.maximum(1.0, np.abs(x0)),
  413. h)
  414. if method == '2-point':
  415. h, use_one_sided = _adjust_scheme_to_bounds(
  416. x0, h, 1, '1-sided', lb, ub)
  417. elif method == '3-point':
  418. h, use_one_sided = _adjust_scheme_to_bounds(
  419. x0, h, 1, '2-sided', lb, ub)
  420. elif method == 'cs':
  421. use_one_sided = False
  422. if sparsity is None:
  423. return _dense_difference(fun_wrapped, x0, f0, h,
  424. use_one_sided, method)
  425. else:
  426. if not issparse(sparsity) and len(sparsity) == 2:
  427. structure, groups = sparsity
  428. else:
  429. structure = sparsity
  430. groups = group_columns(sparsity)
  431. if issparse(structure):
  432. structure = csc_matrix(structure)
  433. else:
  434. structure = np.atleast_2d(structure)
  435. groups = np.atleast_1d(groups)
  436. return _sparse_difference(fun_wrapped, x0, f0, h,
  437. use_one_sided, structure,
  438. groups, method)
  439. def _linear_operator_difference(fun, x0, f0, h, method):
  440. m = f0.size
  441. n = x0.size
  442. if method == '2-point':
  443. def matvec(p):
  444. if np.array_equal(p, np.zeros_like(p)):
  445. return np.zeros(m)
  446. dx = h / norm(p)
  447. x = x0 + dx*p
  448. df = fun(x) - f0
  449. return df / dx
  450. elif method == '3-point':
  451. def matvec(p):
  452. if np.array_equal(p, np.zeros_like(p)):
  453. return np.zeros(m)
  454. dx = 2*h / norm(p)
  455. x1 = x0 - (dx/2)*p
  456. x2 = x0 + (dx/2)*p
  457. f1 = fun(x1)
  458. f2 = fun(x2)
  459. df = f2 - f1
  460. return df / dx
  461. elif method == 'cs':
  462. def matvec(p):
  463. if np.array_equal(p, np.zeros_like(p)):
  464. return np.zeros(m)
  465. dx = h / norm(p)
  466. x = x0 + dx*p*1.j
  467. f1 = fun(x)
  468. df = f1.imag
  469. return df / dx
  470. else:
  471. raise RuntimeError("Never be here.")
  472. return LinearOperator((m, n), matvec)
  473. def _dense_difference(fun, x0, f0, h, use_one_sided, method):
  474. m = f0.size
  475. n = x0.size
  476. J_transposed = np.empty((n, m))
  477. h_vecs = np.diag(h)
  478. for i in range(h.size):
  479. if method == '2-point':
  480. x = x0 + h_vecs[i]
  481. dx = x[i] - x0[i] # Recompute dx as exactly representable number.
  482. df = fun(x) - f0
  483. elif method == '3-point' and use_one_sided[i]:
  484. x1 = x0 + h_vecs[i]
  485. x2 = x0 + 2 * h_vecs[i]
  486. dx = x2[i] - x0[i]
  487. f1 = fun(x1)
  488. f2 = fun(x2)
  489. df = -3.0 * f0 + 4 * f1 - f2
  490. elif method == '3-point' and not use_one_sided[i]:
  491. x1 = x0 - h_vecs[i]
  492. x2 = x0 + h_vecs[i]
  493. dx = x2[i] - x1[i]
  494. f1 = fun(x1)
  495. f2 = fun(x2)
  496. df = f2 - f1
  497. elif method == 'cs':
  498. f1 = fun(x0 + h_vecs[i]*1.j)
  499. df = f1.imag
  500. dx = h_vecs[i, i]
  501. else:
  502. raise RuntimeError("Never be here.")
  503. J_transposed[i] = df / dx
  504. if m == 1:
  505. J_transposed = np.ravel(J_transposed)
  506. return J_transposed.T
  507. def _sparse_difference(fun, x0, f0, h, use_one_sided,
  508. structure, groups, method):
  509. m = f0.size
  510. n = x0.size
  511. row_indices = []
  512. col_indices = []
  513. fractions = []
  514. n_groups = np.max(groups) + 1
  515. for group in range(n_groups):
  516. # Perturb variables which are in the same group simultaneously.
  517. e = np.equal(group, groups)
  518. h_vec = h * e
  519. if method == '2-point':
  520. x = x0 + h_vec
  521. dx = x - x0
  522. df = fun(x) - f0
  523. # The result is written to columns which correspond to perturbed
  524. # variables.
  525. cols, = np.nonzero(e)
  526. # Find all non-zero elements in selected columns of Jacobian.
  527. i, j, _ = find(structure[:, cols])
  528. # Restore column indices in the full array.
  529. j = cols[j]
  530. elif method == '3-point':
  531. # Here we do conceptually the same but separate one-sided
  532. # and two-sided schemes.
  533. x1 = x0.copy()
  534. x2 = x0.copy()
  535. mask_1 = use_one_sided & e
  536. x1[mask_1] += h_vec[mask_1]
  537. x2[mask_1] += 2 * h_vec[mask_1]
  538. mask_2 = ~use_one_sided & e
  539. x1[mask_2] -= h_vec[mask_2]
  540. x2[mask_2] += h_vec[mask_2]
  541. dx = np.zeros(n)
  542. dx[mask_1] = x2[mask_1] - x0[mask_1]
  543. dx[mask_2] = x2[mask_2] - x1[mask_2]
  544. f1 = fun(x1)
  545. f2 = fun(x2)
  546. cols, = np.nonzero(e)
  547. i, j, _ = find(structure[:, cols])
  548. j = cols[j]
  549. mask = use_one_sided[j]
  550. df = np.empty(m)
  551. rows = i[mask]
  552. df[rows] = -3 * f0[rows] + 4 * f1[rows] - f2[rows]
  553. rows = i[~mask]
  554. df[rows] = f2[rows] - f1[rows]
  555. elif method == 'cs':
  556. f1 = fun(x0 + h_vec*1.j)
  557. df = f1.imag
  558. dx = h_vec
  559. cols, = np.nonzero(e)
  560. i, j, _ = find(structure[:, cols])
  561. j = cols[j]
  562. else:
  563. raise ValueError("Never be here.")
  564. # All that's left is to compute the fraction. We store i, j and
  565. # fractions as separate arrays and later construct coo_matrix.
  566. row_indices.append(i)
  567. col_indices.append(j)
  568. fractions.append(df[i] / dx[j])
  569. row_indices = np.hstack(row_indices)
  570. col_indices = np.hstack(col_indices)
  571. fractions = np.hstack(fractions)
  572. J = coo_matrix((fractions, (row_indices, col_indices)), shape=(m, n))
  573. return csr_matrix(J)
  574. def check_derivative(fun, jac, x0, bounds=(-np.inf, np.inf), args=(),
  575. kwargs={}):
  576. """Check correctness of a function computing derivatives (Jacobian or
  577. gradient) by comparison with a finite difference approximation.
  578. Parameters
  579. ----------
  580. fun : callable
  581. Function of which to estimate the derivatives. The argument x
  582. passed to this function is ndarray of shape (n,) (never a scalar
  583. even if n=1). It must return 1-D array_like of shape (m,) or a scalar.
  584. jac : callable
  585. Function which computes Jacobian matrix of `fun`. It must work with
  586. argument x the same way as `fun`. The return value must be array_like
  587. or sparse matrix with an appropriate shape.
  588. x0 : array_like of shape (n,) or float
  589. Point at which to estimate the derivatives. Float will be converted
  590. to 1-D array.
  591. bounds : 2-tuple of array_like, optional
  592. Lower and upper bounds on independent variables. Defaults to no bounds.
  593. Each bound must match the size of `x0` or be a scalar, in the latter
  594. case the bound will be the same for all variables. Use it to limit the
  595. range of function evaluation.
  596. args, kwargs : tuple and dict, optional
  597. Additional arguments passed to `fun` and `jac`. Both empty by default.
  598. The calling signature is ``fun(x, *args, **kwargs)`` and the same
  599. for `jac`.
  600. Returns
  601. -------
  602. accuracy : float
  603. The maximum among all relative errors for elements with absolute values
  604. higher than 1 and absolute errors for elements with absolute values
  605. less or equal than 1. If `accuracy` is on the order of 1e-6 or lower,
  606. then it is likely that your `jac` implementation is correct.
  607. See Also
  608. --------
  609. approx_derivative : Compute finite difference approximation of derivative.
  610. Examples
  611. --------
  612. >>> import numpy as np
  613. >>> from scipy.optimize._numdiff import check_derivative
  614. >>>
  615. >>>
  616. >>> def f(x, c1, c2):
  617. ... return np.array([x[0] * np.sin(c1 * x[1]),
  618. ... x[0] * np.cos(c2 * x[1])])
  619. ...
  620. >>> def jac(x, c1, c2):
  621. ... return np.array([
  622. ... [np.sin(c1 * x[1]), c1 * x[0] * np.cos(c1 * x[1])],
  623. ... [np.cos(c2 * x[1]), -c2 * x[0] * np.sin(c2 * x[1])]
  624. ... ])
  625. ...
  626. >>>
  627. >>> x0 = np.array([1.0, 0.5 * np.pi])
  628. >>> check_derivative(f, jac, x0, args=(1, 2))
  629. 2.4492935982947064e-16
  630. """
  631. J_to_test = jac(x0, *args, **kwargs)
  632. if issparse(J_to_test):
  633. J_diff = approx_derivative(fun, x0, bounds=bounds, sparsity=J_to_test,
  634. args=args, kwargs=kwargs)
  635. J_to_test = csr_matrix(J_to_test)
  636. abs_err = J_to_test - J_diff
  637. i, j, abs_err_data = find(abs_err)
  638. J_diff_data = np.asarray(J_diff[i, j]).ravel()
  639. return np.max(np.abs(abs_err_data) /
  640. np.maximum(1, np.abs(J_diff_data)))
  641. else:
  642. J_diff = approx_derivative(fun, x0, bounds=bounds,
  643. args=args, kwargs=kwargs)
  644. abs_err = np.abs(J_to_test - J_diff)
  645. return np.max(abs_err / np.maximum(1, np.abs(J_diff)))