Actual source code: gamg.c

petsc-master 2018-07-20
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>
  6: #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */

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

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

 22: /* #define GAMG_STAGES */
 23: #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
 24: static PetscLogStage gamg_stages[PETSC_GAMG_MAXLEVELS];
 25: #endif

 27: static PetscFunctionList GAMGList = 0;
 28: static PetscBool PCGAMGPackageInitialized;

 30: /* ----------------------------------------------------------------------------- */
 31: PetscErrorCode PCReset_GAMG(PC pc)
 32: {
 34:   PC_MG          *mg      = (PC_MG*)pc->data;
 35:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

 38:   if (pc_gamg->data) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen, cleaned up in SetUp\n");
 39:   pc_gamg->data_sz = 0;
 40:   PetscFree(pc_gamg->orig_data);
 41:   return(0);
 42: }

 44: /* -------------------------------------------------------------------------- */
 45: /*
 46:    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
 47:      of active processors.

 49:    Input Parameter:
 50:    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
 51:           'pc_gamg->data_sz' are changed via repartitioning/reduction.
 52:    . Amat_fine - matrix on this fine (k) level
 53:    . cr_bs - coarse block size
 54:    In/Output Parameter:
 55:    . a_P_inout - prolongation operator to the next level (k-->k-1)
 56:    . a_nactive_proc - number of active procs
 57:    Output Parameter:
 58:    . a_Amat_crs - coarse matrix that is created (k-1)
 59: */

 61: static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,IS * Pcolumnperm, PetscBool is_last)
 62: {
 63:   PetscErrorCode  ierr;
 64:   PC_MG           *mg         = (PC_MG*)pc->data;
 65:   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
 66:   Mat             Cmat,Pold=*a_P_inout;
 67:   MPI_Comm        comm;
 68:   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
 69:   PetscInt        ncrs_eq,ncrs,f_bs;

 72:   PetscObjectGetComm((PetscObject)Amat_fine,&comm);
 73:   MPI_Comm_rank(comm, &rank);
 74:   MPI_Comm_size(comm, &size);
 75:   MatGetBlockSize(Amat_fine, &f_bs);
 76:   MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);

 78:   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
 79:   MatGetLocalSize(Cmat, &ncrs_eq, NULL);
 80:   if (pc_gamg->data_cell_rows>0) {
 81:     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
 82:   } else {
 83:     PetscInt  bs;
 84:     MatGetBlockSize(Cmat, &bs);
 85:     ncrs = ncrs_eq/bs;
 86:   }

 88:   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
 89:   if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) new_size = 1;
 90:   else {
 91:     PetscInt ncrs_eq_glob;
 92:     MatGetSize(Cmat, &ncrs_eq_glob, NULL);
 93:     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
 94:     if (!new_size) new_size = 1; /* not likely, posible? */
 95:     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
 96:   }

 98:   if (Pcolumnperm) *Pcolumnperm = NULL;

100:   if (!pc_gamg->repart && new_size==nactive) {
101:     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
102:     /* we know that the grid structure can be reused in MatPtAP */
103:   } else {
104:     /* we know that the grid structure can NOT be reused in MatPtAP */
105:     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old;
106:     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;

108:     nloc_old = ncrs_eq/cr_bs;
109:     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);
110: #if defined PETSC_GAMG_USE_LOG
111:     PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);
112: #endif
113:     /* make 'is_eq_newproc' */
114:     PetscMalloc1(size, &counts);
115:     if (pc_gamg->repart) {
116:       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
117:       Mat adj;

119:      PetscInfo3(pc,"Repartition: size (active): %D --> %D, %D local equations\n",*a_nactive_proc,new_size,ncrs_eq);

121:       /* get 'adj' */
122:       if (cr_bs == 1) {
123:         MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);
124:       } else {
125:         /* make a scalar matrix to partition (no Stokes here) */
126:         Mat               tMat;
127:         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
128:         const PetscScalar *vals;
129:         const PetscInt    *idx;
130:         PetscInt          *d_nnz, *o_nnz, M, N;
131:         static PetscInt   llev = 0;
132:         MatType           mtype;

134:         PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);
135:         MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);
136:         MatGetSize(Cmat, &M, &N);
137:         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
138:           MatGetRow(Cmat,Ii,&ncols,0,0);
139:           d_nnz[jj] = ncols/cr_bs;
140:           o_nnz[jj] = ncols/cr_bs;
141:           MatRestoreRow(Cmat,Ii,&ncols,0,0);
142:           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
143:           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
144:         }

