xref: /petsc/src/mat/graphops/order/fnroot.c (revision 53673ba54f5aaba04b9d49ab22cf56c7a7461fe9) !
1 /* fnroot.f -- translated by f2c (version 19931217).*/
2 
3 #include <petscsys.h>
4 #include <petsc/private/matorderimpl.h>
5 
6 /*****************************************************************/
7 /********     FNROOT ..... FIND PSEUDO-PERIPHERAL NODE    ********/
8 /*****************************************************************/
9 /*   PURPOSE - FNROOT IMPLEMENTS A MODIFIED VERSION OF THE       */
10 /*      SCHEME BY GIBBS, POOLE, AND STOCKMEYER TO FIND PSEUDO-   */
11 /*      PERIPHERAL NODES.  IT DETERMINES SUCH A NODE FOR THE     */
12 /*      SECTION SUBGRAPH SPECIFIED BY MASK AND ROOT.             */
13 /*   INPUT PARAMETERS -                                          */
14 /*      (XADJ, ADJNCY) - ADJACENCY STRUCTURE PAIR FOR THE GRAPH. */
15 /*      MASK - SPECIFIES A SECTION SUBGRAPH. NODES FOR WHICH     */
16 /*             MASK IS ZERO ARE IGNORED BY FNROOT.              */
17 /*   UPDATED PARAMETER -                                        */
18 /*      ROOT - ON INPUT, IT (ALONG WITH MASK) DEFINES THE       */
19 /*             COMPONENT FOR WHICH A PSEUDO-PERIPHERAL NODE IS  */
20 /*             TO BE FOUND. ON OUTPUT, IT IS THE NODE OBTAINED. */
21 /*                                                              */
22 /*   OUTPUT PARAMETERS -                                        */
23 /*      NLVL - IS THE NUMBER OF LEVELS IN THE LEVEL STRUCTURE   */
24 /*             ROOTED AT THE NODE ROOT.                         */
25 /*      (XLS,LS) - THE LEVEL STRUCTURE ARRAY PAIR CONTAINING    */
26 /*                 THE LEVEL STRUCTURE FOUND.                   */
27 /*                                                              */
28 /*   PROGRAM SUBROUTINES -                                      */
29 /*      ROOTLS.                                                 */
30 /*                                                              */
31 /****************************************************************/
SPARSEPACKfnroot(PetscInt * root,const PetscInt * xadj,const PetscInt * adjncy,PetscInt * mask,PetscInt * nlvl,PetscInt * xls,PetscInt * ls)32 PetscErrorCode SPARSEPACKfnroot(PetscInt *root, const PetscInt *xadj, const PetscInt *adjncy, PetscInt *mask, PetscInt *nlvl, PetscInt *xls, PetscInt *ls)
33 {
34   /* System generated locals */
35   PetscInt i__1, i__2;
36 
37   /* Local variables */
38   PetscInt ndeg, node, j, k, nabor, kstop, jstrt, kstrt, mindeg, ccsize, nunlvl;
39   /*       DETERMINE THE LEVEL STRUCTURE ROOTED AT ROOT. */
40 
41   PetscFunctionBegin;
42   /* Parameter adjustments */
43   --ls;
44   --xls;
45   --mask;
46   --adjncy;
47   --xadj;
48 
49   PetscCall(SPARSEPACKrootls(root, &xadj[1], &adjncy[1], &mask[1], nlvl, &xls[1], &ls[1]));
50   ccsize = xls[*nlvl + 1] - 1;
51   if (*nlvl == 1 || *nlvl == ccsize) PetscFunctionReturn(PETSC_SUCCESS);
52 
53 /*       PICK A NODE WITH MINIMUM DEGREE FROM THE LAST LEVEL.*/
54 L100:
55   jstrt  = xls[*nlvl];
56   mindeg = ccsize;
57   *root  = ls[jstrt];
58   if (ccsize == jstrt) goto L400;
59   i__1 = ccsize;
60   for (j = jstrt; j <= i__1; ++j) {
61     node  = ls[j];
62     ndeg  = 0;
63     kstrt = xadj[node];
64     kstop = xadj[node + 1] - 1;
65     i__2  = kstop;
66     for (k = kstrt; k <= i__2; ++k) {
67       nabor = adjncy[k];
68       if (mask[nabor] > 0) ++ndeg;
69     }
70     if (ndeg >= mindeg) goto L300;
71     *root  = node;
72     mindeg = ndeg;
73   L300:;
74   }
75 /*       AND GENERATE ITS ROOTED LEVEL STRUCTURE.*/
76 L400:
77   PetscCall(SPARSEPACKrootls(root, &xadj[1], &adjncy[1], &mask[1], &nunlvl, &xls[1], &ls[1]));
78   if (nunlvl <= *nlvl) PetscFunctionReturn(PETSC_SUCCESS);
79   *nlvl = nunlvl;
80   if (*nlvl < ccsize) goto L100;
81   PetscFunctionReturn(PETSC_SUCCESS);
82 }
83