Actual source code: geo.c

petsc-master 2019-08-19
Report Typos and Errors
  1: /*
  2:  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
  3:  */

  5:  #include <../src/ksp/pc/impls/gamg/gamg.h>

  7: #if defined(PETSC_HAVE_TRIANGLE)
  8: #if !defined(ANSI_DECLARATORS)
  9: #define ANSI_DECLARATORS
 10: #endif
 11: #include <triangle.h>
 12: #endif

 14:  #include <petscblaslapack.h>

 16: /* Private context for the GAMG preconditioner */
 17: typedef struct {
 18:   PetscInt lid;            /* local vertex index */
 19:   PetscInt degree;         /* vertex degree */
 20: } GAMGNode;

 22: PETSC_STATIC_INLINE int petsc_geo_mg_compare(const void *a, const void *b)
 23: {
 24:   return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree);
 25: }

 27: /* -------------------------------------------------------------------------- */
 28: /*
 29:    PCSetCoordinates_GEO

 31:    Input Parameter:
 32:    .  pc - the preconditioner context
 33: */
 34: PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
 35: {
 36:   PC_MG          *mg      = (PC_MG*)pc->data;
 37:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
 39:   PetscInt       arrsz,bs,my0,kk,ii,nloc,Iend,aloc;
 40:   Mat            Amat = pc->pmat;

 44:   MatGetBlockSize(Amat, &bs);
 45:   MatGetOwnershipRange(Amat, &my0, &Iend);
 46:   aloc = (Iend-my0);
 47:   nloc = (Iend-my0)/bs;

 49:   if (nloc!=a_nloc && aloc!=a_nloc) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of local blocks %D must be %D or %D.",a_nloc,nloc,aloc);

 51:   pc_gamg->data_cell_rows = 1;
 52:   if (!coords && nloc > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
 53:   pc_gamg->data_cell_cols = ndm; /* coordinates */

 55:   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;

 57:   /* create data - syntactic sugar that should be refactored at some point */
 58:   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
 59:     PetscFree(pc_gamg->data);
 60:     PetscMalloc1(arrsz+1, &pc_gamg->data);
 61:   }
 62:   for (kk=0; kk<arrsz; kk++) pc_gamg->data[kk] = -999.;
 63:   pc_gamg->data[arrsz] = -99.;
 64:   /* copy data in - column oriented */
 65:   if (nloc == a_nloc) {
 66:     for (kk = 0; kk < nloc; kk++) {
 67:       for (ii = 0; ii < ndm; ii++) {
 68:         pc_gamg->data[ii*nloc + kk] =  coords[kk*ndm + ii];
 69:       }
 70:     }
 71:   } else { /* assumes the coordinates are blocked */
 72:     for (kk = 0; kk < nloc; kk++) {
 73:       for (ii = 0; ii < ndm; ii++) {
 74:         pc_gamg->data[ii*nloc + kk] =  coords[bs*kk*ndm + ii];
 75:       }
 76:     }
 77:   }
 78:   if (pc_gamg->data[arrsz] != -99.) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data[arrsz %D] %g != -99.",arrsz,pc_gamg->data[arrsz]);
 79:   pc_gamg->data_sz = arrsz;
 80:   return(0);
 81: }

 83: /* -------------------------------------------------------------------------- */
 84: /*
 85:    PCSetData_GEO

 87:   Input Parameter:
 88:    . pc -
 89: */
 90: PetscErrorCode PCSetData_GEO(PC pc, Mat m)
 91: {
 93:   SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"GEO MG needs coordinates");
 94: }

 96: /* -------------------------------------------------------------------------- */
 97: /*
 98:    PCSetFromOptions_GEO

100:   Input Parameter:
101:    . pc -
102: */
103: PetscErrorCode PCSetFromOptions_GEO(PetscOptionItems *PetscOptionsObject,PC pc)
104: {

108:   PetscOptionsHead(PetscOptionsObject,"GAMG-GEO options");
109:   {
110:     /* -pc_gamg_sa_nsmooths */
111:     /* pc_gamg_sa->smooths = 0; */
112:     /* PetscOptionsInt("-pc_gamg_agg_nsmooths", */
113:     /*                        "smoothing steps for smoothed aggregation, usually 1 (0)", */
114:     /*                        "PCGAMGSetNSmooths_AGG", */
115:     /*                        pc_gamg_sa->smooths, */
116:     /*                        &pc_gamg_sa->smooths, */
117:     /*                        &flag);  */
118:     /*  */
119:   }
120:   PetscOptionsTail();
121:   return(0);
122: }

