Actual source code: gamg.c

petsc-master 2016-06-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: #if defined PETSC_GAMG_USE_LOG
243:         PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);
244: #endif
245:         return(0);
246:       }

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

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

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

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

328:     pc_gamg->data_sz = node_data_sz*ncrs_new;
329:     strideNew        = ncrs_new*ndata_rows;

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

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

374: #if defined PETSC_GAMG_USE_LOG
375:     PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);
376: #endif
377:     /* prolongator */
378:     {
379:       IS       findices;
380:       PetscInt Istart,Iend;
381:       Mat      Pnew;
382: 
383:       MatGetOwnershipRange(Pold, &Istart, &Iend);
384: #if defined PETSC_GAMG_USE_LOG
385:       PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);
386: #endif
387:       ISCreateStride(comm,Iend-Istart,Istart,1,&findices);
388:       ISSetBlockSize(findices,f_bs);
389:       MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);
390:       ISDestroy(&findices);

392: #if defined PETSC_GAMG_USE_LOG
393:       PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);
394: #endif
395:       MatDestroy(a_P_inout);

397:       /* output - repartitioned */
398:       *a_P_inout = Pnew;
399:     }
400:     ISDestroy(&new_eq_indices);

402:     *a_nactive_proc = new_size; /* output */
403:   }
404:   return(0);
405: }

407: /* -------------------------------------------------------------------------- */
408: /*
409:    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
410:                     by setting data structures and options.

412:    Input Parameter:
413: .  pc - the preconditioner context

415: */
418: PetscErrorCode PCSetUp_GAMG(PC pc)
419: {
421:   PC_MG          *mg      = (PC_MG*)pc->data;
422:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
423:   Mat            Pmat     = pc->pmat;
424:   PetscInt       fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS];
425:   MPI_Comm       comm;
426:   PetscMPIInt    rank,size,nactivepe;
427:   Mat            Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS];
428:   IS             *ASMLocalIDsArr[GAMG_MAXLEVELS];
429:   PetscLogDouble nnz0=0.,nnztot=0.;
430:   MatInfo        info;

433:   PetscObjectGetComm((PetscObject)pc,&comm);
434:   MPI_Comm_rank(comm,&rank);
435:   MPI_Comm_size(comm,&size);

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

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

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

460:             mglevels[level]->A = B;
461:           } else {
462:             KSPGetOperators(mglevels[level]->smoothd,NULL,&B);
463:             MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);
464:           }
465:           KSPSetOperators(mglevels[level]->smoothd,B,B);
466:           dB   = B;
467:         }
468:       }

470:       PCSetUp_MG(pc);
471:       return(0);
472:     }
473:   }

475:   if (!pc_gamg->data) {
476:     if (pc_gamg->orig_data) {
477:       MatGetBlockSize(Pmat, &bs);
478:       MatGetLocalSize(Pmat, &qq, NULL);

480:       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
481:       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
482:       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;

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

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

500:   /* get basic dims */
501:   MatGetBlockSize(Pmat, &bs);
502:   MatGetSize(Pmat, &M, &qq);

504:   MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info); /* global reduction */
505:   nnz0   = info.nz_used;
506:   nnztot = info.nz_used;
507:   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);

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

524:       pc_gamg->ops->graph(pc,Aarr[level], &Gmat);
525:       pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);
526:       pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);

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

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

538:         Parr[level1] = Prol11;
539:       } else Parr[level1] = NULL;

541:       if (pc_gamg->use_aggs_in_gasm) {
542:         PetscInt bs;
543:         MatGetBlockSizes(Prol11, &bs, NULL);
544:         PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);
545:       }

547:       MatDestroy(&Gmat);
548:       PetscCDDestroy(agg_lists);
549:     } /* construct prolongator scope */
550: #if defined PETSC_GAMG_USE_LOG
551:     PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);
552: #endif
553:     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
554:     if (!Parr[level1]) {
555:        PetscInfo1(pc,"Stop gridding, level %D\n",level);
556: #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
557:       PetscLogStagePop();
558: #endif
559:       break;
560:     }
561: #if defined PETSC_GAMG_USE_LOG
562:     PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);
563: #endif

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

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, &qq);
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 GASM smoother that uses the 'aggregates' */
610:       if (pc_gamg->use_aggs_in_gasm) {
611:         PetscInt sz;
612:         IS       *is;

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

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

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

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

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

706:   PCReset_GAMG(pc);
707:   if (pc_gamg->ops->destroy) {
708:     (*pc_gamg->ops->destroy)(pc);
709:   }
710:   PetscRandomDestroy(&pc_gamg->random);
711:   PetscFree(pc_gamg->ops);
712:   PetscFree(pc_gamg->gamg_type_name);
713:   PetscFree(pc_gamg);
714:   PCDestroy_MG(pc);
715:   return(0);
716: }

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

723:    Logically Collective on PC

725:    Input Parameters:
726: +  pc - the preconditioner context
727: -  n - the number of equations


730:    Options Database Key:
731: .  -pc_gamg_process_eq_limit <limit>

733:    Level: intermediate

735:    Concepts: Unstructured multigrid preconditioner

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

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

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

757:   if (n>0) pc_gamg->min_eq_proc = n;
758:   return(0);
759: }

763: /*@
764:    PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids.

766:  Collective on PC

768:    Input Parameters:
769: +  pc - the preconditioner context
770: -  n - maximum number of equations to aim for

772:    Options Database Key:
773: .  -pc_gamg_coarse_eq_limit <limit>

775:    Level: intermediate

777:    Concepts: Unstructured multigrid preconditioner

779: .seealso: PCGAMGSetProcEqLim()
780: @*/
781: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
782: {

787:   PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
788:   return(0);
789: }

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

