Actual source code: gamg.c

petsc-master 2016-05-24
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) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen, cleaned up in SetUp\n");
 44:   pc_gamg->data_sz = 0;
 45:   PetscFree(pc_gamg->orig_data);
 46:   return(0);
 47: }

 49: /* -------------------------------------------------------------------------- */
 50: /*
 51:    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
 52:      of active processors.

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

 68: 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)
 69: {
 70:   PetscErrorCode  ierr;
 71:   PC_MG           *mg         = (PC_MG*)pc->data;
 72:   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
 73:   Mat             Cmat,Pold=*a_P_inout;
 74:   MPI_Comm        comm;
 75:   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
 76:   PetscInt        ncrs_eq,ncrs,f_bs;

 79:   PetscObjectGetComm((PetscObject)Amat_fine,&comm);
 80:   MPI_Comm_rank(comm, &rank);
 81:   MPI_Comm_size(comm, &size);
 82:   MatGetBlockSize(Amat_fine, &f_bs);
 83:   MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);

 85:   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
 86:   MatGetLocalSize(Cmat, &ncrs_eq, NULL);
 87:   if (pc_gamg->data_cell_rows>0) {
 88:     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
 89:   } else {
 90:     PetscInt  bs;
 91:     MatGetBlockSize(Cmat, &bs);
 92:     ncrs = ncrs_eq/bs;
 93:   }

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

104:   if (Pcolumnperm) *Pcolumnperm = NULL;

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

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

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

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

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

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

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

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

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

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

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

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

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

238:       if (new_size==nactive) {
239:         *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */
240:         PetscFree(counts);
241:         PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);
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 MatGetSubMatrix
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:       MatGetSubMatrix(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;
379: 
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:       MatGetSubMatrix(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: */
415: PetscErrorCode PCSetUp_GAMG(PC pc)
416: {
418:   PC_MG          *mg      = (PC_MG*)pc->data;
419:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
420:   Mat            Pmat     = pc->pmat;
421:   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
422:   MPI_Comm       comm;
423:   PetscMPIInt    rank,size,nactivepe;
424:   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
425:   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
426:   PetscLogDouble nnz0=0.,nnztot=0.;
427:   MatInfo        info;

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

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

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

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

457:             mglevels[level]->A = B;
458:           } else {
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, &qq);

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;

538:       if (pc_gamg->use_aggs_in_gasm) {
539:         PetscInt bs;
540:         MatGetBlockSizes(Prol11, &bs, NULL);
541:         PetscCDGetASMBlocks(agg_lists, bs, &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]) {
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

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

564: #if defined PETSC_GAMG_USE_LOG
565:     PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);
566: #endif
567:     MatGetSize(Aarr[level1], &M, &qq);
568:     MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);
569:     nnztot += info.nz_used;
570:     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);

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

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

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

595:       PCMGGetSmoother(pc, lidx, &smoother);
596:       KSPGetPC(smoother, &subpc);

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

603:       /* set defaults */
604:       KSPSetType(smoother, KSPCHEBYSHEV);

606:       /* set blocks for GASM smoother that uses the 'aggregates' */
607:       if (pc_gamg->use_aggs_in_gasm) {
608:         PetscInt sz;
609:         IS       *is;

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

660:     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
661:     PetscObjectOptionsBegin((PetscObject)pc);
662:     PCSetFromOptions_MG(PetscOptionsObject,pc);
663:     PetscOptionsEnd();
664:     if (!mg->galerkin) SETERRQ(comm,PETSC_ERR_USER,"PCGAMG must use Galerkin for coarse operators.");
665:     if (mg->galerkin == 1) mg->galerkin = 2;

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

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

689:    Input Parameter:
690: .  pc - the preconditioner context

692:    Application Interface Routine: PCDestroy()
693: */
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:   PetscRandomDestroy(&pc_gamg->random);
708:   PetscFree(pc_gamg->ops);
709:   PetscFree(pc_gamg->gamg_type_name);
710:   PetscFree(pc_gamg);
711:   PCDestroy_MG(pc);
712:   return(0);
713: }

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

720:    Logically Collective on PC

722:    Input Parameters:
723: +  pc - the preconditioner context
724: -  n - the number of equations


