Actual source code: gamg.c

petsc-3.7.5 2017-01-01
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: #if defined PETSC_GAMG_USE_LOG
250:         PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);
251: #endif
252:         return(0);
253:       }

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

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

273:     PetscFree(counts);
274: #if defined PETSC_GAMG_USE_LOG
275:     PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);
276: #endif
277:     /* 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 */
278:     {
279:     Vec            src_crd, dest_crd;
280:     const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
281:     VecScatter     vecscat;
282:     PetscScalar    *array;
283:     IS isscat;

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

335:     pc_gamg->data_sz = node_data_sz*ncrs_new;
336:     strideNew        = ncrs_new*ndata_rows;

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

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

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

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

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

418:       /* output - repartitioned */
419:       *a_P_inout = Pnew;
420:     }
421:     ISDestroy(&new_eq_indices);

423:     *a_nactive_proc = new_size; /* output */
424:   }

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

447: /* -------------------------------------------------------------------------- */
448: /*
449:    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
450:                     by setting data structures and options.

452:    Input Parameter:
453: .  pc - the preconditioner context

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

473:   PetscObjectGetComm((PetscObject)pc,&comm);
474:   MPI_Comm_rank(comm,&rank);
475:   MPI_Comm_size(comm,&size);

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

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

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

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

510:       PCSetUp_MG(pc);

512:       return(0);
513:     }
514:   }

516:   if (!pc_gamg->data) {
517:     if (pc_gamg->orig_data) {
518:       MatGetBlockSize(Pmat, &bs);
519:       MatGetLocalSize(Pmat, &qq, NULL);

521:       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
522:       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
523:       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;

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

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

541:   /* get basic dims */
542:   MatGetBlockSize(Pmat, &bs);
543:   MatGetSize(Pmat, &M, &qq);

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

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

570:       pc_gamg->ops->graph(pc,Aarr[level], &Gmat);
571:       pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);
572:       pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);

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

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

584:         Parr[level1] = Prol11;
585:       } else Parr[level1] = NULL;

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

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

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

613: #if defined PETSC_GAMG_USE_LOG
614:     PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);
615: #endif
616:     MatGetSize(Aarr[level1], &M, &qq);

618:     MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);
619:     nnztot += info.nz_used;
620:     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);

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

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

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

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

661:       PCMGGetSmoother(pc, lidx, &smoother);
662:       KSPGetPC(smoother, &subpc);

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

669:       /* set defaults */
670:       KSPSetType(smoother, KSPCHEBYSHEV);

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

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

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

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

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

751: /* ------------------------------------------------------------------------- */
752: /*
753:  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
754:    that was created with PCCreate_GAMG().

756:    Input Parameter:
757: .  pc - the preconditioner context

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

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


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

788:    Logically Collective on PC

790:    Input Parameters:
791: +  pc - the preconditioner context
792: -  n - the number of equations


795:    Options Database Key:
796: .  -pc_gamg_process_eq_limit <limit>

798:    Level: intermediate

800:    Concepts: Unstructured multigrid preconditioner

802: .seealso: PCGAMGSetCoarseEqLim()
803: @*/
804: PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
805: {

810:   PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));
811:   return(0);
812: }

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

822:   if (n>0) pc_gamg->min_eq_proc = n;
823:   return(0);
824: }

828: /*@
829:    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.

831:  Collective on PC

833:    Input Parameters:
834: +  pc - the preconditioner context
835: -  n - maximum number of equations to aim for

837:    Options Database Key:
838: .  -pc_gamg_coarse_eq_limit <limit>

840:    Level: intermediate

842:    Concepts: Unstructured multigrid preconditioner

844: .seealso: PCGAMGSetProcEqLim()
845: @*/
846: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
847: {

852:   PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
853:   return(0);
854: }

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

864:   if (n>0) pc_gamg->coarse_eq_limit = n;
865:   return(0);
866: }

870: /*@
871:    PCGAMGSetRepartitioning - Repartition the coarse grids

873:    Collective on PC

875:    Input Parameters:
876: +  pc - the preconditioner context
877: -  n - PETSC_TRUE or PETSC_FALSE

879:    Options Database Key:
880: .  -pc_gamg_repartition <true,false>

882:    Level: intermediate

884:    Concepts: Unstructured multigrid preconditioner

886: .seealso: ()
887: @*/
888: PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
889: {

894:   PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));
895:   return(0);
896: }

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

906:   pc_gamg->repart = n;
907:   return(0);
908: }

912: /*@
913:    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner

915:    Collective on PC

917:    Input Parameters:
918: +  pc - the preconditioner context
919: -  n - PETSC_TRUE or PETSC_FALSE

921:    Options Database Key:
922: .  -pc_gamg_reuse_interpolation <true,false>

924:    Level: intermediate

926:    Concepts: Unstructured multigrid preconditioner

928: .seealso: ()
929: @*/
930: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
931: {

936:   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
937:   return(0);
938: }

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

