Actual source code: gamg.c

petsc-master 2016-12-08
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: /* ----------------------------------------------------------------------------- */
 33: PetscErrorCode PCReset_GAMG(PC pc)
 34: {
 36:   PC_MG          *mg      = (PC_MG*)pc->data;
 37:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

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

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

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

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

 76:   PetscObjectGetComm((PetscObject)Amat_fine,&comm);
 77:   MPI_Comm_rank(comm, &rank);
 78:   MPI_Comm_size(comm, &size);
 79:   MatGetBlockSize(Amat_fine, &f_bs);
 80:   MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);

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

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

102:   if (Pcolumnperm) *Pcolumnperm = NULL;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

615:         sz   = nASMBlocksArr[level];
616:         iss   = ASMLocalIDsArr[level];
617:         PCSetType(subpc, PCASM);
618:         PCASMSetOverlap(subpc, 0);
619:         PCASMSetType(subpc,PC_ASM_BASIC);
620:         if (!sz) {
621:           IS       is;
622:           ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);
623:           PCASMSetLocalSubdomains(subpc, 1, NULL, &is);
624:           ISDestroy(&is);
625:         } else {
626:           PetscInt kk;
627:           PCASMSetLocalSubdomains(subpc, sz, NULL, iss);
628:           for (kk=0; kk<sz; kk++) {
629:             ISDestroy(&iss[kk]);
630:           }
631:           PetscFree(iss);
632:         }
633:         ASMLocalIDsArr[level] = NULL;
634:         nASMBlocksArr[level]  = 0;
635:       } else {
636:         PCSetType(subpc, PCSOR);
637:       }
638:     }
639:     {
640:       /* coarse grid */
641:       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
642:       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
643:       PCMGGetSmoother(pc, lidx, &smoother);
644:       KSPSetOperators(smoother, Lmat, Lmat);
645:       if (!pc_gamg->use_parallel_coarse_grid_solver) {
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:       }
662:     }

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

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:   PetscFree(pc_gamg->ops);
711:   PetscFree(pc_gamg->gamg_type_name);
712:   PetscFree(pc_gamg);
713:   PCDestroy_MG(pc);
714:   return(0);
715: }

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

722:    Logically Collective on PC

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


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

732:    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 
733:           that has degrees of freedom

735:    Level: intermediate

737:    Concepts: Unstructured multigrid preconditioner

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

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

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

759:   if (n>0) pc_gamg->min_eq_proc = n;
760:   return(0);
761: }

765: /*@
766:    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.

768:  Collective on PC

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

774:    Options Database Key:
775: .  -pc_gamg_coarse_eq_limit <limit>

777:    Level: intermediate

779:    Concepts: Unstructured multigrid preconditioner

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

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

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

801:   if (n>0) pc_gamg->coarse_eq_limit = n;
802:   return(0);
803: }

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

810:    Collective on PC

812:    Input Parameters:
813: +  pc - the preconditioner context
814: -  n - PETSC_TRUE or PETSC_FALSE

816:    Options Database Key:
817: .  -pc_gamg_repartition <true,false>

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

821:    Level: intermediate

823:    Concepts: Unstructured multigrid preconditioner

825: .seealso: ()
826: @*/
827: PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
828: {

833:   PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));
834:   return(0);
835: }

839: static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
840: {
841:   PC_MG   *mg      = (PC_MG*)pc->data;
842:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

845:   pc_gamg->repart = n;
846:   return(0);
847: }

851: /*@
852:    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner

854:    Collective on PC

856:    Input Parameters:
857: +  pc - the preconditioner context
858: -  n - PETSC_TRUE or PETSC_FALSE

860:    Options Database Key:
861: .  -pc_gamg_reuse_interpolation <true,false>

863:    Level: intermediate

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

868:    Concepts: Unstructured multigrid preconditioner

870: .seealso: ()
871: @*/
872: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
873: {

878:   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
879:   return(0);
880: }

884: static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
885: {
886:   PC_MG   *mg      = (PC_MG*)pc->data;
887:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

890:   pc_gamg->reuse_prol = n;
891:   return(0);
892: }

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

899:    Collective on PC

901:    Input Parameters:
902: +  pc - the preconditioner context
903: -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not