124: /* -------------------------------------------------------------------------- */
125: /*
126:  triangulateAndFormProl

128:    Input Parameter:
129:    . selected_2 - list of selected local ID, includes selected ghosts
130:    . data_stride -
131:    . coords[2*data_stride] - column vector of local coordinates w/ ghosts
132:    . nselected_1 - selected IDs that go with base (1) graph includes selected ghosts
133:    . clid_lid_1[nselected_1] - lids of selected (c) nodes   ???????????
134:    . agg_lists_1 - list of aggregates selected_1 vertices of aggregate unselected vertices
135:    . crsGID[selected.size()] - global index for prolongation operator
136:    . bs - block size
137:   Output Parameter:
138:    . a_Prol - prolongation operator
139:    . a_worst_best - measure of worst missed fine vertex, 0 is no misses
140: */
141: static PetscErrorCode triangulateAndFormProl(IS selected_2,PetscInt data_stride,PetscReal coords[],PetscInt nselected_1,const PetscInt clid_lid_1[],const PetscCoarsenData *agg_lists_1,
142:                                              const PetscInt crsGID[],PetscInt bs,Mat a_Prol,PetscReal *a_worst_best)
143: {
144: #if defined(PETSC_HAVE_TRIANGLE)
145:   PetscErrorCode       ierr;
146:   PetscInt             jj,tid,tt,idx,nselected_2;
147:   struct triangulateio in,mid;
148:   const PetscInt       *selected_idx_2;
149:   PetscMPIInt          rank;
150:   PetscInt             Istart,Iend,nFineLoc,myFine0;
151:   int                  kk,nPlotPts,sid;
152:   MPI_Comm             comm;
153:   PetscReal            tm;

156:   PetscObjectGetComm((PetscObject)a_Prol,&comm);
157:   MPI_Comm_rank(comm,&rank);
158:   ISGetSize(selected_2, &nselected_2);
159:   if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */
160:     *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */
161:   } else *a_worst_best = 0.0;
162:   MPIU_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, comm);
163:   if (tm > 0.0) {
164:     *a_worst_best = 100.0;
165:     return(0);
166:   }
167:   MatGetOwnershipRange(a_Prol, &Istart, &Iend);
168:   nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs;
169:   nPlotPts = nFineLoc; /* locals */
170:   /* triangle */
171:   /* Define input points - in*/
172:   in.numberofpoints          = nselected_2;
173:   in.numberofpointattributes = 0;
174:   /* get nselected points */
175:   PetscMalloc1(2*nselected_2, &in.pointlist);
176:   ISGetIndices(selected_2, &selected_idx_2);