727:    Options Database Key:
728: .  -pc_gamg_process_eq_limit <limit>

730:    Level: intermediate

732:    Concepts: Unstructured multigrid preconditioner

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

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

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

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

760: /*@
761:    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.

763:  Collective on PC

765:    Input Parameters:
766: +  pc - the preconditioner context
767: -  n - maximum number of equations to aim for

769:    Options Database Key:
770: .  -pc_gamg_coarse_eq_limit <limit>

772:    Level: intermediate

774:    Concepts: Unstructured multigrid preconditioner

776: .seealso: PCGAMGSetProcEqLim()
777: @*/
778: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
779: {

784:   PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
785:   return(0);
786: }

790: static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
791: {
792:   PC_MG   *mg      = (PC_MG*)pc->data;
793:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

796:   if (n>0) pc_gamg->coarse_eq_limit = n;
797:   return(0);
798: }

802: /*@
803:    PCGAMGSetRepartitioning - Repartition the coarse grids

805:    Collective on PC

807:    Input Parameters:
808: +  pc - the preconditioner context
809: -  n - PETSC_TRUE or PETSC_FALSE

811:    Options Database Key:
812: .  -pc_gamg_repartition <true,false>

814:    Level: intermediate

816:    Concepts: Unstructured multigrid preconditioner

818: .seealso: ()
819: @*/
820: PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
821: {

826:   PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));
827:   return(0);
828: }

832: static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n)
833: {
834:   PC_MG   *mg      = (PC_MG*)pc->data;
835:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

838:   pc_gamg->repart = n;
839:   return(0);
840: }

844: /*@
845:    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner

847:    Collective on PC

849:    Input Parameters:
850: +  pc - the preconditioner context
851: -  n - PETSC_TRUE or PETSC_FALSE

853:    Options Database Key:
854: .  -pc_gamg_reuse_interpolation <true,false>

856:    Level: intermediate

858:    Concepts: Unstructured multigrid preconditioner

860: .seealso: ()
861: @*/
862: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
863: {

868:   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
869:   return(0);
870: }

874: static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
875: {
876:   PC_MG   *mg      = (PC_MG*)pc->data;
877:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

880:   pc_gamg->reuse_prol = n;
881:   return(0);
882: }

886: /*@
887:    PCGAMGSetUseASMAggs -

889:    Collective on PC

891:    Input Parameters:
892: .  pc - the preconditioner context


895:    Options Database Key:
896: .  -pc_gamg_use_agg_gasm

898:    Level: intermediate

900:    Concepts: Unstructured multigrid preconditioner

902: .seealso: ()
903: @*/
904: PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
905: {

910:   PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));
911:   return(0);
912: }

916: static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n)
917: {
918:   PC_MG   *mg      = (PC_MG*)pc->data;
919:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

922:   pc_gamg->use_aggs_in_gasm = n;
923:   return(0);
924: }

928: /*@
929:    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use

931:    Not collective on PC

933:    Input Parameters:
934: +  pc - the preconditioner
935: -  n - the maximum number of levels to use

937:    Options Database Key:
938: .  -pc_mg_levels

940:    Level: intermediate

942:    Concepts: Unstructured multigrid preconditioner

944: .seealso: ()
945: @*/
946: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
947: {

952:   PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
953:   return(0);
954: }

958: static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
959: {
960:   PC_MG   *mg      = (PC_MG*)pc->data;
961:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

964:   pc_gamg->Nlevels = n;
965:   return(0);
966: }

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

973:    Not collective on PC

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

979:    Options Database Key:
980: .  -pc_gamg_threshold <threshold>

982:    Level: intermediate

984:    Concepts: Unstructured multigrid preconditioner

986: .seealso: ()
987: @*/
988: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
989: {

994:   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));
995:   return(0);
996: }

1000: static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n)
1001: {
1002:   PC_MG   *mg      = (PC_MG*)pc->data;
1003:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1006:   pc_gamg->threshold = n;
1007:   return(0);
1008: }

