projections.py 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406
  1. """Basic linear factorizations needed by the solver."""
  2. from __future__ import division, print_function, absolute_import
  3. from scipy.sparse import (bmat, csc_matrix, eye, issparse)
  4. from scipy.sparse.linalg import LinearOperator
  5. import scipy.linalg
  6. import scipy.sparse.linalg
  7. try:
  8. from sksparse.cholmod import cholesky_AAt
  9. sksparse_available = True
  10. except ImportError:
  11. import warnings
  12. sksparse_available = False
  13. import numpy as np
  14. from warnings import warn
  15. __all__ = [
  16. 'orthogonality',
  17. 'projections',
  18. ]
  19. def orthogonality(A, g):
  20. """Measure orthogonality between a vector and the null space of a matrix.
  21. Compute a measure of orthogonality between the null space
  22. of the (possibly sparse) matrix ``A`` and a given vector ``g``.
  23. The formula is a simplified (and cheaper) version of formula (3.13)
  24. from [1]_.
  25. ``orth = norm(A g, ord=2)/(norm(A, ord='fro')*norm(g, ord=2))``.
  26. References
  27. ----------
  28. .. [1] Gould, Nicholas IM, Mary E. Hribar, and Jorge Nocedal.
  29. "On the solution of equality constrained quadratic
  30. programming problems arising in optimization."
  31. SIAM Journal on Scientific Computing 23.4 (2001): 1376-1395.
  32. """
  33. # Compute vector norms
  34. norm_g = np.linalg.norm(g)
  35. # Compute Frobenius norm of the matrix A
  36. if issparse(A):
  37. norm_A = scipy.sparse.linalg.norm(A, ord='fro')
  38. else:
  39. norm_A = np.linalg.norm(A, ord='fro')
  40. # Check if norms are zero
  41. if norm_g == 0 or norm_A == 0:
  42. return 0
  43. norm_A_g = np.linalg.norm(A.dot(g))
  44. # Orthogonality measure
  45. orth = norm_A_g / (norm_A*norm_g)
  46. return orth
  47. def normal_equation_projections(A, m, n, orth_tol, max_refin, tol):
  48. """Return linear operators for matrix A using ``NormalEquation`` approach.
  49. """
  50. # Cholesky factorization
  51. factor = cholesky_AAt(A)
  52. # z = x - A.T inv(A A.T) A x
  53. def null_space(x):
  54. v = factor(A.dot(x))
  55. z = x - A.T.dot(v)
  56. # Iterative refinement to improve roundoff
  57. # errors described in [2]_, algorithm 5.1.
  58. k = 0
  59. while orthogonality(A, z) > orth_tol:
  60. if k >= max_refin:
  61. break
  62. # z_next = z - A.T inv(A A.T) A z
  63. v = factor(A.dot(z))
  64. z = z - A.T.dot(v)
  65. k += 1
  66. return z
  67. # z = inv(A A.T) A x
  68. def least_squares(x):
  69. return factor(A.dot(x))
  70. # z = A.T inv(A A.T) x
  71. def row_space(x):
  72. return A.T.dot(factor(x))
  73. return null_space, least_squares, row_space
  74. def augmented_system_projections(A, m, n, orth_tol, max_refin, tol):
  75. """Return linear operators for matrix A - ``AugmentedSystem``."""
  76. # Form augmented system
  77. K = csc_matrix(bmat([[eye(n), A.T], [A, None]]))
  78. # LU factorization
  79. # TODO: Use a symmetric indefinite factorization
  80. # to solve the system twice as fast (because
  81. # of the symmetry).
  82. try:
  83. solve = scipy.sparse.linalg.factorized(K)
  84. except RuntimeError:
  85. warn("Singular Jacobian matrix. Using dense SVD decomposition to "
  86. "perform the factorizations.")
  87. return svd_factorization_projections(A.toarray(),
  88. m, n, orth_tol,
  89. max_refin, tol)
  90. # z = x - A.T inv(A A.T) A x
  91. # is computed solving the extended system:
  92. # [I A.T] * [ z ] = [x]
  93. # [A O ] [aux] [0]
  94. def null_space(x):
  95. # v = [x]
  96. # [0]
  97. v = np.hstack([x, np.zeros(m)])
  98. # lu_sol = [ z ]
  99. # [aux]
  100. lu_sol = solve(v)
  101. z = lu_sol[:n]
  102. # Iterative refinement to improve roundoff
  103. # errors described in [2]_, algorithm 5.2.
  104. k = 0
  105. while orthogonality(A, z) > orth_tol:
  106. if k >= max_refin:
  107. break
  108. # new_v = [x] - [I A.T] * [ z ]
  109. # [0] [A O ] [aux]
  110. new_v = v - K.dot(lu_sol)
  111. # [I A.T] * [delta z ] = new_v
  112. # [A O ] [delta aux]
  113. lu_update = solve(new_v)
  114. # [ z ] += [delta z ]
  115. # [aux] [delta aux]
  116. lu_sol += lu_update
  117. z = lu_sol[:n]
  118. k += 1
  119. # return z = x - A.T inv(A A.T) A x
  120. return z
  121. # z = inv(A A.T) A x
  122. # is computed solving the extended system:
  123. # [I A.T] * [aux] = [x]
  124. # [A O ] [ z ] [0]
  125. def least_squares(x):
  126. # v = [x]
  127. # [0]
  128. v = np.hstack([x, np.zeros(m)])
  129. # lu_sol = [aux]
  130. # [ z ]
  131. lu_sol = solve(v)
  132. # return z = inv(A A.T) A x
  133. return lu_sol[n:m+n]
  134. # z = A.T inv(A A.T) x
  135. # is computed solving the extended system:
  136. # [I A.T] * [ z ] = [0]
  137. # [A O ] [aux] [x]
  138. def row_space(x):
  139. # v = [0]
  140. # [x]
  141. v = np.hstack([np.zeros(n), x])
  142. # lu_sol = [ z ]
  143. # [aux]
  144. lu_sol = solve(v)
  145. # return z = A.T inv(A A.T) x
  146. return lu_sol[:n]
  147. return null_space, least_squares, row_space
  148. def qr_factorization_projections(A, m, n, orth_tol, max_refin, tol):
  149. """Return linear operators for matrix A using ``QRFactorization`` approach.
  150. """
  151. # QRFactorization
  152. Q, R, P = scipy.linalg.qr(A.T, pivoting=True, mode='economic')
  153. if np.linalg.norm(R[-1, :], np.inf) < tol:
  154. warn('Singular Jacobian matrix. Using SVD decomposition to ' +
  155. 'perform the factorizations.')
  156. return svd_factorization_projections(A, m, n,
  157. orth_tol,
  158. max_refin,
  159. tol)
  160. # z = x - A.T inv(A A.T) A x
  161. def null_space(x):
  162. # v = P inv(R) Q.T x
  163. aux1 = Q.T.dot(x)
  164. aux2 = scipy.linalg.solve_triangular(R, aux1, lower=False)
  165. v = np.zeros(m)
  166. v[P] = aux2
  167. z = x - A.T.dot(v)
  168. # Iterative refinement to improve roundoff
  169. # errors described in [2]_, algorithm 5.1.
  170. k = 0
  171. while orthogonality(A, z) > orth_tol:
  172. if k >= max_refin:
  173. break
  174. # v = P inv(R) Q.T x
  175. aux1 = Q.T.dot(z)
  176. aux2 = scipy.linalg.solve_triangular(R, aux1, lower=False)
  177. v[P] = aux2
  178. # z_next = z - A.T v
  179. z = z - A.T.dot(v)
  180. k += 1
  181. return z
  182. # z = inv(A A.T) A x
  183. def least_squares(x):
  184. # z = P inv(R) Q.T x
  185. aux1 = Q.T.dot(x)
  186. aux2 = scipy.linalg.solve_triangular(R, aux1, lower=False)
  187. z = np.zeros(m)
  188. z[P] = aux2
  189. return z
  190. # z = A.T inv(A A.T) x
  191. def row_space(x):
  192. # z = Q inv(R.T) P.T x
  193. aux1 = x[P]
  194. aux2 = scipy.linalg.solve_triangular(R, aux1,
  195. lower=False,
  196. trans='T')
  197. z = Q.dot(aux2)
  198. return z
  199. return null_space, least_squares, row_space
  200. def svd_factorization_projections(A, m, n, orth_tol, max_refin, tol):
  201. """Return linear operators for matrix A using ``SVDFactorization`` approach.
  202. """
  203. # SVD Factorization
  204. U, s, Vt = scipy.linalg.svd(A, full_matrices=False)
  205. # Remove dimensions related with very small singular values
  206. U = U[:, s > tol]
  207. Vt = Vt[s > tol, :]
  208. s = s[s > tol]
  209. # z = x - A.T inv(A A.T) A x
  210. def null_space(x):
  211. # v = U 1/s V.T x = inv(A A.T) A x
  212. aux1 = Vt.dot(x)
  213. aux2 = 1/s*aux1
  214. v = U.dot(aux2)
  215. z = x - A.T.dot(v)
  216. # Iterative refinement to improve roundoff
  217. # errors described in [2]_, algorithm 5.1.
  218. k = 0
  219. while orthogonality(A, z) > orth_tol:
  220. if k >= max_refin:
  221. break
  222. # v = U 1/s V.T x = inv(A A.T) A x
  223. aux1 = Vt.dot(z)
  224. aux2 = 1/s*aux1
  225. v = U.dot(aux2)
  226. # z_next = z - A.T v
  227. z = z - A.T.dot(v)
  228. k += 1
  229. return z
  230. # z = inv(A A.T) A x
  231. def least_squares(x):
  232. # z = U 1/s V.T x = inv(A A.T) A x
  233. aux1 = Vt.dot(x)
  234. aux2 = 1/s*aux1
  235. z = U.dot(aux2)
  236. return z
  237. # z = A.T inv(A A.T) x
  238. def row_space(x):
  239. # z = V 1/s U.T x
  240. aux1 = U.T.dot(x)
  241. aux2 = 1/s*aux1
  242. z = Vt.T.dot(aux2)
  243. return z
  244. return null_space, least_squares, row_space
  245. def projections(A, method=None, orth_tol=1e-12, max_refin=3, tol=1e-15):
  246. """Return three linear operators related with a given matrix A.
  247. Parameters
  248. ----------
  249. A : sparse matrix (or ndarray), shape (m, n)
  250. Matrix ``A`` used in the projection.
  251. method : string, optional
  252. Method used for compute the given linear
  253. operators. Should be one of:
  254. - 'NormalEquation': The operators
  255. will be computed using the
  256. so-called normal equation approach
  257. explained in [1]_. In order to do
  258. so the Cholesky factorization of
  259. ``(A A.T)`` is computed. Exclusive
  260. for sparse matrices.
  261. - 'AugmentedSystem': The operators
  262. will be computed using the
  263. so-called augmented system approach
  264. explained in [1]_. Exclusive
  265. for sparse matrices.
  266. - 'QRFactorization': Compute projections
  267. using QR factorization. Exclusive for
  268. dense matrices.
  269. - 'SVDFactorization': Compute projections
  270. using SVD factorization. Exclusive for
  271. dense matrices.
  272. orth_tol : float, optional
  273. Tolerance for iterative refinements.
  274. max_refin : int, optional
  275. Maximum number of iterative refinements
  276. tol : float, optional
  277. Tolerance for singular values
  278. Returns
  279. -------
  280. Z : LinearOperator, shape (n, n)
  281. Null-space operator. For a given vector ``x``,
  282. the null space operator is equivalent to apply
  283. a projection matrix ``P = I - A.T inv(A A.T) A``
  284. to the vector. It can be shown that this is
  285. equivalent to project ``x`` into the null space
  286. of A.
  287. LS : LinearOperator, shape (m, n)
  288. Least-Square operator. For a given vector ``x``,
  289. the least-square operator is equivalent to apply a
  290. pseudoinverse matrix ``pinv(A.T) = inv(A A.T) A``
  291. to the vector. It can be shown that this vector
  292. ``pinv(A.T) x`` is the least_square solution to
  293. ``A.T y = x``.
  294. Y : LinearOperator, shape (n, m)
  295. Row-space operator. For a given vector ``x``,
  296. the row-space operator is equivalent to apply a
  297. projection matrix ``Q = A.T inv(A A.T)``
  298. to the vector. It can be shown that this
  299. vector ``y = Q x`` the minimum norm solution
  300. of ``A y = x``.
  301. Notes
  302. -----
  303. Uses iterative refinements described in [1]
  304. during the computation of ``Z`` in order to
  305. cope with the possibility of large roundoff errors.
  306. References
  307. ----------
  308. .. [1] Gould, Nicholas IM, Mary E. Hribar, and Jorge Nocedal.
  309. "On the solution of equality constrained quadratic
  310. programming problems arising in optimization."
  311. SIAM Journal on Scientific Computing 23.4 (2001): 1376-1395.
  312. """
  313. m, n = np.shape(A)
  314. # The factorization of an empty matrix
  315. # only works for the sparse representation.
  316. if m*n == 0:
  317. A = csc_matrix(A)
  318. # Check Argument
  319. if issparse(A):
  320. if method is None:
  321. method = "AugmentedSystem"
  322. if method not in ("NormalEquation", "AugmentedSystem"):
  323. raise ValueError("Method not allowed for sparse matrix.")
  324. if method == "NormalEquation" and not sksparse_available:
  325. warnings.warn(("Only accepts 'NormalEquation' option when"
  326. " scikit-sparse is available. Using "
  327. "'AugmentedSystem' option instead."),
  328. ImportWarning)
  329. method = 'AugmentedSystem'
  330. else:
  331. if method is None:
  332. method = "QRFactorization"
  333. if method not in ("QRFactorization", "SVDFactorization"):
  334. raise ValueError("Method not allowed for dense array.")
  335. if method == 'NormalEquation':
  336. null_space, least_squares, row_space \
  337. = normal_equation_projections(A, m, n, orth_tol, max_refin, tol)
  338. elif method == 'AugmentedSystem':
  339. null_space, least_squares, row_space \
  340. = augmented_system_projections(A, m, n, orth_tol, max_refin, tol)
  341. elif method == "QRFactorization":
  342. null_space, least_squares, row_space \
  343. = qr_factorization_projections(A, m, n, orth_tol, max_refin, tol)
  344. elif method == "SVDFactorization":
  345. null_space, least_squares, row_space \
  346. = svd_factorization_projections(A, m, n, orth_tol, max_refin, tol)
  347. Z = LinearOperator((n, n), null_space)
  348. LS = LinearOperator((m, n), least_squares)
  349. Y = LinearOperator((n, m), row_space)
  350. return Z, LS, Y