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,
32: PetscInt *deg, PetscInt *marker, PetscInt *rchsze, PetscInt *rchset,
33: PetscInt *nhdsze, PetscInt *nbrhd)
34: {
35: /* System generated locals */
36: PetscInt i__1, i__2;
38: /* Local variables */
39: PetscInt node, i, j, nabor, istop, jstop, istrt, jstrt;
41: /* LOOP THROUGH THE NEIGHBORS OF ../../.. IN THE*/
42: /* QUOTIENT GRAPH.*/
45: /* Parameter adjustments */
46: --nbrhd;
47: --rchset;
48: --marker;
49: --deg;
50: --adjncy;
51: --xadj;
53: *nhdsze = 0;
54: *rchsze = 0;
55: istrt = xadj[*root];
56: istop = xadj[*root + 1] - 1;
57: if (istop < istrt) return(0);
58: i__1 = istop;
59: for (i = istrt; i <= i__1; ++i) {
60: nabor = adjncy[i];
61: if (!nabor) return(0);
62: if (marker[nabor] != 0) goto L600;
63: if (deg[nabor] < 0) goto L200;
65: /* INCLUDE NABOR INTO THE REACHABLE SET.*/
66: ++(*rchsze);
67: rchset[*rchsze] = nabor;
68: marker[nabor] = 1;
69: goto L600;
70: /* NABOR HAS BEEN ELIMINATED. FIND NODES*/
71: /* REACHABLE FROM IT.*/
72: L200:
73: marker[nabor] = -1;
74: ++(*nhdsze);
75: nbrhd[*nhdsze] = nabor;
76: L300:
77: jstrt = xadj[nabor];
78: jstop = xadj[nabor + 1] - 1;
79: i__2 = jstop;
80: for (j = jstrt; j <= i__2; ++j) {
81: node = adjncy[j];
82: nabor = -node;
83: if (node < 0) goto L300;
84: else if (!node) goto L600;
85: else goto L400;
86: L400:
87: if (marker[node] != 0) goto L500;
88: ++(*rchsze);
89: rchset[*rchsze] = node;
90: marker[node] = 1;
91: L500:
92: ;
93: }
94: L600:
95: ;
96: }
97: return(0);
98: }