178:   for (kk=0,sid=0; kk<nselected_2; kk++,sid += 2) {
179:     PetscInt lid = selected_idx_2[kk];
180:     in.pointlist[sid]   = coords[lid];
181:     in.pointlist[sid+1] = coords[data_stride + lid];
182:     if (lid>=nFineLoc) nPlotPts++;
183:   }
184:   if (sid != 2*nselected_2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != 2*nselected_2 %D",sid,nselected_2);

186:   in.numberofsegments      = 0;
187:   in.numberofedges         = 0;
188:   in.numberofholes         = 0;
189:   in.numberofregions       = 0;
190:   in.trianglelist          = 0;
191:   in.segmentmarkerlist     = 0;
192:   in.pointattributelist    = 0;
193:   in.pointmarkerlist       = 0;
194:   in.triangleattributelist = 0;
195:   in.trianglearealist      = 0;
196:   in.segmentlist           = 0;
197:   in.holelist              = 0;
198:   in.regionlist            = 0;
199:   in.edgelist              = 0;
200:   in.edgemarkerlist        = 0;
201:   in.normlist              = 0;

203:   /* triangulate */
204:   mid.pointlist = 0;            /* Not needed if -N switch used. */
205:   /* Not needed if -N switch used or number of point attributes is zero: */
206:   mid.pointattributelist = 0;
207:   mid.pointmarkerlist    = 0; /* Not needed if -N or -B switch used. */
208:   mid.trianglelist       = 0;    /* Not needed if -E switch used. */
209:   /* Not needed if -E switch used or number of triangle attributes is zero: */
210:   mid.triangleattributelist = 0;
211:   mid.neighborlist          = 0; /* Needed only if -n switch used. */
212:   /* Needed only if segments are output (-p or -c) and -P not used: */
213:   mid.segmentlist = 0;
214:   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
215:   mid.segmentmarkerlist = 0;
216:   mid.edgelist          = 0;    /* Needed only if -e switch used. */
217:   mid.edgemarkerlist    = 0; /* Needed if -e used and -B not used. */
218:   mid.numberoftriangles = 0;

220:   /* Triangulate the points.  Switches are chosen to read and write a  */
221:   /*   PSLG (p), preserve the convex hull (c), number everything from  */
222:   /*   zero (z), assign a regional attribute to each element (A), and  */
223:   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
224:   /*   neighbor list (n).                                            */
225:   if (nselected_2 != 0) { /* inactive processor */
226:     char args[] = "npczQ"; /* c is needed ? */
227:     triangulate(args, &in, &mid, (struct triangulateio*) NULL);
228:     /* output .poly files for 'showme' */
229:     if (!PETSC_TRUE) {
230:       static int level = 1;
231:       FILE       *file; char fname[32];

233:       sprintf(fname,"C%d_%d.poly",level,rank); file = fopen(fname, "w");
234:       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
235:       fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0);
236:       /*Following lines: <vertex #> <x> <y> */
237:       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid += 2) {
238:         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
239:       }
240:       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
241:       fprintf(file, "%d  %d\n",0,0);
242:       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
243:       /* One line: <# of holes> */
244:       fprintf(file, "%d\n",0);
245:       /* Following lines: <hole #> <x> <y> */
246:       /* Optional line: <# of regional attributes and/or area constraints> */
247:       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
248:       fclose(file);

250:       /* elems */
251:       sprintf(fname,"C%d_%d.ele",level,rank); file = fopen(fname, "w");
252:       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
253:       fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0);
254:       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
255:       for (kk=0,sid=0; kk<mid.numberoftriangles; kk++,sid += 3) {
256:         fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]);
257:       }
258:       fclose(file);

260:       sprintf(fname,"C%d_%d.node",level,rank); file = fopen(fname, "w");
261:       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
262:       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
263:       fprintf(file, "%d  %d  %d  %d\n",nPlotPts,2,0,0);
264:       /*Following lines: <vertex #> <x> <y> */
265:       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid+=2) {
266:         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
267:       }

