Actual source code: degree.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: /* degree.f -- translated by f2c (version 19931217).*/

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

  6: /*****************************************************************/
  7: /*********     DEGREE ..... DEGREE IN MASKED COMPONENT   *********/
  8: /*****************************************************************/

 10: /*    PURPOSE - THIS ROUTINE COMPUTES THE DEGREES OF THE NODES*/
 11: /*       IN THE CONNECTED COMPONENT SPECIFIED BY MASK AND ../../..*/
 12: /*       NODES FOR WHICH MASK IS ZERO ARE IGNORED.*/

 14: /*    INPUT PARAMETER -*/
 15: /*       ../../.. - IS THE INPUT NODE THAT DEFINES THE COMPONENT.*/
 16: /*       (XADJ, ADJNCY) - ADJACENCY STRUCTURE PAIR.*/
 17: /*       MASK - SPECIFIES A SECTION SUBGRAPH.*/

 19: /*    OUTPUT PARAMETERS -*/
 20: /*       DEG - ARRAY CONTAINING THE DEGREES OF THE NODES IN*/
 21: /*             THE COMPONENT.*/
 22: /*       CCSIZE-SIZE OF THE COMPONENT SPECIFED BY MASK AND ../../..*/
 23: /*    WORKING PARAMETER -*/
 24: /*       LS - A TEMPORARY VECTOR USED TO STORE THE NODES OF THE*/
 25: /*              COMPONENT LEVEL BY LEVEL.*/
 26: /*****************************************************************/
 29: PetscErrorCode SPARSEPACKdegree(const PetscInt *root,const PetscInt *inxadj,const PetscInt *adjncy,PetscInt *mask,PetscInt *deg,PetscInt *ccsize,PetscInt *ls)
 30: {
 31:   PetscInt *xadj = (PetscInt*)inxadj; /* Used as temporary and reset within this function */
 32:   /* System generated locals */
 33:   PetscInt i__1,i__2;

 35:   /* Local variables */
 36:   PetscInt ideg,node,i,j,jstop,jstrt,lbegin,lvlend,lvsize,nbr;
 37: /*       INITIALIZATION ...*/
 38: /*       THE ARRAY XADJ IS USED AS A TEMPORARY MARKER TO*/
 39: /*       INDICATE WHICH NODES HAVE BEEN CONSIDERED SO FAR.*/

 42:   /* Parameter adjustments */
 43:   --ls;
 44:   --deg;
 45:   --mask;
 46:   --adjncy;
 47:   --xadj;

 49:   ls[1]       = *root;
 50:   xadj[*root] = -xadj[*root];
 51:   lvlend      = 0;
 52:   *ccsize     = 1;
 53: /*       LBEGIN IS THE POINTER TO THE BEGINNING OF THE CURRENT*/
 54: /*       LEVEL, AND LVLEND POINTS TO THE END OF THIS LEVEL.*/
 55: L100:
 56:   lbegin = lvlend + 1;
 57:   lvlend = *ccsize;
 58: /*       FIND THE DEGREES OF NODES IN THE CURRENT LEVEL,*/
 59: /*       AND AT THE SAME TIME, GENERATE THE NEXT LEVEL.*/
 60:   i__1 = lvlend;
 61:   for (i = lbegin; i <= i__1; ++i) {
 62:     node  = ls[i];
 63:     jstrt = -xadj[node];
 64:     jstop = (i__2 = xadj[node + 1], (PetscInt)PetscAbsInt(i__2)) - 1;
 65:     ideg  = 0;
 66:     if (jstop < jstrt) goto L300;
 67:     i__2 = jstop;
 68:     for (j = jstrt; j <= i__2; ++j) {
 69:       nbr = adjncy[j];
 70:       if (!mask[nbr]) goto L200;
 71:       ++ideg;
 72:       if (xadj[nbr] < 0) goto L200;
 73:       xadj[nbr] = -xadj[nbr];
 74:       ++(*ccsize);
 75:       ls[*ccsize] = nbr;
 76: L200:
 77:       ;
 78:     }
 79: L300:
 80:     deg[node] = ideg;
 81:   }
 82: /*       COMPUTE THE CURRENT LEVEL WIDTH. */
 83: /*       IF IT IS NONZERO, GENERATE ANOTHER LEVEL.*/
 84:   lvsize = *ccsize - lvlend;
 85:   if (lvsize > 0) goto L100;
 86: /*       RESET XADJ TO ITS CORRECT SIGN AND RETURN. */
 87: /*       ------------------------------------------*/
 88:   i__1 = *ccsize;
 89:   for (i = 1; i <= i__1; ++i) {
 90:     node       = ls[i];
 91:     xadj[node] = -xadj[node];
 92:   }
 93:   return(0);
 94: }