xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision 90fbc344e0a65d91564e1d6edfc837c02f3543ea)
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)
82e68589bSMark F. Adams #define REAL PetscReal
92e68589bSMark F. Adams #include <triangle.h>
102e68589bSMark F. Adams #endif
112e68589bSMark F. Adams 
122e68589bSMark F. Adams #include <petscblaslapack.h>
132e68589bSMark F. Adams 
14c8b0795cSMark F. Adams /* Private context for the GAMG preconditioner */
15c8b0795cSMark F. Adams typedef struct {
16c8b0795cSMark F. Adams   PetscInt lid;            /* local vertex index */
17c8b0795cSMark F. Adams   PetscInt degree;         /* vertex degree */
18c8b0795cSMark F. Adams } GAMGNode;
192fa5cd67SKarl Rupp 
20b817416eSBarry Smith PETSC_STATIC_INLINE int petsc_geo_mg_compare(const void *a, const void *b)
21c8b0795cSMark F. Adams {
22c8b0795cSMark F. Adams   return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree);
23c8b0795cSMark F. Adams }
24c8b0795cSMark F. Adams 
252e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
262e68589bSMark F. Adams /*
272e68589bSMark F. Adams    PCSetCoordinates_GEO
282e68589bSMark F. Adams 
292e68589bSMark F. Adams    Input Parameter:
302e68589bSMark F. Adams    .  pc - the preconditioner context
312e68589bSMark F. Adams */
32302f38e8SMark F. Adams PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
332e68589bSMark F. Adams {
342e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
352e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
362e68589bSMark F. Adams   PetscErrorCode ierr;
37*90fbc344SStefano Zampini   PetscInt       arrsz,bs,my0,kk,ii,nloc,Iend,aloc;
382e68589bSMark F. Adams   Mat            Amat = pc->pmat;
392e68589bSMark F. Adams 
402e68589bSMark F. Adams   PetscFunctionBegin;
412e68589bSMark F. Adams   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
422e68589bSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
432e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr);
44*90fbc344SStefano Zampini   aloc = (Iend-my0);
452e68589bSMark F. Adams   nloc = (Iend-my0)/bs;
46a2f3521dSMark F. Adams 
47*90fbc344SStefano Zampini   if (nloc!=a_nloc && aloc!=a_nloc) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of local blocks %D must be %D or %D.",a_nloc,nloc,aloc);
482e68589bSMark F. Adams 
49c8b0795cSMark F. Adams   pc_gamg->data_cell_rows = 1;
50*90fbc344SStefano Zampini   if (!coords && nloc > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
51c8b0795cSMark F. Adams   pc_gamg->data_cell_cols = ndm; /* coordinates */
522e68589bSMark F. Adams 
53c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
542e68589bSMark F. Adams 
552e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
567f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
572e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
58854ce69bSBarry Smith     ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr);
592e68589bSMark F. Adams   }
602e68589bSMark F. Adams   for (kk=0; kk<arrsz; kk++) pc_gamg->data[kk] = -999.;
612e68589bSMark F. Adams   pc_gamg->data[arrsz] = -99.;
622e68589bSMark F. Adams   /* copy data in - column oriented */
63*90fbc344SStefano Zampini   if (nloc == a_nloc) {
642e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) {
652e68589bSMark F. Adams       for (ii = 0; ii < ndm; ii++) {
662e68589bSMark F. Adams         pc_gamg->data[ii*nloc + kk] =  coords[kk*ndm + ii];
672e68589bSMark F. Adams       }
682e68589bSMark F. Adams     }
69*90fbc344SStefano Zampini   } else { /* assumes the coordinates are blocked */
70*90fbc344SStefano Zampini     for (kk = 0; kk < nloc; kk++) {
71*90fbc344SStefano Zampini       for (ii = 0; ii < ndm; ii++) {
72*90fbc344SStefano Zampini         pc_gamg->data[ii*nloc + kk] =  coords[bs*kk*ndm + ii];
73*90fbc344SStefano Zampini       }
74*90fbc344SStefano Zampini     }
75*90fbc344SStefano Zampini   }
7671959b99SBarry Smith   if (pc_gamg->data[arrsz] != -99.) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data[arrsz %D] %g != -99.",arrsz,pc_gamg->data[arrsz]);
772e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
782e68589bSMark F. Adams   PetscFunctionReturn(0);
792e68589bSMark F. Adams }
802e68589bSMark F. Adams 
812e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
822e68589bSMark F. Adams /*
832e68589bSMark F. Adams    PCSetData_GEO
842e68589bSMark F. Adams 
852e68589bSMark F. Adams   Input Parameter:
862e68589bSMark F. Adams    . pc -
872e68589bSMark F. Adams */
88b8cd405aSMark F. Adams PetscErrorCode PCSetData_GEO(PC pc, Mat m)
892e68589bSMark F. Adams {
902e68589bSMark F. Adams   PetscFunctionBegin;
91ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"GEO MG needs coordinates");
922e68589bSMark F. Adams }
932e68589bSMark F. Adams 
942e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
952e68589bSMark F. Adams /*
962e68589bSMark F. Adams    PCSetFromOptions_GEO
972e68589bSMark F. Adams 
982e68589bSMark F. Adams   Input Parameter:
992e68589bSMark F. Adams    . pc -
1002e68589bSMark F. Adams */
1014416b707SBarry Smith PetscErrorCode PCSetFromOptions_GEO(PetscOptionItems *PetscOptionsObject,PC pc)
1022e68589bSMark F. Adams {
1032e68589bSMark F. Adams   PetscErrorCode ierr;
1042e68589bSMark F. Adams 
1052e68589bSMark F. Adams   PetscFunctionBegin;
106e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-GEO options");CHKERRQ(ierr);
1072e68589bSMark F. Adams   {
1082e68589bSMark F. Adams     /* -pc_gamg_sa_nsmooths */
1092e68589bSMark F. Adams     /* pc_gamg_sa->smooths = 0; */
1102e68589bSMark F. Adams     /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */
1112e68589bSMark F. Adams     /*                        "smoothing steps for smoothed aggregation, usually 1 (0)", */
1122e68589bSMark F. Adams     /*                        "PCGAMGSetNSmooths_AGG", */
1132e68589bSMark F. Adams     /*                        pc_gamg_sa->smooths, */
1142e68589bSMark F. Adams     /*                        &pc_gamg_sa->smooths, */
1152e68589bSMark F. Adams     /*                        &flag);  */
1162e68589bSMark F. Adams     /* CHKERRQ(ierr); */
1172e68589bSMark F. Adams   }
1182e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1192e68589bSMark F. Adams   PetscFunctionReturn(0);
1202e68589bSMark F. Adams }
1212e68589bSMark F. Adams 
1222e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1232e68589bSMark F. Adams /*
1242e68589bSMark F. Adams  triangulateAndFormProl
1252e68589bSMark F. Adams 
1262e68589bSMark F. Adams    Input Parameter:
1272e68589bSMark F. Adams    . selected_2 - list of selected local ID, includes selected ghosts
128a2f3521dSMark F. Adams    . data_stride -
129a2f3521dSMark F. Adams    . coords[2*data_stride] - column vector of local coordinates w/ ghosts
1302adfe9d3SBarry Smith    . nselected_1 - selected IDs that go with base (1) graph includes selected ghosts
1310cbbd2e1SMark F. Adams    . clid_lid_1[nselected_1] - lids of selected (c) nodes   ???????????
1322adfe9d3SBarry Smith    . agg_lists_1 - list of aggregates selected_1 vertices of aggregate unselected vertices
1332e68589bSMark F. Adams    . crsGID[selected.size()] - global index for prolongation operator
1342e68589bSMark F. Adams    . bs - block size
1352e68589bSMark F. Adams   Output Parameter:
1362e68589bSMark F. Adams    . a_Prol - prolongation operator
1372e68589bSMark F. Adams    . a_worst_best - measure of worst missed fine vertex, 0 is no misses
1382e68589bSMark F. Adams */
1392adfe9d3SBarry 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,
1402adfe9d3SBarry Smith                                              const PetscInt crsGID[],PetscInt bs,Mat a_Prol,PetscReal *a_worst_best)
1412e68589bSMark F. Adams {
1422e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE)
1432e68589bSMark F. Adams   PetscErrorCode       ierr;
144b43b03e9SMark F. Adams   PetscInt             jj,tid,tt,idx,nselected_2;
1452e68589bSMark F. Adams   struct triangulateio in,mid;
1460cbbd2e1SMark F. Adams   const PetscInt       *selected_idx_2;
14773911c69SBarry Smith   PetscMPIInt          rank;
1482e68589bSMark F. Adams   PetscInt             Istart,Iend,nFineLoc,myFine0;
1492e68589bSMark F. Adams   int                  kk,nPlotPts,sid;
1503b4367a7SBarry Smith   MPI_Comm             comm;
151c8b0795cSMark F. Adams   PetscReal            tm;
152c8b0795cSMark F. Adams 
1536e111a19SKarl Rupp   PetscFunctionBegin;
1543b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
1553b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
156c8b0795cSMark F. Adams   ierr = ISGetSize(selected_2, &nselected_2);CHKERRQ(ierr);
1572e68589bSMark F. Adams   if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */
1582e68589bSMark F. Adams     *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */
159806fa848SBarry Smith   } else *a_worst_best = 0.0;
160b2566f29SBarry Smith   ierr = MPIU_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm);CHKERRQ(ierr);
161c8b0795cSMark F. Adams   if (tm > 0.0) {
162c8b0795cSMark F. Adams     *a_worst_best = 100.0;
1632e68589bSMark F. Adams     PetscFunctionReturn(0);
1642e68589bSMark F. Adams   }
1652e68589bSMark F. Adams   ierr     = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
1662e68589bSMark F. Adams   nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs;
1672e68589bSMark F. Adams   nPlotPts = nFineLoc; /* locals */
1682e68589bSMark F. Adams   /* traingle */
1692e68589bSMark F. Adams   /* Define input points - in*/
1702e68589bSMark F. Adams   in.numberofpoints          = nselected_2;
1712e68589bSMark F. Adams   in.numberofpointattributes = 0;
1722e68589bSMark F. Adams   /* get nselected points */
173e632b94dSBarry Smith   ierr = PetscMalloc1(2*nselected_2, &in.pointlist);CHKERRQ(ierr);
1742e68589bSMark F. Adams   ierr = ISGetIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
1752e68589bSMark F. Adams 
1762e68589bSMark F. Adams   for (kk=0,sid=0; kk<nselected_2; kk++,sid += 2) {
1772e68589bSMark F. Adams     PetscInt lid = selected_idx_2[kk];
1782e68589bSMark F. Adams     in.pointlist[sid]   = coords[lid];
179a2f3521dSMark F. Adams     in.pointlist[sid+1] = coords[data_stride + lid];
1802e68589bSMark F. Adams     if (lid>=nFineLoc) nPlotPts++;
1812e68589bSMark F. Adams   }
18271959b99SBarry Smith   if (sid != 2*nselected_2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != 2*nselected_2 %D",sid,nselected_2);
1832e68589bSMark F. Adams 
1842e68589bSMark F. Adams   in.numberofsegments      = 0;
1852e68589bSMark F. Adams   in.numberofedges         = 0;
1862e68589bSMark F. Adams   in.numberofholes         = 0;
1872e68589bSMark F. Adams   in.numberofregions       = 0;
1882e68589bSMark F. Adams   in.trianglelist          = 0;
1892e68589bSMark F. Adams   in.segmentmarkerlist     = 0;
1902e68589bSMark F. Adams   in.pointattributelist    = 0;
1912e68589bSMark F. Adams   in.pointmarkerlist       = 0;
1922e68589bSMark F. Adams   in.triangleattributelist = 0;
1932e68589bSMark F. Adams   in.trianglearealist      = 0;
1942e68589bSMark F. Adams   in.segmentlist           = 0;
1952e68589bSMark F. Adams   in.holelist              = 0;
1962e68589bSMark F. Adams   in.regionlist            = 0;
1972e68589bSMark F. Adams   in.edgelist              = 0;
1982e68589bSMark F. Adams   in.edgemarkerlist        = 0;
1992e68589bSMark F. Adams   in.normlist              = 0;
2002fa5cd67SKarl Rupp 
2012e68589bSMark F. Adams   /* triangulate */
2022e68589bSMark F. Adams   mid.pointlist = 0;            /* Not needed if -N switch used. */
2032e68589bSMark F. Adams   /* Not needed if -N switch used or number of point attributes is zero: */
2042e68589bSMark F. Adams   mid.pointattributelist = 0;
2052e68589bSMark F. Adams   mid.pointmarkerlist    = 0; /* Not needed if -N or -B switch used. */
2062e68589bSMark F. Adams   mid.trianglelist       = 0;    /* Not needed if -E switch used. */
2072e68589bSMark F. Adams   /* Not needed if -E switch used or number of triangle attributes is zero: */
2082e68589bSMark F. Adams   mid.triangleattributelist = 0;
2092e68589bSMark F. Adams   mid.neighborlist          = 0; /* Needed only if -n switch used. */
2102e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P not used: */
2112e68589bSMark F. Adams   mid.segmentlist = 0;
2122e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
2132e68589bSMark F. Adams   mid.segmentmarkerlist = 0;
2142e68589bSMark F. Adams   mid.edgelist          = 0;    /* Needed only if -e switch used. */
2152e68589bSMark F. Adams   mid.edgemarkerlist    = 0; /* Needed if -e used and -B not used. */
2162e68589bSMark F. Adams   mid.numberoftriangles = 0;
2172e68589bSMark F. Adams 
2182e68589bSMark F. Adams   /* Triangulate the points.  Switches are chosen to read and write a  */
2192e68589bSMark F. Adams   /*   PSLG (p), preserve the convex hull (c), number everything from  */
2202e68589bSMark F. Adams   /*   zero (z), assign a regional attribute to each element (A), and  */
2212e68589bSMark F. Adams   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
2222e68589bSMark F. Adams   /*   neighbor list (n).                                            */
2232e68589bSMark F. Adams   if (nselected_2 != 0) { /* inactive processor */
2242e68589bSMark F. Adams     char args[] = "npczQ"; /* c is needed ? */
2252e68589bSMark F. Adams     triangulate(args, &in, &mid, (struct triangulateio*) NULL);
2262e68589bSMark F. Adams     /* output .poly files for 'showme' */
2272e68589bSMark F. Adams     if (!PETSC_TRUE) {
2282e68589bSMark F. Adams       static int level = 1;
2292e68589bSMark F. Adams       FILE       *file; char fname[32];
2302e68589bSMark F. Adams 
231c5df96a5SBarry Smith       sprintf(fname,"C%d_%d.poly",level,rank); file = fopen(fname, "w");
2322e68589bSMark F. Adams       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
2332e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0);
2342e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2352e68589bSMark F. Adams       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid += 2) {
2362e68589bSMark F. Adams         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
2372e68589bSMark F. Adams       }
2382e68589bSMark F. Adams       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
2392e68589bSMark F. Adams       fprintf(file, "%d  %d\n",0,0);
2402e68589bSMark F. Adams       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
2412e68589bSMark F. Adams       /* One line: <# of holes> */
2422e68589bSMark F. Adams       fprintf(file, "%d\n",0);
2432e68589bSMark F. Adams       /* Following lines: <hole #> <x> <y> */
2442e68589bSMark F. Adams       /* Optional line: <# of regional attributes and/or area constraints> */
2452e68589bSMark F. Adams       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
2462e68589bSMark F. Adams       fclose(file);
2472e68589bSMark F. Adams 
2482e68589bSMark F. Adams       /* elems */
249c5df96a5SBarry Smith       sprintf(fname,"C%d_%d.ele",level,rank); file = fopen(fname, "w");
2502e68589bSMark F. Adams       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
2512e68589bSMark F. Adams       fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0);
2522e68589bSMark F. Adams       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
2532e68589bSMark F. Adams       for (kk=0,sid=0; kk<mid.numberoftriangles; kk++,sid += 3) {
2542e68589bSMark F. Adams         fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]);
2552e68589bSMark F. Adams       }
2562e68589bSMark F. Adams       fclose(file);
2572e68589bSMark F. Adams 
258c5df96a5SBarry Smith       sprintf(fname,"C%d_%d.node",level,rank); file = fopen(fname, "w");
2592e68589bSMark F. Adams       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
2602e68589bSMark F. Adams       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
2612e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n",nPlotPts,2,0,0);
2622e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2632e68589bSMark F. Adams       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid+=2) {
2642e68589bSMark F. Adams         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
2652e68589bSMark F. Adams       }
2662e68589bSMark F. Adams 
2672e68589bSMark F. Adams       sid /= 2;
2682e68589bSMark F. Adams       for (jj=0; jj<nFineLoc; jj++) {
2692e68589bSMark F. Adams         PetscBool sel = PETSC_TRUE;
2702e68589bSMark F. Adams         for (kk=0; kk<nselected_2 && sel; kk++) {
2712e68589bSMark F. Adams           PetscInt lid = selected_idx_2[kk];
2722e68589bSMark F. Adams           if (lid == jj) sel = PETSC_FALSE;
2732e68589bSMark F. Adams         }
2742fa5cd67SKarl Rupp         if (sel) fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]);
2752e68589bSMark F. Adams       }
2762e68589bSMark F. Adams       fclose(file);
27771959b99SBarry Smith       if (sid != nPlotPts) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != nPlotPts %D",sid,nPlotPts);
2782e68589bSMark F. Adams       level++;
2792e68589bSMark F. Adams     }
2802e68589bSMark F. Adams   }
2810cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2820cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
2832e68589bSMark F. Adams #endif
2842e68589bSMark F. Adams   { /* form P - setup some maps */
2850cbbd2e1SMark F. Adams     PetscInt clid,mm,*nTri,*node_tri;
2862e68589bSMark F. Adams 
287e632b94dSBarry Smith     ierr = PetscMalloc2(nselected_2, &node_tri,nselected_2, &nTri);CHKERRQ(ierr);
2882e68589bSMark F. Adams 
2892e68589bSMark F. Adams     /* need list of triangles on node */
2902e68589bSMark F. Adams     for (kk=0; kk<nselected_2; kk++) nTri[kk] = 0;
2912e68589bSMark F. Adams     for (tid=0,kk=0; tid<mid.numberoftriangles; tid++) {
2922e68589bSMark F. Adams       for (jj=0; jj<3; jj++) {
2932e68589bSMark F. Adams         PetscInt cid = mid.trianglelist[kk++];
2942e68589bSMark F. Adams         if (nTri[cid] == 0) node_tri[cid] = tid;
2952e68589bSMark F. Adams         nTri[cid]++;
2962e68589bSMark F. Adams       }
2972e68589bSMark F. Adams     }
2982e68589bSMark F. Adams #define EPS 1.e-12
2992e68589bSMark F. Adams     /* find points and set prolongation */
3000cbbd2e1SMark F. Adams     for (mm = clid = 0; mm < nFineLoc; mm++) {
301e78576d6SMark F. Adams       PetscBool ise;
302e78576d6SMark F. Adams       ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise);CHKERRQ(ierr);
303e78576d6SMark F. Adams       if (!ise) {
3040cbbd2e1SMark F. Adams         const PetscInt lid = mm;
3052e68589bSMark F. Adams         PetscScalar    AA[3][3];
3062e68589bSMark F. Adams         PetscBLASInt   N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO;
307539c167fSBarry Smith         PetscCDIntNd   *pos;
308539c167fSBarry Smith 
309e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
310e78576d6SMark F. Adams         while (pos) {
311e78576d6SMark F. Adams           PetscInt flid;
312484f0a72SBarry Smith           ierr = PetscCDIntNdGetID(pos, &flid);CHKERRQ(ierr);
313e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
3140cbbd2e1SMark F. Adams 
3152e68589bSMark F. Adams           if (flid < nFineLoc) {  /* could be a ghost */
3162e68589bSMark F. Adams             PetscInt       bestTID = -1; PetscReal best_alpha = 1.e10;
3172e68589bSMark F. Adams             const PetscInt fgid    = flid + myFine0;
3182e68589bSMark F. Adams             /* compute shape function for gid */
319a2f3521dSMark F. Adams             const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0};
3202e68589bSMark F. Adams             PetscBool       haveit    =PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3];
3212fa5cd67SKarl Rupp 
3222e68589bSMark F. Adams             /* look for it */
3230cbbd2e1SMark F. Adams             for (tid = node_tri[clid], jj=0;
3242e68589bSMark F. Adams                  jj < 5 && !haveit && tid != -1;
3252e68589bSMark F. Adams                  jj++) {
3262e68589bSMark F. Adams               for (tt=0; tt<3; tt++) {
3272e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3*tid + tt];
3282e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
329a2f3521dSMark F. Adams                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
3302e68589bSMark F. Adams                 clids[tt] = cid2; /* store for interp */
3312e68589bSMark F. Adams               }
3322e68589bSMark F. Adams 
3332e68589bSMark F. Adams               for (tt=0; tt<3; tt++) alpha[tt] = (PetscScalar)fcoord[tt];
3342e68589bSMark F. Adams 
3352e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
3368b83055fSJed Brown               PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3372e68589bSMark F. Adams               {
3382e68589bSMark F. Adams                 PetscBool have=PETSC_TRUE;  PetscReal lowest=1.e10;
3392e68589bSMark F. Adams                 for (tt = 0, idx = 0; tt < 3; tt++) {
3402e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
3412e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) < lowest) {
3422e68589bSMark F. Adams                     lowest = PetscRealPart(alpha[tt]);
3432e68589bSMark F. Adams                     idx    = tt;
3442e68589bSMark F. Adams                   }
3452e68589bSMark F. Adams                 }
3462e68589bSMark F. Adams                 haveit = have;
3472e68589bSMark F. Adams               }
3482e68589bSMark F. Adams               tid = mid.neighborlist[3*tid + idx];
3492e68589bSMark F. Adams             }
3502e68589bSMark F. Adams 
3512e68589bSMark F. Adams             if (!haveit) {
3522e68589bSMark F. Adams               /* brute force */
3532e68589bSMark F. Adams               for (tid=0; tid<mid.numberoftriangles && !haveit; tid++) {
3542e68589bSMark F. Adams                 for (tt=0; tt<3; tt++) {
3552e68589bSMark F. Adams                   PetscInt cid2 = mid.trianglelist[3*tid + tt];
3562e68589bSMark F. Adams                   PetscInt lid2 = selected_idx_2[cid2];
357a2f3521dSMark F. Adams                   AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
3582e68589bSMark F. Adams                   clids[tt] = cid2; /* store for interp */
3592e68589bSMark F. Adams                 }
3602e68589bSMark F. Adams                 for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
3612e68589bSMark F. Adams                 /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
3628b83055fSJed Brown                 PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3632e68589bSMark F. Adams                 {
3642e68589bSMark F. Adams                   PetscBool have=PETSC_TRUE;  PetscReal worst=0.0, v;
3652e68589bSMark F. Adams                   for (tt=0; tt<3 && have; tt++) {
3662e68589bSMark F. Adams                     if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE;
3672e68589bSMark F. Adams                     if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v;
3682e68589bSMark F. Adams                   }
3692e68589bSMark F. Adams                   if (worst < best_alpha) {
3702e68589bSMark F. Adams                     best_alpha = worst; bestTID = tid;
3712e68589bSMark F. Adams                   }
3722e68589bSMark F. Adams                   haveit = have;
3732e68589bSMark F. Adams                 }
3742e68589bSMark F. Adams               }
3752e68589bSMark F. Adams             }
3762e68589bSMark F. Adams             if (!haveit) {
3772e68589bSMark F. Adams               if (best_alpha > *a_worst_best) *a_worst_best = best_alpha;
3782e68589bSMark F. Adams               /* use best one */
3792e68589bSMark F. Adams               for (tt=0; tt<3; tt++) {
3802e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3*bestTID + tt];
3812e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
382a2f3521dSMark F. Adams                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
3832e68589bSMark F. Adams                 clids[tt] = cid2; /* store for interp */
3842e68589bSMark F. Adams               }
3852e68589bSMark F. Adams               for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
3862e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
3878b83055fSJed Brown               PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3882e68589bSMark F. Adams             }
3892e68589bSMark F. Adams 
3902e68589bSMark F. Adams             /* put in row of P */
3912e68589bSMark F. Adams             for (idx=0; idx<3; idx++) {
3922e68589bSMark F. Adams               PetscScalar shp = alpha[idx];
3932e68589bSMark F. Adams               if (PetscAbs(PetscRealPart(shp)) > 1.e-6) {
3942e68589bSMark F. Adams                 PetscInt cgid = crsGID[clids[idx]];
3952e68589bSMark F. Adams                 PetscInt jj   = cgid*bs, ii = fgid*bs; /* need to gloalize */
3962e68589bSMark F. Adams                 for (tt=0; tt < bs; tt++, ii++, jj++) {
3972e68589bSMark F. Adams                   ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES);CHKERRQ(ierr);
3982e68589bSMark F. Adams                 }
3992e68589bSMark F. Adams               }
4002e68589bSMark F. Adams             }
4012e68589bSMark F. Adams           }
4020cbbd2e1SMark F. Adams         } /* aggregates iterations */
4030cbbd2e1SMark F. Adams         clid++;
4040cbbd2e1SMark F. Adams       } /* a coarse agg */
4050cbbd2e1SMark F. Adams     } /* for all fine nodes */
4060cbbd2e1SMark F. Adams 
4072e68589bSMark F. Adams     ierr = ISRestoreIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
4082e68589bSMark F. Adams     ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4092e68589bSMark F. Adams     ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4102e68589bSMark F. Adams 
411e632b94dSBarry Smith     ierr = PetscFree2(node_tri,nTri);CHKERRQ(ierr);
4122e68589bSMark F. Adams   }
4130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4140cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
4152e68589bSMark F. Adams #endif
4162e68589bSMark F. Adams   free(mid.trianglelist);
4172e68589bSMark F. Adams   free(mid.neighborlist);
4182e68589bSMark F. Adams   ierr = PetscFree(in.pointlist);CHKERRQ(ierr);
4192e68589bSMark F. Adams   PetscFunctionReturn(0);
4202e68589bSMark F. Adams #else
4213b4367a7SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)a_Prol),PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG");
4222e68589bSMark F. Adams #endif
4232e68589bSMark F. Adams }
4242e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
4252e68589bSMark F. Adams /*
4262e68589bSMark F. Adams    getGIDsOnSquareGraph - square graph, get
4272e68589bSMark F. Adams 
4282e68589bSMark F. Adams    Input Parameter:
4290cbbd2e1SMark F. Adams    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
430b43b03e9SMark F. Adams    . clid_lid_1 - [nselected_1] lids of selected nodes
4312e68589bSMark F. Adams    . Gmat1 - graph that goes with 'selected_1'
4322e68589bSMark F. Adams    Output Parameter:
4332e68589bSMark F. Adams    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
4342e68589bSMark F. Adams    . a_Gmat_2 - graph that is squared of 'Gmat_1'
4352e68589bSMark F. Adams    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
4362e68589bSMark F. Adams */
4372adfe9d3SBarry Smith static PetscErrorCode getGIDsOnSquareGraph(PetscInt nselected_1,const PetscInt clid_lid_1[],const Mat Gmat1,IS *a_selected_2,Mat *a_Gmat_2,PetscInt **a_crsGID)
4382e68589bSMark F. Adams {
4392e68589bSMark F. Adams   PetscErrorCode ierr;
44073911c69SBarry Smith   PetscMPIInt    size;
441b43b03e9SMark F. Adams   PetscInt       *crsGID, kk,my0,Iend,nloc;
4423b4367a7SBarry Smith   MPI_Comm       comm;
4432e68589bSMark F. Adams 
4442e68589bSMark F. Adams   PetscFunctionBegin;
4453b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
4463b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
4472e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend);CHKERRQ(ierr); /* AIJ */
4482e68589bSMark F. Adams   nloc = Iend - my0; /* this does not change */
4492e68589bSMark F. Adams 
450c5df96a5SBarry Smith   if (size == 1) { /* not much to do in serial */
451785e854fSJed Brown     ierr = PetscMalloc1(nselected_1, &crsGID);CHKERRQ(ierr);
452b43b03e9SMark F. Adams     for (kk=0; kk<nselected_1; kk++) crsGID[kk] = kk;
4532e68589bSMark F. Adams     *a_Gmat_2 = 0;
454806fa848SBarry Smith     ierr      = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);CHKERRQ(ierr);
455806fa848SBarry Smith   } else {
456b43b03e9SMark F. Adams     PetscInt    idx,num_fine_ghosts,num_crs_ghost,myCrs0;
4572e68589bSMark F. Adams     Mat_MPIAIJ  *mpimat2;
4582e68589bSMark F. Adams     Mat         Gmat2;
4592e68589bSMark F. Adams     Vec         locState;
4602e68589bSMark F. Adams     PetscScalar *cpcol_state;
4612e68589bSMark F. Adams 
4622e68589bSMark F. Adams     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
463b43b03e9SMark F. Adams     kk = nselected_1;
464dbd2ff41SBarry Smith     ierr = MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
465b43b03e9SMark F. Adams     myCrs0 -= nselected_1;
4662e68589bSMark F. Adams 
467b43b03e9SMark F. Adams     if (a_Gmat_2) { /* output */
4682e68589bSMark F. Adams       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
4692e68589bSMark F. Adams       ierr      = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
4702e68589bSMark F. Adams       *a_Gmat_2 = Gmat2; /* output */
471806fa848SBarry Smith     } else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
4722e68589bSMark F. Adams     /* get coarse grid GIDS for selected (locals and ghosts) */
4732e68589bSMark F. Adams     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
4742a7a6963SBarry Smith     ierr    = MatCreateVecs(Gmat2, &locState, 0);CHKERRQ(ierr);
4752e68589bSMark F. Adams     ierr    = VecSet(locState, (PetscScalar)(PetscReal)(-1));CHKERRQ(ierr); /* set with UNKNOWN state */
476b43b03e9SMark F. Adams     for (kk=0; kk<nselected_1; kk++) {
477b43b03e9SMark F. Adams       PetscInt    fgid = clid_lid_1[kk] + my0;
4782e68589bSMark F. Adams       PetscScalar v    = (PetscScalar)(kk+myCrs0);
4792e68589bSMark F. Adams       ierr = VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES);CHKERRQ(ierr); /* set with PID */
4802e68589bSMark F. Adams     }
4812e68589bSMark F. Adams     ierr = VecAssemblyBegin(locState);CHKERRQ(ierr);
4822e68589bSMark F. Adams     ierr = VecAssemblyEnd(locState);CHKERRQ(ierr);
4832e68589bSMark F. Adams     ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4842e68589bSMark F. Adams     ierr =   VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4852e68589bSMark F. Adams     ierr = VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts);CHKERRQ(ierr);
4862e68589bSMark F. Adams     ierr = VecGetArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
4872e68589bSMark F. Adams     for (kk=0,num_crs_ghost=0; kk<num_fine_ghosts; kk++) {
4882e68589bSMark F. Adams       if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++;
4892e68589bSMark F. Adams     }
490854ce69bSBarry Smith     ierr = PetscMalloc1(nselected_1+num_crs_ghost, &crsGID);CHKERRQ(ierr); /* output */
4912e68589bSMark F. Adams     {
4922e68589bSMark F. Adams       PetscInt *selected_set;
493854ce69bSBarry Smith       ierr = PetscMalloc1(nselected_1+num_crs_ghost, &selected_set);CHKERRQ(ierr);
4942e68589bSMark F. Adams       /* do ghost of 'crsGID' */
495b43b03e9SMark F. Adams       for (kk=0,idx=nselected_1; kk<num_fine_ghosts; kk++) {
4962e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
4972e68589bSMark F. Adams           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
4982e68589bSMark F. Adams           selected_set[idx] = nloc + kk;
4992e68589bSMark F. Adams           crsGID[idx++]     = cgid;
5002e68589bSMark F. Adams         }
5012e68589bSMark F. Adams       }
50271959b99SBarry Smith       if (idx != (nselected_1+num_crs_ghost)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != (nselected_1 %D + num_crs_ghost %D)",idx,nselected_1,num_crs_ghost);
5032e68589bSMark F. Adams       ierr = VecRestoreArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
5042e68589bSMark F. Adams       /* do locals in 'crsGID' */
5052e68589bSMark F. Adams       ierr = VecGetArray(locState, &cpcol_state);CHKERRQ(ierr);
5062e68589bSMark F. Adams       for (kk=0,idx=0; kk<nloc; kk++) {
5072e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
5082e68589bSMark F. Adams           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
5092e68589bSMark F. Adams           selected_set[idx] = kk;
5102e68589bSMark F. Adams           crsGID[idx++]     = cgid;
5112e68589bSMark F. Adams         }
5122e68589bSMark F. Adams       }
51371959b99SBarry Smith       if (idx != nselected_1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != nselected_1 %D",idx,nselected_1);
5142e68589bSMark F. Adams       ierr = VecRestoreArray(locState, &cpcol_state);CHKERRQ(ierr);
5152e68589bSMark F. Adams 
5162e68589bSMark F. Adams       if (a_selected_2 != 0) { /* output */
517806fa848SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);CHKERRQ(ierr);
518806fa848SBarry Smith       } else {
5192e68589bSMark F. Adams         ierr = PetscFree(selected_set);CHKERRQ(ierr);
5202e68589bSMark F. Adams       }
5210cbbd2e1SMark F. Adams     }
5222e68589bSMark F. Adams     ierr = VecDestroy(&locState);CHKERRQ(ierr);
5232e68589bSMark F. Adams   }
5242e68589bSMark F. Adams   *a_crsGID = crsGID; /* output */
5252e68589bSMark F. Adams   PetscFunctionReturn(0);
5262e68589bSMark F. Adams }
5272e68589bSMark F. Adams 
5282e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
5292e68589bSMark F. Adams /*
530fd1112cbSBarry Smith    PCGAMGGraph_GEO
5312e68589bSMark F. Adams 
5322e68589bSMark F. Adams   Input Parameter:
5332e68589bSMark F. Adams    . pc - this
5342e68589bSMark F. Adams    . Amat - matrix on this fine level
5352e68589bSMark F. Adams   Output Parameter:
536c8b0795cSMark F. Adams    . a_Gmat
5372e68589bSMark F. Adams */
5382adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_GEO(PC pc,Mat Amat,Mat *a_Gmat)
539c8b0795cSMark F. Adams {
540c8b0795cSMark F. Adams   PetscErrorCode  ierr;
541c8b0795cSMark F. Adams   PC_MG           *mg      = (PC_MG*)pc->data;
542c8b0795cSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
543c1eae691SMark Adams   const PetscReal vfilter  = pc_gamg->threshold[0];
5443b4367a7SBarry Smith   MPI_Comm        comm;
545c8b0795cSMark F. Adams   Mat             Gmat;
5460cbbd2e1SMark F. Adams   PetscBool       set,flg,symm;
5476e111a19SKarl Rupp 
548c8b0795cSMark F. Adams   PetscFunctionBegin;
5493b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
550fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGGraph_GEO,0,0,0,0);CHKERRQ(ierr);
551c8b0795cSMark F. Adams 
5520cbbd2e1SMark F. Adams   ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr);
553263489e9SJed Brown   symm = (PetscBool)!(set && flg);
5540cbbd2e1SMark F. Adams 
5552d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
556bf4339c2SBarry Smith   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr);
557c8b0795cSMark F. Adams 
558c8b0795cSMark F. Adams   *a_Gmat = Gmat;
559fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGGraph_GEO,0,0,0,0);CHKERRQ(ierr);
560c8b0795cSMark F. Adams   PetscFunctionReturn(0);
561c8b0795cSMark F. Adams }
562c8b0795cSMark F. Adams 
563c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
564c8b0795cSMark F. Adams /*
565fd1112cbSBarry Smith    PCGAMGCoarsen_GEO
566c8b0795cSMark F. Adams 
567c8b0795cSMark F. Adams   Input Parameter:
568e0940f08SMark F. Adams    . a_pc - this
569e0940f08SMark F. Adams    . a_Gmat - graph
570c8b0795cSMark F. Adams   Output Parameter:
571c8b0795cSMark F. Adams    . a_llist_parent - linked list from selected indices for data locality only
572c8b0795cSMark F. Adams */
573fd1112cbSBarry Smith PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc,Mat *a_Gmat,PetscCoarsenData **a_llist_parent)
574c8b0795cSMark F. Adams {
575c8b0795cSMark F. Adams   PetscErrorCode ierr;
576c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,kk,Ii,ncols;
5770cbbd2e1SMark F. Adams   IS             perm;
578c8b0795cSMark F. Adams   GAMGNode       *gnodes;
579c8b0795cSMark F. Adams   PetscInt       *permute;
580e0940f08SMark F. Adams   Mat            Gmat  = *a_Gmat;
5813b4367a7SBarry Smith   MPI_Comm       comm;
582b43b03e9SMark F. Adams   MatCoarsen     crs;
583c8b0795cSMark F. Adams 
584c8b0795cSMark F. Adams   PetscFunctionBegin;
5853b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_pc,&comm);CHKERRQ(ierr);
5860cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
587c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr);
588c8b0795cSMark F. Adams   nloc = (Iend-Istart);
589c8b0795cSMark F. Adams 
590c8b0795cSMark F. Adams   /* create random permutation with sort for geo-mg */
591785e854fSJed Brown   ierr = PetscMalloc1(nloc, &gnodes);CHKERRQ(ierr);
592785e854fSJed Brown   ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr);
593c8b0795cSMark F. Adams 
594c8b0795cSMark F. Adams   for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */
595c8b0795cSMark F. Adams     ierr = MatGetRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr);
596c8b0795cSMark F. Adams     {
597c8b0795cSMark F. Adams       PetscInt lid = Ii - Istart;
598c8b0795cSMark F. Adams       gnodes[lid].lid    = lid;
599c8b0795cSMark F. Adams       gnodes[lid].degree = ncols;
600c8b0795cSMark F. Adams     }
601c8b0795cSMark F. Adams     ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr);
602c8b0795cSMark F. Adams   }
603c8b0795cSMark F. Adams   if (PETSC_TRUE) {
604e2c00dcbSBarry Smith     PetscRandom  rand;
605c8b0795cSMark F. Adams     PetscBool    *bIndexSet;
606e2c00dcbSBarry Smith     PetscReal    rr;
607e2c00dcbSBarry Smith     PetscInt     iSwapIndex;
608e2c00dcbSBarry Smith 
609e2c00dcbSBarry Smith     ierr = PetscRandomCreate(comm,&rand);CHKERRQ(ierr);
610e2c00dcbSBarry Smith     ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr);
6112fa5cd67SKarl Rupp     for (Ii = 0; Ii < nloc; Ii++) {
612e2c00dcbSBarry Smith       ierr = PetscRandomGetValueReal(rand,&rr);CHKERRQ(ierr);
613e2c00dcbSBarry Smith       iSwapIndex = (PetscInt) (rr*nloc);
6142fa5cd67SKarl Rupp       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
615c8b0795cSMark F. Adams         GAMGNode iTemp = gnodes[iSwapIndex];
616c8b0795cSMark F. Adams         gnodes[iSwapIndex]    = gnodes[Ii];
617c8b0795cSMark F. Adams         gnodes[Ii]            = iTemp;
618c8b0795cSMark F. Adams         bIndexSet[Ii]         = PETSC_TRUE;
619c8b0795cSMark F. Adams         bIndexSet[iSwapIndex] = PETSC_TRUE;
620c8b0795cSMark F. Adams       }
621c8b0795cSMark F. Adams     }
622e2c00dcbSBarry Smith     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
623c8b0795cSMark F. Adams     ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
624c8b0795cSMark F. Adams   }
625c8b0795cSMark F. Adams   /* only sort locals */
6260cbbd2e1SMark F. Adams   qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare);
627c8b0795cSMark F. Adams   /* create IS of permutation */
6282fa5cd67SKarl Rupp   for (kk=0; kk<nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */
629c8b0795cSMark F. Adams   ierr = PetscFree(gnodes);CHKERRQ(ierr);
630b817416eSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr);
631c8b0795cSMark F. Adams 
632c8b0795cSMark F. Adams   /* get MIS aggs */
633b43b03e9SMark F. Adams 
6343b4367a7SBarry Smith   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
635b43b03e9SMark F. Adams   ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr);
636b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
637b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr);
638b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_FALSE);CHKERRQ(ierr);
639b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
6400cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, a_llist_parent);CHKERRQ(ierr);
641b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
642c8b0795cSMark F. Adams 
643c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
6440cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
645c8b0795cSMark F. Adams   PetscFunctionReturn(0);
646c8b0795cSMark F. Adams }
647c8b0795cSMark F. Adams 
648c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
649c8b0795cSMark F. Adams /*
6500cbbd2e1SMark F. Adams  PCGAMGProlongator_GEO
651c8b0795cSMark F. Adams 
652c8b0795cSMark F. Adams  Input Parameter:
653c8b0795cSMark F. Adams  . pc - this
654c8b0795cSMark F. Adams  . Amat - matrix on this fine level
655c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
6560cbbd2e1SMark F. Adams  . selected_1 - [nselected]
6570cbbd2e1SMark F. Adams  . agg_lists - [nselected]
658c8b0795cSMark F. Adams  Output Parameter:
659c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
660c8b0795cSMark F. Adams  */
6612adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_GEO(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
6622e68589bSMark F. Adams {
6632e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
6642e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
665c8b0795cSMark F. Adams   const PetscInt dim      = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
6662e68589bSMark F. Adams   PetscErrorCode ierr;
667b43b03e9SMark F. Adams   PetscInt       Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid;
668c8b0795cSMark F. Adams   Mat            Prol;
669c5df96a5SBarry Smith   PetscMPIInt    rank, size;
6703b4367a7SBarry Smith   MPI_Comm       comm;
6710cbbd2e1SMark F. Adams   IS             selected_2,selected_1;
6722e68589bSMark F. Adams   const PetscInt *selected_idx;
673d9558ea9SBarry Smith   MatType        mtype;
6742e68589bSMark F. Adams 
6752e68589bSMark F. Adams   PetscFunctionBegin;
6763b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
6770cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
6783b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
6793b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
6802e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
6812e68589bSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
68271959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
68371959b99SBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) % bs %D",Iend,Istart,bs);
6842e68589bSMark F. Adams 
6852e68589bSMark F. Adams   /* get 'nLocalSelected' */
68641b27cdeSMark F. Adams   ierr = PetscCDGetMIS(agg_lists, &selected_1);CHKERRQ(ierr);
687b43b03e9SMark F. Adams   ierr = ISGetSize(selected_1, &jj);CHKERRQ(ierr);
688785e854fSJed Brown   ierr = PetscMalloc1(jj, &clid_flid);CHKERRQ(ierr);
6892e68589bSMark F. Adams   ierr = ISGetIndices(selected_1, &selected_idx);CHKERRQ(ierr);
690c8b0795cSMark F. Adams   for (kk=0,nLocalSelected=0; kk<jj; kk++) {
6912e68589bSMark F. Adams     PetscInt lid = selected_idx[kk];
692b43b03e9SMark F. Adams     if (lid<nloc) {
6930cbbd2e1SMark F. Adams       ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr);
6942fa5cd67SKarl Rupp       if (ncols>1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */
6950cbbd2e1SMark F. Adams       ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr);
696b43b03e9SMark F. Adams     }
6972e68589bSMark F. Adams   }
6982e68589bSMark F. Adams   ierr = ISRestoreIndices(selected_1, &selected_idx);CHKERRQ(ierr);
699a2f3521dSMark F. Adams   ierr = ISDestroy(&selected_1);CHKERRQ(ierr); /* this is selected_1 in serial */
7002e68589bSMark F. Adams 
701d9558ea9SBarry Smith   /* create prolongator  matrix */
702d9558ea9SBarry Smith   ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr);
7033b4367a7SBarry Smith   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
704806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
705a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, bs);CHKERRQ(ierr);
706d9558ea9SBarry Smith   ierr = MatSetType(Prol, mtype);CHKERRQ(ierr);
7070298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(Prol,3*data_cols,NULL);CHKERRQ(ierr);
7080298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(Prol,3*data_cols,NULL,3*data_cols,NULL);CHKERRQ(ierr);
7092e68589bSMark F. Adams 
710c8b0795cSMark F. Adams   /* can get all points "removed" - but not on geomg */
7112e68589bSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &jj);CHKERRQ(ierr);
7127f66b68fSMark Adams   if (!jj) {
713bf4339c2SBarry Smith     ierr = PetscInfo(pc,"ERROE: no selected points on coarse grid\n");CHKERRQ(ierr);
714b43b03e9SMark F. Adams     ierr = PetscFree(clid_flid);CHKERRQ(ierr);
7152e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
7160298fd71SBarry Smith     *a_P_out = NULL;  /* out */
7172e68589bSMark F. Adams     PetscFunctionReturn(0);
7182e68589bSMark F. Adams   }
7192e68589bSMark F. Adams 
7202e68589bSMark F. Adams   {
7212e68589bSMark F. Adams     PetscReal *coords;
722a2f3521dSMark F. Adams     PetscInt  data_stride;
7230298fd71SBarry Smith     PetscInt  *crsGID = NULL;
7242e68589bSMark F. Adams     Mat       Gmat2;
7252e68589bSMark F. Adams 
72671959b99SBarry Smith     if (dim != data_cols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim %D != data_cols %D",dim,data_cols);
7272e68589bSMark F. Adams     /* grow ghost data for better coarse grid cover of fine grid */
7280cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7290cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
7302e68589bSMark F. Adams #endif
731a2f3521dSMark F. Adams     /* messy method, squares graph and gets some data */
732806fa848SBarry Smith     ierr = getGIDsOnSquareGraph(nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID);CHKERRQ(ierr);
7332e68589bSMark F. Adams     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
7340cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7350cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
7362e68589bSMark F. Adams #endif
7372e68589bSMark F. Adams     /* create global vector of coorindates in 'coords' */
738c5df96a5SBarry Smith     if (size > 1) {
739806fa848SBarry Smith       ierr = PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords);CHKERRQ(ierr);
740806fa848SBarry Smith     } else {
741c8b0795cSMark F. Adams       coords      = (PetscReal*)pc_gamg->data;
742a2f3521dSMark F. Adams       data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols;
7432e68589bSMark F. Adams     }
7442e68589bSMark F. Adams     ierr = MatDestroy(&Gmat2);CHKERRQ(ierr);
7452e68589bSMark F. Adams 
7462e68589bSMark F. Adams     /* triangulate */
7472e68589bSMark F. Adams     if (dim == 2) {
748c8b0795cSMark F. Adams       PetscReal metric,tm;
7490cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7500cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
7512e68589bSMark F. Adams #endif
752806fa848SBarry Smith       ierr = triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric);CHKERRQ(ierr);
7530cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7540cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
7552e68589bSMark F. Adams #endif
7562e68589bSMark F. Adams       ierr = PetscFree(crsGID);CHKERRQ(ierr);
7572e68589bSMark F. Adams 
7582e68589bSMark F. Adams       /* clean up and create coordinates for coarse grid (output) */
759c5df96a5SBarry Smith       if (size > 1) ierr = PetscFree(coords);CHKERRQ(ierr);
7602e68589bSMark F. Adams 
761b2566f29SBarry Smith       ierr = MPIU_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm);CHKERRQ(ierr);
762c8b0795cSMark F. Adams       if (tm > 1.) { /* needs to be globalized - should not happen */
763bf4339c2SBarry Smith         ierr = PetscInfo1(pc," failed metric for coarse grid %e\n",(double)tm);CHKERRQ(ierr);
7642e68589bSMark F. Adams         ierr = MatDestroy(&Prol);CHKERRQ(ierr);
765806fa848SBarry Smith       } else if (metric > .0) {
766bf4339c2SBarry Smith         ierr = PetscInfo1(pc,"worst metric for coarse grid = %e\n",(double)metric);CHKERRQ(ierr);
7672e68589bSMark F. Adams       }
7683b4367a7SBarry Smith     } else SETERRQ(comm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG");
7692e68589bSMark F. Adams     { /* create next coords - output */
7702e68589bSMark F. Adams       PetscReal *crs_crds;
771785e854fSJed Brown       ierr = PetscMalloc1(dim*nLocalSelected, &crs_crds);CHKERRQ(ierr);
7722e68589bSMark F. Adams       for (kk=0; kk<nLocalSelected; kk++) { /* grab local select nodes to promote - output */
773b43b03e9SMark F. Adams         PetscInt lid = clid_flid[kk];
774c8b0795cSMark F. Adams         for (jj=0; jj<dim; jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
7752e68589bSMark F. Adams       }
776c8b0795cSMark F. Adams 
777c8b0795cSMark F. Adams       ierr             = PetscFree(pc_gamg->data);CHKERRQ(ierr);
778c8b0795cSMark F. Adams       pc_gamg->data    = crs_crds; /* out */
779c8b0795cSMark F. Adams       pc_gamg->data_sz = dim*nLocalSelected;
7802e68589bSMark F. Adams     }
781a2f3521dSMark F. Adams     ierr = ISDestroy(&selected_2);CHKERRQ(ierr);
7822e68589bSMark F. Adams   }
783a2f3521dSMark F. Adams 
7842e68589bSMark F. Adams   *a_P_out = Prol;  /* out */
785b43b03e9SMark F. Adams   ierr     = PetscFree(clid_flid);CHKERRQ(ierr);
7860cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
787c8b0795cSMark F. Adams   PetscFunctionReturn(0);
788c8b0795cSMark F. Adams }
789c8b0795cSMark F. Adams 
7909b8ffb57SJed Brown static PetscErrorCode PCDestroy_GAMG_GEO(PC pc)
7919b8ffb57SJed Brown {
7929b8ffb57SJed Brown   PetscErrorCode ierr;
7939b8ffb57SJed Brown 
7949b8ffb57SJed Brown   PetscFunctionBegin;
7959b8ffb57SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
7969b8ffb57SJed Brown   PetscFunctionReturn(0);
7979b8ffb57SJed Brown }
7989b8ffb57SJed Brown 
799c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
800c8b0795cSMark F. Adams /*
801c8b0795cSMark F. Adams  PCCreateGAMG_GEO
802c8b0795cSMark F. Adams 
803c8b0795cSMark F. Adams   Input Parameter:
804c8b0795cSMark F. Adams    . pc -
805c8b0795cSMark F. Adams */
806c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_GEO(PC pc)
807c8b0795cSMark F. Adams {
808c8b0795cSMark F. Adams   PetscErrorCode ierr;
809c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
810c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
811c8b0795cSMark F. Adams 
812c8b0795cSMark F. Adams   PetscFunctionBegin;
8131ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO;
8149b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_GEO;
815c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
816c8b0795cSMark F. Adams 
817c8b0795cSMark F. Adams   /* set internal function pointers */
818fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_GEO;
819fd1112cbSBarry Smith   pc_gamg->ops->coarsen           = PCGAMGCoarsen_GEO;
8201ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_GEO;
821fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = NULL;
8221ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_GEO;
823c8b0795cSMark F. Adams 
824bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_GEO);CHKERRQ(ierr);
8252e68589bSMark F. Adams   PetscFunctionReturn(0);
8262e68589bSMark F. Adams }
827