Actual source code: degree.c
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 SPECIFIED BY MASK AND ../../..*/
23: /* WORKING PARAMETER -*/
24: /* LS - A TEMPORARY VECTOR USED TO STORE THE NODES OF THE*/
25: /* COMPONENT LEVEL BY LEVEL.*/
26: /*****************************************************************/
27: PetscErrorCode SPARSEPACKdegree(const PetscInt *root, const PetscInt *inxadj, const PetscInt *adjncy, PetscInt *mask, PetscInt *deg, PetscInt *ccsize, PetscInt *ls)
28: {
29: PetscInt *xadj = (PetscInt *)inxadj; /* Used as temporary and reset within this function */
30: /* System generated locals */
31: PetscInt i__1, i__2;
33: /* Local variables */
34: PetscInt ideg, node, i, j, jstop, jstrt, lbegin, lvlend, lvsize, nbr;
35: /* INITIALIZATION ...*/
36: /* THE ARRAY XADJ IS USED AS A TEMPORARY MARKER TO*/
37: /* INDICATE WHICH NODES HAVE BEEN CONSIDERED SO FAR.*/
39: /* Parameter adjustments */
40: --ls;
41: --deg;
42: --mask;
43: --adjncy;
44: --xadj;
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 = (PetscInt)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: /* ------------------------------------------*/
85: i__1 = *ccsize;
86: for (i = 1; i <= i__1; ++i) {
87: node = ls[i];
88: xadj[node] = -xadj[node];
89: }
90: return 0;
91: }