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