146:         MatGetType(Amat_fine,&mtype);
147:         MatCreate(comm, &tMat);
148:         MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);
149:         MatSetType(tMat,mtype);
150:         MatSeqAIJSetPreallocation(tMat,0,d_nnz);
151:         MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);
152:         PetscFree2(d_nnz,o_nnz);

154:         for (ii = Istart_crs; ii < Iend_crs; ii++) {
155:           PetscInt dest_row = ii/cr_bs;
156:           MatGetRow(Cmat,ii,&ncols,&idx,&vals);
157:           for (jj = 0; jj < ncols; jj++) {
158:             PetscInt    dest_col = idx[jj]/cr_bs;
159:             PetscScalar v        = 1.0;
160:             MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);
161:           }
162:           MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);
163:         }
164:         MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);
165:         MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);

167:         if (llev++ == -1) {
168:           PetscViewer viewer; char fname[32];
169:           PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);
170:           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
171:           MatView(tMat, viewer);
172:           PetscViewerDestroy(&viewer);
173:         }
174:         MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);
175:         MatDestroy(&tMat);
176:       } /* create 'adj' */

178:       { /* partition: get newproc_idx */
179:         char            prefix[256];
180:         const char      *pcpre;
181:         const PetscInt  *is_idx;
182:         MatPartitioning mpart;
183:         IS              proc_is;
184:         PetscInt        targetPE;

186:         MatPartitioningCreate(comm, &mpart);
187:         MatPartitioningSetAdjacency(mpart, adj);
188:         PCGetOptionsPrefix(pc, &pcpre);
189:         PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
190:         PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);
191:         MatPartitioningSetFromOptions(mpart);
192:         MatPartitioningSetNParts(mpart, new_size);
193:         MatPartitioningApply(mpart, &proc_is);
194:         MatPartitioningDestroy(&mpart);

196:         /* collect IS info */
197:         PetscMalloc1(ncrs_eq, &newproc_idx);
198:         ISGetIndices(proc_is, &is_idx);
199:         targetPE = 1; /* bring to "front" of machine */
200:         /*targetPE = size/new_size;*/ /* spread partitioning across machine */
201:         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
202:           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
203:             newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */
204:           }
205:         }
206:         ISRestoreIndices(proc_is, &is_idx);
207:         ISDestroy(&proc_is);
208:       }
209:       MatDestroy(&adj);

211:       ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);
212:       PetscFree(newproc_idx);
213:     } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */
214:       PetscInt rfactor,targetPE;

216:       /* find factor */
217:       if (new_size == 1) rfactor = size; /* easy */
218:       else {
219:         PetscReal best_fact = 0.;
220:         jj = -1;
221:         for (kk = 1 ; kk <= size ; kk++) {
222:           if (!(size%kk)) { /* a candidate */
223:             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
224:             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
225:             if (fact > best_fact) {
226:               best_fact = fact; jj = kk;
227:             }
228:           }
229:         }
230:         if (jj != -1) rfactor = jj;
231:         else rfactor = 1; /* does this happen .. a prime */
232:       }
233:       new_size = size/rfactor;

235:       if (new_size==nactive) {
236:         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
237:         PetscFree(counts);
238:         PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);
239: #if defined PETSC_GAMG_USE_LOG
240:         PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);
241: #endif
242:         return(0);
243:       }

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

250:     /*
251:      Create an index set from the is_eq_newproc index set to indicate the mapping TO
252:      */
253:     ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);
254:     is_eq_num_prim = is_eq_num;
255:     /*
256:       Determine how many equations/vertices are assigned to each processor
257:      */
258:     ISPartitioningCount(is_eq_newproc, size, counts);
259:     ncrs_eq_new = counts[rank];
260:     ISDestroy(&is_eq_newproc);
261:     ncrs_new = ncrs_eq_new/cr_bs; /* eqs */

