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.*/

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

 52:   *nhdsze = 0;
 53:   *rchsze = 0;
 54:   istrt   = xadj[*root];
 55:   istop   = xadj[*root + 1] - 1;
 56:   if (istop < istrt) return 0;
 57:   i__1 = istop;
 58:   for (i = istrt; i <= i__1; ++i) {
 59:     nabor = adjncy[i];
 60:     if (!nabor) return 0;
 61:     if (marker[nabor] != 0) goto L600;
 62:     if (deg[nabor] < 0) goto L200;

 64: /*                   INCLUDE NABOR INTO THE REACHABLE SET.*/
 65:     ++(*rchsze);
 66:     rchset[*rchsze] = nabor;
 67:     marker[nabor]   = 1;
 68:     goto L600;
 69: /*                NABOR HAS BEEN ELIMINATED. FIND NODES*/
 70: /*                REACHABLE FROM IT.*/
 71: L200:
 72:     marker[nabor] = -1;
 73:     ++(*nhdsze);
 74:     nbrhd[*nhdsze] = nabor;
 75: L300:
 76:     jstrt = xadj[nabor];
 77:     jstop = xadj[nabor + 1] - 1;
 78:     i__2  = jstop;
 79:     for (j = jstrt; j <= i__2; ++j) {
 80:       node  = adjncy[j];
 81:       nabor = -node;
 82:       if (node < 0) goto L300;
 83:       else if (!node) goto L600;
 84:       else goto L400;
 85: L400:
 86:       if (marker[node] != 0) goto L500;
 87:       ++(*rchsze);
 88:       rchset[*rchsze] = node;
 89:       marker[node]    = 1;
 90: L500:
 91:       ;
 92:     }
 93: L600:
 94:     ;
 95:   }
 96:   return 0;
 97: }