Actual source code: gamg.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: /*
  2:  GAMG geometric-algebric multigrid PC - Mark Adams 2011
  3:  */
  4:  #include <petsc/private/matimpl.h>
  5:  #include <../src/ksp/pc/impls/gamg/gamg.h>
  6: #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */

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

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

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

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

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

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

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

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

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

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

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

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

 98:   if (Pcolumnperm) *Pcolumnperm = NULL;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

322:     pc_gamg->data_sz = node_data_sz*ncrs_new;
323:     strideNew        = ncrs_new*ndata_rows;

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

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

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

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

386: #if defined PETSC_GAMG_USE_LOG
387:       PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);
388: #endif
389:       MatDestroy(a_P_inout);

391:       /* output - repartitioned */
392:       *a_P_inout = Pnew;
393:     }
394:     ISDestroy(&new_eq_indices);

396:     *a_nactive_proc = new_size; /* output */
397:   }
398:   return(0);
399: }

401: /* -------------------------------------------------------------------------- */
402: /*
403:    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
404:                     by setting data structures and options.

406:    Input Parameter:
407: .  pc - the preconditioner context

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

426:   PetscObjectGetComm((PetscObject)pc,&comm);
427:   MPI_Comm_rank(comm,&rank);
428:   MPI_Comm_size(comm,&size);

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

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

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

453:             mglevels[level]->A = B;
454:           } else {
455:             KSPGetOperators(mglevels[level]->smoothd,NULL,&B);
456:             MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);
457:           }
458:           KSPSetOperators(mglevels[level]->smoothd,B,B);
459:           dB   = B;
460:         }
461:       }

463:       PCSetUp_MG(pc);
464:       return(0);
465:     }
466:   }

468:   if (!pc_gamg->data) {
469:     if (pc_gamg->orig_data) {
470:       MatGetBlockSize(Pmat, &bs);
471:       MatGetLocalSize(Pmat, &qq, NULL);

473:       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
474:       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
475:       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;

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

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

493:   /* get basic dims */
494:   MatGetBlockSize(Pmat, &bs);
495:   MatGetSize(Pmat, &M, &N);

497:   MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info); /* global reduction */
498:   nnz0   = info.nz_used;
499:   nnztot = info.nz_used;
500:   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);

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

517:       pc_gamg->ops->graph(pc,Aarr[level], &Gmat);
518:       pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);
519:       pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);

521:       /* could have failed to create new level */
522:       if (Prol11) {
523:         /* get new block size of coarse matrices */
524:         MatGetBlockSizes(Prol11, NULL, &bs);

526:         if (pc_gamg->ops->optprolongator) {
527:           /* smooth */
528:           pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);
529:         }

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

534:       if (pc_gamg->use_aggs_in_asm) {
535:         PetscInt bs;
536:         MatGetBlockSizes(Prol11, &bs, NULL);
537:         PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);
538:       }

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

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

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

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

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

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

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

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

605:       /* set blocks for ASM smoother that uses the 'aggregates' */
606:       if (pc_gamg->use_aggs_in_asm) {
607:         PetscInt sz;
608:         IS       *iss;

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

659:     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
660:     PetscObjectOptionsBegin((PetscObject)pc);
661:     PCSetFromOptions_MG(PetscOptionsObject,pc);
662:     PetscOptionsEnd();
663:     PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);

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

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

687:    Input Parameter:
688: .  pc - the preconditioner context

690:    Application Interface Routine: PCDestroy()
691: */
692: PetscErrorCode PCDestroy_GAMG(PC pc)
693: {
695:   PC_MG          *mg     = (PC_MG*)pc->data;
696:   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;

699:   PCReset_GAMG(pc);
700:   if (pc_gamg->ops->destroy) {
701:     (*pc_gamg->ops->destroy)(pc);
702:   }
703:   PetscFree(pc_gamg->ops);
704:   PetscFree(pc_gamg->gamg_type_name);
705:   PetscFree(pc_gamg);
706:   PCDestroy_MG(pc);
707:   return(0);
708: }

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

713:    Logically Collective on PC

715:    Input Parameters:
716: +  pc - the preconditioner context
717: -  n - the number of equations


720:    Options Database Key:
721: .  -pc_gamg_process_eq_limit <limit>

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

726:    Level: intermediate

728:    Concepts: Unstructured multigrid preconditioner

730: .seealso: PCGAMGSetCoarseEqLim()
731: @*/
732: PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
733: {

738:   PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));
739:   return(0);
740: }

742: static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
743: {
744:   PC_MG   *mg      = (PC_MG*)pc->data;
745:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

748:   if (n>0) pc_gamg->min_eq_proc = n;
749:   return(0);
750: }

