xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision b2b5984a49e8f33c239e9c22df36548e528714e5)
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:
31a077d33dSBarry Smith . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use
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:
74a077d33dSBarry Smith . -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels on which to square the graph on before aggregating it
75af3c827dSMark Adams 
76ef4ad70eSMark F. Adams   Level: intermediate
77ef4ad70eSMark F. Adams 
78a077d33dSBarry 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()`
79ef4ad70eSMark F. Adams @*/
80d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n)
81d71ae5a4SJacob Faibussowitsch {
82ef4ad70eSMark F. Adams   PetscFunctionBegin;
83ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
84d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
85bae903cbSmarkadams4   PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n));
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87ef4ad70eSMark F. Adams }
88ef4ad70eSMark F. Adams 
89d529f056Smarkadams4 /*@
90d529f056Smarkadams4   PCGAMGMISkSetAggressive - Number (k) distance in MIS coarsening (>2 is 'aggressive')
91d529f056Smarkadams4 
92d529f056Smarkadams4   Logically Collective
93d529f056Smarkadams4 
94d529f056Smarkadams4   Input Parameters:
95d529f056Smarkadams4 + pc - the preconditioner context
96d529f056Smarkadams4 - n  - 1 or more (default = 2)
97d529f056Smarkadams4 
98d529f056Smarkadams4   Options Database Key:
99d529f056Smarkadams4 . -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive')
100d529f056Smarkadams4 
101d529f056Smarkadams4   Level: intermediate
102d529f056Smarkadams4 
103a077d33dSBarry 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()`
104d529f056Smarkadams4 @*/
105d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetAggressive(PC pc, PetscInt n)
106d529f056Smarkadams4 {
107d529f056Smarkadams4   PetscFunctionBegin;
108d529f056Smarkadams4   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
109d529f056Smarkadams4   PetscValidLogicalCollectiveInt(pc, n, 2);
110d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGMISkSetAggressive_C", (PC, PetscInt), (pc, n));
111d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
112d529f056Smarkadams4 }
113d529f056Smarkadams4 
114d529f056Smarkadams4 /*@
115d529f056Smarkadams4   PCGAMGSetAggressiveSquareGraph - Use graph square A'A for aggressive coarsening, old method
116d529f056Smarkadams4 
117d529f056Smarkadams4   Logically Collective
118d529f056Smarkadams4 
119d529f056Smarkadams4   Input Parameters:
120d529f056Smarkadams4 + pc - the preconditioner context
121d529f056Smarkadams4 - b  - default false - MIS-k is faster
122d529f056Smarkadams4 
123d529f056Smarkadams4   Options Database Key:
124d529f056Smarkadams4 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening
125d529f056Smarkadams4 
126d529f056Smarkadams4   Level: intermediate
127d529f056Smarkadams4 
128a077d33dSBarry 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()`
129d529f056Smarkadams4 @*/
130d529f056Smarkadams4 PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC pc, PetscBool b)
131d529f056Smarkadams4 {
132d529f056Smarkadams4   PetscFunctionBegin;
133d529f056Smarkadams4   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
134d529f056Smarkadams4   PetscValidLogicalCollectiveBool(pc, b, 2);
135d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGSetAggressiveSquareGraph_C", (PC, PetscBool), (pc, b));
136d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
137d529f056Smarkadams4 }
138d529f056Smarkadams4 
139d529f056Smarkadams4 /*@
140d529f056Smarkadams4   PCGAMGMISkSetMinDegreeOrdering - Use minimum degree ordering in greedy MIS algorithm
141d529f056Smarkadams4 
142d529f056Smarkadams4   Logically Collective
143d529f056Smarkadams4 
144d529f056Smarkadams4   Input Parameters:
145d529f056Smarkadams4 + pc - the preconditioner context
146d529f056Smarkadams4 - b  - default true
147d529f056Smarkadams4 
148d529f056Smarkadams4   Options Database Key:
149d529f056Smarkadams4 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm
150d529f056Smarkadams4 
151d529f056Smarkadams4   Level: intermediate
152d529f056Smarkadams4 
153a077d33dSBarry 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()`
154d529f056Smarkadams4 @*/
155d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC pc, PetscBool b)
156d529f056Smarkadams4 {
157d529f056Smarkadams4   PetscFunctionBegin;
158d529f056Smarkadams4   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
159d529f056Smarkadams4   PetscValidLogicalCollectiveBool(pc, b, 2);
160d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGMISkSetMinDegreeOrdering_C", (PC, PetscBool), (pc, b));
161d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
162d529f056Smarkadams4 }
163d529f056Smarkadams4 
164a9ccf005SMark Adams /*@
165a9ccf005SMark Adams   PCGAMGSetLowMemoryFilter - Use low memory graph/matrix filter
166a9ccf005SMark Adams 
167a9ccf005SMark Adams   Logically Collective
168a9ccf005SMark Adams 
169a9ccf005SMark Adams   Input Parameters:
170a9ccf005SMark Adams + pc - the preconditioner context
171a9ccf005SMark Adams - b  - default false
172a9ccf005SMark Adams 
173a9ccf005SMark Adams   Options Database Key:
174a9ccf005SMark Adams . -pc_gamg_low_memory_threshold_filter <bool,default=false> - Use low memory graph/matrix filter
175a9ccf005SMark Adams 
176a9ccf005SMark Adams   Level: intermediate
177a9ccf005SMark Adams 
178a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`,
179a077d33dSBarry Smith   `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`
180a9ccf005SMark Adams @*/
181a9ccf005SMark Adams PetscErrorCode PCGAMGSetLowMemoryFilter(PC pc, PetscBool b)
182a9ccf005SMark Adams {
183a9ccf005SMark Adams   PetscFunctionBegin;
184a9ccf005SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
185a9ccf005SMark Adams   PetscValidLogicalCollectiveBool(pc, b, 2);
186a9ccf005SMark Adams   PetscTryMethod(pc, "PCGAMGSetLowMemoryFilter_C", (PC, PetscBool), (pc, b));
187a9ccf005SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
188a9ccf005SMark Adams }
189a9ccf005SMark Adams 
1906c34c54dSStefano Zampini /*@
1916c34c54dSStefano Zampini   PCGAMGSetGraphSymmetrize - Set the flag to symmetrize the graph used in coarsening
1926c34c54dSStefano Zampini 
1936c34c54dSStefano Zampini   Logically Collective
1946c34c54dSStefano Zampini 
1956c34c54dSStefano Zampini   Input Parameters:
1966c34c54dSStefano Zampini + pc - the preconditioner context
1976c34c54dSStefano Zampini - b  - default false
1986c34c54dSStefano Zampini 
1996c34c54dSStefano Zampini   Options Database Key:
2006c34c54dSStefano Zampini . -pc_gamg_graph_symmetrize <bool,default=false> - Symmetrize the graph
2016c34c54dSStefano Zampini 
2026c34c54dSStefano Zampini   Level: intermediate
2036c34c54dSStefano Zampini 
2046c34c54dSStefano Zampini .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`,
2056c34c54dSStefano Zampini   `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`
2066c34c54dSStefano Zampini @*/
2076c34c54dSStefano Zampini PetscErrorCode PCGAMGSetGraphSymmetrize(PC pc, PetscBool b)
2086c34c54dSStefano Zampini {
2096c34c54dSStefano Zampini   PetscFunctionBegin;
2106c34c54dSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2116c34c54dSStefano Zampini   PetscValidLogicalCollectiveBool(pc, b, 2);
2126c34c54dSStefano Zampini   PetscTryMethod(pc, "PCGAMGSetGraphSymmetrize_C", (PC, PetscBool), (pc, b));
2136c34c54dSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
2146c34c54dSStefano Zampini }
2156c34c54dSStefano Zampini 
216d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n)
217d71ae5a4SJacob Faibussowitsch {
218ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
219ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
220ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
221ef4ad70eSMark F. Adams 
222ef4ad70eSMark F. Adams   PetscFunctionBegin;
223bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = n;
2243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
225ef4ad70eSMark F. Adams }
226ef4ad70eSMark F. Adams 
227d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetAggressive_AGG(PC pc, PetscInt n)
228d529f056Smarkadams4 {
229d529f056Smarkadams4   PC_MG       *mg          = (PC_MG *)pc->data;
230d529f056Smarkadams4   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
231d529f056Smarkadams4   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
232d529f056Smarkadams4 
233d529f056Smarkadams4   PetscFunctionBegin;
234d529f056Smarkadams4   pc_gamg_agg->aggressive_mis_k = n;
235d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
236d529f056Smarkadams4 }
237d529f056Smarkadams4 
238d529f056Smarkadams4 static PetscErrorCode PCGAMGSetAggressiveSquareGraph_AGG(PC pc, PetscBool b)
239d529f056Smarkadams4 {
240d529f056Smarkadams4   PC_MG       *mg          = (PC_MG *)pc->data;
241d529f056Smarkadams4   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
242d529f056Smarkadams4   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
243d529f056Smarkadams4 
244d529f056Smarkadams4   PetscFunctionBegin;
245d529f056Smarkadams4   pc_gamg_agg->use_aggressive_square_graph = b;
246d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
247d529f056Smarkadams4 }
248d529f056Smarkadams4 
249a9ccf005SMark Adams static PetscErrorCode PCGAMGSetLowMemoryFilter_AGG(PC pc, PetscBool b)
250a9ccf005SMark Adams {
251a9ccf005SMark Adams   PC_MG       *mg          = (PC_MG *)pc->data;
252a9ccf005SMark Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
253a9ccf005SMark Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
254a9ccf005SMark Adams 
255a9ccf005SMark Adams   PetscFunctionBegin;
256a9ccf005SMark Adams   pc_gamg_agg->use_low_mem_filter = b;
257a9ccf005SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
258a9ccf005SMark Adams }
259a9ccf005SMark Adams 
2606c34c54dSStefano Zampini static PetscErrorCode PCGAMGSetGraphSymmetrize_AGG(PC pc, PetscBool b)
2616c34c54dSStefano Zampini {
2626c34c54dSStefano Zampini   PC_MG       *mg          = (PC_MG *)pc->data;
2636c34c54dSStefano Zampini   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
2646c34c54dSStefano Zampini   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
2656c34c54dSStefano Zampini 
2666c34c54dSStefano Zampini   PetscFunctionBegin;
2676c34c54dSStefano Zampini   pc_gamg_agg->graph_symmetrize = b;
2686c34c54dSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
2696c34c54dSStefano Zampini }
2706c34c54dSStefano Zampini 
271d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetMinDegreeOrdering_AGG(PC pc, PetscBool b)
272d529f056Smarkadams4 {
273d529f056Smarkadams4   PC_MG       *mg          = (PC_MG *)pc->data;
274d529f056Smarkadams4   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
275d529f056Smarkadams4   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
276d529f056Smarkadams4 
277d529f056Smarkadams4   PetscFunctionBegin;
278d529f056Smarkadams4   pc_gamg_agg->use_minimum_degree_ordering = b;
279d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
280d529f056Smarkadams4 }
281d529f056Smarkadams4 
282ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems PetscOptionsObject)
283d71ae5a4SJacob Faibussowitsch {
2842e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
2852e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
286c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
287d4adc10fSMark 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;
288d4adc10fSMark Adams   PetscInt     nsq_graph_old = 0;
2892e68589bSMark F. Adams 
2902e68589bSMark F. Adams   PetscFunctionBegin;
291d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options");
292a077d33dSBarry 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));
293d4adc10fSMark Adams   // aggressive coarsening logic with deprecated -pc_gamg_square_graph
294d4adc10fSMark 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));
295d4adc10fSMark Adams   if (!n_aggressive_flg)
296d4adc10fSMark 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));
297d4adc10fSMark Adams   PetscCall(PetscOptionsBool("-pc_gamg_aggressive_square_graph", "Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening", "PCGAMGSetAggressiveSquareGraph", new_sqr_graph, &pc_gamg_agg->use_aggressive_square_graph, &new_sq_provided));
298d4adc10fSMark Adams   if (!new_sq_provided && old_sq_provided) {
299d4adc10fSMark Adams     pc_gamg_agg->aggressive_coarsening_levels = nsq_graph_old; // could be zero
300d4adc10fSMark Adams     pc_gamg_agg->use_aggressive_square_graph  = PETSC_TRUE;
301bae903cbSmarkadams4   }
302d4adc10fSMark Adams   if (new_sq_provided && old_sq_provided)
303835f2295SStefano 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));
304d529f056Smarkadams4   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));
305a9ccf005SMark 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));
306d529f056Smarkadams4   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));
3076c34c54dSStefano Zampini   PetscCall(PetscOptionsBool("-pc_gamg_graph_symmetrize", "Symmetrize graph for coarsening", "PCGAMGSetGraphSymmetrize", pc_gamg_agg->graph_symmetrize, &pc_gamg_agg->graph_symmetrize, NULL));
308d0609cedSBarry Smith   PetscOptionsHeadEnd();
3093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3102e68589bSMark F. Adams }
3112e68589bSMark F. Adams 
312d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
313d71ae5a4SJacob Faibussowitsch {
3142e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
3152e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
316e0b7e82fSBarry Smith   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
3172e68589bSMark F. Adams 
3182e68589bSMark F. Adams   PetscFunctionBegin;
319e0b7e82fSBarry Smith   PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs));
3209566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
3212e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL));
322bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL));
323d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", NULL));
324d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", NULL));
325a9ccf005SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", NULL));
326d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", NULL));
3276c34c54dSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetGraphSymmetrize_C", NULL));
328bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
3293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3302e68589bSMark F. Adams }
3312e68589bSMark F. Adams 
3322e68589bSMark F. Adams /*
3332e68589bSMark F. Adams    PCSetCoordinates_AGG
33420f4b53cSBarry Smith 
33520f4b53cSBarry Smith    Collective
3362e68589bSMark F. Adams 
3372e68589bSMark F. Adams    Input Parameter:
3382e68589bSMark F. Adams    . pc - the preconditioner context
339145b44c9SPierre Jolivet    . ndm - dimension of data (used for dof/vertex for Stokes)
340302f38e8SMark F. Adams    . a_nloc - number of vertices local
341302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
3422e68589bSMark F. Adams */
3431e6b0712SBarry Smith 
344d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
345d71ae5a4SJacob Faibussowitsch {
3462e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
3472e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
34869344418SMark F. Adams   PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf;
349a2f3521dSMark F. Adams   Mat      mat = pc->pmat;
3502e68589bSMark F. Adams 
3512e68589bSMark F. Adams   PetscFunctionBegin;
352a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
353a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
354302f38e8SMark F. Adams   nloc = a_nloc;
3552e68589bSMark F. Adams 
3562e68589bSMark F. Adams   /* SA: null space vectors */
3579566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &ndf));               /* this does not work for Stokes */
35869344418SMark F. Adams   if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
359a2f3521dSMark F. Adams   else if (coords) {
36063a3b9bcSJacob Faibussowitsch     PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf);
36169344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */
362ad540459SPierre 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);
363806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
36471959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
36563a3b9bcSJacob 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);
366c8b0795cSMark F. Adams   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
3672e68589bSMark F. Adams 
3687f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
3699566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
3709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
3712e68589bSMark F. Adams   }
3725e116b59SBarry Smith   /* copy data in - column-oriented */
3732e68589bSMark F. Adams   for (kk = 0; kk < nloc; kk++) {
37469344418SMark F. Adams     const PetscInt M    = nloc * pc_gamg->data_cell_rows; /* stride into data */
37569344418SMark F. Adams     PetscReal     *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */
376e0b7e82fSBarry Smith 
377c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols == 1) *data = 1.0;
3782e68589bSMark F. Adams     else {
37969344418SMark F. Adams       /* translational modes */
3802fa5cd67SKarl Rupp       for (ii = 0; ii < ndatarows; ii++) {
3812fa5cd67SKarl Rupp         for (jj = 0; jj < ndatarows; jj++) {
38269344418SMark F. Adams           if (ii == jj) data[ii * M + jj] = 1.0;
3832e68589bSMark F. Adams           else data[ii * M + jj] = 0.0;
3842fa5cd67SKarl Rupp         }
3852fa5cd67SKarl Rupp       }
38669344418SMark F. Adams 
38769344418SMark F. Adams       /* rotational modes */
3882e68589bSMark F. Adams       if (coords) {
38969344418SMark F. Adams         if (ndm == 2) {
3902e68589bSMark F. Adams           data += 2 * M;
3912e68589bSMark F. Adams           data[0] = -coords[2 * kk + 1];
3922e68589bSMark F. Adams           data[1] = coords[2 * kk];
393806fa848SBarry Smith         } else {
3942e68589bSMark F. Adams           data += 3 * M;
3959371c9d4SSatish Balay           data[0]         = 0.0;
3969371c9d4SSatish Balay           data[M + 0]     = coords[3 * kk + 2];
3979371c9d4SSatish Balay           data[2 * M + 0] = -coords[3 * kk + 1];
3989371c9d4SSatish Balay           data[1]         = -coords[3 * kk + 2];
3999371c9d4SSatish Balay           data[M + 1]     = 0.0;
4009371c9d4SSatish Balay           data[2 * M + 1] = coords[3 * kk];
4019371c9d4SSatish Balay           data[2]         = coords[3 * kk + 1];
4029371c9d4SSatish Balay           data[M + 2]     = -coords[3 * kk];
4039371c9d4SSatish Balay           data[2 * M + 2] = 0.0;
4042e68589bSMark F. Adams         }
4052e68589bSMark F. Adams       }
4062e68589bSMark F. Adams     }
4072e68589bSMark F. Adams   }
4082e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
4093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4102e68589bSMark F. Adams }
4112e68589bSMark F. Adams 
4122e68589bSMark F. Adams /*
413a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
414a2f3521dSMark F. Adams       Looks in Mat for near null space.
415a2f3521dSMark F. Adams       Does not work for Stokes
4162e68589bSMark F. Adams 
4172e68589bSMark F. Adams   Input Parameter:
4182e68589bSMark F. Adams    . pc -
419a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
4202e68589bSMark F. Adams */
421d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
422d71ae5a4SJacob Faibussowitsch {
423b8cd405aSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
424b8cd405aSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
425b8cd405aSMark F. Adams   MatNullSpace mnull;
42666f2ef4bSMark Adams 
427b817416eSBarry Smith   PetscFunctionBegin;
4289566063dSJacob Faibussowitsch   PetscCall(MatGetNearNullSpace(a_A, &mnull));
429b8cd405aSMark F. Adams   if (!mnull) {
43066f2ef4bSMark Adams     DM dm;
431e0b7e82fSBarry Smith 
4329566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
43348a46eb9SPierre Jolivet     if (!dm) PetscCall(MatGetDM(a_A, &dm));
43466f2ef4bSMark Adams     if (dm) {
43566f2ef4bSMark Adams       PetscObject deformation;
436b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
437b0d30dd6SMatthew G. Knepley 
4389566063dSJacob Faibussowitsch       PetscCall(DMGetNumFields(dm, &Nf));
439b0d30dd6SMatthew G. Knepley       if (Nf) {
4409566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, 0, NULL, &deformation));
4414c4edb6fSBarry Smith         if (deformation) {
442835f2295SStefano Zampini           PetscCall(PetscObjectQuery(deformation, "nearnullspace", (PetscObject *)&mnull));
443835f2295SStefano Zampini           if (!mnull) PetscCall(PetscObjectQuery(deformation, "nullspace", (PetscObject *)&mnull));
44466f2ef4bSMark Adams         }
44566f2ef4bSMark Adams       }
446b0d30dd6SMatthew G. Knepley     }
4474c4edb6fSBarry Smith   }
44866f2ef4bSMark Adams 
44966f2ef4bSMark Adams   if (!mnull) {
450a2f3521dSMark F. Adams     PetscInt bs, NN, MM;
451e0b7e82fSBarry Smith 
4529566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
4539566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
45463a3b9bcSJacob Faibussowitsch     PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs);
4559566063dSJacob Faibussowitsch     PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL));
456806fa848SBarry Smith   } else {
457b8cd405aSMark F. Adams     PetscReal         *nullvec;
458b8cd405aSMark F. Adams     PetscBool          has_const;
459b8cd405aSMark F. Adams     PetscInt           i, j, mlocal, nvec, bs;
460d5d25401SBarry Smith     const Vec         *vecs;
461d5d25401SBarry Smith     const PetscScalar *v;
462b817416eSBarry Smith 
4639566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &mlocal, NULL));
4649566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
465d19c4ebbSmarkadams4     for (i = 0; i < nvec; i++) {
466d19c4ebbSmarkadams4       PetscCall(VecGetLocalSize(vecs[i], &j));
467d19c4ebbSmarkadams4       PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal);
468d19c4ebbSmarkadams4     }
469a0dea87dSMark Adams     pc_gamg->data_sz = (nvec + !!has_const) * mlocal;
4709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
4719371c9d4SSatish Balay     if (has_const)
4729371c9d4SSatish Balay       for (i = 0; i < mlocal; i++) nullvec[i] = 1.0;
473b8cd405aSMark F. Adams     for (i = 0; i < nvec; i++) {
4749566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vecs[i], &v));
475b8cd405aSMark F. Adams       for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]);
4769566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vecs[i], &v));
477b8cd405aSMark F. Adams     }
478b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
479b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec + !!has_const);
4809566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
481b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
482b8cd405aSMark F. Adams   }
4833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4842e68589bSMark F. Adams }
4852e68589bSMark F. Adams 
4862e68589bSMark F. Adams /*
487bae903cbSmarkadams4   formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0
4882e68589bSMark F. Adams 
4892e68589bSMark F. Adams   Input Parameter:
4902adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
4912adfe9d3SBarry Smith    . bs - row block size
4920cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
4930cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
4942adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
4952e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
4962e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
497bae903cbSmarkadams4 
4982e68589bSMark F. Adams   Output Parameter:
4992e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
5002e68589bSMark F. Adams    . a_Prol - prolongation operator
5012e68589bSMark F. Adams */
502d71ae5a4SJacob 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)
503d71ae5a4SJacob Faibussowitsch {
504bd026e97SJed Brown   PetscInt      Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride;
5053b4367a7SBarry Smith   MPI_Comm      comm;
5062e68589bSMark F. Adams   PetscReal    *out_data;
507539c167fSBarry Smith   PetscCDIntNd *pos;
508*b2b5984aSZach Atkins   PetscHMapI    fgid_flid;
5090cbbd2e1SMark F. Adams 
5102e68589bSMark F. Adams   PetscFunctionBegin;
5119566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
5129566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
5139371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
5149371c9d4SSatish Balay   my0  = Istart / bs;
51563a3b9bcSJacob 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);
5160cbbd2e1SMark F. Adams   Iend /= bs;
5170cbbd2e1SMark F. Adams   nghosts = data_stride / bs - nloc;
5182e68589bSMark F. Adams 
519*b2b5984aSZach Atkins   PetscCall(PetscHMapICreateWithSize(2 * nghosts + 1, &fgid_flid));
520*b2b5984aSZach Atkins 
521*b2b5984aSZach Atkins   for (kk = 0; kk < nghosts; kk++) PetscCall(PetscHMapISet(fgid_flid, flid_fgid[nloc + kk], nloc + kk));
5222e68589bSMark F. Adams 
5230cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
5240cbbd2e1SMark F. Adams   for (nSelected = mm = 0; mm < nloc; mm++) {
525e78576d6SMark F. Adams     PetscBool ise;
526e0b7e82fSBarry Smith 
5275e62d202SMark Adams     PetscCall(PetscCDIsEmptyAt(agg_llists, mm, &ise));
528e78576d6SMark F. Adams     if (!ise) nSelected++;
5290cbbd2e1SMark F. Adams   }
5309566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj));
53163a3b9bcSJacob Faibussowitsch   PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT "  != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs);
53263a3b9bcSJacob 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);
5330cbbd2e1SMark F. Adams 
5342e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
5350cbbd2e1SMark F. Adams   out_data_stride = nSelected * nSAvec;
5362fa5cd67SKarl Rupp 
5379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data));
53833812677SSatish Balay   for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL;
5390cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
5402e68589bSMark F. Adams 
5412e68589bSMark F. Adams   /* find points and set prolongation */
542c8b0795cSMark F. Adams   minsz = 100;
5430cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
5445e62d202SMark Adams     PetscCall(PetscCDCountAt(agg_llists, mm, &jj));
545e78576d6SMark F. Adams     if (jj > 0) {
5460cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
5470cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
5486497c311SBarry Smith       PetscBLASInt   asz, M, N, INFO;
5496497c311SBarry Smith       PetscBLASInt   Mdata, LDA, LWORK;
5502e68589bSMark F. Adams       PetscScalar   *qqc, *qqr, *TAU, *WORK;
5512e68589bSMark F. Adams       PetscInt      *fids;
55265d7b583SSatish Balay       PetscReal     *data;
553b817416eSBarry Smith 
5546497c311SBarry Smith       PetscCall(PetscBLASIntCast(jj, &asz));
5556497c311SBarry Smith       PetscCall(PetscBLASIntCast(asz * bs, &M));
5566497c311SBarry Smith       PetscCall(PetscBLASIntCast(nSAvec, &N));
5576497c311SBarry Smith       PetscCall(PetscBLASIntCast(M + ((N - M > 0) ? N - M : 0), &Mdata));
5586497c311SBarry Smith       PetscCall(PetscBLASIntCast(Mdata, &LDA));
5596497c311SBarry Smith       PetscCall(PetscBLASIntCast(N * bs, &LWORK));
5600cbbd2e1SMark F. Adams       /* count agg */
5610cbbd2e1SMark F. Adams       if (asz < minsz) minsz = asz;
5620cbbd2e1SMark F. Adams 
5630cbbd2e1SMark F. Adams       /* get block */
5649566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids));
5652e68589bSMark F. Adams 
5662e68589bSMark F. Adams       aggID = 0;
5679566063dSJacob Faibussowitsch       PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
568e78576d6SMark F. Adams       while (pos) {
569e78576d6SMark F. Adams         PetscInt gid1;
570e0b7e82fSBarry Smith 
5719566063dSJacob Faibussowitsch         PetscCall(PetscCDIntNdGetID(pos, &gid1));
5729566063dSJacob Faibussowitsch         PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
573e78576d6SMark F. Adams 
5740cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
5750cbbd2e1SMark F. Adams         else {
576*b2b5984aSZach Atkins           PetscCall(PetscHMapIGet(fgid_flid, gid1, &flid));
57708401ef6SPierre Jolivet           PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table");
5780cbbd2e1SMark F. Adams         }
5795e116b59SBarry Smith         /* copy in B_i matrix - column-oriented */
58065d7b583SSatish Balay         data = &data_in[flid * bs];
5813d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
5822e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
5830cbbd2e1SMark F. Adams             PetscReal d = data[jj * data_stride + ii];
584e0b7e82fSBarry Smith 
5850cbbd2e1SMark F. Adams             qqc[jj * Mdata + aggID * bs + ii] = d;
5862e68589bSMark F. Adams           }
5872e68589bSMark F. Adams         }
5882e68589bSMark F. Adams         /* set fine IDs */
5892e68589bSMark F. Adams         for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk;
5902e68589bSMark F. Adams         aggID++;
5910cbbd2e1SMark F. Adams       }
5922e68589bSMark F. Adams 
5932e68589bSMark F. Adams       /* pad with zeros */
5942e68589bSMark F. Adams       for (ii = asz * bs; ii < Mdata; ii++) {
595ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0;
5962e68589bSMark F. Adams       }
5972e68589bSMark F. Adams 
5982e68589bSMark F. Adams       /* QR */
5999566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
600792fecdfSBarry Smith       PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
6019566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
60208401ef6SPierre Jolivet       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error");
6035e116b59SBarry Smith       /* get R - column-oriented - output B_{i+1} */
6042e68589bSMark F. Adams       {
6052e68589bSMark F. Adams         PetscReal *data = &out_data[clid * nSAvec];
606e0b7e82fSBarry Smith 
6072e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
6082e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
60908401ef6SPierre 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);
6100cbbd2e1SMark F. Adams             if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]);
6110cbbd2e1SMark F. Adams             else data[jj * out_data_stride + ii] = 0.;
6122e68589bSMark F. Adams           }
6132e68589bSMark F. Adams         }
6142e68589bSMark F. Adams       }
6152e68589bSMark F. Adams 
6165e116b59SBarry Smith       /* get Q - row-oriented */
617792fecdfSBarry Smith       PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
61863a3b9bcSJacob Faibussowitsch       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO);
6192e68589bSMark F. Adams 
6202e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
621ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii];
6222e68589bSMark F. Adams       }
6232e68589bSMark F. Adams 
6242e68589bSMark F. Adams       /* add diagonal block of P0 */
6259371c9d4SSatish Balay       for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ }
6269566063dSJacob Faibussowitsch       PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES));
6279566063dSJacob Faibussowitsch       PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids));
628b43b03e9SMark F. Adams       clid++;
6290cbbd2e1SMark F. Adams     } /* coarse agg */
6300cbbd2e1SMark F. Adams   } /* for all fine nodes */
6319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
6329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
633*b2b5984aSZach Atkins   PetscCall(PetscHMapIDestroy(&fgid_flid));
6343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6352e68589bSMark F. Adams }
6362e68589bSMark F. Adams 
637d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer)
638d71ae5a4SJacob Faibussowitsch {
6395adeb434SBarry Smith   PC_MG       *mg          = (PC_MG *)pc->data;
6405adeb434SBarry Smith   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
6415adeb434SBarry Smith   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
6425adeb434SBarry Smith 
6435adeb434SBarry Smith   PetscFunctionBegin;
6449566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      AGG specific options\n"));
645835f2295SStefano Zampini   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number of levels of aggressive coarsening %" PetscInt_FMT "\n", pc_gamg_agg->aggressive_coarsening_levels));
646d529f056Smarkadams4   if (pc_gamg_agg->aggressive_coarsening_levels > 0) {
647d529f056Smarkadams4     PetscCall(PetscViewerASCIIPrintf(viewer, "        %s aggressive coarsening\n", !pc_gamg_agg->use_aggressive_square_graph ? "MIS-k" : "Square graph"));
648835f2295SStefano 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));
649d529f056Smarkadams4   }
650e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
651e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
652e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
653e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
6540acf423eSBarry Smith   if (pc_gamg_agg->crs) PetscCall(MatCoarsenView(pc_gamg_agg->crs, viewer));
6550acf423eSBarry Smith   else PetscCall(PetscViewerASCIIPrintf(viewer, "Coarsening algorithm not yet selected\n"));
656e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
657e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
658e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
659e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
660835f2295SStefano Zampini   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number smoothing steps to construct prolongation %" PetscInt_FMT "\n", pc_gamg_agg->nsmooths));
6613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6625adeb434SBarry Smith }
6635adeb434SBarry Smith 
6642d776b49SBarry Smith static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat)
665d71ae5a4SJacob Faibussowitsch {
666c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG *)pc->data;
667c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG *)mg->innerctx;
668c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
6692d776b49SBarry Smith   const PetscReal vfilter     = pc_gamg->threshold[pc_gamg->current_level];
6709c15c1aeSMark Adams   PetscBool       ishem, ismis;
6712d776b49SBarry Smith   const char     *prefix;
672d529f056Smarkadams4   MatInfo         info0, info1;
673d529f056Smarkadams4   PetscInt        bs;
674c8b0795cSMark F. Adams 
675c8b0795cSMark F. Adams   PetscFunctionBegin;
676a9ccf005SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
6772d776b49SBarry 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 */
6782d776b49SBarry Smith   /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */
679e0b7e82fSBarry Smith   PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs));
6802d776b49SBarry Smith   PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs));
6812d776b49SBarry Smith   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix));
6822d776b49SBarry Smith   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix));
683e0b7e82fSBarry Smith   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc_gamg_agg->crs, "pc_gamg_"));
6842d776b49SBarry Smith   PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs));
685e02fb3cdSMark Adams   PetscCall(MatGetBlockSize(Amat, &bs));
686e02fb3cdSMark Adams   // check for valid indices wrt bs
687e02fb3cdSMark Adams   for (int ii = 0; ii < pc_gamg_agg->crs->strength_index_size; ii++) {
688835f2295SStefano 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",
689835f2295SStefano Zampini                pc_gamg_agg->crs->strength_index[ii], bs);
690e02fb3cdSMark Adams   }
6912d776b49SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem));
6925e62d202SMark Adams   if (ishem) {
693835f2295SStefano 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));
6945e62d202SMark Adams     pc_gamg_agg->aggressive_coarsening_levels = 0;                                         // aggressive and HEM does not make sense
6955e62d202SMark Adams     PetscCall(MatCoarsenSetMaximumIterations(pc_gamg_agg->crs, pc_gamg_agg->crs->max_it)); // for code coverage
6965e62d202SMark Adams     PetscCall(MatCoarsenSetThreshold(pc_gamg_agg->crs, vfilter));                          // for code coverage
6979c15c1aeSMark Adams   } else {
6989c15c1aeSMark Adams     PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENMIS, &ismis));
6999c15c1aeSMark Adams     if (ismis && pc_gamg_agg->aggressive_coarsening_levels && !pc_gamg_agg->use_aggressive_square_graph) {
7009c15c1aeSMark Adams       PetscCall(PetscInfo(pc, "MIS and aggressive coarsening and no square graph: force square graph\n"));
7019c15c1aeSMark Adams       pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE;
7029c15c1aeSMark Adams     }
7035e62d202SMark Adams   }
704a9ccf005SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
705a9ccf005SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
706d529f056Smarkadams4   PetscCall(MatGetInfo(Amat, MAT_LOCAL, &info0)); /* global reduction */
707a9ccf005SMark Adams 
708a9ccf005SMark Adams   if (ishem || pc_gamg_agg->use_low_mem_filter) {
7096c34c54dSStefano 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));
710a9ccf005SMark Adams   } else {
7116c34c54dSStefano Zampini     // make scalar graph, symmetrize if not known to be symmetric, scale, but do not filter (expensive)
7126c34c54dSStefano 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));
713a9ccf005SMark Adams     if (vfilter >= 0) {
714a9ccf005SMark Adams       PetscInt           Istart, Iend, ncols, nnz0, nnz1, NN, MM, nloc;
715a9ccf005SMark Adams       Mat                tGmat, Gmat = *a_Gmat;
716a9ccf005SMark Adams       MPI_Comm           comm;
717a9ccf005SMark Adams       const PetscScalar *vals;
718a9ccf005SMark Adams       const PetscInt    *idx;
719a9ccf005SMark Adams       PetscInt          *d_nnz, *o_nnz, kk, *garray = NULL, *AJ, maxcols = 0;
720a9ccf005SMark Adams       MatScalar         *AA; // this is checked in graph
721a9ccf005SMark Adams       PetscBool          isseqaij;
722a9ccf005SMark Adams       Mat                a, b, c;
723a9ccf005SMark Adams       MatType            jtype;
724a9ccf005SMark Adams 
725a9ccf005SMark Adams       PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm));
726a9ccf005SMark Adams       PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isseqaij));
727a9ccf005SMark Adams       PetscCall(MatGetType(Gmat, &jtype));
728a9ccf005SMark Adams       PetscCall(MatCreate(comm, &tGmat));
729a9ccf005SMark Adams       PetscCall(MatSetType(tGmat, jtype));
730a9ccf005SMark Adams 
731a9ccf005SMark Adams       /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold?
732a9ccf005SMark Adams         Also, if the matrix is symmetric, can we skip this
733a9ccf005SMark Adams         operation? It can be very expensive on large matrices. */
734a9ccf005SMark Adams 
735a9ccf005SMark Adams       // global sizes
736a9ccf005SMark Adams       PetscCall(MatGetSize(Gmat, &MM, &NN));
737a9ccf005SMark Adams       PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend));
738a9ccf005SMark Adams       nloc = Iend - Istart;
739a9ccf005SMark Adams       PetscCall(PetscMalloc2(nloc, &d_nnz, nloc, &o_nnz));
740a9ccf005SMark Adams       if (isseqaij) {
741a9ccf005SMark Adams         a = Gmat;
742a9ccf005SMark Adams         b = NULL;
743a9ccf005SMark Adams       } else {
744a9ccf005SMark Adams         Mat_MPIAIJ *d = (Mat_MPIAIJ *)Gmat->data;
745e0b7e82fSBarry Smith 
746a9ccf005SMark Adams         a      = d->A;
747a9ccf005SMark Adams         b      = d->B;
748a9ccf005SMark Adams         garray = d->garray;
749a9ccf005SMark Adams       }
750a9ccf005SMark Adams       /* Determine upper bound on non-zeros needed in new filtered matrix */
751a9ccf005SMark Adams       for (PetscInt row = 0; row < nloc; row++) {
752a9ccf005SMark Adams         PetscCall(MatGetRow(a, row, &ncols, NULL, NULL));
753a9ccf005SMark Adams         d_nnz[row] = ncols;
754a9ccf005SMark Adams         if (ncols > maxcols) maxcols = ncols;
755a9ccf005SMark Adams         PetscCall(MatRestoreRow(a, row, &ncols, NULL, NULL));
756a9ccf005SMark Adams       }
757a9ccf005SMark Adams       if (b) {
758a9ccf005SMark Adams         for (PetscInt row = 0; row < nloc; row++) {
759a9ccf005SMark Adams           PetscCall(MatGetRow(b, row, &ncols, NULL, NULL));
760a9ccf005SMark Adams           o_nnz[row] = ncols;
761a9ccf005SMark Adams           if (ncols > maxcols) maxcols = ncols;
762a9ccf005SMark Adams           PetscCall(MatRestoreRow(b, row, &ncols, NULL, NULL));
763a9ccf005SMark Adams         }
764a9ccf005SMark Adams       }
765a9ccf005SMark Adams       PetscCall(MatSetSizes(tGmat, nloc, nloc, MM, MM));
766a9ccf005SMark Adams       PetscCall(MatSetBlockSizes(tGmat, 1, 1));
767a9ccf005SMark Adams       PetscCall(MatSeqAIJSetPreallocation(tGmat, 0, d_nnz));
768a9ccf005SMark Adams       PetscCall(MatMPIAIJSetPreallocation(tGmat, 0, d_nnz, 0, o_nnz));
769a9ccf005SMark Adams       PetscCall(MatSetOption(tGmat, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
770a9ccf005SMark Adams       PetscCall(PetscFree2(d_nnz, o_nnz));
771a9ccf005SMark Adams       PetscCall(PetscMalloc2(maxcols, &AA, maxcols, &AJ));
772a9ccf005SMark Adams       nnz0 = nnz1 = 0;
773a9ccf005SMark Adams       for (c = a, kk = 0; c && kk < 2; c = b, kk++) {
774a9ccf005SMark Adams         for (PetscInt row = 0, grow = Istart, ncol_row, jj; row < nloc; row++, grow++) {
775a9ccf005SMark Adams           PetscCall(MatGetRow(c, row, &ncols, &idx, &vals));
776a9ccf005SMark Adams           for (ncol_row = jj = 0; jj < ncols; jj++, nnz0++) {
777a9ccf005SMark Adams             PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
778a9ccf005SMark Adams             if (PetscRealPart(sv) > vfilter) {
779a9ccf005SMark Adams               PetscInt cid = idx[jj] + Istart; //diag
780e0b7e82fSBarry Smith 
781a9ccf005SMark Adams               nnz1++;
782a9ccf005SMark Adams               if (c != a) cid = garray[idx[jj]];
783a9ccf005SMark Adams               AA[ncol_row] = vals[jj];
784a9ccf005SMark Adams               AJ[ncol_row] = cid;
785a9ccf005SMark Adams               ncol_row++;
786a9ccf005SMark Adams             }
787a9ccf005SMark Adams           }
788a9ccf005SMark Adams           PetscCall(MatRestoreRow(c, row, &ncols, &idx, &vals));
789a9ccf005SMark Adams           PetscCall(MatSetValues(tGmat, 1, &grow, ncol_row, AJ, AA, INSERT_VALUES));
790a9ccf005SMark Adams         }
791a9ccf005SMark Adams       }
792a9ccf005SMark Adams       PetscCall(PetscFree2(AA, AJ));
793a9ccf005SMark Adams       PetscCall(MatAssemblyBegin(tGmat, MAT_FINAL_ASSEMBLY));
794a9ccf005SMark Adams       PetscCall(MatAssemblyEnd(tGmat, MAT_FINAL_ASSEMBLY));
795a9ccf005SMark Adams       PetscCall(MatPropagateSymmetryOptions(Gmat, tGmat)); /* Normal Mat options are not relevant ? */
796a9ccf005SMark 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));
797a9ccf005SMark Adams       PetscCall(MatViewFromOptions(tGmat, NULL, "-mat_filter_graph_view"));
798a9ccf005SMark Adams       PetscCall(MatDestroy(&Gmat));
799a9ccf005SMark Adams       *a_Gmat = tGmat;
800a9ccf005SMark Adams     }
801a9ccf005SMark Adams   }
802a9ccf005SMark Adams 
803d529f056Smarkadams4   PetscCall(MatGetInfo(*a_Gmat, MAT_LOCAL, &info1)); /* global reduction */
804d529f056Smarkadams4   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));
805849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
8063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
807c8b0795cSMark F. Adams }
808c8b0795cSMark F. Adams 
809d529f056Smarkadams4 typedef PetscInt    NState;
810d529f056Smarkadams4 static const NState NOT_DONE = -2;
811d529f056Smarkadams4 static const NState DELETED  = -1;
812d529f056Smarkadams4 static const NState REMOVED  = -3;
813d529f056Smarkadams4 #define IS_SELECTED(s) (s != DELETED && s != NOT_DONE && s != REMOVED)
814d529f056Smarkadams4 
815d529f056Smarkadams4 /*
816d529f056Smarkadams4    fixAggregatesWithSquare - greedy grab of with G1 (unsquared graph) -- AIJ specific -- change to fixAggregatesWithSquare -- TODD
817d529f056Smarkadams4      - AGG-MG specific: clears singletons out of 'selected_2'
818d529f056Smarkadams4 
819d529f056Smarkadams4    Input Parameter:
820d529f056Smarkadams4    . Gmat_2 - global matrix of squared graph (data not defined)
821d529f056Smarkadams4    . Gmat_1 - base graph to grab with base graph
822d529f056Smarkadams4    Input/Output Parameter:
823d529f056Smarkadams4    . aggs_2 - linked list of aggs with gids)
824d529f056Smarkadams4 */
825d529f056Smarkadams4 static PetscErrorCode fixAggregatesWithSquare(PC pc, Mat Gmat_2, Mat Gmat_1, PetscCoarsenData *aggs_2)
826d529f056Smarkadams4 {
827d529f056Smarkadams4   PetscBool      isMPI;
828d529f056Smarkadams4   Mat_SeqAIJ    *matA_1, *matB_1 = NULL;
829d529f056Smarkadams4   MPI_Comm       comm;
830d529f056Smarkadams4   PetscInt       lid, *ii, *idx, ix, Iend, my0, kk, n, j;
831d529f056Smarkadams4   Mat_MPIAIJ    *mpimat_2 = NULL, *mpimat_1 = NULL;
832d529f056Smarkadams4   const PetscInt nloc = Gmat_2->rmap->n;
833d529f056Smarkadams4   PetscScalar   *cpcol_1_state, *cpcol_2_state, *cpcol_2_par_orig, *lid_parent_gid;
834d529f056Smarkadams4   PetscInt      *lid_cprowID_1 = NULL;
835d529f056Smarkadams4   NState        *lid_state;
836d529f056Smarkadams4   Vec            ghost_par_orig2;
837d529f056Smarkadams4   PetscMPIInt    rank;
838d529f056Smarkadams4 
839d529f056Smarkadams4   PetscFunctionBegin;
840d529f056Smarkadams4   PetscCall(PetscObjectGetComm((PetscObject)Gmat_2, &comm));
841d529f056Smarkadams4   PetscCallMPI(MPI_Comm_rank(comm, &rank));
842d529f056Smarkadams4   PetscCall(MatGetOwnershipRange(Gmat_1, &my0, &Iend));
843d529f056Smarkadams4 
844d529f056Smarkadams4   /* get submatrices */
845d529f056Smarkadams4   PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATMPIAIJ, &isMPI));
846d529f056Smarkadams4   PetscCall(PetscInfo(pc, "isMPI = %s\n", isMPI ? "yes" : "no"));
847d529f056Smarkadams4   PetscCall(PetscMalloc3(nloc, &lid_state, nloc, &lid_parent_gid, nloc, &lid_cprowID_1));
848d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
849d529f056Smarkadams4   if (isMPI) {
850d529f056Smarkadams4     /* grab matrix objects */
851d529f056Smarkadams4     mpimat_2 = (Mat_MPIAIJ *)Gmat_2->data;
852d529f056Smarkadams4     mpimat_1 = (Mat_MPIAIJ *)Gmat_1->data;
853d529f056Smarkadams4     matA_1   = (Mat_SeqAIJ *)mpimat_1->A->data;
854d529f056Smarkadams4     matB_1   = (Mat_SeqAIJ *)mpimat_1->B->data;
855d529f056Smarkadams4 
856d529f056Smarkadams4     /* force compressed row storage for B matrix in AuxMat */
857d529f056Smarkadams4     PetscCall(MatCheckCompressedRow(mpimat_1->B, matB_1->nonzerorowcnt, &matB_1->compressedrow, matB_1->i, Gmat_1->rmap->n, -1.0));
858d529f056Smarkadams4     for (ix = 0; ix < matB_1->compressedrow.nrows; ix++) {
859d529f056Smarkadams4       PetscInt lid = matB_1->compressedrow.rindex[ix];
860e0b7e82fSBarry Smith 
861835f2295SStefano Zampini       PetscCheck(lid <= nloc && lid >= -1, PETSC_COMM_SELF, PETSC_ERR_USER, "lid %" PetscInt_FMT " out of range. nloc = %" PetscInt_FMT, lid, nloc);
862d529f056Smarkadams4       if (lid != -1) lid_cprowID_1[lid] = ix;
863d529f056Smarkadams4     }
864d529f056Smarkadams4   } else {
865d529f056Smarkadams4     PetscBool isAIJ;
866e0b7e82fSBarry Smith 
867d529f056Smarkadams4     PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATSEQAIJ, &isAIJ));
868d529f056Smarkadams4     PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix.");
869d529f056Smarkadams4     matA_1 = (Mat_SeqAIJ *)Gmat_1->data;
870d529f056Smarkadams4   }
871d529f056Smarkadams4   if (nloc > 0) { PetscCheck(!matB_1 || matB_1->compressedrow.use, PETSC_COMM_SELF, PETSC_ERR_PLIB, "matB_1 && !matB_1->compressedrow.use: PETSc bug???"); }
872d529f056Smarkadams4   /* get state of locals and selected gid for deleted */
873d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) {
874d529f056Smarkadams4     lid_parent_gid[lid] = -1.0;
875d529f056Smarkadams4     lid_state[lid]      = DELETED;
876d529f056Smarkadams4   }
877d529f056Smarkadams4 
878d529f056Smarkadams4   /* set lid_state */
879d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) {
880d529f056Smarkadams4     PetscCDIntNd *pos;
881e0b7e82fSBarry Smith 
882d529f056Smarkadams4     PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
883d529f056Smarkadams4     if (pos) {
884d529f056Smarkadams4       PetscInt gid1;
885d529f056Smarkadams4 
886d529f056Smarkadams4       PetscCall(PetscCDIntNdGetID(pos, &gid1));
887835f2295SStefano Zampini       PetscCheck(gid1 == lid + my0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gid1 %" PetscInt_FMT " != lid %" PetscInt_FMT " + my0 %" PetscInt_FMT, gid1, lid, my0);
888d529f056Smarkadams4       lid_state[lid] = gid1;
889d529f056Smarkadams4     }
890d529f056Smarkadams4   }
891d529f056Smarkadams4 
892d529f056Smarkadams4   /* map local to selected local, DELETED means a ghost owns it */
89363bfac88SBarry Smith   for (lid = 0; lid < nloc; lid++) {
894d529f056Smarkadams4     NState state = lid_state[lid];
89563bfac88SBarry Smith 
896d529f056Smarkadams4     if (IS_SELECTED(state)) {
897d529f056Smarkadams4       PetscCDIntNd *pos;
898e0b7e82fSBarry Smith 
899d529f056Smarkadams4       PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
900d529f056Smarkadams4       while (pos) {
901d529f056Smarkadams4         PetscInt gid1;
902e0b7e82fSBarry Smith 
903d529f056Smarkadams4         PetscCall(PetscCDIntNdGetID(pos, &gid1));
904d529f056Smarkadams4         PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos));
905d529f056Smarkadams4         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1 - my0] = (PetscScalar)(lid + my0);
906d529f056Smarkadams4       }
907d529f056Smarkadams4     }
908d529f056Smarkadams4   }
909d529f056Smarkadams4   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
910d529f056Smarkadams4   if (isMPI) {
911d529f056Smarkadams4     Vec tempVec;
91263bfac88SBarry Smith 
913d529f056Smarkadams4     /* get 'cpcol_1_state' */
914d529f056Smarkadams4     PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL));
915d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
916d529f056Smarkadams4       PetscScalar v = (PetscScalar)lid_state[kk];
917e0b7e82fSBarry Smith 
918d529f056Smarkadams4       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
919d529f056Smarkadams4     }
920d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
921d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
922d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD));
923d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD));
924d529f056Smarkadams4     PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state));
925d529f056Smarkadams4     /* get 'cpcol_2_state' */
926d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD));
927d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD));
928d529f056Smarkadams4     PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state));
929d529f056Smarkadams4     /* get 'cpcol_2_par_orig' */
930d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
931835f2295SStefano Zampini       PetscScalar v = lid_parent_gid[kk];
932e0b7e82fSBarry Smith 
933d529f056Smarkadams4       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
934d529f056Smarkadams4     }
935d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
936d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
937d529f056Smarkadams4     PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2));
938d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD));
939d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD));
940d529f056Smarkadams4     PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig));
941d529f056Smarkadams4 
942d529f056Smarkadams4     PetscCall(VecDestroy(&tempVec));
943d529f056Smarkadams4   } /* ismpi */
944d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) {
945d529f056Smarkadams4     NState state = lid_state[lid];
946e0b7e82fSBarry Smith 
947d529f056Smarkadams4     if (IS_SELECTED(state)) {
948d529f056Smarkadams4       /* steal locals */
949d529f056Smarkadams4       ii  = matA_1->i;
950d529f056Smarkadams4       n   = ii[lid + 1] - ii[lid];
951d529f056Smarkadams4       idx = matA_1->j + ii[lid];
952d529f056Smarkadams4       for (j = 0; j < n; j++) {
953d529f056Smarkadams4         PetscInt lidj   = idx[j], sgid;
954d529f056Smarkadams4         NState   statej = lid_state[lidj];
955e0b7e82fSBarry Smith 
956d529f056Smarkadams4         if (statej == DELETED && (sgid = (PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid + my0) { /* steal local */
957d529f056Smarkadams4           lid_parent_gid[lidj] = (PetscScalar)(lid + my0);                                              /* send this if sgid is not local */
958d529f056Smarkadams4           if (sgid >= my0 && sgid < Iend) {                                                             /* I'm stealing this local from a local sgid */
959d529f056Smarkadams4             PetscInt      hav = 0, slid = sgid - my0, gidj = lidj + my0;
960d529f056Smarkadams4             PetscCDIntNd *pos, *last = NULL;
961e0b7e82fSBarry Smith 
962d529f056Smarkadams4             /* looking for local from local so id_llist_2 works */
963d529f056Smarkadams4             PetscCall(PetscCDGetHeadPos(aggs_2, slid, &pos));
964d529f056Smarkadams4             while (pos) {
965d529f056Smarkadams4               PetscInt gid;
96663bfac88SBarry Smith 
967d529f056Smarkadams4               PetscCall(PetscCDIntNdGetID(pos, &gid));
968d529f056Smarkadams4               if (gid == gidj) {
969d529f056Smarkadams4                 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null");
970d529f056Smarkadams4                 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last));
971d529f056Smarkadams4                 PetscCall(PetscCDAppendNode(aggs_2, lid, pos));
972d529f056Smarkadams4                 hav = 1;
973d529f056Smarkadams4                 break;
974d529f056Smarkadams4               } else last = pos;
975d529f056Smarkadams4               PetscCall(PetscCDGetNextPos(aggs_2, slid, &pos));
976d529f056Smarkadams4             }
977d529f056Smarkadams4             if (hav != 1) {
978d529f056Smarkadams4               PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find adj in 'selected' lists - structurally unsymmetric matrix");
979835f2295SStefano Zampini               SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %" PetscInt_FMT " times???", hav);
980d529f056Smarkadams4             }
981d529f056Smarkadams4           } else { /* I'm stealing this local, owned by a ghost */
982d529f056Smarkadams4             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.",
983d529f056Smarkadams4                        ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "");
984d529f056Smarkadams4             PetscCall(PetscCDAppendID(aggs_2, lid, lidj + my0));
985d529f056Smarkadams4           }
986d529f056Smarkadams4         }
987d529f056Smarkadams4       } /* local neighbors */
988d529f056Smarkadams4     } else if (state == DELETED /* && lid_cprowID_1 */) {
989d529f056Smarkadams4       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
990e0b7e82fSBarry Smith 
991d529f056Smarkadams4       /* see if I have a selected ghost neighbor that will steal me */
992d529f056Smarkadams4       if ((ix = lid_cprowID_1[lid]) != -1) {
993d529f056Smarkadams4         ii  = matB_1->compressedrow.i;
994d529f056Smarkadams4         n   = ii[ix + 1] - ii[ix];
995d529f056Smarkadams4         idx = matB_1->j + ii[ix];
996d529f056Smarkadams4         for (j = 0; j < n; j++) {
997d529f056Smarkadams4           PetscInt cpid   = idx[j];
998d529f056Smarkadams4           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
999e0b7e82fSBarry Smith 
1000835f2295SStefano Zampini           if (IS_SELECTED(statej) && sgidold != statej) { /* ghost will steal this, remove from my list */
1001d529f056Smarkadams4             lid_parent_gid[lid] = (PetscScalar)statej;    /* send who selected */
1002d529f056Smarkadams4             if (sgidold >= my0 && sgidold < Iend) {       /* this was mine */
1003d529f056Smarkadams4               PetscInt      hav = 0, oldslidj = sgidold - my0;
1004d529f056Smarkadams4               PetscCDIntNd *pos, *last        = NULL;
1005e0b7e82fSBarry Smith 
1006d529f056Smarkadams4               /* remove from 'oldslidj' list */
1007d529f056Smarkadams4               PetscCall(PetscCDGetHeadPos(aggs_2, oldslidj, &pos));
1008d529f056Smarkadams4               while (pos) {
1009d529f056Smarkadams4                 PetscInt gid;
1010e0b7e82fSBarry Smith 
1011d529f056Smarkadams4                 PetscCall(PetscCDIntNdGetID(pos, &gid));
1012d529f056Smarkadams4                 if (lid + my0 == gid) {
1013d529f056Smarkadams4                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
1014d529f056Smarkadams4                   PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null");
1015d529f056Smarkadams4                   PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last));
1016d529f056Smarkadams4                   /* ghost (PetscScalar)statej will add this later */
1017d529f056Smarkadams4                   hav = 1;
1018d529f056Smarkadams4                   break;
1019d529f056Smarkadams4                 } else last = pos;
1020d529f056Smarkadams4                 PetscCall(PetscCDGetNextPos(aggs_2, oldslidj, &pos));
1021d529f056Smarkadams4               }
1022d529f056Smarkadams4               if (hav != 1) {
1023835f2295SStefano Zampini                 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find (hav=%" PetscInt_FMT ") adj in 'selected' lists - structurally unsymmetric matrix", hav);
1024835f2295SStefano Zampini                 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %" PetscInt_FMT " times???", hav);
1025d529f056Smarkadams4               }
1026d529f056Smarkadams4             } else {
1027d529f056Smarkadams4               /* TODO: ghosts remove this later */
1028d529f056Smarkadams4             }
1029d529f056Smarkadams4           }
1030d529f056Smarkadams4         }
1031d529f056Smarkadams4       }
1032d529f056Smarkadams4     } /* selected/deleted */
1033d529f056Smarkadams4   } /* node loop */
1034d529f056Smarkadams4 
1035d529f056Smarkadams4   if (isMPI) {
1036d529f056Smarkadams4     PetscScalar *cpcol_2_parent, *cpcol_2_gid;
1037d529f056Smarkadams4     Vec          tempVec, ghostgids2, ghostparents2;
1038d529f056Smarkadams4     PetscInt     cpid, nghost_2;
1039*b2b5984aSZach Atkins     PetscHMapI   gid_cpid;
1040d529f056Smarkadams4 
1041d529f056Smarkadams4     PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2));
1042d529f056Smarkadams4     PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL));
1043d529f056Smarkadams4 
1044d529f056Smarkadams4     /* get 'cpcol_2_parent' */
1045d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) { PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES)); }
1046d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
1047d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
1048d529f056Smarkadams4     PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2));
1049d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD));
1050d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD));
1051d529f056Smarkadams4     PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent));
1052d529f056Smarkadams4 
1053d529f056Smarkadams4     /* get 'cpcol_2_gid' */
1054d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
1055d529f056Smarkadams4       PetscScalar v = (PetscScalar)j;
1056e0b7e82fSBarry Smith 
1057d529f056Smarkadams4       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
1058d529f056Smarkadams4     }
1059d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
1060d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
1061d529f056Smarkadams4     PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2));
1062d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD));
1063d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD));
1064d529f056Smarkadams4     PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid));
1065d529f056Smarkadams4     PetscCall(VecDestroy(&tempVec));
1066d529f056Smarkadams4 
1067d529f056Smarkadams4     /* look for deleted ghosts and add to table */
1068*b2b5984aSZach Atkins     PetscCall(PetscHMapICreateWithSize(2 * nghost_2 + 1, &gid_cpid));
1069d529f056Smarkadams4     for (cpid = 0; cpid < nghost_2; cpid++) {
1070d529f056Smarkadams4       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
1071e0b7e82fSBarry Smith 
1072d529f056Smarkadams4       if (state == DELETED) {
1073d529f056Smarkadams4         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
1074d529f056Smarkadams4         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
1075e0b7e82fSBarry Smith 
1076d529f056Smarkadams4         if (sgid_old == -1 && sgid_new != -1) {
1077d529f056Smarkadams4           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
1078e0b7e82fSBarry Smith 
1079*b2b5984aSZach Atkins           PetscCall(PetscHMapISet(gid_cpid, gid, cpid));
1080d529f056Smarkadams4         }
1081d529f056Smarkadams4       }
1082d529f056Smarkadams4     }
1083d529f056Smarkadams4 
1084d529f056Smarkadams4     /* look for deleted ghosts and see if they moved - remove it */
1085d529f056Smarkadams4     for (lid = 0; lid < nloc; lid++) {
1086d529f056Smarkadams4       NState state = lid_state[lid];
108763bfac88SBarry Smith 
1088d529f056Smarkadams4       if (IS_SELECTED(state)) {
1089d529f056Smarkadams4         PetscCDIntNd *pos, *last = NULL;
109063bfac88SBarry Smith 
1091d529f056Smarkadams4         /* look for deleted ghosts and see if they moved */
1092d529f056Smarkadams4         PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
1093d529f056Smarkadams4         while (pos) {
1094d529f056Smarkadams4           PetscInt gid;
1095d529f056Smarkadams4 
109663bfac88SBarry Smith           PetscCall(PetscCDIntNdGetID(pos, &gid));
1097d529f056Smarkadams4           if (gid < my0 || gid >= Iend) {
1098*b2b5984aSZach Atkins             PetscCall(PetscHMapIGet(gid_cpid, gid, &cpid));
1099d529f056Smarkadams4             if (cpid != -1) {
1100d529f056Smarkadams4               /* a moved ghost - */
1101d529f056Smarkadams4               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
1102d529f056Smarkadams4               PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last));
1103d529f056Smarkadams4             } else last = pos;
1104d529f056Smarkadams4           } else last = pos;
1105d529f056Smarkadams4 
1106d529f056Smarkadams4           PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos));
1107d529f056Smarkadams4         } /* loop over list of deleted */
1108d529f056Smarkadams4       } /* selected */
1109d529f056Smarkadams4     }
1110*b2b5984aSZach Atkins     PetscCall(PetscHMapIDestroy(&gid_cpid));
1111d529f056Smarkadams4 
1112d529f056Smarkadams4     /* look at ghosts, see if they changed - and it */
1113d529f056Smarkadams4     for (cpid = 0; cpid < nghost_2; cpid++) {
1114d529f056Smarkadams4       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
1115e0b7e82fSBarry Smith 
1116d529f056Smarkadams4       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
1117d529f056Smarkadams4         PetscInt      gid      = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
1118d529f056Smarkadams4         PetscInt      slid_new = sgid_new - my0, hav = 0;
1119d529f056Smarkadams4         PetscCDIntNd *pos;
1120d529f056Smarkadams4 
1121d529f056Smarkadams4         /* search for this gid to see if I have it */
1122d529f056Smarkadams4         PetscCall(PetscCDGetHeadPos(aggs_2, slid_new, &pos));
1123d529f056Smarkadams4         while (pos) {
1124d529f056Smarkadams4           PetscInt gidj;
1125e0b7e82fSBarry Smith 
1126d529f056Smarkadams4           PetscCall(PetscCDIntNdGetID(pos, &gidj));
1127d529f056Smarkadams4           PetscCall(PetscCDGetNextPos(aggs_2, slid_new, &pos));
1128d529f056Smarkadams4 
1129d529f056Smarkadams4           if (gidj == gid) {
1130d529f056Smarkadams4             hav = 1;
1131d529f056Smarkadams4             break;
1132d529f056Smarkadams4           }
1133d529f056Smarkadams4         }
1134d529f056Smarkadams4         if (hav != 1) {
1135d529f056Smarkadams4           /* insert 'flidj' into head of llist */
1136d529f056Smarkadams4           PetscCall(PetscCDAppendID(aggs_2, slid_new, gid));
1137d529f056Smarkadams4         }
1138d529f056Smarkadams4       }
1139d529f056Smarkadams4     }
1140d529f056Smarkadams4     PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state));
1141d529f056Smarkadams4     PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state));
1142d529f056Smarkadams4     PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent));
1143d529f056Smarkadams4     PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid));
1144d529f056Smarkadams4     PetscCall(VecDestroy(&ghostgids2));
1145d529f056Smarkadams4     PetscCall(VecDestroy(&ghostparents2));
1146d529f056Smarkadams4     PetscCall(VecDestroy(&ghost_par_orig2));
1147d529f056Smarkadams4   }
1148d529f056Smarkadams4   PetscCall(PetscFree3(lid_state, lid_parent_gid, lid_cprowID_1));
1149d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
1150d529f056Smarkadams4 }
1151d529f056Smarkadams4 
1152c8b0795cSMark F. Adams /*
1153bae903cbSmarkadams4    PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for
1154bae903cbSmarkadams4      communication of QR data used with HEM and MISk coarsening
1155c8b0795cSMark F. Adams 
1156c8b0795cSMark F. Adams   Input Parameter:
1157e0940f08SMark F. Adams    . a_pc - this
1158bae903cbSmarkadams4 
1159e0940f08SMark F. Adams   Input/Output Parameter:
1160bae903cbSmarkadams4    . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out)
1161bae903cbSmarkadams4 
1162c8b0795cSMark F. Adams   Output Parameter:
11630cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
1164bae903cbSmarkadams4 
1165c8b0795cSMark F. Adams */
1166d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists)
1167d71ae5a4SJacob Faibussowitsch {
1168e0940f08SMark F. Adams   PC_MG       *mg          = (PC_MG *)a_pc->data;
1169c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
1170c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
11718926f930SMark Adams   Mat          Gmat2, Gmat1 = *a_Gmat1; /* aggressive graph */
11720cbbd2e1SMark F. Adams   IS           perm;
1173bae903cbSmarkadams4   PetscInt     Istart, Iend, Ii, nloc, bs, nn;
1174bae903cbSmarkadams4   PetscInt    *permute, *degree;
1175c8b0795cSMark F. Adams   PetscBool   *bIndexSet;
1176aea53286SMark Adams   PetscReal    hashfact;
1177e2c00dcbSBarry Smith   PetscInt     iSwapIndex;
11783b50add6SMark Adams   PetscRandom  random;
1179d529f056Smarkadams4   MPI_Comm     comm;
1180c8b0795cSMark F. Adams 
1181c8b0795cSMark F. Adams   PetscFunctionBegin;
1182d529f056Smarkadams4   PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm));
1183849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
1184bae903cbSmarkadams4   PetscCall(MatGetLocalSize(Gmat1, &nn, NULL));
11859566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Gmat1, &bs));
118663a3b9bcSJacob Faibussowitsch   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs);
1187bae903cbSmarkadams4   nloc = nn / bs;
11885cfd4bd9SMark Adams   /* get MIS aggs - randomize */
1189bae903cbSmarkadams4   PetscCall(PetscMalloc2(nloc, &permute, nloc, &degree));
11909566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nloc, &bIndexSet));
11916e09b0e3SStefano Zampini   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
11929566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random));
11939566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
1194c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
1195bae903cbSmarkadams4     PetscInt nc;
1196e0b7e82fSBarry Smith 
1197bae903cbSmarkadams4     PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
1198bae903cbSmarkadams4     degree[Ii] = nc;
1199bae903cbSmarkadams4     PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
1200bae903cbSmarkadams4   }
1201bae903cbSmarkadams4   for (Ii = 0; Ii < nloc; Ii++) {
12029566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(random, &hashfact));
1203aea53286SMark Adams     iSwapIndex = (PetscInt)(hashfact * nloc) % nloc;
1204c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1205c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
1206e0b7e82fSBarry Smith 
1207c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
1208c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
1209bae903cbSmarkadams4       iTemp                 = degree[iSwapIndex];
1210bae903cbSmarkadams4       degree[iSwapIndex]    = degree[Ii];
1211bae903cbSmarkadams4       degree[Ii]            = iTemp;
1212c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
1213c8b0795cSMark F. Adams     }
1214c8b0795cSMark F. Adams   }
1215d529f056Smarkadams4   // apply minimum degree ordering -- NEW
1216d529f056Smarkadams4   if (pc_gamg_agg->use_minimum_degree_ordering) { PetscCall(PetscSortIntWithArray(nloc, degree, permute)); }
12179566063dSJacob Faibussowitsch   PetscCall(PetscFree(bIndexSet));
12189566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&random));
12199566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
1220849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
1221d529f056Smarkadams4   // square graph
1222d529f056Smarkadams4   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels && pc_gamg_agg->use_aggressive_square_graph) {
1223d529f056Smarkadams4     PetscCall(PCGAMGSquareGraph_GAMG(a_pc, Gmat1, &Gmat2));
1224d529f056Smarkadams4   } else Gmat2 = Gmat1;
1225d529f056Smarkadams4   // switch to old MIS-1 for square graph
1226d529f056Smarkadams4   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) {
1227d529f056Smarkadams4     if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, pc_gamg_agg->aggressive_mis_k)); // hardwire to MIS-2
1228d529f056Smarkadams4     else PetscCall(MatCoarsenSetType(pc_gamg_agg->crs, MATCOARSENMIS));                                                                   // old MIS -- side effect
1229d529f056Smarkadams4   } else if (pc_gamg_agg->use_aggressive_square_graph && pc_gamg_agg->aggressive_coarsening_levels > 0) {                                 // we reset the MIS
1230d529f056Smarkadams4     const char *prefix;
1231e0b7e82fSBarry Smith 
1232d529f056Smarkadams4     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)a_pc, &prefix));
1233d529f056Smarkadams4     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix));
1234d529f056Smarkadams4     PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); // get the default back on non-aggressive levels when square graph switched to old MIS
1235d529f056Smarkadams4   }
1236d529f056Smarkadams4   PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat2));
12372d776b49SBarry Smith   PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE));
1238d529f056Smarkadams4   PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm));
12392d776b49SBarry Smith   PetscCall(MatCoarsenApply(pc_gamg_agg->crs));
12402d776b49SBarry Smith   PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */
1241b43b03e9SMark F. Adams 
12429566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1243bae903cbSmarkadams4   PetscCall(PetscFree2(permute, degree));
1244849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
12459f3f12c8SMark F. Adams 
12469c15c1aeSMark Adams   if (Gmat2 != Gmat1) { // square graph, we need ghosts for selected
1247d529f056Smarkadams4     PetscCoarsenData *llist = *agg_lists;
1248e0b7e82fSBarry Smith 
1249d529f056Smarkadams4     PetscCall(fixAggregatesWithSquare(a_pc, Gmat2, Gmat1, *agg_lists));
1250d529f056Smarkadams4     PetscCall(MatDestroy(&Gmat1));
1251d529f056Smarkadams4     *a_Gmat1 = Gmat2;                          /* output */
12528926f930SMark Adams     PetscCall(PetscCDSetMat(llist, *a_Gmat1)); /* Need a graph with ghosts here */
12530cbbd2e1SMark F. Adams   }
1254849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
12553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1256c8b0795cSMark F. Adams }
1257c8b0795cSMark F. Adams 
1258c8b0795cSMark F. Adams /*
1259e0b7e82fSBarry Smith  PCGAMGConstructProlongator_AGG
1260c8b0795cSMark F. Adams 
1261c8b0795cSMark F. Adams  Input Parameter:
1262c8b0795cSMark F. Adams  . pc - this
1263c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1264c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
12650cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1266c8b0795cSMark F. Adams  Output Parameter:
1267c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1268c8b0795cSMark F. Adams  */
1269e0b7e82fSBarry Smith static PetscErrorCode PCGAMGConstructProlongator_AGG(PC pc, Mat Amat, PetscCoarsenData *agg_lists, Mat *a_P_out)
1270d71ae5a4SJacob Faibussowitsch {
12712e68589bSMark F. Adams   PC_MG         *mg      = (PC_MG *)pc->data;
12722e68589bSMark F. Adams   PC_GAMG       *pc_gamg = (PC_GAMG *)mg->innerctx;
12734a2f8832SBarry Smith   const PetscInt col_bs  = pc_gamg->data_cell_cols;
1274c8b0795cSMark F. Adams   PetscInt       Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs;
12758926f930SMark Adams   Mat            Gmat, Prol;
1276d5d25401SBarry Smith   PetscMPIInt    size;
12773b4367a7SBarry Smith   MPI_Comm       comm;
1278c8b0795cSMark F. Adams   PetscReal     *data_w_ghost;
1279c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes = 0, *flid_fgid;
1280d9558ea9SBarry Smith   MatType        mtype;
12812e68589bSMark F. Adams 
12822e68589bSMark F. Adams   PetscFunctionBegin;
12839566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
128408401ef6SPierre Jolivet   PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1");
1285849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
12869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
12879566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
12889566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
12899371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
12909371c9d4SSatish Balay   my0  = Istart / bs;
129163a3b9bcSJacob 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);
1292d8b4a066SPierre Jolivet   PetscCall(PetscCDGetMat(agg_lists, &Gmat)); // get auxiliary matrix for ghost edges for size > 1
12932e68589bSMark F. Adams 
12942e68589bSMark F. Adams   /* get 'nLocalSelected' */
12950cbbd2e1SMark F. Adams   for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) {
1296e78576d6SMark F. Adams     PetscBool ise;
1297e0b7e82fSBarry Smith 
12980cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
12995e62d202SMark Adams     PetscCall(PetscCDIsEmptyAt(agg_lists, ii, &ise));
1300e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
13012e68589bSMark F. Adams   }
13022e68589bSMark F. Adams 
13032e68589bSMark F. Adams   /* create prolongator, create P matrix */
13049566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat, &mtype));
13059566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
13069566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE));
130744eff04eSMark Adams   PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); // should this be before MatSetSizes?
13089566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
13093742c8caSstefanozampini #if PetscDefined(HAVE_DEVICE)
13103742c8caSstefanozampini   PetscBool flg;
13113742c8caSstefanozampini   PetscCall(MatBoundToCPU(Amat, &flg));
13123742c8caSstefanozampini   PetscCall(MatBindToCPU(Prol, flg));
13133742c8caSstefanozampini   if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE));
13143742c8caSstefanozampini #endif
13159566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL));
13169566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL));
13172e68589bSMark F. Adams 
13182e68589bSMark F. Adams   /* can get all points "removed" */
13199566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &ii));
13207f66b68fSMark Adams   if (!ii) {
132163a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix));
13229566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
13230298fd71SBarry Smith     *a_P_out = NULL; /* out */
1324849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
13253ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
13262e68589bSMark F. Adams   }
132763a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs));
13289566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
13290cbbd2e1SMark F. Adams 
133063a3b9bcSJacob 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);
1331c8b0795cSMark F. Adams   myCrs0 = myCrs0 / col_bs;
133263a3b9bcSJacob 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);
13332e68589bSMark F. Adams 
13342e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
1335849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
1336bae903cbSmarkadams4   if (size > 1) { /* get ghost null space data */
13372e68589bSMark F. Adams     PetscReal *tmp_gdata, *tmp_ldata, *tp2;
1338e0b7e82fSBarry Smith 
13399566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
13404a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
1341c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1342a2f3521dSMark F. Adams         PetscInt         ii, stride;
13438e3a54c0SPierre Jolivet         const PetscReal *tp = PetscSafePointerPlusOffset(pc_gamg->data, jj * bs * nloc + kk);
1344e0b7e82fSBarry Smith 
13452fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
13462fa5cd67SKarl Rupp 
13479566063dSJacob Faibussowitsch         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
1348a2f3521dSMark F. Adams 
1349bae903cbSmarkadams4         if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */
13509566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost));
1351a2f3521dSMark F. Adams           nbnodes = bs * stride;
13522e68589bSMark F. Adams         }
13538e3a54c0SPierre Jolivet         tp2 = PetscSafePointerPlusOffset(data_w_ghost, jj * bs * stride + kk);
13542fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
13559566063dSJacob Faibussowitsch         PetscCall(PetscFree(tmp_gdata));
13562e68589bSMark F. Adams       }
13572e68589bSMark F. Adams     }
13589566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_ldata));
1359806fa848SBarry Smith   } else {
1360c8b0795cSMark F. Adams     nbnodes      = bs * nloc;
1361835f2295SStefano Zampini     data_w_ghost = pc_gamg->data;
13622e68589bSMark F. Adams   }
13632e68589bSMark F. Adams 
1364bae903cbSmarkadams4   /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */
1365c5df96a5SBarry Smith   if (size > 1) {
13662e68589bSMark F. Adams     PetscReal *fid_glid_loc, *fiddata;
1367a2f3521dSMark F. Adams     PetscInt   stride;
13682e68589bSMark F. Adams 
13699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
13702e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk);
13719566063dSJacob Faibussowitsch     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
1372bae903cbSmarkadams4     PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */
1373a2f3521dSMark F. Adams     for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
13749566063dSJacob Faibussowitsch     PetscCall(PetscFree(fiddata));
1375a2f3521dSMark F. Adams 
137663a3b9bcSJacob Faibussowitsch     PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs);
13779566063dSJacob Faibussowitsch     PetscCall(PetscFree(fid_glid_loc));
1378806fa848SBarry Smith   } else {
13799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &flid_fgid));
13802e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk;
13812e68589bSMark F. Adams   }
1382849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
1383bae903cbSmarkadams4   /* get P0 */
1384849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
1385c8b0795cSMark F. Adams   {
13860298fd71SBarry Smith     PetscReal *data_out = NULL;
1387e0b7e82fSBarry Smith 
13889566063dSJacob Faibussowitsch     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol));
13899566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
13902fa5cd67SKarl Rupp 
1391c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
13924a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
13934a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs * col_bs * nLocalSelected;
1394c8b0795cSMark F. Adams   }
1395849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
139648a46eb9SPierre Jolivet   if (size > 1) PetscCall(PetscFree(data_w_ghost));
13979566063dSJacob Faibussowitsch   PetscCall(PetscFree(flid_fgid));
13982e68589bSMark F. Adams 
1399c8b0795cSMark F. Adams   *a_P_out = Prol; /* out */
1400e0b7e82fSBarry Smith   PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_initial_prolongation"));
1401ed4e82eaSPeter Brune 
1402849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
14033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1404c8b0795cSMark F. Adams }
1405c8b0795cSMark F. Adams 
1406c8b0795cSMark F. Adams /*
1407e0b7e82fSBarry Smith    PCGAMGOptimizeProlongator_AGG - given the initial prolongator optimizes it by smoothed aggregation pc_gamg_agg->nsmooths times
1408c8b0795cSMark F. Adams 
1409c8b0795cSMark F. Adams   Input Parameter:
1410c8b0795cSMark F. Adams    . pc - this
1411c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1412c8b0795cSMark F. Adams  In/Output Parameter:
14132a808120SBarry Smith    . a_P - prolongation operator to the next level
1414c8b0795cSMark F. Adams */
1415e0b7e82fSBarry Smith static PetscErrorCode PCGAMGOptimizeProlongator_AGG(PC pc, Mat Amat, Mat *a_P)
1416d71ae5a4SJacob Faibussowitsch {
1417c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
1418c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
1419c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
1420c8b0795cSMark F. Adams   PetscInt     jj;
1421c8b0795cSMark F. Adams   Mat          Prol = *a_P;
14223b4367a7SBarry Smith   MPI_Comm     comm;
14232a808120SBarry Smith   KSP          eksp;
14242a808120SBarry Smith   Vec          bb, xx;
14252a808120SBarry Smith   PC           epc;
14262a808120SBarry Smith   PetscReal    alpha, emax, emin;
1427c8b0795cSMark F. Adams 
1428c8b0795cSMark F. Adams   PetscFunctionBegin;
14299566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
1430849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
1431c8b0795cSMark F. Adams 
1432d5d25401SBarry Smith   /* compute maximum singular value of operator to be used in smoother */
14332a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
143418c3aa7eSMark     /* get eigen estimates */
143518c3aa7eSMark     if (pc_gamg->emax > 0) {
143618c3aa7eSMark       emin = pc_gamg->emin;
143718c3aa7eSMark       emax = pc_gamg->emax;
143818c3aa7eSMark     } else {
14390ed2132dSStefano Zampini       const char *prefix;
14400ed2132dSStefano Zampini 
14419566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &bb, NULL));
14429566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &xx, NULL));
144320654ebbSStefano Zampini       PetscCall(KSPSetNoisy_Private(Amat, bb));
1444e696ed0bSMark F. Adams 
14459566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm, &eksp));
14463821be0aSBarry Smith       PetscCall(KSPSetNestLevel(eksp, pc->kspnestlevel));
14479566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
14489566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(eksp, prefix));
14499566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_"));
145073f7197eSJed Brown       {
1451b94d7dedSBarry Smith         PetscBool isset, sflg;
1452e0b7e82fSBarry Smith 
1453b94d7dedSBarry Smith         PetscCall(MatIsSPDKnown(Amat, &isset, &sflg));
1454b94d7dedSBarry Smith         if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG));
1455d24ecf33SMark       }
14569566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure));
14579566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
1458c2436cacSMark F. Adams 
14599566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
14609566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(eksp, Amat, Amat));
14612e68589bSMark F. Adams 
14629566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(eksp, &epc));
14639566063dSJacob Faibussowitsch       PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */
1464c2436cacSMark F. Adams 
1465fb842aefSJose 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
1466074aec55SMark Adams 
14679566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(eksp));
14689566063dSJacob Faibussowitsch       PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE));
14699566063dSJacob Faibussowitsch       PetscCall(KSPSolve(eksp, bb, xx));
14709566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(eksp, pc, xx));
14712e68589bSMark F. Adams 
14729566063dSJacob Faibussowitsch       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
14739566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI));
14749566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&xx));
14759566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bb));
14769566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&eksp));
14772e68589bSMark F. Adams     }
147818c3aa7eSMark     if (pc_gamg->use_sa_esteig) {
147918c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
148018c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
148163a3b9bcSJacob 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));
148218c3aa7eSMark     } else {
148318c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
148418c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
148518c3aa7eSMark     }
148618c3aa7eSMark   } else {
148718c3aa7eSMark     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
148818c3aa7eSMark     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
148918c3aa7eSMark   }
14902e68589bSMark F. Adams 
14912a808120SBarry Smith   /* smooth P0 */
1492e0b7e82fSBarry Smith   if (pc_gamg_agg->nsmooths > 0) {
14932a808120SBarry Smith     Vec diag;
14942a808120SBarry Smith 
1495e0b7e82fSBarry Smith     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
1496e0b7e82fSBarry Smith     PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero");
14972a808120SBarry Smith 
1498e0b7e82fSBarry Smith     PetscCall(MatCreateVecs(Amat, &diag, NULL));
1499e0b7e82fSBarry Smith     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
1500e0b7e82fSBarry Smith     PetscCall(VecReciprocal(diag));
1501e0b7e82fSBarry Smith 
1502e0b7e82fSBarry Smith     for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
1503e0b7e82fSBarry Smith       Mat tMat;
1504e0b7e82fSBarry Smith 
1505e0b7e82fSBarry Smith       PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
1506e0b7e82fSBarry Smith       /*
1507e0b7e82fSBarry Smith         Smooth aggregation on the prolongator
1508e0b7e82fSBarry Smith 
1509e0b7e82fSBarry Smith         P_{i} := (I - 1.4/emax D^{-1}A) P_i\{i-1}
1510e0b7e82fSBarry Smith       */
15119566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
1512fb842aefSJose E. Roman       PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_CURRENT, &tMat));
15139566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
15149566063dSJacob Faibussowitsch       PetscCall(MatProductClear(tMat));
15159566063dSJacob Faibussowitsch       PetscCall(MatDiagonalScale(tMat, diag, NULL));
1516b4da3a1bSStefano Zampini 
1517d5d25401SBarry Smith       /* TODO: Document the 1.4 and don't hardwire it in this routine */
1518b4ec6429SMark F. Adams       alpha = -1.4 / emax;
15199566063dSJacob Faibussowitsch       PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
15209566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Prol));
15212e68589bSMark F. Adams       Prol = tMat;
1522849bee69SMark Adams       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
15232e68589bSMark F. Adams     }
1524e0b7e82fSBarry Smith     PetscCall(VecDestroy(&diag));
1525e0b7e82fSBarry Smith   }
1526849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
1527e0b7e82fSBarry Smith   PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_prolongation"));
1528c8b0795cSMark F. Adams   *a_P = Prol;
15293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15302e68589bSMark F. Adams }
15312e68589bSMark F. Adams 
1532a077d33dSBarry Smith /*MC
1533a077d33dSBarry Smith   PCGAMGAGG - Smooth aggregation, {cite}`vanek1996algebraic`, {cite}`vanek2001convergence`, variant of PETSc's algebraic multigrid (`PCGAMG`) preconditioner
15342e68589bSMark F. Adams 
1535a077d33dSBarry Smith   Options Database Keys:
1536a077d33dSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation to construct prolongation
1537a077d33dSBarry Smith . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest.
1538a077d33dSBarry Smith . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening
1539a077d33dSBarry Smith . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm
1540a077d33dSBarry Smith . -pc_gamg_pc_gamg_asm_hem_aggs <n,default=0> - Number of HEM aggregation steps for ASM smoother
1541a077d33dSBarry Smith - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive')
1542a077d33dSBarry Smith 
1543a077d33dSBarry Smith   Level: intermediate
1544a077d33dSBarry Smith 
1545a077d33dSBarry Smith   Notes:
1546a077d33dSBarry Smith   To obtain good performance for `PCGAMG` for vector valued problems you must
1547a077d33dSBarry Smith   call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point.
1548a077d33dSBarry Smith   Call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator
1549a077d33dSBarry Smith 
1550a077d33dSBarry Smith   The many options for `PCMG` and `PCGAMG` such as controlling the smoothers on each level etc. also work for `PCGAMGAGG`
1551a077d33dSBarry Smith 
1552a077d33dSBarry Smith .seealso: `PCGAMG`, [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCCreate()`, `PCSetType()`,
1553a077d33dSBarry Smith           `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`,
1554a077d33dSBarry Smith           `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`,
1555a077d33dSBarry Smith           `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`,
1556a077d33dSBarry Smith           `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()`
1557a077d33dSBarry Smith M*/
1558d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc)
1559d71ae5a4SJacob Faibussowitsch {
1560c8b0795cSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
1561c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
1562c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg;
15632e68589bSMark F. Adams 
1564c8b0795cSMark F. Adams   PetscFunctionBegin;
1565c8b0795cSMark F. Adams   /* create sub context for SA */
15664dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg_agg));
1567c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1568c8b0795cSMark F. Adams 
15691ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
15709b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1571c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1572c8b0795cSMark F. Adams 
1573c8b0795cSMark F. Adams   /* set internal function pointers */
15742d776b49SBarry Smith   pc_gamg->ops->creategraph       = PCGAMGCreateGraph_AGG;
15751ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
1576e0b7e82fSBarry Smith   pc_gamg->ops->prolongator       = PCGAMGConstructProlongator_AGG;
1577e0b7e82fSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptimizeProlongator_AGG;
15781ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
15795adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
1580c8b0795cSMark F. Adams 
1581cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths                     = 1;
1582d529f056Smarkadams4   pc_gamg_agg->aggressive_coarsening_levels = 1;
1583d4adc10fSMark Adams   pc_gamg_agg->use_aggressive_square_graph  = PETSC_TRUE;
1584d529f056Smarkadams4   pc_gamg_agg->use_minimum_degree_ordering  = PETSC_FALSE;
1585a9ccf005SMark Adams   pc_gamg_agg->use_low_mem_filter           = PETSC_FALSE;
1586d529f056Smarkadams4   pc_gamg_agg->aggressive_mis_k             = 2;
15876c34c54dSStefano Zampini   pc_gamg_agg->graph_symmetrize             = PETSC_TRUE;
1588cab9ed1eSBarry Smith 
15899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG));
1590bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG));
1591d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", PCGAMGSetAggressiveSquareGraph_AGG));
1592d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", PCGAMGMISkSetMinDegreeOrdering_AGG));
1593a9ccf005SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", PCGAMGSetLowMemoryFilter_AGG));
1594d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", PCGAMGMISkSetAggressive_AGG));
15956c34c54dSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetGraphSymmetrize_C", PCGAMGSetGraphSymmetrize_AGG));
15969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG));
15973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15982e68589bSMark F. Adams }
1599