Actual source code: gamg.c

petsc-3.14.6 2021-03-30
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>
  7: #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h>

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

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

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

 28: static PetscFunctionList GAMGList = NULL;
 29: static PetscBool PCGAMGPackageInitialized;

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

 39:   PetscFree(pc_gamg->data);
 40:   pc_gamg->data_sz = 0;
 41:   PetscFree(pc_gamg->orig_data);
 42:   for (level = 0; level < PETSC_MG_MAXLEVELS ; level++) {
 43:     mg->min_eigen_DinvA[level] = 0;
 44:     mg->max_eigen_DinvA[level] = 0;
 45:   }
 46:   pc_gamg->emin = 0;
 47:   pc_gamg->emax = 0;
 48:   return(0);
 49: }

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

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

 68: static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,IS * Pcolumnperm, PetscBool is_last)
 69: {
 70:   PetscErrorCode  ierr;
 71:   PC_MG           *mg         = (PC_MG*)pc->data;
 72:   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
 73:   Mat             Cmat,Pold=*a_P_inout;
 74:   MPI_Comm        comm;
 75:   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
 76:   PetscInt        ncrs_eq,ncrs,f_bs;

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

 85:   if (Pcolumnperm) *Pcolumnperm = NULL;

 87:   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
 88:   MatGetLocalSize(Cmat, &ncrs_eq, NULL);
 89:   if (pc_gamg->data_cell_rows>0) {
 90:     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
 91:   } else {
 92:     PetscInt  bs;
 93:     MatGetBlockSize(Cmat, &bs);
 94:     ncrs = ncrs_eq/bs;
 95:   }
 96:   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
 97:   if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) new_size = 1;
 98:   else {
 99:     PetscInt ncrs_eq_glob;
100:     MatGetSize(Cmat, &ncrs_eq_glob, NULL);
101:     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
102:     if (!new_size) new_size = 1; /* not likely, posible? */
103:     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
104:   }

106:   if (new_size==nactive) {
107:     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
108:     if (new_size < size) {
109:       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
110:       PetscInfo1(pc,"reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",nactive);
111:       if (pc_gamg->cpu_pin_coarse_grids) {
112:         MatBindToCPU(*a_Amat_crs,PETSC_TRUE);
113:         MatBindToCPU(*a_P_inout,PETSC_TRUE);
114:       }
115:     }
116:     /* we know that the grid structure can be reused in MatPtAP */
117:   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
118:     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1;
119:     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
120:     nloc_old = ncrs_eq/cr_bs;
121:     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);
122:     /* get new_size and rfactor */
123:     if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
124:       /* find factor */
125:       if (new_size == 1) rfactor = size; /* don't modify */
126:       else {
127:         PetscReal best_fact = 0.;
128:         jj = -1;
129:         for (kk = 1 ; kk <= size ; kk++) {
130:           if (!(size%kk)) { /* a candidate */
131:             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
132:             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
133:             if (fact > best_fact) {
134:               best_fact = fact; jj = kk;
135:             }
136:           }
137:         }
138:         if (jj != -1) rfactor = jj;
139:         else rfactor = 1; /* a prime */
140:         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
141:         else expand_factor = rfactor;
142:       }
143:       new_size = size/rfactor; /* make new size one that is factor */
144:       if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
145:         *a_Amat_crs = Cmat;
146:         PetscInfo2(pc,"Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",new_size,ncrs_eq);
147:         return(0);
148:       }
149:     }
150: #if defined PETSC_GAMG_USE_LOG
151:     PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);
152: #endif
153:     /* make 'is_eq_newproc' */
154:     PetscMalloc1(size, &counts);
155:     if (pc_gamg->repart) {
156:       /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */
157:       Mat      adj;
158:       PetscInfo4(pc,"Repartition: size (active): %d --> %d, %D local equations, using %s process layout\n",*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread");
159:       /* get 'adj' */
160:       if (cr_bs == 1) {
161:         MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);
162:       } else {
163:         /* make a scalar matrix to partition (no Stokes here) */
164:         Mat               tMat;
165:         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
166:         const PetscScalar *vals;
167:         const PetscInt    *idx;
168:         PetscInt          *d_nnz, *o_nnz, M, N;
169:         static PetscInt   llev = 0; /* ugly but just used for debugging */
170:         MatType           mtype;

172:         PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);
173:         MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);
174:         MatGetSize(Cmat, &M, &N);
175:         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
176:           MatGetRow(Cmat,Ii,&ncols,NULL,NULL);
177:           d_nnz[jj] = ncols/cr_bs;
178:           o_nnz[jj] = ncols/cr_bs;
179:           MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL);
180:           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
181:           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
182:         }

184:         MatGetType(Amat_fine,&mtype);
185:         MatCreate(comm, &tMat);
186:         MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);
187:         MatSetType(tMat,mtype);
188:         MatSeqAIJSetPreallocation(tMat,0,d_nnz);
189:         MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);
190:         PetscFree2(d_nnz,o_nnz);

192:         for (ii = Istart_crs; ii < Iend_crs; ii++) {
193:           PetscInt dest_row = ii/cr_bs;
194:           MatGetRow(Cmat,ii,&ncols,&idx,&vals);
195:           for (jj = 0; jj < ncols; jj++) {
196:             PetscInt    dest_col = idx[jj]/cr_bs;
197:             PetscScalar v        = 1.0;
198:             MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);
199:           }
200:           MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);
201:         }
202:         MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);
203:         MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);

205:         if (llev++ == -1) {
206:           PetscViewer viewer; char fname[32];
207:           PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);
208:           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
209:           MatView(tMat, viewer);
210:           PetscViewerDestroy(&viewer);
211:         }
212:         MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);
213:         MatDestroy(&tMat);
214:       } /* create 'adj' */

