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