xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision e0b7e82fd3cf27fce84cc3e37e8d70a5c36a2d4e)
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;
172d776b49SBarry Smith   MatCoarsen crs;
182e68589bSMark F. Adams } PC_GAMG_AGG;
192e68589bSMark F. Adams 
202e68589bSMark F. Adams /*@
21a077d33dSBarry Smith   PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used to construct the prolongation operator
222e68589bSMark F. Adams 
23c3339decSBarry Smith   Logically Collective
242e68589bSMark F. Adams 
252e68589bSMark F. Adams   Input Parameters:
2620f4b53cSBarry Smith + pc - the preconditioner context
2720f4b53cSBarry Smith - n  - the number of smooths
282e68589bSMark F. Adams 
292e68589bSMark F. Adams   Options Database Key:
30a077d33dSBarry Smith . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use
312e68589bSMark F. Adams 
322e68589bSMark F. Adams   Level: intermediate
332e68589bSMark F. Adams 
34a077d33dSBarry Smith   Note:
35a077d33dSBarry Smith   This is a different concept from the number smoothing steps used during the linear solution process which
36a077d33dSBarry Smith   can be set with `-mg_levels_ksp_max_it`
37a077d33dSBarry Smith 
38a077d33dSBarry Smith   Developer Note:
39a077d33dSBarry Smith   This should be named `PCGAMGAGGSetNSmooths()`.
40a077d33dSBarry Smith 
41a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCMG`, `PCGAMG`
422e68589bSMark F. Adams @*/
43d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
44d71ae5a4SJacob Faibussowitsch {
452e68589bSMark F. Adams   PetscFunctionBegin;
462e68589bSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
47d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
48cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n));
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
502e68589bSMark F. Adams }
512e68589bSMark F. Adams 
52d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n)
53d71ae5a4SJacob Faibussowitsch {
542e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
552e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
56c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
572e68589bSMark F. Adams 
582e68589bSMark F. Adams   PetscFunctionBegin;
59c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61c8b0795cSMark F. Adams }
62c8b0795cSMark F. Adams 
63c8b0795cSMark F. Adams /*@
64f1580f4eSBarry Smith   PCGAMGSetAggressiveLevels -  Use aggressive coarsening on first n levels
65ef4ad70eSMark F. Adams 
66c3339decSBarry Smith   Logically Collective
67ef4ad70eSMark F. Adams 
68ef4ad70eSMark F. Adams   Input Parameters:
69cab9ed1eSBarry Smith + pc - the preconditioner context
70d5d25401SBarry Smith - n  - 0, 1 or more
71ef4ad70eSMark F. Adams 
72ef4ad70eSMark F. Adams   Options Database Key:
73a077d33dSBarry Smith . -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels on which to square the graph on before aggregating it
74af3c827dSMark Adams 
75ef4ad70eSMark F. Adams   Level: intermediate
76ef4ad70eSMark F. Adams 
77a077d33dSBarry 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()`
78ef4ad70eSMark F. Adams @*/
79d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n)
80d71ae5a4SJacob Faibussowitsch {
81ef4ad70eSMark F. Adams   PetscFunctionBegin;
82ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
83d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
84bae903cbSmarkadams4   PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86ef4ad70eSMark F. Adams }
87ef4ad70eSMark F. Adams 
88d529f056Smarkadams4 /*@
89d529f056Smarkadams4   PCGAMGMISkSetAggressive - Number (k) distance in MIS coarsening (>2 is 'aggressive')
90d529f056Smarkadams4 
91d529f056Smarkadams4   Logically Collective
92d529f056Smarkadams4 
93d529f056Smarkadams4   Input Parameters:
94d529f056Smarkadams4 + pc - the preconditioner context
95d529f056Smarkadams4 - n  - 1 or more (default = 2)
96d529f056Smarkadams4 
97d529f056Smarkadams4   Options Database Key:
98d529f056Smarkadams4 . -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive')
99d529f056Smarkadams4 
100d529f056Smarkadams4   Level: intermediate
101d529f056Smarkadams4 
102a077d33dSBarry 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()`
103d529f056Smarkadams4 @*/
104d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetAggressive(PC pc, PetscInt n)
105d529f056Smarkadams4 {
106d529f056Smarkadams4   PetscFunctionBegin;
107d529f056Smarkadams4   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
108d529f056Smarkadams4   PetscValidLogicalCollectiveInt(pc, n, 2);
109d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGMISkSetAggressive_C", (PC, PetscInt), (pc, n));
110d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
111d529f056Smarkadams4 }
112d529f056Smarkadams4 
113d529f056Smarkadams4 /*@
114d529f056Smarkadams4   PCGAMGSetAggressiveSquareGraph - Use graph square A'A for aggressive coarsening, old method
115d529f056Smarkadams4 
116d529f056Smarkadams4   Logically Collective
117d529f056Smarkadams4 
118d529f056Smarkadams4   Input Parameters:
119d529f056Smarkadams4 + pc - the preconditioner context
120d529f056Smarkadams4 - b  - default false - MIS-k is faster
121d529f056Smarkadams4 
122d529f056Smarkadams4   Options Database Key:
123d529f056Smarkadams4 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening
124d529f056Smarkadams4 
125d529f056Smarkadams4   Level: intermediate
126d529f056Smarkadams4 
127a077d33dSBarry 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()`
128d529f056Smarkadams4 @*/
129d529f056Smarkadams4 PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC pc, PetscBool b)
130d529f056Smarkadams4 {
131d529f056Smarkadams4   PetscFunctionBegin;
132d529f056Smarkadams4   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
133d529f056Smarkadams4   PetscValidLogicalCollectiveBool(pc, b, 2);
134d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGSetAggressiveSquareGraph_C", (PC, PetscBool), (pc, b));
135d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
136d529f056Smarkadams4 }
137d529f056Smarkadams4 
138d529f056Smarkadams4 /*@
139d529f056Smarkadams4   PCGAMGMISkSetMinDegreeOrdering - Use minimum degree ordering in greedy MIS algorithm
140d529f056Smarkadams4 
141d529f056Smarkadams4   Logically Collective
142d529f056Smarkadams4 
143d529f056Smarkadams4   Input Parameters:
144d529f056Smarkadams4 + pc - the preconditioner context
145d529f056Smarkadams4 - b  - default true
146d529f056Smarkadams4 
147d529f056Smarkadams4   Options Database Key:
148d529f056Smarkadams4 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm
149d529f056Smarkadams4 
150d529f056Smarkadams4   Level: intermediate
151d529f056Smarkadams4 
152a077d33dSBarry 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()`
153d529f056Smarkadams4 @*/
154d529f056Smarkadams4 PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC pc, PetscBool b)
155d529f056Smarkadams4 {
156d529f056Smarkadams4   PetscFunctionBegin;
157d529f056Smarkadams4   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
158d529f056Smarkadams4   PetscValidLogicalCollectiveBool(pc, b, 2);
159d529f056Smarkadams4   PetscTryMethod(pc, "PCGAMGMISkSetMinDegreeOrdering_C", (PC, PetscBool), (pc, b));
160d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
161d529f056Smarkadams4 }
162d529f056Smarkadams4 
163a9ccf005SMark Adams /*@
164a9ccf005SMark Adams   PCGAMGSetLowMemoryFilter - Use low memory graph/matrix filter
165a9ccf005SMark Adams 
166a9ccf005SMark Adams   Logically Collective
167a9ccf005SMark Adams 
168a9ccf005SMark Adams   Input Parameters:
169a9ccf005SMark Adams + pc - the preconditioner context
170a9ccf005SMark Adams - b  - default false
171a9ccf005SMark Adams 
172a9ccf005SMark Adams   Options Database Key:
173a9ccf005SMark Adams . -pc_gamg_low_memory_threshold_filter <bool,default=false> - Use low memory graph/matrix filter
174a9ccf005SMark Adams 
175a9ccf005SMark Adams   Level: intermediate
176a9ccf005SMark Adams 
177a077d33dSBarry Smith .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`,
178a077d33dSBarry Smith   `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`
179a9ccf005SMark Adams @*/
180a9ccf005SMark Adams PetscErrorCode PCGAMGSetLowMemoryFilter(PC pc, PetscBool b)
181a9ccf005SMark Adams {
182a9ccf005SMark Adams   PetscFunctionBegin;
183a9ccf005SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
184a9ccf005SMark Adams   PetscValidLogicalCollectiveBool(pc, b, 2);
185a9ccf005SMark Adams   PetscTryMethod(pc, "PCGAMGSetLowMemoryFilter_C", (PC, PetscBool), (pc, b));
186a9ccf005SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
187a9ccf005SMark Adams }
188a9ccf005SMark Adams 
189d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n)
190d71ae5a4SJacob Faibussowitsch {
191ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
192ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
193ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
194ef4ad70eSMark F. Adams 
195ef4ad70eSMark F. Adams   PetscFunctionBegin;
196bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = n;
1973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
198ef4ad70eSMark F. Adams }
199ef4ad70eSMark F. Adams 
200d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetAggressive_AGG(PC pc, PetscInt n)
201d529f056Smarkadams4 {
202d529f056Smarkadams4   PC_MG       *mg          = (PC_MG *)pc->data;
203d529f056Smarkadams4   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
204d529f056Smarkadams4   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
205d529f056Smarkadams4 
206d529f056Smarkadams4   PetscFunctionBegin;
207d529f056Smarkadams4   pc_gamg_agg->aggressive_mis_k = n;
208d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
209d529f056Smarkadams4 }
210d529f056Smarkadams4 
211d529f056Smarkadams4 static PetscErrorCode PCGAMGSetAggressiveSquareGraph_AGG(PC pc, PetscBool b)
212d529f056Smarkadams4 {
213d529f056Smarkadams4   PC_MG       *mg          = (PC_MG *)pc->data;
214d529f056Smarkadams4   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
215d529f056Smarkadams4   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
216d529f056Smarkadams4 
217d529f056Smarkadams4   PetscFunctionBegin;
218d529f056Smarkadams4   pc_gamg_agg->use_aggressive_square_graph = b;
219d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
220d529f056Smarkadams4 }
221d529f056Smarkadams4 
222a9ccf005SMark Adams static PetscErrorCode PCGAMGSetLowMemoryFilter_AGG(PC pc, PetscBool b)
223a9ccf005SMark Adams {
224a9ccf005SMark Adams   PC_MG       *mg          = (PC_MG *)pc->data;
225a9ccf005SMark Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
226a9ccf005SMark Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
227a9ccf005SMark Adams 
228a9ccf005SMark Adams   PetscFunctionBegin;
229a9ccf005SMark Adams   pc_gamg_agg->use_low_mem_filter = b;
230a9ccf005SMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
231a9ccf005SMark Adams }
232a9ccf005SMark Adams 
233d529f056Smarkadams4 static PetscErrorCode PCGAMGMISkSetMinDegreeOrdering_AGG(PC pc, PetscBool b)
234d529f056Smarkadams4 {
235d529f056Smarkadams4   PC_MG       *mg          = (PC_MG *)pc->data;
236d529f056Smarkadams4   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
237d529f056Smarkadams4   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
238d529f056Smarkadams4 
239d529f056Smarkadams4   PetscFunctionBegin;
240d529f056Smarkadams4   pc_gamg_agg->use_minimum_degree_ordering = b;
241d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
242d529f056Smarkadams4 }
243d529f056Smarkadams4 
244d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject)
245d71ae5a4SJacob Faibussowitsch {
2462e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
2472e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
248c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
249d4adc10fSMark 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;
250d4adc10fSMark Adams   PetscInt     nsq_graph_old = 0;
2512e68589bSMark F. Adams 
2522e68589bSMark F. Adams   PetscFunctionBegin;
253d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options");
254a077d33dSBarry 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));
255d4adc10fSMark Adams   // aggressive coarsening logic with deprecated -pc_gamg_square_graph
256d4adc10fSMark 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));
257d4adc10fSMark Adams   if (!n_aggressive_flg)
258d4adc10fSMark 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));
259d4adc10fSMark 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));
260d4adc10fSMark Adams   if (!new_sq_provided && old_sq_provided) {
261d4adc10fSMark Adams     pc_gamg_agg->aggressive_coarsening_levels = nsq_graph_old; // could be zero
262d4adc10fSMark Adams     pc_gamg_agg->use_aggressive_square_graph  = PETSC_TRUE;
263bae903cbSmarkadams4   }
264d4adc10fSMark Adams   if (new_sq_provided && old_sq_provided)
265d4adc10fSMark Adams     PetscCall(PetscInfo(pc, "Warning: both -pc_gamg_square_graph and -pc_gamg_aggressive_coarsening are used. -pc_gamg_square_graph is deprecated, Number of aggressive levels is %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels));
266d529f056Smarkadams4   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));
267a9ccf005SMark 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));
268d529f056Smarkadams4   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));
269d0609cedSBarry Smith   PetscOptionsHeadEnd();
2703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2712e68589bSMark F. Adams }
2722e68589bSMark F. Adams 
273d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
274d71ae5a4SJacob Faibussowitsch {
2752e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
2762e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
277*e0b7e82fSBarry Smith   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
2782e68589bSMark F. Adams 
2792e68589bSMark F. Adams   PetscFunctionBegin;
280*e0b7e82fSBarry Smith   PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs));
2819566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
2822e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL));
283bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL));
284d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", NULL));
285d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", NULL));
286a9ccf005SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", NULL));
287d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", NULL));
288bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
2893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2902e68589bSMark F. Adams }
2912e68589bSMark F. Adams 
2922e68589bSMark F. Adams /*
2932e68589bSMark F. Adams    PCSetCoordinates_AGG
29420f4b53cSBarry Smith 
29520f4b53cSBarry Smith    Collective
2962e68589bSMark F. Adams 
2972e68589bSMark F. Adams    Input Parameter:
2982e68589bSMark F. Adams    . pc - the preconditioner context
299145b44c9SPierre Jolivet    . ndm - dimension of data (used for dof/vertex for Stokes)
300302f38e8SMark F. Adams    . a_nloc - number of vertices local
301302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
3022e68589bSMark F. Adams */
3031e6b0712SBarry Smith 
304d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
305d71ae5a4SJacob Faibussowitsch {
3062e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
3072e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
30869344418SMark F. Adams   PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf;
309a2f3521dSMark F. Adams   Mat      mat = pc->pmat;
3102e68589bSMark F. Adams 
3112e68589bSMark F. Adams   PetscFunctionBegin;
312a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
313a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
314302f38e8SMark F. Adams   nloc = a_nloc;
3152e68589bSMark F. Adams 
3162e68589bSMark F. Adams   /* SA: null space vectors */
3179566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &ndf));               /* this does not work for Stokes */
31869344418SMark F. Adams   if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
319a2f3521dSMark F. Adams   else if (coords) {
32063a3b9bcSJacob Faibussowitsch     PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf);
32169344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */
322ad540459SPierre 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);
323806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
32471959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
32563a3b9bcSJacob 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);
326c8b0795cSMark F. Adams   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
3272e68589bSMark F. Adams 
3287f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
3299566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
3309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
3312e68589bSMark F. Adams   }
3325e116b59SBarry Smith   /* copy data in - column-oriented */
3332e68589bSMark F. Adams   for (kk = 0; kk < nloc; kk++) {
33469344418SMark F. Adams     const PetscInt M    = nloc * pc_gamg->data_cell_rows; /* stride into data */
33569344418SMark F. Adams     PetscReal     *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */
336*e0b7e82fSBarry Smith 
337c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols == 1) *data = 1.0;
3382e68589bSMark F. Adams     else {
33969344418SMark F. Adams       /* translational modes */
3402fa5cd67SKarl Rupp       for (ii = 0; ii < ndatarows; ii++) {
3412fa5cd67SKarl Rupp         for (jj = 0; jj < ndatarows; jj++) {
34269344418SMark F. Adams           if (ii == jj) data[ii * M + jj] = 1.0;
3432e68589bSMark F. Adams           else data[ii * M + jj] = 0.0;
3442fa5cd67SKarl Rupp         }
3452fa5cd67SKarl Rupp       }
34669344418SMark F. Adams 
34769344418SMark F. Adams       /* rotational modes */
3482e68589bSMark F. Adams       if (coords) {
34969344418SMark F. Adams         if (ndm == 2) {
3502e68589bSMark F. Adams           data += 2 * M;
3512e68589bSMark F. Adams           data[0] = -coords[2 * kk + 1];
3522e68589bSMark F. Adams           data[1] = coords[2 * kk];
353806fa848SBarry Smith         } else {
3542e68589bSMark F. Adams           data += 3 * M;
3559371c9d4SSatish Balay           data[0]         = 0.0;
3569371c9d4SSatish Balay           data[M + 0]     = coords[3 * kk + 2];
3579371c9d4SSatish Balay           data[2 * M + 0] = -coords[3 * kk + 1];
3589371c9d4SSatish Balay           data[1]         = -coords[3 * kk + 2];
3599371c9d4SSatish Balay           data[M + 1]     = 0.0;
3609371c9d4SSatish Balay           data[2 * M + 1] = coords[3 * kk];
3619371c9d4SSatish Balay           data[2]         = coords[3 * kk + 1];
3629371c9d4SSatish Balay           data[M + 2]     = -coords[3 * kk];
3639371c9d4SSatish Balay           data[2 * M + 2] = 0.0;
3642e68589bSMark F. Adams         }
3652e68589bSMark F. Adams       }
3662e68589bSMark F. Adams     }
3672e68589bSMark F. Adams   }
3682e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
3693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3702e68589bSMark F. Adams }
3712e68589bSMark F. Adams 
3722e68589bSMark F. Adams /*
373a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
374a2f3521dSMark F. Adams       Looks in Mat for near null space.
375a2f3521dSMark F. Adams       Does not work for Stokes
3762e68589bSMark F. Adams 
3772e68589bSMark F. Adams   Input Parameter:
3782e68589bSMark F. Adams    . pc -
379a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
3802e68589bSMark F. Adams */
381d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
382d71ae5a4SJacob Faibussowitsch {
383b8cd405aSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
384b8cd405aSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
385b8cd405aSMark F. Adams   MatNullSpace mnull;
38666f2ef4bSMark Adams 
387b817416eSBarry Smith   PetscFunctionBegin;
3889566063dSJacob Faibussowitsch   PetscCall(MatGetNearNullSpace(a_A, &mnull));
389b8cd405aSMark F. Adams   if (!mnull) {
39066f2ef4bSMark Adams     DM dm;
391*e0b7e82fSBarry Smith 
3929566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
39348a46eb9SPierre Jolivet     if (!dm) PetscCall(MatGetDM(a_A, &dm));
39466f2ef4bSMark Adams     if (dm) {
39566f2ef4bSMark Adams       PetscObject deformation;
396b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
397b0d30dd6SMatthew G. Knepley 
3989566063dSJacob Faibussowitsch       PetscCall(DMGetNumFields(dm, &Nf));
399b0d30dd6SMatthew G. Knepley       if (Nf) {
4009566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, 0, NULL, &deformation));
4019566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull));
40248a46eb9SPierre Jolivet         if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull));
40366f2ef4bSMark Adams       }
40466f2ef4bSMark Adams     }
405b0d30dd6SMatthew G. Knepley   }
40666f2ef4bSMark Adams 
40766f2ef4bSMark Adams   if (!mnull) {
408a2f3521dSMark F. Adams     PetscInt bs, NN, MM;
409*e0b7e82fSBarry Smith 
4109566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
4119566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
41263a3b9bcSJacob Faibussowitsch     PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs);
4139566063dSJacob Faibussowitsch     PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL));
414806fa848SBarry Smith   } else {
415b8cd405aSMark F. Adams     PetscReal         *nullvec;
416b8cd405aSMark F. Adams     PetscBool          has_const;
417b8cd405aSMark F. Adams     PetscInt           i, j, mlocal, nvec, bs;
418d5d25401SBarry Smith     const Vec         *vecs;
419d5d25401SBarry Smith     const PetscScalar *v;
420b817416eSBarry Smith 
4219566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &mlocal, NULL));
4229566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
423d19c4ebbSmarkadams4     for (i = 0; i < nvec; i++) {
424d19c4ebbSmarkadams4       PetscCall(VecGetLocalSize(vecs[i], &j));
425d19c4ebbSmarkadams4       PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal);
426d19c4ebbSmarkadams4     }
427a0dea87dSMark Adams     pc_gamg->data_sz = (nvec + !!has_const) * mlocal;
4289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
4299371c9d4SSatish Balay     if (has_const)
4309371c9d4SSatish Balay       for (i = 0; i < mlocal; i++) nullvec[i] = 1.0;
431b8cd405aSMark F. Adams     for (i = 0; i < nvec; i++) {
4329566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vecs[i], &v));
433b8cd405aSMark F. Adams       for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]);
4349566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vecs[i], &v));
435b8cd405aSMark F. Adams     }
436b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
437b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec + !!has_const);
4389566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
439b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
440b8cd405aSMark F. Adams   }
4413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4422e68589bSMark F. Adams }
4432e68589bSMark F. Adams 
4442e68589bSMark F. Adams /*
445bae903cbSmarkadams4   formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0
4462e68589bSMark F. Adams 
4472e68589bSMark F. Adams   Input Parameter:
4482adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
4492adfe9d3SBarry Smith    . bs - row block size
4500cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
4510cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
4522adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
4532e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
4542e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
455bae903cbSmarkadams4 
4562e68589bSMark F. Adams   Output Parameter:
4572e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
4582e68589bSMark F. Adams    . a_Prol - prolongation operator
4592e68589bSMark F. Adams */
460d71ae5a4SJacob 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)
461d71ae5a4SJacob Faibussowitsch {
462bd026e97SJed Brown   PetscInt        Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride;
4633b4367a7SBarry Smith   MPI_Comm        comm;
4642e68589bSMark F. Adams   PetscReal      *out_data;
465539c167fSBarry Smith   PetscCDIntNd   *pos;
4661943db53SBarry Smith   PCGAMGHashTable fgid_flid;
4670cbbd2e1SMark F. Adams 
4682e68589bSMark F. Adams   PetscFunctionBegin;
4699566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
4709566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
4719371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
4729371c9d4SSatish Balay   my0  = Istart / bs;
47363a3b9bcSJacob 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);
4740cbbd2e1SMark F. Adams   Iend /= bs;
4750cbbd2e1SMark F. Adams   nghosts = data_stride / bs - nloc;
4762e68589bSMark F. Adams 
4779566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid));
47848a46eb9SPierre Jolivet   for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk));
4792e68589bSMark F. Adams 
4800cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
4810cbbd2e1SMark F. Adams   for (nSelected = mm = 0; mm < nloc; mm++) {
482e78576d6SMark F. Adams     PetscBool ise;
483*e0b7e82fSBarry Smith 
4845e62d202SMark Adams     PetscCall(PetscCDIsEmptyAt(agg_llists, mm, &ise));
485e78576d6SMark F. Adams     if (!ise) nSelected++;
4860cbbd2e1SMark F. Adams   }
4879566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj));
48863a3b9bcSJacob Faibussowitsch   PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT "  != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs);
48963a3b9bcSJacob 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);
4900cbbd2e1SMark F. Adams 
4912e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
4920cbbd2e1SMark F. Adams   out_data_stride = nSelected * nSAvec;
4932fa5cd67SKarl Rupp 
4949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data));
49533812677SSatish Balay   for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL;
4960cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
4972e68589bSMark F. Adams 
4982e68589bSMark F. Adams   /* find points and set prolongation */
499c8b0795cSMark F. Adams   minsz = 100;
5000cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
5015e62d202SMark Adams     PetscCall(PetscCDCountAt(agg_llists, mm, &jj));
502e78576d6SMark F. Adams     if (jj > 0) {
5030cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
5040cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
5050cbbd2e1SMark F. Adams       PetscBLASInt   asz = jj, M = asz * bs, N = nSAvec, INFO;
5062e68589bSMark F. Adams       PetscBLASInt   Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs;
5072e68589bSMark F. Adams       PetscScalar   *qqc, *qqr, *TAU, *WORK;
5082e68589bSMark F. Adams       PetscInt      *fids;
50965d7b583SSatish Balay       PetscReal     *data;
510b817416eSBarry Smith 
5110cbbd2e1SMark F. Adams       /* count agg */
5120cbbd2e1SMark F. Adams       if (asz < minsz) minsz = asz;
5130cbbd2e1SMark F. Adams 
5140cbbd2e1SMark F. Adams       /* get block */
5159566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids));
5162e68589bSMark F. Adams 
5172e68589bSMark F. Adams       aggID = 0;
5189566063dSJacob Faibussowitsch       PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
519e78576d6SMark F. Adams       while (pos) {
520e78576d6SMark F. Adams         PetscInt gid1;
521*e0b7e82fSBarry Smith 
5229566063dSJacob Faibussowitsch         PetscCall(PetscCDIntNdGetID(pos, &gid1));
5239566063dSJacob Faibussowitsch         PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
524e78576d6SMark F. Adams 
5250cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
5260cbbd2e1SMark F. Adams         else {
5279566063dSJacob Faibussowitsch           PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid));
52808401ef6SPierre Jolivet           PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table");
5290cbbd2e1SMark F. Adams         }
5305e116b59SBarry Smith         /* copy in B_i matrix - column-oriented */
53165d7b583SSatish Balay         data = &data_in[flid * bs];
5323d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
5332e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
5340cbbd2e1SMark F. Adams             PetscReal d = data[jj * data_stride + ii];
535*e0b7e82fSBarry Smith 
5360cbbd2e1SMark F. Adams             qqc[jj * Mdata + aggID * bs + ii] = d;
5372e68589bSMark F. Adams           }
5382e68589bSMark F. Adams         }
5392e68589bSMark F. Adams         /* set fine IDs */
5402e68589bSMark F. Adams         for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk;
5412e68589bSMark F. Adams         aggID++;
5420cbbd2e1SMark F. Adams       }
5432e68589bSMark F. Adams 
5442e68589bSMark F. Adams       /* pad with zeros */
5452e68589bSMark F. Adams       for (ii = asz * bs; ii < Mdata; ii++) {
546ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0;
5472e68589bSMark F. Adams       }
5482e68589bSMark F. Adams 
5492e68589bSMark F. Adams       /* QR */
5509566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
551792fecdfSBarry Smith       PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
5529566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
55308401ef6SPierre Jolivet       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error");
5545e116b59SBarry Smith       /* get R - column-oriented - output B_{i+1} */
5552e68589bSMark F. Adams       {
5562e68589bSMark F. Adams         PetscReal *data = &out_data[clid * nSAvec];
557*e0b7e82fSBarry Smith 
5582e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
5592e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
56008401ef6SPierre 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);
5610cbbd2e1SMark F. Adams             if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]);
5620cbbd2e1SMark F. Adams             else data[jj * out_data_stride + ii] = 0.;
5632e68589bSMark F. Adams           }
5642e68589bSMark F. Adams         }
5652e68589bSMark F. Adams       }
5662e68589bSMark F. Adams 
5675e116b59SBarry Smith       /* get Q - row-oriented */
568792fecdfSBarry Smith       PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
56963a3b9bcSJacob Faibussowitsch       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO);
5702e68589bSMark F. Adams 
5712e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
572ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii];
5732e68589bSMark F. Adams       }
5742e68589bSMark F. Adams 
5752e68589bSMark F. Adams       /* add diagonal block of P0 */
5769371c9d4SSatish Balay       for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ }
5779566063dSJacob Faibussowitsch       PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES));
5789566063dSJacob Faibussowitsch       PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids));
579b43b03e9SMark F. Adams       clid++;
5800cbbd2e1SMark F. Adams     } /* coarse agg */
5810cbbd2e1SMark F. Adams   } /* for all fine nodes */
5829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
5839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
5849566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableDestroy(&fgid_flid));
5853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5862e68589bSMark F. Adams }
5872e68589bSMark F. Adams 
588d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer)
589d71ae5a4SJacob Faibussowitsch {
5905adeb434SBarry Smith   PC_MG       *mg          = (PC_MG *)pc->data;
5915adeb434SBarry Smith   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
5925adeb434SBarry Smith   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
5935adeb434SBarry Smith 
5945adeb434SBarry Smith   PetscFunctionBegin;
5959566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      AGG specific options\n"));
596d529f056Smarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number of levels of aggressive coarsening %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels));
597d529f056Smarkadams4   if (pc_gamg_agg->aggressive_coarsening_levels > 0) {
598d529f056Smarkadams4     PetscCall(PetscViewerASCIIPrintf(viewer, "        %s aggressive coarsening\n", !pc_gamg_agg->use_aggressive_square_graph ? "MIS-k" : "Square graph"));
599d529f056Smarkadams4     if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(PetscViewerASCIIPrintf(viewer, "        MIS-%d coarsening on aggressive levels\n", (int)pc_gamg_agg->aggressive_mis_k));
600d529f056Smarkadams4   }
601*e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
602*e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
603*e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
604*e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPushTab(viewer));
605*e0b7e82fSBarry Smith   PetscCall(MatCoarsenView(pc_gamg_agg->crs, viewer));
606*e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
607*e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
608*e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
609*e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPopTab(viewer));
610*e0b7e82fSBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number smoothing steps to construct prolongation %d\n", (int)pc_gamg_agg->nsmooths));
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6125adeb434SBarry Smith }
6135adeb434SBarry Smith 
6142d776b49SBarry Smith static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat)
615d71ae5a4SJacob Faibussowitsch {
616c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG *)pc->data;
617c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG *)mg->innerctx;
618c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
6192d776b49SBarry Smith   const PetscReal vfilter     = pc_gamg->threshold[pc_gamg->current_level];
6209c15c1aeSMark Adams   PetscBool       ishem, ismis;
6212d776b49SBarry Smith   const char     *prefix;
622d529f056Smarkadams4   MatInfo         info0, info1;
623d529f056Smarkadams4   PetscInt        bs;
624c8b0795cSMark F. Adams 
625c8b0795cSMark F. Adams   PetscFunctionBegin;
626a9ccf005SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
6272d776b49SBarry 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 */
6282d776b49SBarry Smith   /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */
629*e0b7e82fSBarry Smith   PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs));
6302d776b49SBarry Smith   PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs));
6312d776b49SBarry Smith   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix));
6322d776b49SBarry Smith   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix));
633*e0b7e82fSBarry Smith   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc_gamg_agg->crs, "pc_gamg_"));
6342d776b49SBarry Smith   PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs));
635e02fb3cdSMark Adams   PetscCall(MatGetBlockSize(Amat, &bs));
636e02fb3cdSMark Adams   // check for valid indices wrt bs
637e02fb3cdSMark Adams   for (int ii = 0; ii < pc_gamg_agg->crs->strength_index_size; ii++) {
638b65aec2dSMark Adams     PetscCheck(pc_gamg_agg->crs->strength_index[ii] >= 0 && pc_gamg_agg->crs->strength_index[ii] < bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Indices (%d) must be non-negative and < block size (%d), NB, can not use -mat_coarsen_strength_index with -mat_coarsen_strength_index",
639b65aec2dSMark Adams                (int)pc_gamg_agg->crs->strength_index[ii], (int)bs);
640e02fb3cdSMark Adams   }
6412d776b49SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem));
6425e62d202SMark Adams   if (ishem) {
6439c15c1aeSMark Adams     if (pc_gamg_agg->aggressive_coarsening_levels) PetscCall(PetscInfo(pc, "HEM and aggressive coarsening ignored: HEM using %d iterations\n", (int)pc_gamg_agg->crs->max_it));
6445e62d202SMark Adams     pc_gamg_agg->aggressive_coarsening_levels = 0;                                         // aggressive and HEM does not make sense
6455e62d202SMark Adams     PetscCall(MatCoarsenSetMaximumIterations(pc_gamg_agg->crs, pc_gamg_agg->crs->max_it)); // for code coverage
6465e62d202SMark Adams     PetscCall(MatCoarsenSetThreshold(pc_gamg_agg->crs, vfilter));                          // for code coverage
6479c15c1aeSMark Adams   } else {
6489c15c1aeSMark Adams     PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENMIS, &ismis));
6499c15c1aeSMark Adams     if (ismis && pc_gamg_agg->aggressive_coarsening_levels && !pc_gamg_agg->use_aggressive_square_graph) {
6509c15c1aeSMark Adams       PetscCall(PetscInfo(pc, "MIS and aggressive coarsening and no square graph: force square graph\n"));
6519c15c1aeSMark Adams       pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE;
6529c15c1aeSMark Adams     }
6535e62d202SMark Adams   }
654a9ccf005SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
655a9ccf005SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
656d529f056Smarkadams4   PetscCall(MatGetInfo(Amat, MAT_LOCAL, &info0)); /* global reduction */
657a9ccf005SMark Adams 
658a9ccf005SMark Adams   if (ishem || pc_gamg_agg->use_low_mem_filter) {
659e02fb3cdSMark Adams     PetscCall(MatCreateGraph(Amat, PETSC_TRUE, (vfilter >= 0 || ishem) ? PETSC_TRUE : PETSC_FALSE, vfilter, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat));
660a9ccf005SMark Adams   } else {
661d8b4a066SPierre Jolivet     // make scalar graph, symetrize if not know to be symmetric, scale, but do not filter (expensive)
662e02fb3cdSMark Adams     PetscCall(MatCreateGraph(Amat, PETSC_TRUE, PETSC_TRUE, -1, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat));
663a9ccf005SMark Adams     if (vfilter >= 0) {
664a9ccf005SMark Adams       PetscInt           Istart, Iend, ncols, nnz0, nnz1, NN, MM, nloc;
665a9ccf005SMark Adams       Mat                tGmat, Gmat = *a_Gmat;
666a9ccf005SMark Adams       MPI_Comm           comm;
667a9ccf005SMark Adams       const PetscScalar *vals;
668a9ccf005SMark Adams       const PetscInt    *idx;
669a9ccf005SMark Adams       PetscInt          *d_nnz, *o_nnz, kk, *garray = NULL, *AJ, maxcols = 0;
670a9ccf005SMark Adams       MatScalar         *AA; // this is checked in graph
671a9ccf005SMark Adams       PetscBool          isseqaij;
672a9ccf005SMark Adams       Mat                a, b, c;
673a9ccf005SMark Adams       MatType            jtype;
674a9ccf005SMark Adams 
675a9ccf005SMark Adams       PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm));
676a9ccf005SMark Adams       PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isseqaij));
677a9ccf005SMark Adams       PetscCall(MatGetType(Gmat, &jtype));
678a9ccf005SMark Adams       PetscCall(MatCreate(comm, &tGmat));
679a9ccf005SMark Adams       PetscCall(MatSetType(tGmat, jtype));
680a9ccf005SMark Adams 
681a9ccf005SMark Adams       /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold?
682a9ccf005SMark Adams         Also, if the matrix is symmetric, can we skip this
683a9ccf005SMark Adams         operation? It can be very expensive on large matrices. */
684a9ccf005SMark Adams 
685a9ccf005SMark Adams       // global sizes
686a9ccf005SMark Adams       PetscCall(MatGetSize(Gmat, &MM, &NN));
687a9ccf005SMark Adams       PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend));
688a9ccf005SMark Adams       nloc = Iend - Istart;
689a9ccf005SMark Adams       PetscCall(PetscMalloc2(nloc, &d_nnz, nloc, &o_nnz));
690a9ccf005SMark Adams       if (isseqaij) {
691a9ccf005SMark Adams         a = Gmat;
692a9ccf005SMark Adams         b = NULL;
693a9ccf005SMark Adams       } else {
694a9ccf005SMark Adams         Mat_MPIAIJ *d = (Mat_MPIAIJ *)Gmat->data;
695*e0b7e82fSBarry Smith 
696a9ccf005SMark Adams         a      = d->A;
697a9ccf005SMark Adams         b      = d->B;
698a9ccf005SMark Adams         garray = d->garray;
699a9ccf005SMark Adams       }
700a9ccf005SMark Adams       /* Determine upper bound on non-zeros needed in new filtered matrix */
701a9ccf005SMark Adams       for (PetscInt row = 0; row < nloc; row++) {
702a9ccf005SMark Adams         PetscCall(MatGetRow(a, row, &ncols, NULL, NULL));
703a9ccf005SMark Adams         d_nnz[row] = ncols;
704a9ccf005SMark Adams         if (ncols > maxcols) maxcols = ncols;
705a9ccf005SMark Adams         PetscCall(MatRestoreRow(a, row, &ncols, NULL, NULL));
706a9ccf005SMark Adams       }
707a9ccf005SMark Adams       if (b) {
708a9ccf005SMark Adams         for (PetscInt row = 0; row < nloc; row++) {
709a9ccf005SMark Adams           PetscCall(MatGetRow(b, row, &ncols, NULL, NULL));
710a9ccf005SMark Adams           o_nnz[row] = ncols;
711a9ccf005SMark Adams           if (ncols > maxcols) maxcols = ncols;
712a9ccf005SMark Adams           PetscCall(MatRestoreRow(b, row, &ncols, NULL, NULL));
713a9ccf005SMark Adams         }
714a9ccf005SMark Adams       }
715a9ccf005SMark Adams       PetscCall(MatSetSizes(tGmat, nloc, nloc, MM, MM));
716a9ccf005SMark Adams       PetscCall(MatSetBlockSizes(tGmat, 1, 1));
717a9ccf005SMark Adams       PetscCall(MatSeqAIJSetPreallocation(tGmat, 0, d_nnz));
718a9ccf005SMark Adams       PetscCall(MatMPIAIJSetPreallocation(tGmat, 0, d_nnz, 0, o_nnz));
719a9ccf005SMark Adams       PetscCall(MatSetOption(tGmat, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
720a9ccf005SMark Adams       PetscCall(PetscFree2(d_nnz, o_nnz));
721a9ccf005SMark Adams       PetscCall(PetscMalloc2(maxcols, &AA, maxcols, &AJ));
722a9ccf005SMark Adams       nnz0 = nnz1 = 0;
723a9ccf005SMark Adams       for (c = a, kk = 0; c && kk < 2; c = b, kk++) {
724a9ccf005SMark Adams         for (PetscInt row = 0, grow = Istart, ncol_row, jj; row < nloc; row++, grow++) {
725a9ccf005SMark Adams           PetscCall(MatGetRow(c, row, &ncols, &idx, &vals));
726a9ccf005SMark Adams           for (ncol_row = jj = 0; jj < ncols; jj++, nnz0++) {
727a9ccf005SMark Adams             PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
728a9ccf005SMark Adams             if (PetscRealPart(sv) > vfilter) {
729a9ccf005SMark Adams               PetscInt cid = idx[jj] + Istart; //diag
730*e0b7e82fSBarry Smith 
731a9ccf005SMark Adams               nnz1++;
732a9ccf005SMark Adams               if (c != a) cid = garray[idx[jj]];
733a9ccf005SMark Adams               AA[ncol_row] = vals[jj];
734a9ccf005SMark Adams               AJ[ncol_row] = cid;
735a9ccf005SMark Adams               ncol_row++;
736a9ccf005SMark Adams             }
737a9ccf005SMark Adams           }
738a9ccf005SMark Adams           PetscCall(MatRestoreRow(c, row, &ncols, &idx, &vals));
739a9ccf005SMark Adams           PetscCall(MatSetValues(tGmat, 1, &grow, ncol_row, AJ, AA, INSERT_VALUES));
740a9ccf005SMark Adams         }
741a9ccf005SMark Adams       }
742a9ccf005SMark Adams       PetscCall(PetscFree2(AA, AJ));
743a9ccf005SMark Adams       PetscCall(MatAssemblyBegin(tGmat, MAT_FINAL_ASSEMBLY));
744a9ccf005SMark Adams       PetscCall(MatAssemblyEnd(tGmat, MAT_FINAL_ASSEMBLY));
745a9ccf005SMark Adams       PetscCall(MatPropagateSymmetryOptions(Gmat, tGmat)); /* Normal Mat options are not relevant ? */
746a9ccf005SMark 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));
747a9ccf005SMark Adams       PetscCall(MatViewFromOptions(tGmat, NULL, "-mat_filter_graph_view"));
748a9ccf005SMark Adams       PetscCall(MatDestroy(&Gmat));
749a9ccf005SMark Adams       *a_Gmat = tGmat;
750a9ccf005SMark Adams     }
751a9ccf005SMark Adams   }
752a9ccf005SMark Adams 
753d529f056Smarkadams4   PetscCall(MatGetInfo(*a_Gmat, MAT_LOCAL, &info1)); /* global reduction */
754d529f056Smarkadams4   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));
755849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
7563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
757c8b0795cSMark F. Adams }
758c8b0795cSMark F. Adams 
759d529f056Smarkadams4 typedef PetscInt    NState;
760d529f056Smarkadams4 static const NState NOT_DONE = -2;
761d529f056Smarkadams4 static const NState DELETED  = -1;
762d529f056Smarkadams4 static const NState REMOVED  = -3;
763d529f056Smarkadams4 #define IS_SELECTED(s) (s != DELETED && s != NOT_DONE && s != REMOVED)
764d529f056Smarkadams4 
765d529f056Smarkadams4 /*
766d529f056Smarkadams4    fixAggregatesWithSquare - greedy grab of with G1 (unsquared graph) -- AIJ specific -- change to fixAggregatesWithSquare -- TODD
767d529f056Smarkadams4      - AGG-MG specific: clears singletons out of 'selected_2'
768d529f056Smarkadams4 
769d529f056Smarkadams4    Input Parameter:
770d529f056Smarkadams4    . Gmat_2 - global matrix of squared graph (data not defined)
771d529f056Smarkadams4    . Gmat_1 - base graph to grab with base graph
772d529f056Smarkadams4    Input/Output Parameter:
773d529f056Smarkadams4    . aggs_2 - linked list of aggs with gids)
774d529f056Smarkadams4 */
775d529f056Smarkadams4 static PetscErrorCode fixAggregatesWithSquare(PC pc, Mat Gmat_2, Mat Gmat_1, PetscCoarsenData *aggs_2)
776d529f056Smarkadams4 {
777d529f056Smarkadams4   PetscBool      isMPI;
778d529f056Smarkadams4   Mat_SeqAIJ    *matA_1, *matB_1 = NULL;
779d529f056Smarkadams4   MPI_Comm       comm;
780d529f056Smarkadams4   PetscInt       lid, *ii, *idx, ix, Iend, my0, kk, n, j;
781d529f056Smarkadams4   Mat_MPIAIJ    *mpimat_2 = NULL, *mpimat_1 = NULL;
782d529f056Smarkadams4   const PetscInt nloc = Gmat_2->rmap->n;
783d529f056Smarkadams4   PetscScalar   *cpcol_1_state, *cpcol_2_state, *cpcol_2_par_orig, *lid_parent_gid;
784d529f056Smarkadams4   PetscInt      *lid_cprowID_1 = NULL;
785d529f056Smarkadams4   NState        *lid_state;
786d529f056Smarkadams4   Vec            ghost_par_orig2;
787d529f056Smarkadams4   PetscMPIInt    rank;
788d529f056Smarkadams4 
789d529f056Smarkadams4   PetscFunctionBegin;
790d529f056Smarkadams4   PetscCall(PetscObjectGetComm((PetscObject)Gmat_2, &comm));
791d529f056Smarkadams4   PetscCallMPI(MPI_Comm_rank(comm, &rank));
792d529f056Smarkadams4   PetscCall(MatGetOwnershipRange(Gmat_1, &my0, &Iend));
793d529f056Smarkadams4 
794d529f056Smarkadams4   /* get submatrices */
795d529f056Smarkadams4   PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATMPIAIJ, &isMPI));
796d529f056Smarkadams4   PetscCall(PetscInfo(pc, "isMPI = %s\n", isMPI ? "yes" : "no"));
797d529f056Smarkadams4   PetscCall(PetscMalloc3(nloc, &lid_state, nloc, &lid_parent_gid, nloc, &lid_cprowID_1));
798d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
799d529f056Smarkadams4   if (isMPI) {
800d529f056Smarkadams4     /* grab matrix objects */
801d529f056Smarkadams4     mpimat_2 = (Mat_MPIAIJ *)Gmat_2->data;
802d529f056Smarkadams4     mpimat_1 = (Mat_MPIAIJ *)Gmat_1->data;
803d529f056Smarkadams4     matA_1   = (Mat_SeqAIJ *)mpimat_1->A->data;
804d529f056Smarkadams4     matB_1   = (Mat_SeqAIJ *)mpimat_1->B->data;
805d529f056Smarkadams4 
806d529f056Smarkadams4     /* force compressed row storage for B matrix in AuxMat */
807d529f056Smarkadams4     PetscCall(MatCheckCompressedRow(mpimat_1->B, matB_1->nonzerorowcnt, &matB_1->compressedrow, matB_1->i, Gmat_1->rmap->n, -1.0));
808d529f056Smarkadams4     for (ix = 0; ix < matB_1->compressedrow.nrows; ix++) {
809d529f056Smarkadams4       PetscInt lid = matB_1->compressedrow.rindex[ix];
810*e0b7e82fSBarry Smith 
811d529f056Smarkadams4       PetscCheck(lid <= nloc && lid >= -1, PETSC_COMM_SELF, PETSC_ERR_USER, "lid %d out of range. nloc = %d", (int)lid, (int)nloc);
812d529f056Smarkadams4       if (lid != -1) lid_cprowID_1[lid] = ix;
813d529f056Smarkadams4     }
814d529f056Smarkadams4   } else {
815d529f056Smarkadams4     PetscBool isAIJ;
816*e0b7e82fSBarry Smith 
817d529f056Smarkadams4     PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATSEQAIJ, &isAIJ));
818d529f056Smarkadams4     PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix.");
819d529f056Smarkadams4     matA_1 = (Mat_SeqAIJ *)Gmat_1->data;
820d529f056Smarkadams4   }
821d529f056Smarkadams4   if (nloc > 0) { PetscCheck(!matB_1 || matB_1->compressedrow.use, PETSC_COMM_SELF, PETSC_ERR_PLIB, "matB_1 && !matB_1->compressedrow.use: PETSc bug???"); }
822d529f056Smarkadams4   /* get state of locals and selected gid for deleted */
823d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) {
824d529f056Smarkadams4     lid_parent_gid[lid] = -1.0;
825d529f056Smarkadams4     lid_state[lid]      = DELETED;
826d529f056Smarkadams4   }
827d529f056Smarkadams4 
828d529f056Smarkadams4   /* set lid_state */
829d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) {
830d529f056Smarkadams4     PetscCDIntNd *pos;
831*e0b7e82fSBarry Smith 
832d529f056Smarkadams4     PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
833d529f056Smarkadams4     if (pos) {
834d529f056Smarkadams4       PetscInt gid1;
835d529f056Smarkadams4 
836d529f056Smarkadams4       PetscCall(PetscCDIntNdGetID(pos, &gid1));
837d529f056Smarkadams4       PetscCheck(gid1 == lid + my0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gid1 %d != lid %d + my0 %d", (int)gid1, (int)lid, (int)my0);
838d529f056Smarkadams4       lid_state[lid] = gid1;
839d529f056Smarkadams4     }
840d529f056Smarkadams4   }
841d529f056Smarkadams4 
842d529f056Smarkadams4   /* map local to selected local, DELETED means a ghost owns it */
843d529f056Smarkadams4   for (lid = kk = 0; lid < nloc; lid++) {
844d529f056Smarkadams4     NState state = lid_state[lid];
845d529f056Smarkadams4     if (IS_SELECTED(state)) {
846d529f056Smarkadams4       PetscCDIntNd *pos;
847*e0b7e82fSBarry Smith 
848d529f056Smarkadams4       PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
849d529f056Smarkadams4       while (pos) {
850d529f056Smarkadams4         PetscInt gid1;
851*e0b7e82fSBarry Smith 
852d529f056Smarkadams4         PetscCall(PetscCDIntNdGetID(pos, &gid1));
853d529f056Smarkadams4         PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos));
854d529f056Smarkadams4         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1 - my0] = (PetscScalar)(lid + my0);
855d529f056Smarkadams4       }
856d529f056Smarkadams4     }
857d529f056Smarkadams4   }
858d529f056Smarkadams4   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
859d529f056Smarkadams4   if (isMPI) {
860d529f056Smarkadams4     Vec tempVec;
861d529f056Smarkadams4     /* get 'cpcol_1_state' */
862d529f056Smarkadams4     PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL));
863d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
864d529f056Smarkadams4       PetscScalar v = (PetscScalar)lid_state[kk];
865*e0b7e82fSBarry Smith 
866d529f056Smarkadams4       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
867d529f056Smarkadams4     }
868d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
869d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
870d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD));
871d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD));
872d529f056Smarkadams4     PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state));
873d529f056Smarkadams4     /* get 'cpcol_2_state' */
874d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD));
875d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD));
876d529f056Smarkadams4     PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state));
877d529f056Smarkadams4     /* get 'cpcol_2_par_orig' */
878d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
879d529f056Smarkadams4       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
880*e0b7e82fSBarry Smith 
881d529f056Smarkadams4       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
882d529f056Smarkadams4     }
883d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
884d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
885d529f056Smarkadams4     PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2));
886d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD));
887d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD));
888d529f056Smarkadams4     PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig));
889d529f056Smarkadams4 
890d529f056Smarkadams4     PetscCall(VecDestroy(&tempVec));
891d529f056Smarkadams4   } /* ismpi */
892d529f056Smarkadams4   for (lid = 0; lid < nloc; lid++) {
893d529f056Smarkadams4     NState state = lid_state[lid];
894*e0b7e82fSBarry Smith 
895d529f056Smarkadams4     if (IS_SELECTED(state)) {
896d529f056Smarkadams4       /* steal locals */
897d529f056Smarkadams4       ii  = matA_1->i;
898d529f056Smarkadams4       n   = ii[lid + 1] - ii[lid];
899d529f056Smarkadams4       idx = matA_1->j + ii[lid];
900d529f056Smarkadams4       for (j = 0; j < n; j++) {
901d529f056Smarkadams4         PetscInt lidj   = idx[j], sgid;
902d529f056Smarkadams4         NState   statej = lid_state[lidj];
903*e0b7e82fSBarry Smith 
904d529f056Smarkadams4         if (statej == DELETED && (sgid = (PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid + my0) { /* steal local */
905d529f056Smarkadams4           lid_parent_gid[lidj] = (PetscScalar)(lid + my0);                                              /* send this if sgid is not local */
906d529f056Smarkadams4           if (sgid >= my0 && sgid < Iend) {                                                             /* I'm stealing this local from a local sgid */
907d529f056Smarkadams4             PetscInt      hav = 0, slid = sgid - my0, gidj = lidj + my0;
908d529f056Smarkadams4             PetscCDIntNd *pos, *last = NULL;
909*e0b7e82fSBarry Smith 
910d529f056Smarkadams4             /* looking for local from local so id_llist_2 works */
911d529f056Smarkadams4             PetscCall(PetscCDGetHeadPos(aggs_2, slid, &pos));
912d529f056Smarkadams4             while (pos) {
913d529f056Smarkadams4               PetscInt gid;
914d529f056Smarkadams4               PetscCall(PetscCDIntNdGetID(pos, &gid));
915d529f056Smarkadams4               if (gid == gidj) {
916d529f056Smarkadams4                 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null");
917d529f056Smarkadams4                 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last));
918d529f056Smarkadams4                 PetscCall(PetscCDAppendNode(aggs_2, lid, pos));
919d529f056Smarkadams4                 hav = 1;
920d529f056Smarkadams4                 break;
921d529f056Smarkadams4               } else last = pos;
922d529f056Smarkadams4               PetscCall(PetscCDGetNextPos(aggs_2, slid, &pos));
923d529f056Smarkadams4             }
924d529f056Smarkadams4             if (hav != 1) {
925d529f056Smarkadams4               PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find adj in 'selected' lists - structurally unsymmetric matrix");
926d529f056Smarkadams4               SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav);
927d529f056Smarkadams4             }
928d529f056Smarkadams4           } else { /* I'm stealing this local, owned by a ghost */
929d529f056Smarkadams4             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.",
930d529f056Smarkadams4                        ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "");
931d529f056Smarkadams4             PetscCall(PetscCDAppendID(aggs_2, lid, lidj + my0));
932d529f056Smarkadams4           }
933d529f056Smarkadams4         }
934d529f056Smarkadams4       } /* local neighbors */
935d529f056Smarkadams4     } else if (state == DELETED /* && lid_cprowID_1 */) {
936d529f056Smarkadams4       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
937*e0b7e82fSBarry Smith 
938d529f056Smarkadams4       /* see if I have a selected ghost neighbor that will steal me */
939d529f056Smarkadams4       if ((ix = lid_cprowID_1[lid]) != -1) {
940d529f056Smarkadams4         ii  = matB_1->compressedrow.i;
941d529f056Smarkadams4         n   = ii[ix + 1] - ii[ix];
942d529f056Smarkadams4         idx = matB_1->j + ii[ix];
943d529f056Smarkadams4         for (j = 0; j < n; j++) {
944d529f056Smarkadams4           PetscInt cpid   = idx[j];
945d529f056Smarkadams4           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
946*e0b7e82fSBarry Smith 
947d529f056Smarkadams4           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
948d529f056Smarkadams4             lid_parent_gid[lid] = (PetscScalar)statej;              /* send who selected */
949d529f056Smarkadams4             if (sgidold >= my0 && sgidold < Iend) {                 /* this was mine */
950d529f056Smarkadams4               PetscInt      hav = 0, oldslidj = sgidold - my0;
951d529f056Smarkadams4               PetscCDIntNd *pos, *last        = NULL;
952*e0b7e82fSBarry Smith 
953d529f056Smarkadams4               /* remove from 'oldslidj' list */
954d529f056Smarkadams4               PetscCall(PetscCDGetHeadPos(aggs_2, oldslidj, &pos));
955d529f056Smarkadams4               while (pos) {
956d529f056Smarkadams4                 PetscInt gid;
957*e0b7e82fSBarry Smith 
958d529f056Smarkadams4                 PetscCall(PetscCDIntNdGetID(pos, &gid));
959d529f056Smarkadams4                 if (lid + my0 == gid) {
960d529f056Smarkadams4                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
961d529f056Smarkadams4                   PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null");
962d529f056Smarkadams4                   PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last));
963d529f056Smarkadams4                   /* ghost (PetscScalar)statej will add this later */
964d529f056Smarkadams4                   hav = 1;
965d529f056Smarkadams4                   break;
966d529f056Smarkadams4                 } else last = pos;
967d529f056Smarkadams4                 PetscCall(PetscCDGetNextPos(aggs_2, oldslidj, &pos));
968d529f056Smarkadams4               }
969d529f056Smarkadams4               if (hav != 1) {
970d529f056Smarkadams4                 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find (hav=%d) adj in 'selected' lists - structurally unsymmetric matrix", (int)hav);
971d529f056Smarkadams4                 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav);
972d529f056Smarkadams4               }
973d529f056Smarkadams4             } else {
974d529f056Smarkadams4               /* TODO: ghosts remove this later */
975d529f056Smarkadams4             }
976d529f056Smarkadams4           }
977d529f056Smarkadams4         }
978d529f056Smarkadams4       }
979d529f056Smarkadams4     } /* selected/deleted */
980d529f056Smarkadams4   } /* node loop */
981d529f056Smarkadams4 
982d529f056Smarkadams4   if (isMPI) {
983d529f056Smarkadams4     PetscScalar    *cpcol_2_parent, *cpcol_2_gid;
984d529f056Smarkadams4     Vec             tempVec, ghostgids2, ghostparents2;
985d529f056Smarkadams4     PetscInt        cpid, nghost_2;
986d529f056Smarkadams4     PCGAMGHashTable gid_cpid;
987d529f056Smarkadams4 
988d529f056Smarkadams4     PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2));
989d529f056Smarkadams4     PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL));
990d529f056Smarkadams4 
991d529f056Smarkadams4     /* get 'cpcol_2_parent' */
992d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) { PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES)); }
993d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
994d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
995d529f056Smarkadams4     PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2));
996d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD));
997d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD));
998d529f056Smarkadams4     PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent));
999d529f056Smarkadams4 
1000d529f056Smarkadams4     /* get 'cpcol_2_gid' */
1001d529f056Smarkadams4     for (kk = 0, j = my0; kk < nloc; kk++, j++) {
1002d529f056Smarkadams4       PetscScalar v = (PetscScalar)j;
1003*e0b7e82fSBarry Smith 
1004d529f056Smarkadams4       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
1005d529f056Smarkadams4     }
1006d529f056Smarkadams4     PetscCall(VecAssemblyBegin(tempVec));
1007d529f056Smarkadams4     PetscCall(VecAssemblyEnd(tempVec));
1008d529f056Smarkadams4     PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2));
1009d529f056Smarkadams4     PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD));
1010d529f056Smarkadams4     PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD));
1011d529f056Smarkadams4     PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid));
1012d529f056Smarkadams4     PetscCall(VecDestroy(&tempVec));
1013d529f056Smarkadams4 
1014d529f056Smarkadams4     /* look for deleted ghosts and add to table */
1015d529f056Smarkadams4     PetscCall(PCGAMGHashTableCreate(2 * nghost_2 + 1, &gid_cpid));
1016d529f056Smarkadams4     for (cpid = 0; cpid < nghost_2; cpid++) {
1017d529f056Smarkadams4       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
1018*e0b7e82fSBarry Smith 
1019d529f056Smarkadams4       if (state == DELETED) {
1020d529f056Smarkadams4         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
1021d529f056Smarkadams4         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
1022*e0b7e82fSBarry Smith 
1023d529f056Smarkadams4         if (sgid_old == -1 && sgid_new != -1) {
1024d529f056Smarkadams4           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
1025*e0b7e82fSBarry Smith 
1026d529f056Smarkadams4           PetscCall(PCGAMGHashTableAdd(&gid_cpid, gid, cpid));
1027d529f056Smarkadams4         }
1028d529f056Smarkadams4       }
1029d529f056Smarkadams4     }
1030d529f056Smarkadams4 
1031d529f056Smarkadams4     /* look for deleted ghosts and see if they moved - remove it */
1032d529f056Smarkadams4     for (lid = 0; lid < nloc; lid++) {
1033d529f056Smarkadams4       NState state = lid_state[lid];
1034d529f056Smarkadams4       if (IS_SELECTED(state)) {
1035d529f056Smarkadams4         PetscCDIntNd *pos, *last = NULL;
1036d529f056Smarkadams4         /* look for deleted ghosts and see if they moved */
1037d529f056Smarkadams4         PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos));
1038d529f056Smarkadams4         while (pos) {
1039d529f056Smarkadams4           PetscInt gid;
1040d529f056Smarkadams4           PetscCall(PetscCDIntNdGetID(pos, &gid));
1041d529f056Smarkadams4 
1042d529f056Smarkadams4           if (gid < my0 || gid >= Iend) {
1043d529f056Smarkadams4             PetscCall(PCGAMGHashTableFind(&gid_cpid, gid, &cpid));
1044d529f056Smarkadams4             if (cpid != -1) {
1045d529f056Smarkadams4               /* a moved ghost - */
1046d529f056Smarkadams4               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
1047d529f056Smarkadams4               PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last));
1048d529f056Smarkadams4             } else last = pos;
1049d529f056Smarkadams4           } else last = pos;
1050d529f056Smarkadams4 
1051d529f056Smarkadams4           PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos));
1052d529f056Smarkadams4         } /* loop over list of deleted */
1053d529f056Smarkadams4       } /* selected */
1054d529f056Smarkadams4     }
1055d529f056Smarkadams4     PetscCall(PCGAMGHashTableDestroy(&gid_cpid));
1056d529f056Smarkadams4 
1057d529f056Smarkadams4     /* look at ghosts, see if they changed - and it */
1058d529f056Smarkadams4     for (cpid = 0; cpid < nghost_2; cpid++) {
1059d529f056Smarkadams4       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
1060*e0b7e82fSBarry Smith 
1061d529f056Smarkadams4       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
1062d529f056Smarkadams4         PetscInt      gid      = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
1063d529f056Smarkadams4         PetscInt      slid_new = sgid_new - my0, hav = 0;
1064d529f056Smarkadams4         PetscCDIntNd *pos;
1065d529f056Smarkadams4 
1066d529f056Smarkadams4         /* search for this gid to see if I have it */
1067d529f056Smarkadams4         PetscCall(PetscCDGetHeadPos(aggs_2, slid_new, &pos));
1068d529f056Smarkadams4         while (pos) {
1069d529f056Smarkadams4           PetscInt gidj;
1070*e0b7e82fSBarry Smith 
1071d529f056Smarkadams4           PetscCall(PetscCDIntNdGetID(pos, &gidj));
1072d529f056Smarkadams4           PetscCall(PetscCDGetNextPos(aggs_2, slid_new, &pos));
1073d529f056Smarkadams4 
1074d529f056Smarkadams4           if (gidj == gid) {
1075d529f056Smarkadams4             hav = 1;
1076d529f056Smarkadams4             break;
1077d529f056Smarkadams4           }
1078d529f056Smarkadams4         }
1079d529f056Smarkadams4         if (hav != 1) {
1080d529f056Smarkadams4           /* insert 'flidj' into head of llist */
1081d529f056Smarkadams4           PetscCall(PetscCDAppendID(aggs_2, slid_new, gid));
1082d529f056Smarkadams4         }
1083d529f056Smarkadams4       }
1084d529f056Smarkadams4     }
1085d529f056Smarkadams4     PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state));
1086d529f056Smarkadams4     PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state));
1087d529f056Smarkadams4     PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent));
1088d529f056Smarkadams4     PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid));
1089d529f056Smarkadams4     PetscCall(VecDestroy(&ghostgids2));
1090d529f056Smarkadams4     PetscCall(VecDestroy(&ghostparents2));
1091d529f056Smarkadams4     PetscCall(VecDestroy(&ghost_par_orig2));
1092d529f056Smarkadams4   }
1093d529f056Smarkadams4   PetscCall(PetscFree3(lid_state, lid_parent_gid, lid_cprowID_1));
1094d529f056Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
1095d529f056Smarkadams4 }
1096d529f056Smarkadams4 
1097c8b0795cSMark F. Adams /*
1098bae903cbSmarkadams4    PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for
1099bae903cbSmarkadams4      communication of QR data used with HEM and MISk coarsening
1100c8b0795cSMark F. Adams 
1101c8b0795cSMark F. Adams   Input Parameter:
1102e0940f08SMark F. Adams    . a_pc - this
1103bae903cbSmarkadams4 
1104e0940f08SMark F. Adams   Input/Output Parameter:
1105bae903cbSmarkadams4    . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out)
1106bae903cbSmarkadams4 
1107c8b0795cSMark F. Adams   Output Parameter:
11080cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
1109bae903cbSmarkadams4 
1110c8b0795cSMark F. Adams */
1111d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists)
1112d71ae5a4SJacob Faibussowitsch {
1113e0940f08SMark F. Adams   PC_MG       *mg          = (PC_MG *)a_pc->data;
1114c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
1115c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
11168926f930SMark Adams   Mat          Gmat2, Gmat1 = *a_Gmat1; /* aggressive graph */
11170cbbd2e1SMark F. Adams   IS           perm;
1118bae903cbSmarkadams4   PetscInt     Istart, Iend, Ii, nloc, bs, nn;
1119bae903cbSmarkadams4   PetscInt    *permute, *degree;
1120c8b0795cSMark F. Adams   PetscBool   *bIndexSet;
1121aea53286SMark Adams   PetscReal    hashfact;
1122e2c00dcbSBarry Smith   PetscInt     iSwapIndex;
11233b50add6SMark Adams   PetscRandom  random;
1124d529f056Smarkadams4   MPI_Comm     comm;
1125c8b0795cSMark F. Adams 
1126c8b0795cSMark F. Adams   PetscFunctionBegin;
1127d529f056Smarkadams4   PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm));
1128849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
1129bae903cbSmarkadams4   PetscCall(MatGetLocalSize(Gmat1, &nn, NULL));
11309566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Gmat1, &bs));
113163a3b9bcSJacob Faibussowitsch   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs);
1132bae903cbSmarkadams4   nloc = nn / bs;
11335cfd4bd9SMark Adams   /* get MIS aggs - randomize */
1134bae903cbSmarkadams4   PetscCall(PetscMalloc2(nloc, &permute, nloc, &degree));
11359566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nloc, &bIndexSet));
11366e09b0e3SStefano Zampini   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
11379566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random));
11389566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
1139c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
1140bae903cbSmarkadams4     PetscInt nc;
1141*e0b7e82fSBarry Smith 
1142bae903cbSmarkadams4     PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
1143bae903cbSmarkadams4     degree[Ii] = nc;
1144bae903cbSmarkadams4     PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
1145bae903cbSmarkadams4   }
1146bae903cbSmarkadams4   for (Ii = 0; Ii < nloc; Ii++) {
11479566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(random, &hashfact));
1148aea53286SMark Adams     iSwapIndex = (PetscInt)(hashfact * nloc) % nloc;
1149c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
1150c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
1151*e0b7e82fSBarry Smith 
1152c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
1153c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
1154bae903cbSmarkadams4       iTemp                 = degree[iSwapIndex];
1155bae903cbSmarkadams4       degree[iSwapIndex]    = degree[Ii];
1156bae903cbSmarkadams4       degree[Ii]            = iTemp;
1157c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
1158c8b0795cSMark F. Adams     }
1159c8b0795cSMark F. Adams   }
1160d529f056Smarkadams4   // apply minimum degree ordering -- NEW
1161d529f056Smarkadams4   if (pc_gamg_agg->use_minimum_degree_ordering) { PetscCall(PetscSortIntWithArray(nloc, degree, permute)); }
11629566063dSJacob Faibussowitsch   PetscCall(PetscFree(bIndexSet));
11639566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&random));
11649566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
1165849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
1166d529f056Smarkadams4   // square graph
1167d529f056Smarkadams4   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels && pc_gamg_agg->use_aggressive_square_graph) {
1168d529f056Smarkadams4     PetscCall(PCGAMGSquareGraph_GAMG(a_pc, Gmat1, &Gmat2));
1169d529f056Smarkadams4   } else Gmat2 = Gmat1;
1170d529f056Smarkadams4   // switch to old MIS-1 for square graph
1171d529f056Smarkadams4   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) {
1172d529f056Smarkadams4     if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, pc_gamg_agg->aggressive_mis_k)); // hardwire to MIS-2
1173d529f056Smarkadams4     else PetscCall(MatCoarsenSetType(pc_gamg_agg->crs, MATCOARSENMIS));                                                                   // old MIS -- side effect
1174d529f056Smarkadams4   } else if (pc_gamg_agg->use_aggressive_square_graph && pc_gamg_agg->aggressive_coarsening_levels > 0) {                                 // we reset the MIS
1175d529f056Smarkadams4     const char *prefix;
1176*e0b7e82fSBarry Smith 
1177d529f056Smarkadams4     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)a_pc, &prefix));
1178d529f056Smarkadams4     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix));
1179d529f056Smarkadams4     PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); // get the default back on non-aggressive levels when square graph switched to old MIS
1180d529f056Smarkadams4   }
1181d529f056Smarkadams4   PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat2));
11822d776b49SBarry Smith   PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE));
1183d529f056Smarkadams4   PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm));
11842d776b49SBarry Smith   PetscCall(MatCoarsenApply(pc_gamg_agg->crs));
11852d776b49SBarry Smith   PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */
1186b43b03e9SMark F. Adams 
11879566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1188bae903cbSmarkadams4   PetscCall(PetscFree2(permute, degree));
1189849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
11909f3f12c8SMark F. Adams 
11919c15c1aeSMark Adams   if (Gmat2 != Gmat1) { // square graph, we need ghosts for selected
1192d529f056Smarkadams4     PetscCoarsenData *llist = *agg_lists;
1193*e0b7e82fSBarry Smith 
1194d529f056Smarkadams4     PetscCall(fixAggregatesWithSquare(a_pc, Gmat2, Gmat1, *agg_lists));
1195d529f056Smarkadams4     PetscCall(MatDestroy(&Gmat1));
1196d529f056Smarkadams4     *a_Gmat1 = Gmat2;                          /* output */
11978926f930SMark Adams     PetscCall(PetscCDSetMat(llist, *a_Gmat1)); /* Need a graph with ghosts here */
11980cbbd2e1SMark F. Adams   }
1199849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
12003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1201c8b0795cSMark F. Adams }
1202c8b0795cSMark F. Adams 
1203c8b0795cSMark F. Adams /*
1204*e0b7e82fSBarry Smith  PCGAMGConstructProlongator_AGG
1205c8b0795cSMark F. Adams 
1206c8b0795cSMark F. Adams  Input Parameter:
1207c8b0795cSMark F. Adams  . pc - this
1208c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1209c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
12100cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1211c8b0795cSMark F. Adams  Output Parameter:
1212c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1213c8b0795cSMark F. Adams  */
1214*e0b7e82fSBarry Smith static PetscErrorCode PCGAMGConstructProlongator_AGG(PC pc, Mat Amat, PetscCoarsenData *agg_lists, Mat *a_P_out)
1215d71ae5a4SJacob Faibussowitsch {
12162e68589bSMark F. Adams   PC_MG         *mg      = (PC_MG *)pc->data;
12172e68589bSMark F. Adams   PC_GAMG       *pc_gamg = (PC_GAMG *)mg->innerctx;
12184a2f8832SBarry Smith   const PetscInt col_bs  = pc_gamg->data_cell_cols;
1219c8b0795cSMark F. Adams   PetscInt       Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs;
12208926f930SMark Adams   Mat            Gmat, Prol;
1221d5d25401SBarry Smith   PetscMPIInt    size;
12223b4367a7SBarry Smith   MPI_Comm       comm;
1223c8b0795cSMark F. Adams   PetscReal     *data_w_ghost;
1224c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes = 0, *flid_fgid;
1225d9558ea9SBarry Smith   MatType        mtype;
12262e68589bSMark F. Adams 
12272e68589bSMark F. Adams   PetscFunctionBegin;
12289566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
122908401ef6SPierre Jolivet   PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1");
1230849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
12319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
12329566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
12339566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
12349371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
12359371c9d4SSatish Balay   my0  = Istart / bs;
123663a3b9bcSJacob 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);
1237d8b4a066SPierre Jolivet   PetscCall(PetscCDGetMat(agg_lists, &Gmat)); // get auxiliary matrix for ghost edges for size > 1
12382e68589bSMark F. Adams 
12392e68589bSMark F. Adams   /* get 'nLocalSelected' */
12400cbbd2e1SMark F. Adams   for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) {
1241e78576d6SMark F. Adams     PetscBool ise;
1242*e0b7e82fSBarry Smith 
12430cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
12445e62d202SMark Adams     PetscCall(PetscCDIsEmptyAt(agg_lists, ii, &ise));
1245e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
12462e68589bSMark F. Adams   }
12472e68589bSMark F. Adams 
12482e68589bSMark F. Adams   /* create prolongator, create P matrix */
12499566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat, &mtype));
12509566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
12519566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE));
125244eff04eSMark Adams   PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); // should this be before MatSetSizes?
12539566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
12543742c8caSstefanozampini #if PetscDefined(HAVE_DEVICE)
12553742c8caSstefanozampini   PetscBool flg;
12563742c8caSstefanozampini   PetscCall(MatBoundToCPU(Amat, &flg));
12573742c8caSstefanozampini   PetscCall(MatBindToCPU(Prol, flg));
12583742c8caSstefanozampini   if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE));
12593742c8caSstefanozampini #endif
12609566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL));
12619566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL));
12622e68589bSMark F. Adams 
12632e68589bSMark F. Adams   /* can get all points "removed" */
12649566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &ii));
12657f66b68fSMark Adams   if (!ii) {
126663a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix));
12679566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
12680298fd71SBarry Smith     *a_P_out = NULL; /* out */
1269849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
12703ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
12712e68589bSMark F. Adams   }
127263a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs));
12739566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
12740cbbd2e1SMark F. Adams 
127563a3b9bcSJacob 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);
1276c8b0795cSMark F. Adams   myCrs0 = myCrs0 / col_bs;
127763a3b9bcSJacob 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);
12782e68589bSMark F. Adams 
12792e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
1280849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
1281bae903cbSmarkadams4   if (size > 1) { /* get ghost null space data */
12822e68589bSMark F. Adams     PetscReal *tmp_gdata, *tmp_ldata, *tp2;
1283*e0b7e82fSBarry Smith 
12849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
12854a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
1286c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1287a2f3521dSMark F. Adams         PetscInt         ii, stride;
12888e3a54c0SPierre Jolivet         const PetscReal *tp = PetscSafePointerPlusOffset(pc_gamg->data, jj * bs * nloc + kk);
1289*e0b7e82fSBarry Smith 
12902fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
12912fa5cd67SKarl Rupp 
12929566063dSJacob Faibussowitsch         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
1293a2f3521dSMark F. Adams 
1294bae903cbSmarkadams4         if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */
12959566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost));
1296a2f3521dSMark F. Adams           nbnodes = bs * stride;
12972e68589bSMark F. Adams         }
12988e3a54c0SPierre Jolivet         tp2 = PetscSafePointerPlusOffset(data_w_ghost, jj * bs * stride + kk);
12992fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
13009566063dSJacob Faibussowitsch         PetscCall(PetscFree(tmp_gdata));
13012e68589bSMark F. Adams       }
13022e68589bSMark F. Adams     }
13039566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_ldata));
1304806fa848SBarry Smith   } else {
1305c8b0795cSMark F. Adams     nbnodes      = bs * nloc;
1306c8b0795cSMark F. Adams     data_w_ghost = (PetscReal *)pc_gamg->data;
13072e68589bSMark F. Adams   }
13082e68589bSMark F. Adams 
1309bae903cbSmarkadams4   /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */
1310c5df96a5SBarry Smith   if (size > 1) {
13112e68589bSMark F. Adams     PetscReal *fid_glid_loc, *fiddata;
1312a2f3521dSMark F. Adams     PetscInt   stride;
13132e68589bSMark F. Adams 
13149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
13152e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk);
13169566063dSJacob Faibussowitsch     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
1317bae903cbSmarkadams4     PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */
1318a2f3521dSMark F. Adams     for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
13199566063dSJacob Faibussowitsch     PetscCall(PetscFree(fiddata));
1320a2f3521dSMark F. Adams 
132163a3b9bcSJacob Faibussowitsch     PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs);
13229566063dSJacob Faibussowitsch     PetscCall(PetscFree(fid_glid_loc));
1323806fa848SBarry Smith   } else {
13249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &flid_fgid));
13252e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk;
13262e68589bSMark F. Adams   }
1327849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
1328bae903cbSmarkadams4   /* get P0 */
1329849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
1330c8b0795cSMark F. Adams   {
13310298fd71SBarry Smith     PetscReal *data_out = NULL;
1332*e0b7e82fSBarry Smith 
13339566063dSJacob Faibussowitsch     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol));
13349566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
13352fa5cd67SKarl Rupp 
1336c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
13374a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
13384a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs * col_bs * nLocalSelected;
1339c8b0795cSMark F. Adams   }
1340849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
134148a46eb9SPierre Jolivet   if (size > 1) PetscCall(PetscFree(data_w_ghost));
13429566063dSJacob Faibussowitsch   PetscCall(PetscFree(flid_fgid));
13432e68589bSMark F. Adams 
1344c8b0795cSMark F. Adams   *a_P_out = Prol; /* out */
1345*e0b7e82fSBarry Smith   PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_initial_prolongation"));
1346ed4e82eaSPeter Brune 
1347849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
13483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1349c8b0795cSMark F. Adams }
1350c8b0795cSMark F. Adams 
1351c8b0795cSMark F. Adams /*
1352*e0b7e82fSBarry Smith    PCGAMGOptimizeProlongator_AGG - given the initial prolongator optimizes it by smoothed aggregation pc_gamg_agg->nsmooths times
1353c8b0795cSMark F. Adams 
1354c8b0795cSMark F. Adams   Input Parameter:
1355c8b0795cSMark F. Adams    . pc - this
1356c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1357c8b0795cSMark F. Adams  In/Output Parameter:
13582a808120SBarry Smith    . a_P - prolongation operator to the next level
1359c8b0795cSMark F. Adams */
1360*e0b7e82fSBarry Smith static PetscErrorCode PCGAMGOptimizeProlongator_AGG(PC pc, Mat Amat, Mat *a_P)
1361d71ae5a4SJacob Faibussowitsch {
1362c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
1363c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
1364c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
1365c8b0795cSMark F. Adams   PetscInt     jj;
1366c8b0795cSMark F. Adams   Mat          Prol = *a_P;
13673b4367a7SBarry Smith   MPI_Comm     comm;
13682a808120SBarry Smith   KSP          eksp;
13692a808120SBarry Smith   Vec          bb, xx;
13702a808120SBarry Smith   PC           epc;
13712a808120SBarry Smith   PetscReal    alpha, emax, emin;
1372c8b0795cSMark F. Adams 
1373c8b0795cSMark F. Adams   PetscFunctionBegin;
13749566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
1375849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
1376c8b0795cSMark F. Adams 
1377d5d25401SBarry Smith   /* compute maximum singular value of operator to be used in smoother */
13782a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
137918c3aa7eSMark     /* get eigen estimates */
138018c3aa7eSMark     if (pc_gamg->emax > 0) {
138118c3aa7eSMark       emin = pc_gamg->emin;
138218c3aa7eSMark       emax = pc_gamg->emax;
138318c3aa7eSMark     } else {
13840ed2132dSStefano Zampini       const char *prefix;
13850ed2132dSStefano Zampini 
13869566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &bb, NULL));
13879566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &xx, NULL));
13889566063dSJacob Faibussowitsch       PetscCall(KSPSetNoisy_Private(bb));
1389e696ed0bSMark F. Adams 
13909566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm, &eksp));
13913821be0aSBarry Smith       PetscCall(KSPSetNestLevel(eksp, pc->kspnestlevel));
13929566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
13939566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(eksp, prefix));
13949566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_"));
139573f7197eSJed Brown       {
1396b94d7dedSBarry Smith         PetscBool isset, sflg;
1397*e0b7e82fSBarry Smith 
1398b94d7dedSBarry Smith         PetscCall(MatIsSPDKnown(Amat, &isset, &sflg));
1399b94d7dedSBarry Smith         if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG));
1400d24ecf33SMark       }
14019566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure));
14029566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
1403c2436cacSMark F. Adams 
14049566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
14059566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(eksp, Amat, Amat));
14062e68589bSMark F. Adams 
14079566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(eksp, &epc));
14089566063dSJacob Faibussowitsch       PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */
1409c2436cacSMark F. Adams 
14109566063dSJacob Faibussowitsch       PetscCall(KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10)); // 10 is safer, but 5 is often fine, can override with -pc_gamg_esteig_ksp_max_it -mg_levels_ksp_chebyshev_esteig 0,0.25,0,1.2
1411074aec55SMark Adams 
14129566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(eksp));
14139566063dSJacob Faibussowitsch       PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE));
14149566063dSJacob Faibussowitsch       PetscCall(KSPSolve(eksp, bb, xx));
14159566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(eksp, pc, xx));
14162e68589bSMark F. Adams 
14179566063dSJacob Faibussowitsch       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
14189566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI));
14199566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&xx));
14209566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bb));
14219566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&eksp));
14222e68589bSMark F. Adams     }
142318c3aa7eSMark     if (pc_gamg->use_sa_esteig) {
142418c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
142518c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
142663a3b9bcSJacob 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));
142718c3aa7eSMark     } else {
142818c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
142918c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
143018c3aa7eSMark     }
143118c3aa7eSMark   } else {
143218c3aa7eSMark     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
143318c3aa7eSMark     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
143418c3aa7eSMark   }
14352e68589bSMark F. Adams 
14362a808120SBarry Smith   /* smooth P0 */
1437*e0b7e82fSBarry Smith   if (pc_gamg_agg->nsmooths > 0) {
14382a808120SBarry Smith     Vec diag;
14392a808120SBarry Smith 
1440*e0b7e82fSBarry Smith     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
1441*e0b7e82fSBarry Smith     PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero");
14422a808120SBarry Smith 
1443*e0b7e82fSBarry Smith     PetscCall(MatCreateVecs(Amat, &diag, NULL));
1444*e0b7e82fSBarry Smith     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
1445*e0b7e82fSBarry Smith     PetscCall(VecReciprocal(diag));
1446*e0b7e82fSBarry Smith 
1447*e0b7e82fSBarry Smith     for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
1448*e0b7e82fSBarry Smith       Mat tMat;
1449*e0b7e82fSBarry Smith 
1450*e0b7e82fSBarry Smith       PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
1451*e0b7e82fSBarry Smith       /*
1452*e0b7e82fSBarry Smith         Smooth aggregation on the prolongator
1453*e0b7e82fSBarry Smith 
1454*e0b7e82fSBarry Smith         P_{i} := (I - 1.4/emax D^{-1}A) P_i\{i-1}
1455*e0b7e82fSBarry Smith       */
14569566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
14579566063dSJacob Faibussowitsch       PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat));
14589566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
14599566063dSJacob Faibussowitsch       PetscCall(MatProductClear(tMat));
14609566063dSJacob Faibussowitsch       PetscCall(MatDiagonalScale(tMat, diag, NULL));
1461b4da3a1bSStefano Zampini 
1462d5d25401SBarry Smith       /* TODO: Document the 1.4 and don't hardwire it in this routine */
1463b4ec6429SMark F. Adams       alpha = -1.4 / emax;
14649566063dSJacob Faibussowitsch       PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
14659566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Prol));
14662e68589bSMark F. Adams       Prol = tMat;
1467849bee69SMark Adams       PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
14682e68589bSMark F. Adams     }
1469*e0b7e82fSBarry Smith     PetscCall(VecDestroy(&diag));
1470*e0b7e82fSBarry Smith   }
1471849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
1472*e0b7e82fSBarry Smith   PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_prolongation"));
1473c8b0795cSMark F. Adams   *a_P = Prol;
14743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14752e68589bSMark F. Adams }
14762e68589bSMark F. Adams 
1477a077d33dSBarry Smith /*MC
1478a077d33dSBarry Smith   PCGAMGAGG - Smooth aggregation, {cite}`vanek1996algebraic`, {cite}`vanek2001convergence`, variant of PETSc's algebraic multigrid (`PCGAMG`) preconditioner
14792e68589bSMark F. Adams 
1480a077d33dSBarry Smith   Options Database Keys:
1481a077d33dSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation to construct prolongation
1482a077d33dSBarry Smith . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest.
1483a077d33dSBarry Smith . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening
1484a077d33dSBarry Smith . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm
1485a077d33dSBarry Smith . -pc_gamg_pc_gamg_asm_hem_aggs <n,default=0> - Number of HEM aggregation steps for ASM smoother
1486a077d33dSBarry Smith - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive')
1487a077d33dSBarry Smith 
1488a077d33dSBarry Smith   Level: intermediate
1489a077d33dSBarry Smith 
1490a077d33dSBarry Smith   Notes:
1491a077d33dSBarry Smith   To obtain good performance for `PCGAMG` for vector valued problems you must
1492a077d33dSBarry Smith   call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point.
1493a077d33dSBarry Smith   Call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator
1494a077d33dSBarry Smith 
1495a077d33dSBarry Smith   The many options for `PCMG` and `PCGAMG` such as controlling the smoothers on each level etc. also work for `PCGAMGAGG`
1496a077d33dSBarry Smith 
1497a077d33dSBarry Smith .seealso: `PCGAMG`, [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCCreate()`, `PCSetType()`,
1498a077d33dSBarry Smith           `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`,
1499a077d33dSBarry Smith           `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`,
1500a077d33dSBarry Smith           `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`,
1501a077d33dSBarry Smith           `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()`
1502a077d33dSBarry Smith M*/
1503d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc)
1504d71ae5a4SJacob Faibussowitsch {
1505c8b0795cSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
1506c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
1507c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg;
15082e68589bSMark F. Adams 
1509c8b0795cSMark F. Adams   PetscFunctionBegin;
1510c8b0795cSMark F. Adams   /* create sub context for SA */
15114dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg_agg));
1512c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1513c8b0795cSMark F. Adams 
15141ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
15159b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1516c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1517c8b0795cSMark F. Adams 
1518c8b0795cSMark F. Adams   /* set internal function pointers */
15192d776b49SBarry Smith   pc_gamg->ops->creategraph       = PCGAMGCreateGraph_AGG;
15201ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
1521*e0b7e82fSBarry Smith   pc_gamg->ops->prolongator       = PCGAMGConstructProlongator_AGG;
1522*e0b7e82fSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptimizeProlongator_AGG;
15231ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
15245adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
1525c8b0795cSMark F. Adams 
1526cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths                     = 1;
1527d529f056Smarkadams4   pc_gamg_agg->aggressive_coarsening_levels = 1;
1528d4adc10fSMark Adams   pc_gamg_agg->use_aggressive_square_graph  = PETSC_TRUE;
1529d529f056Smarkadams4   pc_gamg_agg->use_minimum_degree_ordering  = PETSC_FALSE;
1530a9ccf005SMark Adams   pc_gamg_agg->use_low_mem_filter           = PETSC_FALSE;
1531d529f056Smarkadams4   pc_gamg_agg->aggressive_mis_k             = 2;
1532cab9ed1eSBarry Smith 
15339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG));
1534bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG));
1535d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", PCGAMGSetAggressiveSquareGraph_AGG));
1536d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", PCGAMGMISkSetMinDegreeOrdering_AGG));
1537a9ccf005SMark Adams   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", PCGAMGSetLowMemoryFilter_AGG));
1538d529f056Smarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", PCGAMGMISkSetAggressive_AGG));
15399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG));
15403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15412e68589bSMark F. Adams }
1542