216:       { /* partition: get newproc_idx */
217:         char            prefix[256];
218:         const char      *pcpre;
219:         const PetscInt  *is_idx;
220:         MatPartitioning mpart;
221:         IS              proc_is;

223:         MatPartitioningCreate(comm, &mpart);
224:         MatPartitioningSetAdjacency(mpart, adj);
225:         PCGetOptionsPrefix(pc, &pcpre);
226:         PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
227:         PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);
228:         MatPartitioningSetFromOptions(mpart);
229:         MatPartitioningSetNParts(mpart, new_size);
230:         MatPartitioningApply(mpart, &proc_is);
231:         MatPartitioningDestroy(&mpart);

233:         /* collect IS info */
234:         PetscMalloc1(ncrs_eq, &newproc_idx);
235:         ISGetIndices(proc_is, &is_idx);
236:         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
237:           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
238:             newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
239:           }
240:         }
241:         ISRestoreIndices(proc_is, &is_idx);
242:         ISDestroy(&proc_is);
243:       }
244:       MatDestroy(&adj);

246:       ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);
247:       PetscFree(newproc_idx);
248:     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
249:       PetscInt targetPE;
250:       if (new_size==nactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen");
251:       PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);
252:       targetPE = (rank/rfactor)*expand_factor;
253:       ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);
254:     } /* end simple 'is_eq_newproc' */

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

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

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

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

402:     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
403:     if (pc_gamg->cpu_pin_coarse_grids) {
404: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
405:       static PetscInt llev = 2;
406:       PetscInfo1(pc,"Pinning level %D to the CPU\n",llev++);
407: #endif
408:       MatBindToCPU(*a_Amat_crs,PETSC_TRUE);
409:       MatBindToCPU(*a_P_inout,PETSC_TRUE);
410:       if (1) { /* lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
411:         Mat         A = *a_Amat_crs, P = *a_P_inout;
412:         PetscMPIInt size;
413:         MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
414:         if (size > 1) {
415:           Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data;
416:           VecBindToCPU(a->lvec,PETSC_TRUE);
417:           VecBindToCPU(p->lvec,PETSC_TRUE);
418:         }
419:       }
420:     }
421: #if defined PETSC_GAMG_USE_LOG
422:     PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);
423: #endif
424:   }
425:   return(0);
426: }

428: PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2)
429: {
431:   const char     *prefix;
432:   char           addp[32];
433:   PC_MG          *mg      = (PC_MG*)a_pc->data;
434:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

437:   PCGetOptionsPrefix(a_pc,&prefix);
438:   PetscInfo1(a_pc,"Square Graph on level %D\n",pc_gamg->current_level+1);
439:   MatProductCreate(Gmat1,Gmat1,NULL,Gmat2);
440:   MatSetOptionsPrefix(*Gmat2,prefix);
441:   PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level);
442:   MatAppendOptionsPrefix(*Gmat2,addp);
443:   MatProductSetType(*Gmat2,MATPRODUCT_AtB);
444:   MatProductSetFromOptions(*Gmat2);
445:   MatProductSymbolic(*Gmat2);
446:   /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
447:   (*Gmat2)->assembled = PETSC_TRUE;
448:   return(0);
449: }

451: /* -------------------------------------------------------------------------- */
452: /*
453:    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
454:                     by setting data structures and options.

456:    Input Parameter:
457: .  pc - the preconditioner context

459: */
460: PetscErrorCode PCSetUp_GAMG(PC pc)
461: {
463:   PC_MG          *mg      = (PC_MG*)pc->data;
464:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
465:   Mat            Pmat     = pc->pmat;
466:   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS];
467:   MPI_Comm       comm;
468:   PetscMPIInt    rank,size,nactivepe;
469:   Mat            Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS];
470:   IS             *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
471:   PetscLogDouble nnz0=0.,nnztot=0.;
472:   MatInfo        info;
473:   PetscBool      is_last = PETSC_FALSE;

476:   PetscObjectGetComm((PetscObject)pc,&comm);
477:   MPI_Comm_rank(comm,&rank);
478:   MPI_Comm_size(comm,&size);

480:   if (pc_gamg->setup_count++ > 0) {
481:     if ((PetscBool)(!pc_gamg->reuse_prol)) {
482:       /* reset everything */
483:       PCReset_MG(pc);
484:       pc->setupcalled = 0;
485:     } else {
486:       PC_MG_Levels **mglevels = mg->levels;
487:       /* just do Galerkin grids */
488:       Mat          B,dA,dB;

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

497:         for (level=pc_gamg->Nlevels-2; level>=0; level--) {
498:           /* 2nd solve, matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
499:           if (pc_gamg->setup_count==2 /* && pc_gamg->repart||reduction */) {
500:             PetscInfo2(pc,"new RAP after first solve level %D, %D setup\n",level,pc_gamg->setup_count);
501:             MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,2.0,&B);
502:             MatDestroy(&mglevels[level]->A);
503:             mglevels[level]->A = B;
504:           } else {
505:             PetscInfo2(pc,"RAP after first solve reusing matrix level %D, %D setup\n",level,pc_gamg->setup_count);
506:             KSPGetOperators(mglevels[level]->smoothd,NULL,&B);
507:             MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);
508:           }
509:           KSPSetOperators(mglevels[level]->smoothd,B,B);
510:           dB   = B;
511:         }
512:       }

514:       PCSetUp_MG(pc);
515:       return(0);
516:     }
517:   }

519:   if (!pc_gamg->data) {
520:     if (pc_gamg->orig_data) {
521:       MatGetBlockSize(Pmat, &bs);
522:       MatGetLocalSize(Pmat, &qq, NULL);

524:       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
525:       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
526:       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;

528:       PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);
529:       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
530:     } else {
531:       if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
532:       pc_gamg->ops->createdefaultdata(pc,Pmat);
533:     }
534:   }

