_numdiff.py 23 KB

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