Actual source code: qmdupd.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  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: /******************************************************************/
 35: PetscErrorCode SPARSEPACKqmdupd(const PetscInt *xadj,const PetscInt *adjncy,const PetscInt *nlist,const PetscInt *list, PetscInt *deg,
 36:                                 PetscInt *qsize, PetscInt *qlink, PetscInt *marker, PetscInt *rchset, PetscInt *nbrhd)
 37: {
 38:   /* System generated locals */
 39:   PetscInt i__1, i__2;

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

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


 52:   /* Parameter adjustments */
 53:   --nbrhd;
 54:   --rchset;
 55:   --marker;
 56:   --qlink;
 57:   --qsize;
 58:   --deg;
 59:   --list;
 60:   --adjncy;
 61:   --xadj;

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