905:    Options Database Key:
906: .  -pc_gamg_asm_use_agg

908:    Level: intermediate

910:    Concepts: Unstructured multigrid preconditioner

912: .seealso: ()
913: @*/
914: PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
915: {

920:   PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));
921:   return(0);
922: }

926: static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
927: {
928:   PC_MG   *mg      = (PC_MG*)pc->data;
929:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

932:   pc_gamg->use_aggs_in_asm = flg;
933:   return(0);
934: }

938: /*@
939:    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver

941:    Collective on PC

943:    Input Parameters:
944: +  pc - the preconditioner context
945: -  flg - PETSC_TRUE to not force coarse grid onto one processor

947:    Options Database Key:
948: .  -pc_gamg_use_parallel_coarse_grid_solver

950:    Level: intermediate

952:    Concepts: Unstructured multigrid preconditioner

954: .seealso: ()
955: @*/
956: PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
957: {

962:   PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));
963:   return(0);
964: }

968: static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
969: {
970:   PC_MG   *mg      = (PC_MG*)pc->data;
971:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

974:   pc_gamg->use_parallel_coarse_grid_solver = flg;
975:   return(0);
976: }

980: /*@
981:    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use

983:    Not collective on PC

985:    Input Parameters:
986: +  pc - the preconditioner
987: -  n - the maximum number of levels to use

989:    Options Database Key:
990: .  -pc_mg_levels

992:    Level: intermediate

994:    Concepts: Unstructured multigrid preconditioner

996: .seealso: ()
997: @*/
998: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
999: {

1004:   PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
1005:   return(0);
1006: }

1010: static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1011: {
1012:   PC_MG   *mg      = (PC_MG*)pc->data;
1013:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1016:   pc_gamg->Nlevels = n;
1017:   return(0);
1018: }

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

1025:    Not collective on PC

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

1031:    Options Database Key:
1032: .  -pc_gamg_threshold <threshold>

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

1037:    Level: intermediate

1039:    Concepts: Unstructured multigrid preconditioner

1041: .seealso: ()
1042: @*/
1043: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1044: {

1049:   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));
1050:   return(0);
1051: }

1055: static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1056: {
1057:   PC_MG   *mg      = (PC_MG*)pc->data;
1058:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1059:   PetscInt i;
1061:   for (i=0;i<n;i++) pc_gamg->threshold[i] = v[i];
1062:   do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1063:   return(0);
1064: }

1068: /*@
1069:    PCGAMGSetThresholdScale - Relative threshold reduction at each level

1071:    Not collective on PC

1073:    Input Parameters:
1074: +  pc - the preconditioner context
1075: -  scale - the threshold value reduction, ussually < 1.0

1077:    Options Database Key:
1078: .  -pc_gamg_threshold_scale <v>

1080:    Level: advanced

1082:    Concepts: Unstructured multigrid preconditioner

1084: .seealso: ()
1085: @*/
1086: PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1087: {

1092:   PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));
1093:   return(0);
1094: }

1098: static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1099: {
1100:   PC_MG   *mg      = (PC_MG*)pc->data;
1101:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1103:   pc_gamg->threshold_scale = v;
1104:   return(0);
1105: }

1109: /*@C
1110:    PCGAMGSetType - Set solution method

1112:    Collective on PC

1114:    Input Parameters:
1115: +  pc - the preconditioner context
1116: -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL

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

1121:    Level: intermediate

1123:    Concepts: Unstructured multigrid preconditioner

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

1133:   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1134:   return(0);
1135: }

1139: /*@C
1140:    PCGAMGGetType - Get solution method

1142:    Collective on PC

1144:    Input Parameter:
1145: .  pc - the preconditioner context

1147:    Output Parameter:
1148: .  type - the type of algorithm used

1150:    Level: intermediate

1152:    Concepts: Unstructured multigrid preconditioner

1154: .seealso: PCGAMGSetType(), PCGAMGType
1155: @*/
1156: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1157: {

1162:   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1163:   return(0);
1164: }

1168: static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1169: {
1170:   PC_MG          *mg      = (PC_MG*)pc->data;
1171:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1174:   *type = pc_gamg->type;
1175:   return(0);
1176: }

