projections.py 13 KB

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