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[GAMG_NUM_SET];
 17: PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3];

 19: // #define GAMG_STAGES
 20: #if defined(GAMG_STAGES)
 21: static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
 22: #endif

 24: static PetscFunctionList GAMGList = NULL;
 25: static PetscBool         PCGAMGPackageInitialized;

 27: PetscErrorCode PCReset_GAMG(PC pc)
 28: {
 29:   PC_MG   *mg      = (PC_MG *)pc->data;
 30:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

 32:   PetscFree(pc_gamg->data);
 33:   pc_gamg->data_sz = 0;
 34:   PetscFree(pc_gamg->orig_data);
 35:   for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS; level++) {
 36:     mg->min_eigen_DinvA[level] = 0;
 37:     mg->max_eigen_DinvA[level] = 0;
 38:   }
 39:   pc_gamg->emin = 0;
 40:   pc_gamg->emax = 0;
 41:   return 0;
 42: }

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

 48:    Input Parameter:
 49:    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
 50:           'pc_gamg->data_sz' are changed via repartitioning/reduction.
 51:    . Amat_fine - matrix on this fine (k) level
 52:    . cr_bs - coarse block size
 53:    In/Output Parameter:
 54:    . a_P_inout - prolongation operator to the next level (k-->k-1)
 55:    . a_nactive_proc - number of active procs
 56:    Output Parameter:
 57:    . a_Amat_crs - coarse matrix that is created (k-1)
 58: */
 59: 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)
 60: {
 61:   PC_MG      *mg      = (PC_MG *)pc->data;
 62:   PC_GAMG    *pc_gamg = (PC_GAMG *)mg->innerctx;
 63:   Mat         Cmat, Pold = *a_P_inout;
 64:   MPI_Comm    comm;
 65:   PetscMPIInt rank, size, new_size, nactive = *a_nactive_proc;
 66:   PetscInt    ncrs_eq, ncrs, f_bs;

 68:   PetscObjectGetComm((PetscObject)Amat_fine, &comm);
 69:   MPI_Comm_rank(comm, &rank);
 70:   MPI_Comm_size(comm, &size);
 71:   MatGetBlockSize(Amat_fine, &f_bs);
 72:   PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0);
 73:   PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0);
 74:   MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);
 75:   PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0);
 76:   PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0);

 78:   if (Pcolumnperm) *Pcolumnperm = NULL;

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

202:         PetscMalloc2(ncrs, &d_nnz, ncrs, &o_nnz);
203:         MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);
204:         MatGetSize(Cmat, &M, &N);
205:         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
206:           MatGetRow(Cmat, Ii, &ncols, NULL, NULL);
207:           d_nnz[jj] = ncols / cr_bs;
208:           o_nnz[jj] = ncols / cr_bs;
209:           MatRestoreRow(Cmat, Ii, &ncols, NULL, NULL);
210:           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
211:           if (o_nnz[jj] > (M / cr_bs - ncrs)) o_nnz[jj] = M / cr_bs - ncrs;
212:         }

214:         MatGetType(Amat_fine, &mtype);
215:         MatCreate(comm, &tMat);
216:         MatSetSizes(tMat, ncrs, ncrs, PETSC_DETERMINE, PETSC_DETERMINE);
217:         MatSetType(tMat, mtype);
218:         MatSeqAIJSetPreallocation(tMat, 0, d_nnz);
219:         MatMPIAIJSetPreallocation(tMat, 0, d_nnz, 0, o_nnz);
220:         PetscFree2(d_nnz, o_nnz);

222:         for (ii = Istart_crs; ii < Iend_crs; ii++) {
223:           PetscInt dest_row = ii / cr_bs;
224:           MatGetRow(Cmat, ii, &ncols, &idx, &vals);
225:           for (jj = 0; jj < ncols; jj++) {
226:             PetscInt    dest_col = idx[jj] / cr_bs;
227:             PetscScalar v        = 1.0;
228:             MatSetValues(tMat, 1, &dest_row, 1, &dest_col, &v, ADD_VALUES);
229:           }
230:           MatRestoreRow(Cmat, ii, &ncols, &idx, &vals);
231:         }
232:         MatAssemblyBegin(tMat, MAT_FINAL_ASSEMBLY);
233:         MatAssemblyEnd(tMat, MAT_FINAL_ASSEMBLY);

235:         if (llev++ == -1) {
236:           PetscViewer viewer;
237:           char        fname[32];
238:           PetscSNPrintf(fname, sizeof(fname), "part_mat_%" PetscInt_FMT ".mat", llev);
239:           PetscViewerBinaryOpen(comm, fname, FILE_MODE_WRITE, &viewer);
240:           MatView(tMat, viewer);
241:           PetscViewerDestroy(&viewer);
242:         }
243:         MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);
244:         MatDestroy(&tMat);
245:       } /* create 'adj' */

247:       { /* partition: get newproc_idx */
248:         char            prefix[256];
249:         const char     *pcpre;
250:         const PetscInt *is_idx;
251:         MatPartitioning mpart;
252:         IS              proc_is;

254:         MatPartitioningCreate(comm, &mpart);
255:         MatPartitioningSetAdjacency(mpart, adj);
256:         PCGetOptionsPrefix(pc, &pcpre);
257:         PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : "");
258:         PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix);
259:         MatPartitioningSetFromOptions(mpart);
260:         MatPartitioningSetNParts(mpart, new_size);
261:         MatPartitioningApply(mpart, &proc_is);
262:         MatPartitioningDestroy(&mpart);

