xref: /petsc/src/mat/graphops/order/degree.c (revision be37439ebbbdb2f81c3420c175a94aa72e59929c)
1 /* degree.f -- translated by f2c (version 19931217).*/
2 
3 #include <petsc/private/matorderimpl.h>
4 
5 /*****************************************************************/
6 /*********     DEGREE ..... DEGREE IN MASKED COMPONENT   *********/
7 /*****************************************************************/
8 
9 /*    PURPOSE - THIS ROUTINE COMPUTES THE DEGREES OF THE NODES*/
10 /*       IN THE CONNECTED COMPONENT SPECIFIED BY MASK AND ROOT*/
11 /*       NODES FOR WHICH MASK IS ZERO ARE IGNORED.*/
12 
13 /*    INPUT PARAMETER -*/
14 /*       ROOT - IS THE INPUT NODE THAT DEFINES THE COMPONENT.*/
15 /*       (XADJ, ADJNCY) - ADJACENCY STRUCTURE PAIR.*/
16 /*       MASK - SPECIFIES A SECTION SUBGRAPH.*/
17 
18 /*    OUTPUT PARAMETERS -*/
19 /*       DEG - ARRAY CONTAINING THE DEGREES OF THE NODES IN*/
20 /*             THE COMPONENT.*/
21 /*       CCSIZE-SIZE OF THE COMPONENT SPECIFIED BY MASK AND ROOT*/
22 /*    WORKING PARAMETER -*/
23 /*       LS - A TEMPORARY VECTOR USED TO STORE THE NODES OF THE*/
24 /*              COMPONENT LEVEL BY LEVEL.*/
25 /*****************************************************************/
SPARSEPACKdegree(const PetscInt * root,const PetscInt * inxadj,const PetscInt * adjncy,PetscInt * mask,PetscInt * deg,PetscInt * ccsize,PetscInt * ls)26 PetscErrorCode SPARSEPACKdegree(const PetscInt *root, const PetscInt *inxadj, const PetscInt *adjncy, PetscInt *mask, PetscInt *deg, PetscInt *ccsize, PetscInt *ls)
27 {
28   PetscInt *xadj = (PetscInt *)inxadj; /* Used as temporary and reset within this function */
29   /* System generated locals */
30   PetscInt i__1, i__2;
31 
32   /* Local variables */
33   PetscInt ideg, node, i, j, jstop, jstrt, lbegin, lvlend, lvsize, nbr;
34   /*       INITIALIZATION ...*/
35   /*       THE ARRAY XADJ IS USED AS A TEMPORARY MARKER TO*/
36   /*       INDICATE WHICH NODES HAVE BEEN CONSIDERED SO FAR.*/
37 
38   PetscFunctionBegin;
39   /* Parameter adjustments */
40   --ls;
41   --deg;
42   --mask;
43   --adjncy;
44   --xadj;
45 
46   ls[1]       = *root;
47   xadj[*root] = -xadj[*root];
48   lvlend      = 0;
49   *ccsize     = 1;
50 /*       LBEGIN IS THE POINTER TO THE BEGINNING OF THE CURRENT*/
51 /*       LEVEL, AND LVLEND POINTS TO THE END OF THIS LEVEL.*/
52 L100:
53   lbegin = lvlend + 1;
54   lvlend = *ccsize;
55   /*       FIND THE DEGREES OF NODES IN THE CURRENT LEVEL,*/
56   /*       AND AT THE SAME TIME, GENERATE THE NEXT LEVEL.*/
57   i__1 = lvlend;
58   for (i = lbegin; i <= i__1; ++i) {
59     node  = ls[i];
60     jstrt = -xadj[node];
61     i__2  = xadj[node + 1];
62     jstop = PetscAbsInt(i__2) - 1;
63     ideg  = 0;
64     if (jstop < jstrt) goto L300;
65     i__2 = jstop;
66     for (j = jstrt; j <= i__2; ++j) {
67       nbr = adjncy[j];
68       if (!mask[nbr]) goto L200;
69       ++ideg;
70       if (xadj[nbr] < 0) goto L200;
71       xadj[nbr] = -xadj[nbr];
72       ++(*ccsize);
73       ls[*ccsize] = nbr;
74     L200:;
75     }
76   L300:
77     deg[node] = ideg;
78   }
79   /*       COMPUTE THE CURRENT LEVEL WIDTH. */
80   /*       IF IT IS NONZERO, GENERATE ANOTHER LEVEL.*/
81   lvsize = *ccsize - lvlend;
82   if (lvsize > 0) goto L100;
83   /*       RESET XADJ TO ITS CORRECT SIGN AND RETURN. */
84   i__1 = *ccsize;
85   for (i = 1; i <= i__1; ++i) {
86     node       = ls[i];
87     xadj[node] = -xadj[node];
88   }
89   PetscFunctionReturn(PETSC_SUCCESS);
90 }
91