263:     PetscFree(counts);
264: #if defined PETSC_GAMG_USE_LOG
265:     PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);
266: #endif
267:     /* 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 */
268:     {
269:     Vec            src_crd, dest_crd;
270:     const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
271:     VecScatter     vecscat;
272:     PetscScalar    *array;
273:     IS isscat;

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

325:     pc_gamg->data_sz = node_data_sz*ncrs_new;
326:     strideNew        = ncrs_new*ndata_rows;

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

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

371: #if defined PETSC_GAMG_USE_LOG
372:     PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);
373: #endif
374:     /* prolongator */
375:     {
376:       IS       findices;
377:       PetscInt Istart,Iend;
378:       Mat      Pnew;

380:       MatGetOwnershipRange(Pold, &Istart, &Iend);
381: #if defined PETSC_GAMG_USE_LOG
382:       PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);
383: #endif
384:       ISCreateStride(comm,Iend-Istart,Istart,1,&findices);
385:       ISSetBlockSize(findices,f_bs);
386:       MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);
387:       ISDestroy(&findices);

389: #if defined PETSC_GAMG_USE_LOG
390:       PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);
391: #endif
392:       MatDestroy(a_P_inout);

394:       /* output - repartitioned */
395:       *a_P_inout = Pnew;
396:     }
397:     ISDestroy(&new_eq_indices);

399:     *a_nactive_proc = new_size; /* output */
400:   }
401:   return(0);
402: }

404: /* -------------------------------------------------------------------------- */
405: /*
406:    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
407:                     by setting data structures and options.

409:    Input Parameter:
410: .  pc - the preconditioner context

412: */
413: PetscErrorCode PCSetUp_GAMG(PC pc)
414: {
416:   PC_MG          *mg      = (PC_MG*)pc->data;
417:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
418:   Mat            Pmat     = pc->pmat;
419:   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_GAMG_MAXLEVELS];
420:   MPI_Comm       comm;
421:   PetscMPIInt    rank,size,nactivepe;
422:   Mat            Aarr[PETSC_GAMG_MAXLEVELS],Parr[PETSC_GAMG_MAXLEVELS];
423:   IS             *ASMLocalIDsArr[PETSC_GAMG_MAXLEVELS];
424:   PetscLogDouble nnz0=0.,nnztot=0.;
425:   MatInfo        info;
426:   PetscBool      is_last = PETSC_FALSE;

429:   PetscObjectGetComm((PetscObject)pc,&comm);
430:   MPI_Comm_rank(comm,&rank);
431:   MPI_Comm_size(comm,&size);

433:   if (pc_gamg->setup_count++ > 0) {
434:     if ((PetscBool)(!pc_gamg->reuse_prol)) {
435:       /* reset everything */
436:       PCReset_MG(pc);
437:       pc->setupcalled = 0;
438:     } else {
439:       PC_MG_Levels **mglevels = mg->levels;
440:       /* just do Galerkin grids */
441:       Mat          B,dA,dB;

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

450:         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
451:           /* 2nd solve, matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
452:           if (pc_gamg->setup_count==2 /* && pc_gamg->repart||reduction */) {
453:             PetscInfo2(pc,"new RAP after first solve level %D, %D setup\n",level,pc_gamg->setup_count);
454:             MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,2.0,&B);
455:             MatDestroy(&mglevels[level]->A);
456:             mglevels[level]->A = B;
457:           } else {
458:             PetscInfo2(pc,"RAP after first solve reusing matrix level %D, %D setup\n",level,pc_gamg->setup_count);
459:             KSPGetOperators(mglevels[level]->smoothd,NULL,&B);
460:             MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);
461:           }
462:           KSPSetOperators(mglevels[level]->smoothd,B,B);
463:           dB   = B;
464:         }
465:       }

467:       PCSetUp_MG(pc);
468:       return(0);
469:     }
470:   }

472:   if (!pc_gamg->data) {
473:     if (pc_gamg->orig_data) {
474:       MatGetBlockSize(Pmat, &bs);
475:       MatGetLocalSize(Pmat, &qq, NULL);

477:       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
478:       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
479:       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;

481:       PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);
482:       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
483:     } else {
484:       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
485:       pc_gamg->ops->createdefaultdata(pc,Pmat);
486:     }
487:   }

489:   /* cache original data for reuse */
490:   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
491:     PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);
492:     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
493:     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
494:     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
495:   }

497:   /* get basic dims */
498:   MatGetBlockSize(Pmat, &bs);
499:   MatGetSize(Pmat, &M, &N);

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

506:   /* Get A_i and R_i */
507:   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
508:     pc_gamg->current_level = level;
509:     level1 = level + 1;
510: #if defined PETSC_GAMG_USE_LOG
511:     PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);
512: #if (defined GAMG_STAGES)
513:     PetscLogStagePush(gamg_stages[level]);
514: #endif
515: #endif
516:     { /* construct prolongator */
517:       Mat              Gmat;
518:       PetscCoarsenData *agg_lists;
519:       Mat              Prol11;

521:       pc_gamg->ops->graph(pc,Aarr[level], &Gmat);
522:       pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);
523:       pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);

