Actual source code: genrcm.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: /* genrcm.f -- translated by f2c (version 19931217).*/

  4: #include <petscsys.h>
  5: #include <petsc/private/matorderimpl.h>

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

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

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

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

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

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

 44:   /* Local variables */
 45:   PetscInt nlvl,root,i,ccsize;
 46:   PetscInt num;

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

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