264:         /* collect IS info */
265:         PetscMalloc1(ncrs_eq, &newproc_idx);
266:         ISGetIndices(proc_is, &is_idx);
267:         for (kk = jj = 0; kk < nloc_old; kk++) {
268:           for (ii = 0; ii < cr_bs; ii++, jj++) { newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ }
269:         }
270:         ISRestoreIndices(proc_is, &is_idx);
271:         ISDestroy(&proc_is);
272:       }
273:       MatDestroy(&adj);

275:       ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);
276:       PetscFree(newproc_idx);
277:       PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0);
278:     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
279:       PetscInt targetPE;
281:       PetscInfo(pc, "%s: Number of equations (loc) %" PetscInt_FMT " with simple aggregation\n", ((PetscObject)pc)->prefix, ncrs_eq);
282:       targetPE = (rank / rfactor) * expand_factor;
283:       ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);
284:     } /* end simple 'is_eq_newproc' */

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

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

357:       pc_gamg->data_sz = node_data_sz * ncrs_new;
358:       strideNew        = ncrs_new * ndata_rows;

360:       VecGetArray(dest_crd, &array);
361:       for (jj = 0; jj < ndata_cols; jj++) {
362:         for (ii = 0; ii < ncrs_new; ii++) {
363:           for (kk = 0; kk < ndata_rows; kk++) {
364:             PetscInt ix = ii * ndata_rows + kk + jj * strideNew, jx = ii * node_data_sz + kk * ndata_cols + jj;
365:             pc_gamg->data[ix] = PetscRealPart(array[jx]);
366:           }
367:         }
368:       }
369:       VecRestoreArray(dest_crd, &array);
370:       VecDestroy(&dest_crd);
371:     }
372:     /* move A and P (columns) with new layout */
373:     /* PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0); */
374:     /*
375:       Invert for MatCreateSubMatrix
376:     */
377:     ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);
378:     ISSort(new_eq_indices); /* is this needed? */
379:     ISSetBlockSize(new_eq_indices, cr_bs);
380:     if (is_eq_num != is_eq_num_prim) { ISDestroy(&is_eq_num_prim); /* could be same as 'is_eq_num' */ }
381:     if (Pcolumnperm) {
382:       PetscObjectReference((PetscObject)new_eq_indices);
383:       *Pcolumnperm = new_eq_indices;
384:     }
385:     ISDestroy(&is_eq_num);
386:     /* PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0); */
387:     /* PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0); */

389:     /* 'a_Amat_crs' output */
390:     {
391:       Mat       mat;
392:       PetscBool isset, isspd, isher;
393: #if !defined(PETSC_USE_COMPLEX)
394:       PetscBool issym;
395: #endif

397:       MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);
398:       MatIsSPDKnown(Cmat, &isset, &isspd); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ?
399:       if (isset) MatSetOption(mat, MAT_SPD, isspd);
400:       else {
401:         MatIsHermitianKnown(Cmat, &isset, &isher);
402:         if (isset) MatSetOption(mat, MAT_HERMITIAN, isher);
403:         else {
404: #if !defined(PETSC_USE_COMPLEX)
405:           MatIsSymmetricKnown(Cmat, &isset, &issym);
406:           if (isset) MatSetOption(mat, MAT_SYMMETRIC, issym);
407: #endif
408:         }
409:       }
410:       *a_Amat_crs = mat;
411:     }
412:     MatDestroy(&Cmat);

414:     /* PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0); */
415:     /* prolongator */
416:     {
417:       IS       findices;
418:       PetscInt Istart, Iend;
419:       Mat      Pnew;

421:       MatGetOwnershipRange(Pold, &Istart, &Iend);
422:       /* PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0); */
423:       ISCreateStride(comm, Iend - Istart, Istart, 1, &findices);
424:       ISSetBlockSize(findices, f_bs);
425:       MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);
426:       ISDestroy(&findices);
427:       MatSetOption(Pnew, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE);

429:       /* PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0); */
430:       MatDestroy(a_P_inout);

432:       /* output - repartitioned */
433:       *a_P_inout = Pnew;
434:     }
435:     ISDestroy(&new_eq_indices);

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

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

