xref: /petsc/src/mat/graphops/order/genrcm.c (revision 8be712e46db5d855f641c6bd97b4543e0efe65bd)
1*8be712e4SBarry Smith /* genrcm.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 /*****************************************************************/
7*8be712e4SBarry Smith /*****************************************************************/
8*8be712e4SBarry Smith /*********   GENRCM ..... GENERAL REVERSE CUTHILL MCKEE   ********/
9*8be712e4SBarry Smith /*****************************************************************/
10*8be712e4SBarry Smith 
11*8be712e4SBarry Smith /*    PURPOSE - GENRCM FINDS THE REVERSE CUTHILL-MCKEE*/
12*8be712e4SBarry Smith /*       ORDERING FOR A GENERAL GRAPH. FOR EACH CONNECTED*/
13*8be712e4SBarry Smith /*       COMPONENT IN THE GRAPH, GENRCM OBTAINS THE ORDERING*/
14*8be712e4SBarry Smith /*       BY CALLING THE SUBROUTINE RCM.*/
15*8be712e4SBarry Smith 
16*8be712e4SBarry Smith /*    INPUT PARAMETERS -*/
17*8be712e4SBarry Smith /*       NEQNS - NUMBER OF EQUATIONS*/
18*8be712e4SBarry Smith /*       (XADJ, ADJNCY) - ARRAY PAIR CONTAINING THE ADJACENCY*/
19*8be712e4SBarry Smith /*              STRUCTURE OF THE GRAPH OF THE MATRIX.*/
20*8be712e4SBarry Smith 
21*8be712e4SBarry Smith /*    OUTPUT PARAMETER -*/
22*8be712e4SBarry Smith /*       PERM - VECTOR THAT CONTAINS THE RCM ORDERING.*/
23*8be712e4SBarry Smith 
24*8be712e4SBarry Smith /*    WORKING PARAMETERS -*/
25*8be712e4SBarry Smith /*       MASK - IS USED TO MARK VARIABLES THAT HAVE BEEN*/
26*8be712e4SBarry Smith /*              NUMBERED DURING THE ORDERING PROCESS. IT IS*/
27*8be712e4SBarry Smith /*              INITIALIZED TO 1, AND SET TO ZERO AS EACH NODE*/
28*8be712e4SBarry Smith /*              IS NUMBERED.*/
29*8be712e4SBarry Smith /*       XLS - THE INDEX VECTOR FOR A LEVEL STRUCTURE.  THE*/
30*8be712e4SBarry Smith /*              LEVEL STRUCTURE IS STORED IN THE CURRENTLY*/
31*8be712e4SBarry Smith /*              UNUSED SPACES IN THE PERMUTATION VECTOR PERM.*/
32*8be712e4SBarry Smith 
33*8be712e4SBarry Smith /*    PROGRAM SUBROUTINES -*/
34*8be712e4SBarry Smith /*       FNROOT, RCM.*/
35*8be712e4SBarry Smith /*****************************************************************/
36*8be712e4SBarry Smith PetscErrorCode SPARSEPACKgenrcm(const PetscInt *neqns, const PetscInt *xadj, const PetscInt *adjncy, PetscInt *perm, PetscInt *mask, PetscInt *xls)
37*8be712e4SBarry Smith {
38*8be712e4SBarry Smith   /* System generated locals */
39*8be712e4SBarry Smith   PetscInt i__1;
40*8be712e4SBarry Smith 
41*8be712e4SBarry Smith   /* Local variables */
42*8be712e4SBarry Smith   PetscInt nlvl, root, i, ccsize;
43*8be712e4SBarry Smith   PetscInt num;
44*8be712e4SBarry Smith 
45*8be712e4SBarry Smith   PetscFunctionBegin;
46*8be712e4SBarry Smith   /* Parameter adjustments */
47*8be712e4SBarry Smith   --xls;
48*8be712e4SBarry Smith   --mask;
49*8be712e4SBarry Smith   --perm;
50*8be712e4SBarry Smith   --adjncy;
51*8be712e4SBarry Smith   --xadj;
52*8be712e4SBarry Smith 
53*8be712e4SBarry Smith   i__1 = *neqns;
54*8be712e4SBarry Smith   for (i = 1; i <= i__1; ++i) mask[i] = 1;
55*8be712e4SBarry Smith   num  = 1;
56*8be712e4SBarry Smith   i__1 = *neqns;
57*8be712e4SBarry Smith   for (i = 1; i <= i__1; ++i) {
58*8be712e4SBarry Smith     /*          FOR EACH MASKED CONNECTED COMPONENT ...*/
59*8be712e4SBarry Smith     if (!mask[i]) goto L200;
60*8be712e4SBarry Smith     root = i;
61*8be712e4SBarry Smith     /*             FIRST FIND A PSEUDO-PERIPHERAL NODE ROOT.*/
62*8be712e4SBarry Smith     /*             NOTE THAT THE LEVEL STRUCTURE FOUND BY*/
63*8be712e4SBarry Smith     /*             FNROOT IS STORED STARTING AT PERM(NUM).*/
64*8be712e4SBarry Smith     /*             THEN RCM IS CALLED TO ORDER THE COMPONENT*/
65*8be712e4SBarry Smith     /*             USING ROOT AS THE STARTING NODE.*/
66*8be712e4SBarry Smith     PetscCall(SPARSEPACKfnroot(&root, &xadj[1], &adjncy[1], &mask[1], &nlvl, &xls[1], &perm[num]));
67*8be712e4SBarry Smith     PetscCall(SPARSEPACKrcm(&root, &xadj[1], &adjncy[1], &mask[1], &perm[num], &ccsize, &xls[1]));
68*8be712e4SBarry Smith     num += ccsize;
69*8be712e4SBarry Smith     if (num > *neqns) PetscFunctionReturn(PETSC_SUCCESS);
70*8be712e4SBarry Smith   L200:;
71*8be712e4SBarry Smith   }
72*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
73*8be712e4SBarry Smith }
74