Actual source code: gamg.c

petsc-3.7.0 2016-04-25
Report Typos and Errors
  1: /*
  2:  GAMG geometric-algebric multigrid PC - Mark Adams 2011
  3:  */
  4: #include <petsc/private/matimpl.h>
  5: #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
  6: #include <petsc/private/kspimpl.h>
  7: #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */

  9: #if defined PETSC_GAMG_USE_LOG
 10: PetscLogEvent petsc_gamg_setup_events[NUM_SET];
 11: #endif

 13: #if defined PETSC_USE_LOG
 14: PetscLogEvent PC_GAMGGraph_AGG;
 15: PetscLogEvent PC_GAMGGraph_GEO;
 16: PetscLogEvent PC_GAMGCoarsen_AGG;
 17: PetscLogEvent PC_GAMGCoarsen_GEO;
 18: PetscLogEvent PC_GAMGProlongator_AGG;
 19: PetscLogEvent PC_GAMGProlongator_GEO;
 20: PetscLogEvent PC_GAMGOptProlongator_AGG;
 21: #endif

 23: #define GAMG_MAXLEVELS 30

 25: /* #define GAMG_STAGES */
 26: #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
 27: static PetscLogStage gamg_stages[GAMG_MAXLEVELS];
 28: #endif

 30: static PetscFunctionList GAMGList = 0;
 31: static PetscBool PCGAMGPackageInitialized;

 33: /* ----------------------------------------------------------------------------- */
 36: PetscErrorCode PCReset_GAMG(PC pc)
 37: {
 39:   PC_MG          *mg      = (PC_MG*)pc->data;
 40:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

 43:   if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */
 44:     PetscPrintf(PetscObjectComm((PetscObject)pc),"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__);
 45:     PetscFree(pc_gamg->data);
 46:   }
 47:   pc_gamg->data_sz = 0;
 48:   PetscFree(pc_gamg->orig_data);
 49:   return(0);
 50: }

 52: /* -------------------------------------------------------------------------- */
 53: /*
 54:    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
 55:      of active processors.

 57:    Input Parameter:
 58:    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
 59:           'pc_gamg->data_sz' are changed via repartitioning/reduction.
 60:    . Amat_fine - matrix on this fine (k) level
 61:    . cr_bs - coarse block size
 62:    In/Output Parameter:
 63:    . a_P_inout - prolongation operator to the next level (k-->k-1)
 64:    . a_nactive_proc - number of active procs
 65:    Output Parameter:
 66:    . a_Amat_crs - coarse matrix that is created (k-1)
 67: */

 71: static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,
 72:                                   Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,
 73:                                   IS * Pcolumnperm)
 74: {
 75:   PetscErrorCode  ierr;
 76:   PC_MG           *mg         = (PC_MG*)pc->data;
 77:   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
 78:   Mat             Cmat,Pold=*a_P_inout;
 79:   MPI_Comm        comm;
 80:   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
 81:   PetscInt        ncrs_eq,ncrs,f_bs;

 84:   PetscObjectGetComm((PetscObject)Amat_fine,&comm);
 85:   MPI_Comm_rank(comm, &rank);
 86:   MPI_Comm_size(comm, &size);
 87:   MatGetBlockSize(Amat_fine, &f_bs);
 88:   MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);

 90:   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
 91:   MatGetLocalSize(Cmat, &ncrs_eq, NULL);
 92:   if (pc_gamg->data_cell_rows>0) {
 93:     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
 94:   } else {
 95:     PetscInt  bs;
 96:     MatGetBlockSize(Cmat, &bs);
 97:     ncrs = ncrs_eq/bs;
 98:   }

100:   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
101:   {
102:     PetscInt ncrs_eq_glob;
103:     MatGetSize(Cmat, &ncrs_eq_glob, NULL);
104:     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
105:     if (!new_size) new_size = 1; /* not likely, posible? */
106:     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
107:   }

109:   if (Pcolumnperm) *Pcolumnperm = NULL;

111:   if (!pc_gamg->repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
112:   else {
113:     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
114:     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;

116:     nloc_old = ncrs_eq/cr_bs;
117:     if (ncrs_eq % cr_bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs);
118: #if defined PETSC_GAMG_USE_LOG
119:     PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);
120: #endif
121:     /* make 'is_eq_newproc' */
122:     PetscMalloc1(size, &counts);
123:     if (pc_gamg->repart) {
124:       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
125:       Mat adj;

127:      PetscInfo3(pc,"Repartition: size (active): %D --> %D, neq = %D\n",*a_nactive_proc,new_size,ncrs_eq);

129:       /* get 'adj' */
130:       if (cr_bs == 1) {
131:         MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);
132:       } else {
133:         /* make a scalar matrix to partition (no Stokes here) */
134:         Mat               tMat;
135:         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
136:         const PetscScalar *vals;
137:         const PetscInt    *idx;
138:         PetscInt          *d_nnz, *o_nnz, M, N;
139:         static PetscInt   llev = 0;
140:         MatType           mtype;

142:         PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);
143:         MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);
144:         MatGetSize(Cmat, &M, &N);
145:         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
146:           MatGetRow(Cmat,Ii,&ncols,0,0);
147:           d_nnz[jj] = ncols/cr_bs;
148:           o_nnz[jj] = ncols/cr_bs;
149:           MatRestoreRow(Cmat,Ii,&ncols,0,0);
150:           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
151:           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
152:         }

