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:   PetscFunctionBegin;
 43:   /* Parameter adjustments */
 44:   --nbrhd;
 45:   --rchset;
 46:   --marker;
 47:   --deg;
 48:   --adjncy;
 49:   --xadj;

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

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