xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision aea0ef1414093d6846bfe8b0234593cb4cd71dcc)
12e68589bSMark F. Adams /*
2b817416eSBarry Smith  GAMG geometric-algebric multigrid PC - Mark Adams 2011
32e68589bSMark F. Adams  */
42e68589bSMark F. Adams 
52e68589bSMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/
62e68589bSMark F. Adams #include <petscblaslapack.h>
7539c167fSBarry Smith #include <petscdm.h>
873f7197eSJed Brown #include <petsc/private/kspimpl.h>
92e68589bSMark F. Adams 
102e68589bSMark F. Adams typedef struct {
11a077d33dSBarry Smith   PetscInt   nsmooths;                     // number of smoothing steps to construct prolongation
12bae903cbSmarkadams4   PetscInt   aggressive_coarsening_levels; // number of aggressive coarsening levels (square or MISk)
13d529f056Smarkadams4   PetscInt   aggressive_mis_k;             // the k in MIS-k
14d529f056Smarkadams4   PetscBool  use_aggressive_square_graph;
15d529f056Smarkadams4   PetscBool  use_minimum_degree_ordering;
16a9ccf005SMark Adams   PetscBool  use_low_mem_filter;
176c34c54dSStefano Zampini   PetscBool  graph_symmetrize;
182d776b49SBarry Smith   MatCoarsen crs;
192e68589bSMark F. Adams } PC_GAMG_AGG;
202e68589bSMark F. Adams 
212e68589bSMark F. Adams /*@
22a077d33dSBarry Smith   PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used to construct the prolongation operator
232e68589bSMark F. Adams 
24c3339decSBarry Smith   Logically Collective
252e68589bSMark F. Adams 
262e68589bSMark F. Adams   Input Parameters:
2720f4b53cSBarry Smith + pc - the preconditioner context
2820f4b53cSBarry Smith - n  - the number of smooths
292e68589bSMark F. Adams 
302e68589bSMark F. Adams   Options Database Key:
31*aea0ef14SMark Adams . -pc_gamg_agg_nsmooths <nsmooth, default=1> - the flag
322e68589bSMark F. Adams 
332e68589bSMark F. Adams   Level: intermediate
342e68589bSMark F. Adams 
35a077d33dSBarry Smith   Note:
36a077d33dSBarry Smith   This is a different concept from the number smoothing steps used during the linear solution process which
37a077d33dSBarry Smith   can be set with `-mg_levels_ksp_max_it`
38a077d33dSBarry Smith 
39a077d33dSBarry Smith   Developer Note:
40a077d33dSBarry Smith   This should be named `PCGAMGAGGSetNSmooths()`.
41a077d33dSBarry Smith 
42a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCMG`, `PCGAMG`
432e68589bSMark F. Adams @*/
44d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
45d71ae5a4SJacob Faibussowitsch {
462e68589bSMark F. Adams   PetscFunctionBegin;
472e68589bSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
48d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
49cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n));
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
512e68589bSMark F. Adams }
522e68589bSMark F. Adams 
53d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n)
54d71ae5a4SJacob Faibussowitsch {
552e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
562e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
57c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
582e68589bSMark F. Adams 
592e68589bSMark F. Adams   PetscFunctionBegin;
60c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62c8b0795cSMark F. Adams }
63c8b0795cSMark F. Adams 
64c8b0795cSMark F. Adams /*@
65f1580f4eSBarry Smith   PCGAMGSetAggressiveLevels -  Use aggressive coarsening on first n levels
66ef4ad70eSMark F. Adams 
67c3339decSBarry Smith   Logically Collective
68ef4ad70eSMark F. Adams 
69ef4ad70eSMark F. Adams   Input Parameters:
70cab9ed1eSBarry Smith + pc - the preconditioner context
71d5d25401SBarry Smith - n  - 0, 1 or more
72ef4ad70eSMark F. Adams 
73ef4ad70eSMark F. Adams   Options Database Key:
74*aea0ef14SMark Adams . -pc_gamg_aggressive_coarsening <n,default = 1> - the flag
75af3c827dSMark Adams 
76ef4ad70eSMark F. Adams   Level: intermediate
77ef4ad70eSMark F. Adams 
78*aea0ef14SMark Adams   Note:
79*aea0ef14SMark Adams   By default, aggressive coarsening squares the matrix (computes $ A^T A$) before coarsening. Calling `PCGAMGSetAggressiveSquareGraph()` with a value of `PETSC_FALSE` changes the aggressive coarsening strategy to use MIS-k, see `PCGAMGMISkSetAggressive()`.
80*aea0ef14SMark Adams 
81a077d33dSBarry Smith .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()`
82ef4ad70eSMark F. Adams @*/
83d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n)
84d71ae5a4SJacob Faibussowitsch {
85ef4ad70eSMark F. Adams   PetscFunctionBegin;
86ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
87d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
88bae903cbSmarkadams4   PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n));
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90ef4ad70eSMark F. Adams }
91ef4ad70eSMark F. Adams 
92d529f056Smarkadams4 /*@
93d529f056Smarkadams4   PCGAMGMISkSetAggressive - Number (k) distance in MIS coarsening (>2 is 'aggressive')
94d529f056Smarkadams4 
95d529f056Smarkadams4   Logically Collective
96d529f056Smarkadams4 
97d529f056Smarkadams4   Input Parameters:
98d529f056Smarkadams4 + pc - the preconditioner context
99d529f056Smarkadams4 - n  - 1 or more (default = 2)
100d529f056Smarkadams4 
101d529f056Smarkadams4   Options Database Key:
102*aea0ef14SMark Adams . -pc_gamg_aggressive_mis_k <n,default=2> - the flag
103d529f056Smarkadams4 
104d529f056Smarkadams4   Level: intermediate
105d529f056Smarkadams4 
106a077d33dSBarry Smith .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()`
107d529f056Smarkadams4 @*/
108d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetAggressive(PC pc, PetscInt n)
109d529f056Smarkadams4 {
110d529f056Smarkadams4   PetscFunctionBegin;
111d529f056Smarkadams4   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
112d529f056Smarkadams4   PetscValidLogicalCollectiveInt(pc, n, 2);
113d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGMISkSetAggressive_C", (PC, PetscInt), (pc, n));
114d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
115d529f056Smarkadams4 }
116d529f056Smarkadams4 
117d529f056Smarkadams4 /*@
118*aea0ef14SMark Adams   PCGAMGSetAggressiveSquareGraph - Use graph square, $A^T A$, for aggressive coarsening. Coarsening is slower than the alternative (MIS-2), which is faster and uses less memory
119d529f056Smarkadams4 
120d529f056Smarkadams4   Logically Collective
121d529f056Smarkadams4 
122d529f056Smarkadams4   Input Parameters:
123d529f056Smarkadams4 + pc - the preconditioner context
124*aea0ef14SMark Adams - b  - default true
125d529f056Smarkadams4 
126d529f056Smarkadams4   Options Database Key:
127*aea0ef14SMark Adams . -pc_gamg_aggressive_square_graph <bool,default=true> - the flag
128d529f056Smarkadams4 
129d529f056Smarkadams4   Level: intermediate
130d529f056Smarkadams4 
131*aea0ef14SMark Adams   Notes:
132*aea0ef14SMark Adams   If `b` is `PETSC_FALSE` then MIS-k is used for aggressive coarsening, see `PCGAMGMISkSetAggressive()`
133*aea0ef14SMark Adams 
134*aea0ef14SMark Adams   Squaring the matrix to perform the aggressive coarsening is slower and requires more memory than using MIS-k, but may result in a better preconditioner
135*aea0ef14SMark Adams   that converges faster.
136*aea0ef14SMark Adams 
137a077d33dSBarry Smith .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()`
138d529f056Smarkadams4 @*/
139d529f056Smarkadams4 PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC pc, PetscBool b)
140d529f056Smarkadams4 {
141d529f056Smarkadams4   PetscFunctionBegin;
142d529f056Smarkadams4   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
143d529f056Smarkadams4   PetscValidLogicalCollectiveBool(pc, b, 2);
144d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGSetAggressiveSquareGraph_C", (PC, PetscBool), (pc, b));
145d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
146d529f056Smarkadams4 }
147d529f056Smarkadams4 
148d529f056Smarkadams4 /*@
149d529f056Smarkadams4   PCGAMGMISkSetMinDegreeOrdering - Use minimum degree ordering in greedy MIS algorithm
150d529f056Smarkadams4 
151d529f056Smarkadams4   Logically Collective
152d529f056Smarkadams4 
153d529f056Smarkadams4   Input Parameters:
154d529f056Smarkadams4 + pc - the preconditioner context
155*aea0ef14SMark Adams - b  - default false
156d529f056Smarkadams4 
157d529f056Smarkadams4   Options Database Key:
158*aea0ef14SMark Adams . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=false> - the flag
159d529f056Smarkadams4 
160d529f056Smarkadams4   Level: intermediate
161d529f056Smarkadams4 
162a077d33dSBarry Smith .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()`
163d529f056Smarkadams4 @*/
164d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC pc, PetscBool b)
165d529f056Smarkadams4 {
166d529f056Smarkadams4   PetscFunctionBegin;
167d529f056Smarkadams4   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
168d529f056Smarkadams4   PetscValidLogicalCollectiveBool(pc, b, 2);
169d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGMISkSetMinDegreeOrdering_C", (PC, PetscBool), (pc, b));
170d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
171d529f056Smarkadams4 }
172d529f056Smarkadams4 
173a9ccf005SMark Adams /*@
174a9ccf005SMark Adams   PCGAMGSetLowMemoryFilter - Use low memory graph/matrix filter
175a9ccf005SMark Adams 
176a9ccf005SMark Adams   Logically Collective
177a9ccf005SMark Adams 
178a9ccf005SMark Adams   Input Parameters:
179a9ccf005SMark Adams + pc - the preconditioner context
180a9ccf005SMark Adams - b  - default false
181a9ccf005SMark Adams 
182a9ccf005SMark Adams   Options Database Key:
183*aea0ef14SMark Adams . -pc_gamg_low_memory_threshold_filter <bool,default=false> - the flag
184a9ccf005SMark Adams 
185a9ccf005SMark Adams   Level: intermediate
186a9ccf005SMark Adams 
187a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`,
188a077d33dSBarry Smith   `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`
189a9ccf005SMark Adams @*/
190a9ccf005SMark Adams PetscErrorCode PCGAMGSetLowMemoryFilter(PC pc, PetscBool b)
191a9ccf005SMark Adams {
192a9ccf005SMark Adams   PetscFunctionBegin;
193a9ccf005SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
194a9ccf005SMark Adams   PetscValidLogicalCollectiveBool(pc, b, 2);
195a9ccf005SMark Adams   PetscTryMethod(pc, "PCGAMGSetLowMemoryFilter_C", (PC, PetscBool), (pc, b));
196a9ccf005SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
197a9ccf005SMark Adams }
198a9ccf005SMark Adams 
1996c34c54dSStefano Zampini /*@
200*aea0ef14SMark Adams   PCGAMGSetGraphSymmetrize - Symmetrize graph used for coarsening. Defaults to true, but if matrix has symmetric attribute, then not needed since the graph is already known to be symmetric
2016c34c54dSStefano Zampini 
2026c34c54dSStefano Zampini   Logically Collective
2036c34c54dSStefano Zampini 
2046c34c54dSStefano Zampini   Input Parameters:
2056c34c54dSStefano Zampini + pc - the preconditioner context
206*aea0ef14SMark Adams - b  - default true
2076c34c54dSStefano Zampini 
2086c34c54dSStefano Zampini   Options Database Key:
209*aea0ef14SMark Adams . -pc_gamg_graph_symmetrize <bool,default=true> - the flag
2106c34c54dSStefano Zampini 
2116c34c54dSStefano Zampini   Level: intermediate
2126c34c54dSStefano Zampini 
213*aea0ef14SMark Adams .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `MatCreateGraph()`,
2146c34c54dSStefano Zampini   `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`
2156c34c54dSStefano Zampini @*/
2166c34c54dSStefano Zampini PetscErrorCode PCGAMGSetGraphSymmetrize(PC pc, PetscBool b)
2176c34c54dSStefano Zampini {
2186c34c54dSStefano Zampini   PetscFunctionBegin;
2196c34c54dSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2206c34c54dSStefano Zampini   PetscValidLogicalCollectiveBool(pc, b, 2);
2216c34c54dSStefano Zampini   PetscTryMethod(pc, "PCGAMGSetGraphSymmetrize_C", (PC, PetscBool), (pc, b));
2226c34c54dSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
2236c34c54dSStefano Zampini }
2246c34c54dSStefano Zampini 
225d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n)
226d71ae5a4SJacob Faibussowitsch {
227ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
228ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
229ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
230ef4ad70eSMark F. Adams 
231ef4ad70eSMark F. Adams   PetscFunctionBegin;
232bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = n;
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
234ef4ad70eSMark F. Adams }
235ef4ad70eSMark F. Adams 
236d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetAggressive_AGG(PC pc, PetscInt n)
237d529f056Smarkadams4 {
238d529f056Smarkadams4   PC_MG       *mg          = (PC_MG *)pc->data;
239d529f056Smarkadams4   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
240d529f056Smarkadams4   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
241d529f056Smarkadams4 
242d529f056Smarkadams4   PetscFunctionBegin;
243d529f056Smarkadams4   pc_gamg_agg->aggressive_mis_k = n;
244d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
245d529f056Smarkadams4 }
246d529f056Smarkadams4 
247d529f056Smarkadams4 static PetscErrorCode PCGAMGSetAggressiveSquareGraph_AGG(PC pc, PetscBool b)
248d529f056Smarkadams4 {
249d529f056Smarkadams4   PC_MG       *mg          = (PC_MG *)pc->data;
250d529f056Smarkadams4   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
251d529f056Smarkadams4   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
252d529f056Smarkadams4 
253d529f056Smarkadams4   PetscFunctionBegin;
254d529f056Smarkadams4   pc_gamg_agg->use_aggressive_square_graph = b;
255d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
256d529f056Smarkadams4 }
257d529f056Smarkadams4 
258a9ccf005SMark Adams static PetscErrorCode PCGAMGSetLowMemoryFilter_AGG(PC pc, PetscBool b)
259a9ccf005SMark Adams {
260a9ccf005SMark Adams   PC_MG       *mg          = (PC_MG *)pc->data;
261a9ccf005SMark Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
262a9ccf005SMark Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
263a9ccf005SMark Adams 
264a9ccf005SMark Adams   PetscFunctionBegin;
265a9ccf005SMark Adams   pc_gamg_agg->use_low_mem_filter = b;
266a9ccf005SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
267a9ccf005SMark Adams }
268a9ccf005SMark Adams 
2696c34c54dSStefano Zampini static PetscErrorCode PCGAMGSetGraphSymmetrize_AGG(PC pc, PetscBool b)
2706c34c54dSStefano Zampini {
2716c34c54dSStefano Zampini   PC_MG       *mg          = (PC_MG *)pc->data;
2726c34c54dSStefano Zampini   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
2736c34c54dSStefano Zampini   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
2746c34c54dSStefano Zampini 
2756c34c54dSStefano Zampini   PetscFunctionBegin;
2766c34c54dSStefano Zampini   pc_gamg_agg->graph_symmetrize = b;
2776c34c54dSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
2786c34c54dSStefano Zampini }
2796c34c54dSStefano Zampini 
280d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetMinDegreeOrdering_AGG(PC pc, PetscBool b)
281d529f056Smarkadams4 {
282d529f056Smarkadams4   PC_MG       *mg          = (PC_MG *)pc->data;
283d529f056Smarkadams4   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
284d529f056Smarkadams4   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
285d529f056Smarkadams4 
286d529f056Smarkadams4   PetscFunctionBegin;
287d529f056Smarkadams4   pc_gamg_agg->use_minimum_degree_ordering = b;
288d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
289d529f056Smarkadams4 }
290d529f056Smarkadams4 
291ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems PetscOptionsObject)
292d71ae5a4SJacob Faibussowitsch {
2932e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
2942e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
295c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
296d4adc10fSMark Adams   PetscBool    n_aggressive_flg, old_sq_provided = PETSC_FALSE, new_sq_provided = PETSC_FALSE, new_sqr_graph = pc_gamg_agg->use_aggressive_square_graph;
297d4adc10fSMark Adams   PetscInt     nsq_graph_old = 0;
2982e68589bSMark F. Adams 
2992e68589bSMark F. Adams   PetscFunctionBegin;
300d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options");
301a077d33dSBarry Smith   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));
302d4adc10fSMark Adams   // aggressive coarsening logic with deprecated -pc_gamg_square_graph
303d4adc10fSMark Adams   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));
304d4adc10fSMark Adams   if (!n_aggressive_flg)
305d4adc10fSMark Adams     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));
306*aea0ef14SMark Adams   PetscCall(PetscOptionsBool("-pc_gamg_aggressive_square_graph", "Use square graph $ (A^T A)$ for aggressive coarsening, if false, MIS-k (k=2) is used, see PCGAMGMISkSetAggressive()", "PCGAMGSetAggressiveSquareGraph", new_sqr_graph, &pc_gamg_agg->use_aggressive_square_graph, &new_sq_provided));
307d4adc10fSMark Adams   if (!new_sq_provided && old_sq_provided) {
308d4adc10fSMark Adams     pc_gamg_agg->aggressive_coarsening_levels = nsq_graph_old; // could be zero
309d4adc10fSMark Adams     pc_gamg_agg->use_aggressive_square_graph  = PETSC_TRUE;
310bae903cbSmarkadams4   }
311d4adc10fSMark Adams   if (new_sq_provided && old_sq_provided)
312835f2295SStefano Zampini     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));
313d529f056Smarkadams4   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));
314a9ccf005SMark Adams   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));
315d529f056Smarkadams4   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));
3166c34c54dSStefano Zampini   PetscCall(PetscOptionsBool("-pc_gamg_graph_symmetrize", "Symmetrize graph for coarsening", "PCGAMGSetGraphSymmetrize", pc_gamg_agg->graph_symmetrize, &pc_gamg_agg->graph_symmetrize, NULL));
317d0609cedSBarry Smith   PetscOptionsHeadEnd();
3183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3192e68589bSMark F. Adams }
3202e68589bSMark F. Adams 
321d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
322d71ae5a4SJacob Faibussowitsch {
3232e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
3242e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
325e0b7e82fSBarry Smith   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
3262e68589bSMark F. Adams 
3272e68589bSMark F. Adams   PetscFunctionBegin;
328e0b7e82fSBarry Smith   PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs));
3299566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
3302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL));
331bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL));
332d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", NULL));
333d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", NULL));
334a9ccf005SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", NULL));
335d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", NULL));
3366c34c54dSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetGraphSymmetrize_C", NULL));
337bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3392e68589bSMark F. Adams }
3402e68589bSMark F. Adams 
3412e68589bSMark F. Adams /*
3422e68589bSMark F. Adams    PCSetCoordinates_AGG
34320f4b53cSBarry Smith 
34420f4b53cSBarry Smith    Collective
3452e68589bSMark F. Adams 
3462e68589bSMark F. Adams    Input Parameter:
3472e68589bSMark F. Adams    . pc - the preconditioner context
348145b44c9SPierre Jolivet    . ndm - dimension of data (used for dof/vertex for Stokes)
349302f38e8SMark F. Adams    . a_nloc - number of vertices local
350302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
3512e68589bSMark F. Adams */
3521e6b0712SBarry Smith 
353d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
354d71ae5a4SJacob Faibussowitsch {
3552e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
3562e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
35769344418SMark F. Adams   PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf;
358a2f3521dSMark F. Adams   Mat      mat = pc->pmat;
3592e68589bSMark F. Adams 
3602e68589bSMark F. Adams   PetscFunctionBegin;
361a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
362a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
363302f38e8SMark F. Adams   nloc = a_nloc;
3642e68589bSMark F. Adams 
3652e68589bSMark F. Adams   /* SA: null space vectors */
3669566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &ndf));               /* this does not work for Stokes */
36769344418SMark F. Adams   if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
368a2f3521dSMark F. Adams   else if (coords) {
36963a3b9bcSJacob Faibussowitsch     PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf);
37069344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */
371ad540459SPierre Jolivet     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);
372806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
37371959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
37463a3b9bcSJacob Faibussowitsch   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);
375c8b0795cSMark F. Adams   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
3762e68589bSMark F. Adams 
3777f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
3789566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
3799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
3802e68589bSMark F. Adams   }
3815e116b59SBarry Smith   /* copy data in - column-oriented */
3822e68589bSMark F. Adams   for (kk = 0; kk < nloc; kk++) {
38369344418SMark F. Adams     const PetscInt M    = nloc * pc_gamg->data_cell_rows; /* stride into data */
38469344418SMark F. Adams     PetscReal     *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */
385e0b7e82fSBarry Smith 
386c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols == 1) *data = 1.0;
3872e68589bSMark F. Adams     else {
38869344418SMark F. Adams       /* translational modes */
3892fa5cd67SKarl Rupp       for (ii = 0; ii < ndatarows; ii++) {
3902fa5cd67SKarl Rupp         for (jj = 0; jj < ndatarows; jj++) {
39169344418SMark F. Adams           if (ii == jj) data[ii * M + jj] = 1.0;
3922e68589bSMark F. Adams           else data[ii * M + jj] = 0.0;
3932fa5cd67SKarl Rupp         }
3942fa5cd67SKarl Rupp       }
39569344418SMark F. Adams 
39669344418SMark F. Adams       /* rotational modes */
3972e68589bSMark F. Adams       if (coords) {
39869344418SMark F. Adams         if (ndm == 2) {
3992e68589bSMark F. Adams           data += 2 * M;
4002e68589bSMark F. Adams           data[0] = -coords[2 * kk + 1];
4012e68589bSMark F. Adams           data[1] = coords[2 * kk];
402806fa848SBarry Smith         } else {
4032e68589bSMark F. Adams           data += 3 * M;
4049371c9d4SSatish Balay           data[0]         = 0.0;
4059371c9d4SSatish Balay           data[M + 0]     = coords[3 * kk + 2];
4069371c9d4SSatish Balay           data[2 * M + 0] = -coords[3 * kk + 1];
4079371c9d4SSatish Balay           data[1]         = -coords[3 * kk + 2];
4089371c9d4SSatish Balay           data[M + 1]     = 0.0;
4099371c9d4SSatish Balay           data[2 * M + 1] = coords[3 * kk];
4109371c9d4SSatish Balay           data[2]         = coords[3 * kk + 1];
4119371c9d4SSatish Balay           data[M + 2]     = -coords[3 * kk];
4129371c9d4SSatish Balay           data[2 * M + 2] = 0.0;
4132e68589bSMark F. Adams         }
4142e68589bSMark F. Adams       }
4152e68589bSMark F. Adams     }
4162e68589bSMark F. Adams   }
4172e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
4183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4192e68589bSMark F. Adams }
4202e68589bSMark F. Adams 
4212e68589bSMark F. Adams /*
422a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
423a2f3521dSMark F. Adams       Looks in Mat for near null space.
424a2f3521dSMark F. Adams       Does not work for Stokes
4252e68589bSMark F. Adams 
4262e68589bSMark F. Adams   Input Parameter:
4272e68589bSMark F. Adams    . pc -
428a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
4292e68589bSMark F. Adams */
430d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
431d71ae5a4SJacob Faibussowitsch {
432b8cd405aSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
433b8cd405aSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
434b8cd405aSMark F. Adams   MatNullSpace mnull;
43566f2ef4bSMark Adams 
436b817416eSBarry Smith   PetscFunctionBegin;
4379566063dSJacob Faibussowitsch   PetscCall(MatGetNearNullSpace(a_A, &mnull));
438b8cd405aSMark F. Adams   if (!mnull) {
43966f2ef4bSMark Adams     DM dm;
440e0b7e82fSBarry Smith 
4419566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
44248a46eb9SPierre Jolivet     if (!dm) PetscCall(MatGetDM(a_A, &dm));
44366f2ef4bSMark Adams     if (dm) {
44466f2ef4bSMark Adams       PetscObject deformation;
445b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
446b0d30dd6SMatthew G. Knepley 
4479566063dSJacob Faibussowitsch       PetscCall(DMGetNumFields(dm, &Nf));
448b0d30dd6SMatthew G. Knepley       if (Nf) {
4499566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, 0, NULL, &deformation));
4504c4edb6fSBarry Smith         if (deformation) {
451835f2295SStefano Zampini           PetscCall(PetscObjectQuery(deformation, "nearnullspace", (PetscObject *)&mnull));
452835f2295SStefano Zampini           if (!mnull) PetscCall(PetscObjectQuery(deformation, "nullspace", (PetscObject *)&mnull));
45366f2ef4bSMark Adams         }
45466f2ef4bSMark Adams       }
455b0d30dd6SMatthew G. Knepley     }
4564c4edb6fSBarry Smith   }
45766f2ef4bSMark Adams 
45866f2ef4bSMark Adams   if (!mnull) {
459a2f3521dSMark F. Adams     PetscInt bs, NN, MM;
460e0b7e82fSBarry Smith 
4619566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
4629566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
46363a3b9bcSJacob Faibussowitsch     PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs);
4649566063dSJacob Faibussowitsch     PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL));
465806fa848SBarry Smith   } else {
466b8cd405aSMark F. Adams     PetscReal         *nullvec;
467b8cd405aSMark F. Adams     PetscBool          has_const;
468b8cd405aSMark F. Adams     PetscInt           i, j, mlocal, nvec, bs;
469d5d25401SBarry Smith     const Vec         *vecs;
470d5d25401SBarry Smith     const PetscScalar *v;
471b817416eSBarry Smith 
4729566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &mlocal, NULL));
4739566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
474d19c4ebbSmarkadams4     for (i = 0; i < nvec; i++) {
475d19c4ebbSmarkadams4       PetscCall(VecGetLocalSize(vecs[i], &j));
476d19c4ebbSmarkadams4       PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal);
477d19c4ebbSmarkadams4     }
478a0dea87dSMark Adams     pc_gamg->data_sz = (nvec + !!has_const) * mlocal;
4799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
4809371c9d4SSatish Balay     if (has_const)
4819371c9d4SSatish Balay       for (i = 0; i < mlocal; i++) nullvec[i] = 1.0;
482b8cd405aSMark F. Adams     for (i = 0; i < nvec; i++) {
4839566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vecs[i], &v));
484b8cd405aSMark F. Adams       for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]);
4859566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vecs[i], &v));
486b8cd405aSMark F. Adams     }
487b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
488b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec + !!has_const);
4899566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
490b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
491b8cd405aSMark F. Adams   }
4923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4932e68589bSMark F. Adams }
4942e68589bSMark F. Adams 
4952e68589bSMark F. Adams /*
496bae903cbSmarkadams4   formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0
4972e68589bSMark F. Adams 
4982e68589bSMark F. Adams   Input Parameter:
4992adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
5002adfe9d3SBarry Smith    . bs - row block size
5010cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
5020cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
5032adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
5042e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
5052e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
506bae903cbSmarkadams4 
5072e68589bSMark F. Adams   Output Parameter:
5082e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
5092e68589bSMark F. Adams    . a_Prol - prolongation operator
5102e68589bSMark F. Adams */
511d71ae5a4SJacob Faibussowitsch 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)
512d71ae5a4SJacob Faibussowitsch {
513bd026e97SJed Brown   PetscInt        Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride;
5143b4367a7SBarry Smith   MPI_Comm        comm;
5152e68589bSMark F. Adams   PetscReal      *out_data;
516539c167fSBarry Smith   PetscCDIntNd   *pos;
5171943db53SBarry Smith   PCGAMGHashTable fgid_flid;
5180cbbd2e1SMark F. Adams 
5192e68589bSMark F. Adams   PetscFunctionBegin;
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
5219566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
5229371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
5239371c9d4SSatish Balay   my0  = Istart / bs;
52463a3b9bcSJacob Faibussowitsch   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);
5250cbbd2e1SMark F. Adams   Iend /= bs;
5260cbbd2e1SMark F. Adams   nghosts = data_stride / bs - nloc;
5272e68589bSMark F. Adams 
5289566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid));
52948a46eb9SPierre Jolivet   for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk));
5302e68589bSMark F. Adams 
5310cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
5320cbbd2e1SMark F. Adams   for (nSelected = mm = 0; mm < nloc; mm++) {
533e78576d6SMark F. Adams     PetscBool ise;
534e0b7e82fSBarry Smith 
5355e62d202SMark Adams     PetscCall(PetscCDIsEmptyAt(agg_llists, mm, &ise));
536e78576d6SMark F. Adams     if (!ise) nSelected++;
5370cbbd2e1SMark F. Adams   }
5389566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj));
53963a3b9bcSJacob Faibussowitsch   PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT "  != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs);
54063a3b9bcSJacob Faibussowitsch   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);
5410cbbd2e1SMark F. Adams 
5422e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
5430cbbd2e1SMark F. Adams   out_data_stride = nSelected * nSAvec;
5442fa5cd67SKarl Rupp 
5459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data));
54633812677SSatish Balay   for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL;
5470cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
5482e68589bSMark F. Adams 
5492e68589bSMark F. Adams   /* find points and set prolongation */
550c8b0795cSMark F. Adams   minsz = 100;
5510cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
5525e62d202SMark Adams     PetscCall(PetscCDCountAt(agg_llists, mm, &jj));
553e78576d6SMark F. Adams     if (jj > 0) {
5540cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
5550cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
5566497c311SBarry Smith       PetscBLASInt   asz, M, N, INFO;
5576497c311SBarry Smith       PetscBLASInt   Mdata, LDA, LWORK;
5582e68589bSMark F. Adams       PetscScalar   *qqc, *qqr, *TAU, *WORK;
5592e68589bSMark F. Adams       PetscInt      *fids;
56065d7b583SSatish Balay       PetscReal     *data;
561b817416eSBarry Smith 
5626497c311SBarry Smith       PetscCall(PetscBLASIntCast(jj, &asz));
5636497c311SBarry Smith       PetscCall(PetscBLASIntCast(asz * bs, &M));
5646497c311SBarry Smith       PetscCall(PetscBLASIntCast(nSAvec, &N));
5656497c311SBarry Smith       PetscCall(PetscBLASIntCast(M + ((N - M > 0) ? N - M : 0), &Mdata));
5666497c311SBarry Smith       PetscCall(PetscBLASIntCast(Mdata, &LDA));
5676497c311SBarry Smith       PetscCall(PetscBLASIntCast(N * bs, &LWORK));
5680cbbd2e1SMark F. Adams       /* count agg */
5690cbbd2e1SMark F. Adams       if (asz < minsz) minsz = asz;
5700cbbd2e1SMark F. Adams 
5710cbbd2e1SMark F. Adams       /* get block */
5729566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids));
5732e68589bSMark F. Adams 
5742e68589bSMark F. Adams       aggID = 0;
5759566063dSJacob Faibussowitsch       PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
576e78576d6SMark F. Adams       while (pos) {
577e78576d6SMark F. Adams         PetscInt gid1;
578e0b7e82fSBarry Smith 
5799566063dSJacob Faibussowitsch         PetscCall(PetscCDIntNdGetID(pos, &gid1));
5809566063dSJacob Faibussowitsch         PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
581e78576d6SMark F. Adams 
5820cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
5830cbbd2e1SMark F. Adams         else {
5849566063dSJacob Faibussowitsch           PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid));
58508401ef6SPierre Jolivet           PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table");
5860cbbd2e1SMark F. Adams         }
5875e116b59SBarry Smith         /* copy in B_i matrix - column-oriented */
58865d7b583SSatish Balay         data = &data_in[flid * bs];
5893d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
5902e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
5910cbbd2e1SMark F. Adams             PetscReal d = data[jj * data_stride + ii];
592e0b7e82fSBarry Smith 
5930cbbd2e1SMark F. Adams             qqc[jj * Mdata + aggID * bs + ii] = d;
5942e68589bSMark F. Adams           }
5952e68589bSMark F. Adams         }
5962e68589bSMark F. Adams         /* set fine IDs */
5972e68589bSMark F. Adams         for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk;
5982e68589bSMark F. Adams         aggID++;
5990cbbd2e1SMark F. Adams       }
6002e68589bSMark F. Adams 
6012e68589bSMark F. Adams       /* pad with zeros */
6022e68589bSMark F. Adams       for (ii = asz * bs; ii < Mdata; ii++) {
603ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0;
6042e68589bSMark F. Adams       }
6052e68589bSMark F. Adams 
6062e68589bSMark F. Adams       /* QR */
6079566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
608792fecdfSBarry Smith       PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
6099566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
61008401ef6SPierre Jolivet       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error");
6115e116b59SBarry Smith       /* get R - column-oriented - output B_{i+1} */
6122e68589bSMark F. Adams       {
6132e68589bSMark F. Adams         PetscReal *data = &out_data[clid * nSAvec];
614e0b7e82fSBarry Smith 
6152e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
6162e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
61708401ef6SPierre Jolivet             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);
6180cbbd2e1SMark F. Adams             if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]);
6190cbbd2e1SMark F. Adams             else data[jj * out_data_stride + ii] = 0.;
6202e68589bSMark F. Adams           }
6212e68589bSMark F. Adams         }
6222e68589bSMark F. Adams       }
6232e68589bSMark F. Adams 
6245e116b59SBarry Smith       /* get Q - row-oriented */
625792fecdfSBarry Smith       PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
62663a3b9bcSJacob Faibussowitsch       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO);
6272e68589bSMark F. Adams 
6282e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
629ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii];
6302e68589bSMark F. Adams       }
6312e68589bSMark F. Adams 
6322e68589bSMark F. Adams       /* add diagonal block of P0 */
6339371c9d4SSatish Balay       for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ }
6349566063dSJacob Faibussowitsch       PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES));
6359566063dSJacob Faibussowitsch       PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids));
636b43b03e9SMark F. Adams       clid++;
6370cbbd2e1SMark F. Adams     } /* coarse agg */
6380cbbd2e1SMark F. Adams   } /* for all fine nodes */
6399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
6409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
6419566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableDestroy(&fgid_flid));
6423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6432e68589bSMark F. Adams }
6442e68589bSMark F. Adams 
645d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer)
646d71ae5a4SJacob Faibussowitsch {
6475adeb434SBarry Smith   PC_MG       *mg          = (PC_MG *)pc->data;
6485adeb434SBarry Smith   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
6495adeb434SBarry Smith   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
6505adeb434SBarry Smith 
6515adeb434SBarry Smith   PetscFunctionBegin;
6529566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      AGG specific options\n"));
653835f2295SStefano Zampini   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number of levels of aggressive coarsening %" PetscInt_FMT "\n", pc_gamg_agg->aggressive_coarsening_levels));
654d529f056Smarkadams4   if (pc_gamg_agg->aggressive_coarsening_levels > 0) {
655d529f056Smarkadams4     PetscCall(PetscViewerASCIIPrintf(viewer, "        %s aggressive coarsening\n", !pc_gamg_agg->use_aggressive_square_graph ? "MIS-k" : "Square graph"));
656835f2295SStefano Zampini     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));
657d529f056Smarkadams4   }
658e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
659e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
660e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
661e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
6620acf423eSBarry Smith   if (pc_gamg_agg->crs) PetscCall(MatCoarsenView(pc_gamg_agg->crs, viewer));
6630acf423eSBarry Smith   else PetscCall(PetscViewerASCIIPrintf(viewer, "Coarsening algorithm not yet selected\n"));
664e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
665e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
666e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
667e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
668835f2295SStefano Zampini   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number smoothing steps to construct prolongation %" PetscInt_FMT "\n", pc_gamg_agg->nsmooths));
6693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6705adeb434SBarry Smith }
6715adeb434SBarry Smith 
6722d776b49SBarry Smith static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat)
673d71ae5a4SJacob Faibussowitsch {
674c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG *)pc->data;
675c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG *)mg->innerctx;
676c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
6772d776b49SBarry Smith   const PetscReal vfilter     = pc_gamg->threshold[pc_gamg->current_level];
6789c15c1aeSMark Adams   PetscBool       ishem, ismis;
6792d776b49SBarry Smith   const char     *prefix;
680d529f056Smarkadams4   MatInfo         info0, info1;
681d529f056Smarkadams4   PetscInt        bs;
682c8b0795cSMark F. Adams 
683c8b0795cSMark F. Adams   PetscFunctionBegin;
684a9ccf005SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
6852d776b49SBarry Smith   /* 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 */
6862d776b49SBarry Smith   /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */
687e0b7e82fSBarry Smith   PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs));
6882d776b49SBarry Smith   PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs));
6892d776b49SBarry Smith   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix));
6902d776b49SBarry Smith   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix));
691e0b7e82fSBarry Smith   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc_gamg_agg->crs, "pc_gamg_"));
6922d776b49SBarry Smith   PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs));
693e02fb3cdSMark Adams   PetscCall(MatGetBlockSize(Amat, &bs));
694e02fb3cdSMark Adams   // check for valid indices wrt bs
695e02fb3cdSMark Adams   for (int ii = 0; ii < pc_gamg_agg->crs->strength_index_size; ii++) {
696835f2295SStefano Zampini     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",
697835f2295SStefano Zampini                pc_gamg_agg->crs->strength_index[ii], bs);
698e02fb3cdSMark Adams   }
6992d776b49SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem));
7005e62d202SMark Adams   if (ishem) {
701835f2295SStefano Zampini     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));
7025e62d202SMark Adams     pc_gamg_agg->aggressive_coarsening_levels = 0;                                         // aggressive and HEM does not make sense
7035e62d202SMark Adams     PetscCall(MatCoarsenSetMaximumIterations(pc_gamg_agg->crs, pc_gamg_agg->crs->max_it)); // for code coverage
7045e62d202SMark Adams     PetscCall(MatCoarsenSetThreshold(pc_gamg_agg->crs, vfilter));                          // for code coverage
7059c15c1aeSMark Adams   } else {
7069c15c1aeSMark Adams     PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENMIS, &ismis));
7079c15c1aeSMark Adams     if (ismis && pc_gamg_agg->aggressive_coarsening_levels && !pc_gamg_agg->use_aggressive_square_graph) {
7089c15c1aeSMark Adams       PetscCall(PetscInfo(pc, "MIS and aggressive coarsening and no square graph: force square graph\n"));
7099c15c1aeSMark Adams       pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE;
7109c15c1aeSMark Adams     }
7115e62d202SMark Adams   }
712a9ccf005SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
713a9ccf005SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
714d529f056Smarkadams4   PetscCall(MatGetInfo(Amat, MAT_LOCAL, &info0)); /* global reduction */
715a9ccf005SMark Adams 
716a9ccf005SMark Adams   if (ishem || pc_gamg_agg->use_low_mem_filter) {
7176c34c54dSStefano Zampini     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));
718a9ccf005SMark Adams   } else {
7196c34c54dSStefano Zampini     // make scalar graph, symmetrize if not known to be symmetric, scale, but do not filter (expensive)
7206c34c54dSStefano Zampini     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));
721a9ccf005SMark Adams     if (vfilter >= 0) {
722a9ccf005SMark Adams       PetscInt           Istart, Iend, ncols, nnz0, nnz1, NN, MM, nloc;
723a9ccf005SMark Adams       Mat                tGmat, Gmat = *a_Gmat;
724a9ccf005SMark Adams       MPI_Comm           comm;
725a9ccf005SMark Adams       const PetscScalar *vals;
726a9ccf005SMark Adams       const PetscInt    *idx;
727a9ccf005SMark Adams       PetscInt          *d_nnz, *o_nnz, kk, *garray = NULL, *AJ, maxcols = 0;
728a9ccf005SMark Adams       MatScalar         *AA; // this is checked in graph
729a9ccf005SMark Adams       PetscBool          isseqaij;
730a9ccf005SMark Adams       Mat                a, b, c;
731a9ccf005SMark Adams       MatType            jtype;
732a9ccf005SMark Adams 
733a9ccf005SMark Adams       PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm));
734a9ccf005SMark Adams       PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isseqaij));
735a9ccf005SMark Adams       PetscCall(MatGetType(Gmat, &jtype));
736a9ccf005SMark Adams       PetscCall(MatCreate(comm, &tGmat));
737a9ccf005SMark Adams       PetscCall(MatSetType(tGmat, jtype));
738a9ccf005SMark Adams 
739a9ccf005SMark Adams       /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold?
740a9ccf005SMark Adams         Also, if the matrix is symmetric, can we skip this
741a9ccf005SMark Adams         operation? It can be very expensive on large matrices. */
742a9ccf005SMark Adams 
743a9ccf005SMark Adams       // global sizes
744a9ccf005SMark Adams       PetscCall(MatGetSize(Gmat, &MM, &NN));
745a9ccf005SMark Adams       PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend));
746a9ccf005SMark Adams       nloc = Iend - Istart;
747a9ccf005SMark Adams       PetscCall(PetscMalloc2(nloc, &d_nnz, nloc, &o_nnz));
748a9ccf005SMark Adams       if (isseqaij) {
749a9ccf005SMark Adams         a = Gmat;
750a9ccf005SMark Adams         b = NULL;
751a9ccf005SMark Adams       } else {
752a9ccf005SMark Adams         Mat_MPIAIJ *d = (Mat_MPIAIJ *)Gmat->data;
753e0b7e82fSBarry Smith 
754a9ccf005SMark Adams         a      = d->A;
755a9ccf005SMark Adams         b      = d->B;
756a9ccf005SMark Adams         garray = d->garray;
757a9ccf005SMark Adams       }
758a9ccf005SMark Adams       /* Determine upper bound on non-zeros needed in new filtered matrix */
759a9ccf005SMark Adams       for (PetscInt row = 0; row < nloc; row++) {
760a9ccf005SMark Adams         PetscCall(MatGetRow(a, row, &ncols, NULL, NULL));
761a9ccf005SMark Adams         d_nnz[row] = ncols;
762a9ccf005SMark Adams         if (ncols > maxcols) maxcols = ncols;
763a9ccf005SMark Adams         PetscCall(MatRestoreRow(a, row, &ncols, NULL, NULL));
764a9ccf005SMark Adams       }
765a9ccf005SMark Adams       if (b) {
766a9ccf005SMark Adams         for (PetscInt row = 0; row < nloc; row++) {
767a9ccf005SMark Adams           PetscCall(MatGetRow(b, row, &ncols, NULL, NULL));
768a9ccf005SMark Adams           o_nnz[row] = ncols;
769a9ccf005SMark Adams           if (ncols > maxcols) maxcols = ncols;
770a9ccf005SMark Adams           PetscCall(MatRestoreRow(b, row, &ncols, NULL, NULL));
771a9ccf005SMark Adams         }
772a9ccf005SMark Adams       }
773a9ccf005SMark Adams       PetscCall(MatSetSizes(tGmat, nloc, nloc, MM, MM));
774a9ccf005SMark Adams       PetscCall(MatSetBlockSizes(tGmat, 1, 1));
775a9ccf005SMark Adams       PetscCall(MatSeqAIJSetPreallocation(tGmat, 0, d_nnz));
776a9ccf005SMark Adams       PetscCall(MatMPIAIJSetPreallocation(tGmat, 0, d_nnz, 0, o_nnz));
777a9ccf005SMark Adams       PetscCall(MatSetOption(tGmat, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
778a9ccf005SMark Adams       PetscCall(PetscFree2(d_nnz, o_nnz));
779a9ccf005SMark Adams       PetscCall(PetscMalloc2(maxcols, &AA, maxcols, &AJ));
780a9ccf005SMark Adams       nnz0 = nnz1 = 0;
781a9ccf005SMark Adams       for (c = a, kk = 0; c && kk < 2; c = b, kk++) {
782a9ccf005SMark Adams         for (PetscInt row = 0, grow = Istart, ncol_row, jj; row < nloc; row++, grow++) {
783a9ccf005SMark Adams           PetscCall(MatGetRow(c, row, &ncols, &idx, &vals));
784a9ccf005SMark Adams           for (ncol_row = jj = 0; jj < ncols; jj++, nnz0++) {
785a9ccf005SMark Adams             PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
786a9ccf005SMark Adams             if (PetscRealPart(sv) > vfilter) {
787a9ccf005SMark Adams               PetscInt cid = idx[jj] + Istart; //diag
788e0b7e82fSBarry Smith 
789a9ccf005SMark Adams               nnz1++;
790a9ccf005SMark Adams               if (c != a) cid = garray[idx[jj]];
791a9ccf005SMark Adams               AA[ncol_row] = vals[jj];
792a9ccf005SMark Adams               AJ[ncol_row] = cid;
793a9ccf005SMark Adams               ncol_row++;
794a9ccf005SMark Adams             }
795a9ccf005SMark Adams           }
796a9ccf005SMark Adams           PetscCall(MatRestoreRow(c, row, &ncols, &idx, &vals));
797a9ccf005SMark Adams           PetscCall(MatSetValues(tGmat, 1, &grow, ncol_row, AJ, AA, INSERT_VALUES));
798a9ccf005SMark Adams         }
799a9ccf005SMark Adams       }
800a9ccf005SMark Adams       PetscCall(PetscFree2(AA, AJ));
801a9ccf005SMark Adams       PetscCall(MatAssemblyBegin(tGmat, MAT_FINAL_ASSEMBLY));
802a9ccf005SMark Adams       PetscCall(MatAssemblyEnd(tGmat, MAT_FINAL_ASSEMBLY));
803a9ccf005SMark Adams       PetscCall(MatPropagateSymmetryOptions(Gmat, tGmat)); /* Normal Mat options are not relevant ? */
804a9ccf005SMark Adams       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));
805a9ccf005SMark Adams       PetscCall(MatViewFromOptions(tGmat, NULL, "-mat_filter_graph_view"));
806a9ccf005SMark Adams       PetscCall(MatDestroy(&Gmat));
807a9ccf005SMark Adams       *a_Gmat = tGmat;
808a9ccf005SMark Adams     }
809a9ccf005SMark Adams   }
810a9ccf005SMark Adams 
811d529f056Smarkadams4   PetscCall(MatGetInfo(*a_Gmat, MAT_LOCAL, &info1)); /* global reduction */
812d529f056Smarkadams4   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));
813849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
8143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
815c8b0795cSMark F. Adams }
816c8b0795cSMark F. Adams 
817d529f056Smarkadams4 typedef PetscInt    NState;
818d529f056Smarkadams4 static const NState NOT_DONE = -2;
819d529f056Smarkadams4 static const NState DELETED  = -1;
820d529f056Smarkadams4 static const NState REMOVED  = -3;
821d529f056Smarkadams4 #define IS_SELECTED(s) (s != DELETED && s != NOT_DONE && s != REMOVED)
822d529f056Smarkadams4 
823d529f056Smarkadams4 /*
824d529f056Smarkadams4    fixAggregatesWithSquare - greedy grab of with G1 (unsquared graph) -- AIJ specific -- change to fixAggregatesWithSquare -- TODD
825d529f056Smarkadams4      - AGG-MG specific: clears singletons out of 'selected_2'
826d529f056Smarkadams4 
827d529f056Smarkadams4    Input Parameter:
828d529f056Smarkadams4    . Gmat_2 - global matrix of squared graph (data not defined)
829d529f056Smarkadams4    . Gmat_1 - base graph to grab with base graph
830d529f056Smarkadams4    Input/Output Parameter:
831d529f056Smarkadams4    . aggs_2 - linked list of aggs with gids)
832d529f056Smarkadams4 */
833d529f056Smarkadams4 static PetscErrorCode fixAggregatesWithSquare(PC pc, Mat Gmat_2, Mat Gmat_1, PetscCoarsenData *aggs_2)
834d529f056Smarkadams4 {
835d529f056Smarkadams4   PetscBool      isMPI;
836d529f056Smarkadams4   Mat_SeqAIJ    *matA_1, *matB_1 = NULL;
837d529f056Smarkadams4   MPI_Comm       comm;
838d529f056Smarkadams4   PetscInt       lid, *ii, *idx, ix, Iend, my0, kk, n, j;
839d529f056Smarkadams4   Mat_MPIAIJ    *mpimat_2 = NULL, *mpimat_1 = NULL;
840d529f056Smarkadams4   const PetscInt nloc = Gmat_2->rmap->n;
841d529f056Smarkadams4   PetscScalar   *cpcol_1_state, *cpcol_2_state, *cpcol_2_par_orig, *lid_parent_gid;
842d529f056Smarkadams4   PetscInt      *lid_cprowID_1 = NULL;
843d529f056Smarkadams4   NState        *lid_state;
844d529f056Smarkadams4   Vec            ghost_par_orig2;
845d529f056Smarkadams4   PetscMPIInt    rank;
846d529f056Smarkadams4 
847d529f056Smarkadams4   PetscFunctionBegin;
848d529f056Smarkadams4   PetscCall(PetscObjectGetComm((PetscObject)Gmat_2, &comm));
849d529f056Smarkadams4   PetscCallMPI(MPI_Comm_rank(comm, &rank));
850d529f056Smarkadams4   PetscCall(MatGetOwnershipRange(Gmat_1, &my0, &Iend));
851d529f056Smarkadams4 
852d529f056Smarkadams4   /* get submatrices */
853d529f056Smarkadams4   PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATMPIAIJ, &isMPI));
854d529f056Smarkadams4   PetscCall(PetscInfo(pc, "isMPI = %s\n", isMPI ? "yes" : "no"));
855d529f056Smarkadams4   PetscCall(PetscMalloc3(nloc, &lid_state, nloc, &lid_parent_gid, nloc, &lid_cprowID_1));
856d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
857d529f056Smarkadams4   if (isMPI) {
858d529f056Smarkadams4     /* grab matrix objects */
859d529f056Smarkadams4     mpimat_2 = (Mat_MPIAIJ *)Gmat_2->data;
860d529f056Smarkadams4     mpimat_1 = (Mat_MPIAIJ *)Gmat_1->data;
861d529f056Smarkadams4     matA_1   = (Mat_SeqAIJ *)mpimat_1->A->data;
862d529f056Smarkadams4     matB_1   = (Mat_SeqAIJ *)mpimat_1->B->data;
863d529f056Smarkadams4 
864d529f056Smarkadams4     /* force compressed row storage for B matrix in AuxMat */
865d529f056Smarkadams4     PetscCall(MatCheckCompressedRow(mpimat_1->B, matB_1->nonzerorowcnt, &matB_1->compressedrow, matB_1->i, Gmat_1->rmap->n, -1.0));
866d529f056Smarkadams4     for (ix = 0; ix < matB_1->compressedrow.nrows; ix++) {
867d529f056Smarkadams4       PetscInt lid = matB_1->compressedrow.rindex[ix];
868e0b7e82fSBarry Smith 
869835f2295SStefano Zampini       PetscCheck(lid <= nloc && lid >= -1, PETSC_COMM_SELF, PETSC_ERR_USER, "lid %" PetscInt_FMT " out of range. nloc = %" PetscInt_FMT, lid, nloc);
870d529f056Smarkadams4       if (lid != -1) lid_cprowID_1[lid] = ix;
871d529f056Smarkadams4     }
872d529f056Smarkadams4   } else {
873d529f056Smarkadams4     PetscBool isAIJ;
874e0b7e82fSBarry Smith 
875d529f056Smarkadams4     PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATSEQAIJ, &isAIJ));
876d529f056Smarkadams4     PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix.");
877d529f056Smarkadams4     matA_1 = (Mat_SeqAIJ *)Gmat_1->data;
878d529f056Smarkadams4   }
879d529f056Smarkadams4   if (nloc > 0) { PetscCheck(!matB_1 || matB_1->compressedrow.use, PETSC_COMM_SELF, PETSC_ERR_PLIB, "matB_1 && !matB_1->compressedrow.use: PETSc bug???"); }
880d529f056Smarkadams4   /* get state of locals and selected gid for deleted */
881d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) {
882d529f056Smarkadams4     lid_parent_gid[lid] = -1.0;
883d529f056Smarkadams4     lid_state[lid]      = DELETED;
884d529f056Smarkadams4   }
885d529f056Smarkadams4 
886d529f056Smarkadams4   /* set lid_state */
887d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) {
888d529f056Smarkadams4     PetscCDIntNd *pos;
889e0b7e82fSBarry Smith 
890d529f056Smarkadams4     PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
891d529f056Smarkadams4     if (pos) {
892d529f056Smarkadams4       PetscInt gid1;
893d529f056Smarkadams4 
894d529f056Smarkadams4       PetscCall(PetscCDIntNdGetID(pos, &gid1));
895835f2295SStefano Zampini       PetscCheck(gid1 == lid + my0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gid1 %" PetscInt_FMT " != lid %" PetscInt_FMT " + my0 %" PetscInt_FMT, gid1, lid, my0);
896d529f056Smarkadams4       lid_state[lid] = gid1;
897d529f056Smarkadams4     }
898d529f056Smarkadams4   }
899d529f056Smarkadams4 
900d529f056Smarkadams4   /* map local to selected local, DELETED means a ghost owns it */
90163bfac88SBarry Smith   for (lid = 0; lid < nloc; lid++) {
902d529f056Smarkadams4     NState state = lid_state[lid];
90363bfac88SBarry Smith 
904d529f056Smarkadams4     if (IS_SELECTED(state)) {
905d529f056Smarkadams4       PetscCDIntNd *pos;
906e0b7e82fSBarry Smith 
907d529f056Smarkadams4       PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
908d529f056Smarkadams4       while (pos) {
909d529f056Smarkadams4         PetscInt gid1;
910e0b7e82fSBarry Smith 
911d529f056Smarkadams4         PetscCall(PetscCDIntNdGetID(pos, &gid1));
912d529f056Smarkadams4         PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos));
913d529f056Smarkadams4         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1 - my0] = (PetscScalar)(lid + my0);
914d529f056Smarkadams4       }
915d529f056Smarkadams4     }
916d529f056Smarkadams4   }
917d529f056Smarkadams4   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
918d529f056Smarkadams4   if (isMPI) {
919d529f056Smarkadams4     Vec tempVec;
92063bfac88SBarry Smith 
921d529f056Smarkadams4     /* get 'cpcol_1_state' */
922d529f056Smarkadams4     PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL));
923d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
924d529f056Smarkadams4       PetscScalar v = (PetscScalar)lid_state[kk];
925e0b7e82fSBarry Smith 
926d529f056Smarkadams4       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
927d529f056Smarkadams4     }
928d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
929d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
930d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD));
931d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD));
932d529f056Smarkadams4     PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state));
933d529f056Smarkadams4     /* get 'cpcol_2_state' */
934d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD));
935d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD));
936d529f056Smarkadams4     PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state));
937d529f056Smarkadams4     /* get 'cpcol_2_par_orig' */
938d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
939835f2295SStefano Zampini       PetscScalar v = lid_parent_gid[kk];
940e0b7e82fSBarry Smith 
941d529f056Smarkadams4       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
942d529f056Smarkadams4     }
943d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
944d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
945d529f056Smarkadams4     PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2));
946d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD));
947d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD));
948d529f056Smarkadams4     PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig));
949d529f056Smarkadams4 
950d529f056Smarkadams4     PetscCall(VecDestroy(&tempVec));
951d529f056Smarkadams4   } /* ismpi */
952d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) {
953d529f056Smarkadams4     NState state = lid_state[lid];
954e0b7e82fSBarry Smith 
955d529f056Smarkadams4     if (IS_SELECTED(state)) {
956d529f056Smarkadams4       /* steal locals */
957d529f056Smarkadams4       ii  = matA_1->i;
958d529f056Smarkadams4       n   = ii[lid + 1] - ii[lid];
959d529f056Smarkadams4       idx = matA_1->j + ii[lid];
960d529f056Smarkadams4       for (j = 0; j < n; j++) {
961d529f056Smarkadams4         PetscInt lidj   = idx[j], sgid;
962d529f056Smarkadams4         NState   statej = lid_state[lidj];
963e0b7e82fSBarry Smith 
964d529f056Smarkadams4         if (statej == DELETED && (sgid = (PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid + my0) { /* steal local */
965d529f056Smarkadams4           lid_parent_gid[lidj] = (PetscScalar)(lid + my0);                                              /* send this if sgid is not local */
966d529f056Smarkadams4           if (sgid >= my0 && sgid < Iend) {                                                             /* I'm stealing this local from a local sgid */
967d529f056Smarkadams4             PetscInt      hav = 0, slid = sgid - my0, gidj = lidj + my0;
968d529f056Smarkadams4             PetscCDIntNd *pos, *last = NULL;
969e0b7e82fSBarry Smith 
970d529f056Smarkadams4             /* looking for local from local so id_llist_2 works */
971d529f056Smarkadams4             PetscCall(PetscCDGetHeadPos(aggs_2, slid, &pos));
972d529f056Smarkadams4             while (pos) {
973d529f056Smarkadams4               PetscInt gid;
97463bfac88SBarry Smith 
975d529f056Smarkadams4               PetscCall(PetscCDIntNdGetID(pos, &gid));
976d529f056Smarkadams4               if (gid == gidj) {
977d529f056Smarkadams4                 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null");
978d529f056Smarkadams4                 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last));
979d529f056Smarkadams4                 PetscCall(PetscCDAppendNode(aggs_2, lid, pos));
980d529f056Smarkadams4                 hav = 1;
981d529f056Smarkadams4                 break;
982d529f056Smarkadams4               } else last = pos;
983d529f056Smarkadams4               PetscCall(PetscCDGetNextPos(aggs_2, slid, &pos));
984d529f056Smarkadams4             }
985d529f056Smarkadams4             if (hav != 1) {
986d529f056Smarkadams4               PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find adj in 'selected' lists - structurally unsymmetric matrix");
987835f2295SStefano Zampini               SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %" PetscInt_FMT " times???", hav);
988d529f056Smarkadams4             }
989d529f056Smarkadams4           } else { /* I'm stealing this local, owned by a ghost */
990d529f056Smarkadams4             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.",
991d529f056Smarkadams4                        ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "");
992d529f056Smarkadams4             PetscCall(PetscCDAppendID(aggs_2, lid, lidj + my0));
993d529f056Smarkadams4           }
994d529f056Smarkadams4         }
995d529f056Smarkadams4       } /* local neighbors */
996d529f056Smarkadams4     } else if (state == DELETED /* && lid_cprowID_1 */) {
997d529f056Smarkadams4       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
998e0b7e82fSBarry Smith 
999d529f056Smarkadams4       /* see if I have a selected ghost neighbor that will steal me */
1000d529f056Smarkadams4       if ((ix = lid_cprowID_1[lid]) != -1) {
1001d529f056Smarkadams4         ii  = matB_1->compressedrow.i;
1002d529f056Smarkadams4         n   = ii[ix + 1] - ii[ix];
1003d529f056Smarkadams4         idx = matB_1->j + ii[ix];
1004d529f056Smarkadams4         for (j = 0; j < n; j++) {
1005d529f056Smarkadams4           PetscInt cpid   = idx[j];
1006d529f056Smarkadams4           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
1007e0b7e82fSBarry Smith 
1008835f2295SStefano Zampini           if (IS_SELECTED(statej) && sgidold != statej) { /* ghost will steal this, remove from my list */
1009d529f056Smarkadams4             lid_parent_gid[lid] = (PetscScalar)statej;    /* send who selected */
1010d529f056Smarkadams4             if (sgidold >= my0 && sgidold < Iend) {       /* this was mine */
1011d529f056Smarkadams4               PetscInt      hav = 0, oldslidj = sgidold - my0;
1012d529f056Smarkadams4               PetscCDIntNd *pos, *last        = NULL;
1013e0b7e82fSBarry Smith 
1014d529f056Smarkadams4               /* remove from 'oldslidj' list */
1015d529f056Smarkadams4               PetscCall(PetscCDGetHeadPos(aggs_2, oldslidj, &pos));
1016d529f056Smarkadams4               while (pos) {
1017d529f056Smarkadams4                 PetscInt gid;
1018e0b7e82fSBarry Smith 
1019d529f056Smarkadams4                 PetscCall(PetscCDIntNdGetID(pos, &gid));
1020d529f056Smarkadams4                 if (lid + my0 == gid) {
1021d529f056Smarkadams4                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
1022d529f056Smarkadams4                   PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null");
1023d529f056Smarkadams4                   PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last));
1024d529f056Smarkadams4                   /* ghost (PetscScalar)statej will add this later */
1025d529f056Smarkadams4                   hav = 1;
1026d529f056Smarkadams4                   break;
1027d529f056Smarkadams4                 } else last = pos;
1028d529f056Smarkadams4                 PetscCall(PetscCDGetNextPos(aggs_2, oldslidj, &pos));
1029d529f056Smarkadams4               }
1030d529f056Smarkadams4               if (hav != 1) {
1031835f2295SStefano Zampini                 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find (hav=%" PetscInt_FMT ") adj in 'selected' lists - structurally unsymmetric matrix", hav);
1032835f2295SStefano Zampini                 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %" PetscInt_FMT " times???", hav);
1033d529f056Smarkadams4               }
1034d529f056Smarkadams4             } else {
1035d529f056Smarkadams4               /* TODO: ghosts remove this later */
1036d529f056Smarkadams4             }
1037d529f056Smarkadams4           }
1038d529f056Smarkadams4         }
1039d529f056Smarkadams4       }
1040d529f056Smarkadams4     } /* selected/deleted */
1041d529f056Smarkadams4   } /* node loop */
1042d529f056Smarkadams4 
1043d529f056Smarkadams4   if (isMPI) {
1044d529f056Smarkadams4     PetscScalar    *cpcol_2_parent, *cpcol_2_gid;
1045d529f056Smarkadams4     Vec             tempVec, ghostgids2, ghostparents2;
1046d529f056Smarkadams4     PetscInt        cpid, nghost_2;
1047d529f056Smarkadams4     PCGAMGHashTable gid_cpid;
1048d529f056Smarkadams4 
1049d529f056Smarkadams4     PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2));
1050d529f056Smarkadams4     PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL));
1051d529f056Smarkadams4 
1052d529f056Smarkadams4     /* get 'cpcol_2_parent' */
1053d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) { PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES)); }
1054d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
1055d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
1056d529f056Smarkadams4     PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2));
1057d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD));
1058d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD));
1059d529f056Smarkadams4     PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent));
1060d529f056Smarkadams4 
1061d529f056Smarkadams4     /* get 'cpcol_2_gid' */
1062d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
1063d529f056Smarkadams4       PetscScalar v = (PetscScalar)j;
1064e0b7e82fSBarry Smith 
1065d529f056Smarkadams4       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
1066d529f056Smarkadams4     }
1067d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
1068d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
1069d529f056Smarkadams4     PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2));
1070d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD));
1071d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD));
1072d529f056Smarkadams4     PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid));
1073d529f056Smarkadams4     PetscCall(VecDestroy(&tempVec));
1074d529f056Smarkadams4 
1075d529f056Smarkadams4     /* look for deleted ghosts and add to table */
1076d529f056Smarkadams4     PetscCall(PCGAMGHashTableCreate(2 * nghost_2 + 1, &gid_cpid));
1077d529f056Smarkadams4     for (cpid = 0; cpid < nghost_2; cpid++) {
1078d529f056Smarkadams4       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
1079e0b7e82fSBarry Smith 
1080d529f056Smarkadams4       if (state == DELETED) {
1081d529f056Smarkadams4         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
1082d529f056Smarkadams4         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
1083e0b7e82fSBarry Smith 
1084d529f056Smarkadams4         if (sgid_old == -1 && sgid_new != -1) {
1085d529f056Smarkadams4           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
1086e0b7e82fSBarry Smith 
1087d529f056Smarkadams4           PetscCall(PCGAMGHashTableAdd(&gid_cpid, gid, cpid));
1088d529f056Smarkadams4         }
1089d529f056Smarkadams4       }
1090d529f056Smarkadams4     }
1091d529f056Smarkadams4 
1092d529f056Smarkadams4     /* look for deleted ghosts and see if they moved - remove it */
1093d529f056Smarkadams4     for (lid = 0; lid < nloc; lid++) {
1094d529f056Smarkadams4       NState state = lid_state[lid];
109563bfac88SBarry Smith 
1096d529f056Smarkadams4       if (IS_SELECTED(state)) {
1097d529f056Smarkadams4         PetscCDIntNd *pos, *last = NULL;
109863bfac88SBarry Smith 
1099d529f056Smarkadams4         /* look for deleted ghosts and see if they moved */
1100d529f056Smarkadams4         PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
1101d529f056Smarkadams4         while (pos) {
1102d529f056Smarkadams4           PetscInt gid;
1103d529f056Smarkadams4 
110463bfac88SBarry Smith           PetscCall(PetscCDIntNdGetID(pos, &gid));
1105d529f056Smarkadams4           if (gid < my0 || gid >= Iend) {
1106d529f056Smarkadams4             PetscCall(PCGAMGHashTableFind(&gid_cpid, gid, &cpid));
1107d529f056Smarkadams4             if (cpid != -1) {
1108d529f056Smarkadams4               /* a moved ghost - */
1109d529f056Smarkadams4               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
1110d529f056Smarkadams4               PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last));
1111d529f056Smarkadams4             } else last = pos;
1112d529f056Smarkadams4           } else last = pos;
1113d529f056Smarkadams4 
1114d529f056Smarkadams4           PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos));
1115d529f056Smarkadams4         } /* loop over list of deleted */
1116d529f056Smarkadams4       } /* selected */
1117d529f056Smarkadams4     }
1118d529f056Smarkadams4     PetscCall(PCGAMGHashTableDestroy(&gid_cpid));
1119d529f056Smarkadams4 
1120d529f056Smarkadams4     /* look at ghosts, see if they changed - and it */
1121d529f056Smarkadams4     for (cpid = 0; cpid < nghost_2; cpid++) {
1122d529f056Smarkadams4       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
1123e0b7e82fSBarry Smith 
1124d529f056Smarkadams4       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
1125d529f056Smarkadams4         PetscInt      gid      = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
1126d529f056Smarkadams4         PetscInt      slid_new = sgid_new - my0, hav = 0;
1127d529f056Smarkadams4         PetscCDIntNd *pos;
1128d529f056Smarkadams4 
1129d529f056Smarkadams4         /* search for this gid to see if I have it */
1130d529f056Smarkadams4         PetscCall(PetscCDGetHeadPos(aggs_2, slid_new, &pos));
1131d529f056Smarkadams4         while (pos) {
1132d529f056Smarkadams4           PetscInt gidj;
1133e0b7e82fSBarry Smith 
1134d529f056Smarkadams4           PetscCall(PetscCDIntNdGetID(pos, &gidj));
1135d529f056Smarkadams4           PetscCall(PetscCDGetNextPos(aggs_2, slid_new, &pos));
1136d529f056Smarkadams4 
1137d529f056Smarkadams4           if (gidj == gid) {
1138d529f056Smarkadams4             hav = 1;
1139d529f056Smarkadams4             break;
1140d529f056Smarkadams4           }
1141d529f056Smarkadams4         }
1142d529f056Smarkadams4         if (hav != 1) {
1143d529f056Smarkadams4           /* insert 'flidj' into head of llist */
1144d529f056Smarkadams4           PetscCall(PetscCDAppendID(aggs_2, slid_new, gid));
1145d529f056Smarkadams4         }
1146d529f056Smarkadams4       }
1147d529f056Smarkadams4     }
1148d529f056Smarkadams4     PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state));
1149d529f056Smarkadams4     PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state));
1150d529f056Smarkadams4     PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent));
1151d529f056Smarkadams4     PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid));
1152d529f056Smarkadams4     PetscCall(VecDestroy(&ghostgids2));
1153d529f056Smarkadams4     PetscCall(VecDestroy(&ghostparents2));
1154d529f056Smarkadams4     PetscCall(VecDestroy(&ghost_par_orig2));
1155d529f056Smarkadams4   }
1156d529f056Smarkadams4   PetscCall(PetscFree3(lid_state, lid_parent_gid, lid_cprowID_1));
1157d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
1158d529f056Smarkadams4 }
1159d529f056Smarkadams4 
1160c8b0795cSMark F. Adams /*
1161bae903cbSmarkadams4    PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for
1162bae903cbSmarkadams4      communication of QR data used with HEM and MISk coarsening
1163c8b0795cSMark F. Adams 
1164c8b0795cSMark F. Adams   Input Parameter:
1165e0940f08SMark F. Adams    . a_pc - this
1166bae903cbSmarkadams4 
1167e0940f08SMark F. Adams   Input/Output Parameter:
1168bae903cbSmarkadams4    . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out)
1169bae903cbSmarkadams4 
1170c8b0795cSMark F. Adams   Output Parameter:
11710cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
1172bae903cbSmarkadams4 
1173c8b0795cSMark F. Adams */
1174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists)
1175d71ae5a4SJacob Faibussowitsch {
1176e0940f08SMark F. Adams   PC_MG       *mg          = (PC_MG *)a_pc->data;
1177c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
1178c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
11798926f930SMark Adams   Mat          Gmat2, Gmat1 = *a_Gmat1; /* aggressive graph */
11800cbbd2e1SMark F. Adams   IS           perm;
1181bae903cbSmarkadams4   PetscInt     Istart, Iend, Ii, nloc, bs, nn;
1182bae903cbSmarkadams4   PetscInt    *permute, *degree;
1183c8b0795cSMark F. Adams   PetscBool   *bIndexSet;
1184aea53286SMark Adams   PetscReal    hashfact;
1185e2c00dcbSBarry Smith   PetscInt     iSwapIndex;
11863b50add6SMark Adams   PetscRandom  random;
1187d529f056Smarkadams4   MPI_Comm     comm;
1188c8b0795cSMark F. Adams 
1189c8b0795cSMark F. Adams   PetscFunctionBegin;
1190d529f056Smarkadams4   PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm));
1191849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
1192bae903cbSmarkadams4   PetscCall(MatGetLocalSize(Gmat1, &nn, NULL));
11939566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Gmat1, &bs));
119463a3b9bcSJacob Faibussowitsch   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs);
1195bae903cbSmarkadams4   nloc = nn / bs;
11965cfd4bd9SMark Adams   /* get MIS aggs - randomize */
1197bae903cbSmarkadams4   PetscCall(PetscMalloc2(nloc, &permute, nloc, &degree));
11989566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nloc, &bIndexSet));
11996e09b0e3SStefano Zampini   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
12009566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random));
12019566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
1202c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
1203bae903cbSmarkadams4     PetscInt nc;
1204e0b7e82fSBarry Smith 
1205bae903cbSmarkadams4     PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
1206bae903cbSmarkadams4     degree[Ii] = nc;
1207bae903cbSmarkadams4     PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
1208bae903cbSmarkadams4   }
1209bae903cbSmarkadams4   for (Ii = 0; Ii < nloc; Ii++) {
12109566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(random, &hashfact));
1211aea53286SMark Adams     iSwapIndex = (PetscInt)(hashfact * nloc) % nloc;
1212c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1213c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
1214e0b7e82fSBarry Smith 
1215c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
1216c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
1217bae903cbSmarkadams4       iTemp                 = degree[iSwapIndex];
1218bae903cbSmarkadams4       degree[iSwapIndex]    = degree[Ii];
1219bae903cbSmarkadams4       degree[Ii]            = iTemp;
1220c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
1221c8b0795cSMark F. Adams     }
1222c8b0795cSMark F. Adams   }
1223d529f056Smarkadams4   // apply minimum degree ordering -- NEW
1224d529f056Smarkadams4   if (pc_gamg_agg->use_minimum_degree_ordering) { PetscCall(PetscSortIntWithArray(nloc, degree, permute)); }
12259566063dSJacob Faibussowitsch   PetscCall(PetscFree(bIndexSet));
12269566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&random));
12279566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
1228849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
1229d529f056Smarkadams4   // square graph
1230d529f056Smarkadams4   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels && pc_gamg_agg->use_aggressive_square_graph) {
1231d529f056Smarkadams4     PetscCall(PCGAMGSquareGraph_GAMG(a_pc, Gmat1, &Gmat2));
1232d529f056Smarkadams4   } else Gmat2 = Gmat1;
1233d529f056Smarkadams4   // switch to old MIS-1 for square graph
1234d529f056Smarkadams4   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) {
1235d529f056Smarkadams4     if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, pc_gamg_agg->aggressive_mis_k)); // hardwire to MIS-2
1236d529f056Smarkadams4     else PetscCall(MatCoarsenSetType(pc_gamg_agg->crs, MATCOARSENMIS));                                                                   // old MIS -- side effect
1237d529f056Smarkadams4   } else if (pc_gamg_agg->use_aggressive_square_graph && pc_gamg_agg->aggressive_coarsening_levels > 0) {                                 // we reset the MIS
1238d529f056Smarkadams4     const char *prefix;
1239e0b7e82fSBarry Smith 
1240d529f056Smarkadams4     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)a_pc, &prefix));
1241d529f056Smarkadams4     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix));
1242d529f056Smarkadams4     PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); // get the default back on non-aggressive levels when square graph switched to old MIS
1243d529f056Smarkadams4   }
1244d529f056Smarkadams4   PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat2));
12452d776b49SBarry Smith   PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE));
1246d529f056Smarkadams4   PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm));
12472d776b49SBarry Smith   PetscCall(MatCoarsenApply(pc_gamg_agg->crs));
12482d776b49SBarry Smith   PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */
1249b43b03e9SMark F. Adams 
12509566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1251bae903cbSmarkadams4   PetscCall(PetscFree2(permute, degree));
1252849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
12539f3f12c8SMark F. Adams 
12549c15c1aeSMark Adams   if (Gmat2 != Gmat1) { // square graph, we need ghosts for selected
1255d529f056Smarkadams4     PetscCoarsenData *llist = *agg_lists;
1256e0b7e82fSBarry Smith 
1257d529f056Smarkadams4     PetscCall(fixAggregatesWithSquare(a_pc, Gmat2, Gmat1, *agg_lists));
1258d529f056Smarkadams4     PetscCall(MatDestroy(&Gmat1));
1259d529f056Smarkadams4     *a_Gmat1 = Gmat2;                          /* output */
12608926f930SMark Adams     PetscCall(PetscCDSetMat(llist, *a_Gmat1)); /* Need a graph with ghosts here */
12610cbbd2e1SMark F. Adams   }
1262849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
12633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1264c8b0795cSMark F. Adams }
1265c8b0795cSMark F. Adams 
1266c8b0795cSMark F. Adams /*
1267e0b7e82fSBarry Smith  PCGAMGConstructProlongator_AGG
1268c8b0795cSMark F. Adams 
1269c8b0795cSMark F. Adams  Input Parameter:
1270c8b0795cSMark F. Adams  . pc - this
1271c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1272c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
12730cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1274c8b0795cSMark F. Adams  Output Parameter:
1275c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1276c8b0795cSMark F. Adams  */
1277e0b7e82fSBarry Smith static PetscErrorCode PCGAMGConstructProlongator_AGG(PC pc, Mat Amat, PetscCoarsenData *agg_lists, Mat *a_P_out)
1278d71ae5a4SJacob Faibussowitsch {
12792e68589bSMark F. Adams   PC_MG         *mg      = (PC_MG *)pc->data;
12802e68589bSMark F. Adams   PC_GAMG       *pc_gamg = (PC_GAMG *)mg->innerctx;
12814a2f8832SBarry Smith   const PetscInt col_bs  = pc_gamg->data_cell_cols;
1282c8b0795cSMark F. Adams   PetscInt       Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs;
12838926f930SMark Adams   Mat            Gmat, Prol;
1284d5d25401SBarry Smith   PetscMPIInt    size;
12853b4367a7SBarry Smith   MPI_Comm       comm;
1286c8b0795cSMark F. Adams   PetscReal     *data_w_ghost;
1287c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes = 0, *flid_fgid;
1288d9558ea9SBarry Smith   MatType        mtype;
12892e68589bSMark F. Adams 
12902e68589bSMark F. Adams   PetscFunctionBegin;
12919566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
129208401ef6SPierre Jolivet   PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1");
1293849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
12949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
12959566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
12969566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
12979371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
12989371c9d4SSatish Balay   my0  = Istart / bs;
129963a3b9bcSJacob Faibussowitsch   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);
1300d8b4a066SPierre Jolivet   PetscCall(PetscCDGetMat(agg_lists, &Gmat)); // get auxiliary matrix for ghost edges for size > 1
13012e68589bSMark F. Adams 
13022e68589bSMark F. Adams   /* get 'nLocalSelected' */
13030cbbd2e1SMark F. Adams   for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) {
1304e78576d6SMark F. Adams     PetscBool ise;
1305e0b7e82fSBarry Smith 
13060cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
13075e62d202SMark Adams     PetscCall(PetscCDIsEmptyAt(agg_lists, ii, &ise));
1308e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
13092e68589bSMark F. Adams   }
13102e68589bSMark F. Adams 
13112e68589bSMark F. Adams   /* create prolongator, create P matrix */
13129566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat, &mtype));
13139566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
13149566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE));
131544eff04eSMark Adams   PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); // should this be before MatSetSizes?
13169566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
13173742c8caSstefanozampini #if PetscDefined(HAVE_DEVICE)
13183742c8caSstefanozampini   PetscBool flg;
13193742c8caSstefanozampini   PetscCall(MatBoundToCPU(Amat, &flg));
13203742c8caSstefanozampini   PetscCall(MatBindToCPU(Prol, flg));
13213742c8caSstefanozampini   if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE));
13223742c8caSstefanozampini #endif
13239566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL));
13249566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL));
13252e68589bSMark F. Adams 
13262e68589bSMark F. Adams   /* can get all points "removed" */
13279566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &ii));
13287f66b68fSMark Adams   if (!ii) {
132963a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix));
13309566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
13310298fd71SBarry Smith     *a_P_out = NULL; /* out */
1332849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
13333ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
13342e68589bSMark F. Adams   }
133563a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs));
13369566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
13370cbbd2e1SMark F. Adams 
133863a3b9bcSJacob Faibussowitsch   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);
1339c8b0795cSMark F. Adams   myCrs0 = myCrs0 / col_bs;
134063a3b9bcSJacob Faibussowitsch   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);
13412e68589bSMark F. Adams 
13422e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
1343849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
1344bae903cbSmarkadams4   if (size > 1) { /* get ghost null space data */
13452e68589bSMark F. Adams     PetscReal *tmp_gdata, *tmp_ldata, *tp2;
1346e0b7e82fSBarry Smith 
13479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
13484a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
1349c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1350a2f3521dSMark F. Adams         PetscInt         ii, stride;
13518e3a54c0SPierre Jolivet         const PetscReal *tp = PetscSafePointerPlusOffset(pc_gamg->data, jj * bs * nloc + kk);
1352e0b7e82fSBarry Smith 
13532fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
13542fa5cd67SKarl Rupp 
13559566063dSJacob Faibussowitsch         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
1356a2f3521dSMark F. Adams 
1357bae903cbSmarkadams4         if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */
13589566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost));
1359a2f3521dSMark F. Adams           nbnodes = bs * stride;
13602e68589bSMark F. Adams         }
13618e3a54c0SPierre Jolivet         tp2 = PetscSafePointerPlusOffset(data_w_ghost, jj * bs * stride + kk);
13622fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
13639566063dSJacob Faibussowitsch         PetscCall(PetscFree(tmp_gdata));
13642e68589bSMark F. Adams       }
13652e68589bSMark F. Adams     }
13669566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_ldata));
1367806fa848SBarry Smith   } else {
1368c8b0795cSMark F. Adams     nbnodes      = bs * nloc;
1369835f2295SStefano Zampini     data_w_ghost = pc_gamg->data;
13702e68589bSMark F. Adams   }
13712e68589bSMark F. Adams 
1372bae903cbSmarkadams4   /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */
1373c5df96a5SBarry Smith   if (size > 1) {
13742e68589bSMark F. Adams     PetscReal *fid_glid_loc, *fiddata;
1375a2f3521dSMark F. Adams     PetscInt   stride;
13762e68589bSMark F. Adams 
13779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
13782e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk);
13799566063dSJacob Faibussowitsch     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
1380bae903cbSmarkadams4     PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */
1381a2f3521dSMark F. Adams     for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
13829566063dSJacob Faibussowitsch     PetscCall(PetscFree(fiddata));
1383a2f3521dSMark F. Adams 
138463a3b9bcSJacob Faibussowitsch     PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs);
13859566063dSJacob Faibussowitsch     PetscCall(PetscFree(fid_glid_loc));
1386806fa848SBarry Smith   } else {
13879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &flid_fgid));
13882e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk;
13892e68589bSMark F. Adams   }
1390849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
1391bae903cbSmarkadams4   /* get P0 */
1392849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
1393c8b0795cSMark F. Adams   {
13940298fd71SBarry Smith     PetscReal *data_out = NULL;
1395e0b7e82fSBarry Smith 
13969566063dSJacob Faibussowitsch     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol));
13979566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
13982fa5cd67SKarl Rupp 
1399c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
14004a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
14014a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs * col_bs * nLocalSelected;
1402c8b0795cSMark F. Adams   }
1403849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
140448a46eb9SPierre Jolivet   if (size > 1) PetscCall(PetscFree(data_w_ghost));
14059566063dSJacob Faibussowitsch   PetscCall(PetscFree(flid_fgid));
14062e68589bSMark F. Adams 
1407c8b0795cSMark F. Adams   *a_P_out = Prol; /* out */
1408e0b7e82fSBarry Smith   PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_initial_prolongation"));
1409ed4e82eaSPeter Brune 
1410849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
14113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1412c8b0795cSMark F. Adams }
1413c8b0795cSMark F. Adams 
1414c8b0795cSMark F. Adams /*
1415e0b7e82fSBarry Smith    PCGAMGOptimizeProlongator_AGG - given the initial prolongator optimizes it by smoothed aggregation pc_gamg_agg->nsmooths times
1416c8b0795cSMark F. Adams 
1417c8b0795cSMark F. Adams   Input Parameter:
1418c8b0795cSMark F. Adams    . pc - this
1419c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1420c8b0795cSMark F. Adams  In/Output Parameter:
14212a808120SBarry Smith    . a_P - prolongation operator to the next level
1422c8b0795cSMark F. Adams */
1423e0b7e82fSBarry Smith static PetscErrorCode PCGAMGOptimizeProlongator_AGG(PC pc, Mat Amat, Mat *a_P)
1424d71ae5a4SJacob Faibussowitsch {
1425c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
1426c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
1427c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
1428c8b0795cSMark F. Adams   PetscInt     jj;
1429c8b0795cSMark F. Adams   Mat          Prol = *a_P;
14303b4367a7SBarry Smith   MPI_Comm     comm;
14312a808120SBarry Smith   KSP          eksp;
14322a808120SBarry Smith   Vec          bb, xx;
14332a808120SBarry Smith   PC           epc;
14342a808120SBarry Smith   PetscReal    alpha, emax, emin;
1435c8b0795cSMark F. Adams 
1436c8b0795cSMark F. Adams   PetscFunctionBegin;
14379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
1438849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
1439c8b0795cSMark F. Adams 
1440d5d25401SBarry Smith   /* compute maximum singular value of operator to be used in smoother */
14412a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
144218c3aa7eSMark     /* get eigen estimates */
144318c3aa7eSMark     if (pc_gamg->emax > 0) {
144418c3aa7eSMark       emin = pc_gamg->emin;
144518c3aa7eSMark       emax = pc_gamg->emax;
144618c3aa7eSMark     } else {
14470ed2132dSStefano Zampini       const char *prefix;
14480ed2132dSStefano Zampini 
14499566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &bb, NULL));
14509566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &xx, NULL));
145120654ebbSStefano Zampini       PetscCall(KSPSetNoisy_Private(Amat, bb));
1452e696ed0bSMark F. Adams 
14539566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm, &eksp));
14543821be0aSBarry Smith       PetscCall(KSPSetNestLevel(eksp, pc->kspnestlevel));
14559566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
14569566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(eksp, prefix));
14579566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_"));
145873f7197eSJed Brown       {
1459b94d7dedSBarry Smith         PetscBool isset, sflg;
1460e0b7e82fSBarry Smith 
1461b94d7dedSBarry Smith         PetscCall(MatIsSPDKnown(Amat, &isset, &sflg));
1462b94d7dedSBarry Smith         if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG));
1463d24ecf33SMark       }
14649566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure));
14659566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
1466c2436cacSMark F. Adams 
14679566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
14689566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(eksp, Amat, Amat));
14692e68589bSMark F. Adams 
14709566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(eksp, &epc));
14719566063dSJacob Faibussowitsch       PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */
1472c2436cacSMark F. Adams 
1473fb842aefSJose E. Roman       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
1474074aec55SMark Adams 
14759566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(eksp));
14769566063dSJacob Faibussowitsch       PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE));
14779566063dSJacob Faibussowitsch       PetscCall(KSPSolve(eksp, bb, xx));
14789566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(eksp, pc, xx));
14792e68589bSMark F. Adams 
14809566063dSJacob Faibussowitsch       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
14819566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI));
14829566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&xx));
14839566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bb));
14849566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&eksp));
14852e68589bSMark F. Adams     }
148618c3aa7eSMark     if (pc_gamg->use_sa_esteig) {
148718c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
148818c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
148963a3b9bcSJacob Faibussowitsch       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));
149018c3aa7eSMark     } else {
149118c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
149218c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
149318c3aa7eSMark     }
149418c3aa7eSMark   } else {
149518c3aa7eSMark     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
149618c3aa7eSMark     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
149718c3aa7eSMark   }
14982e68589bSMark F. Adams 
14992a808120SBarry Smith   /* smooth P0 */
1500e0b7e82fSBarry Smith   if (pc_gamg_agg->nsmooths > 0) {
15012a808120SBarry Smith     Vec diag;
15022a808120SBarry Smith 
1503e0b7e82fSBarry Smith     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
1504e0b7e82fSBarry Smith     PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero");
15052a808120SBarry Smith 
1506e0b7e82fSBarry Smith     PetscCall(MatCreateVecs(Amat, &diag, NULL));
1507e0b7e82fSBarry Smith     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
1508e0b7e82fSBarry Smith     PetscCall(VecReciprocal(diag));
1509e0b7e82fSBarry Smith 
1510e0b7e82fSBarry Smith     for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
1511e0b7e82fSBarry Smith       Mat tMat;
1512e0b7e82fSBarry Smith 
1513e0b7e82fSBarry Smith       PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
1514e0b7e82fSBarry Smith       /*
1515e0b7e82fSBarry Smith         Smooth aggregation on the prolongator
1516e0b7e82fSBarry Smith 
1517e0b7e82fSBarry Smith         P_{i} := (I - 1.4/emax D^{-1}A) P_i\{i-1}
1518e0b7e82fSBarry Smith       */
15199566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
1520fb842aefSJose E. Roman       PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_CURRENT, &tMat));
15219566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
15229566063dSJacob Faibussowitsch       PetscCall(MatProductClear(tMat));
15239566063dSJacob Faibussowitsch       PetscCall(MatDiagonalScale(tMat, diag, NULL));
1524b4da3a1bSStefano Zampini 
1525d5d25401SBarry Smith       /* TODO: Document the 1.4 and don't hardwire it in this routine */
1526b4ec6429SMark F. Adams       alpha = -1.4 / emax;
15279566063dSJacob Faibussowitsch       PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
15289566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Prol));
15292e68589bSMark F. Adams       Prol = tMat;
1530849bee69SMark Adams       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
15312e68589bSMark F. Adams     }
1532e0b7e82fSBarry Smith     PetscCall(VecDestroy(&diag));
1533e0b7e82fSBarry Smith   }
1534849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
1535e0b7e82fSBarry Smith   PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_prolongation"));
1536c8b0795cSMark F. Adams   *a_P = Prol;
15373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15382e68589bSMark F. Adams }
15392e68589bSMark F. Adams 
1540a077d33dSBarry Smith /*MC
1541a077d33dSBarry Smith   PCGAMGAGG - Smooth aggregation, {cite}`vanek1996algebraic`, {cite}`vanek2001convergence`, variant of PETSc's algebraic multigrid (`PCGAMG`) preconditioner
15422e68589bSMark F. Adams 
1543a077d33dSBarry Smith   Options Database Keys:
1544a077d33dSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation to construct prolongation
1545a077d33dSBarry Smith . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest.
1546*aea0ef14SMark Adams . -pc_gamg_aggressive_square_graph <bool,default=true> - Use square graph (A'A), alternative is MIS-k (k=2), for aggressive coarsening
1547*aea0ef14SMark Adams . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=false> - Use minimum degree ordering in greedy MIS algorithm
1548a077d33dSBarry Smith . -pc_gamg_pc_gamg_asm_hem_aggs <n,default=0> - Number of HEM aggregation steps for ASM smoother
1549a077d33dSBarry Smith - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive')
1550a077d33dSBarry Smith 
1551a077d33dSBarry Smith   Level: intermediate
1552a077d33dSBarry Smith 
1553a077d33dSBarry Smith   Notes:
1554a077d33dSBarry Smith   To obtain good performance for `PCGAMG` for vector valued problems you must
1555a077d33dSBarry Smith   call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point.
1556a077d33dSBarry Smith   Call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator
1557a077d33dSBarry Smith 
1558a077d33dSBarry Smith   The many options for `PCMG` and `PCGAMG` such as controlling the smoothers on each level etc. also work for `PCGAMGAGG`
1559a077d33dSBarry Smith 
1560a077d33dSBarry Smith .seealso: `PCGAMG`, [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCCreate()`, `PCSetType()`,
1561a077d33dSBarry Smith           `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`,
1562a077d33dSBarry Smith           `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`,
1563a077d33dSBarry Smith           `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`,
1564a077d33dSBarry Smith           `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()`
1565a077d33dSBarry Smith M*/
1566d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc)
1567d71ae5a4SJacob Faibussowitsch {
1568c8b0795cSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
1569c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
1570c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg;
15712e68589bSMark F. Adams 
1572c8b0795cSMark F. Adams   PetscFunctionBegin;
1573c8b0795cSMark F. Adams   /* create sub context for SA */
15744dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg_agg));
1575c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1576c8b0795cSMark F. Adams 
15771ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
15789b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1579c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1580c8b0795cSMark F. Adams 
1581c8b0795cSMark F. Adams   /* set internal function pointers */
15822d776b49SBarry Smith   pc_gamg->ops->creategraph       = PCGAMGCreateGraph_AGG;
15831ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
1584e0b7e82fSBarry Smith   pc_gamg->ops->prolongator       = PCGAMGConstructProlongator_AGG;
1585e0b7e82fSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptimizeProlongator_AGG;
15861ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
15875adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
1588c8b0795cSMark F. Adams 
1589cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths                     = 1;
1590d529f056Smarkadams4   pc_gamg_agg->aggressive_coarsening_levels = 1;
1591d4adc10fSMark Adams   pc_gamg_agg->use_aggressive_square_graph  = PETSC_TRUE;
1592d529f056Smarkadams4   pc_gamg_agg->use_minimum_degree_ordering  = PETSC_FALSE;
1593a9ccf005SMark Adams   pc_gamg_agg->use_low_mem_filter           = PETSC_FALSE;
1594d529f056Smarkadams4   pc_gamg_agg->aggressive_mis_k             = 2;
15956c34c54dSStefano Zampini   pc_gamg_agg->graph_symmetrize             = PETSC_TRUE;
1596cab9ed1eSBarry Smith 
15979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG));
1598bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG));
1599d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", PCGAMGSetAggressiveSquareGraph_AGG));
1600d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", PCGAMGMISkSetMinDegreeOrdering_AGG));
1601a9ccf005SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", PCGAMGSetLowMemoryFilter_AGG));
1602d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", PCGAMGMISkSetAggressive_AGG));
16036c34c54dSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetGraphSymmetrize_C", PCGAMGSetGraphSymmetrize_AGG));
16049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG));
16053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16062e68589bSMark F. Adams }
1607