Actual source code: lognested.c
1: #include <petscviewer.h>
2: #include "lognested.h"
3: #include "xmlviewer.h"
5: PETSC_INTERN PetscErrorCode PetscLogHandlerNestedSetThreshold(PetscLogHandler h, PetscLogDouble newThresh, PetscLogDouble *oldThresh)
6: {
7: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
9: PetscFunctionBegin;
10: if (oldThresh) *oldThresh = nested->threshold;
11: if (newThresh == (PetscLogDouble)PETSC_DECIDE) newThresh = 0.01;
12: if (newThresh == (PetscLogDouble)PETSC_DEFAULT) newThresh = 0.01;
13: nested->threshold = PetscMax(newThresh, 0.0);
14: PetscFunctionReturn(PETSC_SUCCESS);
15: }
17: static PetscErrorCode PetscLogEventGetNestedEvent(PetscLogHandler h, PetscLogEvent e, PetscLogEvent *nested_event)
18: {
19: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
20: NestedIdPair key;
21: PetscHashIter iter;
22: PetscBool missing;
23: PetscLogState state;
25: PetscFunctionBegin;
26: PetscCall(PetscLogHandlerGetState(h, &state));
27: PetscCall(PetscIntStackTop(nested->nested_stack, &key.root));
28: key.leaf = NestedIdFromEvent(e);
29: PetscCall(PetscNestedHashPut(nested->pair_map, key, &iter, &missing));
30: if (missing) {
31: // register a new nested event
32: char name[BUFSIZ];
33: PetscLogEventInfo event_info;
34: PetscLogEventInfo nested_event_info;
36: PetscCall(PetscLogStateEventGetInfo(state, e, &event_info));
37: PetscCall(PetscLogStateEventGetInfo(nested->state, key.root, &nested_event_info));
38: PetscCall(PetscSNPrintf(name, sizeof(name) - 1, "%s;%s", nested_event_info.name, event_info.name));
39: PetscCall(PetscLogStateEventRegister(nested->state, name, event_info.classid, nested_event));
40: PetscCall(PetscNestedHashIterSet(nested->pair_map, iter, *nested_event));
41: } else {
42: PetscCall(PetscNestedHashIterGet(nested->pair_map, iter, nested_event));
43: }
44: PetscFunctionReturn(PETSC_SUCCESS);
45: }
47: static PetscErrorCode PetscLogStageGetNestedEvent(PetscLogHandler h, PetscLogStage stage, PetscLogEvent *nested_event)
48: {
49: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
50: NestedIdPair key;
51: PetscHashIter iter;
52: PetscBool missing;
53: PetscLogState state;
55: PetscFunctionBegin;
56: PetscCall(PetscLogHandlerGetState(h, &state));
57: PetscCall(PetscIntStackTop(nested->nested_stack, &key.root));
58: key.leaf = NestedIdFromStage(stage);
59: PetscCall(PetscNestedHashPut(nested->pair_map, key, &iter, &missing));
60: if (missing) {
61: PetscLogStageInfo stage_info;
62: char name[BUFSIZ];
64: PetscCall(PetscLogStateStageGetInfo(state, stage, &stage_info));
65: if (key.root >= 0) {
66: PetscLogEventInfo nested_event_info;
68: PetscCall(PetscLogStateEventGetInfo(nested->state, key.root, &nested_event_info));
69: PetscCall(PetscSNPrintf(name, sizeof(name) - 1, "%s;%s", nested_event_info.name, stage_info.name));
70: } else {
71: PetscCall(PetscSNPrintf(name, sizeof(name) - 1, "%s", stage_info.name));
72: }
73: PetscCall(PetscLogStateEventRegister(nested->state, name, nested->nested_stage_id, nested_event));
74: PetscCall(PetscNestedHashIterSet(nested->pair_map, iter, *nested_event));
75: } else {
76: PetscCall(PetscNestedHashIterGet(nested->pair_map, iter, nested_event));
77: }
78: PetscFunctionReturn(PETSC_SUCCESS);
79: }
81: static PetscErrorCode PetscLogNestedFindNestedId(PetscLogHandler h, NestedId orig_id, PetscInt *pop_count)
82: {
83: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
84: PetscInt count, i;
86: PetscFunctionBegin;
87: // stop before zero cause there is a null event at the bottom of the stack
88: for (i = nested->orig_stack->top, count = 0; i > 0; i--) {
89: count++;
90: if (nested->orig_stack->stack[i] == orig_id) break;
91: }
92: *pop_count = count;
93: if (count == 1) PetscFunctionReturn(PETSC_SUCCESS); // Normal function, just the top of the stack is being popped.
94: if (orig_id > 0) {
95: PetscLogEvent event_id = NestedIdToEvent(orig_id);
96: PetscLogState state;
97: PetscLogEventInfo event_info;
99: PetscCall(PetscLogHandlerGetState(h, &state));
100: PetscCall(PetscLogStateEventGetInfo(state, event_id, &event_info));
101: PetscCheck(i > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Tried to end event %s, but it is not in the event stack", event_info.name);
102: } else {
103: PetscLogStage stage_id = NestedIdToStage(orig_id);
104: PetscLogState state;
105: PetscLogStageInfo stage_info;
107: PetscCall(PetscLogHandlerGetState(h, &state));
108: PetscCall(PetscLogStateStageGetInfo(state, stage_id, &stage_info));
109: PetscCheck(i > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Tried to pop stage %s, but it is not in the stage stack", stage_info.name);
110: }
111: PetscFunctionReturn(PETSC_SUCCESS);
112: }
114: static PetscErrorCode PetscLogNestedCheckNested(PetscLogHandler h, NestedId leaf, PetscLogEvent nested_event)
115: {
116: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
117: NestedIdPair key;
118: NestedId val;
120: PetscFunctionBegin;
121: PetscCall(PetscIntStackTop(nested->nested_stack, &key.root));
122: key.leaf = leaf;
123: PetscCall(PetscNestedHashGet(nested->pair_map, key, &val));
124: PetscCheck(val == nested_event, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Logging events and stages are not nested, nested logging cannot be used");
125: PetscFunctionReturn(PETSC_SUCCESS);
126: }
128: static PetscErrorCode PetscLogHandlerEventBegin_Nested(PetscLogHandler h, PetscLogEvent e, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
129: {
130: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
131: PetscLogEvent nested_event;
133: PetscFunctionBegin;
134: PetscCall(PetscLogEventGetNestedEvent(h, e, &nested_event));
135: PetscCall(PetscLogHandlerEventBegin(nested->handler, nested_event, o1, o2, o3, o4));
136: PetscCall(PetscIntStackPush(nested->nested_stack, nested_event));
137: PetscCall(PetscIntStackPush(nested->orig_stack, NestedIdFromEvent(e)));
138: PetscFunctionReturn(PETSC_SUCCESS);
139: }
141: static PetscErrorCode PetscLogHandlerNestedEventEnd(PetscLogHandler h, NestedId id, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
142: {
143: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
144: PetscInt pop_count;
146: PetscFunctionBegin;
147: PetscCall(PetscLogNestedFindNestedId(h, id, &pop_count));
148: for (PetscInt c = 0; c < pop_count; c++) {
149: PetscLogEvent nested_event;
150: PetscLogEvent nested_id;
152: PetscCall(PetscIntStackPop(nested->nested_stack, &nested_event));
153: PetscCall(PetscIntStackPop(nested->orig_stack, &nested_id));
154: if (PetscDefined(USE_DEBUG)) PetscCall(PetscLogNestedCheckNested(h, nested_id, nested_event));
155: if ((pop_count > 1) && (c + 1 < pop_count)) {
156: if (nested_id > 0) {
157: PetscLogEvent event_id = NestedIdToEvent(nested_id);
158: PetscLogState state;
159: PetscLogEventInfo event_info;
161: PetscCall(PetscLogHandlerGetState(h, &state));
162: PetscCall(PetscLogStateEventGetInfo(state, event_id, &event_info));
163: PetscCall(PetscInfo(h, "Log event %s wasn't ended, ending it to maintain stack property for nested log handler\n", event_info.name));
164: }
165: }
166: PetscCall(PetscLogHandlerEventEnd(nested->handler, nested_event, o1, o2, o3, o4));
167: }
168: PetscFunctionReturn(PETSC_SUCCESS);
169: }
171: static PetscErrorCode PetscLogHandlerEventEnd_Nested(PetscLogHandler h, PetscLogEvent e, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
172: {
173: PetscFunctionBegin;
174: PetscCall(PetscLogHandlerNestedEventEnd(h, NestedIdFromEvent(e), o1, o2, o3, o4));
175: PetscFunctionReturn(PETSC_SUCCESS);
176: }
178: static PetscErrorCode PetscLogHandlerEventSync_Nested(PetscLogHandler h, PetscLogEvent e, MPI_Comm comm)
179: {
180: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
181: PetscLogEvent nested_event;
183: PetscFunctionBegin;
184: PetscCall(PetscLogEventGetNestedEvent(h, e, &nested_event));
185: PetscCall(PetscLogHandlerEventSync(nested->handler, nested_event, comm));
186: PetscFunctionReturn(PETSC_SUCCESS);
187: }
189: static PetscErrorCode PetscLogHandlerStagePush_Nested(PetscLogHandler h, PetscLogStage stage)
190: {
191: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
192: PetscLogEvent nested_event;
194: PetscFunctionBegin;
195: if (nested->nested_stage_id == -1) PetscCall(PetscClassIdRegister("LogNestedStage", &nested->nested_stage_id));
196: PetscCall(PetscLogStageGetNestedEvent(h, stage, &nested_event));
197: PetscCall(PetscLogHandlerEventBegin(nested->handler, nested_event, NULL, NULL, NULL, NULL));
198: PetscCall(PetscIntStackPush(nested->nested_stack, nested_event));
199: PetscCall(PetscIntStackPush(nested->orig_stack, NestedIdFromStage(stage)));
200: PetscFunctionReturn(PETSC_SUCCESS);
201: }
203: static PetscErrorCode PetscLogHandlerStagePop_Nested(PetscLogHandler h, PetscLogStage stage)
204: {
205: PetscFunctionBegin;
206: PetscCall(PetscLogHandlerNestedEventEnd(h, NestedIdFromStage(stage), NULL, NULL, NULL, NULL));
207: PetscFunctionReturn(PETSC_SUCCESS);
208: }
210: static PetscErrorCode PetscLogHandlerContextCreate_Nested(MPI_Comm comm, PetscLogHandler_Nested *nested_p)
211: {
212: PetscLogStage root_stage;
213: PetscLogHandler_Nested nested;
215: PetscFunctionBegin;
216: PetscCall(PetscNew(nested_p));
217: nested = *nested_p;
218: PetscCall(PetscLogStateCreate(&nested->state));
219: PetscCall(PetscIntStackCreate(&nested->nested_stack));
220: PetscCall(PetscIntStackCreate(&nested->orig_stack));
221: nested->nested_stage_id = -1;
222: nested->threshold = 0.01;
223: PetscCall(PetscNestedHashCreate(&nested->pair_map));
224: PetscCall(PetscLogHandlerCreate(comm, &nested->handler));
225: PetscCall(PetscLogHandlerSetType(nested->handler, PETSCLOGHANDLERDEFAULT));
226: PetscCall(PetscLogHandlerSetState(nested->handler, nested->state));
227: PetscCall(PetscLogStateStageRegister(nested->state, "", &root_stage));
228: PetscAssert(root_stage == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "root stage not zero");
229: PetscCall(PetscLogHandlerStagePush(nested->handler, root_stage));
230: PetscCall(PetscLogStateStagePush(nested->state, root_stage));
231: PetscCall(PetscIntStackPush(nested->nested_stack, -1));
232: PetscCall(PetscIntStackPush(nested->orig_stack, -1));
233: PetscFunctionReturn(PETSC_SUCCESS);
234: }
236: static PetscErrorCode PetscLogHandlerObjectCreate_Nested(PetscLogHandler h, PetscObject obj)
237: {
238: PetscClassId classid;
239: PetscInt num_registered, num_nested_registered;
240: PetscLogState state;
241: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
243: PetscFunctionBegin;
244: // register missing objects
245: PetscCall(PetscObjectGetClassId(obj, &classid));
246: PetscCall(PetscLogHandlerGetState(h, &state));
247: PetscCall(PetscLogStateGetNumClasses(nested->state, &num_nested_registered));
248: PetscCall(PetscLogStateGetNumClasses(state, &num_registered));
249: for (PetscLogClass c = num_nested_registered; c < num_registered; c++) {
250: PetscLogClassInfo class_info;
251: PetscLogClass nested_c;
253: PetscCall(PetscLogStateClassGetInfo(state, c, &class_info));
254: PetscCall(PetscLogStateClassRegister(nested->state, class_info.name, class_info.classid, &nested_c));
255: }
256: PetscCall(PetscLogHandlerObjectCreate(nested->handler, obj));
257: PetscFunctionReturn(PETSC_SUCCESS);
258: }
260: static PetscErrorCode PetscLogHandlerObjectDestroy_Nested(PetscLogHandler h, PetscObject obj)
261: {
262: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
264: PetscFunctionBegin;
265: PetscCall(PetscLogHandlerObjectDestroy(nested->handler, obj));
266: PetscFunctionReturn(PETSC_SUCCESS);
267: }
269: static PetscErrorCode PetscLogHandlerDestroy_Nested(PetscLogHandler h)
270: {
271: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
273: PetscFunctionBegin;
274: PetscCall(PetscLogStateStagePop(nested->state));
275: PetscCall(PetscLogHandlerStagePop(nested->handler, 0));
276: PetscCall(PetscLogStateDestroy(&nested->state));
277: PetscCall(PetscIntStackDestroy(nested->nested_stack));
278: PetscCall(PetscIntStackDestroy(nested->orig_stack));
279: PetscCall(PetscNestedHashDestroy(&nested->pair_map));
280: PetscCall(PetscLogHandlerDestroy(&nested->handler));
281: PetscCall(PetscFree(nested));
282: PetscFunctionReturn(PETSC_SUCCESS);
283: }
285: static PetscErrorCode PetscLogNestedEventNodesOrderDepthFirst(PetscInt num_nodes, PetscInt parent, PetscNestedEventNode tree[], PetscInt *num_descendants)
286: {
287: PetscInt node, start_loc;
289: PetscFunctionBegin;
290: node = 0;
291: start_loc = 0;
292: while (node < num_nodes) {
293: if (tree[node].parent == parent) {
294: PetscInt num_this_descendants = 0;
295: PetscNestedEventNode tmp = tree[start_loc];
296: tree[start_loc] = tree[node];
297: tree[node] = tmp;
298: PetscCall(PetscLogNestedEventNodesOrderDepthFirst(num_nodes - start_loc - 1, tree[start_loc].id, &tree[start_loc + 1], &num_this_descendants));
299: tree[start_loc].num_descendants = num_this_descendants;
300: *num_descendants += 1 + num_this_descendants;
301: start_loc += 1 + num_this_descendants;
302: node = start_loc;
303: } else {
304: node++;
305: }
306: }
307: PetscFunctionReturn(PETSC_SUCCESS);
308: }
310: static PetscErrorCode PetscLogNestedCreatePerfNodes(MPI_Comm comm, PetscLogHandler_Nested nested, PetscLogGlobalNames global_events, PetscNestedEventNode **tree_p, PetscEventPerfInfo **perf_p)
311: {
312: PetscMPIInt size;
313: PetscInt num_nodes;
314: PetscInt num_map_entries;
315: PetscEventPerfInfo *perf;
316: NestedIdPair *keys;
317: NestedId *vals;
318: PetscInt offset;
319: PetscInt num_descendants;
320: PetscNestedEventNode *tree;
322: PetscFunctionBegin;
323: PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &num_nodes));
324: PetscCall(PetscCalloc1(num_nodes, &tree));
325: for (PetscInt node = 0; node < num_nodes; node++) {
326: tree[node].id = node;
327: PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, node, &tree[node].name));
328: tree[node].parent = -1;
329: }
330: PetscCall(PetscNestedHashGetSize(nested->pair_map, &num_map_entries));
331: PetscCall(PetscMalloc2(num_map_entries, &keys, num_map_entries, &vals));
332: offset = 0;
333: PetscCall(PetscNestedHashGetPairs(nested->pair_map, &offset, keys, vals));
334: for (PetscInt k = 0; k < num_map_entries; k++) {
335: NestedId root_local = keys[k].root;
336: NestedId leaf_local = vals[k];
337: PetscInt root_global;
338: PetscInt leaf_global;
340: PetscCall(PetscLogGlobalNamesLocalGetGlobal(global_events, leaf_local, &leaf_global));
341: if (root_local >= 0) {
342: PetscCall(PetscLogGlobalNamesLocalGetGlobal(global_events, root_local, &root_global));
343: tree[leaf_global].parent = root_global;
344: }
345: }
346: PetscCall(PetscFree2(keys, vals));
347: PetscCallMPI(MPI_Comm_size(comm, &size));
348: if (size > 1) { // get missing parents from other processes
349: PetscInt *parents;
351: PetscCall(PetscMalloc1(num_nodes, &parents));
352: for (PetscInt node = 0; node < num_nodes; node++) parents[node] = tree[node].parent;
353: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, parents, num_nodes, MPIU_INT, MPI_MAX, comm));
354: for (PetscInt node = 0; node < num_nodes; node++) tree[node].parent = parents[node];
355: PetscCall(PetscFree(parents));
356: }
358: num_descendants = 0;
359: PetscCall(PetscLogNestedEventNodesOrderDepthFirst(num_nodes, -1, tree, &num_descendants));
360: PetscAssert(num_descendants == num_nodes, comm, PETSC_ERR_PLIB, "Failed tree ordering invariant");
362: PetscCall(PetscCalloc1(num_nodes, &perf));
363: for (PetscInt tree_node = 0; tree_node < num_nodes; tree_node++) {
364: PetscInt global_id = tree[tree_node].id;
365: PetscInt event_id;
367: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, global_id, &event_id));
368: if (event_id >= 0) {
369: PetscEventPerfInfo *event_info;
371: PetscCall(PetscLogHandlerGetEventPerfInfo(nested->handler, 0, event_id, &event_info));
372: perf[tree_node] = *event_info;
373: } else {
374: PetscCall(PetscArrayzero(&perf[tree_node], 1));
375: }
376: }
377: *tree_p = tree;
378: *perf_p = perf;
379: PetscFunctionReturn(PETSC_SUCCESS);
380: }
382: static PetscErrorCode PetscLogHandlerView_Nested(PetscLogHandler handler, PetscViewer viewer)
383: {
384: PetscLogHandler_Nested nested = (PetscLogHandler_Nested)handler->data;
385: PetscNestedEventNode *nodes;
386: PetscEventPerfInfo *perf;
387: PetscLogGlobalNames global_events;
388: PetscNestedEventTree tree;
389: PetscViewerFormat format;
390: MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
392: PetscFunctionBegin;
393: PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, nested->state->registry, &global_events));
394: PetscCall(PetscLogNestedCreatePerfNodes(comm, nested, global_events, &nodes, &perf));
395: tree.comm = comm;
396: tree.global_events = global_events;
397: tree.perf = perf;
398: tree.nodes = nodes;
399: PetscCall(PetscViewerGetFormat(viewer, &format));
400: if (format == PETSC_VIEWER_ASCII_XML) {
401: PetscCall(PetscLogHandlerView_Nested_XML(nested, &tree, viewer));
402: } else if (format == PETSC_VIEWER_ASCII_FLAMEGRAPH) {
403: PetscCall(PetscLogHandlerView_Nested_Flamegraph(nested, &tree, viewer));
404: } else SETERRQ(comm, PETSC_ERR_ARG_INCOMP, "No nested viewer for this format");
405: PetscCall(PetscLogGlobalNamesDestroy(&global_events));
406: PetscCall(PetscFree(tree.nodes));
407: PetscCall(PetscFree(tree.perf));
408: PetscFunctionReturn(PETSC_SUCCESS);
409: }
411: /*MC
412: PETSCLOGHANDLERNESTED - PETSCLOGHANDLERNESTED = "nested" - A `PetscLogHandler` that collects data for PETSc's
413: XML and flamegraph profiling log viewers. A log handler of this type is created and started by
414: by `PetscLogNestedBegin()`.
416: Level: developer
418: .seealso: [](ch_profiling), `PetscLogHandler`
419: M*/
421: PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Nested(PetscLogHandler handler)
422: {
423: PetscFunctionBegin;
424: PetscCall(PetscLogHandlerContextCreate_Nested(PetscObjectComm((PetscObject)handler), (PetscLogHandler_Nested *)&handler->data));
425: handler->ops->destroy = PetscLogHandlerDestroy_Nested;
426: handler->ops->stagepush = PetscLogHandlerStagePush_Nested;
427: handler->ops->stagepop = PetscLogHandlerStagePop_Nested;
428: handler->ops->eventbegin = PetscLogHandlerEventBegin_Nested;
429: handler->ops->eventend = PetscLogHandlerEventEnd_Nested;
430: handler->ops->eventsync = PetscLogHandlerEventSync_Nested;
431: handler->ops->objectcreate = PetscLogHandlerObjectCreate_Nested;
432: handler->ops->objectdestroy = PetscLogHandlerObjectDestroy_Nested;
433: handler->ops->view = PetscLogHandlerView_Nested;
434: PetscFunctionReturn(PETSC_SUCCESS);
435: }