536:   /* cache original data for reuse */
537:   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
538:     PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);
539:     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
540:     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
541:     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
542:   }

544:   /* get basic dims */
545:   MatGetBlockSize(Pmat, &bs);
546:   MatGetSize(Pmat, &M, &N);

548:   MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info); /* global reduction */
549:   nnz0   = info.nz_used;
550:   nnztot = info.nz_used;
551:   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);

553:   /* Get A_i and R_i */
554:   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
555:     pc_gamg->current_level = level;
556:     if (level >= PETSC_MG_MAXLEVELS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level);
557:     level1 = level + 1;
558: #if defined PETSC_GAMG_USE_LOG
559:     PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);
560: #if (defined GAMG_STAGES)
561:     PetscLogStagePush(gamg_stages[level]);
562: #endif
563: #endif
564:     { /* construct prolongator */
565:       Mat              Gmat;
566:       PetscCoarsenData *agg_lists;
567:       Mat              Prol11;

569:       pc_gamg->ops->graph(pc,Aarr[level], &Gmat);
570:       pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);
571:       pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);

573:       /* could have failed to create new level */
574:       if (Prol11) {
575:         const char *prefix;
576:         char       addp[32];

578:         /* get new block size of coarse matrices */
579:         MatGetBlockSizes(Prol11, NULL, &bs);

581:         if (pc_gamg->ops->optprolongator) {
582:           /* smooth */
583:           pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);
584:         }

586:         if (pc_gamg->use_aggs_in_asm) {
587:           PetscInt bs;
588:           MatGetBlockSizes(Prol11, &bs, NULL);
589:           PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);
590:         }

592:         PCGetOptionsPrefix(pc,&prefix);
593:         MatSetOptionsPrefix(Prol11,prefix);
594:         PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",level);
595:         MatAppendOptionsPrefix(Prol11,addp);
596:         /* if we reuse the prolongator, then it is better to generate the transpose by default with CUDA
597:            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
598: #if defined(PETSC_HAVE_CUDA)
599:         {
600:           PetscBool ismpiaij;

602:           PetscObjectBaseTypeCompare((PetscObject)Prol11,MATMPIAIJ,&ismpiaij);
603:           if (ismpiaij) {
604:             Mat Prol_d,Prol_o;

606:             MatMPIAIJGetSeqAIJ(Prol11,&Prol_d,&Prol_o,NULL);
607:             MatSeqAIJCUSPARSESetGenerateTranspose(Prol_d,pc_gamg->reuse_prol);
608:             MatSeqAIJCUSPARSESetGenerateTranspose(Prol_o,pc_gamg->reuse_prol);
609:           } else {
610:             MatSeqAIJCUSPARSESetGenerateTranspose(Prol11,pc_gamg->reuse_prol);
611:           }
612:         }
613: #endif
614:         MatSetFromOptions(Prol11);
615:         Parr[level1] = Prol11;
616:       } else Parr[level1] = NULL; /* failed to coarsen */

618:       MatDestroy(&Gmat);
619:       PetscCDDestroy(agg_lists);
620:     } /* construct prolongator scope */
621: #if defined PETSC_GAMG_USE_LOG
622:     PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);
623: #endif
624:     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
625:     if (!Parr[level1]) { /* failed to coarsen */
626:        PetscInfo1(pc,"Stop gridding, level %D\n",level);
627: #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES
628:       PetscLogStagePop();
629: #endif
630:       break;
631:     }
632: #if defined PETSC_GAMG_USE_LOG
633:     PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);
634: #endif
635:     MatGetSize(Parr[level1], &M, &N); /* N is next M, a loop test variables */
636:     if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ????????");
637:     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
638:     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
639:     pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);

641: #if defined PETSC_GAMG_USE_LOG
642:     PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);
643: #endif
644:     MatGetSize(Aarr[level1], &M, &N); /* M is loop test variables */
645:     MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);
646:     nnztot += info.nz_used;
647:     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);

649: #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES)
650:     PetscLogStagePop();
651: #endif
652:     /* stop if one node or one proc -- could pull back for singular problems */
653:     if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) {
654:        PetscInfo2(pc,"HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",level,M/bs);
655:       level++;
656:       break;
657:     }
658:   } /* levels */
659:   PetscFree(pc_gamg->data);

661:   PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);
662:   pc_gamg->Nlevels = level + 1;
663:   fine_level       = level;
664:   PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);

666:   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
667:     PetscErrorCode (*savesetfromoptions[PETSC_MG_MAXLEVELS])(PetscOptionItems*,KSP);

669:     /* set default smoothers & set operators */
670:     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
671:       KSP smoother;
672:       PC  subpc;

674:       PCMGGetSmoother(pc, lidx, &smoother);
675:       KSPGetPC(smoother, &subpc);

677:       KSPSetNormType(smoother, KSP_NORM_NONE);
678:       /* set ops */
679:       KSPSetOperators(smoother, Aarr[level], Aarr[level]);
680:       PCMGSetInterpolation(pc, lidx, Parr[level+1]);

682:       /* set defaults */
683:       KSPSetType(smoother, KSPCHEBYSHEV);

