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
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: /*@C
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: Use `MatCoarsenGetData()` to access the results of the coarsening
79: The user can define additional coarsens; see `MatCoarsenRegister()`.
81: .seealso: `MatCoarsen`, `MatCoarsenSetFromOptions()`, `MatCoarsenSetType()`, `MatCoarsenRegister()`, `MatCoarsenCreate()`,
82: `MatCoarsenDestroy()`, `MatCoarsenSetAdjacency()`
83: `MatCoarsenGetData()`
84: @*/
85: PetscErrorCode MatCoarsenApply(MatCoarsen coarser)
86: {
87: PetscFunctionBegin;
89: PetscAssertPointer(coarser, 1);
90: PetscCheck(coarser->graph->assembled, PetscObjectComm((PetscObject)coarser), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
91: PetscCheck(!coarser->graph->factortype, PetscObjectComm((PetscObject)coarser), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
92: PetscCall(PetscLogEventBegin(MAT_Coarsen, coarser, 0, 0, 0));
93: PetscUseTypeMethod(coarser, apply);
94: PetscCall(PetscLogEventEnd(MAT_Coarsen, coarser, 0, 0, 0));
95: PetscFunctionReturn(PETSC_SUCCESS);
96: }
98: /*@
99: MatCoarsenSetAdjacency - Sets the adjacency graph (matrix) of the thing to be coarsened.
101: Collective
103: Input Parameters:
104: + agg - the coarsen context
105: - adj - the adjacency matrix
107: Level: advanced
109: .seealso: `MatCoarsen`, `MatCoarsenSetFromOptions()`, `Mat`, `MatCoarsenCreate()`, `MatCoarsenApply()`
110: @*/
111: PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen agg, Mat adj)
112: {
113: PetscFunctionBegin;
116: agg->graph = adj;
117: PetscFunctionReturn(PETSC_SUCCESS);
118: }
120: /*@
121: MatCoarsenSetStrictAggs - Set whether to keep strict (non overlapping) aggregates in the linked list of aggregates for a coarsen context
123: Logically Collective
125: Input Parameters:
126: + agg - the coarsen context
127: - str - `PETSC_TRUE` keep strict aggregates, `PETSC_FALSE` allow overlap
129: Level: advanced
131: .seealso: `MatCoarsen`, `MatCoarsenCreate()`, `MatCoarsenSetFromOptions()`
132: @*/
133: PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen agg, PetscBool str)
134: {
135: PetscFunctionBegin;
137: agg->strict_aggs = str;
138: PetscFunctionReturn(PETSC_SUCCESS);
139: }
141: /*@
142: MatCoarsenDestroy - Destroys the coarsen context.
144: Collective
146: Input Parameter:
147: . agg - the coarsen context
149: Level: advanced
151: .seealso: `MatCoarsen`, `MatCoarsenCreate()`
152: @*/
153: PetscErrorCode MatCoarsenDestroy(MatCoarsen *agg)
154: {
155: PetscFunctionBegin;
156: if (!*agg) PetscFunctionReturn(PETSC_SUCCESS);
158: if (--((PetscObject)*agg)->refct > 0) {
159: *agg = NULL;
160: PetscFunctionReturn(PETSC_SUCCESS);
161: }
163: PetscTryTypeMethod(*agg, destroy);
164: if ((*agg)->agg_lists) PetscCall(PetscCDDestroy((*agg)->agg_lists));
165: PetscCall(PetscObjectComposeFunction((PetscObject)*agg, "MatCoarsenSetMaximumIterations_C", NULL));
166: PetscCall(PetscObjectComposeFunction((PetscObject)*agg, "MatCoarsenSetThreshold_C", NULL));
167: PetscCall(PetscObjectComposeFunction((PetscObject)*agg, "MatCoarsenSetStrengthIndex_C", NULL));
169: PetscCall(PetscHeaderDestroy(agg));
170: PetscFunctionReturn(PETSC_SUCCESS);
171: }
173: /*@C
174: MatCoarsenViewFromOptions - View the coarsener from the options database
176: Collective
178: Input Parameters:
179: + A - the coarsen context
180: . obj - Optional object that provides the prefix for the option name
181: - name - command line option (usually `-mat_coarsen_view`)
183: Options Database Key:
184: . -mat_coarsen_view [viewertype]:... - the viewer and its options
186: Note:
187: .vb
188: If no value is provided ascii:stdout is used
189: ascii[:[filename][:[format][:append]]] defaults to stdout - format can be one of ascii_info, ascii_info_detail, or ascii_matlab,
190: for example ascii::ascii_info prints just the information about the object not all details
191: unless :append is given filename opens in write mode, overwriting what was already there
192: binary[:[filename][:[format][:append]]] defaults to the file binaryoutput
193: draw[:drawtype[:filename]] for example, draw:tikz, draw:tikz:figure.tex or draw:x
194: socket[:port] defaults to the standard output port
195: saws[:communicatorname] publishes object to the Scientific Application Webserver (SAWs)
196: .ve
198: Level: intermediate
200: .seealso: `MatCoarsen`, `MatCoarsenView`, `PetscObjectViewFromOptions()`, `MatCoarsenCreate()`
201: @*/
202: PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A, PetscObject obj, const char name[])
203: {
204: PetscFunctionBegin;
206: PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
207: PetscFunctionReturn(PETSC_SUCCESS);
208: }
210: /*@C
211: MatCoarsenView - Prints the coarsen data structure.
213: Collective
215: Input Parameters:
216: + agg - the coarsen context
217: - viewer - optional visualization context
219: For viewing the options database see `MatCoarsenViewFromOptions()`
221: Level: advanced
223: .seealso: `MatCoarsen`, `PetscViewer`, `PetscViewerASCIIOpen()`, `MatCoarsenViewFromOptions`
224: @*/
225: PetscErrorCode MatCoarsenView(MatCoarsen agg, PetscViewer viewer)
226: {
227: PetscBool iascii;
229: PetscFunctionBegin;
231: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)agg), &viewer));
233: PetscCheckSameComm(agg, 1, viewer, 2);
235: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
236: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)agg, viewer));
237: if (agg->ops->view) {
238: PetscCall(PetscViewerASCIIPushTab(viewer));
239: PetscUseTypeMethod(agg, view, viewer);
240: PetscCall(PetscViewerASCIIPopTab(viewer));
241: }
242: 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]));
243: PetscFunctionReturn(PETSC_SUCCESS);
244: }
246: /*@C
247: MatCoarsenSetType - Sets the type of aggregator to use
249: Collective
251: Input Parameters:
252: + coarser - the coarsen context.
253: - type - a known coarsening method
255: Options Database Key:
256: . -mat_coarsen_type <type> - maximal independent set based; distance k MIS; heavy edge matching
258: Level: advanced
260: .seealso: `MatCoarsen`, `MatCoarsenCreate()`, `MatCoarsenApply()`, `MatCoarsenType`, `MatCoarsenGetType()`
261: @*/
262: PetscErrorCode MatCoarsenSetType(MatCoarsen coarser, MatCoarsenType type)
263: {
264: PetscBool match;
265: PetscErrorCode (*r)(MatCoarsen);
267: PetscFunctionBegin;
269: PetscAssertPointer(type, 2);
271: PetscCall(PetscObjectTypeCompare((PetscObject)coarser, type, &match));
272: if (match) PetscFunctionReturn(PETSC_SUCCESS);
274: PetscTryTypeMethod(coarser, destroy);
275: coarser->ops->destroy = NULL;
276: PetscCall(PetscMemzero(coarser->ops, sizeof(struct _MatCoarsenOps)));
278: PetscCall(PetscFunctionListFind(MatCoarsenList, type, &r));
279: PetscCheck(r, PetscObjectComm((PetscObject)coarser), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown coarsen type %s", type);
280: PetscCall((*r)(coarser));
282: PetscCall(PetscFree(((PetscObject)coarser)->type_name));
283: PetscCall(PetscStrallocpy(type, &((PetscObject)coarser)->type_name));
284: PetscFunctionReturn(PETSC_SUCCESS);
285: }
287: /*@C
288: MatCoarsenSetGreedyOrdering - Sets the ordering of the vertices to use with a greedy coarsening method
290: Logically Collective
292: Input Parameters:
293: + coarser - the coarsen context
294: - perm - vertex ordering of (greedy) algorithm
296: Level: advanced
298: Note:
299: The `IS` weights is freed by PETSc, the user should not destroy it or change it after this call
301: .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenCreate()`, `MatCoarsenSetType()`
302: @*/
303: PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen coarser, const IS perm)
304: {
305: PetscFunctionBegin;
307: coarser->perm = perm;
308: PetscFunctionReturn(PETSC_SUCCESS);
309: }
311: /*@C
312: MatCoarsenGetData - Gets the weights for vertices for a coarsener.
314: Logically Collective
316: Input Parameter:
317: . coarser - the coarsen context
319: Output Parameter:
320: . llist - linked list of aggregates
322: Level: advanced
324: Note:
325: This passes ownership to the caller and nullifies the value of weights (`PetscCoarsenData`) within the `MatCoarsen`
327: .seealso: `MatCoarsen`, `MatCoarsenApply()`, `MatCoarsenCreate()`, `MatCoarsenSetType()`, `PetscCoarsenData`
328: @*/
329: PetscErrorCode MatCoarsenGetData(MatCoarsen coarser, PetscCoarsenData **llist)
330: {
331: PetscFunctionBegin;
333: PetscCheck(coarser->agg_lists, PetscObjectComm((PetscObject)coarser), PETSC_ERR_ARG_WRONGSTATE, "No linked list - generate it or call ApplyCoarsen");
334: *llist = coarser->agg_lists;
335: coarser->agg_lists = NULL; /* giving up ownership */
336: PetscFunctionReturn(PETSC_SUCCESS);
337: }
339: /*@
340: MatCoarsenSetFromOptions - Sets various coarsen options from the options database.
342: Collective
344: Input Parameter:
345: . coarser - the coarsen context.
347: Options Database Key:
348: . -mat_coarsen_type <type> - mis: maximal independent set based; misk: distance k MIS; hem: heavy edge matching
350: Level: advanced
352: Note:
353: Sets the `MatCoarsenType` to `MATCOARSENMISK` if has not been set previously
355: .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenApply()`, `MatCoarsenCreate()`, `MatCoarsenSetType()`
356: @*/
357: PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen coarser)
358: {
359: PetscBool flag;
360: char type[256];
361: const char *def;
363: PetscFunctionBegin;
364: PetscObjectOptionsBegin((PetscObject)coarser);
365: if (!((PetscObject)coarser)->type_name) {
366: def = MATCOARSENMISK;
367: } else {
368: def = ((PetscObject)coarser)->type_name;
369: }
370: PetscCall(PetscOptionsFList("-mat_coarsen_type", "Type of aggregator", "MatCoarsenSetType", MatCoarsenList, def, type, 256, &flag));
371: if (flag) PetscCall(MatCoarsenSetType(coarser, type));
373: PetscCall(PetscOptionsInt("-mat_coarsen_max_it", "Number of iterations (for HEM)", "MatCoarsenSetMaximumIterations", coarser->max_it, &coarser->max_it, NULL));
374: PetscCall(PetscOptionsInt("-mat_coarsen_threshold", "Threshold (for HEM)", "MatCoarsenSetThreshold", coarser->max_it, &coarser->max_it, NULL));
375: coarser->strength_index_size = MAT_COARSEN_STRENGTH_INDEX_SIZE;
376: 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));
377: /*
378: Set the type if it was never set.
379: */
380: if (!((PetscObject)coarser)->type_name) PetscCall(MatCoarsenSetType(coarser, def));
382: PetscTryTypeMethod(coarser, setfromoptions, PetscOptionsObject);
383: PetscOptionsEnd();
384: PetscFunctionReturn(PETSC_SUCCESS);
385: }
387: /*@
388: MatCoarsenSetMaximumIterations - Maximum HEM iterations to use
390: Logically Collective
392: Input Parameters:
393: + coarse - the coarsen context
394: - n - number of HEM iterations
396: Options Database Key:
397: . -mat_coarsen_max_it <default=4> - Max HEM iterations
399: Level: intermediate
401: .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenApply()`, `MatCoarsenCreate()`, `MatCoarsenSetType()`
402: @*/
403: PetscErrorCode MatCoarsenSetMaximumIterations(MatCoarsen coarse, PetscInt n)
404: {
405: PetscFunctionBegin;
408: PetscTryMethod(coarse, "MatCoarsenSetMaximumIterations_C", (MatCoarsen, PetscInt), (coarse, n));
409: PetscFunctionReturn(PETSC_SUCCESS);
410: }
412: static PetscErrorCode MatCoarsenSetMaximumIterations_MATCOARSEN(MatCoarsen coarse, PetscInt b)
413: {
414: PetscFunctionBegin;
415: coarse->max_it = b;
416: PetscFunctionReturn(PETSC_SUCCESS);
417: }
419: /*@
420: MatCoarsenSetStrengthIndex - Index array to use for index to use for strength of connection
422: Logically Collective
424: Input Parameters:
425: + coarse - the coarsen context
426: . n - number of indices
427: - idx - array of indices
429: Options Database Key:
430: . -mat_coarsen_strength_index - array of subset of variables per vertex to use for strength norm, -1 for using all (default)
432: Level: intermediate
434: .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenApply()`, `MatCoarsenCreate()`, `MatCoarsenSetType()`
435: @*/
436: PetscErrorCode MatCoarsenSetStrengthIndex(MatCoarsen coarse, PetscInt n, PetscInt idx[])
437: {
438: PetscFunctionBegin;
441: PetscTryMethod(coarse, "MatCoarsenSetStrengthIndex_C", (MatCoarsen, PetscInt, PetscInt[]), (coarse, n, idx));
442: PetscFunctionReturn(PETSC_SUCCESS);
443: }
445: static PetscErrorCode MatCoarsenSetStrengthIndex_MATCOARSEN(MatCoarsen coarse, PetscInt n, PetscInt idx[])
446: {
447: PetscFunctionBegin;
448: coarse->strength_index_size = n;
449: for (int iii = 0; iii < n; iii++) coarse->strength_index[iii] = idx[iii];
450: PetscFunctionReturn(PETSC_SUCCESS);
451: }
453: /*@
454: MatCoarsenSetThreshold - Set the threshold for HEM
456: Logically Collective
458: Input Parameters:
459: + coarse - the coarsen context
460: - b - threshold value
462: Options Database Key:
463: . -mat_coarsen_threshold <-1> - threshold
465: Level: intermediate
467: Note:
468: It is not documented how this threshold is used
470: .seealso: `MatCoarsen`, `MatCoarsenType`, `MatCoarsenApply()`, `MatCoarsenCreate()`, `MatCoarsenSetType()`
471: @*/
472: PetscErrorCode MatCoarsenSetThreshold(MatCoarsen coarse, PetscReal b)
473: {
474: PetscFunctionBegin;
477: PetscTryMethod(coarse, "MatCoarsenSetThreshold_C", (MatCoarsen, PetscReal), (coarse, b));
478: PetscFunctionReturn(PETSC_SUCCESS);
479: }
481: static PetscErrorCode MatCoarsenSetThreshold_MATCOARSEN(MatCoarsen coarse, PetscReal b)
482: {
483: PetscFunctionBegin;
484: coarse->threshold = b;
485: PetscFunctionReturn(PETSC_SUCCESS);
486: }
488: /*@
489: MatCoarsenCreate - Creates a coarsen context.
491: Collective
493: Input Parameter:
494: . comm - MPI communicator
496: Output Parameter:
497: . newcrs - location to put the context
499: Level: advanced
501: .seealso: `MatCoarsen`, `MatCoarsenSetType()`, `MatCoarsenApply()`, `MatCoarsenDestroy()`,
502: `MatCoarsenSetAdjacency()`, `MatCoarsenGetData()`
503: @*/
504: PetscErrorCode MatCoarsenCreate(MPI_Comm comm, MatCoarsen *newcrs)
505: {
506: MatCoarsen agg;
508: PetscFunctionBegin;
509: *newcrs = NULL;
511: PetscCall(MatInitializePackage());
512: PetscCall(PetscHeaderCreate(agg, MAT_COARSEN_CLASSID, "MatCoarsen", "Matrix/graph coarsen", "MatCoarsen", comm, MatCoarsenDestroy, MatCoarsenView));
513: PetscCall(PetscObjectComposeFunction((PetscObject)agg, "MatCoarsenSetMaximumIterations_C", MatCoarsenSetMaximumIterations_MATCOARSEN));
514: PetscCall(PetscObjectComposeFunction((PetscObject)agg, "MatCoarsenSetThreshold_C", MatCoarsenSetThreshold_MATCOARSEN));
515: PetscCall(PetscObjectComposeFunction((PetscObject)agg, "MatCoarsenSetStrengthIndex_C", MatCoarsenSetStrengthIndex_MATCOARSEN));
517: agg->strength_index_size = 0;
519: *newcrs = agg;
520: PetscFunctionReturn(PETSC_SUCCESS);
521: }