xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
12e68589bSMark F. Adams /*
22e68589bSMark F. Adams  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
32e68589bSMark F. Adams  */
42e68589bSMark F. Adams 
52e68589bSMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/
62e68589bSMark F. Adams 
72e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE)
85a449d36SBarry Smith #if !defined(ANSI_DECLARATORS)
95a449d36SBarry Smith #define ANSI_DECLARATORS
105a449d36SBarry Smith #endif
112e68589bSMark F. Adams #include <triangle.h>
122e68589bSMark F. Adams #endif
132e68589bSMark F. Adams 
142e68589bSMark F. Adams #include <petscblaslapack.h>
152e68589bSMark F. Adams 
16c8b0795cSMark F. Adams /* Private context for the GAMG preconditioner */
17c8b0795cSMark F. Adams typedef struct {
18c8b0795cSMark F. Adams   PetscInt lid;    /* local vertex index */
19c8b0795cSMark F. Adams   PetscInt degree; /* vertex degree */
20c8b0795cSMark F. Adams } GAMGNode;
212fa5cd67SKarl Rupp 
229371c9d4SSatish Balay static inline int petsc_geo_mg_compare(const void *a, const void *b) {
23c8b0795cSMark F. Adams   return (((GAMGNode *)a)->degree - ((GAMGNode *)b)->degree);
24c8b0795cSMark F. Adams }
25c8b0795cSMark F. Adams 
262e68589bSMark F. Adams /*
272e68589bSMark F. Adams    PCSetCoordinates_GEO
282e68589bSMark F. Adams 
292e68589bSMark F. Adams    Input Parameter:
302e68589bSMark F. Adams    .  pc - the preconditioner context
312e68589bSMark F. Adams */
329371c9d4SSatish Balay PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) {
332e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
342e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
3590fbc344SStefano Zampini   PetscInt arrsz, bs, my0, kk, ii, nloc, Iend, aloc;
362e68589bSMark F. Adams   Mat      Amat = pc->pmat;
372e68589bSMark F. Adams 
382e68589bSMark F. Adams   PetscFunctionBegin;
392e68589bSMark F. Adams   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
409566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
419566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &my0, &Iend));
4290fbc344SStefano Zampini   aloc = (Iend - my0);
432e68589bSMark F. Adams   nloc = (Iend - my0) / bs;
44a2f3521dSMark F. Adams 
452472a847SBarry Smith   PetscCheck(nloc == a_nloc || aloc == a_nloc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of local blocks %" PetscInt_FMT " must be %" PetscInt_FMT " or %" PetscInt_FMT ".", a_nloc, nloc, aloc);
462e68589bSMark F. Adams 
47c8b0795cSMark F. Adams   pc_gamg->data_cell_rows = 1;
487827d75bSBarry Smith   PetscCheck(coords || (nloc <= 0), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
49c8b0795cSMark F. Adams   pc_gamg->data_cell_cols = ndm; /* coordinates */
502e68589bSMark F. Adams 
51c8b0795cSMark F. Adams   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
522e68589bSMark F. Adams 
532e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
547f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
559566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
569566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
572e68589bSMark F. Adams   }
582e68589bSMark F. Adams   for (kk = 0; kk < arrsz; kk++) pc_gamg->data[kk] = -999.;
592e68589bSMark F. Adams   pc_gamg->data[arrsz] = -99.;
602e68589bSMark F. Adams   /* copy data in - column oriented */
6190fbc344SStefano Zampini   if (nloc == a_nloc) {
622e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) {
63ad540459SPierre Jolivet       for (ii = 0; ii < ndm; ii++) pc_gamg->data[ii * nloc + kk] = coords[kk * ndm + ii];
642e68589bSMark F. Adams     }
6590fbc344SStefano Zampini   } else { /* assumes the coordinates are blocked */
6690fbc344SStefano Zampini     for (kk = 0; kk < nloc; kk++) {
67ad540459SPierre Jolivet       for (ii = 0; ii < ndm; ii++) pc_gamg->data[ii * nloc + kk] = coords[bs * kk * ndm + ii];
6890fbc344SStefano Zampini     }
6990fbc344SStefano Zampini   }
7063a3b9bcSJacob Faibussowitsch   PetscCheck(pc_gamg->data[arrsz] == -99., PETSC_COMM_SELF, PETSC_ERR_PLIB, "pc_gamg->data[arrsz %" PetscInt_FMT "] %g != -99.", arrsz, (double)pc_gamg->data[arrsz]);
712e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
722e68589bSMark F. Adams   PetscFunctionReturn(0);
732e68589bSMark F. Adams }
742e68589bSMark F. Adams 
752e68589bSMark F. Adams /*
762e68589bSMark F. Adams    PCSetData_GEO
772e68589bSMark F. Adams 
782e68589bSMark F. Adams   Input Parameter:
792e68589bSMark F. Adams    . pc -
802e68589bSMark F. Adams */
819371c9d4SSatish Balay PetscErrorCode PCSetData_GEO(PC pc, Mat m) {
822e68589bSMark F. Adams   PetscFunctionBegin;
83ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "GEO MG needs coordinates");
842e68589bSMark F. Adams }
852e68589bSMark F. Adams 
862e68589bSMark F. Adams /*
872e68589bSMark F. Adams    PCSetFromOptions_GEO
882e68589bSMark F. Adams 
892e68589bSMark F. Adams   Input Parameter:
902e68589bSMark F. Adams    . pc -
912e68589bSMark F. Adams */
929371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_GEO(PC pc, PetscOptionItems *PetscOptionsObject) {
932e68589bSMark F. Adams   PetscFunctionBegin;
94d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-GEO options");
952e68589bSMark F. Adams   {
962e68589bSMark F. Adams     /* -pc_gamg_sa_nsmooths */
972e68589bSMark F. Adams     /* pc_gamg_sa->smooths = 0; */
982e68589bSMark F. Adams     /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */
992e68589bSMark F. Adams     /*                        "smoothing steps for smoothed aggregation, usually 1 (0)", */
1002e68589bSMark F. Adams     /*                        "PCGAMGSetNSmooths_AGG", */
1012e68589bSMark F. Adams     /*                        pc_gamg_sa->smooths, */
1022e68589bSMark F. Adams     /*                        &pc_gamg_sa->smooths, */
1032e68589bSMark F. Adams     /*                        &flag);  */
1042e68589bSMark F. Adams   }
105d0609cedSBarry Smith   PetscOptionsHeadEnd();
1062e68589bSMark F. Adams   PetscFunctionReturn(0);
1072e68589bSMark F. Adams }
1082e68589bSMark F. Adams 
1092e68589bSMark F. Adams /*
1102e68589bSMark F. Adams  triangulateAndFormProl
1112e68589bSMark F. Adams 
1122e68589bSMark F. Adams    Input Parameter:
1132e68589bSMark F. Adams    . selected_2 - list of selected local ID, includes selected ghosts
114a2f3521dSMark F. Adams    . data_stride -
115a2f3521dSMark F. Adams    . coords[2*data_stride] - column vector of local coordinates w/ ghosts
1162adfe9d3SBarry Smith    . nselected_1 - selected IDs that go with base (1) graph includes selected ghosts
1170cbbd2e1SMark F. Adams    . clid_lid_1[nselected_1] - lids of selected (c) nodes   ???????????
1182adfe9d3SBarry Smith    . agg_lists_1 - list of aggregates selected_1 vertices of aggregate unselected vertices
1192e68589bSMark F. Adams    . crsGID[selected.size()] - global index for prolongation operator
1202e68589bSMark F. Adams    . bs - block size
1212e68589bSMark F. Adams   Output Parameter:
1222e68589bSMark F. Adams    . a_Prol - prolongation operator
1232e68589bSMark F. Adams    . a_worst_best - measure of worst missed fine vertex, 0 is no misses
1242e68589bSMark F. Adams */
1259371c9d4SSatish Balay static PetscErrorCode triangulateAndFormProl(IS selected_2, PetscInt data_stride, PetscReal coords[], PetscInt nselected_1, const PetscInt clid_lid_1[], const PetscCoarsenData *agg_lists_1, const PetscInt crsGID[], PetscInt bs, Mat a_Prol, PetscReal *a_worst_best) {
1262e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE)
127b43b03e9SMark F. Adams   PetscInt             jj, tid, tt, idx, nselected_2;
1282e68589bSMark F. Adams   struct triangulateio in, mid;
1290cbbd2e1SMark F. Adams   const PetscInt      *selected_idx_2;
13073911c69SBarry Smith   PetscMPIInt          rank;
1312e68589bSMark F. Adams   PetscInt             Istart, Iend, nFineLoc, myFine0;
1322e68589bSMark F. Adams   int                  kk, nPlotPts, sid;
1333b4367a7SBarry Smith   MPI_Comm             comm;
134c8b0795cSMark F. Adams   PetscReal            tm;
135c8b0795cSMark F. Adams 
1366e111a19SKarl Rupp   PetscFunctionBegin;
1379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
1389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1399566063dSJacob Faibussowitsch   PetscCall(ISGetSize(selected_2, &nselected_2));
1402e68589bSMark F. Adams   if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */
1412e68589bSMark F. Adams     *a_worst_best = 100.0;                    /* this will cause a stop, but not globalized (should not happen) */
142806fa848SBarry Smith   } else *a_worst_best = 0.0;
1431c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm));
144c8b0795cSMark F. Adams   if (tm > 0.0) {
145c8b0795cSMark F. Adams     *a_worst_best = 100.0;
1462e68589bSMark F. Adams     PetscFunctionReturn(0);
1472e68589bSMark F. Adams   }
1489566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
1499371c9d4SSatish Balay   nFineLoc                   = (Iend - Istart) / bs;
1509371c9d4SSatish Balay   myFine0                    = Istart / bs;
1512e68589bSMark F. Adams   nPlotPts                   = nFineLoc; /* locals */
1524c500f23SPierre Jolivet   /* triangle */
1532e68589bSMark F. Adams   /* Define input points - in*/
1542e68589bSMark F. Adams   in.numberofpoints          = nselected_2;
1552e68589bSMark F. Adams   in.numberofpointattributes = 0;
1562e68589bSMark F. Adams   /* get nselected points */
1579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * nselected_2, &in.pointlist));
1589566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(selected_2, &selected_idx_2));
1592e68589bSMark F. Adams 
1602e68589bSMark F. Adams   for (kk = 0, sid = 0; kk < nselected_2; kk++, sid += 2) {
1612e68589bSMark F. Adams     PetscInt lid          = selected_idx_2[kk];
1622e68589bSMark F. Adams     in.pointlist[sid]     = coords[lid];
163a2f3521dSMark F. Adams     in.pointlist[sid + 1] = coords[data_stride + lid];
1642e68589bSMark F. Adams     if (lid >= nFineLoc) nPlotPts++;
1652e68589bSMark F. Adams   }
166c05e6a5bSJed Brown   PetscCheck(sid == 2 * nselected_2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "sid %d != 2*nselected_2 %" PetscInt_FMT, sid, nselected_2);
1672e68589bSMark F. Adams 
1682e68589bSMark F. Adams   in.numberofsegments      = 0;
1692e68589bSMark F. Adams   in.numberofedges         = 0;
1702e68589bSMark F. Adams   in.numberofholes         = 0;
1712e68589bSMark F. Adams   in.numberofregions       = 0;
1720a545947SLisandro Dalcin   in.trianglelist          = NULL;
1730a545947SLisandro Dalcin   in.segmentmarkerlist     = NULL;
1740a545947SLisandro Dalcin   in.pointattributelist    = NULL;
1750a545947SLisandro Dalcin   in.pointmarkerlist       = NULL;
1760a545947SLisandro Dalcin   in.triangleattributelist = NULL;
1770a545947SLisandro Dalcin   in.trianglearealist      = NULL;
1780a545947SLisandro Dalcin   in.segmentlist           = NULL;
1790a545947SLisandro Dalcin   in.holelist              = NULL;
1800a545947SLisandro Dalcin   in.regionlist            = NULL;
1810a545947SLisandro Dalcin   in.edgelist              = NULL;
1820a545947SLisandro Dalcin   in.edgemarkerlist        = NULL;
1830a545947SLisandro Dalcin   in.normlist              = NULL;
1842fa5cd67SKarl Rupp 
1852e68589bSMark F. Adams   /* triangulate */
1860a545947SLisandro Dalcin   mid.pointlist             = NULL; /* Not needed if -N switch used. */
1872e68589bSMark F. Adams   /* Not needed if -N switch used or number of point attributes is zero: */
1880a545947SLisandro Dalcin   mid.pointattributelist    = NULL;
1890a545947SLisandro Dalcin   mid.pointmarkerlist       = NULL; /* Not needed if -N or -B switch used. */
1900a545947SLisandro Dalcin   mid.trianglelist          = NULL; /* Not needed if -E switch used. */
1912e68589bSMark F. Adams   /* Not needed if -E switch used or number of triangle attributes is zero: */
1920a545947SLisandro Dalcin   mid.triangleattributelist = NULL;
1930a545947SLisandro Dalcin   mid.neighborlist          = NULL; /* Needed only if -n switch used. */
1942e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P not used: */
1950a545947SLisandro Dalcin   mid.segmentlist           = NULL;
1962e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
1970a545947SLisandro Dalcin   mid.segmentmarkerlist     = NULL;
1980a545947SLisandro Dalcin   mid.edgelist              = NULL; /* Needed only if -e switch used. */
1990a545947SLisandro Dalcin   mid.edgemarkerlist        = NULL; /* Needed if -e used and -B not used. */
2002e68589bSMark F. Adams   mid.numberoftriangles     = 0;
2012e68589bSMark F. Adams 
2022e68589bSMark F. Adams   /* Triangulate the points.  Switches are chosen to read and write a  */
2032e68589bSMark F. Adams   /*   PSLG (p), preserve the convex hull (c), number everything from  */
2042e68589bSMark F. Adams   /*   zero (z), assign a regional attribute to each element (A), and  */
2052e68589bSMark F. Adams   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
2062e68589bSMark F. Adams   /*   neighbor list (n).                                            */
2072e68589bSMark F. Adams   if (nselected_2 != 0) {  /* inactive processor */
2082e68589bSMark F. Adams     char args[] = "npczQ"; /* c is needed ? */
2092e68589bSMark F. Adams     triangulate(args, &in, &mid, (struct triangulateio *)NULL);
2102e68589bSMark F. Adams     /* output .poly files for 'showme' */
2112e68589bSMark F. Adams     if (!PETSC_TRUE) {
2122e68589bSMark F. Adams       static int level = 1;
2139371c9d4SSatish Balay       FILE      *file;
2149371c9d4SSatish Balay       char       fname[32];
2152e68589bSMark F. Adams 
2169371c9d4SSatish Balay       sprintf(fname, "C%d_%d.poly", level, rank);
2179371c9d4SSatish Balay       file = fopen(fname, "w");
2182e68589bSMark F. Adams       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
2192e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n", in.numberofpoints, 2, 0, 0);
2202e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
221ad540459SPierre Jolivet       for (kk = 0, sid = 0; kk < in.numberofpoints; kk++, sid += 2) fprintf(file, "%d %e %e\n", kk, in.pointlist[sid], in.pointlist[sid + 1]);
2222e68589bSMark F. Adams       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
2232e68589bSMark F. Adams       fprintf(file, "%d  %d\n", 0, 0);
2242e68589bSMark F. Adams       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
2252e68589bSMark F. Adams       /* One line: <# of holes> */
2262e68589bSMark F. Adams       fprintf(file, "%d\n", 0);
2272e68589bSMark F. Adams       /* Following lines: <hole #> <x> <y> */
2282e68589bSMark F. Adams       /* Optional line: <# of regional attributes and/or area constraints> */
2292e68589bSMark F. Adams       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
2302e68589bSMark F. Adams       fclose(file);
2312e68589bSMark F. Adams 
2322e68589bSMark F. Adams       /* elems */
2339371c9d4SSatish Balay       sprintf(fname, "C%d_%d.ele", level, rank);
2349371c9d4SSatish Balay       file = fopen(fname, "w");
2352e68589bSMark F. Adams       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
2362e68589bSMark F. Adams       fprintf(file, "%d %d %d\n", mid.numberoftriangles, 3, 0);
2372e68589bSMark F. Adams       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
238ad540459SPierre Jolivet       for (kk = 0, sid = 0; kk < mid.numberoftriangles; kk++, sid += 3) fprintf(file, "%d %d %d %d\n", kk, mid.trianglelist[sid], mid.trianglelist[sid + 1], mid.trianglelist[sid + 2]);
2392e68589bSMark F. Adams       fclose(file);
2402e68589bSMark F. Adams 
2419371c9d4SSatish Balay       sprintf(fname, "C%d_%d.node", level, rank);
2429371c9d4SSatish Balay       file = fopen(fname, "w");
2432e68589bSMark F. Adams       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
2442e68589bSMark F. Adams       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
2452e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n", nPlotPts, 2, 0, 0);
2462e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
247ad540459SPierre Jolivet       for (kk = 0, sid = 0; kk < in.numberofpoints; kk++, sid += 2) fprintf(file, "%d %e %e\n", kk, in.pointlist[sid], in.pointlist[sid + 1]);
2482e68589bSMark F. Adams 
2492e68589bSMark F. Adams       sid /= 2;
2502e68589bSMark F. Adams       for (jj = 0; jj < nFineLoc; jj++) {
2512e68589bSMark F. Adams         PetscBool sel = PETSC_TRUE;
2522e68589bSMark F. Adams         for (kk = 0; kk < nselected_2 && sel; kk++) {
2532e68589bSMark F. Adams           PetscInt lid = selected_idx_2[kk];
2542e68589bSMark F. Adams           if (lid == jj) sel = PETSC_FALSE;
2552e68589bSMark F. Adams         }
2562fa5cd67SKarl Rupp         if (sel) fprintf(file, "%d %e %e\n", sid++, coords[jj], coords[data_stride + jj]);
2572e68589bSMark F. Adams       }
2582e68589bSMark F. Adams       fclose(file);
259c05e6a5bSJed Brown       PetscCheck(sid == nPlotPts, PETSC_COMM_SELF, PETSC_ERR_PLIB, "sid %d != nPlotPts %d", sid, nPlotPts);
2602e68589bSMark F. Adams       level++;
2612e68589bSMark F. Adams     }
2622e68589bSMark F. Adams   }
2632e68589bSMark F. Adams   { /* form P - setup some maps */
2640cbbd2e1SMark F. Adams     PetscInt clid, mm, *nTri, *node_tri;
2652e68589bSMark F. Adams 
2669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nselected_2, &node_tri, nselected_2, &nTri));
2672e68589bSMark F. Adams 
2682e68589bSMark F. Adams     /* need list of triangles on node */
2692e68589bSMark F. Adams     for (kk = 0; kk < nselected_2; kk++) nTri[kk] = 0;
2702e68589bSMark F. Adams     for (tid = 0, kk = 0; tid < mid.numberoftriangles; tid++) {
2712e68589bSMark F. Adams       for (jj = 0; jj < 3; jj++) {
2722e68589bSMark F. Adams         PetscInt cid = mid.trianglelist[kk++];
2732e68589bSMark F. Adams         if (nTri[cid] == 0) node_tri[cid] = tid;
2742e68589bSMark F. Adams         nTri[cid]++;
2752e68589bSMark F. Adams       }
2762e68589bSMark F. Adams     }
2772e68589bSMark F. Adams #define EPS 1.e-12
2782e68589bSMark F. Adams     /* find points and set prolongation */
2790cbbd2e1SMark F. Adams     for (mm = clid = 0; mm < nFineLoc; mm++) {
280e78576d6SMark F. Adams       PetscBool ise;
2819566063dSJacob Faibussowitsch       PetscCall(PetscCDEmptyAt(agg_lists_1, mm, &ise));
282e78576d6SMark F. Adams       if (!ise) {
2830cbbd2e1SMark F. Adams         const PetscInt lid = mm;
2842e68589bSMark F. Adams         PetscScalar    AA[3][3];
2852e68589bSMark F. Adams         PetscBLASInt   N = 3, NRHS = 1, LDA = 3, IPIV[3], LDB = 3, INFO;
286539c167fSBarry Smith         PetscCDIntNd  *pos;
287539c167fSBarry Smith 
2889566063dSJacob Faibussowitsch         PetscCall(PetscCDGetHeadPos(agg_lists_1, lid, &pos));
289e78576d6SMark F. Adams         while (pos) {
290e78576d6SMark F. Adams           PetscInt flid;
2919566063dSJacob Faibussowitsch           PetscCall(PetscCDIntNdGetID(pos, &flid));
2929566063dSJacob Faibussowitsch           PetscCall(PetscCDGetNextPos(agg_lists_1, lid, &pos));
2930cbbd2e1SMark F. Adams 
2942e68589bSMark F. Adams           if (flid < nFineLoc) { /* could be a ghost */
2959371c9d4SSatish Balay             PetscInt        bestTID    = -1;
2969371c9d4SSatish Balay             PetscReal       best_alpha = 1.e10;
2972e68589bSMark F. Adams             const PetscInt  fgid       = flid + myFine0;
2982e68589bSMark F. Adams             /* compute shape function for gid */
299a2f3521dSMark F. Adams             const PetscReal fcoord[3]  = {coords[flid], coords[data_stride + flid], 1.0};
3009371c9d4SSatish Balay             PetscBool       haveit     = PETSC_FALSE;
3019371c9d4SSatish Balay             PetscScalar     alpha[3];
3029371c9d4SSatish Balay             PetscInt        clids[3];
3032fa5cd67SKarl Rupp 
3042e68589bSMark F. Adams             /* look for it */
3059371c9d4SSatish Balay             for (tid = node_tri[clid], jj = 0; jj < 5 && !haveit && tid != -1; jj++) {
3062e68589bSMark F. Adams               for (tt = 0; tt < 3; tt++) {
3072e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3 * tid + tt];
3082e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
3099371c9d4SSatish Balay                 AA[tt][0]     = coords[lid2];
3109371c9d4SSatish Balay                 AA[tt][1]     = coords[data_stride + lid2];
3119371c9d4SSatish Balay                 AA[tt][2]     = 1.0;
3122e68589bSMark F. Adams                 clids[tt]     = cid2; /* store for interp */
3132e68589bSMark F. Adams               }
3142e68589bSMark F. Adams 
3152e68589bSMark F. Adams               for (tt = 0; tt < 3; tt++) alpha[tt] = (PetscScalar)fcoord[tt];
3162e68589bSMark F. Adams 
3172e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
318792fecdfSBarry Smith               PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&N, &NRHS, (PetscScalar *)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3192e68589bSMark F. Adams               {
3209371c9d4SSatish Balay                 PetscBool have   = PETSC_TRUE;
3219371c9d4SSatish Balay                 PetscReal lowest = 1.e10;
3222e68589bSMark F. Adams                 for (tt = 0, idx = 0; tt < 3; tt++) {
3232e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) > (1.0 + EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
3242e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) < lowest) {
3252e68589bSMark F. Adams                     lowest = PetscRealPart(alpha[tt]);
3262e68589bSMark F. Adams                     idx    = tt;
3272e68589bSMark F. Adams                   }
3282e68589bSMark F. Adams                 }
3292e68589bSMark F. Adams                 haveit = have;
3302e68589bSMark F. Adams               }
3312e68589bSMark F. Adams               tid = mid.neighborlist[3 * tid + idx];
3322e68589bSMark F. Adams             }
3332e68589bSMark F. Adams 
3342e68589bSMark F. Adams             if (!haveit) {
3352e68589bSMark F. Adams               /* brute force */
3362e68589bSMark F. Adams               for (tid = 0; tid < mid.numberoftriangles && !haveit; tid++) {
3372e68589bSMark F. Adams                 for (tt = 0; tt < 3; tt++) {
3382e68589bSMark F. Adams                   PetscInt cid2 = mid.trianglelist[3 * tid + tt];
3392e68589bSMark F. Adams                   PetscInt lid2 = selected_idx_2[cid2];
3409371c9d4SSatish Balay                   AA[tt][0]     = coords[lid2];
3419371c9d4SSatish Balay                   AA[tt][1]     = coords[data_stride + lid2];
3429371c9d4SSatish Balay                   AA[tt][2]     = 1.0;
3432e68589bSMark F. Adams                   clids[tt]     = cid2; /* store for interp */
3442e68589bSMark F. Adams                 }
3452e68589bSMark F. Adams                 for (tt = 0; tt < 3; tt++) alpha[tt] = fcoord[tt];
3462e68589bSMark F. Adams                 /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
347792fecdfSBarry Smith                 PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&N, &NRHS, (PetscScalar *)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3482e68589bSMark F. Adams                 {
3499371c9d4SSatish Balay                   PetscBool have  = PETSC_TRUE;
3509371c9d4SSatish Balay                   PetscReal worst = 0.0, v;
3512e68589bSMark F. Adams                   for (tt = 0; tt < 3 && have; tt++) {
3522e68589bSMark F. Adams                     if (PetscRealPart(alpha[tt]) > 1.0 + EPS || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
3532e68589bSMark F. Adams                     if ((v = PetscAbs(PetscRealPart(alpha[tt]) - 0.5)) > worst) worst = v;
3542e68589bSMark F. Adams                   }
3552e68589bSMark F. Adams                   if (worst < best_alpha) {
3569371c9d4SSatish Balay                     best_alpha = worst;
3579371c9d4SSatish Balay                     bestTID    = tid;
3582e68589bSMark F. Adams                   }
3592e68589bSMark F. Adams                   haveit = have;
3602e68589bSMark F. Adams                 }
3612e68589bSMark F. Adams               }
3622e68589bSMark F. Adams             }
3632e68589bSMark F. Adams             if (!haveit) {
3642e68589bSMark F. Adams               if (best_alpha > *a_worst_best) *a_worst_best = best_alpha;
3652e68589bSMark F. Adams               /* use best one */
3662e68589bSMark F. Adams               for (tt = 0; tt < 3; tt++) {
3672e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3 * bestTID + tt];
3682e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
3699371c9d4SSatish Balay                 AA[tt][0]     = coords[lid2];
3709371c9d4SSatish Balay                 AA[tt][1]     = coords[data_stride + lid2];
3719371c9d4SSatish Balay                 AA[tt][2]     = 1.0;
3722e68589bSMark F. Adams                 clids[tt]     = cid2; /* store for interp */
3732e68589bSMark F. Adams               }
3742e68589bSMark F. Adams               for (tt = 0; tt < 3; tt++) alpha[tt] = fcoord[tt];
3752e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
376792fecdfSBarry Smith               PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&N, &NRHS, (PetscScalar *)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3772e68589bSMark F. Adams             }
3782e68589bSMark F. Adams 
3792e68589bSMark F. Adams             /* put in row of P */
3802e68589bSMark F. Adams             for (idx = 0; idx < 3; idx++) {
3812e68589bSMark F. Adams               PetscScalar shp = alpha[idx];
3822e68589bSMark F. Adams               if (PetscAbs(PetscRealPart(shp)) > 1.e-6) {
3832e68589bSMark F. Adams                 PetscInt cgid = crsGID[clids[idx]];
3842e68589bSMark F. Adams                 PetscInt jj = cgid * bs, ii = fgid * bs; /* need to gloalize */
38548a46eb9SPierre Jolivet                 for (tt = 0; tt < bs; tt++, ii++, jj++) PetscCall(MatSetValues(a_Prol, 1, &ii, 1, &jj, &shp, INSERT_VALUES));
3862e68589bSMark F. Adams               }
3872e68589bSMark F. Adams             }
3882e68589bSMark F. Adams           }
3890cbbd2e1SMark F. Adams         } /* aggregates iterations */
3900cbbd2e1SMark F. Adams         clid++;
3910cbbd2e1SMark F. Adams       } /* a coarse agg */
3920cbbd2e1SMark F. Adams     }   /* for all fine nodes */
3930cbbd2e1SMark F. Adams 
3949566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(selected_2, &selected_idx_2));
3959566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
3969566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
3972e68589bSMark F. Adams 
3989566063dSJacob Faibussowitsch     PetscCall(PetscFree2(node_tri, nTri));
3992e68589bSMark F. Adams   }
4002e68589bSMark F. Adams   free(mid.trianglelist);
4012e68589bSMark F. Adams   free(mid.neighborlist);
402c6b50ea1SMark   free(mid.segmentlist);
403c6b50ea1SMark   free(mid.segmentmarkerlist);
404c6b50ea1SMark   free(mid.pointlist);
405c6b50ea1SMark   free(mid.pointmarkerlist);
4069566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointlist));
4072e68589bSMark F. Adams   PetscFunctionReturn(0);
4082e68589bSMark F. Adams #else
4093b4367a7SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)a_Prol), PETSC_ERR_PLIB, "configure with TRIANGLE to use geometric MG");
4102e68589bSMark F. Adams #endif
4112e68589bSMark F. Adams }
412*f1580f4eSBarry Smith 
4132e68589bSMark F. Adams /*
4142e68589bSMark F. Adams    getGIDsOnSquareGraph - square graph, get
4152e68589bSMark F. Adams 
4162e68589bSMark F. Adams    Input Parameter:
4170cbbd2e1SMark F. Adams    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
418b43b03e9SMark F. Adams    . clid_lid_1 - [nselected_1] lids of selected nodes
4192e68589bSMark F. Adams    . Gmat1 - graph that goes with 'selected_1'
4202e68589bSMark F. Adams    Output Parameter:
4212e68589bSMark F. Adams    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
4222e68589bSMark F. Adams    . a_Gmat_2 - graph that is squared of 'Gmat_1'
4232e68589bSMark F. Adams    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
4242e68589bSMark F. Adams */
4259371c9d4SSatish Balay static PetscErrorCode getGIDsOnSquareGraph(PC pc, PetscInt nselected_1, const PetscInt clid_lid_1[], const Mat Gmat1, IS *a_selected_2, Mat *a_Gmat_2, PetscInt **a_crsGID) {
42673911c69SBarry Smith   PetscMPIInt size;
427b43b03e9SMark F. Adams   PetscInt   *crsGID, kk, my0, Iend, nloc;
4283b4367a7SBarry Smith   MPI_Comm    comm;
4292e68589bSMark F. Adams 
4302e68589bSMark F. Adams   PetscFunctionBegin;
4319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm));
4329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
4339566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &my0, &Iend)); /* AIJ */
4342e68589bSMark F. Adams   nloc = Iend - my0;                                   /* this does not change */
4352e68589bSMark F. Adams 
436c5df96a5SBarry Smith   if (size == 1) { /* not much to do in serial */
4379566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected_1, &crsGID));
438b43b03e9SMark F. Adams     for (kk = 0; kk < nselected_1; kk++) crsGID[kk] = kk;
4390a545947SLisandro Dalcin     *a_Gmat_2 = NULL;
4409566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nselected_1, clid_lid_1, PETSC_COPY_VALUES, a_selected_2));
441806fa848SBarry Smith   } else {
442b43b03e9SMark F. Adams     PetscInt     idx, num_fine_ghosts, num_crs_ghost, myCrs0;
4432e68589bSMark F. Adams     Mat_MPIAIJ  *mpimat2;
4442e68589bSMark F. Adams     Mat          Gmat2;
4452e68589bSMark F. Adams     Vec          locState;
4462e68589bSMark F. Adams     PetscScalar *cpcol_state;
4472e68589bSMark F. Adams 
4482e68589bSMark F. Adams     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
449b43b03e9SMark F. Adams     kk = nselected_1;
4509566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPI_SUM, comm));
451b43b03e9SMark F. Adams     myCrs0 -= nselected_1;
4522e68589bSMark F. Adams 
453b43b03e9SMark F. Adams     if (a_Gmat_2) { /* output */
4542e68589bSMark F. Adams       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
4559566063dSJacob Faibussowitsch       PetscCall(PCGAMGSquareGraph_GAMG(pc, Gmat1, &Gmat2));
4562e68589bSMark F. Adams       *a_Gmat_2 = Gmat2;  /* output */
457806fa848SBarry Smith     } else Gmat2 = Gmat1; /* use local to get crsGIDs at least */
4582e68589bSMark F. Adams     /* get coarse grid GIDS for selected (locals and ghosts) */
4592e68589bSMark F. Adams     mpimat2 = (Mat_MPIAIJ *)Gmat2->data;
4609566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Gmat2, &locState, NULL));
4619566063dSJacob Faibussowitsch     PetscCall(VecSet(locState, (PetscScalar)(PetscReal)(-1))); /* set with UNKNOWN state */
462b43b03e9SMark F. Adams     for (kk = 0; kk < nselected_1; kk++) {
463b43b03e9SMark F. Adams       PetscInt    fgid = clid_lid_1[kk] + my0;
4642e68589bSMark F. Adams       PetscScalar v    = (PetscScalar)(kk + myCrs0);
4659566063dSJacob Faibussowitsch       PetscCall(VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES)); /* set with PID */
4662e68589bSMark F. Adams     }
4679566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(locState));
4689566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(locState));
4699566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mpimat2->Mvctx, locState, mpimat2->lvec, INSERT_VALUES, SCATTER_FORWARD));
4709566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mpimat2->Mvctx, locState, mpimat2->lvec, INSERT_VALUES, SCATTER_FORWARD));
4719566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts));
4729566063dSJacob Faibussowitsch     PetscCall(VecGetArray(mpimat2->lvec, &cpcol_state));
4732e68589bSMark F. Adams     for (kk = 0, num_crs_ghost = 0; kk < num_fine_ghosts; kk++) {
4742e68589bSMark F. Adams       if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++;
4752e68589bSMark F. Adams     }
4769566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected_1 + num_crs_ghost, &crsGID)); /* output */
4772e68589bSMark F. Adams     {
4782e68589bSMark F. Adams       PetscInt *selected_set;
4799566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nselected_1 + num_crs_ghost, &selected_set));
4802e68589bSMark F. Adams       /* do ghost of 'crsGID' */
481b43b03e9SMark F. Adams       for (kk = 0, idx = nselected_1; kk < num_fine_ghosts; kk++) {
4822e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
4832e68589bSMark F. Adams           PetscInt cgid     = (PetscInt)PetscRealPart(cpcol_state[kk]);
4842e68589bSMark F. Adams           selected_set[idx] = nloc + kk;
4852e68589bSMark F. Adams           crsGID[idx++]     = cgid;
4862e68589bSMark F. Adams         }
4872e68589bSMark F. Adams       }
48863a3b9bcSJacob Faibussowitsch       PetscCheck(idx == (nselected_1 + num_crs_ghost), PETSC_COMM_SELF, PETSC_ERR_PLIB, "idx %" PetscInt_FMT " != (nselected_1 %" PetscInt_FMT " + num_crs_ghost %" PetscInt_FMT ")", idx, nselected_1, num_crs_ghost);
4899566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(mpimat2->lvec, &cpcol_state));
4902e68589bSMark F. Adams       /* do locals in 'crsGID' */
4919566063dSJacob Faibussowitsch       PetscCall(VecGetArray(locState, &cpcol_state));
4922e68589bSMark F. Adams       for (kk = 0, idx = 0; kk < nloc; kk++) {
4932e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
4942e68589bSMark F. Adams           PetscInt cgid     = (PetscInt)PetscRealPart(cpcol_state[kk]);
4952e68589bSMark F. Adams           selected_set[idx] = kk;
4962e68589bSMark F. Adams           crsGID[idx++]     = cgid;
4972e68589bSMark F. Adams         }
4982e68589bSMark F. Adams       }
49963a3b9bcSJacob Faibussowitsch       PetscCheck(idx == nselected_1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "idx %" PetscInt_FMT " != nselected_1 %" PetscInt_FMT, idx, nselected_1);
5009566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(locState, &cpcol_state));
5012e68589bSMark F. Adams 
5020a545947SLisandro Dalcin       if (a_selected_2 != NULL) { /* output */
5039566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, (nselected_1 + num_crs_ghost), selected_set, PETSC_OWN_POINTER, a_selected_2));
504806fa848SBarry Smith       } else {
5059566063dSJacob Faibussowitsch         PetscCall(PetscFree(selected_set));
5062e68589bSMark F. Adams       }
5070cbbd2e1SMark F. Adams     }
5089566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&locState));
5092e68589bSMark F. Adams   }
5102e68589bSMark F. Adams   *a_crsGID = crsGID; /* output */
5112e68589bSMark F. Adams   PetscFunctionReturn(0);
5122e68589bSMark F. Adams }
5132e68589bSMark F. Adams 
5142e68589bSMark F. Adams /*
515fd1112cbSBarry Smith    PCGAMGGraph_GEO
5162e68589bSMark F. Adams 
5172e68589bSMark F. Adams   Input Parameter:
5182e68589bSMark F. Adams    . pc - this
5192e68589bSMark F. Adams    . Amat - matrix on this fine level
5202e68589bSMark F. Adams   Output Parameter:
521c8b0795cSMark F. Adams    . a_Gmat
5222e68589bSMark F. Adams */
5239371c9d4SSatish Balay PetscErrorCode PCGAMGGraph_GEO(PC pc, Mat Amat, Mat *a_Gmat) {
524c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG *)pc->data;
525c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG *)mg->innerctx;
526c1eae691SMark Adams   const PetscReal vfilter = pc_gamg->threshold[0];
5273b4367a7SBarry Smith   MPI_Comm        comm;
52872833a62Smarkadams4   Mat             Gmat, F = NULL;
5290cbbd2e1SMark F. Adams   PetscBool       set, flg, symm;
5306e111a19SKarl Rupp 
531c8b0795cSMark F. Adams   PetscFunctionBegin;
5329566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
533c8b0795cSMark F. Adams 
5349566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetricKnown(Amat, &set, &flg));
535263489e9SJed Brown   symm = (PetscBool) !(set && flg);
5360cbbd2e1SMark F. Adams 
53772833a62Smarkadams4   PetscCall(MatCreateGraph(Amat, symm, PETSC_TRUE, &Gmat));
53872833a62Smarkadams4   PetscCall(MatFilter(Gmat, vfilter, &F));
53972833a62Smarkadams4   if (F) {
54072833a62Smarkadams4     PetscCall(MatDestroy(&Gmat));
54172833a62Smarkadams4     Gmat = F;
54272833a62Smarkadams4   }
543c8b0795cSMark F. Adams   *a_Gmat = Gmat;
544849bee69SMark Adams 
545c8b0795cSMark F. Adams   PetscFunctionReturn(0);
546c8b0795cSMark F. Adams }
547c8b0795cSMark F. Adams 
548c8b0795cSMark F. Adams /*
549fd1112cbSBarry Smith    PCGAMGCoarsen_GEO
550c8b0795cSMark F. Adams 
551c8b0795cSMark F. Adams   Input Parameter:
552e0940f08SMark F. Adams    . a_pc - this
553e0940f08SMark F. Adams    . a_Gmat - graph
554c8b0795cSMark F. Adams   Output Parameter:
555c8b0795cSMark F. Adams    . a_llist_parent - linked list from selected indices for data locality only
556c8b0795cSMark F. Adams */
5579371c9d4SSatish Balay PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc, Mat *a_Gmat, PetscCoarsenData **a_llist_parent) {
558c8b0795cSMark F. Adams   PetscInt   Istart, Iend, nloc, kk, Ii, ncols;
5590cbbd2e1SMark F. Adams   IS         perm;
560c8b0795cSMark F. Adams   GAMGNode  *gnodes;
561c8b0795cSMark F. Adams   PetscInt  *permute;
562e0940f08SMark F. Adams   Mat        Gmat = *a_Gmat;
5633b4367a7SBarry Smith   MPI_Comm   comm;
564b43b03e9SMark F. Adams   MatCoarsen crs;
565c8b0795cSMark F. Adams 
566c8b0795cSMark F. Adams   PetscFunctionBegin;
5679566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_pc, &comm));
568849bee69SMark Adams 
5699566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend));
570c8b0795cSMark F. Adams   nloc = (Iend - Istart);
571c8b0795cSMark F. Adams 
572c8b0795cSMark F. Adams   /* create random permutation with sort for geo-mg */
5739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc, &gnodes));
5749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc, &permute));
575c8b0795cSMark F. Adams 
576c8b0795cSMark F. Adams   for (Ii = Istart; Ii < Iend; Ii++) { /* locals only? */
5779566063dSJacob Faibussowitsch     PetscCall(MatGetRow(Gmat, Ii, &ncols, NULL, NULL));
578c8b0795cSMark F. Adams     {
579c8b0795cSMark F. Adams       PetscInt lid       = Ii - Istart;
580c8b0795cSMark F. Adams       gnodes[lid].lid    = lid;
581c8b0795cSMark F. Adams       gnodes[lid].degree = ncols;
582c8b0795cSMark F. Adams     }
5839566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(Gmat, Ii, &ncols, NULL, NULL));
584c8b0795cSMark F. Adams   }
585c8b0795cSMark F. Adams   if (PETSC_TRUE) {
586e2c00dcbSBarry Smith     PetscRandom rand;
587c8b0795cSMark F. Adams     PetscBool  *bIndexSet;
588e2c00dcbSBarry Smith     PetscReal   rr;
589e2c00dcbSBarry Smith     PetscInt    iSwapIndex;
590e2c00dcbSBarry Smith 
5919566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
5929566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nloc, &bIndexSet));
5932fa5cd67SKarl Rupp     for (Ii = 0; Ii < nloc; Ii++) {
5949566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(rand, &rr));
595e2c00dcbSBarry Smith       iSwapIndex = (PetscInt)(rr * nloc);
5962fa5cd67SKarl Rupp       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
597c8b0795cSMark F. Adams         GAMGNode iTemp        = gnodes[iSwapIndex];
598c8b0795cSMark F. Adams         gnodes[iSwapIndex]    = gnodes[Ii];
599c8b0795cSMark F. Adams         gnodes[Ii]            = iTemp;
600c8b0795cSMark F. Adams         bIndexSet[Ii]         = PETSC_TRUE;
601c8b0795cSMark F. Adams         bIndexSet[iSwapIndex] = PETSC_TRUE;
602c8b0795cSMark F. Adams       }
603c8b0795cSMark F. Adams     }
6049566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
6059566063dSJacob Faibussowitsch     PetscCall(PetscFree(bIndexSet));
606c8b0795cSMark F. Adams   }
607c8b0795cSMark F. Adams   /* only sort locals */
6080cbbd2e1SMark F. Adams   qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare);
609c8b0795cSMark F. Adams   /* create IS of permutation */
6102fa5cd67SKarl Rupp   for (kk = 0; kk < nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */
6119566063dSJacob Faibussowitsch   PetscCall(PetscFree(gnodes));
6129566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm));
613c8b0795cSMark F. Adams 
614c8b0795cSMark F. Adams   /* get MIS aggs */
615b43b03e9SMark F. Adams 
6169566063dSJacob Faibussowitsch   PetscCall(MatCoarsenCreate(comm, &crs));
6179566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetType(crs, MATCOARSENMIS));
6189566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetGreedyOrdering(crs, perm));
6199566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetAdjacency(crs, Gmat));
6209566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_FALSE));
6219566063dSJacob Faibussowitsch   PetscCall(MatCoarsenApply(crs));
6229566063dSJacob Faibussowitsch   PetscCall(MatCoarsenGetData(crs, a_llist_parent));
6239566063dSJacob Faibussowitsch   PetscCall(MatCoarsenDestroy(&crs));
624c8b0795cSMark F. Adams 
6259566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
626849bee69SMark Adams 
627c8b0795cSMark F. Adams   PetscFunctionReturn(0);
628c8b0795cSMark F. Adams }
629c8b0795cSMark F. Adams 
630c8b0795cSMark F. Adams /*
6310cbbd2e1SMark F. Adams  PCGAMGProlongator_GEO
632c8b0795cSMark F. Adams 
633c8b0795cSMark F. Adams  Input Parameter:
634c8b0795cSMark F. Adams  . pc - this
635c8b0795cSMark F. Adams  . Amat - matrix on this fine level
636c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
6370cbbd2e1SMark F. Adams  . selected_1 - [nselected]
6380cbbd2e1SMark F. Adams  . agg_lists - [nselected]
639c8b0795cSMark F. Adams  Output Parameter:
640c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
641c8b0795cSMark F. Adams  */
6429371c9d4SSatish Balay PetscErrorCode PCGAMGProlongator_GEO(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out) {
6432e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG *)pc->data;
6442e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG *)mg->innerctx;
645c8b0795cSMark F. Adams   const PetscInt  dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
646b43b03e9SMark F. Adams   PetscInt        Istart, Iend, nloc, my0, jj, kk, ncols, nLocalSelected, bs, *clid_flid;
647c8b0795cSMark F. Adams   Mat             Prol;
648c5df96a5SBarry Smith   PetscMPIInt     rank, size;
6493b4367a7SBarry Smith   MPI_Comm        comm;
6500cbbd2e1SMark F. Adams   IS              selected_2, selected_1;
6512e68589bSMark F. Adams   const PetscInt *selected_idx;
652d9558ea9SBarry Smith   MatType         mtype;
6532e68589bSMark F. Adams 
6542e68589bSMark F. Adams   PetscFunctionBegin;
6559566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
656849bee69SMark Adams 
6579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
6589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6599566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
6609566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
6619371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
6629371c9d4SSatish Balay   my0  = Istart / bs;
66363a3b9bcSJacob Faibussowitsch   PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT ") %% bs %" PetscInt_FMT, Iend, Istart, bs);
6642e68589bSMark F. Adams 
6652e68589bSMark F. Adams   /* get 'nLocalSelected' */
6669566063dSJacob Faibussowitsch   PetscCall(PetscCDGetMIS(agg_lists, &selected_1));
6679566063dSJacob Faibussowitsch   PetscCall(ISGetSize(selected_1, &jj));
6689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(jj, &clid_flid));
6699566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(selected_1, &selected_idx));
670c8b0795cSMark F. Adams   for (kk = 0, nLocalSelected = 0; kk < jj; kk++) {
6712e68589bSMark F. Adams     PetscInt lid = selected_idx[kk];
672b43b03e9SMark F. Adams     if (lid < nloc) {
6739566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Gmat, lid + my0, &ncols, NULL, NULL));
6742fa5cd67SKarl Rupp       if (ncols > 1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */
6759566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Gmat, lid + my0, &ncols, NULL, NULL));
676b43b03e9SMark F. Adams     }
6772e68589bSMark F. Adams   }
6789566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(selected_1, &selected_idx));
6799566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&selected_1)); /* this is selected_1 in serial */
6802e68589bSMark F. Adams 
681d9558ea9SBarry Smith   /* create prolongator  matrix */
6829566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat, &mtype));
6839566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
6849566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * bs, PETSC_DETERMINE, PETSC_DETERMINE));
6859566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(Prol, bs, bs));
6869566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
6879566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol, 3 * data_cols, NULL));
6889566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol, 3 * data_cols, NULL, 3 * data_cols, NULL));
6892e68589bSMark F. Adams 
690c8b0795cSMark F. Adams   /* can get all points "removed" - but not on geomg */
6919566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &jj));
6927f66b68fSMark Adams   if (!jj) {
6939566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "ERROE: no selected points on coarse grid\n"));
6949566063dSJacob Faibussowitsch     PetscCall(PetscFree(clid_flid));
6959566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
6960298fd71SBarry Smith     *a_P_out = NULL; /* out */
6972e68589bSMark F. Adams     PetscFunctionReturn(0);
6982e68589bSMark F. Adams   }
6992e68589bSMark F. Adams 
7002e68589bSMark F. Adams   {
7012e68589bSMark F. Adams     PetscReal *coords;
702a2f3521dSMark F. Adams     PetscInt   data_stride;
7030298fd71SBarry Smith     PetscInt  *crsGID = NULL;
7042e68589bSMark F. Adams     Mat        Gmat2;
7052e68589bSMark F. Adams 
70663a3b9bcSJacob Faibussowitsch     PetscCheck(dim == data_cols, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " != data_cols %" PetscInt_FMT, dim, data_cols);
7072e68589bSMark F. Adams     /* grow ghost data for better coarse grid cover of fine grid */
708a2f3521dSMark F. Adams     /* messy method, squares graph and gets some data */
7099566063dSJacob Faibussowitsch     PetscCall(getGIDsOnSquareGraph(pc, nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID));
7102e68589bSMark F. Adams     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
7112e68589bSMark F. Adams     /* create global vector of coorindates in 'coords' */
712c5df96a5SBarry Smith     if (size > 1) {
7139566063dSJacob Faibussowitsch       PetscCall(PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords));
714806fa848SBarry Smith     } else {
715c8b0795cSMark F. Adams       coords      = (PetscReal *)pc_gamg->data;
716a2f3521dSMark F. Adams       data_stride = pc_gamg->data_sz / pc_gamg->data_cell_cols;
7172e68589bSMark F. Adams     }
7189566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Gmat2));
7192e68589bSMark F. Adams 
7202e68589bSMark F. Adams     /* triangulate */
7212e68589bSMark F. Adams     if (dim == 2) {
722c8b0795cSMark F. Adams       PetscReal metric, tm;
7239566063dSJacob Faibussowitsch       PetscCall(triangulateAndFormProl(selected_2, data_stride, coords, nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric));
7249566063dSJacob Faibussowitsch       PetscCall(PetscFree(crsGID));
7252e68589bSMark F. Adams 
7262e68589bSMark F. Adams       /* clean up and create coordinates for coarse grid (output) */
7279566063dSJacob Faibussowitsch       if (size > 1) PetscCall(PetscFree(coords));
7282e68589bSMark F. Adams 
7291c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm));
730c8b0795cSMark F. Adams       if (tm > 1.) { /* needs to be globalized - should not happen */
7319566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, " failed metric for coarse grid %e\n", (double)tm));
7329566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&Prol));
733806fa848SBarry Smith       } else if (metric > .0) {
7349566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "worst metric for coarse grid = %e\n", (double)metric));
7352e68589bSMark F. Adams       }
7363b4367a7SBarry Smith     } else SETERRQ(comm, PETSC_ERR_PLIB, "3D not implemented for 'geo' AMG");
7372e68589bSMark F. Adams     { /* create next coords - output */
7382e68589bSMark F. Adams       PetscReal *crs_crds;
7399566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dim * nLocalSelected, &crs_crds));
7402e68589bSMark F. Adams       for (kk = 0; kk < nLocalSelected; kk++) { /* grab local select nodes to promote - output */
741b43b03e9SMark F. Adams         PetscInt lid = clid_flid[kk];
742c8b0795cSMark F. Adams         for (jj = 0; jj < dim; jj++) crs_crds[jj * nLocalSelected + kk] = pc_gamg->data[jj * nloc + lid];
7432e68589bSMark F. Adams       }
744c8b0795cSMark F. Adams 
7459566063dSJacob Faibussowitsch       PetscCall(PetscFree(pc_gamg->data));
746c8b0795cSMark F. Adams       pc_gamg->data    = crs_crds; /* out */
747c8b0795cSMark F. Adams       pc_gamg->data_sz = dim * nLocalSelected;
7482e68589bSMark F. Adams     }
7499566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&selected_2));
7502e68589bSMark F. Adams   }
751a2f3521dSMark F. Adams 
7522e68589bSMark F. Adams   *a_P_out = Prol; /* out */
7539566063dSJacob Faibussowitsch   PetscCall(PetscFree(clid_flid));
754849bee69SMark Adams 
755c8b0795cSMark F. Adams   PetscFunctionReturn(0);
756c8b0795cSMark F. Adams }
757c8b0795cSMark F. Adams 
7589371c9d4SSatish Balay static PetscErrorCode PCDestroy_GAMG_GEO(PC pc) {
7599b8ffb57SJed Brown   PetscFunctionBegin;
7609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
7619b8ffb57SJed Brown   PetscFunctionReturn(0);
7629b8ffb57SJed Brown }
7639b8ffb57SJed Brown 
764c8b0795cSMark F. Adams /*
765c8b0795cSMark F. Adams  PCCreateGAMG_GEO
766c8b0795cSMark F. Adams 
767c8b0795cSMark F. Adams   Input Parameter:
768c8b0795cSMark F. Adams    . pc -
769c8b0795cSMark F. Adams */
7709371c9d4SSatish Balay PetscErrorCode PCCreateGAMG_GEO(PC pc) {
771c8b0795cSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
772c8b0795cSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
773c8b0795cSMark F. Adams 
774c8b0795cSMark F. Adams   PetscFunctionBegin;
7751ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO;
7769b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_GEO;
777c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
778c8b0795cSMark F. Adams 
779c8b0795cSMark F. Adams   /* set internal function pointers */
780fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_GEO;
781fd1112cbSBarry Smith   pc_gamg->ops->coarsen           = PCGAMGCoarsen_GEO;
7821ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_GEO;
783fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = NULL;
7841ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_GEO;
785c8b0795cSMark F. Adams 
7869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_GEO));
7872e68589bSMark F. Adams   PetscFunctionReturn(0);
7882e68589bSMark F. Adams }
789