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