463: // used in GEO
464: PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat *Gmat2)
465: {
466:   const char *prefix;
467:   char        addp[32];
468:   PC_MG      *mg      = (PC_MG *)a_pc->data;
469:   PC_GAMG    *pc_gamg = (PC_GAMG *)mg->innerctx;

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

493: /*
494:    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
495:                     by setting data structures and options.

497:    Input Parameter:
498: .  pc - the preconditioner context

500: */
501: PetscErrorCode PCSetUp_GAMG(PC pc)
502: {
503:   PC_MG         *mg      = (PC_MG *)pc->data;
504:   PC_GAMG       *pc_gamg = (PC_GAMG *)mg->innerctx;
505:   Mat            Pmat    = pc->pmat;
506:   PetscInt       fine_level, level, level1, bs, M, N, qq, lidx, nASMBlocksArr[PETSC_MG_MAXLEVELS];
507:   MPI_Comm       comm;
508:   PetscMPIInt    rank, size, nactivepe;
509:   Mat            Aarr[PETSC_MG_MAXLEVELS], Parr[PETSC_MG_MAXLEVELS];
510:   IS            *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
511:   PetscLogDouble nnz0 = 0., nnztot = 0.;
512:   MatInfo        info;
513:   PetscBool      is_last = PETSC_FALSE;

515:   PetscObjectGetComm((PetscObject)pc, &comm);
516:   MPI_Comm_rank(comm, &rank);
517:   MPI_Comm_size(comm, &size);
518:   PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0);
519:   if (pc->setupcalled) {
520:     if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) {
521:       /* reset everything */
522:       PCReset_MG(pc);
523:       pc->setupcalled = 0;
524:     } else {
525:       PC_MG_Levels **mglevels = mg->levels;
526:       /* just do Galerkin grids */
527:       Mat B, dA, dB;
528:       if (pc_gamg->Nlevels > 1) {
529:         PetscInt gl;
530:         /* currently only handle case where mat and pmat are the same on coarser levels */
531:         KSPGetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, &dA, &dB);
532:         /* (re)set to get dirty flag */
533:         KSPSetOperators(mglevels[pc_gamg->Nlevels - 1]->smoothd, dA, dB);

535:         for (level = pc_gamg->Nlevels - 2, gl = 0; level >= 0; level--, gl++) {
536:           MatReuse reuse = MAT_INITIAL_MATRIX;
537: #if defined(GAMG_STAGES)
538:           PetscLogStagePush(gamg_stages[gl]);
539: #endif
540:           /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
541:           KSPGetOperators(mglevels[level]->smoothd, NULL, &B);
542:           if (B->product) {
543:             if (B->product->A == dB && B->product->B == mglevels[level + 1]->interpolate) reuse = MAT_REUSE_MATRIX;
544:           }
545:           if (reuse == MAT_INITIAL_MATRIX) MatDestroy(&mglevels[level]->A);
546:           if (reuse == MAT_REUSE_MATRIX) {
547:             PetscInfo(pc, "%s: RAP after first solve, reuse matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level);
548:           } else {
549:             PetscInfo(pc, "%s: RAP after first solve, new matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level);
550:           }
551:           PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0);
552:           MatPtAP(dB, mglevels[level + 1]->interpolate, reuse, PETSC_DEFAULT, &B);
553:           PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1], 0, 0, 0, 0);
554:           if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B;
555:           KSPSetOperators(mglevels[level]->smoothd, B, B);
556:           dB = B;
557: #if defined(GAMG_STAGES)
558:           PetscLogStagePop();
559: #endif
560:         }
561:       }

563:       PCSetUp_MG(pc);
564:       PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0);
565:       return 0;
566:     }
567:   }

569:   if (!pc_gamg->data) {
570:     if (pc_gamg->orig_data) {
571:       MatGetBlockSize(Pmat, &bs);
572:       MatGetLocalSize(Pmat, &qq, NULL);

574:       pc_gamg->data_sz        = (qq / bs) * pc_gamg->orig_data_cell_rows * pc_gamg->orig_data_cell_cols;
575:       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
576:       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;

578:       PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);
579:       for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
580:     } else {
582:       pc_gamg->ops->createdefaultdata(pc, Pmat);
583:     }
584:   }

586:   /* cache original data for reuse */
587:   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
588:     PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);
589:     for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
590:     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
591:     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
592:   }

594:   /* get basic dims */
595:   MatGetBlockSize(Pmat, &bs);
596:   MatGetSize(Pmat, &M, &N);

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