154:         MatGetType(Amat_fine,&mtype);
155:         MatCreate(comm, &tMat);
156:         MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);
157:         MatSetType(tMat,mtype);
158:         MatSeqAIJSetPreallocation(tMat,0,d_nnz);
159:         MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);
160:         PetscFree2(d_nnz,o_nnz);

162:         for (ii = Istart_crs; ii < Iend_crs; ii++) {
163:           PetscInt dest_row = ii/cr_bs;
164:           MatGetRow(Cmat,ii,&ncols,&idx,&vals);
165:           for (jj = 0; jj < ncols; jj++) {
166:             PetscInt    dest_col = idx[jj]/cr_bs;
167:             PetscScalar v        = 1.0;
168:             MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);
169:           }
170:           MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);
171:         }
172:         MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);
173:         MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);

175:         if (llev++ == -1) {
176:           PetscViewer viewer; char fname[32];
177:           PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);
178:           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
179:           MatView(tMat, viewer);
180:           PetscViewerDestroy(&viewer);
181:         }

183:         MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);

185:         MatDestroy(&tMat);
186:       } /* create 'adj' */

188:       { /* partition: get newproc_idx */
189:         char            prefix[256];
190:         const char      *pcpre;
191:         const PetscInt  *is_idx;
192:         MatPartitioning mpart;
193:         IS              proc_is;
194:         PetscInt        targetPE;

196:         MatPartitioningCreate(comm, &mpart);
197:         MatPartitioningSetAdjacency(mpart, adj);
198:         PCGetOptionsPrefix(pc, &pcpre);
199:         PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
200:         PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);
201:         MatPartitioningSetFromOptions(mpart);
202:         MatPartitioningSetNParts(mpart, new_size);
203:         MatPartitioningApply(mpart, &proc_is);
204:         MatPartitioningDestroy(&mpart);

206:         /* collect IS info */
207:         PetscMalloc1(ncrs_eq, &newproc_idx);
208:         ISGetIndices(proc_is, &is_idx);
209:         targetPE = 1; /* bring to "front" of machine */
210:         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
211:         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
212:           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
213:             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
214:           }
215:         }
216:         ISRestoreIndices(proc_is, &is_idx);
217:         ISDestroy(&proc_is);
218:       }
219:       MatDestroy(&adj);

221:       ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);
222:       PetscFree(newproc_idx);
223:     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */

225:       PetscInt rfactor,targetPE;
226:       /* find factor */
227:       if (new_size == 1) rfactor = size; /* easy */
228:       else {
229:         PetscReal best_fact = 0.;
230:         jj = -1;
231:         for (kk = 1 ; kk <= size ; kk++) {
232:           if (!(size%kk)) { /* a candidate */
233:             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
234:             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
235:             if (fact > best_fact) {
236:               best_fact = fact; jj = kk;
237:             }
238:           }
239:         }
240:         if (jj != -1) rfactor = jj;
241:         else rfactor = 1; /* does this happen .. a prime */
242:       }
243:       new_size = size/rfactor;

245:       if (new_size==nactive) {
246:         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
247:         PetscFree(counts);
248:         PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);
249:         return(0);
250:       }

252:       PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);
253:       targetPE = rank/rfactor;
254:       ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);
255:     } /* end simple 'is_eq_newproc' */

257:     /*
258:      Create an index set from the is_eq_newproc index set to indicate the mapping TO
259:      */
260:     ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);
261:     is_eq_num_prim = is_eq_num;
262:     /*
263:       Determine how many equations/vertices are assigned to each processor
264:      */
265:     ISPartitioningCount(is_eq_newproc, size, counts);
266:     ncrs_eq_new = counts[rank];
267:     ISDestroy(&is_eq_newproc);
268:     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */

