Actual source code: fndsep.c
petsc-3.6.1 2015-08-06
2: /* fndsep.f -- translated by f2c (version 19931217).
3: */
5: #include <petsc/private/matorderimpl.h>
7: /*****************************************************************/
8: /************* FNDSEP ..... FIND SEPARATOR *************/
9: /*****************************************************************/
10: /* PURPOSE - THIS ROUTINE IS USED TO FIND A SMALL */
11: /* SEPARATOR FOR A CONNECTED COMPONENT SPECIFIED */
12: /* BY MASK IN THE GIVEN GRAPH. */
13: /* */
14: /* INPUT PARAMETERS - */
15: /* ../../.. - IS THE NODE THAT DETERMINES THE MASKED */
16: /* COMPONENT. */
17: /* (XADJ, ADJNCY) - THE ADJACENCY STRUCTURE PAIR. */
18: /* */
19: /* OUTPUT PARAMETERS - */
20: /* NSEP - NUMBER OF VARIABLES IN THE SEPARATOR. */
21: /* SEP - VECTOR CONTAINING THE SEPARATOR NODES. */
22: /* */
23: /* UPDATED PARAMETER - */
24: /* MASK - NODES IN THE SEPARATOR HAVE THEIR MASK */
25: /* VALUES SET TO ZERO. */
26: /* */
27: /* WORKING PARAMETERS - */
28: /* (XLS, LS) - LEVEL STRUCTURE PAIR FOR LEVEL STRUCTURE */
29: /* FOUND BY FN../../... */
30: /* */
31: /* PROGRAM SUBROUTINES - */
32: /* FN../../... */
33: /* */
34: /*****************************************************************/
37: PetscErrorCode SPARSEPACKfndsep(PetscInt *root,const PetscInt *inxadj,const PetscInt *adjncy,PetscInt *mask, PetscInt *nsep, PetscInt *sep, PetscInt *xls, PetscInt *ls)
38: {
39: PetscInt *xadj = (PetscInt*)inxadj; /* Used as temporary and reset within this function */
41: /* System generated locals */
42: PetscInt i__1, i__2;
44: /* Local variables */
45: PetscInt node, nlvl, i, j, jstop, jstrt, mp1beg, mp1end, midbeg, midend, midlvl;
46: PetscInt nbr;
49: /* Parameter adjustments */
50: --ls;
51: --xls;
52: --sep;
53: --mask;
54: --adjncy;
55: --xadj;
57: SPARSEPACKfnroot(root, &xadj[1], &adjncy[1], &mask[1], &nlvl, &xls[1], &ls[1]);
58: /* IF THE NUMBER OF LEVELS IS LESS THAN 3, RETURN */
59: /* THE WHOLE COMPONENT AS THE SEPARATOR.*/
60: if (nlvl >= 3) goto L200;
61: *nsep = xls[nlvl + 1] - 1;
62: i__1 = *nsep;
63: for (i = 1; i <= i__1; ++i) {
64: node = ls[i];
65: sep[i] = node;
66: mask[node] = 0;
67: }
68: return(0);
69: /* FIND THE MIDDLE LEVEL OF THE ../../..ED LEVEL STRUCTURE.*/
70: L200:
71: midlvl = (nlvl + 2) / 2;
72: midbeg = xls[midlvl];
73: mp1beg = xls[midlvl + 1];
74: midend = mp1beg - 1;
75: mp1end = xls[midlvl + 2] - 1;
76: /* THE SEPARATOR IS OBTAINED BY INCLUDING ONLY THOSE*/
77: /* MIDDLE-LEVEL NODES WITH NEIGHBORS IN THE MIDDLE+1*/
78: /* LEVEL. XADJ IS USED TEMPORARILY TO MARK THOSE*/
79: /* NODES IN THE MIDDLE+1 LEVEL.*/
80: i__1 = mp1end;
81: for (i = mp1beg; i <= i__1; ++i) {
82: node = ls[i];
83: xadj[node] = -xadj[node];
84: }
85: *nsep = 0;
86: i__1 = midend;
87: for (i = midbeg; i <= i__1; ++i) {
88: node = ls[i];
89: jstrt = xadj[node];
90: jstop = (i__2 = xadj[node + 1], (PetscInt)PetscAbsInt(i__2)) - 1;
91: i__2 = jstop;
92: for (j = jstrt; j <= i__2; ++j) {
93: nbr = adjncy[j];
94: if (xadj[nbr] > 0) goto L400;
95: ++(*nsep);
96: sep[*nsep] = node;
97: mask[node] = 0;
98: goto L500;
99: L400:
100: ;
101: }
102: L500:
103: ;
104: }
105: /* RESET XADJ TO ITS CORRECT SIGN.*/
106: i__1 = mp1end;
107: for (i = mp1beg; i <= i__1; ++i) {
108: node = ls[i];
109: xadj[node] = -xadj[node];
110: }
111: return(0);
112: }