799:   if (n>0) pc_gamg->coarse_eq_limit = n;
800:   return(0);
801: }

805: /*@
806:    PCGAMGSetRepartitioning - Repartition the coarse grids

808:    Collective on PC

810:    Input Parameters:
811: +  pc - the preconditioner context
812: -  n - PETSC_TRUE or PETSC_FALSE

814:    Options Database Key:
815: .  -pc_gamg_repartition <true,false>

817:    Level: intermediate

819:    Concepts: Unstructured multigrid preconditioner

821: .seealso: ()
822: @*/
823: PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n)
824: {

829:   PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));
830:   return(0);
831: }

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

841:   pc_gamg->repart = n;
842:   return(0);
843: }

847: /*@
848:    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner

850:    Collective on PC

852:    Input Parameters:
853: +  pc - the preconditioner context
854: -  n - PETSC_TRUE or PETSC_FALSE

856:    Options Database Key:
857: .  -pc_gamg_reuse_interpolation <true,false>

859:    Level: intermediate

861:    Concepts: Unstructured multigrid preconditioner

863: .seealso: ()
864: @*/
865: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
866: {

871:   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
872:   return(0);
873: }

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

883:   pc_gamg->reuse_prol = n;
884:   return(0);
885: }

889: /*@
890:    PCGAMGSetUseASMAggs -

892:    Collective on PC

894:    Input Parameters:
895: .  pc - the preconditioner context


898:    Options Database Key:
899: .  -pc_gamg_use_agg_gasm

901:    Level: intermediate

903:    Concepts: Unstructured multigrid preconditioner

905: .seealso: ()
906: @*/
907: PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n)
908: {

913:   PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));
914:   return(0);
915: }

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

925:   pc_gamg->use_aggs_in_gasm = n;
926:   return(0);
927: }

931: /*@
932:    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use

934:    Not collective on PC

936:    Input Parameters:
937: +  pc - the preconditioner
938: -  n - the maximum number of levels to use

940:    Options Database Key:
941: .  -pc_mg_levels

943:    Level: intermediate

945:    Concepts: Unstructured multigrid preconditioner

947: .seealso: ()
948: @*/
949: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
950: {

955:   PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
956:   return(0);
957: }

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

967:   pc_gamg->Nlevels = n;
968:   return(0);
969: }

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

976:    Not collective on PC

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

982:    Options Database Key:
983: .  -pc_gamg_threshold <threshold>

985:    Level: intermediate

987:    Concepts: Unstructured multigrid preconditioner

989: .seealso: ()
990: @*/
991: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n)
992: {

997:   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));
998:   return(0);
999: }

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

1009:   pc_gamg->threshold = n;
1010:   return(0);
1011: }

1015: /*@
1016:    PCGAMGSetType - Set solution method

1018:    Collective on PC

1020:    Input Parameters:
1021: +  pc - the preconditioner context
1022: -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL

1024:    Options Database Key:
1025: .  -pc_gamg_type <agg,geo,classical>

1027:    Level: intermediate

1029:    Concepts: Unstructured multigrid preconditioner

1031: .seealso: PCGAMGGetType(), PCGAMG
1032: @*/
1033: PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1034: {

1039:   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1040:   return(0);
1041: }

1045: /*@
1046:    PCGAMGGetType - Get solution method

1048:    Collective on PC

1050:    Input Parameter:
1051: .  pc - the preconditioner context

1053:    Output Parameter:
1054: .  type - the type of algorithm used

1056:    Level: intermediate

1058:    Concepts: Unstructured multigrid preconditioner

1060: .seealso: PCGAMGSetType(), PCGAMGType
1061: @*/
1062: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1063: {

1068:   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1069:   return(0);
1070: }

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

1080:   *type = pc_gamg->type;
1081:   return(0);
1082: }

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

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

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

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

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

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

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

1171: /* -------------------------------------------------------------------------- */
1172: /*MC
1173:      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner

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


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

1188:   Level: intermediate

1190:   Concepts: algebraic multigrid

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

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

1205:   /* register AMG type */
1206:   PCGAMGInitializePackage();

1208:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1209:   PCSetType(pc, PCMG);
1210:   PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);

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

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

1220:   pc_gamg->setup_count = 0;
1221:   /* these should be in subctx but repartitioning needs simple arrays */
1222:   pc_gamg->data_sz = 0;
1223:   pc_gamg->data    = 0;

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

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

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

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

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

1265:  Level: developer

1267:  .keywords: PC, PCGAMG, initialize, package
1268:  .seealso: PetscInitialize()
1269: @*/
1270: PetscErrorCode PCGAMGInitializePackage(void)
1271: {

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

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

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

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

1332: /*@C
1333:  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1334:     called from PetscFinalize() automatically.

1336:  Level: developer

1338:  .keywords: Petsc, destroy, package
1339:  .seealso: PetscFinalize()
1340: @*/
1341: PetscErrorCode PCGAMGFinalizePackage(void)
1342: {

1346:   PCGAMGPackageInitialized = PETSC_FALSE;
1347:   PetscFunctionListDestroy(&GAMGList);
1348:   return(0);
1349: }

1353: /*@C
1354:  PCGAMGRegister - Register a PCGAMG implementation.

1356:  Input Parameters:
1357:  + type - string that will be used as the name of the GAMG type.
1358:  - create - function for creating the gamg context.

1360:   Level: advanced

1362:  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1363: @*/
1364: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1365: {

1369:   PCGAMGInitializePackage();
1370:   PetscFunctionListAdd(&GAMGList,type,create);
1371:   return(0);
1372: }