Actual source code: genrcm.c


  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: /*****************************************************************/
 37: PetscErrorCode SPARSEPACKgenrcm(const PetscInt *neqns, const PetscInt *xadj, const PetscInt *adjncy, PetscInt *perm, PetscInt *mask, PetscInt *xls)
 38: {
 39:   /* System generated locals */
 40:   PetscInt i__1;

 42:   /* Local variables */
 43:   PetscInt nlvl, root, i, ccsize;
 44:   PetscInt num;

 46:   PetscFunctionBegin;
 47:   /* Parameter adjustments */
 48:   --xls;
 49:   --mask;
 50:   --perm;
 51:   --adjncy;
 52:   --xadj;

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