xref: /petsc/src/mat/graphops/order/qmdrch.c (revision 53673ba54f5aaba04b9d49ab22cf56c7a7461fe9)
1*8be712e4SBarry Smith /* qmdrch.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 /**********     QMDRCH ..... QUOT MIN DEG REACH SET    ***********/
8*8be712e4SBarry Smith /*****************************************************************/
9*8be712e4SBarry Smith 
10*8be712e4SBarry Smith /*    PURPOSE - THIS SUBROUTINE DETERMINES THE REACHABLE SET OF*/
11*8be712e4SBarry Smith /*       A NODE THROUGH A GIVEN SUBSET.  THE ADJACENCY STRUCTURE*/
12*8be712e4SBarry Smith /*       IS ASSUMED TO BE STORED IN A QUOTIENT GRAPH FORMAT.*/
13*8be712e4SBarry Smith 
14*8be712e4SBarry Smith /*    INPUT PARAMETERS -*/
15*8be712e4SBarry Smith /*       ROOT - THE GIVEN NODE NOT IN THE SUBSET.*/
16*8be712e4SBarry Smith /*       (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE PAIR.*/
17*8be712e4SBarry Smith /*       DEG - THE DEGREE VECTOR.  DEG(I) LT 0 MEANS THE NODE*/
18*8be712e4SBarry Smith /*              BELONGS TO THE GIVEN SUBSET.*/
19*8be712e4SBarry Smith 
20*8be712e4SBarry Smith /*    OUTPUT PARAMETERS -*/
21*8be712e4SBarry Smith /*       (RCHSZE, RCHSET) - THE REACHABLE SET.*/
22*8be712e4SBarry Smith /*       (NHDSZE, NBRHD) - THE NEIGHBORHOOD SET.*/
23*8be712e4SBarry Smith 
24*8be712e4SBarry Smith /*    UPDATED PARAMETERS -*/
25*8be712e4SBarry Smith /*       MARKER - THE MARKER VECTOR FOR REACH AND NBRHD SETS.*/
26*8be712e4SBarry Smith /*              GT 0 MEANS THE NODE IS IN REACH SET.*/
27*8be712e4SBarry Smith /*              LT 0 MEANS THE NODE HAS BEEN MERGED WITH*/
28*8be712e4SBarry Smith /*              OTHERS IN THE QUOTIENT OR IT IS IN NBRHD SET.*/
29*8be712e4SBarry Smith /*****************************************************************/
SPARSEPACKqmdrch(const PetscInt * root,const PetscInt * xadj,const PetscInt * adjncy,PetscInt * deg,PetscInt * marker,PetscInt * rchsze,PetscInt * rchset,PetscInt * nhdsze,PetscInt * nbrhd)30*8be712e4SBarry Smith PetscErrorCode SPARSEPACKqmdrch(const PetscInt *root, const PetscInt *xadj, const PetscInt *adjncy, PetscInt *deg, PetscInt *marker, PetscInt *rchsze, PetscInt *rchset, PetscInt *nhdsze, PetscInt *nbrhd)
31*8be712e4SBarry Smith {
32*8be712e4SBarry Smith   /* System generated locals */
33*8be712e4SBarry Smith   PetscInt i__1, i__2;
34*8be712e4SBarry Smith 
35*8be712e4SBarry Smith   /* Local variables */
36*8be712e4SBarry Smith   PetscInt node, i, j, nabor, istop, jstop, istrt, jstrt;
37*8be712e4SBarry Smith 
38*8be712e4SBarry Smith   /*       LOOP THROUGH THE NEIGHBORS OF ROOT IN THE*/
39*8be712e4SBarry Smith   /*       QUOTIENT GRAPH.*/
40*8be712e4SBarry Smith 
41*8be712e4SBarry Smith   PetscFunctionBegin;
42*8be712e4SBarry Smith   /* Parameter adjustments */
43*8be712e4SBarry Smith   --nbrhd;
44*8be712e4SBarry Smith   --rchset;
45*8be712e4SBarry Smith   --marker;
46*8be712e4SBarry Smith   --deg;
47*8be712e4SBarry Smith   --adjncy;
48*8be712e4SBarry Smith   --xadj;
49*8be712e4SBarry Smith 
50*8be712e4SBarry Smith   *nhdsze = 0;
51*8be712e4SBarry Smith   *rchsze = 0;
52*8be712e4SBarry Smith   istrt   = xadj[*root];
53*8be712e4SBarry Smith   istop   = xadj[*root + 1] - 1;
54*8be712e4SBarry Smith   if (istop < istrt) PetscFunctionReturn(PETSC_SUCCESS);
55*8be712e4SBarry Smith   i__1 = istop;
56*8be712e4SBarry Smith   for (i = istrt; i <= i__1; ++i) {
57*8be712e4SBarry Smith     nabor = adjncy[i];
58*8be712e4SBarry Smith     if (!nabor) PetscFunctionReturn(PETSC_SUCCESS);
59*8be712e4SBarry Smith     if (marker[nabor] != 0) goto L600;
60*8be712e4SBarry Smith     if (deg[nabor] < 0) goto L200;
61*8be712e4SBarry Smith 
62*8be712e4SBarry Smith     /*                   INCLUDE NABOR INTO THE REACHABLE SET.*/
63*8be712e4SBarry Smith     ++(*rchsze);
64*8be712e4SBarry Smith     rchset[*rchsze] = nabor;
65*8be712e4SBarry Smith     marker[nabor]   = 1;
66*8be712e4SBarry Smith     goto L600;
67*8be712e4SBarry Smith   /*                NABOR HAS BEEN ELIMINATED. FIND NODES*/
68*8be712e4SBarry Smith   /*                REACHABLE FROM IT.*/
69*8be712e4SBarry Smith   L200:
70*8be712e4SBarry Smith     marker[nabor] = -1;
71*8be712e4SBarry Smith     ++(*nhdsze);
72*8be712e4SBarry Smith     nbrhd[*nhdsze] = nabor;
73*8be712e4SBarry Smith   L300:
74*8be712e4SBarry Smith     jstrt = xadj[nabor];
75*8be712e4SBarry Smith     jstop = xadj[nabor + 1] - 1;
76*8be712e4SBarry Smith     i__2  = jstop;
77*8be712e4SBarry Smith     for (j = jstrt; j <= i__2; ++j) {
78*8be712e4SBarry Smith       node  = adjncy[j];
79*8be712e4SBarry Smith       nabor = -node;
80*8be712e4SBarry Smith       if (node < 0) goto L300;
81*8be712e4SBarry Smith       else if (!node) goto L600;
82*8be712e4SBarry Smith       else goto L400;
83*8be712e4SBarry Smith     L400:
84*8be712e4SBarry Smith       if (marker[node] != 0) goto L500;
85*8be712e4SBarry Smith       ++(*rchsze);
86*8be712e4SBarry Smith       rchset[*rchsze] = node;
87*8be712e4SBarry Smith       marker[node]    = 1;
88*8be712e4SBarry Smith     L500:;
89*8be712e4SBarry Smith     }
90*8be712e4SBarry Smith   L600:;
91*8be712e4SBarry Smith   }
92*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
93*8be712e4SBarry Smith }
94