1*8be712e4SBarry Smith /* genqmd.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 /*********** GENQMD ..... QUOT MIN DEGREE ORDERING **********/
8*8be712e4SBarry Smith /******************************************************************/
9*8be712e4SBarry Smith /* PURPOSE - THIS ROUTINE IMPLEMENTS THE MINIMUM DEGREE */
10*8be712e4SBarry Smith /* ALGORITHM. IT MAKES USE OF THE IMPLICIT REPRESENT- */
11*8be712e4SBarry Smith /* ATION OF THE ELIMINATION GRAPHS BY QUOTIENT GRAPHS, */
12*8be712e4SBarry Smith /* AND THE NOTION OF INDISTINGUISHABLE NODES. */
13*8be712e4SBarry Smith /* CAUTION - THE ADJACENCY VECTOR ADJNCY WILL BE */
14*8be712e4SBarry Smith /* DESTROYED. */
15*8be712e4SBarry Smith /* */
16*8be712e4SBarry Smith /* INPUT PARAMETERS - */
17*8be712e4SBarry Smith /* NEQNS - NUMBER OF EQUATIONS. */
18*8be712e4SBarry Smith /* (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE. */
19*8be712e4SBarry Smith /* */
20*8be712e4SBarry Smith /* OUTPUT PARAMETERS - */
21*8be712e4SBarry Smith /* PERM - THE MINIMUM DEGREE ORDERING. */
22*8be712e4SBarry Smith /* INVP - THE INVERSE OF PERM. */
23*8be712e4SBarry Smith /* */
24*8be712e4SBarry Smith /* WORKING PARAMETERS - */
25*8be712e4SBarry Smith /* DEG - THE DEGREE VECTOR. DEG(I) IS NEGATIVE MEANS */
26*8be712e4SBarry Smith /* NODE I HAS BEEN NUMBERED. */
27*8be712e4SBarry Smith /* MARKER - A MARKER VECTOR, WHERE MARKER(I) IS */
28*8be712e4SBarry Smith /* NEGATIVE MEANS NODE I HAS BEEN MERGED WITH */
29*8be712e4SBarry Smith /* ANOTHER NODE AND THUS CAN BE IGNORED. */
30*8be712e4SBarry Smith /* RCHSET - VECTOR USED FOR THE REACHABLE SET. */
31*8be712e4SBarry Smith /* NBRHD - VECTOR USED FOR THE NEIGHBORHOOD SET. */
32*8be712e4SBarry Smith /* QSIZE - VECTOR USED TO STORE THE SIZE OF */
33*8be712e4SBarry Smith /* INDISTINGUISHABLE SUPERNODES. */
34*8be712e4SBarry Smith /* QLINK - VECTOR TO STORE INDISTINGUISHABLE NODES, */
35*8be712e4SBarry Smith /* I, QLINK(I), QLINK(QLINK(I)) ... ARE THE */
36*8be712e4SBarry Smith /* MEMBERS OF THE SUPERNODE REPRESENTED BY I. */
37*8be712e4SBarry Smith /* */
38*8be712e4SBarry Smith /* PROGRAM SUBROUTINES - */
39*8be712e4SBarry Smith /* QMDRCH, QMDQT, QMDUPD. */
40*8be712e4SBarry Smith /* */
41*8be712e4SBarry Smith /******************************************************************/
42*8be712e4SBarry Smith /* */
43*8be712e4SBarry Smith /* */
SPARSEPACKgenqmd(const PetscInt * neqns,const PetscInt * xadj,const PetscInt * adjncy,PetscInt * perm,PetscInt * invp,PetscInt * deg,PetscInt * marker,PetscInt * rchset,PetscInt * nbrhd,PetscInt * qsize,PetscInt * qlink,PetscInt * nofsub)44*8be712e4SBarry Smith PetscErrorCode SPARSEPACKgenqmd(const PetscInt *neqns, const PetscInt *xadj, const PetscInt *adjncy, PetscInt *perm, PetscInt *invp, PetscInt *deg, PetscInt *marker, PetscInt *rchset, PetscInt *nbrhd, PetscInt *qsize, PetscInt *qlink, PetscInt *nofsub)
45*8be712e4SBarry Smith {
46*8be712e4SBarry Smith /* System generated locals */
47*8be712e4SBarry Smith PetscInt i__1;
48*8be712e4SBarry Smith
49*8be712e4SBarry Smith /* Local variables */
50*8be712e4SBarry Smith PetscInt ndeg, irch, node, nump1, j, inode;
51*8be712e4SBarry Smith PetscInt ip, np, mindeg, search;
52*8be712e4SBarry Smith PetscInt nhdsze, nxnode, rchsze, thresh, num;
53*8be712e4SBarry Smith
54*8be712e4SBarry Smith /* INITIALIZE DEGREE VECTOR AND OTHER WORKING VARIABLES. */
55*8be712e4SBarry Smith
56*8be712e4SBarry Smith PetscFunctionBegin;
57*8be712e4SBarry Smith /* Parameter adjustments */
58*8be712e4SBarry Smith --qlink;
59*8be712e4SBarry Smith --qsize;
60*8be712e4SBarry Smith --nbrhd;
61*8be712e4SBarry Smith --rchset;
62*8be712e4SBarry Smith --marker;
63*8be712e4SBarry Smith --deg;
64*8be712e4SBarry Smith --invp;
65*8be712e4SBarry Smith --perm;
66*8be712e4SBarry Smith --adjncy;
67*8be712e4SBarry Smith --xadj;
68*8be712e4SBarry Smith
69*8be712e4SBarry Smith mindeg = *neqns;
70*8be712e4SBarry Smith *nofsub = 0;
71*8be712e4SBarry Smith i__1 = *neqns;
72*8be712e4SBarry Smith for (node = 1; node <= i__1; ++node) {
73*8be712e4SBarry Smith perm[node] = node;
74*8be712e4SBarry Smith invp[node] = node;
75*8be712e4SBarry Smith marker[node] = 0;
76*8be712e4SBarry Smith qsize[node] = 1;
77*8be712e4SBarry Smith qlink[node] = 0;
78*8be712e4SBarry Smith ndeg = xadj[node + 1] - xadj[node];
79*8be712e4SBarry Smith deg[node] = ndeg;
80*8be712e4SBarry Smith if (ndeg < mindeg) mindeg = ndeg;
81*8be712e4SBarry Smith }
82*8be712e4SBarry Smith num = 0;
83*8be712e4SBarry Smith /* PERFORM THRESHOLD SEARCH TO GET A NODE OF MIN DEGREE. */
84*8be712e4SBarry Smith /* VARIABLE SEARCH POINTS TO WHERE SEARCH SHOULD START. */
85*8be712e4SBarry Smith L200:
86*8be712e4SBarry Smith search = 1;
87*8be712e4SBarry Smith thresh = mindeg;
88*8be712e4SBarry Smith mindeg = *neqns;
89*8be712e4SBarry Smith L300:
90*8be712e4SBarry Smith nump1 = num + 1;
91*8be712e4SBarry Smith if (nump1 > search) search = nump1;
92*8be712e4SBarry Smith i__1 = *neqns;
93*8be712e4SBarry Smith for (j = search; j <= i__1; ++j) {
94*8be712e4SBarry Smith node = perm[j];
95*8be712e4SBarry Smith if (marker[node] < 0) goto L400;
96*8be712e4SBarry Smith ndeg = deg[node];
97*8be712e4SBarry Smith if (ndeg <= thresh) goto L500;
98*8be712e4SBarry Smith if (ndeg < mindeg) mindeg = ndeg;
99*8be712e4SBarry Smith L400:;
100*8be712e4SBarry Smith }
101*8be712e4SBarry Smith goto L200;
102*8be712e4SBarry Smith /* NODE HAS MINIMUM DEGREE. FIND ITS REACHABLE SETS BY */
103*8be712e4SBarry Smith /* CALLING QMDRCH. */
104*8be712e4SBarry Smith L500:
105*8be712e4SBarry Smith search = j;
106*8be712e4SBarry Smith *nofsub += deg[node];
107*8be712e4SBarry Smith marker[node] = 1;
108*8be712e4SBarry Smith PetscCall(SPARSEPACKqmdrch(&node, &xadj[1], &adjncy[1], °[1], &marker[1], &rchsze, &rchset[1], &nhdsze, &nbrhd[1]));
109*8be712e4SBarry Smith /* ELIMINATE ALL NODES INDISTINGUISHABLE FROM NODE. */
110*8be712e4SBarry Smith /* THEY ARE GIVEN BY NODE, QLINK(NODE), .... */
111*8be712e4SBarry Smith nxnode = node;
112*8be712e4SBarry Smith L600:
113*8be712e4SBarry Smith ++num;
114*8be712e4SBarry Smith np = invp[nxnode];
115*8be712e4SBarry Smith ip = perm[num];
116*8be712e4SBarry Smith perm[np] = ip;
117*8be712e4SBarry Smith invp[ip] = np;
118*8be712e4SBarry Smith perm[num] = nxnode;
119*8be712e4SBarry Smith invp[nxnode] = num;
120*8be712e4SBarry Smith deg[nxnode] = -1;
121*8be712e4SBarry Smith nxnode = qlink[nxnode];
122*8be712e4SBarry Smith if (nxnode > 0) goto L600;
123*8be712e4SBarry Smith if (rchsze <= 0) goto L800;
124*8be712e4SBarry Smith
125*8be712e4SBarry Smith /* UPDATE THE DEGREES OF THE NODES IN THE REACHABLE */
126*8be712e4SBarry Smith /* SET AND IDENTIFY INDISTINGUISHABLE NODES. */
127*8be712e4SBarry Smith PetscCall(SPARSEPACKqmdupd(&xadj[1], &adjncy[1], &rchsze, &rchset[1], °[1], &qsize[1], &qlink[1], &marker[1], &rchset[rchsze + 1], &nbrhd[nhdsze + 1]));
128*8be712e4SBarry Smith
129*8be712e4SBarry Smith /* RESET MARKER VALUE OF NODES IN REACH SET. */
130*8be712e4SBarry Smith /* UPDATE THRESHOLD VALUE FOR CYCLIC SEARCH. */
131*8be712e4SBarry Smith /* ALSO CALL QMDQT TO FORM NEW QUOTIENT GRAPH. */
132*8be712e4SBarry Smith marker[node] = 0;
133*8be712e4SBarry Smith i__1 = rchsze;
134*8be712e4SBarry Smith for (irch = 1; irch <= i__1; ++irch) {
135*8be712e4SBarry Smith inode = rchset[irch];
136*8be712e4SBarry Smith if (marker[inode] < 0) goto L700;
137*8be712e4SBarry Smith
138*8be712e4SBarry Smith marker[inode] = 0;
139*8be712e4SBarry Smith ndeg = deg[inode];
140*8be712e4SBarry Smith if (ndeg < mindeg) mindeg = ndeg;
141*8be712e4SBarry Smith if (ndeg > thresh) goto L700;
142*8be712e4SBarry Smith mindeg = thresh;
143*8be712e4SBarry Smith thresh = ndeg;
144*8be712e4SBarry Smith search = invp[inode];
145*8be712e4SBarry Smith L700:;
146*8be712e4SBarry Smith }
147*8be712e4SBarry Smith if (nhdsze > 0) PetscCall(SPARSEPACKqmdqt(&node, &xadj[1], &adjncy[1], &marker[1], &rchsze, &rchset[1], &nbrhd[1]));
148*8be712e4SBarry Smith L800:
149*8be712e4SBarry Smith if (num < *neqns) goto L300;
150*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
151*8be712e4SBarry Smith }
152