603:   /* Get A_i and R_i */
604:   for (level = 0, Aarr[0] = Pmat, nactivepe = size; level < (pc_gamg->Nlevels - 1) && (!level || M > pc_gamg->coarse_eq_limit); level++) {
605:     pc_gamg->current_level = level;
607:     level1 = level + 1;
608: #if defined(GAMG_STAGES)
609:     if (!gamg_stages[level]) {
610:       char str[32];
611:       sprintf(str, "GAMG Level %d", (int)level);
612:       PetscLogStageRegister(str, &gamg_stages[level]);
613:     }
614:     PetscLogStagePush(gamg_stages[level]);
615: #endif
616:     { /* construct prolongator */
617:       Mat               Gmat;
618:       PetscCoarsenData *agg_lists;
619:       Mat               Prol11;

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

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

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

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

638:         if (pc_gamg->use_aggs_in_asm) {
639:           PetscInt bs;
640:           MatGetBlockSizes(Prol11, &bs, NULL); // not timed directly, ugly, could remove, but good ASM method
641:           PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);
642:         }

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

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

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

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

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

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

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

703:       PCMGGetSmoother(pc, lidx, &smoother);
704:       KSPGetPC(smoother, &subpc);

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

711:       /* set defaults */
712:       KSPSetType(smoother, KSPCHEBYSHEV);

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

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

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

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

772:     /* set cheby eigen estimates from SA to use in the solver */
773:     if (pc_gamg->use_sa_esteig) {
774:       for (lidx = 1, level = pc_gamg->Nlevels - 2; level >= 0; lidx++, level--) {
775:         KSP       smoother;
776:         PetscBool ischeb;

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

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

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

799:     PCSetUp_MG(pc);

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

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

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

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

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

833:   PCReset_GAMG(pc);
834:   if (pc_gamg->ops->destroy) (*pc_gamg->ops->destroy)(pc);
835:   PetscFree(pc_gamg->ops);
836:   PetscFree(pc_gamg->gamg_type_name);
837:   PetscFree(pc_gamg);
838:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", NULL);
839:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", NULL);
840:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", NULL);
841:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", NULL);
842:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", NULL);
843:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", NULL);
844:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", NULL);
845:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", NULL);
846:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", NULL);
847:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", NULL);
848:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", NULL);
849:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", NULL);
850:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", NULL);
851:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", NULL);
852:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", NULL);
853:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", NULL);
854:   PCDestroy_MG(pc);
855:   return 0;
856: }

858: /*@
859:    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction in `PCGAMG`

861:    Logically Collective on pc

863:    Input Parameters:
864: +  pc - the preconditioner context
865: -  n - the number of equations

867:    Options Database Key:
868: .  -pc_gamg_process_eq_limit <limit> - set the limit

870:    Note:
871:    `PCGAMG` will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
872:    that has degrees of freedom

874:    Level: intermediate

876: .seealso: `PCGAMG`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`
877: @*/
878: PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n)
879: {
881:   PetscTryMethod(pc, "PCGAMGSetProcEqLim_C", (PC, PetscInt), (pc, n));
882:   return 0;
883: }

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

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

894: /*@
895:    PCGAMGSetCoarseEqLim - Set maximum number of equations on the coarsest grid of `PCGAMG`

897:    Collective on pc

899:    Input Parameters:
900: +  pc - the preconditioner context
901: -  n - maximum number of equations to aim for

903:    Options Database Key:
904: .  -pc_gamg_coarse_eq_limit <limit> - set the limit

906:    Note:
907:    For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
908:    has less than 1000 unknowns.

910:    Level: intermediate

912: .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`
913: @*/
914: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
915: {
917:   PetscTryMethod(pc, "PCGAMGSetCoarseEqLim_C", (PC, PetscInt), (pc, n));
918:   return 0;
919: }

921: static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
922: {
923:   PC_MG   *mg      = (PC_MG *)pc->data;
924:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

926:   if (n > 0) pc_gamg->coarse_eq_limit = n;
927:   return 0;
928: }

930: /*@
931:    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids when reducing the number of MPI ranks to use

933:    Collective on pc

935:    Input Parameters:
936: +  pc - the preconditioner context
937: -  n - `PETSC_TRUE` or `PETSC_FALSE`

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

942:    Note:
943:    This will generally improve the loading balancing of the work on each level

945:    Level: intermediate

947: .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`
948: @*/
949: PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
950: {
952:   PetscTryMethod(pc, "PCGAMGSetRepartition_C", (PC, PetscBool), (pc, n));
953:   return 0;
954: }

956: static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
957: {
958:   PC_MG   *mg      = (PC_MG *)pc->data;
959:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

961:   pc_gamg->repart = n;
962:   return 0;
963: }

965: /*@
966:    PCGAMGSetUseSAEstEig - Use the eigen estimate from smoothed aggregation for the Chebyshev smoother during the solution process

968:    Collective on pc

970:    Input Parameters:
971: +  pc - the preconditioner context
972: -  n - number of its

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

977:    Notes:
978:    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$.
979:    Eigenvalue estimates (based on a few `PCCG` or `PCGMRES` iterations) are computed to choose $\omega$ so that this is a stable smoothing operation.
980:    If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused during the solution process
981:    This option is only used when the smoother uses Jacobi, and should be turned off if a different `PCJacobiType` is used.
982:    It became default in PETSc 3.17.

984:    Level: advanced

986: .seealso: `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`
987: @*/
988: PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
989: {
991:   PetscTryMethod(pc, "PCGAMGSetUseSAEstEig_C", (PC, PetscBool), (pc, n));
992:   return 0;
993: }

995: static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
996: {
997:   PC_MG   *mg      = (PC_MG *)pc->data;
998:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1000:   pc_gamg->use_sa_esteig = n;
1001:   return 0;
1002: }

1004: /*@
1005:    PCGAMGSetEigenvalues - Set WHAT eigenvalues WHY?

1007:    Collective on pc

1009:    Input Parameters:
1010: +  pc - the preconditioner context
1011: -  emax - max eigenvalue
1012: .  emin - min eigenvalue

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

1017:    Level: intermediate

1019: .seealso: `PCGAMG`, `PCGAMGSetUseSAEstEig()`
1020: @*/
1021: PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax, PetscReal emin)
1022: {
1024:   PetscTryMethod(pc, "PCGAMGSetEigenvalues_C", (PC, PetscReal, PetscReal), (pc, emax, emin));
1025:   return 0;
1026: }

1028: static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc, PetscReal emax, PetscReal emin)
1029: {
1030:   PC_MG   *mg      = (PC_MG *)pc->data;
1031:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1035:   pc_gamg->emax = emax;
1036:   pc_gamg->emin = emin;
1037:   return 0;
1038: }

1040: /*@
1041:    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding a `PCGAMG` algebraic multigrid preconditioner

1043:    Collective on pc

1045:    Input Parameters:
1046: +  pc - the preconditioner context
1047: -  n - `PETSC_TRUE` or `PETSC_FALSE`

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

1052:    Level: intermediate

1054:    Note:
1055:    May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1056:    rebuilding the preconditioner quicker.

1058: .seealso: `PCGAMG`
1059: @*/
1060: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1061: {
1063:   PetscTryMethod(pc, "PCGAMGSetReuseInterpolation_C", (PC, PetscBool), (pc, n));
1064:   return 0;
1065: }

1067: static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1068: {
1069:   PC_MG   *mg      = (PC_MG *)pc->data;
1070:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1072:   pc_gamg->reuse_prol = n;
1073:   return 0;
1074: }

1076: /*@
1077:    PCGAMGASMSetUseAggs - Have the `PCGAMG` smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner
1078:    used as the smoother

1080:    Collective on pc

1082:    Input Parameters:
1083: +  pc - the preconditioner context
1084: -  flg - `PETSC_TRUE` to use aggregates, `PETSC_FALSE` to not

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

1089:    Level: intermediate

1091: .seealso: `PCGAMG`, `PCASM`, `PCSetType`
1092: @*/
1093: PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1094: {
1096:   PetscTryMethod(pc, "PCGAMGASMSetUseAggs_C", (PC, PetscBool), (pc, flg));
1097:   return 0;
1098: }

1100: static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1101: {
1102:   PC_MG   *mg      = (PC_MG *)pc->data;
1103:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1105:   pc_gamg->use_aggs_in_asm = flg;
1106:   return 0;
1107: }

