|
| 1 | +import numpy as np |
| 2 | +from scipy.linalg import get_lapack_funcs |
| 3 | + |
| 4 | +from pytensor.graph import Apply, Op |
| 5 | +from pytensor.tensor.basic import as_tensor, diagonal |
| 6 | +from pytensor.tensor.blockwise import Blockwise |
| 7 | +from pytensor.tensor.type import tensor, vector |
| 8 | + |
| 9 | + |
| 10 | +class LUFactorTridiagonal(Op): |
| 11 | + """Compute LU factorization of a tridiagonal matrix (lapack gttrf)""" |
| 12 | + |
| 13 | + __props__ = ( |
| 14 | + "overwrite_dl", |
| 15 | + "overwrite_d", |
| 16 | + "overwrite_du", |
| 17 | + ) |
| 18 | + gufunc_signature = "(dl),(d),(dl)->(dl),(d),(dl),(du2),(d)" |
| 19 | + |
| 20 | + def __init__(self, overwrite_dl=False, overwrite_d=False, overwrite_du=False): |
| 21 | + self.destroy_map = dm = {} |
| 22 | + if overwrite_dl: |
| 23 | + dm[0] = [0] |
| 24 | + if overwrite_d: |
| 25 | + dm[1] = [1] |
| 26 | + if overwrite_du: |
| 27 | + dm[2] = [2] |
| 28 | + self.overwrite_dl = overwrite_dl |
| 29 | + self.overwrite_d = overwrite_d |
| 30 | + self.overwrite_du = overwrite_du |
| 31 | + super().__init__() |
| 32 | + |
| 33 | + def make_node(self, dl, d, du): |
| 34 | + dl, d, du = map(as_tensor, (dl, d, du)) |
| 35 | + |
| 36 | + if not all(inp.type.ndim == 1 for inp in (dl, d, du)): |
| 37 | + raise ValueError("Diagonals must be vectors") |
| 38 | + |
| 39 | + ndl, nd, ndu = (inp.type.shape[-1] for inp in (dl, d, du)) |
| 40 | + n = ( |
| 41 | + ndl + 1 |
| 42 | + if ndl is not None |
| 43 | + else (nd if nd is not None else (ndu + 1 if ndu is not None else None)) |
| 44 | + ) |
| 45 | + dummy_arrays = [np.zeros((), dtype=inp.type.dtype) for inp in (dl, d, du)] |
| 46 | + out_dtype = get_lapack_funcs("gttrf", dummy_arrays).dtype |
| 47 | + outputs = [ |
| 48 | + vector(shape=(None if n is None else (n - 1),), dtype=out_dtype), |
| 49 | + vector(shape=(n,), dtype=out_dtype), |
| 50 | + vector(shape=(None if n is None else n - 1,), dtype=out_dtype), |
| 51 | + vector(shape=(None if n is None else n - 2,), dtype=out_dtype), |
| 52 | + vector(shape=(n,), dtype=np.int32), |
| 53 | + ] |
| 54 | + return Apply(self, [dl, d, du], outputs) |
| 55 | + |
| 56 | + def perform(self, node, inputs, output_storage): |
| 57 | + gttrf = get_lapack_funcs("gttrf", dtype=node.outputs[0].type.dtype) |
| 58 | + dl, d, du, du2, ipiv, _ = gttrf( |
| 59 | + *inputs, |
| 60 | + overwrite_dl=self.overwrite_dl, |
| 61 | + overwrite_d=self.overwrite_d, |
| 62 | + overwrite_du=self.overwrite_du, |
| 63 | + ) |
| 64 | + output_storage[0][0] = dl |
| 65 | + output_storage[1][0] = d |
| 66 | + output_storage[2][0] = du |
| 67 | + output_storage[3][0] = du2 |
| 68 | + output_storage[4][0] = ipiv |
| 69 | + |
| 70 | + |
| 71 | +class SolveLUFactorTridiagonal(Op): |
| 72 | + """Solve a system of linear equations with a tridiagonal coefficient matrix (lapack gttrs).""" |
| 73 | + |
| 74 | + __props__ = ("b_ndim", "overwrite_b", "transposed") |
| 75 | + |
| 76 | + def __init__(self, b_ndim: int, transposed: bool, overwrite_b=False): |
| 77 | + if b_ndim not in (1, 2): |
| 78 | + raise ValueError("b_ndim must be 1 or 2") |
| 79 | + if b_ndim == 1: |
| 80 | + self.gufunc_signature = "(dl),(d),(dl),(du2),(d),(d)->(d)" |
| 81 | + else: |
| 82 | + self.gufunc_signature = "(dl),(d),(dl),(du2),(d),(d,rhs)->(d,rhs)" |
| 83 | + if overwrite_b: |
| 84 | + self.destroy_map = {0: [5]} |
| 85 | + self.b_ndim = b_ndim |
| 86 | + self.transposed = transposed |
| 87 | + self.overwrite_b = overwrite_b |
| 88 | + super().__init__() |
| 89 | + |
| 90 | + def make_node(self, dl, d, du, du2, ipiv, b): |
| 91 | + dl, d, du, du2, ipiv, b = map(as_tensor, (dl, d, du, du2, ipiv, b)) |
| 92 | + |
| 93 | + if b.type.ndim != self.b_ndim: |
| 94 | + raise ValueError("Wrang number of dimensions for input b.") |
| 95 | + |
| 96 | + if not all(inp.type.ndim == 1 for inp in (dl, d, du, du2, ipiv)): |
| 97 | + raise ValueError("Inputs must be vectors") |
| 98 | + |
| 99 | + ndl, nd, ndu, ndu2, nipiv = ( |
| 100 | + inp.type.shape[-1] for inp in (dl, d, du, du2, ipiv) |
| 101 | + ) |
| 102 | + nb = b.type.shape[0] |
| 103 | + n = ( |
| 104 | + ndl + 1 |
| 105 | + if ndl is not None |
| 106 | + else ( |
| 107 | + nd |
| 108 | + if nd is not None |
| 109 | + else ( |
| 110 | + ndu + 1 |
| 111 | + if ndu is not None |
| 112 | + else ( |
| 113 | + ndu2 + 2 |
| 114 | + if ndu2 is not None |
| 115 | + else (nipiv if nipiv is not None else nb) |
| 116 | + ) |
| 117 | + ) |
| 118 | + ) |
| 119 | + ) |
| 120 | + dummy_arrays = [ |
| 121 | + np.zeros((), dtype=inp.type.dtype) for inp in (dl, d, du, du2, ipiv) |
| 122 | + ] |
| 123 | + # Seems to always be float64? |
| 124 | + out_dtype = get_lapack_funcs("gttrs", dummy_arrays).dtype |
| 125 | + if self.b_ndim == 1: |
| 126 | + output_shape = (n,) |
| 127 | + else: |
| 128 | + output_shape = (n, b.type.shape[-1]) |
| 129 | + |
| 130 | + outputs = [tensor(shape=output_shape, dtype=out_dtype)] |
| 131 | + return Apply(self, [dl, d, du, du2, ipiv, b], outputs) |
| 132 | + |
| 133 | + def perform(self, node, inputs, output_storage): |
| 134 | + gttrs = get_lapack_funcs("gttrs", dtype=node.outputs[0].type.dtype) |
| 135 | + x, _ = gttrs( |
| 136 | + *inputs, |
| 137 | + overwrite_b=self.overwrite_b, |
| 138 | + trans="N" if not self.transposed else "T", |
| 139 | + ) |
| 140 | + output_storage[0][0] = x |
| 141 | + |
| 142 | + |
| 143 | +def tridiagonal_lu_factor(a): |
| 144 | + # Return the decomposition of A implied by a solve tridiagonal |
| 145 | + dl, d, du = (diagonal(a, offset=o, axis1=-2, axis2=-1) for o in (-1, 0, 1)) |
| 146 | + dl, d, du, du2, ipiv = Blockwise(LUFactorTridiagonal())(dl, d, du) |
| 147 | + return dl, d, du, du2, ipiv |
| 148 | + |
| 149 | + |
| 150 | +def tridiagonal_lu_solve(a_diagonals, b, *, b_ndim: int, transposed: bool = False): |
| 151 | + dl, d, du, du2, ipiv = a_diagonals |
| 152 | + return Blockwise(SolveLUFactorTridiagonal(b_ndim=b_ndim, transposed=transposed))( |
| 153 | + dl, d, du, du2, ipiv, b |
| 154 | + ) |
0 commit comments