projections.py 13 KB

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