270:     PetscFree(counts);
271: #if defined PETSC_GAMG_USE_LOG
272:     PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);
273: #endif
274:     /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxilary data into some complex abstracted thing */
275:     {
276:     Vec            src_crd, dest_crd;
277:     const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
278:     VecScatter     vecscat;
279:     PetscScalar    *array;
280:     IS isscat;

282:     /* move data (for primal equations only) */
283:     /* Create a vector to contain the newly ordered element information */
284:     VecCreate(comm, &dest_crd);
285:     VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);
286:     VecSetType(dest_crd,VECSTANDARD); /* this is needed! */
287:     /*
288:      There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
289:      a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
290:      */
291:     PetscMalloc1(ncrs*node_data_sz, &tidx);
292:     ISGetIndices(is_eq_num_prim, &idx);
293:     for (ii=0,jj=0; ii<ncrs; ii++) {
294:       PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
295:       for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
296:     }
297:     ISRestoreIndices(is_eq_num_prim, &idx);
298:     ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);
299:     PetscFree(tidx);
300:     /*
301:      Create a vector to contain the original vertex information for each element
302:      */
303:     VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);
304:     for (jj=0; jj<ndata_cols; jj++) {
305:       const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
306:       for (ii=0; ii<ncrs; ii++) {
307:         for (kk=0; kk<ndata_rows; kk++) {
308:           PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
309:           PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
310:           VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);
311:         }
312:       }
313:     }
314:     VecAssemblyBegin(src_crd);
315:     VecAssemblyEnd(src_crd);
316:     /*
317:       Scatter the element vertex information (still in the original vertex ordering)
318:       to the correct processor
319:     */
320:     VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);
321:     ISDestroy(&isscat);
322:     VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);
323:     VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);
324:     VecScatterDestroy(&vecscat);
325:     VecDestroy(&src_crd);
326:     /*
327:       Put the element vertex data into a new allocation of the gdata->ele
328:     */
329:     PetscFree(pc_gamg->data);
330:     PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);

332:     pc_gamg->data_sz = node_data_sz*ncrs_new;
333:     strideNew        = ncrs_new*ndata_rows;

335:     VecGetArray(dest_crd, &array);
336:     for (jj=0; jj<ndata_cols; jj++) {
337:       for (ii=0; ii<ncrs_new; ii++) {
338:         for (kk=0; kk<ndata_rows; kk++) {
339:           PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
340:           pc_gamg->data[ix] = PetscRealPart(array[jx]);
341:         }
342:       }
343:     }
344:     VecRestoreArray(dest_crd, &array);
345:     VecDestroy(&dest_crd);
346:     }
347:     /* move A and P (columns) with new layout */
348: #if defined PETSC_GAMG_USE_LOG
349:     PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);
350: #endif

352:     /*
353:       Invert for MatGetSubMatrix
354:     */
355:     ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);
356:     ISSort(new_eq_indices); /* is this needed? */
357:     ISSetBlockSize(new_eq_indices, cr_bs);
358:     if (is_eq_num != is_eq_num_prim) {
359:       ISDestroy(&is_eq_num_prim); /* could be same as 'is_eq_num' */
360:     }
361:     if (Pcolumnperm) {
362:       PetscObjectReference((PetscObject)new_eq_indices);
363:       *Pcolumnperm = new_eq_indices;
364:     }
365:     ISDestroy(&is_eq_num);
366: #if defined PETSC_GAMG_USE_LOG
367:     PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);
368:     PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);
369: #endif
370:     /* 'a_Amat_crs' output */
371:     {
372:       Mat mat;
373:       MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);
374:       *a_Amat_crs = mat;

376:       if (!PETSC_TRUE) {
377:         PetscInt cbs, rbs;
378:         MatGetBlockSizes(Cmat, &rbs, &cbs);
379:         PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);
380:         MatGetBlockSizes(mat, &rbs, &cbs);
381:         PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",rank,__FUNCT__,rbs,cbs,cr_bs);
382:       }
383:     }
384:     MatDestroy(&Cmat);

386: #if defined PETSC_GAMG_USE_LOG
387:     PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);
388: #endif
389:     /* prolongator */
390:     {
391:       IS       findices;
392:       PetscInt Istart,Iend;
393:       Mat      Pnew;
394:       MatGetOwnershipRange(Pold, &Istart, &Iend);
395: #if defined PETSC_GAMG_USE_LOG
396:       PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);
397: #endif
398:       ISCreateStride(comm,Iend-Istart,Istart,1,&findices);
399:       ISSetBlockSize(findices,f_bs);
400:       MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);
401:       ISDestroy(&findices);

403:       if (!PETSC_TRUE) {
404:         PetscInt cbs, rbs;
405:         MatGetBlockSizes(Pold, &rbs, &cbs);
406:         PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);
407:         MatGetBlockSizes(Pnew, &rbs, &cbs);
408:         PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);
409:       }
410: #if defined PETSC_GAMG_USE_LOG
411:       PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);
412: #endif
413:       MatDestroy(a_P_inout);

415:       /* output - repartitioned */
416:       *a_P_inout = Pnew;
417:     }
418:     ISDestroy(&new_eq_indices);

420:     *a_nactive_proc = new_size; /* output */
421:   }

423:   /* outout matrix data */
424:   if (!PETSC_TRUE) {
425:     PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs;
426:     if (!llev) {
427:       sprintf(fname,"Cmat_%d.m",llev++);
428:       PetscViewerASCIIOpen(comm,fname,&viewer);
429:       PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
430:       MatView(Amat_fine, viewer);
431:       PetscViewerPopFormat(viewer);
432:       PetscViewerDestroy(&viewer);
433:     }
434:     sprintf(fname,"Cmat_%d.m",llev++);
435:     PetscViewerASCIIOpen(comm,fname,&viewer);
436:     PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
437:     MatView(Cmat, viewer);
438:     PetscViewerPopFormat(viewer);
439:     PetscViewerDestroy(&viewer);
440:   }
441:   return(0);
442: }

