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