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