Actual source code: gamg.c

  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/ksp/impls/cheby/chebyshevimpl.h>

  8: #if defined(PETSC_HAVE_CUDA)
  9:   #include <cuda_runtime.h>
 10: #endif

 12: #if defined(PETSC_HAVE_HIP)
 13:   #include <hip/hip_runtime.h>
 14: #endif

 16: PetscLogEvent petsc_gamg_setup_events[NUM_SET];
 17: PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3];
 18: PetscLogEvent PC_GAMGGraph_AGG;
 19: PetscLogEvent PC_GAMGGraph_GEO;
 20: PetscLogEvent PC_GAMGCoarsen_AGG;
 21: PetscLogEvent PC_GAMGCoarsen_GEO;
 22: PetscLogEvent PC_GAMGProlongator_AGG;
 23: PetscLogEvent PC_GAMGProlongator_GEO;
 24: PetscLogEvent PC_GAMGOptProlongator_AGG;

 26: /* #define GAMG_STAGES */
 27: #if defined(GAMG_STAGES)
 28: static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
 29: #endif

 31: static PetscFunctionList GAMGList = NULL;
 32: static PetscBool PCGAMGPackageInitialized;

 34: /* ----------------------------------------------------------------------------- */
 35: PetscErrorCode PCReset_GAMG(PC pc)
 36: {
 37:   PC_MG          *mg      = (PC_MG*)pc->data;
 38:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

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

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

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

 69: 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)
 70: {
 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;

 78:   PetscObjectGetComm((PetscObject)Amat_fine,&comm);
 79:   MPI_Comm_rank(comm, &rank);
 80:   MPI_Comm_size(comm, &size);
 81:   MatGetBlockSize(Amat_fine, &f_bs);
 82:   PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0);
 83:   MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);
 84:   PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0);

 86:   if (Pcolumnperm) *Pcolumnperm = NULL;

 88:   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
 89:   MatGetLocalSize(Cmat, &ncrs_eq, NULL);
 90:   if (pc_gamg->data_cell_rows>0) {
 91:     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
 92:   } else {
 93:     PetscInt  bs;
 94:     MatGetBlockSize(Cmat, &bs);
 95:     ncrs = ncrs_eq/bs;
 96:   }
 97:   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
 98:   if (pc_gamg->level_reduction_factors[pc_gamg->current_level] == 0 && PetscDefined(HAVE_CUDA) && pc_gamg->current_level==0) { /* 0 turns reducing to 1 process/device on; do for HIP, etc. */
 99: #if defined(PETSC_HAVE_CUDA)
100:     PetscShmComm pshmcomm;
101:     PetscMPIInt  locrank;
102:     MPI_Comm     loccomm;
103:     PetscInt     s_nnodes,r_nnodes, new_new_size;
104:     cudaError_t  cerr;
105:     int          devCount;
106:     PetscShmCommGet(comm,&pshmcomm);
107:     PetscShmCommGetMpiShmComm(pshmcomm,&loccomm);
108:     MPI_Comm_rank(loccomm, &locrank);
109:     s_nnodes = !locrank;
110:     MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm);
112:     devCount = 0;
113:     cerr = cudaGetDeviceCount(&devCount);
114:     cudaGetLastError(); /* Reset the last error */
115:     if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */
116:       new_new_size = r_nnodes * devCount;
117:       new_size = new_new_size;
118:       PetscInfo(pc,"%s: Fine grid with Cuda. %D nodes. Change new active set size %d --> %d (devCount=%d #nodes=%D)\n",((PetscObject)pc)->prefix,r_nnodes,nactive,new_size,devCount,r_nnodes);
119:     } else {
120:       PetscInfo(pc,"%s: With Cuda but no device. Use heuristics.");
121:       goto HEURISTIC;
122:     }
123: #else
124:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here");
125: #endif
126:   } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) {
128:     new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level];
129:     PetscInfo(pc,"%s: Manually setting reduction to %d active processes (%d/%D)\n",((PetscObject)pc)->prefix,new_size,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]);
130:   } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) {
131:     new_size = 1;
132:     PetscInfo(pc,"%s: Force coarsest grid reduction to %d active processes\n",((PetscObject)pc)->prefix,new_size);
133:   } else {
134:     PetscInt ncrs_eq_glob;
135: #if defined(PETSC_HAVE_CUDA)
136:     HEURISTIC:
137: #endif
138:     MatGetSize(Cmat, &ncrs_eq_glob, NULL);
139:     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
140:     if (!new_size) new_size = 1; /* not likely, posible? */
141:     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
142:     PetscInfo(pc,"%s: Coarse grid reduction from %d to %d active processes\n",((PetscObject)pc)->prefix,nactive,new_size);
143:   }
144:   if (new_size==nactive) {
145:     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
146:     if (new_size < size) {
147:       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
148:       PetscInfo(pc,"%s: reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",((PetscObject)pc)->prefix,nactive);
149:       if (pc_gamg->cpu_pin_coarse_grids) {
150:         MatBindToCPU(*a_Amat_crs,PETSC_TRUE);
151:         MatBindToCPU(*a_P_inout,PETSC_TRUE);
152:       }
153:     }
154:     /* we know that the grid structure can be reused in MatPtAP */
155:   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
156:     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1;
157:     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
158:     nloc_old = ncrs_eq/cr_bs;
160:     /* get new_size and rfactor */
161:     if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
162:       /* find factor */
163:       if (new_size == 1) rfactor = size; /* don't modify */
164:       else {
165:         PetscReal best_fact = 0.;
166:         jj = -1;
167:         for (kk = 1 ; kk <= size ; kk++) {
168:           if (!(size%kk)) { /* a candidate */
169:             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
170:             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
171:             if (fact > best_fact) {
172:               best_fact = fact; jj = kk;
173:             }
174:           }
175:         }
176:         if (jj != -1) rfactor = jj;
177:         else rfactor = 1; /* a prime */
178:         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
179:         else expand_factor = rfactor;
180:       }
181:       new_size = size/rfactor; /* make new size one that is factor */
182:       if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
183:         *a_Amat_crs = Cmat;
184:         PetscInfo(pc,"%s: Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",((PetscObject)pc)->prefix,new_size,ncrs_eq);
185:         return 0;
186:       }
187:     }
188:     PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);
189:     /* make 'is_eq_newproc' */
190:     PetscMalloc1(size, &counts);
191:     if (pc_gamg->repart) {
192:       /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */
193:       Mat      adj;
194:       PetscInfo(pc,"%s: Repartition: size (active): %d --> %d, %D local equations, using %s process layout\n",((PetscObject)pc)->prefix,*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread");
195:       /* get 'adj' */
196:       if (cr_bs == 1) {
197:         MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);
198:       } else {
199:         /* make a scalar matrix to partition (no Stokes here) */
200:         Mat               tMat;
201:         PetscInt          Istart_crs,Iend_crs,ncols,jj,Ii;
202:         const PetscScalar *vals;
203:         const PetscInt    *idx;
204:         PetscInt          *d_nnz, *o_nnz, M, N;
205:         static PetscInt   llev = 0; /* ugly but just used for debugging */
206:         MatType           mtype;

208:         PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);
209:         MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);
210:         MatGetSize(Cmat, &M, &N);
211:         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
212:           MatGetRow(Cmat,Ii,&ncols,NULL,NULL);
213:           d_nnz[jj] = ncols/cr_bs;
214:           o_nnz[jj] = ncols/cr_bs;
215:           MatRestoreRow(Cmat,Ii,&ncols,NULL,NULL);
216:           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
217:           if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs;
218:         }