444: /* -------------------------------------------------------------------------- */
445: /*
446:    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
447:                     by setting data structures and options.

449:    Input Parameter:
450: .  pc - the preconditioner context

452: */
455: PetscErrorCode PCSetUp_GAMG(PC pc)
456: {
458:   PC_MG          *mg      = (PC_MG*)pc->data;
459:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
460:   Mat            Pmat     = pc->pmat;
461:   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
462:   MPI_Comm       comm;
463:   PetscMPIInt    rank,size,nactivepe;
464:   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
465:   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
466:   PetscLogDouble nnz0=0.,nnztot=0.;
467:   MatInfo        info;

470:   PetscObjectGetComm((PetscObject)pc,&comm);
471:   MPI_Comm_rank(comm,&rank);
472:   MPI_Comm_size(comm,&size);

474:   if (pc_gamg->setup_count++ > 0) {
475:     if ((PetscBool)(!pc_gamg->reuse_prol)) {
476:       /* reset everything */
477:       PCReset_MG(pc);
478:       pc->setupcalled = 0;
479:     } else {
480:       PC_MG_Levels **mglevels = mg->levels;
481:       /* just do Galerkin grids */
482:       Mat          B,dA,dB;

484:      if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
485:       if (pc_gamg->Nlevels > 1) {
486:         /* currently only handle case where mat and pmat are the same on coarser levels */
487:         KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);
488:         /* (re)set to get dirty flag */
489:         KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);

491:         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
492:           /* the first time through the matrix structure has changed from repartitioning */
493:           if (pc_gamg->setup_count==2) {
494:             MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);
495:             MatDestroy(&mglevels[level]->A);

497:             mglevels[level]->A = B;
498:           } else {
499:             KSPGetOperators(mglevels[level]->smoothd,NULL,&B);
500:             MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);
501:           }
502:           KSPSetOperators(mglevels[level]->smoothd,B,B);
503:           dB   = B;
504:         }
505:       }

507:       PCSetUp_MG(pc);

509:       return(0);
510:     }
511:   }

513:   if (!pc_gamg->data) {
514:     if (pc_gamg->orig_data) {
515:       MatGetBlockSize(Pmat, &bs);
516:       MatGetLocalSize(Pmat, &qq, NULL);

518:       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
519:       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
520:       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;

522:       PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);
523:       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
524:     } else {
525:       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
526:       pc_gamg->ops->createdefaultdata(pc,Pmat);
527:     }
528:   }

530:   /* cache original data for reuse */
531:   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
532:     PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);
533:     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
534:     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
535:     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
536:   }

538:   /* get basic dims */
539:   MatGetBlockSize(Pmat, &bs);
540:   MatGetSize(Pmat, &M, &qq);

542:   MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info); /* global reduction */
543:   nnz0   = info.nz_used;
544:   nnztot = info.nz_used;
545:   PetscInfo6(pc,"level %d) N=%D, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",
546:                     0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,
547:                     (int)(nnz0/(PetscReal)M+0.5),size);
548: 

550:   /* Get A_i and R_i */
551:   for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */
552:        level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit);
553:        level++) {
554:     pc_gamg->current_level = level;
555:     level1 = level + 1;
556: #if defined PETSC_GAMG_USE_LOG
557:     PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);
558: #if (defined GAMG_STAGES)
559:     PetscLogStagePush(gamg_stages[level]);
560: #endif
561: #endif
562:     { /* construct prolongator */
563:       Mat              Gmat;
564:       PetscCoarsenData *agg_lists;
565:       Mat              Prol11;

567:       pc_gamg->ops->graph(pc,Aarr[level], &Gmat);
568:       pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);
569:       pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);

571:       /* could have failed to create new level */
572:       if (Prol11) {
573:         /* get new block size of coarse matrices */
574:         MatGetBlockSizes(Prol11, NULL, &bs);

576:         if (pc_gamg->ops->optprolongator) {
577:           /* smooth */
578:           pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);
579:         }

581:         Parr[level1] = Prol11;
582:       } else Parr[level1] = NULL;

584:       if (pc_gamg->use_aggs_in_gasm) {
585:         PetscInt bs;
586:         MatGetBlockSizes(Prol11, &bs, NULL);
587:         PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);
588:       }

590:       MatDestroy(&Gmat);
591:       PetscCDDestroy(agg_lists);
592:     } /* construct prolongator scope */
593: #if defined PETSC_GAMG_USE_LOG
594:     PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);
595: #endif
596:     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
597:     if (!Parr[level1]) {
598:        PetscInfo1(pc,"Stop gridding, level %D\n",level);
599: #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
600:       PetscLogStagePop();
601: #endif
602:       break;
603:     }
604: #if defined PETSC_GAMG_USE_LOG
605:     PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);
606: #endif

608:     pc_gamg->ops->createlevel(pc, Aarr[level], bs,&Parr[level1], &Aarr[level1], &nactivepe, NULL);

610: #if defined PETSC_GAMG_USE_LOG
611:     PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);
612: #endif
613:     MatGetSize(Aarr[level1], &M, &qq);

