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:   PetscFunctionBegin;
 33:   PetscCall(PetscFree(pc_gamg->data));
 34:   pc_gamg->data_sz = 0;
 35:   PetscCall(PetscFree(pc_gamg->orig_data));
 36:   for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS; level++) {
 37:     mg->min_eigen_DinvA[level] = 0;
 38:     mg->max_eigen_DinvA[level] = 0;
 39:   }
 40:   pc_gamg->emin = 0;
 41:   pc_gamg->emax = 0;
 42:   PetscFunctionReturn(PETSC_SUCCESS);
 43: }

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

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

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

 80:   if (Pcolumnperm) *Pcolumnperm = NULL;

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

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

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

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

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

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

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

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

277:       PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc));
278:       PetscCall(PetscFree(newproc_idx));
279:       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0));
280:     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
281:       PetscInt targetPE;
282:       PetscCheck(new_size != nactive, PETSC_COMM_SELF, PETSC_ERR_PLIB, "new_size==nactive. Should not happen");
283:       PetscCall(PetscInfo(pc, "%s: Number of equations (loc) %" PetscInt_FMT " with simple aggregation\n", ((PetscObject)pc)->prefix, ncrs_eq));
284:       targetPE = (rank / rfactor) * expand_factor;
285:       PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc));
286:     } /* end simple 'is_eq_newproc' */

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

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

359:       pc_gamg->data_sz = node_data_sz * ncrs_new;
360:       strideNew        = ncrs_new * ndata_rows;

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

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

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

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

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

431:       /* PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0)); */
432:       PetscCall(MatDestroy(a_P_inout));

434:       /* output - repartitioned */
435:       *a_P_inout = Pnew;
436:     }
437:     PetscCall(ISDestroy(&new_eq_indices));

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

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

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

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

496: /*
497:    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
498:                     by setting data structures and options.

500:    Input Parameter:
501: .  pc - the preconditioner context

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

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

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

569:       PetscCall(PCSetUp_MG(pc));
570:       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0));
571:       PetscFunctionReturn(PETSC_SUCCESS);
572:     }
573:   }

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

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

584:       PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data));
585:       for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
586:     } else {
587:       PetscCheck(pc_gamg->ops->createdefaultdata, comm, PETSC_ERR_PLIB, "'createdefaultdata' not set(?) need to support NULL data");
588:       PetscCall(pc_gamg->ops->createdefaultdata(pc, Pmat));
589:     }
590:   }

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

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

604: #if defined(PETSC_USE_INFO)
605:   PetscCall(MatGetInfo(Pmat, MAT_GLOBAL_SUM, &info)); /* global reduction */
606:   nnz0   = info.nz_used;
607:   nnztot = info.nz_used;
608: #endif
609:   PetscCall(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));

611:   /* Get A_i and R_i */
612:   for (level = 0, Aarr[0] = Pmat, nactivepe = size; level < (pc_gamg->Nlevels - 1) && (!level || M > pc_gamg->coarse_eq_limit); level++) {
613:     pc_gamg->current_level = level;
614:     PetscCheck(level < PETSC_MG_MAXLEVELS, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Too many levels %" PetscInt_FMT, level);
615:     level1 = level + 1;
616: #if defined(GAMG_STAGES)
617:     if (!gamg_stages[level]) {
618:       char str[32];
619:       PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", (int)level));
620:       PetscCall(PetscLogStageRegister(str, &gamg_stages[level]));
621:     }
622:     PetscCall(PetscLogStagePush(gamg_stages[level]));
623: #endif
624:     { /* construct prolongator */
625:       Mat               Gmat;
626:       PetscCoarsenData *agg_lists;
627:       Mat               Prol11;

629:       PetscCall(PCGAMGCreateGraph(pc, Aarr[level], &Gmat));
630:       PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists));
631:       PetscCall(pc_gamg->ops->prolongator(pc, Aarr[level], Gmat, agg_lists, &Prol11));

633:       /* could have failed to create new level */
634:       if (Prol11) {
635:         const char *prefix;
636:         char        addp[32];

638:         /* get new block size of coarse matrices */
639:         PetscCall(MatGetBlockSizes(Prol11, NULL, &bs));

641:         if (pc_gamg->ops->optprolongator) {
642:           /* smooth */
643:           PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11));
644:         }

646:         if (pc_gamg->use_aggs_in_asm) {
647:           PetscInt bs;
648:           PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // not timed directly, ugly, could remove, but good ASM method
649:           PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]));
650:         }

652:         PetscCall(PCGetOptionsPrefix(pc, &prefix));
653:         PetscCall(MatSetOptionsPrefix(Prol11, prefix));
654:         PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_prolongator_%d_", (int)level));
655:         PetscCall(MatAppendOptionsPrefix(Prol11, addp));
656:         /* Always generate the transpose with CUDA
657:            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
658:         PetscCall(MatSetOption(Prol11, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE));
659:         PetscCall(MatSetFromOptions(Prol11));
660:         Parr[level1] = Prol11;
661:       } else Parr[level1] = NULL; /* failed to coarsen */

