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