Actual source code: genqmd.c
2: /* genqmd.f -- translated by f2c (version 19931217).*/
4: #include <petscsys.h>
5: #include <petsc/private/matorderimpl.h>
7: /******************************************************************/
8: /*********** GENQMD ..... QUOT MIN DEGREE ORDERING **********/
9: /******************************************************************/
10: /* PURPOSE - THIS ROUTINE IMPLEMENTS THE MINIMUM DEGREE */
11: /* ALGORITHM. IT MAKES USE OF THE IMPLICIT REPRESENT- */
12: /* ATION OF THE ELIMINATION GRAPHS BY QUOTIENT GRAPHS, */
13: /* AND THE NOTION OF INDISTINGUISHABLE NODES. */
14: /* CAUTION - THE ADJACENCY VECTOR ADJNCY WILL BE */
15: /* DESTROYED. */
16: /* */
17: /* INPUT PARAMETERS - */
18: /* NEQNS - NUMBER OF EQUATIONS. */
19: /* (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE. */
20: /* */
21: /* OUTPUT PARAMETERS - */
22: /* PERM - THE MINIMUM DEGREE ORDERING. */
23: /* INVP - THE INVERSE OF PERM. */
24: /* */
25: /* WORKING PARAMETERS - */
26: /* DEG - THE DEGREE VECTOR. DEG(I) IS NEGATIVE MEANS */
27: /* NODE I HAS BEEN NUMBERED. */
28: /* MARKER - A MARKER VECTOR, WHERE MARKER(I) IS */
29: /* NEGATIVE MEANS NODE I HAS BEEN MERGED WITH */
30: /* ANOTHER NODE AND THUS CAN BE IGNORED. */
31: /* RCHSET - VECTOR USED FOR THE REACHABLE SET. */
32: /* NBRHD - VECTOR USED FOR THE NEIGHBORHOOD SET. */
33: /* QSIZE - VECTOR USED TO STORE THE SIZE OF */
34: /* INDISTINGUISHABLE SUPERNODES. */
35: /* QLINK - VECTOR TO STORE INDISTINGUISHABLE NODES, */
36: /* I, QLINK(I), QLINK(QLINK(I)) ... ARE THE */
37: /* MEMBERS OF THE SUPERNODE REPRESENTED BY I. */
38: /* */
39: /* PROGRAM SUBROUTINES - */
40: /* QMDRCH, QMDQT, QMDUPD. */
41: /* */
42: /******************************************************************/
43: /* */
44: /* */
45: PetscErrorCode SPARSEPACKgenqmd(const PetscInt *neqns,const PetscInt *xadj,const PetscInt *adjncy,
46: PetscInt *perm, PetscInt *invp, PetscInt *deg, PetscInt *marker,
47: PetscInt *rchset, PetscInt *nbrhd, PetscInt *qsize, PetscInt *qlink, PetscInt *nofsub)
48: {
49: /* System generated locals */
50: PetscInt i__1;
52: /* Local variables */
53: PetscInt ndeg, irch, node, nump1, j, inode;
54: PetscInt ip, np, mindeg, search;
55: PetscInt nhdsze, nxnode, rchsze, thresh, num;
57: /* INITIALIZE DEGREE VECTOR AND OTHER WORKING VARIABLES. */
59: /* Parameter adjustments */
60: --qlink;
61: --qsize;
62: --nbrhd;
63: --rchset;
64: --marker;
65: --deg;
66: --invp;
67: --perm;
68: --adjncy;
69: --xadj;
71: mindeg = *neqns;
72: *nofsub = 0;
73: i__1 = *neqns;
74: for (node = 1; node <= i__1; ++node) {
75: perm[node] = node;
76: invp[node] = node;
77: marker[node] = 0;
78: qsize[node] = 1;
79: qlink[node] = 0;
80: ndeg = xadj[node + 1] - xadj[node];
81: deg[node] = ndeg;
82: if (ndeg < mindeg) mindeg = ndeg;
83: }
84: num = 0;
85: /* PERFORM THRESHOLD SEARCH TO GET A NODE OF MIN DEGREE. */
86: /* VARIABLE SEARCH POINTS TO WHERE SEARCH SHOULD START. */
87: L200:
88: search = 1;
89: thresh = mindeg;
90: mindeg = *neqns;
91: L300:
92: nump1 = num + 1;
93: if (nump1 > search) search = nump1;
94: i__1 = *neqns;
95: for (j = search; j <= i__1; ++j) {
96: node = perm[j];
97: if (marker[node] < 0) goto L400;
98: ndeg = deg[node];
99: if (ndeg <= thresh) goto L500;
100: if (ndeg < mindeg) mindeg = ndeg;
101: L400:
102: ;
103: }
104: goto L200;
105: /* NODE HAS MINIMUM DEGREE. FIND ITS REACHABLE SETS BY */
106: /* CALLING QMDRCH. */
107: L500:
108: search = j;
109: *nofsub += deg[node];
110: marker[node] = 1;
111: SPARSEPACKqmdrch(&node, &xadj[1], &adjncy[1], °[1], &marker[1], &rchsze, &
112: rchset[1], &nhdsze, &nbrhd[1]);
113: /* ELIMINATE ALL NODES INDISTINGUISHABLE FROM NODE. */
114: /* THEY ARE GIVEN BY NODE, QLINK(NODE), .... */
115: nxnode = node;
116: L600:
117: ++num;
118: np = invp[nxnode];
119: ip = perm[num];
120: perm[np] = ip;
121: invp[ip] = np;
122: perm[num] = nxnode;
123: invp[nxnode] = num;
124: deg[nxnode] = -1;
125: nxnode = qlink[nxnode];
126: if (nxnode > 0) goto L600;
127: if (rchsze <= 0) goto L800;
129: /* UPDATE THE DEGREES OF THE NODES IN THE REACHABLE */
130: /* SET AND IDENTIFY INDISTINGUISHABLE NODES. */
131: SPARSEPACKqmdupd(&xadj[1], &adjncy[1], &rchsze, &rchset[1], °[1], &qsize[1], &
132: qlink[1], &marker[1], &rchset[rchsze + 1], &nbrhd[nhdsze + 1]);
134: /* RESET MARKER VALUE OF NODES IN REACH SET. */
135: /* UPDATE THRESHOLD VALUE FOR CYCLIC SEARCH. */
136: /* ALSO CALL QMDQT TO FORM NEW QUOTIENT GRAPH. */
137: marker[node] = 0;
138: i__1 = rchsze;
139: for (irch = 1; irch <= i__1; ++irch) {
140: inode = rchset[irch];
141: if (marker[inode] < 0) goto L700;
143: marker[inode] = 0;
144: ndeg = deg[inode];
145: if (ndeg < mindeg) mindeg = ndeg;
146: if (ndeg > thresh) goto L700;
147: mindeg = thresh;
148: thresh = ndeg;
149: search = invp[inode];
150: L700:
151: ;
152: }
153: if (nhdsze > 0) {
154: SPARSEPACKqmdqt(&node, &xadj[1], &adjncy[1], &marker[1], &rchsze, &rchset[1], &
155: nbrhd[1]);
156: }
157: L800:
158: if (num < *neqns) goto L300;
159: return 0;
160: }