685:       /* set blocks for ASM smoother that uses the 'aggregates' */
686:       if (pc_gamg->use_aggs_in_asm) {
687:         PetscInt sz;
688:         IS       *iss;

690:         sz   = nASMBlocksArr[level];
691:         iss  = ASMLocalIDsArr[level];
692:         PCSetType(subpc, PCASM);
693:         PCASMSetOverlap(subpc, 0);
694:         PCASMSetType(subpc,PC_ASM_BASIC);
695:         if (!sz) {
696:           IS       is;
697:           ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);
698:           PCASMSetLocalSubdomains(subpc, 1, NULL, &is);
699:           ISDestroy(&is);
700:         } else {
701:           PetscInt kk;
702:           PCASMSetLocalSubdomains(subpc, sz, NULL, iss);
703:           for (kk=0; kk<sz; kk++) {
704:             ISDestroy(&iss[kk]);
705:           }
706:           PetscFree(iss);
707:         }
708:         ASMLocalIDsArr[level] = NULL;
709:         nASMBlocksArr[level]  = 0;
710:       } else {
711:         PCSetType(subpc, PCSOR);
712:       }
713:     }
714:     {
715:       /* coarse grid */
716:       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
717:       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;

719:       PCMGGetSmoother(pc, lidx, &smoother);
720:       KSPSetOperators(smoother, Lmat, Lmat);
721:       if (!pc_gamg->use_parallel_coarse_grid_solver) {
722:         KSPSetNormType(smoother, KSP_NORM_NONE);
723:         KSPGetPC(smoother, &subpc);
724:         PCSetType(subpc, PCBJACOBI);
725:         PCSetUp(subpc);
726:         PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);
727:         if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
728:         KSPGetPC(k2[0],&pc2);
729:         PCSetType(pc2, PCLU);
730:         PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);
731:         KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);
732:         KSPSetType(k2[0], KSPPREONLY);
733:         /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in
734:          * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to
735:          * view every subdomain as though they were different. */
736:         ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE;
737:       }
738:     }

740:     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
741:     PetscObjectOptionsBegin((PetscObject)pc);
742:     PCSetFromOptions_MG(PetscOptionsObject,pc);
743:     PetscOptionsEnd();
744:     PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);

746:     /* setup cheby eigen estimates from SA */
747:     if (pc_gamg->use_sa_esteig==1) {
748:       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
749:         KSP       smoother;
750:         PetscBool ischeb;

752:         savesetfromoptions[level] = NULL;
753:         PCMGGetSmoother(pc, lidx, &smoother);
754:         PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb);
755:         if (ischeb) {
756:           KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data;

758:           KSPSetFromOptions(smoother); /* let command line emax override using SA's eigenvalues */
759:           if (mg->max_eigen_DinvA[level] > 0 && cheb->emax == 0.) {
760:             PC        subpc;
761:             PetscBool isjac;
762:             KSPGetPC(smoother, &subpc);
763:             PetscObjectTypeCompare((PetscObject)subpc,PCJACOBI,&isjac);
764:             if (isjac && pc_gamg->use_sa_esteig==1) {
765:               PetscReal emax,emin;

767:               emin = mg->min_eigen_DinvA[level];
768:               emax = mg->max_eigen_DinvA[level];
769:               PetscInfo4(pc,"PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %D (N=%D) with emax = %g emin = %g\n",level,Aarr[level]->rmap->N,(double)emax,(double)emin);
770:               cheb->emin_computed = emin;
771:               cheb->emax_computed = emax;
772:               KSPChebyshevSetEigenvalues(smoother, cheb->tform[2]*emin + cheb->tform[3]*emax, cheb->tform[0]*emin + cheb->tform[1]*emax);

774:               /* We have set the eigenvalues and consumed the transformation values
775:                  prevent from flagging the recomputation of the eigenvalues again in PCSetUp_MG
776:                  below when setfromoptions will be called again */
777:               savesetfromoptions[level] = smoother->ops->setfromoptions;
778:               smoother->ops->setfromoptions = NULL;
779:             }
780:           }
781:         }
782:       }
783:     }

785:     PCSetUp_MG(pc);

787:     /* restore Chebyshev smoother for next calls */
788:     if (pc_gamg->use_sa_esteig==1) {
789:       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
790:         if (savesetfromoptions[level]) {
791:           KSP smoother;
792:           PCMGGetSmoother(pc, lidx, &smoother);
793:           smoother->ops->setfromoptions = savesetfromoptions[level];
794:         }
795:       }
796:     }

798:     /* clean up */
799:     for (level=1; level<pc_gamg->Nlevels; level++) {
800:       MatDestroy(&Parr[level]);
801:       MatDestroy(&Aarr[level]);
802:     }
803:   } else {
804:     KSP smoother;

806:     PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");
807:     PCMGGetSmoother(pc, 0, &smoother);
808:     KSPSetOperators(smoother, Aarr[0], Aarr[0]);
809:     KSPSetType(smoother, KSPPREONLY);
810:     PCSetUp_MG(pc);
811:   }
812:   return(0);
813: }

815: /* ------------------------------------------------------------------------- */
816: /*
817:  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
818:    that was created with PCCreate_GAMG().

820:    Input Parameter:
821: .  pc - the preconditioner context

823:    Application Interface Routine: PCDestroy()
824: */
825: PetscErrorCode PCDestroy_GAMG(PC pc)
826: {
828:   PC_MG          *mg     = (PC_MG*)pc->data;
829:   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;

832:   PCReset_GAMG(pc);
833:   if (pc_gamg->ops->destroy) {
834:     (*pc_gamg->ops->destroy)(pc);
835:   }
836:   PetscFree(pc_gamg->ops);
837:   PetscFree(pc_gamg->gamg_type_name);
838:   PetscFree(pc_gamg);
839:   PCDestroy_MG(pc);
840:   return(0);
841: }

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

846:    Logically Collective on PC

848:    Input Parameters:
849: +  pc - the preconditioner context
850: -  n - the number of equations


853:    Options Database Key:
854: .  -pc_gamg_process_eq_limit <limit>

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

860:    Level: intermediate

862: .seealso: PCGAMGSetCoarseEqLim()
863: @*/
864: PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
865: {

870:   PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));
871:   return(0);
872: }

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

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

884: /*@
885:    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.

887:  Collective on PC

889:    Input Parameters:
890: +  pc - the preconditioner context
891: -  n - maximum number of equations to aim for

893:    Options Database Key:
894: .  -pc_gamg_coarse_eq_limit <limit>

896:    Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
897:      has less than 1000 unknowns.

899:    Level: intermediate

901: .seealso: PCGAMGSetProcEqLim()
902: @*/
903: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
904: {

909:   PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
910:   return(0);
911: }