269:       sid /= 2;
270:       for (jj=0; jj<nFineLoc; jj++) {
271:         PetscBool sel = PETSC_TRUE;
272:         for (kk=0; kk<nselected_2 && sel; kk++) {
273:           PetscInt lid = selected_idx_2[kk];
274:           if (lid == jj) sel = PETSC_FALSE;
275:         }
276:         if (sel) fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]);
277:       }
278:       fclose(file);
279:       if (sid != nPlotPts) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != nPlotPts %D",sid,nPlotPts);
280:       level++;
281:     }
282:   }
283: #if defined PETSC_GAMG_USE_LOG
284:   PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);
285: #endif
286:   { /* form P - setup some maps */
287:     PetscInt clid,mm,*nTri,*node_tri;

289:     PetscMalloc2(nselected_2, &node_tri,nselected_2, &nTri);

291:     /* need list of triangles on node */
292:     for (kk=0; kk<nselected_2; kk++) nTri[kk] = 0;
293:     for (tid=0,kk=0; tid<mid.numberoftriangles; tid++) {
294:       for (jj=0; jj<3; jj++) {
295:         PetscInt cid = mid.trianglelist[kk++];
296:         if (nTri[cid] == 0) node_tri[cid] = tid;
297:         nTri[cid]++;
298:       }
299:     }
300: #define EPS 1.e-12
301:     /* find points and set prolongation */
302:     for (mm = clid = 0; mm < nFineLoc; mm++) {
303:       PetscBool ise;
304:       PetscCDEmptyAt(agg_lists_1,mm,&ise);
305:       if (!ise) {
306:         const PetscInt lid = mm;
307:         PetscScalar    AA[3][3];
308:         PetscBLASInt   N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO;
309:         PetscCDIntNd   *pos;

311:         PetscCDGetHeadPos(agg_lists_1,lid,&pos);
312:         while (pos) {
313:           PetscInt flid;
314:           PetscCDIntNdGetID(pos, &flid);
315:           PetscCDGetNextPos(agg_lists_1,lid,&pos);

317:           if (flid < nFineLoc) {  /* could be a ghost */
318:             PetscInt       bestTID = -1; PetscReal best_alpha = 1.e10;
319:             const PetscInt fgid    = flid + myFine0;
320:             /* compute shape function for gid */
321:             const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0};
322:             PetscBool       haveit    =PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3];

324:             /* look for it */
325:             for (tid = node_tri[clid], jj=0;
326:                  jj < 5 && !haveit && tid != -1;
327:                  jj++) {
328:               for (tt=0; tt<3; tt++) {
329:                 PetscInt cid2 = mid.trianglelist[3*tid + tt];
330:                 PetscInt lid2 = selected_idx_2[cid2];
331:                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
332:                 clids[tt] = cid2; /* store for interp */
333:               }

335:               for (tt=0; tt<3; tt++) alpha[tt] = (PetscScalar)fcoord[tt];

337:               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
338:               PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
339:               {
340:                 PetscBool have=PETSC_TRUE;  PetscReal lowest=1.e10;
341:                 for (tt = 0, idx = 0; tt < 3; tt++) {
342:                   if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
343:                   if (PetscRealPart(alpha[tt]) < lowest) {
344:                     lowest = PetscRealPart(alpha[tt]);
345:                     idx    = tt;
346:                   }
347:                 }
348:                 haveit = have;
349:               }
350:               tid = mid.neighborlist[3*tid + idx];
351:             }

353:             if (!haveit) {
354:               /* brute force */
355:               for (tid=0; tid<mid.numberoftriangles && !haveit; tid++) {
356:                 for (tt=0; tt<3; tt++) {
357:                   PetscInt cid2 = mid.trianglelist[3*tid + tt];
358:                   PetscInt lid2 = selected_idx_2[cid2];
359:                   AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
360:                   clids[tt] = cid2; /* store for interp */
361:                 }
362:                 for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
363:                 /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
364:                 PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
365:                 {
366:                   PetscBool have=PETSC_TRUE;  PetscReal worst=0.0, v;
367:                   for (tt=0; tt<3 && have; tt++) {
368:                     if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE;
369:                     if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v;
370:                   }
371:                   if (worst < best_alpha) {
372:                     best_alpha = worst; bestTID = tid;
373:                   }
374:                   haveit = have;
375:                 }
376:               }
377:             }
378:             if (!haveit) {
379:               if (best_alpha > *a_worst_best) *a_worst_best = best_alpha;
380:               /* use best one */
381:               for (tt=0; tt<3; tt++) {
382:                 PetscInt cid2 = mid.trianglelist[3*bestTID + tt];
383:                 PetscInt lid2 = selected_idx_2[cid2];
384:                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
385:                 clids[tt] = cid2; /* store for interp */
386:               }
387:               for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
388:               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
389:               PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
390:             }

392:             /* put in row of P */
393:             for (idx=0; idx<3; idx++) {
394:               PetscScalar shp = alpha[idx];
395:               if (PetscAbs(PetscRealPart(shp)) > 1.e-6) {
396:                 PetscInt cgid = crsGID[clids[idx]];
397:                 PetscInt jj   = cgid*bs, ii = fgid*bs; /* need to gloalize */
398:                 for (tt=0; tt < bs; tt++, ii++, jj++) {
399:                   MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES);
400:                 }
401:               }
402:             }
403:           }
404:         } /* aggregates iterations */
405:         clid++;
406:       } /* a coarse agg */
407:     } /* for all fine nodes */