220:         MatGetType(Amat_fine,&mtype);
221:         MatCreate(comm, &tMat);
222:         MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);
223:         MatSetType(tMat,mtype);
224:         MatSeqAIJSetPreallocation(tMat,0,d_nnz);
225:         MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);
226:         PetscFree2(d_nnz,o_nnz);

228:         for (ii = Istart_crs; ii < Iend_crs; ii++) {
229:           PetscInt dest_row = ii/cr_bs;
230:           MatGetRow(Cmat,ii,&ncols,&idx,&vals);
231:           for (jj = 0; jj < ncols; jj++) {
232:             PetscInt    dest_col = idx[jj]/cr_bs;
233:             PetscScalar v        = 1.0;
234:             MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);
235:           }
236:           MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);
237:         }
238:         MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);
239:         MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);

241:         if (llev++ == -1) {
242:           PetscViewer viewer; char fname[32];
243:           PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);
244:           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
245:           MatView(tMat, viewer);
246:           PetscViewerDestroy(&viewer);
247:         }
248:         MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);
249:         MatDestroy(&tMat);
250:       } /* create 'adj' */

252:       { /* partition: get newproc_idx */
253:         char            prefix[256];
254:         const char      *pcpre;
255:         const PetscInt  *is_idx;
256:         MatPartitioning mpart;
257:         IS              proc_is;

259:         MatPartitioningCreate(comm, &mpart);
260:         MatPartitioningSetAdjacency(mpart, adj);
261:         PCGetOptionsPrefix(pc, &pcpre);
262:         PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
263:         PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);
264:         MatPartitioningSetFromOptions(mpart);
265:         MatPartitioningSetNParts(mpart, new_size);
266:         MatPartitioningApply(mpart, &proc_is);
267:         MatPartitioningDestroy(&mpart);

269:         /* collect IS info */
270:         PetscMalloc1(ncrs_eq, &newproc_idx);
271:         ISGetIndices(proc_is, &is_idx);
272:         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
273:           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
274:             newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
275:           }
276:         }
277:         ISRestoreIndices(proc_is, &is_idx);
278:         ISDestroy(&proc_is);
279:       }
280:       MatDestroy(&adj);

282:       ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);
283:       PetscFree(newproc_idx);
284:     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
285:       PetscInt targetPE;
287:       PetscInfo(pc,"%s: Number of equations (loc) %D with simple aggregation\n",((PetscObject)pc)->prefix,ncrs_eq);
288:       targetPE = (rank/rfactor)*expand_factor;
289:       ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);
290:     } /* end simple 'is_eq_newproc' */

292:     /*
293:       Create an index set from the is_eq_newproc index set to indicate the mapping TO
294:     */
295:     ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);
296:     is_eq_num_prim = is_eq_num;
297:     /*
298:       Determine how many equations/vertices are assigned to each processor
299:     */
300:     ISPartitioningCount(is_eq_newproc, size, counts);
301:     ncrs_eq_new = counts[rank];
302:     ISDestroy(&is_eq_newproc);
303:     ncrs_new = ncrs_eq_new/cr_bs;

305:     PetscFree(counts);
306:     /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxiliary data into some complex abstracted thing */
307:     {
308:       Vec            src_crd, dest_crd;
309:       const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols;
310:       VecScatter     vecscat;
311:       PetscScalar    *array;
312:       IS isscat;
313:       /* move data (for primal equations only) */
314:       /* Create a vector to contain the newly ordered element information */
315:       VecCreate(comm, &dest_crd);
316:       VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);
317:       VecSetType(dest_crd,VECSTANDARD); /* this is needed! */
318:       /*
319:         There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
320:         a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
321:       */
322:       PetscMalloc1(ncrs*node_data_sz, &tidx);
323:       ISGetIndices(is_eq_num_prim, &idx);
324:       for (ii=0,jj=0; ii<ncrs; ii++) {
325:         PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */
326:         for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk;
327:       }
328:       ISRestoreIndices(is_eq_num_prim, &idx);
329:       ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);
330:       PetscFree(tidx);
331:       /*
332:         Create a vector to contain the original vertex information for each element
333:       */
334:       VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);
335:       for (jj=0; jj<ndata_cols; jj++) {
336:         const PetscInt stride0=ncrs*pc_gamg->data_cell_rows;
337:         for (ii=0; ii<ncrs; ii++) {
338:           for (kk=0; kk<ndata_rows; kk++) {
339:             PetscInt    ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj;
340:             PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
341:             VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);
342:           }
343:         }
344:       }
345:       VecAssemblyBegin(src_crd);
346:       VecAssemblyEnd(src_crd);
347:       /*
348:         Scatter the element vertex information (still in the original vertex ordering)
349:         to the correct processor
350:       */
351:       VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);
352:       ISDestroy(&isscat);
353:       VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);
354:       VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);
355:       VecScatterDestroy(&vecscat);
356:       VecDestroy(&src_crd);
357:       /*
358:         Put the element vertex data into a new allocation of the gdata->ele
359:       */
360:       PetscFree(pc_gamg->data);
361:       PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);