525:       /* could have failed to create new level */
526:       if (Prol11) {
527:         /* get new block size of coarse matrices */
528:         MatGetBlockSizes(Prol11, NULL, &bs);

530:         if (pc_gamg->ops->optprolongator) {
531:           /* smooth */
532:           pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);
533:         }

535:         Parr[level1] = Prol11;
536:       } else Parr[level1] = NULL; /* failed to coarsen */

538:       if (pc_gamg->use_aggs_in_asm) {
539:         PetscInt bs;
540:         MatGetBlockSizes(Prol11, &bs, NULL);
541:         PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);
542:       }

544:       MatDestroy(&Gmat);
545:       PetscCDDestroy(agg_lists);
546:     } /* construct prolongator scope */
547: #if defined PETSC_GAMG_USE_LOG
548:     PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);
549: #endif
550:     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
551:     if (!Parr[level1]) { /* failed to coarsen */
552:        PetscInfo1(pc,"Stop gridding, level %D\n",level);
553: #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
554:       PetscLogStagePop();
555: #endif
556:       break;
557:     }
558: #if defined PETSC_GAMG_USE_LOG
559:     PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);
560: #endif
561:     MatGetSize(Parr[level1], &M, &N); /* N is next M, a loop test variables */
562:     if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ????????");
563:     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
564:     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
565:     pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);

567: #if defined PETSC_GAMG_USE_LOG
568:     PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);
569: #endif
570:     MatGetSize(Aarr[level1], &M, &N); /* M is loop test variables */
571:     MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);
572:     nnztot += info.nz_used;
573:     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);

575: #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
576:     PetscLogStagePop();
577: #endif
578:     /* stop if one node or one proc -- could pull back for singular problems */
579:     if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) {
580:        PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);
581:       level++;
582:       break;
583:     }
584:   } /* levels */
585:   PetscFree(pc_gamg->data);

587:   PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);
588:   pc_gamg->Nlevels = level + 1;
589:   fine_level       = level;
590:   PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);

592:   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
593:     /* set default smoothers & set operators */
594:     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
595:       KSP smoother;
596:       PC  subpc;

598:       PCMGGetSmoother(pc, lidx, &smoother);
599:       KSPGetPC(smoother, &subpc);

601:       KSPSetNormType(smoother, KSP_NORM_NONE);
602:       /* set ops */
603:       KSPSetOperators(smoother, Aarr[level], Aarr[level]);
604:       PCMGSetInterpolation(pc, lidx, Parr[level+1]);

606:       /* set defaults */
607:       KSPSetType(smoother, KSPCHEBYSHEV);

609:       /* set blocks for ASM smoother that uses the 'aggregates' */
610:       if (pc_gamg->use_aggs_in_asm) {
611:         PetscInt sz;
612:         IS       *iss;

614:         sz   = nASMBlocksArr[level];
615:         iss   = ASMLocalIDsArr[level];
616:         PCSetType(subpc, PCASM);
617:         PCASMSetOverlap(subpc, 0);
618:         PCASMSetType(subpc,PC_ASM_BASIC);
619:         if (!sz) {
620:           IS       is;
621:           ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);
622:           PCASMSetLocalSubdomains(subpc, 1, NULL, &is);
623:           ISDestroy(&is);
624:         } else {
625:           PetscInt kk;
626:           PCASMSetLocalSubdomains(subpc, sz, NULL, iss);
627:           for (kk=0; kk<sz; kk++) {
628:             ISDestroy(&iss[kk]);
629:           }
630:           PetscFree(iss);
631:         }
632:         ASMLocalIDsArr[level] = NULL;
633:         nASMBlocksArr[level]  = 0;
634:       } else {
635:         PCSetType(subpc, PCSOR);
636:       }
637:     }
638:     {
639:       /* coarse grid */
640:       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
641:       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
642:       PCMGGetSmoother(pc, lidx, &smoother);
643:       KSPSetOperators(smoother, Lmat, Lmat);
644:       if (!pc_gamg->use_parallel_coarse_grid_solver) {
645:         KSPSetNormType(smoother, KSP_NORM_NONE);
646:         KSPGetPC(smoother, &subpc);
647:         PCSetType(subpc, PCBJACOBI);
648:         PCSetUp(subpc);
649:         PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);
650:         if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
651:         KSPGetPC(k2[0],&pc2);
652:         PCSetType(pc2, PCLU);
653:         PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);
654:         KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);
655:         KSPSetType(k2[0], KSPPREONLY);
656:         /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
657:          * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
658:          * view every subdomain as though they were different. */
659:         ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
660:       }
661:     }