409:     ISRestoreIndices(selected_2, &selected_idx_2);
410:     MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);
411:     MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);

413:     PetscFree2(node_tri,nTri);
414:   }
415: #if defined PETSC_GAMG_USE_LOG
416:   PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);
417: #endif
418:   free(mid.trianglelist);
419:   free(mid.neighborlist);
420:   PetscFree(in.pointlist);
421:   return(0);
422: #else
423:   SETERRQ(PetscObjectComm((PetscObject)a_Prol),PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG");
424: #endif
425: }
426: /* -------------------------------------------------------------------------- */
427: /*
428:    getGIDsOnSquareGraph - square graph, get

430:    Input Parameter:
431:    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
432:    . clid_lid_1 - [nselected_1] lids of selected nodes
433:    . Gmat1 - graph that goes with 'selected_1'
434:    Output Parameter:
435:    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
436:    . a_Gmat_2 - graph that is squared of 'Gmat_1'
437:    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
438: */
439: static PetscErrorCode getGIDsOnSquareGraph(PetscInt nselected_1,const PetscInt clid_lid_1[],const Mat Gmat1,IS *a_selected_2,Mat *a_Gmat_2,PetscInt **a_crsGID)
440: {
442:   PetscMPIInt    size;
443:   PetscInt       *crsGID, kk,my0,Iend,nloc;
444:   MPI_Comm       comm;

447:   PetscObjectGetComm((PetscObject)Gmat1,&comm);
448:   MPI_Comm_size(comm,&size);
449:   MatGetOwnershipRange(Gmat1,&my0,&Iend); /* AIJ */
450:   nloc = Iend - my0; /* this does not change */

452:   if (size == 1) { /* not much to do in serial */
453:     PetscMalloc1(nselected_1, &crsGID);
454:     for (kk=0; kk<nselected_1; kk++) crsGID[kk] = kk;
455:     *a_Gmat_2 = 0;
456:     ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);
457:   } else {
458:     PetscInt    idx,num_fine_ghosts,num_crs_ghost,myCrs0;
459:     Mat_MPIAIJ  *mpimat2;
460:     Mat         Gmat2;
461:     Vec         locState;
462:     PetscScalar *cpcol_state;

464:     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
465:     kk = nselected_1;
466:     MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPI_SUM, comm);
467:     myCrs0 -= nselected_1;

469:     if (a_Gmat_2) { /* output */
470:       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
471:       MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);
472:       *a_Gmat_2 = Gmat2; /* output */
473:     } else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
474:     /* get coarse grid GIDS for selected (locals and ghosts) */
475:     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
476:     MatCreateVecs(Gmat2, &locState, 0);
477:     VecSet(locState, (PetscScalar)(PetscReal)(-1)); /* set with UNKNOWN state */
478:     for (kk=0; kk<nselected_1; kk++) {
479:       PetscInt    fgid = clid_lid_1[kk] + my0;
480:       PetscScalar v    = (PetscScalar)(kk+myCrs0);
481:       VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES); /* set with PID */
482:     }
483:     VecAssemblyBegin(locState);
484:     VecAssemblyEnd(locState);
485:     VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);
486:       VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);
487:     VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts);
488:     VecGetArray(mpimat2->lvec, &cpcol_state);
489:     for (kk=0,num_crs_ghost=0; kk<num_fine_ghosts; kk++) {
490:       if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++;
491:     }
492:     PetscMalloc1(nselected_1+num_crs_ghost, &crsGID); /* output */
493:     {
494:       PetscInt *selected_set;
495:       PetscMalloc1(nselected_1+num_crs_ghost, &selected_set);
496:       /* do ghost of 'crsGID' */
497:       for (kk=0,idx=nselected_1; kk<num_fine_ghosts; kk++) {
498:         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
499:           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
500:           selected_set[idx] = nloc + kk;
501:           crsGID[idx++]     = cgid;
502:         }
503:       }
504:       if (idx != (nselected_1+num_crs_ghost)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != (nselected_1 %D + num_crs_ghost %D)",idx,nselected_1,num_crs_ghost);
505:       VecRestoreArray(mpimat2->lvec, &cpcol_state);
506:       /* do locals in 'crsGID' */
507:       VecGetArray(locState, &cpcol_state);
508:       for (kk=0,idx=0; kk<nloc; kk++) {
509:         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
510:           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
511:           selected_set[idx] = kk;
512:           crsGID[idx++]     = cgid;
513:         }
514:       }
515:       if (idx != nselected_1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != nselected_1 %D",idx,nselected_1);
516:       VecRestoreArray(locState, &cpcol_state);

518:       if (a_selected_2 != 0) { /* output */
519:         ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);
520:       } else {
521:         PetscFree(selected_set);
522:       }
523:     }
524:     VecDestroy(&locState);
525:   }
526:   *a_crsGID = crsGID; /* output */
527:   return(0);
528: }