363:       pc_gamg->data_sz = node_data_sz*ncrs_new;
364:       strideNew        = ncrs_new*ndata_rows;

366:       VecGetArray(dest_crd, &array);
367:       for (jj=0; jj<ndata_cols; jj++) {
368:         for (ii=0; ii<ncrs_new; ii++) {
369:           for (kk=0; kk<ndata_rows; kk++) {
370:             PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj;
371:             pc_gamg->data[ix] = PetscRealPart(array[jx]);
372:           }
373:         }
374:       }
375:       VecRestoreArray(dest_crd, &array);
376:       VecDestroy(&dest_crd);
377:     }
378:     /* move A and P (columns) with new layout */
379:     PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);
380:     /*
381:       Invert for MatCreateSubMatrix
382:     */
383:     ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);
384:     ISSort(new_eq_indices); /* is this needed? */
385:     ISSetBlockSize(new_eq_indices, cr_bs);
386:     if (is_eq_num != is_eq_num_prim) {
387:       ISDestroy(&is_eq_num_prim); /* could be same as 'is_eq_num' */
388:     }
389:     if (Pcolumnperm) {
390:       PetscObjectReference((PetscObject)new_eq_indices);
391:       *Pcolumnperm = new_eq_indices;
392:     }
393:     ISDestroy(&is_eq_num);
394:     PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);
395:     PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);
396:     /* 'a_Amat_crs' output */
397:     {
398:       Mat       mat;
399:       PetscBool flg;
400:       MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);
401:       MatGetOption(Cmat, MAT_SPD, &flg);
402:       if (flg) {
403:         MatSetOption(mat, MAT_SPD,PETSC_TRUE);
404:       } else {
405:         MatGetOption(Cmat, MAT_HERMITIAN, &flg);
406:         if (flg) {
407:           MatSetOption(mat, MAT_HERMITIAN,PETSC_TRUE);
408:         } else {
409: #if !defined(PETSC_USE_COMPLEX)
410:           MatGetOption(Cmat, MAT_SYMMETRIC, &flg);
411:           if (flg) {
412:             MatSetOption(mat, MAT_SYMMETRIC,PETSC_TRUE);
413:           }
414: #endif
415:         }
416:       }
417:       *a_Amat_crs = mat;
418:     }
419:     MatDestroy(&Cmat);

421:     PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);
422:     /* prolongator */
423:     {
424:       IS       findices;
425:       PetscInt Istart,Iend;
426:       Mat      Pnew;

428:       MatGetOwnershipRange(Pold, &Istart, &Iend);
429:       PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);
430:       ISCreateStride(comm,Iend-Istart,Istart,1,&findices);
431:       ISSetBlockSize(findices,f_bs);
432:       MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);
433:       ISDestroy(&findices);
434:       MatSetOption(Pnew,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);

436:       PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);
437:       MatDestroy(a_P_inout);

439:       /* output - repartitioned */
440:       *a_P_inout = Pnew;
441:     }
442:     ISDestroy(&new_eq_indices);

444:     *a_nactive_proc = new_size; /* output */

446:     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
447:     if (pc_gamg->cpu_pin_coarse_grids) {
448: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
449:       static PetscInt llev = 2;
450:       PetscInfo(pc,"%s: Pinning level %D to the CPU\n",((PetscObject)pc)->prefix,llev++);
451: #endif
452:       MatBindToCPU(*a_Amat_crs,PETSC_TRUE);
453:       MatBindToCPU(*a_P_inout,PETSC_TRUE);
454:       if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
455:         Mat         A = *a_Amat_crs, P = *a_P_inout;
456:         PetscMPIInt size;
457:         MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
458:         if (size > 1) {
459:           Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data;
460:           VecBindToCPU(a->lvec,PETSC_TRUE);
461:           VecBindToCPU(p->lvec,PETSC_TRUE);
462:         }
463:       }
464:     }
465:     PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);
466:   }
467:   return 0;
468: }

470: PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2)
471: {
472:   const char     *prefix;
473:   char           addp[32];
474:   PC_MG          *mg      = (PC_MG*)a_pc->data;
475:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

477:   PCGetOptionsPrefix(a_pc,&prefix);
478:   PetscInfo(a_pc,"%s: Square Graph on level %D\n",((PetscObject)a_pc)->prefix,pc_gamg->current_level+1);
479:   MatProductCreate(Gmat1,Gmat1,NULL,Gmat2);
480:   MatSetOptionsPrefix(*Gmat2,prefix);
481:   PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level);
482:   MatAppendOptionsPrefix(*Gmat2,addp);
483:   if ((*Gmat2)->structurally_symmetric) {
484:     MatProductSetType(*Gmat2,MATPRODUCT_AB);
485:   } else {
486:     MatSetOption(Gmat1,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);
487:     MatProductSetType(*Gmat2,MATPRODUCT_AtB);
488:   }
489:   MatProductSetFromOptions(*Gmat2);
490:   PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0);
491:   MatProductSymbolic(*Gmat2);
492:   PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0);
493:   MatProductClear(*Gmat2);
494:   /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
495:   (*Gmat2)->assembled = PETSC_TRUE;
496:   return 0;
497: }

499: /* -------------------------------------------------------------------------- */
500: /*
501:    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
502:                     by setting data structures and options.

504:    Input Parameter:
505: .  pc - the preconditioner context

507: */
508: PetscErrorCode PCSetUp_GAMG(PC pc)
509: {
511:   PC_MG          *mg      = (PC_MG*)pc->data;
512:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
513:   Mat            Pmat     = pc->pmat;
514:   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS];
515:   MPI_Comm       comm;
516:   PetscMPIInt    rank,size,nactivepe;
517:   Mat            Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS];
518:   IS             *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
519:   PetscLogDouble nnz0=0.,nnztot=0.;
520:   MatInfo        info;
521:   PetscBool      is_last = PETSC_FALSE;

