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