xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision e0b7e82fd3cf27fce84cc3e37e8d70a5c36a2d4e) !
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;                     // number of smoothing steps to construct prolongation
12   PetscInt   aggressive_coarsening_levels; // number of aggressive coarsening levels (square or MISk)
13   PetscInt   aggressive_mis_k;             // the k in MIS-k
14   PetscBool  use_aggressive_square_graph;
15   PetscBool  use_minimum_degree_ordering;
16   PetscBool  use_low_mem_filter;
17   MatCoarsen crs;
18 } PC_GAMG_AGG;
19 
20 /*@
21   PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used to construct the prolongation operator
22 
23   Logically Collective
24 
25   Input Parameters:
26 + pc - the preconditioner context
27 - n  - the number of smooths
28 
29   Options Database Key:
30 . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use
31 
32   Level: intermediate
33 
34   Note:
35   This is a different concept from the number smoothing steps used during the linear solution process which
36   can be set with `-mg_levels_ksp_max_it`
37 
38   Developer Note:
39   This should be named `PCGAMGAGGSetNSmooths()`.
40 
41 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCMG`, `PCGAMG`
42 @*/
43 PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
44 {
45   PetscFunctionBegin;
46   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
47   PetscValidLogicalCollectiveInt(pc, n, 2);
48   PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n));
49   PetscFunctionReturn(PETSC_SUCCESS);
50 }
51 
52 static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n)
53 {
54   PC_MG       *mg          = (PC_MG *)pc->data;
55   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
56   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
57 
58   PetscFunctionBegin;
59   pc_gamg_agg->nsmooths = n;
60   PetscFunctionReturn(PETSC_SUCCESS);
61 }
62 
63 /*@
64   PCGAMGSetAggressiveLevels -  Use aggressive coarsening on first n levels
65 
66   Logically Collective
67 
68   Input Parameters:
69 + pc - the preconditioner context
70 - n  - 0, 1 or more
71 
72   Options Database Key:
73 . -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels on which to square the graph on before aggregating it
74 
75   Level: intermediate
76 
77 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()`
78 @*/
79 PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n)
80 {
81   PetscFunctionBegin;
82   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
83   PetscValidLogicalCollectiveInt(pc, n, 2);
84   PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n));
85   PetscFunctionReturn(PETSC_SUCCESS);
86 }
87 
88 /*@
89   PCGAMGMISkSetAggressive - Number (k) distance in MIS coarsening (>2 is 'aggressive')
90 
91   Logically Collective
92 
93   Input Parameters:
94 + pc - the preconditioner context
95 - n  - 1 or more (default = 2)
96 
97   Options Database Key:
98 . -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive')
99 
100   Level: intermediate
101 
102 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()`
103 @*/
104 PetscErrorCode PCGAMGMISkSetAggressive(PC pc, PetscInt n)
105 {
106   PetscFunctionBegin;
107   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
108   PetscValidLogicalCollectiveInt(pc, n, 2);
109   PetscTryMethod(pc, "PCGAMGMISkSetAggressive_C", (PC, PetscInt), (pc, n));
110   PetscFunctionReturn(PETSC_SUCCESS);
111 }
112 
113 /*@
114   PCGAMGSetAggressiveSquareGraph - Use graph square A'A for aggressive coarsening, old method
115 
116   Logically Collective
117 
118   Input Parameters:
119 + pc - the preconditioner context
120 - b  - default false - MIS-k is faster
121 
122   Options Database Key:
123 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening
124 
125   Level: intermediate
126 
127 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()`
128 @*/
129 PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC pc, PetscBool b)
130 {
131   PetscFunctionBegin;
132   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
133   PetscValidLogicalCollectiveBool(pc, b, 2);
134   PetscTryMethod(pc, "PCGAMGSetAggressiveSquareGraph_C", (PC, PetscBool), (pc, b));
135   PetscFunctionReturn(PETSC_SUCCESS);
136 }
137 
138 /*@
139   PCGAMGMISkSetMinDegreeOrdering - Use minimum degree ordering in greedy MIS algorithm
140 
141   Logically Collective
142 
143   Input Parameters:
144 + pc - the preconditioner context
145 - b  - default true
146 
147   Options Database Key:
148 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm
149 
150   Level: intermediate
151 
152 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGSetLowMemoryFilter()`
153 @*/
154 PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC pc, PetscBool b)
155 {
156   PetscFunctionBegin;
157   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
158   PetscValidLogicalCollectiveBool(pc, b, 2);
159   PetscTryMethod(pc, "PCGAMGMISkSetMinDegreeOrdering_C", (PC, PetscBool), (pc, b));
160   PetscFunctionReturn(PETSC_SUCCESS);
161 }
162 
163 /*@
164   PCGAMGSetLowMemoryFilter - Use low memory graph/matrix filter
165 
166   Logically Collective
167 
168   Input Parameters:
169 + pc - the preconditioner context
170 - b  - default false
171 
172   Options Database Key:
173 . -pc_gamg_low_memory_threshold_filter <bool,default=false> - Use low memory graph/matrix filter
174 
175   Level: intermediate
176 
177 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`,
178   `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`
179 @*/
180 PetscErrorCode PCGAMGSetLowMemoryFilter(PC pc, PetscBool b)
181 {
182   PetscFunctionBegin;
183   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
184   PetscValidLogicalCollectiveBool(pc, b, 2);
185   PetscTryMethod(pc, "PCGAMGSetLowMemoryFilter_C", (PC, PetscBool), (pc, b));
186   PetscFunctionReturn(PETSC_SUCCESS);
187 }
188 
189 static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n)
190 {
191   PC_MG       *mg          = (PC_MG *)pc->data;
192   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
193   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
194 
195   PetscFunctionBegin;
196   pc_gamg_agg->aggressive_coarsening_levels = n;
197   PetscFunctionReturn(PETSC_SUCCESS);
198 }
199 
200 static PetscErrorCode PCGAMGMISkSetAggressive_AGG(PC pc, PetscInt n)
201 {
202   PC_MG       *mg          = (PC_MG *)pc->data;
203   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
204   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
205 
206   PetscFunctionBegin;
207   pc_gamg_agg->aggressive_mis_k = n;
208   PetscFunctionReturn(PETSC_SUCCESS);
209 }
210 
211 static PetscErrorCode PCGAMGSetAggressiveSquareGraph_AGG(PC pc, PetscBool b)
212 {
213   PC_MG       *mg          = (PC_MG *)pc->data;
214   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
215   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
216 
217   PetscFunctionBegin;
218   pc_gamg_agg->use_aggressive_square_graph = b;
219   PetscFunctionReturn(PETSC_SUCCESS);
220 }
221 
222 static PetscErrorCode PCGAMGSetLowMemoryFilter_AGG(PC pc, PetscBool b)
223 {
224   PC_MG       *mg          = (PC_MG *)pc->data;
225   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
226   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
227 
228   PetscFunctionBegin;
229   pc_gamg_agg->use_low_mem_filter = b;
230   PetscFunctionReturn(PETSC_SUCCESS);
231 }
232 
233 static PetscErrorCode PCGAMGMISkSetMinDegreeOrdering_AGG(PC pc, PetscBool b)
234 {
235   PC_MG       *mg          = (PC_MG *)pc->data;
236   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
237   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
238 
239   PetscFunctionBegin;
240   pc_gamg_agg->use_minimum_degree_ordering = b;
241   PetscFunctionReturn(PETSC_SUCCESS);
242 }
243 
244 static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject)
245 {
246   PC_MG       *mg          = (PC_MG *)pc->data;
247   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
248   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
249   PetscBool    n_aggressive_flg, old_sq_provided = PETSC_FALSE, new_sq_provided = PETSC_FALSE, new_sqr_graph = pc_gamg_agg->use_aggressive_square_graph;
250   PetscInt     nsq_graph_old = 0;
251 
252   PetscFunctionBegin;
253   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options");
254   PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "number of smoothing steps to construct prolongation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL));
255   // aggressive coarsening logic with deprecated -pc_gamg_square_graph
256   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, &n_aggressive_flg));
257   if (!n_aggressive_flg)
258     PetscCall(PetscOptionsInt("-pc_gamg_square_graph", "Number of aggressive coarsening (MIS-2) levels from finest (deprecated alias for -pc_gamg_aggressive_coarsening)", "PCGAMGSetAggressiveLevels", nsq_graph_old, &nsq_graph_old, &old_sq_provided));
259   PetscCall(PetscOptionsBool("-pc_gamg_aggressive_square_graph", "Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening", "PCGAMGSetAggressiveSquareGraph", new_sqr_graph, &pc_gamg_agg->use_aggressive_square_graph, &new_sq_provided));
260   if (!new_sq_provided && old_sq_provided) {
261     pc_gamg_agg->aggressive_coarsening_levels = nsq_graph_old; // could be zero
262     pc_gamg_agg->use_aggressive_square_graph  = PETSC_TRUE;
263   }
264   if (new_sq_provided && old_sq_provided)
265     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));
266   PetscCall(PetscOptionsBool("-pc_gamg_mis_k_minimum_degree_ordering", "Use minimum degree ordering for greedy MIS", "PCGAMGMISkSetMinDegreeOrdering", pc_gamg_agg->use_minimum_degree_ordering, &pc_gamg_agg->use_minimum_degree_ordering, NULL));
267   PetscCall(PetscOptionsBool("-pc_gamg_low_memory_threshold_filter", "Use the (built-in) low memory graph/matrix filter", "PCGAMGSetLowMemoryFilter", pc_gamg_agg->use_low_mem_filter, &pc_gamg_agg->use_low_mem_filter, NULL));
268   PetscCall(PetscOptionsInt("-pc_gamg_aggressive_mis_k", "Number of levels of multigrid to use.", "PCGAMGMISkSetAggressive", pc_gamg_agg->aggressive_mis_k, &pc_gamg_agg->aggressive_mis_k, NULL));
269   PetscOptionsHeadEnd();
270   PetscFunctionReturn(PETSC_SUCCESS);
271 }
272 
273 static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
274 {
275   PC_MG       *mg          = (PC_MG *)pc->data;
276   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
277   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
278 
279   PetscFunctionBegin;
280   PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs));
281   PetscCall(PetscFree(pc_gamg->subctx));
282   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL));
283   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL));
284   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", NULL));
285   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", NULL));
286   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", NULL));
287   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", NULL));
288   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
289   PetscFunctionReturn(PETSC_SUCCESS);
290 }
291 
292 /*
293    PCSetCoordinates_AGG
294 
295    Collective
296 
297    Input Parameter:
298    . pc - the preconditioner context
299    . ndm - dimension of data (used for dof/vertex for Stokes)
300    . a_nloc - number of vertices local
301    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
302 */
303 
304 static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
305 {
306   PC_MG   *mg      = (PC_MG *)pc->data;
307   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
308   PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf;
309   Mat      mat = pc->pmat;
310 
311   PetscFunctionBegin;
312   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
313   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
314   nloc = a_nloc;
315 
316   /* SA: null space vectors */
317   PetscCall(MatGetBlockSize(mat, &ndf));               /* this does not work for Stokes */
318   if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
319   else if (coords) {
320     PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf);
321     pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */
322     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);
323   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
324   pc_gamg->data_cell_rows = ndatarows = ndf;
325   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);
326   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
327 
328   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
329     PetscCall(PetscFree(pc_gamg->data));
330     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
331   }
332   /* copy data in - column-oriented */
333   for (kk = 0; kk < nloc; kk++) {
334     const PetscInt M    = nloc * pc_gamg->data_cell_rows; /* stride into data */
335     PetscReal     *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */
336 
337     if (pc_gamg->data_cell_cols == 1) *data = 1.0;
338     else {
339       /* translational modes */
340       for (ii = 0; ii < ndatarows; ii++) {
341         for (jj = 0; jj < ndatarows; jj++) {
342           if (ii == jj) data[ii * M + jj] = 1.0;
343           else data[ii * M + jj] = 0.0;
344         }
345       }
346 
347       /* rotational modes */
348       if (coords) {
349         if (ndm == 2) {
350           data += 2 * M;
351           data[0] = -coords[2 * kk + 1];
352           data[1] = coords[2 * kk];
353         } else {
354           data += 3 * M;
355           data[0]         = 0.0;
356           data[M + 0]     = coords[3 * kk + 2];
357           data[2 * M + 0] = -coords[3 * kk + 1];
358           data[1]         = -coords[3 * kk + 2];
359           data[M + 1]     = 0.0;
360           data[2 * M + 1] = coords[3 * kk];
361           data[2]         = coords[3 * kk + 1];
362           data[M + 2]     = -coords[3 * kk];
363           data[2 * M + 2] = 0.0;
364         }
365       }
366     }
367   }
368   pc_gamg->data_sz = arrsz;
369   PetscFunctionReturn(PETSC_SUCCESS);
370 }
371 
372 /*
373    PCSetData_AGG - called if data is not set with PCSetCoordinates.
374       Looks in Mat for near null space.
375       Does not work for Stokes
376 
377   Input Parameter:
378    . pc -
379    . a_A - matrix to get (near) null space out of.
380 */
381 static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
382 {
383   PC_MG       *mg      = (PC_MG *)pc->data;
384   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
385   MatNullSpace mnull;
386 
387   PetscFunctionBegin;
388   PetscCall(MatGetNearNullSpace(a_A, &mnull));
389   if (!mnull) {
390     DM dm;
391 
392     PetscCall(PCGetDM(pc, &dm));
393     if (!dm) PetscCall(MatGetDM(a_A, &dm));
394     if (dm) {
395       PetscObject deformation;
396       PetscInt    Nf;
397 
398       PetscCall(DMGetNumFields(dm, &Nf));
399       if (Nf) {
400         PetscCall(DMGetField(dm, 0, NULL, &deformation));
401         PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull));
402         if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull));
403       }
404     }
405   }
406 
407   if (!mnull) {
408     PetscInt bs, NN, MM;
409 
410     PetscCall(MatGetBlockSize(a_A, &bs));
411     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
412     PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs);
413     PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL));
414   } else {
415     PetscReal         *nullvec;
416     PetscBool          has_const;
417     PetscInt           i, j, mlocal, nvec, bs;
418     const Vec         *vecs;
419     const PetscScalar *v;
420 
421     PetscCall(MatGetLocalSize(a_A, &mlocal, NULL));
422     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
423     for (i = 0; i < nvec; i++) {
424       PetscCall(VecGetLocalSize(vecs[i], &j));
425       PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal);
426     }
427     pc_gamg->data_sz = (nvec + !!has_const) * mlocal;
428     PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
429     if (has_const)
430       for (i = 0; i < mlocal; i++) nullvec[i] = 1.0;
431     for (i = 0; i < nvec; i++) {
432       PetscCall(VecGetArrayRead(vecs[i], &v));
433       for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]);
434       PetscCall(VecRestoreArrayRead(vecs[i], &v));
435     }
436     pc_gamg->data           = nullvec;
437     pc_gamg->data_cell_cols = (nvec + !!has_const);
438     PetscCall(MatGetBlockSize(a_A, &bs));
439     pc_gamg->data_cell_rows = bs;
440   }
441   PetscFunctionReturn(PETSC_SUCCESS);
442 }
443 
444 /*
445   formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0
446 
447   Input Parameter:
448    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
449    . bs - row block size
450    . nSAvec - column bs of new P
451    . my0crs - global index of start of locals
452    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
453    . data_in[data_stride*nSAvec] - local data on fine grid
454    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
455 
456   Output Parameter:
457    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
458    . a_Prol - prolongation operator
459 */
460 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)
461 {
462   PetscInt        Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride;
463   MPI_Comm        comm;
464   PetscReal      *out_data;
465   PetscCDIntNd   *pos;
466   PCGAMGHashTable fgid_flid;
467 
468   PetscFunctionBegin;
469   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
470   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
471   nloc = (Iend - Istart) / bs;
472   my0  = Istart / bs;
473   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);
474   Iend /= bs;
475   nghosts = data_stride / bs - nloc;
476 
477   PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid));
478   for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk));
479 
480   /* count selected -- same as number of cols of P */
481   for (nSelected = mm = 0; mm < nloc; mm++) {
482     PetscBool ise;
483 
484     PetscCall(PetscCDIsEmptyAt(agg_llists, mm, &ise));
485     if (!ise) nSelected++;
486   }
487   PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj));
488   PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT "  != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs);
489   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);
490 
491   /* aloc space for coarse point data (output) */
492   out_data_stride = nSelected * nSAvec;
493 
494   PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data));
495   for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL;
496   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
497 
498   /* find points and set prolongation */
499   minsz = 100;
500   for (mm = clid = 0; mm < nloc; mm++) {
501     PetscCall(PetscCDCountAt(agg_llists, mm, &jj));
502     if (jj > 0) {
503       const PetscInt lid = mm, cgid = my0crs + clid;
504       PetscInt       cids[100]; /* max bs */
505       PetscBLASInt   asz = jj, M = asz * bs, N = nSAvec, INFO;
506       PetscBLASInt   Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs;
507       PetscScalar   *qqc, *qqr, *TAU, *WORK;
508       PetscInt      *fids;
509       PetscReal     *data;
510 
511       /* count agg */
512       if (asz < minsz) minsz = asz;
513 
514       /* get block */
515       PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids));
516 
517       aggID = 0;
518       PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
519       while (pos) {
520         PetscInt gid1;
521 
522         PetscCall(PetscCDIntNdGetID(pos, &gid1));
523         PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
524 
525         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
526         else {
527           PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid));
528           PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table");
529         }
530         /* copy in B_i matrix - column-oriented */
531         data = &data_in[flid * bs];
532         for (ii = 0; ii < bs; ii++) {
533           for (jj = 0; jj < N; jj++) {
534             PetscReal d = data[jj * data_stride + ii];
535 
536             qqc[jj * Mdata + aggID * bs + ii] = d;
537           }
538         }
539         /* set fine IDs */
540         for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk;
541         aggID++;
542       }
543 
544       /* pad with zeros */
545       for (ii = asz * bs; ii < Mdata; ii++) {
546         for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0;
547       }
548 
549       /* QR */
550       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
551       PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
552       PetscCall(PetscFPTrapPop());
553       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error");
554       /* get R - column-oriented - output B_{i+1} */
555       {
556         PetscReal *data = &out_data[clid * nSAvec];
557 
558         for (jj = 0; jj < nSAvec; jj++) {
559           for (ii = 0; ii < nSAvec; ii++) {
560             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);
561             if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]);
562             else data[jj * out_data_stride + ii] = 0.;
563           }
564         }
565       }
566 
567       /* get Q - row-oriented */
568       PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
569       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO);
570 
571       for (ii = 0; ii < M; ii++) {
572         for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii];
573       }
574 
575       /* add diagonal block of P0 */
576       for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ }
577       PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES));
578       PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids));
579       clid++;
580     } /* coarse agg */
581   } /* for all fine nodes */
582   PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
583   PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
584   PetscCall(PCGAMGHashTableDestroy(&fgid_flid));
585   PetscFunctionReturn(PETSC_SUCCESS);
586 }
587 
588 static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer)
589 {
590   PC_MG       *mg          = (PC_MG *)pc->data;
591   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
592   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
593 
594   PetscFunctionBegin;
595   PetscCall(PetscViewerASCIIPrintf(viewer, "      AGG specific options\n"));
596   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number of levels of aggressive coarsening %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels));
597   if (pc_gamg_agg->aggressive_coarsening_levels > 0) {
598     PetscCall(PetscViewerASCIIPrintf(viewer, "        %s aggressive coarsening\n", !pc_gamg_agg->use_aggressive_square_graph ? "MIS-k" : "Square graph"));
599     if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(PetscViewerASCIIPrintf(viewer, "        MIS-%d coarsening on aggressive levels\n", (int)pc_gamg_agg->aggressive_mis_k));
600   }
601   PetscCall(PetscViewerASCIIPushTab(viewer));
602   PetscCall(PetscViewerASCIIPushTab(viewer));
603   PetscCall(PetscViewerASCIIPushTab(viewer));
604   PetscCall(PetscViewerASCIIPushTab(viewer));
605   PetscCall(MatCoarsenView(pc_gamg_agg->crs, viewer));
606   PetscCall(PetscViewerASCIIPopTab(viewer));
607   PetscCall(PetscViewerASCIIPopTab(viewer));
608   PetscCall(PetscViewerASCIIPopTab(viewer));
609   PetscCall(PetscViewerASCIIPopTab(viewer));
610   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number smoothing steps to construct prolongation %d\n", (int)pc_gamg_agg->nsmooths));
611   PetscFunctionReturn(PETSC_SUCCESS);
612 }
613 
614 static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat)
615 {
616   PC_MG          *mg          = (PC_MG *)pc->data;
617   PC_GAMG        *pc_gamg     = (PC_GAMG *)mg->innerctx;
618   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
619   const PetscReal vfilter     = pc_gamg->threshold[pc_gamg->current_level];
620   PetscBool       ishem, ismis;
621   const char     *prefix;
622   MatInfo         info0, info1;
623   PetscInt        bs;
624 
625   PetscFunctionBegin;
626   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
627   /* Note: depending on the algorithm that will be used for computing the coarse grid points this should pass PETSC_TRUE or PETSC_FALSE as the first argument */
628   /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */
629   PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs));
630   PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs));
631   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix));
632   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix));
633   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc_gamg_agg->crs, "pc_gamg_"));
634   PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs));
635   PetscCall(MatGetBlockSize(Amat, &bs));
636   // check for valid indices wrt bs
637   for (int ii = 0; ii < pc_gamg_agg->crs->strength_index_size; ii++) {
638     PetscCheck(pc_gamg_agg->crs->strength_index[ii] >= 0 && pc_gamg_agg->crs->strength_index[ii] < bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Indices (%d) must be non-negative and < block size (%d), NB, can not use -mat_coarsen_strength_index with -mat_coarsen_strength_index",
639                (int)pc_gamg_agg->crs->strength_index[ii], (int)bs);
640   }
641   PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem));
642   if (ishem) {
643     if (pc_gamg_agg->aggressive_coarsening_levels) PetscCall(PetscInfo(pc, "HEM and aggressive coarsening ignored: HEM using %d iterations\n", (int)pc_gamg_agg->crs->max_it));
644     pc_gamg_agg->aggressive_coarsening_levels = 0;                                         // aggressive and HEM does not make sense
645     PetscCall(MatCoarsenSetMaximumIterations(pc_gamg_agg->crs, pc_gamg_agg->crs->max_it)); // for code coverage
646     PetscCall(MatCoarsenSetThreshold(pc_gamg_agg->crs, vfilter));                          // for code coverage
647   } else {
648     PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENMIS, &ismis));
649     if (ismis && pc_gamg_agg->aggressive_coarsening_levels && !pc_gamg_agg->use_aggressive_square_graph) {
650       PetscCall(PetscInfo(pc, "MIS and aggressive coarsening and no square graph: force square graph\n"));
651       pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE;
652     }
653   }
654   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
655   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
656   PetscCall(MatGetInfo(Amat, MAT_LOCAL, &info0)); /* global reduction */
657 
658   if (ishem || pc_gamg_agg->use_low_mem_filter) {
659     PetscCall(MatCreateGraph(Amat, PETSC_TRUE, (vfilter >= 0 || ishem) ? PETSC_TRUE : PETSC_FALSE, vfilter, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat));
660   } else {
661     // make scalar graph, symetrize if not know to be symmetric, scale, but do not filter (expensive)
662     PetscCall(MatCreateGraph(Amat, PETSC_TRUE, PETSC_TRUE, -1, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat));
663     if (vfilter >= 0) {
664       PetscInt           Istart, Iend, ncols, nnz0, nnz1, NN, MM, nloc;
665       Mat                tGmat, Gmat = *a_Gmat;
666       MPI_Comm           comm;
667       const PetscScalar *vals;
668       const PetscInt    *idx;
669       PetscInt          *d_nnz, *o_nnz, kk, *garray = NULL, *AJ, maxcols = 0;
670       MatScalar         *AA; // this is checked in graph
671       PetscBool          isseqaij;
672       Mat                a, b, c;
673       MatType            jtype;
674 
675       PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm));
676       PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isseqaij));
677       PetscCall(MatGetType(Gmat, &jtype));
678       PetscCall(MatCreate(comm, &tGmat));
679       PetscCall(MatSetType(tGmat, jtype));
680 
681       /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold?
682         Also, if the matrix is symmetric, can we skip this
683         operation? It can be very expensive on large matrices. */
684 
685       // global sizes
686       PetscCall(MatGetSize(Gmat, &MM, &NN));
687       PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend));
688       nloc = Iend - Istart;
689       PetscCall(PetscMalloc2(nloc, &d_nnz, nloc, &o_nnz));
690       if (isseqaij) {
691         a = Gmat;
692         b = NULL;
693       } else {
694         Mat_MPIAIJ *d = (Mat_MPIAIJ *)Gmat->data;
695 
696         a      = d->A;
697         b      = d->B;
698         garray = d->garray;
699       }
700       /* Determine upper bound on non-zeros needed in new filtered matrix */
701       for (PetscInt row = 0; row < nloc; row++) {
702         PetscCall(MatGetRow(a, row, &ncols, NULL, NULL));
703         d_nnz[row] = ncols;
704         if (ncols > maxcols) maxcols = ncols;
705         PetscCall(MatRestoreRow(a, row, &ncols, NULL, NULL));
706       }
707       if (b) {
708         for (PetscInt row = 0; row < nloc; row++) {
709           PetscCall(MatGetRow(b, row, &ncols, NULL, NULL));
710           o_nnz[row] = ncols;
711           if (ncols > maxcols) maxcols = ncols;
712           PetscCall(MatRestoreRow(b, row, &ncols, NULL, NULL));
713         }
714       }
715       PetscCall(MatSetSizes(tGmat, nloc, nloc, MM, MM));
716       PetscCall(MatSetBlockSizes(tGmat, 1, 1));
717       PetscCall(MatSeqAIJSetPreallocation(tGmat, 0, d_nnz));
718       PetscCall(MatMPIAIJSetPreallocation(tGmat, 0, d_nnz, 0, o_nnz));
719       PetscCall(MatSetOption(tGmat, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
720       PetscCall(PetscFree2(d_nnz, o_nnz));
721       PetscCall(PetscMalloc2(maxcols, &AA, maxcols, &AJ));
722       nnz0 = nnz1 = 0;
723       for (c = a, kk = 0; c && kk < 2; c = b, kk++) {
724         for (PetscInt row = 0, grow = Istart, ncol_row, jj; row < nloc; row++, grow++) {
725           PetscCall(MatGetRow(c, row, &ncols, &idx, &vals));
726           for (ncol_row = jj = 0; jj < ncols; jj++, nnz0++) {
727             PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
728             if (PetscRealPart(sv) > vfilter) {
729               PetscInt cid = idx[jj] + Istart; //diag
730 
731               nnz1++;
732               if (c != a) cid = garray[idx[jj]];
733               AA[ncol_row] = vals[jj];
734               AJ[ncol_row] = cid;
735               ncol_row++;
736             }
737           }
738           PetscCall(MatRestoreRow(c, row, &ncols, &idx, &vals));
739           PetscCall(MatSetValues(tGmat, 1, &grow, ncol_row, AJ, AA, INSERT_VALUES));
740         }
741       }
742       PetscCall(PetscFree2(AA, AJ));
743       PetscCall(MatAssemblyBegin(tGmat, MAT_FINAL_ASSEMBLY));
744       PetscCall(MatAssemblyEnd(tGmat, MAT_FINAL_ASSEMBLY));
745       PetscCall(MatPropagateSymmetryOptions(Gmat, tGmat)); /* Normal Mat options are not relevant ? */
746       PetscCall(PetscInfo(pc, "\t %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%" PetscInt_FMT ", max row size %" PetscInt_FMT "\n", (!nnz0) ? 1. : 100. * (double)nnz1 / (double)nnz0, (double)vfilter, (!nloc) ? 1. : (double)nnz0 / (double)nloc, MM, maxcols));
747       PetscCall(MatViewFromOptions(tGmat, NULL, "-mat_filter_graph_view"));
748       PetscCall(MatDestroy(&Gmat));
749       *a_Gmat = tGmat;
750     }
751   }
752 
753   PetscCall(MatGetInfo(*a_Gmat, MAT_LOCAL, &info1)); /* global reduction */
754   if (info0.nz_used > 0) PetscCall(PetscInfo(pc, "Filtering left %g %% edges in graph (%e %e)\n", 100.0 * info1.nz_used * (double)(bs * bs) / info0.nz_used, info0.nz_used, info1.nz_used));
755   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
756   PetscFunctionReturn(PETSC_SUCCESS);
757 }
758 
759 typedef PetscInt    NState;
760 static const NState NOT_DONE = -2;
761 static const NState DELETED  = -1;
762 static const NState REMOVED  = -3;
763 #define IS_SELECTED(s) (s != DELETED && s != NOT_DONE && s != REMOVED)
764 
765 /*
766    fixAggregatesWithSquare - greedy grab of with G1 (unsquared graph) -- AIJ specific -- change to fixAggregatesWithSquare -- TODD
767      - AGG-MG specific: clears singletons out of 'selected_2'
768 
769    Input Parameter:
770    . Gmat_2 - global matrix of squared graph (data not defined)
771    . Gmat_1 - base graph to grab with base graph
772    Input/Output Parameter:
773    . aggs_2 - linked list of aggs with gids)
774 */
775 static PetscErrorCode fixAggregatesWithSquare(PC pc, Mat Gmat_2, Mat Gmat_1, PetscCoarsenData *aggs_2)
776 {
777   PetscBool      isMPI;
778   Mat_SeqAIJ    *matA_1, *matB_1 = NULL;
779   MPI_Comm       comm;
780   PetscInt       lid, *ii, *idx, ix, Iend, my0, kk, n, j;
781   Mat_MPIAIJ    *mpimat_2 = NULL, *mpimat_1 = NULL;
782   const PetscInt nloc = Gmat_2->rmap->n;
783   PetscScalar   *cpcol_1_state, *cpcol_2_state, *cpcol_2_par_orig, *lid_parent_gid;
784   PetscInt      *lid_cprowID_1 = NULL;
785   NState        *lid_state;
786   Vec            ghost_par_orig2;
787   PetscMPIInt    rank;
788 
789   PetscFunctionBegin;
790   PetscCall(PetscObjectGetComm((PetscObject)Gmat_2, &comm));
791   PetscCallMPI(MPI_Comm_rank(comm, &rank));
792   PetscCall(MatGetOwnershipRange(Gmat_1, &my0, &Iend));
793 
794   /* get submatrices */
795   PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATMPIAIJ, &isMPI));
796   PetscCall(PetscInfo(pc, "isMPI = %s\n", isMPI ? "yes" : "no"));
797   PetscCall(PetscMalloc3(nloc, &lid_state, nloc, &lid_parent_gid, nloc, &lid_cprowID_1));
798   for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
799   if (isMPI) {
800     /* grab matrix objects */
801     mpimat_2 = (Mat_MPIAIJ *)Gmat_2->data;
802     mpimat_1 = (Mat_MPIAIJ *)Gmat_1->data;
803     matA_1   = (Mat_SeqAIJ *)mpimat_1->A->data;
804     matB_1   = (Mat_SeqAIJ *)mpimat_1->B->data;
805 
806     /* force compressed row storage for B matrix in AuxMat */
807     PetscCall(MatCheckCompressedRow(mpimat_1->B, matB_1->nonzerorowcnt, &matB_1->compressedrow, matB_1->i, Gmat_1->rmap->n, -1.0));
808     for (ix = 0; ix < matB_1->compressedrow.nrows; ix++) {
809       PetscInt lid = matB_1->compressedrow.rindex[ix];
810 
811       PetscCheck(lid <= nloc && lid >= -1, PETSC_COMM_SELF, PETSC_ERR_USER, "lid %d out of range. nloc = %d", (int)lid, (int)nloc);
812       if (lid != -1) lid_cprowID_1[lid] = ix;
813     }
814   } else {
815     PetscBool isAIJ;
816 
817     PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATSEQAIJ, &isAIJ));
818     PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix.");
819     matA_1 = (Mat_SeqAIJ *)Gmat_1->data;
820   }
821   if (nloc > 0) { PetscCheck(!matB_1 || matB_1->compressedrow.use, PETSC_COMM_SELF, PETSC_ERR_PLIB, "matB_1 && !matB_1->compressedrow.use: PETSc bug???"); }
822   /* get state of locals and selected gid for deleted */
823   for (lid = 0; lid < nloc; lid++) {
824     lid_parent_gid[lid] = -1.0;
825     lid_state[lid]      = DELETED;
826   }
827 
828   /* set lid_state */
829   for (lid = 0; lid < nloc; lid++) {
830     PetscCDIntNd *pos;
831 
832     PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
833     if (pos) {
834       PetscInt gid1;
835 
836       PetscCall(PetscCDIntNdGetID(pos, &gid1));
837       PetscCheck(gid1 == lid + my0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gid1 %d != lid %d + my0 %d", (int)gid1, (int)lid, (int)my0);
838       lid_state[lid] = gid1;
839     }
840   }
841 
842   /* map local to selected local, DELETED means a ghost owns it */
843   for (lid = kk = 0; lid < nloc; lid++) {
844     NState state = lid_state[lid];
845     if (IS_SELECTED(state)) {
846       PetscCDIntNd *pos;
847 
848       PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
849       while (pos) {
850         PetscInt gid1;
851 
852         PetscCall(PetscCDIntNdGetID(pos, &gid1));
853         PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos));
854         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1 - my0] = (PetscScalar)(lid + my0);
855       }
856     }
857   }
858   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
859   if (isMPI) {
860     Vec tempVec;
861     /* get 'cpcol_1_state' */
862     PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL));
863     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
864       PetscScalar v = (PetscScalar)lid_state[kk];
865 
866       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
867     }
868     PetscCall(VecAssemblyBegin(tempVec));
869     PetscCall(VecAssemblyEnd(tempVec));
870     PetscCall(VecScatterBegin(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD));
871     PetscCall(VecScatterEnd(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD));
872     PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state));
873     /* get 'cpcol_2_state' */
874     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD));
875     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD));
876     PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state));
877     /* get 'cpcol_2_par_orig' */
878     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
879       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
880 
881       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
882     }
883     PetscCall(VecAssemblyBegin(tempVec));
884     PetscCall(VecAssemblyEnd(tempVec));
885     PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2));
886     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD));
887     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD));
888     PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig));
889 
890     PetscCall(VecDestroy(&tempVec));
891   } /* ismpi */
892   for (lid = 0; lid < nloc; lid++) {
893     NState state = lid_state[lid];
894 
895     if (IS_SELECTED(state)) {
896       /* steal locals */
897       ii  = matA_1->i;
898       n   = ii[lid + 1] - ii[lid];
899       idx = matA_1->j + ii[lid];
900       for (j = 0; j < n; j++) {
901         PetscInt lidj   = idx[j], sgid;
902         NState   statej = lid_state[lidj];
903 
904         if (statej == DELETED && (sgid = (PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid + my0) { /* steal local */
905           lid_parent_gid[lidj] = (PetscScalar)(lid + my0);                                              /* send this if sgid is not local */
906           if (sgid >= my0 && sgid < Iend) {                                                             /* I'm stealing this local from a local sgid */
907             PetscInt      hav = 0, slid = sgid - my0, gidj = lidj + my0;
908             PetscCDIntNd *pos, *last = NULL;
909 
910             /* looking for local from local so id_llist_2 works */
911             PetscCall(PetscCDGetHeadPos(aggs_2, slid, &pos));
912             while (pos) {
913               PetscInt gid;
914               PetscCall(PetscCDIntNdGetID(pos, &gid));
915               if (gid == gidj) {
916                 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null");
917                 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last));
918                 PetscCall(PetscCDAppendNode(aggs_2, lid, pos));
919                 hav = 1;
920                 break;
921               } else last = pos;
922               PetscCall(PetscCDGetNextPos(aggs_2, slid, &pos));
923             }
924             if (hav != 1) {
925               PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find adj in 'selected' lists - structurally unsymmetric matrix");
926               SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav);
927             }
928           } else { /* I'm stealing this local, owned by a ghost */
929             PetscCheck(sgid == -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mat has an un-symmetric graph. Use '-%spc_gamg_sym_graph true' to symmetrize the graph or '-%spc_gamg_threshold -1' if the matrix is structurally symmetric.",
930                        ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "");
931             PetscCall(PetscCDAppendID(aggs_2, lid, lidj + my0));
932           }
933         }
934       } /* local neighbors */
935     } else if (state == DELETED /* && lid_cprowID_1 */) {
936       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
937 
938       /* see if I have a selected ghost neighbor that will steal me */
939       if ((ix = lid_cprowID_1[lid]) != -1) {
940         ii  = matB_1->compressedrow.i;
941         n   = ii[ix + 1] - ii[ix];
942         idx = matB_1->j + ii[ix];
943         for (j = 0; j < n; j++) {
944           PetscInt cpid   = idx[j];
945           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
946 
947           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
948             lid_parent_gid[lid] = (PetscScalar)statej;              /* send who selected */
949             if (sgidold >= my0 && sgidold < Iend) {                 /* this was mine */
950               PetscInt      hav = 0, oldslidj = sgidold - my0;
951               PetscCDIntNd *pos, *last        = NULL;
952 
953               /* remove from 'oldslidj' list */
954               PetscCall(PetscCDGetHeadPos(aggs_2, oldslidj, &pos));
955               while (pos) {
956                 PetscInt gid;
957 
958                 PetscCall(PetscCDIntNdGetID(pos, &gid));
959                 if (lid + my0 == gid) {
960                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
961                   PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null");
962                   PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last));
963                   /* ghost (PetscScalar)statej will add this later */
964                   hav = 1;
965                   break;
966                 } else last = pos;
967                 PetscCall(PetscCDGetNextPos(aggs_2, oldslidj, &pos));
968               }
969               if (hav != 1) {
970                 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find (hav=%d) adj in 'selected' lists - structurally unsymmetric matrix", (int)hav);
971                 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav);
972               }
973             } else {
974               /* TODO: ghosts remove this later */
975             }
976           }
977         }
978       }
979     } /* selected/deleted */
980   } /* node loop */
981 
982   if (isMPI) {
983     PetscScalar    *cpcol_2_parent, *cpcol_2_gid;
984     Vec             tempVec, ghostgids2, ghostparents2;
985     PetscInt        cpid, nghost_2;
986     PCGAMGHashTable gid_cpid;
987 
988     PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2));
989     PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL));
990 
991     /* get 'cpcol_2_parent' */
992     for (kk = 0, j = my0; kk < nloc; kk++, j++) { PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES)); }
993     PetscCall(VecAssemblyBegin(tempVec));
994     PetscCall(VecAssemblyEnd(tempVec));
995     PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2));
996     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD));
997     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD));
998     PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent));
999 
1000     /* get 'cpcol_2_gid' */
1001     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
1002       PetscScalar v = (PetscScalar)j;
1003 
1004       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
1005     }
1006     PetscCall(VecAssemblyBegin(tempVec));
1007     PetscCall(VecAssemblyEnd(tempVec));
1008     PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2));
1009     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD));
1010     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD));
1011     PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid));
1012     PetscCall(VecDestroy(&tempVec));
1013 
1014     /* look for deleted ghosts and add to table */
1015     PetscCall(PCGAMGHashTableCreate(2 * nghost_2 + 1, &gid_cpid));
1016     for (cpid = 0; cpid < nghost_2; cpid++) {
1017       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
1018 
1019       if (state == DELETED) {
1020         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
1021         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
1022 
1023         if (sgid_old == -1 && sgid_new != -1) {
1024           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
1025 
1026           PetscCall(PCGAMGHashTableAdd(&gid_cpid, gid, cpid));
1027         }
1028       }
1029     }
1030 
1031     /* look for deleted ghosts and see if they moved - remove it */
1032     for (lid = 0; lid < nloc; lid++) {
1033       NState state = lid_state[lid];
1034       if (IS_SELECTED(state)) {
1035         PetscCDIntNd *pos, *last = NULL;
1036         /* look for deleted ghosts and see if they moved */
1037         PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
1038         while (pos) {
1039           PetscInt gid;
1040           PetscCall(PetscCDIntNdGetID(pos, &gid));
1041 
1042           if (gid < my0 || gid >= Iend) {
1043             PetscCall(PCGAMGHashTableFind(&gid_cpid, gid, &cpid));
1044             if (cpid != -1) {
1045               /* a moved ghost - */
1046               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
1047               PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last));
1048             } else last = pos;
1049           } else last = pos;
1050 
1051           PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos));
1052         } /* loop over list of deleted */
1053       } /* selected */
1054     }
1055     PetscCall(PCGAMGHashTableDestroy(&gid_cpid));
1056 
1057     /* look at ghosts, see if they changed - and it */
1058     for (cpid = 0; cpid < nghost_2; cpid++) {
1059       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
1060 
1061       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
1062         PetscInt      gid      = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
1063         PetscInt      slid_new = sgid_new - my0, hav = 0;
1064         PetscCDIntNd *pos;
1065 
1066         /* search for this gid to see if I have it */
1067         PetscCall(PetscCDGetHeadPos(aggs_2, slid_new, &pos));
1068         while (pos) {
1069           PetscInt gidj;
1070 
1071           PetscCall(PetscCDIntNdGetID(pos, &gidj));
1072           PetscCall(PetscCDGetNextPos(aggs_2, slid_new, &pos));
1073 
1074           if (gidj == gid) {
1075             hav = 1;
1076             break;
1077           }
1078         }
1079         if (hav != 1) {
1080           /* insert 'flidj' into head of llist */
1081           PetscCall(PetscCDAppendID(aggs_2, slid_new, gid));
1082         }
1083       }
1084     }
1085     PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state));
1086     PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state));
1087     PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent));
1088     PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid));
1089     PetscCall(VecDestroy(&ghostgids2));
1090     PetscCall(VecDestroy(&ghostparents2));
1091     PetscCall(VecDestroy(&ghost_par_orig2));
1092   }
1093   PetscCall(PetscFree3(lid_state, lid_parent_gid, lid_cprowID_1));
1094   PetscFunctionReturn(PETSC_SUCCESS);
1095 }
1096 
1097 /*
1098    PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for
1099      communication of QR data used with HEM and MISk coarsening
1100 
1101   Input Parameter:
1102    . a_pc - this
1103 
1104   Input/Output Parameter:
1105    . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out)
1106 
1107   Output Parameter:
1108    . agg_lists - list of aggregates
1109 
1110 */
1111 static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists)
1112 {
1113   PC_MG       *mg          = (PC_MG *)a_pc->data;
1114   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
1115   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
1116   Mat          Gmat2, Gmat1 = *a_Gmat1; /* aggressive graph */
1117   IS           perm;
1118   PetscInt     Istart, Iend, Ii, nloc, bs, nn;
1119   PetscInt    *permute, *degree;
1120   PetscBool   *bIndexSet;
1121   PetscReal    hashfact;
1122   PetscInt     iSwapIndex;
1123   PetscRandom  random;
1124   MPI_Comm     comm;
1125 
1126   PetscFunctionBegin;
1127   PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm));
1128   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
1129   PetscCall(MatGetLocalSize(Gmat1, &nn, NULL));
1130   PetscCall(MatGetBlockSize(Gmat1, &bs));
1131   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs);
1132   nloc = nn / bs;
1133   /* get MIS aggs - randomize */
1134   PetscCall(PetscMalloc2(nloc, &permute, nloc, &degree));
1135   PetscCall(PetscCalloc1(nloc, &bIndexSet));
1136   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
1137   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random));
1138   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
1139   for (Ii = 0; Ii < nloc; Ii++) {
1140     PetscInt nc;
1141 
1142     PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
1143     degree[Ii] = nc;
1144     PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
1145   }
1146   for (Ii = 0; Ii < nloc; Ii++) {
1147     PetscCall(PetscRandomGetValueReal(random, &hashfact));
1148     iSwapIndex = (PetscInt)(hashfact * nloc) % nloc;
1149     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1150       PetscInt iTemp = permute[iSwapIndex];
1151 
1152       permute[iSwapIndex]   = permute[Ii];
1153       permute[Ii]           = iTemp;
1154       iTemp                 = degree[iSwapIndex];
1155       degree[iSwapIndex]    = degree[Ii];
1156       degree[Ii]            = iTemp;
1157       bIndexSet[iSwapIndex] = PETSC_TRUE;
1158     }
1159   }
1160   // apply minimum degree ordering -- NEW
1161   if (pc_gamg_agg->use_minimum_degree_ordering) { PetscCall(PetscSortIntWithArray(nloc, degree, permute)); }
1162   PetscCall(PetscFree(bIndexSet));
1163   PetscCall(PetscRandomDestroy(&random));
1164   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
1165   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
1166   // square graph
1167   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels && pc_gamg_agg->use_aggressive_square_graph) {
1168     PetscCall(PCGAMGSquareGraph_GAMG(a_pc, Gmat1, &Gmat2));
1169   } else Gmat2 = Gmat1;
1170   // switch to old MIS-1 for square graph
1171   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) {
1172     if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, pc_gamg_agg->aggressive_mis_k)); // hardwire to MIS-2
1173     else PetscCall(MatCoarsenSetType(pc_gamg_agg->crs, MATCOARSENMIS));                                                                   // old MIS -- side effect
1174   } else if (pc_gamg_agg->use_aggressive_square_graph && pc_gamg_agg->aggressive_coarsening_levels > 0) {                                 // we reset the MIS
1175     const char *prefix;
1176 
1177     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)a_pc, &prefix));
1178     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix));
1179     PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); // get the default back on non-aggressive levels when square graph switched to old MIS
1180   }
1181   PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat2));
1182   PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE));
1183   PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm));
1184   PetscCall(MatCoarsenApply(pc_gamg_agg->crs));
1185   PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */
1186 
1187   PetscCall(ISDestroy(&perm));
1188   PetscCall(PetscFree2(permute, degree));
1189   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
1190 
1191   if (Gmat2 != Gmat1) { // square graph, we need ghosts for selected
1192     PetscCoarsenData *llist = *agg_lists;
1193 
1194     PetscCall(fixAggregatesWithSquare(a_pc, Gmat2, Gmat1, *agg_lists));
1195     PetscCall(MatDestroy(&Gmat1));
1196     *a_Gmat1 = Gmat2;                          /* output */
1197     PetscCall(PetscCDSetMat(llist, *a_Gmat1)); /* Need a graph with ghosts here */
1198   }
1199   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
1200   PetscFunctionReturn(PETSC_SUCCESS);
1201 }
1202 
1203 /*
1204  PCGAMGConstructProlongator_AGG
1205 
1206  Input Parameter:
1207  . pc - this
1208  . Amat - matrix on this fine level
1209  . Graph - used to get ghost data for nodes in
1210  . agg_lists - list of aggregates
1211  Output Parameter:
1212  . a_P_out - prolongation operator to the next level
1213  */
1214 static PetscErrorCode PCGAMGConstructProlongator_AGG(PC pc, Mat Amat, PetscCoarsenData *agg_lists, Mat *a_P_out)
1215 {
1216   PC_MG         *mg      = (PC_MG *)pc->data;
1217   PC_GAMG       *pc_gamg = (PC_GAMG *)mg->innerctx;
1218   const PetscInt col_bs  = pc_gamg->data_cell_cols;
1219   PetscInt       Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs;
1220   Mat            Gmat, Prol;
1221   PetscMPIInt    size;
1222   MPI_Comm       comm;
1223   PetscReal     *data_w_ghost;
1224   PetscInt       myCrs0, nbnodes = 0, *flid_fgid;
1225   MatType        mtype;
1226 
1227   PetscFunctionBegin;
1228   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
1229   PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1");
1230   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
1231   PetscCallMPI(MPI_Comm_size(comm, &size));
1232   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
1233   PetscCall(MatGetBlockSize(Amat, &bs));
1234   nloc = (Iend - Istart) / bs;
1235   my0  = Istart / bs;
1236   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);
1237   PetscCall(PetscCDGetMat(agg_lists, &Gmat)); // get auxiliary matrix for ghost edges for size > 1
1238 
1239   /* get 'nLocalSelected' */
1240   for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) {
1241     PetscBool ise;
1242 
1243     /* filter out singletons 0 or 1? */
1244     PetscCall(PetscCDIsEmptyAt(agg_lists, ii, &ise));
1245     if (!ise) nLocalSelected++;
1246   }
1247 
1248   /* create prolongator, create P matrix */
1249   PetscCall(MatGetType(Amat, &mtype));
1250   PetscCall(MatCreate(comm, &Prol));
1251   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE));
1252   PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); // should this be before MatSetSizes?
1253   PetscCall(MatSetType(Prol, mtype));
1254 #if PetscDefined(HAVE_DEVICE)
1255   PetscBool flg;
1256   PetscCall(MatBoundToCPU(Amat, &flg));
1257   PetscCall(MatBindToCPU(Prol, flg));
1258   if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE));
1259 #endif
1260   PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL));
1261   PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL));
1262 
1263   /* can get all points "removed" */
1264   PetscCall(MatGetSize(Prol, &kk, &ii));
1265   if (!ii) {
1266     PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix));
1267     PetscCall(MatDestroy(&Prol));
1268     *a_P_out = NULL; /* out */
1269     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
1270     PetscFunctionReturn(PETSC_SUCCESS);
1271   }
1272   PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs));
1273   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
1274 
1275   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);
1276   myCrs0 = myCrs0 / col_bs;
1277   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);
1278 
1279   /* create global vector of data in 'data_w_ghost' */
1280   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
1281   if (size > 1) { /* get ghost null space data */
1282     PetscReal *tmp_gdata, *tmp_ldata, *tp2;
1283 
1284     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
1285     for (jj = 0; jj < col_bs; jj++) {
1286       for (kk = 0; kk < bs; kk++) {
1287         PetscInt         ii, stride;
1288         const PetscReal *tp = PetscSafePointerPlusOffset(pc_gamg->data, jj * bs * nloc + kk);
1289 
1290         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
1291 
1292         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
1293 
1294         if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */
1295           PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost));
1296           nbnodes = bs * stride;
1297         }
1298         tp2 = PetscSafePointerPlusOffset(data_w_ghost, jj * bs * stride + kk);
1299         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
1300         PetscCall(PetscFree(tmp_gdata));
1301       }
1302     }
1303     PetscCall(PetscFree(tmp_ldata));
1304   } else {
1305     nbnodes      = bs * nloc;
1306     data_w_ghost = (PetscReal *)pc_gamg->data;
1307   }
1308 
1309   /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */
1310   if (size > 1) {
1311     PetscReal *fid_glid_loc, *fiddata;
1312     PetscInt   stride;
1313 
1314     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
1315     for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk);
1316     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
1317     PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */
1318     for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
1319     PetscCall(PetscFree(fiddata));
1320 
1321     PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs);
1322     PetscCall(PetscFree(fid_glid_loc));
1323   } else {
1324     PetscCall(PetscMalloc1(nloc, &flid_fgid));
1325     for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk;
1326   }
1327   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
1328   /* get P0 */
1329   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
1330   {
1331     PetscReal *data_out = NULL;
1332 
1333     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol));
1334     PetscCall(PetscFree(pc_gamg->data));
1335 
1336     pc_gamg->data           = data_out;
1337     pc_gamg->data_cell_rows = col_bs;
1338     pc_gamg->data_sz        = col_bs * col_bs * nLocalSelected;
1339   }
1340   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
1341   if (size > 1) PetscCall(PetscFree(data_w_ghost));
1342   PetscCall(PetscFree(flid_fgid));
1343 
1344   *a_P_out = Prol; /* out */
1345   PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_initial_prolongation"));
1346 
1347   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
1348   PetscFunctionReturn(PETSC_SUCCESS);
1349 }
1350 
1351 /*
1352    PCGAMGOptimizeProlongator_AGG - given the initial prolongator optimizes it by smoothed aggregation pc_gamg_agg->nsmooths times
1353 
1354   Input Parameter:
1355    . pc - this
1356    . Amat - matrix on this fine level
1357  In/Output Parameter:
1358    . a_P - prolongation operator to the next level
1359 */
1360 static PetscErrorCode PCGAMGOptimizeProlongator_AGG(PC pc, Mat Amat, Mat *a_P)
1361 {
1362   PC_MG       *mg          = (PC_MG *)pc->data;
1363   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
1364   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
1365   PetscInt     jj;
1366   Mat          Prol = *a_P;
1367   MPI_Comm     comm;
1368   KSP          eksp;
1369   Vec          bb, xx;
1370   PC           epc;
1371   PetscReal    alpha, emax, emin;
1372 
1373   PetscFunctionBegin;
1374   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
1375   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
1376 
1377   /* compute maximum singular value of operator to be used in smoother */
1378   if (0 < pc_gamg_agg->nsmooths) {
1379     /* get eigen estimates */
1380     if (pc_gamg->emax > 0) {
1381       emin = pc_gamg->emin;
1382       emax = pc_gamg->emax;
1383     } else {
1384       const char *prefix;
1385 
1386       PetscCall(MatCreateVecs(Amat, &bb, NULL));
1387       PetscCall(MatCreateVecs(Amat, &xx, NULL));
1388       PetscCall(KSPSetNoisy_Private(bb));
1389 
1390       PetscCall(KSPCreate(comm, &eksp));
1391       PetscCall(KSPSetNestLevel(eksp, pc->kspnestlevel));
1392       PetscCall(PCGetOptionsPrefix(pc, &prefix));
1393       PetscCall(KSPSetOptionsPrefix(eksp, prefix));
1394       PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_"));
1395       {
1396         PetscBool isset, sflg;
1397 
1398         PetscCall(MatIsSPDKnown(Amat, &isset, &sflg));
1399         if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG));
1400       }
1401       PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure));
1402       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
1403 
1404       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
1405       PetscCall(KSPSetOperators(eksp, Amat, Amat));
1406 
1407       PetscCall(KSPGetPC(eksp, &epc));
1408       PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */
1409 
1410       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
1411 
1412       PetscCall(KSPSetFromOptions(eksp));
1413       PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE));
1414       PetscCall(KSPSolve(eksp, bb, xx));
1415       PetscCall(KSPCheckSolve(eksp, pc, xx));
1416 
1417       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
1418       PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI));
1419       PetscCall(VecDestroy(&xx));
1420       PetscCall(VecDestroy(&bb));
1421       PetscCall(KSPDestroy(&eksp));
1422     }
1423     if (pc_gamg->use_sa_esteig) {
1424       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
1425       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
1426       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));
1427     } else {
1428       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
1429       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
1430     }
1431   } else {
1432     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
1433     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
1434   }
1435 
1436   /* smooth P0 */
1437   if (pc_gamg_agg->nsmooths > 0) {
1438     Vec diag;
1439 
1440     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
1441     PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero");
1442 
1443     PetscCall(MatCreateVecs(Amat, &diag, NULL));
1444     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
1445     PetscCall(VecReciprocal(diag));
1446 
1447     for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
1448       Mat tMat;
1449 
1450       PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
1451       /*
1452         Smooth aggregation on the prolongator
1453 
1454         P_{i} := (I - 1.4/emax D^{-1}A) P_i\{i-1}
1455       */
1456       PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
1457       PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat));
1458       PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
1459       PetscCall(MatProductClear(tMat));
1460       PetscCall(MatDiagonalScale(tMat, diag, NULL));
1461 
1462       /* TODO: Document the 1.4 and don't hardwire it in this routine */
1463       alpha = -1.4 / emax;
1464       PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
1465       PetscCall(MatDestroy(&Prol));
1466       Prol = tMat;
1467       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
1468     }
1469     PetscCall(VecDestroy(&diag));
1470   }
1471   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
1472   PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_prolongation"));
1473   *a_P = Prol;
1474   PetscFunctionReturn(PETSC_SUCCESS);
1475 }
1476 
1477 /*MC
1478   PCGAMGAGG - Smooth aggregation, {cite}`vanek1996algebraic`, {cite}`vanek2001convergence`, variant of PETSc's algebraic multigrid (`PCGAMG`) preconditioner
1479 
1480   Options Database Keys:
1481 + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation to construct prolongation
1482 . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest.
1483 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening
1484 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm
1485 . -pc_gamg_pc_gamg_asm_hem_aggs <n,default=0> - Number of HEM aggregation steps for ASM smoother
1486 - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive')
1487 
1488   Level: intermediate
1489 
1490   Notes:
1491   To obtain good performance for `PCGAMG` for vector valued problems you must
1492   call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point.
1493   Call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator
1494 
1495   The many options for `PCMG` and `PCGAMG` such as controlling the smoothers on each level etc. also work for `PCGAMGAGG`
1496 
1497 .seealso: `PCGAMG`, [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCCreate()`, `PCSetType()`,
1498           `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`,
1499           `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`,
1500           `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`,
1501           `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()`
1502 M*/
1503 PetscErrorCode PCCreateGAMG_AGG(PC pc)
1504 {
1505   PC_MG       *mg      = (PC_MG *)pc->data;
1506   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
1507   PC_GAMG_AGG *pc_gamg_agg;
1508 
1509   PetscFunctionBegin;
1510   /* create sub context for SA */
1511   PetscCall(PetscNew(&pc_gamg_agg));
1512   pc_gamg->subctx = pc_gamg_agg;
1513 
1514   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
1515   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1516   /* reset does not do anything; setup not virtual */
1517 
1518   /* set internal function pointers */
1519   pc_gamg->ops->creategraph       = PCGAMGCreateGraph_AGG;
1520   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
1521   pc_gamg->ops->prolongator       = PCGAMGConstructProlongator_AGG;
1522   pc_gamg->ops->optprolongator    = PCGAMGOptimizeProlongator_AGG;
1523   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
1524   pc_gamg->ops->view              = PCView_GAMG_AGG;
1525 
1526   pc_gamg_agg->nsmooths                     = 1;
1527   pc_gamg_agg->aggressive_coarsening_levels = 1;
1528   pc_gamg_agg->use_aggressive_square_graph  = PETSC_TRUE;
1529   pc_gamg_agg->use_minimum_degree_ordering  = PETSC_FALSE;
1530   pc_gamg_agg->use_low_mem_filter           = PETSC_FALSE;
1531   pc_gamg_agg->aggressive_mis_k             = 2;
1532 
1533   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG));
1534   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG));
1535   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", PCGAMGSetAggressiveSquareGraph_AGG));
1536   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", PCGAMGMISkSetMinDegreeOrdering_AGG));
1537   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", PCGAMGSetLowMemoryFilter_AGG));
1538   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", PCGAMGMISkSetAggressive_AGG));
1539   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG));
1540   PetscFunctionReturn(PETSC_SUCCESS);
1541 }
1542