948:   pc_gamg->reuse_prol = n;
949:   return(0);
950: }

954: /*@
955:    PCGAMGSetUseASMAggs -

957:    Collective on PC

959:    Input Parameters:
960: .  pc - the preconditioner context


963:    Options Database Key:
964: .  -pc_gamg_use_agg_gasm

966:    Level: intermediate

968:    Concepts: Unstructured multigrid preconditioner

970: .seealso: ()
971: @*/
972: PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
973: {

978:   PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));
979:   return(0);
980: }

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

990:   pc_gamg->use_aggs_in_gasm = n;
991:   return(0);
992: }

996: /*@
997:    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use

999:    Not collective on PC

1001:    Input Parameters:
1002: +  pc - the preconditioner
1003: -  n - the maximum number of levels to use

1005:    Options Database Key:
1006: .  -pc_mg_levels

1008:    Level: intermediate

1010:    Concepts: Unstructured multigrid preconditioner

1012: .seealso: ()
1013: @*/
1014: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1015: {

1020:   PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
1021:   return(0);
1022: }

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

1032:   pc_gamg->Nlevels = n;
1033:   return(0);
1034: }

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

1041:    Not collective on PC

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

1047:    Options Database Key:
1048: .  -pc_gamg_threshold <threshold>

1050:    Level: intermediate

1052:    Concepts: Unstructured multigrid preconditioner

1054: .seealso: ()
1055: @*/
1056: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
1057: {

1062:   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));
1063:   return(0);
1064: }

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

1074:   pc_gamg->threshold = n;
1075:   return(0);
1076: }

1080: /*@C
1081:    PCGAMGSetType - Set solution method

1083:    Collective on PC

1085:    Input Parameters:
1086: +  pc - the preconditioner context
1087: -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL

1089:    Options Database Key:
1090: .  -pc_gamg_type <agg,geo,classical>

1092:    Level: intermediate

1094:    Concepts: Unstructured multigrid preconditioner

1096: .seealso: PCGAMGGetType(), PCGAMG
1097: @*/
1098: PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1099: {

1104:   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1105:   return(0);
1106: }

1110: /*@C
1111:    PCGAMGGetType - Get solution method

1113:    Collective on PC

1115:    Input Parameter:
1116: .  pc - the preconditioner context

1118:    Output Parameter:
1119: .  type - the type of algorithm used

1121:    Level: intermediate

1123:    Concepts: Unstructured multigrid preconditioner

1125: .seealso: PCGAMGSetType(), PCGAMGType
1126: @*/
1127: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1128: {

1133:   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1134:   return(0);
1135: }

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

1145:   *type = pc_gamg->type;
1146:   return(0);
1147: }

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

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

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

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

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

1209:   PetscObjectGetComm((PetscObject)pc,&comm);
1210:   PetscOptionsHead(PetscOptionsObject,"GAMG options");
1211:   {
1212:     char tname[256];
1213:     PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);
1214:     if (flag) {
1215:       PCGAMGSetType(pc,tname);
1216:     }
1217:     PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);
1218:     PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);
1219:     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);
1220:     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);
1221:     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);
1222:     PetscOptionsReal("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&pc_gamg->threshold,&flag);
1223:     PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);

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

1236: /* -------------------------------------------------------------------------- */
1237: /*MC
1238:      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner

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


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

1253:   Level: intermediate

1255:   Concepts: algebraic multigrid

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

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

1270:   /* register AMG type */
1271:   PCGAMGInitializePackage();

1273:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1274:   PCSetType(pc, PCMG);
1275:   PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);

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

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

1285:   pc_gamg->setup_count = 0;
1286:   /* these should be in subctx but repartitioning needs simple arrays */
1287:   pc_gamg->data_sz = 0;
1288:   pc_gamg->data    = 0;

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

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

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

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

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

1330:  Level: developer

1332:  .keywords: PC, PCGAMG, initialize, package
1333:  .seealso: PetscInitialize()
1334: @*/
1335: PetscErrorCode PCGAMGInitializePackage(void)
1336: {

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

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

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

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

1397: /*@C
1398:  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1399:     called from PetscFinalize() automatically.

1401:  Level: developer

1403:  .keywords: Petsc, destroy, package
1404:  .seealso: PetscFinalize()
1405: @*/
1406: PetscErrorCode PCGAMGFinalizePackage(void)
1407: {

1411:   PCGAMGPackageInitialized = PETSC_FALSE;
1412:   PetscFunctionListDestroy(&GAMGList);
1413:   return(0);
1414: }

1418: /*@C
1419:  PCGAMGRegister - Register a PCGAMG implementation.

1421:  Input Parameters:
1422:  + type - string that will be used as the name of the GAMG type.
1423:  - create - function for creating the gamg context.

1425:   Level: advanced

1427:  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1428: @*/
1429: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1430: {

1434:   PCGAMGInitializePackage();
1435:   PetscFunctionListAdd(&GAMGList,type,create);
1436:   return(0);
1437: }