752: /*@
753:    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.

755:  Collective on PC

757:    Input Parameters:
758: +  pc - the preconditioner context
759: -  n - maximum number of equations to aim for

761:    Options Database Key:
762: .  -pc_gamg_coarse_eq_limit <limit>

764:    Level: intermediate

766:    Concepts: Unstructured multigrid preconditioner

768: .seealso: PCGAMGSetProcEqLim()
769: @*/
770: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
771: {

776:   PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
777:   return(0);
778: }

780: static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
781: {
782:   PC_MG   *mg      = (PC_MG*)pc->data;
783:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

786:   if (n>0) pc_gamg->coarse_eq_limit = n;
787:   return(0);
788: }

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

793:    Collective on PC

795:    Input Parameters:
796: +  pc - the preconditioner context
797: -  n - PETSC_TRUE or PETSC_FALSE

799:    Options Database Key:
800: .  -pc_gamg_repartition <true,false>

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

804:    Level: intermediate

806:    Concepts: Unstructured multigrid preconditioner

808: .seealso: ()
809: @*/
810: PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
811: {

816:   PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));
817:   return(0);
818: }

820: static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
821: {
822:   PC_MG   *mg      = (PC_MG*)pc->data;
823:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

826:   pc_gamg->repart = n;
827:   return(0);
828: }

830: /*@
831:    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner

833:    Collective on PC

835:    Input Parameters:
836: +  pc - the preconditioner context
837: -  n - PETSC_TRUE or PETSC_FALSE

839:    Options Database Key:
840: .  -pc_gamg_reuse_interpolation <true,false>

842:    Level: intermediate

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

847:    Concepts: Unstructured multigrid preconditioner

849: .seealso: ()
850: @*/
851: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
852: {

857:   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
858:   return(0);
859: }

861: static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
862: {
863:   PC_MG   *mg      = (PC_MG*)pc->data;
864:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

867:   pc_gamg->reuse_prol = n;
868:   return(0);
869: }

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

874:    Collective on PC

876:    Input Parameters:
877: +  pc - the preconditioner context
878: -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not

880:    Options Database Key:
881: .  -pc_gamg_asm_use_agg

883:    Level: intermediate

885:    Concepts: Unstructured multigrid preconditioner

887: .seealso: ()
888: @*/
889: PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
890: {

895:   PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));
896:   return(0);
897: }

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

905:   pc_gamg->use_aggs_in_asm = flg;
906:   return(0);
907: }

909: /*@
910:    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver

912:    Collective on PC

914:    Input Parameters:
915: +  pc - the preconditioner context
916: -  flg - PETSC_TRUE to not force coarse grid onto one processor

918:    Options Database Key:
919: .  -pc_gamg_use_parallel_coarse_grid_solver

921:    Level: intermediate

923:    Concepts: Unstructured multigrid preconditioner

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

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

937: static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
938: {
939:   PC_MG   *mg      = (PC_MG*)pc->data;
940:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

943:   pc_gamg->use_parallel_coarse_grid_solver = flg;
944:   return(0);
945: }

947: /*@
948:    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use

950:    Not collective on PC

952:    Input Parameters:
953: +  pc - the preconditioner
954: -  n - the maximum number of levels to use

956:    Options Database Key:
957: .  -pc_mg_levels

959:    Level: intermediate

961:    Concepts: Unstructured multigrid preconditioner

963: .seealso: ()
964: @*/
965: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
966: {

971:   PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
972:   return(0);
973: }

975: static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
976: {
977:   PC_MG   *mg      = (PC_MG*)pc->data;
978:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

981:   pc_gamg->Nlevels = n;
982:   return(0);
983: }

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

988:    Not collective on PC

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

994:    Options Database Key:
995: .  -pc_gamg_threshold <threshold>

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

1000:    Level: intermediate

1002:    Concepts: Unstructured multigrid preconditioner

1004: .seealso: PCGAMGFilterGraph()
1005: @*/
1006: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1007: {

1012:   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));
1013:   return(0);
1014: }

1016: static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1017: {
1018:   PC_MG   *mg      = (PC_MG*)pc->data;
1019:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1020:   PetscInt i;
1022:   for (i=0;i<n;i++) pc_gamg->threshold[i] = v[i];
1023:   do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1024:   return(0);
1025: }

1027: /*@
1028:    PCGAMGSetThresholdScale - Relative threshold reduction at each level

1030:    Not collective on PC

1032:    Input Parameters:
1033: +  pc - the preconditioner context
1034: -  scale - the threshold value reduction, ussually < 1.0

1036:    Options Database Key:
1037: .  -pc_gamg_threshold_scale <v>

1039:    Level: advanced

1041:    Concepts: Unstructured multigrid preconditioner

1043: .seealso: ()
1044: @*/
1045: PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1046: {

1051:   PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));
1052:   return(0);
1053: }

