Actual source code: qmdupd.c

  1: #define PETSCMAT_DLL

  3: /* qmdupd.f -- translated by f2c (version 19931217).*/

 5:  #include petsc.h

  7: /******************************************************************/
  8: /***********     QMDUPD ..... QUOT MIN DEG UPDATE      ************/
  9: /******************************************************************/
 10: /******************************************************************/

 12: /*    PURPOSE - THIS ROUTINE PERFORMS DEGREE UPDATE FOR A SET*/
 13: /*       OF NODES IN THE MINIMUM DEGREE ALGORITHM.*/

 15: /*    INPUT PARAMETERS -*/
 16: /*       (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE.*/
 17: /*       (NLIST, LIST) - THE LIST OF NODES WHOSE DEGREE HAS TO*/
 18: /*              BE UPDATED.*/

 20: /*    UPDATED PARAMETERS -*/
 21: /*       DEG - THE DEGREE VECTOR.*/
 22: /*       QSIZE - SIZE OF INDISTINGUISHABLE SUPERNODES.*/
 23: /*       QLINK - LINKED LIST FOR INDISTINGUISHABLE NODES.*/
 24: /*       MARKER - USED TO MARK THOSE NODES IN REACH/NBRHD SETS.*/

 26: /*    WORKING PARAMETERS -*/
 27: /*       RCHSET - THE REACHABLE SET.*/
 28: /*       NBRHD -  THE NEIGHBORHOOD SET.*/

 30: /*    PROGRAM SUBROUTINES -*/
 31: /*       QMDMRG.*/
 32: /******************************************************************/
 35: PetscErrorCode SPARSEPACKqmdupd(PetscInt *xadj, PetscInt *adjncy, PetscInt *nlist, 
 36:         PetscInt *list, PetscInt *deg, PetscInt *qsize, PetscInt *qlink, PetscInt *
 37:         marker, PetscInt *rchset, PetscInt *nbrhd)
 38: {
 39:     /* System generated locals */
 40:     PetscInt i__1, i__2;

 42:     /* Local variables */
 43:     PetscInt inhd, irch, node, mark, j, inode, nabor, jstop, jstrt, il;
 44:     EXTERN PetscErrorCode SPARSEPACKqmdrch(PetscInt*, PetscInt *, PetscInt *,
 45:             PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *),
 46:              SPARSEPACKqmdmrg(PetscInt*, PetscInt *, PetscInt *, PetscInt *, PetscInt *,
 47:             PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
 48:     PetscInt nhdsze, rchsze, deg0, deg1;

 50: /*       FIND ALL ELIMINATED SUPERNODES THAT ARE ADJACENT*/
 51: /*       TO SOME NODES IN THE GIVEN LIST. PUT THEM INTO.*/
 52: /*       (NHDSZE, NBRHD). DEG0 CONTAINS THE NUMBER OF*/
 53: /*       NODES IN THE LIST.*/


 57:     /* Parameter adjustments */
 58:     --nbrhd;
 59:     --rchset;
 60:     --marker;
 61:     --qlink;
 62:     --qsize;
 63:     --deg;
 64:     --list;
 65:     --adjncy;
 66:     --xadj;

 68:     if (*nlist <= 0) {
 69:         return(0);
 70:     }
 71:     deg0 = 0;
 72:     nhdsze = 0;
 73:     i__1 = *nlist;
 74:     for (il = 1; il <= i__1; ++il) {
 75:         node = list[il];
 76:         deg0 += qsize[node];
 77:         jstrt = xadj[node];
 78:         jstop = xadj[node + 1] - 1;
 79:         i__2 = jstop;
 80:         for (j = jstrt; j <= i__2; ++j) {
 81:             nabor = adjncy[j];
 82:             if (marker[nabor] != 0 || deg[nabor] >= 0) {
 83:                 goto L100;
 84:             }
 85:             marker[nabor] = -1;
 86:             ++nhdsze;
 87:             nbrhd[nhdsze] = nabor;
 88: L100:
 89:             ;
 90:         }
 91:     }
 92: /*       MERGE INDISTINGUISHABLE NODES IN THE LIST BY*/
 93: /*       CALLING THE SUBROUTINE QMDMRG.*/
 94:     if (nhdsze > 0) {
 95:         SPARSEPACKqmdmrg(&xadj[1], &adjncy[1], &deg[1], &qsize[1], &qlink[1], &marker[
 96:                 1], &deg0, &nhdsze, &nbrhd[1], &rchset[1], &nbrhd[nhdsze + 1]);
 97:     }
 98: /*       FIND THE NEW DEGREES OF THE NODES THAT HAVE NOT BEEN*/
 99: /*       MERGED.*/
100:     i__1 = *nlist;
101:     for (il = 1; il <= i__1; ++il) {
102:         node = list[il];
103:         mark = marker[node];
104:         if (mark > 1 || mark < 0) {
105:             goto L600;
106:         }
107:         marker[node] = 2;
108:         SPARSEPACKqmdrch(&node, &xadj[1], &adjncy[1], &deg[1], &marker[1], &rchsze, &
109:                 rchset[1], &nhdsze, &nbrhd[1]);
110:         deg1 = deg0;
111:         if (rchsze <= 0) {
112:             goto L400;
113:         }
114:         i__2 = rchsze;
115:         for (irch = 1; irch <= i__2; ++irch) {
116:             inode = rchset[irch];
117:             deg1 += qsize[inode];
118:             marker[inode] = 0;
119:         }
120: L400:
121:         deg[node] = deg1 - 1;
122:         if (nhdsze <= 0) {
123:             goto L600;
124:         }
125:         i__2 = nhdsze;
126:         for (inhd = 1; inhd <= i__2; ++inhd) {
127:             inode = nbrhd[inhd];
128:             marker[inode] = 0;
129:         }
130: L600:
131:         ;
132:     }
133:     return(0);
134: }