xref: /petsc/src/mat/graphops/order/qmdupd.c (revision 53673ba54f5aaba04b9d49ab22cf56c7a7461fe9)
1*8be712e4SBarry Smith /* qmdupd.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 /***********     QMDUPD ..... QUOT MIN DEG UPDATE      ************/
8*8be712e4SBarry Smith /******************************************************************/
9*8be712e4SBarry Smith /******************************************************************/
10*8be712e4SBarry Smith 
11*8be712e4SBarry Smith /*    PURPOSE - THIS ROUTINE PERFORMS DEGREE UPDATE FOR A SET*/
12*8be712e4SBarry Smith /*       OF NODES IN THE MINIMUM DEGREE ALGORITHM.*/
13*8be712e4SBarry Smith 
14*8be712e4SBarry Smith /*    INPUT PARAMETERS -*/
15*8be712e4SBarry Smith /*       (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE.*/
16*8be712e4SBarry Smith /*       (NLIST, LIST) - THE LIST OF NODES WHOSE DEGREE HAS TO*/
17*8be712e4SBarry Smith /*              BE UPDATED.*/
18*8be712e4SBarry Smith 
19*8be712e4SBarry Smith /*    UPDATED PARAMETERS -*/
20*8be712e4SBarry Smith /*       DEG - THE DEGREE VECTOR.*/
21*8be712e4SBarry Smith /*       QSIZE - SIZE OF INDISTINGUISHABLE SUPERNODES.*/
22*8be712e4SBarry Smith /*       QLINK - LINKED LIST FOR INDISTINGUISHABLE NODES.*/
23*8be712e4SBarry Smith /*       MARKER - USED TO MARK THOSE NODES IN REACH/NBRHD SETS.*/
24*8be712e4SBarry Smith 
25*8be712e4SBarry Smith /*    WORKING PARAMETERS -*/
26*8be712e4SBarry Smith /*       RCHSET - THE REACHABLE SET.*/
27*8be712e4SBarry Smith /*       NBRHD -  THE NEIGHBORHOOD SET.*/
28*8be712e4SBarry Smith 
29*8be712e4SBarry Smith /*    PROGRAM SUBROUTINES -*/
30*8be712e4SBarry Smith /*       QMDMRG.*/
31*8be712e4SBarry Smith /******************************************************************/
SPARSEPACKqmdupd(const PetscInt * xadj,const PetscInt * adjncy,const PetscInt * nlist,const PetscInt * list,PetscInt * deg,PetscInt * qsize,PetscInt * qlink,PetscInt * marker,PetscInt * rchset,PetscInt * nbrhd)32*8be712e4SBarry Smith PetscErrorCode SPARSEPACKqmdupd(const PetscInt *xadj, const PetscInt *adjncy, const PetscInt *nlist, const PetscInt *list, PetscInt *deg, PetscInt *qsize, PetscInt *qlink, PetscInt *marker, PetscInt *rchset, PetscInt *nbrhd)
33*8be712e4SBarry Smith {
34*8be712e4SBarry Smith   /* System generated locals */
35*8be712e4SBarry Smith   PetscInt i__1, i__2;
36*8be712e4SBarry Smith 
37*8be712e4SBarry Smith   /* Local variables */
38*8be712e4SBarry Smith   PetscInt inhd, irch, node, mark, j, inode, nabor, jstop, jstrt, il;
39*8be712e4SBarry Smith   PetscInt nhdsze, rchsze, deg0, deg1;
40*8be712e4SBarry Smith 
41*8be712e4SBarry Smith   /*       FIND ALL ELIMINATED SUPERNODES THAT ARE ADJACENT*/
42*8be712e4SBarry Smith   /*       TO SOME NODES IN THE GIVEN LIST. PUT THEM INTO.*/
43*8be712e4SBarry Smith   /*       (NHDSZE, NBRHD). DEG0 CONTAINS THE NUMBER OF*/
44*8be712e4SBarry Smith   /*       NODES IN THE LIST.*/
45*8be712e4SBarry Smith 
46*8be712e4SBarry Smith   PetscFunctionBegin;
47*8be712e4SBarry Smith   /* Parameter adjustments */
48*8be712e4SBarry Smith   --nbrhd;
49*8be712e4SBarry Smith   --rchset;
50*8be712e4SBarry Smith   --marker;
51*8be712e4SBarry Smith   --qlink;
52*8be712e4SBarry Smith   --qsize;
53*8be712e4SBarry Smith   --deg;
54*8be712e4SBarry Smith   --list;
55*8be712e4SBarry Smith   --adjncy;
56*8be712e4SBarry Smith   --xadj;
57*8be712e4SBarry Smith 
58*8be712e4SBarry Smith   if (*nlist <= 0) PetscFunctionReturn(PETSC_SUCCESS);
59*8be712e4SBarry Smith   deg0   = 0;
60*8be712e4SBarry Smith   nhdsze = 0;
61*8be712e4SBarry Smith   i__1   = *nlist;
62*8be712e4SBarry Smith   for (il = 1; il <= i__1; ++il) {
63*8be712e4SBarry Smith     node = list[il];
64*8be712e4SBarry Smith     deg0 += qsize[node];
65*8be712e4SBarry Smith     jstrt = xadj[node];
66*8be712e4SBarry Smith     jstop = xadj[node + 1] - 1;
67*8be712e4SBarry Smith     i__2  = jstop;
68*8be712e4SBarry Smith     for (j = jstrt; j <= i__2; ++j) {
69*8be712e4SBarry Smith       nabor = adjncy[j];
70*8be712e4SBarry Smith       if (marker[nabor] != 0 || deg[nabor] >= 0) goto L100;
71*8be712e4SBarry Smith       marker[nabor] = -1;
72*8be712e4SBarry Smith       ++nhdsze;
73*8be712e4SBarry Smith       nbrhd[nhdsze] = nabor;
74*8be712e4SBarry Smith     L100:;
75*8be712e4SBarry Smith     }
76*8be712e4SBarry Smith   }
77*8be712e4SBarry Smith   /*       MERGE INDISTINGUISHABLE NODES IN THE LIST BY*/
78*8be712e4SBarry Smith   /*       CALLING THE SUBROUTINE QMDMRG.*/
79*8be712e4SBarry Smith   if (nhdsze > 0) PetscCall(SPARSEPACKqmdmrg(&xadj[1], &adjncy[1], &deg[1], &qsize[1], &qlink[1], &marker[1], &deg0, &nhdsze, &nbrhd[1], &rchset[1], &nbrhd[nhdsze + 1]));
80*8be712e4SBarry Smith   /*       FIND THE NEW DEGREES OF THE NODES THAT HAVE NOT BEEN*/
81*8be712e4SBarry Smith   /*       MERGED.*/
82*8be712e4SBarry Smith   i__1 = *nlist;
83*8be712e4SBarry Smith   for (il = 1; il <= i__1; ++il) {
84*8be712e4SBarry Smith     node = list[il];
85*8be712e4SBarry Smith     mark = marker[node];
86*8be712e4SBarry Smith     if (mark > 1 || mark < 0) goto L600;
87*8be712e4SBarry Smith     marker[node] = 2;
88*8be712e4SBarry Smith     PetscCall(SPARSEPACKqmdrch(&node, &xadj[1], &adjncy[1], &deg[1], &marker[1], &rchsze, &rchset[1], &nhdsze, &nbrhd[1]));
89*8be712e4SBarry Smith     deg1 = deg0;
90*8be712e4SBarry Smith     if (rchsze <= 0) goto L400;
91*8be712e4SBarry Smith     i__2 = rchsze;
92*8be712e4SBarry Smith     for (irch = 1; irch <= i__2; ++irch) {
93*8be712e4SBarry Smith       inode = rchset[irch];
94*8be712e4SBarry Smith       deg1 += qsize[inode];
95*8be712e4SBarry Smith       marker[inode] = 0;
96*8be712e4SBarry Smith     }
97*8be712e4SBarry Smith   L400:
98*8be712e4SBarry Smith     deg[node] = deg1 - 1;
99*8be712e4SBarry Smith     if (nhdsze <= 0) goto L600;
100*8be712e4SBarry Smith     i__2 = nhdsze;
101*8be712e4SBarry Smith     for (inhd = 1; inhd <= i__2; ++inhd) {
102*8be712e4SBarry Smith       inode         = nbrhd[inhd];
103*8be712e4SBarry Smith       marker[inode] = 0;
104*8be712e4SBarry Smith     }
105*8be712e4SBarry Smith   L600:;
106*8be712e4SBarry Smith   }
107*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
108*8be712e4SBarry Smith }
109