Actual source code: fn1wd.c
petsc-3.6.1 2015-08-06
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: /*****************************************************************/
33: PetscErrorCode SPARSEPACKfn1wd(PetscInt *root,const PetscInt *inxadj,const PetscInt *adjncy,
34: PetscInt *mask, PetscInt *nsep, PetscInt *sep, PetscInt *nlvl, PetscInt *
35: xls, PetscInt *ls)
36: {
37: PetscInt *xadj = (PetscInt*)inxadj; /* Used as temporary and reset */
38: /* System generated locals */
39: PetscInt i__1, i__2;
41: /* Local variables */
42: PetscInt node, i, j, k;
43: PetscReal width, fnlvl;
44: PetscInt kstop, kstrt, lp1beg, lp1end;
45: PetscReal deltp1;
46: PetscInt lvlbeg, lvlend;
47: PetscInt nbr, lvl;
50: /* Parameter adjustments */
51: --ls;
52: --xls;
53: --sep;
54: --mask;
55: --adjncy;
56: --xadj;
58: SPARSEPACKfnroot(root, &xadj[1], &adjncy[1], &mask[1], nlvl, &xls[1], &ls[1]);
59: fnlvl = (PetscReal) (*nlvl);
60: *nsep = xls[*nlvl + 1] - 1;
61: width = (PetscReal) (*nsep) / fnlvl;
62: deltp1 = PetscSqrtReal((width * 3. + 13.) / 2.) + 1.;
63: if (*nsep >= 50 && deltp1 <= fnlvl * .5f) goto L300;
65: /* THE COMPONENT IS TOO SMALL, OR THE LEVEL STRUCTURE */
66: /* IS VERY LONG AND NARROW. RETURN THE WHOLE COMPONENT.*/
67: i__1 = *nsep;
68: for (i = 1; i <= i__1; ++i) {
69: node = ls[i];
70: sep[i] = node;
71: mask[node] = 0;
72: }
73: return(0);
74: /* FIND THE PARALLEL DISSECTORS.*/
75: L300:
76: *nsep = 0;
77: i = 0;
78: L400:
79: ++i;
80: lvl = (PetscInt)((PetscReal) i * deltp1 + .5f);
81: if (lvl >= *nlvl) return(0);
82: lvlbeg = xls[lvl];
83: lp1beg = xls[lvl + 1];
84: lvlend = lp1beg - 1;
85: lp1end = xls[lvl + 2] - 1;
86: i__1 = lp1end;
87: for (j = lp1beg; j <= i__1; ++j) {
88: node = ls[j];
89: xadj[node] = -xadj[node];
90: }
91: /* NODES IN LEVEL LVL ARE CHOSEN TO FORM DISSECTOR. */
92: /* INCLUDE ONLY THOSE WITH NEIGHBORS IN LVL+1 LEVEL. */
93: /* XADJ IS USED TEMPORARILY TO MARK NODES IN LVL+1. */
94: i__1 = lvlend;
95: for (j = lvlbeg; j <= i__1; ++j) {
96: node = ls[j];
97: kstrt = xadj[node];
98: kstop = (i__2 = xadj[node + 1], (PetscInt)PetscAbsInt(i__2)) - 1;
99: i__2 = kstop;
100: for (k = kstrt; k <= i__2; ++k) {
101: nbr = adjncy[k];
102: if (xadj[nbr] > 0) goto L600;
103: ++(*nsep);
104: sep[*nsep] = node;
105: mask[node] = 0;
106: goto L700;
107: L600:
108: ;
109: }
110: L700:
111: ;
112: }
113: i__1 = lp1end;
114: for (j = lp1beg; j <= i__1; ++j) {
115: node = ls[j];
116: xadj[node] = -xadj[node];
117: }
118: goto L400;
119: }