xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision 77c8aea5b73b50700a790636621daa7beaee9d6e)
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, 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   ierr  = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr);
48   nloc = (Iend-my0)/bs;
49   if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
50 
51   pc_gamg->data_cell_rows = 1;
52   if( coords==0 && nloc > 0 ) {
53     SETERRQ(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
54   }
55   pc_gamg->data_cell_cols = ndm; /* coordinates */
56 
57   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
58 
59   /* create data - syntactic sugar that should be refactored at some point */
60   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
61     ierr = PetscFree( pc_gamg->data );  CHKERRQ(ierr);
62     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr);
63   }
64   for(kk=0;kk<arrsz;kk++)pc_gamg->data[kk] = -999.;
65   pc_gamg->data[arrsz] = -99.;
66   /* copy data in - column oriented */
67   for( kk = 0 ; kk < nloc ; kk++ ){
68     for( ii = 0 ; ii < ndm ; ii++ ) {
69       pc_gamg->data[ii*nloc + kk] =  coords[kk*ndm + ii];
70     }
71   }
72   assert(pc_gamg->data[arrsz] == -99.);
73 
74   pc_gamg->data_sz = arrsz;
75 
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 )
90 {
91   PetscFunctionBegin;
92   SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"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     PetscPrintf(PETSC_COMM_WORLD,"[%d]%s done\n",0,__FUNCT__);
130   }
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    . nnodes -
142    . coords[2*nnodes] - 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 nnodes,
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 {
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   PetscFunctionBegin;
177 
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   }
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[nnodes + 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[nnodes + 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       if( (jj=AILSizeAt(agg_lists_1,mm)) > 0 ) {
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         LLNPos         pos;
335 
336         for( pos=AILGetHeadPos(agg_lists_1,lid) ;
337              pos ;
338              pos=AILGetNextPos(agg_lists_1,lid,pos) ){
339           PetscInt flid = LLNGetID(pos);
340           if( flid < nFineLoc ) {  /* could be a ghost */
341             PetscInt bestTID = -1; PetscReal best_alpha = 1.e10;
342             const PetscInt fgid = flid + myFine0;
343             /* compute shape function for gid */
344             const PetscReal fcoord[3] = {coords[flid],coords[nnodes+flid],1.0};
345             PetscBool haveit=PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3];
346             /* look for it */
347             for( tid = node_tri[clid], jj=0;
348                  jj < 5 && !haveit && tid != -1;
349                  jj++ ){
350               for(tt=0;tt<3;tt++){
351                 PetscInt cid2 = mid.trianglelist[3*tid + tt];
352                 PetscInt lid2 = selected_idx_2[cid2];
353                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0;
354                 clids[tt] = cid2; /* store for interp */
355               }
356 
357               for(tt=0;tt<3;tt++) alpha[tt] = (PetscScalar)fcoord[tt];
358 
359               /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */
360               LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
361               {
362                 PetscBool have=PETSC_TRUE;  PetscReal lowest=1.e10;
363                 for( tt = 0, idx = 0 ; tt < 3 ; tt++ ) {
364                   if( PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS ) have = PETSC_FALSE;
365                   if( PetscRealPart(alpha[tt]) < lowest ){
366                     lowest = PetscRealPart(alpha[tt]);
367                     idx = tt;
368                   }
369                 }
370                 haveit = have;
371               }
372               tid = mid.neighborlist[3*tid + idx];
373             }
374 
375             if( !haveit ) {
376               /* brute force */
377               for(tid=0 ; tid<mid.numberoftriangles && !haveit ; tid++ ){
378                 for(tt=0;tt<3;tt++){
379                   PetscInt cid2 = mid.trianglelist[3*tid + tt];
380                   PetscInt lid2 = selected_idx_2[cid2];
381                   AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0;
382                   clids[tt] = cid2; /* store for interp */
383                 }
384                 for(tt=0;tt<3;tt++) alpha[tt] = fcoord[tt];
385                 /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */
386                 LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
387                 {
388                   PetscBool have=PETSC_TRUE;  PetscReal worst=0.0, v;
389                   for(tt=0; tt<3 && have ;tt++) {
390                     if( PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS ) have=PETSC_FALSE;
391                     if( (v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst ) worst = v;
392                   }
393                   if( worst < best_alpha ) {
394                     best_alpha = worst; bestTID = tid;
395                   }
396                   haveit = have;
397                 }
398               }
399             }
400             if( !haveit ) {
401               if( best_alpha > *a_worst_best ) *a_worst_best = best_alpha;
402               /* use best one */
403               for(tt=0;tt<3;tt++){
404                 PetscInt cid2 = mid.trianglelist[3*bestTID + tt];
405                 PetscInt lid2 = selected_idx_2[cid2];
406                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0;
407                 clids[tt] = cid2; /* store for interp */
408               }
409               for(tt=0;tt<3;tt++) alpha[tt] = fcoord[tt];
410               /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */
411               LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO);
412             }
413 
414             /* put in row of P */
415             for(idx=0;idx<3;idx++){
416               PetscScalar shp = alpha[idx];
417               if( PetscAbs(PetscRealPart(shp)) > 1.e-6 ) {
418                 PetscInt cgid = crsGID[clids[idx]];
419                 PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */
420                 for(tt=0 ; tt < bs ; tt++, ii++, jj++ ){
421                   ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES); CHKERRQ(ierr);
422                 }
423               }
424             }
425           }
426         } /* aggregates iterations */
427         clid++;
428       } /* a coarse agg */
429     } /* for all fine nodes */
430 
431     ierr = ISRestoreIndices( selected_2, &selected_idx_2 );     CHKERRQ(ierr);
432     ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
433     ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
434 
435     ierr = PetscFree( node_tri );  CHKERRQ(ierr);
436     ierr = PetscFree( nTri );  CHKERRQ(ierr);
437   }
438 #if defined PETSC_GAMG_USE_LOG
439   ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
440 #endif
441   free( mid.trianglelist );
442   free( mid.neighborlist );
443   ierr = PetscFree( in.pointlist );  CHKERRQ(ierr);
444 
445   PetscFunctionReturn(0);
446 #else
447   SETERRQ(((PetscObject)a_Prol)->comm,PETSC_ERR_LIB,"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 {
473   PetscErrorCode ierr;
474   PetscMPIInt    mype,npe;
475   PetscInt       *crsGID, kk,my0,Iend,nloc;
476   MPI_Comm       wcomm = ((PetscObject)Gmat1)->comm;
477 
478   PetscFunctionBegin;
479   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
480   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
481   ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend); CHKERRQ(ierr); /* AIJ */
482   nloc = Iend - my0; /* this does not change */
483 
484   if (npe == 1) { /* not much to do in serial */
485     ierr = PetscMalloc( nselected_1*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr);
486     for(kk=0;kk<nselected_1;kk++) crsGID[kk] = kk;
487     *a_Gmat_2 = 0;
488     ierr = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);
489     CHKERRQ(ierr);
490   }
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     }
508     else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
509     /* get coarse grid GIDS for selected (locals and ghosts) */
510     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
511     ierr = MatGetVecs( Gmat2, &locState, 0 );         CHKERRQ(ierr);
512     ierr = VecSet( locState, (PetscScalar)(PetscReal)(-1) );  CHKERRQ(ierr); /* set with UNKNOWN state */
513     for(kk=0;kk<nselected_1;kk++){
514       PetscInt fgid = clid_lid_1[kk] + my0;
515       PetscScalar v = (PetscScalar)(kk+myCrs0);
516       ierr = VecSetValues( locState, 1, &fgid, &v, INSERT_VALUES );  CHKERRQ(ierr); /* set with PID */
517     }
518     ierr = VecAssemblyBegin( locState ); CHKERRQ(ierr);
519     ierr = VecAssemblyEnd( locState ); CHKERRQ(ierr);
520     ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
521     ierr =   VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
522     ierr = VecGetLocalSize( mpimat2->lvec, &num_fine_ghosts ); CHKERRQ(ierr);
523     ierr = VecGetArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr);
524     for(kk=0,num_crs_ghost=0;kk<num_fine_ghosts;kk++){
525       if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ) num_crs_ghost++;
526     }
527     ierr = PetscMalloc( (nselected_1+num_crs_ghost)*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr); /* output */
528     {
529       PetscInt *selected_set;
530       ierr = PetscMalloc( (nselected_1+num_crs_ghost)*sizeof(PetscInt), &selected_set ); CHKERRQ(ierr);
531       /* do ghost of 'crsGID' */
532       for(kk=0,idx=nselected_1;kk<num_fine_ghosts;kk++){
533         if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){
534           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
535           selected_set[idx] = nloc + kk;
536           crsGID[idx++] = cgid;
537         }
538       }
539       assert(idx==(nselected_1+num_crs_ghost));
540       ierr = VecRestoreArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr);
541       /* do locals in 'crsGID' */
542       ierr = VecGetArray( locState, &cpcol_state ); CHKERRQ(ierr);
543       for(kk=0,idx=0;kk<nloc;kk++){
544         if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){
545           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
546           selected_set[idx] = kk;
547           crsGID[idx++] = cgid;
548         }
549       }
550       assert(idx==nselected_1);
551       ierr = VecRestoreArray( locState, &cpcol_state ); CHKERRQ(ierr);
552 
553       if( a_selected_2 != 0 ) { /* output */
554         ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);
555         CHKERRQ(ierr);
556       }
557       else {
558         ierr = PetscFree( selected_set );  CHKERRQ(ierr);
559       }
560     }
561     ierr = VecDestroy( &locState );                    CHKERRQ(ierr);
562   }
563   *a_crsGID = crsGID; /* output */
564 
565   PetscFunctionReturn(0);
566 }
567 
568 /* -------------------------------------------------------------------------- */
569 /*
570    PCGAMGgraph_GEO
571 
572   Input Parameter:
573    . pc - this
574    . Amat - matrix on this fine level
575   Output Parameter:
576    . a_Gmat
577 */
578 #undef __FUNCT__
579 #define __FUNCT__ "PCGAMGgraph_GEO"
580 PetscErrorCode PCGAMGgraph_GEO( PC pc,
581                                 const Mat Amat,
582                                 Mat *a_Gmat
583                                 )
584 {
585   PetscErrorCode ierr;
586   PC_MG          *mg = (PC_MG*)pc->data;
587   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
588   const PetscInt verbose = pc_gamg->verbose;
589   const PetscReal vfilter = pc_gamg->threshold;
590   PetscMPIInt    mype,npe;
591   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
592   Mat            Gmat;
593   PetscBool  set,flg,symm;
594   PetscFunctionBegin;
595 #if defined PETSC_USE_LOG
596   ierr = PetscLogEventBegin(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr);
597 #endif
598   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
599   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
600 
601   ierr = MatIsSymmetricKnown(Amat, &set, &flg);        CHKERRQ(ierr);
602   symm = !(set && flg);
603 
604   ierr  = PCGAMGCreateSimpleGraph( Amat, &Gmat ); CHKERRQ( ierr );
605   ierr  = PCGAMGScaleFilterGraph( &Gmat, vfilter, symm, verbose ); CHKERRQ( ierr );
606 
607   *a_Gmat = Gmat;
608 #if defined PETSC_USE_LOG
609   ierr = PetscLogEventEnd(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr);
610 #endif
611   PetscFunctionReturn(0);
612 }
613 
614 /* -------------------------------------------------------------------------- */
615 /*
616    PCGAMGcoarsen_GEO
617 
618   Input Parameter:
619    . a_pc - this
620    . a_Gmat - graph
621   Output Parameter:
622    . a_llist_parent - linked list from selected indices for data locality only
623 */
624 #undef __FUNCT__
625 #define __FUNCT__ "PCGAMGcoarsen_GEO"
626 PetscErrorCode PCGAMGcoarsen_GEO( PC a_pc,
627                                   Mat *a_Gmat,
628                                   PetscCoarsenData **a_llist_parent
629                                   )
630 {
631   PetscErrorCode ierr;
632   PC_MG          *mg = (PC_MG*)a_pc->data;
633   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
634   PetscInt       Istart,Iend,nloc,kk,Ii,ncols;
635   PetscMPIInt    mype,npe;
636   IS             perm;
637   GAMGNode *gnodes;
638   PetscInt *permute;
639   Mat       Gmat = *a_Gmat;
640   MPI_Comm  wcomm = ((PetscObject)Gmat)->comm;
641   MatCoarsen crs;
642 
643   PetscFunctionBegin;
644 #if defined PETSC_USE_LOG
645   ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
646 #endif
647   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
648   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
649   ierr = MatGetOwnershipRange( Gmat, &Istart, &Iend ); CHKERRQ(ierr);
650   nloc = (Iend-Istart);
651 
652   /* create random permutation with sort for geo-mg */
653   ierr = PetscMalloc( nloc*sizeof(GAMGNode), &gnodes ); CHKERRQ(ierr);
654   ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr);
655 
656   for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */
657     ierr = MatGetRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr);
658     {
659       PetscInt lid = Ii - Istart;
660       gnodes[lid].lid = lid;
661       gnodes[lid].degree = ncols;
662     }
663     ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr);
664   }
665   /* randomize */
666   srand(1); /* make deterministic */
667   if( PETSC_TRUE ) {
668     PetscBool *bIndexSet;
669     ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr);
670     for ( Ii = 0; Ii < nloc ; Ii++) bIndexSet[Ii] = PETSC_FALSE;
671     for ( Ii = 0; Ii < nloc ; Ii++)
672     {
673       PetscInt iSwapIndex = rand()%nloc;
674       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii)
675       {
676         GAMGNode iTemp = gnodes[iSwapIndex];
677         gnodes[iSwapIndex] = gnodes[Ii];
678         gnodes[Ii] = iTemp;
679         bIndexSet[Ii] = PETSC_TRUE;
680         bIndexSet[iSwapIndex] = PETSC_TRUE;
681       }
682     }
683     ierr = PetscFree( bIndexSet );  CHKERRQ(ierr);
684   }
685   /* only sort locals */
686   qsort( gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare );
687   /* create IS of permutation */
688   for(kk=0;kk<nloc;kk++) { /* locals only */
689     permute[kk] = gnodes[kk].lid;
690   }
691   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);
692   CHKERRQ(ierr);
693 
694   ierr = PetscFree( gnodes );  CHKERRQ(ierr);
695 
696   /* get MIS aggs */
697 
698   ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr);
699   ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr);
700   ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr);
701   ierr = MatCoarsenSetAdjacency( crs, Gmat ); CHKERRQ(ierr);
702   ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr);
703   ierr = MatCoarsenSetStrictAggs( crs, PETSC_FALSE ); CHKERRQ(ierr);
704   ierr = MatCoarsenApply( crs ); CHKERRQ(ierr);
705   ierr = MatCoarsenGetData( crs, a_llist_parent ); CHKERRQ(ierr);
706   ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr);
707 
708   ierr = ISDestroy( &perm );                    CHKERRQ(ierr);
709 #if defined PETSC_USE_LOG
710   ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
711 #endif
712   PetscFunctionReturn(0);
713 }
714 
715 /* -------------------------------------------------------------------------- */
716 /*
717  PCGAMGProlongator_GEO
718 
719  Input Parameter:
720  . pc - this
721  . Amat - matrix on this fine level
722  . Graph - used to get ghost data for nodes in
723  . selected_1 - [nselected]
724  . agg_lists - [nselected]
725  Output Parameter:
726  . a_P_out - prolongation operator to the next level
727  */
728 #undef __FUNCT__
729 #define __FUNCT__ "PCGAMGProlongator_GEO"
730 PetscErrorCode PCGAMGProlongator_GEO( PC pc,
731                                       const Mat Amat,
732                                       const Mat Gmat,
733                                       PetscCoarsenData *agg_lists,
734                                       Mat *a_P_out
735                                       )
736 {
737   PC_MG          *mg = (PC_MG*)pc->data;
738   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
739   const PetscInt  verbose = pc_gamg->verbose;
740   const PetscInt  dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
741   PetscErrorCode ierr;
742   PetscInt       Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid;
743   Mat            Prol;
744   PetscMPIInt    mype, npe;
745   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
746   IS             selected_2,selected_1;
747   const PetscInt *selected_idx;
748 
749   PetscFunctionBegin;
750 #if defined PETSC_USE_LOG
751   ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
752 #endif
753   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
754   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
755   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
756   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
757   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
758 
759   /* get 'nLocalSelected' */
760   ierr = AILGetMIS( agg_lists, &selected_1 );        CHKERRQ(ierr);
761   ierr = ISGetSize( selected_1, &jj );               CHKERRQ(ierr);
762   ierr = PetscMalloc( jj*sizeof(PetscInt), &clid_flid ); CHKERRQ(ierr);
763   ierr = ISGetIndices( selected_1, &selected_idx );     CHKERRQ(ierr);
764   for(kk=0,nLocalSelected=0;kk<jj;kk++) {
765     PetscInt lid = selected_idx[kk];
766     if( lid<nloc ) {
767       ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0); CHKERRQ(ierr);
768       if( ncols>1 ) { /* fiter out singletons */
769         clid_flid[nLocalSelected++] = lid;
770       }
771       else assert(0); /* filtered in coarsening */
772       ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0); CHKERRQ(ierr);
773     }
774   }
775   ierr = ISRestoreIndices( selected_1, &selected_idx );     CHKERRQ(ierr);
776 
777   /* create prolongator, create P matrix */
778   ierr = MatCreateAIJ(wcomm,
779                       nloc*bs, nLocalSelected*bs,
780                       PETSC_DETERMINE, PETSC_DETERMINE,
781                       3*data_cols, PETSC_NULL, /* don't have a good way to set this!!! */
782                       3*data_cols, PETSC_NULL,
783                       &Prol );
784   CHKERRQ(ierr);
785 
786   /* can get all points "removed" - but not on geomg */
787   ierr =  MatGetSize( Prol, &kk, &jj ); CHKERRQ(ierr);
788   if( jj==0 ) {
789     if( verbose ) {
790       PetscPrintf(PETSC_COMM_WORLD,"[%d]%s ERROE: no selected points on coarse grid\n",mype,__FUNCT__);
791     }
792     ierr = PetscFree( clid_flid );  CHKERRQ(ierr);
793     ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
794     *a_P_out = PETSC_NULL;  /* out */
795     PetscFunctionReturn(0);
796   }
797 
798   {
799     PetscReal *coords;
800     PetscInt nnodes;
801     PetscInt  *crsGID;
802     Mat        Gmat2;
803 
804     assert(dim==data_cols);
805     /* grow ghost data for better coarse grid cover of fine grid */
806 #if defined PETSC_GAMG_USE_LOG
807     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
808 #endif
809     ierr = getGIDsOnSquareGraph( nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID );
810     CHKERRQ(ierr);
811     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
812 #if defined PETSC_GAMG_USE_LOG
813     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
814 #endif
815     /* create global vector of coorindates in 'coords' */
816     if (npe > 1) {
817       ierr = PCGAMGGetDataWithGhosts( Gmat2, dim, pc_gamg->data, &nnodes, &coords );
818       CHKERRQ(ierr);
819     }
820     else {
821       coords = (PetscReal*)pc_gamg->data;
822       nnodes = nloc;
823     }
824     ierr = MatDestroy( &Gmat2 );  CHKERRQ(ierr);
825 
826     /* triangulate */
827     if( dim == 2 ) {
828       PetscReal metric,tm;
829 #if defined PETSC_GAMG_USE_LOG
830       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
831 #endif
832       ierr = triangulateAndFormProl( selected_2, nnodes, coords,
833                                      nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric );
834       CHKERRQ(ierr);
835 #if defined PETSC_GAMG_USE_LOG
836       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0); CHKERRQ(ierr);
837 #endif
838       ierr = PetscFree( crsGID );  CHKERRQ(ierr);
839 
840       /* clean up and create coordinates for coarse grid (output) */
841       if (npe > 1) ierr = PetscFree( coords ); CHKERRQ(ierr);
842 
843       ierr = MPI_Allreduce( &metric, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm );  CHKERRQ(ierr);
844       if( tm > 1. ) { /* needs to be globalized - should not happen */
845         if( verbose ) {
846           PetscPrintf(PETSC_COMM_WORLD,"[%d]%s failed metric for coarse grid %e\n",mype,__FUNCT__,tm);
847         }
848         ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
849         Prol = PETSC_NULL;
850       }
851       else if( metric > .0 ) {
852         if( verbose ) {
853           PetscPrintf(PETSC_COMM_WORLD,"[%d]%s worst metric for coarse grid = %e\n",mype,__FUNCT__,metric);
854         }
855       }
856     } else {
857       SETERRQ(wcomm,PETSC_ERR_LIB,"3D not implemented for 'geo' AMG");
858     }
859     { /* create next coords - output */
860       PetscReal *crs_crds;
861       ierr = PetscMalloc( dim*nLocalSelected*sizeof(PetscReal), &crs_crds );
862       CHKERRQ(ierr);
863       for(kk=0;kk<nLocalSelected;kk++){/* grab local select nodes to promote - output */
864         PetscInt lid = clid_flid[kk];
865         for(jj=0;jj<dim;jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
866       }
867 
868       ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
869       pc_gamg->data = crs_crds; /* out */
870       pc_gamg->data_sz = dim*nLocalSelected;
871     }
872     ierr = ISDestroy( &selected_2 ); CHKERRQ(ierr); /* this is selected_1 in serial */
873   }
874   *a_P_out = Prol;  /* out */
875   ierr = PetscFree( clid_flid );  CHKERRQ(ierr);
876 #if defined PETSC_USE_LOG
877   ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
878 #endif
879   PetscFunctionReturn(0);
880 }
881 
882 /* -------------------------------------------------------------------------- */
883 /*
884  PCCreateGAMG_GEO
885 
886   Input Parameter:
887    . pc -
888 */
889 #undef __FUNCT__
890 #define __FUNCT__ "PCCreateGAMG_GEO"
891 PetscErrorCode  PCCreateGAMG_GEO( PC pc )
892 {
893   PetscErrorCode  ierr;
894   PC_MG           *mg = (PC_MG*)pc->data;
895   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
896 
897   PetscFunctionBegin;
898   pc->ops->setfromoptions = PCSetFromOptions_GEO;
899   /* pc->ops->destroy        = PCDestroy_GEO; */
900   /* reset does not do anything; setup not virtual */
901 
902   /* set internal function pointers */
903   pc_gamg->graph = PCGAMGgraph_GEO;
904   pc_gamg->coarsen = PCGAMGcoarsen_GEO;
905   pc_gamg->prolongator = PCGAMGProlongator_GEO;
906   pc_gamg->optprol = 0;
907 
908   pc_gamg->createdefaultdata = PCSetData_GEO;
909 
910   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
911                                             "PCSetCoordinates_C",
912                                             "PCSetCoordinates_GEO",
913                                             PCSetCoordinates_GEO);CHKERRQ(ierr);
914 
915   PetscFunctionReturn(0);
916 }
917