1109: /*@
1110:    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver

1112:    Collective on pc

1114:    Input Parameters:
1115: +  pc - the preconditioner context
1116: -  flg - `PETSC_TRUE` to not force coarse grid onto one processor

1118:    Options Database Key:
1119: .  -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver

1121:    Level: intermediate

1123: .seealso: `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()`
1124: @*/
1125: PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1126: {
1128:   PetscTryMethod(pc, "PCGAMGSetUseParallelCoarseGridSolve_C", (PC, PetscBool), (pc, flg));
1129:   return 0;
1130: }

1132: static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1133: {
1134:   PC_MG   *mg      = (PC_MG *)pc->data;
1135:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1137:   pc_gamg->use_parallel_coarse_grid_solver = flg;
1138:   return 0;
1139: }

1141: /*@
1142:    PCGAMGSetCpuPinCoarseGrids - pin the coarse grids created in `PCGAMG` to run only on the CPU since the problems may be too small to run efficiently on the GPUs

1144:    Collective on pc

1146:    Input Parameters:
1147: +  pc - the preconditioner context
1148: -  flg - `PETSC_TRUE` to pin coarse grids to the CPU

1150:    Options Database Key:
1151: .  -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU

1153:    Level: intermediate

1155: .seealso: `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetUseParallelCoarseGridSolve()`
1156: @*/
1157: PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1158: {
1160:   PetscTryMethod(pc, "PCGAMGSetCpuPinCoarseGrids_C", (PC, PetscBool), (pc, flg));
1161:   return 0;
1162: }

1164: static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1165: {
1166:   PC_MG   *mg      = (PC_MG *)pc->data;
1167:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1169:   pc_gamg->cpu_pin_coarse_grids = flg;
1170:   return 0;
1171: }

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

1176:    Collective on pc

1178:    Input Parameters:
1179: +  pc - the preconditioner context
1180: -  flg - `PCGAMGLayoutType` type, either `PCGAMG_LAYOUT_COMPACT` or `PCGAMG_LAYOUT_SPREAD`

1182:    Options Database Key:
1183: .  -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering

1185:    Level: intermediate

1187: .seealso: `PCGAMG`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGLayoutType`, `PCGAMG_LAYOUT_COMPACT`, `PCGAMG_LAYOUT_SPREAD`
1188: @*/
1189: PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1190: {
1192:   PetscTryMethod(pc, "PCGAMGSetCoarseGridLayoutType_C", (PC, PCGAMGLayoutType), (pc, flg));
1193:   return 0;
1194: }

1196: static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1197: {
1198:   PC_MG   *mg      = (PC_MG *)pc->data;
1199:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1201:   pc_gamg->layout_type = flg;
1202:   return 0;
1203: }

1205: /*@
1206:    PCGAMGSetNlevels -  Sets the maximum number of levels `PCGAMG` will use

1208:    Not collective on pc

1210:    Input Parameters:
1211: +  pc - the preconditioner
1212: -  n - the maximum number of levels to use

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

1217:    Level: intermediate

1219:    Developer Note:
1220:    Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG`

1222: .seealso: `PCGAMG`
1223: @*/
1224: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1225: {
1227:   PetscTryMethod(pc, "PCGAMGSetNlevels_C", (PC, PetscInt), (pc, n));
1228:   return 0;
1229: }

1231: static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1232: {
1233:   PC_MG   *mg      = (PC_MG *)pc->data;
1234:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1236:   pc_gamg->Nlevels = n;
1237:   return 0;
1238: }

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

1243:    Not collective on pc

1245:    Input Parameters:
1246: +  pc - the preconditioner context
1247: .  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
1248: -  n - number of threshold values provided in array

1250:    Options Database Key:
1251: .  -pc_gamg_threshold <threshold> - the threshold to drop edges

1253:    Notes:
1254:     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.
1255:     Before coarsening or aggregating the graph, `PCGAMG` 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.

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

1261:    Level: intermediate

1263: .seealso: `PCGAMG`, `PCGAMGFilterGraph()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGSetThresholdScale()`
1264: @*/
1265: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1266: {
1269:   PetscTryMethod(pc, "PCGAMGSetThreshold_C", (PC, PetscReal[], PetscInt), (pc, v, n));
1270:   return 0;
1271: }

1273: static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1274: {
1275:   PC_MG   *mg      = (PC_MG *)pc->data;
1276:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1277:   PetscInt i;
1278:   for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1279:   for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale;
1280:   return 0;
1281: }

1283: /*@
1284:    PCGAMGSetRankReductionFactors - Set a manual schedule for MPI rank reduction on coarse grids

1286:    Collective on pc

1288:    Input Parameters:
1289: +  pc - the preconditioner context
1290: .  v - array of reduction factors. 0 for first value forces a reduction to one process/device on first level in CUDA
1291: -  n - number of values provided in array

1293:    Options Database Key:
1294: .  -pc_gamg_rank_reduction_factors <factors> - provide the schedule

1296:    Level: intermediate

1298: .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`
1299: @*/
1300: PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1301: {
1304:   PetscTryMethod(pc, "PCGAMGSetRankReductionFactors_C", (PC, PetscInt[], PetscInt), (pc, v, n));
1305:   return 0;
1306: }

1308: static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1309: {
1310:   PC_MG   *mg      = (PC_MG *)pc->data;
1311:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1312:   PetscInt i;
1313:   for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1314:   for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
1315:   return 0;
1316: }

