Actual source code: qmdqt.c

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

  4: #include <petscsys.h>

  6: /***************************************************************/
  7: /********     QMDQT  ..... QUOT MIN DEG QUOT TRANSFORM  ********/
  8: /***************************************************************/

 10: /*    PURPOSE - THIS SUBROUTINE PERFORMS THE QUOTIENT GRAPH  */
 11: /*       TRANSFORMATION AFTER A NODE HAS BEEN ELIMINATED.*/

 13: /*    INPUT PARAMETERS -*/
 14: /*       ../../.. - THE NODE JUST ELIMINATED. IT BECOMES THE*/
 15: /*              REPRESENTATIVE OF THE NEW SUPERNODE.*/
 16: /*       (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE.*/
 17: /*       (RCHSZE, RCHSET) - THE REACHABLE SET OF ../../.. IN THE*/
 18: /*              OLD QUOTIENT GRAPH.*/
 19: /*       NBRHD - THE NEIGHBORHOOD SET WHICH WILL BE MERGED*/
 20: /*              WITH ../../.. TO FORM THE NEW SUPERNODE.*/
 21: /*       MARKER - THE MARKER VECTOR.*/

 23: /*    UPDATED PARAMETER -*/
 24: /*       ADJNCY - BECOMES THE ADJNCY OF THE QUOTIENT GRAPH.*/
 25: /***************************************************************/
 28: PetscErrorCode SPARSEPACKqmdqt(PetscInt *root, PetscInt *xadj, PetscInt *adjncy, 
 29:         PetscInt *marker, PetscInt *rchsze, PetscInt *rchset, PetscInt *nbrhd)
 30: {
 31:     /* System generated locals */
 32:     PetscInt i__1, i__2;

 34:     /* Local variables */
 35:     PetscInt inhd, irch, node, ilink, j, nabor, jstop, jstrt;

 38:     /* Parameter adjustments */
 39:     --nbrhd;
 40:     --rchset;
 41:     --marker;
 42:     --adjncy;
 43:     --xadj;

 45:     irch = 0;
 46:     inhd = 0;
 47:     node = *root;
 48: L100:
 49:     jstrt = xadj[node];
 50:     jstop = xadj[node + 1] - 2;
 51:     if (jstop < jstrt) {
 52:         goto L300;
 53:     }
 54: /*          PLACE REACH NODES INTO THE ADJACENT LIST OF NODE*/
 55:     i__1 = jstop;
 56:     for (j = jstrt; j <= i__1; ++j) {
 57:         ++irch;
 58:         adjncy[j] = rchset[irch];
 59:         if (irch >= *rchsze) {
 60:             goto L400;
 61:         }
 62:     }
 63: /*       LINK TO OTHER SPACE PROVIDED BY THE NBRHD SET.*/
 64: L300:
 65:     ilink = adjncy[jstop + 1];
 66:     node = -ilink;
 67:     if (ilink < 0) {
 68:         goto L100;
 69:     }
 70:     ++inhd;
 71:     node = nbrhd[inhd];
 72:     adjncy[jstop + 1] = -node;
 73:     goto L100;
 74: /*       ALL REACHABLE NODES HAVE BEEN SAVED.  END THE ADJ LIST.*/
 75: /*       ADD ../../.. TO THE NBR LIST OF EACH NODE IN THE REACH SET.*/
 76: L400:
 77:     adjncy[j + 1] = 0;
 78:     i__1 = *rchsze;
 79:     for (irch = 1; irch <= i__1; ++irch) {
 80:         node = rchset[irch];
 81:         if (marker[node] < 0) {
 82:             goto L600;
 83:         }
 84:         jstrt = xadj[node];
 85:         jstop = xadj[node + 1] - 1;
 86:         i__2 = jstop;
 87:         for (j = jstrt; j <= i__2; ++j) {
 88:             nabor = adjncy[j];
 89:             if (marker[nabor] >= 0) {
 90:                 goto L500;
 91:             }
 92:             adjncy[j] = *root;
 93:             goto L600;
 94: L500:
 95:             ;
 96:         }
 97: L600:
 98:         ;
 99:     }
100:     return(0);
101: }