Actual source code: fn1wd.c
1: /*$Id: fn1wd.c,v 1.23 2001/03/23 23:22:51 balay Exp $*/
2: /* fn1wd.f -- translated by f2c (version 19931217).*/
4: #include petsc.h
5: #include src/mat/order/order.h
7: /*****************************************************************/
8: /******** FN1WD ..... FIND ONE-WAY DISSECTORS *********/
9: /*****************************************************************/
10: /* PURPOSE - THIS SUBROUTINE FINDS ONE-WAY DISSECTORS OF */
11: /* A CONNECTED COMPONENT SPECIFIED BY MASK AND ../../... */
12: /* */
13: /* INPUT PARAMETERS - */
14: /* ../../.. - A NODE THAT DEFINES (ALONG WITH MASK) THE */
15: /* COMPONENT TO BE PROCESSED. */
16: /* (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE. */
17: /* */
18: /* OUTPUT PARAMETERS - */
19: /* NSEP - NUMBER OF NODES IN THE ONE-WAY DISSECTORS. */
20: /* SEP - VECTOR CONTAINING THE DISSECTOR NODES. */
21: /* */
22: /* UPDATED PARAMETER - */
23: /* MASK - NODES IN THE DISSECTOR HAVE THEIR MASK VALUES */
24: /* SET TO ZERO. */
25: /* */
26: /* WORKING PARAMETERS- */
27: /* (XLS, LS) - LEVEL STRUCTURE USED BY THE ROUTINE FN../../... */
28: /* */
29: /* PROGRAM SUBROUTINE - */
30: /* FN../../... */
31: /*****************************************************************/
34: int SPARSEPACKfn1wd(int *root, int *xadj, int *adjncy,
35: int *mask, int *nsep, int *sep, int *nlvl, int *
36: xls, int *ls)
37: {
38: /* System generated locals */
39: int i__1, i__2;
41: /* Local variables */
42: int node, i, j, k;
43: PetscReal width, fnlvl;
44: int kstop, kstrt, lp1beg, lp1end;
45: PetscReal deltp1;
46: int lvlbeg, lvlend;
47: EXTERN int SPARSEPACKfnroot(int *, int *, int *,
48: int *, int *, int *, int *);
49: int nbr, lvl;
52: /* Parameter adjustments */
53: --ls;
54: --xls;
55: --sep;
56: --mask;
57: --adjncy;
58: --xadj;
60: SPARSEPACKfnroot(root, &xadj[1], &adjncy[1], &mask[1], nlvl, &xls[1], &ls[1]);
61: fnlvl = (PetscReal) (*nlvl);
62: *nsep = xls[*nlvl + 1] - 1;
63: width = (PetscReal) (*nsep) / fnlvl;
64: deltp1 = sqrt((width * 3.f + 13.f) / 2.f) + 1.f;
65: if (*nsep >= 50 && deltp1 <= fnlvl * .5f) {
66: goto L300;
67: }
68: /* THE COMPONENT IS TOO SMALL, OR THE LEVEL STRUCTURE */
69: /* IS VERY LONG AND NARROW. RETURN THE WHOLE COMPONENT.*/
70: i__1 = *nsep;
71: for (i = 1; i <= i__1; ++i) {
72: node = ls[i];
73: sep[i] = node;
74: mask[node] = 0;
75: }
76: return(0);
77: /* FIND THE PARALLEL DISSECTORS.*/
78: L300:
79: *nsep = 0;
80: i = 0;
81: L400:
82: ++i;
83: lvl = (int)((PetscReal) i * deltp1 + .5f);
84: if (lvl >= *nlvl) {
85: return(0);
86: }
87: lvlbeg = xls[lvl];
88: lp1beg = xls[lvl + 1];
89: lvlend = lp1beg - 1;
90: lp1end = xls[lvl + 2] - 1;
91: i__1 = lp1end;
92: for (j = lp1beg; j <= i__1; ++j) {
93: node = ls[j];
94: xadj[node] = -xadj[node];
95: }
96: /* NODES IN LEVEL LVL ARE CHOSEN TO FORM DISSECTOR. */
97: /* INCLUDE ONLY THOSE WITH NEIGHBORS IN LVL+1 LEVEL. */
98: /* XADJ IS USED TEMPORARILY TO MARK NODES IN LVL+1. */
99: i__1 = lvlend;
100: for (j = lvlbeg; j <= i__1; ++j) {
101: node = ls[j];
102: kstrt = xadj[node];
103: kstop = (i__2 = xadj[node + 1], (int)PetscAbsInt(i__2)) - 1;
104: i__2 = kstop;
105: for (k = kstrt; k <= i__2; ++k) {
106: nbr = adjncy[k];
107: if (xadj[nbr] > 0) {
108: goto L600;
109: }
110: ++(*nsep);
111: sep[*nsep] = node;
112: mask[node] = 0;
113: goto L700;
114: L600:
115: ;
116: }
117: L700:
118: ;
119: }
120: i__1 = lp1end;
121: for (j = lp1beg; j <= i__1; ++j) {
122: node = ls[j];
123: xadj[node] = -xadj[node];
124: }
125: goto L400;
126: }