xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision 5f042095cd6794de959e41580753d6ff4dec299c)
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     MPI_Comm  wcomm = ((PetscObject)pc)->comm;
130     PetscPrintf(wcomm,"[%d]%s done\n",0,__FUNCT__);
131   }
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    . nnodes -
143    . coords[2*nnodes] - 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 nnodes,
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 {
167 #if defined(PETSC_HAVE_TRIANGLE)
168   PetscErrorCode       ierr;
169   PetscInt             jj,tid,tt,idx,nselected_2;
170   struct triangulateio in,mid;
171   const PetscInt      *selected_idx_2;
172   PetscMPIInt          mype,npe;
173   PetscInt             Istart,Iend,nFineLoc,myFine0;
174   int                  kk,nPlotPts,sid;
175   MPI_Comm             wcomm = ((PetscObject)a_Prol)->comm;
176   PetscReal            tm;
177   PetscFunctionBegin;
178 
179   ierr = MPI_Comm_rank(wcomm,&mype);    CHKERRQ(ierr);
180   ierr = MPI_Comm_size(wcomm,&npe);     CHKERRQ(ierr);
181   ierr = ISGetSize( selected_2, &nselected_2 );        CHKERRQ(ierr);
182   if(nselected_2 == 1 || nselected_2 == 2 ){ /* 0 happens on idle processors */
183     *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */
184   }
185   else *a_worst_best = 0.0;
186   ierr = MPI_Allreduce( a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm );  CHKERRQ(ierr);
187   if( tm > 0.0 ) {
188     *a_worst_best = 100.0;
189     PetscFunctionReturn(0);
190   }
191   ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend );  CHKERRQ(ierr);
192   nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs;
193   nPlotPts = nFineLoc; /* locals */
194   /* traingle */
195   /* Define input points - in*/
196   in.numberofpoints = nselected_2;
197   in.numberofpointattributes = 0;
198   /* get nselected points */
199   ierr = PetscMalloc( 2*(nselected_2)*sizeof(REAL), &in.pointlist ); CHKERRQ(ierr);
200   ierr = ISGetIndices( selected_2, &selected_idx_2 );     CHKERRQ(ierr);
201 
202   for(kk=0,sid=0;kk<nselected_2;kk++,sid += 2){
203     PetscInt lid = selected_idx_2[kk];
204     in.pointlist[sid] = coords[lid];
205     in.pointlist[sid+1] = coords[nnodes + lid];
206     if(lid>=nFineLoc) nPlotPts++;
207   }
208   assert(sid==2*nselected_2);
209 
210   in.numberofsegments = 0;
211   in.numberofedges = 0;
212   in.numberofholes = 0;
213   in.numberofregions = 0;
214   in.trianglelist = 0;
215   in.segmentmarkerlist = 0;
216   in.pointattributelist = 0;
217   in.pointmarkerlist = 0;
218   in.triangleattributelist = 0;
219   in.trianglearealist = 0;
220   in.segmentlist = 0;
221   in.holelist = 0;
222   in.regionlist = 0;
223   in.edgelist = 0;
224   in.edgemarkerlist = 0;
225   in.normlist = 0;
226   /* triangulate */
227   mid.pointlist = 0;            /* Not needed if -N switch used. */
228   /* Not needed if -N switch used or number of point attributes is zero: */
229   mid.pointattributelist = 0;
230   mid.pointmarkerlist = 0; /* Not needed if -N or -B switch used. */
231   mid.trianglelist = 0;          /* Not needed if -E switch used. */
232   /* Not needed if -E switch used or number of triangle attributes is zero: */
233   mid.triangleattributelist = 0;
234   mid.neighborlist = 0;         /* Needed only if -n switch used. */
235   /* Needed only if segments are output (-p or -c) and -P not used: */
236   mid.segmentlist = 0;
237   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
238   mid.segmentmarkerlist = 0;
239   mid.edgelist = 0;             /* Needed only if -e switch used. */
240   mid.edgemarkerlist = 0;   /* Needed if -e used and -B not used. */
241   mid.numberoftriangles = 0;
242 
243   /* Triangulate the points.  Switches are chosen to read and write a  */
244   /*   PSLG (p), preserve the convex hull (c), number everything from  */
245   /*   zero (z), assign a regional attribute to each element (A), and  */
246   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
247   /*   neighbor list (n).                                            */
248   if(nselected_2 != 0){ /* inactive processor */
249     char args[] = "npczQ"; /* c is needed ? */
250     triangulate(args, &in, &mid, (struct triangulateio *) NULL );
251     /* output .poly files for 'showme' */
252     if( !PETSC_TRUE ) {
253       static int level = 1;
254       FILE *file; char fname[32];
255 
256       sprintf(fname,"C%d_%d.poly",level,mype); file = fopen(fname, "w");
257       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
258       fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0);
259       /*Following lines: <vertex #> <x> <y> */
260       for(kk=0,sid=0;kk<in.numberofpoints;kk++,sid += 2){
261         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
262       }
263       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
264       fprintf(file, "%d  %d\n",0,0);
265       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
266       /* One line: <# of holes> */
267       fprintf(file, "%d\n",0);
268       /* Following lines: <hole #> <x> <y> */
269       /* Optional line: <# of regional attributes and/or area constraints> */
270       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
271       fclose(file);
272 
273       /* elems */
274       sprintf(fname,"C%d_%d.ele",level,mype); file = fopen(fname, "w");
275       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
276       fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0);
277       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
278       for(kk=0,sid=0;kk<mid.numberoftriangles;kk++,sid += 3){
279         fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]);
280       }
281       fclose(file);
282 
283       sprintf(fname,"C%d_%d.node",level,mype); file = fopen(fname, "w");
284       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
285       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
286       fprintf(file, "%d  %d  %d  %d\n",nPlotPts,2,0,0);
287       /*Following lines: <vertex #> <x> <y> */
288       for(kk=0,sid=0;kk<in.numberofpoints;kk++,sid+=2){
289         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
290       }
291 
292       sid /= 2;
293       for(jj=0;jj<nFineLoc;jj++){
294         PetscBool sel = PETSC_TRUE;
295         for( kk=0 ; kk<nselected_2 && sel ; kk++ ){
296           PetscInt lid = selected_idx_2[kk];
297           if( lid == jj ) sel = PETSC_FALSE;
298         }
299         if( sel ) {
300           fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[nnodes + jj]);
301         }
302       }
303       fclose(file);
304       assert(sid==nPlotPts);
305       level++;
306     }
307   }
308 #if defined PETSC_GAMG_USE_LOG
309   ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
310 #endif
311   { /* form P - setup some maps */
312     PetscInt clid,mm,*nTri,*node_tri;
313 
314     ierr = PetscMalloc( nselected_2*sizeof(PetscInt), &node_tri ); CHKERRQ(ierr);
315     ierr = PetscMalloc( nselected_2*sizeof(PetscInt), &nTri ); CHKERRQ(ierr);
316 
317     /* need list of triangles on node */
318     for(kk=0;kk<nselected_2;kk++) nTri[kk] = 0;
319     for(tid=0,kk=0;tid<mid.numberoftriangles;tid++){
320       for(jj=0;jj<3;jj++) {
321         PetscInt cid = mid.trianglelist[kk++];
322         if( nTri[cid] == 0 ) node_tri[cid] = tid;
323         nTri[cid]++;
324       }
325     }
326 #define EPS 1.e-12
327     /* find points and set prolongation */
328     for( mm = clid = 0 ; mm < nFineLoc ; mm++ ){
329       if( (jj=PetscCDSizeAt(agg_lists_1,mm)) > 0 ) {
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 
337         for( pos=PetscCDGetHeadPos(agg_lists_1,lid) ;
338              pos ;
339              pos=PetscCDGetNextPos(agg_lists_1,lid,pos) ){
340           PetscInt flid = LLNGetID(pos);
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[nnodes+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[nnodes + 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[nnodes + 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[nnodes + 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 
446   PetscFunctionReturn(0);
447 #else
448   SETERRQ(((PetscObject)a_Prol)->comm,PETSC_ERR_LIB,"configure with TRIANGLE to use geometric MG");
449 #endif
450 }
451 /* -------------------------------------------------------------------------- */
452 /*
453    getGIDsOnSquareGraph - square graph, get
454 
455    Input Parameter:
456    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
457    . clid_lid_1 - [nselected_1] lids of selected nodes
458    . Gmat1 - graph that goes with 'selected_1'
459    Output Parameter:
460    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
461    . a_Gmat_2 - graph that is squared of 'Gmat_1'
462    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
463 */
464 #undef __FUNCT__
465 #define __FUNCT__ "getGIDsOnSquareGraph"
466 static PetscErrorCode getGIDsOnSquareGraph( const PetscInt nselected_1,
467                                             const PetscInt clid_lid_1[],
468                                             const Mat Gmat1,
469                                             IS *a_selected_2,
470                                             Mat *a_Gmat_2,
471                                             PetscInt **a_crsGID
472                                             )
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);
490     CHKERRQ(ierr);
491   }
492   else {
493     PetscInt      idx,num_fine_ghosts,num_crs_ghost,myCrs0;
494     Mat_MPIAIJ   *mpimat2;
495     Mat           Gmat2;
496     Vec           locState;
497     PetscScalar   *cpcol_state;
498 
499     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
500     kk = nselected_1;
501     MPI_Scan( &kk, &myCrs0, 1, MPIU_INT, MPIU_SUM, wcomm );
502     myCrs0 -= nselected_1;
503 
504     if( a_Gmat_2 ) { /* output */
505       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
506       ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 );   CHKERRQ(ierr);
507       *a_Gmat_2 = Gmat2; /* output */
508     }
509     else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
510     /* get coarse grid GIDS for selected (locals and ghosts) */
511     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
512     ierr = MatGetVecs( Gmat2, &locState, 0 );         CHKERRQ(ierr);
513     ierr = VecSet( locState, (PetscScalar)(PetscReal)(-1) );  CHKERRQ(ierr); /* set with UNKNOWN state */
514     for(kk=0;kk<nselected_1;kk++){
515       PetscInt fgid = clid_lid_1[kk] + my0;
516       PetscScalar v = (PetscScalar)(kk+myCrs0);
517       ierr = VecSetValues( locState, 1, &fgid, &v, INSERT_VALUES );  CHKERRQ(ierr); /* set with PID */
518     }
519     ierr = VecAssemblyBegin( locState ); CHKERRQ(ierr);
520     ierr = VecAssemblyEnd( locState ); CHKERRQ(ierr);
521     ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
522     ierr =   VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
523     ierr = VecGetLocalSize( mpimat2->lvec, &num_fine_ghosts ); CHKERRQ(ierr);
524     ierr = VecGetArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr);
525     for(kk=0,num_crs_ghost=0;kk<num_fine_ghosts;kk++){
526       if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ) num_crs_ghost++;
527     }
528     ierr = PetscMalloc( (nselected_1+num_crs_ghost)*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr); /* output */
529     {
530       PetscInt *selected_set;
531       ierr = PetscMalloc( (nselected_1+num_crs_ghost)*sizeof(PetscInt), &selected_set ); CHKERRQ(ierr);
532       /* do ghost of 'crsGID' */
533       for(kk=0,idx=nselected_1;kk<num_fine_ghosts;kk++){
534         if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){
535           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
536           selected_set[idx] = nloc + kk;
537           crsGID[idx++] = cgid;
538         }
539       }
540       assert(idx==(nselected_1+num_crs_ghost));
541       ierr = VecRestoreArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr);
542       /* do locals in 'crsGID' */
543       ierr = VecGetArray( locState, &cpcol_state ); CHKERRQ(ierr);
544       for(kk=0,idx=0;kk<nloc;kk++){
545         if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){
546           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
547           selected_set[idx] = kk;
548           crsGID[idx++] = cgid;
549         }
550       }
551       assert(idx==nselected_1);
552       ierr = VecRestoreArray( locState, &cpcol_state ); CHKERRQ(ierr);
553 
554       if( a_selected_2 != 0 ) { /* output */
555         ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);
556         CHKERRQ(ierr);
557       }
558       else {
559         ierr = PetscFree( selected_set );  CHKERRQ(ierr);
560       }
561     }
562     ierr = VecDestroy( &locState );                    CHKERRQ(ierr);
563   }
564   *a_crsGID = crsGID; /* output */
565 
566   PetscFunctionReturn(0);
567 }
568 
569 /* -------------------------------------------------------------------------- */
570 /*
571    PCGAMGgraph_GEO
572 
573   Input Parameter:
574    . pc - this
575    . Amat - matrix on this fine level
576   Output Parameter:
577    . a_Gmat
578 */
579 #undef __FUNCT__
580 #define __FUNCT__ "PCGAMGgraph_GEO"
581 PetscErrorCode PCGAMGgraph_GEO( PC pc,
582                                 const Mat Amat,
583                                 Mat *a_Gmat
584                                 )
585 {
586   PetscErrorCode ierr;
587   PC_MG          *mg = (PC_MG*)pc->data;
588   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
589   const PetscInt verbose = pc_gamg->verbose;
590   const PetscReal vfilter = pc_gamg->threshold;
591   PetscMPIInt    mype,npe;
592   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
593   Mat            Gmat;
594   PetscBool  set,flg,symm;
595   PetscFunctionBegin;
596 #if defined PETSC_USE_LOG
597   ierr = PetscLogEventBegin(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr);
598 #endif
599   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
600   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
601 
602   ierr = MatIsSymmetricKnown(Amat, &set, &flg);        CHKERRQ(ierr);
603   symm = (PetscBool)!(set && flg);
604 
605   ierr  = PCGAMGCreateSimpleGraph( Amat, &Gmat ); CHKERRQ( ierr );
606   ierr  = PCGAMGScaleFilterGraph( &Gmat, vfilter, symm, verbose ); CHKERRQ( ierr );
607 
608   *a_Gmat = Gmat;
609 #if defined PETSC_USE_LOG
610   ierr = PetscLogEventEnd(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr);
611 #endif
612   PetscFunctionReturn(0);
613 }
614 
615 /* -------------------------------------------------------------------------- */
616 /*
617    PCGAMGcoarsen_GEO
618 
619   Input Parameter:
620    . a_pc - this
621    . a_Gmat - graph
622   Output Parameter:
623    . a_llist_parent - linked list from selected indices for data locality only
624 */
625 #undef __FUNCT__
626 #define __FUNCT__ "PCGAMGcoarsen_GEO"
627 PetscErrorCode PCGAMGcoarsen_GEO( PC a_pc,
628                                   Mat *a_Gmat,
629                                   PetscCoarsenData **a_llist_parent
630                                   )
631 {
632   PetscErrorCode ierr;
633   PC_MG          *mg = (PC_MG*)a_pc->data;
634   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
635   PetscInt       Istart,Iend,nloc,kk,Ii,ncols;
636   PetscMPIInt    mype,npe;
637   IS             perm;
638   GAMGNode *gnodes;
639   PetscInt *permute;
640   Mat       Gmat = *a_Gmat;
641   MPI_Comm  wcomm = ((PetscObject)Gmat)->comm;
642   MatCoarsen crs;
643 
644   PetscFunctionBegin;
645 #if defined PETSC_USE_LOG
646   ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
647 #endif
648   ierr = MPI_Comm_rank( wcomm, &mype);  CHKERRQ(ierr);
649   ierr = MPI_Comm_size( wcomm, &npe);   CHKERRQ(ierr);
650   ierr = MatGetOwnershipRange( Gmat, &Istart, &Iend ); CHKERRQ(ierr);
651   nloc = (Iend-Istart);
652 
653   /* create random permutation with sort for geo-mg */
654   ierr = PetscMalloc( nloc*sizeof(GAMGNode), &gnodes ); CHKERRQ(ierr);
655   ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr);
656 
657   for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */
658     ierr = MatGetRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr);
659     {
660       PetscInt lid = Ii - Istart;
661       gnodes[lid].lid = lid;
662       gnodes[lid].degree = ncols;
663     }
664     ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr);
665   }
666   /* randomize */
667   srand(1); /* make deterministic */
668   if( PETSC_TRUE ) {
669     PetscBool *bIndexSet;
670     ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr);
671     for ( Ii = 0; Ii < nloc ; Ii++) bIndexSet[Ii] = PETSC_FALSE;
672     for ( Ii = 0; Ii < nloc ; Ii++)
673     {
674       PetscInt iSwapIndex = rand()%nloc;
675       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii)
676       {
677         GAMGNode iTemp = gnodes[iSwapIndex];
678         gnodes[iSwapIndex] = gnodes[Ii];
679         gnodes[Ii] = iTemp;
680         bIndexSet[Ii] = PETSC_TRUE;
681         bIndexSet[iSwapIndex] = PETSC_TRUE;
682       }
683     }
684     ierr = PetscFree( bIndexSet );  CHKERRQ(ierr);
685   }
686   /* only sort locals */
687   qsort( gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare );
688   /* create IS of permutation */
689   for(kk=0;kk<nloc;kk++) { /* locals only */
690     permute[kk] = gnodes[kk].lid;
691   }
692   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);
693   CHKERRQ(ierr);
694 
695   ierr = PetscFree( gnodes );  CHKERRQ(ierr);
696 
697   /* get MIS aggs */
698 
699   ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr);
700   ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr);
701   ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr);
702   ierr = MatCoarsenSetAdjacency( crs, Gmat ); CHKERRQ(ierr);
703   ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr);
704   ierr = MatCoarsenSetStrictAggs( crs, PETSC_FALSE ); CHKERRQ(ierr);
705   ierr = MatCoarsenApply( crs ); CHKERRQ(ierr);
706   ierr = MatCoarsenGetData( crs, a_llist_parent ); CHKERRQ(ierr);
707   ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr);
708 
709   ierr = ISDestroy( &perm );                    CHKERRQ(ierr);
710 #if defined PETSC_USE_LOG
711   ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
712 #endif
713   PetscFunctionReturn(0);
714 }
715 
716 /* -------------------------------------------------------------------------- */
717 /*
718  PCGAMGProlongator_GEO
719 
720  Input Parameter:
721  . pc - this
722  . Amat - matrix on this fine level
723  . Graph - used to get ghost data for nodes in
724  . selected_1 - [nselected]
725  . agg_lists - [nselected]
726  Output Parameter:
727  . a_P_out - prolongation operator to the next level
728  */
729 #undef __FUNCT__
730 #define __FUNCT__ "PCGAMGProlongator_GEO"
731 PetscErrorCode PCGAMGProlongator_GEO( PC pc,
732                                       const Mat Amat,
733                                       const Mat Gmat,
734                                       PetscCoarsenData *agg_lists,
735                                       Mat *a_P_out
736                                       )
737 {
738   PC_MG          *mg = (PC_MG*)pc->data;
739   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
740   const PetscInt  verbose = pc_gamg->verbose;
741   const PetscInt  dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
742   PetscErrorCode ierr;
743   PetscInt       Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid;
744   Mat            Prol;
745   PetscMPIInt    mype, npe;
746   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
747   IS             selected_2,selected_1;
748   const PetscInt *selected_idx;
749 
750   PetscFunctionBegin;
751 #if defined PETSC_USE_LOG
752   ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
753 #endif
754   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
755   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
756   ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr);
757   ierr  = MatGetBlockSize( Amat, &bs );               CHKERRQ( ierr );
758   nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0);
759 
760   /* get 'nLocalSelected' */
761   ierr = PetscCDGetMIS( agg_lists, &selected_1 );        CHKERRQ(ierr);
762   ierr = ISGetSize( selected_1, &jj );               CHKERRQ(ierr);
763   ierr = PetscMalloc( jj*sizeof(PetscInt), &clid_flid ); CHKERRQ(ierr);
764   ierr = ISGetIndices( selected_1, &selected_idx );     CHKERRQ(ierr);
765   for(kk=0,nLocalSelected=0;kk<jj;kk++) {
766     PetscInt lid = selected_idx[kk];
767     if( lid<nloc ) {
768       ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0); CHKERRQ(ierr);
769       if( ncols>1 ) { /* fiter out singletons */
770         clid_flid[nLocalSelected++] = lid;
771       }
772       else assert(0); /* filtered in coarsening */
773       ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0); CHKERRQ(ierr);
774     }
775   }
776   ierr = ISRestoreIndices( selected_1, &selected_idx );     CHKERRQ(ierr);
777 
778   /* create prolongator, create P matrix */
779   ierr = MatCreateAIJ(wcomm,
780                       nloc*bs, nLocalSelected*bs,
781                       PETSC_DETERMINE, PETSC_DETERMINE,
782                       3*data_cols, PETSC_NULL, /* don't have a good way to set this!!! */
783                       3*data_cols, PETSC_NULL,
784                       &Prol );
785   CHKERRQ(ierr);
786 
787   /* can get all points "removed" - but not on geomg */
788   ierr =  MatGetSize( Prol, &kk, &jj ); CHKERRQ(ierr);
789   if( jj==0 ) {
790     if( verbose ) {
791       PetscPrintf(wcomm,"[%d]%s ERROE: no selected points on coarse grid\n",mype,__FUNCT__);
792     }
793     ierr = PetscFree( clid_flid );  CHKERRQ(ierr);
794     ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
795     *a_P_out = PETSC_NULL;  /* out */
796     PetscFunctionReturn(0);
797   }
798 
799   {
800     PetscReal *coords;
801     PetscInt nnodes;
802     PetscInt  *crsGID;
803     Mat        Gmat2;
804 
805     assert(dim==data_cols);
806     /* grow ghost data for better coarse grid cover of fine grid */
807 #if defined PETSC_GAMG_USE_LOG
808     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
809 #endif
810     ierr = getGIDsOnSquareGraph( nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID );
811     CHKERRQ(ierr);
812     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
813 #if defined PETSC_GAMG_USE_LOG
814     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
815 #endif
816     /* create global vector of coorindates in 'coords' */
817     if (npe > 1) {
818       ierr = PCGAMGGetDataWithGhosts( Gmat2, dim, pc_gamg->data, &nnodes, &coords );
819       CHKERRQ(ierr);
820     }
821     else {
822       coords = (PetscReal*)pc_gamg->data;
823       nnodes = nloc;
824     }
825     ierr = MatDestroy( &Gmat2 );  CHKERRQ(ierr);
826 
827     /* triangulate */
828     if( dim == 2 ) {
829       PetscReal metric,tm;
830 #if defined PETSC_GAMG_USE_LOG
831       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
832 #endif
833       ierr = triangulateAndFormProl( selected_2, nnodes, coords,
834                                      nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric );
835       CHKERRQ(ierr);
836 #if defined PETSC_GAMG_USE_LOG
837       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0); CHKERRQ(ierr);
838 #endif
839       ierr = PetscFree( crsGID );  CHKERRQ(ierr);
840 
841       /* clean up and create coordinates for coarse grid (output) */
842       if (npe > 1) ierr = PetscFree( coords ); CHKERRQ(ierr);
843 
844       ierr = MPI_Allreduce( &metric, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm );  CHKERRQ(ierr);
845       if( tm > 1. ) { /* needs to be globalized - should not happen */
846         if( verbose ) {
847           PetscPrintf(wcomm,"[%d]%s failed metric for coarse grid %e\n",mype,__FUNCT__,tm);
848         }
849         ierr = MatDestroy( &Prol );  CHKERRQ(ierr);
850         Prol = PETSC_NULL;
851       }
852       else if( metric > .0 ) {
853         if( verbose ) {
854           PetscPrintf(wcomm,"[%d]%s worst metric for coarse grid = %e\n",mype,__FUNCT__,metric);
855         }
856       }
857     } else {
858       SETERRQ(wcomm,PETSC_ERR_LIB,"3D not implemented for 'geo' AMG");
859     }
860     { /* create next coords - output */
861       PetscReal *crs_crds;
862       ierr = PetscMalloc( dim*nLocalSelected*sizeof(PetscReal), &crs_crds );
863       CHKERRQ(ierr);
864       for(kk=0;kk<nLocalSelected;kk++){/* grab local select nodes to promote - output */
865         PetscInt lid = clid_flid[kk];
866         for(jj=0;jj<dim;jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
867       }
868 
869       ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr );
870       pc_gamg->data = crs_crds; /* out */
871       pc_gamg->data_sz = dim*nLocalSelected;
872     }
873     ierr = ISDestroy( &selected_2 ); CHKERRQ(ierr); /* this is selected_1 in serial */
874   }
875   *a_P_out = Prol;  /* out */
876   ierr = PetscFree( clid_flid );  CHKERRQ(ierr);
877 #if defined PETSC_USE_LOG
878   ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
879 #endif
880   PetscFunctionReturn(0);
881 }
882 
883 /* -------------------------------------------------------------------------- */
884 /*
885  PCCreateGAMG_GEO
886 
887   Input Parameter:
888    . pc -
889 */
890 #undef __FUNCT__
891 #define __FUNCT__ "PCCreateGAMG_GEO"
892 PetscErrorCode  PCCreateGAMG_GEO( PC pc )
893 {
894   PetscErrorCode  ierr;
895   PC_MG           *mg = (PC_MG*)pc->data;
896   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
897 
898   PetscFunctionBegin;
899   pc->ops->setfromoptions = PCSetFromOptions_GEO;
900   /* pc->ops->destroy        = PCDestroy_GEO; */
901   /* reset does not do anything; setup not virtual */
902 
903   /* set internal function pointers */
904   pc_gamg->graph = PCGAMGgraph_GEO;
905   pc_gamg->coarsen = PCGAMGcoarsen_GEO;
906   pc_gamg->prolongator = PCGAMGProlongator_GEO;
907   pc_gamg->optprol = 0;
908 
909   pc_gamg->createdefaultdata = PCSetData_GEO;
910 
911   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
912                                             "PCSetCoordinates_C",
913                                             "PCSetCoordinates_GEO",
914                                             PCSetCoordinates_GEO);CHKERRQ(ierr);
915 
916   PetscFunctionReturn(0);
917 }
918