530: /* -------------------------------------------------------------------------- */
531: /*
532:    PCGAMGGraph_GEO

534:   Input Parameter:
535:    . pc - this
536:    . Amat - matrix on this fine level
537:   Output Parameter:
538:    . a_Gmat
539: */
540: PetscErrorCode PCGAMGGraph_GEO(PC pc,Mat Amat,Mat *a_Gmat)
541: {
542:   PetscErrorCode  ierr;
543:   PC_MG           *mg      = (PC_MG*)pc->data;
544:   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
545:   const PetscReal vfilter  = pc_gamg->threshold[0];
546:   MPI_Comm        comm;
547:   Mat             Gmat;
548:   PetscBool       set,flg,symm;

551:   PetscObjectGetComm((PetscObject)Amat,&comm);
552:   PetscLogEventBegin(PC_GAMGGraph_GEO,0,0,0,0);

554:   MatIsSymmetricKnown(Amat, &set, &flg);
555:   symm = (PetscBool)!(set && flg);

557:   PCGAMGCreateGraph(Amat, &Gmat);
558:   PCGAMGFilterGraph(&Gmat, vfilter, symm);

560:   *a_Gmat = Gmat;
561:   PetscLogEventEnd(PC_GAMGGraph_GEO,0,0,0,0);
562:   return(0);
563: }

565: /* -------------------------------------------------------------------------- */
566: /*
567:    PCGAMGCoarsen_GEO

569:   Input Parameter:
570:    . a_pc - this
571:    . a_Gmat - graph
572:   Output Parameter:
573:    . a_llist_parent - linked list from selected indices for data locality only
574: */
575: PetscErrorCode PCGAMGCoarsen_GEO(PC a_pc,Mat *a_Gmat,PetscCoarsenData **a_llist_parent)
576: {
578:   PetscInt       Istart,Iend,nloc,kk,Ii,ncols;
579:   IS             perm;
580:   GAMGNode       *gnodes;
581:   PetscInt       *permute;
582:   Mat            Gmat  = *a_Gmat;
583:   MPI_Comm       comm;
584:   MatCoarsen     crs;

587:   PetscObjectGetComm((PetscObject)a_pc,&comm);
588:   PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);
589:   MatGetOwnershipRange(Gmat, &Istart, &Iend);
590:   nloc = (Iend-Istart);

592:   /* create random permutation with sort for geo-mg */
593:   PetscMalloc1(nloc, &gnodes);
594:   PetscMalloc1(nloc, &permute);

