Actual source code: gen1wd.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

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

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

  7: /*****************************************************************/
  8: /***********     GEN1WD ..... GENERAL ONE-WAY DISSECTION  ********/
  9: /*****************************************************************/

 11: /*    PURPOSE - GEN1WD FINDS A ONE-WAY DISSECTION PARTITIONING*/
 12: /*       FOR A GENERAL GRAPH.  FN1WD IS USED FOR EACH CONNECTED*/
 13: /*       COMPONENT.*/

 15: /*    INPUT PARAMETERS -*/
 16: /*       NEQNS - NUMBER OF EQUATIONS.*/
 17: /*       (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE PAIR.*/

 19: /*    OUTPUT PARAMETERS -*/
 20: /*       (NBLKS, XBLK) - THE PARTITIONING FOUND.*/
 21: /*       PERM - THE ONE-WAY DISSECTION ORDERING.*/

 23: /*    WORKING VECTORS -*/
 24: /*       MASK - IS USED TO MARK VARIABLES THAT HAVE*/
 25: /*              BEEN NUMBERED DURING THE ORDERING PROCESS.*/
 26: /*       (XLS, LS) - LEVEL STRUCTURE USED BY ../../..LS.*/

 28: /*    PROGRAM SUBROUTINES -*/
 29: /*       FN1WD, REVRSE, ../../..LS.*/
 30: /****************************************************************/
 31: PetscErrorCode SPARSEPACKgen1wd(const PetscInt *neqns,const PetscInt *xadj,const PetscInt *adjncy,PetscInt *mask, PetscInt *nblks, PetscInt *xblk, PetscInt *perm, PetscInt *xls, PetscInt *ls)
 32: {
 33:   /* System generated locals */
 34:   PetscInt i__1, i__2, i__3;

 36:   /* Local variables */
 37:   PetscInt node, nsep, lnum, nlvl, root;
 38:   PetscInt i, j, k, ccsize;
 39:   PetscInt num;

 42:   /* Parameter adjustments */
 43:   --ls;
 44:   --xls;
 45:   --perm;
 46:   --xblk;
 47:   --mask;
 48:   --xadj;
 49:   --adjncy;

 51:   i__1 = *neqns;
 52:   for (i = 1; i <= i__1; ++i) mask[i] = 1;
 53:   *nblks = 0;
 54:   num    = 0;
 55:   i__1   = *neqns;
 56:   for (i = 1; i <= i__1; ++i) {
 57:     if (!mask[i]) goto L400;
 58: /*             FIND A ONE-WAY DISSECTOR FOR EACH COMPONENT.*/
 59:     root = i;
 60:     SPARSEPACKfn1wd(&root, &xadj[1], &adjncy[1], &mask[1], &nsep, &perm[num + 1], &nlvl, &xls[1], &ls[1]);
 61:     num += nsep;
 62:     ++(*nblks);
 63:     xblk[*nblks] = *neqns - num + 1;
 64:     ccsize       = xls[nlvl + 1] - 1;
 65: /*             NUMBER THE REMAINING NODES IN THE COMPONENT.*/
 66: /*             EACH COMPONENT IN THE REMAINING SUBGRAPH FORMS*/
 67: /*             A NEW BLOCK IN THE PARTITIONING.*/
 68:     i__2 = ccsize;
 69:     for (j = 1; j <= i__2; ++j) {
 70:       node = ls[j];
 71:       if (!mask[node]) goto L300;
 72:       SPARSEPACKrootls(&node, &xadj[1], &adjncy[1], &mask[1], &nlvl, &xls[1], &perm[num + 1]);
 73:       lnum = num + 1;
 74:       num  = num + xls[nlvl + 1] - 1;
 75:       ++(*nblks);
 76:       xblk[*nblks] = *neqns - num + 1;
 77:       i__3         = num;
 78:       for (k = lnum; k <= i__3; ++k) {
 79:         node       = perm[k];
 80:         mask[node] = 0;
 81:       }
 82:       if (num > *neqns) goto L500;
 83: L300:
 84:       ;
 85:     }
 86: L400:
 87:     ;
 88:   }
 89: /*       SINCE DISSECTORS FOUND FIRST SHOULD BE ORDERED LAST,*/
 90: /*       ROUTINE REVRSE IS CALLED TO ADJUST THE ORDERING*/
 91: /*       VECTOR, AND THE BLOCK INDEX VECTOR.*/
 92: L500:
 93:   SPARSEPACKrevrse(neqns, &perm[1]);
 94:   SPARSEPACKrevrse(nblks, &xblk[1]);
 95:   xblk[*nblks + 1] = *neqns + 1;
 96:   return(0);
 97: }