xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision 5c8f6a953e7ed1c81f507d64aebddb11080b60e9)
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   PetscFunctionReturn(0);
134 }
135 
136 /* -------------------------------------------------------------------------- */
137 /*
138  triangulateAndFormProl
139 
140    Input Parameter:
141    . selected_2 - list of selected local ID, includes selected ghosts
142    . data_stride -
143    . coords[2*data_stride] - column vector of local coordinates w/ ghosts
144    . nselected_1 - selected IDs that go with base (1) graph
145    . clid_lid_1[nselected_1] - lids of selected (c) nodes   ???????????
146    . agg_lists_1 - list of aggregates
147    . crsGID[selected.size()] - global index for prolongation operator
148    . bs - block size
149   Output Parameter:
150    . a_Prol - prolongation operator
151    . a_worst_best - measure of worst missed fine vertex, 0 is no misses
152 */
153 #undef __FUNCT__
154 #define __FUNCT__ "triangulateAndFormProl"
155 static PetscErrorCode triangulateAndFormProl(IS  selected_2, /* list of selected local ID, includes selected ghosts */
156                                               const PetscInt data_stride,
157                                               const PetscReal coords[], /* column vector of local coordinates w/ ghosts */
158                                               const PetscInt nselected_1, /* list of selected local ID, includes selected ghosts */
159                                               const PetscInt clid_lid_1[],
160                                               const PetscCoarsenData *agg_lists_1, /* selected_1 vertices of aggregate unselected vertices */
161                                               const PetscInt crsGID[],
162                                               const PetscInt bs,
163                                               Mat a_Prol, /* prolongation operator (output) */
164                                               PetscReal *a_worst_best) /* measure of worst missed fine vertex, 0 is no misses */
165 {
166 #if defined(PETSC_HAVE_TRIANGLE)
167   PetscErrorCode       ierr;
168   PetscInt             jj,tid,tt,idx,nselected_2;
169   struct triangulateio in,mid;
170   const PetscInt      *selected_idx_2;
171   PetscMPIInt          mype,npe;
172   PetscInt             Istart,Iend,nFineLoc,myFine0;
173   int                  kk,nPlotPts,sid;
174   MPI_Comm             wcomm = ((PetscObject)a_Prol)->comm;
175   PetscReal            tm;
176 
177   PetscFunctionBegin;
178   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
179   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
180   ierr = ISGetSize(selected_2, &nselected_2);CHKERRQ(ierr);
181   if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */
182     *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */
183   } else *a_worst_best = 0.0;
184   ierr = MPI_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm);CHKERRQ(ierr);
185   if (tm > 0.0) {
186     *a_worst_best = 100.0;
187     PetscFunctionReturn(0);
188   }
189   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
190   nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs;
191   nPlotPts = nFineLoc; /* locals */
192   /* traingle */
193   /* Define input points - in*/
194   in.numberofpoints = nselected_2;
195   in.numberofpointattributes = 0;
196   /* get nselected points */
197   ierr = PetscMalloc(2*(nselected_2)*sizeof(REAL), &in.pointlist);CHKERRQ(ierr);
198   ierr = ISGetIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
199 
200   for (kk=0,sid=0;kk<nselected_2;kk++,sid += 2) {
201     PetscInt lid = selected_idx_2[kk];
202     in.pointlist[sid] = coords[lid];
203     in.pointlist[sid+1] = coords[data_stride + lid];
204     if (lid>=nFineLoc) nPlotPts++;
205   }
206   assert(sid==2*nselected_2);
207 
208   in.numberofsegments = 0;
209   in.numberofedges = 0;
210   in.numberofholes = 0;
211   in.numberofregions = 0;
212   in.trianglelist = 0;
213   in.segmentmarkerlist = 0;
214   in.pointattributelist = 0;
215   in.pointmarkerlist = 0;
216   in.triangleattributelist = 0;
217   in.trianglearealist = 0;
218   in.segmentlist = 0;
219   in.holelist = 0;
220   in.regionlist = 0;
221   in.edgelist = 0;
222   in.edgemarkerlist = 0;
223   in.normlist = 0;
224   /* triangulate */
225   mid.pointlist = 0;            /* Not needed if -N switch used. */
226   /* Not needed if -N switch used or number of point attributes is zero: */
227   mid.pointattributelist = 0;
228   mid.pointmarkerlist = 0; /* Not needed if -N or -B switch used. */
229   mid.trianglelist = 0;          /* Not needed if -E switch used. */
230   /* Not needed if -E switch used or number of triangle attributes is zero: */
231   mid.triangleattributelist = 0;
232   mid.neighborlist = 0;         /* Needed only if -n switch used. */
233   /* Needed only if segments are output (-p or -c) and -P not used: */
234   mid.segmentlist = 0;
235   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
236   mid.segmentmarkerlist = 0;
237   mid.edgelist = 0;             /* Needed only if -e switch used. */
238   mid.edgemarkerlist = 0;   /* Needed if -e used and -B not used. */
239   mid.numberoftriangles = 0;
240 
241   /* Triangulate the points.  Switches are chosen to read and write a  */
242   /*   PSLG (p), preserve the convex hull (c), number everything from  */
243   /*   zero (z), assign a regional attribute to each element (A), and  */
244   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
245   /*   neighbor list (n).                                            */
246   if (nselected_2 != 0) { /* inactive processor */
247     char args[] = "npczQ"; /* c is needed ? */
248     triangulate(args, &in, &mid, (struct triangulateio *) NULL);
249     /* output .poly files for 'showme' */
250     if (!PETSC_TRUE) {
251       static int level = 1;
252       FILE *file; char fname[32];
253 
254       sprintf(fname,"C%d_%d.poly",level,mype); file = fopen(fname, "w");
255       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
256       fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0);
257       /*Following lines: <vertex #> <x> <y> */
258       for (kk=0,sid=0;kk<in.numberofpoints;kk++,sid += 2) {
259         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
260       }
261       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
262       fprintf(file, "%d  %d\n",0,0);
263       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
264       /* One line: <# of holes> */
265       fprintf(file, "%d\n",0);
266       /* Following lines: <hole #> <x> <y> */
267       /* Optional line: <# of regional attributes and/or area constraints> */
268       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
269       fclose(file);
270 
271       /* elems */
272       sprintf(fname,"C%d_%d.ele",level,mype); file = fopen(fname, "w");
273       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
274       fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0);
275       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
276       for (kk=0,sid=0;kk<mid.numberoftriangles;kk++,sid += 3) {
277         fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]);
278       }
279       fclose(file);
280 
281       sprintf(fname,"C%d_%d.node",level,mype); file = fopen(fname, "w");
282       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
283       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
284       fprintf(file, "%d  %d  %d  %d\n",nPlotPts,2,0,0);
285       /*Following lines: <vertex #> <x> <y> */
286       for (kk=0,sid=0;kk<in.numberofpoints;kk++,sid+=2) {
287         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
288       }
289 
290       sid /= 2;
291       for (jj=0;jj<nFineLoc;jj++) {
292         PetscBool sel = PETSC_TRUE;
293         for (kk=0 ; kk<nselected_2 && sel ; kk++) {
294           PetscInt lid = selected_idx_2[kk];
295           if (lid == jj) sel = PETSC_FALSE;
296         }
297         if (sel) {
298           fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]);
299         }
300       }
301       fclose(file);
302       assert(sid==nPlotPts);
303       level++;
304     }
305   }
306 #if defined PETSC_GAMG_USE_LOG
307   ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
308 #endif
309   { /* form P - setup some maps */
310     PetscInt clid,mm,*nTri,*node_tri;
311 
312     ierr = PetscMalloc(nselected_2*sizeof(PetscInt), &node_tri);CHKERRQ(ierr);
313     ierr = PetscMalloc(nselected_2*sizeof(PetscInt), &nTri);CHKERRQ(ierr);
314 
315     /* need list of triangles on node */
316     for (kk=0;kk<nselected_2;kk++) nTri[kk] = 0;
317     for (tid=0,kk=0;tid<mid.numberoftriangles;tid++) {
318       for (jj=0;jj<3;jj++) {
319         PetscInt cid = mid.trianglelist[kk++];
320         if (nTri[cid] == 0) node_tri[cid] = tid;
321         nTri[cid]++;
322       }
323     }
324 #define EPS 1.e-12
325     /* find points and set prolongation */
326     for (mm = clid = 0 ; mm < nFineLoc ; mm++) {
327       PetscBool ise;
328       ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise);CHKERRQ(ierr);
329       if (!ise) {
330         const PetscInt lid = mm;
331         /* for (clid_iterator=0;clid_iterator<nselected_1;clid_iterator++) { */
332         /* PetscInt flid = clid_lid_1[clid_iterator]; assert(flid != -1); */
333         PetscScalar AA[3][3];
334         PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO;
335         PetscCDPos         pos;
336         ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
337         while (pos) {
338           PetscInt flid;
339           ierr = PetscLLNGetID(pos, &flid);CHKERRQ(ierr);
340           ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
341 
342           if (flid < nFineLoc) {  /* could be a ghost */
343             PetscInt bestTID = -1; PetscReal best_alpha = 1.e10;
344             const PetscInt fgid = flid + myFine0;
345             /* compute shape function for gid */
346             const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0};
347             PetscBool haveit=PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3];
348             /* look for it */
349             for (tid = node_tri[clid], jj=0;
350                  jj < 5 && !haveit && tid != -1;
351                  jj++) {
352               for (tt=0;tt<3;tt++) {
353                 PetscInt cid2 = mid.trianglelist[3*tid + tt];
354                 PetscInt lid2 = selected_idx_2[cid2];
355                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
356                 clids[tt] = cid2; /* store for interp */
357               }
358 
359               for (tt=0;tt<3;tt++) alpha[tt] = (PetscScalar)fcoord[tt];
360 
361               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
362               LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
363               {
364                 PetscBool have=PETSC_TRUE;  PetscReal lowest=1.e10;
365                 for (tt = 0, idx = 0 ; tt < 3 ; tt++) {
366                   if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
367                   if (PetscRealPart(alpha[tt]) < lowest) {
368                     lowest = PetscRealPart(alpha[tt]);
369                     idx = tt;
370                   }
371                 }
372                 haveit = have;
373               }
374               tid = mid.neighborlist[3*tid + idx];
375             }
376 
377             if (!haveit) {
378               /* brute force */
379               for (tid=0 ; tid<mid.numberoftriangles && !haveit ; tid++) {
380                 for (tt=0;tt<3;tt++) {
381                   PetscInt cid2 = mid.trianglelist[3*tid + 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                 LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
389                 {
390                   PetscBool have=PETSC_TRUE;  PetscReal worst=0.0, v;
391                   for (tt=0; tt<3 && have ;tt++) {
392                     if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE;
393                     if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v;
394                   }
395                   if (worst < best_alpha) {
396                     best_alpha = worst; bestTID = tid;
397                   }
398                   haveit = have;
399                 }
400               }
401             }
402             if (!haveit) {
403               if (best_alpha > *a_worst_best) *a_worst_best = best_alpha;
404               /* use best one */
405               for (tt=0;tt<3;tt++) {
406                 PetscInt cid2 = mid.trianglelist[3*bestTID + tt];
407                 PetscInt lid2 = selected_idx_2[cid2];
408                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
409                 clids[tt] = cid2; /* store for interp */
410               }
411               for (tt=0;tt<3;tt++) alpha[tt] = fcoord[tt];
412               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
413               LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
414             }
415 
416             /* put in row of P */
417             for (idx=0;idx<3;idx++) {
418               PetscScalar shp = alpha[idx];
419               if (PetscAbs(PetscRealPart(shp)) > 1.e-6) {
420                 PetscInt cgid = crsGID[clids[idx]];
421                 PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */
422                 for (tt=0 ; tt < bs ; tt++, ii++, jj++) {
423                   ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES);CHKERRQ(ierr);
424                 }
425               }
426             }
427           }
428         } /* aggregates iterations */
429         clid++;
430       } /* a coarse agg */
431     } /* for all fine nodes */
432 
433     ierr = ISRestoreIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
434     ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
435     ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
436 
437     ierr = PetscFree(node_tri);CHKERRQ(ierr);
438     ierr = PetscFree(nTri);CHKERRQ(ierr);
439   }
440 #if defined PETSC_GAMG_USE_LOG
441   ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
442 #endif
443   free(mid.trianglelist);
444   free(mid.neighborlist);
445   ierr = PetscFree(in.pointlist);CHKERRQ(ierr);
446 
447   PetscFunctionReturn(0);
448 #else
449   SETERRQ(((PetscObject)a_Prol)->comm,PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG");
450 #endif
451 }
452 /* -------------------------------------------------------------------------- */
453 /*
454    getGIDsOnSquareGraph - square graph, get
455 
456    Input Parameter:
457    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
458    . clid_lid_1 - [nselected_1] lids of selected nodes
459    . Gmat1 - graph that goes with 'selected_1'
460    Output Parameter:
461    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
462    . a_Gmat_2 - graph that is squared of 'Gmat_1'
463    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
464 */
465 #undef __FUNCT__
466 #define __FUNCT__ "getGIDsOnSquareGraph"
467 static PetscErrorCode getGIDsOnSquareGraph(const PetscInt nselected_1,
468                                             const PetscInt clid_lid_1[],
469                                             const Mat Gmat1,
470                                             IS *a_selected_2,
471                                             Mat *a_Gmat_2,
472                                             PetscInt **a_crsGID)
473 {
474   PetscErrorCode ierr;
475   PetscMPIInt    mype,npe;
476   PetscInt       *crsGID, kk,my0,Iend,nloc;
477   MPI_Comm       wcomm = ((PetscObject)Gmat1)->comm;
478 
479   PetscFunctionBegin;
480   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
481   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
482   ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend);CHKERRQ(ierr); /* AIJ */
483   nloc = Iend - my0; /* this does not change */
484 
485   if (npe == 1) { /* not much to do in serial */
486     ierr = PetscMalloc(nselected_1*sizeof(PetscInt), &crsGID);CHKERRQ(ierr);
487     for (kk=0;kk<nselected_1;kk++) crsGID[kk] = kk;
488     *a_Gmat_2 = 0;
489     ierr = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);CHKERRQ(ierr);
490   } else {
491     PetscInt      idx,num_fine_ghosts,num_crs_ghost,myCrs0;
492     Mat_MPIAIJ   *mpimat2;
493     Mat           Gmat2;
494     Vec           locState;
495     PetscScalar   *cpcol_state;
496 
497     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
498     kk = nselected_1;
499     MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPIU_SUM, wcomm);
500     myCrs0 -= nselected_1;
501 
502     if (a_Gmat_2) { /* output */
503       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
504       ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
505       *a_Gmat_2 = Gmat2; /* output */
506     } else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
507     /* get coarse grid GIDS for selected (locals and ghosts) */
508     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
509     ierr = MatGetVecs(Gmat2, &locState, 0);CHKERRQ(ierr);
510     ierr = VecSet(locState, (PetscScalar)(PetscReal)(-1));CHKERRQ(ierr); /* set with UNKNOWN state */
511     for (kk=0;kk<nselected_1;kk++) {
512       PetscInt fgid = clid_lid_1[kk] + my0;
513       PetscScalar v = (PetscScalar)(kk+myCrs0);
514       ierr = VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES);CHKERRQ(ierr); /* set with PID */
515     }
516     ierr = VecAssemblyBegin(locState);CHKERRQ(ierr);
517     ierr = VecAssemblyEnd(locState);CHKERRQ(ierr);
518     ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
519     ierr =   VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
520     ierr = VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts);CHKERRQ(ierr);
521     ierr = VecGetArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
522     for (kk=0,num_crs_ghost=0;kk<num_fine_ghosts;kk++) {
523       if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++;
524     }
525     ierr = PetscMalloc((nselected_1+num_crs_ghost)*sizeof(PetscInt), &crsGID);CHKERRQ(ierr); /* output */
526     {
527       PetscInt *selected_set;
528       ierr = PetscMalloc((nselected_1+num_crs_ghost)*sizeof(PetscInt), &selected_set);CHKERRQ(ierr);
529       /* do ghost of 'crsGID' */
530       for (kk=0,idx=nselected_1;kk<num_fine_ghosts;kk++) {
531         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
532           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
533           selected_set[idx] = nloc + kk;
534           crsGID[idx++] = cgid;
535         }
536       }
537       assert(idx==(nselected_1+num_crs_ghost));
538       ierr = VecRestoreArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
539       /* do locals in 'crsGID' */
540       ierr = VecGetArray(locState, &cpcol_state);CHKERRQ(ierr);
541       for (kk=0,idx=0;kk<nloc;kk++) {
542         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
543           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
544           selected_set[idx] = kk;
545           crsGID[idx++] = cgid;
546         }
547       }
548       assert(idx==nselected_1);
549       ierr = VecRestoreArray(locState, &cpcol_state);CHKERRQ(ierr);
550 
551       if (a_selected_2 != 0) { /* output */
552         ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);CHKERRQ(ierr);
553       } else {
554         ierr = PetscFree(selected_set);CHKERRQ(ierr);
555       }
556     }
557     ierr = VecDestroy(&locState);CHKERRQ(ierr);
558   }
559   *a_crsGID = crsGID; /* output */
560 
561   PetscFunctionReturn(0);
562 }
563 
564 /* -------------------------------------------------------------------------- */
565 /*
566    PCGAMGgraph_GEO
567 
568   Input Parameter:
569    . pc - this
570    . Amat - matrix on this fine level
571   Output Parameter:
572    . a_Gmat
573 */
574 #undef __FUNCT__
575 #define __FUNCT__ "PCGAMGgraph_GEO"
576 PetscErrorCode PCGAMGgraph_GEO(PC pc,
577                                 const Mat Amat,
578                                 Mat *a_Gmat)
579 {
580   PetscErrorCode ierr;
581   PC_MG          *mg = (PC_MG*)pc->data;
582   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
583   const PetscInt verbose = pc_gamg->verbose;
584   const PetscReal vfilter = pc_gamg->threshold;
585   PetscMPIInt    mype,npe;
586   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
587   Mat            Gmat;
588   PetscBool  set,flg,symm;
589 
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 SETERRQ(wcomm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG");
851     { /* create next coords - output */
852       PetscReal *crs_crds;
853       ierr = PetscMalloc(dim*nLocalSelected*sizeof(PetscReal), &crs_crds);CHKERRQ(ierr);
854       for (kk=0;kk<nLocalSelected;kk++) {/* grab local select nodes to promote - output */
855         PetscInt lid = clid_flid[kk];
856         for (jj=0;jj<dim;jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
857       }
858 
859       ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
860       pc_gamg->data = crs_crds; /* out */
861       pc_gamg->data_sz = dim*nLocalSelected;
862     }
863     ierr = ISDestroy(&selected_2);CHKERRQ(ierr);
864   }
865 
866   *a_P_out = Prol;  /* out */
867   ierr = PetscFree(clid_flid);CHKERRQ(ierr);
868 #if defined PETSC_USE_LOG
869   ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
870 #endif
871   PetscFunctionReturn(0);
872 }
873 
874 /* -------------------------------------------------------------------------- */
875 /*
876  PCCreateGAMG_GEO
877 
878   Input Parameter:
879    . pc -
880 */
881 #undef __FUNCT__
882 #define __FUNCT__ "PCCreateGAMG_GEO"
883 PetscErrorCode  PCCreateGAMG_GEO(PC pc)
884 {
885   PetscErrorCode  ierr;
886   PC_MG           *mg = (PC_MG*)pc->data;
887   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
888 
889   PetscFunctionBegin;
890   pc->ops->setfromoptions = PCSetFromOptions_GEO;
891   /* pc->ops->destroy        = PCDestroy_GEO; */
892   /* reset does not do anything; setup not virtual */
893 
894   /* set internal function pointers */
895   pc_gamg->graph = PCGAMGgraph_GEO;
896   pc_gamg->coarsen = PCGAMGcoarsen_GEO;
897   pc_gamg->prolongator = PCGAMGProlongator_GEO;
898   pc_gamg->optprol = 0;
899   pc_gamg->formkktprol = 0;
900 
901   pc_gamg->createdefaultdata = PCSetData_GEO;
902 
903   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
904                                             "PCSetCoordinates_C",
905                                             "PCSetCoordinates_GEO",
906                                             PCSetCoordinates_GEO);CHKERRQ(ierr);
907 
908   PetscFunctionReturn(0);
909 }
910