596:   for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */
597:     MatGetRow(Gmat,Ii,&ncols,0,0);
598:     {
599:       PetscInt lid = Ii - Istart;
600:       gnodes[lid].lid    = lid;
601:       gnodes[lid].degree = ncols;
602:     }
603:     MatRestoreRow(Gmat,Ii,&ncols,0,0);
604:   }
605:   if (PETSC_TRUE) {
606:     PetscRandom  rand;
607:     PetscBool    *bIndexSet;
608:     PetscReal    rr;
609:     PetscInt     iSwapIndex;

611:     PetscRandomCreate(comm,&rand);
612:     PetscCalloc1(nloc, &bIndexSet);
613:     for (Ii = 0; Ii < nloc; Ii++) {
614:       PetscRandomGetValueReal(rand,&rr);
615:       iSwapIndex = (PetscInt) (rr*nloc);
616:       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
617:         GAMGNode iTemp = gnodes[iSwapIndex];
618:         gnodes[iSwapIndex]    = gnodes[Ii];
619:         gnodes[Ii]            = iTemp;
620:         bIndexSet[Ii]         = PETSC_TRUE;
621:         bIndexSet[iSwapIndex] = PETSC_TRUE;
622:       }
623:     }
624:     PetscRandomDestroy(&rand);
625:     PetscFree(bIndexSet);
626:   }
627:   /* only sort locals */
628:   qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare);
629:   /* create IS of permutation */
630:   for (kk=0; kk<nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */
631:   PetscFree(gnodes);
632:   ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);

634:   /* get MIS aggs */

636:   MatCoarsenCreate(comm, &crs);
637:   MatCoarsenSetType(crs, MATCOARSENMIS);
638:   MatCoarsenSetGreedyOrdering(crs, perm);
639:   MatCoarsenSetAdjacency(crs, Gmat);
640:   MatCoarsenSetStrictAggs(crs, PETSC_FALSE);
641:   MatCoarsenApply(crs);
642:   MatCoarsenGetData(crs, a_llist_parent);
643:   MatCoarsenDestroy(&crs);

645:   ISDestroy(&perm);
646:   PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);
647:   return(0);
648: }

650: /* -------------------------------------------------------------------------- */
651: /*
652:  PCGAMGProlongator_GEO

654:  Input Parameter:
655:  . pc - this
656:  . Amat - matrix on this fine level
657:  . Graph - used to get ghost data for nodes in
658:  . selected_1 - [nselected]
659:  . agg_lists - [nselected]
660:  Output Parameter:
661:  . a_P_out - prolongation operator to the next level
662:  */
663: PetscErrorCode PCGAMGProlongator_GEO(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
664: {
665:   PC_MG          *mg      = (PC_MG*)pc->data;
666:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
667:   const PetscInt dim      = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
669:   PetscInt       Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid;
670:   Mat            Prol;
671:   PetscMPIInt    rank, size;
672:   MPI_Comm       comm;
673:   IS             selected_2,selected_1;
674:   const PetscInt *selected_idx;
675:   MatType        mtype;

678:   PetscObjectGetComm((PetscObject)Amat,&comm);
679:   PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);
680:   MPI_Comm_rank(comm,&rank);
681:   MPI_Comm_size(comm,&size);
682:   MatGetOwnershipRange(Amat, &Istart, &Iend);
683:   MatGetBlockSize(Amat, &bs);
684:   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
685:   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) % bs %D",Iend,Istart,bs);

687:   /* get 'nLocalSelected' */
688:   PetscCDGetMIS(agg_lists, &selected_1);
689:   ISGetSize(selected_1, &jj);
690:   PetscMalloc1(jj, &clid_flid);
691:   ISGetIndices(selected_1, &selected_idx);
692:   for (kk=0,nLocalSelected=0; kk<jj; kk++) {
693:     PetscInt lid = selected_idx[kk];
694:     if (lid<nloc) {
695:       MatGetRow(Gmat,lid+my0,&ncols,0,0);
696:       if (ncols>1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */
697:       MatRestoreRow(Gmat,lid+my0,&ncols,0,0);
698:     }
699:   }
700:   ISRestoreIndices(selected_1, &selected_idx);
701:   ISDestroy(&selected_1); /* this is selected_1 in serial */

703:   /* create prolongator  matrix */
704:   MatGetType(Amat,&mtype);
705:   MatCreate(comm, &Prol);
706:   MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE);
707:   MatSetBlockSizes(Prol, bs, bs);
708:   MatSetType(Prol, mtype);
709:   MatSeqAIJSetPreallocation(Prol,3*data_cols,NULL);
710:   MatMPIAIJSetPreallocation(Prol,3*data_cols,NULL,3*data_cols,NULL);

