Actual source code: qmdrch.c

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

 3:  #include petsc.h

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

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

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

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

 23: /*    UPDATED PARAMETERS -*/
 24: /*       MARKER - THE MARKER VECTOR FOR REACH AND NBRHD SETS.*/
 25: /*              GT 0 MEANS THE NODE IS IN REACH SET.*/
 26: /*              LT 0 MEANS THE NODE HAS BEEN MERGED WITH*/
 27: /*              OTHERS IN THE QUOTIENT OR IT IS IN NBRHD SET.*/
 28: /*****************************************************************/
 31: PetscErrorCode SPARSEPACKqmdrch(PetscInt *root, PetscInt *xadj, 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.*/


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

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