Actual source code: rootls.c

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

  4: #include <petscsys.h>

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

 32:     /* Local variables */
 33:     PetscInt node, i, j, jstop, jstrt, lbegin, ccsize, lvlend, lvsize,
 34:             nbr;

 36: /*       INITIALIZATION ...*/


 40:     /* Parameter adjustments */
 41:     --ls;
 42:     --xls;
 43:     --mask;
 44:     --adjncy;
 45:     --xadj;

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