Actual source code: rootls.c


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

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

  7: /*****************************************************************/
  8: /*********     ../../..LS ..... ../../..ED LEVEL STRUCTURE      **********/
  9: /*****************************************************************/
 10: /*    PURPOSE - ../../..LS GENERATES THE LEVEL STRUCTURE ../../..ED */
 11: /*       AT THE INPUT NODE CALLED ../../... ONLY THOSE NODES FOR*/
 12: /*       WHICH MASK IS NONZERO WILL BE CONSIDERED.*/
 13: /*                                                */
 14: /*    INPUT PARAMETERS -                          */
 15: /*       ../../.. - THE NODE AT WHICH THE LEVEL STRUCTURE IS TO*/
 16: /*              BE ../../..ED.*/
 17: /*       (XADJ, ADJNCY) - ADJACENCY STRUCTURE PAIR FOR THE*/
 18: /*              GIVEN GRAPH.*/
 19: /*       MASK - IS USED TO SPECIFY A SECTION SUBGRAPH. NODES*/
 20: /*              WITH MASK(I)=0 ARE IGNORED.*/
 21: /*    OUTPUT PARAMETERS -*/
 22: /*       NLVL - IS THE NUMBER OF LEVELS IN THE LEVEL STRUCTURE.*/
 23: /*       (XLS, LS) - ARRAY PAIR FOR THE ../../..ED LEVEL STRUCTURE.*/
 24: /*****************************************************************/
 25: PetscErrorCode SPARSEPACKrootls(const PetscInt *root, const PetscInt *xadj, const PetscInt *adjncy, PetscInt *mask, PetscInt *nlvl, PetscInt *xls, PetscInt *ls)
 26: {
 27:   /* System generated locals */
 28:   PetscInt i__1, i__2;

 30:   /* Local variables */
 31:   PetscInt node, i, j, jstop, jstrt, lbegin, ccsize, lvlend, lvsize, nbr;

 33:   /*       INITIALIZATION ...*/

 35:   PetscFunctionBegin;
 36:   /* Parameter adjustments */
 37:   --ls;
 38:   --xls;
 39:   --mask;
 40:   --adjncy;
 41:   --xadj;

 43:   mask[*root] = 0;
 44:   ls[1]       = *root;
 45:   *nlvl       = 0;
 46:   lvlend      = 0;
 47:   ccsize      = 1;
 48: /*       LBEGIN IS THE POINTER TO THE BEGINNING OF THE CURRENT*/
 49: /*       LEVEL, AND LVLEND POINTS TO THE END OF THIS LEVEL.*/
 50: L200:
 51:   lbegin = lvlend + 1;
 52:   lvlend = ccsize;
 53:   ++(*nlvl);
 54:   xls[*nlvl] = lbegin;
 55:   /*       GENERATE THE NEXT LEVEL BY FINDING ALL THE MASKED */
 56:   /*       NEIGHBORS OF NODES IN THE CURRENT LEVEL.*/
 57:   i__1 = lvlend;
 58:   for (i = lbegin; i <= i__1; ++i) {
 59:     node  = ls[i];
 60:     jstrt = xadj[node];
 61:     jstop = xadj[node + 1] - 1;
 62:     if (jstop < jstrt) goto L400;
 63:     i__2 = jstop;
 64:     for (j = jstrt; j <= i__2; ++j) {
 65:       nbr = adjncy[j];
 66:       if (!mask[nbr]) goto L300;
 67:       ++ccsize;
 68:       ls[ccsize] = nbr;
 69:       mask[nbr]  = 0;
 70:     L300:;
 71:     }
 72:   L400:;
 73:   }
 74:   /*       COMPUTE THE CURRENT LEVEL WIDTH.*/
 75:   /*       IF IT IS NONZERO, GENERATE THE NEXT LEVEL.*/
 76:   lvsize = ccsize - lvlend;
 77:   if (lvsize > 0) goto L200;
 78:   /*       RESET MASK TO ONE FOR THE NODES IN THE LEVEL STRUCTURE.*/
 79:   xls[*nlvl + 1] = lvlend + 1;
 80:   i__1           = ccsize;
 81:   for (i = 1; i <= i__1; ++i) {
 82:     node       = ls[i];
 83:     mask[node] = 1;
 84:   }
 85:   PetscFunctionReturn(PETSC_SUCCESS);
 86: }