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