1318: /*@
1319:    PCGAMGSetThresholdScale - Relative threshold reduction at each level

1321:    Not collective on pc

1323:    Input Parameters:
1324: +  pc - the preconditioner context
1325: -  scale - the threshold value reduction, usually < 1.0

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

1330:    Note:
1331:    The initial threshold (for an arbitrary number of levels starting from the finest) can be set with `PCGAMGSetThreshold()`.
1332:    This scaling is used for each subsequent coarsening, but must be called before `PCGAMGSetThreshold()`.

1334:    Level: advanced

1336: .seealso: `PCGAMG`, `PCGAMGSetThreshold()`
1337: @*/
1338: PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1339: {
1341:   PetscTryMethod(pc, "PCGAMGSetThresholdScale_C", (PC, PetscReal), (pc, v));
1342:   return 0;
1343: }

1345: static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1346: {
1347:   PC_MG   *mg      = (PC_MG *)pc->data;
1348:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1349:   pc_gamg->threshold_scale = v;
1350:   return 0;
1351: }

1353: /*@C
1354:    PCGAMGSetType - Set the type of algorithm `PCGAMG` should use

1356:    Collective on pc

1358:    Input Parameters:
1359: +  pc - the preconditioner context
1360: -  type - `PCGAMGAGG`, `PCGAMGGEO`, or `PCGAMGCLASSICAL`

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

1365:    Level: intermediate

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

1376: /*@C
1377:    PCGAMGGetType - Get the type of algorithm `PCGAMG` will use

1379:    Collective on pc

1381:    Input Parameter:
1382: .  pc - the preconditioner context

1384:    Output Parameter:
1385: .  type - the type of algorithm used

1387:    Level: intermediate

1389: .seealso: `PCGAMG`, `PCGAMGSetType()`, `PCGAMGType`
1390: @*/
1391: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1392: {
1394:   PetscUseMethod(pc, "PCGAMGGetType_C", (PC, PCGAMGType *), (pc, type));
1395:   return 0;
1396: }

1398: static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1399: {
1400:   PC_MG   *mg      = (PC_MG *)pc->data;
1401:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1403:   *type = pc_gamg->type;
1404:   return 0;
1405: }

1407: static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1408: {
1409:   PC_MG   *mg      = (PC_MG *)pc->data;
1410:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1411:   PetscErrorCode (*r)(PC);

1413:   pc_gamg->type = type;
1414:   PetscFunctionListFind(GAMGList, type, &r);
1416:   if (pc_gamg->ops->destroy) {
1417:     (*pc_gamg->ops->destroy)(pc);
1418:     PetscMemzero(pc_gamg->ops, sizeof(struct _PCGAMGOps));
1419:     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1420:     /* cleaning up common data in pc_gamg - this should disapear someday */
1421:     pc_gamg->data_cell_cols      = 0;
1422:     pc_gamg->data_cell_rows      = 0;
1423:     pc_gamg->orig_data_cell_cols = 0;
1424:     pc_gamg->orig_data_cell_rows = 0;
1425:     PetscFree(pc_gamg->data);
1426:     pc_gamg->data_sz = 0;
1427:   }
1428:   PetscFree(pc_gamg->gamg_type_name);
1429:   PetscStrallocpy(type, &pc_gamg->gamg_type_name);
1430:   (*r)(pc);
1431:   return 0;
1432: }

1434: static PetscErrorCode PCView_GAMG(PC pc, PetscViewer viewer)
1435: {
1436:   PC_MG    *mg      = (PC_MG *)pc->data;
1437:   PC_GAMG  *pc_gamg = (PC_GAMG *)mg->innerctx;
1438:   PetscReal gc = 0, oc = 0;

1440:   PetscViewerASCIIPrintf(viewer, "    GAMG specific options\n");
1441:   PetscViewerASCIIPrintf(viewer, "      Threshold for dropping small values in graph on each level =");
1442:   for (PetscInt i = 0; i < mg->nlevels; i++) PetscViewerASCIIPrintf(viewer, " %g", (double)pc_gamg->threshold[i]);
1443:   PetscViewerASCIIPrintf(viewer, "\n");
1444:   PetscViewerASCIIPrintf(viewer, "      Threshold scaling factor for each level not specified = %g\n", (double)pc_gamg->threshold_scale);
1445:   if (pc_gamg->use_aggs_in_asm) PetscViewerASCIIPrintf(viewer, "      Using aggregates from coarsening process to define subdomains for PCASM\n");
1446:   if (pc_gamg->use_parallel_coarse_grid_solver) PetscViewerASCIIPrintf(viewer, "      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");
1447:   if (pc_gamg->ops->view) (*pc_gamg->ops->view)(pc, viewer);
1448:   PCMGGetGridComplexity(pc, &gc, &oc);
1449:   PetscViewerASCIIPrintf(viewer, "      Complexity:    grid = %g    operator = %g\n", (double)gc, (double)oc);
1450:   return 0;
1451: }

