xref: /petsc/src/mat/graphops/order/qmdmrg.c (revision 53673ba54f5aaba04b9d49ab22cf56c7a7461fe9)
1*8be712e4SBarry Smith /* qmdmrg.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 /***********     QMDMRG ..... QUOT MIN DEG MERGE       ************/
8*8be712e4SBarry Smith /******************************************************************/
9*8be712e4SBarry Smith /*    PURPOSE - THIS ROUTINE MERGES INDISTINGUISHABLE NODES IN   */
10*8be712e4SBarry Smith /*              THE MINIMUM DEGREE ORDERING ALGORITHM.           */
11*8be712e4SBarry Smith /*              IT ALSO COMPUTES THE NEW DEGREES OF THESE        */
12*8be712e4SBarry Smith /*              NEW SUPERNODES.                                  */
13*8be712e4SBarry Smith /*                                                               */
14*8be712e4SBarry Smith /*    INPUT PARAMETERS -                                         */
15*8be712e4SBarry Smith /*       (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE.               */
16*8be712e4SBarry Smith /*       DEG0 - THE NUMBER OF NODES IN THE GIVEN SET.            */
17*8be712e4SBarry Smith /*       (NHDSZE, NBRHD) - THE SET OF ELIMINATED SUPERNODES      */
18*8be712e4SBarry Smith /*              ADJACENT TO SOME NODES IN THE SET.               */
19*8be712e4SBarry Smith /*                                                               */
20*8be712e4SBarry Smith /*    UPDATED PARAMETERS -                                       */
21*8be712e4SBarry Smith /*       DEG - THE DEGREE VECTOR.                                */
22*8be712e4SBarry Smith /*       QSIZE - SIZE OF INDISTINGUISHABLE NODES.                */
23*8be712e4SBarry Smith /*       QLINK - LINKED LIST FOR INDISTINGUISHABLE NODES.        */
24*8be712e4SBarry Smith /*       MARKER - THE GIVEN SET IS GIVEN BY THOSE NODES WITH     */
25*8be712e4SBarry Smith /*              MARKER VALUE SET TO 1.  THOSE NODES WITH DEGREE  */
26*8be712e4SBarry Smith /*              UPDATED WILL HAVE MARKER VALUE SET TO 2.         */
27*8be712e4SBarry Smith /*                                                               */
28*8be712e4SBarry Smith /*    WORKING PARAMETERS -                                       */
29*8be712e4SBarry Smith /*       RCHSET - THE REACHABLE SET.                             */
30*8be712e4SBarry Smith /*       OVRLP -  TEMP VECTOR TO STORE THE INTERSECTION OF TWO   */
31*8be712e4SBarry Smith /*              REACHABLE SETS.                                  */
32*8be712e4SBarry Smith /*                                                               */
33*8be712e4SBarry Smith /*****************************************************************/
SPARSEPACKqmdmrg(const PetscInt * xadj,const PetscInt * adjncy,PetscInt * deg,PetscInt * qsize,PetscInt * qlink,PetscInt * marker,PetscInt * deg0,PetscInt * nhdsze,PetscInt * nbrhd,PetscInt * rchset,PetscInt * ovrlp)34*8be712e4SBarry Smith PetscErrorCode SPARSEPACKqmdmrg(const PetscInt *xadj, const PetscInt *adjncy, PetscInt *deg, PetscInt *qsize, PetscInt *qlink, PetscInt *marker, PetscInt *deg0, PetscInt *nhdsze, PetscInt *nbrhd, PetscInt *rchset, PetscInt *ovrlp)
35*8be712e4SBarry Smith {
36*8be712e4SBarry Smith   /* System generated locals */
37*8be712e4SBarry Smith   PetscInt i__1, i__2, i__3;
38*8be712e4SBarry Smith 
39*8be712e4SBarry Smith   /* Local variables */
40*8be712e4SBarry Smith   PetscInt head, inhd, irch, node, mark, ilink, root, j, lnode, nabor, jstop, jstrt, rchsze, mrgsze, novrlp, iov, deg1;
41*8be712e4SBarry Smith 
42*8be712e4SBarry Smith   PetscFunctionBegin;
43*8be712e4SBarry Smith   /* Parameter adjustments */
44*8be712e4SBarry Smith   --ovrlp;
45*8be712e4SBarry Smith   --rchset;
46*8be712e4SBarry Smith   --nbrhd;
47*8be712e4SBarry Smith   --marker;
48*8be712e4SBarry Smith   --qlink;
49*8be712e4SBarry Smith   --qsize;
50*8be712e4SBarry Smith   --deg;
51*8be712e4SBarry Smith   --adjncy;
52*8be712e4SBarry Smith   --xadj;
53*8be712e4SBarry Smith 
54*8be712e4SBarry Smith   if (*nhdsze <= 0) PetscFunctionReturn(PETSC_SUCCESS);
55*8be712e4SBarry Smith   i__1 = *nhdsze;
56*8be712e4SBarry Smith   for (inhd = 1; inhd <= i__1; ++inhd) {
57*8be712e4SBarry Smith     root         = nbrhd[inhd];
58*8be712e4SBarry Smith     marker[root] = 0;
59*8be712e4SBarry Smith   }
60*8be712e4SBarry Smith   /*       LOOP THROUGH EACH ELIMINATED SUPERNODE IN THE SET     */
61*8be712e4SBarry Smith   /*       (NHDSZE, NBRHD).                                      */
62*8be712e4SBarry Smith   i__1 = *nhdsze;
63*8be712e4SBarry Smith   for (inhd = 1; inhd <= i__1; ++inhd) {
64*8be712e4SBarry Smith     root         = nbrhd[inhd];
65*8be712e4SBarry Smith     marker[root] = -1;
66*8be712e4SBarry Smith     rchsze       = 0;
67*8be712e4SBarry Smith     novrlp       = 0;
68*8be712e4SBarry Smith     deg1         = 0;
69*8be712e4SBarry Smith   L200:
70*8be712e4SBarry Smith     jstrt = xadj[root];
71*8be712e4SBarry Smith     jstop = xadj[root + 1] - 1;
72*8be712e4SBarry Smith     /*          DETERMINE THE REACHABLE SET AND ITS PETSCINTERSECT-     */
73*8be712e4SBarry Smith     /*          ION WITH THE INPUT REACHABLE SET.                  */
74*8be712e4SBarry Smith     i__2 = jstop;
75*8be712e4SBarry Smith     for (j = jstrt; j <= i__2; ++j) {
76*8be712e4SBarry Smith       nabor = adjncy[j];
77*8be712e4SBarry Smith       root  = -nabor;
78*8be712e4SBarry Smith       if (nabor < 0) goto L200;
79*8be712e4SBarry Smith       else if (!nabor) goto L700;
80*8be712e4SBarry Smith       else goto L300;
81*8be712e4SBarry Smith     L300:
82*8be712e4SBarry Smith       mark = marker[nabor];
83*8be712e4SBarry Smith       if (mark < 0) goto L600;
84*8be712e4SBarry Smith       else if (!mark) goto L400;
85*8be712e4SBarry Smith       else goto L500;
86*8be712e4SBarry Smith     L400:
87*8be712e4SBarry Smith       ++rchsze;
88*8be712e4SBarry Smith       rchset[rchsze] = nabor;
89*8be712e4SBarry Smith       deg1 += qsize[nabor];
90*8be712e4SBarry Smith       marker[nabor] = 1;
91*8be712e4SBarry Smith       goto L600;
92*8be712e4SBarry Smith     L500:
93*8be712e4SBarry Smith       if (mark > 1) goto L600;
94*8be712e4SBarry Smith       ++novrlp;
95*8be712e4SBarry Smith       ovrlp[novrlp] = nabor;
96*8be712e4SBarry Smith       marker[nabor] = 2;
97*8be712e4SBarry Smith     L600:;
98*8be712e4SBarry Smith     }
99*8be712e4SBarry Smith   /*          FROM THE OVERLAPPED SET, DETERMINE THE NODES        */
100*8be712e4SBarry Smith   /*          THAT CAN BE MERGED TOGETHER.                        */
101*8be712e4SBarry Smith   L700:
102*8be712e4SBarry Smith     head   = 0;
103*8be712e4SBarry Smith     mrgsze = 0;
104*8be712e4SBarry Smith     i__2   = novrlp;
105*8be712e4SBarry Smith     for (iov = 1; iov <= i__2; ++iov) {
106*8be712e4SBarry Smith       node  = ovrlp[iov];
107*8be712e4SBarry Smith       jstrt = xadj[node];
108*8be712e4SBarry Smith       jstop = xadj[node + 1] - 1;
109*8be712e4SBarry Smith       i__3  = jstop;
110*8be712e4SBarry Smith       for (j = jstrt; j <= i__3; ++j) {
111*8be712e4SBarry Smith         nabor = adjncy[j];
112*8be712e4SBarry Smith         if (marker[nabor] != 0) goto L800;
113*8be712e4SBarry Smith         marker[node] = 1;
114*8be712e4SBarry Smith         goto L1100;
115*8be712e4SBarry Smith       L800:;
116*8be712e4SBarry Smith       }
117*8be712e4SBarry Smith       /*             NODE BELONGS TO THE NEW MERGED SUPERNODE.      */
118*8be712e4SBarry Smith       /*             UPDATE THE VECTORS QLINK AND QSIZE.            */
119*8be712e4SBarry Smith       mrgsze += qsize[node];
120*8be712e4SBarry Smith       marker[node] = -1;
121*8be712e4SBarry Smith       lnode        = node;
122*8be712e4SBarry Smith     L900:
123*8be712e4SBarry Smith       ilink = qlink[lnode];
124*8be712e4SBarry Smith       if (ilink <= 0) goto L1000;
125*8be712e4SBarry Smith       lnode = ilink;
126*8be712e4SBarry Smith       goto L900;
127*8be712e4SBarry Smith     L1000:
128*8be712e4SBarry Smith       qlink[lnode] = head;
129*8be712e4SBarry Smith       head         = node;
130*8be712e4SBarry Smith     L1100:;
131*8be712e4SBarry Smith     }
132*8be712e4SBarry Smith     if (head <= 0) goto L1200;
133*8be712e4SBarry Smith     qsize[head]  = mrgsze;
134*8be712e4SBarry Smith     deg[head]    = *deg0 + deg1 - 1;
135*8be712e4SBarry Smith     marker[head] = 2;
136*8be712e4SBarry Smith   /*          RESET MARKER VALUES.          */
137*8be712e4SBarry Smith   L1200:
138*8be712e4SBarry Smith     root         = nbrhd[inhd];
139*8be712e4SBarry Smith     marker[root] = 0;
140*8be712e4SBarry Smith     if (rchsze <= 0) goto L1400;
141*8be712e4SBarry Smith     i__2 = rchsze;
142*8be712e4SBarry Smith     for (irch = 1; irch <= i__2; ++irch) {
143*8be712e4SBarry Smith       node         = rchset[irch];
144*8be712e4SBarry Smith       marker[node] = 0;
145*8be712e4SBarry Smith     }
146*8be712e4SBarry Smith   L1400:;
147*8be712e4SBarry Smith   }
148*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
149*8be712e4SBarry Smith }
150