663:       PetscCall(MatDestroy(&Gmat));
664:       PetscCall(PetscCDDestroy(agg_lists));
665:     }                           /* construct prolongator scope */
666:     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
667:     if (!Parr[level1]) {        /* failed to coarsen */
668:       PetscCall(PetscInfo(pc, "%s: Stop gridding, level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level));
669: #if defined(GAMG_STAGES)
670:       PetscCall(PetscLogStagePop());
671: #endif
672:       break;
673:     }
674:     PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */
675:     PetscCheck(!is_last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Is last ?");
676:     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
677:     if (level1 == pc_gamg->Nlevels - 1) is_last = PETSC_TRUE;
678:     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0));
679:     PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last));
680:     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0));

682:     PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */
683: #if defined(PETSC_USE_INFO)
684:     PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info));
685:     nnztot += info.nz_used;
686: #endif
687:     PetscCall(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));

689: #if defined(GAMG_STAGES)
690:     PetscCall(PetscLogStagePop());
691: #endif
692:     /* stop if one node or one proc -- could pull back for singular problems */
693:     if ((pc_gamg->data_cell_cols && M / pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M / bs < 2)) {
694:       PetscCall(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));
695:       level++;
696:       break;
697:     }
698:   } /* levels */
699:   PetscCall(PetscFree(pc_gamg->data));

701:   PetscCall(PetscInfo(pc, "%s: %" PetscInt_FMT " levels, operator complexity = %g\n", ((PetscObject)pc)->prefix, level + 1, nnztot / nnz0));
702:   pc_gamg->Nlevels = level + 1;
703:   fine_level       = level;
704:   PetscCall(PCMGSetLevels(pc, pc_gamg->Nlevels, NULL));

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

708:     /* set default smoothers & set operators */
709:     for (lidx = 1, level = pc_gamg->Nlevels - 2; lidx <= fine_level; lidx++, level--) {
710:       KSP smoother;
711:       PC  subpc;

713:       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
714:       PetscCall(KSPGetPC(smoother, &subpc));

716:       PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
717:       /* set ops */
718:       PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level]));
719:       PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level + 1]));

721:       /* set defaults */
722:       PetscCall(KSPSetType(smoother, KSPCHEBYSHEV));

724:       /* set blocks for ASM smoother that uses the 'aggregates' */
725:       if (pc_gamg->use_aggs_in_asm) {
726:         PetscInt sz;
727:         IS      *iss;

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

759:       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
760:       PetscCall(KSPSetOperators(smoother, Lmat, Lmat));
761:       if (!pc_gamg->use_parallel_coarse_grid_solver) {
762:         PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
763:         PetscCall(KSPGetPC(smoother, &subpc));
764:         PetscCall(PCSetType(subpc, PCBJACOBI));
765:         PetscCall(PCSetUp(subpc));
766:         PetscCall(PCBJacobiGetSubKSP(subpc, &ii, &first, &k2));
767:         PetscCheck(ii == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " is not one", ii);
768:         PetscCall(KSPGetPC(k2[0], &pc2));
769:         PetscCall(PCSetType(pc2, PCLU));
770:         PetscCall(PCFactorSetShiftType(pc2, MAT_SHIFT_INBLOCKS));
771:         PetscCall(KSPSetTolerances(k2[0], PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1));
772:         PetscCall(KSPSetType(k2[0], KSPPREONLY));
773:       }
774:     }

776:     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
777:     PetscObjectOptionsBegin((PetscObject)pc);
778:     PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject));
779:     PetscOptionsEnd();
780:     PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));

782:     /* set cheby eigen estimates from SA to use in the solver */
783:     if (pc_gamg->use_sa_esteig) {
784:       for (lidx = 1, level = pc_gamg->Nlevels - 2; level >= 0; lidx++, level--) {
785:         KSP       smoother;
786:         PetscBool ischeb;

788:         PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
789:         PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb));
790:         if (ischeb) {
791:           KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data;

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

799:             emin = mg->min_eigen_DinvA[level];
800:             emax = mg->max_eigen_DinvA[level];
801:             PetscCall(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));
802:             cheb->emin_provided = emin;
803:             cheb->emax_provided = emax;
804:           }
805:         }
806:       }
807:     }

809:     PetscCall(PCSetUp_MG(pc));

811:     /* clean up */
812:     for (level = 1; level < pc_gamg->Nlevels; level++) {
813:       PetscCall(MatDestroy(&Parr[level]));
814:       PetscCall(MatDestroy(&Aarr[level]));
815:     }
816:   } else {
817:     KSP smoother;

819:     PetscCall(PetscInfo(pc, "%s: One level solver used (system is seen as DD). Using default solver.\n", ((PetscObject)pc)->prefix));
820:     PetscCall(PCMGGetSmoother(pc, 0, &smoother));
821:     PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0]));
822:     PetscCall(KSPSetType(smoother, KSPPREONLY));
823:     PetscCall(PCSetUp_MG(pc));
824:   }
825:   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0));
826:   PetscFunctionReturn(PETSC_SUCCESS);
827: }

