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