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