913: static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
914: {
915:   PC_MG   *mg      = (PC_MG*)pc->data;
916:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

919:   if (n>0) pc_gamg->coarse_eq_limit = n;
920:   return(0);
921: }

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

926:    Collective on PC

928:    Input Parameters:
929: +  pc - the preconditioner context
930: -  n - PETSC_TRUE or PETSC_FALSE

932:    Options Database Key:
933: .  -pc_gamg_repartition <true,false>

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

938:    Level: intermediate

940: .seealso: ()
941: @*/
942: PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
943: {

948:   PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));
949:   return(0);
950: }

952: static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
953: {
954:   PC_MG   *mg      = (PC_MG*)pc->data;
955:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

958:   pc_gamg->repart = n;
959:   return(0);
960: }

962: /*@
963:    PCGAMGSetEstEigKSPMaxIt - Set number of KSP iterations in eigen estimator (for Cheby)

965:    Collective on PC

967:    Input Parameters:
968: +  pc - the preconditioner context
969: -  n - number of its

971:    Options Database Key:
972: .  -pc_gamg_esteig_ksp_max_it <its>

974:    Notes:

976:    Level: intermediate

978: .seealso: ()
979: @*/
980: PetscErrorCode PCGAMGSetEstEigKSPMaxIt(PC pc, PetscInt n)
981: {

986:   PetscTryMethod(pc,"PCGAMGSetEstEigKSPMaxIt_C",(PC,PetscInt),(pc,n));
987:   return(0);
988: }

990: static PetscErrorCode PCGAMGSetEstEigKSPMaxIt_GAMG(PC pc, PetscInt n)
991: {
992:   PC_MG   *mg      = (PC_MG*)pc->data;
993:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

996:   pc_gamg->esteig_max_it = n;
997:   return(0);
998: }

1000: /*@
1001:    PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Cheby smoother

1003:    Collective on PC

1005:    Input Parameters:
1006: +  pc - the preconditioner context
1007: -  n - number of its

1009:    Options Database Key:
1010: .  -pc_gamg_use_sa_esteig <true,false>

1012:    Notes:

1014:    Level: intermediate

1016: .seealso: ()
1017: @*/
1018: PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
1019: {

1024:   PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));
1025:   return(0);
1026: }

1028: static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
1029: {
1030:   PC_MG   *mg      = (PC_MG*)pc->data;
1031:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1034:   pc_gamg->use_sa_esteig = n ? 1 : 0;
1035:   return(0);
1036: }

1038: /*@C
1039:    PCGAMGSetEstEigKSPType - Set type of KSP in eigen estimator (for Cheby)

1041:    Collective on PC

1043:    Input Parameters:
1044: +  pc - the preconditioner context
1045: -  t - ksp type

1047:    Options Database Key:
1048: .  -pc_gamg_esteig_ksp_type <type>

1050:    Notes:

1052:    Level: intermediate

1054: .seealso: ()
1055: @*/
1056: PetscErrorCode PCGAMGSetEstEigKSPType(PC pc, char t[])
1057: {

1062:   PetscTryMethod(pc,"PCGAMGSetEstEigKSPType_C",(PC,char[]),(pc,t));
1063:   return(0);
1064: }

1066: static PetscErrorCode PCGAMGSetEstEigKSPType_GAMG(PC pc, char t[])
1067: {
1069:   PC_MG   *mg      = (PC_MG*)pc->data;
1070:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1073:   PetscStrcpy(pc_gamg->esteig_type,t);
1074:   return(0);
1075: }

1077: /*@
1078:    PCGAMGSetEigenvalues - Set eigenvalues

1080:    Collective on PC

1082:    Input Parameters:
1083: +  pc - the preconditioner context
1084: -  emax - max eigenvalue
1085: .  emin - min eigenvalue

1087:    Options Database Key:
1088: .  -pc_gamg_eigenvalues

1090:    Level: intermediate

1092: .seealso: PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPType()
1093: @*/
1094: PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin)
1095: {

1100:   PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));
1101:   return(0);
1102: }
1103: static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin)
1104: {
1105:   PC_MG          *mg      = (PC_MG*)pc->data;
1106:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1109:   if (emax <= emin) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin);
1110:   if (emax*emin <= 0.0) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin);
1111:   pc_gamg->emax = emax;
1112:   pc_gamg->emin = emin;

1114:   return(0);
1115: }

1117: /*@
1118:    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner

1120:    Collective on PC

1122:    Input Parameters:
1123: +  pc - the preconditioner context
1124: -  n - PETSC_TRUE or PETSC_FALSE

1126:    Options Database Key:
1127: .  -pc_gamg_reuse_interpolation <true,false>

1129:    Level: intermediate

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

1135: .seealso: ()
1136: @*/
1137: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1138: {

1143:   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
1144:   return(0);
1145: }

1147: static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1148: {
1149:   PC_MG   *mg      = (PC_MG*)pc->data;
1150:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1153:   pc_gamg->reuse_prol = n;
1154:   return(0);
1155: }

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

1160:    Collective on PC

1162:    Input Parameters:
1163: +  pc - the preconditioner context
1164: -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not

1166:    Options Database Key:
1167: .  -pc_gamg_asm_use_agg

1169:    Level: intermediate

1171: .seealso: ()
1172: @*/
1173: PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1174: {

1179:   PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));
1180:   return(0);
1181: }

1183: static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1184: {
1185:   PC_MG   *mg      = (PC_MG*)pc->data;
1186:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1189:   pc_gamg->use_aggs_in_asm = flg;
1190:   return(0);
1191: }