829: /*
830:  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
831:    that was created with PCCreate_GAMG().

833:    Input Parameter:
834: .  pc - the preconditioner context

836:    Application Interface Routine: PCDestroy()
837: */
838: PetscErrorCode PCDestroy_GAMG(PC pc)
839: {
840:   PC_MG   *mg      = (PC_MG *)pc->data;
841:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

843:   PetscFunctionBegin;
844:   PetscCall(PCReset_GAMG(pc));
845:   if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc));
846:   PetscCall(PetscFree(pc_gamg->ops));
847:   PetscCall(PetscFree(pc_gamg->gamg_type_name));
848:   PetscCall(PetscFree(pc_gamg));
849:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", NULL));
850:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", NULL));
851:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", NULL));
852:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", NULL));
853:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", NULL));
854:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", NULL));
855:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", NULL));
856:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", NULL));
857:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", NULL));
858:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", NULL));
859:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", NULL));
860:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", NULL));
861:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", NULL));
862:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", NULL));
863:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", NULL));
864:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", NULL));
865:   PetscCall(PCDestroy_MG(pc));
866:   PetscFunctionReturn(PETSC_SUCCESS);
867: }

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

872:    Logically Collective

874:    Input Parameters:
875: +  pc - the preconditioner context
876: -  n - the number of equations

878:    Options Database Key:
879: .  -pc_gamg_process_eq_limit <limit> - set the limit

881:    Level: intermediate

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

887: .seealso: `PCGAMG`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`
888: @*/
889: PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n)
890: {
891:   PetscFunctionBegin;
893:   PetscTryMethod(pc, "PCGAMGSetProcEqLim_C", (PC, PetscInt), (pc, n));
894:   PetscFunctionReturn(PETSC_SUCCESS);
895: }

897: static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
898: {
899:   PC_MG   *mg      = (PC_MG *)pc->data;
900:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

902:   PetscFunctionBegin;
903:   if (n > 0) pc_gamg->min_eq_proc = n;
904:   PetscFunctionReturn(PETSC_SUCCESS);
905: }

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

910:    Collective

912:    Input Parameters:
913: +  pc - the preconditioner context
914: -  n - maximum number of equations to aim for

916:    Options Database Key:
917: .  -pc_gamg_coarse_eq_limit <limit> - set the limit

919:    Level: intermediate

921:    Note:
922:    For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
923:    has less than 1000 unknowns.

925: .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`
926: @*/
927: PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
928: {
929:   PetscFunctionBegin;
931:   PetscTryMethod(pc, "PCGAMGSetCoarseEqLim_C", (PC, PetscInt), (pc, n));
932:   PetscFunctionReturn(PETSC_SUCCESS);
933: }

935: static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
936: {
937:   PC_MG   *mg      = (PC_MG *)pc->data;
938:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

940:   PetscFunctionBegin;
941:   if (n > 0) pc_gamg->coarse_eq_limit = n;
942:   PetscFunctionReturn(PETSC_SUCCESS);
943: }

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

948:    Collective

950:    Input Parameters:
951: +  pc - the preconditioner context
952: -  n - `PETSC_TRUE` or `PETSC_FALSE`

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

957:    Level: intermediate

959:    Note:
960:    This will generally improve the loading balancing of the work on each level