1453: PetscErrorCode PCSetFromOptions_GAMG(PC pc, PetscOptionItems *PetscOptionsObject)
1454: {
1455:   PC_MG             *mg      = (PC_MG *)pc->data;
1456:   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
1457:   PetscBool          flag;
1458:   MPI_Comm           comm;
1459:   char               prefix[256], tname[32];
1460:   PetscInt           i, n;
1461:   const char        *pcpre;
1462:   static const char *LayoutTypes[] = {"compact", "spread", "PCGAMGLayoutType", "PC_GAMG_LAYOUT", NULL};
1463:   PetscObjectGetComm((PetscObject)pc, &comm);
1464:   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG options");
1465:   PetscOptionsFList("-pc_gamg_type", "Type of AMG method", "PCGAMGSetType", GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);
1466:   if (flag) PCGAMGSetType(pc, tname);
1467:   PetscOptionsBool("-pc_gamg_repartition", "Repartion coarse grids", "PCGAMGSetRepartition", pc_gamg->repart, &pc_gamg->repart, NULL);
1468:   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);
1469:   PetscOptionsBool("-pc_gamg_reuse_interpolation", "Reuse prolongation operator", "PCGAMGReuseInterpolation", pc_gamg->reuse_prol, &pc_gamg->reuse_prol, NULL);
1470:   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);
1471:   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);
1472:   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);
1473:   PetscCall(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,
1474:                              (PetscEnum)pc_gamg->layout_type, (PetscEnum *)&pc_gamg->layout_type, NULL));
1475:   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);
1476:   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);
1477:   PetscOptionsReal("-pc_gamg_threshold_scale", "Scaling of threshold for each level not specified", "PCGAMGSetThresholdScale", pc_gamg->threshold_scale, &pc_gamg->threshold_scale, NULL);
1478:   n = PETSC_MG_MAXLEVELS;
1479:   PetscOptionsRealArray("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", pc_gamg->threshold, &n, &flag);
1480:   if (!flag || n < PETSC_MG_MAXLEVELS) {
1481:     if (!flag) n = 1;
1482:     i = n;
1483:     do {
1484:       pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale;
1485:     } while (++i < PETSC_MG_MAXLEVELS);
1486:   }
1487:   n = PETSC_MG_MAXLEVELS;
1488:   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);
1489:   if (!flag) i = 0;
1490:   else i = n;
1491:   do {
1492:     pc_gamg->level_reduction_factors[i] = -1;
1493:   } while (++i < PETSC_MG_MAXLEVELS);
1494:   PetscOptionsInt("-pc_mg_levels", "Set number of MG levels", "PCGAMGSetNlevels", pc_gamg->Nlevels, &pc_gamg->Nlevels, NULL);
1495:   {
1496:     PetscReal eminmax[2] = {0., 0.};
1497:     n                    = 2;
1498:     PetscOptionsRealArray("-pc_gamg_eigenvalues", "extreme eigenvalues for smoothed aggregation", "PCGAMGSetEigenvalues", eminmax, &n, &flag);
1499:     if (flag) {
1501:       PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);
1502:     }
1503:   }
1504:   /* set options for subtype */
1505:   (*pc_gamg->ops->setfromoptions)(pc, PetscOptionsObject);

1507:   PCGetOptionsPrefix(pc, &pcpre);
1508:   PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : "");
1509:   PetscOptionsHeadEnd();
1510:   return 0;
1511: }

1513: /*MC
1514:      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner

1516:    Options Database Keys:
1517: +   -pc_gamg_type <type,default=agg> - one of agg, geo, or classical
1518: .   -pc_gamg_repartition  <bool,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1519: .   -pc_gamg_asm_use_agg <bool,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother
1520: .   -pc_gamg_process_eq_limit <limit, default=50> - `PCGAMG` will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit>
1521:                                         equations on each process that has degrees of freedom
1522: .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
1523: .   -pc_gamg_reuse_interpolation <bool,default=true> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations (should always be true)
1524: .   -pc_gamg_threshold[] <thresh,default=[-1,...]> - Before aggregating the graph `PCGAMG` will remove small values from the graph on each level (< 0 does no filtering)
1525: -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified

1527:    Options Database Keys for Aggregation:
1528: +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1529: .  -pc_gamg_square_graph <n,default=1> - alias for pc_gamg_aggressive_coarsening (deprecated)
1530: -  -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest.

1532:    Options Database Keys for Multigrid:
1533: +  -pc_mg_cycles <v> - v or w, see `PCMGSetCycleType()`
1534: .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1535: .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1536: -  -pc_mg_levels <levels> - Number of levels of multigrid to use. GAMG has a heuristic so pc_mg_levels is not usually used with GAMG

1538:   Notes:
1539:   To obtain good performance for `PCGAMG` for vector valued problems you must
1540:   call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point
1541:   call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator

1543:   See [the Users Manual section on PCGAMG](sec_amg) for more details.

1545:   Level: intermediate

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

1551: PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1552: {
1553:   PC_GAMG *pc_gamg;
1554:   PC_MG   *mg;

1556:   /* register AMG type */
1557:   PCGAMGInitializePackage();

1559:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1560:   PCSetType(pc, PCMG);
1561:   PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);

1563:   /* create a supporting struct and attach it to pc */
1564:   PetscNew(&pc_gamg);
1565:   PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL);
1566:   mg           = (PC_MG *)pc->data;
1567:   mg->innerctx = pc_gamg;

1569:   PetscNew(&pc_gamg->ops);

1571:   /* these should be in subctx but repartitioning needs simple arrays */
1572:   pc_gamg->data_sz = 0;
1573:   pc_gamg->data    = NULL;

1575:   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1576:   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1577:   pc->ops->setup          = PCSetUp_GAMG;
1578:   pc->ops->reset          = PCReset_GAMG;
1579:   pc->ops->destroy        = PCDestroy_GAMG;
1580:   mg->view                = PCView_GAMG;