523:   PetscObjectGetComm((PetscObject)pc,&comm);
524:   MPI_Comm_rank(comm,&rank);
525:   MPI_Comm_size(comm,&size);

527:   if (pc->setupcalled) {
528:     if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) {
529:       /* reset everything */
530:       PCReset_MG(pc);
531:       pc->setupcalled = 0;
532:     } else {
533:       PC_MG_Levels **mglevels = mg->levels;
534:       /* just do Galerkin grids */
535:       Mat          B,dA,dB;

537:       if (pc_gamg->Nlevels > 1) {
538:         PetscInt gl;
539:         /* currently only handle case where mat and pmat are the same on coarser levels */
540:         KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);
541:         /* (re)set to get dirty flag */
542:         KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);

544:         for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) {
545:           MatReuse reuse = MAT_INITIAL_MATRIX ;

547:           /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
548:           KSPGetOperators(mglevels[level]->smoothd,NULL,&B);
549:           if (B->product) {
550:             if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) {
551:               reuse = MAT_REUSE_MATRIX;
552:             }
553:           }
554:           if (reuse == MAT_INITIAL_MATRIX) MatDestroy(&mglevels[level]->A);
555:           if (reuse == MAT_REUSE_MATRIX) {
556:             PetscInfo(pc,"%s: RAP after first solve, reuse matrix level %D\n",((PetscObject)pc)->prefix,level);
557:           } else {
558:             PetscInfo(pc,"%s: RAP after first solve, new matrix level %D\n",((PetscObject)pc)->prefix,level);
559:           }
560:           PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0);
561:           MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B);
562:           PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0);
563:           if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B;
564:           KSPSetOperators(mglevels[level]->smoothd,B,B);
565:           dB   = B;
566:         }
567:       }

569:       PCSetUp_MG(pc);
570:       return 0;
571:     }
572:   }

574:   if (!pc_gamg->data) {
575:     if (pc_gamg->orig_data) {
576:       MatGetBlockSize(Pmat, &bs);
577:       MatGetLocalSize(Pmat, &qq, NULL);

579:       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
580:       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
581:       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;

583:       PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);
584:       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
585:     } else {
587:       pc_gamg->ops->createdefaultdata(pc,Pmat);
588:     }
589:   }

591:   /* cache original data for reuse */
592:   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
593:     PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);
594:     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
595:     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
596:     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
597:   }

599:   /* get basic dims */
600:   MatGetBlockSize(Pmat, &bs);
601:   MatGetSize(Pmat, &M, &N);

603:   MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info); /* global reduction */
604:   nnz0   = info.nz_used;
605:   nnztot = info.nz_used;
606:   PetscInfo(pc,"%s: level %D) N=%D, n data rows=%D, n data cols=%D, nnz/row (ave)=%d, np=%D\n",((PetscObject)pc)->prefix,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(int)(nnz0/(PetscReal)M+0.5),size);

608:   /* Get A_i and R_i */
609:   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
610:     pc_gamg->current_level = level;
612:     level1 = level + 1;
613:     PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);
614: #if defined(GAMG_STAGES)
615:     PetscLogStagePush(gamg_stages[level]);
616: #endif
617:     { /* construct prolongator */
618:       Mat              Gmat;
619:       PetscCoarsenData *agg_lists;
620:       Mat              Prol11;

622:       pc_gamg->ops->graph(pc,Aarr[level], &Gmat);
623:       pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);
624:       pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);

626:       /* could have failed to create new level */
627:       if (Prol11) {
628:         const char *prefix;
629:         char       addp[32];

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

634:         if (pc_gamg->ops->optprolongator) {
635:           /* smooth */
636:           pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);
637:         }

639:         if (pc_gamg->use_aggs_in_asm) {
640:           PetscInt bs;
641:           MatGetBlockSizes(Prol11, &bs, NULL);
642:           PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);
643:         }

645:         PCGetOptionsPrefix(pc,&prefix);
646:         MatSetOptionsPrefix(Prol11,prefix);
647:         PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level);
648:         MatAppendOptionsPrefix(Prol11,addp);
649:         /* Always generate the transpose with CUDA
650:            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
651:         MatSetOption(Prol11,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE);
652:         MatSetFromOptions(Prol11);
653:         Parr[level1] = Prol11;
654:       } else Parr[level1] = NULL; /* failed to coarsen */

656:       MatDestroy(&Gmat);
657:       PetscCDDestroy(agg_lists);
658:     } /* construct prolongator scope */
659:     PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);
660:     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
661:     if (!Parr[level1]) { /* failed to coarsen */
662:       PetscInfo(pc,"%s: Stop gridding, level %D\n",((PetscObject)pc)->prefix,level);
663: #if defined(GAMG_STAGES)
664:       PetscLogStagePop();
665: #endif
666:       break;
667:     }
668:     PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);
669:     MatGetSize(Parr[level1], &M, &N); /* N is next M, a loop test variables */
671:     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
672:     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
673:     pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);

675:     PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);
676:     MatGetSize(Aarr[level1], &M, &N); /* M is loop test variables */
677:     MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);
678:     nnztot += info.nz_used;
679:     PetscInfo(pc,"%s: %D) N=%D, n data cols=%D, nnz/row (ave)=%d, %D active pes\n",((PetscObject)pc)->prefix,level1,M,pc_gamg->data_cell_cols,(int)(info.nz_used/(PetscReal)M),nactivepe);

681: #if defined(GAMG_STAGES)
682:     PetscLogStagePop();
683: #endif
684:     /* stop if one node or one proc -- could pull back for singular problems */
685:     if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) {
686:       PetscInfo(pc,"%s: HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",((PetscObject)pc)->prefix,level,M/bs);
687:       level++;
688:       break;
689:     }
690:   } /* levels */
691:   PetscFree(pc_gamg->data);

693:   PetscInfo(pc,"%s: %D levels, grid complexity = %g\n",((PetscObject)pc)->prefix,level+1,nnztot/nnz0);
694:   pc_gamg->Nlevels = level + 1;
695:   fine_level       = level;
696:   PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);