962: .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`
963: @*/
964: PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
965: {
966:   PetscFunctionBegin;
968:   PetscTryMethod(pc, "PCGAMGSetRepartition_C", (PC, PetscBool), (pc, n));
969:   PetscFunctionReturn(PETSC_SUCCESS);
970: }

972: static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
973: {
974:   PC_MG   *mg      = (PC_MG *)pc->data;
975:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

977:   PetscFunctionBegin;
978:   pc_gamg->repart = n;
979:   PetscFunctionReturn(PETSC_SUCCESS);
980: }

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

985:    Collective

987:    Input Parameters:
988: +  pc - the preconditioner context
989: -  n - number of its

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

994:    Level: advanced

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

1003: .seealso: `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`
1004: @*/
1005: PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
1006: {
1007:   PetscFunctionBegin;
1009:   PetscTryMethod(pc, "PCGAMGSetUseSAEstEig_C", (PC, PetscBool), (pc, n));
1010:   PetscFunctionReturn(PETSC_SUCCESS);
1011: }

1013: static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
1014: {
1015:   PC_MG   *mg      = (PC_MG *)pc->data;
1016:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1018:   PetscFunctionBegin;
1019:   pc_gamg->use_sa_esteig = n;
1020:   PetscFunctionReturn(PETSC_SUCCESS);
1021: }

1023: /*@
1024:    PCGAMGSetEigenvalues - Set WHAT eigenvalues WHY?

1026:    Collective

1028:    Input Parameters:
1029: +  pc - the preconditioner context
1030: -  emax - max eigenvalue
1031: .  emin - min eigenvalue

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

1036:    Level: intermediate

1038: .seealso: `PCGAMG`, `PCGAMGSetUseSAEstEig()`
1039: @*/
1040: PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax, PetscReal emin)
1041: {
1042:   PetscFunctionBegin;
1044:   PetscTryMethod(pc, "PCGAMGSetEigenvalues_C", (PC, PetscReal, PetscReal), (pc, emax, emin));
1045:   PetscFunctionReturn(PETSC_SUCCESS);
1046: }

1048: static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc, PetscReal emax, PetscReal emin)
1049: {
1050:   PC_MG   *mg      = (PC_MG *)pc->data;
1051:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1053:   PetscFunctionBegin;
1054:   PetscCheck(emax > emin, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Maximum eigenvalue must be larger than minimum: max %g min %g", (double)emax, (double)emin);
1055:   PetscCheck(emax * emin > 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Both eigenvalues must be of the same sign: max %g min %g", (double)emax, (double)emin);
1056:   pc_gamg->emax = emax;
1057:   pc_gamg->emin = emin;
1058:   PetscFunctionReturn(PETSC_SUCCESS);
1059: }

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

1064:    Collective

1066:    Input Parameters:
1067: +  pc - the preconditioner context
1068: -  n - `PETSC_TRUE` or `PETSC_FALSE`

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

1073:    Level: intermediate

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

1079: .seealso: `PCGAMG`
1080: @*/
1081: PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1082: {
1083:   PetscFunctionBegin;
1085:   PetscTryMethod(pc, "PCGAMGSetReuseInterpolation_C", (PC, PetscBool), (pc, n));
1086:   PetscFunctionReturn(PETSC_SUCCESS);
1087: }

1089: static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1090: {
1091:   PC_MG   *mg      = (PC_MG *)pc->data;
1092:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1094:   PetscFunctionBegin;
1095:   pc_gamg->reuse_prol = n;
1096:   PetscFunctionReturn(PETSC_SUCCESS);
1097: }

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

1103:    Collective

1105:    Input Parameters:
1106: +  pc - the preconditioner context
1107: -  flg - `PETSC_TRUE` to use aggregates, `PETSC_FALSE` to not

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

1112:    Level: intermediate

1114: .seealso: `PCGAMG`, `PCASM`, `PCSetType`
1115: @*/
1116: PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1117: {
1118:   PetscFunctionBegin;
1120:   PetscTryMethod(pc, "PCGAMGASMSetUseAggs_C", (PC, PetscBool), (pc, flg));
1121:   PetscFunctionReturn(PETSC_SUCCESS);
1122: }

1124: static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1125: {
1126:   PC_MG   *mg      = (PC_MG *)pc->data;
1127:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1129:   PetscFunctionBegin;
1130:   pc_gamg->use_aggs_in_asm = flg;
1131:   PetscFunctionReturn(PETSC_SUCCESS);
1132: }

1134: /*@
1135:    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver

1137:    Collective

1139:    Input Parameters:
1140: +  pc - the preconditioner context
1141: -  flg - `PETSC_TRUE` to not force coarse grid onto one processor

1143:    Options Database Key:
1144: .  -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver

1146:    Level: intermediate

1148: .seealso: `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()`
1149: @*/
1150: PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1151: {
1152:   PetscFunctionBegin;
1154:   PetscTryMethod(pc, "PCGAMGSetUseParallelCoarseGridSolve_C", (PC, PetscBool), (pc, flg));
1155:   PetscFunctionReturn(PETSC_SUCCESS);
1156: }

1158: static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1159: {
1160:   PC_MG   *mg      = (PC_MG *)pc->data;
1161:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1163:   PetscFunctionBegin;
1164:   pc_gamg->use_parallel_coarse_grid_solver = flg;
1165:   PetscFunctionReturn(PETSC_SUCCESS);
1166: }

1168: /*@
1169:    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

1171:    Collective

1173:    Input Parameters:
1174: +  pc - the preconditioner context
1175: -  flg - `PETSC_TRUE` to pin coarse grids to the CPU

1177:    Options Database Key:
1178: .  -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU

1180:    Level: intermediate

1182: .seealso: `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetUseParallelCoarseGridSolve()`
1183: @*/
1184: PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1185: {
1186:   PetscFunctionBegin;
1188:   PetscTryMethod(pc, "PCGAMGSetCpuPinCoarseGrids_C", (PC, PetscBool), (pc, flg));
1189:   PetscFunctionReturn(PETSC_SUCCESS);
1190: }

1192: static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1193: {
1194:   PC_MG   *mg      = (PC_MG *)pc->data;
1195:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1197:   PetscFunctionBegin;
1198:   pc_gamg->cpu_pin_coarse_grids = flg;
1199:   PetscFunctionReturn(PETSC_SUCCESS);
1200: }

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

1205:    Collective

1207:    Input Parameters:
1208: +  pc - the preconditioner context
1209: -  flg - `PCGAMGLayoutType` type, either `PCGAMG_LAYOUT_COMPACT` or `PCGAMG_LAYOUT_SPREAD`

1211:    Options Database Key:
1212: .  -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering

1214:    Level: intermediate

1216: .seealso: `PCGAMG`, `PCGAMGSetUseParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGLayoutType`, `PCGAMG_LAYOUT_COMPACT`, `PCGAMG_LAYOUT_SPREAD`
1217: @*/
1218: PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1219: {
1220:   PetscFunctionBegin;
1222:   PetscTryMethod(pc, "PCGAMGSetCoarseGridLayoutType_C", (PC, PCGAMGLayoutType), (pc, flg));
1223:   PetscFunctionReturn(PETSC_SUCCESS);
1224: }

1226: static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1227: {
1228:   PC_MG   *mg      = (PC_MG *)pc->data;
1229:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1231:   PetscFunctionBegin;
1232:   pc_gamg->layout_type = flg;
1233:   PetscFunctionReturn(PETSC_SUCCESS);
1234: }

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

1239:    Not Collective

1241:    Input Parameters:
1242: +  pc - the preconditioner
1243: -  n - the maximum number of levels to use

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

1248:    Level: intermediate

1250:    Developer Note:
1251:    Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG`

