Actual source code: rcm.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  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: /****************************************************************/
 43: PetscErrorCode SPARSEPACKrcm(const PetscInt *root,const PetscInt *xadj,const PetscInt *adjncy,PetscInt *mask, PetscInt *perm, PetscInt *ccsize, PetscInt *deg)
 44: {
 45:   /* System generated locals */
 46:   PetscInt i__1, i__2;

 48:   /* Local variables */
 49:   PetscInt node, fnbr, lnbr, i, j, k, l, lperm, jstop, jstrt;
 50:   PetscInt lbegin, lvlend, nbr;

 52: /*       FIND THE DEGREES OF THE NODES IN THE                  */
 53: /*       COMPONENT SPECIFIED BY MASK AND ../../...                 */
 54: /*       -------------------------------------                 */


 58:   /* Parameter adjustments */
 59:   --deg;
 60:   --perm;
 61:   --mask;
 62:   --adjncy;
 63:   --xadj;


 66:   SPARSEPACKdegree(root, &xadj[1], &adjncy[1], &mask[1], &deg[1], ccsize, &perm[1]);
 67:   mask[*root] = 0;
 68:   if (*ccsize <= 1) return(0);
 69:   lvlend = 0;
 70:   lnbr   = 1;
 71: /*       LBEGIN AND LVLEND POINT TO THE BEGINNING AND */
 72: /*       THE END OF THE CURRENT LEVEL RESPECTIVELY.  */
 73: L100:
 74:   lbegin = lvlend + 1;
 75:   lvlend = lnbr;
 76:   i__1   = lvlend;
 77:   for (i = lbegin; i <= i__1; ++i) {
 78: /*          FOR EACH NODE IN CURRENT LEVEL ...     */
 79:     node  = perm[i];
 80:     jstrt = xadj[node];
 81:     jstop = xadj[node + 1] - 1;

 83: /*          FIND THE UNNUMBERED NEIGHBORS OF NODE.   */
 84: /*          FNBR AND LNBR POINT TO THE FIRST AND LAST  */
 85: /*          UNNUMBERED NEIGHBORS RESPECTIVELY OF THE CURRENT  */
 86: /*          NODE IN PERM. */
 87:     fnbr = lnbr + 1;
 88:     i__2 = jstop;
 89:     for (j = jstrt; j <= i__2; ++j) {
 90:       nbr = adjncy[j];
 91:       if (!mask[nbr]) goto L200;
 92:       ++lnbr;
 93:       mask[nbr]  = 0;
 94:       perm[lnbr] = nbr;
 95: L200:
 96:       ;
 97:     }
 98:     if (fnbr >= lnbr) goto L600;

100: /*             SORT THE NEIGHBORS OF NODE IN INCREASING    */
101: /*             ORDER BY DEGREE. LINEAR INSERTION IS USED.*/
102:     k = fnbr;
103: L300:
104:     l = k;
105:     ++k;
106:     nbr = perm[k];
107: L400:
108:     if (l < fnbr) goto L500;
109:     lperm = perm[l];
110:     if (deg[lperm] <= deg[nbr]) goto L500;
111:     perm[l + 1] = lperm;
112:     --l;
113:     goto L400;
114: L500:
115:     perm[l + 1] = nbr;
116:     if (k < lnbr) goto L300;
117: L600:
118:     ;
119:   }
120:   if (lnbr > lvlend) goto L100;

122: /*       WE NOW HAVE THE CUTHILL MCKEE ORDERING.*/
123: /*       REVERSE IT BELOW ...*/
124:   k    = *ccsize / 2;
125:   l    = *ccsize;
126:   i__1 = k;
127:   for (i = 1; i <= i__1; ++i) {
128:     lperm   = perm[l];
129:     perm[l] = perm[i];
130:     perm[i] = lperm;
131:     --l;
132:   }
133:   return(0);
134: }