Actual source code: qmdmrg.c
petsc-3.6.1 2015-08-06
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: /*****************************************************************/
37: PetscErrorCode SPARSEPACKqmdmrg(const PetscInt *xadj,const PetscInt *adjncy, PetscInt *deg,
38: PetscInt *qsize, PetscInt *qlink, PetscInt *marker, PetscInt *deg0,
39: PetscInt *nhdsze, PetscInt *nbrhd, PetscInt *rchset, PetscInt *ovrlp)
40: {
41: /* System generated locals */
42: PetscInt i__1, i__2, i__3;
44: /* Local variables */
45: PetscInt head, inhd, irch, node, mark, ilink, root, j, lnode, nabor,
46: jstop, jstrt, rchsze, mrgsze, novrlp, iov, deg1;
49: /* Parameter adjustments */
50: --ovrlp;
51: --rchset;
52: --nbrhd;
53: --marker;
54: --qlink;
55: --qsize;
56: --deg;
57: --adjncy;
58: --xadj;
60: if (*nhdsze <= 0) return(0);
61: i__1 = *nhdsze;
62: for (inhd = 1; inhd <= i__1; ++inhd) {
63: root = nbrhd[inhd];
64: marker[root] = 0;
65: }
66: /* LOOP THROUGH EACH ELIMINATED SUPERNODE IN THE SET */
67: /* (NHDSZE, NBRHD). */
68: i__1 = *nhdsze;
69: for (inhd = 1; inhd <= i__1; ++inhd) {
70: root = nbrhd[inhd];
71: marker[root] = -1;
72: rchsze = 0;
73: novrlp = 0;
74: deg1 = 0;
75: L200:
76: jstrt = xadj[root];
77: jstop = xadj[root + 1] - 1;
78: /* DETERMINE THE REACHABLE SET AND ITS PETSCINTERSECT- */
79: /* ION WITH THE INPUT REACHABLE SET. */
80: i__2 = jstop;
81: for (j = jstrt; j <= i__2; ++j) {
82: nabor = adjncy[j];
83: root = -nabor;
84: if (nabor < 0) goto L200;
85: else if (!nabor) goto L700;
86: else goto L300;
87: L300:
88: mark = marker[nabor];
89: if (mark < 0) goto L600;
90: else if (!mark) goto L400;
91: else goto L500;
92: L400:
93: ++rchsze;
94: rchset[rchsze] = nabor;
95: deg1 += qsize[nabor];
96: marker[nabor] = 1;
97: goto L600;
98: L500:
99: if (mark > 1) goto L600;
100: ++novrlp;
101: ovrlp[novrlp] = nabor;
102: marker[nabor] = 2;
103: L600:
104: ;
105: }
106: /* FROM THE OVERLAPPED SET, DETERMINE THE NODES */
107: /* THAT CAN BE MERGED TOGETHER. */
108: L700:
109: head = 0;
110: mrgsze = 0;
111: i__2 = novrlp;
112: for (iov = 1; iov <= i__2; ++iov) {
113: node = ovrlp[iov];
114: jstrt = xadj[node];
115: jstop = xadj[node + 1] - 1;
116: i__3 = jstop;
117: for (j = jstrt; j <= i__3; ++j) {
118: nabor = adjncy[j];
119: if (marker[nabor] != 0) goto L800;
120: marker[node] = 1;
121: goto L1100;
122: L800:
123: ;
124: }
125: /* NODE BELONGS TO THE NEW MERGED SUPERNODE. */
126: /* UPDATE THE VECTORS QLINK AND QSIZE. */
127: mrgsze += qsize[node];
128: marker[node] = -1;
129: lnode = node;
130: L900:
131: ilink = qlink[lnode];
132: if (ilink <= 0) goto L1000;
133: lnode = ilink;
134: goto L900;
135: L1000:
136: qlink[lnode] = head;
137: head = node;
138: L1100:
139: ;
140: }
141: if (head <= 0) goto L1200;
142: qsize[head] = mrgsze;
143: deg[head] = *deg0 + deg1 - 1;
144: marker[head] = 2;
145: /* RESET MARKER VALUES. */
146: L1200:
147: root = nbrhd[inhd];
148: marker[root] = 0;
149: if (rchsze <= 0) goto L1400;
150: i__2 = rchsze;
151: for (irch = 1; irch <= i__2; ++irch) {
152: node = rchset[irch];
153: marker[node] = 0;
154: }
155: L1400:
156: ;
157: }
158: return(0);
159: }