1253: .seealso: `PCGAMG`
1254: @*/
1255: PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1256: {
1257:   PetscFunctionBegin;
1259:   PetscTryMethod(pc, "PCGAMGSetNlevels_C", (PC, PetscInt), (pc, n));
1260:   PetscFunctionReturn(PETSC_SUCCESS);
1261: }

1263: static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1264: {
1265:   PC_MG   *mg      = (PC_MG *)pc->data;
1266:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1268:   PetscFunctionBegin;
1269:   pc_gamg->Nlevels = n;
1270:   PetscFunctionReturn(PETSC_SUCCESS);
1271: }

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

1276:    Not Collective

1278:    Input Parameters:
1279: +  pc - the preconditioner context
1280: .  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
1281: -  n - number of threshold values provided in array

1283:    Options Database Key:
1284: .  -pc_gamg_threshold <threshold> - the threshold to drop edges

1286:    Level: intermediate

1288:    Notes:
1289:     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.
1290:     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.

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

1296: .seealso: `PCGAMG`, `PCGAMGSetAggressiveLevels()`, `PCGAMGSetThresholdScale()`
1297: @*/
1298: PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1299: {
1300:   PetscFunctionBegin;
1303:   PetscTryMethod(pc, "PCGAMGSetThreshold_C", (PC, PetscReal[], PetscInt), (pc, v, n));
1304:   PetscFunctionReturn(PETSC_SUCCESS);
1305: }

1307: static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1308: {
1309:   PC_MG   *mg      = (PC_MG *)pc->data;
1310:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1311:   PetscInt i;
1312:   PetscFunctionBegin;
1313:   for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1314:   for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale;
1315:   PetscFunctionReturn(PETSC_SUCCESS);
1316: }

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

1321:    Collective

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

1328:    Options Database Key:
1329: .  -pc_gamg_rank_reduction_factors <factors> - provide the schedule

1331:    Level: intermediate

1333: .seealso: `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`
1334: @*/
1335: PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1336: {
1337:   PetscFunctionBegin;
1340:   PetscTryMethod(pc, "PCGAMGSetRankReductionFactors_C", (PC, PetscInt[], PetscInt), (pc, v, n));
1341:   PetscFunctionReturn(PETSC_SUCCESS);
1342: }

1344: static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1345: {
1346:   PC_MG   *mg      = (PC_MG *)pc->data;
1347:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1348:   PetscInt i;
1349:   PetscFunctionBegin;
1350:   for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1351:   for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
1352:   PetscFunctionReturn(PETSC_SUCCESS);
1353: }

1355: /*@
1356:    PCGAMGSetThresholdScale - Relative threshold reduction at each level

1358:    Not Collective

1360:    Input Parameters:
1361: +  pc - the preconditioner context
1362: -  scale - the threshold value reduction, usually < 1.0

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

1367:    Level: advanced

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

1373: .seealso: `PCGAMG`, `PCGAMGSetThreshold()`
1374: @*/
1375: PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1376: {
1377:   PetscFunctionBegin;
1379:   PetscTryMethod(pc, "PCGAMGSetThresholdScale_C", (PC, PetscReal), (pc, v));
1380:   PetscFunctionReturn(PETSC_SUCCESS);
1381: }

1383: static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1384: {
1385:   PC_MG   *mg      = (PC_MG *)pc->data;
1386:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1387:   PetscFunctionBegin;
1388:   pc_gamg->threshold_scale = v;
1389:   PetscFunctionReturn(PETSC_SUCCESS);
1390: }

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

1395:    Collective

1397:    Input Parameters:
1398: +  pc - the preconditioner context
1399: -  type - `PCGAMGAGG`, `PCGAMGGEO`, or `PCGAMGCLASSICAL`

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

1404:    Level: intermediate

1406: .seealso: `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType`
1407: @*/
1408: PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1409: {
1410:   PetscFunctionBegin;
1412:   PetscTryMethod(pc, "PCGAMGSetType_C", (PC, PCGAMGType), (pc, type));
1413:   PetscFunctionReturn(PETSC_SUCCESS);
1414: }

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

1419:    Collective

1421:    Input Parameter:
1422: .  pc - the preconditioner context

1424:    Output Parameter:
1425: .  type - the type of algorithm used

1427:    Level: intermediate

