For any \(n \times m\) (typically) sparse matrix x compute the Dulmage-Mendelsohn row and columns permutations which at first splits the \(n\) rows and m columns into coarse partitions each; and then a finer one, reordering rows and columns such that the permutated matrix is “as upper triangular” as possible.

dmperm(x, nAns = 6L, seed = 0L)

Arguments

x

a typically sparse matrix; internally coerced to either "dgCMatrix" or "dtCMatrix".

nAns

an integer specifying the length of the resulting list. Must be 2, 4, or 6.

seed

an integer code in -1,0,1; determining the (initial) permutation; by default, seed = 0, no (or the identity) permutation; seed = -1 uses the “reverse” permutation k:1; for seed = 1, it is a random permutation (using R's RNG, seed, etc).

Details

See the book section by Tim Davis; page 122–127, in the References.

Value

a named list with (by default) 6 components,

p

integer vector with the permutation p, of length nrow(x).

q

integer vector with the permutation q, of length ncol(x).

r

integer vector of length nb+1, where block k is rows r[k] to r[k+1]-1 in A[p,q].

s

integer vector of length nb+1, where block k is cols s[k] to s[k+1]-1 in A[p,q].

rr5

integer vector of length 5, defining the coarse row decomposition.

cc5

integer vector of length 5, defining the coarse column decomposition.

References

Section 7.4 Dulmage-Mendelsohn decomposition, pp. 122 ff of
Timothy A. Davis (2006) Direct Methods for Sparse Linear Systems, SIAM Series “Fundamentals of Algorithms”.

Author

Martin Maechler, with a lot of “encouragement” by Mauricio Vargas.

See also

Schur, the class of permutation matrices; "pMatrix".

Examples

set.seed(17)
(S9 <- rsparsematrix(9, 9, nnz = 10, symmetric=TRUE)) # dsCMatrix
#> 9 x 9 sparse Matrix of class "dsCMatrix"
#>                                                   
#>  [1,] -0.055 .    .    .    .   .    1.6 .    .   
#>  [2,]  .     .    .    .    .   0.84 .   0.64 .   
#>  [3,]  .     .    0.63 .    .   .    .   .    1.30
#>  [4,]  .     .    .    0.16 .   .    .   .    0.37
#>  [5,]  .     .    .    .    .   .    .   1.20 .   
#>  [6,]  .     0.84 .    .    .   .    .   .    0.19
#>  [7,]  1.600 .    .    .    .   .    .   .    .   
#>  [8,]  .     0.64 .    .    1.2 .    .   .    .   
#>  [9,]  .     .    1.30 0.37 .   0.19 .   .    .   
str( dm9 <- dmperm(S9) )
#> List of 6
#>  $ p  : int [1:9] 1 7 8 6 9 4 3 2 5
#>  $ q  : int [1:9] 7 1 5 2 3 4 9 6 8
#>  $ r  : int [1:8] 0 1 2 3 4 7 8 9
#>  $ s  : int [1:8] 0 1 2 3 4 7 8 9
#>  $ rr5: int [1:5] 0 0 9 9 9
#>  $ cc5: int [1:5] 0 0 0 9 9
(S9p <- with(dm9, S9[p, q]))
#> 9 x 9 sparse Matrix of class "dgCMatrix"
#>                                                   
#>  [1,] 1.6 -0.055 .   .    .    .    .    .    .   
#>  [2,] .    1.600 .   .    .    .    .    .    .   
#>  [3,] .    .     1.2 0.64 .    .    .    .    .   
#>  [4,] .    .     .   0.84 .    .    0.19 .    .   
#>  [5,] .    .     .   .    1.30 0.37 .    0.19 .   
#>  [6,] .    .     .   .    .    0.16 0.37 .    .   
#>  [7,] .    .     .   .    0.63 .    1.30 .    .   
#>  [8,] .    .     .   .    .    .    .    0.84 0.64
#>  [9,] .    .     .   .    .    .    .    .    1.20
## looks good, but *not* quite upper triangular; these, too:
str( dm9.0 <- dmperm(S9, seed=-1)) # non-random too.
#> List of 6
#>  $ p  : int [1:9] 1 7 8 6 9 4 3 2 5
#>  $ q  : int [1:9] 7 1 5 2 3 4 9 6 8
#>  $ r  : int [1:8] 0 1 2 3 4 7 8 9
#>  $ s  : int [1:8] 0 1 2 3 4 7 8 9
#>  $ rr5: int [1:5] 0 0 9 9 9
#>  $ cc5: int [1:5] 0 0 0 9 9
str( dm9_1 <- dmperm(S9, seed= 1)) # a random one
#> List of 6
#>  $ p  : int [1:9] 1 7 8 6 3 9 4 2 5
#>  $ q  : int [1:9] 7 1 5 2 3 4 9 6 8
#>  $ r  : int [1:8] 0 1 2 3 4 7 8 9
#>  $ s  : int [1:8] 0 1 2 3 4 7 8 9
#>  $ rr5: int [1:5] 0 0 9 9 9
#>  $ cc5: int [1:5] 0 0 0 9 9
## The last two permutations differ, but have the same effect!
(S9p0 <- with(dm9.0, S9[p, q])) # .. hmm ..
#> 9 x 9 sparse Matrix of class "dgCMatrix"
#>                                                   
#>  [1,] 1.6 -0.055 .   .    .    .    .    .    .   
#>  [2,] .    1.600 .   .    .    .    .    .    .   
#>  [3,] .    .     1.2 0.64 .    .    .    .    .   
#>  [4,] .    .     .   0.84 .    .    0.19 .    .   
#>  [5,] .    .     .   .    1.30 0.37 .    0.19 .   
#>  [6,] .    .     .   .    .    0.16 0.37 .    .   
#>  [7,] .    .     .   .    0.63 .    1.30 .    .   
#>  [8,] .    .     .   .    .    .    .    0.84 0.64
#>  [9,] .    .     .   .    .    .    .    .    1.20
stopifnot(all.equal(S9p0, S9p))# same as as default, but different from the random one


set.seed(11)
(M <- triu(rsparsematrix(9,11, 1/4)))
#> 9 x 11 sparse Matrix of class "dgCMatrix"
#>                                                         
#>  [1,] .  0.89 . .  0.0072  .     .    .    .  .     .   
#>  [2,] . -0.34 . .  .      -0.61 -0.44 .    .  .    -0.94
#>  [3,] .  .    . . -1.6000  .     .    .    .  .     0.89
#>  [4,] .  .    . .  1.5000 -0.19  .    0.02 . -0.22  .   
#>  [5,] .  .    . .  .       .     .    .    . -0.98  .   
#>  [6,] .  .    . .  .       .    -0.19 .    . -0.35  .   
#>  [7,] .  .    . .  .       .     .    .    .  .     .   
#>  [8,] .  .    . .  .       .     .    .    .  .     .   
#>  [9,] .  .    . .  .       .     .    .    .  .     0.48
dM <- dmperm(M); with(dM, M[p, q])
#> 9 x 11 sparse Matrix of class "dgCMatrix"
#>                                                         
#>  [1,] . . . . 0.02 -0.19  .     1.5000  .    -0.22  .   
#>  [2,] . . . . .    -0.61 -0.34  .      -0.44  .    -0.94
#>  [3,] . . . . .     .     0.89  0.0072  .     .     .   
#>  [4,] . . . . .     .     .    -1.6000  .     .     0.89
#>  [5,] . . . . .     .     .     .      -0.19 -0.35  .   
#>  [6,] . . . . .     .     .     .       .    -0.98  .   
#>  [7,] . . . . .     .     .     .       .     .     0.48
#>  [8,] . . . . .     .     .     .       .     .     .   
#>  [9,] . . . . .     .     .     .       .     .     .   
(Mp <- M[sample.int(nrow(M)), sample.int(ncol(M))])
#> 9 x 11 sparse Matrix of class "dgCMatrix"
#>                                                         
#>  [1,]  0.89 -1.6000 . . .     .    . .  .     .     .   
#>  [2,]  0.48  .      . . .     .    . .  .     .     .   
#>  [3,]  .     1.5000 . . 0.02  .    . . -0.22 -0.19  .   
#>  [4,]  .     .      . . .     .    . .  .     .     .   
#>  [5,]  .     .      . . .     .    . . -0.35  .    -0.19
#>  [6,]  .     .      . . .     .    . . -0.98  .     .   
#>  [7,]  .     .      . . .     .    . .  .     .     .   
#>  [8,] -0.94  .      . . .    -0.34 . .  .    -0.61 -0.44
#>  [9,]  .     0.0072 . . .     0.89 . .  .     .     .   
dMp <- dmperm(Mp); with(dMp, Mp[p, q])
#> 9 x 11 sparse Matrix of class "dgCMatrix"
#>                                                         
#>  [1,] . . . . 0.02 -0.19  .     1.5000  .     .    -0.22
#>  [2,] . . . . .    -0.61 -0.34  .      -0.94 -0.44  .   
#>  [3,] . . . . .     .     0.89  0.0072  .     .     .   
#>  [4,] . . . . .     .     .    -1.6000  0.89  .     .   
#>  [5,] . . . . .     .     .     .       0.48  .     .   
#>  [6,] . . . . .     .     .     .       .    -0.19 -0.35
#>  [7,] . . . . .     .     .     .       .     .    -0.98
#>  [8,] . . . . .     .     .     .       .     .     .   
#>  [9,] . . . . .     .     .     .       .     .     .   


set.seed(7)
(n7 <- rsparsematrix(5, 12, nnz = 10, rand.x = NULL))
#> 5 x 12 sparse Matrix of class "ngCMatrix"
#>                             
#> [1,] . . . . . | | . . . . .
#> [2,] | . . . . . . . | . . .
#> [3,] | | . . . | . . . . . .
#> [4,] . . . | . . . | . . . .
#> [5,] . . | . . . . . . . . .
str( dm.7 <- dmperm(n7) )
#> List of 6
#>  $ p  : int [1:5] 2 3 4 1 5
#>  $ q  : int [1:12] 5 7 8 9 10 11 12 1 2 4 ...
#>  $ r  : int [1:3] 0 4 5
#>  $ s  : int [1:3] 0 11 12
#>  $ rr5: int [1:5] 0 4 5 5 5
#>  $ cc5: int [1:5] 0 7 11 12 12
stopifnot(exprs = {
  lengths(dm.7[1:2]) == dim(n7)
  identical(dm.7,      dmperm(as(n7, "dMatrix")))
  identical(dm.7[1:4], dmperm(n7, nAns=4))
  identical(dm.7[1:2], dmperm(n7, nAns=2))
})