615:     MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);
616:     nnztot += info.nz_used;
617:     PetscInfo5(pc,"%d) N=%D, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",level1,M,pc_gamg->data_cell_cols,(int)(info.nz_used/(PetscReal)M),nactivepe);

619: #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
620:     PetscLogStagePop();
621: #endif
622:     /* stop if one node or one proc -- could pull back for singular problems */
623:     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) {
624:        PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);
625:       level++;
626:       break;
627:     }
628:   } /* levels */
629:   PetscFree(pc_gamg->data);

631:   PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);
632:   pc_gamg->Nlevels = level + 1;
633:   fine_level       = level;
634:   PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);

636:   /* simple setup */
637:   if (!PETSC_TRUE) {
638:     PC_MG_Levels **mglevels = mg->levels;
639:     for (lidx=0,level=pc_gamg->Nlevels-1;
640:          lidx<fine_level;
641:          lidx++, level--) {
642:       PCMGSetInterpolation(pc, lidx+1, Parr[level]);
643:       KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level]);
644:       MatDestroy(&Parr[level]);
645:       MatDestroy(&Aarr[level]);
646:     }
647:     KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0]);

649:     PCSetUp_MG(pc);
650:   } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
651:     /* set default smoothers & set operators */
652:     for (lidx = 1, level = pc_gamg->Nlevels-2;
653:          lidx <= fine_level;
654:          lidx++, level--) {
655:       KSP smoother;
656:       PC  subpc;

658:       PCMGGetSmoother(pc, lidx, &smoother);
659:       KSPGetPC(smoother, &subpc);

661:       KSPSetNormType(smoother, KSP_NORM_NONE);
662:       /* set ops */
663:       KSPSetOperators(smoother, Aarr[level], Aarr[level]);
664:       PCMGSetInterpolation(pc, lidx, Parr[level+1]);

666:       /* set defaults */
667:       KSPSetType(smoother, KSPCHEBYSHEV);

669:       /* set blocks for GASM smoother that uses the 'aggregates' */
670:       if (pc_gamg->use_aggs_in_gasm) {
671:         PetscInt sz;
672:         IS       *is;

674:         sz   = nASMBlocksArr[level];
675:         is   = ASMLocalIDsArr[level];
676:         PCSetType(subpc, PCGASM);
677:         PCGASMSetOverlap(subpc, 0);
678:         if (!sz) {
679:           IS       is;
680:           PetscInt my0,kk;
681:           MatGetOwnershipRange(Aarr[level], &my0, &kk);
682:           ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);
683:           PCGASMSetSubdomains(subpc, 1, &is, NULL);
684:           ISDestroy(&is);
685:         } else {
686:           PetscInt kk;
687:           PCGASMSetSubdomains(subpc, sz, is, NULL);
688:           for (kk=0; kk<sz; kk++) {
689:             ISDestroy(&is[kk]);
690:           }
691:           PetscFree(is);
692:         }
693:         ASMLocalIDsArr[level] = NULL;
694:         nASMBlocksArr[level]  = 0;
695:         PCGASMSetType(subpc, PC_GASM_BASIC);
696:       } else {
697:         PCSetType(subpc, PCSOR);
698:       }
699:     }
700:     {
701:       /* coarse grid */
702:       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
703:       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
704:       PCMGGetSmoother(pc, lidx, &smoother);
705:       KSPSetOperators(smoother, Lmat, Lmat);
706:       KSPSetNormType(smoother, KSP_NORM_NONE);
707:       KSPGetPC(smoother, &subpc);
708:       PCSetType(subpc, PCBJACOBI);
709:       PCSetUp(subpc);
710:       PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);
711:       if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
712:       KSPGetPC(k2[0],&pc2);
713:       PCSetType(pc2, PCLU);
714:       PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);
715:       KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);
716:       KSPSetType(k2[0], KSPPREONLY);
717:       /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
718:        * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
719:        * view every subdomain as though they were different. */
720:       ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
721:     }

723:     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
724:     PetscObjectOptionsBegin((PetscObject)pc);
725:     PCSetFromOptions_MG(PetscOptionsObject,pc);
726:     PetscOptionsEnd();
727:     if (!mg->galerkin) SETERRQ(comm,PETSC_ERR_USER,"PCGAMG must use Galerkin for coarse operators.");
728:     if (mg->galerkin == 1) mg->galerkin = 2;

730:     /* clean up */
731:     for (level=1; level<pc_gamg->Nlevels; level++) {
732:       MatDestroy(&Parr[level]);
733:       MatDestroy(&Aarr[level]);
734:     }

736:     PCSetUp_MG(pc);
737:   } else {
738:     KSP smoother;
739:     PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");
740:     PCMGGetSmoother(pc, 0, &smoother);
741:     KSPSetOperators(smoother, Aarr[0], Aarr[0]);
742:     KSPSetType(smoother, KSPPREONLY);
743:     PCSetUp_MG(pc);
744:   }
745:   return(0);
746: }

