xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision e0f5bfbec699682fa3e8b8532b1176849ea4e12a)
1 /*
2  GAMG geometric-algebric multigrid PC - Mark Adams 2011
3  */
4 
5 #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/
6 #include <petscblaslapack.h>
7 #include <petscdm.h>
8 #include <petsc/private/kspimpl.h>
9 
10 typedef struct {
11   PetscInt  nsmooths;
12   PetscBool symmetrize_graph;
13   PetscInt  aggressive_coarsening_levels; // number of aggressive coarsening levels (square or MISk)
14 } PC_GAMG_AGG;
15 
16 /*@
17    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used for multigrid on all the levels
18 
19    Logically Collective on pc
20 
21    Input Parameters:
22 .  pc - the preconditioner context
23 
24    Options Database Key:
25 .  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
26 
27    Level: intermediate
28 
29 .seealso: `PCMG`, `PCGAMG`
30 @*/
31 PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) {
32   PetscFunctionBegin;
33   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
34   PetscValidLogicalCollectiveInt(pc, n, 2);
35   PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n));
36   PetscFunctionReturn(0);
37 }
38 
39 static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) {
40   PC_MG       *mg          = (PC_MG *)pc->data;
41   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
42   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
43 
44   PetscFunctionBegin;
45   pc_gamg_agg->nsmooths = n;
46   PetscFunctionReturn(0);
47 }
48 
49 /*@
50    PCGAMGSetSymmetrizeGraph - Symmetrize the graph before computing the aggregation. Some algorithms require the graph be symmetric
51 
52    Logically Collective on pc
53 
54    Input Parameters:
55 +  pc - the preconditioner context
56 -  n - `PETSC_TRUE` or `PETSC_FALSE`
57 
58    Options Database Key:
59 .  -pc_gamg_symmetrize_graph <true,default=false> - symmetrize the graph before computing the aggregation
60 
61    Level: intermediate
62 
63    Developer Note:
64    If the aggregation can hang with a nonsymmetric graph then the graph should always be symmetrized before calling the aggregation,
65    it should not be up to the user.
66 
67 .seealso: `PCGAMG`, `PCGAMGSetAggressiveLevels()`
68 @*/
69 PetscErrorCode PCGAMGSetSymmetrizeGraph(PC pc, PetscBool n) {
70   PetscFunctionBegin;
71   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
72   PetscValidLogicalCollectiveBool(pc, n, 2);
73   PetscTryMethod(pc, "PCGAMGSetSymmetrizeGraph_C", (PC, PetscBool), (pc, n));
74   PetscFunctionReturn(0);
75 }
76 
77 static PetscErrorCode PCGAMGSetSymmetrizeGraph_AGG(PC pc, PetscBool n) {
78   PC_MG       *mg          = (PC_MG *)pc->data;
79   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
80   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
81 
82   PetscFunctionBegin;
83   pc_gamg_agg->symmetrize_graph = n;
84   PetscFunctionReturn(0);
85 }
86 
87 /*@
88    PCGAMGSetAggressiveLevels -  Use aggressive coarsening on first n levels
89 
90    Logically Collective on pc
91 
92    Input Parameters:
93 +  pc - the preconditioner context
94 -  n - 0, 1 or more
95 
96    Options Database Key:
97 .  -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels to square the graph on before aggregating it
98 
99    Level: intermediate
100 
101 .seealso: `PCGAMG`, `PCGAMGSetSymmetrizeGraph()`, `PCGAMGSetThreshold()`
102 @*/
103 PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n) {
104   PetscFunctionBegin;
105   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
106   PetscValidLogicalCollectiveInt(pc, n, 2);
107   PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n));
108   PetscFunctionReturn(0);
109 }
110 
111 static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n) {
112   PC_MG       *mg          = (PC_MG *)pc->data;
113   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
114   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
115 
116   PetscFunctionBegin;
117   pc_gamg_agg->aggressive_coarsening_levels = n;
118   PetscFunctionReturn(0);
119 }
120 
121 static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject) {
122   PC_MG       *mg          = (PC_MG *)pc->data;
123   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
124   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
125 
126   PetscFunctionBegin;
127   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options");
128   {
129     PetscBool flg;
130     PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "smoothing steps for smoothed aggregation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL));
131     PetscCall(PetscOptionsBool("-pc_gamg_symmetrize_graph", "Set for asymmetric matrices", "PCGAMGSetSymmetrizeGraph", pc_gamg_agg->symmetrize_graph, &pc_gamg_agg->symmetrize_graph, NULL));
132     pc_gamg_agg->aggressive_coarsening_levels = 1;
133     PetscCall(
134       PetscOptionsInt("-pc_gamg_square_graph", "Number of aggressive coarsening (MIS-2) levels from finest (alias for -pc_gamg_aggressive_coarsening, deprecated)", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, &flg));
135     if (!flg) {
136       PetscCall(PetscOptionsInt("-pc_gamg_aggressive_coarsening", "Number of aggressive coarsening (MIS-2) levels from finest", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, NULL));
137     } else {
138       PetscCall(PetscOptionsInt("-pc_gamg_aggressive_coarsening", "Number of aggressive coarsening (MIS-2) levels from finest", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, &flg));
139       if (flg) PetscCall(PetscInfo(pc, "Warning: both -pc_gamg_square_graph and -pc_gamg_aggressive_coarsening are used. -pc_gamg_square_graph is deprecated, Number of aggressive levels is %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels));
140     }
141   }
142   PetscOptionsHeadEnd();
143   PetscFunctionReturn(0);
144 }
145 
146 static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) {
147   PC_MG   *mg      = (PC_MG *)pc->data;
148   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
149 
150   PetscFunctionBegin;
151   PetscCall(PetscFree(pc_gamg->subctx));
152   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL));
153   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetSymmetrizeGraph_C", NULL));
154   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL));
155   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
156   PetscFunctionReturn(0);
157 }
158 
159 /*
160    PCSetCoordinates_AGG
161      - collective
162 
163    Input Parameter:
164    . pc - the preconditioner context
165    . ndm - dimesion of data (used for dof/vertex for Stokes)
166    . a_nloc - number of vertices local
167    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
168 */
169 
170 static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) {
171   PC_MG   *mg      = (PC_MG *)pc->data;
172   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
173   PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf;
174   Mat      mat = pc->pmat;
175 
176   PetscFunctionBegin;
177   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
178   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
179   nloc = a_nloc;
180 
181   /* SA: null space vectors */
182   PetscCall(MatGetBlockSize(mat, &ndf));               /* this does not work for Stokes */
183   if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
184   else if (coords) {
185     PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf);
186     pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */
187     if (ndm != ndf) PetscCheck(pc_gamg->data_cell_cols == ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Don't know how to create null space for ndm=%" PetscInt_FMT ", ndf=%" PetscInt_FMT ".  Use MatSetNearNullSpace().", ndm, ndf);
188   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
189   pc_gamg->data_cell_rows = ndatarows = ndf;
190   PetscCheck(pc_gamg->data_cell_cols > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "pc_gamg->data_cell_cols %" PetscInt_FMT " <= 0", pc_gamg->data_cell_cols);
191   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
192 
193   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
194     PetscCall(PetscFree(pc_gamg->data));
195     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
196   }
197   /* copy data in - column oriented */
198   for (kk = 0; kk < nloc; kk++) {
199     const PetscInt M    = nloc * pc_gamg->data_cell_rows; /* stride into data */
200     PetscReal     *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */
201     if (pc_gamg->data_cell_cols == 1) *data = 1.0;
202     else {
203       /* translational modes */
204       for (ii = 0; ii < ndatarows; ii++) {
205         for (jj = 0; jj < ndatarows; jj++) {
206           if (ii == jj) data[ii * M + jj] = 1.0;
207           else data[ii * M + jj] = 0.0;
208         }
209       }
210 
211       /* rotational modes */
212       if (coords) {
213         if (ndm == 2) {
214           data += 2 * M;
215           data[0] = -coords[2 * kk + 1];
216           data[1] = coords[2 * kk];
217         } else {
218           data += 3 * M;
219           data[0]         = 0.0;
220           data[M + 0]     = coords[3 * kk + 2];
221           data[2 * M + 0] = -coords[3 * kk + 1];
222           data[1]         = -coords[3 * kk + 2];
223           data[M + 1]     = 0.0;
224           data[2 * M + 1] = coords[3 * kk];
225           data[2]         = coords[3 * kk + 1];
226           data[M + 2]     = -coords[3 * kk];
227           data[2 * M + 2] = 0.0;
228         }
229       }
230     }
231   }
232   pc_gamg->data_sz = arrsz;
233   PetscFunctionReturn(0);
234 }
235 
236 /*
237    PCSetData_AGG - called if data is not set with PCSetCoordinates.
238       Looks in Mat for near null space.
239       Does not work for Stokes
240 
241   Input Parameter:
242    . pc -
243    . a_A - matrix to get (near) null space out of.
244 */
245 static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) {
246   PC_MG       *mg      = (PC_MG *)pc->data;
247   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
248   MatNullSpace mnull;
249 
250   PetscFunctionBegin;
251   PetscCall(MatGetNearNullSpace(a_A, &mnull));
252   if (!mnull) {
253     DM dm;
254     PetscCall(PCGetDM(pc, &dm));
255     if (!dm) PetscCall(MatGetDM(a_A, &dm));
256     if (dm) {
257       PetscObject deformation;
258       PetscInt    Nf;
259 
260       PetscCall(DMGetNumFields(dm, &Nf));
261       if (Nf) {
262         PetscCall(DMGetField(dm, 0, NULL, &deformation));
263         PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull));
264         if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull));
265       }
266     }
267   }
268 
269   if (!mnull) {
270     PetscInt bs, NN, MM;
271     PetscCall(MatGetBlockSize(a_A, &bs));
272     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
273     PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs);
274     PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL));
275   } else {
276     PetscReal         *nullvec;
277     PetscBool          has_const;
278     PetscInt           i, j, mlocal, nvec, bs;
279     const Vec         *vecs;
280     const PetscScalar *v;
281 
282     PetscCall(MatGetLocalSize(a_A, &mlocal, NULL));
283     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
284     pc_gamg->data_sz = (nvec + !!has_const) * mlocal;
285     PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
286     if (has_const)
287       for (i = 0; i < mlocal; i++) nullvec[i] = 1.0;
288     for (i = 0; i < nvec; i++) {
289       PetscCall(VecGetArrayRead(vecs[i], &v));
290       for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]);
291       PetscCall(VecRestoreArrayRead(vecs[i], &v));
292     }
293     pc_gamg->data           = nullvec;
294     pc_gamg->data_cell_cols = (nvec + !!has_const);
295     PetscCall(MatGetBlockSize(a_A, &bs));
296     pc_gamg->data_cell_rows = bs;
297   }
298   PetscFunctionReturn(0);
299 }
300 
301 /*
302   formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0
303 
304   Input Parameter:
305    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
306    . bs - row block size
307    . nSAvec - column bs of new P
308    . my0crs - global index of start of locals
309    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
310    . data_in[data_stride*nSAvec] - local data on fine grid
311    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
312 
313   Output Parameter:
314    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
315    . a_Prol - prolongation operator
316 */
317 static PetscErrorCode formProl0(PetscCoarsenData *agg_llists, PetscInt bs, PetscInt nSAvec, PetscInt my0crs, PetscInt data_stride, PetscReal data_in[], const PetscInt flid_fgid[], PetscReal **a_data_out, Mat a_Prol) {
318   PetscInt        Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride;
319   MPI_Comm        comm;
320   PetscReal      *out_data;
321   PetscCDIntNd   *pos;
322   PCGAMGHashTable fgid_flid;
323 
324   PetscFunctionBegin;
325   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
326   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
327   nloc = (Iend - Istart) / bs;
328   my0  = Istart / bs;
329   PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, Iend, Istart, bs);
330   Iend /= bs;
331   nghosts = data_stride / bs - nloc;
332 
333   PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid));
334   for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk));
335 
336   /* count selected -- same as number of cols of P */
337   for (nSelected = mm = 0; mm < nloc; mm++) {
338     PetscBool ise;
339     PetscCall(PetscCDEmptyAt(agg_llists, mm, &ise));
340     if (!ise) nSelected++;
341   }
342   PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj));
343   PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT "  != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs);
344   PetscCheck(nSelected == (jj - ii) / nSAvec, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nSelected %" PetscInt_FMT " != (jj %" PetscInt_FMT " - ii %" PetscInt_FMT ")/nSAvec %" PetscInt_FMT, nSelected, jj, ii, nSAvec);
345 
346   /* aloc space for coarse point data (output) */
347   out_data_stride = nSelected * nSAvec;
348 
349   PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data));
350   for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL;
351   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
352 
353   /* find points and set prolongation */
354   minsz = 100;
355   for (mm = clid = 0; mm < nloc; mm++) {
356     PetscCall(PetscCDSizeAt(agg_llists, mm, &jj));
357     if (jj > 0) {
358       const PetscInt lid = mm, cgid = my0crs + clid;
359       PetscInt       cids[100]; /* max bs */
360       PetscBLASInt   asz = jj, M = asz * bs, N = nSAvec, INFO;
361       PetscBLASInt   Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs;
362       PetscScalar   *qqc, *qqr, *TAU, *WORK;
363       PetscInt      *fids;
364       PetscReal     *data;
365 
366       /* count agg */
367       if (asz < minsz) minsz = asz;
368 
369       /* get block */
370       PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids));
371 
372       aggID = 0;
373       PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
374       while (pos) {
375         PetscInt gid1;
376         PetscCall(PetscCDIntNdGetID(pos, &gid1));
377         PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
378 
379         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
380         else {
381           PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid));
382           PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table");
383         }
384         /* copy in B_i matrix - column oriented */
385         data = &data_in[flid * bs];
386         for (ii = 0; ii < bs; ii++) {
387           for (jj = 0; jj < N; jj++) {
388             PetscReal d                       = data[jj * data_stride + ii];
389             qqc[jj * Mdata + aggID * bs + ii] = d;
390           }
391         }
392         /* set fine IDs */
393         for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk;
394         aggID++;
395       }
396 
397       /* pad with zeros */
398       for (ii = asz * bs; ii < Mdata; ii++) {
399         for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0;
400       }
401 
402       /* QR */
403       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
404       PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
405       PetscCall(PetscFPTrapPop());
406       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error");
407       /* get R - column oriented - output B_{i+1} */
408       {
409         PetscReal *data = &out_data[clid * nSAvec];
410         for (jj = 0; jj < nSAvec; jj++) {
411           for (ii = 0; ii < nSAvec; ii++) {
412             PetscCheck(data[jj * out_data_stride + ii] == PETSC_MAX_REAL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "data[jj*out_data_stride + ii] != %e", (double)PETSC_MAX_REAL);
413             if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]);
414             else data[jj * out_data_stride + ii] = 0.;
415           }
416         }
417       }
418 
419       /* get Q - row oriented */
420       PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
421       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO);
422 
423       for (ii = 0; ii < M; ii++) {
424         for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii];
425       }
426 
427       /* add diagonal block of P0 */
428       for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ }
429       PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES));
430       PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids));
431       clid++;
432     } /* coarse agg */
433   }   /* for all fine nodes */
434   PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
435   PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
436   PetscCall(PCGAMGHashTableDestroy(&fgid_flid));
437   PetscFunctionReturn(0);
438 }
439 
440 static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) {
441   PC_MG       *mg          = (PC_MG *)pc->data;
442   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
443   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
444 
445   PetscFunctionBegin;
446   PetscCall(PetscViewerASCIIPrintf(viewer, "      AGG specific options\n"));
447   PetscCall(PetscViewerASCIIPrintf(viewer, "        Symmetric graph %s\n", pc_gamg_agg->symmetrize_graph ? "true" : "false"));
448   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number of levels to square graph %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels));
449   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths));
450   PetscFunctionReturn(0);
451 }
452 
453 /*
454    PCGAMGGraph_AGG
455 
456   Input Parameter:
457    . pc - this
458    . Amat - matrix on this fine level
459 
460   Output Parameter:
461    . a_Gmat -
462 */
463 static PetscErrorCode PCGAMGGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) {
464   PC_MG                    *mg          = (PC_MG *)pc->data;
465   PC_GAMG                  *pc_gamg     = (PC_GAMG *)mg->innerctx;
466   const PetscReal           vfilter     = pc_gamg->threshold[pc_gamg->current_level];
467   PC_GAMG_AGG              *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
468   Mat                       Gmat, F = NULL;
469   MPI_Comm                  comm;
470   PetscBool /* set,flg , */ symm;
471 
472   PetscFunctionBegin;
473   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
474 
475   /* PetscCall(MatIsSymmetricKnown(Amat, &set, &flg)); || !(set && flg) -- this causes lot of symm calls */
476   symm = (PetscBool)(pc_gamg_agg->symmetrize_graph); /* && !pc_gamg_agg->square_graph; */
477 
478   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
479   PetscCall(MatCreateGraph(Amat, symm, PETSC_TRUE, &Gmat));
480   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
481 
482   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_FILTER], 0, 0, 0, 0));
483   PetscCall(MatFilter(Gmat, vfilter, &F));
484   if (F) {
485     PetscCall(MatDestroy(&Gmat));
486     Gmat = F;
487   }
488   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_FILTER], 0, 0, 0, 0));
489 
490   *a_Gmat = Gmat;
491 
492   PetscFunctionReturn(0);
493 }
494 
495 /*
496    PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for
497      communication of QR data used with HEM and MISk coarsening
498 
499   Input Parameter:
500    . a_pc - this
501 
502   Input/Output Parameter:
503    . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out)
504 
505   Output Parameter:
506    . agg_lists - list of aggregates
507 
508 */
509 static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) {
510   PC_MG       *mg          = (PC_MG *)a_pc->data;
511   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
512   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
513   Mat          mat, Gmat1 = *a_Gmat1; /* aggressive graph */
514   IS           perm;
515   PetscInt     Istart, Iend, Ii, nloc, bs, nn;
516   PetscInt    *permute, *degree;
517   PetscBool   *bIndexSet;
518   MatCoarsen   crs;
519   MPI_Comm     comm;
520   PetscReal    hashfact;
521   PetscInt     iSwapIndex;
522   PetscRandom  random;
523 
524   PetscFunctionBegin;
525   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
526   PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm));
527   PetscCall(MatGetLocalSize(Gmat1, &nn, NULL));
528   PetscCall(MatGetBlockSize(Gmat1, &bs));
529   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs);
530   nloc = nn / bs;
531 
532   /* get MIS aggs - randomize */
533   PetscCall(PetscMalloc2(nloc, &permute, nloc, &degree));
534   PetscCall(PetscCalloc1(nloc, &bIndexSet));
535   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
536   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random));
537   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
538   for (Ii = 0; Ii < nloc; Ii++) {
539     PetscInt nc;
540     PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
541     degree[Ii] = nc;
542     PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
543   }
544   for (Ii = 0; Ii < nloc; Ii++) {
545     PetscCall(PetscRandomGetValueReal(random, &hashfact));
546     iSwapIndex = (PetscInt)(hashfact * nloc) % nloc;
547     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
548       PetscInt iTemp        = permute[iSwapIndex];
549       permute[iSwapIndex]   = permute[Ii];
550       permute[Ii]           = iTemp;
551       iTemp                 = degree[iSwapIndex];
552       degree[iSwapIndex]    = degree[Ii];
553       degree[Ii]            = iTemp;
554       bIndexSet[iSwapIndex] = PETSC_TRUE;
555     }
556   }
557   // create minimum degree ordering
558   PetscCall(PetscSortIntWithArray(nloc, degree, permute));
559 
560   PetscCall(PetscFree(bIndexSet));
561   PetscCall(PetscRandomDestroy(&random));
562   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
563   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
564   PetscCall(MatCoarsenCreate(comm, &crs));
565   PetscCall(MatCoarsenSetFromOptions(crs));
566   PetscCall(MatCoarsenSetGreedyOrdering(crs, perm));
567   PetscCall(MatCoarsenSetAdjacency(crs, Gmat1));
568   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE));
569   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) PetscCall(MatCoarsenMISKSetDistance(crs, 2)); // hardwire to MIS-2
570   else PetscCall(MatCoarsenMISKSetDistance(crs, 1));                                                                    // MIS
571   PetscCall(MatCoarsenApply(crs));
572   PetscCall(MatCoarsenViewFromOptions(crs, NULL, "-mat_coarsen_view"));
573   PetscCall(MatCoarsenGetData(crs, agg_lists)); /* output */
574   PetscCall(MatCoarsenDestroy(&crs));
575 
576   PetscCall(ISDestroy(&perm));
577   PetscCall(PetscFree2(permute, degree));
578   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
579 
580   {
581     PetscCoarsenData *llist = *agg_lists;
582     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
583     PetscCall(PetscCDGetMat(llist, &mat));
584     if (mat) {
585       PetscCall(MatDestroy(&Gmat1));
586       *a_Gmat1 = mat; /* output */
587     }
588   }
589   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
590   PetscFunctionReturn(0);
591 }
592 
593 /*
594  PCGAMGProlongator_AGG
595 
596  Input Parameter:
597  . pc - this
598  . Amat - matrix on this fine level
599  . Graph - used to get ghost data for nodes in
600  . agg_lists - list of aggregates
601  Output Parameter:
602  . a_P_out - prolongation operator to the next level
603  */
604 static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out) {
605   PC_MG         *mg      = (PC_MG *)pc->data;
606   PC_GAMG       *pc_gamg = (PC_GAMG *)mg->innerctx;
607   const PetscInt col_bs  = pc_gamg->data_cell_cols;
608   PetscInt       Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs;
609   Mat            Prol;
610   PetscMPIInt    size;
611   MPI_Comm       comm;
612   PetscReal     *data_w_ghost;
613   PetscInt       myCrs0, nbnodes = 0, *flid_fgid;
614   MatType        mtype;
615 
616   PetscFunctionBegin;
617   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
618   PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1");
619   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
620   PetscCallMPI(MPI_Comm_size(comm, &size));
621   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
622   PetscCall(MatGetBlockSize(Amat, &bs));
623   nloc = (Iend - Istart) / bs;
624   my0  = Istart / bs;
625   PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT ") not divisible by bs %" PetscInt_FMT, Iend, Istart, bs);
626 
627   /* get 'nLocalSelected' */
628   for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) {
629     PetscBool ise;
630     /* filter out singletons 0 or 1? */
631     PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise));
632     if (!ise) nLocalSelected++;
633   }
634 
635   /* create prolongator, create P matrix */
636   PetscCall(MatGetType(Amat, &mtype));
637   PetscCall(MatCreate(comm, &Prol));
638   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE));
639   PetscCall(MatSetBlockSizes(Prol, bs, col_bs));
640   PetscCall(MatSetType(Prol, mtype));
641   PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL));
642   PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL));
643 
644   /* can get all points "removed" */
645   PetscCall(MatGetSize(Prol, &kk, &ii));
646   if (!ii) {
647     PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix));
648     PetscCall(MatDestroy(&Prol));
649     *a_P_out = NULL; /* out */
650     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
651     PetscFunctionReturn(0);
652   }
653   PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs));
654   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
655 
656   PetscCheck((kk - myCrs0) % col_bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(kk %" PetscInt_FMT " -myCrs0 %" PetscInt_FMT ") not divisible by col_bs %" PetscInt_FMT, kk, myCrs0, col_bs);
657   myCrs0 = myCrs0 / col_bs;
658   PetscCheck((kk / col_bs - myCrs0) == nLocalSelected, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(kk %" PetscInt_FMT "/col_bs %" PetscInt_FMT " - myCrs0 %" PetscInt_FMT ") != nLocalSelected %" PetscInt_FMT ")", kk, col_bs, myCrs0, nLocalSelected);
659 
660   /* create global vector of data in 'data_w_ghost' */
661   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
662   if (size > 1) { /* get ghost null space data */
663     PetscReal *tmp_gdata, *tmp_ldata, *tp2;
664     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
665     for (jj = 0; jj < col_bs; jj++) {
666       for (kk = 0; kk < bs; kk++) {
667         PetscInt         ii, stride;
668         const PetscReal *tp = pc_gamg->data + jj * bs * nloc + kk;
669         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
670 
671         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
672 
673         if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */
674           PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost));
675           nbnodes = bs * stride;
676         }
677         tp2 = data_w_ghost + jj * bs * stride + kk;
678         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
679         PetscCall(PetscFree(tmp_gdata));
680       }
681     }
682     PetscCall(PetscFree(tmp_ldata));
683   } else {
684     nbnodes      = bs * nloc;
685     data_w_ghost = (PetscReal *)pc_gamg->data;
686   }
687 
688   /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */
689   if (size > 1) {
690     PetscReal *fid_glid_loc, *fiddata;
691     PetscInt   stride;
692 
693     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
694     for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk);
695     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
696     PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */
697     for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
698     PetscCall(PetscFree(fiddata));
699 
700     PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs);
701     PetscCall(PetscFree(fid_glid_loc));
702   } else {
703     PetscCall(PetscMalloc1(nloc, &flid_fgid));
704     for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk;
705   }
706   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
707   /* get P0 */
708   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
709   {
710     PetscReal *data_out = NULL;
711     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol));
712     PetscCall(PetscFree(pc_gamg->data));
713 
714     pc_gamg->data           = data_out;
715     pc_gamg->data_cell_rows = col_bs;
716     pc_gamg->data_sz        = col_bs * col_bs * nLocalSelected;
717   }
718   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
719   if (size > 1) PetscCall(PetscFree(data_w_ghost));
720   PetscCall(PetscFree(flid_fgid));
721 
722   *a_P_out = Prol; /* out */
723 
724   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
725   PetscFunctionReturn(0);
726 }
727 
728 /*
729    PCGAMGOptProlongator_AGG
730 
731   Input Parameter:
732    . pc - this
733    . Amat - matrix on this fine level
734  In/Output Parameter:
735    . a_P - prolongation operator to the next level
736 */
737 static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P) {
738   PC_MG       *mg          = (PC_MG *)pc->data;
739   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
740   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
741   PetscInt     jj;
742   Mat          Prol = *a_P;
743   MPI_Comm     comm;
744   KSP          eksp;
745   Vec          bb, xx;
746   PC           epc;
747   PetscReal    alpha, emax, emin;
748 
749   PetscFunctionBegin;
750   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
751   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
752 
753   /* compute maximum singular value of operator to be used in smoother */
754   if (0 < pc_gamg_agg->nsmooths) {
755     /* get eigen estimates */
756     if (pc_gamg->emax > 0) {
757       emin = pc_gamg->emin;
758       emax = pc_gamg->emax;
759     } else {
760       const char *prefix;
761 
762       PetscCall(MatCreateVecs(Amat, &bb, NULL));
763       PetscCall(MatCreateVecs(Amat, &xx, NULL));
764       PetscCall(KSPSetNoisy_Private(bb));
765 
766       PetscCall(KSPCreate(comm, &eksp));
767       PetscCall(PCGetOptionsPrefix(pc, &prefix));
768       PetscCall(KSPSetOptionsPrefix(eksp, prefix));
769       PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_"));
770       {
771         PetscBool isset, sflg;
772         PetscCall(MatIsSPDKnown(Amat, &isset, &sflg));
773         if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG));
774       }
775       PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure));
776       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
777 
778       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
779       PetscCall(KSPSetOperators(eksp, Amat, Amat));
780 
781       PetscCall(KSPGetPC(eksp, &epc));
782       PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */
783 
784       PetscCall(KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10)); // 10 is safer, but 5 is often fine, can override with -pc_gamg_esteig_ksp_max_it -mg_levels_ksp_chebyshev_esteig 0,0.25,0,1.2
785 
786       PetscCall(KSPSetFromOptions(eksp));
787       PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE));
788       PetscCall(KSPSolve(eksp, bb, xx));
789       PetscCall(KSPCheckSolve(eksp, pc, xx));
790 
791       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
792       PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI));
793       PetscCall(VecDestroy(&xx));
794       PetscCall(VecDestroy(&bb));
795       PetscCall(KSPDestroy(&eksp));
796     }
797     if (pc_gamg->use_sa_esteig) {
798       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
799       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
800       PetscCall(PetscInfo(pc, "%s: Smooth P0: level %" PetscInt_FMT ", cache spectra %g %g\n", ((PetscObject)pc)->prefix, pc_gamg->current_level, (double)emin, (double)emax));
801     } else {
802       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
803       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
804     }
805   } else {
806     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
807     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
808   }
809 
810   /* smooth P0 */
811   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
812     Mat tMat;
813     Vec diag;
814 
815     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
816 
817     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
818     PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
819     PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat));
820     PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
821     PetscCall(MatProductClear(tMat));
822     PetscCall(MatCreateVecs(Amat, &diag, NULL));
823     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
824     PetscCall(VecReciprocal(diag));
825     PetscCall(MatDiagonalScale(tMat, diag, NULL));
826     PetscCall(VecDestroy(&diag));
827 
828     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
829     PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero");
830     /* TODO: Document the 1.4 and don't hardwire it in this routine */
831     alpha = -1.4 / emax;
832 
833     PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
834     PetscCall(MatDestroy(&Prol));
835     Prol = tMat;
836     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
837   }
838   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
839   *a_P = Prol;
840   PetscFunctionReturn(0);
841 }
842 
843 /*
844    PCCreateGAMG_AGG
845 
846   Input Parameter:
847    . pc -
848 */
849 PetscErrorCode PCCreateGAMG_AGG(PC pc) {
850   PC_MG       *mg      = (PC_MG *)pc->data;
851   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
852   PC_GAMG_AGG *pc_gamg_agg;
853 
854   PetscFunctionBegin;
855   /* create sub context for SA */
856   PetscCall(PetscNew(&pc_gamg_agg));
857   pc_gamg->subctx = pc_gamg_agg;
858 
859   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
860   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
861   /* reset does not do anything; setup not virtual */
862 
863   /* set internal function pointers */
864   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
865   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
866   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
867   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
868   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
869   pc_gamg->ops->view              = PCView_GAMG_AGG;
870 
871   pc_gamg_agg->aggressive_coarsening_levels = 0;
872   pc_gamg_agg->symmetrize_graph             = PETSC_FALSE;
873   pc_gamg_agg->nsmooths                     = 1;
874 
875   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG));
876   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetSymmetrizeGraph_C", PCGAMGSetSymmetrizeGraph_AGG));
877   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG));
878   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG));
879   PetscFunctionReturn(0);
880 }
881