698:   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */

700:     /* set default smoothers & set operators */
701:     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
702:       KSP smoother;
703:       PC  subpc;

705:       PCMGGetSmoother(pc, lidx, &smoother);
706:       KSPGetPC(smoother, &subpc);

708:       KSPSetNormType(smoother, KSP_NORM_NONE);
709:       /* set ops */
710:       KSPSetOperators(smoother, Aarr[level], Aarr[level]);
711:       PCMGSetInterpolation(pc, lidx, Parr[level+1]);

713:       /* set defaults */
714:       KSPSetType(smoother, KSPCHEBYSHEV);

716:       /* set blocks for ASM smoother that uses the 'aggregates' */
717:       if (pc_gamg->use_aggs_in_asm) {
718:         PetscInt sz;
719:         IS       *iss;

721:         sz   = nASMBlocksArr[level];
722:         iss  = ASMLocalIDsArr[level];
723:         PCSetType(subpc, PCASM);
724:         PCASMSetOverlap(subpc, 0);
725:         PCASMSetType(subpc,PC_ASM_BASIC);
726:         if (!sz) {
727:           IS       is;
728:           ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);
729:           PCASMSetLocalSubdomains(subpc, 1, NULL, &is);
730:           ISDestroy(&is);
731:         } else {
732:           PetscInt kk;
733:           PCASMSetLocalSubdomains(subpc, sz, NULL, iss);
734:           for (kk=0; kk<sz; kk++) {
735:             ISDestroy(&iss[kk]);
736:           }
737:           PetscFree(iss);
738:         }
739:         ASMLocalIDsArr[level] = NULL;
740:         nASMBlocksArr[level]  = 0;
741:       } else {
742:         PCSetType(subpc, PCJACOBI);
743:       }
744:     }
745:     {
746:       /* coarse grid */
747:       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
748:       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;

750:       PCMGGetSmoother(pc, lidx, &smoother);
751:       KSPSetOperators(smoother, Lmat, Lmat);
752:       if (!pc_gamg->use_parallel_coarse_grid_solver) {
753:         KSPSetNormType(smoother, KSP_NORM_NONE);
754:         KSPGetPC(smoother, &subpc);
755:         PCSetType(subpc, PCBJACOBI);
756:         PCSetUp(subpc);
757:         PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);
759:         KSPGetPC(k2[0],&pc2);
760:         PCSetType(pc2, PCLU);
761:         PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);
762:         KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);
763:         KSPSetType(k2[0], KSPPREONLY);
764:       }
765:     }

767:     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
768:     PetscObjectOptionsBegin((PetscObject)pc);
769:     PCSetFromOptions_MG(PetscOptionsObject,pc);
770:     PetscOptionsEnd();
771:     PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);

773:     /* setup cheby eigen estimates from SA */
774:     if (pc_gamg->use_sa_esteig) {
775:       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
776:         KSP       smoother;
777:         PetscBool ischeb;

779:         PCMGGetSmoother(pc, lidx, &smoother);
780:         PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb);
781:         if (ischeb) {
782:           KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data;

784:           // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG
785:           if (mg->max_eigen_DinvA[level] > 0) {
786:             // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it.
787:             // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.)
788:             PetscReal emax,emin;

790:             emin = mg->min_eigen_DinvA[level];
791:             emax = mg->max_eigen_DinvA[level];
792:             PetscInfo(pc,"%s: PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %D (N=%D) with emax = %g emin = %g\n",((PetscObject)pc)->prefix,level,Aarr[level]->rmap->N,(double)emax,(double)emin);
793:             cheb->emin_provided = emin;
794:             cheb->emax_provided = emax;
795:           }
796:         }
797:       }
798:     }

800:     PCSetUp_MG(pc);

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

810:     PetscInfo(pc,"%s: One level solver used (system is seen as DD). Using default solver.\n");
811:     PCMGGetSmoother(pc, 0, &smoother);
812:     KSPSetOperators(smoother, Aarr[0], Aarr[0]);
813:     KSPSetType(smoother, KSPPREONLY);
814:     PCSetUp_MG(pc);
815:   }
816:   return 0;
817: }

819: /* ------------------------------------------------------------------------- */
820: /*
821:  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
822:    that was created with PCCreate_GAMG().

824:    Input Parameter:
825: .  pc - the preconditioner context

827:    Application Interface Routine: PCDestroy()
828: */
829: PetscErrorCode PCDestroy_GAMG(PC pc)
830: {
831:   PC_MG          *mg     = (PC_MG*)pc->data;
832:   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;

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

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

848:    Logically Collective on PC

850:    Input Parameters:
851: +  pc - the preconditioner context
852: -  n - the number of equations

854:    Options Database Key:
855: .  -pc_gamg_process_eq_limit <limit> - set the limit

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

861:    Level: intermediate

863: .seealso: PCGAMGSetCoarseEqLim(), PCGAMGSetRankReductionFactors()
864: @*/
865: PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
866: {
868:   PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));
869:   return 0;
870: }

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

877:   if (n>0) pc_gamg->min_eq_proc = n;
878:   return 0;
879: }

881: /*@
882:    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.

884:  Collective on PC

886:    Input Parameters:
887: +  pc - the preconditioner context
888: -  n - maximum number of equations to aim for

890:    Options Database Key:
891: .  -pc_gamg_coarse_eq_limit <limit> - set the limit

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

896:    Level: intermediate

898: .seealso: PCGAMGSetProcEqLim(), PCGAMGSetRankReductionFactors()
899: @*/
900: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
901: {
903:   PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
904:   return 0;
905: }

907: static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
908: {
909:   PC_MG   *mg      = (PC_MG*)pc->data;
910:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

912:   if (n>0) pc_gamg->coarse_eq_limit = n;
913:   return 0;
914: }

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

919:    Collective on PC

921:    Input Parameters:
922: +  pc - the preconditioner context
923: -  n - PETSC_TRUE or PETSC_FALSE

925:    Options Database Key:
926: .  -pc_gamg_repartition <true,false> - turn on the repartitioning

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

931:    Level: intermediate

