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