Actual source code: rcm.c
1: /* rcm.f -- translated by f2c (version 19931217).*/
3: #include <petscsys.h>
4: #include <petsc/private/matorderimpl.h>
6: /********* RCM ..... REVERSE CUTHILL-MCKEE ORDERING *******/
7: /* PURPOSE - RCM NUMBERS A CONNECTED COMPONENT SPECIFIED BY */
8: /* MASK AND ROOT, USING THE RCM ALGORITHM. */
9: /* THE NUMBERING IS TO BE STARTED AT THE NODE ROOT. */
10: /* */
11: /* INPUT PARAMETERS - */
12: /* ROOT - IS THE NODE THAT DEFINES THE CONNECTED */
13: /* COMPONENT AND IT IS USED AS THE STARTING */
14: /* NODE FOR THE RCM ORDERING. */
15: /* (XADJ, ADJNCY) - ADJACENCY STRUCTURE PAIR FOR */
16: /* THE GRAPH. */
17: /* */
18: /* UPDATED PARAMETERS - */
19: /* MASK - ONLY THOSE NODES WITH NONZERO INPUT MASK */
20: /* VALUES ARE CONSIDERED BY THE ROUTINE. THE */
21: /* NODES NUMBERED BY RCM WILL HAVE THEIR */
22: /* MASK VALUES SET TO ZERO. */
23: /* */
24: /* OUTPUT PARAMETERS - */
25: /* PERM - WILL CONTAIN THE RCM ORDERING. */
26: /* CCSIZE - IS THE SIZE OF THE CONNECTED COMPONENT */
27: /* THAT HAS BEEN NUMBERED BY RCM. */
28: /* */
29: /* WORKING PARAMETER - */
30: /* DEG - IS A TEMPORARY VECTOR USED TO HOLD THE DEGREE */
31: /* OF THE NODES IN THE SECTION GRAPH SPECIFIED */
32: /* BY MASK AND ROOT. */
33: /* */
34: /* PROGRAM SUBROUTINES - */
35: /* DEGREE. */
36: /* */
37: PetscErrorCode SPARSEPACKrcm(const PetscInt *root, const PetscInt *xadj, const PetscInt *adjncy, PetscInt *mask, PetscInt *perm, PetscInt *ccsize, PetscInt *deg)
38: {
39: /* System generated locals */
40: PetscInt i__1, i__2;
42: /* Local variables */
43: PetscInt node, fnbr, lnbr, i, j, k, l, lperm, jstop, jstrt;
44: PetscInt lbegin, lvlend, nbr;
46: /* FIND THE DEGREES OF THE NODES IN THE */
47: /* COMPONENT SPECIFIED BY MASK AND ROOT. */
49: PetscFunctionBegin;
50: /* Parameter adjustments */
51: --deg;
52: --perm;
53: --mask;
54: --adjncy;
55: --xadj;
57: PetscCall(SPARSEPACKdegree(root, &xadj[1], &adjncy[1], &mask[1], °[1], ccsize, &perm[1]));
58: mask[*root] = 0;
59: if (*ccsize <= 1) PetscFunctionReturn(PETSC_SUCCESS);
60: lvlend = 0;
61: lnbr = 1;
62: /* LBEGIN AND LVLEND POINT TO THE BEGINNING AND */
63: /* THE END OF THE CURRENT LEVEL RESPECTIVELY. */
64: L100:
65: lbegin = lvlend + 1;
66: lvlend = lnbr;
67: i__1 = lvlend;
68: for (i = lbegin; i <= i__1; ++i) {
69: /* FOR EACH NODE IN CURRENT LEVEL ... */
70: node = perm[i];
71: jstrt = xadj[node];
72: jstop = xadj[node + 1] - 1;
74: /* FIND THE UNNUMBERED NEIGHBORS OF NODE. */
75: /* FNBR AND LNBR POINT TO THE FIRST AND LAST */
76: /* UNNUMBERED NEIGHBORS RESPECTIVELY OF THE CURRENT */
77: /* NODE IN PERM. */
78: fnbr = lnbr + 1;
79: i__2 = jstop;
80: for (j = jstrt; j <= i__2; ++j) {
81: nbr = adjncy[j];
82: if (!mask[nbr]) goto L200;
83: ++lnbr;
84: mask[nbr] = 0;
85: perm[lnbr] = nbr;
86: L200:;
87: }
88: if (fnbr >= lnbr) goto L600;
90: /* SORT THE NEIGHBORS OF NODE IN INCREASING */
91: /* ORDER BY DEGREE. LINEAR INSERTION IS USED.*/
92: k = fnbr;
93: L300:
94: l = k;
95: ++k;
96: nbr = perm[k];
97: L400:
98: if (l < fnbr) goto L500;
99: lperm = perm[l];
100: if (deg[lperm] <= deg[nbr]) goto L500;
101: perm[l + 1] = lperm;
102: --l;
103: goto L400;
104: L500:
105: perm[l + 1] = nbr;
106: if (k < lnbr) goto L300;
107: L600:;
108: }
109: if (lnbr > lvlend) goto L100;
111: /* WE NOW HAVE THE CUTHILL MCKEE ORDERING.*/
112: /* REVERSE IT BELOW ...*/
113: k = *ccsize / 2;
114: l = *ccsize;
115: i__1 = k;
116: for (i = 1; i <= i__1; ++i) {
117: lperm = perm[l];
118: perm[l] = perm[i];
119: perm[i] = lperm;
120: --l;
121: }
122: PetscFunctionReturn(PETSC_SUCCESS);
123: }