933: @*/
934: PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
935: {
937:   PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));
938:   return 0;
939: }

941: static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
942: {
943:   PC_MG   *mg      = (PC_MG*)pc->data;
944:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

946:   pc_gamg->repart = n;
947:   return 0;
948: }

950: /*@
951:    PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Chebyshev smoother

953:    Collective on PC

955:    Input Parameters:
956: +  pc - the preconditioner context
957: -  n - number of its

959:    Options Database Key:
960: .  -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate

962:    Notes:
963:    Smoothed aggregation constructs the smoothed prolongator $P = (I - \omega D^{-1} A) T$ where $T$ is the tentative prolongator and $D$ is the diagonal of $A$.
964:    Eigenvalue estimates (based on a few CG or GMRES iterations) are computed to choose $\omega$ so that this is a stable smoothing operation.
965:    If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused.
966:    This option is only used when the smoother uses Jacobi, and should be turned off if a different PCJacobiType is used.
967:    It became default in PETSc 3.17.

969:    Level: advanced

971: .seealso: KSPChebyshevSetEigenvalues(), KSPChebyshevEstEigSet()
972: @*/
973: PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
974: {
976:   PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));
977:   return 0;
978: }

980: static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
981: {
982:   PC_MG   *mg      = (PC_MG*)pc->data;
983:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

985:   pc_gamg->use_sa_esteig = n;
986:   return 0;
987: }

989: /*@
990:    PCGAMGSetEigenvalues - Set eigenvalues

992:    Collective on PC

994:    Input Parameters:
995: +  pc - the preconditioner context
996: -  emax - max eigenvalue
997: .  emin - min eigenvalue

999:    Options Database Key:
1000: .  -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues

1002:    Level: intermediate

1004: .seealso: PCGAMGSetUseSAEstEig()
1005: @*/
1006: PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin)
1007: {
1009:   PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));
1010:   return 0;
1011: }

1013: static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin)
1014: {
1015:   PC_MG          *mg      = (PC_MG*)pc->data;
1016:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1020:   pc_gamg->emax = emax;
1021:   pc_gamg->emin = emin;
1022:   return 0;
1023: }

1025: /*@
1026:    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner

1028:    Collective on PC

1030:    Input Parameters:
1031: +  pc - the preconditioner context
1032: -  n - PETSC_TRUE or PETSC_FALSE

1034:    Options Database Key:
1035: .  -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation

1037:    Level: intermediate

1039:    Notes:
1040:     May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1041:     rebuilding the preconditioner quicker.

1043: @*/
1044: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1045: {
1047:   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
1048:   return 0;
1049: }

1051: static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1052: {
1053:   PC_MG   *mg      = (PC_MG*)pc->data;
1054:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1056:   pc_gamg->reuse_prol = n;
1057:   return 0;
1058: }

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

1063:    Collective on PC

1065:    Input Parameters:
1066: +  pc - the preconditioner context
1067: -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not

1069:    Options Database Key:
1070: .  -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains

1072:    Level: intermediate

1074: @*/
1075: PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1076: {
1078:   PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));
1079:   return 0;
1080: }

1082: static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1083: {
1084:   PC_MG   *mg      = (PC_MG*)pc->data;
1085:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1087:   pc_gamg->use_aggs_in_asm = flg;
1088:   return 0;
1089: }

1091: /*@
1092:    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver

1094:    Collective on PC

1096:    Input Parameters:
1097: +  pc - the preconditioner context
1098: -  flg - PETSC_TRUE to not force coarse grid onto one processor

1100:    Options Database Key:
1101: .  -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver

1103:    Level: intermediate

1105: .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids()
1106: @*/
1107: PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1108: {
1110:   PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));
1111:   return 0;
1112: }

1114: static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1115: {
1116:   PC_MG   *mg      = (PC_MG*)pc->data;
1117:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1119:   pc_gamg->use_parallel_coarse_grid_solver = flg;
1120:   return 0;
1121: }

1123: /*@
1124:    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU

1126:    Collective on PC

1128:    Input Parameters:
1129: +  pc - the preconditioner context
1130: -  flg - PETSC_TRUE to pin coarse grids to CPU

1132:    Options Database Key:
1133: .  -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU

1135:    Level: intermediate

1137: .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve()
1138: @*/
1139: PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1140: {
1142:   PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));
1143:   return 0;
1144: }

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

1151:   pc_gamg->cpu_pin_coarse_grids = flg;
1152:   return 0;
1153: }

1155: /*@
1156:    PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type)

1158:    Collective on PC

1160:    Input Parameters:
1161: +  pc - the preconditioner context
1162: -  flg - Layout type

1164:    Options Database Key:
1165: .  -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering

1167:    Level: intermediate

1169: .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids()
1170: @*/
1171: PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1172: {
1174:   PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));
1175:   return 0;
1176: }

1178: static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1179: {
1180:   PC_MG   *mg      = (PC_MG*)pc->data;
1181:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1183:   pc_gamg->layout_type = flg;
1184:   return 0;
1185: }

1187: /*@
1188:    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use

1190:    Not collective on PC

1192:    Input Parameters:
1193: +  pc - the preconditioner
1194: -  n - the maximum number of levels to use

1196:    Options Database Key:
1197: .  -pc_mg_levels <n> - set the maximum number of levels to allow

1199:    Level: intermediate

1201: @*/
1202: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1203: {
1205:   PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
1206:   return 0;
1207: }

1209: static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1210: {
1211:   PC_MG   *mg      = (PC_MG*)pc->data;
1212:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;

1214:   pc_gamg->Nlevels = n;
1215:   return 0;
1216: }

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

1221:    Not collective on PC

1223:    Input Parameters:
1224: +  pc - the preconditioner context
1225: .  v - 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
1226: -  n - number of threshold values provided in array

1228:    Options Database Key:
1229: .  -pc_gamg_threshold <threshold> - the threshold to drop edges

1231:    Notes:
1232:     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.
1233:     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.

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

1239:    Level: intermediate

1241: .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
1242: @*/
1243: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1244: {
1247:   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));
1248:   return 0;
1249: }

