Quick startΒΆ
The core functionality of CHOMPACK is contained in two types of
objects: the symbolic object and the cspmatrix
(chordal sparse matrix) object. A symbolic object
represents a symbolic factorization of a sparse symmetric matrix
, and it can be created as follows:
from cvxopt import spmatrix, amd
import chompack as cp
# generate sparse matrix
I = [0, 1, 3, 1, 5, 2, 6, 3, 4, 5, 4, 5, 6, 5, 6]
J = [0, 0, 0, 1, 1, 2, 2, 3, 3, 3, 4, 4, 4, 5, 6]
A = spmatrix(1.0, I, J, (7,7))
# compute symbolic factorization using AMD ordering
symb = cp.symbolic(A, p=amd.order)
The argument is a so-called elimination order, and it can
be either an ordering routine or a permutation vector. In the above
example we used the “approximate minimum degree” (AMD) ordering
routine. Note that
is a lower-triangular sparse matrix that represents a
symmetric matrix; upper-triangular entries in
are ignored in the
symbolic factorization.
Now let’s inspect the sparsity pattern of A and its chordal embedding (i.e., the filled pattern):
>>> print(A)
[ 1.00e+00 0 0 0 0 0 0 ]
[ 1.00e+00 1.00e+00 0 0 0 0 0 ]
[ 0 0 1.00e+00 0 0 0 0 ]
[ 1.00e+00 0 0 1.00e+00 0 0 0 ]
[ 0 0 0 1.00e+00 1.00e+00 0 0 ]
[ 0 1.00e+00 0 1.00e+00 1.00e+00 1.00e+00 0 ]
[ 0 0 1.00e+00 0 1.00e+00 0 1.00e+00]
>>> print(symb.sparsity_pattern(reordered=False, symmetric=False))
[ 1.00e+00 0 0 0 0 0 0 ]
[ 1.00e+00 1.00e+00 0 0 0 0 0 ]
[ 0 0 1.00e+00 0 0 0 0 ]
[ 1.00e+00 0 0 1.00e+00 0 0 0 ]
[ 0 0 0 1.00e+00 1.00e+00 0 0 ]
[ 1.00e+00 1.00e+00 0 1.00e+00 1.00e+00 1.00e+00 0 ]
[ 0 0 1.00e+00 0 1.00e+00 0 1.00e+00]
The reordered pattern and its cliques can be inspected using the following commands:
>>> print(symb)
[X X ]
[X X X ]
[ X X X X ]
[ X X X]
[ X X X X]
[ X X X X X]
[ X X X X]
>>> print(symb.cliques())
[[0, 1], [1, 2], [2, 4, 5], [3, 5, 6], [4, 5, 6]]
Similarly, the clique tree, the supernodes, and the separator sets are:
>>> print(symb.parent())
[1, 2, 4, 4, 4]
>>> print(symb.supernodes())
[[0], [1], [2], [3], [4, 5, 6]]
>>> print(symb.separators())
[[1], [2], [4, 5], [5, 6], []]
The cspmatrix object represents a chordal sparse matrix,
and it contains lower-triangular numerical values as well as a
reference to a symbolic factorization that defines the sparsity
pattern. Given a symbolic object symb and a sparse
matrix , we can create a cspmatrix as follows:
from cvxopt import spmatrix, amd, printing
import chompack as cp
printing.options['dformat'] = '%3.1f'
# generate sparse matrix and compute symbolic factorization
I = [0, 1, 3, 1, 5, 2, 6, 3, 4, 5, 4, 5, 6, 5, 6]
J = [0, 0, 0, 1, 1, 2, 2, 3, 3, 3, 4, 4, 4, 5, 6]
A = spmatrix([1.0*i for i in range(1,15+1)], I, J, (7,7))
symb = cp.symbolic(A, p=amd.order)
L = cp.cspmatrix(symb)
L += A
Now let us take a look at and
:
>>> print(A)
[ 1.0 0 0 0 0 0 0 ]
[ 2.0 4.0 0 0 0 0 0 ]
[ 0 0 6.0 0 0 0 0 ]
[ 3.0 0 0 8.0 0 0 0 ]
[ 0 0 0 9.0 11.0 0 0 ]
[ 0 5.0 0 10.0 12.0 14.0 0 ]
[ 0 0 7.0 0 13.0 0 15.0]
>>> print(L)
[ 6.0 0 0 0 0 0 0 ]
[ 7.0 15.0 0 0 0 0 0 ]
[ 0 13.0 11.0 0 0 0 0 ]
[ 0 0 0 4.0 0 0 0 ]
[ 0 0 9.0 0 8.0 0 0 ]
[ 0 0 12.0 5.0 10.0 14.0 0 ]
[ 0 0 0 2.0 3.0 0.0 1.0]
Notice that is a reordered lower-triangular representation
of
. We can convert
to an spmatrix using
the spmatrix() method:
>>> print(L.spmatrix(reordered = False))
[ 1.0 0 0 0 0 0 0 ]
[ 2.0 4.0 0 0 0 0 0 ]
[ 0 0 6.0 0 0 0 0 ]
[ 3.0 0 0 8.0 0 0 0 ]
[ 0 0 0 9.0 11.0 0 0 ]
[ 0.0 5.0 0 10.0 12.0 14.0 0 ]
[ 0 0 7.0 0 13.0 0 15.0]
This returns an spmatrix with the same ordering
as , i.e., the inverse permutation is applied to
.
The following example illustrates how to use the Cholesky routine:
from cvxopt import spmatrix, amd, normal
from chompack import symbolic, cspmatrix, cholesky
# generate sparse matrix and compute symbolic factorization
I = [0, 1, 3, 1, 5, 2, 6, 3, 4, 5, 4, 5, 6, 5, 6]
J = [0, 0, 0, 1, 1, 2, 2, 3, 3, 3, 4, 4, 4, 5, 6]
A = spmatrix([0.1*(i+1) for i in range(15)], I, J, (7,7)) + spmatrix(10.0,range(7),range(7))
symb = symbolic(A, p=amd.order)
# create cspmatrix
L = cspmatrix(symb)
L += A
# compute numeric factorization
cholesky(L)
>>> print(L)
[ 3.3 0 0 0 0 0 0 ]
[ 0.2 3.4 0 0 0 0 0 ]
[ 0 0.4 3.3 0 0 0 0 ]
[ 0 0 0 3.2 0 0 0 ]
[ 0 0 0.3 0 3.3 0 0 ]
[ 0 0 0.4 0.2 0.3 3.3 0 ]
[ 0 0 0 0.1 0.1 -0.0 3.2]
Given a sparse matrix , we can check if it is chordal by
checking whether the permutation
returned by maximum cardinality
search is a perfect elimination ordering:
from cvxopt import spmatrix, printing
printing.options['width'] = -1
import chompack as cp
# Define chordal sparse matrix
I = range(17)+[2,2,3,3,4,14,4,14,8,14,15,8,15,7,8,14,8,14,14,\
15,10,12,13,16,12,13,16,12,13,15,16,13,15,16,15,16,15,16,16]
J = range(17)+[0,1,1,2,2,2,3,3,4,4,4,5,5,6,6,6,7,7,8,\
8,9,9,9,9,10,10,10,11,11,11,11,12,12,12,13,13,14,14,15]
A = spmatrix(1.0,I,J,(17,17))
# Compute maximum cardinality search
p = cp.maxcardsearch(A)
Is a perfect elimination ordering?
>>> cp.peo(A,p)
True
Let’s verify that no fill is generated by the symbolic factorization:
>>> symb = cp.symbolic(A,p)
>>> print(symb.fill)
(0, 0)