748: /* ------------------------------------------------------------------------- */
749: /*
750:  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
751:    that was created with PCCreate_GAMG().

753:    Input Parameter:
754: .  pc - the preconditioner context

756:    Application Interface Routine: PCDestroy()
757: */
760: PetscErrorCode PCDestroy_GAMG(PC pc)
761: {
763:   PC_MG          *mg     = (PC_MG*)pc->data;
764:   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;

767:   PCReset_GAMG(pc);
768:   if (pc_gamg->ops->destroy) {
769:     (*pc_gamg->ops->destroy)(pc);
770:   }
771:   PetscRandomDestroy(&pc_gamg->random);
772:   PetscFree(pc_gamg->ops);
773:   PetscFree(pc_gamg->gamg_type_name);
774:   PetscFree(pc_gamg);
775:   PCDestroy_MG(pc);
776:   return(0);
777: }


782: /*@
783:    PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via processor reduction.

785:    Logically Collective on PC

787:    Input Parameters:
788: +  pc - the preconditioner context
789: -  n - the number of equations


792:    Options Database Key:
793: .  -pc_gamg_process_eq_limit <limit>

795:    Level: intermediate

797:    Concepts: Unstructured multigrid preconditioner

799: .seealso: PCGAMGSetCoarseEqLim()
800: @*/
801: PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
802: {

807:   PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));
808:   return(0);
809: }

813: static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
814: {
815:   PC_MG   *mg      = (PC_MG*)pc->data;
816:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

819:   if (n>0) pc_gamg->min_eq_proc = n;
820:   return(0);
821: }

825: /*@
826:    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.

828:  Collective on PC

830:    Input Parameters:
831: +  pc - the preconditioner context
832: -  n - maximum number of equations to aim for

834:    Options Database Key:
835: .  -pc_gamg_coarse_eq_limit <limit>

837:    Level: intermediate

839:    Concepts: Unstructured multigrid preconditioner

841: .seealso: PCGAMGSetProcEqLim()
842: @*/
843: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
844: {

849:   PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
850:   return(0);
851: }

855: static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
856: {
857:   PC_MG   *mg      = (PC_MG*)pc->data;
858:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

861:   if (n>0) pc_gamg->coarse_eq_limit = n;
862:   return(0);
863: }

867: /*@
868:    PCGAMGSetRepartitioning - Repartition the coarse grids

870:    Collective on PC

872:    Input Parameters:
873: +  pc - the preconditioner context
874: -  n - PETSC_TRUE or PETSC_FALSE

876:    Options Database Key:
877: .  -pc_gamg_repartition <true,false>

879:    Level: intermediate

881:    Concepts: Unstructured multigrid preconditioner

883: .seealso: ()
884: @*/
885: PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
886: {

891:   PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));
892:   return(0);
893: }

897: static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
898: {
899:   PC_MG   *mg      = (PC_MG*)pc->data;
900:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

903:   pc_gamg->repart = n;
904:   return(0);
905: }

909: /*@
910:    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner

912:    Collective on PC

914:    Input Parameters:
915: +  pc - the preconditioner context
916: -  n - PETSC_TRUE or PETSC_FALSE

918:    Options Database Key:
919: .  -pc_gamg_reuse_interpolation <true,false>

921:    Level: intermediate

923:    Concepts: Unstructured multigrid preconditioner

925: .seealso: ()
926: @*/
927: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
928: {

933:   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
934:   return(0);
935: }

939: static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
940: {
941:   PC_MG   *mg      = (PC_MG*)pc->data;
942:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

945:   pc_gamg->reuse_prol = n;
946:   return(0);
947: }

951: /*@
952:    PCGAMGSetUseASMAggs -

954:    Collective on PC

956:    Input Parameters:
957: .  pc - the preconditioner context


960:    Options Database Key:
961: .  -pc_gamg_use_agg_gasm

963:    Level: intermediate

965:    Concepts: Unstructured multigrid preconditioner

967: .seealso: ()
968: @*/
969: PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
970: {

975:   PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));
976:   return(0);
977: }

981: static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
982: {
983:   PC_MG   *mg      = (PC_MG*)pc->data;
984:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

987:   pc_gamg->use_aggs_in_gasm = n;
988:   return(0);
989: }

993: /*@
994:    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use

996:    Not collective on PC

998:    Input Parameters:
999: +  pc - the preconditioner
1000: -  n - the maximum number of levels to use

1002:    Options Database Key:
1003: .  -pc_mg_levels

1005:    Level: intermediate

1007:    Concepts: Unstructured multigrid preconditioner

1009: .seealso: ()
1010: @*/
1011: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1012: {

1017:   PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
1018:   return(0);
1019: }

1023: static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1024: {
1025:   PC_MG   *mg      = (PC_MG*)pc->data;
1026:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1029:   pc_gamg->Nlevels = n;
1030:   return(0);
1031: }