1251: static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1252: {
1253:   PC_MG   *mg      = (PC_MG*)pc->data;
1254:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1255:   PetscInt i;
1256:   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1257:   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;
1258:   return 0;
1259: }

1261: /*@
1262:    PCGAMGSetRankReductionFactors - Set manual schedule for process reduction on coarse grids

1264:    Collective on PC

1266:    Input Parameters:
1267: +  pc - the preconditioner context
1268: .  v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda
1269: -  n - number of values provided in array

1271:    Options Database Key:
1272: .  -pc_gamg_rank_reduction_factors <factors> - provide the schedule

1274:    Level: intermediate

1276: .seealso: PCGAMGSetProcEqLim(), PCGAMGSetCoarseEqLim()
1277: @*/
1278: PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1279: {
1282:   PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n));
1283:   return 0;
1284: }

1286: static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1287: {
1288:   PC_MG   *mg      = (PC_MG*)pc->data;
1289:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1290:   PetscInt i;
1291:   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1292:   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
1293:   return 0;
1294: }

1296: /*@
1297:    PCGAMGSetThresholdScale - Relative threshold reduction at each level

1299:    Not collective on PC

1301:    Input Parameters:
1302: +  pc - the preconditioner context
1303: -  scale - the threshold value reduction, usually < 1.0

1305:    Options Database Key:
1306: .  -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level

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

1312:    Level: advanced

1314: .seealso: PCGAMGSetThreshold()
1315: @*/
1316: PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1317: {
1319:   PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));
1320:   return 0;
1321: }

1323: static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1324: {
1325:   PC_MG   *mg      = (PC_MG*)pc->data;
1326:   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1327:   pc_gamg->threshold_scale = v;
1328:   return 0;
1329: }

1331: /*@C
1332:    PCGAMGSetType - Set solution method

1334:    Collective on PC

1336:    Input Parameters:
1337: +  pc - the preconditioner context
1338: -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL

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

1343:    Level: intermediate

1345: .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1346: @*/
1347: PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1348: {
1350:   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1351:   return 0;
1352: }

1354: /*@C
1355:    PCGAMGGetType - Get solution method

1357:    Collective on PC

1359:    Input Parameter:
1360: .  pc - the preconditioner context

1362:    Output Parameter:
1363: .  type - the type of algorithm used

1365:    Level: intermediate

1367: .seealso: PCGAMGSetType(), PCGAMGType
1368: @*/
1369: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1370: {
1372:   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1373:   return 0;
1374: }

1376: static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1377: {
1378:   PC_MG          *mg      = (PC_MG*)pc->data;
1379:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;

1381:   *type = pc_gamg->type;
1382:   return 0;
1383: }

1385: static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1386: {
1387:   PC_MG          *mg      = (PC_MG*)pc->data;
1388:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1389:   PetscErrorCode (*r)(PC);

1391:   pc_gamg->type = type;
1392:   PetscFunctionListFind(GAMGList,type,&r);
1394:   if (pc_gamg->ops->destroy) {
1395:     (*pc_gamg->ops->destroy)(pc);
1396:     PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));
1397:     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1398:     /* cleaning up common data in pc_gamg - this should disapear someday */
1399:     pc_gamg->data_cell_cols = 0;
1400:     pc_gamg->data_cell_rows = 0;
1401:     pc_gamg->orig_data_cell_cols = 0;
1402:     pc_gamg->orig_data_cell_rows = 0;
1403:     PetscFree(pc_gamg->data);
1404:     pc_gamg->data_sz = 0;
1405:   }
1406:   PetscFree(pc_gamg->gamg_type_name);
1407:   PetscStrallocpy(type,&pc_gamg->gamg_type_name);
1408:   (*r)(pc);
1409:   return 0;
1410: }

1412: static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
1413: {
1414:   PC_MG          *mg      = (PC_MG*)pc->data;
1415:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1416:   PetscReal       gc=0, oc=0;

1418:   PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n");
1419:   PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level =");
1420:   for (PetscInt i=0;i<mg->nlevels; i++) PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);
1421:   PetscViewerASCIIPrintf(viewer,"\n");
1422:   PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);
1423:   if (pc_gamg->use_aggs_in_asm) {
1424:     PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n");
1425:   }
1426:   if (pc_gamg->use_parallel_coarse_grid_solver) {
1427:     PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");
1428:   }
1429: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1430:   if (pc_gamg->cpu_pin_coarse_grids) {
1431:     /* PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n"); */
1432:   }
1433: #endif
1434:   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
1435:   /*   PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n"); */
1436:   /* } else { */
1437:   /*   PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n"); */
1438:   /* } */
1439:   if (pc_gamg->ops->view) {
1440:     (*pc_gamg->ops->view)(pc,viewer);
1441:   }
1442:   PCMGGetGridComplexity(pc,&gc,&oc);
1443:   PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g    operator = %g\n",gc,oc);
1444:   return 0;
1445: }

1447: PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
1448: {
1449:   PC_MG          *mg      = (PC_MG*)pc->data;
1450:   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1451:   PetscBool      flag;
1452:   MPI_Comm       comm;
1453:   char           prefix[256],tname[32];
1454:   PetscInt       i,n;
1455:   const char     *pcpre;
1456:   static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL};
1457:   PetscObjectGetComm((PetscObject)pc,&comm);
1458:   PetscOptionsHead(PetscOptionsObject,"GAMG options");
1459:   PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);
1460:   if (flag) {
1461:     PCGAMGSetType(pc,tname);
1462:   }
1463:   PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);
1464:   PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",pc_gamg->use_sa_esteig,&pc_gamg->use_sa_esteig,NULL);
1465:   PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);
1466:   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);
1467:   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);
1468:   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);
1469:   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);
1470:   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);
1471:   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);
1472:   PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL);
1473:   n = PETSC_MG_MAXLEVELS;
1474:   PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);
1475:   if (!flag || n < PETSC_MG_MAXLEVELS) {
1476:     if (!flag) n = 1;
1477:     i = n;
1478:     do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS);
1479:   }
1480:   n = PETSC_MG_MAXLEVELS;
1481:   PetscOptionsIntArray("-pc_gamg_rank_reduction_factors","Manual schedule of coarse grid reduction factors that overrides internal heuristics (0 for first reduction puts one process/device)","PCGAMGSetRankReductionFactors",pc_gamg->level_reduction_factors,&n,&flag);
1482:   if (!flag) i = 0;
1483:   else i = n;
1484:   do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS);
1485:   PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);
1486:   {
1487:     PetscReal eminmax[2] = {0., 0.};
1488:     n = 2;
1489:     PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag);
1490:     if (flag) {
1492:       PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);
1493:     }
1494:   }
1495:   /* set options for subtype */
1496:   if (pc_gamg->ops->setfromoptions) (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);

