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