Actual source code: qmdrch.c

  2: /* qmdrch.f -- translated by f2c (version 19931217).*/

  4: #include <petscsys.h>

  6: /*****************************************************************/
  7: /**********     QMDRCH ..... QUOT MIN DEG REACH SET    ***********/
  8: /*****************************************************************/

 10: /*    PURPOSE - THIS SUBROUTINE DETERMINES THE REACHABLE SET OF*/
 11: /*       A NODE THROUGH A GIVEN SUBSET.  THE ADJACENCY STRUCTURE*/
 12: /*       IS ASSUMED TO BE STORED IN A QUOTIENT GRAPH FORMAT.*/

 14: /*    INPUT PARAMETERS -*/
 15: /*       ../../.. - THE GIVEN NODE NOT IN THE SUBSET.*/
 16: /*       (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE PAIR.*/
 17: /*       DEG - THE DEGREE VECTOR.  DEG(I) LT 0 MEANS THE NODE*/
 18: /*              BELONGS TO THE GIVEN SUBSET.*/

 20: /*    OUTPUT PARAMETERS -*/
 21: /*       (RCHSZE, RCHSET) - THE REACHABLE SET.*/
 22: /*       (NHDSZE, NBRHD) - THE NEIGHBORHOOD SET.*/

 24: /*    UPDATED PARAMETERS -*/
 25: /*       MARKER - THE MARKER VECTOR FOR REACH AND NBRHD SETS.*/
 26: /*              GT 0 MEANS THE NODE IS IN REACH SET.*/
 27: /*              LT 0 MEANS THE NODE HAS BEEN MERGED WITH*/
 28: /*              OTHERS IN THE QUOTIENT OR IT IS IN NBRHD SET.*/
 29: /*****************************************************************/
 32: PetscErrorCode SPARSEPACKqmdrch(PetscInt *root, PetscInt *xadj, PetscInt *adjncy, 
 33:         PetscInt *deg, PetscInt *marker, PetscInt *rchsze, PetscInt *rchset, 
 34:         PetscInt *nhdsze, PetscInt *nbrhd)
 35: {
 36:     /* System generated locals */
 37:     PetscInt i__1, i__2;

 39:     /* Local variables */
 40:     PetscInt node, i, j, nabor, istop, jstop, istrt, jstrt;

 42: /*       LOOP THROUGH THE NEIGHBORS OF ../../.. IN THE*/
 43: /*       QUOTIENT GRAPH.*/


 47:     /* Parameter adjustments */
 48:     --nbrhd;
 49:     --rchset;
 50:     --marker;
 51:     --deg;
 52:     --adjncy;
 53:     --xadj;

 55:     *nhdsze = 0;
 56:     *rchsze = 0;
 57:     istrt = xadj[*root];
 58:     istop = xadj[*root + 1] - 1;
 59:     if (istop < istrt) {
 60:         return(0);
 61:     }
 62:     i__1 = istop;
 63:     for (i = istrt; i <= i__1; ++i) {
 64:         nabor = adjncy[i];
 65:         if (!nabor) {
 66:             return(0);
 67:         }
 68:         if (marker[nabor] != 0) {
 69:             goto L600;
 70:         }
 71:         if (deg[nabor] < 0) {
 72:             goto L200;
 73:         }
 74: /*                   INCLUDE NABOR INTO THE REACHABLE SET.*/
 75:         ++(*rchsze);
 76:         rchset[*rchsze] = nabor;
 77:         marker[nabor] = 1;
 78:         goto L600;
 79: /*                NABOR HAS BEEN ELIMINATED. FIND NODES*/
 80: /*                REACHABLE FROM IT.*/
 81: L200:
 82:         marker[nabor] = -1;
 83:         ++(*nhdsze);
 84:         nbrhd[*nhdsze] = nabor;
 85: L300:
 86:         jstrt = xadj[nabor];
 87:         jstop = xadj[nabor + 1] - 1;
 88:         i__2 = jstop;
 89:         for (j = jstrt; j <= i__2; ++j) {
 90:             node = adjncy[j];
 91:             nabor = -node;
 92:             if (node < 0) {
 93:                 goto L300;
 94:             } else if (!node) {
 95:                 goto L600;
 96:             } else {
 97:                 goto L400;
 98:             }
 99: L400:
100:             if (marker[node] != 0) {
101:                 goto L500;
102:             }
103:             ++(*rchsze);
104:             rchset[*rchsze] = node;
105:             marker[node] = 1;
106: L500:
107:             ;
108:         }
109: L600:
110:         ;
111:     }
112:     return(0);
113: }