xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision e0b7e82fd3cf27fce84cc3e37e8d70a5c36a2d4e)
1 /*
2  GAMG geometric-algebric multigrid PC - Mark Adams 2011
3  */
4 #include <../src/ksp/pc/impls/gamg/gamg.h>            /*I "petscpc.h" I*/
5 #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h> /*I "petscksp.h" I*/
6 
7 #if defined(PETSC_HAVE_CUDA)
8   #include <cuda_runtime.h>
9 #endif
10 
11 #if defined(PETSC_HAVE_HIP)
12   #include <hip/hip_runtime.h>
13 #endif
14 
15 PetscLogEvent petsc_gamg_setup_events[GAMG_NUM_SET];
16 PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3];
17 
18 // #define GAMG_STAGES
19 #if defined(GAMG_STAGES)
20 static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
21 #endif
22 
23 static PetscFunctionList GAMGList = NULL;
24 static PetscBool         PCGAMGPackageInitialized;
25 
26 static PetscErrorCode PCReset_GAMG(PC pc)
27 {
28   PC_MG   *mg      = (PC_MG *)pc->data;
29   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
30 
31   PetscFunctionBegin;
32   PetscCall(PetscFree(pc_gamg->data));
33   pc_gamg->data_sz = 0;
34   PetscCall(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   PetscCall(PCReset_MG(pc));
42   PetscCall(MatCoarsenDestroy(&pc_gamg->asm_crs));
43   PetscFunctionReturn(PETSC_SUCCESS);
44 }
45 
46 /*
47    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
48      of active processors.
49 
50    Input Parameter:
51    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
52           'pc_gamg->data_sz' are changed via repartitioning/reduction.
53    . Amat_fine - matrix on this fine (k) level
54    . cr_bs - coarse block size
55    In/Output Parameter:
56    . a_P_inout - prolongation operator to the next level (k-->k-1)
57    . a_nactive_proc - number of active procs
58    Output Parameter:
59    . a_Amat_crs - coarse matrix that is created (k-1)
60 */
61 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)
62 {
63   PC_MG      *mg      = (PC_MG *)pc->data;
64   PC_GAMG    *pc_gamg = (PC_GAMG *)mg->innerctx;
65   Mat         Cmat, Pold = *a_P_inout;
66   MPI_Comm    comm;
67   PetscMPIInt rank, size, new_size, nactive = *a_nactive_proc;
68   PetscInt    ncrs_eq, ncrs, f_bs;
69 
70   PetscFunctionBegin;
71   PetscCall(PetscObjectGetComm((PetscObject)Amat_fine, &comm));
72   PetscCallMPI(MPI_Comm_rank(comm, &rank));
73   PetscCallMPI(MPI_Comm_size(comm, &size));
74   PetscCall(MatGetBlockSize(Amat_fine, &f_bs));
75   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0));
76   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0));
77   PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat));
78   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1], 0, 0, 0, 0));
79   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PTAP], 0, 0, 0, 0));
80 
81   if (Pcolumnperm) *Pcolumnperm = NULL;
82 
83   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
84   PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL));
85   if (pc_gamg->data_cell_rows > 0) {
86     ncrs = pc_gamg->data_sz / pc_gamg->data_cell_cols / pc_gamg->data_cell_rows;
87   } else {
88     PetscInt bs;
89     PetscCall(MatGetBlockSize(Cmat, &bs));
90     ncrs = ncrs_eq / bs;
91   }
92   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
93   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. */
94 #if defined(PETSC_HAVE_CUDA)
95     PetscShmComm pshmcomm;
96     PetscMPIInt  locrank;
97     MPI_Comm     loccomm;
98     PetscInt     s_nnodes, r_nnodes, new_new_size;
99     cudaError_t  cerr;
100     int          devCount;
101     PetscCall(PetscShmCommGet(comm, &pshmcomm));
102     PetscCall(PetscShmCommGetMpiShmComm(pshmcomm, &loccomm));
103     PetscCallMPI(MPI_Comm_rank(loccomm, &locrank));
104     s_nnodes = !locrank;
105     PetscCall(MPIU_Allreduce(&s_nnodes, &r_nnodes, 1, MPIU_INT, MPI_SUM, comm));
106     PetscCheck((size % r_nnodes) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "odd number of nodes np=%d nnodes%" PetscInt_FMT, size, r_nnodes);
107     devCount = 0;
108     cerr     = cudaGetDeviceCount(&devCount);
109     cudaGetLastError();                         /* Reset the last error */
110     if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */
111       new_new_size = r_nnodes * devCount;
112       new_size     = new_new_size;
113       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));
114     } else {
115       PetscCall(PetscInfo(pc, "%s: With Cuda but no device. Use heuristics.\n", ((PetscObject)pc)->prefix));
116       goto HEURISTIC;
117     }
118 #else
119     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "should not be here");
120 #endif
121   } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) {
122     if (nactive < pc_gamg->level_reduction_factors[pc_gamg->current_level]) {
123       new_size = 1;
124       PetscCall(PetscInfo(pc, "%s: reduction factor too small for %d active processes: reduce to one process\n", ((PetscObject)pc)->prefix, new_size));
125     } else {
126       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]);
127       new_size = nactive / pc_gamg->level_reduction_factors[pc_gamg->current_level];
128       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]));
129     }
130   } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) {
131     new_size = 1;
132     PetscCall(PetscInfo(pc, "%s: Force coarsest grid reduction to %d active processes\n", ((PetscObject)pc)->prefix, new_size));
133   } else {
134     PetscInt ncrs_eq_glob;
135 #if defined(PETSC_HAVE_CUDA)
136   HEURISTIC:
137 #endif
138     PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL));
139     new_size = (PetscMPIInt)((float)ncrs_eq_glob / (float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
140     if (!new_size) new_size = 1;                                                       /* not likely, possible? */
141     else if (new_size >= nactive) new_size = nactive;                                  /* no change, rare */
142     PetscCall(PetscInfo(pc, "%s: Coarse grid reduction from %d to %d active processes\n", ((PetscObject)pc)->prefix, nactive, new_size));
143   }
144   if (new_size == nactive) {
145     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
146     if (new_size < size) {
147       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
148       PetscCall(PetscInfo(pc, "%s: reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n", ((PetscObject)pc)->prefix, nactive));
149       if (pc_gamg->cpu_pin_coarse_grids) {
150         PetscCall(MatBindToCPU(*a_Amat_crs, PETSC_TRUE));
151         PetscCall(MatBindToCPU(*a_P_inout, PETSC_TRUE));
152       }
153     }
154     /* we know that the grid structure can be reused in MatPtAP */
155   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
156     PetscInt *counts, *newproc_idx, ii, jj, kk, strideNew, *tidx, ncrs_new, ncrs_eq_new, nloc_old, expand_factor = 1, rfactor = 1;
157     IS        is_eq_newproc, is_eq_num, is_eq_num_prim, new_eq_indices;
158     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0));
159     nloc_old = ncrs_eq / cr_bs;
160     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);
161     /* get new_size and rfactor */
162     if (pc_gamg->layout_type == PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
163       /* find factor */
164       if (new_size == 1) rfactor = size; /* don't modify */
165       else {
166         PetscReal best_fact = 0.;
167         jj                  = -1;
168         for (kk = 1; kk <= size; kk++) {
169           if (!(size % kk)) { /* a candidate */
170             PetscReal nactpe = (PetscReal)size / (PetscReal)kk, fact = nactpe / (PetscReal)new_size;
171             if (fact > 1.0) fact = 1. / fact; /* keep fact < 1 */
172             if (fact > best_fact) {
173               best_fact = fact;
174               jj        = kk;
175             }
176           }
177         }
178         if (jj != -1) rfactor = jj;
179         else rfactor = 1; /* a prime */
180         if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1;
181         else expand_factor = rfactor;
182       }
183       new_size = size / rfactor; /* make new size one that is factor */
184       if (new_size == nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */
185         *a_Amat_crs = Cmat;
186         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));
187         PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REDUCE], 0, 0, 0, 0));
188         PetscFunctionReturn(PETSC_SUCCESS);
189       }
190     }
191     /* make 'is_eq_newproc' */
192     PetscCall(PetscMalloc1(size, &counts));
193     if (pc_gamg->repart) { /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */
194       Mat adj;
195       PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0));
196       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"));
197       /* get 'adj' */
198       if (cr_bs == 1) {
199         PetscCall(MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
200       } else {
201         /* make a scalar matrix to partition (no Stokes here) */
202         Mat                tMat;
203         PetscInt           Istart_crs, Iend_crs, ncols, jj, Ii;
204         const PetscScalar *vals;
205         const PetscInt    *idx;
206         PetscInt          *d_nnz, *o_nnz, M, N;
207         static PetscInt    llev = 0; /* ugly but just used for debugging */
208         MatType            mtype;
209 
210         PetscCall(PetscMalloc2(ncrs, &d_nnz, ncrs, &o_nnz));
211         PetscCall(MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs));
212         PetscCall(MatGetSize(Cmat, &M, &N));
213         for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) {
214           PetscCall(MatGetRow(Cmat, Ii, &ncols, NULL, NULL));
215           d_nnz[jj] = ncols / cr_bs;
216           o_nnz[jj] = ncols / cr_bs;
217           PetscCall(MatRestoreRow(Cmat, Ii, &ncols, NULL, NULL));
218           if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs;
219           if (o_nnz[jj] > (M / cr_bs - ncrs)) o_nnz[jj] = M / cr_bs - ncrs;
220         }
221 
222         PetscCall(MatGetType(Amat_fine, &mtype));
223         PetscCall(MatCreate(comm, &tMat));
224         PetscCall(MatSetSizes(tMat, ncrs, ncrs, PETSC_DETERMINE, PETSC_DETERMINE));
225         PetscCall(MatSetType(tMat, mtype));
226         PetscCall(MatSeqAIJSetPreallocation(tMat, 0, d_nnz));
227         PetscCall(MatMPIAIJSetPreallocation(tMat, 0, d_nnz, 0, o_nnz));
228         PetscCall(PetscFree2(d_nnz, o_nnz));
229 
230         for (ii = Istart_crs; ii < Iend_crs; ii++) {
231           PetscInt dest_row = ii / cr_bs;
232           PetscCall(MatGetRow(Cmat, ii, &ncols, &idx, &vals));
233           for (jj = 0; jj < ncols; jj++) {
234             PetscInt    dest_col = idx[jj] / cr_bs;
235             PetscScalar v        = 1.0;
236             PetscCall(MatSetValues(tMat, 1, &dest_row, 1, &dest_col, &v, ADD_VALUES));
237           }
238           PetscCall(MatRestoreRow(Cmat, ii, &ncols, &idx, &vals));
239         }
240         PetscCall(MatAssemblyBegin(tMat, MAT_FINAL_ASSEMBLY));
241         PetscCall(MatAssemblyEnd(tMat, MAT_FINAL_ASSEMBLY));
242 
243         if (llev++ == -1) {
244           PetscViewer viewer;
245           char        fname[32];
246           PetscCall(PetscSNPrintf(fname, sizeof(fname), "part_mat_%" PetscInt_FMT ".mat", llev));
247           PetscCall(PetscViewerBinaryOpen(comm, fname, FILE_MODE_WRITE, &viewer));
248           PetscCall(MatView(tMat, viewer));
249           PetscCall(PetscViewerDestroy(&viewer));
250         }
251         PetscCall(MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
252         PetscCall(MatDestroy(&tMat));
253       } /* create 'adj' */
254 
255       { /* partition: get newproc_idx */
256         char            prefix[256];
257         const char     *pcpre;
258         const PetscInt *is_idx;
259         MatPartitioning mpart;
260         IS              proc_is;
261 
262         PetscCall(MatPartitioningCreate(comm, &mpart));
263         PetscCall(MatPartitioningSetAdjacency(mpart, adj));
264         PetscCall(PCGetOptionsPrefix(pc, &pcpre));
265         PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : ""));
266         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix));
267         PetscCall(MatPartitioningSetFromOptions(mpart));
268         PetscCall(MatPartitioningSetNParts(mpart, new_size));
269         PetscCall(MatPartitioningApply(mpart, &proc_is));
270         PetscCall(MatPartitioningDestroy(&mpart));
271 
272         /* collect IS info */
273         PetscCall(PetscMalloc1(ncrs_eq, &newproc_idx));
274         PetscCall(ISGetIndices(proc_is, &is_idx));
275         for (kk = jj = 0; kk < nloc_old; kk++) {
276           for (ii = 0; ii < cr_bs; ii++, jj++) { newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ }
277         }
278         PetscCall(ISRestoreIndices(proc_is, &is_idx));
279         PetscCall(ISDestroy(&proc_is));
280       }
281       PetscCall(MatDestroy(&adj));
282 
283       PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc));
284       PetscCall(PetscFree(newproc_idx));
285       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_REPART], 0, 0, 0, 0));
286     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
287       PetscInt targetPE;
288       PetscCheck(new_size != nactive, PETSC_COMM_SELF, PETSC_ERR_PLIB, "new_size==nactive. Should not happen");
289       PetscCall(PetscInfo(pc, "%s: Number of equations (loc) %" PetscInt_FMT " with simple aggregation\n", ((PetscObject)pc)->prefix, ncrs_eq));
290       targetPE = (rank / rfactor) * expand_factor;
291       PetscCall(ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc));
292     } /* end simple 'is_eq_newproc' */
293 
294     /*
295       Create an index set from the is_eq_newproc index set to indicate the mapping TO
296     */
297     PetscCall(ISPartitioningToNumbering(is_eq_newproc, &is_eq_num));
298     is_eq_num_prim = is_eq_num;
299     /*
300       Determine how many equations/vertices are assigned to each processor
301     */
302     PetscCall(ISPartitioningCount(is_eq_newproc, size, counts));
303     ncrs_eq_new = counts[rank];
304     PetscCall(ISDestroy(&is_eq_newproc));
305     ncrs_new = ncrs_eq_new / cr_bs;
306 
307     PetscCall(PetscFree(counts));
308     /* 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 */
309     {
310       Vec             src_crd, dest_crd;
311       const PetscInt *idx, ndata_rows = pc_gamg->data_cell_rows, ndata_cols = pc_gamg->data_cell_cols, node_data_sz = ndata_rows * ndata_cols;
312       VecScatter      vecscat;
313       PetscScalar    *array;
314       IS              isscat;
315       /* move data (for primal equations only) */
316       /* Create a vector to contain the newly ordered element information */
317       PetscCall(VecCreate(comm, &dest_crd));
318       PetscCall(VecSetSizes(dest_crd, node_data_sz * ncrs_new, PETSC_DECIDE));
319       PetscCall(VecSetType(dest_crd, VECSTANDARD)); /* this is needed! */
320       /*
321         There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having
322         a block size of ...).  Note, ISs are expanded into equation space by 'cr_bs'.
323       */
324       PetscCall(PetscMalloc1(ncrs * node_data_sz, &tidx));
325       PetscCall(ISGetIndices(is_eq_num_prim, &idx));
326       for (ii = 0, jj = 0; ii < ncrs; ii++) {
327         PetscInt id = idx[ii * cr_bs] / cr_bs; /* get node back */
328         for (kk = 0; kk < node_data_sz; kk++, jj++) tidx[jj] = id * node_data_sz + kk;
329       }
330       PetscCall(ISRestoreIndices(is_eq_num_prim, &idx));
331       PetscCall(ISCreateGeneral(comm, node_data_sz * ncrs, tidx, PETSC_COPY_VALUES, &isscat));
332       PetscCall(PetscFree(tidx));
333       /*
334         Create a vector to contain the original vertex information for each element
335       */
336       PetscCall(VecCreateSeq(PETSC_COMM_SELF, node_data_sz * ncrs, &src_crd));
337       for (jj = 0; jj < ndata_cols; jj++) {
338         const PetscInt stride0 = ncrs * pc_gamg->data_cell_rows;
339         for (ii = 0; ii < ncrs; ii++) {
340           for (kk = 0; kk < ndata_rows; kk++) {
341             PetscInt    ix = ii * ndata_rows + kk + jj * stride0, jx = ii * node_data_sz + kk * ndata_cols + jj;
342             PetscScalar tt = (PetscScalar)pc_gamg->data[ix];
343             PetscCall(VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES));
344           }
345         }
346       }
347       PetscCall(VecAssemblyBegin(src_crd));
348       PetscCall(VecAssemblyEnd(src_crd));
349       /*
350         Scatter the element vertex information (still in the original vertex ordering)
351         to the correct processor
352       */
353       PetscCall(VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat));
354       PetscCall(ISDestroy(&isscat));
355       PetscCall(VecScatterBegin(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD));
356       PetscCall(VecScatterEnd(vecscat, src_crd, dest_crd, INSERT_VALUES, SCATTER_FORWARD));
357       PetscCall(VecScatterDestroy(&vecscat));
358       PetscCall(VecDestroy(&src_crd));
359       /*
360         Put the element vertex data into a new allocation of the gdata->ele
361       */
362       PetscCall(PetscFree(pc_gamg->data));
363       PetscCall(PetscMalloc1(node_data_sz * ncrs_new, &pc_gamg->data));
364 
365       pc_gamg->data_sz = node_data_sz * ncrs_new;
366       strideNew        = ncrs_new * ndata_rows;
367 
368       PetscCall(VecGetArray(dest_crd, &array));
369       for (jj = 0; jj < ndata_cols; jj++) {
370         for (ii = 0; ii < ncrs_new; ii++) {
371           for (kk = 0; kk < ndata_rows; kk++) {
372             PetscInt ix = ii * ndata_rows + kk + jj * strideNew, jx = ii * node_data_sz + kk * ndata_cols + jj;
373             pc_gamg->data[ix] = PetscRealPart(array[jx]);
374           }
375         }
376       }
377       PetscCall(VecRestoreArray(dest_crd, &array));
378       PetscCall(VecDestroy(&dest_crd));
379     }
380     /* move A and P (columns) with new layout */
381     /*
382       Invert for MatCreateSubMatrix
383     */
384     PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices));
385     PetscCall(ISSort(new_eq_indices)); /* is this needed? */
386     PetscCall(ISSetBlockSize(new_eq_indices, cr_bs));
387     if (is_eq_num != is_eq_num_prim) { PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */ }
388     if (Pcolumnperm) {
389       PetscCall(PetscObjectReference((PetscObject)new_eq_indices));
390       *Pcolumnperm = new_eq_indices;
391     }
392     PetscCall(ISDestroy(&is_eq_num));
393 
394     /* 'a_Amat_crs' output */
395     {
396       Mat       mat;
397       PetscBool isset, isspd, isher;
398 #if !defined(PETSC_USE_COMPLEX)
399       PetscBool issym;
400 #endif
401 
402       PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat));
403       PetscCall(MatIsSPDKnown(Cmat, &isset, &isspd)); // like MatPropagateSymmetryOptions, but should set MAT_STRUCTURALLY_SYMMETRIC ?
404       if (isset) PetscCall(MatSetOption(mat, MAT_SPD, isspd));
405       else {
406         PetscCall(MatIsHermitianKnown(Cmat, &isset, &isher));
407         if (isset) PetscCall(MatSetOption(mat, MAT_HERMITIAN, isher));
408         else {
409 #if !defined(PETSC_USE_COMPLEX)
410           PetscCall(MatIsSymmetricKnown(Cmat, &isset, &issym));
411           if (isset) PetscCall(MatSetOption(mat, MAT_SYMMETRIC, issym));
412 #endif
413         }
414       }
415       *a_Amat_crs = mat;
416     }
417     PetscCall(MatDestroy(&Cmat));
418 
419     /* prolongator */
420     {
421       IS       findices;
422       PetscInt Istart, Iend;
423       Mat      Pnew;
424 
425       PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend));
426       PetscCall(ISCreateStride(comm, Iend - Istart, Istart, 1, &findices));
427       PetscCall(ISSetBlockSize(findices, f_bs));
428       PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew));
429       PetscCall(ISDestroy(&findices));
430       PetscCall(MatSetOption(Pnew, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE));
431 
432       PetscCall(MatDestroy(a_P_inout));
433 
434       /* output - repartitioned */
435       *a_P_inout = Pnew;
436     }
437     PetscCall(ISDestroy(&new_eq_indices));
438 
439     *a_nactive_proc = new_size; /* output */
440 
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 }
464 
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;
472 
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 }
495 
496 /*
497    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
498                     by setting data structures and options.
499 
500    Input Parameter:
501 .  pc - the preconditioner context
502 
503 */
504 static 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], cr_bs;
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
519 
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));
540 
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 initial setup, reuse matrix level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level));
554           } else {
555             PetscCall(PetscInfo(pc, "%s: RAP after initial setup, with repartitioning (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           // check for redoing eigen estimates
563           if (pc_gamg->recompute_esteig) {
564             PetscBool ischeb;
565             KSP       smoother;
566             PetscCall(PCMGGetSmoother(pc, level + 1, &smoother));
567             PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb));
568             if (ischeb) {
569               KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data;
570               cheb->emin_provided = 0;
571               cheb->emax_provided = 0;
572             }
573             /* we could call PetscCall(KSPChebyshevSetEigenvalues(smoother, 0, 0)); but the logic does not work properly */
574           }
575           // inc
576           dB = B;
577 #if defined(GAMG_STAGES)
578           PetscCall(PetscLogStagePop());
579 #endif
580         }
581       }
582 
583       PetscCall(PCSetUp_MG(pc));
584       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0));
585       PetscFunctionReturn(PETSC_SUCCESS);
586     }
587   }
588 
589   if (!pc_gamg->data) {
590     if (pc_gamg->orig_data) {
591       PetscCall(MatGetBlockSize(Pmat, &bs));
592       PetscCall(MatGetLocalSize(Pmat, &qq, NULL));
593 
594       pc_gamg->data_sz        = (qq / bs) * pc_gamg->orig_data_cell_rows * pc_gamg->orig_data_cell_cols;
595       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
596       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
597 
598       PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data));
599       for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
600     } else {
601       PetscCheck(pc_gamg->ops->createdefaultdata, comm, PETSC_ERR_PLIB, "'createdefaultdata' not set(?) need to support NULL data");
602       PetscCall(pc_gamg->ops->createdefaultdata(pc, Pmat));
603     }
604   }
605 
606   /* cache original data for reuse */
607   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
608     PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data));
609     for (qq = 0; qq < pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
610     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
611     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
612   }
613 
614   /* get basic dims */
615   PetscCall(MatGetBlockSize(Pmat, &bs));
616   PetscCall(MatGetSize(Pmat, &M, NULL));
617 
618 #if defined(PETSC_USE_INFO)
619   PetscCall(MatGetInfo(Pmat, MAT_GLOBAL_SUM, &info)); /* global reduction */
620   nnz0   = info.nz_used;
621   nnztot = info.nz_used;
622 #endif
623   PetscCall(PetscInfo(pc, "%s: level %d) N=%" PetscInt_FMT ", n data rows=%" PetscInt_FMT ", n data cols=%" PetscInt_FMT ", nnz/row (ave)=%" PetscInt_FMT ", block size %d, np=%d\n", ((PetscObject)pc)->prefix, 0, M, pc_gamg->data_cell_rows, pc_gamg->data_cell_cols, (PetscInt)(nnz0 / (PetscReal)M + 0.5), (int)bs, size));
624 
625   /* Get A_i and R_i */
626   for (level = 0, Aarr[0] = Pmat, nactivepe = size; level < (pc_gamg->Nlevels - 1) && (level == 0 || M > pc_gamg->coarse_eq_limit); level++) {
627     pc_gamg->current_level = level;
628     level1                 = level + 1;
629 #if defined(GAMG_STAGES)
630     if (!gamg_stages[level]) {
631       char str[32];
632       PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", (int)level));
633       PetscCall(PetscLogStageRegister(str, &gamg_stages[level]));
634     }
635     PetscCall(PetscLogStagePush(gamg_stages[level]));
636 #endif
637     /* construct prolongator - Parr[level1] */
638     if (level == 0 && pc_gamg->injection_index_size > 0) {
639       Mat      Prol;
640       MatType  mtype;
641       PetscInt prol_m, prol_n, Prol_N = (M / bs) * pc_gamg->injection_index_size, Istart, Iend, nn, row;
642       PetscCall(PetscInfo(pc, "Create fine grid injection space prolongation %" PetscInt_FMT " x %" PetscInt_FMT ". %s\n", M, Prol_N, pc_gamg->data ? "delete null space data" : ""));
643       PetscCall(MatGetOwnershipRange(Pmat, &Istart, &Iend));
644       PetscCall(MatGetLocalSize(Pmat, &prol_m, NULL)); // rows m x n
645       prol_n = (prol_m / bs) * pc_gamg->injection_index_size;
646       PetscCheck(pc_gamg->injection_index_size < bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Injection size %d must be less that block size %d", (int)pc_gamg->injection_index_size, (int)bs);
647       PetscCall(MatGetType(Pmat, &mtype));
648       PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &Prol));
649       PetscCall(MatSetBlockSizes(Prol, bs, pc_gamg->injection_index_size));
650       PetscCall(MatSetSizes(Prol, prol_m, prol_n, M, Prol_N));
651       PetscCall(MatSetType(Prol, mtype));
652 #if PetscDefined(HAVE_DEVICE)
653       PetscBool flg;
654       PetscCall(MatBoundToCPU(Pmat, &flg));
655       PetscCall(MatBindToCPU(Prol, flg));
656       if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE));
657 #endif
658       PetscCall(MatSeqAIJSetPreallocation(Prol, 1, NULL));
659       PetscCall(MatMPIAIJSetPreallocation(Prol, 1, NULL, 0, NULL));
660       // set I \kron [1, 1, ... ]^T
661       for (PetscInt ii = Istart, col = (Istart / bs) * pc_gamg->injection_index_size; ii < Iend; ii += bs) {
662         const PetscScalar one = 1;
663         for (PetscInt jj = 0; jj < pc_gamg->injection_index_size; jj++, col++) {
664           PetscInt row = ii + pc_gamg->injection_index[jj];
665           PetscCall(MatSetValues(Prol, 1, &row, 1, &col, &one, INSERT_VALUES));
666         }
667       }
668       PetscCall(MatAssemblyBegin(Prol, MAT_FINAL_ASSEMBLY));
669       PetscCall(MatAssemblyEnd(Prol, MAT_FINAL_ASSEMBLY));
670       PetscCall(MatViewFromOptions(Prol, NULL, "-mat_view_injection"));
671       PetscCall(MatGetBlockSizes(Prol, NULL, &cr_bs)); // column size
672       Parr[level1] = Prol;
673       // can not deal with null space -- with array of 'injection cols' we could take 'injection rows and 'injection cols' to 'data'
674       if (pc_gamg->data) {
675         pc_gamg->data_cell_cols      = pc_gamg->injection_index_size;
676         pc_gamg->data_cell_rows      = pc_gamg->injection_index_size;
677         pc_gamg->orig_data_cell_cols = 0;
678         pc_gamg->orig_data_cell_rows = 0;
679         PetscCall(PetscFree(pc_gamg->data));
680         pc_gamg->data_sz = pc_gamg->injection_index_size * prol_n;
681         PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data));
682         for (row = nn = 0; row < prol_n; row += pc_gamg->injection_index_size) {
683           for (int jj = 0; jj < pc_gamg->injection_index_size; jj++) {
684             int idx = row * pc_gamg->injection_index_size + jj * pc_gamg->injection_index_size;
685             for (int kk = 0; kk < pc_gamg->injection_index_size; kk++, nn++) { pc_gamg->data[idx + kk] = (jj == kk) ? 1 : 0; }
686           }
687         }
688         PetscCheck(nn == pc_gamg->data_sz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nn != pc_gamg->data_sz %" PetscInt_FMT " %" PetscInt_FMT, pc_gamg->data_sz, nn);
689       }
690     } else {
691       Mat               Gmat, mat;
692       PetscCoarsenData *agg_lists;
693       Mat               Prol11;
694 
695       PetscCall(PCGAMGCreateGraph(pc, Aarr[level], &Gmat));
696       PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists)); // Gmat may have ghosts for QR aggregates not in matrix
697       PetscCall(PetscCDGetMat(agg_lists, &mat));
698       if (!mat) PetscCall(PetscCDSetMat(agg_lists, Gmat));
699       PetscCall(pc_gamg->ops->prolongator(pc, Aarr[level], agg_lists, &Prol11));
700       /* could have failed to create new level */
701       if (Prol11) {
702         const char *prefix;
703         char        addp[32];
704 
705         /* get new block size of coarse matrices */
706         PetscCall(MatGetBlockSizes(Prol11, NULL, &cr_bs)); // column size
707 
708         if (pc_gamg->ops->optprolongator) {
709           /* smooth */
710           PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11));
711         }
712 
713         if (pc_gamg->use_aggs_in_asm) {
714           PetscInt bs;
715           PetscCall(MatGetBlockSizes(Prol11, &bs, NULL)); // row block size
716           PetscCall(PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]));
717           PetscCall(PetscInfo(pc, "%d: %" PetscInt_FMT " ASM local domains,  bs = %d\n", (int)level, nASMBlocksArr[level], (int)bs));
718         } else if (pc_gamg->asm_hem_aggs) {
719           const char *prefix;
720           PetscInt    bs;
721 
722           /*
723              Do not use aggs created for defining coarser problems, instead create aggs specifically to use
724              to define PCASM blocks.
725           */
726           PetscCall(PetscCDGetMat(agg_lists, &mat));
727           if (mat == Gmat) PetscCall(PetscCDClearMat(agg_lists)); // take the Mat away from the list (yuck)
728           PetscCall(PetscCDDestroy(agg_lists));
729           PetscCall(PetscInfo(pc, "HEM ASM passes = %d\n", (int)pc_gamg->asm_hem_aggs));
730           PetscCall(MatCoarsenDestroy(&pc_gamg->asm_crs));
731           PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg->asm_crs));
732           PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix));
733           PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg->asm_crs, prefix));
734           PetscCall(MatCoarsenSetFromOptions(pc_gamg->asm_crs)); // get strength args
735           PetscCall(MatCoarsenSetType(pc_gamg->asm_crs, MATCOARSENHEM));
736           PetscCall(MatCoarsenSetMaximumIterations(pc_gamg->asm_crs, pc_gamg->asm_hem_aggs));
737           PetscCall(MatCoarsenSetAdjacency(pc_gamg->asm_crs, Gmat));
738           PetscCall(MatCoarsenSetStrictAggs(pc_gamg->asm_crs, PETSC_TRUE));
739           PetscCall(MatCoarsenApply(pc_gamg->asm_crs));
740           PetscCall(MatCoarsenGetData(pc_gamg->asm_crs, &agg_lists)); /* output */
741           // create aggregates
742           PetscCall(MatGetBlockSizes(Aarr[level], &bs, NULL)); // row block size
743           PetscCall(PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]));
744         }
745         PetscCall(PCGetOptionsPrefix(pc, &prefix));
746         PetscCall(MatSetOptionsPrefix(Prol11, prefix));
747         PetscCall(PetscSNPrintf(addp, sizeof(addp), "pc_gamg_prolongator_%d_", (int)level));
748         PetscCall(MatAppendOptionsPrefix(Prol11, addp));
749         /* Always generate the transpose with CUDA
750            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
751         PetscCall(MatSetOption(Prol11, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE));
752         PetscCall(MatSetFromOptions(Prol11));
753         Parr[level1] = Prol11;
754       } else Parr[level1] = NULL; /* failed to coarsen */
755       PetscCall(PetscCDGetMat(agg_lists, &mat));
756       if (mat == Gmat) PetscCall(PetscCDClearMat(agg_lists)); // take the Mat away from the list (yuck)
757       PetscCall(MatDestroy(&Gmat));
758       PetscCall(PetscCDDestroy(agg_lists));
759     } /* construct prolongator scope */
760     if (level == 0) Aarr[0] = Pmat; /* use Pmat for finest level setup */
761     if (!Parr[level1]) {            /* failed to coarsen */
762       PetscCall(PetscInfo(pc, "%s: Stop gridding, level %" PetscInt_FMT "\n", ((PetscObject)pc)->prefix, level));
763 #if defined(GAMG_STAGES)
764       PetscCall(PetscLogStagePop());
765 #endif
766       break;
767     }
768     PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */
769     PetscCheck(!is_last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Is last ?");
770     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
771     if (level1 == pc_gamg->Nlevels - 1) is_last = PETSC_TRUE;
772     if (level == PETSC_MG_MAXLEVELS - 2) is_last = PETSC_TRUE;
773     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0));
774     PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], cr_bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last));
775     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_LEVEL], 0, 0, 0, 0));
776 
777     PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */
778 #if defined(PETSC_USE_INFO)
779     PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info));
780     nnztot += info.nz_used;
781 #endif
782     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));
783 
784 #if defined(GAMG_STAGES)
785     PetscCall(PetscLogStagePop());
786 #endif
787     /* stop if one node or one proc -- could pull back for singular problems */
788     if ((pc_gamg->data_cell_cols && M / pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M / bs < 2)) {
789       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));
790       level++;
791       break;
792     } else if (level == PETSC_MG_MAXLEVELS - 2) { /* stop if we are limited by PC_MG_MAXLEVELS */
793       PetscCall(PetscInfo(pc, "%s: HARD stop of coarsening on level %" PetscInt_FMT ".  PC_MG_MAXLEVELS reached\n", ((PetscObject)pc)->prefix, level));
794       level++;
795       break;
796     }
797   } /* levels */
798   PetscCall(PetscFree(pc_gamg->data));
799 
800   PetscCall(PetscInfo(pc, "%s: %" PetscInt_FMT " levels, operator complexity = %g\n", ((PetscObject)pc)->prefix, level + 1, nnztot / nnz0));
801   pc_gamg->Nlevels = level + 1;
802   fine_level       = level;
803   PetscCall(PCMGSetLevels(pc, pc_gamg->Nlevels, NULL));
804 
805   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
806 
807     /* set default smoothers & set operators */
808     for (lidx = 1, level = pc_gamg->Nlevels - 2; lidx <= fine_level; lidx++, level--) {
809       KSP smoother;
810       PC  subpc;
811 
812       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
813       PetscCall(KSPGetPC(smoother, &subpc));
814 
815       PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
816       /* set ops */
817       PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level]));
818       PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level + 1]));
819 
820       /* set defaults */
821       PetscCall(KSPSetType(smoother, KSPCHEBYSHEV));
822 
823       /* set blocks for ASM smoother that uses the 'aggregates' */
824       if (pc_gamg->use_aggs_in_asm || pc_gamg->asm_hem_aggs) {
825         PetscInt sz;
826         IS      *iss;
827 
828         sz  = nASMBlocksArr[level];
829         iss = ASMLocalIDsArr[level];
830         PetscCall(PCSetType(subpc, PCASM));
831         PetscCall(PCASMSetOverlap(subpc, 0));
832         PetscCall(PCASMSetType(subpc, PC_ASM_BASIC));
833         if (!sz) {
834           IS is;
835           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is));
836           PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is));
837           PetscCall(ISDestroy(&is));
838         } else {
839           PetscInt kk;
840           PetscCall(PCASMSetLocalSubdomains(subpc, sz, iss, NULL));
841           for (kk = 0; kk < sz; kk++) PetscCall(ISDestroy(&iss[kk]));
842           PetscCall(PetscFree(iss));
843         }
844         ASMLocalIDsArr[level] = NULL;
845         nASMBlocksArr[level]  = 0;
846       } else {
847         PetscCall(PCSetType(subpc, PCJACOBI));
848       }
849     }
850     {
851       /* coarse grid */
852       KSP      smoother, *k2;
853       PC       subpc, pc2;
854       PetscInt ii, first;
855       Mat      Lmat = Aarr[(level = pc_gamg->Nlevels - 1)];
856       lidx          = 0;
857 
858       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
859       PetscCall(KSPSetOperators(smoother, Lmat, Lmat));
860       if (!pc_gamg->use_parallel_coarse_grid_solver) {
861         PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
862         PetscCall(KSPGetPC(smoother, &subpc));
863         PetscCall(PCSetType(subpc, PCBJACOBI));
864         PetscCall(PCSetUp(subpc));
865         PetscCall(PCBJacobiGetSubKSP(subpc, &ii, &first, &k2));
866         PetscCheck(ii == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " is not one", ii);
867         PetscCall(KSPGetPC(k2[0], &pc2));
868         PetscCall(PCSetType(pc2, PCLU));
869         PetscCall(PCFactorSetShiftType(pc2, MAT_SHIFT_INBLOCKS));
870         PetscCall(KSPSetTolerances(k2[0], PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1));
871         PetscCall(KSPSetType(k2[0], KSPPREONLY));
872       }
873     }
874 
875     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
876     PetscObjectOptionsBegin((PetscObject)pc);
877     PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject));
878     PetscOptionsEnd();
879     PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
880 
881     /* set cheby eigen estimates from SA to use in the solver */
882     if (pc_gamg->use_sa_esteig) {
883       for (lidx = 1, level = pc_gamg->Nlevels - 2; level >= 0; lidx++, level--) {
884         KSP       smoother;
885         PetscBool ischeb;
886 
887         PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
888         PetscCall(PetscObjectTypeCompare((PetscObject)smoother, KSPCHEBYSHEV, &ischeb));
889         if (ischeb) {
890           KSP_Chebyshev *cheb = (KSP_Chebyshev *)smoother->data;
891 
892           // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG
893           if (mg->max_eigen_DinvA[level] > 0) {
894             // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it.
895             // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.)
896             PetscReal emax, emin;
897 
898             emin = mg->min_eigen_DinvA[level];
899             emax = mg->max_eigen_DinvA[level];
900             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));
901             cheb->emin_provided = emin;
902             cheb->emax_provided = emax;
903           }
904         }
905       }
906     }
907 
908     PetscCall(PCSetUp_MG(pc));
909 
910     /* clean up */
911     for (level = 1; level < pc_gamg->Nlevels; level++) {
912       PetscCall(MatDestroy(&Parr[level]));
913       PetscCall(MatDestroy(&Aarr[level]));
914     }
915   } else {
916     KSP smoother;
917 
918     PetscCall(PetscInfo(pc, "%s: One level solver used (system is seen as DD). Using default solver.\n", ((PetscObject)pc)->prefix));
919     PetscCall(PCMGGetSmoother(pc, 0, &smoother));
920     PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0]));
921     PetscCall(KSPSetType(smoother, KSPPREONLY));
922     PetscCall(PCSetUp_MG(pc));
923   }
924   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SETUP], 0, 0, 0, 0));
925   PetscFunctionReturn(PETSC_SUCCESS);
926 }
927 
928 /*
929  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
930    that was created with PCCreate_GAMG().
931 
932    Input Parameter:
933 .  pc - the preconditioner context
934 
935    Application Interface Routine: PCDestroy()
936 */
937 PetscErrorCode PCDestroy_GAMG(PC pc)
938 {
939   PC_MG   *mg      = (PC_MG *)pc->data;
940   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
941 
942   PetscFunctionBegin;
943   PetscCall(PCReset_GAMG(pc));
944   if (pc_gamg->ops->destroy) PetscCall((*pc_gamg->ops->destroy)(pc));
945   PetscCall(PetscFree(pc_gamg->ops));
946   PetscCall(PetscFree(pc_gamg->gamg_type_name));
947   PetscCall(PetscFree(pc_gamg));
948   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", NULL));
949   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", NULL));
950   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", NULL));
951   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", NULL));
952   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", NULL));
953   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRecomputeEstEig_C", NULL));
954   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", NULL));
955   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", NULL));
956   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetParallelCoarseGridSolve_C", NULL));
957   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", NULL));
958   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", NULL));
959   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", NULL));
960   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", NULL));
961   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", NULL));
962   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", NULL));
963   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", NULL));
964   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", NULL));
965   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetHEM_C", NULL));
966   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetInjectionIndices_C", NULL));
967   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetInjectionIndex_C", NULL));
968   PetscCall(PCDestroy_MG(pc));
969   PetscFunctionReturn(PETSC_SUCCESS);
970 }
971 
972 /*@
973   PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction in `PCGAMG`
974 
975   Logically Collective
976 
977   Input Parameters:
978 + pc - the preconditioner context
979 - n  - the number of equations
980 
981   Options Database Key:
982 . -pc_gamg_process_eq_limit <limit> - set the limit
983 
984   Level: intermediate
985 
986   Note:
987   `PCGAMG` will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
988   that has degrees of freedom
989 
990   Developer Note:
991   Should be named `PCGAMGSetProcessEquationLimit()`.
992 
993 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`
994 @*/
995 PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n)
996 {
997   PetscFunctionBegin;
998   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
999   PetscTryMethod(pc, "PCGAMGSetProcEqLim_C", (PC, PetscInt), (pc, n));
1000   PetscFunctionReturn(PETSC_SUCCESS);
1001 }
1002 
1003 static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
1004 {
1005   PC_MG   *mg      = (PC_MG *)pc->data;
1006   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1007 
1008   PetscFunctionBegin;
1009   if (n > 0) pc_gamg->min_eq_proc = n;
1010   PetscFunctionReturn(PETSC_SUCCESS);
1011 }
1012 
1013 /*@
1014   PCGAMGSetCoarseEqLim - Set maximum number of equations on the coarsest grid of `PCGAMG`
1015 
1016   Collective
1017 
1018   Input Parameters:
1019 + pc - the preconditioner context
1020 - n  - maximum number of equations to aim for
1021 
1022   Options Database Key:
1023 . -pc_gamg_coarse_eq_limit <limit> - set the limit
1024 
1025   Level: intermediate
1026 
1027   Note:
1028   For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
1029   has less than 1000 unknowns.
1030 
1031 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`, `PCGAMGSetRepartition()`,
1032           `PCGAMGSetParallelCoarseGridSolve()`
1033 @*/
1034 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
1035 {
1036   PetscFunctionBegin;
1037   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1038   PetscTryMethod(pc, "PCGAMGSetCoarseEqLim_C", (PC, PetscInt), (pc, n));
1039   PetscFunctionReturn(PETSC_SUCCESS);
1040 }
1041 
1042 static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
1043 {
1044   PC_MG   *mg      = (PC_MG *)pc->data;
1045   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1046 
1047   PetscFunctionBegin;
1048   if (n > 0) pc_gamg->coarse_eq_limit = n;
1049   PetscFunctionReturn(PETSC_SUCCESS);
1050 }
1051 
1052 /*@
1053   PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids when reducing the number of MPI ranks to use
1054 
1055   Collective
1056 
1057   Input Parameters:
1058 + pc - the preconditioner context
1059 - n  - `PETSC_TRUE` or `PETSC_FALSE`
1060 
1061   Options Database Key:
1062 . -pc_gamg_repartition <true,false> - turn on the repartitioning
1063 
1064   Level: intermediate
1065 
1066   Note:
1067   This will generally improve the loading balancing of the work on each level so the solves will be faster but it adds to the
1068   preconditioner setup time.
1069 
1070 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetRankReductionFactors()`
1071 @*/
1072 PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
1073 {
1074   PetscFunctionBegin;
1075   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1076   PetscTryMethod(pc, "PCGAMGSetRepartition_C", (PC, PetscBool), (pc, n));
1077   PetscFunctionReturn(PETSC_SUCCESS);
1078 }
1079 
1080 static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
1081 {
1082   PC_MG   *mg      = (PC_MG *)pc->data;
1083   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1084 
1085   PetscFunctionBegin;
1086   pc_gamg->repart = n;
1087   PetscFunctionReturn(PETSC_SUCCESS);
1088 }
1089 
1090 /*@
1091   PCGAMGSetUseSAEstEig - Use the eigen estimate from smoothed aggregation for the Chebyshev smoother during the solution process
1092 
1093   Collective
1094 
1095   Input Parameters:
1096 + pc - the preconditioner context
1097 - b  - flag
1098 
1099   Options Database Key:
1100 . -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate
1101 
1102   Level: advanced
1103 
1104   Notes:
1105   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$.
1106   Eigenvalue estimates (based on a few `PCCG` or `PCGMRES` iterations) are computed to choose $\omega$ so that this is a stable smoothing operation.
1107   If `KSPCHEBYSHEV` with `PCJACOBI` (diagonal) preconditioning is used for smoothing on the finest level, then the eigenvalue estimates
1108   can be reused during the solution process.
1109   This option is only used when the smoother uses `PCJACOBI`, and should be turned off when a different `PCJacobiType` is used.
1110   It became default in PETSc 3.17.
1111 
1112 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`, `PCGAMGSetRecomputeEstEig()`
1113 @*/
1114 PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool b)
1115 {
1116   PetscFunctionBegin;
1117   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1118   PetscTryMethod(pc, "PCGAMGSetUseSAEstEig_C", (PC, PetscBool), (pc, b));
1119   PetscFunctionReturn(PETSC_SUCCESS);
1120 }
1121 
1122 static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool b)
1123 {
1124   PC_MG   *mg      = (PC_MG *)pc->data;
1125   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1126 
1127   PetscFunctionBegin;
1128   pc_gamg->use_sa_esteig = b;
1129   PetscFunctionReturn(PETSC_SUCCESS);
1130 }
1131 
1132 /*@
1133   PCGAMGSetRecomputeEstEig - Set flag for Chebyshev smoothers to recompute the eigen estimates when a new matrix is used
1134 
1135   Collective
1136 
1137   Input Parameters:
1138 + pc - the preconditioner context
1139 - b  - flag, default is `PETSC_TRUE`
1140 
1141   Options Database Key:
1142 . -pc_gamg_recompute_esteig <true> - use the eigen estimate
1143 
1144   Level: advanced
1145 
1146   Note:
1147   If the matrix changes only slightly in a new solve using ``PETSC_FALSE`` will save time in the setting up of the preconditioner
1148   and may not affect the solution time much.
1149 
1150 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`
1151 @*/
1152 PetscErrorCode PCGAMGSetRecomputeEstEig(PC pc, PetscBool b)
1153 {
1154   PetscFunctionBegin;
1155   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1156   PetscTryMethod(pc, "PCGAMGSetRecomputeEstEig_C", (PC, PetscBool), (pc, b));
1157   PetscFunctionReturn(PETSC_SUCCESS);
1158 }
1159 
1160 static PetscErrorCode PCGAMGSetRecomputeEstEig_GAMG(PC pc, PetscBool b)
1161 {
1162   PC_MG   *mg      = (PC_MG *)pc->data;
1163   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1164 
1165   PetscFunctionBegin;
1166   pc_gamg->recompute_esteig = b;
1167   PetscFunctionReturn(PETSC_SUCCESS);
1168 }
1169 
1170 /*@
1171   PCGAMGSetEigenvalues - Set WHAT eigenvalues WHY?
1172 
1173   Collective
1174 
1175   Input Parameters:
1176 + pc   - the preconditioner context
1177 . emax - max eigenvalue
1178 - emin - min eigenvalue
1179 
1180   Options Database Key:
1181 . -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues
1182 
1183   Level: intermediate
1184 
1185 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetUseSAEstEig()`
1186 @*/
1187 PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax, PetscReal emin)
1188 {
1189   PetscFunctionBegin;
1190   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1191   PetscTryMethod(pc, "PCGAMGSetEigenvalues_C", (PC, PetscReal, PetscReal), (pc, emax, emin));
1192   PetscFunctionReturn(PETSC_SUCCESS);
1193 }
1194 
1195 static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc, PetscReal emax, PetscReal emin)
1196 {
1197   PC_MG   *mg      = (PC_MG *)pc->data;
1198   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1199 
1200   PetscFunctionBegin;
1201   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);
1202   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);
1203   pc_gamg->emax = emax;
1204   pc_gamg->emin = emin;
1205   PetscFunctionReturn(PETSC_SUCCESS);
1206 }
1207 
1208 /*@
1209   PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding a `PCGAMG` algebraic multigrid preconditioner
1210 
1211   Collective
1212 
1213   Input Parameters:
1214 + pc - the preconditioner context
1215 - n  - `PETSC_TRUE` or `PETSC_FALSE`
1216 
1217   Options Database Key:
1218 . -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation
1219 
1220   Level: intermediate
1221 
1222   Note:
1223   May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1224   rebuilding the preconditioner quicker.
1225 
1226 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`
1227 @*/
1228 PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1229 {
1230   PetscFunctionBegin;
1231   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1232   PetscTryMethod(pc, "PCGAMGSetReuseInterpolation_C", (PC, PetscBool), (pc, n));
1233   PetscFunctionReturn(PETSC_SUCCESS);
1234 }
1235 
1236 static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1237 {
1238   PC_MG   *mg      = (PC_MG *)pc->data;
1239   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1240 
1241   PetscFunctionBegin;
1242   pc_gamg->reuse_prol = n;
1243   PetscFunctionReturn(PETSC_SUCCESS);
1244 }
1245 
1246 /*@
1247   PCGAMGASMSetUseAggs - Have the `PCGAMG` smoother on each level use `PCASM` where the aggregates defined by the coarsening process are
1248   the subdomains for the additive Schwarz preconditioner used as the smoother
1249 
1250   Collective
1251 
1252   Input Parameters:
1253 + pc  - the preconditioner context
1254 - flg - `PETSC_TRUE` to use aggregates, `PETSC_FALSE` to not
1255 
1256   Options Database Key:
1257 . -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains
1258 
1259   Level: intermediate
1260 
1261   Note:
1262   This option automatically sets the preconditioner on the levels to be `PCASM`.
1263 
1264 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCASM`, `PCSetType`
1265 @*/
1266 PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1267 {
1268   PetscFunctionBegin;
1269   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1270   PetscTryMethod(pc, "PCGAMGASMSetUseAggs_C", (PC, PetscBool), (pc, flg));
1271   PetscFunctionReturn(PETSC_SUCCESS);
1272 }
1273 
1274 static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1275 {
1276   PC_MG   *mg      = (PC_MG *)pc->data;
1277   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1278 
1279   PetscFunctionBegin;
1280   pc_gamg->use_aggs_in_asm = flg;
1281   PetscFunctionReturn(PETSC_SUCCESS);
1282 }
1283 
1284 /*@
1285   PCGAMGSetParallelCoarseGridSolve - allow a parallel coarse grid solver
1286 
1287   Collective
1288 
1289   Input Parameters:
1290 + pc  - the preconditioner context
1291 - flg - `PETSC_TRUE` to not force coarse grid onto one processor
1292 
1293   Options Database Key:
1294 . -pc_gamg_parallel_coarse_grid_solver - use a parallel coarse grid solver
1295 
1296   Level: intermediate
1297 
1298 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGSetRankReductionFactors()`
1299 @*/
1300 PetscErrorCode PCGAMGSetParallelCoarseGridSolve(PC pc, PetscBool flg)
1301 {
1302   PetscFunctionBegin;
1303   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1304   PetscTryMethod(pc, "PCGAMGSetParallelCoarseGridSolve_C", (PC, PetscBool), (pc, flg));
1305   PetscFunctionReturn(PETSC_SUCCESS);
1306 }
1307 
1308 static PetscErrorCode PCGAMGSetParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1309 {
1310   PC_MG   *mg      = (PC_MG *)pc->data;
1311   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1312 
1313   PetscFunctionBegin;
1314   pc_gamg->use_parallel_coarse_grid_solver = flg;
1315   PetscFunctionReturn(PETSC_SUCCESS);
1316 }
1317 
1318 /*@
1319   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
1320 
1321   Collective
1322 
1323   Input Parameters:
1324 + pc  - the preconditioner context
1325 - flg - `PETSC_TRUE` to pin coarse grids to the CPU
1326 
1327   Options Database Key:
1328 . -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU
1329 
1330   Level: intermediate
1331 
1332 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetCoarseGridLayoutType()`, `PCGAMGSetParallelCoarseGridSolve()`
1333 @*/
1334 PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1335 {
1336   PetscFunctionBegin;
1337   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1338   PetscTryMethod(pc, "PCGAMGSetCpuPinCoarseGrids_C", (PC, PetscBool), (pc, flg));
1339   PetscFunctionReturn(PETSC_SUCCESS);
1340 }
1341 
1342 static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1343 {
1344   PC_MG   *mg      = (PC_MG *)pc->data;
1345   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1346 
1347   PetscFunctionBegin;
1348   pc_gamg->cpu_pin_coarse_grids = flg;
1349   PetscFunctionReturn(PETSC_SUCCESS);
1350 }
1351 
1352 /*@
1353   PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type)
1354 
1355   Collective
1356 
1357   Input Parameters:
1358 + pc  - the preconditioner context
1359 - flg - `PCGAMGLayoutType` type, either `PCGAMG_LAYOUT_COMPACT` or `PCGAMG_LAYOUT_SPREAD`
1360 
1361   Options Database Key:
1362 . -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering
1363 
1364   Level: intermediate
1365 
1366 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetCpuPinCoarseGrids()`, `PCGAMGLayoutType`, `PCGAMG_LAYOUT_COMPACT`, `PCGAMG_LAYOUT_SPREAD`
1367 @*/
1368 PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1369 {
1370   PetscFunctionBegin;
1371   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1372   PetscTryMethod(pc, "PCGAMGSetCoarseGridLayoutType_C", (PC, PCGAMGLayoutType), (pc, flg));
1373   PetscFunctionReturn(PETSC_SUCCESS);
1374 }
1375 
1376 static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1377 {
1378   PC_MG   *mg      = (PC_MG *)pc->data;
1379   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1380 
1381   PetscFunctionBegin;
1382   pc_gamg->layout_type = flg;
1383   PetscFunctionReturn(PETSC_SUCCESS);
1384 }
1385 
1386 /*@
1387   PCGAMGSetNlevels -  Sets the maximum number of levels `PCGAMG` will use
1388 
1389   Collective
1390 
1391   Input Parameters:
1392 + pc - the preconditioner
1393 - n  - the maximum number of levels to use
1394 
1395   Options Database Key:
1396 . -pc_mg_levels <n> - set the maximum number of levels to allow
1397 
1398   Level: intermediate
1399 
1400   Developer Notes:
1401   Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG`
1402 
1403 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`
1404 @*/
1405 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1406 {
1407   PetscFunctionBegin;
1408   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1409   PetscTryMethod(pc, "PCGAMGSetNlevels_C", (PC, PetscInt), (pc, n));
1410   PetscFunctionReturn(PETSC_SUCCESS);
1411 }
1412 
1413 static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1414 {
1415   PC_MG   *mg      = (PC_MG *)pc->data;
1416   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1417 
1418   PetscFunctionBegin;
1419   pc_gamg->Nlevels = n;
1420   PetscFunctionReturn(PETSC_SUCCESS);
1421 }
1422 
1423 /*@
1424   PCGAMGASMSetHEM -  Sets the number of HEM matching passed
1425 
1426   Collective
1427 
1428   Input Parameters:
1429 + pc - the preconditioner
1430 - n  - number of HEM matching passed to construct ASM subdomains
1431 
1432   Options Database Key:
1433 . -pc_gamg_asm_hem <n> - set the number of HEM matching passed
1434 
1435   Level: intermediate
1436 
1437   Developer Notes:
1438   Should be called `PCGAMGSetMaximumNumberlevels()` and possible be shared with `PCMG`
1439 
1440 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`
1441 @*/
1442 PetscErrorCode PCGAMGASMSetHEM(PC pc, PetscInt n)
1443 {
1444   PetscFunctionBegin;
1445   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1446   PetscTryMethod(pc, "PCGAMGASMSetHEM_C", (PC, PetscInt), (pc, n));
1447   PetscFunctionReturn(PETSC_SUCCESS);
1448 }
1449 
1450 static PetscErrorCode PCGAMGASMSetHEM_GAMG(PC pc, PetscInt n)
1451 {
1452   PC_MG   *mg      = (PC_MG *)pc->data;
1453   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1454 
1455   PetscFunctionBegin;
1456   pc_gamg->asm_hem_aggs = n;
1457   PetscFunctionReturn(PETSC_SUCCESS);
1458 }
1459 
1460 /*@
1461   PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
1462 
1463   Not Collective
1464 
1465   Input Parameters:
1466 + pc - the preconditioner context
1467 . 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
1468 - n  - number of threshold values provided in array
1469 
1470   Options Database Key:
1471 . -pc_gamg_threshold <threshold> - the threshold to drop edges
1472 
1473   Level: intermediate
1474 
1475   Notes:
1476   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.
1477   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.
1478 
1479   If `n` is less than the total number of coarsenings (see `PCGAMGSetNlevels()`), then threshold scaling (see `PCGAMGSetThresholdScale()`) is used for each successive coarsening.
1480   In this case, `PCGAMGSetThresholdScale()` must be called before `PCGAMGSetThreshold()`.
1481   If `n` is greater than the total number of levels, the excess entries in threshold will not be used.
1482 
1483 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetMinDegreeOrderingMISk()`, `PCGAMGSetThresholdScale()`
1484 @*/
1485 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1486 {
1487   PetscFunctionBegin;
1488   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1489   if (n) PetscAssertPointer(v, 2);
1490   PetscTryMethod(pc, "PCGAMGSetThreshold_C", (PC, PetscReal[], PetscInt), (pc, v, n));
1491   PetscFunctionReturn(PETSC_SUCCESS);
1492 }
1493 
1494 static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1495 {
1496   PC_MG   *mg      = (PC_MG *)pc->data;
1497   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1498   PetscInt i;
1499 
1500   PetscFunctionBegin;
1501   for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1502   for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale;
1503   PetscFunctionReturn(PETSC_SUCCESS);
1504 }
1505 
1506 /*@
1507   PCGAMGSetRankReductionFactors - Set a manual schedule for MPI rank reduction on coarse grids
1508 
1509   Collective
1510 
1511   Input Parameters:
1512 + pc - the preconditioner context
1513 . v  - array of reduction factors. 0 for first value forces a reduction to one process/device on first level in CUDA
1514 - n  - number of values provided in array
1515 
1516   Options Database Key:
1517 . -pc_gamg_rank_reduction_factors <factors> - provide the schedule
1518 
1519   Level: intermediate
1520 
1521 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetParallelCoarseGridSolve()`
1522 @*/
1523 PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1524 {
1525   PetscFunctionBegin;
1526   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1527   if (n) PetscAssertPointer(v, 2);
1528   PetscTryMethod(pc, "PCGAMGSetRankReductionFactors_C", (PC, PetscInt[], PetscInt), (pc, v, n));
1529   PetscFunctionReturn(PETSC_SUCCESS);
1530 }
1531 
1532 static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1533 {
1534   PC_MG   *mg      = (PC_MG *)pc->data;
1535   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1536   PetscInt i;
1537 
1538   PetscFunctionBegin;
1539   for (i = 0; i < PetscMin(n, PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1540   for (; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
1541   PetscFunctionReturn(PETSC_SUCCESS);
1542 }
1543 
1544 /*@
1545   PCGAMGSetThresholdScale - Relative threshold reduction at each level
1546 
1547   Not Collective
1548 
1549   Input Parameters:
1550 + pc - the preconditioner context
1551 - v  - the threshold value reduction, usually < 1.0
1552 
1553   Options Database Key:
1554 . -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level
1555 
1556   Level: advanced
1557 
1558   Note:
1559   The initial threshold (for an arbitrary number of levels starting from the finest) can be set with `PCGAMGSetThreshold()`.
1560   This scaling is used for each subsequent coarsening, but must be called before `PCGAMGSetThreshold()`.
1561 
1562 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`
1563 @*/
1564 PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1565 {
1566   PetscFunctionBegin;
1567   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1568   PetscTryMethod(pc, "PCGAMGSetThresholdScale_C", (PC, PetscReal), (pc, v));
1569   PetscFunctionReturn(PETSC_SUCCESS);
1570 }
1571 
1572 static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1573 {
1574   PC_MG   *mg      = (PC_MG *)pc->data;
1575   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1576 
1577   PetscFunctionBegin;
1578   pc_gamg->threshold_scale = v;
1579   PetscFunctionReturn(PETSC_SUCCESS);
1580 }
1581 
1582 /*@
1583   PCGAMGSetType - Set the type of algorithm `PCGAMG` should use
1584 
1585   Collective
1586 
1587   Input Parameters:
1588 + pc   - the preconditioner context
1589 - type - `PCGAMGAGG`, `PCGAMGGEO`, or `PCGAMGCLASSICAL`
1590 
1591   Options Database Key:
1592 . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply - only agg is supported
1593 
1594   Level: intermediate
1595 
1596 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMGGetType()`, `PCGAMG`, `PCGAMGType`
1597 @*/
1598 PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1599 {
1600   PetscFunctionBegin;
1601   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1602   PetscTryMethod(pc, "PCGAMGSetType_C", (PC, PCGAMGType), (pc, type));
1603   PetscFunctionReturn(PETSC_SUCCESS);
1604 }
1605 
1606 /*@
1607   PCGAMGGetType - Get the type of algorithm `PCGAMG` will use
1608 
1609   Collective
1610 
1611   Input Parameter:
1612 . pc - the preconditioner context
1613 
1614   Output Parameter:
1615 . type - the type of algorithm used
1616 
1617   Level: intermediate
1618 
1619 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetType()`, `PCGAMGType`
1620 @*/
1621 PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1622 {
1623   PetscFunctionBegin;
1624   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1625   PetscUseMethod(pc, "PCGAMGGetType_C", (PC, PCGAMGType *), (pc, type));
1626   PetscFunctionReturn(PETSC_SUCCESS);
1627 }
1628 
1629 static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1630 {
1631   PC_MG   *mg      = (PC_MG *)pc->data;
1632   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1633 
1634   PetscFunctionBegin;
1635   *type = pc_gamg->type;
1636   PetscFunctionReturn(PETSC_SUCCESS);
1637 }
1638 
1639 static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1640 {
1641   PC_MG   *mg      = (PC_MG *)pc->data;
1642   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1643   PetscErrorCode (*r)(PC);
1644 
1645   PetscFunctionBegin;
1646   pc_gamg->type = type;
1647   PetscCall(PetscFunctionListFind(GAMGList, type, &r));
1648   PetscCheck(r, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown GAMG type %s given", type);
1649   if (pc_gamg->ops->destroy) {
1650     PetscCall((*pc_gamg->ops->destroy)(pc));
1651     PetscCall(PetscMemzero(pc_gamg->ops, sizeof(struct _PCGAMGOps)));
1652     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1653     /* cleaning up common data in pc_gamg - this should disappear someday */
1654     pc_gamg->data_cell_cols      = 0;
1655     pc_gamg->data_cell_rows      = 0;
1656     pc_gamg->orig_data_cell_cols = 0;
1657     pc_gamg->orig_data_cell_rows = 0;
1658     PetscCall(PetscFree(pc_gamg->data));
1659     pc_gamg->data_sz = 0;
1660   }
1661   PetscCall(PetscFree(pc_gamg->gamg_type_name));
1662   PetscCall(PetscStrallocpy(type, &pc_gamg->gamg_type_name));
1663   PetscCall((*r)(pc));
1664   PetscFunctionReturn(PETSC_SUCCESS);
1665 }
1666 
1667 static PetscErrorCode PCView_GAMG(PC pc, PetscViewer viewer)
1668 {
1669   PC_MG    *mg      = (PC_MG *)pc->data;
1670   PC_GAMG  *pc_gamg = (PC_GAMG *)mg->innerctx;
1671   PetscReal gc = 0, oc = 0;
1672 
1673   PetscFunctionBegin;
1674   PetscCall(PetscViewerASCIIPrintf(viewer, "    GAMG specific options\n"));
1675   PetscCall(PetscViewerASCIIPrintf(viewer, "      Threshold for dropping small values in graph on each level ="));
1676   for (PetscInt i = 0; i < mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)pc_gamg->threshold[i]));
1677   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1678   PetscCall(PetscViewerASCIIPrintf(viewer, "      Threshold scaling factor for each level not specified = %g\n", (double)pc_gamg->threshold_scale));
1679   if (pc_gamg->use_aggs_in_asm) PetscCall(PetscViewerASCIIPrintf(viewer, "      Using aggregates from GAMG coarsening process to define subdomains for PCASM\n")); // this take presedence
1680   else if (pc_gamg->asm_hem_aggs) {
1681     PetscCall(PetscViewerASCIIPrintf(viewer, "      Using aggregates made with %d applications of heavy edge matching (HEM) to define subdomains for PCASM\n", (int)pc_gamg->asm_hem_aggs));
1682     PetscCall(PetscViewerASCIIPushTab(viewer));
1683     PetscCall(PetscViewerASCIIPushTab(viewer));
1684     PetscCall(PetscViewerASCIIPushTab(viewer));
1685     PetscCall(PetscViewerASCIIPushTab(viewer));
1686     PetscCall(MatCoarsenView(pc_gamg->asm_crs, viewer));
1687     PetscCall(PetscViewerASCIIPopTab(viewer));
1688     PetscCall(PetscViewerASCIIPopTab(viewer));
1689     PetscCall(PetscViewerASCIIPopTab(viewer));
1690   }
1691   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"));
1692   if (pc_gamg->injection_index_size) {
1693     PetscCall(PetscViewerASCIIPrintf(viewer, "      Using injection restriction/prolongation on first level, dofs:"));
1694     for (int i = 0; i < pc_gamg->injection_index_size; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %d", (int)pc_gamg->injection_index[i]));
1695     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1696   }
1697   if (pc_gamg->ops->view) PetscCall((*pc_gamg->ops->view)(pc, viewer));
1698   PetscCall(PCMGGetGridComplexity(pc, &gc, &oc));
1699   PetscCall(PetscViewerASCIIPrintf(viewer, "      Complexity:    grid = %g    operator = %g\n", (double)gc, (double)oc));
1700   PetscFunctionReturn(PETSC_SUCCESS);
1701 }
1702 
1703 /*@
1704   PCGAMGSetInjectionIndex - Array of subset of variables per vertex to inject into coarse grid space
1705 
1706   Logically Collective
1707 
1708   Input Parameters:
1709 + pc  - the coarsen context
1710 . n   - number of indices
1711 - idx - array of indices
1712 
1713   Options Database Key:
1714 . -pc_gamg_injection_index - array of subset of variables per vertex to use for injection coarse grid space
1715 
1716   Level: intermediate
1717 
1718 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`
1719 @*/
1720 PetscErrorCode PCGAMGSetInjectionIndex(PC pc, PetscInt n, PetscInt idx[])
1721 {
1722   PetscFunctionBegin;
1723   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1724   PetscValidLogicalCollectiveInt(pc, n, 2);
1725   PetscTryMethod(pc, "PCGAMGSetInjectionIndex_C", (PC, PetscInt, PetscInt[]), (pc, n, idx));
1726   PetscFunctionReturn(PETSC_SUCCESS);
1727 }
1728 
1729 static PetscErrorCode PCGAMGSetInjectionIndex_GAMG(PC pc, PetscInt n, PetscInt idx[])
1730 {
1731   PC_MG   *mg      = (PC_MG *)pc->data;
1732   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1733 
1734   PetscFunctionBegin;
1735   pc_gamg->injection_index_size = n;
1736   PetscCheck(n < MAT_COARSEN_STRENGTH_INDEX_SIZE, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "array size %d larger than max %d", (int)n, MAT_COARSEN_STRENGTH_INDEX_SIZE);
1737   for (PetscInt i = 0; i < n; i++) pc_gamg->injection_index[i] = idx[i];
1738   PetscFunctionReturn(PETSC_SUCCESS);
1739 }
1740 
1741 static PetscErrorCode PCSetFromOptions_GAMG(PC pc, PetscOptionItems *PetscOptionsObject)
1742 {
1743   PC_MG             *mg      = (PC_MG *)pc->data;
1744   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
1745   PetscBool          flag;
1746   MPI_Comm           comm;
1747   char               prefix[256], tname[32];
1748   PetscInt           i, n;
1749   const char        *pcpre;
1750   static const char *LayoutTypes[] = {"compact", "spread", "PCGAMGLayoutType", "PC_GAMG_LAYOUT", NULL};
1751 
1752   PetscFunctionBegin;
1753   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1754   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG options");
1755   PetscCall(PetscOptionsFList("-pc_gamg_type", "Type of AMG method (only 'agg' supported and useful)", "PCGAMGSetType", GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag));
1756   if (flag) PetscCall(PCGAMGSetType(pc, tname));
1757   PetscCall(PetscOptionsBool("-pc_gamg_repartition", "Repartion coarse grids", "PCGAMGSetRepartition", pc_gamg->repart, &pc_gamg->repart, NULL));
1758   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));
1759   PetscCall(PetscOptionsBool("-pc_gamg_recompute_esteig", "Set flag to recompute eigen estimates for Chebyshev when matrix changes", "PCGAMGSetRecomputeEstEig", pc_gamg->recompute_esteig, &pc_gamg->recompute_esteig, NULL));
1760   PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation", "Reuse prolongation operator", "PCGAMGReuseInterpolation", pc_gamg->reuse_prol, &pc_gamg->reuse_prol, NULL));
1761   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));
1762   PetscCall(PetscOptionsBool("-pc_gamg_parallel_coarse_grid_solver", "Use parallel coarse grid solver (otherwise put last grid on one process)", "PCGAMGSetParallelCoarseGridSolve", pc_gamg->use_parallel_coarse_grid_solver, &pc_gamg->use_parallel_coarse_grid_solver, NULL));
1763   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));
1764   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,
1765                              (PetscEnum)pc_gamg->layout_type, (PetscEnum *)&pc_gamg->layout_type, NULL));
1766   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));
1767   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));
1768   PetscCall(PetscOptionsInt("-pc_gamg_asm_hem_aggs", "Number of HEM matching passed in aggregates for ASM smoother", "PCGAMGASMSetHEM", pc_gamg->asm_hem_aggs, &pc_gamg->asm_hem_aggs, NULL));
1769   PetscCall(PetscOptionsReal("-pc_gamg_threshold_scale", "Scaling of threshold for each level not specified", "PCGAMGSetThresholdScale", pc_gamg->threshold_scale, &pc_gamg->threshold_scale, NULL));
1770   n = PETSC_MG_MAXLEVELS;
1771   PetscCall(PetscOptionsRealArray("-pc_gamg_threshold", "Relative threshold to use for dropping edges in aggregation graph", "PCGAMGSetThreshold", pc_gamg->threshold, &n, &flag));
1772   if (!flag || n < PETSC_MG_MAXLEVELS) {
1773     if (!flag) n = 1;
1774     i = n;
1775     do {
1776       pc_gamg->threshold[i] = pc_gamg->threshold[i - 1] * pc_gamg->threshold_scale;
1777     } while (++i < PETSC_MG_MAXLEVELS);
1778   }
1779   PetscCall(PetscOptionsInt("-pc_mg_levels", "Set number of MG levels (should get from base class)", "PCGAMGSetNlevels", pc_gamg->Nlevels, &pc_gamg->Nlevels, NULL));
1780   PetscCheck(pc_gamg->Nlevels <= PETSC_MG_MAXLEVELS, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_mg_levels (%d) >= PETSC_MG_MAXLEVELS (%d)", (int)pc_gamg->Nlevels, (int)PETSC_MG_MAXLEVELS);
1781   n = PETSC_MG_MAXLEVELS;
1782   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));
1783   if (!flag) i = 0;
1784   else i = n;
1785   do {
1786     pc_gamg->level_reduction_factors[i] = -1;
1787   } while (++i < PETSC_MG_MAXLEVELS);
1788   {
1789     PetscReal eminmax[2] = {0., 0.};
1790     n                    = 2;
1791     PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues", "extreme eigenvalues for smoothed aggregation", "PCGAMGSetEigenvalues", eminmax, &n, &flag));
1792     if (flag) {
1793       PetscCheck(n == 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
1794       PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]));
1795     }
1796   }
1797   pc_gamg->injection_index_size = MAT_COARSEN_STRENGTH_INDEX_SIZE;
1798   PetscCall(PetscOptionsIntArray("-pc_gamg_injection_index", "Array of indices to use to use injection coarse grid space", "PCGAMGSetInjectionIndex", pc_gamg->injection_index, &pc_gamg->injection_index_size, NULL));
1799   /* set options for subtype */
1800   PetscCall((*pc_gamg->ops->setfromoptions)(pc, PetscOptionsObject));
1801 
1802   PetscCall(PCGetOptionsPrefix(pc, &pcpre));
1803   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%spc_gamg_", pcpre ? pcpre : ""));
1804   PetscOptionsHeadEnd();
1805   PetscFunctionReturn(PETSC_SUCCESS);
1806 }
1807 
1808 /*MC
1809   PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
1810 
1811   Options Database Keys:
1812 + -pc_gamg_type <type,default=agg> - one of agg, geo, or classical (only smoothed aggregation, agg, supported)
1813 . -pc_gamg_repartition  <bool,default=false> - repartition the degrees of freedom across the coarse grids as they are determined
1814 . -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
1815 . -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>
1816                                         equations on each process that has degrees of freedom
1817 . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
1818 . -pc_gamg_reuse_interpolation <bool,default=true> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations (should always be true)
1819 . -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)
1820 - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1821 
1822   Options Database Keys for Aggregation:
1823 + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation to construct prolongation
1824 . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest.
1825 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening
1826 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm
1827 . -pc_gamg_pc_gamg_asm_hem_aggs <n,default=0> - Number of HEM aggregation steps for ASM smoother
1828 - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive')
1829 
1830   Options Database Keys for Multigrid:
1831 + -pc_mg_cycle_type <v> - v or w, see `PCMGSetCycleType()`
1832 . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1833 . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1834 - -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
1835 
1836   Level: intermediate
1837 
1838   Notes:
1839   To obtain good performance for `PCGAMG` for vector valued problems you must
1840   call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point
1841   call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator
1842 
1843   The many options for `PCMG` also work directly for `PCGAMG` such as controlling the smoothers on each level etc.
1844 
1845 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCCreate()`, `PCSetType()`,
1846           `MatSetBlockSize()`,
1847           `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`,
1848           `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`,
1849           `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()`
1850 M*/
1851 PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1852 {
1853   PC_GAMG *pc_gamg;
1854   PC_MG   *mg;
1855 
1856   PetscFunctionBegin;
1857   /* register AMG type */
1858   PetscCall(PCGAMGInitializePackage());
1859 
1860   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1861   PetscCall(PCSetType(pc, PCMG));
1862   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG));
1863 
1864   /* create a supporting struct and attach it to pc */
1865   PetscCall(PetscNew(&pc_gamg));
1866   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
1867   mg           = (PC_MG *)pc->data;
1868   mg->innerctx = pc_gamg;
1869 
1870   PetscCall(PetscNew(&pc_gamg->ops));
1871 
1872   /* these should be in subctx but repartitioning needs simple arrays */
1873   pc_gamg->data_sz = 0;
1874   pc_gamg->data    = NULL;
1875 
1876   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1877   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1878   pc->ops->setup          = PCSetUp_GAMG;
1879   pc->ops->reset          = PCReset_GAMG;
1880   pc->ops->destroy        = PCDestroy_GAMG;
1881   mg->view                = PCView_GAMG;
1882 
1883   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetProcEqLim_C", PCGAMGSetProcEqLim_GAMG));
1884   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseEqLim_C", PCGAMGSetCoarseEqLim_GAMG));
1885   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRepartition_C", PCGAMGSetRepartition_GAMG));
1886   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetEigenvalues_C", PCGAMGSetEigenvalues_GAMG));
1887   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetUseSAEstEig_C", PCGAMGSetUseSAEstEig_GAMG));
1888   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRecomputeEstEig_C", PCGAMGSetRecomputeEstEig_GAMG));
1889   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetReuseInterpolation_C", PCGAMGSetReuseInterpolation_GAMG));
1890   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetUseAggs_C", PCGAMGASMSetUseAggs_GAMG));
1891   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetParallelCoarseGridSolve_C", PCGAMGSetParallelCoarseGridSolve_GAMG));
1892   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCpuPinCoarseGrids_C", PCGAMGSetCpuPinCoarseGrids_GAMG));
1893   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetCoarseGridLayoutType_C", PCGAMGSetCoarseGridLayoutType_GAMG));
1894   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThreshold_C", PCGAMGSetThreshold_GAMG));
1895   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetRankReductionFactors_C", PCGAMGSetRankReductionFactors_GAMG));
1896   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetThresholdScale_C", PCGAMGSetThresholdScale_GAMG));
1897   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetType_C", PCGAMGSetType_GAMG));
1898   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGGetType_C", PCGAMGGetType_GAMG));
1899   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNlevels_C", PCGAMGSetNlevels_GAMG));
1900   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGASMSetHEM_C", PCGAMGASMSetHEM_GAMG));
1901   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetInjectionIndex_C", PCGAMGSetInjectionIndex_GAMG));
1902   pc_gamg->repart                          = PETSC_FALSE;
1903   pc_gamg->reuse_prol                      = PETSC_TRUE;
1904   pc_gamg->use_aggs_in_asm                 = PETSC_FALSE;
1905   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1906   pc_gamg->cpu_pin_coarse_grids            = PETSC_FALSE;
1907   pc_gamg->layout_type                     = PCGAMG_LAYOUT_SPREAD;
1908   pc_gamg->min_eq_proc                     = 50;
1909   pc_gamg->asm_hem_aggs                    = 0;
1910   pc_gamg->coarse_eq_limit                 = 50;
1911   for (int i = 0; i < PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = -1;
1912   pc_gamg->threshold_scale  = 1.;
1913   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
1914   pc_gamg->current_level    = 0; /* don't need to init really */
1915   pc_gamg->use_sa_esteig    = PETSC_TRUE;
1916   pc_gamg->recompute_esteig = PETSC_TRUE;
1917   pc_gamg->emin             = 0;
1918   pc_gamg->emax             = 0;
1919 
1920   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1921 
1922   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1923   PetscCall(PCGAMGSetType(pc, PCGAMGAGG));
1924   PetscFunctionReturn(PETSC_SUCCESS);
1925 }
1926 
1927 /*@C
1928   PCGAMGInitializePackage - This function initializes everything in the `PCGAMG` package. It is called
1929   from `PCInitializePackage()`.
1930 
1931   Level: developer
1932 
1933 .seealso: [](ch_ksp), `PetscInitialize()`
1934 @*/
1935 PetscErrorCode PCGAMGInitializePackage(void)
1936 {
1937   PetscInt l;
1938 
1939   PetscFunctionBegin;
1940   if (PCGAMGPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
1941   PCGAMGPackageInitialized = PETSC_TRUE;
1942   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGGEO, PCCreateGAMG_GEO));
1943   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGAGG, PCCreateGAMG_AGG));
1944   PetscCall(PetscFunctionListAdd(&GAMGList, PCGAMGCLASSICAL, PCCreateGAMG_Classical));
1945   PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage));
1946 
1947   /* general events */
1948   PetscCall(PetscLogEventRegister("PCSetUp_GAMG+", PC_CLASSID, &petsc_gamg_setup_events[GAMG_SETUP]));
1949   PetscCall(PetscLogEventRegister(" PCGAMGCreateG", PC_CLASSID, &petsc_gamg_setup_events[GAMG_GRAPH]));
1950   PetscCall(PetscLogEventRegister(" GAMG Coarsen", PC_CLASSID, &petsc_gamg_setup_events[GAMG_COARSEN]));
1951   PetscCall(PetscLogEventRegister("  GAMG MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[GAMG_MIS]));
1952   PetscCall(PetscLogEventRegister(" PCGAMGProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROL]));
1953   PetscCall(PetscLogEventRegister("  GAMG Prol-col", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLA]));
1954   PetscCall(PetscLogEventRegister("  GAMG Prol-lift", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PROLB]));
1955   PetscCall(PetscLogEventRegister(" PCGAMGOptProl", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPT]));
1956   PetscCall(PetscLogEventRegister("  GAMG smooth", PC_CLASSID, &petsc_gamg_setup_events[GAMG_OPTSM]));
1957   PetscCall(PetscLogEventRegister(" PCGAMGCreateL", PC_CLASSID, &petsc_gamg_setup_events[GAMG_LEVEL]));
1958   PetscCall(PetscLogEventRegister("  GAMG PtAP", PC_CLASSID, &petsc_gamg_setup_events[GAMG_PTAP]));
1959   PetscCall(PetscLogEventRegister("  GAMG Reduce", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REDUCE]));
1960   PetscCall(PetscLogEventRegister("   GAMG Repart", PC_CLASSID, &petsc_gamg_setup_events[GAMG_REPART]));
1961   /* PetscCall(PetscLogEventRegister("   GAMG Inv-Srt", PC_CLASSID, &petsc_gamg_setup_events[SET13])); */
1962   /* PetscCall(PetscLogEventRegister("   GAMG Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14])); */
1963   /* PetscCall(PetscLogEventRegister("   GAMG Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15])); */
1964   for (l = 0; l < PETSC_MG_MAXLEVELS; l++) {
1965     char ename[32];
1966 
1967     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Squ l%02" PetscInt_FMT, l));
1968     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]));
1969     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Gal l%02" PetscInt_FMT, l));
1970     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]));
1971     PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCGAMG Opt l%02" PetscInt_FMT, l));
1972     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]));
1973   }
1974 #if defined(GAMG_STAGES)
1975   { /* create timer stages */
1976     char str[32];
1977     PetscCall(PetscSNPrintf(str, PETSC_STATIC_ARRAY_LENGTH(str), "GAMG Level %d", 0));
1978     PetscCall(PetscLogStageRegister(str, &gamg_stages[0]));
1979   }
1980 #endif
1981   PetscFunctionReturn(PETSC_SUCCESS);
1982 }
1983 
1984 /*@C
1985   PCGAMGFinalizePackage - This function frees everything from the `PCGAMG` package. It is
1986   called from `PetscFinalize()` automatically.
1987 
1988   Level: developer
1989 
1990 .seealso: [](ch_ksp), `PetscFinalize()`
1991 @*/
1992 PetscErrorCode PCGAMGFinalizePackage(void)
1993 {
1994   PetscFunctionBegin;
1995   PCGAMGPackageInitialized = PETSC_FALSE;
1996   PetscCall(PetscFunctionListDestroy(&GAMGList));
1997   PetscFunctionReturn(PETSC_SUCCESS);
1998 }
1999 
2000 /*@C
2001   PCGAMGRegister - Register a `PCGAMG` implementation.
2002 
2003   Input Parameters:
2004 + type   - string that will be used as the name of the `PCGAMG` type.
2005 - create - function for creating the gamg context.
2006 
2007   Level: developer
2008 
2009 .seealso: [](ch_ksp), `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
2010 @*/
2011 PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
2012 {
2013   PetscFunctionBegin;
2014   PetscCall(PCGAMGInitializePackage());
2015   PetscCall(PetscFunctionListAdd(&GAMGList, type, create));
2016   PetscFunctionReturn(PETSC_SUCCESS);
2017 }
2018 
2019 /*@
2020   PCGAMGCreateGraph - Creates a graph that is used by the ``PCGAMGType`` in the coarsening process
2021 
2022   Input Parameters:
2023 + pc - the `PCGAMG`
2024 - A  - the matrix, for any level
2025 
2026   Output Parameter:
2027 . G - the graph
2028 
2029   Level: advanced
2030 
2031 .seealso: [](ch_ksp), `PCGAMGType`, `PCGAMG`, `PCGAMGSetType()`
2032 @*/
2033 PetscErrorCode PCGAMGCreateGraph(PC pc, Mat A, Mat *G)
2034 {
2035   PC_MG   *mg      = (PC_MG *)pc->data;
2036   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
2037 
2038   PetscFunctionBegin;
2039   PetscCall(pc_gamg->ops->creategraph(pc, A, G));
2040   PetscFunctionReturn(PETSC_SUCCESS);
2041 }
2042