Actual source code: qmdmrg.c


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

  4: #include <petscsys.h>
  5: #include <petsc/private/matorderimpl.h>

  7: /******************************************************************/
  8: /***********     QMDMRG ..... QUOT MIN DEG MERGE       ************/
  9: /******************************************************************/
 10: /*    PURPOSE - THIS ROUTINE MERGES INDISTINGUISHABLE NODES IN   */
 11: /*              THE MINIMUM DEGREE ORDERING ALGORITHM.           */
 12: /*              IT ALSO COMPUTES THE NEW DEGREES OF THESE        */
 13: /*              NEW SUPERNODES.                                  */
 14: /*                                                               */
 15: /*    INPUT PARAMETERS -                                         */
 16: /*       (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE.               */
 17: /*       DEG0 - THE NUMBER OF NODES IN THE GIVEN SET.            */
 18: /*       (NHDSZE, NBRHD) - THE SET OF ELIMINATED SUPERNODES      */
 19: /*              ADJACENT TO SOME NODES IN THE SET.               */
 20: /*                                                               */
 21: /*    UPDATED PARAMETERS -                                       */
 22: /*       DEG - THE DEGREE VECTOR.                                */
 23: /*       QSIZE - SIZE OF INDISTINGUISHABLE NODES.                */
 24: /*       QLINK - LINKED LIST FOR INDISTINGUISHABLE NODES.        */
 25: /*       MARKER - THE GIVEN SET IS GIVEN BY THOSE NODES WITH     */
 26: /*              MARKER VALUE SET TO 1.  THOSE NODES WITH DEGREE  */
 27: /*              UPDATED WILL HAVE MARKER VALUE SET TO 2.         */
 28: /*                                                               */
 29: /*    WORKING PARAMETERS -                                       */
 30: /*       RCHSET - THE REACHABLE SET.                             */
 31: /*       OVRLP -  TEMP VECTOR TO STORE THE INTERSECTION OF TWO   */
 32: /*              REACHABLE SETS.                                  */
 33: /*                                                               */
 34: /*****************************************************************/
 35: PetscErrorCode SPARSEPACKqmdmrg(const PetscInt *xadj,const PetscInt *adjncy, PetscInt *deg,
 36:                                 PetscInt *qsize, PetscInt *qlink, PetscInt *marker, PetscInt *deg0,
 37:                                 PetscInt *nhdsze, PetscInt *nbrhd, PetscInt *rchset, PetscInt *ovrlp)
 38: {
 39:   /* System generated locals */
 40:   PetscInt i__1, i__2, i__3;

 42:   /* Local variables */
 43:   PetscInt head, inhd, irch, node, mark, ilink, root, j, lnode, nabor,
 44:            jstop, jstrt, rchsze, mrgsze, novrlp, iov, deg1;

 46:   /* Parameter adjustments */
 47:   --ovrlp;
 48:   --rchset;
 49:   --nbrhd;
 50:   --marker;
 51:   --qlink;
 52:   --qsize;
 53:   --deg;
 54:   --adjncy;
 55:   --xadj;

 57:   if (*nhdsze <= 0) return 0;
 58:   i__1 = *nhdsze;
 59:   for (inhd = 1; inhd <= i__1; ++inhd) {
 60:     root         = nbrhd[inhd];
 61:     marker[root] = 0;
 62:   }
 63: /*       LOOP THROUGH EACH ELIMINATED SUPERNODE IN THE SET     */
 64: /*       (NHDSZE, NBRHD).                                      */
 65:   i__1 = *nhdsze;
 66:   for (inhd = 1; inhd <= i__1; ++inhd) {
 67:     root         = nbrhd[inhd];
 68:     marker[root] = -1;
 69:     rchsze       = 0;
 70:     novrlp       = 0;
 71:     deg1         = 0;
 72: L200:
 73:     jstrt = xadj[root];
 74:     jstop = xadj[root + 1] - 1;
 75: /*          DETERMINE THE REACHABLE SET AND ITS PETSCINTERSECT-     */
 76: /*          ION WITH THE INPUT REACHABLE SET.                  */
 77:     i__2 = jstop;
 78:     for (j = jstrt; j <= i__2; ++j) {
 79:       nabor = adjncy[j];
 80:       root  = -nabor;
 81:       if (nabor < 0) goto L200;
 82:       else if (!nabor) goto L700;
 83:       else goto L300;
 84: L300:
 85:       mark = marker[nabor];
 86:       if (mark < 0) goto L600;
 87:       else if (!mark) goto L400;
 88:       else goto L500;
 89: L400:
 90:       ++rchsze;
 91:       rchset[rchsze] = nabor;
 92:       deg1          += qsize[nabor];
 93:       marker[nabor]  = 1;
 94:       goto L600;
 95: L500:
 96:       if (mark > 1) goto L600;
 97:       ++novrlp;
 98:       ovrlp[novrlp] = nabor;
 99:       marker[nabor] = 2;
100: L600:
101:       ;
102:     }
103: /*          FROM THE OVERLAPPED SET, DETERMINE THE NODES        */
104: /*          THAT CAN BE MERGED TOGETHER.                        */
105: L700:
106:     head   = 0;
107:     mrgsze = 0;
108:     i__2   = novrlp;
109:     for (iov = 1; iov <= i__2; ++iov) {
110:       node  = ovrlp[iov];
111:       jstrt = xadj[node];
112:       jstop = xadj[node + 1] - 1;
113:       i__3  = jstop;
114:       for (j = jstrt; j <= i__3; ++j) {
115:         nabor = adjncy[j];
116:         if (marker[nabor] != 0) goto L800;
117:         marker[node] = 1;
118:         goto L1100;
119: L800:
120:         ;
121:       }
122: /*             NODE BELONGS TO THE NEW MERGED SUPERNODE.      */
123: /*             UPDATE THE VECTORS QLINK AND QSIZE.            */
124:       mrgsze      += qsize[node];
125:       marker[node] = -1;
126:       lnode        = node;
127: L900:
128:       ilink = qlink[lnode];
129:       if (ilink <= 0) goto L1000;
130:       lnode = ilink;
131:       goto L900;
132: L1000:
133:       qlink[lnode] = head;
134:       head         = node;
135: L1100:
136:       ;
137:     }
138:     if (head <= 0) goto L1200;
139:     qsize[head]  = mrgsze;
140:     deg[head]    = *deg0 + deg1 - 1;
141:     marker[head] = 2;
142: /*          RESET MARKER VALUES.          */
143: L1200:
144:     root         = nbrhd[inhd];
145:     marker[root] = 0;
146:     if (rchsze <= 0) goto L1400;
147:     i__2 = rchsze;
148:     for (irch = 1; irch <= i__2; ++irch) {
149:       node         = rchset[irch];
150:       marker[node] = 0;
151:     }
152: L1400:
153:     ;
154:   }
155:   return 0;
156: }