xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision e2c00dcb78d252d07ecb8cd54e33a12cb62c5f18)
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*/
6af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
72e68589bSMark F. Adams 
82e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE)
92e68589bSMark F. Adams #define REAL PetscReal
102e68589bSMark F. Adams #include <triangle.h>
112e68589bSMark F. Adams #endif
122e68589bSMark F. Adams 
132e68589bSMark F. Adams #include <petscblaslapack.h>
142e68589bSMark F. Adams 
15c8b0795cSMark F. Adams /* Private context for the GAMG preconditioner */
16c8b0795cSMark F. Adams typedef struct {
17c8b0795cSMark F. Adams   PetscInt lid;            /* local vertex index */
18c8b0795cSMark F. Adams   PetscInt degree;         /* vertex degree */
19c8b0795cSMark F. Adams } GAMGNode;
202fa5cd67SKarl Rupp 
210cbbd2e1SMark F. Adams int petsc_geo_mg_compare(const void *a, const void *b)
22c8b0795cSMark F. Adams {
23c8b0795cSMark F. Adams   return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree);
24c8b0795cSMark F. Adams }
25c8b0795cSMark F. Adams 
262e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
272e68589bSMark F. Adams /*
282e68589bSMark F. Adams    PCSetCoordinates_GEO
292e68589bSMark F. Adams 
302e68589bSMark F. Adams    Input Parameter:
312e68589bSMark F. Adams    .  pc - the preconditioner context
322e68589bSMark F. Adams */
332e68589bSMark F. Adams #undef __FUNCT__
342e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_GEO"
35302f38e8SMark F. Adams PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
362e68589bSMark F. Adams {
372e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
382e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
392e68589bSMark F. Adams   PetscErrorCode ierr;
402e68589bSMark F. Adams   PetscInt       arrsz,bs,my0,kk,ii,nloc,Iend;
412e68589bSMark F. Adams   Mat            Amat = pc->pmat;
422e68589bSMark F. Adams 
432e68589bSMark F. Adams   PetscFunctionBegin;
442e68589bSMark F. Adams   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
452e68589bSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
46a2f3521dSMark F. Adams 
472e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr);
482e68589bSMark F. Adams   nloc = (Iend-my0)/bs;
49a2f3521dSMark F. Adams 
50ce94432eSBarry Smith   if (nloc!=a_nloc) SETERRQ2(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Stokes not supported nloc = %d %d.",a_nloc,nloc);
51ce94432eSBarry Smith   if ((Iend-my0)%bs!=0) SETERRQ1(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
522e68589bSMark F. Adams 
53c8b0795cSMark F. Adams   pc_gamg->data_cell_rows = 1;
547f66b68fSMark Adams   if (!coords && nloc > 0) SETERRQ(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
55c8b0795cSMark F. Adams   pc_gamg->data_cell_cols = ndm; /* coordinates */
562e68589bSMark F. Adams 
57c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
582e68589bSMark F. Adams 
592e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
607f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
612e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
62854ce69bSBarry Smith     ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr);
632e68589bSMark F. Adams   }
642e68589bSMark F. Adams   for (kk=0; kk<arrsz; kk++) pc_gamg->data[kk] = -999.;
652e68589bSMark F. Adams   pc_gamg->data[arrsz] = -99.;
662e68589bSMark F. Adams   /* copy data in - column oriented */
672e68589bSMark F. Adams   for (kk = 0; kk < nloc; kk++) {
682e68589bSMark F. Adams     for (ii = 0; ii < ndm; ii++) {
692e68589bSMark F. Adams       pc_gamg->data[ii*nloc + kk] =  coords[kk*ndm + ii];
702e68589bSMark F. Adams     }
712e68589bSMark F. Adams   }
7271959b99SBarry 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]);
732e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
742e68589bSMark F. Adams   PetscFunctionReturn(0);
752e68589bSMark F. Adams }
762e68589bSMark F. Adams 
772e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
782e68589bSMark F. Adams /*
792e68589bSMark F. Adams    PCSetData_GEO
802e68589bSMark F. Adams 
812e68589bSMark F. Adams   Input Parameter:
822e68589bSMark F. Adams    . pc -
832e68589bSMark F. Adams */
842e68589bSMark F. Adams #undef __FUNCT__
852e68589bSMark F. Adams #define __FUNCT__ "PCSetData_GEO"
86b8cd405aSMark F. Adams PetscErrorCode PCSetData_GEO(PC pc, Mat m)
872e68589bSMark F. Adams {
882e68589bSMark F. Adams   PetscFunctionBegin;
89ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"GEO MG needs coordinates");
902e68589bSMark F. Adams }
912e68589bSMark F. Adams 
922e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
932e68589bSMark F. Adams /*
942e68589bSMark F. Adams    PCSetFromOptions_GEO
952e68589bSMark F. Adams 
962e68589bSMark F. Adams   Input Parameter:
972e68589bSMark F. Adams    . pc -
982e68589bSMark F. Adams */
992e68589bSMark F. Adams #undef __FUNCT__
1002e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GEO"
1018c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_GEO(PetscOptions *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 */
1392e68589bSMark F. Adams #undef __FUNCT__
1402e68589bSMark F. Adams #define __FUNCT__ "triangulateAndFormProl"
1412adfe9d3SBarry 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,
1422adfe9d3SBarry Smith                                              const PetscInt crsGID[],PetscInt bs,Mat a_Prol,PetscReal *a_worst_best)
1432e68589bSMark F. Adams {
1442e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE)
1452e68589bSMark F. Adams   PetscErrorCode       ierr;
146b43b03e9SMark F. Adams   PetscInt             jj,tid,tt,idx,nselected_2;
1472e68589bSMark F. Adams   struct triangulateio in,mid;
1480cbbd2e1SMark F. Adams   const PetscInt       *selected_idx_2;
14973911c69SBarry Smith   PetscMPIInt          rank;
1502e68589bSMark F. Adams   PetscInt             Istart,Iend,nFineLoc,myFine0;
1512e68589bSMark F. Adams   int                  kk,nPlotPts,sid;
1523b4367a7SBarry Smith   MPI_Comm             comm;
153c8b0795cSMark F. Adams   PetscReal            tm;
154c8b0795cSMark F. Adams 
1556e111a19SKarl Rupp   PetscFunctionBegin;
1563b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
1573b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
158c8b0795cSMark F. Adams   ierr = ISGetSize(selected_2, &nselected_2);CHKERRQ(ierr);
1592e68589bSMark F. Adams   if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */
1602e68589bSMark F. Adams     *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */
161806fa848SBarry Smith   } else *a_worst_best = 0.0;
1623b4367a7SBarry Smith   ierr = MPI_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm);CHKERRQ(ierr);
163c8b0795cSMark F. Adams   if (tm > 0.0) {
164c8b0795cSMark F. Adams     *a_worst_best = 100.0;
1652e68589bSMark F. Adams     PetscFunctionReturn(0);
1662e68589bSMark F. Adams   }
1672e68589bSMark F. Adams   ierr     = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
1682e68589bSMark F. Adams   nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs;
1692e68589bSMark F. Adams   nPlotPts = nFineLoc; /* locals */
1702e68589bSMark F. Adams   /* traingle */
1712e68589bSMark F. Adams   /* Define input points - in*/
1722e68589bSMark F. Adams   in.numberofpoints          = nselected_2;
1732e68589bSMark F. Adams   in.numberofpointattributes = 0;
1742e68589bSMark F. Adams   /* get nselected points */
175e632b94dSBarry Smith   ierr = PetscMalloc1(2*nselected_2, &in.pointlist);CHKERRQ(ierr);
1762e68589bSMark F. Adams   ierr = ISGetIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
1772e68589bSMark F. Adams 
1782e68589bSMark F. Adams   for (kk=0,sid=0; kk<nselected_2; kk++,sid += 2) {
1792e68589bSMark F. Adams     PetscInt lid = selected_idx_2[kk];
1802e68589bSMark F. Adams     in.pointlist[sid]   = coords[lid];
181a2f3521dSMark F. Adams     in.pointlist[sid+1] = coords[data_stride + lid];
1822e68589bSMark F. Adams     if (lid>=nFineLoc) nPlotPts++;
1832e68589bSMark F. Adams   }
18471959b99SBarry Smith   if (sid != 2*nselected_2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != 2*nselected_2 %D",sid,nselected_2);
1852e68589bSMark F. Adams 
1862e68589bSMark F. Adams   in.numberofsegments      = 0;
1872e68589bSMark F. Adams   in.numberofedges         = 0;
1882e68589bSMark F. Adams   in.numberofholes         = 0;
1892e68589bSMark F. Adams   in.numberofregions       = 0;
1902e68589bSMark F. Adams   in.trianglelist          = 0;
1912e68589bSMark F. Adams   in.segmentmarkerlist     = 0;
1922e68589bSMark F. Adams   in.pointattributelist    = 0;
1932e68589bSMark F. Adams   in.pointmarkerlist       = 0;
1942e68589bSMark F. Adams   in.triangleattributelist = 0;
1952e68589bSMark F. Adams   in.trianglearealist      = 0;
1962e68589bSMark F. Adams   in.segmentlist           = 0;
1972e68589bSMark F. Adams   in.holelist              = 0;
1982e68589bSMark F. Adams   in.regionlist            = 0;
1992e68589bSMark F. Adams   in.edgelist              = 0;
2002e68589bSMark F. Adams   in.edgemarkerlist        = 0;
2012e68589bSMark F. Adams   in.normlist              = 0;
2022fa5cd67SKarl Rupp 
2032e68589bSMark F. Adams   /* triangulate */
2042e68589bSMark F. Adams   mid.pointlist = 0;            /* Not needed if -N switch used. */
2052e68589bSMark F. Adams   /* Not needed if -N switch used or number of point attributes is zero: */
2062e68589bSMark F. Adams   mid.pointattributelist = 0;
2072e68589bSMark F. Adams   mid.pointmarkerlist    = 0; /* Not needed if -N or -B switch used. */
2082e68589bSMark F. Adams   mid.trianglelist       = 0;    /* Not needed if -E switch used. */
2092e68589bSMark F. Adams   /* Not needed if -E switch used or number of triangle attributes is zero: */
2102e68589bSMark F. Adams   mid.triangleattributelist = 0;
2112e68589bSMark F. Adams   mid.neighborlist          = 0; /* Needed only if -n switch used. */
2122e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P not used: */
2132e68589bSMark F. Adams   mid.segmentlist = 0;
2142e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
2152e68589bSMark F. Adams   mid.segmentmarkerlist = 0;
2162e68589bSMark F. Adams   mid.edgelist          = 0;    /* Needed only if -e switch used. */
2172e68589bSMark F. Adams   mid.edgemarkerlist    = 0; /* Needed if -e used and -B not used. */
2182e68589bSMark F. Adams   mid.numberoftriangles = 0;
2192e68589bSMark F. Adams 
2202e68589bSMark F. Adams   /* Triangulate the points.  Switches are chosen to read and write a  */
2212e68589bSMark F. Adams   /*   PSLG (p), preserve the convex hull (c), number everything from  */
2222e68589bSMark F. Adams   /*   zero (z), assign a regional attribute to each element (A), and  */
2232e68589bSMark F. Adams   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
2242e68589bSMark F. Adams   /*   neighbor list (n).                                            */
2252e68589bSMark F. Adams   if (nselected_2 != 0) { /* inactive processor */
2262e68589bSMark F. Adams     char args[] = "npczQ"; /* c is needed ? */
2272e68589bSMark F. Adams     triangulate(args, &in, &mid, (struct triangulateio*) NULL);
2282e68589bSMark F. Adams     /* output .poly files for 'showme' */
2292e68589bSMark F. Adams     if (!PETSC_TRUE) {
2302e68589bSMark F. Adams       static int level = 1;
2312e68589bSMark F. Adams       FILE       *file; char fname[32];
2322e68589bSMark F. Adams 
233c5df96a5SBarry Smith       sprintf(fname,"C%d_%d.poly",level,rank); file = fopen(fname, "w");
2342e68589bSMark F. Adams       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
2352e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0);
2362e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2372e68589bSMark F. Adams       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid += 2) {
2382e68589bSMark F. Adams         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
2392e68589bSMark F. Adams       }
2402e68589bSMark F. Adams       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
2412e68589bSMark F. Adams       fprintf(file, "%d  %d\n",0,0);
2422e68589bSMark F. Adams       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
2432e68589bSMark F. Adams       /* One line: <# of holes> */
2442e68589bSMark F. Adams       fprintf(file, "%d\n",0);
2452e68589bSMark F. Adams       /* Following lines: <hole #> <x> <y> */
2462e68589bSMark F. Adams       /* Optional line: <# of regional attributes and/or area constraints> */
2472e68589bSMark F. Adams       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
2482e68589bSMark F. Adams       fclose(file);
2492e68589bSMark F. Adams 
2502e68589bSMark F. Adams       /* elems */
251c5df96a5SBarry Smith       sprintf(fname,"C%d_%d.ele",level,rank); file = fopen(fname, "w");
2522e68589bSMark F. Adams       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
2532e68589bSMark F. Adams       fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0);
2542e68589bSMark F. Adams       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
2552e68589bSMark F. Adams       for (kk=0,sid=0; kk<mid.numberoftriangles; kk++,sid += 3) {
2562e68589bSMark F. Adams         fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]);
2572e68589bSMark F. Adams       }
2582e68589bSMark F. Adams       fclose(file);
2592e68589bSMark F. Adams 
260c5df96a5SBarry Smith       sprintf(fname,"C%d_%d.node",level,rank); file = fopen(fname, "w");
2612e68589bSMark F. Adams       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
2622e68589bSMark F. Adams       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
2632e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n",nPlotPts,2,0,0);
2642e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2652e68589bSMark F. Adams       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid+=2) {
2662e68589bSMark F. Adams         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
2672e68589bSMark F. Adams       }
2682e68589bSMark F. Adams 
2692e68589bSMark F. Adams       sid /= 2;
2702e68589bSMark F. Adams       for (jj=0; jj<nFineLoc; jj++) {
2712e68589bSMark F. Adams         PetscBool sel = PETSC_TRUE;
2722e68589bSMark F. Adams         for (kk=0; kk<nselected_2 && sel; kk++) {
2732e68589bSMark F. Adams           PetscInt lid = selected_idx_2[kk];
2742e68589bSMark F. Adams           if (lid == jj) sel = PETSC_FALSE;
2752e68589bSMark F. Adams         }
2762fa5cd67SKarl Rupp         if (sel) fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]);
2772e68589bSMark F. Adams       }
2782e68589bSMark F. Adams       fclose(file);
27971959b99SBarry Smith       if (sid != nPlotPts) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != nPlotPts %D",sid,nPlotPts);
2802e68589bSMark F. Adams       level++;
2812e68589bSMark F. Adams     }
2822e68589bSMark F. Adams   }
2830cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
2840cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
2852e68589bSMark F. Adams #endif
2862e68589bSMark F. Adams   { /* form P - setup some maps */
2870cbbd2e1SMark F. Adams     PetscInt clid,mm,*nTri,*node_tri;
2882e68589bSMark F. Adams 
289e632b94dSBarry Smith     ierr = PetscMalloc2(nselected_2, &node_tri,nselected_2, &nTri);CHKERRQ(ierr);
2902e68589bSMark F. Adams 
2912e68589bSMark F. Adams     /* need list of triangles on node */
2922e68589bSMark F. Adams     for (kk=0; kk<nselected_2; kk++) nTri[kk] = 0;
2932e68589bSMark F. Adams     for (tid=0,kk=0; tid<mid.numberoftriangles; tid++) {
2942e68589bSMark F. Adams       for (jj=0; jj<3; jj++) {
2952e68589bSMark F. Adams         PetscInt cid = mid.trianglelist[kk++];
2962e68589bSMark F. Adams         if (nTri[cid] == 0) node_tri[cid] = tid;
2972e68589bSMark F. Adams         nTri[cid]++;
2982e68589bSMark F. Adams       }
2992e68589bSMark F. Adams     }
3002e68589bSMark F. Adams #define EPS 1.e-12
3012e68589bSMark F. Adams     /* find points and set prolongation */
3020cbbd2e1SMark F. Adams     for (mm = clid = 0; mm < nFineLoc; mm++) {
303e78576d6SMark F. Adams       PetscBool ise;
304e78576d6SMark F. Adams       ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise);CHKERRQ(ierr);
305e78576d6SMark F. Adams       if (!ise) {
3060cbbd2e1SMark F. Adams         const PetscInt lid = mm;
307c41eeb4cSKarl Rupp         /* for (clid_iterator=0;clid_iterator<nselected_1;clid_iterator++) { */
3082e68589bSMark F. Adams         PetscScalar  AA[3][3];
3092e68589bSMark F. Adams         PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO;
31041b27cdeSMark F. Adams         PetscCDPos   pos;
311e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
312e78576d6SMark F. Adams         while (pos) {
313e78576d6SMark F. Adams           PetscInt flid;
314ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &flid);CHKERRQ(ierr);
315e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
3160cbbd2e1SMark F. Adams 
3172e68589bSMark F. Adams           if (flid < nFineLoc) {  /* could be a ghost */
3182e68589bSMark F. Adams             PetscInt       bestTID = -1; PetscReal best_alpha = 1.e10;
3192e68589bSMark F. Adams             const PetscInt fgid    = flid + myFine0;
3202e68589bSMark F. Adams             /* compute shape function for gid */
321a2f3521dSMark F. Adams             const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0};
3222e68589bSMark F. Adams             PetscBool       haveit    =PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3];
3232fa5cd67SKarl Rupp 
3242e68589bSMark F. Adams             /* look for it */
3250cbbd2e1SMark F. Adams             for (tid = node_tri[clid], jj=0;
3262e68589bSMark F. Adams                  jj < 5 && !haveit && tid != -1;
3272e68589bSMark F. Adams                  jj++) {
3282e68589bSMark F. Adams               for (tt=0; tt<3; tt++) {
3292e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3*tid + tt];
3302e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
331a2f3521dSMark F. Adams                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
3322e68589bSMark F. Adams                 clids[tt] = cid2; /* store for interp */
3332e68589bSMark F. Adams               }
3342e68589bSMark F. Adams 
3352e68589bSMark F. Adams               for (tt=0; tt<3; tt++) alpha[tt] = (PetscScalar)fcoord[tt];
3362e68589bSMark F. Adams 
3372e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
3388b83055fSJed Brown               PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3392e68589bSMark F. Adams               {
3402e68589bSMark F. Adams                 PetscBool have=PETSC_TRUE;  PetscReal lowest=1.e10;
3412e68589bSMark F. Adams                 for (tt = 0, idx = 0; tt < 3; tt++) {
3422e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
3432e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) < lowest) {
3442e68589bSMark F. Adams                     lowest = PetscRealPart(alpha[tt]);
3452e68589bSMark F. Adams                     idx    = tt;
3462e68589bSMark F. Adams                   }
3472e68589bSMark F. Adams                 }
3482e68589bSMark F. Adams                 haveit = have;
3492e68589bSMark F. Adams               }
3502e68589bSMark F. Adams               tid = mid.neighborlist[3*tid + idx];
3512e68589bSMark F. Adams             }
3522e68589bSMark F. Adams 
3532e68589bSMark F. Adams             if (!haveit) {
3542e68589bSMark F. Adams               /* brute force */
3552e68589bSMark F. Adams               for (tid=0; tid<mid.numberoftriangles && !haveit; tid++) {
3562e68589bSMark F. Adams                 for (tt=0; tt<3; tt++) {
3572e68589bSMark F. Adams                   PetscInt cid2 = mid.trianglelist[3*tid + tt];
3582e68589bSMark F. Adams                   PetscInt lid2 = selected_idx_2[cid2];
359a2f3521dSMark F. Adams                   AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
3602e68589bSMark F. Adams                   clids[tt] = cid2; /* store for interp */
3612e68589bSMark F. Adams                 }
3622e68589bSMark F. Adams                 for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
3632e68589bSMark F. Adams                 /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
3648b83055fSJed Brown                 PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3652e68589bSMark F. Adams                 {
3662e68589bSMark F. Adams                   PetscBool have=PETSC_TRUE;  PetscReal worst=0.0, v;
3672e68589bSMark F. Adams                   for (tt=0; tt<3 && have; tt++) {
3682e68589bSMark F. Adams                     if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE;
3692e68589bSMark F. Adams                     if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v;
3702e68589bSMark F. Adams                   }
3712e68589bSMark F. Adams                   if (worst < best_alpha) {
3722e68589bSMark F. Adams                     best_alpha = worst; bestTID = tid;
3732e68589bSMark F. Adams                   }
3742e68589bSMark F. Adams                   haveit = have;
3752e68589bSMark F. Adams                 }
3762e68589bSMark F. Adams               }
3772e68589bSMark F. Adams             }
3782e68589bSMark F. Adams             if (!haveit) {
3792e68589bSMark F. Adams               if (best_alpha > *a_worst_best) *a_worst_best = best_alpha;
3802e68589bSMark F. Adams               /* use best one */
3812e68589bSMark F. Adams               for (tt=0; tt<3; tt++) {
3822e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3*bestTID + tt];
3832e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
384a2f3521dSMark F. Adams                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
3852e68589bSMark F. Adams                 clids[tt] = cid2; /* store for interp */
3862e68589bSMark F. Adams               }
3872e68589bSMark F. Adams               for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
3882e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
3898b83055fSJed Brown               PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
3902e68589bSMark F. Adams             }
3912e68589bSMark F. Adams 
3922e68589bSMark F. Adams             /* put in row of P */
3932e68589bSMark F. Adams             for (idx=0; idx<3; idx++) {
3942e68589bSMark F. Adams               PetscScalar shp = alpha[idx];
3952e68589bSMark F. Adams               if (PetscAbs(PetscRealPart(shp)) > 1.e-6) {
3962e68589bSMark F. Adams                 PetscInt cgid = crsGID[clids[idx]];
3972e68589bSMark F. Adams                 PetscInt jj   = cgid*bs, ii = fgid*bs; /* need to gloalize */
3982e68589bSMark F. Adams                 for (tt=0; tt < bs; tt++, ii++, jj++) {
3992e68589bSMark F. Adams                   ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES);CHKERRQ(ierr);
4002e68589bSMark F. Adams                 }
4012e68589bSMark F. Adams               }
4022e68589bSMark F. Adams             }
4032e68589bSMark F. Adams           }
4040cbbd2e1SMark F. Adams         } /* aggregates iterations */
4050cbbd2e1SMark F. Adams         clid++;
4060cbbd2e1SMark F. Adams       } /* a coarse agg */
4070cbbd2e1SMark F. Adams     } /* for all fine nodes */
4080cbbd2e1SMark F. Adams 
4092e68589bSMark F. Adams     ierr = ISRestoreIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
4102e68589bSMark F. Adams     ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4112e68589bSMark F. Adams     ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4122e68589bSMark F. Adams 
413e632b94dSBarry Smith     ierr = PetscFree2(node_tri,nTri);CHKERRQ(ierr);
4142e68589bSMark F. Adams   }
4150cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4160cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
4172e68589bSMark F. Adams #endif
4182e68589bSMark F. Adams   free(mid.trianglelist);
4192e68589bSMark F. Adams   free(mid.neighborlist);
4202e68589bSMark F. Adams   ierr = PetscFree(in.pointlist);CHKERRQ(ierr);
4212e68589bSMark F. Adams   PetscFunctionReturn(0);
4222e68589bSMark F. Adams #else
4233b4367a7SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)a_Prol),PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG");
4242e68589bSMark F. Adams #endif
4252e68589bSMark F. Adams }
4262e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
4272e68589bSMark F. Adams /*
4282e68589bSMark F. Adams    getGIDsOnSquareGraph - square graph, get
4292e68589bSMark F. Adams 
4302e68589bSMark F. Adams    Input Parameter:
4310cbbd2e1SMark F. Adams    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
432b43b03e9SMark F. Adams    . clid_lid_1 - [nselected_1] lids of selected nodes
4332e68589bSMark F. Adams    . Gmat1 - graph that goes with 'selected_1'
4342e68589bSMark F. Adams    Output Parameter:
4352e68589bSMark F. Adams    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
4362e68589bSMark F. Adams    . a_Gmat_2 - graph that is squared of 'Gmat_1'
4372e68589bSMark F. Adams    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
4382e68589bSMark F. Adams */
4392e68589bSMark F. Adams #undef __FUNCT__
4402e68589bSMark F. Adams #define __FUNCT__ "getGIDsOnSquareGraph"
4412adfe9d3SBarry 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)
4422e68589bSMark F. Adams {
4432e68589bSMark F. Adams   PetscErrorCode ierr;
44473911c69SBarry Smith   PetscMPIInt    size;
445b43b03e9SMark F. Adams   PetscInt       *crsGID, kk,my0,Iend,nloc;
4463b4367a7SBarry Smith   MPI_Comm       comm;
4472e68589bSMark F. Adams 
4482e68589bSMark F. Adams   PetscFunctionBegin;
4493b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
4503b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
4512e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend);CHKERRQ(ierr); /* AIJ */
4522e68589bSMark F. Adams   nloc = Iend - my0; /* this does not change */
4532e68589bSMark F. Adams 
454c5df96a5SBarry Smith   if (size == 1) { /* not much to do in serial */
455785e854fSJed Brown     ierr = PetscMalloc1(nselected_1, &crsGID);CHKERRQ(ierr);
456b43b03e9SMark F. Adams     for (kk=0; kk<nselected_1; kk++) crsGID[kk] = kk;
4572e68589bSMark F. Adams     *a_Gmat_2 = 0;
458806fa848SBarry Smith     ierr      = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);CHKERRQ(ierr);
459806fa848SBarry Smith   } else {
460b43b03e9SMark F. Adams     PetscInt    idx,num_fine_ghosts,num_crs_ghost,myCrs0;
4612e68589bSMark F. Adams     Mat_MPIAIJ  *mpimat2;
4622e68589bSMark F. Adams     Mat         Gmat2;
4632e68589bSMark F. Adams     Vec         locState;
4642e68589bSMark F. Adams     PetscScalar *cpcol_state;
4652e68589bSMark F. Adams 
4662e68589bSMark F. Adams     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
467b43b03e9SMark F. Adams     kk = nselected_1;
468dbd2ff41SBarry Smith     ierr = MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
469b43b03e9SMark F. Adams     myCrs0 -= nselected_1;
4702e68589bSMark F. Adams 
471b43b03e9SMark F. Adams     if (a_Gmat_2) { /* output */
4722e68589bSMark F. Adams       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
4732e68589bSMark F. Adams       ierr      = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
4742e68589bSMark F. Adams       *a_Gmat_2 = Gmat2; /* output */
475806fa848SBarry Smith     } else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
4762e68589bSMark F. Adams     /* get coarse grid GIDS for selected (locals and ghosts) */
4772e68589bSMark F. Adams     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
4782a7a6963SBarry Smith     ierr    = MatCreateVecs(Gmat2, &locState, 0);CHKERRQ(ierr);
4792e68589bSMark F. Adams     ierr    = VecSet(locState, (PetscScalar)(PetscReal)(-1));CHKERRQ(ierr); /* set with UNKNOWN state */
480b43b03e9SMark F. Adams     for (kk=0; kk<nselected_1; kk++) {
481b43b03e9SMark F. Adams       PetscInt    fgid = clid_lid_1[kk] + my0;
4822e68589bSMark F. Adams       PetscScalar v    = (PetscScalar)(kk+myCrs0);
4832e68589bSMark F. Adams       ierr = VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES);CHKERRQ(ierr); /* set with PID */
4842e68589bSMark F. Adams     }
4852e68589bSMark F. Adams     ierr = VecAssemblyBegin(locState);CHKERRQ(ierr);
4862e68589bSMark F. Adams     ierr = VecAssemblyEnd(locState);CHKERRQ(ierr);
4872e68589bSMark F. Adams     ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4882e68589bSMark F. Adams     ierr =   VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4892e68589bSMark F. Adams     ierr = VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts);CHKERRQ(ierr);
4902e68589bSMark F. Adams     ierr = VecGetArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
4912e68589bSMark F. Adams     for (kk=0,num_crs_ghost=0; kk<num_fine_ghosts; kk++) {
4922e68589bSMark F. Adams       if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++;
4932e68589bSMark F. Adams     }
494854ce69bSBarry Smith     ierr = PetscMalloc1(nselected_1+num_crs_ghost, &crsGID);CHKERRQ(ierr); /* output */
4952e68589bSMark F. Adams     {
4962e68589bSMark F. Adams       PetscInt *selected_set;
497854ce69bSBarry Smith       ierr = PetscMalloc1(nselected_1+num_crs_ghost, &selected_set);CHKERRQ(ierr);
4982e68589bSMark F. Adams       /* do ghost of 'crsGID' */
499b43b03e9SMark F. Adams       for (kk=0,idx=nselected_1; kk<num_fine_ghosts; kk++) {
5002e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
5012e68589bSMark F. Adams           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
5022e68589bSMark F. Adams           selected_set[idx] = nloc + kk;
5032e68589bSMark F. Adams           crsGID[idx++]     = cgid;
5042e68589bSMark F. Adams         }
5052e68589bSMark F. Adams       }
50671959b99SBarry 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);
5072e68589bSMark F. Adams       ierr = VecRestoreArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
5082e68589bSMark F. Adams       /* do locals in 'crsGID' */
5092e68589bSMark F. Adams       ierr = VecGetArray(locState, &cpcol_state);CHKERRQ(ierr);
5102e68589bSMark F. Adams       for (kk=0,idx=0; kk<nloc; kk++) {
5112e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
5122e68589bSMark F. Adams           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
5132e68589bSMark F. Adams           selected_set[idx] = kk;
5142e68589bSMark F. Adams           crsGID[idx++]     = cgid;
5152e68589bSMark F. Adams         }
5162e68589bSMark F. Adams       }
51771959b99SBarry Smith       if (idx != nselected_1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != nselected_1 %D",idx,nselected_1);
5182e68589bSMark F. Adams       ierr = VecRestoreArray(locState, &cpcol_state);CHKERRQ(ierr);
5192e68589bSMark F. Adams 
5202e68589bSMark F. Adams       if (a_selected_2 != 0) { /* output */
521806fa848SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);CHKERRQ(ierr);
522806fa848SBarry Smith       } else {
5232e68589bSMark F. Adams         ierr = PetscFree(selected_set);CHKERRQ(ierr);
5242e68589bSMark F. Adams       }
5250cbbd2e1SMark F. Adams     }
5262e68589bSMark F. Adams     ierr = VecDestroy(&locState);CHKERRQ(ierr);
5272e68589bSMark F. Adams   }
5282e68589bSMark F. Adams   *a_crsGID = crsGID; /* output */
5292e68589bSMark F. Adams   PetscFunctionReturn(0);
5302e68589bSMark F. Adams }
5312e68589bSMark F. Adams 
5322e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
5332e68589bSMark F. Adams /*
534fd1112cbSBarry Smith    PCGAMGGraph_GEO
5352e68589bSMark F. Adams 
5362e68589bSMark F. Adams   Input Parameter:
5372e68589bSMark F. Adams    . pc - this
5382e68589bSMark F. Adams    . Amat - matrix on this fine level
5392e68589bSMark F. Adams   Output Parameter:
540c8b0795cSMark F. Adams    . a_Gmat
5412e68589bSMark F. Adams */
5422e68589bSMark F. Adams #undef __FUNCT__
543fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGGraph_GEO"
5442adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_GEO(PC pc,Mat Amat,Mat *a_Gmat)
545c8b0795cSMark F. Adams {
546c8b0795cSMark F. Adams   PetscErrorCode  ierr;
547c8b0795cSMark F. Adams   PC_MG           *mg      = (PC_MG*)pc->data;
548c8b0795cSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
549c8b0795cSMark F. Adams   const PetscReal vfilter  = pc_gamg->threshold;
5503b4367a7SBarry Smith   MPI_Comm        comm;
551c8b0795cSMark F. Adams   Mat             Gmat;
5520cbbd2e1SMark F. Adams   PetscBool       set,flg,symm;
5536e111a19SKarl Rupp 
554c8b0795cSMark F. Adams   PetscFunctionBegin;
5553b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
556fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGGraph_GEO,0,0,0,0);CHKERRQ(ierr);
557c8b0795cSMark F. Adams 
5580cbbd2e1SMark F. Adams   ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr);
559263489e9SJed Brown   symm = (PetscBool)!(set && flg);
5600cbbd2e1SMark F. Adams 
5612d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
562bf4339c2SBarry Smith   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr);
563c8b0795cSMark F. Adams 
564c8b0795cSMark F. Adams   *a_Gmat = Gmat;
565fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGGraph_GEO,0,0,0,0);CHKERRQ(ierr);
566c8b0795cSMark F. Adams   PetscFunctionReturn(0);
567c8b0795cSMark F. Adams }
568c8b0795cSMark F. Adams 
569c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
570c8b0795cSMark F. Adams /*
571fd1112cbSBarry Smith    PCGAMGCoarsen_GEO
572c8b0795cSMark F. Adams 
573c8b0795cSMark F. Adams   Input Parameter:
574e0940f08SMark F. Adams    . a_pc - this
575e0940f08SMark F. Adams    . a_Gmat - graph
576c8b0795cSMark F. Adams   Output Parameter:
577c8b0795cSMark F. Adams    . a_llist_parent - linked list from selected indices for data locality only
578c8b0795cSMark F. Adams */
579c8b0795cSMark F. Adams #undef __FUNCT__
580fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGCoarsen_GEO"
581fd1112cbSBarry Smith PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc,Mat *a_Gmat,PetscCoarsenData **a_llist_parent)
582c8b0795cSMark F. Adams {
583c8b0795cSMark F. Adams   PetscErrorCode ierr;
584c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,kk,Ii,ncols;
5850cbbd2e1SMark F. Adams   IS             perm;
586c8b0795cSMark F. Adams   GAMGNode       *gnodes;
587c8b0795cSMark F. Adams   PetscInt       *permute;
588e0940f08SMark F. Adams   Mat            Gmat  = *a_Gmat;
5893b4367a7SBarry Smith   MPI_Comm       comm;
590b43b03e9SMark F. Adams   MatCoarsen     crs;
591c8b0795cSMark F. Adams 
592c8b0795cSMark F. Adams   PetscFunctionBegin;
5933b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_pc,&comm);CHKERRQ(ierr);
5940cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
595c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr);
596c8b0795cSMark F. Adams   nloc = (Iend-Istart);
597c8b0795cSMark F. Adams 
598c8b0795cSMark F. Adams   /* create random permutation with sort for geo-mg */
599785e854fSJed Brown   ierr = PetscMalloc1(nloc, &gnodes);CHKERRQ(ierr);
600785e854fSJed Brown   ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr);
601c8b0795cSMark F. Adams 
602c8b0795cSMark F. Adams   for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */
603c8b0795cSMark F. Adams     ierr = MatGetRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr);
604c8b0795cSMark F. Adams     {
605c8b0795cSMark F. Adams       PetscInt lid = Ii - Istart;
606c8b0795cSMark F. Adams       gnodes[lid].lid    = lid;
607c8b0795cSMark F. Adams       gnodes[lid].degree = ncols;
608c8b0795cSMark F. Adams     }
609c8b0795cSMark F. Adams     ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr);
610c8b0795cSMark F. Adams   }
611c8b0795cSMark F. Adams   if (PETSC_TRUE) {
612*e2c00dcbSBarry Smith     PetscRandom  rand;
613c8b0795cSMark F. Adams     PetscBool    *bIndexSet;
614*e2c00dcbSBarry Smith     PetscReal    rr;
615*e2c00dcbSBarry Smith     PetscInt     iSwapIndex;
616*e2c00dcbSBarry Smith 
617*e2c00dcbSBarry Smith     ierr = PetscRandomCreate(comm,&rand);CHKERRQ(ierr);
618*e2c00dcbSBarry Smith     ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr);
6192fa5cd67SKarl Rupp     for (Ii = 0; Ii < nloc; Ii++) {
620*e2c00dcbSBarry Smith       ierr = PetscRandomGetValueReal(rand,&rr);CHKERRQ(ierr);
621*e2c00dcbSBarry Smith       iSwapIndex = (PetscInt) (rr*nloc);
6222fa5cd67SKarl Rupp       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
623c8b0795cSMark F. Adams         GAMGNode iTemp = gnodes[iSwapIndex];
624c8b0795cSMark F. Adams         gnodes[iSwapIndex]    = gnodes[Ii];
625c8b0795cSMark F. Adams         gnodes[Ii]            = iTemp;
626c8b0795cSMark F. Adams         bIndexSet[Ii]         = PETSC_TRUE;
627c8b0795cSMark F. Adams         bIndexSet[iSwapIndex] = PETSC_TRUE;
628c8b0795cSMark F. Adams       }
629c8b0795cSMark F. Adams     }
630*e2c00dcbSBarry Smith     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
631c8b0795cSMark F. Adams     ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
632c8b0795cSMark F. Adams   }
633c8b0795cSMark F. Adams   /* only sort locals */
6340cbbd2e1SMark F. Adams   qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare);
635c8b0795cSMark F. Adams   /* create IS of permutation */
6362fa5cd67SKarl Rupp   for (kk=0; kk<nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */
637806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr);
638c8b0795cSMark F. Adams 
639c8b0795cSMark F. Adams   ierr = PetscFree(gnodes);CHKERRQ(ierr);
640c8b0795cSMark F. Adams 
641c8b0795cSMark F. Adams   /* get MIS aggs */
642b43b03e9SMark F. Adams 
6433b4367a7SBarry Smith   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
644b43b03e9SMark F. Adams   ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr);
645b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
646b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr);
647b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_FALSE);CHKERRQ(ierr);
648b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
6490cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, a_llist_parent);CHKERRQ(ierr);
650b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
651c8b0795cSMark F. Adams 
652c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
6530cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
654c8b0795cSMark F. Adams   PetscFunctionReturn(0);
655c8b0795cSMark F. Adams }
656c8b0795cSMark F. Adams 
657c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
658c8b0795cSMark F. Adams /*
6590cbbd2e1SMark F. Adams  PCGAMGProlongator_GEO
660c8b0795cSMark F. Adams 
661c8b0795cSMark F. Adams  Input Parameter:
662c8b0795cSMark F. Adams  . pc - this
663c8b0795cSMark F. Adams  . Amat - matrix on this fine level
664c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
6650cbbd2e1SMark F. Adams  . selected_1 - [nselected]
6660cbbd2e1SMark F. Adams  . agg_lists - [nselected]
667c8b0795cSMark F. Adams  Output Parameter:
668c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
669c8b0795cSMark F. Adams  */
670c8b0795cSMark F. Adams #undef __FUNCT__
6710cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_GEO"
6722adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_GEO(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
6732e68589bSMark F. Adams {
6742e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
6752e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
676c8b0795cSMark F. Adams   const PetscInt dim      = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
6772e68589bSMark F. Adams   PetscErrorCode ierr;
678b43b03e9SMark F. Adams   PetscInt       Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid;
679c8b0795cSMark F. Adams   Mat            Prol;
680c5df96a5SBarry Smith   PetscMPIInt    rank, size;
6813b4367a7SBarry Smith   MPI_Comm       comm;
6820cbbd2e1SMark F. Adams   IS             selected_2,selected_1;
6832e68589bSMark F. Adams   const PetscInt *selected_idx;
684d9558ea9SBarry Smith   MatType        mtype;
6852e68589bSMark F. Adams 
6862e68589bSMark F. Adams   PetscFunctionBegin;
6873b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
6880cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
6893b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
6903b4367a7SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
6912e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
6922e68589bSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
69371959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
69471959b99SBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) % bs %D",Iend,Istart,bs);
6952e68589bSMark F. Adams 
6962e68589bSMark F. Adams   /* get 'nLocalSelected' */
69741b27cdeSMark F. Adams   ierr = PetscCDGetMIS(agg_lists, &selected_1);CHKERRQ(ierr);
698b43b03e9SMark F. Adams   ierr = ISGetSize(selected_1, &jj);CHKERRQ(ierr);
699785e854fSJed Brown   ierr = PetscMalloc1(jj, &clid_flid);CHKERRQ(ierr);
7002e68589bSMark F. Adams   ierr = ISGetIndices(selected_1, &selected_idx);CHKERRQ(ierr);
701c8b0795cSMark F. Adams   for (kk=0,nLocalSelected=0; kk<jj; kk++) {
7022e68589bSMark F. Adams     PetscInt lid = selected_idx[kk];
703b43b03e9SMark F. Adams     if (lid<nloc) {
7040cbbd2e1SMark F. Adams       ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr);
7052fa5cd67SKarl Rupp       if (ncols>1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */
7060cbbd2e1SMark F. Adams       ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr);
707b43b03e9SMark F. Adams     }
7082e68589bSMark F. Adams   }
7092e68589bSMark F. Adams   ierr = ISRestoreIndices(selected_1, &selected_idx);CHKERRQ(ierr);
710a2f3521dSMark F. Adams   ierr = ISDestroy(&selected_1);CHKERRQ(ierr); /* this is selected_1 in serial */
7112e68589bSMark F. Adams 
712d9558ea9SBarry Smith   /* create prolongator  matrix */
713d9558ea9SBarry Smith   ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr);
7143b4367a7SBarry Smith   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
715806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
716a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, bs);CHKERRQ(ierr);
717d9558ea9SBarry Smith   ierr = MatSetType(Prol, mtype);CHKERRQ(ierr);
7180298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(Prol,3*data_cols,NULL);CHKERRQ(ierr);
7190298fd71SBarry Smith   ierr = MatMPIAIJSetPreallocation(Prol,3*data_cols,NULL,3*data_cols,NULL);CHKERRQ(ierr);
7202e68589bSMark F. Adams 
721c8b0795cSMark F. Adams   /* can get all points "removed" - but not on geomg */
7222e68589bSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &jj);CHKERRQ(ierr);
7237f66b68fSMark Adams   if (!jj) {
724bf4339c2SBarry Smith     ierr = PetscInfo(pc,"ERROE: no selected points on coarse grid\n");CHKERRQ(ierr);
725b43b03e9SMark F. Adams     ierr = PetscFree(clid_flid);CHKERRQ(ierr);
7262e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
7270298fd71SBarry Smith     *a_P_out = NULL;  /* out */
7282e68589bSMark F. Adams     PetscFunctionReturn(0);
7292e68589bSMark F. Adams   }
7302e68589bSMark F. Adams 
7312e68589bSMark F. Adams   {
7322e68589bSMark F. Adams     PetscReal *coords;
733a2f3521dSMark F. Adams     PetscInt  data_stride;
7340298fd71SBarry Smith     PetscInt  *crsGID = NULL;
7352e68589bSMark F. Adams     Mat       Gmat2;
7362e68589bSMark F. Adams 
73771959b99SBarry Smith     if (dim != data_cols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim %D != data_cols %D",dim,data_cols);
7382e68589bSMark F. Adams     /* grow ghost data for better coarse grid cover of fine grid */
7390cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7400cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
7412e68589bSMark F. Adams #endif
742a2f3521dSMark F. Adams     /* messy method, squares graph and gets some data */
743806fa848SBarry Smith     ierr = getGIDsOnSquareGraph(nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID);CHKERRQ(ierr);
7442e68589bSMark F. Adams     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
7450cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7460cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
7472e68589bSMark F. Adams #endif
7482e68589bSMark F. Adams     /* create global vector of coorindates in 'coords' */
749c5df96a5SBarry Smith     if (size > 1) {
750806fa848SBarry Smith       ierr = PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords);CHKERRQ(ierr);
751806fa848SBarry Smith     } else {
752c8b0795cSMark F. Adams       coords      = (PetscReal*)pc_gamg->data;
753a2f3521dSMark F. Adams       data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols;
7542e68589bSMark F. Adams     }
7552e68589bSMark F. Adams     ierr = MatDestroy(&Gmat2);CHKERRQ(ierr);
7562e68589bSMark F. Adams 
7572e68589bSMark F. Adams     /* triangulate */
7582e68589bSMark F. Adams     if (dim == 2) {
759c8b0795cSMark F. Adams       PetscReal metric,tm;
7600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7610cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
7622e68589bSMark F. Adams #endif
763806fa848SBarry Smith       ierr = triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric);CHKERRQ(ierr);
7640cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
7650cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
7662e68589bSMark F. Adams #endif
7672e68589bSMark F. Adams       ierr = PetscFree(crsGID);CHKERRQ(ierr);
7682e68589bSMark F. Adams 
7692e68589bSMark F. Adams       /* clean up and create coordinates for coarse grid (output) */
770c5df96a5SBarry Smith       if (size > 1) ierr = PetscFree(coords);CHKERRQ(ierr);
7712e68589bSMark F. Adams 
7723b4367a7SBarry Smith       ierr = MPI_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm);CHKERRQ(ierr);
773c8b0795cSMark F. Adams       if (tm > 1.) { /* needs to be globalized - should not happen */
774bf4339c2SBarry Smith         ierr = PetscInfo1(pc," failed metric for coarse grid %e\n",(double)tm);CHKERRQ(ierr);
7752e68589bSMark F. Adams         ierr = MatDestroy(&Prol);CHKERRQ(ierr);
776806fa848SBarry Smith       } else if (metric > .0) {
777bf4339c2SBarry Smith         ierr = PetscInfo1(pc,"worst metric for coarse grid = %e\n",(double)metric);CHKERRQ(ierr);
7782e68589bSMark F. Adams       }
7793b4367a7SBarry Smith     } else SETERRQ(comm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG");
7802e68589bSMark F. Adams     { /* create next coords - output */
7812e68589bSMark F. Adams       PetscReal *crs_crds;
782785e854fSJed Brown       ierr = PetscMalloc1(dim*nLocalSelected, &crs_crds);CHKERRQ(ierr);
7832e68589bSMark F. Adams       for (kk=0; kk<nLocalSelected; kk++) { /* grab local select nodes to promote - output */
784b43b03e9SMark F. Adams         PetscInt lid = clid_flid[kk];
785c8b0795cSMark F. Adams         for (jj=0; jj<dim; jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
7862e68589bSMark F. Adams       }
787c8b0795cSMark F. Adams 
788c8b0795cSMark F. Adams       ierr             = PetscFree(pc_gamg->data);CHKERRQ(ierr);
789c8b0795cSMark F. Adams       pc_gamg->data    = crs_crds; /* out */
790c8b0795cSMark F. Adams       pc_gamg->data_sz = dim*nLocalSelected;
7912e68589bSMark F. Adams     }
792a2f3521dSMark F. Adams     ierr = ISDestroy(&selected_2);CHKERRQ(ierr);
7932e68589bSMark F. Adams   }
794a2f3521dSMark F. Adams 
7952e68589bSMark F. Adams   *a_P_out = Prol;  /* out */
796b43b03e9SMark F. Adams   ierr     = PetscFree(clid_flid);CHKERRQ(ierr);
7970cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
798c8b0795cSMark F. Adams   PetscFunctionReturn(0);
799c8b0795cSMark F. Adams }
800c8b0795cSMark F. Adams 
8019b8ffb57SJed Brown #undef __FUNCT__
8029b8ffb57SJed Brown #define __FUNCT__ "PCDestroy_GAMG_GEO"
8039b8ffb57SJed Brown static PetscErrorCode PCDestroy_GAMG_GEO(PC pc)
8049b8ffb57SJed Brown {
8059b8ffb57SJed Brown   PetscErrorCode ierr;
8069b8ffb57SJed Brown 
8079b8ffb57SJed Brown   PetscFunctionBegin;
8089b8ffb57SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
8099b8ffb57SJed Brown   PetscFunctionReturn(0);
8109b8ffb57SJed Brown }
8119b8ffb57SJed Brown 
812c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
813c8b0795cSMark F. Adams /*
814c8b0795cSMark F. Adams  PCCreateGAMG_GEO
815c8b0795cSMark F. Adams 
816c8b0795cSMark F. Adams   Input Parameter:
817c8b0795cSMark F. Adams    . pc -
818c8b0795cSMark F. Adams */
819c8b0795cSMark F. Adams #undef __FUNCT__
820c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_GEO"
821c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_GEO(PC pc)
822c8b0795cSMark F. Adams {
823c8b0795cSMark F. Adams   PetscErrorCode ierr;
824c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
825c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
826c8b0795cSMark F. Adams 
827c8b0795cSMark F. Adams   PetscFunctionBegin;
8281ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO;
8299b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_GEO;
830c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
831c8b0795cSMark F. Adams 
832c8b0795cSMark F. Adams   /* set internal function pointers */
833fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_GEO;
834fd1112cbSBarry Smith   pc_gamg->ops->coarsen           = PCGAMGCoarsen_GEO;
8351ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_GEO;
836fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = NULL;
8371ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_GEO;
838c8b0795cSMark F. Adams 
839bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_GEO);CHKERRQ(ierr);
8402e68589bSMark F. Adams   PetscFunctionReturn(0);
8412e68589bSMark F. Adams }
842