Skip to content

quatica.decomp.chol

quatica.decomp.chol

Quaternion Cholesky decompositions.

This module adds
  • Native dense quaternion Cholesky for Hermitian PD matrices (no embedding)
  • Sparse Cholesky via complex embedding (CHOLMOD via scikit-sparse, optional)

QuatSparseCholeskyFactor dataclass

logdet()

Return (\log\det(A)) for quaternion Hermitian PD (A) (when supported).

For Hermitian quaternion PD matrices, the complex embedding satisfies (\det(\chi(A)) = \det(A)^2) (Moore determinant relationship), hence: (\log\det(A) = \tfrac12 \log\det(\chi(A))).

solve(b)

Solve (Ax=b) for quaternion b using the embedded complex factorization.

chol_quat_dense(A, *, tol=1e-12, hermitianize=False, jitter=0.0)

Compute a dense quaternion Cholesky factorization (A = L L^*).

This is a native quaternion implementation for Hermitian positive definite (HPD) matrices, without any real/complex embedding.

The algorithm is the quaternion analogue of classical left-looking Cholesky, with the key property that each pivot is real:

  • Pivot (real): (s_k = a_{kk} - \sum_{j<k} L_{kj}\overline{L_{kj}} = a_{kk} - \sum_{j0}).
  • Column update: (t_{ik} = a_{ik} - \sum_{j<k} L_{ij}\overline{L_{kj}}), then (L_{ik} = t_{ik}/L_{kk}) (safe since (L_{kk}\) is real).

Complexity is (O(n^3)). The inner accumulation for the column update is implemented with NumPy vectorization over the row index (i) (and reduction over (j)), avoiding Python loops over (i,j).

Parameters:

Name Type Description Default
A ndarray

Dense quaternion Hermitian matrix of shape (n, n) with dtype np.quaternion.

required
tol float

Pivot tolerance. If any pivot (s_k \le tol), the matrix is treated as not HPD and a LinAlgError is raised.

1e-12
hermitianize bool

If True, symmetrize input as 0.5*(A + A^*) before factorization (useful when small numerical asymmetry is present).

False
jitter float

Optional diagonal shift added to the real diagonal of A before factorization (common stabilization trick).

0.0

Returns:

Name Type Description
L ndarray

Dense quaternion lower-triangular matrix (n, n) such that

ndarray

(A \approx L L^*), with real positive diagonal.

Raises:

Type Description
ValueError

If A is not a square 2D quaternion array, or if diagonal entries are not approximately real.

LinAlgError

If a non-positive pivot is encountered.

chol_quat_sparse(Aq, *, tol=1e-12, jitter=0.0, ordering='cholmod')

Sparse quaternion Cholesky via complex embedding (CHOLMOD backend).

This factors a sparse quaternion Hermitian PD matrix (A\) by embedding it into its (2n\times 2n) complex adjoint/symplectic representation:

[ \chi(A)=\begin{bmatrix}X & Y\\-\overline{Y} & \overline{X}\end{bmatrix}, \quad X = A_w + i A_x,\; Y = A_y + i A_z. ]

Then it calls CHOLMOD on the complex sparse matrix (\chi(A)). Solves are performed by packing quaternion right-hand sides (b) into a complex vector ([u; v]) (with (u=b_w + i b_x), (v=b_y + i b_z)).

Notes
  • This requires the optional dependency scikit-sparse (CHOLMOD).
  • The ordering="paired" option is currently a placeholder hook; it still uses CHOLMOD's internal ordering.

Parameters:

Name Type Description Default
Aq SparseQuaternionMatrix

Sparse quaternion matrix in component form (SparseQuaternionMatrix) with shape (n, n).

required
tol float

Reserved for future numerical checks (kept for API stability).

1e-12
jitter float

Optional diagonal shift added to (\chi(A)) before factoring.

0.0
ordering Literal['cholmod', 'paired']

"cholmod" (default) or "paired" (hook for lifted ordering).

'cholmod'

Returns:

Type Description
QuatSparseCholeskyFactor

A QuatSparseCholeskyFactor object exposing .solve(b) and (optionally)

QuatSparseCholeskyFactor

.logdet().

Raises:

Type Description
ValueError

If Aq is not square or not a SparseQuaternionMatrix.

ImportError

If scikit-sparse / CHOLMOD is not available.

solve_chol_quat_dense(L, b)

Solve a quaternion linear system using a Cholesky factor.

Solves (Ax=b) given (A = L L^*) where L is the output of chol_quat_dense.

This uses forward/back substitution: - Forward: solve (Ly=b) - Backward: solve (L^x=y) where ((L^){ij} = \overline{L{ji}})

Parameters:

Name Type Description Default
L ndarray

Lower-triangular quaternion matrix of shape (n, n), with real positive diagonal, typically returned by chol_quat_dense.

required
b ndarray

Quaternion right-hand side vector (n,) or matrix (n, nrhs).

required

Returns:

Name Type Description
x ndarray

Quaternion solution with same shape as b.

Raises:

Type Description
ValueError

If inputs have incompatible shapes/dtypes.

LinAlgError

If L has a non-positive diagonal.