1180: static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1181: {
1182:   PetscErrorCode ierr,(*r)(PC);
1183:   PC_MG          *mg      = (PC_MG*)pc->data;
1184:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1187:   pc_gamg->type = type;
1188:   PetscFunctionListFind(GAMGList,type,&r);
1189:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1190:   if (pc_gamg->ops->destroy) {
1191:     (*pc_gamg->ops->destroy)(pc);
1192:     PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));
1193:     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1194:     /* cleaning up common data in pc_gamg - this should disapear someday */
1195:     pc_gamg->data_cell_cols = 0;
1196:     pc_gamg->data_cell_rows = 0;
1197:     pc_gamg->orig_data_cell_cols = 0;
1198:     pc_gamg->orig_data_cell_rows = 0;
1199:     PetscFree(pc_gamg->data);
1200:     pc_gamg->data_sz = 0;
1201:   }
1202:   PetscFree(pc_gamg->gamg_type_name);
1203:   PetscStrallocpy(type,&pc_gamg->gamg_type_name);
1204:   (*r)(pc);
1205:   return(0);
1206: }

1210: static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
1211: {
1212:   PetscErrorCode ierr,i;
1213:   PC_MG          *mg      = (PC_MG*)pc->data;
1214:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1217:   PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");
1218:   PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");
1219:   for (i=0;i<pc_gamg->current_level;i++) {
1220:     PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);
1221:   }
1222:   PetscViewerASCIIPrintf(viewer,"\n");
1223:   PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);
1224:   if (pc_gamg->use_aggs_in_asm) {
1225:     PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");
1226:   }
1227:   if (pc_gamg->use_parallel_coarse_grid_solver) {
1228:     PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");
1229:   }
1230:   if (pc_gamg->ops->view) {
1231:     (*pc_gamg->ops->view)(pc,viewer);
1232:   }
1233:   return(0);
1234: }

1238: PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
1239: {
1241:   PC_MG          *mg      = (PC_MG*)pc->data;
1242:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1243:   PetscBool      flag;
1244:   MPI_Comm       comm;
1245:   char           prefix[256];
1246:   PetscInt       i,n;
1247:   const char     *pcpre;

1250:   PetscObjectGetComm((PetscObject)pc,&comm);
1251:   PetscOptionsHead(PetscOptionsObject,"GAMG options");
1252:   {
1253:     char tname[256];
1254:     PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);
1255:     if (flag) {
1256:       PCGAMGSetType(pc,tname);
1257:     }
1258:     PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);
1259:     PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);
1260:     PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation agragates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL);
1261:     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);
1262:     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);
1263:     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);
1264:     PetscOptionsReal("-pc_gamg_threshold_scale","scaleing of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL);
1265:     n = PETSC_GAMG_MAXLEVELS;
1266:     PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);
1267:     if (!flag || n < PETSC_GAMG_MAXLEVELS) {
1268:       if (!flag) n = 1;
1269:       i = n;
1270:       do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_GAMG_MAXLEVELS);
1271:     }
1272:     PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);

1274:     /* set options for subtype */
1275:     if (pc_gamg->ops->setfromoptions) {(*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);}
1276:   }
1277:   PCGetOptionsPrefix(pc, &pcpre);
1278:   PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
1279:   PetscOptionsTail();
1280:   return(0);
1281: }

1283: /* -------------------------------------------------------------------------- */
1284: /*MC
1285:      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner

1287:    Options Database Keys:
1288: +   -pc_gamg_type <type> - one of agg, geo, or classical
1289: .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1290: .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1291: .   -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
1292: .   -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>
1293:                                         equations on each process that has degrees of freedom
1294: .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
1295: -   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1296: -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified

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

1303:    Multigrid options(inherited):
1304: +  -pc_mg_cycles <v>: v or w (PCMGSetCycleType())
1305: .  -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp)
1306: .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown)
1307: .  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
1308: -  -pc_mg_levels <levels> - Number of levels of multigrid to use.


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

1316:   Level: intermediate

1318:   Concepts: algebraic multigrid

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

1326: PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1327: {
1328:   PetscErrorCode ierr,i;
1329:   PC_GAMG        *pc_gamg;
1330:   PC_MG          *mg;

1333:    /* register AMG type */
1334:   PCGAMGInitializePackage();

