xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision bf4339c22951fe2ed4cd378ee6bed5e3925f33c6)
1 /*
2  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
3  */
4 
5 #include <../src/ksp/pc/impls/gamg/gamg.h>        /*I "petscpc.h" I*/
6 #include <petsc-private/kspimpl.h>
7 
8 #if defined(PETSC_HAVE_TRIANGLE)
9 #define REAL PetscReal
10 #include <triangle.h>
11 #endif
12 
13 #include <petscblaslapack.h>
14 
15 /* Private context for the GAMG preconditioner */
16 typedef struct {
17   PetscInt lid;            /* local vertex index */
18   PetscInt degree;         /* vertex degree */
19 } GAMGNode;
20 
21 int petsc_geo_mg_compare(const void *a, const void *b)
22 {
23   return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree);
24 }
25 
26 /* -------------------------------------------------------------------------- */
27 /*
28    PCSetCoordinates_GEO
29 
30    Input Parameter:
31    .  pc - the preconditioner context
32 */
33 #undef __FUNCT__
34 #define __FUNCT__ "PCSetCoordinates_GEO"
35 PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
36 {
37   PC_MG          *mg      = (PC_MG*)pc->data;
38   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
39   PetscErrorCode ierr;
40   PetscInt       arrsz,bs,my0,kk,ii,nloc,Iend;
41   Mat            Amat = pc->pmat;
42 
43   PetscFunctionBegin;
44   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
45   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
46 
47   ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr);
48   nloc = (Iend-my0)/bs;
49 
50   if (nloc!=a_nloc) SETERRQ2(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Stokes not supported nloc = %d %d.",a_nloc,nloc);
51   if ((Iend-my0)%bs!=0) SETERRQ1(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
52 
53   pc_gamg->data_cell_rows = 1;
54   if (coords==0 && nloc > 0) SETERRQ(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
55   pc_gamg->data_cell_cols = ndm; /* coordinates */
56 
57   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
58 
59   /* create data - syntactic sugar that should be refactored at some point */
60   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
61     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
62     ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr);
63   }
64   for (kk=0; kk<arrsz; kk++) pc_gamg->data[kk] = -999.;
65   pc_gamg->data[arrsz] = -99.;
66   /* copy data in - column oriented */
67   for (kk = 0; kk < nloc; kk++) {
68     for (ii = 0; ii < ndm; ii++) {
69       pc_gamg->data[ii*nloc + kk] =  coords[kk*ndm + ii];
70     }
71   }
72   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]);
73   pc_gamg->data_sz = arrsz;
74   PetscFunctionReturn(0);
75 }
76 
77 /* -------------------------------------------------------------------------- */
78 /*
79    PCSetData_GEO
80 
81   Input Parameter:
82    . pc -
83 */
84 #undef __FUNCT__
85 #define __FUNCT__ "PCSetData_GEO"
86 PetscErrorCode PCSetData_GEO(PC pc, Mat m)
87 {
88   PetscFunctionBegin;
89   SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"GEO MG needs coordinates");
90 }
91 
92 /* -------------------------------------------------------------------------- */
93 /*
94    PCSetFromOptions_GEO
95 
96   Input Parameter:
97    . pc -
98 */
99 #undef __FUNCT__
100 #define __FUNCT__ "PCSetFromOptions_GEO"
101 PetscErrorCode PCSetFromOptions_GEO(PetscOptions *PetscOptionsObject,PC pc)
102 {
103   PetscErrorCode ierr;
104 
105   PetscFunctionBegin;
106   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-GEO options");CHKERRQ(ierr);
107   {
108     /* -pc_gamg_sa_nsmooths */
109     /* pc_gamg_sa->smooths = 0; */
110     /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */
111     /*                        "smoothing steps for smoothed aggregation, usually 1 (0)", */
112     /*                        "PCGAMGSetNSmooths_AGG", */
113     /*                        pc_gamg_sa->smooths, */
114     /*                        &pc_gamg_sa->smooths, */
115     /*                        &flag);  */
116     /* CHKERRQ(ierr); */
117   }
118   ierr = PetscOptionsTail();CHKERRQ(ierr);
119   PetscFunctionReturn(0);
120 }
121 
122 /* -------------------------------------------------------------------------- */
123 /*
124  triangulateAndFormProl
125 
126    Input Parameter:
127    . selected_2 - list of selected local ID, includes selected ghosts
128    . data_stride -
129    . coords[2*data_stride] - column vector of local coordinates w/ ghosts
130    . nselected_1 - selected IDs that go with base (1) graph
131    . clid_lid_1[nselected_1] - lids of selected (c) nodes   ???????????
132    . agg_lists_1 - list of aggregates
133    . crsGID[selected.size()] - global index for prolongation operator
134    . bs - block size
135   Output Parameter:
136    . a_Prol - prolongation operator
137    . a_worst_best - measure of worst missed fine vertex, 0 is no misses
138 */
139 #undef __FUNCT__
140 #define __FUNCT__ "triangulateAndFormProl"
141 static PetscErrorCode triangulateAndFormProl(IS  selected_2, /* list of selected local ID, includes selected ghosts */
142                                              const PetscInt data_stride,
143                                              const PetscReal coords[], /* column vector of local coordinates w/ ghosts */
144                                              const PetscInt nselected_1, /* list of selected local ID, includes selected ghosts */
145                                              const PetscInt clid_lid_1[],
146                                              const PetscCoarsenData *agg_lists_1, /* selected_1 vertices of aggregate unselected vertices */
147                                              const PetscInt crsGID[],
148                                              const PetscInt bs,
149                                              Mat a_Prol, /* prolongation operator (output) */
150                                              PetscReal *a_worst_best) /* measure of worst missed fine vertex, 0 is no misses */
151 {
152 #if defined(PETSC_HAVE_TRIANGLE)
153   PetscErrorCode       ierr;
154   PetscInt             jj,tid,tt,idx,nselected_2;
155   struct triangulateio in,mid;
156   const PetscInt       *selected_idx_2;
157   PetscMPIInt          rank;
158   PetscInt             Istart,Iend,nFineLoc,myFine0;
159   int                  kk,nPlotPts,sid;
160   MPI_Comm             comm;
161   PetscReal            tm;
162 
163   PetscFunctionBegin;
164   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
165   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
166   ierr = ISGetSize(selected_2, &nselected_2);CHKERRQ(ierr);
167   if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */
168     *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */
169   } else *a_worst_best = 0.0;
170   ierr = MPI_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm);CHKERRQ(ierr);
171   if (tm > 0.0) {
172     *a_worst_best = 100.0;
173     PetscFunctionReturn(0);
174   }
175   ierr     = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
176   nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs;
177   nPlotPts = nFineLoc; /* locals */
178   /* traingle */
179   /* Define input points - in*/
180   in.numberofpoints          = nselected_2;
181   in.numberofpointattributes = 0;
182   /* get nselected points */
183   ierr = PetscMalloc1(2*(nselected_2), &in.pointlist);CHKERRQ(ierr);
184   ierr = ISGetIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
185 
186   for (kk=0,sid=0; kk<nselected_2; kk++,sid += 2) {
187     PetscInt lid = selected_idx_2[kk];
188     in.pointlist[sid]   = coords[lid];
189     in.pointlist[sid+1] = coords[data_stride + lid];
190     if (lid>=nFineLoc) nPlotPts++;
191   }
192   if (sid != 2*nselected_2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != 2*nselected_2 %D",sid,nselected_2);
193 
194   in.numberofsegments      = 0;
195   in.numberofedges         = 0;
196   in.numberofholes         = 0;
197   in.numberofregions       = 0;
198   in.trianglelist          = 0;
199   in.segmentmarkerlist     = 0;
200   in.pointattributelist    = 0;
201   in.pointmarkerlist       = 0;
202   in.triangleattributelist = 0;
203   in.trianglearealist      = 0;
204   in.segmentlist           = 0;
205   in.holelist              = 0;
206   in.regionlist            = 0;
207   in.edgelist              = 0;
208   in.edgemarkerlist        = 0;
209   in.normlist              = 0;
210 
211   /* triangulate */
212   mid.pointlist = 0;            /* Not needed if -N switch used. */
213   /* Not needed if -N switch used or number of point attributes is zero: */
214   mid.pointattributelist = 0;
215   mid.pointmarkerlist    = 0; /* Not needed if -N or -B switch used. */
216   mid.trianglelist       = 0;    /* Not needed if -E switch used. */
217   /* Not needed if -E switch used or number of triangle attributes is zero: */
218   mid.triangleattributelist = 0;
219   mid.neighborlist          = 0; /* Needed only if -n switch used. */
220   /* Needed only if segments are output (-p or -c) and -P not used: */
221   mid.segmentlist = 0;
222   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
223   mid.segmentmarkerlist = 0;
224   mid.edgelist          = 0;    /* Needed only if -e switch used. */
225   mid.edgemarkerlist    = 0; /* Needed if -e used and -B not used. */
226   mid.numberoftriangles = 0;
227 
228   /* Triangulate the points.  Switches are chosen to read and write a  */
229   /*   PSLG (p), preserve the convex hull (c), number everything from  */
230   /*   zero (z), assign a regional attribute to each element (A), and  */
231   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
232   /*   neighbor list (n).                                            */
233   if (nselected_2 != 0) { /* inactive processor */
234     char args[] = "npczQ"; /* c is needed ? */
235     triangulate(args, &in, &mid, (struct triangulateio*) NULL);
236     /* output .poly files for 'showme' */
237     if (!PETSC_TRUE) {
238       static int level = 1;
239       FILE       *file; char fname[32];
240 
241       sprintf(fname,"C%d_%d.poly",level,rank); file = fopen(fname, "w");
242       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
243       fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0);
244       /*Following lines: <vertex #> <x> <y> */
245       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid += 2) {
246         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
247       }
248       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
249       fprintf(file, "%d  %d\n",0,0);
250       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
251       /* One line: <# of holes> */
252       fprintf(file, "%d\n",0);
253       /* Following lines: <hole #> <x> <y> */
254       /* Optional line: <# of regional attributes and/or area constraints> */
255       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
256       fclose(file);
257 
258       /* elems */
259       sprintf(fname,"C%d_%d.ele",level,rank); file = fopen(fname, "w");
260       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
261       fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0);
262       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
263       for (kk=0,sid=0; kk<mid.numberoftriangles; kk++,sid += 3) {
264         fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]);
265       }
266       fclose(file);
267 
268       sprintf(fname,"C%d_%d.node",level,rank); file = fopen(fname, "w");
269       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
270       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
271       fprintf(file, "%d  %d  %d  %d\n",nPlotPts,2,0,0);
272       /*Following lines: <vertex #> <x> <y> */
273       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid+=2) {
274         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
275       }
276 
277       sid /= 2;
278       for (jj=0; jj<nFineLoc; jj++) {
279         PetscBool sel = PETSC_TRUE;
280         for (kk=0; kk<nselected_2 && sel; kk++) {
281           PetscInt lid = selected_idx_2[kk];
282           if (lid == jj) sel = PETSC_FALSE;
283         }
284         if (sel) fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]);
285       }
286       fclose(file);
287       if (sid != nPlotPts) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != nPlotPts %D",sid,nPlotPts);
288       level++;
289     }
290   }
291 #if defined PETSC_GAMG_USE_LOG
292   ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
293 #endif
294   { /* form P - setup some maps */
295     PetscInt clid,mm,*nTri,*node_tri;
296 
297     ierr = PetscMalloc1(nselected_2, &node_tri);CHKERRQ(ierr);
298     ierr = PetscMalloc1(nselected_2, &nTri);CHKERRQ(ierr);
299 
300     /* need list of triangles on node */
301     for (kk=0; kk<nselected_2; kk++) nTri[kk] = 0;
302     for (tid=0,kk=0; tid<mid.numberoftriangles; tid++) {
303       for (jj=0; jj<3; jj++) {
304         PetscInt cid = mid.trianglelist[kk++];
305         if (nTri[cid] == 0) node_tri[cid] = tid;
306         nTri[cid]++;
307       }
308     }
309 #define EPS 1.e-12
310     /* find points and set prolongation */
311     for (mm = clid = 0; mm < nFineLoc; mm++) {
312       PetscBool ise;
313       ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise);CHKERRQ(ierr);
314       if (!ise) {
315         const PetscInt lid = mm;
316         /* for (clid_iterator=0;clid_iterator<nselected_1;clid_iterator++) { */
317         PetscScalar  AA[3][3];
318         PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO;
319         PetscCDPos   pos;
320         ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
321         while (pos) {
322           PetscInt flid;
323           ierr = PetscLLNGetID(pos, &flid);CHKERRQ(ierr);
324           ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
325 
326           if (flid < nFineLoc) {  /* could be a ghost */
327             PetscInt       bestTID = -1; PetscReal best_alpha = 1.e10;
328             const PetscInt fgid    = flid + myFine0;
329             /* compute shape function for gid */
330             const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0};
331             PetscBool       haveit    =PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3];
332 
333             /* look for it */
334             for (tid = node_tri[clid], jj=0;
335                  jj < 5 && !haveit && tid != -1;
336                  jj++) {
337               for (tt=0; tt<3; tt++) {
338                 PetscInt cid2 = mid.trianglelist[3*tid + tt];
339                 PetscInt lid2 = selected_idx_2[cid2];
340                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
341                 clids[tt] = cid2; /* store for interp */
342               }
343 
344               for (tt=0; tt<3; tt++) alpha[tt] = (PetscScalar)fcoord[tt];
345 
346               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
347               PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
348               {
349                 PetscBool have=PETSC_TRUE;  PetscReal lowest=1.e10;
350                 for (tt = 0, idx = 0; tt < 3; tt++) {
351                   if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
352                   if (PetscRealPart(alpha[tt]) < lowest) {
353                     lowest = PetscRealPart(alpha[tt]);
354                     idx    = tt;
355                   }
356                 }
357                 haveit = have;
358               }
359               tid = mid.neighborlist[3*tid + idx];
360             }
361 
362             if (!haveit) {
363               /* brute force */
364               for (tid=0; tid<mid.numberoftriangles && !haveit; tid++) {
365                 for (tt=0; tt<3; tt++) {
366                   PetscInt cid2 = mid.trianglelist[3*tid + tt];
367                   PetscInt lid2 = selected_idx_2[cid2];
368                   AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
369                   clids[tt] = cid2; /* store for interp */
370                 }
371                 for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
372                 /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
373                 PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
374                 {
375                   PetscBool have=PETSC_TRUE;  PetscReal worst=0.0, v;
376                   for (tt=0; tt<3 && have; tt++) {
377                     if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE;
378                     if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v;
379                   }
380                   if (worst < best_alpha) {
381                     best_alpha = worst; bestTID = tid;
382                   }
383                   haveit = have;
384                 }
385               }
386             }
387             if (!haveit) {
388               if (best_alpha > *a_worst_best) *a_worst_best = best_alpha;
389               /* use best one */
390               for (tt=0; tt<3; tt++) {
391                 PetscInt cid2 = mid.trianglelist[3*bestTID + tt];
392                 PetscInt lid2 = selected_idx_2[cid2];
393                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
394                 clids[tt] = cid2; /* store for interp */
395               }
396               for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
397               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
398               PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
399             }
400 
401             /* put in row of P */
402             for (idx=0; idx<3; idx++) {
403               PetscScalar shp = alpha[idx];
404               if (PetscAbs(PetscRealPart(shp)) > 1.e-6) {
405                 PetscInt cgid = crsGID[clids[idx]];
406                 PetscInt jj   = cgid*bs, ii = fgid*bs; /* need to gloalize */
407                 for (tt=0; tt < bs; tt++, ii++, jj++) {
408                   ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES);CHKERRQ(ierr);
409                 }
410               }
411             }
412           }
413         } /* aggregates iterations */
414         clid++;
415       } /* a coarse agg */
416     } /* for all fine nodes */
417 
418     ierr = ISRestoreIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
419     ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
420     ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
421 
422     ierr = PetscFree(node_tri);CHKERRQ(ierr);
423     ierr = PetscFree(nTri);CHKERRQ(ierr);
424   }
425 #if defined PETSC_GAMG_USE_LOG
426   ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
427 #endif
428   free(mid.trianglelist);
429   free(mid.neighborlist);
430   ierr = PetscFree(in.pointlist);CHKERRQ(ierr);
431   PetscFunctionReturn(0);
432 #else
433   SETERRQ(PetscObjectComm((PetscObject)a_Prol),PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG");
434 #endif
435 }
436 /* -------------------------------------------------------------------------- */
437 /*
438    getGIDsOnSquareGraph - square graph, get
439 
440    Input Parameter:
441    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
442    . clid_lid_1 - [nselected_1] lids of selected nodes
443    . Gmat1 - graph that goes with 'selected_1'
444    Output Parameter:
445    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
446    . a_Gmat_2 - graph that is squared of 'Gmat_1'
447    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
448 */
449 #undef __FUNCT__
450 #define __FUNCT__ "getGIDsOnSquareGraph"
451 static PetscErrorCode getGIDsOnSquareGraph(const PetscInt nselected_1,const PetscInt clid_lid_1[],const Mat Gmat1,IS *a_selected_2,Mat *a_Gmat_2,PetscInt **a_crsGID)
452 {
453   PetscErrorCode ierr;
454   PetscMPIInt    size;
455   PetscInt       *crsGID, kk,my0,Iend,nloc;
456   MPI_Comm       comm;
457 
458   PetscFunctionBegin;
459   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
460   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
461   ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend);CHKERRQ(ierr); /* AIJ */
462   nloc = Iend - my0; /* this does not change */
463 
464   if (size == 1) { /* not much to do in serial */
465     ierr = PetscMalloc1(nselected_1, &crsGID);CHKERRQ(ierr);
466     for (kk=0; kk<nselected_1; kk++) crsGID[kk] = kk;
467     *a_Gmat_2 = 0;
468     ierr      = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);CHKERRQ(ierr);
469   } else {
470     PetscInt    idx,num_fine_ghosts,num_crs_ghost,myCrs0;
471     Mat_MPIAIJ  *mpimat2;
472     Mat         Gmat2;
473     Vec         locState;
474     PetscScalar *cpcol_state;
475 
476     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
477     kk = nselected_1;
478     MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPIU_SUM, comm);
479     myCrs0 -= nselected_1;
480 
481     if (a_Gmat_2) { /* output */
482       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
483       ierr      = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
484       *a_Gmat_2 = Gmat2; /* output */
485     } else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
486     /* get coarse grid GIDS for selected (locals and ghosts) */
487     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
488     ierr    = MatCreateVecs(Gmat2, &locState, 0);CHKERRQ(ierr);
489     ierr    = VecSet(locState, (PetscScalar)(PetscReal)(-1));CHKERRQ(ierr); /* set with UNKNOWN state */
490     for (kk=0; kk<nselected_1; kk++) {
491       PetscInt    fgid = clid_lid_1[kk] + my0;
492       PetscScalar v    = (PetscScalar)(kk+myCrs0);
493       ierr = VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES);CHKERRQ(ierr); /* set with PID */
494     }
495     ierr = VecAssemblyBegin(locState);CHKERRQ(ierr);
496     ierr = VecAssemblyEnd(locState);CHKERRQ(ierr);
497     ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
498     ierr =   VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
499     ierr = VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts);CHKERRQ(ierr);
500     ierr = VecGetArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
501     for (kk=0,num_crs_ghost=0; kk<num_fine_ghosts; kk++) {
502       if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++;
503     }
504     ierr = PetscMalloc1(nselected_1+num_crs_ghost, &crsGID);CHKERRQ(ierr); /* output */
505     {
506       PetscInt *selected_set;
507       ierr = PetscMalloc1(nselected_1+num_crs_ghost, &selected_set);CHKERRQ(ierr);
508       /* do ghost of 'crsGID' */
509       for (kk=0,idx=nselected_1; kk<num_fine_ghosts; kk++) {
510         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
511           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
512           selected_set[idx] = nloc + kk;
513           crsGID[idx++]     = cgid;
514         }
515       }
516       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);
517       ierr = VecRestoreArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
518       /* do locals in 'crsGID' */
519       ierr = VecGetArray(locState, &cpcol_state);CHKERRQ(ierr);
520       for (kk=0,idx=0; kk<nloc; kk++) {
521         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
522           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
523           selected_set[idx] = kk;
524           crsGID[idx++]     = cgid;
525         }
526       }
527       if (idx != nselected_1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != nselected_1 %D",idx,nselected_1);
528       ierr = VecRestoreArray(locState, &cpcol_state);CHKERRQ(ierr);
529 
530       if (a_selected_2 != 0) { /* output */
531         ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);CHKERRQ(ierr);
532       } else {
533         ierr = PetscFree(selected_set);CHKERRQ(ierr);
534       }
535     }
536     ierr = VecDestroy(&locState);CHKERRQ(ierr);
537   }
538   *a_crsGID = crsGID; /* output */
539   PetscFunctionReturn(0);
540 }
541 
542 /* -------------------------------------------------------------------------- */
543 /*
544    PCGAMGGraph_GEO
545 
546   Input Parameter:
547    . pc - this
548    . Amat - matrix on this fine level
549   Output Parameter:
550    . a_Gmat
551 */
552 #undef __FUNCT__
553 #define __FUNCT__ "PCGAMGGraph_GEO"
554 PetscErrorCode PCGAMGGraph_GEO(PC pc,const Mat Amat,Mat *a_Gmat)
555 {
556   PetscErrorCode  ierr;
557   PC_MG           *mg      = (PC_MG*)pc->data;
558   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
559   const PetscReal vfilter  = pc_gamg->threshold;
560   MPI_Comm        comm;
561   Mat             Gmat;
562   PetscBool       set,flg,symm;
563 
564   PetscFunctionBegin;
565   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
566   ierr = PetscLogEventBegin(PC_GAMGGraph_GEO,0,0,0,0);CHKERRQ(ierr);
567 
568   ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr);
569   symm = (PetscBool)!(set && flg);
570 
571   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
572   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr);
573 
574   *a_Gmat = Gmat;
575   ierr = PetscLogEventEnd(PC_GAMGGraph_GEO,0,0,0,0);CHKERRQ(ierr);
576   PetscFunctionReturn(0);
577 }
578 
579 /* -------------------------------------------------------------------------- */
580 /*
581    PCGAMGCoarsen_GEO
582 
583   Input Parameter:
584    . a_pc - this
585    . a_Gmat - graph
586   Output Parameter:
587    . a_llist_parent - linked list from selected indices for data locality only
588 */
589 #undef __FUNCT__
590 #define __FUNCT__ "PCGAMGCoarsen_GEO"
591 PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc,Mat *a_Gmat,PetscCoarsenData **a_llist_parent)
592 {
593   PetscErrorCode ierr;
594   PC_MG          *mg      = (PC_MG*)a_pc->data;
595   PetscInt       Istart,Iend,nloc,kk,Ii,ncols;
596   IS             perm;
597   GAMGNode       *gnodes;
598   PetscInt       *permute;
599   Mat            Gmat  = *a_Gmat;
600   MPI_Comm       comm;
601   MatCoarsen     crs;
602 
603   PetscFunctionBegin;
604   ierr = PetscObjectGetComm((PetscObject)a_pc,&comm);CHKERRQ(ierr);
605   ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
606   ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr);
607   nloc = (Iend-Istart);
608 
609   /* create random permutation with sort for geo-mg */
610   ierr = PetscMalloc1(nloc, &gnodes);CHKERRQ(ierr);
611   ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr);
612 
613   for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */
614     ierr = MatGetRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr);
615     {
616       PetscInt lid = Ii - Istart;
617       gnodes[lid].lid    = lid;
618       gnodes[lid].degree = ncols;
619     }
620     ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr);
621   }
622   /* randomize */
623   srand(1); /* make deterministic */
624   if (PETSC_TRUE) {
625     PetscBool *bIndexSet;
626     ierr = PetscMalloc1(nloc, &bIndexSet);CHKERRQ(ierr);
627     for (Ii = 0; Ii < nloc; Ii++) bIndexSet[Ii] = PETSC_FALSE;
628     for (Ii = 0; Ii < nloc; Ii++) {
629       PetscInt iSwapIndex = rand()%nloc;
630       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
631         GAMGNode iTemp = gnodes[iSwapIndex];
632         gnodes[iSwapIndex]    = gnodes[Ii];
633         gnodes[Ii]            = iTemp;
634         bIndexSet[Ii]         = PETSC_TRUE;
635         bIndexSet[iSwapIndex] = PETSC_TRUE;
636       }
637     }
638     ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
639   }
640   /* only sort locals */
641   qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare);
642   /* create IS of permutation */
643   for (kk=0; kk<nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */
644   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr);
645 
646   ierr = PetscFree(gnodes);CHKERRQ(ierr);
647 
648   /* get MIS aggs */
649 
650   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
651   ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr);
652   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
653   ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr);
654   ierr = MatCoarsenSetStrictAggs(crs, PETSC_FALSE);CHKERRQ(ierr);
655   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
656   ierr = MatCoarsenGetData(crs, a_llist_parent);CHKERRQ(ierr);
657   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
658 
659   ierr = ISDestroy(&perm);CHKERRQ(ierr);
660   ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
661   PetscFunctionReturn(0);
662 }
663 
664 /* -------------------------------------------------------------------------- */
665 /*
666  PCGAMGProlongator_GEO
667 
668  Input Parameter:
669  . pc - this
670  . Amat - matrix on this fine level
671  . Graph - used to get ghost data for nodes in
672  . selected_1 - [nselected]
673  . agg_lists - [nselected]
674  Output Parameter:
675  . a_P_out - prolongation operator to the next level
676  */
677 #undef __FUNCT__
678 #define __FUNCT__ "PCGAMGProlongator_GEO"
679 PetscErrorCode PCGAMGProlongator_GEO(PC pc,const Mat Amat,const Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
680 {
681   PC_MG          *mg      = (PC_MG*)pc->data;
682   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
683   const PetscInt dim      = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
684   PetscErrorCode ierr;
685   PetscInt       Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid;
686   Mat            Prol;
687   PetscMPIInt    rank, size;
688   MPI_Comm       comm;
689   IS             selected_2,selected_1;
690   const PetscInt *selected_idx;
691 
692   PetscFunctionBegin;
693   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
694   ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
695   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
696   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
697   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
698   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
699   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
700   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) % bs %D",Iend,Istart,bs);
701 
702   /* get 'nLocalSelected' */
703   ierr = PetscCDGetMIS(agg_lists, &selected_1);CHKERRQ(ierr);
704   ierr = ISGetSize(selected_1, &jj);CHKERRQ(ierr);
705   ierr = PetscMalloc1(jj, &clid_flid);CHKERRQ(ierr);
706   ierr = ISGetIndices(selected_1, &selected_idx);CHKERRQ(ierr);
707   for (kk=0,nLocalSelected=0; kk<jj; kk++) {
708     PetscInt lid = selected_idx[kk];
709     if (lid<nloc) {
710       ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr);
711       if (ncols>1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */
712       ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr);
713     }
714   }
715   ierr = ISRestoreIndices(selected_1, &selected_idx);CHKERRQ(ierr);
716   ierr = ISDestroy(&selected_1);CHKERRQ(ierr); /* this is selected_1 in serial */
717 
718   /* create prolongator, create P matrix */
719   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
720   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
721   ierr = MatSetBlockSizes(Prol, bs, bs);CHKERRQ(ierr);
722   ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr);
723   ierr = MatSeqAIJSetPreallocation(Prol,3*data_cols,NULL);CHKERRQ(ierr);
724   ierr = MatMPIAIJSetPreallocation(Prol,3*data_cols,NULL,3*data_cols,NULL);CHKERRQ(ierr);
725   /* ierr = MatCreateAIJ(comm,  */
726   /*                      nloc*bs, nLocalSelected*bs, */
727   /*                      PETSC_DETERMINE, PETSC_DETERMINE, */
728   /*                      3*data_cols, NULL,  */
729   /*                      3*data_cols, NULL, */
730   /*                      &Prol); */
731   /* CHKERRQ(ierr); */
732 
733   /* can get all points "removed" - but not on geomg */
734   ierr =  MatGetSize(Prol, &kk, &jj);CHKERRQ(ierr);
735   if (jj==0) {
736     ierr = PetscInfo(pc,"ERROE: no selected points on coarse grid\n");CHKERRQ(ierr);
737     ierr = PetscFree(clid_flid);CHKERRQ(ierr);
738     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
739     *a_P_out = NULL;  /* out */
740     PetscFunctionReturn(0);
741   }
742 
743   {
744     PetscReal *coords;
745     PetscInt  data_stride;
746     PetscInt  *crsGID = NULL;
747     Mat       Gmat2;
748 
749     if (dim != data_cols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim %D != data_cols %D",dim,data_cols);
750     /* grow ghost data for better coarse grid cover of fine grid */
751 #if defined PETSC_GAMG_USE_LOG
752     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
753 #endif
754     /* messy method, squares graph and gets some data */
755     ierr = getGIDsOnSquareGraph(nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID);CHKERRQ(ierr);
756     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
757 #if defined PETSC_GAMG_USE_LOG
758     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
759 #endif
760     /* create global vector of coorindates in 'coords' */
761     if (size > 1) {
762       ierr = PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords);CHKERRQ(ierr);
763     } else {
764       coords      = (PetscReal*)pc_gamg->data;
765       data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols;
766     }
767     ierr = MatDestroy(&Gmat2);CHKERRQ(ierr);
768 
769     /* triangulate */
770     if (dim == 2) {
771       PetscReal metric,tm;
772 #if defined PETSC_GAMG_USE_LOG
773       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
774 #endif
775       ierr = triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric);CHKERRQ(ierr);
776 #if defined PETSC_GAMG_USE_LOG
777       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
778 #endif
779       ierr = PetscFree(crsGID);CHKERRQ(ierr);
780 
781       /* clean up and create coordinates for coarse grid (output) */
782       if (size > 1) ierr = PetscFree(coords);CHKERRQ(ierr);
783 
784       ierr = MPI_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm);CHKERRQ(ierr);
785       if (tm > 1.) { /* needs to be globalized - should not happen */
786         ierr = PetscInfo1(pc," failed metric for coarse grid %e\n",(double)tm);CHKERRQ(ierr);
787         ierr = MatDestroy(&Prol);CHKERRQ(ierr);
788       } else if (metric > .0) {
789         ierr = PetscInfo1(pc,"worst metric for coarse grid = %e\n",(double)metric);CHKERRQ(ierr);
790       }
791     } else SETERRQ(comm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG");
792     { /* create next coords - output */
793       PetscReal *crs_crds;
794       ierr = PetscMalloc1(dim*nLocalSelected, &crs_crds);CHKERRQ(ierr);
795       for (kk=0; kk<nLocalSelected; kk++) { /* grab local select nodes to promote - output */
796         PetscInt lid = clid_flid[kk];
797         for (jj=0; jj<dim; jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
798       }
799 
800       ierr             = PetscFree(pc_gamg->data);CHKERRQ(ierr);
801       pc_gamg->data    = crs_crds; /* out */
802       pc_gamg->data_sz = dim*nLocalSelected;
803     }
804     ierr = ISDestroy(&selected_2);CHKERRQ(ierr);
805   }
806 
807   *a_P_out = Prol;  /* out */
808   ierr     = PetscFree(clid_flid);CHKERRQ(ierr);
809   ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
810   PetscFunctionReturn(0);
811 }
812 
813 #undef __FUNCT__
814 #define __FUNCT__ "PCDestroy_GAMG_GEO"
815 static PetscErrorCode PCDestroy_GAMG_GEO(PC pc)
816 {
817   PetscErrorCode ierr;
818 
819   PetscFunctionBegin;
820   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
821   PetscFunctionReturn(0);
822 }
823 
824 /* -------------------------------------------------------------------------- */
825 /*
826  PCCreateGAMG_GEO
827 
828   Input Parameter:
829    . pc -
830 */
831 #undef __FUNCT__
832 #define __FUNCT__ "PCCreateGAMG_GEO"
833 PetscErrorCode  PCCreateGAMG_GEO(PC pc)
834 {
835   PetscErrorCode ierr;
836   PC_MG          *mg      = (PC_MG*)pc->data;
837   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
838 
839   PetscFunctionBegin;
840   pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO;
841   pc_gamg->ops->destroy        = PCDestroy_GAMG_GEO;
842   /* reset does not do anything; setup not virtual */
843 
844   /* set internal function pointers */
845   pc_gamg->ops->graph             = PCGAMGGraph_GEO;
846   pc_gamg->ops->coarsen           = PCGAMGCoarsen_GEO;
847   pc_gamg->ops->prolongator       = PCGAMGProlongator_GEO;
848   pc_gamg->ops->optprolongator    = NULL;
849   pc_gamg->ops->createdefaultdata = PCSetData_GEO;
850 
851   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_GEO);CHKERRQ(ierr);
852   PetscFunctionReturn(0);
853 }
854