712:   /* can get all points "removed" - but not on geomg */
713:    MatGetSize(Prol, &kk, &jj);
714:   if (!jj) {
715:     PetscInfo(pc,"ERROE: no selected points on coarse grid\n");
716:     PetscFree(clid_flid);
717:     MatDestroy(&Prol);
718:     *a_P_out = NULL;  /* out */
719:     return(0);
720:   }

722:   {
723:     PetscReal *coords;
724:     PetscInt  data_stride;
725:     PetscInt  *crsGID = NULL;
726:     Mat       Gmat2;

728:     if (dim != data_cols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim %D != data_cols %D",dim,data_cols);
729:     /* grow ghost data for better coarse grid cover of fine grid */
730: #if defined PETSC_GAMG_USE_LOG
731:     PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);
732: #endif
733:     /* messy method, squares graph and gets some data */
734:     getGIDsOnSquareGraph(nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID);
735:     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
736: #if defined PETSC_GAMG_USE_LOG
737:     PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);
738: #endif
739:     /* create global vector of coorindates in 'coords' */
740:     if (size > 1) {
741:       PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords);
742:     } else {
743:       coords      = (PetscReal*)pc_gamg->data;
744:       data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols;
745:     }
746:     MatDestroy(&Gmat2);

748:     /* triangulate */
749:     if (dim == 2) {
750:       PetscReal metric,tm;
751: #if defined PETSC_GAMG_USE_LOG
752:       PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);
753: #endif
754:       triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric);
755: #if defined PETSC_GAMG_USE_LOG
756:       PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0);
757: #endif
758:       PetscFree(crsGID);

760:       /* clean up and create coordinates for coarse grid (output) */
761:       if (size > 1) PetscFree(coords);

763:       MPIU_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, comm);
764:       if (tm > 1.) { /* needs to be globalized - should not happen */
765:         PetscInfo1(pc," failed metric for coarse grid %e\n",(double)tm);
766:         MatDestroy(&Prol);
767:       } else if (metric > .0) {
768:         PetscInfo1(pc,"worst metric for coarse grid = %e\n",(double)metric);
769:       }
770:     } else SETERRQ(comm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG");
771:     { /* create next coords - output */
772:       PetscReal *crs_crds;
773:       PetscMalloc1(dim*nLocalSelected, &crs_crds);
774:       for (kk=0; kk<nLocalSelected; kk++) { /* grab local select nodes to promote - output */
775:         PetscInt lid = clid_flid[kk];
776:         for (jj=0; jj<dim; jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
777:       }

779:       PetscFree(pc_gamg->data);
780:       pc_gamg->data    = crs_crds; /* out */
781:       pc_gamg->data_sz = dim*nLocalSelected;
782:     }
783:     ISDestroy(&selected_2);
784:   }

786:   *a_P_out = Prol;  /* out */
787:   PetscFree(clid_flid);
788:   PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);
789:   return(0);
790: }

792: static PetscErrorCode PCDestroy_GAMG_GEO(PC pc)
793: {

797:   PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);
798:   return(0);
799: }

801: /* -------------------------------------------------------------------------- */
802: /*
803:  PCCreateGAMG_GEO

805:   Input Parameter:
806:    . pc -
807: */
808: PetscErrorCode  PCCreateGAMG_GEO(PC pc)
809: {
811:   PC_MG          *mg      = (PC_MG*)pc->data;
812:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

815:   pc_gamg->ops->setfromoptions = PCSetFromOptions_GEO;
816:   pc_gamg->ops->destroy        = PCDestroy_GAMG_GEO;
817:   /* reset does not do anything; setup not virtual */

819:   /* set internal function pointers */
820:   pc_gamg->ops->graph             = PCGAMGGraph_GEO;
821:   pc_gamg->ops->coarsen           = PCGAMGCoarsen_GEO;
822:   pc_gamg->ops->prolongator       = PCGAMGProlongator_GEO;
823:   pc_gamg->ops->optprolongator    = NULL;
824:   pc_gamg->ops->createdefaultdata = PCSetData_GEO;

826:   PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_GEO);
827:   return(0);
828: }