Actual source code: rcm.c
petsc-3.10.5 2019-03-28
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: /* ------------------------------------- */
56: /* Parameter adjustments */
57: --deg;
58: --perm;
59: --mask;
60: --adjncy;
61: --xadj;
64: SPARSEPACKdegree(root, &xadj[1], &adjncy[1], &mask[1], °[1], ccsize, &perm[1]);
65: mask[*root] = 0;
66: if (*ccsize <= 1) return(0);
67: lvlend = 0;
68: lnbr = 1;
69: /* LBEGIN AND LVLEND POINT TO THE BEGINNING AND */
70: /* THE END OF THE CURRENT LEVEL RESPECTIVELY. */
71: L100:
72: lbegin = lvlend + 1;
73: lvlend = lnbr;
74: i__1 = lvlend;
75: for (i = lbegin; i <= i__1; ++i) {
76: /* FOR EACH NODE IN CURRENT LEVEL ... */
77: node = perm[i];
78: jstrt = xadj[node];
79: jstop = xadj[node + 1] - 1;
81: /* FIND THE UNNUMBERED NEIGHBORS OF NODE. */
82: /* FNBR AND LNBR POINT TO THE FIRST AND LAST */
83: /* UNNUMBERED NEIGHBORS RESPECTIVELY OF THE CURRENT */
84: /* NODE IN PERM. */
85: fnbr = lnbr + 1;
86: i__2 = jstop;
87: for (j = jstrt; j <= i__2; ++j) {
88: nbr = adjncy[j];
89: if (!mask[nbr]) goto L200;
90: ++lnbr;
91: mask[nbr] = 0;
92: perm[lnbr] = nbr;
93: L200:
94: ;
95: }
96: if (fnbr >= lnbr) goto L600;
98: /* SORT THE NEIGHBORS OF NODE IN INCREASING */
99: /* ORDER BY DEGREE. LINEAR INSERTION IS USED.*/
100: k = fnbr;
101: L300:
102: l = k;
103: ++k;
104: nbr = perm[k];
105: L400:
106: if (l < fnbr) goto L500;
107: lperm = perm[l];
108: if (deg[lperm] <= deg[nbr]) goto L500;
109: perm[l + 1] = lperm;
110: --l;
111: goto L400;
112: L500:
113: perm[l + 1] = nbr;
114: if (k < lnbr) goto L300;
115: L600:
116: ;
117: }
118: if (lnbr > lvlend) goto L100;
120: /* WE NOW HAVE THE CUTHILL MCKEE ORDERING.*/
121: /* REVERSE IT BELOW ...*/
122: k = *ccsize / 2;
123: l = *ccsize;
124: i__1 = k;
125: for (i = 1; i <= i__1; ++i) {
126: lperm = perm[l];
127: perm[l] = perm[i];
128: perm[i] = lperm;
129: --l;
130: }
131: return(0);
132: }