_lobpcg.py 43 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173
  1. # mypy: allow-untyped-defs
  2. """Locally Optimal Block Preconditioned Conjugate Gradient methods."""
  3. # Author: Pearu Peterson
  4. # Created: February 2020
  5. import torch
  6. from torch import _linalg_utils as _utils, Tensor
  7. from torch.overrides import handle_torch_function, has_torch_function
  8. __all__ = ["lobpcg"]
  9. def _symeig_backward_complete_eigenspace(D_grad, U_grad, A, D, U):
  10. # compute F, such that F_ij = (d_j - d_i)^{-1} for i != j, F_ii = 0
  11. F = D.unsqueeze(-2) - D.unsqueeze(-1)
  12. F.diagonal(dim1=-2, dim2=-1).fill_(float("inf"))
  13. F.pow_(-1)
  14. # A.grad = U (D.grad + (U^T U.grad * F)) U^T
  15. Ut = U.mT.contiguous()
  16. res = torch.matmul(
  17. U, torch.matmul(torch.diag_embed(D_grad) + torch.matmul(Ut, U_grad) * F, Ut)
  18. )
  19. return res
  20. def _polynomial_coefficients_given_roots(roots):
  21. """
  22. Given the `roots` of a polynomial, find the polynomial's coefficients.
  23. If roots = (r_1, ..., r_n), then the method returns
  24. coefficients (a_0, a_1, ..., a_n (== 1)) so that
  25. p(x) = (x - r_1) * ... * (x - r_n)
  26. = x^n + a_{n-1} * x^{n-1} + ... a_1 * x_1 + a_0
  27. Note: for better performance requires writing a low-level kernel
  28. """
  29. poly_order = roots.shape[-1]
  30. poly_coeffs_shape = list(roots.shape)
  31. # we assume p(x) = x^n + a_{n-1} * x^{n-1} + ... + a_1 * x + a_0,
  32. # so poly_coeffs = {a_0, ..., a_n, a_{n+1}(== 1)},
  33. # but we insert one extra coefficient to enable better vectorization below
  34. poly_coeffs_shape[-1] += 2
  35. poly_coeffs = roots.new_zeros(poly_coeffs_shape)
  36. poly_coeffs[..., 0] = 1
  37. poly_coeffs[..., -1] = 1
  38. # perform the Horner's rule
  39. for i in range(1, poly_order + 1):
  40. # note that it is computationally hard to compute backward for this method,
  41. # because then given the coefficients it would require finding the roots and/or
  42. # calculating the sensitivity based on the Vieta's theorem.
  43. # So the code below tries to circumvent the explicit root finding by series
  44. # of operations on memory copies imitating the Horner's method.
  45. # The memory copies are required to construct nodes in the computational graph
  46. # by exploiting the explicit (not in-place, separate node for each step)
  47. # recursion of the Horner's method.
  48. # Needs more memory, O(... * k^2), but with only O(... * k^2) complexity.
  49. poly_coeffs_new = poly_coeffs.clone() if roots.requires_grad else poly_coeffs
  50. out = poly_coeffs_new.narrow(-1, poly_order - i, i + 1)
  51. out -= roots.narrow(-1, i - 1, 1) * poly_coeffs.narrow(
  52. -1, poly_order - i + 1, i + 1
  53. )
  54. poly_coeffs = poly_coeffs_new
  55. return poly_coeffs.narrow(-1, 1, poly_order + 1)
  56. def _polynomial_value(poly, x, zero_power, transition):
  57. """
  58. A generic method for computing poly(x) using the Horner's rule.
  59. Args:
  60. poly (Tensor): the (possibly batched) 1D Tensor representing
  61. polynomial coefficients such that
  62. poly[..., i] = (a_{i_0}, ..., a{i_n} (==1)), and
  63. poly(x) = poly[..., 0] * zero_power + ... + poly[..., n] * x^n
  64. x (Tensor): the value (possible batched) to evaluate the polynomial `poly` at.
  65. zero_power (Tensor): the representation of `x^0`. It is application-specific.
  66. transition (Callable): the function that accepts some intermediate result `int_val`,
  67. the `x` and a specific polynomial coefficient
  68. `poly[..., k]` for some iteration `k`.
  69. It basically performs one iteration of the Horner's rule
  70. defined as `x * int_val + poly[..., k] * zero_power`.
  71. Note that `zero_power` is not a parameter,
  72. because the step `+ poly[..., k] * zero_power` depends on `x`,
  73. whether it is a vector, a matrix, or something else, so this
  74. functionality is delegated to the user.
  75. """
  76. res = zero_power.clone()
  77. for k in range(poly.size(-1) - 2, -1, -1):
  78. res = transition(res, x, poly[..., k])
  79. return res
  80. def _matrix_polynomial_value(poly, x, zero_power=None):
  81. """
  82. Evaluates `poly(x)` for the (batched) matrix input `x`.
  83. Check out `_polynomial_value` function for more details.
  84. """
  85. # matrix-aware Horner's rule iteration
  86. def transition(curr_poly_val, x, poly_coeff):
  87. res = x.matmul(curr_poly_val)
  88. res.diagonal(dim1=-2, dim2=-1).add_(poly_coeff.unsqueeze(-1))
  89. return res
  90. if zero_power is None:
  91. zero_power = torch.eye(
  92. x.size(-1), x.size(-1), dtype=x.dtype, device=x.device
  93. ).view(*([1] * len(list(x.shape[:-2]))), x.size(-1), x.size(-1))
  94. return _polynomial_value(poly, x, zero_power, transition)
  95. def _vector_polynomial_value(poly, x, zero_power=None):
  96. """
  97. Evaluates `poly(x)` for the (batched) vector input `x`.
  98. Check out `_polynomial_value` function for more details.
  99. """
  100. # vector-aware Horner's rule iteration
  101. def transition(curr_poly_val, x, poly_coeff):
  102. res = torch.addcmul(poly_coeff.unsqueeze(-1), x, curr_poly_val)
  103. return res
  104. if zero_power is None:
  105. zero_power = x.new_ones(1).expand(x.shape)
  106. return _polynomial_value(poly, x, zero_power, transition)
  107. def _symeig_backward_partial_eigenspace(D_grad, U_grad, A, D, U, largest):
  108. # compute a projection operator onto an orthogonal subspace spanned by the
  109. # columns of U defined as (I - UU^T)
  110. Ut = U.mT.contiguous()
  111. proj_U_ortho = -U.matmul(Ut)
  112. proj_U_ortho.diagonal(dim1=-2, dim2=-1).add_(1)
  113. # compute U_ortho, a basis for the orthogonal complement to the span(U),
  114. # by projecting a random [..., m, m - k] matrix onto the subspace spanned
  115. # by the columns of U.
  116. #
  117. # fix generator for determinism
  118. gen = torch.Generator(A.device)
  119. # orthogonal complement to the span(U)
  120. U_ortho = proj_U_ortho.matmul(
  121. torch.randn(
  122. (*A.shape[:-1], A.size(-1) - D.size(-1)),
  123. dtype=A.dtype,
  124. device=A.device,
  125. generator=gen,
  126. )
  127. )
  128. U_ortho_t = U_ortho.mT.contiguous()
  129. # compute the coefficients of the characteristic polynomial of the tensor D.
  130. # Note that D is diagonal, so the diagonal elements are exactly the roots
  131. # of the characteristic polynomial.
  132. chr_poly_D = _polynomial_coefficients_given_roots(D)
  133. # the code below finds the explicit solution to the Sylvester equation
  134. # U_ortho^T A U_ortho dX - dX D = -U_ortho^T A U
  135. # and incorporates it into the whole gradient stored in the `res` variable.
  136. #
  137. # Equivalent to the following naive implementation:
  138. # res = A.new_zeros(A.shape)
  139. # p_res = A.new_zeros(*A.shape[:-1], D.size(-1))
  140. # for k in range(1, chr_poly_D.size(-1)):
  141. # p_res.zero_()
  142. # for i in range(0, k):
  143. # p_res += (A.matrix_power(k - 1 - i) @ U_grad) * D.pow(i).unsqueeze(-2)
  144. # res -= chr_poly_D[k] * (U_ortho @ poly_D_at_A.inverse() @ U_ortho_t @ p_res @ U.t())
  145. #
  146. # Note that dX is a differential, so the gradient contribution comes from the backward sensitivity
  147. # Tr(f(U_grad, D_grad, A, U, D)^T dX) = Tr(g(U_grad, A, U, D)^T dA) for some functions f and g,
  148. # and we need to compute g(U_grad, A, U, D)
  149. #
  150. # The naive implementation is based on the paper
  151. # Hu, Qingxi, and Daizhan Cheng.
  152. # "The polynomial solution to the Sylvester matrix equation."
  153. # Applied mathematics letters 19.9 (2006): 859-864.
  154. #
  155. # We can modify the computation of `p_res` from above in a more efficient way
  156. # p_res = U_grad * (chr_poly_D[1] * D.pow(0) + ... + chr_poly_D[k] * D.pow(k)).unsqueeze(-2)
  157. # + A U_grad * (chr_poly_D[2] * D.pow(0) + ... + chr_poly_D[k] * D.pow(k - 1)).unsqueeze(-2)
  158. # + ...
  159. # + A.matrix_power(k - 1) U_grad * chr_poly_D[k]
  160. # Note that this saves us from redundant matrix products with A (elimination of matrix_power)
  161. U_grad_projected = U_grad
  162. series_acc = U_grad_projected.new_zeros(U_grad_projected.shape)
  163. for k in range(1, chr_poly_D.size(-1)):
  164. poly_D = _vector_polynomial_value(chr_poly_D[..., k:], D)
  165. series_acc += U_grad_projected * poly_D.unsqueeze(-2)
  166. U_grad_projected = A.matmul(U_grad_projected)
  167. # compute chr_poly_D(A) which essentially is:
  168. #
  169. # chr_poly_D_at_A = A.new_zeros(A.shape)
  170. # for k in range(chr_poly_D.size(-1)):
  171. # chr_poly_D_at_A += chr_poly_D[k] * A.matrix_power(k)
  172. #
  173. # Note, however, for better performance we use the Horner's rule
  174. chr_poly_D_at_A = _matrix_polynomial_value(chr_poly_D, A)
  175. # compute the action of `chr_poly_D_at_A` restricted to U_ortho_t
  176. chr_poly_D_at_A_to_U_ortho = torch.matmul(
  177. U_ortho_t, torch.matmul(chr_poly_D_at_A, U_ortho)
  178. )
  179. # we need to invert 'chr_poly_D_at_A_to_U_ortho`, for that we compute its
  180. # Cholesky decomposition and then use `torch.cholesky_solve` for better stability.
  181. # Cholesky decomposition requires the input to be positive-definite.
  182. # Note that `chr_poly_D_at_A_to_U_ortho` is positive-definite if
  183. # 1. `largest` == False, or
  184. # 2. `largest` == True and `k` is even
  185. # under the assumption that `A` has distinct eigenvalues.
  186. #
  187. # check if `chr_poly_D_at_A_to_U_ortho` is positive-definite or negative-definite
  188. chr_poly_D_at_A_to_U_ortho_sign = -1 if (largest and (k % 2 == 1)) else +1
  189. chr_poly_D_at_A_to_U_ortho_L = torch.linalg.cholesky(
  190. chr_poly_D_at_A_to_U_ortho_sign * chr_poly_D_at_A_to_U_ortho
  191. )
  192. # compute the gradient part in span(U)
  193. res = _symeig_backward_complete_eigenspace(D_grad, U_grad, A, D, U)
  194. # incorporate the Sylvester equation solution into the full gradient
  195. # it resides in span(U_ortho)
  196. res -= U_ortho.matmul(
  197. chr_poly_D_at_A_to_U_ortho_sign
  198. * torch.cholesky_solve(
  199. U_ortho_t.matmul(series_acc), chr_poly_D_at_A_to_U_ortho_L
  200. )
  201. ).matmul(Ut)
  202. return res
  203. def _symeig_backward(D_grad, U_grad, A, D, U, largest):
  204. # if `U` is square, then the columns of `U` is a complete eigenspace
  205. if U.size(-1) == U.size(-2):
  206. return _symeig_backward_complete_eigenspace(D_grad, U_grad, A, D, U)
  207. else:
  208. return _symeig_backward_partial_eigenspace(D_grad, U_grad, A, D, U, largest)
  209. class LOBPCGAutogradFunction(torch.autograd.Function):
  210. @staticmethod
  211. def forward( # type: ignore[override]
  212. ctx,
  213. A: Tensor,
  214. k: int | None = None,
  215. B: Tensor | None = None,
  216. X: Tensor | None = None,
  217. n: int | None = None,
  218. iK: Tensor | None = None,
  219. niter: int | None = None,
  220. tol: float | None = None,
  221. largest: bool | None = None,
  222. method: str | None = None,
  223. tracker: None = None,
  224. ortho_iparams: dict[str, int] | None = None,
  225. ortho_fparams: dict[str, float] | None = None,
  226. ortho_bparams: dict[str, bool] | None = None,
  227. ) -> tuple[Tensor, Tensor]:
  228. # makes sure that input is contiguous for efficiency.
  229. # Note: autograd does not support dense gradients for sparse input yet.
  230. A = A.contiguous() if (not A.is_sparse) else A
  231. if B is not None:
  232. B = B.contiguous() if (not B.is_sparse) else B
  233. D, U = _lobpcg(
  234. A,
  235. k,
  236. B,
  237. X,
  238. n,
  239. iK,
  240. niter,
  241. tol,
  242. largest,
  243. method,
  244. tracker,
  245. ortho_iparams,
  246. ortho_fparams,
  247. ortho_bparams,
  248. )
  249. ctx.save_for_backward(A, B, D, U)
  250. ctx.largest = largest
  251. return D, U
  252. @staticmethod
  253. def backward(ctx, D_grad, U_grad): # pyrefly: ignore # bad-override
  254. A_grad = B_grad = None
  255. grads = [None] * 14
  256. A, B, D, U = ctx.saved_tensors
  257. largest = ctx.largest
  258. # lobpcg.backward has some limitations. Checks for unsupported input
  259. if A.is_sparse or (B is not None and B.is_sparse and ctx.needs_input_grad[2]):
  260. raise ValueError(
  261. "lobpcg.backward does not support sparse input yet."
  262. "Note that lobpcg.forward does though."
  263. )
  264. if (
  265. A.dtype in (torch.complex64, torch.complex128)
  266. or B is not None
  267. and B.dtype in (torch.complex64, torch.complex128)
  268. ):
  269. raise ValueError(
  270. "lobpcg.backward does not support complex input yet."
  271. "Note that lobpcg.forward does though."
  272. )
  273. if B is not None:
  274. raise ValueError(
  275. "lobpcg.backward does not support backward with B != I yet."
  276. )
  277. if largest is None:
  278. largest = True
  279. # symeig backward
  280. if B is None:
  281. A_grad = _symeig_backward(D_grad, U_grad, A, D, U, largest)
  282. # A has index 0
  283. grads[0] = A_grad
  284. # B has index 2
  285. grads[2] = B_grad
  286. return tuple(grads)
  287. def lobpcg(
  288. A: Tensor,
  289. k: int | None = None,
  290. B: Tensor | None = None,
  291. X: Tensor | None = None,
  292. n: int | None = None,
  293. iK: Tensor | None = None,
  294. niter: int | None = None,
  295. tol: float | None = None,
  296. largest: bool | None = None,
  297. method: str | None = None,
  298. tracker: None = None,
  299. ortho_iparams: dict[str, int] | None = None,
  300. ortho_fparams: dict[str, float] | None = None,
  301. ortho_bparams: dict[str, bool] | None = None,
  302. ) -> tuple[Tensor, Tensor]:
  303. """Find the k largest (or smallest) eigenvalues and the corresponding
  304. eigenvectors of a symmetric positive definite generalized
  305. eigenvalue problem using matrix-free LOBPCG methods.
  306. This function is a front-end to the following LOBPCG algorithms
  307. selectable via `method` argument:
  308. `method="basic"` - the LOBPCG method introduced by Andrew
  309. Knyazev, see [Knyazev2001]. A less robust method, may fail when
  310. Cholesky is applied to singular input.
  311. `method="ortho"` - the LOBPCG method with orthogonal basis
  312. selection [StathopoulosEtal2002]. A robust method.
  313. Supported inputs are dense, sparse, and batches of dense matrices.
  314. .. note:: In general, the basic method spends least time per
  315. iteration. However, the robust methods converge much faster and
  316. are more stable. So, the usage of the basic method is generally
  317. not recommended but there exist cases where the usage of the
  318. basic method may be preferred.
  319. .. warning:: The backward method does not support sparse and complex inputs.
  320. It works only when `B` is not provided (i.e. `B == None`).
  321. We are actively working on extensions, and the details of
  322. the algorithms are going to be published promptly.
  323. .. warning:: While it is assumed that `A` is symmetric, `A.grad` is not.
  324. To make sure that `A.grad` is symmetric, so that `A - t * A.grad` is symmetric
  325. in first-order optimization routines, prior to running `lobpcg`
  326. we do the following symmetrization map: `A -> (A + A.t()) / 2`.
  327. The map is performed only when the `A` requires gradients.
  328. .. warning:: LOBPCG algorithm is not applicable when the number of `A`'s rows
  329. is smaller than 3x the number of requested eigenpairs `n`.
  330. Args:
  331. A (Tensor): the input tensor of size :math:`(*, m, m)`
  332. k (integer, optional): the number of requested
  333. eigenpairs. Default is the number of :math:`X`
  334. columns (when specified) or `1`.
  335. B (Tensor, optional): the input tensor of size :math:`(*, m,
  336. m)`. When not specified, `B` is interpreted as
  337. identity matrix.
  338. X (tensor, optional): the input tensor of size :math:`(*, m, n)`
  339. where `k <= n <= m`. When specified, it is used as
  340. initial approximation of eigenvectors. X must be a
  341. dense tensor.
  342. n (integer, optional): if :math:`X` is not specified then `n`
  343. specifies the size of the generated random
  344. approximation of eigenvectors. Default value for `n`
  345. is `k`. If :math:`X` is specified, any provided value of `n` is
  346. ignored and `n` is automatically set to the number of
  347. columns in :math:`X`.
  348. iK (tensor, optional): the input tensor of size :math:`(*, m,
  349. m)`. When specified, it will be used as preconditioner.
  350. niter (int, optional): maximum number of iterations. When
  351. reached, the iteration process is hard-stopped and
  352. the current approximation of eigenpairs is returned.
  353. For infinite iteration but until convergence criteria
  354. is met, use `-1`.
  355. tol (float, optional): residual tolerance for stopping
  356. criterion. Default is `feps ** 0.5` where `feps` is
  357. smallest non-zero floating-point number of the given
  358. input tensor `A` data type.
  359. largest (bool, optional): when True, solve the eigenproblem for
  360. the largest eigenvalues. Otherwise, solve the
  361. eigenproblem for smallest eigenvalues. Default is
  362. `True`.
  363. method (str, optional): select LOBPCG method. See the
  364. description of the function above. Default is
  365. "ortho".
  366. tracker (callable, optional) : a function for tracing the
  367. iteration process. When specified, it is called at
  368. each iteration step with LOBPCG instance as an
  369. argument. The LOBPCG instance holds the full state of
  370. the iteration process in the following attributes:
  371. `iparams`, `fparams`, `bparams` - dictionaries of
  372. integer, float, and boolean valued input
  373. parameters, respectively
  374. `ivars`, `fvars`, `bvars`, `tvars` - dictionaries
  375. of integer, float, boolean, and Tensor valued
  376. iteration variables, respectively.
  377. `A`, `B`, `iK` - input Tensor arguments.
  378. `E`, `X`, `S`, `R` - iteration Tensor variables.
  379. For instance:
  380. `ivars["istep"]` - the current iteration step
  381. `X` - the current approximation of eigenvectors
  382. `E` - the current approximation of eigenvalues
  383. `R` - the current residual
  384. `ivars["converged_count"]` - the current number of converged eigenpairs
  385. `tvars["rerr"]` - the current state of convergence criteria
  386. Note that when `tracker` stores Tensor objects from
  387. the LOBPCG instance, it must make copies of these.
  388. If `tracker` sets `bvars["force_stop"] = True`, the
  389. iteration process will be hard-stopped.
  390. ortho_iparams, ortho_fparams, ortho_bparams (dict, optional):
  391. various parameters to LOBPCG algorithm when using
  392. `method="ortho"`.
  393. Returns:
  394. E (Tensor): tensor of eigenvalues of size :math:`(*, k)`
  395. X (Tensor): tensor of eigenvectors of size :math:`(*, m, k)`
  396. References:
  397. [Knyazev2001] Andrew V. Knyazev. (2001) Toward the Optimal
  398. Preconditioned Eigensolver: Locally Optimal Block Preconditioned
  399. Conjugate Gradient Method. SIAM J. Sci. Comput., 23(2),
  400. 517-541. (25 pages)
  401. https://epubs.siam.org/doi/abs/10.1137/S1064827500366124
  402. [StathopoulosEtal2002] Andreas Stathopoulos and Kesheng
  403. Wu. (2002) A Block Orthogonalization Procedure with Constant
  404. Synchronization Requirements. SIAM J. Sci. Comput., 23(6),
  405. 2165-2182. (18 pages)
  406. https://epubs.siam.org/doi/10.1137/S1064827500370883
  407. [DuerschEtal2018] Jed A. Duersch, Meiyue Shao, Chao Yang, Ming
  408. Gu. (2018) A Robust and Efficient Implementation of LOBPCG.
  409. SIAM J. Sci. Comput., 40(5), C655-C676. (22 pages)
  410. https://arxiv.org/abs/1704.07458
  411. """
  412. if not torch.jit.is_scripting():
  413. tensor_ops = (A, B, X, iK)
  414. if not set(map(type, tensor_ops)).issubset(
  415. (torch.Tensor, type(None))
  416. ) and has_torch_function(tensor_ops):
  417. return handle_torch_function(
  418. lobpcg,
  419. tensor_ops,
  420. A,
  421. k=k,
  422. B=B,
  423. X=X,
  424. n=n,
  425. iK=iK,
  426. niter=niter,
  427. tol=tol,
  428. largest=largest,
  429. method=method,
  430. tracker=tracker,
  431. ortho_iparams=ortho_iparams,
  432. ortho_fparams=ortho_fparams,
  433. ortho_bparams=ortho_bparams,
  434. )
  435. if not torch._jit_internal.is_scripting():
  436. if A.requires_grad or (B is not None and B.requires_grad):
  437. # While it is expected that `A` is symmetric,
  438. # the `A_grad` might be not. Therefore we perform the trick below,
  439. # so that `A_grad` becomes symmetric.
  440. # The symmetrization is important for first-order optimization methods,
  441. # so that (A - alpha * A_grad) is still a symmetric matrix.
  442. # Same holds for `B`.
  443. A_sym = (A + A.mT) / 2
  444. B_sym = (B + B.mT) / 2 if (B is not None) else None
  445. return LOBPCGAutogradFunction.apply(
  446. A_sym,
  447. k,
  448. B_sym,
  449. X,
  450. n,
  451. iK,
  452. niter,
  453. tol,
  454. largest,
  455. method,
  456. tracker,
  457. ortho_iparams,
  458. ortho_fparams,
  459. ortho_bparams,
  460. )
  461. else:
  462. if A.requires_grad or (B is not None and B.requires_grad):
  463. raise RuntimeError(
  464. "Script and require grads is not supported atm."
  465. "If you just want to do the forward, use .detach()"
  466. "on A and B before calling into lobpcg"
  467. )
  468. return _lobpcg(
  469. A,
  470. k,
  471. B,
  472. X,
  473. n,
  474. iK,
  475. niter,
  476. tol,
  477. largest,
  478. method,
  479. tracker,
  480. ortho_iparams,
  481. ortho_fparams,
  482. ortho_bparams,
  483. )
  484. def _lobpcg(
  485. A: Tensor,
  486. k: int | None = None,
  487. B: Tensor | None = None,
  488. X: Tensor | None = None,
  489. n: int | None = None,
  490. iK: Tensor | None = None,
  491. niter: int | None = None,
  492. tol: float | None = None,
  493. largest: bool | None = None,
  494. method: str | None = None,
  495. tracker: None = None,
  496. ortho_iparams: dict[str, int] | None = None,
  497. ortho_fparams: dict[str, float] | None = None,
  498. ortho_bparams: dict[str, bool] | None = None,
  499. ) -> tuple[Tensor, Tensor]:
  500. # A must be square:
  501. if A.shape[-2] != A.shape[-1]:
  502. raise AssertionError(f"A must be square, got shape {A.shape}")
  503. if B is not None:
  504. # A and B must have the same shapes:
  505. if A.shape != B.shape:
  506. raise AssertionError(
  507. f"A and B must have same shape, got {A.shape} vs {B.shape}"
  508. )
  509. dtype = _utils.get_floating_dtype(A)
  510. device = A.device
  511. if tol is None:
  512. feps = {torch.float32: 1.2e-07, torch.float64: 2.23e-16}[dtype]
  513. tol = feps**0.5
  514. m = A.shape[-1]
  515. k = (1 if X is None else X.shape[-1]) if k is None else k
  516. n = (k if n is None else n) if X is None else X.shape[-1]
  517. if m < 3 * n:
  518. raise ValueError(
  519. f"LPBPCG algorithm is not applicable when the number of A rows (={m})"
  520. f" is smaller than 3 x the number of requested eigenpairs (={n})"
  521. )
  522. method = "ortho" if method is None else method
  523. iparams = {
  524. "m": m,
  525. "n": n,
  526. "k": k,
  527. "niter": 1000 if niter is None else niter,
  528. }
  529. fparams = {
  530. "tol": tol,
  531. }
  532. bparams = {"largest": True if largest is None else largest}
  533. if method == "ortho":
  534. if ortho_iparams is not None:
  535. iparams.update(ortho_iparams)
  536. if ortho_fparams is not None:
  537. fparams.update(ortho_fparams)
  538. if ortho_bparams is not None:
  539. bparams.update(ortho_bparams)
  540. iparams["ortho_i_max"] = iparams.get("ortho_i_max", 3)
  541. iparams["ortho_j_max"] = iparams.get("ortho_j_max", 3)
  542. fparams["ortho_tol"] = fparams.get("ortho_tol", tol)
  543. fparams["ortho_tol_drop"] = fparams.get("ortho_tol_drop", tol)
  544. fparams["ortho_tol_replace"] = fparams.get("ortho_tol_replace", tol)
  545. bparams["ortho_use_drop"] = bparams.get("ortho_use_drop", False)
  546. if not torch.jit.is_scripting():
  547. LOBPCG.call_tracker = LOBPCG_call_tracker # type: ignore[method-assign]
  548. if len(A.shape) > 2:
  549. N = int(torch.prod(torch.tensor(A.shape[:-2])))
  550. bA = A.reshape((N,) + A.shape[-2:])
  551. bB = B.reshape((N,) + A.shape[-2:]) if B is not None else None
  552. bX = X.reshape((N,) + X.shape[-2:]) if X is not None else None
  553. bE = torch.empty((N, k), dtype=dtype, device=device)
  554. bXret = torch.empty((N, m, k), dtype=dtype, device=device)
  555. for i in range(N):
  556. A_ = bA[i]
  557. B_ = bB[i] if bB is not None else None
  558. X_ = (
  559. torch.randn((m, n), dtype=dtype, device=device) if bX is None else bX[i]
  560. )
  561. if not (len(X_.shape) == 2 and X_.shape == (m, n)):
  562. raise AssertionError(
  563. f"X_ shape mismatch: got {X_.shape}, expected {(m, n)}"
  564. )
  565. iparams["batch_index"] = i
  566. worker = LOBPCG(A_, B_, X_, iK, iparams, fparams, bparams, method, tracker)
  567. worker.run()
  568. bE[i] = worker.E[:k]
  569. bXret[i] = worker.X[:, :k]
  570. if not torch.jit.is_scripting():
  571. LOBPCG.call_tracker = LOBPCG_call_tracker_orig # type: ignore[method-assign]
  572. return bE.reshape(A.shape[:-2] + (k,)), bXret.reshape(A.shape[:-2] + (m, k))
  573. X = torch.randn((m, n), dtype=dtype, device=device) if X is None else X
  574. if not (len(X.shape) == 2 and X.shape == (m, n)):
  575. raise AssertionError(f"X shape mismatch: got {X.shape}, expected {(m, n)}")
  576. worker = LOBPCG(A, B, X, iK, iparams, fparams, bparams, method, tracker)
  577. worker.run()
  578. if not torch.jit.is_scripting():
  579. LOBPCG.call_tracker = LOBPCG_call_tracker_orig # type: ignore[method-assign]
  580. return worker.E[:k], worker.X[:, :k]
  581. class LOBPCG:
  582. """Worker class of LOBPCG methods."""
  583. def __init__(
  584. self,
  585. A: Tensor | None,
  586. B: Tensor | None,
  587. X: Tensor,
  588. iK: Tensor | None,
  589. iparams: dict[str, int],
  590. fparams: dict[str, float],
  591. bparams: dict[str, bool],
  592. method: str,
  593. tracker: None,
  594. ) -> None:
  595. # constant parameters
  596. self.A = A
  597. self.B = B
  598. self.iK = iK
  599. self.iparams = iparams
  600. self.fparams = fparams
  601. self.bparams = bparams
  602. self.method = method
  603. self.tracker = tracker
  604. m = iparams["m"]
  605. n = iparams["n"]
  606. # variable parameters
  607. self.X = X
  608. self.E = torch.zeros((n,), dtype=X.dtype, device=X.device)
  609. self.R = torch.zeros((m, n), dtype=X.dtype, device=X.device)
  610. self.S = torch.zeros((m, 3 * n), dtype=X.dtype, device=X.device)
  611. self.tvars: dict[str, Tensor] = {}
  612. self.ivars: dict[str, int] = {"istep": 0}
  613. self.fvars: dict[str, float] = {"_": 0.0}
  614. self.bvars: dict[str, bool] = {"_": False}
  615. def __str__(self):
  616. lines = ["LOPBCG:"]
  617. lines += [f" iparams={self.iparams}"]
  618. lines += [f" fparams={self.fparams}"]
  619. lines += [f" bparams={self.bparams}"]
  620. lines += [f" ivars={self.ivars}"]
  621. lines += [f" fvars={self.fvars}"]
  622. lines += [f" bvars={self.bvars}"]
  623. lines += [f" tvars={self.tvars}"]
  624. lines += [f" A={self.A}"]
  625. lines += [f" B={self.B}"]
  626. lines += [f" iK={self.iK}"]
  627. lines += [f" X={self.X}"]
  628. lines += [f" E={self.E}"]
  629. r = ""
  630. for line in lines:
  631. r += line + "\n"
  632. return r
  633. def update(self):
  634. """Set and update iteration variables."""
  635. if self.ivars["istep"] == 0:
  636. X_norm = float(torch.norm(self.X))
  637. iX_norm = X_norm**-1
  638. A_norm = float(torch.norm(_utils.matmul(self.A, self.X))) * iX_norm
  639. B_norm = float(torch.norm(_utils.matmul(self.B, self.X))) * iX_norm
  640. self.fvars["X_norm"] = X_norm
  641. self.fvars["A_norm"] = A_norm
  642. self.fvars["B_norm"] = B_norm
  643. self.ivars["iterations_left"] = self.iparams["niter"]
  644. self.ivars["converged_count"] = 0
  645. self.ivars["converged_end"] = 0
  646. if self.method == "ortho":
  647. self._update_ortho()
  648. else:
  649. self._update_basic()
  650. self.ivars["iterations_left"] = self.ivars["iterations_left"] - 1
  651. self.ivars["istep"] = self.ivars["istep"] + 1
  652. def update_residual(self):
  653. """Update residual R from A, B, X, E."""
  654. mm = _utils.matmul
  655. self.R = mm(self.A, self.X) - mm(self.B, self.X) * self.E
  656. def update_converged_count(self):
  657. """Determine the number of converged eigenpairs using backward stable
  658. convergence criterion, see discussion in Sec 4.3 of [DuerschEtal2018].
  659. Users may redefine this method for custom convergence criteria.
  660. """
  661. # (...) -> int
  662. prev_count = self.ivars["converged_count"]
  663. tol = self.fparams["tol"]
  664. A_norm = self.fvars["A_norm"]
  665. B_norm = self.fvars["B_norm"]
  666. E, X, R = self.E, self.X, self.R
  667. rerr = torch.norm(R, 2, (0,)) / (
  668. torch.norm(X, 2, (0,)) * (A_norm + torch.abs(E[: X.shape[-1]]) * B_norm)
  669. )
  670. converged = rerr < tol
  671. count = 0
  672. for b in converged:
  673. if not b:
  674. # ignore convergence of following pairs to ensure
  675. # strict ordering of eigenpairs
  676. break
  677. count += 1
  678. if count < prev_count:
  679. raise AssertionError(
  680. f"the number of converged eigenpairs (was {prev_count}, got {count}) cannot decrease"
  681. )
  682. self.ivars["converged_count"] = count
  683. self.tvars["rerr"] = rerr
  684. return count
  685. def stop_iteration(self):
  686. """Return True to stop iterations.
  687. Note that tracker (if defined) can force-stop iterations by
  688. setting ``worker.bvars['force_stop'] = True``.
  689. """
  690. return (
  691. self.bvars.get("force_stop", False)
  692. or self.ivars["iterations_left"] == 0
  693. or self.ivars["converged_count"] >= self.iparams["k"]
  694. )
  695. def run(self):
  696. """Run LOBPCG iterations.
  697. Use this method as a template for implementing LOBPCG
  698. iteration scheme with custom tracker that is compatible with
  699. TorchScript.
  700. """
  701. self.update()
  702. if not torch.jit.is_scripting() and self.tracker is not None:
  703. self.call_tracker()
  704. while not self.stop_iteration():
  705. self.update()
  706. if not torch.jit.is_scripting() and self.tracker is not None:
  707. self.call_tracker()
  708. @torch.jit.unused
  709. def call_tracker(self):
  710. """Interface for tracking iteration process in Python mode.
  711. Tracking the iteration process is disabled in TorchScript
  712. mode. In fact, one should specify tracker=None when JIT
  713. compiling functions using lobpcg.
  714. """
  715. # do nothing when in TorchScript mode
  716. # Internal methods
  717. def _update_basic(self):
  718. """
  719. Update or initialize iteration variables when `method == "basic"`.
  720. """
  721. mm = torch.matmul
  722. ns = self.ivars["converged_end"]
  723. nc = self.ivars["converged_count"]
  724. n = self.iparams["n"]
  725. largest = self.bparams["largest"]
  726. if self.ivars["istep"] == 0:
  727. Ri = self._get_rayleigh_ritz_transform(self.X)
  728. M = _utils.qform(_utils.qform(self.A, self.X), Ri)
  729. E, Z = _utils.symeig(M, largest)
  730. self.X[:] = mm(self.X, mm(Ri, Z))
  731. self.E[:] = E
  732. np = 0
  733. self.update_residual()
  734. nc = self.update_converged_count()
  735. self.S[..., :n] = self.X
  736. W = _utils.matmul(self.iK, self.R)
  737. self.ivars["converged_end"] = ns = n + np + W.shape[-1]
  738. self.S[:, n + np : ns] = W
  739. else:
  740. S_ = self.S[:, nc:ns]
  741. Ri = self._get_rayleigh_ritz_transform(S_)
  742. M = _utils.qform(_utils.qform(self.A, S_), Ri)
  743. E_, Z = _utils.symeig(M, largest)
  744. self.X[:, nc:] = mm(S_, mm(Ri, Z[:, : n - nc]))
  745. self.E[nc:] = E_[: n - nc]
  746. P = mm(S_, mm(Ri, Z[:, n : 2 * n - nc]))
  747. np = P.shape[-1]
  748. self.update_residual()
  749. nc = self.update_converged_count()
  750. self.S[..., :n] = self.X
  751. self.S[:, n : n + np] = P
  752. W = _utils.matmul(self.iK, self.R[:, nc:])
  753. self.ivars["converged_end"] = ns = n + np + W.shape[-1]
  754. self.S[:, n + np : ns] = W
  755. def _update_ortho(self):
  756. """
  757. Update or initialize iteration variables when `method == "ortho"`.
  758. """
  759. mm = torch.matmul
  760. ns = self.ivars["converged_end"]
  761. nc = self.ivars["converged_count"]
  762. n = self.iparams["n"]
  763. largest = self.bparams["largest"]
  764. if self.ivars["istep"] == 0:
  765. Ri = self._get_rayleigh_ritz_transform(self.X)
  766. M = _utils.qform(_utils.qform(self.A, self.X), Ri)
  767. _E, Z = _utils.symeig(M, largest)
  768. self.X = mm(self.X, mm(Ri, Z))
  769. self.update_residual()
  770. np = 0
  771. nc = self.update_converged_count()
  772. self.S[:, :n] = self.X
  773. W = self._get_ortho(self.R, self.X)
  774. ns = self.ivars["converged_end"] = n + np + W.shape[-1]
  775. self.S[:, n + np : ns] = W
  776. else:
  777. S_ = self.S[:, nc:ns]
  778. # Rayleigh-Ritz procedure
  779. E_, Z = _utils.symeig(_utils.qform(self.A, S_), largest)
  780. # Update E, X, P
  781. self.X[:, nc:] = mm(S_, Z[:, : n - nc])
  782. self.E[nc:] = E_[: n - nc]
  783. P = mm(S_, mm(Z[:, n - nc :], _utils.basis(Z[: n - nc, n - nc :].mT)))
  784. np = P.shape[-1]
  785. # check convergence
  786. self.update_residual()
  787. nc = self.update_converged_count()
  788. # update S
  789. self.S[:, :n] = self.X
  790. self.S[:, n : n + np] = P
  791. W = self._get_ortho(self.R[:, nc:], self.S[:, : n + np])
  792. ns = self.ivars["converged_end"] = n + np + W.shape[-1]
  793. self.S[:, n + np : ns] = W
  794. def _get_rayleigh_ritz_transform(self, S):
  795. """Return a transformation matrix that is used in Rayleigh-Ritz
  796. procedure for reducing a general eigenvalue problem :math:`(S^TAS)
  797. C = (S^TBS) C E` to a standard eigenvalue problem :math: `(Ri^T
  798. S^TAS Ri) Z = Z E` where `C = Ri Z`.
  799. .. note:: In the original Rayleight-Ritz procedure in
  800. [DuerschEtal2018], the problem is formulated as follows::
  801. SAS = S^T A S
  802. SBS = S^T B S
  803. D = (<diagonal matrix of SBS>) ** -1/2
  804. R^T R = Cholesky(D SBS D)
  805. Ri = D R^-1
  806. solve symeig problem Ri^T SAS Ri Z = Theta Z
  807. C = Ri Z
  808. To reduce the number of matrix products (denoted by empty
  809. space between matrices), here we introduce element-wise
  810. products (denoted by symbol `*`) so that the Rayleight-Ritz
  811. procedure becomes::
  812. SAS = S^T A S
  813. SBS = S^T B S
  814. d = (<diagonal of SBS>) ** -1/2 # this is 1-d column vector
  815. dd = d d^T # this is 2-d matrix
  816. R^T R = Cholesky(dd * SBS)
  817. Ri = R^-1 * d # broadcasting
  818. solve symeig problem Ri^T SAS Ri Z = Theta Z
  819. C = Ri Z
  820. where `dd` is 2-d matrix that replaces matrix products `D M
  821. D` with one element-wise product `M * dd`; and `d` replaces
  822. matrix product `D M` with element-wise product `M *
  823. d`. Also, creating the diagonal matrix `D` is avoided.
  824. Args:
  825. S (Tensor): the matrix basis for the search subspace, size is
  826. :math:`(m, n)`.
  827. Returns:
  828. Ri (tensor): upper-triangular transformation matrix of size
  829. :math:`(n, n)`.
  830. """
  831. B = self.B
  832. SBS = _utils.qform(B, S)
  833. d_row = SBS.diagonal(0, -2, -1) ** -0.5
  834. d_col = d_row.reshape(d_row.shape[0], 1)
  835. # TODO use torch.linalg.cholesky_solve once it is implemented
  836. R = torch.linalg.cholesky((SBS * d_row) * d_col, upper=True)
  837. return torch.linalg.solve_triangular(
  838. R, d_row.diag_embed(), upper=True, left=False
  839. )
  840. def _get_svqb(self, U: Tensor, drop: bool, tau: float) -> Tensor:
  841. """Return B-orthonormal U.
  842. .. note:: When `drop` is `False` then `svqb` is based on the
  843. Algorithm 4 from [DuerschPhD2015] that is a slight
  844. modification of the corresponding algorithm
  845. introduced in [StathopolousWu2002].
  846. Args:
  847. U (Tensor) : initial approximation, size is (m, n)
  848. drop (bool) : when True, drop columns that
  849. contribution to the `span([U])` is small.
  850. tau (float) : positive tolerance
  851. Returns:
  852. U (Tensor) : B-orthonormal columns (:math:`U^T B U = I`), size
  853. is (m, n1), where `n1 = n` if `drop` is `False,
  854. otherwise `n1 <= n`.
  855. """
  856. if torch.numel(U) == 0:
  857. return U
  858. UBU = _utils.qform(self.B, U)
  859. d = UBU.diagonal(0, -2, -1)
  860. # Detect and drop exact zero columns from U. While the test
  861. # `abs(d) == 0` is unlikely to be True for random data, it is
  862. # possible to construct input data to lobpcg where it will be
  863. # True leading to a failure (notice the `d ** -0.5` operation
  864. # in the original algorithm). To prevent the failure, we drop
  865. # the exact zero columns here and then continue with the
  866. # original algorithm below.
  867. nz = torch.where(abs(d) != 0.0)
  868. if len(nz) != 1:
  869. raise AssertionError(f"expected nz to have length 1, got {nz}")
  870. if len(nz[0]) < len(d):
  871. U = U[:, nz[0]]
  872. if torch.numel(U) == 0:
  873. return U
  874. UBU = _utils.qform(self.B, U)
  875. d = UBU.diagonal(0, -2, -1)
  876. nz = torch.where(abs(d) != 0.0)
  877. if len(nz[0]) != len(d):
  878. raise AssertionError(
  879. f"expected nz[0] length {len(d)}, got {len(nz[0])}"
  880. )
  881. # The original algorithm 4 from [DuerschPhD2015].
  882. d_col = (d**-0.5).reshape(d.shape[0], 1)
  883. DUBUD = (UBU * d_col) * d_col.mT
  884. E, Z = _utils.symeig(DUBUD)
  885. t = tau * abs(E).max()
  886. if drop:
  887. keep = torch.where(E > t)
  888. if len(keep) != 1:
  889. raise AssertionError(f"expected keep to have length 1, got {keep}")
  890. E = E[keep[0]]
  891. Z = Z[:, keep[0]]
  892. d_col = d_col[keep[0]]
  893. else:
  894. E[(torch.where(E < t))[0]] = t
  895. return torch.matmul(
  896. U * d_col.mT,
  897. Z * E**-0.5,
  898. )
  899. def _get_ortho(self, U, V):
  900. """Return B-orthonormal U with columns are B-orthogonal to V.
  901. .. note:: When `bparams["ortho_use_drop"] == False` then
  902. `_get_ortho` is based on the Algorithm 3 from
  903. [DuerschPhD2015] that is a slight modification of
  904. the corresponding algorithm introduced in
  905. [StathopolousWu2002]. Otherwise, the method
  906. implements Algorithm 6 from [DuerschPhD2015]
  907. .. note:: If all U columns are B-collinear to V then the
  908. returned tensor U will be empty.
  909. Args:
  910. U (Tensor) : initial approximation, size is (m, n)
  911. V (Tensor) : B-orthogonal external basis, size is (m, k)
  912. Returns:
  913. U (Tensor) : B-orthonormal columns (:math:`U^T B U = I`)
  914. such that :math:`V^T B U=0`, size is (m, n1),
  915. where `n1 = n` if `drop` is `False, otherwise
  916. `n1 <= n`.
  917. """
  918. mm = torch.matmul
  919. mm_B = _utils.matmul
  920. m = self.iparams["m"]
  921. tau_ortho = self.fparams["ortho_tol"]
  922. tau_drop = self.fparams["ortho_tol_drop"]
  923. tau_replace = self.fparams["ortho_tol_replace"]
  924. i_max = self.iparams["ortho_i_max"]
  925. j_max = self.iparams["ortho_j_max"]
  926. # when use_drop==True, enable dropping U columns that have
  927. # small contribution to the `span([U, V])`.
  928. use_drop = self.bparams["ortho_use_drop"]
  929. # clean up variables from the previous call
  930. for vkey in list(self.fvars.keys()):
  931. if vkey.startswith("ortho_") and vkey.endswith("_rerr"):
  932. self.fvars.pop(vkey)
  933. self.ivars.pop("ortho_i", 0)
  934. self.ivars.pop("ortho_j", 0)
  935. BV_norm = torch.norm(mm_B(self.B, V))
  936. BU = mm_B(self.B, U)
  937. VBU = mm(V.mT, BU)
  938. i = j = 0
  939. for i in range(i_max):
  940. U = U - mm(V, VBU)
  941. drop = False
  942. tau_svqb = tau_drop
  943. for j in range(j_max):
  944. if use_drop:
  945. U = self._get_svqb(U, drop, tau_svqb)
  946. drop = True
  947. tau_svqb = tau_replace
  948. else:
  949. U = self._get_svqb(U, False, tau_replace)
  950. if torch.numel(U) == 0:
  951. # all initial U columns are B-collinear to V
  952. self.ivars["ortho_i"] = i
  953. self.ivars["ortho_j"] = j
  954. return U
  955. BU = mm_B(self.B, U)
  956. UBU = mm(U.mT, BU)
  957. U_norm = torch.norm(U)
  958. BU_norm = torch.norm(BU)
  959. R = UBU - torch.eye(UBU.shape[-1], device=UBU.device, dtype=UBU.dtype)
  960. R_norm = torch.norm(R)
  961. # https://github.com/pytorch/pytorch/issues/33810 workaround:
  962. rerr = float(R_norm) * float(BU_norm * U_norm) ** -1
  963. vkey = f"ortho_UBUmI_rerr[{i}, {j}]"
  964. self.fvars[vkey] = rerr
  965. if rerr < tau_ortho:
  966. break
  967. VBU = mm(V.mT, BU)
  968. VBU_norm = torch.norm(VBU)
  969. U_norm = torch.norm(U)
  970. rerr = float(VBU_norm) * float(BV_norm * U_norm) ** -1
  971. vkey = f"ortho_VBU_rerr[{i}]"
  972. self.fvars[vkey] = rerr
  973. if rerr < tau_ortho:
  974. break
  975. if m < U.shape[-1] + V.shape[-1]:
  976. # TorchScript needs the class var to be assigned to a local to
  977. # do optional type refinement
  978. B = self.B
  979. if B is None:
  980. raise AssertionError("B must not be None")
  981. raise ValueError(
  982. "Overdetermined shape of U:"
  983. f" #B-cols(={B.shape[-1]}) >= #U-cols(={U.shape[-1]}) + #V-cols(={V.shape[-1]}) must hold"
  984. )
  985. self.ivars["ortho_i"] = i
  986. self.ivars["ortho_j"] = j
  987. return U
  988. # Calling tracker is separated from LOBPCG definitions because
  989. # TorchScript does not support user-defined callback arguments:
  990. LOBPCG_call_tracker_orig = LOBPCG.call_tracker
  991. def LOBPCG_call_tracker(self):
  992. self.tracker(self)