xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 /*
282e68589bSMark F. Adams    PCSetCoordinates_GEO
292e68589bSMark F. Adams 
302e68589bSMark F. Adams    Input Parameter:
312e68589bSMark F. Adams    .  pc - the preconditioner context
322e68589bSMark F. Adams */
339371c9d4SSatish Balay PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) {
342e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
352e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
3690fbc344SStefano Zampini   PetscInt arrsz, bs, my0, kk, ii, nloc, Iend, aloc;
372e68589bSMark F. Adams   Mat      Amat = pc->pmat;
382e68589bSMark F. Adams 
392e68589bSMark F. Adams   PetscFunctionBegin;
402e68589bSMark F. Adams   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
419566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
429566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &my0, &Iend));
4390fbc344SStefano Zampini   aloc = (Iend - my0);
442e68589bSMark F. Adams   nloc = (Iend - my0) / bs;
45a2f3521dSMark F. Adams 
462472a847SBarry 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);
472e68589bSMark F. Adams 
48c8b0795cSMark F. Adams   pc_gamg->data_cell_rows = 1;
497827d75bSBarry Smith   PetscCheck(coords || (nloc <= 0), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
50c8b0795cSMark F. Adams   pc_gamg->data_cell_cols = ndm; /* coordinates */
512e68589bSMark F. Adams 
52c8b0795cSMark F. Adams   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
532e68589bSMark F. Adams 
542e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
557f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
569566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
582e68589bSMark F. Adams   }
592e68589bSMark F. Adams   for (kk = 0; kk < arrsz; kk++) pc_gamg->data[kk] = -999.;
602e68589bSMark F. Adams   pc_gamg->data[arrsz] = -99.;
612e68589bSMark F. Adams   /* copy data in - column oriented */
6290fbc344SStefano Zampini   if (nloc == a_nloc) {
632e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) {
649371c9d4SSatish Balay       for (ii = 0; ii < ndm; ii++) { pc_gamg->data[ii * nloc + kk] = coords[kk * ndm + ii]; }
652e68589bSMark F. Adams     }
6690fbc344SStefano Zampini   } else { /* assumes the coordinates are blocked */
6790fbc344SStefano Zampini     for (kk = 0; kk < nloc; kk++) {
689371c9d4SSatish Balay       for (ii = 0; ii < ndm; ii++) { pc_gamg->data[ii * nloc + kk] = coords[bs * kk * ndm + ii]; }
6990fbc344SStefano Zampini     }
7090fbc344SStefano Zampini   }
7163a3b9bcSJacob 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]);
722e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
732e68589bSMark F. Adams   PetscFunctionReturn(0);
742e68589bSMark F. Adams }
752e68589bSMark F. Adams 
762e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
772e68589bSMark F. Adams /*
782e68589bSMark F. Adams    PCSetData_GEO
792e68589bSMark F. Adams 
802e68589bSMark F. Adams   Input Parameter:
812e68589bSMark F. Adams    . pc -
822e68589bSMark F. Adams */
839371c9d4SSatish Balay PetscErrorCode PCSetData_GEO(PC pc, Mat m) {
842e68589bSMark F. Adams   PetscFunctionBegin;
85ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "GEO MG needs coordinates");
862e68589bSMark F. Adams }
872e68589bSMark F. Adams 
882e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
892e68589bSMark F. Adams /*
902e68589bSMark F. Adams    PCSetFromOptions_GEO
912e68589bSMark F. Adams 
922e68589bSMark F. Adams   Input Parameter:
932e68589bSMark F. Adams    . pc -
942e68589bSMark F. Adams */
959371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_GEO(PC pc, PetscOptionItems *PetscOptionsObject) {
962e68589bSMark F. Adams   PetscFunctionBegin;
97d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-GEO options");
982e68589bSMark F. Adams   {
992e68589bSMark F. Adams     /* -pc_gamg_sa_nsmooths */
1002e68589bSMark F. Adams     /* pc_gamg_sa->smooths = 0; */
1012e68589bSMark F. Adams     /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */
1022e68589bSMark F. Adams     /*                        "smoothing steps for smoothed aggregation, usually 1 (0)", */
1032e68589bSMark F. Adams     /*                        "PCGAMGSetNSmooths_AGG", */
1042e68589bSMark F. Adams     /*                        pc_gamg_sa->smooths, */
1052e68589bSMark F. Adams     /*                        &pc_gamg_sa->smooths, */
1062e68589bSMark F. Adams     /*                        &flag);  */
1072e68589bSMark F. Adams   }
108d0609cedSBarry Smith   PetscOptionsHeadEnd();
1092e68589bSMark F. Adams   PetscFunctionReturn(0);
1102e68589bSMark F. Adams }
1112e68589bSMark F. Adams 
1122e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1132e68589bSMark F. Adams /*
1142e68589bSMark F. Adams  triangulateAndFormProl
1152e68589bSMark F. Adams 
1162e68589bSMark F. Adams    Input Parameter:
1172e68589bSMark F. Adams    . selected_2 - list of selected local ID, includes selected ghosts
118a2f3521dSMark F. Adams    . data_stride -
119a2f3521dSMark F. Adams    . coords[2*data_stride] - column vector of local coordinates w/ ghosts
1202adfe9d3SBarry Smith    . nselected_1 - selected IDs that go with base (1) graph includes selected ghosts
1210cbbd2e1SMark F. Adams    . clid_lid_1[nselected_1] - lids of selected (c) nodes   ???????????
1222adfe9d3SBarry Smith    . agg_lists_1 - list of aggregates selected_1 vertices of aggregate unselected vertices
1232e68589bSMark F. Adams    . crsGID[selected.size()] - global index for prolongation operator
1242e68589bSMark F. Adams    . bs - block size
1252e68589bSMark F. Adams   Output Parameter:
1262e68589bSMark F. Adams    . a_Prol - prolongation operator
1272e68589bSMark F. Adams    . a_worst_best - measure of worst missed fine vertex, 0 is no misses
1282e68589bSMark F. Adams */
1299371c9d4SSatish 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) {
1302e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE)
131b43b03e9SMark F. Adams   PetscInt             jj, tid, tt, idx, nselected_2;
1322e68589bSMark F. Adams   struct triangulateio in, mid;
1330cbbd2e1SMark F. Adams   const PetscInt      *selected_idx_2;
13473911c69SBarry Smith   PetscMPIInt          rank;
1352e68589bSMark F. Adams   PetscInt             Istart, Iend, nFineLoc, myFine0;
1362e68589bSMark F. Adams   int                  kk, nPlotPts, sid;
1373b4367a7SBarry Smith   MPI_Comm             comm;
138c8b0795cSMark F. Adams   PetscReal            tm;
139c8b0795cSMark F. Adams 
1406e111a19SKarl Rupp   PetscFunctionBegin;
1419566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
1429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1439566063dSJacob Faibussowitsch   PetscCall(ISGetSize(selected_2, &nselected_2));
1442e68589bSMark F. Adams   if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */
1452e68589bSMark F. Adams     *a_worst_best = 100.0;                    /* this will cause a stop, but not globalized (should not happen) */
146806fa848SBarry Smith   } else *a_worst_best = 0.0;
1471c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm));
148c8b0795cSMark F. Adams   if (tm > 0.0) {
149c8b0795cSMark F. Adams     *a_worst_best = 100.0;
1502e68589bSMark F. Adams     PetscFunctionReturn(0);
1512e68589bSMark F. Adams   }
1529566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
1539371c9d4SSatish Balay   nFineLoc                   = (Iend - Istart) / bs;
1549371c9d4SSatish Balay   myFine0                    = Istart / bs;
1552e68589bSMark F. Adams   nPlotPts                   = nFineLoc; /* locals */
1564c500f23SPierre Jolivet   /* triangle */
1572e68589bSMark F. Adams   /* Define input points - in*/
1582e68589bSMark F. Adams   in.numberofpoints          = nselected_2;
1592e68589bSMark F. Adams   in.numberofpointattributes = 0;
1602e68589bSMark F. Adams   /* get nselected points */
1619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * nselected_2, &in.pointlist));
1629566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(selected_2, &selected_idx_2));
1632e68589bSMark F. Adams 
1642e68589bSMark F. Adams   for (kk = 0, sid = 0; kk < nselected_2; kk++, sid += 2) {
1652e68589bSMark F. Adams     PetscInt lid          = selected_idx_2[kk];
1662e68589bSMark F. Adams     in.pointlist[sid]     = coords[lid];
167a2f3521dSMark F. Adams     in.pointlist[sid + 1] = coords[data_stride + lid];
1682e68589bSMark F. Adams     if (lid >= nFineLoc) nPlotPts++;
1692e68589bSMark F. Adams   }
170c05e6a5bSJed Brown   PetscCheck(sid == 2 * nselected_2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "sid %d != 2*nselected_2 %" PetscInt_FMT, sid, nselected_2);
1712e68589bSMark F. Adams 
1722e68589bSMark F. Adams   in.numberofsegments      = 0;
1732e68589bSMark F. Adams   in.numberofedges         = 0;
1742e68589bSMark F. Adams   in.numberofholes         = 0;
1752e68589bSMark F. Adams   in.numberofregions       = 0;
1760a545947SLisandro Dalcin   in.trianglelist          = NULL;
1770a545947SLisandro Dalcin   in.segmentmarkerlist     = NULL;
1780a545947SLisandro Dalcin   in.pointattributelist    = NULL;
1790a545947SLisandro Dalcin   in.pointmarkerlist       = NULL;
1800a545947SLisandro Dalcin   in.triangleattributelist = NULL;
1810a545947SLisandro Dalcin   in.trianglearealist      = NULL;
1820a545947SLisandro Dalcin   in.segmentlist           = NULL;
1830a545947SLisandro Dalcin   in.holelist              = NULL;
1840a545947SLisandro Dalcin   in.regionlist            = NULL;
1850a545947SLisandro Dalcin   in.edgelist              = NULL;
1860a545947SLisandro Dalcin   in.edgemarkerlist        = NULL;
1870a545947SLisandro Dalcin   in.normlist              = NULL;
1882fa5cd67SKarl Rupp 
1892e68589bSMark F. Adams   /* triangulate */
1900a545947SLisandro Dalcin   mid.pointlist             = NULL; /* Not needed if -N switch used. */
1912e68589bSMark F. Adams   /* Not needed if -N switch used or number of point attributes is zero: */
1920a545947SLisandro Dalcin   mid.pointattributelist    = NULL;
1930a545947SLisandro Dalcin   mid.pointmarkerlist       = NULL; /* Not needed if -N or -B switch used. */
1940a545947SLisandro Dalcin   mid.trianglelist          = NULL; /* Not needed if -E switch used. */
1952e68589bSMark F. Adams   /* Not needed if -E switch used or number of triangle attributes is zero: */
1960a545947SLisandro Dalcin   mid.triangleattributelist = NULL;
1970a545947SLisandro Dalcin   mid.neighborlist          = NULL; /* Needed only if -n switch used. */
1982e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P not used: */
1990a545947SLisandro Dalcin   mid.segmentlist           = NULL;
2002e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
2010a545947SLisandro Dalcin   mid.segmentmarkerlist     = NULL;
2020a545947SLisandro Dalcin   mid.edgelist              = NULL; /* Needed only if -e switch used. */
2030a545947SLisandro Dalcin   mid.edgemarkerlist        = NULL; /* Needed if -e used and -B not used. */
2042e68589bSMark F. Adams   mid.numberoftriangles     = 0;
2052e68589bSMark F. Adams 
2062e68589bSMark F. Adams   /* Triangulate the points.  Switches are chosen to read and write a  */
2072e68589bSMark F. Adams   /*   PSLG (p), preserve the convex hull (c), number everything from  */
2082e68589bSMark F. Adams   /*   zero (z), assign a regional attribute to each element (A), and  */
2092e68589bSMark F. Adams   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
2102e68589bSMark F. Adams   /*   neighbor list (n).                                            */
2112e68589bSMark F. Adams   if (nselected_2 != 0) {  /* inactive processor */
2122e68589bSMark F. Adams     char args[] = "npczQ"; /* c is needed ? */
2132e68589bSMark F. Adams     triangulate(args, &in, &mid, (struct triangulateio *)NULL);
2142e68589bSMark F. Adams     /* output .poly files for 'showme' */
2152e68589bSMark F. Adams     if (!PETSC_TRUE) {
2162e68589bSMark F. Adams       static int level = 1;
2179371c9d4SSatish Balay       FILE      *file;
2189371c9d4SSatish Balay       char       fname[32];
2192e68589bSMark F. Adams 
2209371c9d4SSatish Balay       sprintf(fname, "C%d_%d.poly", level, rank);
2219371c9d4SSatish Balay       file = fopen(fname, "w");
2222e68589bSMark F. Adams       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
2232e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n", in.numberofpoints, 2, 0, 0);
2242e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2259371c9d4SSatish Balay       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]); }
2262e68589bSMark F. Adams       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
2272e68589bSMark F. Adams       fprintf(file, "%d  %d\n", 0, 0);
2282e68589bSMark F. Adams       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
2292e68589bSMark F. Adams       /* One line: <# of holes> */
2302e68589bSMark F. Adams       fprintf(file, "%d\n", 0);
2312e68589bSMark F. Adams       /* Following lines: <hole #> <x> <y> */
2322e68589bSMark F. Adams       /* Optional line: <# of regional attributes and/or area constraints> */
2332e68589bSMark F. Adams       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
2342e68589bSMark F. Adams       fclose(file);
2352e68589bSMark F. Adams 
2362e68589bSMark F. Adams       /* elems */
2379371c9d4SSatish Balay       sprintf(fname, "C%d_%d.ele", level, rank);
2389371c9d4SSatish Balay       file = fopen(fname, "w");
2392e68589bSMark F. Adams       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
2402e68589bSMark F. Adams       fprintf(file, "%d %d %d\n", mid.numberoftriangles, 3, 0);
2412e68589bSMark F. Adams       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
2429371c9d4SSatish Balay       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]); }
2432e68589bSMark F. Adams       fclose(file);
2442e68589bSMark F. Adams 
2459371c9d4SSatish Balay       sprintf(fname, "C%d_%d.node", level, rank);
2469371c9d4SSatish Balay       file = fopen(fname, "w");
2472e68589bSMark F. Adams       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
2482e68589bSMark F. Adams       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
2492e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n", nPlotPts, 2, 0, 0);
2502e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2519371c9d4SSatish Balay       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]); }
2522e68589bSMark F. Adams 
2532e68589bSMark F. Adams       sid /= 2;
2542e68589bSMark F. Adams       for (jj = 0; jj < nFineLoc; jj++) {
2552e68589bSMark F. Adams         PetscBool sel = PETSC_TRUE;
2562e68589bSMark F. Adams         for (kk = 0; kk < nselected_2 && sel; kk++) {
2572e68589bSMark F. Adams           PetscInt lid = selected_idx_2[kk];
2582e68589bSMark F. Adams           if (lid == jj) sel = PETSC_FALSE;
2592e68589bSMark F. Adams         }
2602fa5cd67SKarl Rupp         if (sel) fprintf(file, "%d %e %e\n", sid++, coords[jj], coords[data_stride + jj]);
2612e68589bSMark F. Adams       }
2622e68589bSMark F. Adams       fclose(file);
263c05e6a5bSJed Brown       PetscCheck(sid == nPlotPts, PETSC_COMM_SELF, PETSC_ERR_PLIB, "sid %d != nPlotPts %d", sid, nPlotPts);
2642e68589bSMark F. Adams       level++;
2652e68589bSMark F. Adams     }
2662e68589bSMark F. Adams   }
2672e68589bSMark F. Adams   { /* form P - setup some maps */
2680cbbd2e1SMark F. Adams     PetscInt clid, mm, *nTri, *node_tri;
2692e68589bSMark F. Adams 
2709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nselected_2, &node_tri, nselected_2, &nTri));
2712e68589bSMark F. Adams 
2722e68589bSMark F. Adams     /* need list of triangles on node */
2732e68589bSMark F. Adams     for (kk = 0; kk < nselected_2; kk++) nTri[kk] = 0;
2742e68589bSMark F. Adams     for (tid = 0, kk = 0; tid < mid.numberoftriangles; tid++) {
2752e68589bSMark F. Adams       for (jj = 0; jj < 3; jj++) {
2762e68589bSMark F. Adams         PetscInt cid = mid.trianglelist[kk++];
2772e68589bSMark F. Adams         if (nTri[cid] == 0) node_tri[cid] = tid;
2782e68589bSMark F. Adams         nTri[cid]++;
2792e68589bSMark F. Adams       }
2802e68589bSMark F. Adams     }
2812e68589bSMark F. Adams #define EPS 1.e-12
2822e68589bSMark F. Adams     /* find points and set prolongation */
2830cbbd2e1SMark F. Adams     for (mm = clid = 0; mm < nFineLoc; mm++) {
284e78576d6SMark F. Adams       PetscBool ise;
2859566063dSJacob Faibussowitsch       PetscCall(PetscCDEmptyAt(agg_lists_1, mm, &ise));
286e78576d6SMark F. Adams       if (!ise) {
2870cbbd2e1SMark F. Adams         const PetscInt lid = mm;
2882e68589bSMark F. Adams         PetscScalar    AA[3][3];
2892e68589bSMark F. Adams         PetscBLASInt   N = 3, NRHS = 1, LDA = 3, IPIV[3], LDB = 3, INFO;
290539c167fSBarry Smith         PetscCDIntNd  *pos;
291539c167fSBarry Smith 
2929566063dSJacob Faibussowitsch         PetscCall(PetscCDGetHeadPos(agg_lists_1, lid, &pos));
293e78576d6SMark F. Adams         while (pos) {
294e78576d6SMark F. Adams           PetscInt flid;
2959566063dSJacob Faibussowitsch           PetscCall(PetscCDIntNdGetID(pos, &flid));
2969566063dSJacob Faibussowitsch           PetscCall(PetscCDGetNextPos(agg_lists_1, lid, &pos));
2970cbbd2e1SMark F. Adams 
2982e68589bSMark F. Adams           if (flid < nFineLoc) { /* could be a ghost */
2999371c9d4SSatish Balay             PetscInt        bestTID    = -1;
3009371c9d4SSatish Balay             PetscReal       best_alpha = 1.e10;
3012e68589bSMark F. Adams             const PetscInt  fgid       = flid + myFine0;
3022e68589bSMark F. Adams             /* compute shape function for gid */
303a2f3521dSMark F. Adams             const PetscReal fcoord[3]  = {coords[flid], coords[data_stride + flid], 1.0};
3049371c9d4SSatish Balay             PetscBool       haveit     = PETSC_FALSE;
3059371c9d4SSatish Balay             PetscScalar     alpha[3];
3069371c9d4SSatish Balay             PetscInt        clids[3];
3072fa5cd67SKarl Rupp 
3082e68589bSMark F. Adams             /* look for it */
3099371c9d4SSatish Balay             for (tid = node_tri[clid], jj = 0; jj < 5 && !haveit && tid != -1; jj++) {
3102e68589bSMark F. Adams               for (tt = 0; tt < 3; tt++) {
3112e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3 * tid + tt];
3122e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
3139371c9d4SSatish Balay                 AA[tt][0]     = coords[lid2];
3149371c9d4SSatish Balay                 AA[tt][1]     = coords[data_stride + lid2];
3159371c9d4SSatish Balay                 AA[tt][2]     = 1.0;
3162e68589bSMark F. Adams                 clids[tt]     = cid2; /* store for interp */
3172e68589bSMark F. Adams               }
3182e68589bSMark F. Adams 
3192e68589bSMark F. Adams               for (tt = 0; tt < 3; tt++) alpha[tt] = (PetscScalar)fcoord[tt];
3202e68589bSMark F. Adams 
3212e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
322792fecdfSBarry Smith               PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&N, &NRHS, (PetscScalar *)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3232e68589bSMark F. Adams               {
3249371c9d4SSatish Balay                 PetscBool have   = PETSC_TRUE;
3259371c9d4SSatish Balay                 PetscReal lowest = 1.e10;
3262e68589bSMark F. Adams                 for (tt = 0, idx = 0; tt < 3; tt++) {
3272e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) > (1.0 + EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
3282e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) < lowest) {
3292e68589bSMark F. Adams                     lowest = PetscRealPart(alpha[tt]);
3302e68589bSMark F. Adams                     idx    = tt;
3312e68589bSMark F. Adams                   }
3322e68589bSMark F. Adams                 }
3332e68589bSMark F. Adams                 haveit = have;
3342e68589bSMark F. Adams               }
3352e68589bSMark F. Adams               tid = mid.neighborlist[3 * tid + idx];
3362e68589bSMark F. Adams             }
3372e68589bSMark F. Adams 
3382e68589bSMark F. Adams             if (!haveit) {
3392e68589bSMark F. Adams               /* brute force */
3402e68589bSMark F. Adams               for (tid = 0; tid < mid.numberoftriangles && !haveit; tid++) {
3412e68589bSMark F. Adams                 for (tt = 0; tt < 3; tt++) {
3422e68589bSMark F. Adams                   PetscInt cid2 = mid.trianglelist[3 * tid + tt];
3432e68589bSMark F. Adams                   PetscInt lid2 = selected_idx_2[cid2];
3449371c9d4SSatish Balay                   AA[tt][0]     = coords[lid2];
3459371c9d4SSatish Balay                   AA[tt][1]     = coords[data_stride + lid2];
3469371c9d4SSatish Balay                   AA[tt][2]     = 1.0;
3472e68589bSMark F. Adams                   clids[tt]     = cid2; /* store for interp */
3482e68589bSMark F. Adams                 }
3492e68589bSMark F. Adams                 for (tt = 0; tt < 3; tt++) alpha[tt] = fcoord[tt];
3502e68589bSMark F. Adams                 /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
351792fecdfSBarry Smith                 PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&N, &NRHS, (PetscScalar *)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3522e68589bSMark F. Adams                 {
3539371c9d4SSatish Balay                   PetscBool have  = PETSC_TRUE;
3549371c9d4SSatish Balay                   PetscReal worst = 0.0, v;
3552e68589bSMark F. Adams                   for (tt = 0; tt < 3 && have; tt++) {
3562e68589bSMark F. Adams                     if (PetscRealPart(alpha[tt]) > 1.0 + EPS || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
3572e68589bSMark F. Adams                     if ((v = PetscAbs(PetscRealPart(alpha[tt]) - 0.5)) > worst) worst = v;
3582e68589bSMark F. Adams                   }
3592e68589bSMark F. Adams                   if (worst < best_alpha) {
3609371c9d4SSatish Balay                     best_alpha = worst;
3619371c9d4SSatish Balay                     bestTID    = tid;
3622e68589bSMark F. Adams                   }
3632e68589bSMark F. Adams                   haveit = have;
3642e68589bSMark F. Adams                 }
3652e68589bSMark F. Adams               }
3662e68589bSMark F. Adams             }
3672e68589bSMark F. Adams             if (!haveit) {
3682e68589bSMark F. Adams               if (best_alpha > *a_worst_best) *a_worst_best = best_alpha;
3692e68589bSMark F. Adams               /* use best one */
3702e68589bSMark F. Adams               for (tt = 0; tt < 3; tt++) {
3712e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3 * bestTID + tt];
3722e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
3739371c9d4SSatish Balay                 AA[tt][0]     = coords[lid2];
3749371c9d4SSatish Balay                 AA[tt][1]     = coords[data_stride + lid2];
3759371c9d4SSatish Balay                 AA[tt][2]     = 1.0;
3762e68589bSMark F. Adams                 clids[tt]     = cid2; /* store for interp */
3772e68589bSMark F. Adams               }
3782e68589bSMark F. Adams               for (tt = 0; tt < 3; tt++) alpha[tt] = fcoord[tt];
3792e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
380792fecdfSBarry Smith               PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&N, &NRHS, (PetscScalar *)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3812e68589bSMark F. Adams             }
3822e68589bSMark F. Adams 
3832e68589bSMark F. Adams             /* put in row of P */
3842e68589bSMark F. Adams             for (idx = 0; idx < 3; idx++) {
3852e68589bSMark F. Adams               PetscScalar shp = alpha[idx];
3862e68589bSMark F. Adams               if (PetscAbs(PetscRealPart(shp)) > 1.e-6) {
3872e68589bSMark F. Adams                 PetscInt cgid = crsGID[clids[idx]];
3882e68589bSMark F. Adams                 PetscInt jj = cgid * bs, ii = fgid * bs; /* need to gloalize */
389*48a46eb9SPierre Jolivet                 for (tt = 0; tt < bs; tt++, ii++, jj++) PetscCall(MatSetValues(a_Prol, 1, &ii, 1, &jj, &shp, INSERT_VALUES));
3902e68589bSMark F. Adams               }
3912e68589bSMark F. Adams             }
3922e68589bSMark F. Adams           }
3930cbbd2e1SMark F. Adams         } /* aggregates iterations */
3940cbbd2e1SMark F. Adams         clid++;
3950cbbd2e1SMark F. Adams       } /* a coarse agg */
3960cbbd2e1SMark F. Adams     }   /* for all fine nodes */
3970cbbd2e1SMark F. Adams 
3989566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(selected_2, &selected_idx_2));
3999566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
4009566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
4012e68589bSMark F. Adams 
4029566063dSJacob Faibussowitsch     PetscCall(PetscFree2(node_tri, nTri));
4032e68589bSMark F. Adams   }
4042e68589bSMark F. Adams   free(mid.trianglelist);
4052e68589bSMark F. Adams   free(mid.neighborlist);
406c6b50ea1SMark   free(mid.segmentlist);
407c6b50ea1SMark   free(mid.segmentmarkerlist);
408c6b50ea1SMark   free(mid.pointlist);
409c6b50ea1SMark   free(mid.pointmarkerlist);
4109566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointlist));
4112e68589bSMark F. Adams   PetscFunctionReturn(0);
4122e68589bSMark F. Adams #else
4133b4367a7SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)a_Prol), PETSC_ERR_PLIB, "configure with TRIANGLE to use geometric MG");
4142e68589bSMark F. Adams #endif
4152e68589bSMark F. Adams }
4162e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
4172e68589bSMark F. Adams /*
4182e68589bSMark F. Adams    getGIDsOnSquareGraph - square graph, get
4192e68589bSMark F. Adams 
4202e68589bSMark F. Adams    Input Parameter:
4210cbbd2e1SMark F. Adams    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
422b43b03e9SMark F. Adams    . clid_lid_1 - [nselected_1] lids of selected nodes
4232e68589bSMark F. Adams    . Gmat1 - graph that goes with 'selected_1'
4242e68589bSMark F. Adams    Output Parameter:
4252e68589bSMark F. Adams    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
4262e68589bSMark F. Adams    . a_Gmat_2 - graph that is squared of 'Gmat_1'
4272e68589bSMark F. Adams    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
4282e68589bSMark F. Adams */
4299371c9d4SSatish 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) {
43073911c69SBarry Smith   PetscMPIInt size;
431b43b03e9SMark F. Adams   PetscInt   *crsGID, kk, my0, Iend, nloc;
4323b4367a7SBarry Smith   MPI_Comm    comm;
4332e68589bSMark F. Adams 
4342e68589bSMark F. Adams   PetscFunctionBegin;
4359566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm));
4369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
4379566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &my0, &Iend)); /* AIJ */
4382e68589bSMark F. Adams   nloc = Iend - my0;                                   /* this does not change */
4392e68589bSMark F. Adams 
440c5df96a5SBarry Smith   if (size == 1) { /* not much to do in serial */
4419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected_1, &crsGID));
442b43b03e9SMark F. Adams     for (kk = 0; kk < nselected_1; kk++) crsGID[kk] = kk;
4430a545947SLisandro Dalcin     *a_Gmat_2 = NULL;
4449566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nselected_1, clid_lid_1, PETSC_COPY_VALUES, a_selected_2));
445806fa848SBarry Smith   } else {
446b43b03e9SMark F. Adams     PetscInt     idx, num_fine_ghosts, num_crs_ghost, myCrs0;
4472e68589bSMark F. Adams     Mat_MPIAIJ  *mpimat2;
4482e68589bSMark F. Adams     Mat          Gmat2;
4492e68589bSMark F. Adams     Vec          locState;
4502e68589bSMark F. Adams     PetscScalar *cpcol_state;
4512e68589bSMark F. Adams 
4522e68589bSMark F. Adams     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
453b43b03e9SMark F. Adams     kk = nselected_1;
4549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPI_SUM, comm));
455b43b03e9SMark F. Adams     myCrs0 -= nselected_1;
4562e68589bSMark F. Adams 
457b43b03e9SMark F. Adams     if (a_Gmat_2) { /* output */
4582e68589bSMark F. Adams       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
4599566063dSJacob Faibussowitsch       PetscCall(PCGAMGSquareGraph_GAMG(pc, Gmat1, &Gmat2));
4602e68589bSMark F. Adams       *a_Gmat_2 = Gmat2;  /* output */
461806fa848SBarry Smith     } else Gmat2 = Gmat1; /* use local to get crsGIDs at least */
4622e68589bSMark F. Adams     /* get coarse grid GIDS for selected (locals and ghosts) */
4632e68589bSMark F. Adams     mpimat2 = (Mat_MPIAIJ *)Gmat2->data;
4649566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Gmat2, &locState, NULL));
4659566063dSJacob Faibussowitsch     PetscCall(VecSet(locState, (PetscScalar)(PetscReal)(-1))); /* set with UNKNOWN state */
466b43b03e9SMark F. Adams     for (kk = 0; kk < nselected_1; kk++) {
467b43b03e9SMark F. Adams       PetscInt    fgid = clid_lid_1[kk] + my0;
4682e68589bSMark F. Adams       PetscScalar v    = (PetscScalar)(kk + myCrs0);
4699566063dSJacob Faibussowitsch       PetscCall(VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES)); /* set with PID */
4702e68589bSMark F. Adams     }
4719566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(locState));
4729566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(locState));
4739566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mpimat2->Mvctx, locState, mpimat2->lvec, INSERT_VALUES, SCATTER_FORWARD));
4749566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mpimat2->Mvctx, locState, mpimat2->lvec, INSERT_VALUES, SCATTER_FORWARD));
4759566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts));
4769566063dSJacob Faibussowitsch     PetscCall(VecGetArray(mpimat2->lvec, &cpcol_state));
4772e68589bSMark F. Adams     for (kk = 0, num_crs_ghost = 0; kk < num_fine_ghosts; kk++) {
4782e68589bSMark F. Adams       if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++;
4792e68589bSMark F. Adams     }
4809566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected_1 + num_crs_ghost, &crsGID)); /* output */
4812e68589bSMark F. Adams     {
4822e68589bSMark F. Adams       PetscInt *selected_set;
4839566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nselected_1 + num_crs_ghost, &selected_set));
4842e68589bSMark F. Adams       /* do ghost of 'crsGID' */
485b43b03e9SMark F. Adams       for (kk = 0, idx = nselected_1; kk < num_fine_ghosts; kk++) {
4862e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
4872e68589bSMark F. Adams           PetscInt cgid     = (PetscInt)PetscRealPart(cpcol_state[kk]);
4882e68589bSMark F. Adams           selected_set[idx] = nloc + kk;
4892e68589bSMark F. Adams           crsGID[idx++]     = cgid;
4902e68589bSMark F. Adams         }
4912e68589bSMark F. Adams       }
49263a3b9bcSJacob 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);
4939566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(mpimat2->lvec, &cpcol_state));
4942e68589bSMark F. Adams       /* do locals in 'crsGID' */
4959566063dSJacob Faibussowitsch       PetscCall(VecGetArray(locState, &cpcol_state));
4962e68589bSMark F. Adams       for (kk = 0, idx = 0; kk < nloc; kk++) {
4972e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
4982e68589bSMark F. Adams           PetscInt cgid     = (PetscInt)PetscRealPart(cpcol_state[kk]);
4992e68589bSMark F. Adams           selected_set[idx] = kk;
5002e68589bSMark F. Adams           crsGID[idx++]     = cgid;
5012e68589bSMark F. Adams         }
5022e68589bSMark F. Adams       }
50363a3b9bcSJacob Faibussowitsch       PetscCheck(idx == nselected_1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "idx %" PetscInt_FMT " != nselected_1 %" PetscInt_FMT, idx, nselected_1);
5049566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(locState, &cpcol_state));
5052e68589bSMark F. Adams 
5060a545947SLisandro Dalcin       if (a_selected_2 != NULL) { /* output */
5079566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, (nselected_1 + num_crs_ghost), selected_set, PETSC_OWN_POINTER, a_selected_2));
508806fa848SBarry Smith       } else {
5099566063dSJacob Faibussowitsch         PetscCall(PetscFree(selected_set));
5102e68589bSMark F. Adams       }
5110cbbd2e1SMark F. Adams     }
5129566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&locState));
5132e68589bSMark F. Adams   }
5142e68589bSMark F. Adams   *a_crsGID = crsGID; /* output */
5152e68589bSMark F. Adams   PetscFunctionReturn(0);
5162e68589bSMark F. Adams }
5172e68589bSMark F. Adams 
5182e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
5192e68589bSMark F. Adams /*
520fd1112cbSBarry Smith    PCGAMGGraph_GEO
5212e68589bSMark F. Adams 
5222e68589bSMark F. Adams   Input Parameter:
5232e68589bSMark F. Adams    . pc - this
5242e68589bSMark F. Adams    . Amat - matrix on this fine level
5252e68589bSMark F. Adams   Output Parameter:
526c8b0795cSMark F. Adams    . a_Gmat
5272e68589bSMark F. Adams */
5289371c9d4SSatish Balay PetscErrorCode PCGAMGGraph_GEO(PC pc, Mat Amat, Mat *a_Gmat) {
529c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG *)pc->data;
530c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG *)mg->innerctx;
531c1eae691SMark Adams   const PetscReal vfilter = pc_gamg->threshold[0];
5323b4367a7SBarry Smith   MPI_Comm        comm;
53372833a62Smarkadams4   Mat             Gmat, F = NULL;
5340cbbd2e1SMark F. Adams   PetscBool       set, flg, symm;
5356e111a19SKarl Rupp 
536c8b0795cSMark F. Adams   PetscFunctionBegin;
5379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
538c8b0795cSMark F. Adams 
5399566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetricKnown(Amat, &set, &flg));
540263489e9SJed Brown   symm = (PetscBool) !(set && flg);
5410cbbd2e1SMark F. Adams 
54272833a62Smarkadams4   PetscCall(MatCreateGraph(Amat, symm, PETSC_TRUE, &Gmat));
54372833a62Smarkadams4   PetscCall(MatFilter(Gmat, vfilter, &F));
54472833a62Smarkadams4   if (F) {
54572833a62Smarkadams4     PetscCall(MatDestroy(&Gmat));
54672833a62Smarkadams4     Gmat = F;
54772833a62Smarkadams4   }
548c8b0795cSMark F. Adams   *a_Gmat = Gmat;
549849bee69SMark Adams 
550c8b0795cSMark F. Adams   PetscFunctionReturn(0);
551c8b0795cSMark F. Adams }
552c8b0795cSMark F. Adams 
553c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
554c8b0795cSMark F. Adams /*
555fd1112cbSBarry Smith    PCGAMGCoarsen_GEO
556c8b0795cSMark F. Adams 
557c8b0795cSMark F. Adams   Input Parameter:
558e0940f08SMark F. Adams    . a_pc - this
559e0940f08SMark F. Adams    . a_Gmat - graph
560c8b0795cSMark F. Adams   Output Parameter:
561c8b0795cSMark F. Adams    . a_llist_parent - linked list from selected indices for data locality only
562c8b0795cSMark F. Adams */
5639371c9d4SSatish Balay PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc, Mat *a_Gmat, PetscCoarsenData **a_llist_parent) {
564c8b0795cSMark F. Adams   PetscInt   Istart, Iend, nloc, kk, Ii, ncols;
5650cbbd2e1SMark F. Adams   IS         perm;
566c8b0795cSMark F. Adams   GAMGNode  *gnodes;
567c8b0795cSMark F. Adams   PetscInt  *permute;
568e0940f08SMark F. Adams   Mat        Gmat = *a_Gmat;
5693b4367a7SBarry Smith   MPI_Comm   comm;
570b43b03e9SMark F. Adams   MatCoarsen crs;
571c8b0795cSMark F. Adams 
572c8b0795cSMark F. Adams   PetscFunctionBegin;
5739566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_pc, &comm));
574849bee69SMark Adams 
5759566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend));
576c8b0795cSMark F. Adams   nloc = (Iend - Istart);
577c8b0795cSMark F. Adams 
578c8b0795cSMark F. Adams   /* create random permutation with sort for geo-mg */
5799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc, &gnodes));
5809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc, &permute));
581c8b0795cSMark F. Adams 
582c8b0795cSMark F. Adams   for (Ii = Istart; Ii < Iend; Ii++) { /* locals only? */
5839566063dSJacob Faibussowitsch     PetscCall(MatGetRow(Gmat, Ii, &ncols, NULL, NULL));
584c8b0795cSMark F. Adams     {
585c8b0795cSMark F. Adams       PetscInt lid       = Ii - Istart;
586c8b0795cSMark F. Adams       gnodes[lid].lid    = lid;
587c8b0795cSMark F. Adams       gnodes[lid].degree = ncols;
588c8b0795cSMark F. Adams     }
5899566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(Gmat, Ii, &ncols, NULL, NULL));
590c8b0795cSMark F. Adams   }
591c8b0795cSMark F. Adams   if (PETSC_TRUE) {
592e2c00dcbSBarry Smith     PetscRandom rand;
593c8b0795cSMark F. Adams     PetscBool  *bIndexSet;
594e2c00dcbSBarry Smith     PetscReal   rr;
595e2c00dcbSBarry Smith     PetscInt    iSwapIndex;
596e2c00dcbSBarry Smith 
5979566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
5989566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nloc, &bIndexSet));
5992fa5cd67SKarl Rupp     for (Ii = 0; Ii < nloc; Ii++) {
6009566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(rand, &rr));
601e2c00dcbSBarry Smith       iSwapIndex = (PetscInt)(rr * nloc);
6022fa5cd67SKarl Rupp       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
603c8b0795cSMark F. Adams         GAMGNode iTemp        = gnodes[iSwapIndex];
604c8b0795cSMark F. Adams         gnodes[iSwapIndex]    = gnodes[Ii];
605c8b0795cSMark F. Adams         gnodes[Ii]            = iTemp;
606c8b0795cSMark F. Adams         bIndexSet[Ii]         = PETSC_TRUE;
607c8b0795cSMark F. Adams         bIndexSet[iSwapIndex] = PETSC_TRUE;
608c8b0795cSMark F. Adams       }
609c8b0795cSMark F. Adams     }
6109566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
6119566063dSJacob Faibussowitsch     PetscCall(PetscFree(bIndexSet));
612c8b0795cSMark F. Adams   }
613c8b0795cSMark F. Adams   /* only sort locals */
6140cbbd2e1SMark F. Adams   qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare);
615c8b0795cSMark F. Adams   /* create IS of permutation */
6162fa5cd67SKarl Rupp   for (kk = 0; kk < nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */
6179566063dSJacob Faibussowitsch   PetscCall(PetscFree(gnodes));
6189566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm));
619c8b0795cSMark F. Adams 
620c8b0795cSMark F. Adams   /* get MIS aggs */
621b43b03e9SMark F. Adams 
6229566063dSJacob Faibussowitsch   PetscCall(MatCoarsenCreate(comm, &crs));
6239566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetType(crs, MATCOARSENMIS));
6249566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetGreedyOrdering(crs, perm));
6259566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetAdjacency(crs, Gmat));
6269566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_FALSE));
6279566063dSJacob Faibussowitsch   PetscCall(MatCoarsenApply(crs));
6289566063dSJacob Faibussowitsch   PetscCall(MatCoarsenGetData(crs, a_llist_parent));
6299566063dSJacob Faibussowitsch   PetscCall(MatCoarsenDestroy(&crs));
630c8b0795cSMark F. Adams 
6319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
632849bee69SMark Adams 
633c8b0795cSMark F. Adams   PetscFunctionReturn(0);
634c8b0795cSMark F. Adams }
635c8b0795cSMark F. Adams 
636c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
637c8b0795cSMark F. Adams /*
6380cbbd2e1SMark F. Adams  PCGAMGProlongator_GEO
639c8b0795cSMark F. Adams 
640c8b0795cSMark F. Adams  Input Parameter:
641c8b0795cSMark F. Adams  . pc - this
642c8b0795cSMark F. Adams  . Amat - matrix on this fine level
643c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
6440cbbd2e1SMark F. Adams  . selected_1 - [nselected]
6450cbbd2e1SMark F. Adams  . agg_lists - [nselected]
646c8b0795cSMark F. Adams  Output Parameter:
647c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
648c8b0795cSMark F. Adams  */
6499371c9d4SSatish Balay PetscErrorCode PCGAMGProlongator_GEO(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out) {
6502e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG *)pc->data;
6512e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG *)mg->innerctx;
652c8b0795cSMark F. Adams   const PetscInt  dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
653b43b03e9SMark F. Adams   PetscInt        Istart, Iend, nloc, my0, jj, kk, ncols, nLocalSelected, bs, *clid_flid;
654c8b0795cSMark F. Adams   Mat             Prol;
655c5df96a5SBarry Smith   PetscMPIInt     rank, size;
6563b4367a7SBarry Smith   MPI_Comm        comm;
6570cbbd2e1SMark F. Adams   IS              selected_2, selected_1;
6582e68589bSMark F. Adams   const PetscInt *selected_idx;
659d9558ea9SBarry Smith   MatType         mtype;
6602e68589bSMark F. Adams 
6612e68589bSMark F. Adams   PetscFunctionBegin;
6629566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
663849bee69SMark Adams 
6649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
6659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6669566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
6679566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
6689371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
6699371c9d4SSatish Balay   my0  = Istart / bs;
67063a3b9bcSJacob Faibussowitsch   PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT ") %% bs %" PetscInt_FMT, Iend, Istart, bs);
6712e68589bSMark F. Adams 
6722e68589bSMark F. Adams   /* get 'nLocalSelected' */
6739566063dSJacob Faibussowitsch   PetscCall(PetscCDGetMIS(agg_lists, &selected_1));
6749566063dSJacob Faibussowitsch   PetscCall(ISGetSize(selected_1, &jj));
6759566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(jj, &clid_flid));
6769566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(selected_1, &selected_idx));
677c8b0795cSMark F. Adams   for (kk = 0, nLocalSelected = 0; kk < jj; kk++) {
6782e68589bSMark F. Adams     PetscInt lid = selected_idx[kk];
679b43b03e9SMark F. Adams     if (lid < nloc) {
6809566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Gmat, lid + my0, &ncols, NULL, NULL));
6812fa5cd67SKarl Rupp       if (ncols > 1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */
6829566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Gmat, lid + my0, &ncols, NULL, NULL));
683b43b03e9SMark F. Adams     }
6842e68589bSMark F. Adams   }
6859566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(selected_1, &selected_idx));
6869566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&selected_1)); /* this is selected_1 in serial */
6872e68589bSMark F. Adams 
688d9558ea9SBarry Smith   /* create prolongator  matrix */
6899566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat, &mtype));
6909566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
6919566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * bs, PETSC_DETERMINE, PETSC_DETERMINE));
6929566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(Prol, bs, bs));
6939566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
6949566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol, 3 * data_cols, NULL));
6959566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol, 3 * data_cols, NULL, 3 * data_cols, NULL));
6962e68589bSMark F. Adams 
697c8b0795cSMark F. Adams   /* can get all points "removed" - but not on geomg */
6989566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &jj));
6997f66b68fSMark Adams   if (!jj) {
7009566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "ERROE: no selected points on coarse grid\n"));
7019566063dSJacob Faibussowitsch     PetscCall(PetscFree(clid_flid));
7029566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
7030298fd71SBarry Smith     *a_P_out = NULL; /* out */
7042e68589bSMark F. Adams     PetscFunctionReturn(0);
7052e68589bSMark F. Adams   }
7062e68589bSMark F. Adams 
7072e68589bSMark F. Adams   {
7082e68589bSMark F. Adams     PetscReal *coords;
709a2f3521dSMark F. Adams     PetscInt   data_stride;
7100298fd71SBarry Smith     PetscInt  *crsGID = NULL;
7112e68589bSMark F. Adams     Mat        Gmat2;
7122e68589bSMark F. Adams 
71363a3b9bcSJacob Faibussowitsch     PetscCheck(dim == data_cols, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " != data_cols %" PetscInt_FMT, dim, data_cols);
7142e68589bSMark F. Adams     /* grow ghost data for better coarse grid cover of fine grid */
715a2f3521dSMark F. Adams     /* messy method, squares graph and gets some data */
7169566063dSJacob Faibussowitsch     PetscCall(getGIDsOnSquareGraph(pc, nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID));
7172e68589bSMark F. Adams     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
7182e68589bSMark F. Adams     /* create global vector of coorindates in 'coords' */
719c5df96a5SBarry Smith     if (size > 1) {
7209566063dSJacob Faibussowitsch       PetscCall(PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords));
721806fa848SBarry Smith     } else {
722c8b0795cSMark F. Adams       coords      = (PetscReal *)pc_gamg->data;
723a2f3521dSMark F. Adams       data_stride = pc_gamg->data_sz / pc_gamg->data_cell_cols;
7242e68589bSMark F. Adams     }
7259566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Gmat2));
7262e68589bSMark F. Adams 
7272e68589bSMark F. Adams     /* triangulate */
7282e68589bSMark F. Adams     if (dim == 2) {
729c8b0795cSMark F. Adams       PetscReal metric, tm;
7309566063dSJacob Faibussowitsch       PetscCall(triangulateAndFormProl(selected_2, data_stride, coords, nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric));
7319566063dSJacob Faibussowitsch       PetscCall(PetscFree(crsGID));
7322e68589bSMark F. Adams 
7332e68589bSMark F. Adams       /* clean up and create coordinates for coarse grid (output) */
7349566063dSJacob Faibussowitsch       if (size > 1) PetscCall(PetscFree(coords));
7352e68589bSMark F. Adams 
7361c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm));
737c8b0795cSMark F. Adams       if (tm > 1.) { /* needs to be globalized - should not happen */
7389566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, " failed metric for coarse grid %e\n", (double)tm));
7399566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&Prol));
740806fa848SBarry Smith       } else if (metric > .0) {
7419566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "worst metric for coarse grid = %e\n", (double)metric));
7422e68589bSMark F. Adams       }
7433b4367a7SBarry Smith     } else SETERRQ(comm, PETSC_ERR_PLIB, "3D not implemented for 'geo' AMG");
7442e68589bSMark F. Adams     { /* create next coords - output */
7452e68589bSMark F. Adams       PetscReal *crs_crds;
7469566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dim * nLocalSelected, &crs_crds));
7472e68589bSMark F. Adams       for (kk = 0; kk < nLocalSelected; kk++) { /* grab local select nodes to promote - output */
748b43b03e9SMark F. Adams         PetscInt lid = clid_flid[kk];
749c8b0795cSMark F. Adams         for (jj = 0; jj < dim; jj++) crs_crds[jj * nLocalSelected + kk] = pc_gamg->data[jj * nloc + lid];
7502e68589bSMark F. Adams       }
751c8b0795cSMark F. Adams 
7529566063dSJacob Faibussowitsch       PetscCall(PetscFree(pc_gamg->data));
753c8b0795cSMark F. Adams       pc_gamg->data    = crs_crds; /* out */
754c8b0795cSMark F. Adams       pc_gamg->data_sz = dim * nLocalSelected;
7552e68589bSMark F. Adams     }
7569566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&selected_2));
7572e68589bSMark F. Adams   }
758a2f3521dSMark F. Adams 
7592e68589bSMark F. Adams   *a_P_out = Prol; /* out */
7609566063dSJacob Faibussowitsch   PetscCall(PetscFree(clid_flid));
761849bee69SMark Adams 
762c8b0795cSMark F. Adams   PetscFunctionReturn(0);
763c8b0795cSMark F. Adams }
764c8b0795cSMark F. Adams 
7659371c9d4SSatish Balay static PetscErrorCode PCDestroy_GAMG_GEO(PC pc) {
7669b8ffb57SJed Brown   PetscFunctionBegin;
7679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
7689b8ffb57SJed Brown   PetscFunctionReturn(0);
7699b8ffb57SJed Brown }
7709b8ffb57SJed Brown 
771c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
772c8b0795cSMark F. Adams /*
773c8b0795cSMark F. Adams  PCCreateGAMG_GEO
774c8b0795cSMark F. Adams 
775c8b0795cSMark F. Adams   Input Parameter:
776c8b0795cSMark F. Adams    . pc -
777c8b0795cSMark F. Adams */
7789371c9d4SSatish Balay PetscErrorCode PCCreateGAMG_GEO(PC pc) {
779c8b0795cSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
780c8b0795cSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
781c8b0795cSMark F. Adams 
782c8b0795cSMark F. Adams   PetscFunctionBegin;
7831ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO;
7849b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_GEO;
785c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
786c8b0795cSMark F. Adams 
787c8b0795cSMark F. Adams   /* set internal function pointers */
788fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_GEO;
789fd1112cbSBarry Smith   pc_gamg->ops->coarsen           = PCGAMGCoarsen_GEO;
7901ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_GEO;
791fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = NULL;
7921ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_GEO;
793c8b0795cSMark F. Adams 
7949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_GEO));
7952e68589bSMark F. Adams   PetscFunctionReturn(0);
7962e68589bSMark F. Adams }
797