Actual source code: genrcm.c

  2: /* genrcm.f -- translated by f2c (version 19931217).*/

  4: #include <petscsys.h>

  6: /*****************************************************************/
  7: /*****************************************************************/
  8: /*********   GENRCM ..... GENERAL REVERSE CUTHILL MCKEE   ********/
  9: /*****************************************************************/

 11: /*    PURPOSE - GENRCM FINDS THE REVERSE CUTHILL-MCKEE*/
 12: /*       ORDERING FOR A GENERAL GRAPH. FOR EACH CONNECTED*/
 13: /*       COMPONENT IN THE GRAPH, GENRCM OBTAINS THE ORDERING*/
 14: /*       BY CALLING THE SUBROUTINE RCM.*/

 16: /*    INPUT PARAMETERS -*/
 17: /*       NEQNS - NUMBER OF EQUATIONS*/
 18: /*       (XADJ, ADJNCY) - ARRAY PAIR CONTAINING THE ADJACENCY*/
 19: /*              STRUCTURE OF THE GRAPH OF THE MATRIX.*/

 21: /*    OUTPUT PARAMETER -*/
 22: /*       PERM - VECTOR THAT CONTAINS THE RCM ORDERING.*/

 24: /*    WORKING PARAMETERS -*/
 25: /*       MASK - IS USED TO MARK VARIABLES THAT HAVE BEEN*/
 26: /*              NUMBERED DURING THE ORDERING PROCESS. IT IS*/
 27: /*              INITIALIZED TO 1, AND SET TO ZERO AS EACH NODE*/
 28: /*              IS NUMBERED.*/
 29: /*       XLS - THE INDEX VECTOR FOR A LEVEL STRUCTURE.  THE*/
 30: /*              LEVEL STRUCTURE IS STORED IN THE CURRENTLY*/
 31: /*              UNUSED SPACES IN THE PERMUTATION VECTOR PERM.*/

 33: /*    PROGRAM SUBROUTINES -*/
 34: /*       FN../../.., RCM.*/
 35: /*****************************************************************/
 38: PetscErrorCode SPARSEPACKgenrcm(PetscInt *neqns,PetscInt *xadj,PetscInt *adjncy,PetscInt *perm,PetscInt *mask,PetscInt *xls)
 39: {
 40:     /* System generated locals */
 41:     PetscInt i__1;

 43:     /* Local variables */
 44:     PetscInt nlvl,root,i,ccsize;
 46:                SPARSEPACKrcm(PetscInt*,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *);
 47:     PetscInt num;

 50:     /* Parameter adjustments */
 51:     --xls;
 52:     --mask;
 53:     --perm;
 54:     --adjncy;
 55:     --xadj;

 57:     i__1 = *neqns;
 58:     for (i = 1; i <= i__1; ++i) {
 59:         mask[i] = 1;
 60:     }
 61:     num = 1;
 62:     i__1 = *neqns;
 63:     for (i = 1; i <= i__1; ++i) {
 64: /*          FOR EACH MASKED CONNECTED COMPONENT ...*/
 65:         if (!mask[i]) {
 66:             goto L200;
 67:         }
 68:         root = i;
 69: /*             FIRST FIND A PSEUDO-PERIPHERAL NODE ../../...*/
 70: /*             NOTE THAT THE LEVEL STRUCTURE FOUND BY*/
 71: /*             FN../../.. IS STORED STARTING AT PERM(NUM).*/
 72: /*             THEN RCM IS CALLED TO ORDER THE COMPONENT*/
 73: /*             USING ../../.. AS THE STARTING NODE.*/
 74:         SPARSEPACKfnroot(&root,&xadj[1],&adjncy[1],&mask[1],&nlvl,&xls[1],&perm[num]);
 75:         SPARSEPACKrcm(&root,&xadj[1],&adjncy[1],&mask[1],&perm[num],&ccsize,&xls[1]);
 76:         num += ccsize;
 77:         if (num > *neqns) {
 78:             return(0);
 79:         }
 80: L200:
 81:         ;
 82:     }
 83:     return(0);
 84: }