1012: /*@
1013:    PCGAMGSetType - Set solution method

1015:    Collective on PC

1017:    Input Parameters:
1018: +  pc - the preconditioner context
1019: -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL

1021:    Options Database Key:
1022: .  -pc_gamg_type <agg,geo,classical>

1024:    Level: intermediate

1026:    Concepts: Unstructured multigrid preconditioner

1028: .seealso: PCGAMGGetType(), PCGAMG
1029: @*/
1030: PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1031: {

1036:   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1037:   return(0);
1038: }

1042: /*@
1043:    PCGAMGGetType - Get solution method

1045:    Collective on PC

1047:    Input Parameter:
1048: .  pc - the preconditioner context

1050:    Output Parameter:
1051: .  type - the type of algorithm used

1053:    Level: intermediate

1055:    Concepts: Unstructured multigrid preconditioner

1057: .seealso: PCGAMGSetType(), PCGAMGType
1058: @*/
1059: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1060: {

1065:   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1066:   return(0);
1067: }

1071: static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1072: {
1073:   PC_MG          *mg      = (PC_MG*)pc->data;
1074:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1077:   *type = pc_gamg->type;
1078:   return(0);
1079: }

1083: static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1084: {
1085:   PetscErrorCode ierr,(*r)(PC);
1086:   PC_MG          *mg      = (PC_MG*)pc->data;
1087:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1090:   pc_gamg->type = type;
1091:   PetscFunctionListFind(GAMGList,type,&r);
1092:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1093:   if (pc_gamg->ops->destroy) {
1094:     (*pc_gamg->ops->destroy)(pc);
1095:     PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));
1096:     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1097:     /* cleaning up common data in pc_gamg - this should disapear someday */
1098:     pc_gamg->data_cell_cols = 0;
1099:     pc_gamg->data_cell_rows = 0;
1100:     pc_gamg->orig_data_cell_cols = 0;
1101:     pc_gamg->orig_data_cell_rows = 0;
1102:     PetscFree(pc_gamg->data);
1103:     pc_gamg->data_sz = 0;
1104:   }
1105:   PetscFree(pc_gamg->gamg_type_name);
1106:   PetscStrallocpy(type,&pc_gamg->gamg_type_name);
1107:   (*r)(pc);
1108:   return(0);
1109: }

1113: static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
1114: {
1116:   PC_MG          *mg      = (PC_MG*)pc->data;
1117:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1120:   PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");
1121:   PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values from graph %g\n",(double)pc_gamg->threshold);
1122:   if (pc_gamg->ops->view) {
1123:     (*pc_gamg->ops->view)(pc,viewer);
1124:   }
1125:   return(0);
1126: }

1130: PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
1131: {
1133:   PC_MG          *mg      = (PC_MG*)pc->data;
1134:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1135:   PetscBool      flag;
1136:   MPI_Comm       comm;
1137:   char           prefix[256];
1138:   const char     *pcpre;

1141:   PetscObjectGetComm((PetscObject)pc,&comm);
1142:   PetscOptionsHead(PetscOptionsObject,"GAMG options");
1143:   {
1144:     char tname[256];
1145:     PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);
1146:     if (flag) {
1147:       PCGAMGSetType(pc,tname);
1148:     }
1149:     PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);
1150:     PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);
1151:     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);
1152:     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);
1153:     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);
1154:     PetscOptionsReal("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&pc_gamg->threshold,&flag);
1155:     PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);

1157:     /* set options for subtype */
1158:     if (pc_gamg->ops->setfromoptions) {(*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);}
1159:   }
1160:   PCGetOptionsPrefix(pc, &pcpre);
1161:   PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
1162:   PetscObjectSetOptionsPrefix((PetscObject)pc_gamg->random,prefix);
1163:   PetscRandomSetFromOptions(pc_gamg->random);
1164:   PetscOptionsTail();
1165:   return(0);
1166: }

1168: /* -------------------------------------------------------------------------- */
1169: /*MC
1170:      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner

1172:    Options Database Keys:
1173:    Multigrid options(inherited)
1174: +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1175: .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1176: .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1177: -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade


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

1185:   Level: intermediate

1187:   Concepts: algebraic multigrid

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

1195: PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1196: {
1198:   PC_GAMG        *pc_gamg;
1199:   PC_MG          *mg;

1202:   /* register AMG type */
1203:   PCGAMGInitializePackage();

1205:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1206:   PCSetType(pc, PCMG);
1207:   PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);