1193: /*@
1194:    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver

1196:    Collective on PC

1198:    Input Parameters:
1199: +  pc - the preconditioner context
1200: -  flg - PETSC_TRUE to not force coarse grid onto one processor

1202:    Options Database Key:
1203: .  -pc_gamg_use_parallel_coarse_grid_solver

1205:    Level: intermediate

1207: .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids()
1208: @*/
1209: PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1210: {

1215:   PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));
1216:   return(0);
1217: }

1219: static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1220: {
1221:   PC_MG   *mg      = (PC_MG*)pc->data;
1222:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1225:   pc_gamg->use_parallel_coarse_grid_solver = flg;
1226:   return(0);
1227: }

1229: /*@
1230:    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU

1232:    Collective on PC

1234:    Input Parameters:
1235: +  pc - the preconditioner context
1236: -  flg - PETSC_TRUE to pin coarse grids to CPU

1238:    Options Database Key:
1239: .  -pc_gamg_cpu_pin_coarse_grids

1241:    Level: intermediate

1243: .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve()
1244: @*/
1245: PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1246: {

1251:   PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));
1252:   return(0);
1253: }

1255: static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1256: {
1257:   PC_MG   *mg      = (PC_MG*)pc->data;
1258:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1261:   pc_gamg->cpu_pin_coarse_grids = flg;
1262:   return(0);
1263: }

1265: /*@
1266:    PCGAMGSetCoarseGridLayoutType - place reduce grids on processors with natural order (compact type)

1268:    Collective on PC

1270:    Input Parameters:
1271: +  pc - the preconditioner context
1272: -  flg - Layout type

1274:    Options Database Key:
1275: .  -pc_gamg_coarse_grid_layout_type

1277:    Level: intermediate

1279: .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids()
1280: @*/
1281: PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1282: {

1287:   PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));
1288:   return(0);
1289: }

1291: static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1292: {
1293:   PC_MG   *mg      = (PC_MG*)pc->data;
1294:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1297:   pc_gamg->layout_type = flg;
1298:   return(0);
1299: }

1301: /*@
1302:    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use

1304:    Not collective on PC

1306:    Input Parameters:
1307: +  pc - the preconditioner
1308: -  n - the maximum number of levels to use

1310:    Options Database Key:
1311: .  -pc_mg_levels

1313:    Level: intermediate

1315: .seealso: ()
1316: @*/
1317: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1318: {

1323:   PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
1324:   return(0);
1325: }

1327: static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1328: {
1329:   PC_MG   *mg      = (PC_MG*)pc->data;
1330:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1333:   pc_gamg->Nlevels = n;
1334:   return(0);
1335: }

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

1340:    Not collective on PC

1342:    Input Parameters:
1343: +  pc - the preconditioner context
1344: .  threshold - array of threshold values for finest n levels; 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph
1345: -  n - number of threshold values provided in array

1347:    Options Database Key:
1348: .  -pc_gamg_threshold <threshold>

1350:    Notes:
1351:     Increasing the threshold decreases the rate of coarsening. Conversely reducing the threshold increases the rate of coarsening (aggressive coarsening) and thereby reduces the complexity of the coarse grids, and generally results in slower solver converge rates. Reducing coarse grid complexity reduced the complexity of Galerkin coarse grid construction considerably.
1352:     Before coarsening or aggregating the graph, GAMG removes small values from the graph with this threshold, and thus reducing the coupling in the graph and a different (perhaps better) coarser set of points.

1354:     If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening.
1355:     In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold().
1356:     If n is greater than the total number of levels, the excess entries in threshold will not be used.

1358:    Level: intermediate

1360: .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
1361: @*/
1362: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1363: {

1369:   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));
1370:   return(0);
1371: }

1373: static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1374: {
1375:   PC_MG   *mg      = (PC_MG*)pc->data;
1376:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1377:   PetscInt i;
1379:   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1380:   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;
1381:   return(0);
1382: }

1384: /*@
1385:    PCGAMGSetThresholdScale - Relative threshold reduction at each level

1387:    Not collective on PC

1389:    Input Parameters:
1390: +  pc - the preconditioner context
1391: -  scale - the threshold value reduction, ussually < 1.0

1393:    Options Database Key:
1394: .  -pc_gamg_threshold_scale <v>

1396:    Notes:
1397:    The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold().
1398:    This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold().

1400:    Level: advanced

1402: .seealso: PCGAMGSetThreshold()
1403: @*/
1404: PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1405: {

1410:   PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));
1411:   return(0);
1412: }

1414: static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1415: {
1416:   PC_MG   *mg      = (PC_MG*)pc->data;
1417:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1419:   pc_gamg->threshold_scale = v;
1420:   return(0);
1421: }

1423: /*@C
1424:    PCGAMGSetType - Set solution method

1426:    Collective on PC

1428:    Input Parameters:
1429: +  pc - the preconditioner context
1430: -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL

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

1435:    Level: intermediate

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

1445:   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1446:   return(0);
1447: }

1449: /*@C
1450:    PCGAMGGetType - Get solution method

1452:    Collective on PC

1454:    Input Parameter:
1455: .  pc - the preconditioner context

1457:    Output Parameter:
1458: .  type - the type of algorithm used

1460:    Level: intermediate

1462: .seealso: PCGAMGSetType(), PCGAMGType
1463: @*/
1464: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1465: {

1470:   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1471:   return(0);
1472: }

1474: static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1475: {
1476:   PC_MG          *mg      = (PC_MG*)pc->data;
1477:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1480:   *type = pc_gamg->type;
1481:   return(0);
1482: }

