Actual source code: qmdupd.c


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

  4: #include <petscsys.h>
  5: #include <petsc/private/matorderimpl.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: /******************************************************************/
 33: PetscErrorCode SPARSEPACKqmdupd(const PetscInt *xadj,const PetscInt *adjncy,const PetscInt *nlist,const PetscInt *list, PetscInt *deg,
 34:                                 PetscInt *qsize, PetscInt *qlink, PetscInt *marker, PetscInt *rchset, PetscInt *nbrhd)
 35: {
 36:   /* System generated locals */
 37:   PetscInt i__1, i__2;

 39:   /* Local variables */
 40:   PetscInt inhd, irch, node, mark, j, inode, nabor, jstop, jstrt, il;
 41:   PetscInt nhdsze, rchsze, deg0, deg1;

 43: /*       FIND ALL ELIMINATED SUPERNODES THAT ARE ADJACENT*/
 44: /*       TO SOME NODES IN THE GIVEN LIST. PUT THEM INTO.*/
 45: /*       (NHDSZE, NBRHD). DEG0 CONTAINS THE NUMBER OF*/
 46: /*       NODES IN THE LIST.*/

 48:   /* Parameter adjustments */
 49:   --nbrhd;
 50:   --rchset;
 51:   --marker;
 52:   --qlink;
 53:   --qsize;
 54:   --deg;
 55:   --list;
 56:   --adjncy;
 57:   --xadj;

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