Actual source code: dsm.c
1: #define PETSCMAT_DLL
2: /* dsm.f -- translated by f2c (version of 25 March 1992 12:58:56). */
4: #include ../src/mat/color/color.h
6: static PetscInt c_n1 = -1;
10: PetscErrorCode MINPACKdsm(PetscInt *m,PetscInt *n,PetscInt *npairs,PetscInt *indrow,PetscInt *indcol,PetscInt *ngrp,PetscInt *maxgrp,
11: PetscInt *mingrp,PetscInt *info,PetscInt *ipntr,PetscInt *jpntr,PetscInt *iwa,PetscInt *liwa)
12: {
13: /* System generated locals */
14: PetscInt i__1,i__2,i__3;
16: /* Local variables */
17: PetscInt i,j,maxclq,numgrp;
19: /* Given the sparsity pattern of an m by n matrix A, this */
20: /* subroutine determines a partition of the columns of A */
21: /* consistent with the direct determination of A. */
22: /* The sparsity pattern of the matrix A is specified by */
23: /* the arrays indrow and indcol. On input the indices */
24: /* for the non-zero elements of A are */
25: /* indrow(k),indcol(k), k = 1,2,...,npairs. */
26: /* The (indrow,indcol) pairs may be specified in any order. */
27: /* Duplicate input pairs are permitted, but the subroutine */
28: /* eliminates them. */
29: /* The subroutine partitions the columns of A into groups */
30: /* such that columns in the same group do not have a */
31: /* non-zero in the same row position. A partition of the */
32: /* columns of A with this property is consistent with the */
33: /* direct determination of A. */
34: /* The subroutine statement is */
35: /* subroutine dsm(m,n,npairs,indrow,indcol,ngrp,maxgrp,mingrp, */
36: /* info,ipntr,jpntr,iwa,liwa) */
37: /* where */
38: /* m is a positive integer input variable set to the number */
39: /* of rows of A. */
40: /* n is a positive integer input variable set to the number */
41: /* of columns of A. */
42: /* npairs is a positive integer input variable set to the */
43: /* number of (indrow,indcol) pairs used to describe the */
44: /* sparsity pattern of A. */
45: /* indrow is an integer array of length npairs. On input indrow */
46: /* must contain the row indices of the non-zero elements of A. */
47: /* On output indrow is permuted so that the corresponding */
48: /* column indices are in non-decreasing order. The column */
49: /* indices can be recovered from the array jpntr. */
50: /* indcol is an integer array of length npairs. On input indcol */
51: /* must contain the column indices of the non-zero elements of */
52: /* A. On output indcol is permuted so that the corresponding */
53: /* row indices are in non-decreasing order. The row indices */
54: /* can be recovered from the array ipntr. */
55: /* ngrp is an integer output array of length n which specifies */
56: /* the partition of the columns of A. Column jcol belongs */
57: /* to group ngrp(jcol). */
58: /* maxgrp is an integer output variable which specifies the */
59: /* number of groups in the partition of the columns of A. */
60: /* mingrp is an integer output variable which specifies a lower */
61: /* bound for the number of groups in any consistent partition */
62: /* of the columns of A. */
63: /* info is an integer output variable set as follows. For */
64: /* normal termination info = 1. If m, n, or npairs is not */
65: /* positive or liwa is less than max(m,6*n), then info = 0. */
66: /* If the k-th element of indrow is not an integer between */
67: /* 1 and m or the k-th element of indcol is not an integer */
68: /* between 1 and n, then info = -k. */
69: /* ipntr is an integer output array of length m + 1 which */
70: /* specifies the locations of the column indices in indcol. */
71: /* The column indices for row i are */
72: /* indcol(k), k = ipntr(i),...,ipntr(i+1)-1. */
73: /* Note that ipntr(m+1)-1 is then the number of non-zero */
74: /* elements of the matrix A. */
75: /* jpntr is an integer output array of length n + 1 which */
76: /* specifies the locations of the row indices in indrow. */
77: /* The row indices for column j are */
78: /* indrow(k), k = jpntr(j),...,jpntr(j+1)-1. */
79: /* Note that jpntr(n+1)-1 is then the number of non-zero */
80: /* elements of the matrix A. */
81: /* iwa is an integer work array of length liwa. */
82: /* liwa is a positive integer input variable not less than */
83: /* max(m,6*n). */
84: /* Subprograms called */
85: /* MINPACK-supplied ... degr,ido,numsrt,seq,setr,slo,srtdat */
86: /* FORTRAN-supplied ... max */
87: /* Argonne National Laboratory. MINPACK Project. December 1984. */
88: /* Thomas F. Coleman, Burton S. Garbow, Jorge J. More' */
91: /* Parameter adjustments */
92: --iwa;
93: --jpntr;
94: --ipntr;
95: --ngrp;
96: --indcol;
97: --indrow;
99: *info = 0;
101: /* Determine a lower bound for the number of groups. */
103: *mingrp = 0;
104: i__1 = *m;
105: for (i = 1; i <= i__1; ++i) {
106: /* Computing MAX */
107: i__2 = *mingrp,i__3 = ipntr[i + 1] - ipntr[i];
108: *mingrp = PetscMax(i__2,i__3);
109: }
111: /* Determine the degree sequence for the intersection */
112: /* graph of the columns of A. */
114: MINPACKdegr(n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[*n * 5 + 1],&
115: iwa[*n + 1]);
117: /* Color the intersection graph of the columns of A */
118: /* with the smallest-last (SL) ordering. */
120: MINPACKslo(n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[*n * 5 + 1],&
121: iwa[(*n << 2) + 1],&maxclq,&iwa[1],&iwa[*n + 1],&iwa[(*n << 1)
122: + 1],&iwa[*n * 3 + 1]);
123: MINPACKseq(n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[(*n << 2) + 1],
124: &ngrp[1],maxgrp,&iwa[*n + 1]);
125: *mingrp = PetscMax(*mingrp,maxclq);
127: /* Exit if the smallest-last ordering is optimal. */
129: if (*maxgrp == *mingrp) {
130: return(0);
131: }
133: /* Color the intersection graph of the columns of A */
134: /* with the incidence-degree (ID) ordering. */
136: MINPACKido(m,n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[*n * 5 + 1],
137: &iwa[(*n << 2) + 1],&maxclq,&iwa[1],&iwa[*n + 1],&iwa[(*n <<
138: 1) + 1],&iwa[*n * 3 + 1]);
139: MINPACKseq(n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[(*n << 2) + 1],
140: &iwa[1],&numgrp,&iwa[*n + 1]);
141: *mingrp = PetscMax(*mingrp,maxclq);
143: /* Retain the better of the two orderings so far. */
145: if (numgrp < *maxgrp) {
146: *maxgrp = numgrp;
147: i__1 = *n;
148: for (j = 1; j <= i__1; ++j) {
149: ngrp[j] = iwa[j];
150: }
152: /* Exit if the incidence-degree ordering is optimal. */
154: if (*maxgrp == *mingrp) {
155: return(0);
156: }
157: }
159: /* Color the intersection graph of the columns of A */
160: /* with the largest-first (LF) ordering. */
162: i__1 = *n - 1;
163: MINPACKnumsrt(n,&i__1,&iwa[*n * 5 + 1],&c_n1,&iwa[(*n << 2) + 1],&iwa[(*n
164: << 1) + 1],&iwa[*n + 1]);
165: MINPACKseq(n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[(*n << 2) + 1],
166: &iwa[1],&numgrp,&iwa[*n + 1]);
168: /* Retain the best of the three orderings and exit. */
170: if (numgrp < *maxgrp) {
171: *maxgrp = numgrp;
172: i__1 = *n;
173: for (j = 1; j <= i__1; ++j) {
174: ngrp[j] = iwa[j];
175: }
176: }
177: return(0);
178: }