xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision 792fecdfe9134cce4d631112660ddd34f063bc17)
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 
229fbee547SJacob Faibussowitsch static inline int petsc_geo_mg_compare(const void *a, const void *b)
23c8b0795cSMark F. Adams {
24c8b0795cSMark F. Adams   return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree);
25c8b0795cSMark F. Adams }
26c8b0795cSMark F. Adams 
272e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
282e68589bSMark F. Adams /*
292e68589bSMark F. Adams    PCSetCoordinates_GEO
302e68589bSMark F. Adams 
312e68589bSMark F. Adams    Input Parameter:
322e68589bSMark F. Adams    .  pc - the preconditioner context
332e68589bSMark F. Adams */
34302f38e8SMark F. Adams PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
352e68589bSMark F. Adams {
362e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
372e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
3890fbc344SStefano Zampini   PetscInt       arrsz,bs,my0,kk,ii,nloc,Iend,aloc;
392e68589bSMark F. Adams   Mat            Amat = pc->pmat;
402e68589bSMark F. Adams 
412e68589bSMark F. Adams   PetscFunctionBegin;
422e68589bSMark F. Adams   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
439566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
449566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &my0, &Iend));
4590fbc344SStefano Zampini   aloc = (Iend-my0);
462e68589bSMark F. Adams   nloc = (Iend-my0)/bs;
47a2f3521dSMark F. Adams 
482472a847SBarry 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);
492e68589bSMark F. Adams 
50c8b0795cSMark F. Adams   pc_gamg->data_cell_rows = 1;
517827d75bSBarry Smith   PetscCheck(coords || (nloc <= 0),PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
52c8b0795cSMark F. Adams   pc_gamg->data_cell_cols = ndm; /* coordinates */
532e68589bSMark F. Adams 
54c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
552e68589bSMark F. Adams 
562e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
577f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
589566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz+1, &pc_gamg->data));
602e68589bSMark F. Adams   }
612e68589bSMark F. Adams   for (kk=0; kk<arrsz; kk++) pc_gamg->data[kk] = -999.;
622e68589bSMark F. Adams   pc_gamg->data[arrsz] = -99.;
632e68589bSMark F. Adams   /* copy data in - column oriented */
6490fbc344SStefano Zampini   if (nloc == a_nloc) {
652e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) {
662e68589bSMark F. Adams       for (ii = 0; ii < ndm; ii++) {
672e68589bSMark F. Adams         pc_gamg->data[ii*nloc + kk] =  coords[kk*ndm + ii];
682e68589bSMark F. Adams       }
692e68589bSMark F. Adams     }
7090fbc344SStefano Zampini   } else { /* assumes the coordinates are blocked */
7190fbc344SStefano Zampini     for (kk = 0; kk < nloc; kk++) {
7290fbc344SStefano Zampini       for (ii = 0; ii < ndm; ii++) {
7390fbc344SStefano Zampini         pc_gamg->data[ii*nloc + kk] =  coords[bs*kk*ndm + ii];
7490fbc344SStefano Zampini       }
7590fbc344SStefano Zampini     }
7690fbc344SStefano Zampini   }
7763a3b9bcSJacob 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]);
782e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
792e68589bSMark F. Adams   PetscFunctionReturn(0);
802e68589bSMark F. Adams }
812e68589bSMark F. Adams 
822e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
832e68589bSMark F. Adams /*
842e68589bSMark F. Adams    PCSetData_GEO
852e68589bSMark F. Adams 
862e68589bSMark F. Adams   Input Parameter:
872e68589bSMark F. Adams    . pc -
882e68589bSMark F. Adams */
89b8cd405aSMark F. Adams PetscErrorCode PCSetData_GEO(PC pc, Mat m)
902e68589bSMark F. Adams {
912e68589bSMark F. Adams   PetscFunctionBegin;
92ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"GEO MG needs coordinates");
932e68589bSMark F. Adams }
942e68589bSMark F. Adams 
952e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
962e68589bSMark F. Adams /*
972e68589bSMark F. Adams    PCSetFromOptions_GEO
982e68589bSMark F. Adams 
992e68589bSMark F. Adams   Input Parameter:
1002e68589bSMark F. Adams    . pc -
1012e68589bSMark F. Adams */
1024416b707SBarry Smith PetscErrorCode PCSetFromOptions_GEO(PetscOptionItems *PetscOptionsObject,PC pc)
1032e68589bSMark F. Adams {
1042e68589bSMark F. Adams   PetscFunctionBegin;
105d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"GAMG-GEO options");
1062e68589bSMark F. Adams   {
1072e68589bSMark F. Adams     /* -pc_gamg_sa_nsmooths */
1082e68589bSMark F. Adams     /* pc_gamg_sa->smooths = 0; */
1092e68589bSMark F. Adams     /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */
1102e68589bSMark F. Adams     /*                        "smoothing steps for smoothed aggregation, usually 1 (0)", */
1112e68589bSMark F. Adams     /*                        "PCGAMGSetNSmooths_AGG", */
1122e68589bSMark F. Adams     /*                        pc_gamg_sa->smooths, */
1132e68589bSMark F. Adams     /*                        &pc_gamg_sa->smooths, */
1142e68589bSMark F. Adams     /*                        &flag);  */
1152e68589bSMark F. Adams   }
116d0609cedSBarry Smith   PetscOptionsHeadEnd();
1172e68589bSMark F. Adams   PetscFunctionReturn(0);
1182e68589bSMark F. Adams }
1192e68589bSMark F. Adams 
1202e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1212e68589bSMark F. Adams /*
1222e68589bSMark F. Adams  triangulateAndFormProl
1232e68589bSMark F. Adams 
1242e68589bSMark F. Adams    Input Parameter:
1252e68589bSMark F. Adams    . selected_2 - list of selected local ID, includes selected ghosts
126a2f3521dSMark F. Adams    . data_stride -
127a2f3521dSMark F. Adams    . coords[2*data_stride] - column vector of local coordinates w/ ghosts
1282adfe9d3SBarry Smith    . nselected_1 - selected IDs that go with base (1) graph includes selected ghosts
1290cbbd2e1SMark F. Adams    . clid_lid_1[nselected_1] - lids of selected (c) nodes   ???????????
1302adfe9d3SBarry Smith    . agg_lists_1 - list of aggregates selected_1 vertices of aggregate unselected vertices
1312e68589bSMark F. Adams    . crsGID[selected.size()] - global index for prolongation operator
1322e68589bSMark F. Adams    . bs - block size
1332e68589bSMark F. Adams   Output Parameter:
1342e68589bSMark F. Adams    . a_Prol - prolongation operator
1352e68589bSMark F. Adams    . a_worst_best - measure of worst missed fine vertex, 0 is no misses
1362e68589bSMark F. Adams */
1372adfe9d3SBarry Smith static PetscErrorCode triangulateAndFormProl(IS selected_2,PetscInt data_stride,PetscReal coords[],PetscInt nselected_1,const PetscInt clid_lid_1[],const PetscCoarsenData *agg_lists_1,
1382adfe9d3SBarry Smith                                              const PetscInt crsGID[],PetscInt bs,Mat a_Prol,PetscReal *a_worst_best)
1392e68589bSMark F. Adams {
1402e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE)
141b43b03e9SMark F. Adams   PetscInt             jj,tid,tt,idx,nselected_2;
1422e68589bSMark F. Adams   struct triangulateio in,mid;
1430cbbd2e1SMark F. Adams   const PetscInt       *selected_idx_2;
14473911c69SBarry Smith   PetscMPIInt          rank;
1452e68589bSMark F. Adams   PetscInt             Istart,Iend,nFineLoc,myFine0;
1462e68589bSMark F. Adams   int                  kk,nPlotPts,sid;
1473b4367a7SBarry Smith   MPI_Comm             comm;
148c8b0795cSMark F. Adams   PetscReal            tm;
149c8b0795cSMark F. Adams 
1506e111a19SKarl Rupp   PetscFunctionBegin;
1519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol,&comm));
1529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
1539566063dSJacob Faibussowitsch   PetscCall(ISGetSize(selected_2, &nselected_2));
1542e68589bSMark F. Adams   if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */
1552e68589bSMark F. Adams     *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */
156806fa848SBarry Smith   } else *a_worst_best = 0.0;
1571c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm));
158c8b0795cSMark F. Adams   if (tm > 0.0) {
159c8b0795cSMark F. Adams     *a_worst_best = 100.0;
1602e68589bSMark F. Adams     PetscFunctionReturn(0);
1612e68589bSMark F. Adams   }
1629566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
1632e68589bSMark F. Adams   nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs;
1642e68589bSMark F. Adams   nPlotPts = nFineLoc; /* locals */
1654c500f23SPierre Jolivet   /* triangle */
1662e68589bSMark F. Adams   /* Define input points - in*/
1672e68589bSMark F. Adams   in.numberofpoints          = nselected_2;
1682e68589bSMark F. Adams   in.numberofpointattributes = 0;
1692e68589bSMark F. Adams   /* get nselected points */
1709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*nselected_2, &in.pointlist));
1719566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(selected_2, &selected_idx_2));
1722e68589bSMark F. Adams 
1732e68589bSMark F. Adams   for (kk=0,sid=0; kk<nselected_2; kk++,sid += 2) {
1742e68589bSMark F. Adams     PetscInt lid = selected_idx_2[kk];
1752e68589bSMark F. Adams     in.pointlist[sid]   = coords[lid];
176a2f3521dSMark F. Adams     in.pointlist[sid+1] = coords[data_stride + lid];
1772e68589bSMark F. Adams     if (lid>=nFineLoc) nPlotPts++;
1782e68589bSMark F. Adams   }
17963a3b9bcSJacob Faibussowitsch   PetscCheck(sid == 2*nselected_2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %" PetscInt_FMT " != 2*nselected_2 %" PetscInt_FMT,sid,nselected_2);
1802e68589bSMark F. Adams 
1812e68589bSMark F. Adams   in.numberofsegments      = 0;
1822e68589bSMark F. Adams   in.numberofedges         = 0;
1832e68589bSMark F. Adams   in.numberofholes         = 0;
1842e68589bSMark F. Adams   in.numberofregions       = 0;
1850a545947SLisandro Dalcin   in.trianglelist          = NULL;
1860a545947SLisandro Dalcin   in.segmentmarkerlist     = NULL;
1870a545947SLisandro Dalcin   in.pointattributelist    = NULL;
1880a545947SLisandro Dalcin   in.pointmarkerlist       = NULL;
1890a545947SLisandro Dalcin   in.triangleattributelist = NULL;
1900a545947SLisandro Dalcin   in.trianglearealist      = NULL;
1910a545947SLisandro Dalcin   in.segmentlist           = NULL;
1920a545947SLisandro Dalcin   in.holelist              = NULL;
1930a545947SLisandro Dalcin   in.regionlist            = NULL;
1940a545947SLisandro Dalcin   in.edgelist              = NULL;
1950a545947SLisandro Dalcin   in.edgemarkerlist        = NULL;
1960a545947SLisandro Dalcin   in.normlist              = NULL;
1972fa5cd67SKarl Rupp 
1982e68589bSMark F. Adams   /* triangulate */
1990a545947SLisandro Dalcin   mid.pointlist = NULL;          /* Not needed if -N switch used. */
2002e68589bSMark F. Adams   /* Not needed if -N switch used or number of point attributes is zero: */
2010a545947SLisandro Dalcin   mid.pointattributelist = NULL;
2020a545947SLisandro Dalcin   mid.pointmarkerlist    = NULL; /* Not needed if -N or -B switch used. */
2030a545947SLisandro Dalcin   mid.trianglelist       = NULL; /* Not needed if -E switch used. */
2042e68589bSMark F. Adams   /* Not needed if -E switch used or number of triangle attributes is zero: */
2050a545947SLisandro Dalcin   mid.triangleattributelist = NULL;
2060a545947SLisandro Dalcin   mid.neighborlist          = NULL; /* Needed only if -n switch used. */
2072e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P not used: */
2080a545947SLisandro Dalcin   mid.segmentlist = NULL;
2092e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
2100a545947SLisandro Dalcin   mid.segmentmarkerlist = NULL;
2110a545947SLisandro Dalcin   mid.edgelist          = NULL; /* Needed only if -e switch used. */
2120a545947SLisandro Dalcin   mid.edgemarkerlist    = NULL; /* Needed if -e used and -B not used. */
2132e68589bSMark F. Adams   mid.numberoftriangles = 0;
2142e68589bSMark F. Adams 
2152e68589bSMark F. Adams   /* Triangulate the points.  Switches are chosen to read and write a  */
2162e68589bSMark F. Adams   /*   PSLG (p), preserve the convex hull (c), number everything from  */
2172e68589bSMark F. Adams   /*   zero (z), assign a regional attribute to each element (A), and  */
2182e68589bSMark F. Adams   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
2192e68589bSMark F. Adams   /*   neighbor list (n).                                            */
2202e68589bSMark F. Adams   if (nselected_2 != 0) { /* inactive processor */
2212e68589bSMark F. Adams     char args[] = "npczQ"; /* c is needed ? */
2222e68589bSMark F. Adams     triangulate(args, &in, &mid, (struct triangulateio*) NULL);
2232e68589bSMark F. Adams     /* output .poly files for 'showme' */
2242e68589bSMark F. Adams     if (!PETSC_TRUE) {
2252e68589bSMark F. Adams       static int level = 1;
2262e68589bSMark F. Adams       FILE       *file; char fname[32];
2272e68589bSMark F. Adams 
228c5df96a5SBarry Smith       sprintf(fname,"C%d_%d.poly",level,rank); file = fopen(fname, "w");
2292e68589bSMark F. Adams       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
2302e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0);
2312e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2322e68589bSMark F. Adams       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid += 2) {
2332e68589bSMark F. Adams         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
2342e68589bSMark F. Adams       }
2352e68589bSMark F. Adams       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
2362e68589bSMark F. Adams       fprintf(file, "%d  %d\n",0,0);
2372e68589bSMark F. Adams       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
2382e68589bSMark F. Adams       /* One line: <# of holes> */
2392e68589bSMark F. Adams       fprintf(file, "%d\n",0);
2402e68589bSMark F. Adams       /* Following lines: <hole #> <x> <y> */
2412e68589bSMark F. Adams       /* Optional line: <# of regional attributes and/or area constraints> */
2422e68589bSMark F. Adams       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
2432e68589bSMark F. Adams       fclose(file);
2442e68589bSMark F. Adams 
2452e68589bSMark F. Adams       /* elems */
246c5df96a5SBarry Smith       sprintf(fname,"C%d_%d.ele",level,rank); file = fopen(fname, "w");
2472e68589bSMark F. Adams       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
2482e68589bSMark F. Adams       fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0);
2492e68589bSMark F. Adams       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
2502e68589bSMark F. Adams       for (kk=0,sid=0; kk<mid.numberoftriangles; kk++,sid += 3) {
2512e68589bSMark F. Adams         fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]);
2522e68589bSMark F. Adams       }
2532e68589bSMark F. Adams       fclose(file);
2542e68589bSMark F. Adams 
255c5df96a5SBarry Smith       sprintf(fname,"C%d_%d.node",level,rank); file = fopen(fname, "w");
2562e68589bSMark F. Adams       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
2572e68589bSMark F. Adams       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
2582e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n",nPlotPts,2,0,0);
2592e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2602e68589bSMark F. Adams       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid+=2) {
2612e68589bSMark F. Adams         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
2622e68589bSMark F. Adams       }
2632e68589bSMark F. Adams 
2642e68589bSMark F. Adams       sid /= 2;
2652e68589bSMark F. Adams       for (jj=0; jj<nFineLoc; jj++) {
2662e68589bSMark F. Adams         PetscBool sel = PETSC_TRUE;
2672e68589bSMark F. Adams         for (kk=0; kk<nselected_2 && sel; kk++) {
2682e68589bSMark F. Adams           PetscInt lid = selected_idx_2[kk];
2692e68589bSMark F. Adams           if (lid == jj) sel = PETSC_FALSE;
2702e68589bSMark F. Adams         }
2712fa5cd67SKarl Rupp         if (sel) fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]);
2722e68589bSMark F. Adams       }
2732e68589bSMark F. Adams       fclose(file);
27463a3b9bcSJacob Faibussowitsch       PetscCheck(sid == nPlotPts,PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %" PetscInt_FMT " != nPlotPts %" PetscInt_FMT,sid,nPlotPts);
2752e68589bSMark F. Adams       level++;
2762e68589bSMark F. Adams     }
2772e68589bSMark F. Adams   }
2782e68589bSMark F. Adams   { /* form P - setup some maps */
2790cbbd2e1SMark F. Adams     PetscInt clid,mm,*nTri,*node_tri;
2802e68589bSMark F. Adams 
2819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nselected_2, &node_tri,nselected_2, &nTri));
2822e68589bSMark F. Adams 
2832e68589bSMark F. Adams     /* need list of triangles on node */
2842e68589bSMark F. Adams     for (kk=0; kk<nselected_2; kk++) nTri[kk] = 0;
2852e68589bSMark F. Adams     for (tid=0,kk=0; tid<mid.numberoftriangles; tid++) {
2862e68589bSMark F. Adams       for (jj=0; jj<3; jj++) {
2872e68589bSMark F. Adams         PetscInt cid = mid.trianglelist[kk++];
2882e68589bSMark F. Adams         if (nTri[cid] == 0) node_tri[cid] = tid;
2892e68589bSMark F. Adams         nTri[cid]++;
2902e68589bSMark F. Adams       }
2912e68589bSMark F. Adams     }
2922e68589bSMark F. Adams #define EPS 1.e-12
2932e68589bSMark F. Adams     /* find points and set prolongation */
2940cbbd2e1SMark F. Adams     for (mm = clid = 0; mm < nFineLoc; mm++) {
295e78576d6SMark F. Adams       PetscBool ise;
2969566063dSJacob Faibussowitsch       PetscCall(PetscCDEmptyAt(agg_lists_1,mm,&ise));
297e78576d6SMark F. Adams       if (!ise) {
2980cbbd2e1SMark F. Adams         const PetscInt lid = mm;
2992e68589bSMark F. Adams         PetscScalar    AA[3][3];
3002e68589bSMark F. Adams         PetscBLASInt   N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO;
301539c167fSBarry Smith         PetscCDIntNd   *pos;
302539c167fSBarry Smith 
3039566063dSJacob Faibussowitsch         PetscCall(PetscCDGetHeadPos(agg_lists_1,lid,&pos));
304e78576d6SMark F. Adams         while (pos) {
305e78576d6SMark F. Adams           PetscInt flid;
3069566063dSJacob Faibussowitsch           PetscCall(PetscCDIntNdGetID(pos, &flid));
3079566063dSJacob Faibussowitsch           PetscCall(PetscCDGetNextPos(agg_lists_1,lid,&pos));
3080cbbd2e1SMark F. Adams 
3092e68589bSMark F. Adams           if (flid < nFineLoc) {  /* could be a ghost */
3102e68589bSMark F. Adams             PetscInt       bestTID = -1; PetscReal best_alpha = 1.e10;
3112e68589bSMark F. Adams             const PetscInt fgid    = flid + myFine0;
3122e68589bSMark F. Adams             /* compute shape function for gid */
313a2f3521dSMark F. Adams             const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0};
3142e68589bSMark F. Adams             PetscBool       haveit    =PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3];
3152fa5cd67SKarl Rupp 
3162e68589bSMark F. Adams             /* look for it */
3170cbbd2e1SMark F. Adams             for (tid = node_tri[clid], jj=0;
3182e68589bSMark F. Adams                  jj < 5 && !haveit && tid != -1;
3192e68589bSMark F. Adams                  jj++) {
3202e68589bSMark F. Adams               for (tt=0; tt<3; tt++) {
3212e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3*tid + tt];
3222e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
323a2f3521dSMark F. Adams                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
3242e68589bSMark F. Adams                 clids[tt] = cid2; /* store for interp */
3252e68589bSMark F. Adams               }
3262e68589bSMark F. Adams 
3272e68589bSMark F. Adams               for (tt=0; tt<3; tt++) alpha[tt] = (PetscScalar)fcoord[tt];
3282e68589bSMark F. Adams 
3292e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
330*792fecdfSBarry Smith               PetscCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3312e68589bSMark F. Adams               {
3322e68589bSMark F. Adams                 PetscBool have=PETSC_TRUE;  PetscReal lowest=1.e10;
3332e68589bSMark F. Adams                 for (tt = 0, idx = 0; tt < 3; tt++) {
3342e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
3352e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) < lowest) {
3362e68589bSMark F. Adams                     lowest = PetscRealPart(alpha[tt]);
3372e68589bSMark F. Adams                     idx    = tt;
3382e68589bSMark F. Adams                   }
3392e68589bSMark F. Adams                 }
3402e68589bSMark F. Adams                 haveit = have;
3412e68589bSMark F. Adams               }
3422e68589bSMark F. Adams               tid = mid.neighborlist[3*tid + idx];
3432e68589bSMark F. Adams             }
3442e68589bSMark F. Adams 
3452e68589bSMark F. Adams             if (!haveit) {
3462e68589bSMark F. Adams               /* brute force */
3472e68589bSMark F. Adams               for (tid=0; tid<mid.numberoftriangles && !haveit; tid++) {
3482e68589bSMark F. Adams                 for (tt=0; tt<3; tt++) {
3492e68589bSMark F. Adams                   PetscInt cid2 = mid.trianglelist[3*tid + tt];
3502e68589bSMark F. Adams                   PetscInt lid2 = selected_idx_2[cid2];
351a2f3521dSMark F. Adams                   AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
3522e68589bSMark F. Adams                   clids[tt] = cid2; /* store for interp */
3532e68589bSMark F. Adams                 }
3542e68589bSMark F. Adams                 for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
3552e68589bSMark F. Adams                 /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
356*792fecdfSBarry Smith                 PetscCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3572e68589bSMark F. Adams                 {
3582e68589bSMark F. Adams                   PetscBool have=PETSC_TRUE;  PetscReal worst=0.0, v;
3592e68589bSMark F. Adams                   for (tt=0; tt<3 && have; tt++) {
3602e68589bSMark F. Adams                     if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE;
3612e68589bSMark F. Adams                     if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v;
3622e68589bSMark F. Adams                   }
3632e68589bSMark F. Adams                   if (worst < best_alpha) {
3642e68589bSMark F. Adams                     best_alpha = worst; bestTID = tid;
3652e68589bSMark F. Adams                   }
3662e68589bSMark F. Adams                   haveit = have;
3672e68589bSMark F. Adams                 }
3682e68589bSMark F. Adams               }
3692e68589bSMark F. Adams             }
3702e68589bSMark F. Adams             if (!haveit) {
3712e68589bSMark F. Adams               if (best_alpha > *a_worst_best) *a_worst_best = best_alpha;
3722e68589bSMark F. Adams               /* use best one */
3732e68589bSMark F. Adams               for (tt=0; tt<3; tt++) {
3742e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3*bestTID + tt];
3752e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
376a2f3521dSMark F. Adams                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
3772e68589bSMark F. Adams                 clids[tt] = cid2; /* store for interp */
3782e68589bSMark F. Adams               }
3792e68589bSMark F. Adams               for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
3802e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
381*792fecdfSBarry Smith               PetscCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3822e68589bSMark F. Adams             }
3832e68589bSMark F. Adams 
3842e68589bSMark F. Adams             /* put in row of P */
3852e68589bSMark F. Adams             for (idx=0; idx<3; idx++) {
3862e68589bSMark F. Adams               PetscScalar shp = alpha[idx];
3872e68589bSMark F. Adams               if (PetscAbs(PetscRealPart(shp)) > 1.e-6) {
3882e68589bSMark F. Adams                 PetscInt cgid = crsGID[clids[idx]];
3892e68589bSMark F. Adams                 PetscInt jj   = cgid*bs, ii = fgid*bs; /* need to gloalize */
3902e68589bSMark F. Adams                 for (tt=0; tt < bs; tt++, ii++, jj++) {
3919566063dSJacob Faibussowitsch                   PetscCall(MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES));
3922e68589bSMark F. Adams                 }
3932e68589bSMark F. Adams               }
3942e68589bSMark F. Adams             }
3952e68589bSMark F. Adams           }
3960cbbd2e1SMark F. Adams         } /* aggregates iterations */
3970cbbd2e1SMark F. Adams         clid++;
3980cbbd2e1SMark F. Adams       } /* a coarse agg */
3990cbbd2e1SMark F. Adams     } /* for all fine nodes */
4000cbbd2e1SMark F. Adams 
4019566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(selected_2, &selected_idx_2));
4029566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY));
4039566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY));
4042e68589bSMark F. Adams 
4059566063dSJacob Faibussowitsch     PetscCall(PetscFree2(node_tri,nTri));
4062e68589bSMark F. Adams   }
4072e68589bSMark F. Adams   free(mid.trianglelist);
4082e68589bSMark F. Adams   free(mid.neighborlist);
409c6b50ea1SMark   free(mid.segmentlist);
410c6b50ea1SMark   free(mid.segmentmarkerlist);
411c6b50ea1SMark   free(mid.pointlist);
412c6b50ea1SMark   free(mid.pointmarkerlist);
4139566063dSJacob Faibussowitsch   PetscCall(PetscFree(in.pointlist));
4142e68589bSMark F. Adams   PetscFunctionReturn(0);
4152e68589bSMark F. Adams #else
4163b4367a7SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)a_Prol),PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG");
4172e68589bSMark F. Adams #endif
4182e68589bSMark F. Adams }
4192e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
4202e68589bSMark F. Adams /*
4212e68589bSMark F. Adams    getGIDsOnSquareGraph - square graph, get
4222e68589bSMark F. Adams 
4232e68589bSMark F. Adams    Input Parameter:
4240cbbd2e1SMark F. Adams    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
425b43b03e9SMark F. Adams    . clid_lid_1 - [nselected_1] lids of selected nodes
4262e68589bSMark F. Adams    . Gmat1 - graph that goes with 'selected_1'
4272e68589bSMark F. Adams    Output Parameter:
4282e68589bSMark F. Adams    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
4292e68589bSMark F. Adams    . a_Gmat_2 - graph that is squared of 'Gmat_1'
4302e68589bSMark F. Adams    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
4312e68589bSMark F. Adams */
4324b1575e2SStefano Zampini 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)
4332e68589bSMark F. Adams {
43473911c69SBarry Smith   PetscMPIInt    size;
435b43b03e9SMark F. Adams   PetscInt       *crsGID, kk,my0,Iend,nloc;
4363b4367a7SBarry Smith   MPI_Comm       comm;
4372e68589bSMark F. Adams 
4382e68589bSMark F. Adams   PetscFunctionBegin;
4399566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Gmat1,&comm));
4409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
4419566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1,&my0,&Iend)); /* AIJ */
4422e68589bSMark F. Adams   nloc = Iend - my0; /* this does not change */
4432e68589bSMark F. Adams 
444c5df96a5SBarry Smith   if (size == 1) { /* not much to do in serial */
4459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected_1, &crsGID));
446b43b03e9SMark F. Adams     for (kk=0; kk<nselected_1; kk++) crsGID[kk] = kk;
4470a545947SLisandro Dalcin     *a_Gmat_2 = NULL;
4489566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2));
449806fa848SBarry Smith   } else {
450b43b03e9SMark F. Adams     PetscInt    idx,num_fine_ghosts,num_crs_ghost,myCrs0;
4512e68589bSMark F. Adams     Mat_MPIAIJ  *mpimat2;
4522e68589bSMark F. Adams     Mat         Gmat2;
4532e68589bSMark F. Adams     Vec         locState;
4542e68589bSMark F. Adams     PetscScalar *cpcol_state;
4552e68589bSMark F. Adams 
4562e68589bSMark F. Adams     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
457b43b03e9SMark F. Adams     kk = nselected_1;
4589566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPI_SUM, comm));
459b43b03e9SMark F. Adams     myCrs0 -= nselected_1;
4602e68589bSMark F. Adams 
461b43b03e9SMark F. Adams     if (a_Gmat_2) { /* output */
4622e68589bSMark F. Adams       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
4639566063dSJacob Faibussowitsch       PetscCall(PCGAMGSquareGraph_GAMG(pc,Gmat1,&Gmat2));
4642e68589bSMark F. Adams       *a_Gmat_2 = Gmat2; /* output */
465806fa848SBarry Smith     } else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
4662e68589bSMark F. Adams     /* get coarse grid GIDS for selected (locals and ghosts) */
4672e68589bSMark F. Adams     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
4689566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Gmat2, &locState, NULL));
4699566063dSJacob Faibussowitsch     PetscCall(VecSet(locState, (PetscScalar)(PetscReal)(-1))); /* set with UNKNOWN state */
470b43b03e9SMark F. Adams     for (kk=0; kk<nselected_1; kk++) {
471b43b03e9SMark F. Adams       PetscInt    fgid = clid_lid_1[kk] + my0;
4722e68589bSMark F. Adams       PetscScalar v    = (PetscScalar)(kk+myCrs0);
4739566063dSJacob Faibussowitsch       PetscCall(VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES)); /* set with PID */
4742e68589bSMark F. Adams     }
4759566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(locState));
4769566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(locState));
4779566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD));
4789566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD));
4799566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts));
4809566063dSJacob Faibussowitsch     PetscCall(VecGetArray(mpimat2->lvec, &cpcol_state));
4812e68589bSMark F. Adams     for (kk=0,num_crs_ghost=0; kk<num_fine_ghosts; kk++) {
4822e68589bSMark F. Adams       if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++;
4832e68589bSMark F. Adams     }
4849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected_1+num_crs_ghost, &crsGID)); /* output */
4852e68589bSMark F. Adams     {
4862e68589bSMark F. Adams       PetscInt *selected_set;
4879566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nselected_1+num_crs_ghost, &selected_set));
4882e68589bSMark F. Adams       /* do ghost of 'crsGID' */
489b43b03e9SMark F. Adams       for (kk=0,idx=nselected_1; kk<num_fine_ghosts; kk++) {
4902e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
4912e68589bSMark F. Adams           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
4922e68589bSMark F. Adams           selected_set[idx] = nloc + kk;
4932e68589bSMark F. Adams           crsGID[idx++]     = cgid;
4942e68589bSMark F. Adams         }
4952e68589bSMark F. Adams       }
49663a3b9bcSJacob 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);
4979566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(mpimat2->lvec, &cpcol_state));
4982e68589bSMark F. Adams       /* do locals in 'crsGID' */
4999566063dSJacob Faibussowitsch       PetscCall(VecGetArray(locState, &cpcol_state));
5002e68589bSMark F. Adams       for (kk=0,idx=0; kk<nloc; kk++) {
5012e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
5022e68589bSMark F. Adams           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
5032e68589bSMark F. Adams           selected_set[idx] = kk;
5042e68589bSMark F. Adams           crsGID[idx++]     = cgid;
5052e68589bSMark F. Adams         }
5062e68589bSMark F. Adams       }
50763a3b9bcSJacob Faibussowitsch       PetscCheck(idx == nselected_1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %" PetscInt_FMT " != nselected_1 %" PetscInt_FMT,idx,nselected_1);
5089566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(locState, &cpcol_state));
5092e68589bSMark F. Adams 
5100a545947SLisandro Dalcin       if (a_selected_2 != NULL) { /* output */
5119566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2));
512806fa848SBarry Smith       } else {
5139566063dSJacob Faibussowitsch         PetscCall(PetscFree(selected_set));
5142e68589bSMark F. Adams       }
5150cbbd2e1SMark F. Adams     }
5169566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&locState));
5172e68589bSMark F. Adams   }
5182e68589bSMark F. Adams   *a_crsGID = crsGID; /* output */
5192e68589bSMark F. Adams   PetscFunctionReturn(0);
5202e68589bSMark F. Adams }
5212e68589bSMark F. Adams 
5222e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
5232e68589bSMark F. Adams /*
524fd1112cbSBarry Smith    PCGAMGGraph_GEO
5252e68589bSMark F. Adams 
5262e68589bSMark F. Adams   Input Parameter:
5272e68589bSMark F. Adams    . pc - this
5282e68589bSMark F. Adams    . Amat - matrix on this fine level
5292e68589bSMark F. Adams   Output Parameter:
530c8b0795cSMark F. Adams    . a_Gmat
5312e68589bSMark F. Adams */
5322adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_GEO(PC pc,Mat Amat,Mat *a_Gmat)
533c8b0795cSMark F. Adams {
534c8b0795cSMark F. Adams   PC_MG           *mg      = (PC_MG*)pc->data;
535c8b0795cSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
536c1eae691SMark Adams   const PetscReal vfilter  = pc_gamg->threshold[0];
5373b4367a7SBarry Smith   MPI_Comm        comm;
53872833a62Smarkadams4   Mat             Gmat,F=NULL;
5390cbbd2e1SMark F. Adams   PetscBool       set,flg,symm;
5406e111a19SKarl Rupp 
541c8b0795cSMark F. Adams   PetscFunctionBegin;
5429566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm));
543c8b0795cSMark F. Adams 
5449566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetricKnown(Amat, &set, &flg));
545263489e9SJed Brown   symm = (PetscBool)!(set && flg);
5460cbbd2e1SMark F. Adams 
54772833a62Smarkadams4   PetscCall(MatCreateGraph(Amat, symm, PETSC_TRUE, &Gmat));
54872833a62Smarkadams4   PetscCall(MatFilter(Gmat, vfilter, &F));
54972833a62Smarkadams4   if (F) {
55072833a62Smarkadams4     PetscCall(MatDestroy(&Gmat));
55172833a62Smarkadams4     Gmat = F;
55272833a62Smarkadams4   }
553c8b0795cSMark F. Adams   *a_Gmat = Gmat;
554849bee69SMark Adams 
555c8b0795cSMark F. Adams   PetscFunctionReturn(0);
556c8b0795cSMark F. Adams }
557c8b0795cSMark F. Adams 
558c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
559c8b0795cSMark F. Adams /*
560fd1112cbSBarry Smith    PCGAMGCoarsen_GEO
561c8b0795cSMark F. Adams 
562c8b0795cSMark F. Adams   Input Parameter:
563e0940f08SMark F. Adams    . a_pc - this
564e0940f08SMark F. Adams    . a_Gmat - graph
565c8b0795cSMark F. Adams   Output Parameter:
566c8b0795cSMark F. Adams    . a_llist_parent - linked list from selected indices for data locality only
567c8b0795cSMark F. Adams */
568fd1112cbSBarry Smith PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc,Mat *a_Gmat,PetscCoarsenData **a_llist_parent)
569c8b0795cSMark F. Adams {
570c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,kk,Ii,ncols;
5710cbbd2e1SMark F. Adams   IS             perm;
572c8b0795cSMark F. Adams   GAMGNode       *gnodes;
573c8b0795cSMark F. Adams   PetscInt       *permute;
574e0940f08SMark F. Adams   Mat            Gmat  = *a_Gmat;
5753b4367a7SBarry Smith   MPI_Comm       comm;
576b43b03e9SMark F. Adams   MatCoarsen     crs;
577c8b0795cSMark F. Adams 
578c8b0795cSMark F. Adams   PetscFunctionBegin;
5799566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_pc,&comm));
580849bee69SMark Adams 
5819566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend));
582c8b0795cSMark F. Adams   nloc = (Iend-Istart);
583c8b0795cSMark F. Adams 
584c8b0795cSMark F. Adams   /* create random permutation with sort for geo-mg */
5859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc, &gnodes));
5869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc, &permute));
587c8b0795cSMark F. Adams 
588c8b0795cSMark F. Adams   for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */
5899566063dSJacob Faibussowitsch     PetscCall(MatGetRow(Gmat,Ii,&ncols,NULL,NULL));
590c8b0795cSMark F. Adams     {
591c8b0795cSMark F. Adams       PetscInt lid = Ii - Istart;
592c8b0795cSMark F. Adams       gnodes[lid].lid    = lid;
593c8b0795cSMark F. Adams       gnodes[lid].degree = ncols;
594c8b0795cSMark F. Adams     }
5959566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(Gmat,Ii,&ncols,NULL,NULL));
596c8b0795cSMark F. Adams   }
597c8b0795cSMark F. Adams   if (PETSC_TRUE) {
598e2c00dcbSBarry Smith     PetscRandom  rand;
599c8b0795cSMark F. Adams     PetscBool    *bIndexSet;
600e2c00dcbSBarry Smith     PetscReal    rr;
601e2c00dcbSBarry Smith     PetscInt     iSwapIndex;
602e2c00dcbSBarry Smith 
6039566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm,&rand));
6049566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nloc, &bIndexSet));
6052fa5cd67SKarl Rupp     for (Ii = 0; Ii < nloc; Ii++) {
6069566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(rand,&rr));
607e2c00dcbSBarry Smith       iSwapIndex = (PetscInt) (rr*nloc);
6082fa5cd67SKarl Rupp       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
609c8b0795cSMark F. Adams         GAMGNode iTemp = gnodes[iSwapIndex];
610c8b0795cSMark F. Adams         gnodes[iSwapIndex]    = gnodes[Ii];
611c8b0795cSMark F. Adams         gnodes[Ii]            = iTemp;
612c8b0795cSMark F. Adams         bIndexSet[Ii]         = PETSC_TRUE;
613c8b0795cSMark F. Adams         bIndexSet[iSwapIndex] = PETSC_TRUE;
614c8b0795cSMark F. Adams       }
615c8b0795cSMark F. Adams     }
6169566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
6179566063dSJacob Faibussowitsch     PetscCall(PetscFree(bIndexSet));
618c8b0795cSMark F. Adams   }
619c8b0795cSMark F. Adams   /* only sort locals */
6200cbbd2e1SMark F. Adams   qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare);
621c8b0795cSMark F. Adams   /* create IS of permutation */
6222fa5cd67SKarl Rupp   for (kk=0; kk<nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */
6239566063dSJacob Faibussowitsch   PetscCall(PetscFree(gnodes));
6249566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm));
625c8b0795cSMark F. Adams 
626c8b0795cSMark F. Adams   /* get MIS aggs */
627b43b03e9SMark F. Adams 
6289566063dSJacob Faibussowitsch   PetscCall(MatCoarsenCreate(comm, &crs));
6299566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetType(crs, MATCOARSENMIS));
6309566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetGreedyOrdering(crs, perm));
6319566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetAdjacency(crs, Gmat));
6329566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_FALSE));
6339566063dSJacob Faibussowitsch   PetscCall(MatCoarsenApply(crs));
6349566063dSJacob Faibussowitsch   PetscCall(MatCoarsenGetData(crs, a_llist_parent));
6359566063dSJacob Faibussowitsch   PetscCall(MatCoarsenDestroy(&crs));
636c8b0795cSMark F. Adams 
6379566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
638849bee69SMark Adams 
639c8b0795cSMark F. Adams   PetscFunctionReturn(0);
640c8b0795cSMark F. Adams }
641c8b0795cSMark F. Adams 
642c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
643c8b0795cSMark F. Adams /*
6440cbbd2e1SMark F. Adams  PCGAMGProlongator_GEO
645c8b0795cSMark F. Adams 
646c8b0795cSMark F. Adams  Input Parameter:
647c8b0795cSMark F. Adams  . pc - this
648c8b0795cSMark F. Adams  . Amat - matrix on this fine level
649c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
6500cbbd2e1SMark F. Adams  . selected_1 - [nselected]
6510cbbd2e1SMark F. Adams  . agg_lists - [nselected]
652c8b0795cSMark F. Adams  Output Parameter:
653c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
654c8b0795cSMark F. Adams  */
6552adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_GEO(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
6562e68589bSMark F. Adams {
6572e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
6582e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
659c8b0795cSMark F. Adams   const PetscInt dim      = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
660b43b03e9SMark F. Adams   PetscInt       Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid;
661c8b0795cSMark F. Adams   Mat            Prol;
662c5df96a5SBarry Smith   PetscMPIInt    rank, size;
6633b4367a7SBarry Smith   MPI_Comm       comm;
6640cbbd2e1SMark F. Adams   IS             selected_2,selected_1;
6652e68589bSMark F. Adams   const PetscInt *selected_idx;
666d9558ea9SBarry Smith   MatType        mtype;
6672e68589bSMark F. Adams 
6682e68589bSMark F. Adams   PetscFunctionBegin;
6699566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm));
670849bee69SMark Adams 
6719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
6729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
6739566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
6749566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
67571959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
67663a3b9bcSJacob Faibussowitsch   PetscCheck((Iend-Istart) % bs == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT ") %% bs %" PetscInt_FMT,Iend,Istart,bs);
6772e68589bSMark F. Adams 
6782e68589bSMark F. Adams   /* get 'nLocalSelected' */
6799566063dSJacob Faibussowitsch   PetscCall(PetscCDGetMIS(agg_lists, &selected_1));
6809566063dSJacob Faibussowitsch   PetscCall(ISGetSize(selected_1, &jj));
6819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(jj, &clid_flid));
6829566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(selected_1, &selected_idx));
683c8b0795cSMark F. Adams   for (kk=0,nLocalSelected=0; kk<jj; kk++) {
6842e68589bSMark F. Adams     PetscInt lid = selected_idx[kk];
685b43b03e9SMark F. Adams     if (lid<nloc) {
6869566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Gmat,lid+my0,&ncols,NULL,NULL));
6872fa5cd67SKarl Rupp       if (ncols>1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */
6889566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Gmat,lid+my0,&ncols,NULL,NULL));
689b43b03e9SMark F. Adams     }
6902e68589bSMark F. Adams   }
6919566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(selected_1, &selected_idx));
6929566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&selected_1)); /* this is selected_1 in serial */
6932e68589bSMark F. Adams 
694d9558ea9SBarry Smith   /* create prolongator  matrix */
6959566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat,&mtype));
6969566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
6979566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE));
6989566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(Prol, bs, bs));
6999566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
7009566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol,3*data_cols,NULL));
7019566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol,3*data_cols,NULL,3*data_cols,NULL));
7022e68589bSMark F. Adams 
703c8b0795cSMark F. Adams   /* can get all points "removed" - but not on geomg */
7049566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &jj));
7057f66b68fSMark Adams   if (!jj) {
7069566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"ERROE: no selected points on coarse grid\n"));
7079566063dSJacob Faibussowitsch     PetscCall(PetscFree(clid_flid));
7089566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
7090298fd71SBarry Smith     *a_P_out = NULL;  /* out */
7102e68589bSMark F. Adams     PetscFunctionReturn(0);
7112e68589bSMark F. Adams   }
7122e68589bSMark F. Adams 
7132e68589bSMark F. Adams   {
7142e68589bSMark F. Adams     PetscReal *coords;
715a2f3521dSMark F. Adams     PetscInt  data_stride;
7160298fd71SBarry Smith     PetscInt  *crsGID = NULL;
7172e68589bSMark F. Adams     Mat       Gmat2;
7182e68589bSMark F. Adams 
71963a3b9bcSJacob Faibussowitsch     PetscCheck(dim == data_cols,PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim %" PetscInt_FMT " != data_cols %" PetscInt_FMT,dim,data_cols);
7202e68589bSMark F. Adams     /* grow ghost data for better coarse grid cover of fine grid */
721a2f3521dSMark F. Adams     /* messy method, squares graph and gets some data */
7229566063dSJacob Faibussowitsch     PetscCall(getGIDsOnSquareGraph(pc, nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID));
7232e68589bSMark F. Adams     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
7242e68589bSMark F. Adams     /* create global vector of coorindates in 'coords' */
725c5df96a5SBarry Smith     if (size > 1) {
7269566063dSJacob Faibussowitsch       PetscCall(PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords));
727806fa848SBarry Smith     } else {
728c8b0795cSMark F. Adams       coords      = (PetscReal*)pc_gamg->data;
729a2f3521dSMark F. Adams       data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols;
7302e68589bSMark F. Adams     }
7319566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Gmat2));
7322e68589bSMark F. Adams 
7332e68589bSMark F. Adams     /* triangulate */
7342e68589bSMark F. Adams     if (dim == 2) {
735c8b0795cSMark F. Adams       PetscReal metric,tm;
7369566063dSJacob Faibussowitsch       PetscCall(triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric));
7379566063dSJacob Faibussowitsch       PetscCall(PetscFree(crsGID));
7382e68589bSMark F. Adams 
7392e68589bSMark F. Adams       /* clean up and create coordinates for coarse grid (output) */
7409566063dSJacob Faibussowitsch       if (size > 1) PetscCall(PetscFree(coords));
7412e68589bSMark F. Adams 
7421c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm));
743c8b0795cSMark F. Adams       if (tm > 1.) { /* needs to be globalized - should not happen */
7449566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc," failed metric for coarse grid %e\n",(double)tm));
7459566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&Prol));
746806fa848SBarry Smith       } else if (metric > .0) {
7479566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc,"worst metric for coarse grid = %e\n",(double)metric));
7482e68589bSMark F. Adams       }
7493b4367a7SBarry Smith     } else SETERRQ(comm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG");
7502e68589bSMark F. Adams     { /* create next coords - output */
7512e68589bSMark F. Adams       PetscReal *crs_crds;
7529566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dim*nLocalSelected, &crs_crds));
7532e68589bSMark F. Adams       for (kk=0; kk<nLocalSelected; kk++) { /* grab local select nodes to promote - output */
754b43b03e9SMark F. Adams         PetscInt lid = clid_flid[kk];
755c8b0795cSMark F. Adams         for (jj=0; jj<dim; jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
7562e68589bSMark F. Adams       }
757c8b0795cSMark F. Adams 
7589566063dSJacob Faibussowitsch       PetscCall(PetscFree(pc_gamg->data));
759c8b0795cSMark F. Adams       pc_gamg->data    = crs_crds; /* out */
760c8b0795cSMark F. Adams       pc_gamg->data_sz = dim*nLocalSelected;
7612e68589bSMark F. Adams     }
7629566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&selected_2));
7632e68589bSMark F. Adams   }
764a2f3521dSMark F. Adams 
7652e68589bSMark F. Adams   *a_P_out = Prol;  /* out */
7669566063dSJacob Faibussowitsch   PetscCall(PetscFree(clid_flid));
767849bee69SMark Adams 
768c8b0795cSMark F. Adams   PetscFunctionReturn(0);
769c8b0795cSMark F. Adams }
770c8b0795cSMark F. Adams 
7719b8ffb57SJed Brown static PetscErrorCode PCDestroy_GAMG_GEO(PC pc)
7729b8ffb57SJed Brown {
7739b8ffb57SJed Brown   PetscFunctionBegin;
7749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL));
7759b8ffb57SJed Brown   PetscFunctionReturn(0);
7769b8ffb57SJed Brown }
7779b8ffb57SJed Brown 
778c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
779c8b0795cSMark F. Adams /*
780c8b0795cSMark F. Adams  PCCreateGAMG_GEO
781c8b0795cSMark F. Adams 
782c8b0795cSMark F. Adams   Input Parameter:
783c8b0795cSMark F. Adams    . pc -
784c8b0795cSMark F. Adams */
785c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_GEO(PC pc)
786c8b0795cSMark F. Adams {
787c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
788c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
789c8b0795cSMark F. Adams 
790c8b0795cSMark F. Adams   PetscFunctionBegin;
7911ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO;
7929b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_GEO;
793c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
794c8b0795cSMark F. Adams 
795c8b0795cSMark F. Adams   /* set internal function pointers */
796fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_GEO;
797fd1112cbSBarry Smith   pc_gamg->ops->coarsen           = PCGAMGCoarsen_GEO;
7981ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_GEO;
799fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = NULL;
8001ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_GEO;
801c8b0795cSMark F. Adams 
8029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_GEO));
8032e68589bSMark F. Adams   PetscFunctionReturn(0);
8042e68589bSMark F. Adams }
805