1484: static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1485: {
1486:   PetscErrorCode ierr,(*r)(PC);
1487:   PC_MG          *mg      = (PC_MG*)pc->data;
1488:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1491:   pc_gamg->type = type;
1492:   PetscFunctionListFind(GAMGList,type,&r);
1493:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1494:   if (pc_gamg->ops->destroy) {
1495:     (*pc_gamg->ops->destroy)(pc);
1496:     PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));
1497:     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1498:     /* cleaning up common data in pc_gamg - this should disapear someday */
1499:     pc_gamg->data_cell_cols = 0;
1500:     pc_gamg->data_cell_rows = 0;
1501:     pc_gamg->orig_data_cell_cols = 0;
1502:     pc_gamg->orig_data_cell_rows = 0;
1503:     PetscFree(pc_gamg->data);
1504:     pc_gamg->data_sz = 0;
1505:   }
1506:   PetscFree(pc_gamg->gamg_type_name);
1507:   PetscStrallocpy(type,&pc_gamg->gamg_type_name);
1508:   (*r)(pc);
1509:   return(0);
1510: }

1512: /* -------------------------------------------------------------------------- */
1513: /*
1514:    PCMGGetGridComplexity - compute coarse grid complexity of MG hierarchy

1516:    Input Parameter:
1517: .  pc - the preconditioner context

1519:    Output Parameter:
1520: .  gc - grid complexity = sum_i(nnz_i) / nnz_0

1522:    Level: advanced
1523: */
1524: static PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc)
1525: {
1527:   PC_MG          *mg      = (PC_MG*)pc->data;
1528:   PC_MG_Levels   **mglevels = mg->levels;
1529:   PetscInt       lev;
1530:   PetscLogDouble nnz0 = 0, sgc = 0;
1531:   MatInfo        info;

1534:   if (!pc->setupcalled) {
1535:     *gc = 0;
1536:     return(0);
1537:   }
1538:   if (!mg->nlevels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels");
1539:   for (lev=0; lev<mg->nlevels; lev++) {
1540:     Mat dB;
1541:     KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB);
1542:     MatGetInfo(dB,MAT_GLOBAL_SUM,&info); /* global reduction */
1543:     sgc += info.nz_used;
1544:     if (lev==mg->nlevels-1) nnz0 = info.nz_used;
1545:   }
1546:   if (nnz0 > 0) *gc = (PetscReal)(sgc/nnz0);
1547:   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number for grid points on finest level is not available");
1548:   return(0);
1549: }

1551: static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
1552: {
1553:   PetscErrorCode ierr,i;
1554:   PC_MG          *mg      = (PC_MG*)pc->data;
1555:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1556:   PetscReal       gc=0;
1558:   PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");
1559:   PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");
1560:   for (i=0;i<pc_gamg->current_level;i++) {
1561:     PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);
1562:   }
1563:   PetscViewerASCIIPrintf(viewer,"\n");
1564:   PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);
1565:   if (pc_gamg->use_aggs_in_asm) {
1566:     PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");
1567:   }
1568:   if (pc_gamg->use_parallel_coarse_grid_solver) {
1569:     PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");
1570:   }
1571: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1572:   if (pc_gamg->cpu_pin_coarse_grids) {
1573:     /* PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n"); */
1574:   }
1575: #endif
1576:   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
1577:   /*   PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n"); */
1578:   /* } else { */
1579:   /*   PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n"); */
1580:   /* } */
1581:   if (pc_gamg->ops->view) {
1582:     (*pc_gamg->ops->view)(pc,viewer);
1583:   }
1584:   PCMGGetGridComplexity(pc,&gc);
1585:   PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g\n",gc);
1586:   return(0);
1587: }

1589: PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
1590: {
1592:   PC_MG          *mg      = (PC_MG*)pc->data;
1593:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1594:   PetscBool      flag,f2;
1595:   MPI_Comm       comm;
1596:   char           prefix[256],tname[32];
1597:   PetscInt       i,n;
1598:   const char     *pcpre;
1599:   static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL};
1601:   PetscObjectGetComm((PetscObject)pc,&comm);
1602:   PetscOptionsHead(PetscOptionsObject,"GAMG options");
1603:   PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);
1604:   if (flag) {
1605:     PCGAMGSetType(pc,tname);
1606:   }
1607:   PetscOptionsFList("-pc_gamg_esteig_ksp_type","Krylov method for eigen estimator","PCGAMGSetEstEigKSPType",KSPList,pc_gamg->esteig_type,tname,sizeof(tname),&flag);
1608:   if (flag) {
1609:     PCGAMGSetEstEigKSPType(pc,tname);
1610:   }
1611:   PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);
1612:   f2 = PETSC_TRUE;
1613:   PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",f2,&f2,&flag);
1614:   if (flag) pc_gamg->use_sa_esteig = f2 ? 1 : 0;
1615:   PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);
1616:   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);
1617:   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);
1618:   PetscOptionsBool("-pc_gamg_cpu_pin_coarse_grids","Pin coarse grids to the CPU","PCGAMGSetCpuPinCoarseGrids",pc_gamg->cpu_pin_coarse_grids,&pc_gamg->cpu_pin_coarse_grids,NULL);
1619:   PetscOptionsEnum("-pc_gamg_coarse_grid_layout_type","compact: place reduced grids on processes in natural order; spread: distribute to whole machine for more memory bandwidth","PCGAMGSetCoarseGridLayoutType",LayoutTypes,(PetscEnum)pc_gamg->layout_type,(PetscEnum*)&pc_gamg->layout_type,NULL);
1620:   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);
1621:   PetscOptionsInt("-pc_gamg_esteig_ksp_max_it","Number of iterations of eigen estimator","PCGAMGSetEstEigKSPMaxIt",pc_gamg->esteig_max_it,&pc_gamg->esteig_max_it,NULL);
1622:   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);
1623:   PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL);
1624:   n = PETSC_MG_MAXLEVELS;
1625:   PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);
1626:   if (!flag || n < PETSC_MG_MAXLEVELS) {
1627:     if (!flag) n = 1;
1628:     i = n;
1629:     do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS);
1630:   }
1631:   PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);
1632:   {
1633:     PetscReal eminmax[2] = {0., 0.};
1634:     n = 2;
1635:     PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag);
1636:     if (flag) {
1637:       if (n != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
1638:       PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);
1639:     }
1640:   }
1641:   /* set options for subtype */
1642:   if (pc_gamg->ops->setfromoptions) {(*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);}

