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