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