1429: .seealso: `PCGAMG`, `PCGAMGSetType()`, `PCGAMGType`
1430: @*/
1431: PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1432: {
1433:   PetscFunctionBegin;
1435:   PetscUseMethod(pc, "PCGAMGGetType_C", (PC, PCGAMGType *), (pc, type));
1436:   PetscFunctionReturn(PETSC_SUCCESS);
1437: }

1439: static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1440: {
1441:   PC_MG   *mg      = (PC_MG *)pc->data;
1442:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1444:   PetscFunctionBegin;
1445:   *type = pc_gamg->type;
1446:   PetscFunctionReturn(PETSC_SUCCESS);
1447: }

1449: static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1450: {
1451:   PC_MG   *mg      = (PC_MG *)pc->data;
1452:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1453:   PetscErrorCode (*r)(PC);

1455:   PetscFunctionBegin;
1456:   pc_gamg->type = type;
1457:   PetscCall(PetscFunctionListFind(GAMGList, type, &r));
1458:   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown GAMG type %s given", type);
1459:   if (pc_gamg->ops->destroy) {
1460:     PetscCall((*pc_gamg->ops->destroy)(pc));
1461:     PetscCall(PetscMemzero(pc_gamg->ops, sizeof(struct _PCGAMGOps)));
1462:     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1463:     /* cleaning up common data in pc_gamg - this should disappear someday */
1464:     pc_gamg->data_cell_cols      = 0;
1465:     pc_gamg->data_cell_rows      = 0;
1466:     pc_gamg->orig_data_cell_cols = 0;
1467:     pc_gamg->orig_data_cell_rows = 0;
1468:     PetscCall(PetscFree(pc_gamg->data));
1469:     pc_gamg->data_sz = 0;
1470:   }
1471:   PetscCall(PetscFree(pc_gamg->gamg_type_name));
1472:   PetscCall(PetscStrallocpy(type, &pc_gamg->gamg_type_name));
1473:   PetscCall((*r)(pc));
1474:   PetscFunctionReturn(PETSC_SUCCESS);
1475: }

1477: static PetscErrorCode PCView_GAMG(PC pc, PetscViewer viewer)
1478: {
1479:   PC_MG    *mg      = (PC_MG *)pc->data;
1480:   PC_GAMG  *pc_gamg = (PC_GAMG *)mg->innerctx;
1481:   PetscReal gc = 0, oc = 0;

1483:   PetscFunctionBegin;
1484:   PetscCall(PetscViewerASCIIPrintf(viewer, "    GAMG specific options\n"));
1485:   PetscCall(PetscViewerASCIIPrintf(viewer, "      Threshold for dropping small values in graph on each level ="));
1486:   for (PetscInt i = 0; i < mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)pc_gamg->threshold[i]));
1487:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1488:   PetscCall(PetscViewerASCIIPrintf(viewer, "      Threshold scaling factor for each level not specified = %g\n", (double)pc_gamg->threshold_scale));
1489:   if (pc_gamg->use_aggs_in_asm) PetscCall(PetscViewerASCIIPrintf(viewer, "      Using aggregates from coarsening process to define subdomains for PCASM\n"));
1490:   if (pc_gamg->use_parallel_coarse_grid_solver) PetscCall(PetscViewerASCIIPrintf(viewer, "      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n"));
1491:   if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc, viewer));
1492:   PetscCall(PCMGGetGridComplexity(pc, &gc, &oc));
1493:   PetscCall(PetscViewerASCIIPrintf(viewer, "      Complexity:    grid = %g    operator = %g\n", (double)gc, (double)oc));
1494:   PetscFunctionReturn(PETSC_SUCCESS);
1495: }

