Actual source code: qmdmrg.c
2: /* qmdmrg.f -- translated by f2c (version 19931217).*/
4: #include <petscsys.h>
5: #include <petsc/private/matorderimpl.h>
7: /******************************************************************/
8: /*********** QMDMRG ..... QUOT MIN DEG MERGE ************/
9: /******************************************************************/
10: /* PURPOSE - THIS ROUTINE MERGES INDISTINGUISHABLE NODES IN */
11: /* THE MINIMUM DEGREE ORDERING ALGORITHM. */
12: /* IT ALSO COMPUTES THE NEW DEGREES OF THESE */
13: /* NEW SUPERNODES. */
14: /* */
15: /* INPUT PARAMETERS - */
16: /* (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE. */
17: /* DEG0 - THE NUMBER OF NODES IN THE GIVEN SET. */
18: /* (NHDSZE, NBRHD) - THE SET OF ELIMINATED SUPERNODES */
19: /* ADJACENT TO SOME NODES IN THE SET. */
20: /* */
21: /* UPDATED PARAMETERS - */
22: /* DEG - THE DEGREE VECTOR. */
23: /* QSIZE - SIZE OF INDISTINGUISHABLE NODES. */
24: /* QLINK - LINKED LIST FOR INDISTINGUISHABLE NODES. */
25: /* MARKER - THE GIVEN SET IS GIVEN BY THOSE NODES WITH */
26: /* MARKER VALUE SET TO 1. THOSE NODES WITH DEGREE */
27: /* UPDATED WILL HAVE MARKER VALUE SET TO 2. */
28: /* */
29: /* WORKING PARAMETERS - */
30: /* RCHSET - THE REACHABLE SET. */
31: /* OVRLP - TEMP VECTOR TO STORE THE INTERSECTION OF TWO */
32: /* REACHABLE SETS. */
33: /* */
34: /*****************************************************************/
35: PetscErrorCode SPARSEPACKqmdmrg(const PetscInt *xadj, const PetscInt *adjncy, PetscInt *deg, PetscInt *qsize, PetscInt *qlink, PetscInt *marker, PetscInt *deg0, PetscInt *nhdsze, PetscInt *nbrhd, PetscInt *rchset, PetscInt *ovrlp)
36: {
37: /* System generated locals */
38: PetscInt i__1, i__2, i__3;
40: /* Local variables */
41: PetscInt head, inhd, irch, node, mark, ilink, root, j, lnode, nabor, jstop, jstrt, rchsze, mrgsze, novrlp, iov, deg1;
43: PetscFunctionBegin;
44: /* Parameter adjustments */
45: --ovrlp;
46: --rchset;
47: --nbrhd;
48: --marker;
49: --qlink;
50: --qsize;
51: --deg;
52: --adjncy;
53: --xadj;
55: if (*nhdsze <= 0) PetscFunctionReturn(PETSC_SUCCESS);
56: i__1 = *nhdsze;
57: for (inhd = 1; inhd <= i__1; ++inhd) {
58: root = nbrhd[inhd];
59: marker[root] = 0;
60: }
61: /* LOOP THROUGH EACH ELIMINATED SUPERNODE IN THE SET */
62: /* (NHDSZE, NBRHD). */
63: i__1 = *nhdsze;
64: for (inhd = 1; inhd <= i__1; ++inhd) {
65: root = nbrhd[inhd];
66: marker[root] = -1;
67: rchsze = 0;
68: novrlp = 0;
69: deg1 = 0;
70: L200:
71: jstrt = xadj[root];
72: jstop = xadj[root + 1] - 1;
73: /* DETERMINE THE REACHABLE SET AND ITS PETSCINTERSECT- */
74: /* ION WITH THE INPUT REACHABLE SET. */
75: i__2 = jstop;
76: for (j = jstrt; j <= i__2; ++j) {
77: nabor = adjncy[j];
78: root = -nabor;
79: if (nabor < 0) goto L200;
80: else if (!nabor) goto L700;
81: else goto L300;
82: L300:
83: mark = marker[nabor];
84: if (mark < 0) goto L600;
85: else if (!mark) goto L400;
86: else goto L500;
87: L400:
88: ++rchsze;
89: rchset[rchsze] = nabor;
90: deg1 += qsize[nabor];
91: marker[nabor] = 1;
92: goto L600;
93: L500:
94: if (mark > 1) goto L600;
95: ++novrlp;
96: ovrlp[novrlp] = nabor;
97: marker[nabor] = 2;
98: L600:;
99: }
100: /* FROM THE OVERLAPPED SET, DETERMINE THE NODES */
101: /* THAT CAN BE MERGED TOGETHER. */
102: L700:
103: head = 0;
104: mrgsze = 0;
105: i__2 = novrlp;
106: for (iov = 1; iov <= i__2; ++iov) {
107: node = ovrlp[iov];
108: jstrt = xadj[node];
109: jstop = xadj[node + 1] - 1;
110: i__3 = jstop;
111: for (j = jstrt; j <= i__3; ++j) {
112: nabor = adjncy[j];
113: if (marker[nabor] != 0) goto L800;
114: marker[node] = 1;
115: goto L1100;
116: L800:;
117: }
118: /* NODE BELONGS TO THE NEW MERGED SUPERNODE. */
119: /* UPDATE THE VECTORS QLINK AND QSIZE. */
120: mrgsze += qsize[node];
121: marker[node] = -1;
122: lnode = node;
123: L900:
124: ilink = qlink[lnode];
125: if (ilink <= 0) goto L1000;
126: lnode = ilink;
127: goto L900;
128: L1000:
129: qlink[lnode] = head;
130: head = node;
131: L1100:;
132: }
133: if (head <= 0) goto L1200;
134: qsize[head] = mrgsze;
135: deg[head] = *deg0 + deg1 - 1;
136: marker[head] = 2;
137: /* RESET MARKER VALUES. */
138: L1200:
139: root = nbrhd[inhd];
140: marker[root] = 0;
141: if (rchsze <= 0) goto L1400;
142: i__2 = rchsze;
143: for (irch = 1; irch <= i__2; ++irch) {
144: node = rchset[irch];
145: marker[node] = 0;
146: }
147: L1400:;
148: }
149: PetscFunctionReturn(PETSC_SUCCESS);
150: }