xref: /petsc/src/mat/graphops/order/rootls.c (revision 53673ba54f5aaba04b9d49ab22cf56c7a7461fe9)
1*8be712e4SBarry Smith /* rootls.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 /*********     ROOTLS ..... ROOTED LEVEL STRUCTURE      **********/
8*8be712e4SBarry Smith /*****************************************************************/
9*8be712e4SBarry Smith /*    PURPOSE - ROOTLS GENERATES THE LEVEL STRUCTURE ROOTED */
10*8be712e4SBarry Smith /*       AT THE INPUT NODE CALLED ROOT. ONLY THOSE NODES FOR*/
11*8be712e4SBarry Smith /*       WHICH MASK IS NONZERO WILL BE CONSIDERED.*/
12*8be712e4SBarry Smith /*                                                */
13*8be712e4SBarry Smith /*    INPUT PARAMETERS -                          */
14*8be712e4SBarry Smith /*       ROOT - THE NODE AT WHICH THE LEVEL STRUCTURE IS TO*/
15*8be712e4SBarry Smith /*              BE ROOTED.*/
16*8be712e4SBarry Smith /*       (XADJ, ADJNCY) - ADJACENCY STRUCTURE PAIR FOR THE*/
17*8be712e4SBarry Smith /*              GIVEN GRAPH.*/
18*8be712e4SBarry Smith /*       MASK - IS USED TO SPECIFY A SECTION SUBGRAPH. NODES*/
19*8be712e4SBarry Smith /*              WITH MASK(I)=0 ARE IGNORED.*/
20*8be712e4SBarry Smith /*    OUTPUT PARAMETERS -*/
21*8be712e4SBarry Smith /*       NLVL - IS THE NUMBER OF LEVELS IN THE LEVEL STRUCTURE.*/
22*8be712e4SBarry Smith /*       (XLS, LS) - ARRAY PAIR FOR THE ROOTED LEVEL STRUCTURE.*/
23*8be712e4SBarry Smith /*****************************************************************/
SPARSEPACKrootls(const PetscInt * root,const PetscInt * xadj,const PetscInt * adjncy,PetscInt * mask,PetscInt * nlvl,PetscInt * xls,PetscInt * ls)24*8be712e4SBarry Smith PetscErrorCode SPARSEPACKrootls(const PetscInt *root, const PetscInt *xadj, const PetscInt *adjncy, PetscInt *mask, PetscInt *nlvl, PetscInt *xls, PetscInt *ls)
25*8be712e4SBarry Smith {
26*8be712e4SBarry Smith   /* System generated locals */
27*8be712e4SBarry Smith   PetscInt i__1, i__2;
28*8be712e4SBarry Smith 
29*8be712e4SBarry Smith   /* Local variables */
30*8be712e4SBarry Smith   PetscInt node, i, j, jstop, jstrt, lbegin, ccsize, lvlend, lvsize, nbr;
31*8be712e4SBarry Smith 
32*8be712e4SBarry Smith   /*       INITIALIZATION ...*/
33*8be712e4SBarry Smith 
34*8be712e4SBarry Smith   PetscFunctionBegin;
35*8be712e4SBarry Smith   /* Parameter adjustments */
36*8be712e4SBarry Smith   --ls;
37*8be712e4SBarry Smith   --xls;
38*8be712e4SBarry Smith   --mask;
39*8be712e4SBarry Smith   --adjncy;
40*8be712e4SBarry Smith   --xadj;
41*8be712e4SBarry Smith 
42*8be712e4SBarry Smith   mask[*root] = 0;
43*8be712e4SBarry Smith   ls[1]       = *root;
44*8be712e4SBarry Smith   *nlvl       = 0;
45*8be712e4SBarry Smith   lvlend      = 0;
46*8be712e4SBarry Smith   ccsize      = 1;
47*8be712e4SBarry Smith /*       LBEGIN IS THE POINTER TO THE BEGINNING OF THE CURRENT*/
48*8be712e4SBarry Smith /*       LEVEL, AND LVLEND POINTS TO THE END OF THIS LEVEL.*/
49*8be712e4SBarry Smith L200:
50*8be712e4SBarry Smith   lbegin = lvlend + 1;
51*8be712e4SBarry Smith   lvlend = ccsize;
52*8be712e4SBarry Smith   ++(*nlvl);
53*8be712e4SBarry Smith   xls[*nlvl] = lbegin;
54*8be712e4SBarry Smith   /*       GENERATE THE NEXT LEVEL BY FINDING ALL THE MASKED */
55*8be712e4SBarry Smith   /*       NEIGHBORS OF NODES IN THE CURRENT LEVEL.*/
56*8be712e4SBarry Smith   i__1 = lvlend;
57*8be712e4SBarry Smith   for (i = lbegin; i <= i__1; ++i) {
58*8be712e4SBarry Smith     node  = ls[i];
59*8be712e4SBarry Smith     jstrt = xadj[node];
60*8be712e4SBarry Smith     jstop = xadj[node + 1] - 1;
61*8be712e4SBarry Smith     if (jstop < jstrt) goto L400;
62*8be712e4SBarry Smith     i__2 = jstop;
63*8be712e4SBarry Smith     for (j = jstrt; j <= i__2; ++j) {
64*8be712e4SBarry Smith       nbr = adjncy[j];
65*8be712e4SBarry Smith       if (!mask[nbr]) goto L300;
66*8be712e4SBarry Smith       ++ccsize;
67*8be712e4SBarry Smith       ls[ccsize] = nbr;
68*8be712e4SBarry Smith       mask[nbr]  = 0;
69*8be712e4SBarry Smith     L300:;
70*8be712e4SBarry Smith     }
71*8be712e4SBarry Smith   L400:;
72*8be712e4SBarry Smith   }
73*8be712e4SBarry Smith   /*       COMPUTE THE CURRENT LEVEL WIDTH.*/
74*8be712e4SBarry Smith   /*       IF IT IS NONZERO, GENERATE THE NEXT LEVEL.*/
75*8be712e4SBarry Smith   lvsize = ccsize - lvlend;
76*8be712e4SBarry Smith   if (lvsize > 0) goto L200;
77*8be712e4SBarry Smith   /*       RESET MASK TO ONE FOR THE NODES IN THE LEVEL STRUCTURE.*/
78*8be712e4SBarry Smith   xls[*nlvl + 1] = lvlend + 1;
79*8be712e4SBarry Smith   i__1           = ccsize;
80*8be712e4SBarry Smith   for (i = 1; i <= i__1; ++i) {
81*8be712e4SBarry Smith     node       = ls[i];
82*8be712e4SBarry Smith     mask[node] = 1;
83*8be712e4SBarry Smith   }
84*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
85*8be712e4SBarry Smith }
86