1497: PetscErrorCode PCSetFromOptions_GAMG(PC pc, PetscOptionItems *PetscOptionsObject)
1498: {
1499:   PC_MG             *mg      = (PC_MG *)pc->data;
1500:   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
1501:   PetscBool          flag;
1502:   MPI_Comm           comm;
1503:   char               prefix[256], tname[32];
1504:   PetscInt           i, n;
1505:   const char        *pcpre;
1506:   static const char *LayoutTypes[] = {"compact", "spread", "PCGAMGLayoutType", "PC_GAMG_LAYOUT", NULL};
1507:   PetscFunctionBegin;
1508:   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1509:   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG options");
1510:   PetscCall(PetscOptionsFList("-pc_gamg_type", "Type of AMG method", "PCGAMGSetType", GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag));
1511:   if (flag) PetscCall(PCGAMGSetType(pc, tname));
1512:   PetscCall(PetscOptionsBool("-pc_gamg_repartition", "Repartion coarse grids", "PCGAMGSetRepartition", pc_gamg->repart, &pc_gamg->repart, NULL));
1513:   PetscCall(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));
1514:   PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation", "Reuse prolongation operator", "PCGAMGReuseInterpolation", pc_gamg->reuse_prol, &pc_gamg->reuse_prol, NULL));
1515:   PetscCall(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));
1516:   PetscCall(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));
1517:   PetscCall(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));
1518:   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,
1519:                              (PetscEnum)pc_gamg->layout_type, (PetscEnum *)&pc_gamg->layout_type, NULL));
1520:   PetscCall(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));
1521:   PetscCall(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));
1522:   PetscCall(PetscOptionsReal("-pc_gamg_threshold_scale", "Scaling of threshold for each level not specified", "PCGAMGSetThresholdScale", pc_gamg->threshold_scale, &pc_gamg->threshold_scale, NULL));
1523:   n = PETSC_MG_MAXLEVELS;
1524:   PetscCall(PetscOptionsRealArray("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", pc_gamg->threshold, &n, &flag));
1525:   if (!flag || n < PETSC_MG_MAXLEVELS) {
1526:     if (!flag) n = 1;
1527:     i = n;
1528:     do {
1529:       pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale;
1530:     } while (++i < PETSC_MG_MAXLEVELS);
1531:   }
1532:   n = PETSC_MG_MAXLEVELS;
1533:   PetscCall(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));
1534:   if (!flag) i = 0;
1535:   else i = n;
1536:   do {
1537:     pc_gamg->level_reduction_factors[i] = -1;
1538:   } while (++i < PETSC_MG_MAXLEVELS);
1539:   PetscCall(PetscOptionsInt("-pc_mg_levels", "Set number of MG levels", "PCGAMGSetNlevels", pc_gamg->Nlevels, &pc_gamg->Nlevels, NULL));
1540:   {
1541:     PetscReal eminmax[2] = {0., 0.};
1542:     n                    = 2;
1543:     PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues", "extreme eigenvalues for smoothed aggregation", "PCGAMGSetEigenvalues", eminmax, &n, &flag));
1544:     if (flag) {
1545:       PetscCheck(n == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
1546:       PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]));
1547:     }
1548:   }
1549:   /* set options for subtype */
1550:   PetscCall((*pc_gamg->ops->setfromoptions)(pc, PetscOptionsObject));

1552:   PetscCall(PCGetOptionsPrefix(pc, &pcpre));
1553:   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : ""));
1554:   PetscOptionsHeadEnd();
1555:   PetscFunctionReturn(PETSC_SUCCESS);
1556: }

1558: /*MC
1559:      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner

1561:    Options Database Keys:
1562: +   -pc_gamg_type <type,default=agg> - one of agg, geo, or classical
1563: .   -pc_gamg_repartition  <bool,default=false> - repartition the degrees of freedom across the coarse grids as they are determined
1564: .   -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
1565: .   -pc_gamg_process_eq_limit <limit, default=50> - `PCGAMG` will reduce the number of MPI ranks used directly on the coarse grids so that there are around <limit>
1566:                                         equations on each process that has degrees of freedom
1567: .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
1568: .   -pc_gamg_reuse_interpolation <bool,default=true> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations (should always be true)
1569: .   -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)
1570: -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified

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

1577:    Options Database Keys for Multigrid:
1578: +  -pc_mg_cycle_type <v> - v or w, see `PCMGSetCycleType()`
1579: .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1580: .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1581: -  -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

1583:   Level: intermediate

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

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

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

1596: PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1597: {
1598:   PC_GAMG *pc_gamg;
1599:   PC_MG   *mg;

1601:   PetscFunctionBegin;
1602:   /* register AMG type */
1603:   PetscCall(PCGAMGInitializePackage());

1605:   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1606:   PetscCall(PCSetType(pc, PCMG));
1607:   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG));

1609:   /* create a supporting struct and attach it to pc */
1610:   PetscCall(PetscNew(&pc_gamg));
1611:   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
1612:   mg           = (PC_MG *)pc->data;
1613:   mg->innerctx = pc_gamg;

1615:   PetscCall(PetscNew(&pc_gamg->ops));

1617:   /* these should be in subctx but repartitioning needs simple arrays */
1618:   pc_gamg->data_sz = 0;
1619:   pc_gamg->data    = NULL;

1621:   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1622:   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1623:   pc->ops->setup          = PCSetUp_GAMG;
1624:   pc->ops->reset          = PCReset_GAMG;
1625:   pc->ops->destroy        = PCDestroy_GAMG;
1626:   mg->view                = PCView_GAMG;