663:     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
664:     PetscObjectOptionsBegin((PetscObject)pc);
665:     PCSetFromOptions_MG(PetscOptionsObject,pc);
666:     PetscOptionsEnd();
667:     PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);

669:     /* clean up */
670:     for (level=1; level<pc_gamg->Nlevels; level++) {
671:       MatDestroy(&Parr[level]);
672:       MatDestroy(&Aarr[level]);
673:     }
674:     PCSetUp_MG(pc);
675:   } else {
676:     KSP smoother;
677:     PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");
678:     PCMGGetSmoother(pc, 0, &smoother);
679:     KSPSetOperators(smoother, Aarr[0], Aarr[0]);
680:     KSPSetType(smoother, KSPPREONLY);
681:     PCSetUp_MG(pc);
682:   }
683:   return(0);
684: }

686: /* ------------------------------------------------------------------------- */
687: /*
688:  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
689:    that was created with PCCreate_GAMG().

691:    Input Parameter:
692: .  pc - the preconditioner context

694:    Application Interface Routine: PCDestroy()
695: */
696: PetscErrorCode PCDestroy_GAMG(PC pc)
697: {
699:   PC_MG          *mg     = (PC_MG*)pc->data;
700:   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;

703:   PCReset_GAMG(pc);
704:   if (pc_gamg->ops->destroy) {
705:     (*pc_gamg->ops->destroy)(pc);
706:   }
707:   PetscFree(pc_gamg->ops);
708:   PetscFree(pc_gamg->gamg_type_name);
709:   PetscFree(pc_gamg);
710:   PCDestroy_MG(pc);
711:   return(0);
712: }

714: /*@
715:    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.

717:    Logically Collective on PC

719:    Input Parameters:
720: +  pc - the preconditioner context
721: -  n - the number of equations


724:    Options Database Key:
725: .  -pc_gamg_process_eq_limit <limit>

727:    Notes:
728:     GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process 
729:           that has degrees of freedom

731:    Level: intermediate

733:    Concepts: Unstructured multigrid preconditioner

735: .seealso: PCGAMGSetCoarseEqLim()
736: @*/
737: PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
738: {

743:   PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));
744:   return(0);
745: }

747: static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
748: {
749:   PC_MG   *mg      = (PC_MG*)pc->data;
750:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

753:   if (n>0) pc_gamg->min_eq_proc = n;
754:   return(0);
755: }

757: /*@
758:    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.

760:  Collective on PC

762:    Input Parameters:
763: +  pc - the preconditioner context
764: -  n - maximum number of equations to aim for

766:    Options Database Key:
767: .  -pc_gamg_coarse_eq_limit <limit>

769:    Level: intermediate

771:    Concepts: Unstructured multigrid preconditioner

773: .seealso: PCGAMGSetProcEqLim()
774: @*/
775: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
776: {

781:   PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
782:   return(0);
783: }

785: static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
786: {
787:   PC_MG   *mg      = (PC_MG*)pc->data;
788:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

791:   if (n>0) pc_gamg->coarse_eq_limit = n;
792:   return(0);
793: }

795: /*@
796:    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids

798:    Collective on PC

800:    Input Parameters:
801: +  pc - the preconditioner context
802: -  n - PETSC_TRUE or PETSC_FALSE

804:    Options Database Key:
805: .  -pc_gamg_repartition <true,false>

807:    Notes:
808:     this will generally improve the loading balancing of the work on each level

810:    Level: intermediate

812:    Concepts: Unstructured multigrid preconditioner

814: .seealso: ()
815: @*/
816: PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
817: {

822:   PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));
823:   return(0);
824: }

826: static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
827: {
828:   PC_MG   *mg      = (PC_MG*)pc->data;
829:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

832:   pc_gamg->repart = n;
833:   return(0);
834: }