1498:   PCGetOptionsPrefix(pc, &pcpre);
1499:   PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");
1500:   PetscOptionsTail();
1501:   return 0;
1502: }

1504: /* -------------------------------------------------------------------------- */
1505: /*MC
1506:      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner

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

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

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

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

1536:   Level: intermediate

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

1542: PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1543: {
1544:   PC_GAMG *pc_gamg;
1545:   PC_MG   *mg;

1547:    /* register AMG type */
1548:   PCGAMGInitializePackage();

1550:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1551:   PCSetType(pc, PCMG);
1552:   PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);

1554:   /* create a supporting struct and attach it to pc */
1555:   PetscNewLog(pc,&pc_gamg);
1556:   PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);
1557:   mg           = (PC_MG*)pc->data;
1558:   mg->innerctx = pc_gamg;

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

1562:   /* these should be in subctx but repartitioning needs simple arrays */
1563:   pc_gamg->data_sz = 0;
1564:   pc_gamg->data    = NULL;

1566:   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1567:   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1568:   pc->ops->setup          = PCSetUp_GAMG;
1569:   pc->ops->reset          = PCReset_GAMG;
1570:   pc->ops->destroy        = PCDestroy_GAMG;
1571:   mg->view                = PCView_GAMG;

1573:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);
1574:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);
1575:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);
1576:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);
1577:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);
1578:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG);
1579:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG);
1580:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);
1581:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);
1582:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);
1583:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);
1584:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);
1585:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);
1586:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG);
1587:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);
1588:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);
1589:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);
1590:   PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);
1591:   pc_gamg->repart           = PETSC_FALSE;
1592:   pc_gamg->reuse_prol       = PETSC_FALSE;
1593:   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1594:   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1595:   pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1596:   pc_gamg->layout_type      = PCGAMG_LAYOUT_SPREAD;
1597:   pc_gamg->min_eq_proc      = 50;
1598:   pc_gamg->coarse_eq_limit  = 50;
1599:   PetscArrayzero(pc_gamg->threshold,PETSC_MG_MAXLEVELS);
1600:   pc_gamg->threshold_scale = 1.;
1601:   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
1602:   pc_gamg->current_level    = 0; /* don't need to init really */
1603:   pc_gamg->use_sa_esteig    = PETSC_TRUE;
1604:   pc_gamg->emin             = 0;
1605:   pc_gamg->emax             = 0;

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

1609:   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1610:   PCGAMGSetType(pc,PCGAMGAGG);
1611:   return 0;
1612: }

1614: /*@C
1615:  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
1616:     from PCInitializePackage().

1618:  Level: developer

1620:  .seealso: PetscInitialize()
1621: @*/
1622: PetscErrorCode PCGAMGInitializePackage(void)
1623: {
1624:   PetscInt       l;

1626:   if (PCGAMGPackageInitialized) return 0;
1627:   PCGAMGPackageInitialized = PETSC_TRUE;
1628:   PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);
1629:   PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);
1630:   PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);
1631:   PetscRegisterFinalize(PCGAMGFinalizePackage);

1633:   /* general events */
1634:   PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);
1635:   PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);
1636:   PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);
1637:   PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);
1638:   PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);
1639:   PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);
1640:   PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);

1642:   PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);
1643:   PetscLogEventRegister("  Create Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);
1644:   PetscLogEventRegister("  Filter Graph", PC_CLASSID, &petsc_gamg_setup_events[SET16]);
1645:   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1646:   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1647:   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1648:   PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);
1649:   PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);
1650:   PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);
1651:   PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);
1652:   PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);
1653:   PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);
1654:   PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);
1655:   PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);
1656:   PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);
1657:   PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);
1658:   PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);
1659:   PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);
1660:   for (l=0;l<PETSC_MG_MAXLEVELS;l++) {
1661:     char ename[32];

1663:     PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02d",l);
1664:     PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]);
1665:     PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02d",l);
1666:     PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]);
1667:     PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02d",l);
1668:     PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]);
1669:   }
1670:   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1671:   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1672:   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1673:   /* create timer stages */
1674: #if defined(GAMG_STAGES)
1675:   {
1676:     char     str[32];
1677:     PetscInt lidx;
1678:     sprintf(str,"MG Level %d (finest)",0);
1679:     PetscLogStageRegister(str, &gamg_stages[0]);
1680:     for (lidx=1; lidx<9; lidx++) {
1681:       sprintf(str,"MG Level %d",(int)lidx);
1682:       PetscLogStageRegister(str, &gamg_stages[lidx]);
1683:     }
1684:   }
1685: #endif
1686:   return 0;
1687: }

1689: /*@C
1690:  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1691:     called from PetscFinalize() automatically.

1693:  Level: developer

1695:  .seealso: PetscFinalize()
1696: @*/
1697: PetscErrorCode PCGAMGFinalizePackage(void)
1698: {
1699:   PCGAMGPackageInitialized = PETSC_FALSE;
1700:   PetscFunctionListDestroy(&GAMGList);
1701:   return 0;
1702: }

1704: /*@C
1705:  PCGAMGRegister - Register a PCGAMG implementation.

1707:  Input Parameters:
1708:  + type - string that will be used as the name of the GAMG type.
1709:  - create - function for creating the gamg context.

1711:   Level: advanced

1713:  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1714: @*/
1715: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1716: {
1717:   PCGAMGInitializePackage();
1718:   PetscFunctionListAdd(&GAMGList,type,create);
1719:   return 0;
1720: }