1628:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", PCGAMGSetProcEqLim_GAMG));
1629:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", PCGAMGSetCoarseEqLim_GAMG));
1630:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", PCGAMGSetRepartition_GAMG));
1631:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", PCGAMGSetEigenvalues_GAMG));
1632:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", PCGAMGSetUseSAEstEig_GAMG));
1633:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", PCGAMGSetReuseInterpolation_GAMG));
1634:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", PCGAMGASMSetUseAggs_GAMG));
1635:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseParallelCoarseGridSolve_C", PCGAMGSetUseParallelCoarseGridSolve_GAMG));
1636:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", PCGAMGSetCpuPinCoarseGrids_GAMG));
1637:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", PCGAMGSetCoarseGridLayoutType_GAMG));
1638:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", PCGAMGSetThreshold_GAMG));
1639:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", PCGAMGSetRankReductionFactors_GAMG));
1640:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", PCGAMGSetThresholdScale_GAMG));
1641:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", PCGAMGSetType_GAMG));
1642:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", PCGAMGGetType_GAMG));
1643:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", PCGAMGSetNlevels_GAMG));
1644:   pc_gamg->repart                          = PETSC_FALSE;
1645:   pc_gamg->reuse_prol                      = PETSC_TRUE;
1646:   pc_gamg->use_aggs_in_asm                 = PETSC_FALSE;
1647:   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1648:   pc_gamg->cpu_pin_coarse_grids            = PETSC_FALSE;
1649:   pc_gamg->layout_type                     = PCGAMG_LAYOUT_SPREAD;
1650:   pc_gamg->min_eq_proc                     = 50;
1651:   pc_gamg->coarse_eq_limit                 = 50;
1652:   for (int i = 0; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = -1;
1653:   pc_gamg->threshold_scale = 1.;
1654:   pc_gamg->Nlevels         = PETSC_MG_MAXLEVELS;
1655:   pc_gamg->current_level   = 0; /* don't need to init really */
1656:   pc_gamg->use_sa_esteig   = PETSC_TRUE;
1657:   pc_gamg->emin            = 0;
1658:   pc_gamg->emax            = 0;

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

1662:   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1663:   PetscCall(PCGAMGSetType(pc, PCGAMGAGG));
1664:   PetscFunctionReturn(PETSC_SUCCESS);
1665: }

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

1671:  Level: developer

1673:  .seealso: `PetscInitialize()`
1674: @*/
1675: PetscErrorCode PCGAMGInitializePackage(void)
1676: {
1677:   PetscInt l;

1679:   PetscFunctionBegin;
1680:   if (PCGAMGPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
1681:   PCGAMGPackageInitialized = PETSC_TRUE;
1682:   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGGEO, PCCreateGAMG_GEO));
1683:   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGAGG, PCCreateGAMG_AGG));
1684:   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGCLASSICAL, PCCreateGAMG_Classical));
1685:   PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage));

1687:   /* general events */
1688:   PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP]));
1689:   PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH]));
1690:   PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN]));
1691:   PetscCall(PetscLogEventRegister("  GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS]));
1692:   PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL]));
1693:   PetscCall(PetscLogEventRegister("  GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA]));
1694:   PetscCall(PetscLogEventRegister("  GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB]));
1695:   PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT]));
1696:   PetscCall(PetscLogEventRegister("  GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM]));
1697:   PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL]));
1698:   PetscCall(PetscLogEventRegister("  GAMG PtAP", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP]));
1699:   PetscCall(PetscLogEventRegister("  GAMG Reduce", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE]));
1700:   PetscCall(PetscLogEventRegister("   GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART]));
1701:   /* PetscCall(PetscLogEventRegister("   GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */
1702:   /* PetscCall(PetscLogEventRegister("   GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */
1703:   /* PetscCall(PetscLogEventRegister("   GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */
1704:   for (l = 0; l < PETSC_MG_MAXLEVELS; l++) {
1705:     char ename[32];

1707:     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Squ l%02" PetscInt_FMT, l));
1708:     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]));
1709:     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Gal l%02" PetscInt_FMT, l));
1710:     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]));
1711:     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Opt l%02" PetscInt_FMT, l));
1712:     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]));
1713:   }
1714: #if defined(GAMG_STAGES)
1715:   { /* create timer stages */
1716:     char str[32];
1717:     PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", 0));
1718:     PetscCall(PetscLogStageRegister(str, &gamg_stages[0]));
1719:   }
1720: #endif
1721:   PetscFunctionReturn(PETSC_SUCCESS);
1722: }

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

1728:  Level: developer

1730:  .seealso: `PetscFinalize()`
1731: @*/
1732: PetscErrorCode PCGAMGFinalizePackage(void)
1733: {
1734:   PetscFunctionBegin;
1735:   PCGAMGPackageInitialized = PETSC_FALSE;
1736:   PetscCall(PetscFunctionListDestroy(&GAMGList));
1737:   PetscFunctionReturn(PETSC_SUCCESS);
1738: }

1740: /*@C
1741:  PCGAMGRegister - Register a `PCGAMG` implementation.

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

1747:   Level: developer

1749:  .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
1750: @*/
1751: PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1752: {
1753:   PetscFunctionBegin;
1754:   PetscCall(PCGAMGInitializePackage());
1755:   PetscCall(PetscFunctionListAdd(&GAMGList, type, create));
1756:   PetscFunctionReturn(PETSC_SUCCESS);
1757: }

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

1762:    Input Parameters:
1763: +  pc - the `PCGAMG`
1764: -  A - the matrix, for any level

1766:    Output Parameter:
1767: .  G - the graph

1769:   Level: advanced

1771:  .seealso: `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
1772: @*/
1773: PetscErrorCode PCGAMGCreateGraph(PC pc, Mat A, Mat *G)
1774: {
1775:   PC_MG   *mg      = (PC_MG *)pc->data;
1776:   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;

1778:   PetscFunctionBegin;
1779:   PetscCall(pc_gamg->ops->creategraph(pc, A, G));
1780:   PetscFunctionReturn(PETSC_SUCCESS);
1781: }