1644:   PCGetOptionsPrefix(pc, &pcpre);
1645:   PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
1646:   PetscOptionsTail();
1647:   return(0);
1648: }

1650: /* -------------------------------------------------------------------------- */
1651: /*MC
1652:      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner

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

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

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


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

1683:   Level: intermediate

1685: .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1686:            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetEstEigKSPType()
1687: M*/

1689: PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1690: {
1691:   PetscErrorCode ierr,i;
1692:   PC_GAMG        *pc_gamg;
1693:   PC_MG          *mg;

1696:    /* register AMG type */
1697:   PCGAMGInitializePackage();

1699:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1700:   PCSetType(pc, PCMG);
1701:   PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);

1703:   /* create a supporting struct and attach it to pc */
1704:   PetscNewLog(pc,&pc_gamg);
1705:   PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);
1706:   mg           = (PC_MG*)pc->data;
1707:   mg->innerctx = pc_gamg;

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

1711:   pc_gamg->setup_count = 0;
1712:   /* these should be in subctx but repartitioning needs simple arrays */
1713:   pc_gamg->data_sz = 0;
1714:   pc_gamg->data    = NULL;

1716:   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1717:   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1718:   pc->ops->setup          = PCSetUp_GAMG;
1719:   pc->ops->reset          = PCReset_GAMG;
1720:   pc->ops->destroy        = PCDestroy_GAMG;
1721:   mg->view                = PCView_GAMG;

1723:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);
1724:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);
1725:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);
1726:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);
1727:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);
1728:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPType_C",PCGAMGSetEstEigKSPType_GAMG);
1729:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPMaxIt_C",PCGAMGSetEstEigKSPMaxIt_GAMG);
1730:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG);
1731:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG);
1732:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);
1733:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);
1734:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);
1735:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);
1736:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);
1737:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);
1738:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);
1739:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);
1740:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);
1741:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);
1742:   pc_gamg->repart           = PETSC_FALSE;
1743:   pc_gamg->reuse_prol       = PETSC_FALSE;
1744:   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1745:   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1746:   pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1747:   pc_gamg->layout_type      = PCGAMG_LAYOUT_SPREAD;
1748:   pc_gamg->min_eq_proc      = 50;
1749:   pc_gamg->coarse_eq_limit  = 50;
1750:   for (i=0;i<PETSC_MG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.;
1751:   pc_gamg->threshold_scale = 1.;
1752:   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
1753:   pc_gamg->current_level    = 0; /* don't need to init really */
1754:   PetscStrcpy(pc_gamg->esteig_type,NULL);
1755:   pc_gamg->esteig_max_it    = 10;
1756:   pc_gamg->use_sa_esteig    = -1;
1757:   pc_gamg->emin             = 0;
1758:   pc_gamg->emax             = 0;

1760:   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;

1762:   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1763:   PCGAMGSetType(pc,PCGAMGAGG);
1764:   return(0);
1765: }

1767: /*@C
1768:  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
1769:     from PCInitializePackage().

1771:  Level: developer

1773:  .seealso: PetscInitialize()
1774: @*/
1775: PetscErrorCode PCGAMGInitializePackage(void)
1776: {

1780:   if (PCGAMGPackageInitialized) return(0);
1781:   PCGAMGPackageInitialized = PETSC_TRUE;
1782:   PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);
1783:   PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);
1784:   PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);
1785:   PetscRegisterFinalize(PCGAMGFinalizePackage);

1787:   /* general events */
1788:   PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);
1789:   PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);
1790:   PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1791:   PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1792:   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1793:   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1794:   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);

1796: #if defined PETSC_GAMG_USE_LOG
1797:   PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1798:   PetscLogEventRegister("  Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1799:   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1800:   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1801:   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1802:   PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1803:   PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1804:   PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1805:   PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1806:   PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1807:   PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1808:   PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1809:   PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1810:   PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1811:   PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1812:   PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1813:   PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);

1815:   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1816:   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1817:   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1818:   /* create timer stages */
1819: #if defined GAMG_STAGES
1820:   {
1821:     char     str[32];
1822:     PetscInt lidx;
1823:     sprintf(str,"MG Level %d (finest)",0);
1824:     PetscLogStageRegister(str, &gamg_stages[0]);
1825:     for (lidx=1; lidx<9; lidx++) {
1826:       sprintf(str,"MG Level %d",lidx);
1827:       PetscLogStageRegister(str, &gamg_stages[lidx]);
1828:     }
1829:   }
1830: #endif
1831: #endif
1832:   return(0);
1833: }

1835: /*@C
1836:  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1837:     called from PetscFinalize() automatically.

1839:  Level: developer

1841:  .seealso: PetscFinalize()
1842: @*/
1843: PetscErrorCode PCGAMGFinalizePackage(void)
1844: {

1848:   PCGAMGPackageInitialized = PETSC_FALSE;
1849:   PetscFunctionListDestroy(&GAMGList);
1850:   return(0);
1851: }

1853: /*@C
1854:  PCGAMGRegister - Register a PCGAMG implementation.

1856:  Input Parameters:
1857:  + type - string that will be used as the name of the GAMG type.
1858:  - create - function for creating the gamg context.

1860:   Level: advanced

1862:  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1863: @*/
1864: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1865: {

1869:   PCGAMGInitializePackage();
1870:   PetscFunctionListAdd(&GAMGList,type,create);
1871:   return(0);
1872: }