1035: /*@
1036:    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph

1038:    Not collective on PC

1040:    Input Parameters:
1041: +  pc - the preconditioner context
1042: -  threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph

1044:    Options Database Key:
1045: .  -pc_gamg_threshold <threshold>

1047:    Level: intermediate

1049:    Concepts: Unstructured multigrid preconditioner

1051: .seealso: ()
1052: @*/
1053: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
1054: {

1059:   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));
1060:   return(0);
1061: }

1065: static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
1066: {
1067:   PC_MG   *mg      = (PC_MG*)pc->data;
1068:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1071:   pc_gamg->threshold = n;
1072:   return(0);
1073: }

1077: /*@
1078:    PCGAMGSetType - Set solution method

1080:    Collective on PC

1082:    Input Parameters:
1083: +  pc - the preconditioner context
1084: -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL

1086:    Options Database Key:
1087: .  -pc_gamg_type <agg,geo,classical>

1089:    Level: intermediate

1091:    Concepts: Unstructured multigrid preconditioner

1093: .seealso: PCGAMGGetType(), PCGAMG
1094: @*/
1095: PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1096: {

1101:   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1102:   return(0);
1103: }

1107: /*@
1108:    PCGAMGGetType - Get solution method

1110:    Collective on PC

1112:    Input Parameter:
1113: .  pc - the preconditioner context

1115:    Output Parameter:
1116: .  type - the type of algorithm used

1118:    Level: intermediate

1120:    Concepts: Unstructured multigrid preconditioner

1122: .seealso: PCGAMGSetType(), PCGAMGType
1123: @*/
1124: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1125: {

1130:   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1131:   return(0);
1132: }

1136: static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1137: {
1138:   PC_MG          *mg      = (PC_MG*)pc->data;
1139:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1142:   *type = pc_gamg->type;
1143:   return(0);
1144: }

1148: static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1149: {
1150:   PetscErrorCode ierr,(*r)(PC);
1151:   PC_MG          *mg      = (PC_MG*)pc->data;
1152:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1155:   pc_gamg->type = type;
1156:   PetscFunctionListFind(GAMGList,type,&r);
1157:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1158:   if (pc_gamg->ops->destroy) {
1159:     (*pc_gamg->ops->destroy)(pc);
1160:     PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));
1161:     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1162:     /* cleaning up common data in pc_gamg - this should disapear someday */
1163:     pc_gamg->data_cell_cols = 0;
1164:     pc_gamg->data_cell_rows = 0;
1165:     pc_gamg->orig_data_cell_cols = 0;
1166:     pc_gamg->orig_data_cell_rows = 0;
1167:     PetscFree(pc_gamg->data);
1168:     pc_gamg->data_sz = 0;
1169:   }
1170:   PetscFree(pc_gamg->gamg_type_name);
1171:   PetscStrallocpy(type,&pc_gamg->gamg_type_name);
1172:   (*r)(pc);
1173:   return(0);
1174: }

1178: static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
1179: {
1181:   PC_MG          *mg      = (PC_MG*)pc->data;
1182:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1185:   PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");
1186:   PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values from graph %g\n",(double)pc_gamg->threshold);
1187:   if (pc_gamg->ops->view) {
1188:     (*pc_gamg->ops->view)(pc,viewer);
1189:   }
1190:   return(0);
1191: }

1195: PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
1196: {
1198:   PC_MG          *mg      = (PC_MG*)pc->data;
1199:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1200:   PetscBool      flag;
1201:   MPI_Comm       comm;
1202:   char           prefix[256];
1203:   const char     *pcpre;

1206:   PetscObjectGetComm((PetscObject)pc,&comm);
1207:   PetscOptionsHead(PetscOptionsObject,"GAMG options");
1208:   {
1209:     char tname[256];
1210:     PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);
1211:     if (flag) {
1212:       PCGAMGSetType(pc,tname);
1213:     }
1214:     PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);
1215:     PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);
1216:     PetscOptionsBool("-pc_gamg_use_agg_gasm","Use aggregation agragates for GASM smoother","PCGAMGUseASMAggs",pc_gamg->use_aggs_in_gasm,&pc_gamg->use_aggs_in_gasm,NULL);
1217:     PetscOptionsInt("-pc_gamg_process_eq_limit","Limit (goal) on number of equations per process on coarse grids","PCGAMGSetProcEqLim",pc_gamg->min_eq_proc,&pc_gamg->min_eq_proc,NULL);
1218:     PetscOptionsInt("-pc_gamg_coarse_eq_limit","Limit on number of equations for the coarse grid","PCGAMGSetCoarseEqLim",pc_gamg->coarse_eq_limit,&pc_gamg->coarse_eq_limit,NULL);
1219:     PetscOptionsReal("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&pc_gamg->threshold,&flag);
1220:     PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);

1222:     /* set options for subtype */
1223:     if (pc_gamg->ops->setfromoptions) {(*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);}
1224:   }
1225:   PCGetOptionsPrefix(pc, &pcpre);
1226:   PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
1227:   PetscObjectSetOptionsPrefix((PetscObject)pc_gamg->random,prefix);
1228:   PetscRandomSetFromOptions(pc_gamg->random);
1229:   PetscOptionsTail();
1230:   return(0);
1231: }

