xref: /petsc/src/mat/graphops/coarsen/interface/coarsen.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
1a77d671bSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
2a77d671bSBarry Smith 
3a77d671bSBarry Smith /* Logging support */
4a77d671bSBarry Smith PetscClassId MAT_COARSEN_CLASSID;
5a77d671bSBarry Smith 
6a77d671bSBarry Smith PetscFunctionList MatCoarsenList              = NULL;
7a77d671bSBarry Smith PetscBool         MatCoarsenRegisterAllCalled = PETSC_FALSE;
8a77d671bSBarry Smith 
9a77d671bSBarry Smith /*@C
10a77d671bSBarry Smith   MatCoarsenRegister - Adds a new sparse matrix coarsening algorithm to the matrix package.
11a77d671bSBarry Smith 
12a77d671bSBarry Smith   Logically Collective, No Fortran Support
13a77d671bSBarry Smith 
14a77d671bSBarry Smith   Input Parameters:
15a77d671bSBarry Smith + sname    - name of coarsen (for example `MATCOARSENMIS`)
16a77d671bSBarry Smith - function - function pointer that creates the coarsen type
17a77d671bSBarry Smith 
18a77d671bSBarry Smith   Level: developer
19a77d671bSBarry Smith 
20a77d671bSBarry Smith   Example Usage:
21a77d671bSBarry Smith .vb
22a77d671bSBarry Smith    MatCoarsenRegister("my_agg", MyAggCreate);
23a77d671bSBarry Smith .ve
24a77d671bSBarry Smith 
25a77d671bSBarry Smith   Then, your aggregator can be chosen with the procedural interface via `MatCoarsenSetType(agg, "my_agg")` or at runtime via the option `-mat_coarsen_type my_agg`
26a77d671bSBarry Smith 
27a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenSetType()`, `MatCoarsenCreate()`, `MatCoarsenRegisterDestroy()`, `MatCoarsenRegisterAll()`
28a77d671bSBarry Smith @*/
MatCoarsenRegister(const char sname[],PetscErrorCode (* function)(MatCoarsen))29a77d671bSBarry Smith PetscErrorCode MatCoarsenRegister(const char sname[], PetscErrorCode (*function)(MatCoarsen))
30a77d671bSBarry Smith {
31a77d671bSBarry Smith   PetscFunctionBegin;
32a77d671bSBarry Smith   PetscCall(MatInitializePackage());
33a77d671bSBarry Smith   PetscCall(PetscFunctionListAdd(&MatCoarsenList, sname, function));
34a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
35a77d671bSBarry Smith }
36a77d671bSBarry Smith 
37a77d671bSBarry Smith /*@
38a77d671bSBarry Smith   MatCoarsenGetType - Gets the Coarsen method type and name (as a string)
39a77d671bSBarry Smith   from the coarsen context.
40a77d671bSBarry Smith 
41a77d671bSBarry Smith   Not Collective
42a77d671bSBarry Smith 
43a77d671bSBarry Smith   Input Parameter:
44a77d671bSBarry Smith . coarsen - the coarsen context
45a77d671bSBarry Smith 
46a77d671bSBarry Smith   Output Parameter:
47a77d671bSBarry Smith . type - coarsener type
48a77d671bSBarry Smith 
49a77d671bSBarry Smith   Level: advanced
50a77d671bSBarry Smith 
51a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenCreate()`, `MatCoarsenType`, `MatCoarsenSetType()`, `MatCoarsenRegister()`
52a77d671bSBarry Smith @*/
MatCoarsenGetType(MatCoarsen coarsen,MatCoarsenType * type)53a77d671bSBarry Smith PetscErrorCode MatCoarsenGetType(MatCoarsen coarsen, MatCoarsenType *type)
54a77d671bSBarry Smith {
55a77d671bSBarry Smith   PetscFunctionBegin;
56a77d671bSBarry Smith   PetscValidHeaderSpecific(coarsen, MAT_COARSEN_CLASSID, 1);
57a77d671bSBarry Smith   PetscAssertPointer(type, 2);
58a77d671bSBarry Smith   *type = ((PetscObject)coarsen)->type_name;
59a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
60a77d671bSBarry Smith }
61a77d671bSBarry Smith 
62a77d671bSBarry Smith /*@
63a77d671bSBarry Smith   MatCoarsenApply - Gets a coarsen for a matrix.
64a77d671bSBarry Smith 
65a77d671bSBarry Smith   Collective
66a77d671bSBarry Smith 
67a77d671bSBarry Smith   Input Parameter:
68a77d671bSBarry Smith . coarser - the coarsen
69a77d671bSBarry Smith 
70a77d671bSBarry Smith   Options Database Keys:
71a77d671bSBarry Smith + -mat_coarsen_type mis|hem|misk - mis: maximal independent set based; misk: distance k MIS; hem: heavy edge matching
72a77d671bSBarry Smith - -mat_coarsen_view              - view the coarsening object
73a77d671bSBarry Smith 
74a77d671bSBarry Smith   Level: advanced
75a77d671bSBarry Smith 
76a77d671bSBarry Smith   Notes:
77a77d671bSBarry Smith   When the coarsening is used inside `PCGAMG` then the options database keys are prefixed with `-pc_gamg_`
78a77d671bSBarry Smith 
79a77d671bSBarry Smith   Use `MatCoarsenGetData()` to access the results of the coarsening
80a77d671bSBarry Smith 
81a77d671bSBarry Smith   The user can define additional coarsens; see `MatCoarsenRegister()`.
82a77d671bSBarry Smith 
83eca7e54bSSatish Balay .seealso: `MatCoarsen`, `MatCoarsenSetFromOptions()`, `MatCoarsenSetType()`, `MatCoarsenRegister()`, `MatCoarsenCreate()`,
84a77d671bSBarry Smith           `MatCoarsenDestroy()`, `MatCoarsenSetAdjacency()`
85a77d671bSBarry Smith           `MatCoarsenGetData()`
86a77d671bSBarry Smith @*/
MatCoarsenApply(MatCoarsen coarser)87a77d671bSBarry Smith PetscErrorCode MatCoarsenApply(MatCoarsen coarser)
88a77d671bSBarry Smith {
89a77d671bSBarry Smith   PetscFunctionBegin;
90a77d671bSBarry Smith   PetscValidHeaderSpecific(coarser, MAT_COARSEN_CLASSID, 1);
91a77d671bSBarry Smith   PetscAssertPointer(coarser, 1);
92a77d671bSBarry Smith   PetscCheck(coarser->graph->assembled, PetscObjectComm((PetscObject)coarser), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
93a77d671bSBarry Smith   PetscCheck(!coarser->graph->factortype, PetscObjectComm((PetscObject)coarser), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
94a77d671bSBarry Smith   PetscCall(PetscLogEventBegin(MAT_Coarsen, coarser, 0, 0, 0));
95a77d671bSBarry Smith   PetscUseTypeMethod(coarser, apply);
96a77d671bSBarry Smith   PetscCall(PetscLogEventEnd(MAT_Coarsen, coarser, 0, 0, 0));
97a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
98a77d671bSBarry Smith }
99a77d671bSBarry Smith 
100a77d671bSBarry Smith /*@
101a77d671bSBarry Smith   MatCoarsenSetAdjacency - Sets the adjacency graph (matrix) of the thing to be coarsened.
102a77d671bSBarry Smith 
103a77d671bSBarry Smith   Collective
104a77d671bSBarry Smith 
105a77d671bSBarry Smith   Input Parameters:
106a77d671bSBarry Smith + agg - the coarsen context
107a77d671bSBarry Smith - adj - the adjacency matrix
108a77d671bSBarry Smith 
109a77d671bSBarry Smith   Level: advanced
110a77d671bSBarry Smith 
111a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenSetFromOptions()`, `Mat`, `MatCoarsenCreate()`, `MatCoarsenApply()`
112a77d671bSBarry Smith @*/
MatCoarsenSetAdjacency(MatCoarsen agg,Mat adj)113a77d671bSBarry Smith PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen agg, Mat adj)
114a77d671bSBarry Smith {
115a77d671bSBarry Smith   PetscFunctionBegin;
116a77d671bSBarry Smith   PetscValidHeaderSpecific(agg, MAT_COARSEN_CLASSID, 1);
117a77d671bSBarry Smith   PetscValidHeaderSpecific(adj, MAT_CLASSID, 2);
118a77d671bSBarry Smith   agg->graph = adj;
119a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
120a77d671bSBarry Smith }
121a77d671bSBarry Smith 
122a77d671bSBarry Smith /*@
123a77d671bSBarry Smith   MatCoarsenSetStrictAggs - Set whether to keep strict (non overlapping) aggregates in the linked list of aggregates for a coarsen context
124a77d671bSBarry Smith 
125a77d671bSBarry Smith   Logically Collective
126a77d671bSBarry Smith 
127a77d671bSBarry Smith   Input Parameters:
128a77d671bSBarry Smith + agg - the coarsen context
129a77d671bSBarry Smith - str - `PETSC_TRUE` keep strict aggregates, `PETSC_FALSE` allow overlap
130a77d671bSBarry Smith 
131a77d671bSBarry Smith   Level: advanced
132a77d671bSBarry Smith 
133a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenCreate()`, `MatCoarsenSetFromOptions()`
134a77d671bSBarry Smith @*/
MatCoarsenSetStrictAggs(MatCoarsen agg,PetscBool str)135a77d671bSBarry Smith PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen agg, PetscBool str)
136a77d671bSBarry Smith {
137a77d671bSBarry Smith   PetscFunctionBegin;
138a77d671bSBarry Smith   PetscValidHeaderSpecific(agg, MAT_COARSEN_CLASSID, 1);
139a77d671bSBarry Smith   agg->strict_aggs = str;
140a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
141a77d671bSBarry Smith }
142a77d671bSBarry Smith 
143a77d671bSBarry Smith /*@
144a77d671bSBarry Smith   MatCoarsenDestroy - Destroys the coarsen context.
145a77d671bSBarry Smith 
146a77d671bSBarry Smith   Collective
147a77d671bSBarry Smith 
148a77d671bSBarry Smith   Input Parameter:
149a77d671bSBarry Smith . agg - the coarsen context
150a77d671bSBarry Smith 
151a77d671bSBarry Smith   Level: advanced
152a77d671bSBarry Smith 
153a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenCreate()`
154a77d671bSBarry Smith @*/
MatCoarsenDestroy(MatCoarsen * agg)155a77d671bSBarry Smith PetscErrorCode MatCoarsenDestroy(MatCoarsen *agg)
156a77d671bSBarry Smith {
157a77d671bSBarry Smith   PetscFunctionBegin;
158a77d671bSBarry Smith   if (!*agg) PetscFunctionReturn(PETSC_SUCCESS);
159a77d671bSBarry Smith   PetscValidHeaderSpecific(*agg, MAT_COARSEN_CLASSID, 1);
160a77d671bSBarry Smith   if (--((PetscObject)*agg)->refct > 0) {
161a77d671bSBarry Smith     *agg = NULL;
162a77d671bSBarry Smith     PetscFunctionReturn(PETSC_SUCCESS);
163a77d671bSBarry Smith   }
164a77d671bSBarry Smith 
165a77d671bSBarry Smith   PetscTryTypeMethod(*agg, destroy);
166a77d671bSBarry Smith   if ((*agg)->agg_lists) PetscCall(PetscCDDestroy((*agg)->agg_lists));
167a77d671bSBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)*agg, "MatCoarsenSetMaximumIterations_C", NULL));
168a77d671bSBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)*agg, "MatCoarsenSetThreshold_C", NULL));
169a77d671bSBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)*agg, "MatCoarsenSetStrengthIndex_C", NULL));
170a77d671bSBarry Smith 
171a77d671bSBarry Smith   PetscCall(PetscHeaderDestroy(agg));
172a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
173a77d671bSBarry Smith }
174a77d671bSBarry Smith 
175a77d671bSBarry Smith /*@
176a77d671bSBarry Smith   MatCoarsenViewFromOptions - View the coarsener from the options database
177a77d671bSBarry Smith 
178a77d671bSBarry Smith   Collective
179a77d671bSBarry Smith 
180a77d671bSBarry Smith   Input Parameters:
181a77d671bSBarry Smith + A    - the coarsen context
182a77d671bSBarry Smith . obj  - Optional object that provides the prefix for the option name
183a77d671bSBarry Smith - name - command line option (usually `-mat_coarsen_view`)
184a77d671bSBarry Smith 
185a77d671bSBarry Smith   Options Database Key:
186a77d671bSBarry Smith . -mat_coarsen_view [viewertype]:... - the viewer and its options
187a77d671bSBarry Smith 
188a77d671bSBarry Smith   Note:
189a77d671bSBarry Smith .vb
190a77d671bSBarry Smith     If no value is provided ascii:stdout is used
191a77d671bSBarry Smith        ascii[:[filename][:[format][:append]]]    defaults to stdout - format can be one of ascii_info, ascii_info_detail, or ascii_matlab,
192a77d671bSBarry Smith                                                   for example ascii::ascii_info prints just the information about the object not all details
193a77d671bSBarry Smith                                                   unless :append is given filename opens in write mode, overwriting what was already there
194a77d671bSBarry Smith        binary[:[filename][:[format][:append]]]   defaults to the file binaryoutput
195a77d671bSBarry Smith        draw[:drawtype[:filename]]                for example, draw:tikz, draw:tikz:figure.tex  or draw:x
196a77d671bSBarry Smith        socket[:port]                             defaults to the standard output port
197a77d671bSBarry Smith        saws[:communicatorname]                    publishes object to the Scientific Application Webserver (SAWs)
198a77d671bSBarry Smith .ve
199a77d671bSBarry Smith 
200a77d671bSBarry Smith   Level: intermediate
201a77d671bSBarry Smith 
202a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenView`, `PetscObjectViewFromOptions()`, `MatCoarsenCreate()`
203a77d671bSBarry Smith @*/
MatCoarsenViewFromOptions(MatCoarsen A,PetscObject obj,const char name[])204a77d671bSBarry Smith PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A, PetscObject obj, const char name[])
205a77d671bSBarry Smith {
206a77d671bSBarry Smith   PetscFunctionBegin;
207a77d671bSBarry Smith   PetscValidHeaderSpecific(A, MAT_COARSEN_CLASSID, 1);
208a77d671bSBarry Smith   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
209a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
210a77d671bSBarry Smith }
211a77d671bSBarry Smith 
212a77d671bSBarry Smith /*@
213a77d671bSBarry Smith   MatCoarsenView - Prints the coarsen data structure.
214a77d671bSBarry Smith 
215a77d671bSBarry Smith   Collective
216a77d671bSBarry Smith 
217a77d671bSBarry Smith   Input Parameters:
218a77d671bSBarry Smith + agg    - the coarsen context
219a77d671bSBarry Smith - viewer - optional visualization context
220a77d671bSBarry Smith 
221a77d671bSBarry Smith    For viewing the options database see `MatCoarsenViewFromOptions()`
222a77d671bSBarry Smith 
223a77d671bSBarry Smith   Level: advanced
224a77d671bSBarry Smith 
225a77d671bSBarry Smith .seealso: `MatCoarsen`, `PetscViewer`, `PetscViewerASCIIOpen()`, `MatCoarsenViewFromOptions`
226a77d671bSBarry Smith @*/
MatCoarsenView(MatCoarsen agg,PetscViewer viewer)227a77d671bSBarry Smith PetscErrorCode MatCoarsenView(MatCoarsen agg, PetscViewer viewer)
228a77d671bSBarry Smith {
229*9f196a02SMartin Diehl   PetscBool isascii;
230a77d671bSBarry Smith 
231a77d671bSBarry Smith   PetscFunctionBegin;
232a77d671bSBarry Smith   PetscValidHeaderSpecific(agg, MAT_COARSEN_CLASSID, 1);
233a77d671bSBarry Smith   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)agg), &viewer));
234a77d671bSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
235a77d671bSBarry Smith   PetscCheckSameComm(agg, 1, viewer, 2);
236a77d671bSBarry Smith 
237*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
238a77d671bSBarry Smith   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)agg, viewer));
239a77d671bSBarry Smith   if (agg->ops->view) {
240a77d671bSBarry Smith     PetscCall(PetscViewerASCIIPushTab(viewer));
241a77d671bSBarry Smith     PetscUseTypeMethod(agg, view, viewer);
242a77d671bSBarry Smith     PetscCall(PetscViewerASCIIPopTab(viewer));
243a77d671bSBarry Smith   }
244835f2295SStefano Zampini   if (agg->strength_index_size > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " Using scalar strength-of-connection index[%" PetscInt_FMT "] = {%" PetscInt_FMT ", ..}\n", agg->strength_index_size, agg->strength_index[0]));
245a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
246a77d671bSBarry Smith }
247a77d671bSBarry Smith 
248a77d671bSBarry Smith /*@
249a77d671bSBarry Smith   MatCoarsenSetType - Sets the type of aggregator to use
250a77d671bSBarry Smith 
251a77d671bSBarry Smith   Collective
252a77d671bSBarry Smith 
253a77d671bSBarry Smith   Input Parameters:
254a77d671bSBarry Smith + coarser - the coarsen context.
255a77d671bSBarry Smith - type    - a known coarsening method
256a77d671bSBarry Smith 
257a77d671bSBarry Smith   Options Database Key:
258a77d671bSBarry Smith . -mat_coarsen_type  <type> - maximal independent set based; distance k MIS; heavy edge matching
259a77d671bSBarry Smith 
260a77d671bSBarry Smith   Level: advanced
261a77d671bSBarry Smith 
262a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenCreate()`, `MatCoarsenApply()`, `MatCoarsenType`, `MatCoarsenGetType()`
263a77d671bSBarry Smith @*/
MatCoarsenSetType(MatCoarsen coarser,MatCoarsenType type)264a77d671bSBarry Smith PetscErrorCode MatCoarsenSetType(MatCoarsen coarser, MatCoarsenType type)
265a77d671bSBarry Smith {
266a77d671bSBarry Smith   PetscBool match;
267a77d671bSBarry Smith   PetscErrorCode (*r)(MatCoarsen);
268a77d671bSBarry Smith 
269a77d671bSBarry Smith   PetscFunctionBegin;
270a77d671bSBarry Smith   PetscValidHeaderSpecific(coarser, MAT_COARSEN_CLASSID, 1);
271a77d671bSBarry Smith   PetscAssertPointer(type, 2);
272a77d671bSBarry Smith 
273a77d671bSBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)coarser, type, &match));
274a77d671bSBarry Smith   if (match) PetscFunctionReturn(PETSC_SUCCESS);
275a77d671bSBarry Smith 
276a77d671bSBarry Smith   PetscTryTypeMethod(coarser, destroy);
277a77d671bSBarry Smith   coarser->ops->destroy = NULL;
278a77d671bSBarry Smith   PetscCall(PetscMemzero(coarser->ops, sizeof(struct _MatCoarsenOps)));
279a77d671bSBarry Smith 
280a77d671bSBarry Smith   PetscCall(PetscFunctionListFind(MatCoarsenList, type, &r));
281a77d671bSBarry Smith   PetscCheck(r, PetscObjectComm((PetscObject)coarser), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown coarsen type %s", type);
282a77d671bSBarry Smith   PetscCall((*r)(coarser));
283a77d671bSBarry Smith 
284a77d671bSBarry Smith   PetscCall(PetscFree(((PetscObject)coarser)->type_name));
285a77d671bSBarry Smith   PetscCall(PetscStrallocpy(type, &((PetscObject)coarser)->type_name));
286a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
287a77d671bSBarry Smith }
288a77d671bSBarry Smith 
289a77d671bSBarry Smith /*@
290a77d671bSBarry Smith   MatCoarsenSetGreedyOrdering - Sets the ordering of the vertices to use with a greedy coarsening method
291a77d671bSBarry Smith 
292a77d671bSBarry Smith   Logically Collective
293a77d671bSBarry Smith 
294a77d671bSBarry Smith   Input Parameters:
295a77d671bSBarry Smith + coarser - the coarsen context
296a77d671bSBarry Smith - perm    - vertex ordering of (greedy) algorithm
297a77d671bSBarry Smith 
298a77d671bSBarry Smith   Level: advanced
299a77d671bSBarry Smith 
300a77d671bSBarry Smith   Note:
301a77d671bSBarry Smith   The `IS` weights is freed by PETSc, the user should not destroy it or change it after this call
302a77d671bSBarry Smith 
303a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenCreate()`, `MatCoarsenSetType()`
304a77d671bSBarry Smith @*/
MatCoarsenSetGreedyOrdering(MatCoarsen coarser,const IS perm)305a77d671bSBarry Smith PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen coarser, const IS perm)
306a77d671bSBarry Smith {
307a77d671bSBarry Smith   PetscFunctionBegin;
308a77d671bSBarry Smith   PetscValidHeaderSpecific(coarser, MAT_COARSEN_CLASSID, 1);
309a77d671bSBarry Smith   coarser->perm = perm;
310a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
311a77d671bSBarry Smith }
312a77d671bSBarry Smith 
313a77d671bSBarry Smith /*@C
314a77d671bSBarry Smith   MatCoarsenGetData - Gets the weights for vertices for a coarsener.
315a77d671bSBarry Smith 
316a77d671bSBarry Smith   Logically Collective, No Fortran Support
317a77d671bSBarry Smith 
318a77d671bSBarry Smith   Input Parameter:
319a77d671bSBarry Smith . coarser - the coarsen context
320a77d671bSBarry Smith 
321a77d671bSBarry Smith   Output Parameter:
322a77d671bSBarry Smith . llist - linked list of aggregates
323a77d671bSBarry Smith 
324a77d671bSBarry Smith   Level: advanced
325a77d671bSBarry Smith 
326a77d671bSBarry Smith   Note:
327a77d671bSBarry Smith   This passes ownership to the caller and nullifies the value of weights (`PetscCoarsenData`) within the `MatCoarsen`
328a77d671bSBarry Smith 
329a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenApply()`, `MatCoarsenCreate()`, `MatCoarsenSetType()`, `PetscCoarsenData`
330a77d671bSBarry Smith @*/
MatCoarsenGetData(MatCoarsen coarser,PetscCoarsenData ** llist)331a77d671bSBarry Smith PetscErrorCode MatCoarsenGetData(MatCoarsen coarser, PetscCoarsenData **llist)
332a77d671bSBarry Smith {
333a77d671bSBarry Smith   PetscFunctionBegin;
334a77d671bSBarry Smith   PetscValidHeaderSpecific(coarser, MAT_COARSEN_CLASSID, 1);
335a77d671bSBarry Smith   PetscCheck(coarser->agg_lists, PetscObjectComm((PetscObject)coarser), PETSC_ERR_ARG_WRONGSTATE, "No linked list - generate it or call ApplyCoarsen");
336a77d671bSBarry Smith   *llist             = coarser->agg_lists;
337a77d671bSBarry Smith   coarser->agg_lists = NULL; /* giving up ownership */
338a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
339a77d671bSBarry Smith }
340a77d671bSBarry Smith 
341a77d671bSBarry Smith /*@
342a77d671bSBarry Smith   MatCoarsenSetFromOptions - Sets various coarsen options from the options database.
343a77d671bSBarry Smith 
344a77d671bSBarry Smith   Collective
345a77d671bSBarry Smith 
346a77d671bSBarry Smith   Input Parameter:
347a77d671bSBarry Smith . coarser - the coarsen context.
348a77d671bSBarry Smith 
349a77d671bSBarry Smith   Options Database Key:
350a77d671bSBarry Smith + -mat_coarsen_type  <type>                                                       - mis: maximal independent set based; misk: distance k MIS; hem: heavy edge matching
351a77d671bSBarry Smith - -mat_coarsen_max_it <its> number of iterations to use in the coarsening process - see `MatCoarsenSetMaximumIterations()`
352a77d671bSBarry Smith 
353a77d671bSBarry Smith   Level: advanced
354a77d671bSBarry Smith 
355a77d671bSBarry Smith   Notes:
356a77d671bSBarry Smith   When the coarsening is used inside `PCGAMG` then the options database keys are prefixed with `-pc_gamg_`
357a77d671bSBarry Smith 
358a77d671bSBarry Smith   Sets the `MatCoarsenType` to `MATCOARSENMISK` if has not been set previously
359a77d671bSBarry Smith 
360a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenApply()`, `MatCoarsenCreate()`, `MatCoarsenSetType()`,
361a77d671bSBarry Smith           `MatCoarsenSetMaximumIterations()`
362a77d671bSBarry Smith @*/
MatCoarsenSetFromOptions(MatCoarsen coarser)363a77d671bSBarry Smith PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen coarser)
364a77d671bSBarry Smith {
365a77d671bSBarry Smith   PetscBool   flag;
366a77d671bSBarry Smith   char        type[256];
367a77d671bSBarry Smith   const char *def;
368a77d671bSBarry Smith 
369a77d671bSBarry Smith   PetscFunctionBegin;
370a77d671bSBarry Smith   PetscObjectOptionsBegin((PetscObject)coarser);
371a77d671bSBarry Smith   if (!((PetscObject)coarser)->type_name) {
372a77d671bSBarry Smith     def = MATCOARSENMISK;
373a77d671bSBarry Smith   } else {
374a77d671bSBarry Smith     def = ((PetscObject)coarser)->type_name;
375a77d671bSBarry Smith   }
376a77d671bSBarry Smith   PetscCall(PetscOptionsFList("-mat_coarsen_type", "Type of aggregator", "MatCoarsenSetType", MatCoarsenList, def, type, 256, &flag));
377a77d671bSBarry Smith   if (flag) PetscCall(MatCoarsenSetType(coarser, type));
378a77d671bSBarry Smith 
379a77d671bSBarry Smith   PetscCall(PetscOptionsInt("-mat_coarsen_max_it", "Number of iterations (for HEM)", "MatCoarsenSetMaximumIterations", coarser->max_it, &coarser->max_it, NULL));
380a77d671bSBarry Smith   PetscCall(PetscOptionsInt("-mat_coarsen_threshold", "Threshold (for HEM)", "MatCoarsenSetThreshold", coarser->max_it, &coarser->max_it, NULL));
381a77d671bSBarry Smith   coarser->strength_index_size = MAT_COARSEN_STRENGTH_INDEX_SIZE;
382a77d671bSBarry Smith   PetscCall(PetscOptionsIntArray("-mat_coarsen_strength_index", "Array of indices to use strength of connection measure (default is all indices)", "MatCoarsenSetStrengthIndex", coarser->strength_index, &coarser->strength_index_size, NULL));
383a77d671bSBarry Smith   /*
384a77d671bSBarry Smith    Set the type if it was never set.
385a77d671bSBarry Smith    */
386a77d671bSBarry Smith   if (!((PetscObject)coarser)->type_name) PetscCall(MatCoarsenSetType(coarser, def));
387a77d671bSBarry Smith 
388a77d671bSBarry Smith   PetscTryTypeMethod(coarser, setfromoptions, PetscOptionsObject);
389a77d671bSBarry Smith   PetscOptionsEnd();
390a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
391a77d671bSBarry Smith }
392a77d671bSBarry Smith 
393a77d671bSBarry Smith /*@
394a77d671bSBarry Smith   MatCoarsenSetMaximumIterations - Maximum `MATCOARSENHEM` iterations to use
395a77d671bSBarry Smith 
396a77d671bSBarry Smith   Logically Collective
397a77d671bSBarry Smith 
398a77d671bSBarry Smith   Input Parameters:
399a77d671bSBarry Smith + coarse - the coarsen context
400a77d671bSBarry Smith - n      - number of HEM iterations
401a77d671bSBarry Smith 
402a77d671bSBarry Smith   Options Database Key:
403a77d671bSBarry Smith . -mat_coarsen_max_it <default=4> - Maximum `MATCOARSENHEM` iterations to use
404a77d671bSBarry Smith 
405a77d671bSBarry Smith   Level: intermediate
406a77d671bSBarry Smith 
407a77d671bSBarry Smith   Note:
408a77d671bSBarry Smith   When the coarsening is used inside `PCGAMG` then the options database keys are prefixed with `-pc_gamg_`
409a77d671bSBarry Smith 
410a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenApply()`, `MatCoarsenCreate()`, `MatCoarsenSetType()`
411a77d671bSBarry Smith @*/
MatCoarsenSetMaximumIterations(MatCoarsen coarse,PetscInt n)412a77d671bSBarry Smith PetscErrorCode MatCoarsenSetMaximumIterations(MatCoarsen coarse, PetscInt n)
413a77d671bSBarry Smith {
414a77d671bSBarry Smith   PetscFunctionBegin;
415a77d671bSBarry Smith   PetscValidHeaderSpecific(coarse, MAT_COARSEN_CLASSID, 1);
416a77d671bSBarry Smith   PetscValidLogicalCollectiveInt(coarse, n, 2);
417a77d671bSBarry Smith   PetscTryMethod(coarse, "MatCoarsenSetMaximumIterations_C", (MatCoarsen, PetscInt), (coarse, n));
418a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
419a77d671bSBarry Smith }
420a77d671bSBarry Smith 
MatCoarsenSetMaximumIterations_MATCOARSEN(MatCoarsen coarse,PetscInt b)421a77d671bSBarry Smith static PetscErrorCode MatCoarsenSetMaximumIterations_MATCOARSEN(MatCoarsen coarse, PetscInt b)
422a77d671bSBarry Smith {
423a77d671bSBarry Smith   PetscFunctionBegin;
424a77d671bSBarry Smith   coarse->max_it = b;
425a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
426a77d671bSBarry Smith }
427a77d671bSBarry Smith 
428a77d671bSBarry Smith /*@
429a77d671bSBarry Smith   MatCoarsenSetStrengthIndex -  Index array to use for index to use for strength of connection
430a77d671bSBarry Smith 
431a77d671bSBarry Smith   Logically Collective
432a77d671bSBarry Smith 
433a77d671bSBarry Smith   Input Parameters:
434a77d671bSBarry Smith + coarse - the coarsen context
435a77d671bSBarry Smith . n      - number of indices
436a77d671bSBarry Smith - idx    - array of indices
437a77d671bSBarry Smith 
438a77d671bSBarry Smith   Options Database Key:
439a77d671bSBarry Smith . -mat_coarsen_strength_index - array of subset of variables per vertex to use for strength norm, -1 for using all (default)
440a77d671bSBarry Smith 
441a77d671bSBarry Smith   Level: intermediate
442a77d671bSBarry Smith 
443a77d671bSBarry Smith   Note:
444a77d671bSBarry Smith   When the coarsening is used inside `PCGAMG` then the options database keys are prefixed with `-pc_gamg_`
445a77d671bSBarry Smith 
446a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenApply()`, `MatCoarsenCreate()`, `MatCoarsenSetType()`
447a77d671bSBarry Smith @*/
MatCoarsenSetStrengthIndex(MatCoarsen coarse,PetscInt n,PetscInt idx[])448a77d671bSBarry Smith PetscErrorCode MatCoarsenSetStrengthIndex(MatCoarsen coarse, PetscInt n, PetscInt idx[])
449a77d671bSBarry Smith {
450a77d671bSBarry Smith   PetscFunctionBegin;
451a77d671bSBarry Smith   PetscValidHeaderSpecific(coarse, MAT_COARSEN_CLASSID, 1);
452a77d671bSBarry Smith   PetscValidLogicalCollectiveInt(coarse, n, 2);
453a77d671bSBarry Smith   PetscTryMethod(coarse, "MatCoarsenSetStrengthIndex_C", (MatCoarsen, PetscInt, PetscInt[]), (coarse, n, idx));
454a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
455a77d671bSBarry Smith }
456a77d671bSBarry Smith 
MatCoarsenSetStrengthIndex_MATCOARSEN(MatCoarsen coarse,PetscInt n,PetscInt idx[])457a77d671bSBarry Smith static PetscErrorCode MatCoarsenSetStrengthIndex_MATCOARSEN(MatCoarsen coarse, PetscInt n, PetscInt idx[])
458a77d671bSBarry Smith {
459a77d671bSBarry Smith   PetscFunctionBegin;
460a77d671bSBarry Smith   coarse->strength_index_size = n;
461a77d671bSBarry Smith   for (int iii = 0; iii < n; iii++) coarse->strength_index[iii] = idx[iii];
462a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
463a77d671bSBarry Smith }
464a77d671bSBarry Smith 
465a77d671bSBarry Smith /*@
466a77d671bSBarry Smith   MatCoarsenSetThreshold - Set the threshold for HEM
467a77d671bSBarry Smith 
468a77d671bSBarry Smith   Logically Collective
469a77d671bSBarry Smith 
470a77d671bSBarry Smith   Input Parameters:
471a77d671bSBarry Smith + coarse - the coarsen context
472a77d671bSBarry Smith - b      - threshold value
473a77d671bSBarry Smith 
474a77d671bSBarry Smith   Options Database Key:
475a77d671bSBarry Smith . -mat_coarsen_threshold <-1> - threshold
476a77d671bSBarry Smith 
477a77d671bSBarry Smith   Level: intermediate
478a77d671bSBarry Smith 
479a77d671bSBarry Smith   Note:
480a77d671bSBarry Smith   When the coarsening is used inside `PCGAMG` then the options database keys are prefixed with `-pc_gamg_`
481a77d671bSBarry Smith 
482a77d671bSBarry Smith   Developer Note:
483a77d671bSBarry Smith   It is not documented how this threshold is used
484a77d671bSBarry Smith 
485a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenApply()`, `MatCoarsenCreate()`, `MatCoarsenSetType()`
486a77d671bSBarry Smith @*/
MatCoarsenSetThreshold(MatCoarsen coarse,PetscReal b)487a77d671bSBarry Smith PetscErrorCode MatCoarsenSetThreshold(MatCoarsen coarse, PetscReal b)
488a77d671bSBarry Smith {
489a77d671bSBarry Smith   PetscFunctionBegin;
490a77d671bSBarry Smith   PetscValidHeaderSpecific(coarse, MAT_COARSEN_CLASSID, 1);
491a77d671bSBarry Smith   PetscValidLogicalCollectiveReal(coarse, b, 2);
492a77d671bSBarry Smith   PetscTryMethod(coarse, "MatCoarsenSetThreshold_C", (MatCoarsen, PetscReal), (coarse, b));
493a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
494a77d671bSBarry Smith }
495a77d671bSBarry Smith 
MatCoarsenSetThreshold_MATCOARSEN(MatCoarsen coarse,PetscReal b)496a77d671bSBarry Smith static PetscErrorCode MatCoarsenSetThreshold_MATCOARSEN(MatCoarsen coarse, PetscReal b)
497a77d671bSBarry Smith {
498a77d671bSBarry Smith   PetscFunctionBegin;
499a77d671bSBarry Smith   coarse->threshold = b;
500a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
501a77d671bSBarry Smith }
502a77d671bSBarry Smith 
503a77d671bSBarry Smith /*@
504a77d671bSBarry Smith   MatCoarsenCreate - Creates a coarsen context.
505a77d671bSBarry Smith 
506a77d671bSBarry Smith   Collective
507a77d671bSBarry Smith 
508a77d671bSBarry Smith   Input Parameter:
509a77d671bSBarry Smith . comm - MPI communicator
510a77d671bSBarry Smith 
511a77d671bSBarry Smith   Output Parameter:
512a77d671bSBarry Smith . newcrs - location to put the context
513a77d671bSBarry Smith 
514a77d671bSBarry Smith   Level: advanced
515a77d671bSBarry Smith 
516a77d671bSBarry Smith .seealso: `MatCoarsen`, `MatCoarsenSetType()`, `MatCoarsenApply()`, `MatCoarsenDestroy()`,
517a77d671bSBarry Smith           `MatCoarsenSetAdjacency()`, `MatCoarsenGetData()`
518a77d671bSBarry Smith @*/
MatCoarsenCreate(MPI_Comm comm,MatCoarsen * newcrs)519a77d671bSBarry Smith PetscErrorCode MatCoarsenCreate(MPI_Comm comm, MatCoarsen *newcrs)
520a77d671bSBarry Smith {
521a77d671bSBarry Smith   MatCoarsen agg;
522a77d671bSBarry Smith 
523a77d671bSBarry Smith   PetscFunctionBegin;
524377f809aSBarry Smith   PetscAssertPointer(newcrs, 2);
525a77d671bSBarry Smith   PetscCall(MatInitializePackage());
526377f809aSBarry Smith 
527a77d671bSBarry Smith   PetscCall(PetscHeaderCreate(agg, MAT_COARSEN_CLASSID, "MatCoarsen", "Matrix/graph coarsen", "MatCoarsen", comm, MatCoarsenDestroy, MatCoarsenView));
528a77d671bSBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)agg, "MatCoarsenSetMaximumIterations_C", MatCoarsenSetMaximumIterations_MATCOARSEN));
529a77d671bSBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)agg, "MatCoarsenSetThreshold_C", MatCoarsenSetThreshold_MATCOARSEN));
530a77d671bSBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)agg, "MatCoarsenSetStrengthIndex_C", MatCoarsenSetStrengthIndex_MATCOARSEN));
531a77d671bSBarry Smith   agg->strength_index_size = 0;
532a77d671bSBarry Smith   *newcrs                  = agg;
533a77d671bSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
534a77d671bSBarry Smith }
535