1209:   /* create a supporting struct and attach it to pc */
1210:   PetscNewLog(pc,&pc_gamg);
1211:   mg           = (PC_MG*)pc->data;
1212:   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally from PCMG by GAMG code */
1213:   mg->innerctx = pc_gamg;

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

1217:   pc_gamg->setup_count = 0;
1218:   /* these should be in subctx but repartitioning needs simple arrays */
1219:   pc_gamg->data_sz = 0;
1220:   pc_gamg->data    = 0;

1222:   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1223:   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1224:   pc->ops->setup          = PCSetUp_GAMG;
1225:   pc->ops->reset          = PCReset_GAMG;
1226:   pc->ops->destroy        = PCDestroy_GAMG;
1227:   mg->view                = PCView_GAMG;

1229:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);
1230:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);
1231:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);
1232:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);
1233:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);
1234:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);
1235:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);
1236:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);
1237:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);
1238:   pc_gamg->repart           = PETSC_FALSE;
1239:   pc_gamg->reuse_prol       = PETSC_FALSE;
1240:   pc_gamg->use_aggs_in_gasm = PETSC_FALSE;
1241:   pc_gamg->min_eq_proc      = 50;
1242:   pc_gamg->coarse_eq_limit  = 50;
1243:   pc_gamg->threshold        = 0.;
1244:   pc_gamg->Nlevels          = GAMG_MAXLEVELS;
1245:   pc_gamg->current_level    = 0; /* don't need to init really */
1246:   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;

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

1250:   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1251:   PCGAMGSetType(pc,PCGAMGAGG);
1252:   return(0);
1253: }

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

1262:  Level: developer

1264:  .keywords: PC, PCGAMG, initialize, package
1265:  .seealso: PetscInitialize()
1266: @*/
1267: PetscErrorCode PCGAMGInitializePackage(void)
1268: {

1272:   if (PCGAMGPackageInitialized) return(0);
1273:   PCGAMGPackageInitialized = PETSC_TRUE;
1274:   PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);
1275:   PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);
1276:   PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);
1277:   PetscRegisterFinalize(PCGAMGFinalizePackage);

1279:   /* general events */
1280:   PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);
1281:   PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);
1282:   PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1283:   PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1284:   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1285:   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1286:   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);

1288: #if defined PETSC_GAMG_USE_LOG
1289:   PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1290:   PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1291:   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1292:   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1293:   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1294:   PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1295:   PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1296:   PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1297:   PetscLogEventRegister("    search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1298:   PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1299:   PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1300:   PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1301:   PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1302:   PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1303:   PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1304:   PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1305:   PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);

1307:   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1308:   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1309:   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1310:   /* create timer stages */
1311: #if defined GAMG_STAGES
1312:   {
1313:     char     str[32];
1314:     PetscInt lidx;
1315:     sprintf(str,"MG Level %d (finest)",0);
1316:     PetscLogStageRegister(str, &gamg_stages[0]);
1317:     for (lidx=1; lidx<9; lidx++) {
1318:       sprintf(str,"MG Level %d",lidx);
1319:       PetscLogStageRegister(str, &gamg_stages[lidx]);
1320:     }
1321:   }
1322: #endif
1323: #endif
1324:   return(0);
1325: }

1329: /*@C
1330:  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1331:     called from PetscFinalize() automatically.

1333:  Level: developer

1335:  .keywords: Petsc, destroy, package
1336:  .seealso: PetscFinalize()
1337: @*/
1338: PetscErrorCode PCGAMGFinalizePackage(void)
1339: {

1343:   PCGAMGPackageInitialized = PETSC_FALSE;
1344:   PetscFunctionListDestroy(&GAMGList);
1345:   return(0);
1346: }

1350: /*@C
1351:  PCGAMGRegister - Register a PCGAMG implementation.

1353:  Input Parameters:
1354:  + type - string that will be used as the name of the GAMG type.
1355:  - create - function for creating the gamg context.

1357:   Level: advanced

1359:  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1360: @*/
1361: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1362: {

1366:   PCGAMGInitializePackage();
1367:   PetscFunctionListAdd(&GAMGList,type,create);
1368:   return(0);
1369: }