1582:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", PCGAMGSetProcEqLim_GAMG);
1583:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", PCGAMGSetCoarseEqLim_GAMG);
1584:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", PCGAMGSetRepartition_GAMG);
1585:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", PCGAMGSetEigenvalues_GAMG);
1586:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", PCGAMGSetUseSAEstEig_GAMG);
1587:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", PCGAMGSetReuseInterpolation_GAMG);
1588:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", PCGAMGASMSetUseAggs_GAMG);
1589:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", PCGAMGSetUseParallelCoarseGridSolve_GAMG);
1590:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", PCGAMGSetCpuPinCoarseGrids_GAMG);
1591:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", PCGAMGSetCoarseGridLayoutType_GAMG);
1592:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", PCGAMGSetThreshold_GAMG);
1593:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", PCGAMGSetRankReductionFactors_GAMG);
1594:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", PCGAMGSetThresholdScale_GAMG);
1595:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", PCGAMGSetType_GAMG);
1596:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", PCGAMGGetType_GAMG);
1597:   PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", PCGAMGSetNlevels_GAMG);
1598:   pc_gamg->repart                          = PETSC_FALSE;
1599:   pc_gamg->reuse_prol                      = PETSC_TRUE;
1600:   pc_gamg->use_aggs_in_asm                 = PETSC_FALSE;
1601:   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1602:   pc_gamg->cpu_pin_coarse_grids            = PETSC_FALSE;
1603:   pc_gamg->layout_type                     = PCGAMG_LAYOUT_SPREAD;
1604:   pc_gamg->min_eq_proc                     = 50;
1605:   pc_gamg->coarse_eq_limit                 = 50;
1606:   for (int i = 0; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = -1;
1607:   pc_gamg->threshold_scale = 1.;
1608:   pc_gamg->Nlevels         = PETSC_MG_MAXLEVELS;
1609:   pc_gamg->current_level   = 0; /* don't need to init really */
1610:   pc_gamg->use_sa_esteig   = PETSC_TRUE;
1611:   pc_gamg->emin            = 0;
1612:   pc_gamg->emax            = 0;

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

1616:   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1617:   PCGAMGSetType(pc, PCGAMGAGG);
1618:   return 0;
1619: }

1621: /*@C
1622:  PCGAMGInitializePackage - This function initializes everything in the `PCGAMG` package. It is called
1623:     from `PCInitializePackage()`.

1625:  Level: developer

1627:  .seealso: `PetscInitialize()`
1628: @*/
1629: PetscErrorCode PCGAMGInitializePackage(void)
1630: {
1631:   PetscInt l;

1633:   if (PCGAMGPackageInitialized) return 0;
1634:   PCGAMGPackageInitialized = PETSC_TRUE;
1635:   PetscFunctionListAdd(&GAMGList, PCGAMGGEO, PCCreateGAMG_GEO);
1636:   PetscFunctionListAdd(&GAMGList, PCGAMGAGG, PCCreateGAMG_AGG);
1637:   PetscFunctionListAdd(&GAMGList, PCGAMGCLASSICAL, PCCreateGAMG_Classical);
1638:   PetscRegisterFinalize(PCGAMGFinalizePackage);

1640:   /* general events */
1641:   PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP]);
1642:   PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH]);
1643:   PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN]);
1644:   PetscLogEventRegister("  GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS]);
1645:   PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL]);
1646:   PetscLogEventRegister("  GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA]);
1647:   PetscLogEventRegister("  GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB]);
1648:   PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT]);
1649:   PetscLogEventRegister("  GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM]);
1650:   PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL]);
1651:   PetscLogEventRegister("  GAMG PtAP", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP]);
1652:   PetscLogEventRegister("  GAMG Reduce", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE]);
1653:   PetscLogEventRegister("   GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART]);
1654:   /* PetscLogEventRegister("   GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1655:   /* PetscLogEventRegister("   GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]); */
1656:   /* PetscLogEventRegister("   GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]); */
1657:   for (l = 0; l < PETSC_MG_MAXLEVELS; l++) {
1658:     char ename[32];

1660:     PetscSNPrintf(ename, sizeof(ename), "PCGAMG Squ l%02" PetscInt_FMT, l);
1661:     PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]);
1662:     PetscSNPrintf(ename, sizeof(ename), "PCGAMG Gal l%02" PetscInt_FMT, l);
1663:     PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]);
1664:     PetscSNPrintf(ename, sizeof(ename), "PCGAMG Opt l%02" PetscInt_FMT, l);
1665:     PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]);
1666:   }
1667: #if defined(GAMG_STAGES)
1668:   { /* create timer stages */
1669:     char str[32];
1670:     sprintf(str, "GAMG Level %d", 0);
1671:     PetscLogStageRegister(str, &gamg_stages[0]);
1672:   }
1673: #endif
1674:   return 0;
1675: }

1677: /*@C
1678:  PCGAMGFinalizePackage - This function frees everything from the `PCGAMG` package. It is
1679:     called from `PetscFinalize()` automatically.

1681:  Level: developer

1683:  .seealso: `PetscFinalize()`
1684: @*/
1685: PetscErrorCode PCGAMGFinalizePackage(void)
1686: {
1687:   PCGAMGPackageInitialized = PETSC_FALSE;
1688:   PetscFunctionListDestroy(&GAMGList);
1689:   return 0;
1690: }

1692: /*@C
1693:  PCGAMGRegister - Register a `PCGAMG` implementation.

1695:  Input Parameters:
1696:  + type - string that will be used as the name of the `PCGAMG` type.
1697:  - create - function for creating the gamg context.

1699:   Level: developer

1701:  .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
1702: @*/
1703: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1704: {
1705:   PCGAMGInitializePackage();
1706:   PetscFunctionListAdd(&GAMGList, type, create);
1707:   return 0;
1708: }

1710: /*@C
1711:     PCGAMGCreateGraph - Creates a graph that is used by the ``PCGAMGType`` in the coarsening process

1713:    Input Parameters:
1714: +  pc - the `PCGAMG`
1715: -  A - the matrix, for any level

1717:    Output Parameter:
1718: .  G - the graph

1720:   Level: advanced

1722:  .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
1723: @*/
1724: PetscErrorCode PCGAMGCreateGraph(PC pc, Mat A, Mat *G)
1725: {
1726:   PC_MG   *mg      = (PC_MG *)pc->data;
1727:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1729:   pc_gamg->ops->creategraph(pc, A, G);
1730:   return 0;
1731: }