1336:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1337:   PCSetType(pc, PCMG);
1338:   PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);

1340:   /* create a supporting struct and attach it to pc */
1341:   PetscNewLog(pc,&pc_gamg);
1342:   PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);
1343:   mg           = (PC_MG*)pc->data;
1344:   mg->innerctx = pc_gamg;

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

1348:   pc_gamg->setup_count = 0;
1349:   /* these should be in subctx but repartitioning needs simple arrays */
1350:   pc_gamg->data_sz = 0;
1351:   pc_gamg->data    = 0;

1353:   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1354:   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1355:   pc->ops->setup          = PCSetUp_GAMG;
1356:   pc->ops->reset          = PCReset_GAMG;
1357:   pc->ops->destroy        = PCDestroy_GAMG;
1358:   mg->view                = PCView_GAMG;

1360:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);
1361:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);
1362:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);
1363:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);
1364:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);
1365:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);
1366:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);
1367:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);
1368:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);
1369:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);
1370:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);
1371:   pc_gamg->repart           = PETSC_FALSE;
1372:   pc_gamg->reuse_prol       = PETSC_FALSE;
1373:   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1374:   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1375:   pc_gamg->min_eq_proc      = 50;
1376:   pc_gamg->coarse_eq_limit  = 50;
1377:   for (i=0;i<PETSC_GAMG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1378:   pc_gamg->threshold_scale = 1.;
1379:   pc_gamg->Nlevels          = PETSC_GAMG_MAXLEVELS;
1380:   pc_gamg->current_level    = 0; /* don't need to init really */
1381:   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;

1383:   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1384:   PCGAMGSetType(pc,PCGAMGAGG);
1385:   return(0);
1386: }

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

1395:  Level: developer

1397:  .keywords: PC, PCGAMG, initialize, package
1398:  .seealso: PetscInitialize()
1399: @*/
1400: PetscErrorCode PCGAMGInitializePackage(void)
1401: {

1405:   if (PCGAMGPackageInitialized) return(0);
1406:   PCGAMGPackageInitialized = PETSC_TRUE;
1407:   PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);
1408:   PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);
1409:   PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);
1410:   PetscRegisterFinalize(PCGAMGFinalizePackage);

1412:   /* general events */
1413:   PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);
1414:   PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);
1415:   PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1416:   PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1417:   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1418:   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1419:   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);

1421: #if defined PETSC_GAMG_USE_LOG
1422:   PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1423:   PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1424:   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1425:   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1426:   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1427:   PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1428:   PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1429:   PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1430:   PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1431:   PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1432:   PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1433:   PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1434:   PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1435:   PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1436:   PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1437:   PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1438:   PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);

1440:   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1441:   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1442:   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1443:   /* create timer stages */
1444: #if defined GAMG_STAGES
1445:   {
1446:     char     str[32];
1447:     PetscInt lidx;
1448:     sprintf(str,"MG Level %d (finest)",0);
1449:     PetscLogStageRegister(str, &gamg_stages[0]);
1450:     for (lidx=1; lidx<9; lidx++) {
1451:       sprintf(str,"MG Level %d",lidx);
1452:       PetscLogStageRegister(str, &gamg_stages[lidx]);
1453:     }
1454:   }
1455: #endif
1456: #endif
1457:   return(0);
1458: }

1462: /*@C
1463:  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1464:     called from PetscFinalize() automatically.

1466:  Level: developer

1468:  .keywords: Petsc, destroy, package
1469:  .seealso: PetscFinalize()
1470: @*/
1471: PetscErrorCode PCGAMGFinalizePackage(void)
1472: {

1476:   PCGAMGPackageInitialized = PETSC_FALSE;
1477:   PetscFunctionListDestroy(&GAMGList);
1478:   return(0);
1479: }

1483: /*@C
1484:  PCGAMGRegister - Register a PCGAMG implementation.

1486:  Input Parameters:
1487:  + type - string that will be used as the name of the GAMG type.
1488:  - create - function for creating the gamg context.

1490:   Level: advanced

1492:  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1493: @*/
1494: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1495: {

1499:   PCGAMGInitializePackage();
1500:   PetscFunctionListAdd(&GAMGList,type,create);
1501:   return(0);
1502: }