836: /*@
837:    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner

839:    Collective on PC

841:    Input Parameters:
842: +  pc - the preconditioner context
843: -  n - PETSC_TRUE or PETSC_FALSE

845:    Options Database Key:
846: .  -pc_gamg_reuse_interpolation <true,false>

848:    Level: intermediate

850:    Notes:
851:     this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
852:           rebuilding the preconditioner quicker.

854:    Concepts: Unstructured multigrid preconditioner

856: .seealso: ()
857: @*/
858: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
859: {

864:   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
865:   return(0);
866: }

868: static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
869: {
870:   PC_MG   *mg      = (PC_MG*)pc->data;
871:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

874:   pc_gamg->reuse_prol = n;
875:   return(0);
876: }

878: /*@
879:    PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner.

881:    Collective on PC

883:    Input Parameters:
884: +  pc - the preconditioner context
885: -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not

887:    Options Database Key:
888: .  -pc_gamg_asm_use_agg

890:    Level: intermediate

892:    Concepts: Unstructured multigrid preconditioner

894: .seealso: ()
895: @*/
896: PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
897: {

902:   PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));
903:   return(0);
904: }

906: static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
907: {
908:   PC_MG   *mg      = (PC_MG*)pc->data;
909:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

912:   pc_gamg->use_aggs_in_asm = flg;
913:   return(0);
914: }

916: /*@
917:    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver

919:    Collective on PC

921:    Input Parameters:
922: +  pc - the preconditioner context
923: -  flg - PETSC_TRUE to not force coarse grid onto one processor

925:    Options Database Key:
926: .  -pc_gamg_use_parallel_coarse_grid_solver

928:    Level: intermediate

930:    Concepts: Unstructured multigrid preconditioner

932: .seealso: ()
933: @*/
934: PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
935: {

940:   PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));
941:   return(0);
942: }

944: static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
945: {
946:   PC_MG   *mg      = (PC_MG*)pc->data;
947:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

950:   pc_gamg->use_parallel_coarse_grid_solver = flg;
951:   return(0);
952: }

954: /*@
955:    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use

957:    Not collective on PC

959:    Input Parameters:
960: +  pc - the preconditioner
961: -  n - the maximum number of levels to use

963:    Options Database Key:
964: .  -pc_mg_levels

966:    Level: intermediate

968:    Concepts: Unstructured multigrid preconditioner

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

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

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

988:   pc_gamg->Nlevels = n;
989:   return(0);
990: }

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

995:    Not collective on PC

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

1001:    Options Database Key:
1002: .  -pc_gamg_threshold <threshold>

1004:    Notes:
1005:     Before aggregating the graph GAMG will remove small values from the graph thus reducing the coupling in the graph and a different
1006:     (perhaps better) coarser set of points.

1008:    Level: intermediate

1010:    Concepts: Unstructured multigrid preconditioner

1012: .seealso: PCGAMGFilterGraph()
1013: @*/
1014: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1015: {

1020:   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));
1021:   return(0);
1022: }

1024: static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1025: {
1026:   PC_MG   *mg      = (PC_MG*)pc->data;
1027:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1028:   PetscInt i;
1030:   for (i=0;i<n;i++) pc_gamg->threshold[i] = v[i];
1031:   do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1032:   return(0);
1033: }

1035: /*@
1036:    PCGAMGSetThresholdScale - Relative threshold reduction at each level

1038:    Not collective on PC

1040:    Input Parameters:
1041: +  pc - the preconditioner context
1042: -  scale - the threshold value reduction, ussually < 1.0

1044:    Options Database Key:
1045: .  -pc_gamg_threshold_scale <v>

1047:    Level: advanced

1049:    Concepts: Unstructured multigrid preconditioner

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

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

1063: static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1064: {
1065:   PC_MG   *mg      = (PC_MG*)pc->data;
1066:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1068:   pc_gamg->threshold_scale = v;
1069:   return(0);
1070: }

1072: /*@C
1073:    PCGAMGSetType - Set solution method

1075:    Collective on PC

1077:    Input Parameters:
1078: +  pc - the preconditioner context
1079: -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL

1081:    Options Database Key:
1082: .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply

1084:    Level: intermediate

1086:    Concepts: Unstructured multigrid preconditioner

1088: .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1089: @*/
1090: PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1091: {

1096:   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1097:   return(0);
1098: }