1055: static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1056: {
1057:   PC_MG   *mg      = (PC_MG*)pc->data;
1058:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1060:   pc_gamg->threshold_scale = v;
1061:   return(0);
1062: }

1064: /*@C
1065:    PCGAMGSetType - Set solution method

1067:    Collective on PC

1069:    Input Parameters:
1070: +  pc - the preconditioner context
1071: -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL

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

1076:    Level: intermediate

1078:    Concepts: Unstructured multigrid preconditioner

1080: .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1081: @*/
1082: PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1083: {

1088:   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1089:   return(0);
1090: }

1092: /*@C
1093:    PCGAMGGetType - Get solution method

1095:    Collective on PC

1097:    Input Parameter:
1098: .  pc - the preconditioner context

1100:    Output Parameter:
1101: .  type - the type of algorithm used

1103:    Level: intermediate

1105:    Concepts: Unstructured multigrid preconditioner

1107: .seealso: PCGAMGSetType(), PCGAMGType
1108: @*/
1109: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1110: {

1115:   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1116:   return(0);
1117: }

1119: static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1120: {
1121:   PC_MG          *mg      = (PC_MG*)pc->data;
1122:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1125:   *type = pc_gamg->type;
1126:   return(0);
1127: }

1129: static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1130: {
1131:   PetscErrorCode ierr,(*r)(PC);
1132:   PC_MG          *mg      = (PC_MG*)pc->data;
1133:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

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

1157: static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
1158: {
1159:   PetscErrorCode ierr,i;
1160:   PC_MG          *mg      = (PC_MG*)pc->data;
1161:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

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

1183: PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
1184: {
1186:   PC_MG          *mg      = (PC_MG*)pc->data;
1187:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1188:   PetscBool      flag;
1189:   MPI_Comm       comm;
1190:   char           prefix[256];
1191:   PetscInt       i,n;
1192:   const char     *pcpre;

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

1219:     /* set options for subtype */
1220:     if (pc_gamg->ops->setfromoptions) {(*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);}
1221:   }
1222:   PCGetOptionsPrefix(pc, &pcpre);
1223:   PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
1224:   PetscOptionsTail();
1225:   return(0);
1226: }

1228: /* -------------------------------------------------------------------------- */
1229: /*MC
1230:      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner

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

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

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


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

1260:   Level: intermediate

1262:   Concepts: algebraic multigrid

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

1268: PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1269: {
1270:   PetscErrorCode ierr,i;
1271:   PC_GAMG        *pc_gamg;
1272:   PC_MG          *mg;

1275:    /* register AMG type */
1276:   PCGAMGInitializePackage();

1278:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1279:   PCSetType(pc, PCMG);
1280:   PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);

1282:   /* create a supporting struct and attach it to pc */
1283:   PetscNewLog(pc,&pc_gamg);
1284:   PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);
1285:   mg           = (PC_MG*)pc->data;
1286:   mg->innerctx = pc_gamg;

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

1290:   pc_gamg->setup_count = 0;
1291:   /* these should be in subctx but repartitioning needs simple arrays */
1292:   pc_gamg->data_sz = 0;
1293:   pc_gamg->data    = 0;

1295:   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1296:   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1297:   pc->ops->setup          = PCSetUp_GAMG;
1298:   pc->ops->reset          = PCReset_GAMG;
1299:   pc->ops->destroy        = PCDestroy_GAMG;
1300:   mg->view                = PCView_GAMG;

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

1325:   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1326:   PCGAMGSetType(pc,PCGAMGAGG);
1327:   return(0);
1328: }

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

1335:  Level: developer

1337:  .keywords: PC, PCGAMG, initialize, package
1338:  .seealso: PetscInitialize()
1339: @*/
1340: PetscErrorCode PCGAMGInitializePackage(void)
1341: {

1345:   if (PCGAMGPackageInitialized) return(0);
1346:   PCGAMGPackageInitialized = PETSC_TRUE;
1347:   PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);
1348:   PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);
1349:   PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);
1350:   PetscRegisterFinalize(PCGAMGFinalizePackage);

1352:   /* general events */
1353:   PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);
1354:   PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);
1355:   PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1356:   PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1357:   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1358:   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1359:   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);

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

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

1400: /*@C
1401:  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1402:     called from PetscFinalize() automatically.

1404:  Level: developer

1406:  .keywords: Petsc, destroy, package
1407:  .seealso: PetscFinalize()
1408: @*/
1409: PetscErrorCode PCGAMGFinalizePackage(void)
1410: {

1414:   PCGAMGPackageInitialized = PETSC_FALSE;
1415:   PetscFunctionListDestroy(&GAMGList);
1416:   return(0);
1417: }

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

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

1426:   Level: advanced

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

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