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