1100: /*@C
1101:    PCGAMGGetType - Get solution method

1103:    Collective on PC

1105:    Input Parameter:
1106: .  pc - the preconditioner context

1108:    Output Parameter:
1109: .  type - the type of algorithm used

1111:    Level: intermediate

1113:    Concepts: Unstructured multigrid preconditioner

1115: .seealso: PCGAMGSetType(), PCGAMGType
1116: @*/
1117: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1118: {

1123:   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1124:   return(0);
1125: }

1127: static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1128: {
1129:   PC_MG          *mg      = (PC_MG*)pc->data;
1130:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1133:   *type = pc_gamg->type;
1134:   return(0);
1135: }

1137: static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1138: {
1139:   PetscErrorCode ierr,(*r)(PC);
1140:   PC_MG          *mg      = (PC_MG*)pc->data;
1141:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1144:   pc_gamg->type = type;
1145:   PetscFunctionListFind(GAMGList,type,&r);
1146:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1147:   if (pc_gamg->ops->destroy) {
1148:     (*pc_gamg->ops->destroy)(pc);
1149:     PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));
1150:     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1151:     /* cleaning up common data in pc_gamg - this should disapear someday */
1152:     pc_gamg->data_cell_cols = 0;
1153:     pc_gamg->data_cell_rows = 0;
1154:     pc_gamg->orig_data_cell_cols = 0;
1155:     pc_gamg->orig_data_cell_rows = 0;
1156:     PetscFree(pc_gamg->data);
1157:     pc_gamg->data_sz = 0;
1158:   }
1159:   PetscFree(pc_gamg->gamg_type_name);
1160:   PetscStrallocpy(type,&pc_gamg->gamg_type_name);
1161:   (*r)(pc);
1162:   return(0);
1163: }

1165: static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
1166: {
1167:   PetscErrorCode ierr,i;
1168:   PC_MG          *mg      = (PC_MG*)pc->data;
1169:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1172:   PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");
1173:   PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");
1174:   for (i=0;i<pc_gamg->current_level;i++) {
1175:     PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);
1176:   }
1177:   PetscViewerASCIIPrintf(viewer,"\n");
1178:   PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);
1179:   if (pc_gamg->use_aggs_in_asm) {
1180:     PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");
1181:   }
1182:   if (pc_gamg->use_parallel_coarse_grid_solver) {
1183:     PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");
1184:   }
1185:   if (pc_gamg->ops->view) {
1186:     (*pc_gamg->ops->view)(pc,viewer);
1187:   }
1188:   return(0);
1189: }

1191: PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
1192: {
1194:   PC_MG          *mg      = (PC_MG*)pc->data;
1195:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1196:   PetscBool      flag;
1197:   MPI_Comm       comm;
1198:   char           prefix[256];
1199:   PetscInt       i,n;
1200:   const char     *pcpre;

1203:   PetscObjectGetComm((PetscObject)pc,&comm);
1204:   PetscOptionsHead(PetscOptionsObject,"GAMG options");
1205:   {
1206:     char tname[256];
1207:     PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);
1208:     if (flag) {
1209:       PCGAMGSetType(pc,tname);
1210:     }
1211:     PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);
1212:     PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);
1213:     PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation aggregates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL);
1214:     PetscOptionsBool("-pc_gamg_use_parallel_coarse_grid_solver","Use parallel coarse grid solver (otherwise put last grid on one process)","PCGAMGSetUseParallelCoarseGridSolve",pc_gamg->use_parallel_coarse_grid_solver,&pc_gamg->use_parallel_coarse_grid_solver,NULL);
1215:     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);
1216:     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);
1217:     PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL);
1218:     n = PETSC_GAMG_MAXLEVELS;
1219:     PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);
1220:     if (!flag || n < PETSC_GAMG_MAXLEVELS) {
1221:       if (!flag) n = 1;
1222:       i = n;
1223:       do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1224:     }
1225:     PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);

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

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

1240:    Options Database Keys:
1241: +   -pc_gamg_type <type> - one of agg, geo, or classical
1242: .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1243: .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1244: .   -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother
1245: .   -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit>
1246:                                         equations on each process that has degrees of freedom
1247: .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
1248: .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1249: -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified

1251:    Options Database Keys for default Aggregation:
1252: +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1253: .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1254: -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it

1256:    Multigrid options:
1257: +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1258: .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1259: .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1260: -  -pc_mg_levels <levels> - Number of levels of multigrid to use.


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

1269:   Level: intermediate

1271:   Concepts: algebraic multigrid

1273: .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1274:            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation()
1275: M*/