1233: /* -------------------------------------------------------------------------- */
1234: /*MC
1235:      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner

1237:    Options Database Keys:
1238:    Multigrid options(inherited)
1239: +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1240: .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1241: .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1242: -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade


1245:   Notes: In order to obtain good performance for PCGAMG for vector valued problems you must
1246: $       Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1247: $       Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1248: $       See the Users Manual Chapter 4 for more details

1250:   Level: intermediate

1252:   Concepts: algebraic multigrid

1254: .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1255:            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartitioning(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseASMAggs(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType()
1256: M*/

1260: PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1261: {
1263:   PC_GAMG        *pc_gamg;
1264:   PC_MG          *mg;

1267:   /* register AMG type */
1268:   PCGAMGInitializePackage();

1270:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1271:   PCSetType(pc, PCMG);
1272:   PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);

1274:   /* create a supporting struct and attach it to pc */
1275:   PetscNewLog(pc,&pc_gamg);
1276:   mg           = (PC_MG*)pc->data;
1277:   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally from PCMG by GAMG code */
1278:   mg->innerctx = pc_gamg;

1280:   PetscNewLog(pc,&pc_gamg->ops);

1282:   pc_gamg->setup_count = 0;
1283:   /* these should be in subctx but repartitioning needs simple arrays */
1284:   pc_gamg->data_sz = 0;
1285:   pc_gamg->data    = 0;

1287:   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1288:   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1289:   pc->ops->setup          = PCSetUp_GAMG;
1290:   pc->ops->reset          = PCReset_GAMG;
1291:   pc->ops->destroy        = PCDestroy_GAMG;
1292:   mg->view                = PCView_GAMG;

1294:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);
1295:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);
1296:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);
1297:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);
1298:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);
1299:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);
1300:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);
1301:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);
1302:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);
1303:   pc_gamg->repart           = PETSC_FALSE;
1304:   pc_gamg->reuse_prol       = PETSC_FALSE;
1305:   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1306:   pc_gamg->min_eq_proc      = 50;
1307:   pc_gamg->coarse_eq_limit  = 50;
1308:   pc_gamg->threshold        = 0.;
1309:   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
1310:   pc_gamg->current_level    = 0; /* don't need to init really */
1311:   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;

1313:   PetscRandomCreate(PetscObjectComm((PetscObject)pc),&pc_gamg->random);

1315:   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1316:   PCGAMGSetType(pc,PCGAMGAGG);
1317:   return(0);
1318: }

1322: /*@C
1323:  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
1324:     from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG()
1325:     when using static libraries.

1327:  Level: developer

1329:  .keywords: PC, PCGAMG, initialize, package
1330:  .seealso: PetscInitialize()
1331: @*/
1332: PetscErrorCode PCGAMGInitializePackage(void)
1333: {

1337:   if (PCGAMGPackageInitialized) return(0);
1338:   PCGAMGPackageInitialized = PETSC_TRUE;
1339:   PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);
1340:   PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);
1341:   PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);
1342:   PetscRegisterFinalize(PCGAMGFinalizePackage);

1344:   /* general events */
1345:   PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);
1346:   PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);
1347:   PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1348:   PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1349:   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1350:   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1351:   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);

1353: #if defined PETSC_GAMG_USE_LOG
1354:   PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1355:   PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1356:   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1357:   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1358:   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1359:   PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1360:   PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1361:   PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1362:   PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1363:   PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1364:   PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1365:   PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1366:   PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1367:   PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1368:   PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1369:   PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1370:   PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);

1372:   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1373:   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1374:   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1375:   /* create timer stages */
1376: #if defined GAMG_STAGES
1377:   {
1378:     char     str[32];
1379:     PetscInt lidx;
1380:     sprintf(str,"MG Level %d (finest)",0);
1381:     PetscLogStageRegister(str, &gamg_stages[0]);
1382:     for (lidx=1; lidx<9; lidx++) {
1383:       sprintf(str,"MG Level %d",lidx);
1384:       PetscLogStageRegister(str, &gamg_stages[lidx]);
1385:     }
1386:   }
1387: #endif
1388: #endif
1389:   return(0);
1390: }

1394: /*@C
1395:  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1396:     called from PetscFinalize() automatically.

1398:  Level: developer

1400:  .keywords: Petsc, destroy, package
1401:  .seealso: PetscFinalize()
1402: @*/
1403: PetscErrorCode PCGAMGFinalizePackage(void)
1404: {

1408:   PCGAMGPackageInitialized = PETSC_FALSE;
1409:   PetscFunctionListDestroy(&GAMGList);
1410:   return(0);
1411: }

1415: /*@C
1416:  PCGAMGRegister - Register a PCGAMG implementation.

1418:  Input Parameters:
1419:  + type - string that will be used as the name of the GAMG type.
1420:  - create - function for creating the gamg context.

1422:   Level: advanced

1424:  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1425: @*/
1426: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1427: {

1431:   PCGAMGInitializePackage();
1432:   PetscFunctionListAdd(&GAMGList,type,create);
1433:   return(0);
1434: }