imate.logdet(method=’cholesky’)#

imate.logdet(A, gram=False, p=1.0, return_info=False, method='cholesky', cholmod=None)

Log-determinant of non-singular matrix using Cholesky method.

See also

This page describes only the cholesky method. For other methods, see imate.logdet().

Given the matrix \(\mathbf{A}\) and the real exponent \(p\), the following is computed:

\[\mathrm{logdet} \left(\mathbf{A}^p \right) = p \log_e \vert \det (\mathbf{A}) \vert.\]

If gram is True, then \(\mathbf{A}\) in the above is replaced by the Gramian matrix \(\mathbf{A}^{\intercal} \mathbf{A}\). In this case, if the matrix \(\mathbf{A}\) is square, then the following is instead computed:

\[\mathrm{logdet} \left((\mathbf{A}^{\intercal}\mathbf{A})^p \right) = 2p \log_e \vert \det (\mathbf{A}) \vert.\]
Parameters:
Anumpy.ndarray, scipy.sparse

A positive-definite sparse or dense matrix.

Warning

This function does not pre-check whether the input matrix is positive-definite.

Note

In the Cholesky method, the matrix cannot be a type of Matrix or imate.AffineMatrixFunction classes.

grambool, default=False

If True, the log-determinant of the Gramian matrix, \((\mathbf{A}^{\intercal}\mathbf{A})^p\), is computed. The Gramian matrix itself is not directly computed. If False, the log-determinant of \(\mathbf{A}^p\) is computed.

pfloat, default=1.0

The real exponent \(p\) in \(\mathbf{A}^p\).

return_infobool, default=False

If True, this function also returns a dictionary containing information about the inner computation, such as process time, algorithm settings, etc.

cholmodbool, default=None

If set to True, it uses the Cholmod library from scikit-sparse package to compute the Cholesky decomposition. If set to False, it uses scipy.sparse.cholesky method. If set to None, first, it tries to use Cholmod library, but if Cholmod is not available, then it uses scipy.sparse.cholesky method.

Returns:
logdetfloat or numpy.array

Log-determinant of A.

infodict

(Only if return_info is True) A dictionary of information with the following keys.

  • matrix:
    • data_type: str, {float32, float64, float128}. Type of the matrix data.

    • gram: bool, whether the matrix A or its Gramian is considered.

    • exponent: float, the exponent p in \(\mathbf{A}^p\).

    • size: (int, int), the size of matrix A.

    • sparse: bool, whether the matrix A is sparse or dense.

    • nnz: int, if A is sparse, the number of non-zero elements of A.

    • density: float, if A is sparse, the density of A, which is the nnz divided by size squared.

    • num_inquiries: int, for the cholesky method, this is always 1.

  • device:
    • num_cpu_threads: int, number of CPU threads used in shared memory parallel processing.

    • num_gpu_devices: int, for the cholesky method, this is always 0.

    • num_gpu_multiprocessors: int, for the cholesky method, this is always 0.

    • num_gpu_threads_per_multiprocessor: int, for the cholesky method, this is always 0.

  • time:
    • tot_wall_time: float, total elapsed time of computation.

    • alg_wall_time: float, elapsed time of computation during only the algorithm execution.

    • cpu_proc_time: float, CPU processing time of computation.

  • solver:
    • cholmod_used: bool, whether the Cholmod from SparseSuite library was used.

    • version: str, version of imate.

    • method: ‘cholesky’

Notes

Algorithm:

The log-determinant is computed from the Cholesky decomposition \(\mathbf{A} = \mathbf{L} \mathbf{L}^{\intercal}\) as

\[\log | \mathbf{A} | = 2 \mathrm{trace}( \log \mathrm{diag}(\mathbf{L})).\]

The result is exact (no approximation) and could be used as a benchmark to test other methods.

Computational Complexity:

The computational complexity of this method is \(\mathcal{O}(\frac{1}{3}n^3)\) for dense matrices.

For sparse matrices obtained from 1D, 2D, and 3D meshes, the computational complexity of this method is respectively \(\mathcal{O}(n)\), \(\mathcal{O}(n^{\frac{3}{2}})\), and \(\mathcal{O}(n^2)\) (see for instance [1] and [2]).

Implementation:

This function is essentially a wrapper for the Cholesky function of the scipy and scikit-sparse packages and is primarily used for testing and comparison (benchmarking) against the randomized methods that are implemented in this package. If cholmod is set to True, this function uses the Suite Sparse package to compute the Cholesky decomposition.

References

[1]

George, A. and Ng, E. (1988). On the Complexity of Sparse QR and LU Factorization of Finite-Element Matrices. SIAM Journal on Scientific and Statistical Computing, volume 9, number 5, pp. 849-861. doi: 10.1137/0909057.

[2]

Davis, T. (2006). Direct Methods for Sparse Linear Systems. SIAM. doi: 10.1137/1.9780898718881.

Examples

Compute the log-determinant of a sparse positive-definite Toeplitz matrix:

>>> # Import packages
>>> from imate import toeplitz, logdet

>>> # Generate a sample symmetric and positive-definite matrix
>>> A = toeplitz(2, 1, size=100, gram=True)

>>> # Compute log-determinant with Cholesky method (default method)
>>> logdet(A, method='cholesky')
138.62943611198907

Print information about the inner computation:

>>> ld, info = logdet(A, method='cholesky', return_info=True)
>>> print(ld)
138.6294361119891

>>> # Print dictionary neatly using pprint
>>> from pprint import pprint
>>> pprint(info)
{
    'matrix': {
        'data_type': b'float64',
        'density': 0.0298,
        'exponent': 1.0,
        'gram': False,
        'nnz': 298,
        'num_inquiries': 1,
        'size': (100, 100),
        'sparse': True
    },
    'solver': {
        'cholmod_used': True,
        'method': 'cholesky',
        'version': '0.13.0'
    },
    'device': {
        'num_cpu_threads': 8,
        'num_gpu_devices': 0,
        'num_gpu_multiprocessors': 0,
        'num_gpu_threads_per_multiprocessor': 0
    },
    'time': {
        'alg_wall_time': 0.0007234140066429973,
        'cpu_proc_time': 0.0009358710000000325,
        'tot_wall_time': 0.0007234140066429973
    }
}