xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision c3b5f7ba6bc5ce25a01a67bb37ba5d34b02bbbd7)
1 /*
2  GAMG geometric-algebric multigrid PC - Mark Adams 2011
3  */
4 #include <petsc/private/matimpl.h>
5 #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
6 #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h>    /*I "petscksp.h" I*/
7 
8 #if defined(PETSC_HAVE_CUDA)
9   #include <cuda_runtime.h>
10 #endif
11 
12 #if defined(PETSC_HAVE_HIP)
13   #include <hip/hip_runtime.h>
14 #endif
15 
16 PetscLogEvent petsc_gamg_setup_events[NUM_SET];
17 PetscLogEvent petsc_gamg_setup_matmat_events[PETSC_MG_MAXLEVELS][3];
18 PetscLogEvent PC_GAMGGraph_AGG;
19 PetscLogEvent PC_GAMGGraph_GEO;
20 PetscLogEvent PC_GAMGCoarsen_AGG;
21 PetscLogEvent PC_GAMGCoarsen_GEO;
22 PetscLogEvent PC_GAMGProlongator_AGG;
23 PetscLogEvent PC_GAMGProlongator_GEO;
24 PetscLogEvent PC_GAMGOptProlongator_AGG;
25 
26 /* #define GAMG_STAGES */
27 #if defined(GAMG_STAGES)
28 static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS];
29 #endif
30 
31 static PetscFunctionList GAMGList = NULL;
32 static PetscBool PCGAMGPackageInitialized;
33 
34 /* ----------------------------------------------------------------------------- */
35 PetscErrorCode PCReset_GAMG(PC pc)
36 {
37   PC_MG          *mg      = (PC_MG*)pc->data;
38   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
39 
40   PetscFunctionBegin;
41   PetscCall(PetscFree(pc_gamg->data));
42   pc_gamg->data_sz = 0;
43   PetscCall(PetscFree(pc_gamg->orig_data));
44   for (PetscInt level = 0; level < PETSC_MG_MAXLEVELS ; level++) {
45     mg->min_eigen_DinvA[level] = 0;
46     mg->max_eigen_DinvA[level] = 0;
47   }
48   pc_gamg->emin = 0;
49   pc_gamg->emax = 0;
50   PetscFunctionReturn(0);
51 }
52 
53 /* -------------------------------------------------------------------------- */
54 /*
55    PCGAMGCreateLevel_GAMG: create coarse op with RAP.  repartition and/or reduce number
56      of active processors.
57 
58    Input Parameter:
59    . pc - parameters + side effect: coarse data in 'pc_gamg->data' and
60           'pc_gamg->data_sz' are changed via repartitioning/reduction.
61    . Amat_fine - matrix on this fine (k) level
62    . cr_bs - coarse block size
63    In/Output Parameter:
64    . a_P_inout - prolongation operator to the next level (k-->k-1)
65    . a_nactive_proc - number of active procs
66    Output Parameter:
67    . a_Amat_crs - coarse matrix that is created (k-1)
68 */
69 
70 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)
71 {
72   PC_MG           *mg         = (PC_MG*)pc->data;
73   PC_GAMG         *pc_gamg    = (PC_GAMG*)mg->innerctx;
74   Mat             Cmat,Pold=*a_P_inout;
75   MPI_Comm        comm;
76   PetscMPIInt     rank,size,new_size,nactive=*a_nactive_proc;
77   PetscInt        ncrs_eq,ncrs,f_bs;
78 
79   PetscFunctionBegin;
80   PetscCall(PetscObjectGetComm((PetscObject)Amat_fine,&comm));
81   PetscCallMPI(MPI_Comm_rank(comm, &rank));
82   PetscCallMPI(MPI_Comm_size(comm, &size));
83   PetscCall(MatGetBlockSize(Amat_fine, &f_bs));
84   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0));
85   PetscCall(MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat));
86   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][1],0,0,0,0));
87 
88   if (Pcolumnperm) *Pcolumnperm = NULL;
89 
90   /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/
91   PetscCall(MatGetLocalSize(Cmat, &ncrs_eq, NULL));
92   if (pc_gamg->data_cell_rows>0) {
93     ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows;
94   } else {
95     PetscInt  bs;
96     PetscCall(MatGetBlockSize(Cmat, &bs));
97     ncrs = ncrs_eq/bs;
98   }
99   /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */
100   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. */
101 #if defined(PETSC_HAVE_CUDA)
102     PetscShmComm pshmcomm;
103     PetscMPIInt  locrank;
104     MPI_Comm     loccomm;
105     PetscInt     s_nnodes,r_nnodes, new_new_size;
106     cudaError_t  cerr;
107     int          devCount;
108     PetscCall(PetscShmCommGet(comm,&pshmcomm));
109     PetscCall(PetscShmCommGetMpiShmComm(pshmcomm,&loccomm));
110     PetscCallMPI(MPI_Comm_rank(loccomm, &locrank));
111     s_nnodes = !locrank;
112     PetscCallMPI(MPI_Allreduce(&s_nnodes,&r_nnodes,1,MPIU_INT,MPI_SUM,comm));
113     PetscCheck((size%r_nnodes) == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of nodes np=%D nnodes%D",size,r_nnodes);
114     devCount = 0;
115     cerr = cudaGetDeviceCount(&devCount);
116     cudaGetLastError(); /* Reset the last error */
117     if (cerr == cudaSuccess && devCount >= 1) { /* There are devices, else go to heuristic */
118       new_new_size = r_nnodes * devCount;
119       new_size = new_new_size;
120       PetscCall(PetscInfo(pc,"%s: Fine grid with Cuda. %D nodes. Change new active set size %d --> %d (devCount=%d #nodes=%D)\n",((PetscObject)pc)->prefix,r_nnodes,nactive,new_size,devCount,r_nnodes));
121     } else {
122       PetscCall(PetscInfo(pc,"%s: With Cuda but no device. Use heuristics."));
123       goto HEURISTIC;
124     }
125 #else
126     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not be here");
127 #endif
128   } else if (pc_gamg->level_reduction_factors[pc_gamg->current_level] > 0) {
129     PetscCheckFalse(nactive%pc_gamg->level_reduction_factors[pc_gamg->current_level],PETSC_COMM_SELF,PETSC_ERR_PLIB,"odd number of active process %D wrt reduction factor %D",nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]);
130     new_size = nactive/pc_gamg->level_reduction_factors[pc_gamg->current_level];
131     PetscCall(PetscInfo(pc,"%s: Manually setting reduction to %d active processes (%d/%D)\n",((PetscObject)pc)->prefix,new_size,nactive,pc_gamg->level_reduction_factors[pc_gamg->current_level]));
132   } else if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) {
133     new_size = 1;
134     PetscCall(PetscInfo(pc,"%s: Force coarsest grid reduction to %d active processes\n",((PetscObject)pc)->prefix,new_size));
135   } else {
136     PetscInt ncrs_eq_glob;
137 #if defined(PETSC_HAVE_CUDA)
138     HEURISTIC:
139 #endif
140     PetscCall(MatGetSize(Cmat, &ncrs_eq_glob, NULL));
141     new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */
142     if (!new_size) new_size = 1; /* not likely, posible? */
143     else if (new_size >= nactive) new_size = nactive; /* no change, rare */
144     PetscCall(PetscInfo(pc,"%s: Coarse grid reduction from %d to %d active processes\n",((PetscObject)pc)->prefix,nactive,new_size));
145   }
146   if (new_size==nactive) {
147     *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */
148     if (new_size < size) {
149       /* odd case where multiple coarse grids are on one processor or no coarsening ... */
150       PetscCall(PetscInfo(pc,"%s: reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",((PetscObject)pc)->prefix,nactive));
151       if (pc_gamg->cpu_pin_coarse_grids) {
152         PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE));
153         PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE));
154       }
155     }
156     /* we know that the grid structure can be reused in MatPtAP */
157   } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */
158     PetscInt       *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1;
159     IS             is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices;
160     nloc_old = ncrs_eq/cr_bs;
161     PetscCheckFalse(ncrs_eq % cr_bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs);
162     /* get new_size and rfactor */
163     if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) {
164       /* find factor */
165       if (new_size == 1) rfactor = size; /* don't modify */
166       else {
167         PetscReal best_fact = 0.;
168         jj = -1;
169         for (kk = 1 ; kk <= size ; kk++) {
170           if (!(size%kk)) { /* a candidate */
171             PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size;
172             if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */
173             if (fact > best_fact) {
174               best_fact = fact; 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)=%D\n",((PetscObject)pc)->prefix,new_size,ncrs_eq));
187         PetscFunctionReturn(0);
188       }
189     }
190     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0));
191     /* make 'is_eq_newproc' */
192     PetscCall(PetscMalloc1(size, &counts));
193     if (pc_gamg->repart) {
194       /* Repartition Cmat_{k} and move columns of P^{k}_{k-1} and coordinates of primal part accordingly */
195       Mat      adj;
196       PetscCall(PetscInfo(pc,"%s: Repartition: size (active): %d --> %d, %D local equations, using %s process layout\n",((PetscObject)pc)->prefix,*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread"));
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; char fname[32];
245           PetscCall(PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev));
246           PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer);
247           PetscCall(MatView(tMat, viewer));
248           PetscCall(PetscViewerDestroy(&viewer));
249         }
250         PetscCall(MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
251         PetscCall(MatDestroy(&tMat));
252       } /* create 'adj' */
253 
254       { /* partition: get newproc_idx */
255         char            prefix[256];
256         const char      *pcpre;
257         const PetscInt  *is_idx;
258         MatPartitioning mpart;
259         IS              proc_is;
260 
261         PetscCall(MatPartitioningCreate(comm, &mpart));
262         PetscCall(MatPartitioningSetAdjacency(mpart, adj));
263         PetscCall(PCGetOptionsPrefix(pc, &pcpre));
264         PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : ""));
265         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix));
266         PetscCall(MatPartitioningSetFromOptions(mpart));
267         PetscCall(MatPartitioningSetNParts(mpart, new_size));
268         PetscCall(MatPartitioningApply(mpart, &proc_is));
269         PetscCall(MatPartitioningDestroy(&mpart));
270 
271         /* collect IS info */
272         PetscCall(PetscMalloc1(ncrs_eq, &newproc_idx));
273         PetscCall(ISGetIndices(proc_is, &is_idx));
274         for (kk = jj = 0 ; kk < nloc_old ; kk++) {
275           for (ii = 0 ; ii < cr_bs ; ii++, jj++) {
276             newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */
277           }
278         }
279         PetscCall(ISRestoreIndices(proc_is, &is_idx));
280         PetscCall(ISDestroy(&proc_is));
281       }
282       PetscCall(MatDestroy(&adj));
283 
284       PetscCall(ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc));
285       PetscCall(PetscFree(newproc_idx));
286     } else { /* simple aggregation of parts -- 'is_eq_newproc' */
287       PetscInt targetPE;
288       PetscCheckFalse(new_size==nactive,PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen");
289       PetscCall(PetscInfo(pc,"%s: Number of equations (loc) %D 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     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0));
382     /*
383       Invert for MatCreateSubMatrix
384     */
385     PetscCall(ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices));
386     PetscCall(ISSort(new_eq_indices)); /* is this needed? */
387     PetscCall(ISSetBlockSize(new_eq_indices, cr_bs));
388     if (is_eq_num != is_eq_num_prim) {
389       PetscCall(ISDestroy(&is_eq_num_prim)); /* could be same as 'is_eq_num' */
390     }
391     if (Pcolumnperm) {
392       PetscCall(PetscObjectReference((PetscObject)new_eq_indices));
393       *Pcolumnperm = new_eq_indices;
394     }
395     PetscCall(ISDestroy(&is_eq_num));
396     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0));
397     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0));
398     /* 'a_Amat_crs' output */
399     {
400       Mat       mat;
401       PetscBool flg;
402       PetscCall(MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat));
403       PetscCall(MatGetOption(Cmat, MAT_SPD, &flg));
404       if (flg) {
405         PetscCall(MatSetOption(mat, MAT_SPD,PETSC_TRUE));
406       } else {
407         PetscCall(MatGetOption(Cmat, MAT_HERMITIAN, &flg));
408         if (flg) {
409           PetscCall(MatSetOption(mat, MAT_HERMITIAN,PETSC_TRUE));
410         } else {
411 #if !defined(PETSC_USE_COMPLEX)
412           PetscCall(MatGetOption(Cmat, MAT_SYMMETRIC, &flg));
413           if (flg) {
414             PetscCall(MatSetOption(mat, MAT_SYMMETRIC,PETSC_TRUE));
415           }
416 #endif
417         }
418       }
419       *a_Amat_crs = mat;
420     }
421     PetscCall(MatDestroy(&Cmat));
422 
423     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0));
424     /* prolongator */
425     {
426       IS       findices;
427       PetscInt Istart,Iend;
428       Mat      Pnew;
429 
430       PetscCall(MatGetOwnershipRange(Pold, &Istart, &Iend));
431       PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0));
432       PetscCall(ISCreateStride(comm,Iend-Istart,Istart,1,&findices));
433       PetscCall(ISSetBlockSize(findices,f_bs));
434       PetscCall(MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew));
435       PetscCall(ISDestroy(&findices));
436       PetscCall(MatSetOption(Pnew,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE));
437 
438       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0));
439       PetscCall(MatDestroy(a_P_inout));
440 
441       /* output - repartitioned */
442       *a_P_inout = Pnew;
443     }
444     PetscCall(ISDestroy(&new_eq_indices));
445 
446     *a_nactive_proc = new_size; /* output */
447 
448     /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */
449     if (pc_gamg->cpu_pin_coarse_grids) {
450 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
451       static PetscInt llev = 2;
452       PetscCall(PetscInfo(pc,"%s: Pinning level %D to the CPU\n",((PetscObject)pc)->prefix,llev++));
453 #endif
454       PetscCall(MatBindToCPU(*a_Amat_crs,PETSC_TRUE));
455       PetscCall(MatBindToCPU(*a_P_inout,PETSC_TRUE));
456       if (1) { /* HACK: move this to MatBindCPU_MPIAIJXXX; lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */
457         Mat         A = *a_Amat_crs, P = *a_P_inout;
458         PetscMPIInt size;
459         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
460         if (size > 1) {
461           Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data;
462           PetscCall(VecBindToCPU(a->lvec,PETSC_TRUE));
463           PetscCall(VecBindToCPU(p->lvec,PETSC_TRUE));
464         }
465       }
466     }
467     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0));
468   }
469   PetscFunctionReturn(0);
470 }
471 
472 PetscErrorCode PCGAMGSquareGraph_GAMG(PC a_pc, Mat Gmat1, Mat* Gmat2)
473 {
474   const char     *prefix;
475   char           addp[32];
476   PC_MG          *mg      = (PC_MG*)a_pc->data;
477   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
478 
479   PetscFunctionBegin;
480   PetscCall(PCGetOptionsPrefix(a_pc,&prefix));
481   PetscCall(PetscInfo(a_pc,"%s: Square Graph on level %D\n",((PetscObject)a_pc)->prefix,pc_gamg->current_level+1));
482   PetscCall(MatProductCreate(Gmat1,Gmat1,NULL,Gmat2));
483   PetscCall(MatSetOptionsPrefix(*Gmat2,prefix));
484   PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_square_%d_",pc_gamg->current_level));
485   PetscCall(MatAppendOptionsPrefix(*Gmat2,addp));
486   if ((*Gmat2)->structurally_symmetric) {
487     PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AB));
488   } else {
489     PetscCall(MatSetOption(Gmat1,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE));
490     PetscCall(MatProductSetType(*Gmat2,MATPRODUCT_AtB));
491   }
492   PetscCall(MatProductSetFromOptions(*Gmat2));
493   PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0));
494   PetscCall(MatProductSymbolic(*Gmat2));
495   PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][0],0,0,0,0));
496   PetscCall(MatProductClear(*Gmat2));
497   /* we only need the sparsity, cheat and tell PETSc the matrix has been assembled */
498   (*Gmat2)->assembled = PETSC_TRUE;
499   PetscFunctionReturn(0);
500 }
501 
502 /* -------------------------------------------------------------------------- */
503 /*
504    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
505                     by setting data structures and options.
506 
507    Input Parameter:
508 .  pc - the preconditioner context
509 
510 */
511 PetscErrorCode PCSetUp_GAMG(PC pc)
512 {
513   PetscErrorCode ierr;
514   PC_MG          *mg      = (PC_MG*)pc->data;
515   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
516   Mat            Pmat     = pc->pmat;
517   PetscInt       fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS];
518   MPI_Comm       comm;
519   PetscMPIInt    rank,size,nactivepe;
520   Mat            Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS];
521   IS             *ASMLocalIDsArr[PETSC_MG_MAXLEVELS];
522   PetscLogDouble nnz0=0.,nnztot=0.;
523   MatInfo        info;
524   PetscBool      is_last = PETSC_FALSE;
525 
526   PetscFunctionBegin;
527   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
528   PetscCallMPI(MPI_Comm_rank(comm,&rank));
529   PetscCallMPI(MPI_Comm_size(comm,&size));
530 
531   if (pc->setupcalled) {
532     if (!pc_gamg->reuse_prol || pc->flag == DIFFERENT_NONZERO_PATTERN) {
533       /* reset everything */
534       PetscCall(PCReset_MG(pc));
535       pc->setupcalled = 0;
536     } else {
537       PC_MG_Levels **mglevels = mg->levels;
538       /* just do Galerkin grids */
539       Mat          B,dA,dB;
540 
541       if (pc_gamg->Nlevels > 1) {
542         PetscInt gl;
543         /* currently only handle case where mat and pmat are the same on coarser levels */
544         PetscCall(KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB));
545         /* (re)set to get dirty flag */
546         PetscCall(KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB));
547 
548         for (level=pc_gamg->Nlevels-2,gl=0; level>=0; level--,gl++) {
549           MatReuse reuse = MAT_INITIAL_MATRIX ;
550 
551           /* matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */
552           PetscCall(KSPGetOperators(mglevels[level]->smoothd,NULL,&B));
553           if (B->product) {
554             if (B->product->A == dB && B->product->B == mglevels[level+1]->interpolate) {
555               reuse = MAT_REUSE_MATRIX;
556             }
557           }
558           if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDestroy(&mglevels[level]->A));
559           if (reuse == MAT_REUSE_MATRIX) {
560             PetscCall(PetscInfo(pc,"%s: RAP after first solve, reuse matrix level %D\n",((PetscObject)pc)->prefix,level));
561           } else {
562             PetscCall(PetscInfo(pc,"%s: RAP after first solve, new matrix level %D\n",((PetscObject)pc)->prefix,level));
563           }
564           PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0));
565           PetscCall(MatPtAP(dB,mglevels[level+1]->interpolate,reuse,PETSC_DEFAULT,&B));
566           PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[gl][1],0,0,0,0));
567           if (reuse == MAT_INITIAL_MATRIX) mglevels[level]->A = B;
568           PetscCall(KSPSetOperators(mglevels[level]->smoothd,B,B));
569           dB   = B;
570         }
571       }
572 
573       PetscCall(PCSetUp_MG(pc));
574       PetscFunctionReturn(0);
575     }
576   }
577 
578   if (!pc_gamg->data) {
579     if (pc_gamg->orig_data) {
580       PetscCall(MatGetBlockSize(Pmat, &bs));
581       PetscCall(MatGetLocalSize(Pmat, &qq, NULL));
582 
583       pc_gamg->data_sz        = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols;
584       pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows;
585       pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols;
586 
587       PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data));
588       for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq];
589     } else {
590       PetscCheck(pc_gamg->ops->createdefaultdata,comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data");
591       PetscCall(pc_gamg->ops->createdefaultdata(pc,Pmat));
592     }
593   }
594 
595   /* cache original data for reuse */
596   if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) {
597     PetscCall(PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data));
598     for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq];
599     pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows;
600     pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols;
601   }
602 
603   /* get basic dims */
604   PetscCall(MatGetBlockSize(Pmat, &bs));
605   PetscCall(MatGetSize(Pmat, &M, &N));
606 
607   PetscCall(MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info)); /* global reduction */
608   nnz0   = info.nz_used;
609   nnztot = info.nz_used;
610   PetscCall(PetscInfo(pc,"%s: level %D) N=%D, n data rows=%D, n data cols=%D, nnz/row (ave)=%d, np=%D\n",((PetscObject)pc)->prefix,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(int)(nnz0/(PetscReal)M+0.5),size));
611 
612   /* Get A_i and R_i */
613   for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) {
614     pc_gamg->current_level = level;
615     PetscCheck(level < PETSC_MG_MAXLEVELS,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level);
616     level1 = level + 1;
617     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0));
618 #if defined(GAMG_STAGES)
619     PetscCall(PetscLogStagePush(gamg_stages[level]));
620 #endif
621     { /* construct prolongator */
622       Mat              Gmat;
623       PetscCoarsenData *agg_lists;
624       Mat              Prol11;
625 
626       PetscCall(pc_gamg->ops->graph(pc,Aarr[level], &Gmat));
627       PetscCall(pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists));
628       PetscCall(pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11));
629 
630       /* could have failed to create new level */
631       if (Prol11) {
632         const char *prefix;
633         char       addp[32];
634 
635         /* get new block size of coarse matrices */
636         PetscCall(MatGetBlockSizes(Prol11, NULL, &bs));
637 
638         if (pc_gamg->ops->optprolongator) {
639           /* smooth */
640           PetscCall(pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11));
641         }
642 
643         if (pc_gamg->use_aggs_in_asm) {
644           PetscInt bs;
645           PetscCall(MatGetBlockSizes(Prol11, &bs, NULL));
646           PetscCall(PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]));
647         }
648 
649         PetscCall(PCGetOptionsPrefix(pc,&prefix));
650         PetscCall(MatSetOptionsPrefix(Prol11,prefix));
651         PetscCall(PetscSNPrintf(addp,sizeof(addp),"pc_gamg_prolongator_%d_",(int)level));
652         PetscCall(MatAppendOptionsPrefix(Prol11,addp));
653         /* Always generate the transpose with CUDA
654            Such behaviour can be adapted with -pc_gamg_prolongator_ prefixed options */
655         PetscCall(MatSetOption(Prol11,MAT_FORM_EXPLICIT_TRANSPOSE,PETSC_TRUE));
656         PetscCall(MatSetFromOptions(Prol11));
657         Parr[level1] = Prol11;
658       } else Parr[level1] = NULL; /* failed to coarsen */
659 
660       PetscCall(MatDestroy(&Gmat));
661       PetscCall(PetscCDDestroy(agg_lists));
662     } /* construct prolongator scope */
663     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0));
664     if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */
665     if (!Parr[level1]) { /* failed to coarsen */
666       PetscCall(PetscInfo(pc,"%s: Stop gridding, level %D\n",((PetscObject)pc)->prefix,level));
667 #if defined(GAMG_STAGES)
668       PetscCall(PetscLogStagePop());
669 #endif
670       break;
671     }
672     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0));
673     PetscCall(MatGetSize(Parr[level1], &M, &N)); /* N is next M, a loop test variables */
674     PetscCheck(!is_last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ?");
675     if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE;
676     if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE;
677     PetscCall(pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last));
678 
679     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0));
680     PetscCall(MatGetSize(Aarr[level1], &M, &N)); /* M is loop test variables */
681     PetscCall(MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info));
682     nnztot += info.nz_used;
683     PetscCall(PetscInfo(pc,"%s: %D) N=%D, n data cols=%D, nnz/row (ave)=%d, %D active pes\n",((PetscObject)pc)->prefix,level1,M,pc_gamg->data_cell_cols,(int)(info.nz_used/(PetscReal)M),nactivepe));
684 
685 #if defined(GAMG_STAGES)
686     PetscCall(PetscLogStagePop());
687 #endif
688     /* stop if one node or one proc -- could pull back for singular problems */
689     if ((pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2)) {
690       PetscCall(PetscInfo(pc,"%s: HARD stop of coarsening on level %D.  Grid too small: %D block nodes\n",((PetscObject)pc)->prefix,level,M/bs));
691       level++;
692       break;
693     }
694   } /* levels */
695   PetscCall(PetscFree(pc_gamg->data));
696 
697   PetscCall(PetscInfo(pc,"%s: %D levels, grid complexity = %g\n",((PetscObject)pc)->prefix,level+1,nnztot/nnz0));
698   pc_gamg->Nlevels = level + 1;
699   fine_level       = level;
700   PetscCall(PCMGSetLevels(pc,pc_gamg->Nlevels,NULL));
701 
702   if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */
703 
704     /* set default smoothers & set operators */
705     for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) {
706       KSP smoother;
707       PC  subpc;
708 
709       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
710       PetscCall(KSPGetPC(smoother, &subpc));
711 
712       PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
713       /* set ops */
714       PetscCall(KSPSetOperators(smoother, Aarr[level], Aarr[level]));
715       PetscCall(PCMGSetInterpolation(pc, lidx, Parr[level+1]));
716 
717       /* set defaults */
718       PetscCall(KSPSetType(smoother, KSPCHEBYSHEV));
719 
720       /* set blocks for ASM smoother that uses the 'aggregates' */
721       if (pc_gamg->use_aggs_in_asm) {
722         PetscInt sz;
723         IS       *iss;
724 
725         sz   = nASMBlocksArr[level];
726         iss  = ASMLocalIDsArr[level];
727         PetscCall(PCSetType(subpc, PCASM));
728         PetscCall(PCASMSetOverlap(subpc, 0));
729         PetscCall(PCASMSetType(subpc,PC_ASM_BASIC));
730         if (!sz) {
731           IS       is;
732           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is));
733           PetscCall(PCASMSetLocalSubdomains(subpc, 1, NULL, &is));
734           PetscCall(ISDestroy(&is));
735         } else {
736           PetscInt kk;
737           PetscCall(PCASMSetLocalSubdomains(subpc, sz, NULL, iss));
738           for (kk=0; kk<sz; kk++) {
739             PetscCall(ISDestroy(&iss[kk]));
740           }
741           PetscCall(PetscFree(iss));
742         }
743         ASMLocalIDsArr[level] = NULL;
744         nASMBlocksArr[level]  = 0;
745       } else {
746         PetscCall(PCSetType(subpc, PCJACOBI));
747       }
748     }
749     {
750       /* coarse grid */
751       KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first;
752       Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0;
753 
754       PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
755       PetscCall(KSPSetOperators(smoother, Lmat, Lmat));
756       if (!pc_gamg->use_parallel_coarse_grid_solver) {
757         PetscCall(KSPSetNormType(smoother, KSP_NORM_NONE));
758         PetscCall(KSPGetPC(smoother, &subpc));
759         PetscCall(PCSetType(subpc, PCBJACOBI));
760         PetscCall(PCSetUp(subpc));
761         PetscCall(PCBJacobiGetSubKSP(subpc,&ii,&first,&k2));
762         PetscCheck(ii == 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii);
763         PetscCall(KSPGetPC(k2[0],&pc2));
764         PetscCall(PCSetType(pc2, PCLU));
765         PetscCall(PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS));
766         PetscCall(KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1));
767         PetscCall(KSPSetType(k2[0], KSPPREONLY));
768       }
769     }
770 
771     /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
772     ierr = PetscObjectOptionsBegin((PetscObject)pc);PetscCall(ierr);
773     PetscCall(PCSetFromOptions_MG(PetscOptionsObject,pc));
774     ierr = PetscOptionsEnd();PetscCall(ierr);
775     PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL));
776 
777     /* setup cheby eigen estimates from SA */
778     if (pc_gamg->use_sa_esteig) {
779       for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) {
780         KSP       smoother;
781         PetscBool ischeb;
782 
783         PetscCall(PCMGGetSmoother(pc, lidx, &smoother));
784         PetscCall(PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb));
785         if (ischeb) {
786           KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data;
787 
788           // The command line will override these settings because KSPSetFromOptions is called in PCSetUp_MG
789           if (mg->max_eigen_DinvA[level] > 0) {
790             // SA uses Jacobi for P; we use SA estimates if the smoother is also Jacobi or if the user explicitly requested it.
791             // TODO: This should test whether it's the same Jacobi variant (DIAG, ROWSUM, etc.)
792             PetscReal emax,emin;
793 
794             emin = mg->min_eigen_DinvA[level];
795             emax = mg->max_eigen_DinvA[level];
796             PetscCall(PetscInfo(pc,"%s: PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %D (N=%D) with emax = %g emin = %g\n",((PetscObject)pc)->prefix,level,Aarr[level]->rmap->N,(double)emax,(double)emin));
797             cheb->emin_provided = emin;
798             cheb->emax_provided = emax;
799           }
800         }
801       }
802     }
803 
804     PetscCall(PCSetUp_MG(pc));
805 
806     /* clean up */
807     for (level=1; level<pc_gamg->Nlevels; level++) {
808       PetscCall(MatDestroy(&Parr[level]));
809       PetscCall(MatDestroy(&Aarr[level]));
810     }
811   } else {
812     KSP smoother;
813 
814     PetscCall(PetscInfo(pc,"%s: One level solver used (system is seen as DD). Using default solver.\n"));
815     PetscCall(PCMGGetSmoother(pc, 0, &smoother));
816     PetscCall(KSPSetOperators(smoother, Aarr[0], Aarr[0]));
817     PetscCall(KSPSetType(smoother, KSPPREONLY));
818     PetscCall(PCSetUp_MG(pc));
819   }
820   PetscFunctionReturn(0);
821 }
822 
823 /* ------------------------------------------------------------------------- */
824 /*
825  PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
826    that was created with PCCreate_GAMG().
827 
828    Input Parameter:
829 .  pc - the preconditioner context
830 
831    Application Interface Routine: PCDestroy()
832 */
833 PetscErrorCode PCDestroy_GAMG(PC pc)
834 {
835   PC_MG          *mg     = (PC_MG*)pc->data;
836   PC_GAMG        *pc_gamg= (PC_GAMG*)mg->innerctx;
837 
838   PetscFunctionBegin;
839   PetscCall(PCReset_GAMG(pc));
840   if (pc_gamg->ops->destroy) {
841     PetscCall((*pc_gamg->ops->destroy)(pc));
842   }
843   PetscCall(PetscFree(pc_gamg->ops));
844   PetscCall(PetscFree(pc_gamg->gamg_type_name));
845   PetscCall(PetscFree(pc_gamg));
846   PetscCall(PCDestroy_MG(pc));
847   PetscFunctionReturn(0);
848 }
849 
850 /*@
851    PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction.
852 
853    Logically Collective on PC
854 
855    Input Parameters:
856 +  pc - the preconditioner context
857 -  n - the number of equations
858 
859    Options Database Key:
860 .  -pc_gamg_process_eq_limit <limit> - set the limit
861 
862    Notes:
863     GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process
864     that has degrees of freedom
865 
866    Level: intermediate
867 
868 .seealso: PCGAMGSetCoarseEqLim(), PCGAMGSetRankReductionFactors()
869 @*/
870 PetscErrorCode  PCGAMGSetProcEqLim(PC pc, PetscInt n)
871 {
872   PetscFunctionBegin;
873   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
874   PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));
875   PetscFunctionReturn(0);
876 }
877 
878 static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n)
879 {
880   PC_MG   *mg      = (PC_MG*)pc->data;
881   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
882 
883   PetscFunctionBegin;
884   if (n>0) pc_gamg->min_eq_proc = n;
885   PetscFunctionReturn(0);
886 }
887 
888 /*@
889    PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid.
890 
891  Collective on PC
892 
893    Input Parameters:
894 +  pc - the preconditioner context
895 -  n - maximum number of equations to aim for
896 
897    Options Database Key:
898 .  -pc_gamg_coarse_eq_limit <limit> - set the limit
899 
900    Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid
901      has less than 1000 unknowns.
902 
903    Level: intermediate
904 
905 .seealso: PCGAMGSetProcEqLim(), PCGAMGSetRankReductionFactors()
906 @*/
907 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n)
908 {
909   PetscFunctionBegin;
910   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
911   PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));
912   PetscFunctionReturn(0);
913 }
914 
915 static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n)
916 {
917   PC_MG   *mg      = (PC_MG*)pc->data;
918   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
919 
920   PetscFunctionBegin;
921   if (n>0) pc_gamg->coarse_eq_limit = n;
922   PetscFunctionReturn(0);
923 }
924 
925 /*@
926    PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids
927 
928    Collective on PC
929 
930    Input Parameters:
931 +  pc - the preconditioner context
932 -  n - PETSC_TRUE or PETSC_FALSE
933 
934    Options Database Key:
935 .  -pc_gamg_repartition <true,false> - turn on the repartitioning
936 
937    Notes:
938     this will generally improve the loading balancing of the work on each level
939 
940    Level: intermediate
941 
942 @*/
943 PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n)
944 {
945   PetscFunctionBegin;
946   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
947   PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));
948   PetscFunctionReturn(0);
949 }
950 
951 static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n)
952 {
953   PC_MG   *mg      = (PC_MG*)pc->data;
954   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
955 
956   PetscFunctionBegin;
957   pc_gamg->repart = n;
958   PetscFunctionReturn(0);
959 }
960 
961 /*@
962    PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Chebyshev smoother
963 
964    Collective on PC
965 
966    Input Parameters:
967 +  pc - the preconditioner context
968 -  n - number of its
969 
970    Options Database Key:
971 .  -pc_gamg_use_sa_esteig <true,false> - use the eigen estimate
972 
973    Notes:
974    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$.
975    Eigenvalue estimates (based on a few CG or GMRES iterations) are computed to choose $\omega$ so that this is a stable smoothing operation.
976    If Chebyshev with Jacobi (diagonal) preconditioning is used for smoothing, then the eigenvalue estimates can be reused.
977    This option is only used when the smoother uses Jacobi, and should be turned off if a different PCJacobiType is used.
978    It became default in PETSc 3.17.
979 
980    Level: advanced
981 
982 .seealso: KSPChebyshevSetEigenvalues(), KSPChebyshevEstEigSet()
983 @*/
984 PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n)
985 {
986   PetscFunctionBegin;
987   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
988   PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));
989   PetscFunctionReturn(0);
990 }
991 
992 static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscBool n)
993 {
994   PC_MG   *mg      = (PC_MG*)pc->data;
995   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
996 
997   PetscFunctionBegin;
998   pc_gamg->use_sa_esteig = n;
999   PetscFunctionReturn(0);
1000 }
1001 
1002 /*@
1003    PCGAMGSetEigenvalues - Set eigenvalues
1004 
1005    Collective on PC
1006 
1007    Input Parameters:
1008 +  pc - the preconditioner context
1009 -  emax - max eigenvalue
1010 .  emin - min eigenvalue
1011 
1012    Options Database Key:
1013 .  -pc_gamg_eigenvalues <emin,emax> - estimates of the eigenvalues
1014 
1015    Level: intermediate
1016 
1017 .seealso: PCGAMGSetUseSAEstEig()
1018 @*/
1019 PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin)
1020 {
1021   PetscFunctionBegin;
1022   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1023   PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));
1024   PetscFunctionReturn(0);
1025 }
1026 
1027 static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin)
1028 {
1029   PC_MG          *mg      = (PC_MG*)pc->data;
1030   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1031 
1032   PetscFunctionBegin;
1033   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);
1034   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);
1035   pc_gamg->emax = emax;
1036   pc_gamg->emin = emin;
1037   PetscFunctionReturn(0);
1038 }
1039 
1040 /*@
1041    PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner
1042 
1043    Collective on PC
1044 
1045    Input Parameters:
1046 +  pc - the preconditioner context
1047 -  n - PETSC_TRUE or PETSC_FALSE
1048 
1049    Options Database Key:
1050 .  -pc_gamg_reuse_interpolation <true,false> - reuse the previous interpolation
1051 
1052    Level: intermediate
1053 
1054    Notes:
1055     May negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows
1056     rebuilding the preconditioner quicker.
1057 
1058 @*/
1059 PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n)
1060 {
1061   PetscFunctionBegin;
1062   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1063   PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));
1064   PetscFunctionReturn(0);
1065 }
1066 
1067 static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n)
1068 {
1069   PC_MG   *mg      = (PC_MG*)pc->data;
1070   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1071 
1072   PetscFunctionBegin;
1073   pc_gamg->reuse_prol = n;
1074   PetscFunctionReturn(0);
1075 }
1076 
1077 /*@
1078    PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner.
1079 
1080    Collective on PC
1081 
1082    Input Parameters:
1083 +  pc - the preconditioner context
1084 -  flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not
1085 
1086    Options Database Key:
1087 .  -pc_gamg_asm_use_agg <true,false> - use aggregates to define the additive Schwarz subdomains
1088 
1089    Level: intermediate
1090 
1091 @*/
1092 PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg)
1093 {
1094   PetscFunctionBegin;
1095   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1096   PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));
1097   PetscFunctionReturn(0);
1098 }
1099 
1100 static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg)
1101 {
1102   PC_MG   *mg      = (PC_MG*)pc->data;
1103   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1104 
1105   PetscFunctionBegin;
1106   pc_gamg->use_aggs_in_asm = flg;
1107   PetscFunctionReturn(0);
1108 }
1109 
1110 /*@
1111    PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver
1112 
1113    Collective on PC
1114 
1115    Input Parameters:
1116 +  pc - the preconditioner context
1117 -  flg - PETSC_TRUE to not force coarse grid onto one processor
1118 
1119    Options Database Key:
1120 .  -pc_gamg_use_parallel_coarse_grid_solver - use a parallel coarse grid direct solver
1121 
1122    Level: intermediate
1123 
1124 .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids()
1125 @*/
1126 PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg)
1127 {
1128   PetscFunctionBegin;
1129   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1130   PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));
1131   PetscFunctionReturn(0);
1132 }
1133 
1134 static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg)
1135 {
1136   PC_MG   *mg      = (PC_MG*)pc->data;
1137   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1138 
1139   PetscFunctionBegin;
1140   pc_gamg->use_parallel_coarse_grid_solver = flg;
1141   PetscFunctionReturn(0);
1142 }
1143 
1144 /*@
1145    PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU
1146 
1147    Collective on PC
1148 
1149    Input Parameters:
1150 +  pc - the preconditioner context
1151 -  flg - PETSC_TRUE to pin coarse grids to CPU
1152 
1153    Options Database Key:
1154 .  -pc_gamg_cpu_pin_coarse_grids - pin the coarse grids to the CPU
1155 
1156    Level: intermediate
1157 
1158 .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve()
1159 @*/
1160 PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg)
1161 {
1162   PetscFunctionBegin;
1163   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1164   PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));
1165   PetscFunctionReturn(0);
1166 }
1167 
1168 static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg)
1169 {
1170   PC_MG   *mg      = (PC_MG*)pc->data;
1171   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1172 
1173   PetscFunctionBegin;
1174   pc_gamg->cpu_pin_coarse_grids = flg;
1175   PetscFunctionReturn(0);
1176 }
1177 
1178 /*@
1179    PCGAMGSetCoarseGridLayoutType - place coarse grids on processors with natural order (compact type)
1180 
1181    Collective on PC
1182 
1183    Input Parameters:
1184 +  pc - the preconditioner context
1185 -  flg - Layout type
1186 
1187    Options Database Key:
1188 .  -pc_gamg_coarse_grid_layout_type - place the coarse grids with natural ordering
1189 
1190    Level: intermediate
1191 
1192 .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids()
1193 @*/
1194 PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg)
1195 {
1196   PetscFunctionBegin;
1197   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1198   PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));
1199   PetscFunctionReturn(0);
1200 }
1201 
1202 static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg)
1203 {
1204   PC_MG   *mg      = (PC_MG*)pc->data;
1205   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1206 
1207   PetscFunctionBegin;
1208   pc_gamg->layout_type = flg;
1209   PetscFunctionReturn(0);
1210 }
1211 
1212 /*@
1213    PCGAMGSetNlevels -  Sets the maximum number of levels PCGAMG will use
1214 
1215    Not collective on PC
1216 
1217    Input Parameters:
1218 +  pc - the preconditioner
1219 -  n - the maximum number of levels to use
1220 
1221    Options Database Key:
1222 .  -pc_mg_levels <n> - set the maximum number of levels to allow
1223 
1224    Level: intermediate
1225 
1226 @*/
1227 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n)
1228 {
1229   PetscFunctionBegin;
1230   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1231   PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));
1232   PetscFunctionReturn(0);
1233 }
1234 
1235 static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n)
1236 {
1237   PC_MG   *mg      = (PC_MG*)pc->data;
1238   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1239 
1240   PetscFunctionBegin;
1241   pc_gamg->Nlevels = n;
1242   PetscFunctionReturn(0);
1243 }
1244 
1245 /*@
1246    PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph
1247 
1248    Not collective on PC
1249 
1250    Input Parameters:
1251 +  pc - the preconditioner context
1252 .  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
1253 -  n - number of threshold values provided in array
1254 
1255    Options Database Key:
1256 .  -pc_gamg_threshold <threshold> - the threshold to drop edges
1257 
1258    Notes:
1259     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.
1260     Before coarsening or aggregating the graph, GAMG removes small values from the graph with this threshold, and thus reducing the coupling in the graph and a different (perhaps better) coarser set of points.
1261 
1262     If n is less than the total number of coarsenings (see PCGAMGSetNlevels()), then threshold scaling (see PCGAMGSetThresholdScale()) is used for each successive coarsening.
1263     In this case, PCGAMGSetThresholdScale() must be called before PCGAMGSetThreshold().
1264     If n is greater than the total number of levels, the excess entries in threshold will not be used.
1265 
1266    Level: intermediate
1267 
1268 .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph()
1269 @*/
1270 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n)
1271 {
1272   PetscFunctionBegin;
1273   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1274   if (n) PetscValidRealPointer(v,2);
1275   PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));
1276   PetscFunctionReturn(0);
1277 }
1278 
1279 static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n)
1280 {
1281   PC_MG   *mg      = (PC_MG*)pc->data;
1282   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1283   PetscInt i;
1284   PetscFunctionBegin;
1285   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->threshold[i] = v[i];
1286   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;
1287   PetscFunctionReturn(0);
1288 }
1289 
1290 /*@
1291    PCGAMGSetRankReductionFactors - Set manual schedule for process reduction on coarse grids
1292 
1293    Collective on PC
1294 
1295    Input Parameters:
1296 +  pc - the preconditioner context
1297 .  v - array of reduction factors. 0 for fist value forces a reduction to one process/device on first level in Cuda
1298 -  n - number of values provided in array
1299 
1300    Options Database Key:
1301 .  -pc_gamg_rank_reduction_factors <factors> - provide the schedule
1302 
1303    Level: intermediate
1304 
1305 .seealso: PCGAMGSetProcEqLim(), PCGAMGSetCoarseEqLim()
1306 @*/
1307 PetscErrorCode PCGAMGSetRankReductionFactors(PC pc, PetscInt v[], PetscInt n)
1308 {
1309   PetscFunctionBegin;
1310   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1311   if (n) PetscValidIntPointer(v,2);
1312   PetscTryMethod(pc,"PCGAMGSetRankReductionFactors_C",(PC,PetscInt[],PetscInt),(pc,v,n));
1313   PetscFunctionReturn(0);
1314 }
1315 
1316 static PetscErrorCode PCGAMGSetRankReductionFactors_GAMG(PC pc, PetscInt v[], PetscInt n)
1317 {
1318   PC_MG   *mg      = (PC_MG*)pc->data;
1319   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1320   PetscInt i;
1321   PetscFunctionBegin;
1322   for (i=0; i<PetscMin(n,PETSC_MG_MAXLEVELS); i++) pc_gamg->level_reduction_factors[i] = v[i];
1323   for (; i<PETSC_MG_MAXLEVELS; i++) pc_gamg->level_reduction_factors[i] = -1; /* 0 stop putting one process/device on first level */
1324   PetscFunctionReturn(0);
1325 }
1326 
1327 /*@
1328    PCGAMGSetThresholdScale - Relative threshold reduction at each level
1329 
1330    Not collective on PC
1331 
1332    Input Parameters:
1333 +  pc - the preconditioner context
1334 -  scale - the threshold value reduction, usually < 1.0
1335 
1336    Options Database Key:
1337 .  -pc_gamg_threshold_scale <v> - set the relative threshold reduction on each level
1338 
1339    Notes:
1340    The initial threshold (for an arbitrary number of levels starting from the finest) can be set with PCGAMGSetThreshold().
1341    This scaling is used for each subsequent coarsening, but must be called before PCGAMGSetThreshold().
1342 
1343    Level: advanced
1344 
1345 .seealso: PCGAMGSetThreshold()
1346 @*/
1347 PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v)
1348 {
1349   PetscFunctionBegin;
1350   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1351   PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));
1352   PetscFunctionReturn(0);
1353 }
1354 
1355 static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v)
1356 {
1357   PC_MG   *mg      = (PC_MG*)pc->data;
1358   PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx;
1359   PetscFunctionBegin;
1360   pc_gamg->threshold_scale = v;
1361   PetscFunctionReturn(0);
1362 }
1363 
1364 /*@C
1365    PCGAMGSetType - Set solution method
1366 
1367    Collective on PC
1368 
1369    Input Parameters:
1370 +  pc - the preconditioner context
1371 -  type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL
1372 
1373    Options Database Key:
1374 .  -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply
1375 
1376    Level: intermediate
1377 
1378 .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType
1379 @*/
1380 PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type)
1381 {
1382   PetscFunctionBegin;
1383   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1384   PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));
1385   PetscFunctionReturn(0);
1386 }
1387 
1388 /*@C
1389    PCGAMGGetType - Get solution method
1390 
1391    Collective on PC
1392 
1393    Input Parameter:
1394 .  pc - the preconditioner context
1395 
1396    Output Parameter:
1397 .  type - the type of algorithm used
1398 
1399    Level: intermediate
1400 
1401 .seealso: PCGAMGSetType(), PCGAMGType
1402 @*/
1403 PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type)
1404 {
1405   PetscFunctionBegin;
1406   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1407   PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));
1408   PetscFunctionReturn(0);
1409 }
1410 
1411 static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type)
1412 {
1413   PC_MG          *mg      = (PC_MG*)pc->data;
1414   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1415 
1416   PetscFunctionBegin;
1417   *type = pc_gamg->type;
1418   PetscFunctionReturn(0);
1419 }
1420 
1421 static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type)
1422 {
1423   PC_MG          *mg      = (PC_MG*)pc->data;
1424   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1425   PetscErrorCode (*r)(PC);
1426 
1427   PetscFunctionBegin;
1428   pc_gamg->type = type;
1429   PetscCall(PetscFunctionListFind(GAMGList,type,&r));
1430   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type);
1431   if (pc_gamg->ops->destroy) {
1432     PetscCall((*pc_gamg->ops->destroy)(pc));
1433     PetscCall(PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps)));
1434     pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1435     /* cleaning up common data in pc_gamg - this should disapear someday */
1436     pc_gamg->data_cell_cols = 0;
1437     pc_gamg->data_cell_rows = 0;
1438     pc_gamg->orig_data_cell_cols = 0;
1439     pc_gamg->orig_data_cell_rows = 0;
1440     PetscCall(PetscFree(pc_gamg->data));
1441     pc_gamg->data_sz = 0;
1442   }
1443   PetscCall(PetscFree(pc_gamg->gamg_type_name));
1444   PetscCall(PetscStrallocpy(type,&pc_gamg->gamg_type_name));
1445   PetscCall((*r)(pc));
1446   PetscFunctionReturn(0);
1447 }
1448 
1449 static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer)
1450 {
1451   PC_MG          *mg      = (PC_MG*)pc->data;
1452   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1453   PetscReal       gc=0, oc=0;
1454 
1455   PetscFunctionBegin;
1456   PetscCall(PetscViewerASCIIPrintf(viewer,"    GAMG specific options\n"));
1457   PetscCall(PetscViewerASCIIPrintf(viewer,"      Threshold for dropping small values in graph on each level ="));
1458   for (PetscInt i=0;i<mg->nlevels; i++) PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]));
1459   PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
1460   PetscCall(PetscViewerASCIIPrintf(viewer,"      Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale));
1461   if (pc_gamg->use_aggs_in_asm) {
1462     PetscCall(PetscViewerASCIIPrintf(viewer,"      Using aggregates from coarsening process to define subdomains for PCASM\n"));
1463   }
1464   if (pc_gamg->use_parallel_coarse_grid_solver) {
1465     PetscCall(PetscViewerASCIIPrintf(viewer,"      Using parallel coarse grid solver (all coarse grid equations not put on one process)\n"));
1466   }
1467 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1468   if (pc_gamg->cpu_pin_coarse_grids) {
1469     /* PetscCall(PetscViewerASCIIPrintf(viewer,"      Pinning coarse grids to the CPU)\n")); */
1470   }
1471 #endif
1472   /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */
1473   /*   PetscCall(PetscViewerASCIIPrintf(viewer,"      Put reduced grids on processes in natural order (ie, 0,1,2...)\n")); */
1474   /* } else { */
1475   /*   PetscCall(PetscViewerASCIIPrintf(viewer,"      Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n")); */
1476   /* } */
1477   if (pc_gamg->ops->view) {
1478     PetscCall((*pc_gamg->ops->view)(pc,viewer));
1479   }
1480   PetscCall(PCMGGetGridComplexity(pc,&gc,&oc));
1481   PetscCall(PetscViewerASCIIPrintf(viewer,"      Complexity:    grid = %g    operator = %g\n",gc,oc));
1482   PetscFunctionReturn(0);
1483 }
1484 
1485 PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc)
1486 {
1487   PC_MG          *mg      = (PC_MG*)pc->data;
1488   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1489   PetscBool      flag;
1490   MPI_Comm       comm;
1491   char           prefix[256],tname[32];
1492   PetscInt       i,n;
1493   const char     *pcpre;
1494   static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",NULL};
1495   PetscFunctionBegin;
1496   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
1497   PetscCall(PetscOptionsHead(PetscOptionsObject,"GAMG options"));
1498   PetscCall(PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag));
1499   if (flag) {
1500     PetscCall(PCGAMGSetType(pc,tname));
1501   }
1502   PetscCall(PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL));
1503   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));
1504   PetscCall(PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL));
1505   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));
1506   PetscCall(PetscOptionsBool("-pc_gamg_use_parallel_coarse_grid_solver","Use parallel coarse grid solver (otherwise put last grid on one process)","PCGAMGSetUseParallelCoarseGridSolve",pc_gamg->use_parallel_coarse_grid_solver,&pc_gamg->use_parallel_coarse_grid_solver,NULL));
1507   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));
1508   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,(PetscEnum)pc_gamg->layout_type,(PetscEnum*)&pc_gamg->layout_type,NULL));
1509   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));
1510   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));
1511   PetscCall(PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL));
1512   n = PETSC_MG_MAXLEVELS;
1513   PetscCall(PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag));
1514   if (!flag || n < PETSC_MG_MAXLEVELS) {
1515     if (!flag) n = 1;
1516     i = n;
1517     do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS);
1518   }
1519   n = PETSC_MG_MAXLEVELS;
1520   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));
1521   if (!flag) i = 0;
1522   else i = n;
1523   do {pc_gamg->level_reduction_factors[i] = -1;} while (++i<PETSC_MG_MAXLEVELS);
1524   PetscCall(PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL));
1525   {
1526     PetscReal eminmax[2] = {0., 0.};
1527     n = 2;
1528     PetscCall(PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag));
1529     if (flag) {
1530       PetscCheckFalse(n != 2,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues");
1531       PetscCall(PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]));
1532     }
1533   }
1534   /* set options for subtype */
1535   if (pc_gamg->ops->setfromoptions) PetscCall((*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc));
1536 
1537   PetscCall(PCGetOptionsPrefix(pc, &pcpre));
1538   PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : ""));
1539   PetscCall(PetscOptionsTail());
1540   PetscFunctionReturn(0);
1541 }
1542 
1543 /* -------------------------------------------------------------------------- */
1544 /*MC
1545      PCGAMG - Geometric algebraic multigrid (AMG) preconditioner
1546 
1547    Options Database Keys:
1548 +   -pc_gamg_type <type> - one of agg, geo, or classical
1549 .   -pc_gamg_repartition  <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined
1550 .   -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations
1551 .   -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother
1552 .   -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit>
1553                                         equations on each process that has degrees of freedom
1554 .   -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for.
1555 .   -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level
1556 -   -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified
1557 
1558    Options Database Keys for default Aggregation:
1559 +  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
1560 .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
1561 -  -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it
1562 
1563    Multigrid options:
1564 +  -pc_mg_cycles <v> - v or w, see PCMGSetCycleType()
1565 .  -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp()
1566 .  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1567 -  -pc_mg_levels <levels> - Number of levels of multigrid to use.
1568 
1569   Notes:
1570     In order to obtain good performance for PCGAMG for vector valued problems you must
1571        Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point
1572        Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator
1573        See the Users Manual Chapter 4 for more details
1574 
1575   Level: intermediate
1576 
1577 .seealso:  PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(),
1578            PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig()
1579 M*/
1580 
1581 PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc)
1582 {
1583   PC_GAMG *pc_gamg;
1584   PC_MG   *mg;
1585 
1586   PetscFunctionBegin;
1587    /* register AMG type */
1588   PetscCall(PCGAMGInitializePackage());
1589 
1590   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
1591   PetscCall(PCSetType(pc, PCMG));
1592   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCGAMG));
1593 
1594   /* create a supporting struct and attach it to pc */
1595   PetscCall(PetscNewLog(pc,&pc_gamg));
1596   PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL));
1597   mg           = (PC_MG*)pc->data;
1598   mg->innerctx = pc_gamg;
1599 
1600   PetscCall(PetscNewLog(pc,&pc_gamg->ops));
1601 
1602   /* these should be in subctx but repartitioning needs simple arrays */
1603   pc_gamg->data_sz = 0;
1604   pc_gamg->data    = NULL;
1605 
1606   /* overwrite the pointers of PCMG by the functions of base class PCGAMG */
1607   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
1608   pc->ops->setup          = PCSetUp_GAMG;
1609   pc->ops->reset          = PCReset_GAMG;
1610   pc->ops->destroy        = PCDestroy_GAMG;
1611   mg->view                = PCView_GAMG;
1612 
1613   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG));
1614   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG));
1615   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG));
1616   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG));
1617   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG));
1618   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG));
1619   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG));
1620   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG));
1621   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG));
1622   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG));
1623   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG));
1624   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG));
1625   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG));
1626   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRankReductionFactors_C",PCGAMGSetRankReductionFactors_GAMG));
1627   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG));
1628   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG));
1629   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG));
1630   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG));
1631   pc_gamg->repart           = PETSC_FALSE;
1632   pc_gamg->reuse_prol       = PETSC_FALSE;
1633   pc_gamg->use_aggs_in_asm  = PETSC_FALSE;
1634   pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE;
1635   pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE;
1636   pc_gamg->layout_type      = PCGAMG_LAYOUT_SPREAD;
1637   pc_gamg->min_eq_proc      = 50;
1638   pc_gamg->coarse_eq_limit  = 50;
1639   PetscCall(PetscArrayzero(pc_gamg->threshold,PETSC_MG_MAXLEVELS));
1640   pc_gamg->threshold_scale = 1.;
1641   pc_gamg->Nlevels          = PETSC_MG_MAXLEVELS;
1642   pc_gamg->current_level    = 0; /* don't need to init really */
1643   pc_gamg->use_sa_esteig    = PETSC_TRUE;
1644   pc_gamg->emin             = 0;
1645   pc_gamg->emax             = 0;
1646 
1647   pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG;
1648 
1649   /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */
1650   PetscCall(PCGAMGSetType(pc,PCGAMGAGG));
1651   PetscFunctionReturn(0);
1652 }
1653 
1654 /*@C
1655  PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called
1656     from PCInitializePackage().
1657 
1658  Level: developer
1659 
1660  .seealso: PetscInitialize()
1661 @*/
1662 PetscErrorCode PCGAMGInitializePackage(void)
1663 {
1664   PetscInt       l;
1665 
1666   PetscFunctionBegin;
1667   if (PCGAMGPackageInitialized) PetscFunctionReturn(0);
1668   PCGAMGPackageInitialized = PETSC_TRUE;
1669   PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO));
1670   PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG));
1671   PetscCall(PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical));
1672   PetscCall(PetscRegisterFinalize(PCGAMGFinalizePackage));
1673 
1674   /* general events */
1675   PetscCall(PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG));
1676   PetscCall(PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO));
1677   PetscCall(PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG));
1678   PetscCall(PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO));
1679   PetscCall(PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG));
1680   PetscCall(PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO));
1681   PetscCall(PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG));
1682 
1683   PetscCall(PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]));
1684   PetscCall(PetscLogEventRegister("  Create Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]));
1685   PetscCall(PetscLogEventRegister("  Filter Graph", PC_CLASSID, &petsc_gamg_setup_events[SET16]));
1686   /* PetscLogEventRegister("    G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */
1687   /* PetscLogEventRegister("    G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */
1688   /* PetscLogEventRegister("    G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */
1689   PetscCall(PetscLogEventRegister("  MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]));
1690   PetscCall(PetscLogEventRegister("  geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]));
1691   PetscCall(PetscLogEventRegister("  geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]));
1692   PetscCall(PetscLogEventRegister("    search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]));
1693   PetscCall(PetscLogEventRegister("  SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]));
1694   PetscCall(PetscLogEventRegister("  SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]));
1695   PetscCall(PetscLogEventRegister("  SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]));
1696   PetscCall(PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]));
1697   PetscCall(PetscLogEventRegister("  repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]));
1698   PetscCall(PetscLogEventRegister("  Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]));
1699   PetscCall(PetscLogEventRegister("  Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]));
1700   PetscCall(PetscLogEventRegister("  Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]));
1701   for (l=0;l<PETSC_MG_MAXLEVELS;l++) {
1702     char ename[32];
1703 
1704     PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Squ l%02d",l));
1705     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][0]));
1706     PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Gal l%02d",l));
1707     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][1]));
1708     PetscCall(PetscSNPrintf(ename,sizeof(ename),"PCGAMG Opt l%02d",l));
1709     PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &petsc_gamg_setup_matmat_events[l][2]));
1710   }
1711   /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */
1712   /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */
1713   /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */
1714   /* create timer stages */
1715 #if defined(GAMG_STAGES)
1716   {
1717     char     str[32];
1718     PetscInt lidx;
1719     sprintf(str,"MG Level %d (finest)",0);
1720     PetscCall(PetscLogStageRegister(str, &gamg_stages[0]));
1721     for (lidx=1; lidx<9; lidx++) {
1722       sprintf(str,"MG Level %d",(int)lidx);
1723       PetscCall(PetscLogStageRegister(str, &gamg_stages[lidx]));
1724     }
1725   }
1726 #endif
1727   PetscFunctionReturn(0);
1728 }
1729 
1730 /*@C
1731  PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is
1732     called from PetscFinalize() automatically.
1733 
1734  Level: developer
1735 
1736  .seealso: PetscFinalize()
1737 @*/
1738 PetscErrorCode PCGAMGFinalizePackage(void)
1739 {
1740   PetscFunctionBegin;
1741   PCGAMGPackageInitialized = PETSC_FALSE;
1742   PetscCall(PetscFunctionListDestroy(&GAMGList));
1743   PetscFunctionReturn(0);
1744 }
1745 
1746 /*@C
1747  PCGAMGRegister - Register a PCGAMG implementation.
1748 
1749  Input Parameters:
1750  + type - string that will be used as the name of the GAMG type.
1751  - create - function for creating the gamg context.
1752 
1753   Level: advanced
1754 
1755  .seealso: PCGAMGType, PCGAMG, PCGAMGSetType()
1756 @*/
1757 PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC))
1758 {
1759   PetscFunctionBegin;
1760   PetscCall(PCGAMGInitializePackage());
1761   PetscCall(PetscFunctionListAdd(&GAMGList,type,create));
1762   PetscFunctionReturn(0);
1763 }
1764