xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision 1147fc2a6bf7922defbbca15fc6c33f234803ff0)
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*/
6b45d2f2cSJed Brown #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 <assert.h>
142e68589bSMark F. Adams #include <petscblaslapack.h>
152e68589bSMark F. Adams 
16c8b0795cSMark F. Adams /* Private context for the GAMG preconditioner */
17c8b0795cSMark F. Adams typedef struct{
18c8b0795cSMark F. Adams   PetscInt       lid;      /* local vertex index */
19c8b0795cSMark F. Adams   PetscInt       degree;   /* vertex degree */
20c8b0795cSMark F. Adams } GAMGNode;
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 EXTERN_C_BEGIN
342e68589bSMark F. Adams #undef __FUNCT__
352e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_GEO"
36302f38e8SMark F. Adams PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
372e68589bSMark F. Adams {
382e68589bSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
392e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
402e68589bSMark F. Adams   PetscErrorCode ierr;
412e68589bSMark F. Adams   PetscInt       arrsz,bs,my0,kk,ii,nloc,Iend;
422e68589bSMark F. Adams   Mat            Amat = pc->pmat;
432e68589bSMark F. Adams 
442e68589bSMark F. Adams   PetscFunctionBegin;
452e68589bSMark F. Adams   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
462e68589bSMark F. Adams   ierr  = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
47a2f3521dSMark F. Adams 
482e68589bSMark F. Adams   ierr  = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr);
492e68589bSMark F. Adams   nloc = (Iend-my0)/bs;
50a2f3521dSMark F. Adams 
51302f38e8SMark F. Adams   if (nloc!=a_nloc)SETERRQ2(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Stokes not supported nloc = %d %d.",a_nloc,nloc);
522e68589bSMark F. Adams   if ((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
532e68589bSMark F. Adams 
54c8b0795cSMark F. Adams   pc_gamg->data_cell_rows = 1;
55806fa848SBarry Smith   if (coords==0 && nloc > 0) SETERRQ(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
56c8b0795cSMark F. Adams   pc_gamg->data_cell_cols = ndm; /* coordinates */
572e68589bSMark F. Adams 
58c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
592e68589bSMark F. Adams 
602e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
612e68589bSMark F. Adams   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
622e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
632e68589bSMark F. Adams     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
642e68589bSMark F. Adams   }
652e68589bSMark F. Adams   for (kk=0;kk<arrsz;kk++)pc_gamg->data[kk] = -999.;
662e68589bSMark F. Adams   pc_gamg->data[arrsz] = -99.;
672e68589bSMark F. Adams   /* copy data in - column oriented */
682e68589bSMark F. Adams   for (kk = 0 ; kk < nloc ; kk++) {
692e68589bSMark F. Adams     for (ii = 0 ; ii < ndm ; ii++) {
702e68589bSMark F. Adams       pc_gamg->data[ii*nloc + kk] =  coords[kk*ndm + ii];
712e68589bSMark F. Adams     }
722e68589bSMark F. Adams   }
732e68589bSMark F. Adams   assert(pc_gamg->data[arrsz] == -99.);
742e68589bSMark F. Adams 
752e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
762e68589bSMark F. Adams 
772e68589bSMark F. Adams   PetscFunctionReturn(0);
782e68589bSMark F. Adams }
792e68589bSMark F. Adams EXTERN_C_END
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 */
882e68589bSMark F. Adams #undef __FUNCT__
892e68589bSMark F. Adams #define __FUNCT__ "PCSetData_GEO"
90b8cd405aSMark F. Adams PetscErrorCode PCSetData_GEO(PC pc, Mat m)
912e68589bSMark F. Adams {
922e68589bSMark F. Adams   PetscFunctionBegin;
93c666adbfSMark F. Adams   SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_PLIB,"GEO MG needs coordinates");
942e68589bSMark F. Adams }
952e68589bSMark F. Adams 
962e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
972e68589bSMark F. Adams /*
982e68589bSMark F. Adams    PCSetFromOptions_GEO
992e68589bSMark F. Adams 
1002e68589bSMark F. Adams   Input Parameter:
1012e68589bSMark F. Adams    . pc -
1022e68589bSMark F. Adams */
1032e68589bSMark F. Adams #undef __FUNCT__
1042e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GEO"
1052e68589bSMark F. Adams PetscErrorCode PCSetFromOptions_GEO(PC pc)
1062e68589bSMark F. Adams {
1072e68589bSMark F. Adams   PetscErrorCode  ierr;
1082e68589bSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1092e68589bSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1102e68589bSMark F. Adams 
1112e68589bSMark F. Adams   PetscFunctionBegin;
1122e68589bSMark F. Adams   ierr = PetscOptionsHead("GAMG-GEO options");CHKERRQ(ierr);
1132e68589bSMark F. Adams   {
1142e68589bSMark F. Adams     /* -pc_gamg_sa_nsmooths */
1152e68589bSMark F. Adams     /* pc_gamg_sa->smooths = 0; */
1162e68589bSMark F. Adams     /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */
1172e68589bSMark F. Adams     /*                        "smoothing steps for smoothed aggregation, usually 1 (0)", */
1182e68589bSMark F. Adams     /*                        "PCGAMGSetNSmooths_AGG", */
1192e68589bSMark F. Adams     /*                        pc_gamg_sa->smooths, */
1202e68589bSMark F. Adams     /*                        &pc_gamg_sa->smooths, */
1212e68589bSMark F. Adams     /*                        &flag);  */
1222e68589bSMark F. Adams     /* CHKERRQ(ierr); */
1232e68589bSMark F. Adams   }
1242e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1252e68589bSMark F. Adams 
1262e68589bSMark F. Adams   /* call base class */
1272e68589bSMark F. Adams   ierr = PCSetFromOptions_GAMG(pc);CHKERRQ(ierr);
1282e68589bSMark F. Adams 
1292e68589bSMark F. Adams   if (pc_gamg->verbose) {
13041b27cdeSMark F. Adams     MPI_Comm  wcomm = ((PetscObject)pc)->comm;
13141b27cdeSMark F. Adams     PetscPrintf(wcomm,"[%d]%s done\n",0,__FUNCT__);
1322e68589bSMark F. Adams   }
1332e68589bSMark F. Adams 
1342e68589bSMark F. Adams   PetscFunctionReturn(0);
1352e68589bSMark F. Adams }
1362e68589bSMark F. Adams 
1372e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1382e68589bSMark F. Adams /*
1392e68589bSMark F. Adams  triangulateAndFormProl
1402e68589bSMark F. Adams 
1412e68589bSMark F. Adams    Input Parameter:
1422e68589bSMark F. Adams    . selected_2 - list of selected local ID, includes selected ghosts
143a2f3521dSMark F. Adams    . data_stride -
144a2f3521dSMark F. Adams    . coords[2*data_stride] - column vector of local coordinates w/ ghosts
145b43b03e9SMark F. Adams    . nselected_1 - selected IDs that go with base (1) graph
1460cbbd2e1SMark F. Adams    . clid_lid_1[nselected_1] - lids of selected (c) nodes   ???????????
1470cbbd2e1SMark F. Adams    . agg_lists_1 - list of aggregates
1482e68589bSMark F. Adams    . crsGID[selected.size()] - global index for prolongation operator
1492e68589bSMark F. Adams    . bs - block size
1502e68589bSMark F. Adams   Output Parameter:
1512e68589bSMark F. Adams    . a_Prol - prolongation operator
1522e68589bSMark F. Adams    . a_worst_best - measure of worst missed fine vertex, 0 is no misses
1532e68589bSMark F. Adams */
1542e68589bSMark F. Adams #undef __FUNCT__
1552e68589bSMark F. Adams #define __FUNCT__ "triangulateAndFormProl"
1560cbbd2e1SMark F. Adams static PetscErrorCode triangulateAndFormProl(IS  selected_2, /* list of selected local ID, includes selected ghosts */
157a2f3521dSMark F. Adams                                               const PetscInt data_stride,
1582e68589bSMark F. Adams                                               const PetscReal coords[], /* column vector of local coordinates w/ ghosts */
159b43b03e9SMark F. Adams                                               const PetscInt nselected_1, /* list of selected local ID, includes selected ghosts */
160b43b03e9SMark F. Adams                                               const PetscInt clid_lid_1[],
1610cbbd2e1SMark F. Adams                                               const PetscCoarsenData *agg_lists_1, /* selected_1 vertices of aggregate unselected vertices */
1622e68589bSMark F. Adams                                               const PetscInt crsGID[],
1632e68589bSMark F. Adams                                               const PetscInt bs,
1642e68589bSMark F. Adams                                               Mat a_Prol, /* prolongation operator (output) */
165*1147fc2aSKarl Rupp                                               PetscReal *a_worst_best) /* measure of worst missed fine vertex, 0 is no misses */
1662e68589bSMark F. Adams {
1672e68589bSMark F. Adams #if defined(PETSC_HAVE_TRIANGLE)
1682e68589bSMark F. Adams   PetscErrorCode       ierr;
169b43b03e9SMark F. Adams   PetscInt             jj,tid,tt,idx,nselected_2;
1702e68589bSMark F. Adams   struct triangulateio in,mid;
1710cbbd2e1SMark F. Adams   const PetscInt      *selected_idx_2;
1722e68589bSMark F. Adams   PetscMPIInt          mype,npe;
1732e68589bSMark F. Adams   PetscInt             Istart,Iend,nFineLoc,myFine0;
1742e68589bSMark F. Adams   int                  kk,nPlotPts,sid;
175c8b0795cSMark F. Adams   MPI_Comm             wcomm = ((PetscObject)a_Prol)->comm;
176c8b0795cSMark F. Adams   PetscReal            tm;
1772e68589bSMark F. Adams   PetscFunctionBegin;
178c8b0795cSMark F. Adams 
179c8b0795cSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
180c8b0795cSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
181c8b0795cSMark F. Adams   ierr = ISGetSize(selected_2, &nselected_2);CHKERRQ(ierr);
1822e68589bSMark F. Adams   if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */
1832e68589bSMark F. Adams     *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */
184806fa848SBarry Smith   } else *a_worst_best = 0.0;
185c8b0795cSMark F. Adams   ierr = MPI_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm);CHKERRQ(ierr);
186c8b0795cSMark F. Adams   if (tm > 0.0) {
187c8b0795cSMark F. Adams     *a_worst_best = 100.0;
1882e68589bSMark F. Adams     PetscFunctionReturn(0);
1892e68589bSMark F. Adams   }
1902e68589bSMark F. Adams   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
1912e68589bSMark F. Adams   nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs;
1922e68589bSMark F. Adams   nPlotPts = nFineLoc; /* locals */
1932e68589bSMark F. Adams   /* traingle */
1942e68589bSMark F. Adams   /* Define input points - in*/
1952e68589bSMark F. Adams   in.numberofpoints = nselected_2;
1962e68589bSMark F. Adams   in.numberofpointattributes = 0;
1972e68589bSMark F. Adams   /* get nselected points */
1982e68589bSMark F. Adams   ierr = PetscMalloc(2*(nselected_2)*sizeof(REAL), &in.pointlist);CHKERRQ(ierr);
1992e68589bSMark F. Adams   ierr = ISGetIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
2002e68589bSMark F. Adams 
2012e68589bSMark F. Adams   for (kk=0,sid=0;kk<nselected_2;kk++,sid += 2) {
2022e68589bSMark F. Adams     PetscInt lid = selected_idx_2[kk];
2032e68589bSMark F. Adams     in.pointlist[sid] = coords[lid];
204a2f3521dSMark F. Adams     in.pointlist[sid+1] = coords[data_stride + lid];
2052e68589bSMark F. Adams     if (lid>=nFineLoc) nPlotPts++;
2062e68589bSMark F. Adams   }
2072e68589bSMark F. Adams   assert(sid==2*nselected_2);
2082e68589bSMark F. Adams 
2092e68589bSMark F. Adams   in.numberofsegments = 0;
2102e68589bSMark F. Adams   in.numberofedges = 0;
2112e68589bSMark F. Adams   in.numberofholes = 0;
2122e68589bSMark F. Adams   in.numberofregions = 0;
2132e68589bSMark F. Adams   in.trianglelist = 0;
2142e68589bSMark F. Adams   in.segmentmarkerlist = 0;
2152e68589bSMark F. Adams   in.pointattributelist = 0;
2162e68589bSMark F. Adams   in.pointmarkerlist = 0;
2172e68589bSMark F. Adams   in.triangleattributelist = 0;
2182e68589bSMark F. Adams   in.trianglearealist = 0;
2192e68589bSMark F. Adams   in.segmentlist = 0;
2202e68589bSMark F. Adams   in.holelist = 0;
2212e68589bSMark F. Adams   in.regionlist = 0;
2222e68589bSMark F. Adams   in.edgelist = 0;
2232e68589bSMark F. Adams   in.edgemarkerlist = 0;
2242e68589bSMark F. Adams   in.normlist = 0;
2252e68589bSMark F. Adams   /* triangulate */
2262e68589bSMark F. Adams   mid.pointlist = 0;            /* Not needed if -N switch used. */
2272e68589bSMark F. Adams   /* Not needed if -N switch used or number of point attributes is zero: */
2282e68589bSMark F. Adams   mid.pointattributelist = 0;
2292e68589bSMark F. Adams   mid.pointmarkerlist = 0; /* Not needed if -N or -B switch used. */
2302e68589bSMark F. Adams   mid.trianglelist = 0;          /* Not needed if -E switch used. */
2312e68589bSMark F. Adams   /* Not needed if -E switch used or number of triangle attributes is zero: */
2322e68589bSMark F. Adams   mid.triangleattributelist = 0;
2332e68589bSMark F. Adams   mid.neighborlist = 0;         /* Needed only if -n switch used. */
2342e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P not used: */
2352e68589bSMark F. Adams   mid.segmentlist = 0;
2362e68589bSMark F. Adams   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
2372e68589bSMark F. Adams   mid.segmentmarkerlist = 0;
2382e68589bSMark F. Adams   mid.edgelist = 0;             /* Needed only if -e switch used. */
2392e68589bSMark F. Adams   mid.edgemarkerlist = 0;   /* Needed if -e used and -B not used. */
2402e68589bSMark F. Adams   mid.numberoftriangles = 0;
2412e68589bSMark F. Adams 
2422e68589bSMark F. Adams   /* Triangulate the points.  Switches are chosen to read and write a  */
2432e68589bSMark F. Adams   /*   PSLG (p), preserve the convex hull (c), number everything from  */
2442e68589bSMark F. Adams   /*   zero (z), assign a regional attribute to each element (A), and  */
2452e68589bSMark F. Adams   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
2462e68589bSMark F. Adams   /*   neighbor list (n).                                            */
2472e68589bSMark F. Adams   if (nselected_2 != 0) { /* inactive processor */
2482e68589bSMark F. Adams     char args[] = "npczQ"; /* c is needed ? */
2492e68589bSMark F. Adams     triangulate(args, &in, &mid, (struct triangulateio *) NULL);
2502e68589bSMark F. Adams     /* output .poly files for 'showme' */
2512e68589bSMark F. Adams     if (!PETSC_TRUE) {
2522e68589bSMark F. Adams       static int level = 1;
2532e68589bSMark F. Adams       FILE *file; char fname[32];
2542e68589bSMark F. Adams 
2552e68589bSMark F. Adams       sprintf(fname,"C%d_%d.poly",level,mype); file = fopen(fname, "w");
2562e68589bSMark F. Adams       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
2572e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0);
2582e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2592e68589bSMark F. Adams       for (kk=0,sid=0;kk<in.numberofpoints;kk++,sid += 2) {
2602e68589bSMark F. Adams         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
2612e68589bSMark F. Adams       }
2622e68589bSMark F. Adams       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
2632e68589bSMark F. Adams       fprintf(file, "%d  %d\n",0,0);
2642e68589bSMark F. Adams       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
2652e68589bSMark F. Adams       /* One line: <# of holes> */
2662e68589bSMark F. Adams       fprintf(file, "%d\n",0);
2672e68589bSMark F. Adams       /* Following lines: <hole #> <x> <y> */
2682e68589bSMark F. Adams       /* Optional line: <# of regional attributes and/or area constraints> */
2692e68589bSMark F. Adams       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
2702e68589bSMark F. Adams       fclose(file);
2712e68589bSMark F. Adams 
2722e68589bSMark F. Adams       /* elems */
2732e68589bSMark F. Adams       sprintf(fname,"C%d_%d.ele",level,mype); file = fopen(fname, "w");
2742e68589bSMark F. Adams       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
2752e68589bSMark F. Adams       fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0);
2762e68589bSMark F. Adams       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
2772e68589bSMark F. Adams       for (kk=0,sid=0;kk<mid.numberoftriangles;kk++,sid += 3) {
2782e68589bSMark F. Adams         fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]);
2792e68589bSMark F. Adams       }
2802e68589bSMark F. Adams       fclose(file);
2812e68589bSMark F. Adams 
2822e68589bSMark F. Adams       sprintf(fname,"C%d_%d.node",level,mype); file = fopen(fname, "w");
2832e68589bSMark F. Adams       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
2842e68589bSMark F. Adams       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
2852e68589bSMark F. Adams       fprintf(file, "%d  %d  %d  %d\n",nPlotPts,2,0,0);
2862e68589bSMark F. Adams       /*Following lines: <vertex #> <x> <y> */
2872e68589bSMark F. Adams       for (kk=0,sid=0;kk<in.numberofpoints;kk++,sid+=2) {
2882e68589bSMark F. Adams         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
2892e68589bSMark F. Adams       }
2902e68589bSMark F. Adams 
2912e68589bSMark F. Adams       sid /= 2;
2922e68589bSMark F. Adams       for (jj=0;jj<nFineLoc;jj++) {
2932e68589bSMark F. Adams         PetscBool sel = PETSC_TRUE;
2942e68589bSMark F. Adams         for (kk=0 ; kk<nselected_2 && sel ; kk++) {
2952e68589bSMark F. Adams           PetscInt lid = selected_idx_2[kk];
2962e68589bSMark F. Adams           if (lid == jj) sel = PETSC_FALSE;
2972e68589bSMark F. Adams         }
2982e68589bSMark F. Adams         if (sel) {
299a2f3521dSMark F. Adams           fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]);
3002e68589bSMark F. Adams         }
3012e68589bSMark F. Adams       }
3022e68589bSMark F. Adams       fclose(file);
3032e68589bSMark F. Adams       assert(sid==nPlotPts);
3042e68589bSMark F. Adams       level++;
3052e68589bSMark F. Adams     }
3062e68589bSMark F. Adams   }
3070cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
3080cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
3092e68589bSMark F. Adams #endif
3102e68589bSMark F. Adams   { /* form P - setup some maps */
3110cbbd2e1SMark F. Adams     PetscInt clid,mm,*nTri,*node_tri;
3122e68589bSMark F. Adams 
3132e68589bSMark F. Adams     ierr = PetscMalloc(nselected_2*sizeof(PetscInt), &node_tri);CHKERRQ(ierr);
3142e68589bSMark F. Adams     ierr = PetscMalloc(nselected_2*sizeof(PetscInt), &nTri);CHKERRQ(ierr);
3152e68589bSMark F. Adams 
3162e68589bSMark F. Adams     /* need list of triangles on node */
3172e68589bSMark F. Adams     for (kk=0;kk<nselected_2;kk++) nTri[kk] = 0;
3182e68589bSMark F. Adams     for (tid=0,kk=0;tid<mid.numberoftriangles;tid++) {
3192e68589bSMark F. Adams       for (jj=0;jj<3;jj++) {
3202e68589bSMark F. Adams         PetscInt cid = mid.trianglelist[kk++];
3212e68589bSMark F. Adams         if (nTri[cid] == 0) node_tri[cid] = tid;
3222e68589bSMark F. Adams         nTri[cid]++;
3232e68589bSMark F. Adams       }
3242e68589bSMark F. Adams     }
3252e68589bSMark F. Adams #define EPS 1.e-12
3262e68589bSMark F. Adams     /* find points and set prolongation */
3270cbbd2e1SMark F. Adams     for (mm = clid = 0 ; mm < nFineLoc ; mm++) {
328e78576d6SMark F. Adams       PetscBool ise;
329e78576d6SMark F. Adams       ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise);CHKERRQ(ierr);
330e78576d6SMark F. Adams       if (!ise) {
3310cbbd2e1SMark F. Adams         const PetscInt lid = mm;
332c41eeb4cSKarl Rupp         /* for (clid_iterator=0;clid_iterator<nselected_1;clid_iterator++) { */
333c41eeb4cSKarl Rupp         /* PetscInt flid = clid_lid_1[clid_iterator]; assert(flid != -1); */
3342e68589bSMark F. Adams         PetscScalar AA[3][3];
3352e68589bSMark F. Adams         PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO;
33641b27cdeSMark F. Adams         PetscCDPos         pos;
337e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
338e78576d6SMark F. Adams         while (pos) {
339e78576d6SMark F. Adams           PetscInt flid;
340ffc955d6SMark F. Adams           ierr = PetscLLNGetID(pos, &flid);CHKERRQ(ierr);
341e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
3420cbbd2e1SMark F. Adams 
3432e68589bSMark F. Adams           if (flid < nFineLoc) {  /* could be a ghost */
3442e68589bSMark F. Adams             PetscInt bestTID = -1; PetscReal best_alpha = 1.e10;
3452e68589bSMark F. Adams             const PetscInt fgid = flid + myFine0;
3462e68589bSMark F. Adams             /* compute shape function for gid */
347a2f3521dSMark F. Adams             const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0};
3482e68589bSMark F. Adams             PetscBool haveit=PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3];
3492e68589bSMark F. Adams             /* look for it */
3500cbbd2e1SMark F. Adams             for (tid = node_tri[clid], jj=0;
3512e68589bSMark F. Adams                  jj < 5 && !haveit && tid != -1;
3522e68589bSMark F. Adams                  jj++) {
3532e68589bSMark F. Adams               for (tt=0;tt<3;tt++) {
3542e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3*tid + tt];
3552e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
356a2f3521dSMark F. Adams                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
3572e68589bSMark F. Adams                 clids[tt] = cid2; /* store for interp */
3582e68589bSMark F. Adams               }
3592e68589bSMark F. Adams 
3602e68589bSMark F. Adams               for (tt=0;tt<3;tt++) alpha[tt] = (PetscScalar)fcoord[tt];
3612e68589bSMark F. Adams 
3622e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
3632e68589bSMark F. Adams               LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
3642e68589bSMark F. Adams               {
3652e68589bSMark F. Adams                 PetscBool have=PETSC_TRUE;  PetscReal lowest=1.e10;
3662e68589bSMark F. Adams                 for (tt = 0, idx = 0 ; tt < 3 ; tt++) {
3672e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
3682e68589bSMark F. Adams                   if (PetscRealPart(alpha[tt]) < lowest) {
3692e68589bSMark F. Adams                     lowest = PetscRealPart(alpha[tt]);
3702e68589bSMark F. Adams                     idx = tt;
3712e68589bSMark F. Adams                   }
3722e68589bSMark F. Adams                 }
3732e68589bSMark F. Adams                 haveit = have;
3742e68589bSMark F. Adams               }
3752e68589bSMark F. Adams               tid = mid.neighborlist[3*tid + idx];
3762e68589bSMark F. Adams             }
3772e68589bSMark F. Adams 
3782e68589bSMark F. Adams             if (!haveit) {
3792e68589bSMark F. Adams               /* brute force */
3802e68589bSMark F. Adams               for (tid=0 ; tid<mid.numberoftriangles && !haveit ; tid++) {
3812e68589bSMark F. Adams                 for (tt=0;tt<3;tt++) {
3822e68589bSMark F. Adams                   PetscInt cid2 = mid.trianglelist[3*tid + 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) */
3892e68589bSMark F. Adams                 LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
3902e68589bSMark F. Adams                 {
3912e68589bSMark F. Adams                   PetscBool have=PETSC_TRUE;  PetscReal worst=0.0, v;
3922e68589bSMark F. Adams                   for (tt=0; tt<3 && have ;tt++) {
3932e68589bSMark F. Adams                     if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE;
3942e68589bSMark F. Adams                     if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v;
3952e68589bSMark F. Adams                   }
3962e68589bSMark F. Adams                   if (worst < best_alpha) {
3972e68589bSMark F. Adams                     best_alpha = worst; bestTID = tid;
3982e68589bSMark F. Adams                   }
3992e68589bSMark F. Adams                   haveit = have;
4002e68589bSMark F. Adams                 }
4012e68589bSMark F. Adams               }
4022e68589bSMark F. Adams             }
4032e68589bSMark F. Adams             if (!haveit) {
4042e68589bSMark F. Adams               if (best_alpha > *a_worst_best) *a_worst_best = best_alpha;
4052e68589bSMark F. Adams               /* use best one */
4062e68589bSMark F. Adams               for (tt=0;tt<3;tt++) {
4072e68589bSMark F. Adams                 PetscInt cid2 = mid.trianglelist[3*bestTID + tt];
4082e68589bSMark F. Adams                 PetscInt lid2 = selected_idx_2[cid2];
409a2f3521dSMark F. Adams                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
4102e68589bSMark F. Adams                 clids[tt] = cid2; /* store for interp */
4112e68589bSMark F. Adams               }
4122e68589bSMark F. Adams               for (tt=0;tt<3;tt++) alpha[tt] = fcoord[tt];
4132e68589bSMark F. Adams               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
4142e68589bSMark F. Adams               LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
4152e68589bSMark F. Adams             }
4162e68589bSMark F. Adams 
4172e68589bSMark F. Adams             /* put in row of P */
4182e68589bSMark F. Adams             for (idx=0;idx<3;idx++) {
4192e68589bSMark F. Adams               PetscScalar shp = alpha[idx];
4202e68589bSMark F. Adams               if (PetscAbs(PetscRealPart(shp)) > 1.e-6) {
4212e68589bSMark F. Adams                 PetscInt cgid = crsGID[clids[idx]];
4222e68589bSMark F. Adams                 PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */
4232e68589bSMark F. Adams                 for (tt=0 ; tt < bs ; tt++, ii++, jj++) {
4242e68589bSMark F. Adams                   ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES);CHKERRQ(ierr);
4252e68589bSMark F. Adams                 }
4262e68589bSMark F. Adams               }
4272e68589bSMark F. Adams             }
4282e68589bSMark F. Adams           }
4290cbbd2e1SMark F. Adams         } /* aggregates iterations */
4300cbbd2e1SMark F. Adams         clid++;
4310cbbd2e1SMark F. Adams       } /* a coarse agg */
4320cbbd2e1SMark F. Adams     } /* for all fine nodes */
4330cbbd2e1SMark F. Adams 
4342e68589bSMark F. Adams     ierr = ISRestoreIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
4352e68589bSMark F. Adams     ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4362e68589bSMark F. Adams     ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4372e68589bSMark F. Adams 
4382e68589bSMark F. Adams     ierr = PetscFree(node_tri);CHKERRQ(ierr);
4392e68589bSMark F. Adams     ierr = PetscFree(nTri);CHKERRQ(ierr);
4402e68589bSMark F. Adams   }
4410cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
4420cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
4432e68589bSMark F. Adams #endif
4442e68589bSMark F. Adams   free(mid.trianglelist);
4452e68589bSMark F. Adams   free(mid.neighborlist);
4462e68589bSMark F. Adams   ierr = PetscFree(in.pointlist);CHKERRQ(ierr);
4472e68589bSMark F. Adams 
4482e68589bSMark F. Adams   PetscFunctionReturn(0);
4492e68589bSMark F. Adams #else
450c666adbfSMark F. Adams   SETERRQ(((PetscObject)a_Prol)->comm,PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG");
4512e68589bSMark F. Adams #endif
4522e68589bSMark F. Adams }
4532e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
4542e68589bSMark F. Adams /*
4552e68589bSMark F. Adams    getGIDsOnSquareGraph - square graph, get
4562e68589bSMark F. Adams 
4572e68589bSMark F. Adams    Input Parameter:
4580cbbd2e1SMark F. Adams    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
459b43b03e9SMark F. Adams    . clid_lid_1 - [nselected_1] lids of selected nodes
4602e68589bSMark F. Adams    . Gmat1 - graph that goes with 'selected_1'
4612e68589bSMark F. Adams    Output Parameter:
4622e68589bSMark F. Adams    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
4632e68589bSMark F. Adams    . a_Gmat_2 - graph that is squared of 'Gmat_1'
4642e68589bSMark F. Adams    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
4652e68589bSMark F. Adams */
4662e68589bSMark F. Adams #undef __FUNCT__
4672e68589bSMark F. Adams #define __FUNCT__ "getGIDsOnSquareGraph"
4680cbbd2e1SMark F. Adams static PetscErrorCode getGIDsOnSquareGraph(const PetscInt nselected_1,
469b43b03e9SMark F. Adams                                             const PetscInt clid_lid_1[],
4702e68589bSMark F. Adams                                             const Mat Gmat1,
4712e68589bSMark F. Adams                                             IS *a_selected_2,
4722e68589bSMark F. Adams                                             Mat *a_Gmat_2,
473*1147fc2aSKarl Rupp                                             PetscInt **a_crsGID)
4742e68589bSMark F. Adams {
4752e68589bSMark F. Adams   PetscErrorCode ierr;
4762e68589bSMark F. Adams   PetscMPIInt    mype,npe;
477b43b03e9SMark F. Adams   PetscInt       *crsGID, kk,my0,Iend,nloc;
4782e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Gmat1)->comm;
4792e68589bSMark F. Adams 
4802e68589bSMark F. Adams   PetscFunctionBegin;
4812e68589bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
4822e68589bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
4832e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend);CHKERRQ(ierr); /* AIJ */
4842e68589bSMark F. Adams   nloc = Iend - my0; /* this does not change */
4852e68589bSMark F. Adams 
4862e68589bSMark F. Adams   if (npe == 1) { /* not much to do in serial */
487b43b03e9SMark F. Adams     ierr = PetscMalloc(nselected_1*sizeof(PetscInt), &crsGID);CHKERRQ(ierr);
488b43b03e9SMark F. Adams     for (kk=0;kk<nselected_1;kk++) crsGID[kk] = kk;
4892e68589bSMark F. Adams     *a_Gmat_2 = 0;
490806fa848SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);CHKERRQ(ierr);
491806fa848SBarry Smith   } else {
492b43b03e9SMark F. Adams     PetscInt      idx,num_fine_ghosts,num_crs_ghost,myCrs0;
4932e68589bSMark F. Adams     Mat_MPIAIJ   *mpimat2;
4942e68589bSMark F. Adams     Mat           Gmat2;
4952e68589bSMark F. Adams     Vec           locState;
4962e68589bSMark F. Adams     PetscScalar   *cpcol_state;
4972e68589bSMark F. Adams 
4982e68589bSMark F. Adams     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
499b43b03e9SMark F. Adams     kk = nselected_1;
500b43b03e9SMark F. Adams     MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPIU_SUM, wcomm);
501b43b03e9SMark F. Adams     myCrs0 -= nselected_1;
5022e68589bSMark F. Adams 
503b43b03e9SMark F. Adams     if (a_Gmat_2) { /* output */
5042e68589bSMark F. Adams       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
5052e68589bSMark F. Adams       ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
5062e68589bSMark F. Adams       *a_Gmat_2 = Gmat2; /* output */
507806fa848SBarry Smith     } else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
5082e68589bSMark F. Adams     /* get coarse grid GIDS for selected (locals and ghosts) */
5092e68589bSMark F. Adams     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
5102e68589bSMark F. Adams     ierr = MatGetVecs(Gmat2, &locState, 0);CHKERRQ(ierr);
5112e68589bSMark F. Adams     ierr = VecSet(locState, (PetscScalar)(PetscReal)(-1));CHKERRQ(ierr); /* set with UNKNOWN state */
512b43b03e9SMark F. Adams     for (kk=0;kk<nselected_1;kk++) {
513b43b03e9SMark F. Adams       PetscInt fgid = clid_lid_1[kk] + my0;
5142e68589bSMark F. Adams       PetscScalar v = (PetscScalar)(kk+myCrs0);
5152e68589bSMark F. Adams       ierr = VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES);CHKERRQ(ierr); /* set with PID */
5162e68589bSMark F. Adams     }
5172e68589bSMark F. Adams     ierr = VecAssemblyBegin(locState);CHKERRQ(ierr);
5182e68589bSMark F. Adams     ierr = VecAssemblyEnd(locState);CHKERRQ(ierr);
5192e68589bSMark F. Adams     ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5202e68589bSMark F. Adams     ierr =   VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5212e68589bSMark F. Adams     ierr = VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts);CHKERRQ(ierr);
5222e68589bSMark F. Adams     ierr = VecGetArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
5232e68589bSMark F. Adams     for (kk=0,num_crs_ghost=0;kk<num_fine_ghosts;kk++) {
5242e68589bSMark F. Adams       if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++;
5252e68589bSMark F. Adams     }
526b43b03e9SMark F. Adams     ierr = PetscMalloc((nselected_1+num_crs_ghost)*sizeof(PetscInt), &crsGID);CHKERRQ(ierr); /* output */
5272e68589bSMark F. Adams     {
5282e68589bSMark F. Adams       PetscInt *selected_set;
529b43b03e9SMark F. Adams       ierr = PetscMalloc((nselected_1+num_crs_ghost)*sizeof(PetscInt), &selected_set);CHKERRQ(ierr);
5302e68589bSMark F. Adams       /* do ghost of 'crsGID' */
531b43b03e9SMark F. Adams       for (kk=0,idx=nselected_1;kk<num_fine_ghosts;kk++) {
5322e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
5332e68589bSMark F. Adams           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
5342e68589bSMark F. Adams           selected_set[idx] = nloc + kk;
5352e68589bSMark F. Adams           crsGID[idx++] = cgid;
5362e68589bSMark F. Adams         }
5372e68589bSMark F. Adams       }
538b43b03e9SMark F. Adams       assert(idx==(nselected_1+num_crs_ghost));
5392e68589bSMark F. Adams       ierr = VecRestoreArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
5402e68589bSMark F. Adams       /* do locals in 'crsGID' */
5412e68589bSMark F. Adams       ierr = VecGetArray(locState, &cpcol_state);CHKERRQ(ierr);
5422e68589bSMark F. Adams       for (kk=0,idx=0;kk<nloc;kk++) {
5432e68589bSMark F. Adams         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
5442e68589bSMark F. Adams           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
5452e68589bSMark F. Adams           selected_set[idx] = kk;
5462e68589bSMark F. Adams           crsGID[idx++] = cgid;
5472e68589bSMark F. Adams         }
5482e68589bSMark F. Adams       }
549b43b03e9SMark F. Adams       assert(idx==nselected_1);
5502e68589bSMark F. Adams       ierr = VecRestoreArray(locState, &cpcol_state);CHKERRQ(ierr);
5512e68589bSMark F. Adams 
5522e68589bSMark F. Adams       if (a_selected_2 != 0) { /* output */
553806fa848SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);CHKERRQ(ierr);
554806fa848SBarry Smith       } else {
5552e68589bSMark F. Adams         ierr = PetscFree(selected_set);CHKERRQ(ierr);
5562e68589bSMark F. Adams       }
5570cbbd2e1SMark F. Adams     }
5582e68589bSMark F. Adams     ierr = VecDestroy(&locState);CHKERRQ(ierr);
5592e68589bSMark F. Adams   }
5602e68589bSMark F. Adams   *a_crsGID = crsGID; /* output */
5612e68589bSMark F. Adams 
5622e68589bSMark F. Adams   PetscFunctionReturn(0);
5632e68589bSMark F. Adams }
5642e68589bSMark F. Adams 
5652e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
5662e68589bSMark F. Adams /*
567c8b0795cSMark F. Adams    PCGAMGgraph_GEO
5682e68589bSMark F. Adams 
5692e68589bSMark F. Adams   Input Parameter:
5702e68589bSMark F. Adams    . pc - this
5712e68589bSMark F. Adams    . Amat - matrix on this fine level
5722e68589bSMark F. Adams   Output Parameter:
573c8b0795cSMark F. Adams    . a_Gmat
5742e68589bSMark F. Adams */
5752e68589bSMark F. Adams #undef __FUNCT__
576c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_GEO"
577c8b0795cSMark F. Adams PetscErrorCode PCGAMGgraph_GEO(PC pc,
5782e68589bSMark F. Adams                                 const Mat Amat,
579*1147fc2aSKarl Rupp                                 Mat *a_Gmat)
580c8b0795cSMark F. Adams {
581c8b0795cSMark F. Adams   PetscErrorCode ierr;
582c8b0795cSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
583c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
584c8b0795cSMark F. Adams   const PetscInt verbose = pc_gamg->verbose;
585c8b0795cSMark F. Adams   const PetscReal vfilter = pc_gamg->threshold;
586c8b0795cSMark F. Adams   PetscMPIInt    mype,npe;
587c8b0795cSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
588c8b0795cSMark F. Adams   Mat            Gmat;
5890cbbd2e1SMark F. Adams   PetscBool  set,flg,symm;
590c8b0795cSMark F. Adams   PetscFunctionBegin;
5910cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
5920cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr);
5930cbbd2e1SMark F. Adams #endif
594c8b0795cSMark F. Adams   ierr = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr);
595c8b0795cSMark F. Adams   ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr);
596c8b0795cSMark F. Adams 
5970cbbd2e1SMark F. Adams   ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr);
598263489e9SJed Brown   symm = (PetscBool)!(set && flg);
5990cbbd2e1SMark F. Adams 
6002d7fac45SMark F. Adams   ierr  = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
6012d7fac45SMark F. Adams   ierr  = PCGAMGFilterGraph(&Gmat, vfilter, symm, verbose);CHKERRQ(ierr);
602c8b0795cSMark F. Adams 
603c8b0795cSMark F. Adams   *a_Gmat = Gmat;
6040cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
6050cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr);
6060cbbd2e1SMark F. Adams #endif
607c8b0795cSMark F. Adams   PetscFunctionReturn(0);
608c8b0795cSMark F. Adams }
609c8b0795cSMark F. Adams 
610c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
611c8b0795cSMark F. Adams /*
612c8b0795cSMark F. Adams    PCGAMGcoarsen_GEO
613c8b0795cSMark F. Adams 
614c8b0795cSMark F. Adams   Input Parameter:
615e0940f08SMark F. Adams    . a_pc - this
616e0940f08SMark F. Adams    . a_Gmat - graph
617c8b0795cSMark F. Adams   Output Parameter:
618c8b0795cSMark F. Adams    . a_llist_parent - linked list from selected indices for data locality only
619c8b0795cSMark F. Adams */
620c8b0795cSMark F. Adams #undef __FUNCT__
621c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGcoarsen_GEO"
622e0940f08SMark F. Adams PetscErrorCode PCGAMGcoarsen_GEO(PC a_pc,
623e0940f08SMark F. Adams                                   Mat *a_Gmat,
624*1147fc2aSKarl Rupp                                   PetscCoarsenData **a_llist_parent)
625c8b0795cSMark F. Adams {
626c8b0795cSMark F. Adams   PetscErrorCode ierr;
627e0940f08SMark F. Adams   PC_MG          *mg = (PC_MG*)a_pc->data;
628c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
629c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,kk,Ii,ncols;
630c8b0795cSMark F. Adams   PetscMPIInt    mype,npe;
6310cbbd2e1SMark F. Adams   IS             perm;
632c8b0795cSMark F. Adams   GAMGNode *gnodes;
633c8b0795cSMark F. Adams   PetscInt *permute;
634e0940f08SMark F. Adams   Mat       Gmat = *a_Gmat;
635c8b0795cSMark F. Adams   MPI_Comm  wcomm = ((PetscObject)Gmat)->comm;
636b43b03e9SMark F. Adams   MatCoarsen crs;
637c8b0795cSMark F. Adams 
638c8b0795cSMark F. Adams   PetscFunctionBegin;
6390cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
6400cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
6410cbbd2e1SMark F. Adams #endif
642c8b0795cSMark F. Adams   ierr = MPI_Comm_rank(wcomm, &mype);CHKERRQ(ierr);
643c8b0795cSMark F. Adams   ierr = MPI_Comm_size(wcomm, &npe);CHKERRQ(ierr);
644c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr);
645c8b0795cSMark F. Adams   nloc = (Iend-Istart);
646c8b0795cSMark F. Adams 
647c8b0795cSMark F. Adams   /* create random permutation with sort for geo-mg */
648c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(GAMGNode), &gnodes);CHKERRQ(ierr);
649c8b0795cSMark F. Adams   ierr = PetscMalloc(nloc*sizeof(PetscInt), &permute);CHKERRQ(ierr);
650c8b0795cSMark F. Adams 
651c8b0795cSMark F. Adams   for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */
652c8b0795cSMark F. Adams     ierr = MatGetRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr);
653c8b0795cSMark F. Adams     {
654c8b0795cSMark F. Adams       PetscInt lid = Ii - Istart;
655c8b0795cSMark F. Adams       gnodes[lid].lid = lid;
656c8b0795cSMark F. Adams       gnodes[lid].degree = ncols;
657c8b0795cSMark F. Adams     }
658c8b0795cSMark F. Adams     ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr);
659c8b0795cSMark F. Adams   }
660c8b0795cSMark F. Adams   /* randomize */
661c8b0795cSMark F. Adams   srand(1); /* make deterministic */
662c8b0795cSMark F. Adams   if (PETSC_TRUE) {
663c8b0795cSMark F. Adams     PetscBool *bIndexSet;
664c8b0795cSMark F. Adams     ierr = PetscMalloc(nloc*sizeof(PetscBool), &bIndexSet);CHKERRQ(ierr);
665c8b0795cSMark F. Adams     for (Ii = 0; Ii < nloc ; Ii++) bIndexSet[Ii] = PETSC_FALSE;
666c8b0795cSMark F. Adams     for (Ii = 0; Ii < nloc ; Ii++)
667c8b0795cSMark F. Adams     {
668c8b0795cSMark F. Adams       PetscInt iSwapIndex = rand()%nloc;
669c8b0795cSMark F. Adams       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii)
670c8b0795cSMark F. Adams       {
671c8b0795cSMark F. Adams         GAMGNode iTemp = gnodes[iSwapIndex];
672c8b0795cSMark F. Adams         gnodes[iSwapIndex] = gnodes[Ii];
673c8b0795cSMark F. Adams         gnodes[Ii] = iTemp;
674c8b0795cSMark F. Adams         bIndexSet[Ii] = PETSC_TRUE;
675c8b0795cSMark F. Adams         bIndexSet[iSwapIndex] = PETSC_TRUE;
676c8b0795cSMark F. Adams       }
677c8b0795cSMark F. Adams     }
678c8b0795cSMark F. Adams     ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
679c8b0795cSMark F. Adams   }
680c8b0795cSMark F. Adams   /* only sort locals */
6810cbbd2e1SMark F. Adams   qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare);
682c8b0795cSMark F. Adams   /* create IS of permutation */
683c8b0795cSMark F. Adams   for (kk=0;kk<nloc;kk++) { /* locals only */
684c8b0795cSMark F. Adams     permute[kk] = gnodes[kk].lid;
685c8b0795cSMark F. Adams   }
686806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr);
687c8b0795cSMark F. Adams 
688c8b0795cSMark F. Adams   ierr = PetscFree(gnodes);CHKERRQ(ierr);
689c8b0795cSMark F. Adams 
690c8b0795cSMark F. Adams   /* get MIS aggs */
691b43b03e9SMark F. Adams 
692b43b03e9SMark F. Adams   ierr = MatCoarsenCreate(wcomm, &crs);CHKERRQ(ierr);
693b43b03e9SMark F. Adams   ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr);
694b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
695b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr);
696b43b03e9SMark F. Adams   ierr = MatCoarsenSetVerbose(crs, pc_gamg->verbose);CHKERRQ(ierr);
697b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_FALSE);CHKERRQ(ierr);
698b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
6990cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, a_llist_parent);CHKERRQ(ierr);
700b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
701c8b0795cSMark F. Adams 
702c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
7030cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
7040cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
7050cbbd2e1SMark F. Adams #endif
706c8b0795cSMark F. Adams   PetscFunctionReturn(0);
707c8b0795cSMark F. Adams }
708c8b0795cSMark F. Adams 
709c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
710c8b0795cSMark F. Adams /*
7110cbbd2e1SMark F. Adams  PCGAMGProlongator_GEO
712c8b0795cSMark F. Adams 
713c8b0795cSMark F. Adams  Input Parameter:
714c8b0795cSMark F. Adams  . pc - this
715c8b0795cSMark F. Adams  . Amat - matrix on this fine level
716c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
7170cbbd2e1SMark F. Adams  . selected_1 - [nselected]
7180cbbd2e1SMark F. Adams  . agg_lists - [nselected]
719c8b0795cSMark F. Adams  Output Parameter:
720c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
721c8b0795cSMark F. Adams  */
722c8b0795cSMark F. Adams #undef __FUNCT__
7230cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_GEO"
7240cbbd2e1SMark F. Adams PetscErrorCode PCGAMGProlongator_GEO(PC pc,
725c8b0795cSMark F. Adams                                       const Mat Amat,
726c8b0795cSMark F. Adams                                       const Mat Gmat,
7270cbbd2e1SMark F. Adams                                       PetscCoarsenData *agg_lists,
728*1147fc2aSKarl Rupp                                       Mat *a_P_out)
7292e68589bSMark F. Adams {
7302e68589bSMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
7312e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
7322e68589bSMark F. Adams   const PetscInt  verbose = pc_gamg->verbose;
733c8b0795cSMark F. Adams   const PetscInt  dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
7342e68589bSMark F. Adams   PetscErrorCode ierr;
735b43b03e9SMark F. Adams   PetscInt       Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid;
736c8b0795cSMark F. Adams   Mat            Prol;
7372e68589bSMark F. Adams   PetscMPIInt    mype, npe;
7382e68589bSMark F. Adams   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
7390cbbd2e1SMark F. Adams   IS             selected_2,selected_1;
7402e68589bSMark F. Adams   const PetscInt *selected_idx;
7412e68589bSMark F. Adams 
7422e68589bSMark F. Adams   PetscFunctionBegin;
7430cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
7440cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
7450cbbd2e1SMark F. Adams #endif
7462e68589bSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
7472e68589bSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
7482e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
7492e68589bSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
750b43b03e9SMark F. Adams   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
7512e68589bSMark F. Adams 
7522e68589bSMark F. Adams   /* get 'nLocalSelected' */
75341b27cdeSMark F. Adams   ierr = PetscCDGetMIS(agg_lists, &selected_1);CHKERRQ(ierr);
754b43b03e9SMark F. Adams   ierr = ISGetSize(selected_1, &jj);CHKERRQ(ierr);
755b43b03e9SMark F. Adams   ierr = PetscMalloc(jj*sizeof(PetscInt), &clid_flid);CHKERRQ(ierr);
7562e68589bSMark F. Adams   ierr = ISGetIndices(selected_1, &selected_idx);CHKERRQ(ierr);
757c8b0795cSMark F. Adams   for (kk=0,nLocalSelected=0;kk<jj;kk++) {
7582e68589bSMark F. Adams     PetscInt lid = selected_idx[kk];
759b43b03e9SMark F. Adams     if (lid<nloc) {
7600cbbd2e1SMark F. Adams       ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr);
761b43b03e9SMark F. Adams       if (ncols>1) { /* fiter out singletons */
762b43b03e9SMark F. Adams         clid_flid[nLocalSelected++] = lid;
763806fa848SBarry Smith       } else assert(0); /* filtered in coarsening */
7640cbbd2e1SMark F. Adams       ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr);
765b43b03e9SMark F. Adams     }
7662e68589bSMark F. Adams   }
7672e68589bSMark F. Adams   ierr = ISRestoreIndices(selected_1, &selected_idx);CHKERRQ(ierr);
768a2f3521dSMark F. Adams   ierr = ISDestroy(&selected_1);CHKERRQ(ierr); /* this is selected_1 in serial */
7692e68589bSMark F. Adams 
7702e68589bSMark F. Adams   /* create prolongator, create P matrix */
771a2f3521dSMark F. Adams   ierr = MatCreate(wcomm, &Prol);CHKERRQ(ierr);
772806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
773a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, bs);CHKERRQ(ierr);
774a2f3521dSMark F. Adams   ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr);
775a2f3521dSMark F. Adams   ierr = MatSeqAIJSetPreallocation(Prol,3*data_cols,PETSC_NULL);CHKERRQ(ierr);
776a2f3521dSMark F. Adams   ierr = MatMPIAIJSetPreallocation(Prol,3*data_cols,PETSC_NULL,3*data_cols,PETSC_NULL);CHKERRQ(ierr);
777a2f3521dSMark F. Adams   /* ierr = MatCreateAIJ(wcomm,  */
778a2f3521dSMark F. Adams   /*                      nloc*bs, nLocalSelected*bs, */
779a2f3521dSMark F. Adams   /*                      PETSC_DETERMINE, PETSC_DETERMINE, */
780a2f3521dSMark F. Adams   /*                      3*data_cols, PETSC_NULL,  */
781a2f3521dSMark F. Adams   /*                      3*data_cols, PETSC_NULL, */
782a2f3521dSMark F. Adams   /*                      &Prol); */
783a2f3521dSMark F. Adams   /* CHKERRQ(ierr); */
7842e68589bSMark F. Adams 
785c8b0795cSMark F. Adams   /* can get all points "removed" - but not on geomg */
7862e68589bSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &jj);CHKERRQ(ierr);
7872e68589bSMark F. Adams   if (jj==0) {
7882e68589bSMark F. Adams     if (verbose) {
78941b27cdeSMark F. Adams       PetscPrintf(wcomm,"[%d]%s ERROE: no selected points on coarse grid\n",mype,__FUNCT__);
7902e68589bSMark F. Adams     }
791b43b03e9SMark F. Adams     ierr = PetscFree(clid_flid);CHKERRQ(ierr);
7922e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
7932e68589bSMark F. Adams     *a_P_out = PETSC_NULL;  /* out */
7942e68589bSMark F. Adams     PetscFunctionReturn(0);
7952e68589bSMark F. Adams   }
7962e68589bSMark F. Adams 
7972e68589bSMark F. Adams   {
7982e68589bSMark F. Adams     PetscReal *coords;
799a2f3521dSMark F. Adams     PetscInt   data_stride;
800ffc955d6SMark F. Adams     PetscInt  *crsGID = PETSC_NULL;
8012e68589bSMark F. Adams     Mat        Gmat2;
8022e68589bSMark F. Adams 
8032e68589bSMark F. Adams     assert(dim==data_cols);
8042e68589bSMark F. Adams     /* grow ghost data for better coarse grid cover of fine grid */
8050cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
8060cbbd2e1SMark F. Adams     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
8072e68589bSMark F. Adams #endif
808a2f3521dSMark F. Adams     /* messy method, squares graph and gets some data */
809806fa848SBarry Smith     ierr = getGIDsOnSquareGraph(nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID);CHKERRQ(ierr);
8102e68589bSMark F. Adams     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
8110cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
8120cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
8132e68589bSMark F. Adams #endif
8142e68589bSMark F. Adams     /* create global vector of coorindates in 'coords' */
8152e68589bSMark F. Adams     if (npe > 1) {
816806fa848SBarry Smith       ierr = PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords);CHKERRQ(ierr);
817806fa848SBarry Smith     } else {
818c8b0795cSMark F. Adams       coords = (PetscReal*)pc_gamg->data;
819a2f3521dSMark F. Adams       data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols;
8202e68589bSMark F. Adams     }
8212e68589bSMark F. Adams     ierr = MatDestroy(&Gmat2);CHKERRQ(ierr);
8222e68589bSMark F. Adams 
8232e68589bSMark F. Adams     /* triangulate */
8242e68589bSMark F. Adams     if (dim == 2) {
825c8b0795cSMark F. Adams       PetscReal metric,tm;
8260cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
8270cbbd2e1SMark F. Adams       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
8282e68589bSMark F. Adams #endif
829806fa848SBarry Smith       ierr = triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric);CHKERRQ(ierr);
8300cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
8310cbbd2e1SMark F. Adams       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
8322e68589bSMark F. Adams #endif
8332e68589bSMark F. Adams       ierr = PetscFree(crsGID);CHKERRQ(ierr);
8342e68589bSMark F. Adams 
8352e68589bSMark F. Adams       /* clean up and create coordinates for coarse grid (output) */
8362e68589bSMark F. Adams       if (npe > 1) ierr = PetscFree(coords);CHKERRQ(ierr);
8372e68589bSMark F. Adams 
838c8b0795cSMark F. Adams       ierr = MPI_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm);CHKERRQ(ierr);
839c8b0795cSMark F. Adams       if (tm > 1.) { /* needs to be globalized - should not happen */
8402e68589bSMark F. Adams         if (verbose) {
84141b27cdeSMark F. Adams           PetscPrintf(wcomm,"[%d]%s failed metric for coarse grid %e\n",mype,__FUNCT__,tm);
8422e68589bSMark F. Adams         }
8432e68589bSMark F. Adams         ierr = MatDestroy(&Prol);CHKERRQ(ierr);
8442e68589bSMark F. Adams         Prol = PETSC_NULL;
845806fa848SBarry Smith       } else if (metric > .0) {
8462e68589bSMark F. Adams         if (verbose) {
84741b27cdeSMark F. Adams           PetscPrintf(wcomm,"[%d]%s worst metric for coarse grid = %e\n",mype,__FUNCT__,metric);
8482e68589bSMark F. Adams         }
8492e68589bSMark F. Adams       }
8502e68589bSMark F. Adams     } else {
851c666adbfSMark F. Adams       SETERRQ(wcomm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG");
8522e68589bSMark F. Adams     }
8532e68589bSMark F. Adams     { /* create next coords - output */
8542e68589bSMark F. Adams       PetscReal *crs_crds;
855806fa848SBarry Smith       ierr = PetscMalloc(dim*nLocalSelected*sizeof(PetscReal), &crs_crds);CHKERRQ(ierr);
8562e68589bSMark F. Adams       for (kk=0;kk<nLocalSelected;kk++) {/* grab local select nodes to promote - output */
857b43b03e9SMark F. Adams         PetscInt lid = clid_flid[kk];
858c8b0795cSMark F. Adams         for (jj=0;jj<dim;jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
8592e68589bSMark F. Adams       }
860c8b0795cSMark F. Adams 
861c8b0795cSMark F. Adams       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
862c8b0795cSMark F. Adams       pc_gamg->data = crs_crds; /* out */
863c8b0795cSMark F. Adams       pc_gamg->data_sz = dim*nLocalSelected;
8642e68589bSMark F. Adams     }
865a2f3521dSMark F. Adams     ierr = ISDestroy(&selected_2);CHKERRQ(ierr);
8662e68589bSMark F. Adams   }
867a2f3521dSMark F. Adams 
8682e68589bSMark F. Adams   *a_P_out = Prol;  /* out */
869b43b03e9SMark F. Adams   ierr = PetscFree(clid_flid);CHKERRQ(ierr);
8700cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG
8710cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
8720cbbd2e1SMark F. Adams #endif
873c8b0795cSMark F. Adams   PetscFunctionReturn(0);
874c8b0795cSMark F. Adams }
875c8b0795cSMark F. Adams 
876c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
877c8b0795cSMark F. Adams /*
878c8b0795cSMark F. Adams  PCCreateGAMG_GEO
879c8b0795cSMark F. Adams 
880c8b0795cSMark F. Adams   Input Parameter:
881c8b0795cSMark F. Adams    . pc -
882c8b0795cSMark F. Adams */
883c8b0795cSMark F. Adams #undef __FUNCT__
884c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_GEO"
885c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_GEO(PC pc)
886c8b0795cSMark F. Adams {
887c8b0795cSMark F. Adams   PetscErrorCode  ierr;
888c8b0795cSMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
889c8b0795cSMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
890c8b0795cSMark F. Adams 
891c8b0795cSMark F. Adams   PetscFunctionBegin;
892c8b0795cSMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GEO;
893c8b0795cSMark F. Adams   /* pc->ops->destroy        = PCDestroy_GEO; */
894c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
895c8b0795cSMark F. Adams 
896c8b0795cSMark F. Adams   /* set internal function pointers */
897c8b0795cSMark F. Adams   pc_gamg->graph = PCGAMGgraph_GEO;
898c8b0795cSMark F. Adams   pc_gamg->coarsen = PCGAMGcoarsen_GEO;
8990cbbd2e1SMark F. Adams   pc_gamg->prolongator = PCGAMGProlongator_GEO;
900c8b0795cSMark F. Adams   pc_gamg->optprol = 0;
901a2f3521dSMark F. Adams   pc_gamg->formkktprol = 0;
902c8b0795cSMark F. Adams 
903c8b0795cSMark F. Adams   pc_gamg->createdefaultdata = PCSetData_GEO;
904c8b0795cSMark F. Adams 
905c8b0795cSMark F. Adams   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
906c8b0795cSMark F. Adams                                             "PCSetCoordinates_C",
907c8b0795cSMark F. Adams                                             "PCSetCoordinates_GEO",
908c8b0795cSMark F. Adams                                             PCSetCoordinates_GEO);CHKERRQ(ierr);
9092e68589bSMark F. Adams 
9102e68589bSMark F. Adams   PetscFunctionReturn(0);
9112e68589bSMark F. Adams }
912