1b9321188SToby Isaac #include <petscviewer.h>
2b9321188SToby Isaac #include "lognested.h"
3b9321188SToby Isaac #include "xmlviewer.h"
4b9321188SToby Isaac
PetscLogHandlerNestedSetThreshold(PetscLogHandler h,PetscLogDouble newThresh,PetscLogDouble * oldThresh)5b9321188SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerNestedSetThreshold(PetscLogHandler h, PetscLogDouble newThresh, PetscLogDouble *oldThresh)
6b9321188SToby Isaac {
7b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
8b9321188SToby Isaac
9b9321188SToby Isaac PetscFunctionBegin;
10b9321188SToby Isaac if (oldThresh) *oldThresh = nested->threshold;
11b9321188SToby Isaac if (newThresh == (PetscLogDouble)PETSC_DECIDE) newThresh = 0.01;
12b9321188SToby Isaac if (newThresh == (PetscLogDouble)PETSC_DEFAULT) newThresh = 0.01;
13b9321188SToby Isaac nested->threshold = PetscMax(newThresh, 0.0);
14b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
15b9321188SToby Isaac }
16b9321188SToby Isaac
PetscLogEventGetNestedEvent(PetscLogHandler h,PetscLogEvent e,PetscLogEvent * nested_event)17b9321188SToby Isaac static PetscErrorCode PetscLogEventGetNestedEvent(PetscLogHandler h, PetscLogEvent e, PetscLogEvent *nested_event)
18b9321188SToby Isaac {
19b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
20b9321188SToby Isaac NestedIdPair key;
21b9321188SToby Isaac PetscHashIter iter;
22b9321188SToby Isaac PetscBool missing;
23b9321188SToby Isaac PetscLogState state;
24b9321188SToby Isaac
25b9321188SToby Isaac PetscFunctionBegin;
26b9321188SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state));
27f4f49eeaSPierre Jolivet PetscCall(PetscIntStackTop(nested->nested_stack, &key.root));
28b9321188SToby Isaac key.leaf = NestedIdFromEvent(e);
29b9321188SToby Isaac PetscCall(PetscNestedHashPut(nested->pair_map, key, &iter, &missing));
30b9321188SToby Isaac if (missing) {
31b9321188SToby Isaac // register a new nested event
32b9321188SToby Isaac char name[BUFSIZ];
33b9321188SToby Isaac PetscLogEventInfo event_info;
34b9321188SToby Isaac PetscLogEventInfo nested_event_info;
35b9321188SToby Isaac
36b9321188SToby Isaac PetscCall(PetscLogStateEventGetInfo(state, e, &event_info));
37b9321188SToby Isaac PetscCall(PetscLogStateEventGetInfo(nested->state, key.root, &nested_event_info));
38b9321188SToby Isaac PetscCall(PetscSNPrintf(name, sizeof(name) - 1, "%s;%s", nested_event_info.name, event_info.name));
39b9321188SToby Isaac PetscCall(PetscLogStateEventRegister(nested->state, name, event_info.classid, nested_event));
409a6bacb0SStefano Zampini PetscCall(PetscLogStateEventSetCollective(nested->state, *nested_event, event_info.collective));
41b9321188SToby Isaac PetscCall(PetscNestedHashIterSet(nested->pair_map, iter, *nested_event));
42b9321188SToby Isaac } else {
43b9321188SToby Isaac PetscCall(PetscNestedHashIterGet(nested->pair_map, iter, nested_event));
44b9321188SToby Isaac }
45b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
46b9321188SToby Isaac }
47b9321188SToby Isaac
PetscLogStageGetNestedEvent(PetscLogHandler h,PetscLogStage stage,PetscLogEvent * nested_event)48b9321188SToby Isaac static PetscErrorCode PetscLogStageGetNestedEvent(PetscLogHandler h, PetscLogStage stage, PetscLogEvent *nested_event)
49b9321188SToby Isaac {
50b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
51b9321188SToby Isaac NestedIdPair key;
52b9321188SToby Isaac PetscHashIter iter;
53b9321188SToby Isaac PetscBool missing;
54b9321188SToby Isaac PetscLogState state;
55b9321188SToby Isaac
56b9321188SToby Isaac PetscFunctionBegin;
57b9321188SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state));
58f4f49eeaSPierre Jolivet PetscCall(PetscIntStackTop(nested->nested_stack, &key.root));
59b9321188SToby Isaac key.leaf = NestedIdFromStage(stage);
60b9321188SToby Isaac PetscCall(PetscNestedHashPut(nested->pair_map, key, &iter, &missing));
61b9321188SToby Isaac if (missing) {
62b9321188SToby Isaac PetscLogStageInfo stage_info;
63b9321188SToby Isaac char name[BUFSIZ];
649a6bacb0SStefano Zampini PetscBool collective = PETSC_TRUE;
65b9321188SToby Isaac
66b9321188SToby Isaac PetscCall(PetscLogStateStageGetInfo(state, stage, &stage_info));
67b9321188SToby Isaac if (key.root >= 0) {
68b9321188SToby Isaac PetscLogEventInfo nested_event_info;
69b9321188SToby Isaac
70b9321188SToby Isaac PetscCall(PetscLogStateEventGetInfo(nested->state, key.root, &nested_event_info));
71b9321188SToby Isaac PetscCall(PetscSNPrintf(name, sizeof(name) - 1, "%s;%s", nested_event_info.name, stage_info.name));
729a6bacb0SStefano Zampini collective = nested_event_info.collective;
73b9321188SToby Isaac } else {
74b9321188SToby Isaac PetscCall(PetscSNPrintf(name, sizeof(name) - 1, "%s", stage_info.name));
75b9321188SToby Isaac }
76b9321188SToby Isaac PetscCall(PetscLogStateEventRegister(nested->state, name, nested->nested_stage_id, nested_event));
779a6bacb0SStefano Zampini PetscCall(PetscLogStateEventSetCollective(nested->state, *nested_event, collective));
78b9321188SToby Isaac PetscCall(PetscNestedHashIterSet(nested->pair_map, iter, *nested_event));
79b9321188SToby Isaac } else {
80b9321188SToby Isaac PetscCall(PetscNestedHashIterGet(nested->pair_map, iter, nested_event));
81b9321188SToby Isaac }
82b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
83b9321188SToby Isaac }
84b9321188SToby Isaac
PetscLogNestedFindNestedId(PetscLogHandler h,NestedId orig_id,PetscInt * pop_count)856b2a052cSToby Isaac static PetscErrorCode PetscLogNestedFindNestedId(PetscLogHandler h, NestedId orig_id, PetscInt *pop_count)
866b2a052cSToby Isaac {
876b2a052cSToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
886b2a052cSToby Isaac PetscInt count, i;
896b2a052cSToby Isaac
906b2a052cSToby Isaac PetscFunctionBegin;
916b2a052cSToby Isaac // stop before zero cause there is a null event at the bottom of the stack
926b2a052cSToby Isaac for (i = nested->orig_stack->top, count = 0; i > 0; i--) {
936b2a052cSToby Isaac count++;
946b2a052cSToby Isaac if (nested->orig_stack->stack[i] == orig_id) break;
956b2a052cSToby Isaac }
966b2a052cSToby Isaac *pop_count = count;
976b2a052cSToby Isaac if (count == 1) PetscFunctionReturn(PETSC_SUCCESS); // Normal function, just the top of the stack is being popped.
986b2a052cSToby Isaac if (orig_id > 0) {
996b2a052cSToby Isaac PetscLogEvent event_id = NestedIdToEvent(orig_id);
1006b2a052cSToby Isaac PetscLogState state;
1016b2a052cSToby Isaac PetscLogEventInfo event_info;
1026b2a052cSToby Isaac
1036b2a052cSToby Isaac PetscCall(PetscLogHandlerGetState(h, &state));
1046b2a052cSToby Isaac PetscCall(PetscLogStateEventGetInfo(state, event_id, &event_info));
10500045ab3SPierre Jolivet 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);
1066b2a052cSToby Isaac } else {
1076b2a052cSToby Isaac PetscLogStage stage_id = NestedIdToStage(orig_id);
1086b2a052cSToby Isaac PetscLogState state;
1096b2a052cSToby Isaac PetscLogStageInfo stage_info;
1106b2a052cSToby Isaac
1116b2a052cSToby Isaac PetscCall(PetscLogHandlerGetState(h, &state));
1126b2a052cSToby Isaac PetscCall(PetscLogStateStageGetInfo(state, stage_id, &stage_info));
11300045ab3SPierre Jolivet 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);
1146b2a052cSToby Isaac }
1156b2a052cSToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
1166b2a052cSToby Isaac }
1176b2a052cSToby Isaac
PetscLogNestedCheckNested(PetscLogHandler h,NestedId leaf,PetscLogEvent nested_event)118b9321188SToby Isaac static PetscErrorCode PetscLogNestedCheckNested(PetscLogHandler h, NestedId leaf, PetscLogEvent nested_event)
119b9321188SToby Isaac {
120b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
121b9321188SToby Isaac NestedIdPair key;
122b9321188SToby Isaac NestedId val;
123b9321188SToby Isaac
124b9321188SToby Isaac PetscFunctionBegin;
125f4f49eeaSPierre Jolivet PetscCall(PetscIntStackTop(nested->nested_stack, &key.root));
126b9321188SToby Isaac key.leaf = leaf;
127b9321188SToby Isaac PetscCall(PetscNestedHashGet(nested->pair_map, key, &val));
128b9321188SToby Isaac PetscCheck(val == nested_event, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Logging events and stages are not nested, nested logging cannot be used");
129b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
130b9321188SToby Isaac }
131b9321188SToby Isaac
PetscLogHandlerEventBegin_Nested(PetscLogHandler h,PetscLogEvent e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)132b9321188SToby Isaac static PetscErrorCode PetscLogHandlerEventBegin_Nested(PetscLogHandler h, PetscLogEvent e, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
133b9321188SToby Isaac {
134b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
135b9321188SToby Isaac PetscLogEvent nested_event;
136b9321188SToby Isaac
137b9321188SToby Isaac PetscFunctionBegin;
138b9321188SToby Isaac PetscCall(PetscLogEventGetNestedEvent(h, e, &nested_event));
139b9321188SToby Isaac PetscCall(PetscLogHandlerEventBegin(nested->handler, nested_event, o1, o2, o3, o4));
1406b2a052cSToby Isaac PetscCall(PetscIntStackPush(nested->nested_stack, nested_event));
1416b2a052cSToby Isaac PetscCall(PetscIntStackPush(nested->orig_stack, NestedIdFromEvent(e)));
1426b2a052cSToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
1436b2a052cSToby Isaac }
1446b2a052cSToby Isaac
PetscLogHandlerNestedEventEnd(PetscLogHandler h,NestedId id,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)1456b2a052cSToby Isaac static PetscErrorCode PetscLogHandlerNestedEventEnd(PetscLogHandler h, NestedId id, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
1466b2a052cSToby Isaac {
1476b2a052cSToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
1486b2a052cSToby Isaac PetscInt pop_count;
1496b2a052cSToby Isaac
1506b2a052cSToby Isaac PetscFunctionBegin;
1516b2a052cSToby Isaac PetscCall(PetscLogNestedFindNestedId(h, id, &pop_count));
1526b2a052cSToby Isaac for (PetscInt c = 0; c < pop_count; c++) {
1536b2a052cSToby Isaac PetscLogEvent nested_event;
1546b2a052cSToby Isaac PetscLogEvent nested_id;
1556b2a052cSToby Isaac
1566b2a052cSToby Isaac PetscCall(PetscIntStackPop(nested->nested_stack, &nested_event));
1576b2a052cSToby Isaac PetscCall(PetscIntStackPop(nested->orig_stack, &nested_id));
1586b2a052cSToby Isaac if (PetscDefined(USE_DEBUG)) PetscCall(PetscLogNestedCheckNested(h, nested_id, nested_event));
1596b2a052cSToby Isaac if ((pop_count > 1) && (c + 1 < pop_count)) {
1606b2a052cSToby Isaac if (nested_id > 0) {
1616b2a052cSToby Isaac PetscLogEvent event_id = NestedIdToEvent(nested_id);
1626b2a052cSToby Isaac PetscLogState state;
1636b2a052cSToby Isaac PetscLogEventInfo event_info;
1646b2a052cSToby Isaac
1656b2a052cSToby Isaac PetscCall(PetscLogHandlerGetState(h, &state));
1666b2a052cSToby Isaac PetscCall(PetscLogStateEventGetInfo(state, event_id, &event_info));
1676b2a052cSToby Isaac PetscCall(PetscInfo(h, "Log event %s wasn't ended, ending it to maintain stack property for nested log handler\n", event_info.name));
1686b2a052cSToby Isaac }
1696b2a052cSToby Isaac }
1706b2a052cSToby Isaac PetscCall(PetscLogHandlerEventEnd(nested->handler, nested_event, o1, o2, o3, o4));
1716b2a052cSToby Isaac }
172b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
173b9321188SToby Isaac }
174b9321188SToby Isaac
PetscLogHandlerEventEnd_Nested(PetscLogHandler h,PetscLogEvent e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)175b9321188SToby Isaac static PetscErrorCode PetscLogHandlerEventEnd_Nested(PetscLogHandler h, PetscLogEvent e, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
176b9321188SToby Isaac {
177b9321188SToby Isaac PetscFunctionBegin;
1786b2a052cSToby Isaac PetscCall(PetscLogHandlerNestedEventEnd(h, NestedIdFromEvent(e), o1, o2, o3, o4));
179b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
180b9321188SToby Isaac }
181b9321188SToby Isaac
PetscLogHandlerEventSync_Nested(PetscLogHandler h,PetscLogEvent e,MPI_Comm comm)182b9321188SToby Isaac static PetscErrorCode PetscLogHandlerEventSync_Nested(PetscLogHandler h, PetscLogEvent e, MPI_Comm comm)
183b9321188SToby Isaac {
184b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
185b9321188SToby Isaac PetscLogEvent nested_event;
186b9321188SToby Isaac
187b9321188SToby Isaac PetscFunctionBegin;
188b9321188SToby Isaac PetscCall(PetscLogEventGetNestedEvent(h, e, &nested_event));
189b9321188SToby Isaac PetscCall(PetscLogHandlerEventSync(nested->handler, nested_event, comm));
190b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
191b9321188SToby Isaac }
192b9321188SToby Isaac
PetscLogHandlerStagePush_Nested(PetscLogHandler h,PetscLogStage stage)193b9321188SToby Isaac static PetscErrorCode PetscLogHandlerStagePush_Nested(PetscLogHandler h, PetscLogStage stage)
194b9321188SToby Isaac {
195b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
196b9321188SToby Isaac PetscLogEvent nested_event;
197b9321188SToby Isaac
198b9321188SToby Isaac PetscFunctionBegin;
199b9321188SToby Isaac if (nested->nested_stage_id == -1) PetscCall(PetscClassIdRegister("LogNestedStage", &nested->nested_stage_id));
200b9321188SToby Isaac PetscCall(PetscLogStageGetNestedEvent(h, stage, &nested_event));
201b9321188SToby Isaac PetscCall(PetscLogHandlerEventBegin(nested->handler, nested_event, NULL, NULL, NULL, NULL));
2026b2a052cSToby Isaac PetscCall(PetscIntStackPush(nested->nested_stack, nested_event));
2036b2a052cSToby Isaac PetscCall(PetscIntStackPush(nested->orig_stack, NestedIdFromStage(stage)));
204b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
205b9321188SToby Isaac }
206b9321188SToby Isaac
PetscLogHandlerStagePop_Nested(PetscLogHandler h,PetscLogStage stage)207b9321188SToby Isaac static PetscErrorCode PetscLogHandlerStagePop_Nested(PetscLogHandler h, PetscLogStage stage)
208b9321188SToby Isaac {
209b9321188SToby Isaac PetscFunctionBegin;
2106b2a052cSToby Isaac PetscCall(PetscLogHandlerNestedEventEnd(h, NestedIdFromStage(stage), NULL, NULL, NULL, NULL));
211b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
212b9321188SToby Isaac }
213b9321188SToby Isaac
PetscLogHandlerContextCreate_Nested(MPI_Comm comm,PetscLogHandler_Nested * nested_p)214b9321188SToby Isaac static PetscErrorCode PetscLogHandlerContextCreate_Nested(MPI_Comm comm, PetscLogHandler_Nested *nested_p)
215b9321188SToby Isaac {
216b9321188SToby Isaac PetscLogStage root_stage;
217b9321188SToby Isaac PetscLogHandler_Nested nested;
218b9321188SToby Isaac
219b9321188SToby Isaac PetscFunctionBegin;
220b9321188SToby Isaac PetscCall(PetscNew(nested_p));
221b9321188SToby Isaac nested = *nested_p;
222b9321188SToby Isaac PetscCall(PetscLogStateCreate(&nested->state));
2236b2a052cSToby Isaac PetscCall(PetscIntStackCreate(&nested->nested_stack));
2246b2a052cSToby Isaac PetscCall(PetscIntStackCreate(&nested->orig_stack));
225b9321188SToby Isaac nested->nested_stage_id = -1;
226b9321188SToby Isaac nested->threshold = 0.01;
227b9321188SToby Isaac PetscCall(PetscNestedHashCreate(&nested->pair_map));
228b9321188SToby Isaac PetscCall(PetscLogHandlerCreate(comm, &nested->handler));
229294de794SToby Isaac PetscCall(PetscLogHandlerSetType(nested->handler, PETSCLOGHANDLERDEFAULT));
230b9321188SToby Isaac PetscCall(PetscLogHandlerSetState(nested->handler, nested->state));
231b9321188SToby Isaac PetscCall(PetscLogStateStageRegister(nested->state, "", &root_stage));
232b9321188SToby Isaac PetscAssert(root_stage == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "root stage not zero");
233b9321188SToby Isaac PetscCall(PetscLogHandlerStagePush(nested->handler, root_stage));
234b9321188SToby Isaac PetscCall(PetscLogStateStagePush(nested->state, root_stage));
2356b2a052cSToby Isaac PetscCall(PetscIntStackPush(nested->nested_stack, -1));
2366b2a052cSToby Isaac PetscCall(PetscIntStackPush(nested->orig_stack, -1));
237b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
238b9321188SToby Isaac }
239b9321188SToby Isaac
240*31a765c4SPierre Jolivet PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wconversion")
PetscLogHandlerObjectCreate_Nested(PetscLogHandler h,PetscObject obj)241b9321188SToby Isaac static PetscErrorCode PetscLogHandlerObjectCreate_Nested(PetscLogHandler h, PetscObject obj)
242b9321188SToby Isaac {
243b9321188SToby Isaac PetscClassId classid;
244b9321188SToby Isaac PetscInt num_registered, num_nested_registered;
245b9321188SToby Isaac PetscLogState state;
246b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
247b9321188SToby Isaac
248b9321188SToby Isaac PetscFunctionBegin;
249b9321188SToby Isaac // register missing objects
250b9321188SToby Isaac PetscCall(PetscObjectGetClassId(obj, &classid));
251b9321188SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state));
252b9321188SToby Isaac PetscCall(PetscLogStateGetNumClasses(nested->state, &num_nested_registered));
253b9321188SToby Isaac PetscCall(PetscLogStateGetNumClasses(state, &num_registered));
254b9321188SToby Isaac for (PetscLogClass c = num_nested_registered; c < num_registered; c++) {
255b9321188SToby Isaac PetscLogClassInfo class_info;
256b9321188SToby Isaac PetscLogClass nested_c;
257b9321188SToby Isaac
258b9321188SToby Isaac PetscCall(PetscLogStateClassGetInfo(state, c, &class_info));
259b9321188SToby Isaac PetscCall(PetscLogStateClassRegister(nested->state, class_info.name, class_info.classid, &nested_c));
260b9321188SToby Isaac }
261b9321188SToby Isaac PetscCall(PetscLogHandlerObjectCreate(nested->handler, obj));
262b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
263b9321188SToby Isaac }
264b9321188SToby Isaac
PetscLogHandlerObjectDestroy_Nested(PetscLogHandler h,PetscObject obj)265b9321188SToby Isaac static PetscErrorCode PetscLogHandlerObjectDestroy_Nested(PetscLogHandler h, PetscObject obj)
266b9321188SToby Isaac {
267b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
268b9321188SToby Isaac
269b9321188SToby Isaac PetscFunctionBegin;
270b9321188SToby Isaac PetscCall(PetscLogHandlerObjectDestroy(nested->handler, obj));
271b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
272b9321188SToby Isaac }
273b9321188SToby Isaac
PetscLogHandlerDestroy_Nested(PetscLogHandler h)274b9321188SToby Isaac static PetscErrorCode PetscLogHandlerDestroy_Nested(PetscLogHandler h)
275b9321188SToby Isaac {
276b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
277b9321188SToby Isaac
278b9321188SToby Isaac PetscFunctionBegin;
279b9321188SToby Isaac PetscCall(PetscLogStateStagePop(nested->state));
280b9321188SToby Isaac PetscCall(PetscLogHandlerStagePop(nested->handler, 0));
281b9321188SToby Isaac PetscCall(PetscLogStateDestroy(&nested->state));
2826b2a052cSToby Isaac PetscCall(PetscIntStackDestroy(nested->nested_stack));
2836b2a052cSToby Isaac PetscCall(PetscIntStackDestroy(nested->orig_stack));
284b9321188SToby Isaac PetscCall(PetscNestedHashDestroy(&nested->pair_map));
285b9321188SToby Isaac PetscCall(PetscLogHandlerDestroy(&nested->handler));
286b9321188SToby Isaac PetscCall(PetscFree(nested));
287b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
288b9321188SToby Isaac }
289b9321188SToby Isaac
PetscLogNestedEventNodesOrderDepthFirst(PetscInt num_nodes,PetscInt parent,PetscNestedEventNode tree[],PetscInt * num_descendants)290b9321188SToby Isaac static PetscErrorCode PetscLogNestedEventNodesOrderDepthFirst(PetscInt num_nodes, PetscInt parent, PetscNestedEventNode tree[], PetscInt *num_descendants)
291b9321188SToby Isaac {
292b9321188SToby Isaac PetscInt node, start_loc;
293b9321188SToby Isaac
2944d86920dSPierre Jolivet PetscFunctionBegin;
295b9321188SToby Isaac node = 0;
296b9321188SToby Isaac start_loc = 0;
297b9321188SToby Isaac while (node < num_nodes) {
298b9321188SToby Isaac if (tree[node].parent == parent) {
299b9321188SToby Isaac PetscInt num_this_descendants = 0;
300b9321188SToby Isaac PetscNestedEventNode tmp = tree[start_loc];
301b9321188SToby Isaac tree[start_loc] = tree[node];
302b9321188SToby Isaac tree[node] = tmp;
303b9321188SToby Isaac PetscCall(PetscLogNestedEventNodesOrderDepthFirst(num_nodes - start_loc - 1, tree[start_loc].id, &tree[start_loc + 1], &num_this_descendants));
304b9321188SToby Isaac tree[start_loc].num_descendants = num_this_descendants;
305b9321188SToby Isaac *num_descendants += 1 + num_this_descendants;
306b9321188SToby Isaac start_loc += 1 + num_this_descendants;
307b9321188SToby Isaac node = start_loc;
308b9321188SToby Isaac } else {
309b9321188SToby Isaac node++;
310b9321188SToby Isaac }
311b9321188SToby Isaac }
312b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
313b9321188SToby Isaac }
314b9321188SToby Isaac
PetscLogNestedCreatePerfNodes(MPI_Comm comm,PetscLogHandler_Nested nested,PetscLogGlobalNames global_events,PetscNestedEventNode ** tree_p,PetscEventPerfInfo ** perf_p)315b9321188SToby Isaac static PetscErrorCode PetscLogNestedCreatePerfNodes(MPI_Comm comm, PetscLogHandler_Nested nested, PetscLogGlobalNames global_events, PetscNestedEventNode **tree_p, PetscEventPerfInfo **perf_p)
316b9321188SToby Isaac {
317b9321188SToby Isaac PetscMPIInt size;
318b9321188SToby Isaac PetscInt num_nodes;
319b9321188SToby Isaac PetscInt num_map_entries;
320b9321188SToby Isaac PetscEventPerfInfo *perf;
321b9321188SToby Isaac NestedIdPair *keys;
322b9321188SToby Isaac NestedId *vals;
323b9321188SToby Isaac PetscInt offset;
324b9321188SToby Isaac PetscInt num_descendants;
325b9321188SToby Isaac PetscNestedEventNode *tree;
326b9321188SToby Isaac
327b9321188SToby Isaac PetscFunctionBegin;
328b9321188SToby Isaac PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &num_nodes));
329b9321188SToby Isaac PetscCall(PetscCalloc1(num_nodes, &tree));
330b9321188SToby Isaac for (PetscInt node = 0; node < num_nodes; node++) {
331b9321188SToby Isaac tree[node].id = node;
332b9321188SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, node, &tree[node].name));
333b9321188SToby Isaac tree[node].parent = -1;
334b9321188SToby Isaac }
335b9321188SToby Isaac PetscCall(PetscNestedHashGetSize(nested->pair_map, &num_map_entries));
336b9321188SToby Isaac PetscCall(PetscMalloc2(num_map_entries, &keys, num_map_entries, &vals));
337b9321188SToby Isaac offset = 0;
338b9321188SToby Isaac PetscCall(PetscNestedHashGetPairs(nested->pair_map, &offset, keys, vals));
339b9321188SToby Isaac for (PetscInt k = 0; k < num_map_entries; k++) {
340b9321188SToby Isaac NestedId root_local = keys[k].root;
341b9321188SToby Isaac NestedId leaf_local = vals[k];
342b9321188SToby Isaac PetscInt root_global;
343b9321188SToby Isaac PetscInt leaf_global;
344b9321188SToby Isaac
345b9321188SToby Isaac PetscCall(PetscLogGlobalNamesLocalGetGlobal(global_events, leaf_local, &leaf_global));
346b9321188SToby Isaac if (root_local >= 0) {
347b9321188SToby Isaac PetscCall(PetscLogGlobalNamesLocalGetGlobal(global_events, root_local, &root_global));
348b9321188SToby Isaac tree[leaf_global].parent = root_global;
349b9321188SToby Isaac }
350b9321188SToby Isaac }
351b9321188SToby Isaac PetscCall(PetscFree2(keys, vals));
352b9321188SToby Isaac PetscCallMPI(MPI_Comm_size(comm, &size));
353b9321188SToby Isaac if (size > 1) { // get missing parents from other processes
354b9321188SToby Isaac PetscInt *parents;
355b9321188SToby Isaac
356b9321188SToby Isaac PetscCall(PetscMalloc1(num_nodes, &parents));
357b9321188SToby Isaac for (PetscInt node = 0; node < num_nodes; node++) parents[node] = tree[node].parent;
358462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, parents, num_nodes, MPIU_INT, MPI_MAX, comm));
359b9321188SToby Isaac for (PetscInt node = 0; node < num_nodes; node++) tree[node].parent = parents[node];
360b9321188SToby Isaac PetscCall(PetscFree(parents));
361b9321188SToby Isaac }
362b9321188SToby Isaac
363b9321188SToby Isaac num_descendants = 0;
364b9321188SToby Isaac PetscCall(PetscLogNestedEventNodesOrderDepthFirst(num_nodes, -1, tree, &num_descendants));
365b9321188SToby Isaac PetscAssert(num_descendants == num_nodes, comm, PETSC_ERR_PLIB, "Failed tree ordering invariant");
366b9321188SToby Isaac
367b9321188SToby Isaac PetscCall(PetscCalloc1(num_nodes, &perf));
368c244758eSToby Isaac for (PetscInt tree_node = 0; tree_node < num_nodes; tree_node++) {
369c244758eSToby Isaac PetscInt global_id = tree[tree_node].id;
370b9321188SToby Isaac PetscInt event_id;
371b9321188SToby Isaac
372c244758eSToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, global_id, &event_id));
373b9321188SToby Isaac if (event_id >= 0) {
374b9321188SToby Isaac PetscEventPerfInfo *event_info;
375b9321188SToby Isaac
376dff009beSToby Isaac PetscCall(PetscLogHandlerGetEventPerfInfo(nested->handler, 0, event_id, &event_info));
377c244758eSToby Isaac perf[tree_node] = *event_info;
378b9321188SToby Isaac } else {
379c244758eSToby Isaac PetscCall(PetscArrayzero(&perf[tree_node], 1));
380b9321188SToby Isaac }
381b9321188SToby Isaac }
382b9321188SToby Isaac *tree_p = tree;
383b9321188SToby Isaac *perf_p = perf;
384b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
385b9321188SToby Isaac }
PETSC_PRAGMA_DIAGNOSTIC_IGNORED_END()386*31a765c4SPierre Jolivet PETSC_PRAGMA_DIAGNOSTIC_IGNORED_END()
387b9321188SToby Isaac
388b9321188SToby Isaac static PetscErrorCode PetscLogHandlerView_Nested(PetscLogHandler handler, PetscViewer viewer)
389b9321188SToby Isaac {
390b9321188SToby Isaac PetscLogHandler_Nested nested = (PetscLogHandler_Nested)handler->data;
391b9321188SToby Isaac PetscNestedEventNode *nodes;
392b9321188SToby Isaac PetscEventPerfInfo *perf;
393b9321188SToby Isaac PetscLogGlobalNames global_events;
394b9321188SToby Isaac PetscNestedEventTree tree;
395b9321188SToby Isaac PetscViewerFormat format;
396b9321188SToby Isaac MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
397b9321188SToby Isaac
398b9321188SToby Isaac PetscFunctionBegin;
399b9321188SToby Isaac PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, nested->state->registry, &global_events));
400b9321188SToby Isaac PetscCall(PetscLogNestedCreatePerfNodes(comm, nested, global_events, &nodes, &perf));
401b9321188SToby Isaac tree.comm = comm;
402b9321188SToby Isaac tree.global_events = global_events;
403b9321188SToby Isaac tree.perf = perf;
404b9321188SToby Isaac tree.nodes = nodes;
405b9321188SToby Isaac PetscCall(PetscViewerGetFormat(viewer, &format));
406b9321188SToby Isaac if (format == PETSC_VIEWER_ASCII_XML) {
407b9321188SToby Isaac PetscCall(PetscLogHandlerView_Nested_XML(nested, &tree, viewer));
408b9321188SToby Isaac } else if (format == PETSC_VIEWER_ASCII_FLAMEGRAPH) {
409b9321188SToby Isaac PetscCall(PetscLogHandlerView_Nested_Flamegraph(nested, &tree, viewer));
410b9321188SToby Isaac } else SETERRQ(comm, PETSC_ERR_ARG_INCOMP, "No nested viewer for this format");
411b9321188SToby Isaac PetscCall(PetscLogGlobalNamesDestroy(&global_events));
412b9321188SToby Isaac PetscCall(PetscFree(tree.nodes));
413b9321188SToby Isaac PetscCall(PetscFree(tree.perf));
414b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
415b9321188SToby Isaac }
416b9321188SToby Isaac
417b9321188SToby Isaac /*MC
418294de794SToby Isaac PETSCLOGHANDLERNESTED - PETSCLOGHANDLERNESTED = "nested" - A `PetscLogHandler` that collects data for PETSc's
419b665b14eSToby Isaac XML and flamegraph profiling log viewers. A log handler of this type is created and started by
420b665b14eSToby Isaac by `PetscLogNestedBegin()`.
421b9321188SToby Isaac
422b9321188SToby Isaac Level: developer
423b9321188SToby Isaac
424b9321188SToby Isaac .seealso: [](ch_profiling), `PetscLogHandler`
425b9321188SToby Isaac M*/
426b9321188SToby Isaac
PetscLogHandlerCreate_Nested(PetscLogHandler handler)427b9321188SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Nested(PetscLogHandler handler)
428b9321188SToby Isaac {
429b9321188SToby Isaac PetscFunctionBegin;
430b9321188SToby Isaac PetscCall(PetscLogHandlerContextCreate_Nested(PetscObjectComm((PetscObject)handler), (PetscLogHandler_Nested *)&handler->data));
431b9321188SToby Isaac handler->ops->destroy = PetscLogHandlerDestroy_Nested;
432b9321188SToby Isaac handler->ops->stagepush = PetscLogHandlerStagePush_Nested;
433b9321188SToby Isaac handler->ops->stagepop = PetscLogHandlerStagePop_Nested;
434b9321188SToby Isaac handler->ops->eventbegin = PetscLogHandlerEventBegin_Nested;
435b9321188SToby Isaac handler->ops->eventend = PetscLogHandlerEventEnd_Nested;
436b9321188SToby Isaac handler->ops->eventsync = PetscLogHandlerEventSync_Nested;
437b9321188SToby Isaac handler->ops->objectcreate = PetscLogHandlerObjectCreate_Nested;
438b9321188SToby Isaac handler->ops->objectdestroy = PetscLogHandlerObjectDestroy_Nested;
439b9321188SToby Isaac handler->ops->view = PetscLogHandlerView_Nested;
440b9321188SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
441b9321188SToby Isaac }
442