Actual source code: fnroot.c
1: /*$Id: fnroot.c,v 1.16 2001/03/23 23:22:51 balay Exp $*/
2: /* fnroot.f -- translated by f2c (version 19931217).*/
4: #include petsc.h
6: /*****************************************************************/
7: /******** FN../../.. ..... FIND PSEUDO-PERIPHERAL NODE ********/
8: /*****************************************************************/
9: /* PURPOSE - FN../../.. IMPLEMENTS A MODIFIED VERSION OF THE */
10: /* SCHEME BY GIBBS, POOLE, AND STOCKMEYER TO FIND PSEUDO- */
11: /* PERIPHERAL NODES. IT DETERMINES SUCH A NODE FOR THE */
12: /* SECTION SUBGRAPH SPECIFIED BY MASK AND ../../... */
13: /* INPUT PARAMETERS - */
14: /* (XADJ, ADJNCY) - ADJACENCY STRUCTURE PAIR FOR THE GRAPH. */
15: /* MASK - SPECIFIES A SECTION SUBGRAPH. NODES FOR WHICH */
16: /* MASK IS ZERO ARE IGNORED BY FN../../... */
17: /* UPDATED PARAMETER - */
18: /* ../../.. - ON INPUT, IT (ALONG WITH MASK) DEFINES THE */
19: /* COMPONENT FOR WHICH A PSEUDO-PERIPHERAL NODE IS */
20: /* TO BE FOUND. ON OUTPUT, IT IS THE NODE OBTAINED. */
21: /* */
22: /* OUTPUT PARAMETERS - */
23: /* NLVL - IS THE NUMBER OF LEVELS IN THE LEVEL STRUCTURE */
24: /* ../../..ED AT THE NODE ../../... */
25: /* (XLS,LS) - THE LEVEL STRUCTURE ARRAY PAIR CONTAINING */
26: /* THE LEVEL STRUCTURE FOUND. */
27: /* */
28: /* PROGRAM SUBROUTINES - */
29: /* ../../..LS. */
30: /* */
31: /****************************************************************/
34: int SPARSEPACKfnroot(int *root, int *xadj, int *adjncy,
35: int *mask, int *nlvl, int *xls, int *ls)
36: {
37: /* System generated locals */
38: int i__1, i__2;
40: /* Local variables */
41: int ndeg, node, j, k, nabor, kstop, jstrt, kstrt, mindeg,
42: ccsize, nunlvl;
43: EXTERN int SPARSEPACKrootls(int *, int *, int *,
44: int *, int *, int *, int *);
45: /* DETERMINE THE LEVEL STRUCTURE ../../..ED AT ../../... */
48: /* Parameter adjustments */
49: --ls;
50: --xls;
51: --mask;
52: --adjncy;
53: --xadj;
55: SPARSEPACKrootls(root, &xadj[1], &adjncy[1], &mask[1], nlvl, &xls[1], &ls[1]);
56: ccsize = xls[*nlvl + 1] - 1;
57: if (*nlvl == 1 || *nlvl == ccsize) {
58: return(0);
59: }
60: /* PICK A NODE WITH MINIMUM DEGREE FROM THE LAST LEVEL.*/
61: L100:
62: jstrt = xls[*nlvl];
63: mindeg = ccsize;
64: *root = ls[jstrt];
65: if (ccsize == jstrt) {
66: goto L400;
67: }
68: i__1 = ccsize;
69: for (j = jstrt; j <= i__1; ++j) {
70: node = ls[j];
71: ndeg = 0;
72: kstrt = xadj[node];
73: kstop = xadj[node + 1] - 1;
74: i__2 = kstop;
75: for (k = kstrt; k <= i__2; ++k) {
76: nabor = adjncy[k];
77: if (mask[nabor] > 0) {
78: ++ndeg;
79: }
80: }
81: if (ndeg >= mindeg) {
82: goto L300;
83: }
84: *root = node;
85: mindeg = ndeg;
86: L300:
87: ;
88: }
89: /* AND GENERATE ITS ../../..ED LEVEL STRUCTURE.*/
90: L400:
91: SPARSEPACKrootls(root, &xadj[1], &adjncy[1], &mask[1], &nunlvl, &xls[1], &ls[1]);
92: if (nunlvl <= *nlvl) {
93: return(0);
94: }
95: *nlvl = nunlvl;
96: if (*nlvl < ccsize) {
97: goto L100;
98: }
99: return(0);
100: }