1277: PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1278: {
1279:   PetscErrorCode ierr,i;
1280:   PC_GAMG        *pc_gamg;
1281:   PC_MG          *mg;

1284:    /* register AMG type */
1285:   PCGAMGInitializePackage();

1287:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1288:   PCSetType(pc, PCMG);
1289:   PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);

1291:   /* create a supporting struct and attach it to pc */
1292:   PetscNewLog(pc,&pc_gamg);
1293:   PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);
1294:   mg           = (PC_MG*)pc->data;
1295:   mg->innerctx = pc_gamg;

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

1299:   pc_gamg->setup_count = 0;
1300:   /* these should be in subctx but repartitioning needs simple arrays */
1301:   pc_gamg->data_sz = 0;
1302:   pc_gamg->data    = 0;

1304:   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1305:   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1306:   pc->ops->setup          = PCSetUp_GAMG;
1307:   pc->ops->reset          = PCReset_GAMG;
1308:   pc->ops->destroy        = PCDestroy_GAMG;
1309:   mg->view                = PCView_GAMG;

1311:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);
1312:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);
1313:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);
1314:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);
1315:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);
1316:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);
1317:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);
1318:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);
1319:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);
1320:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);
1321:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);
1322:   pc_gamg->repart           = PETSC_FALSE;
1323:   pc_gamg->reuse_prol       = PETSC_FALSE;
1324:   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1325:   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1326:   pc_gamg->min_eq_proc      = 50;
1327:   pc_gamg->coarse_eq_limit  = 50;
1328:   for (i=0;i<PETSC_GAMG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1329:   pc_gamg->threshold_scale = 1.;
1330:   pc_gamg->Nlevels          = PETSC_GAMG_MAXLEVELS;
1331:   pc_gamg->current_level    = 0; /* don't need to init really */
1332:   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;

1334:   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1335:   PCGAMGSetType(pc,PCGAMGAGG);
1336:   return(0);
1337: }

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

1344:  Level: developer

1346:  .keywords: PC, PCGAMG, initialize, package
1347:  .seealso: PetscInitialize()
1348: @*/
1349: PetscErrorCode PCGAMGInitializePackage(void)
1350: {

1354:   if (PCGAMGPackageInitialized) return(0);
1355:   PCGAMGPackageInitialized = PETSC_TRUE;
1356:   PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);
1357:   PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);
1358:   PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);
1359:   PetscRegisterFinalize(PCGAMGFinalizePackage);

1361:   /* general events */
1362:   PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);
1363:   PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);
1364:   PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1365:   PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1366:   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1367:   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1368:   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);

1370: #if defined PETSC_GAMG_USE_LOG
1371:   PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1372:   PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1373:   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1374:   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1375:   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1376:   PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1377:   PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1378:   PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1379:   PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1380:   PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1381:   PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1382:   PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1383:   PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1384:   PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1385:   PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1386:   PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1387:   PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);

1389:   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1390:   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1391:   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1392:   /* create timer stages */
1393: #if defined GAMG_STAGES
1394:   {
1395:     char     str[32];
1396:     PetscInt lidx;
1397:     sprintf(str,"MG Level %d (finest)",0);
1398:     PetscLogStageRegister(str, &gamg_stages[0]);
1399:     for (lidx=1; lidx<9; lidx++) {
1400:       sprintf(str,"MG Level %d",lidx);
1401:       PetscLogStageRegister(str, &gamg_stages[lidx]);
1402:     }
1403:   }
1404: #endif
1405: #endif
1406:   return(0);
1407: }

1409: /*@C
1410:  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1411:     called from PetscFinalize() automatically.

1413:  Level: developer

1415:  .keywords: Petsc, destroy, package
1416:  .seealso: PetscFinalize()
1417: @*/
1418: PetscErrorCode PCGAMGFinalizePackage(void)
1419: {

1423:   PCGAMGPackageInitialized = PETSC_FALSE;
1424:   PetscFunctionListDestroy(&GAMGList);
1425:   return(0);
1426: }

1428: /*@C
1429:  PCGAMGRegister - Register a PCGAMG implementation.

1431:  Input Parameters:
1432:  + type - string that will be used as the name of the GAMG type.
1433:  - create - function for creating the gamg context.

1435:   Level: advanced

1437:  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1438: @*/
1439: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1440: {

1444:   PCGAMGInitializePackage();
1445:   PetscFunctionListAdd(&GAMGList,type,create);
1446:   return(0);
1447: }