xref: /petsc/src/sys/logging/handler/impls/nested/lognested.c (revision dff009be1503132eb32dbaedc8512f49b68e5e4c)
1b9321188SToby Isaac 
2b9321188SToby Isaac #include <petscviewer.h>
3b9321188SToby Isaac #include "lognested.h"
4b9321188SToby Isaac #include "xmlviewer.h"
5b9321188SToby Isaac 
6b9321188SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerNestedSetThreshold(PetscLogHandler h, PetscLogDouble newThresh, PetscLogDouble *oldThresh)
7b9321188SToby Isaac {
8b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
9b9321188SToby Isaac 
10b9321188SToby Isaac   PetscFunctionBegin;
11b9321188SToby Isaac   if (oldThresh) *oldThresh = nested->threshold;
12b9321188SToby Isaac   if (newThresh == (PetscLogDouble)PETSC_DECIDE) newThresh = 0.01;
13b9321188SToby Isaac   if (newThresh == (PetscLogDouble)PETSC_DEFAULT) newThresh = 0.01;
14b9321188SToby Isaac   nested->threshold = PetscMax(newThresh, 0.0);
15b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
16b9321188SToby Isaac }
17b9321188SToby Isaac 
18b9321188SToby Isaac static PetscErrorCode PetscLogEventGetNestedEvent(PetscLogHandler h, PetscLogEvent e, PetscLogEvent *nested_event)
19b9321188SToby Isaac {
20b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
21b9321188SToby Isaac   NestedIdPair           key;
22b9321188SToby Isaac   PetscHashIter          iter;
23b9321188SToby Isaac   PetscBool              missing;
24b9321188SToby Isaac   PetscLogState          state;
25b9321188SToby Isaac 
26b9321188SToby Isaac   PetscFunctionBegin;
27b9321188SToby Isaac   PetscCall(PetscLogHandlerGetState(h, &state));
28b9321188SToby Isaac   PetscCall(PetscIntStackTop(nested->stack, &(key.root)));
29b9321188SToby Isaac   key.leaf = NestedIdFromEvent(e);
30b9321188SToby Isaac   PetscCall(PetscNestedHashPut(nested->pair_map, key, &iter, &missing));
31b9321188SToby Isaac   if (missing) {
32b9321188SToby Isaac     // register a new nested event
33b9321188SToby Isaac     char              name[BUFSIZ];
34b9321188SToby Isaac     PetscLogEventInfo event_info;
35b9321188SToby Isaac     PetscLogEventInfo nested_event_info;
36b9321188SToby Isaac 
37b9321188SToby Isaac     PetscCall(PetscLogStateEventGetInfo(state, e, &event_info));
38b9321188SToby Isaac     PetscCall(PetscLogStateEventGetInfo(nested->state, key.root, &nested_event_info));
39b9321188SToby Isaac     PetscCall(PetscSNPrintf(name, sizeof(name) - 1, "%s;%s", nested_event_info.name, event_info.name));
40b9321188SToby Isaac     PetscCall(PetscLogStateEventRegister(nested->state, name, event_info.classid, nested_event));
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 
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));
58b9321188SToby Isaac   PetscCall(PetscIntStackTop(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];
64b9321188SToby Isaac 
65b9321188SToby Isaac     PetscCall(PetscLogStateStageGetInfo(state, stage, &stage_info));
66b9321188SToby Isaac     if (key.root >= 0) {
67b9321188SToby Isaac       PetscLogEventInfo nested_event_info;
68b9321188SToby Isaac 
69b9321188SToby Isaac       PetscCall(PetscLogStateEventGetInfo(nested->state, key.root, &nested_event_info));
70b9321188SToby Isaac       PetscCall(PetscSNPrintf(name, sizeof(name) - 1, "%s;%s", nested_event_info.name, stage_info.name));
71b9321188SToby Isaac     } else {
72b9321188SToby Isaac       PetscCall(PetscSNPrintf(name, sizeof(name) - 1, "%s", stage_info.name));
73b9321188SToby Isaac     }
74b9321188SToby Isaac     PetscCall(PetscLogStateEventRegister(nested->state, name, nested->nested_stage_id, nested_event));
75b9321188SToby Isaac     PetscCall(PetscNestedHashIterSet(nested->pair_map, iter, *nested_event));
76b9321188SToby Isaac   } else {
77b9321188SToby Isaac     PetscCall(PetscNestedHashIterGet(nested->pair_map, iter, nested_event));
78b9321188SToby Isaac   }
79b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
80b9321188SToby Isaac }
81b9321188SToby Isaac 
82b9321188SToby Isaac static PetscErrorCode PetscLogNestedCheckNested(PetscLogHandler h, NestedId leaf, PetscLogEvent nested_event)
83b9321188SToby Isaac {
84b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
85b9321188SToby Isaac   NestedIdPair           key;
86b9321188SToby Isaac   NestedId               val;
87b9321188SToby Isaac 
88b9321188SToby Isaac   PetscFunctionBegin;
89b9321188SToby Isaac   PetscCall(PetscIntStackTop(nested->stack, &(key.root)));
90b9321188SToby Isaac   key.leaf = leaf;
91b9321188SToby Isaac   PetscCall(PetscNestedHashGet(nested->pair_map, key, &val));
92b9321188SToby Isaac   PetscCheck(val == nested_event, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Logging events and stages are not nested, nested logging cannot be used");
93b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
94b9321188SToby Isaac }
95b9321188SToby Isaac 
96b9321188SToby Isaac static PetscErrorCode PetscLogHandlerEventBegin_Nested(PetscLogHandler h, PetscLogEvent e, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
97b9321188SToby Isaac {
98b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
99b9321188SToby Isaac   PetscLogEvent          nested_event;
100b9321188SToby Isaac 
101b9321188SToby Isaac   PetscFunctionBegin;
102b9321188SToby Isaac   PetscCall(PetscLogEventGetNestedEvent(h, e, &nested_event));
103b9321188SToby Isaac   PetscCall(PetscLogHandlerEventBegin(nested->handler, nested_event, o1, o2, o3, o4));
104b9321188SToby Isaac   PetscCall(PetscIntStackPush(nested->stack, nested_event));
105b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
106b9321188SToby Isaac }
107b9321188SToby Isaac 
108b9321188SToby Isaac static PetscErrorCode PetscLogHandlerEventEnd_Nested(PetscLogHandler h, PetscLogEvent e, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
109b9321188SToby Isaac {
110b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
111b9321188SToby Isaac   PetscLogEvent          nested_event;
112b9321188SToby Isaac 
113b9321188SToby Isaac   PetscFunctionBegin;
114b9321188SToby Isaac   PetscCall(PetscIntStackPop(nested->stack, &nested_event));
115b9321188SToby Isaac   if (PetscDefined(USE_DEBUG)) PetscCall(PetscLogNestedCheckNested(h, NestedIdFromEvent(e), nested_event));
116b9321188SToby Isaac   PetscCall(PetscLogHandlerEventEnd(nested->handler, nested_event, o1, o2, o3, o4));
117b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
118b9321188SToby Isaac }
119b9321188SToby Isaac 
120b9321188SToby Isaac static PetscErrorCode PetscLogHandlerEventSync_Nested(PetscLogHandler h, PetscLogEvent e, MPI_Comm comm)
121b9321188SToby Isaac {
122b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
123b9321188SToby Isaac   PetscLogEvent          nested_event;
124b9321188SToby Isaac 
125b9321188SToby Isaac   PetscFunctionBegin;
126b9321188SToby Isaac   PetscCall(PetscLogEventGetNestedEvent(h, e, &nested_event));
127b9321188SToby Isaac   PetscCall(PetscLogHandlerEventSync(nested->handler, nested_event, comm));
128b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
129b9321188SToby Isaac }
130b9321188SToby Isaac 
131b9321188SToby Isaac static PetscErrorCode PetscLogHandlerStagePush_Nested(PetscLogHandler h, PetscLogStage stage)
132b9321188SToby Isaac {
133b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
134b9321188SToby Isaac   PetscLogEvent          nested_event;
135b9321188SToby Isaac 
136b9321188SToby Isaac   PetscFunctionBegin;
137b9321188SToby Isaac   if (nested->nested_stage_id == -1) PetscCall(PetscClassIdRegister("LogNestedStage", &nested->nested_stage_id));
138b9321188SToby Isaac   PetscCall(PetscLogStageGetNestedEvent(h, stage, &nested_event));
139b9321188SToby Isaac   PetscCall(PetscLogHandlerEventBegin(nested->handler, nested_event, NULL, NULL, NULL, NULL));
140b9321188SToby Isaac   PetscCall(PetscIntStackPush(nested->stack, nested_event));
141b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
142b9321188SToby Isaac }
143b9321188SToby Isaac 
144b9321188SToby Isaac static PetscErrorCode PetscLogHandlerStagePop_Nested(PetscLogHandler h, PetscLogStage stage)
145b9321188SToby Isaac {
146b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
147b9321188SToby Isaac   PetscLogEvent          nested_event;
148b9321188SToby Isaac 
149b9321188SToby Isaac   PetscFunctionBegin;
150b9321188SToby Isaac   PetscCall(PetscIntStackPop(nested->stack, &nested_event));
151b9321188SToby Isaac   if (PetscDefined(USE_DEBUG)) PetscCall(PetscLogNestedCheckNested(h, NestedIdFromStage(stage), nested_event));
152b9321188SToby Isaac   PetscCall(PetscLogHandlerEventEnd(nested->handler, nested_event, NULL, NULL, NULL, NULL));
153b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
154b9321188SToby Isaac }
155b9321188SToby Isaac 
156b9321188SToby Isaac static PetscErrorCode PetscLogHandlerContextCreate_Nested(MPI_Comm comm, PetscLogHandler_Nested *nested_p)
157b9321188SToby Isaac {
158b9321188SToby Isaac   PetscLogStage          root_stage;
159b9321188SToby Isaac   PetscLogHandler_Nested nested;
160b9321188SToby Isaac 
161b9321188SToby Isaac   PetscFunctionBegin;
162b9321188SToby Isaac   PetscCall(PetscNew(nested_p));
163b9321188SToby Isaac   nested = *nested_p;
164b9321188SToby Isaac   PetscCall(PetscLogStateCreate(&nested->state));
165b9321188SToby Isaac   PetscCall(PetscIntStackCreate(&nested->stack));
166b9321188SToby Isaac   nested->nested_stage_id = -1;
167b9321188SToby Isaac   nested->threshold       = 0.01;
168b9321188SToby Isaac   PetscCall(PetscNestedHashCreate(&nested->pair_map));
169b9321188SToby Isaac   PetscCall(PetscLogHandlerCreate(comm, &nested->handler));
170294de794SToby Isaac   PetscCall(PetscLogHandlerSetType(nested->handler, PETSCLOGHANDLERDEFAULT));
171b9321188SToby Isaac   PetscCall(PetscLogHandlerSetState(nested->handler, nested->state));
172b9321188SToby Isaac   PetscCall(PetscLogStateStageRegister(nested->state, "", &root_stage));
173b9321188SToby Isaac   PetscAssert(root_stage == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "root stage not zero");
174b9321188SToby Isaac   PetscCall(PetscLogHandlerStagePush(nested->handler, root_stage));
175b9321188SToby Isaac   PetscCall(PetscLogStateStagePush(nested->state, root_stage));
176b9321188SToby Isaac   PetscCall(PetscIntStackPush(nested->stack, -1));
177b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
178b9321188SToby Isaac }
179b9321188SToby Isaac 
180b9321188SToby Isaac static PetscErrorCode PetscLogHandlerObjectCreate_Nested(PetscLogHandler h, PetscObject obj)
181b9321188SToby Isaac {
182b9321188SToby Isaac   PetscClassId           classid;
183b9321188SToby Isaac   PetscInt               num_registered, num_nested_registered;
184b9321188SToby Isaac   PetscLogState          state;
185b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
186b9321188SToby Isaac 
187b9321188SToby Isaac   PetscFunctionBegin;
188b9321188SToby Isaac   // register missing objects
189b9321188SToby Isaac   PetscCall(PetscObjectGetClassId(obj, &classid));
190b9321188SToby Isaac   PetscCall(PetscLogHandlerGetState(h, &state));
191b9321188SToby Isaac   PetscCall(PetscLogStateGetNumClasses(nested->state, &num_nested_registered));
192b9321188SToby Isaac   PetscCall(PetscLogStateGetNumClasses(state, &num_registered));
193b9321188SToby Isaac   for (PetscLogClass c = num_nested_registered; c < num_registered; c++) {
194b9321188SToby Isaac     PetscLogClassInfo class_info;
195b9321188SToby Isaac     PetscLogClass     nested_c;
196b9321188SToby Isaac 
197b9321188SToby Isaac     PetscCall(PetscLogStateClassGetInfo(state, c, &class_info));
198b9321188SToby Isaac     PetscCall(PetscLogStateClassRegister(nested->state, class_info.name, class_info.classid, &nested_c));
199b9321188SToby Isaac   }
200b9321188SToby Isaac   PetscCall(PetscLogHandlerObjectCreate(nested->handler, obj));
201b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
202b9321188SToby Isaac }
203b9321188SToby Isaac 
204b9321188SToby Isaac static PetscErrorCode PetscLogHandlerObjectDestroy_Nested(PetscLogHandler h, PetscObject obj)
205b9321188SToby Isaac {
206b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
207b9321188SToby Isaac 
208b9321188SToby Isaac   PetscFunctionBegin;
209b9321188SToby Isaac   PetscCall(PetscLogHandlerObjectDestroy(nested->handler, obj));
210b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
211b9321188SToby Isaac }
212b9321188SToby Isaac 
213b9321188SToby Isaac static PetscErrorCode PetscLogHandlerDestroy_Nested(PetscLogHandler h)
214b9321188SToby Isaac {
215b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)h->data;
216b9321188SToby Isaac 
217b9321188SToby Isaac   PetscFunctionBegin;
218b9321188SToby Isaac   PetscCall(PetscLogStateStagePop(nested->state));
219b9321188SToby Isaac   PetscCall(PetscLogHandlerStagePop(nested->handler, 0));
220b9321188SToby Isaac   PetscCall(PetscLogStateDestroy(&nested->state));
221b9321188SToby Isaac   PetscCall(PetscIntStackDestroy(nested->stack));
222b9321188SToby Isaac   PetscCall(PetscNestedHashDestroy(&nested->pair_map));
223b9321188SToby Isaac   PetscCall(PetscLogHandlerDestroy(&nested->handler));
224b9321188SToby Isaac   PetscCall(PetscFree(nested));
225b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
226b9321188SToby Isaac }
227b9321188SToby Isaac 
228b9321188SToby Isaac static PetscErrorCode PetscLogNestedEventNodesOrderDepthFirst(PetscInt num_nodes, PetscInt parent, PetscNestedEventNode tree[], PetscInt *num_descendants)
229b9321188SToby Isaac {
230b9321188SToby Isaac   PetscInt node, start_loc;
231b9321188SToby Isaac   PetscFunctionBegin;
232b9321188SToby Isaac 
233b9321188SToby Isaac   node      = 0;
234b9321188SToby Isaac   start_loc = 0;
235b9321188SToby Isaac   while (node < num_nodes) {
236b9321188SToby Isaac     if (tree[node].parent == parent) {
237b9321188SToby Isaac       PetscInt             num_this_descendants = 0;
238b9321188SToby Isaac       PetscNestedEventNode tmp                  = tree[start_loc];
239b9321188SToby Isaac       tree[start_loc]                           = tree[node];
240b9321188SToby Isaac       tree[node]                                = tmp;
241b9321188SToby Isaac       PetscCall(PetscLogNestedEventNodesOrderDepthFirst(num_nodes - start_loc - 1, tree[start_loc].id, &tree[start_loc + 1], &num_this_descendants));
242b9321188SToby Isaac       tree[start_loc].num_descendants = num_this_descendants;
243b9321188SToby Isaac       *num_descendants += 1 + num_this_descendants;
244b9321188SToby Isaac       start_loc += 1 + num_this_descendants;
245b9321188SToby Isaac       node = start_loc;
246b9321188SToby Isaac     } else {
247b9321188SToby Isaac       node++;
248b9321188SToby Isaac     }
249b9321188SToby Isaac   }
250b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
251b9321188SToby Isaac }
252b9321188SToby Isaac 
253b9321188SToby Isaac static PetscErrorCode PetscLogNestedCreatePerfNodes(MPI_Comm comm, PetscLogHandler_Nested nested, PetscLogGlobalNames global_events, PetscNestedEventNode **tree_p, PetscEventPerfInfo **perf_p)
254b9321188SToby Isaac {
255b9321188SToby Isaac   PetscMPIInt           size;
256b9321188SToby Isaac   PetscInt              num_nodes;
257b9321188SToby Isaac   PetscInt              num_map_entries;
258b9321188SToby Isaac   PetscEventPerfInfo   *perf;
259b9321188SToby Isaac   NestedIdPair         *keys;
260b9321188SToby Isaac   NestedId             *vals;
261b9321188SToby Isaac   PetscInt              offset;
262b9321188SToby Isaac   PetscInt              num_descendants;
263b9321188SToby Isaac   PetscNestedEventNode *tree;
264b9321188SToby Isaac 
265b9321188SToby Isaac   PetscFunctionBegin;
266b9321188SToby Isaac   PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &num_nodes));
267b9321188SToby Isaac   PetscCall(PetscCalloc1(num_nodes, &tree));
268b9321188SToby Isaac   for (PetscInt node = 0; node < num_nodes; node++) {
269b9321188SToby Isaac     tree[node].id = node;
270b9321188SToby Isaac     PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, node, &tree[node].name));
271b9321188SToby Isaac     tree[node].parent = -1;
272b9321188SToby Isaac   }
273b9321188SToby Isaac   PetscCall(PetscNestedHashGetSize(nested->pair_map, &num_map_entries));
274b9321188SToby Isaac   PetscCall(PetscMalloc2(num_map_entries, &keys, num_map_entries, &vals));
275b9321188SToby Isaac   offset = 0;
276b9321188SToby Isaac   PetscCall(PetscNestedHashGetPairs(nested->pair_map, &offset, keys, vals));
277b9321188SToby Isaac   for (PetscInt k = 0; k < num_map_entries; k++) {
278b9321188SToby Isaac     NestedId root_local = keys[k].root;
279b9321188SToby Isaac     NestedId leaf_local = vals[k];
280b9321188SToby Isaac     PetscInt root_global;
281b9321188SToby Isaac     PetscInt leaf_global;
282b9321188SToby Isaac 
283b9321188SToby Isaac     PetscCall(PetscLogGlobalNamesLocalGetGlobal(global_events, leaf_local, &leaf_global));
284b9321188SToby Isaac     if (root_local >= 0) {
285b9321188SToby Isaac       PetscCall(PetscLogGlobalNamesLocalGetGlobal(global_events, root_local, &root_global));
286b9321188SToby Isaac       tree[leaf_global].parent = root_global;
287b9321188SToby Isaac     }
288b9321188SToby Isaac   }
289b9321188SToby Isaac   PetscCall(PetscFree2(keys, vals));
290b9321188SToby Isaac   PetscCallMPI(MPI_Comm_size(comm, &size));
291b9321188SToby Isaac   if (size > 1) { // get missing parents from other processes
292b9321188SToby Isaac     PetscInt *parents;
293b9321188SToby Isaac 
294b9321188SToby Isaac     PetscCall(PetscMalloc1(num_nodes, &parents));
295b9321188SToby Isaac     for (PetscInt node = 0; node < num_nodes; node++) parents[node] = tree[node].parent;
296b9321188SToby Isaac     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, parents, num_nodes, MPIU_INT, MPI_MAX, comm));
297b9321188SToby Isaac     for (PetscInt node = 0; node < num_nodes; node++) tree[node].parent = parents[node];
298b9321188SToby Isaac     PetscCall(PetscFree(parents));
299b9321188SToby Isaac   }
300b9321188SToby Isaac 
301b9321188SToby Isaac   num_descendants = 0;
302b9321188SToby Isaac   PetscCall(PetscLogNestedEventNodesOrderDepthFirst(num_nodes, -1, tree, &num_descendants));
303b9321188SToby Isaac   PetscAssert(num_descendants == num_nodes, comm, PETSC_ERR_PLIB, "Failed tree ordering invariant");
304b9321188SToby Isaac 
305b9321188SToby Isaac   PetscCall(PetscCalloc1(num_nodes, &perf));
306b9321188SToby Isaac   for (PetscInt node = 0; node < num_nodes; node++) {
307b9321188SToby Isaac     PetscInt event_id;
308b9321188SToby Isaac 
309b9321188SToby Isaac     PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, node, &event_id));
310b9321188SToby Isaac     if (event_id >= 0) {
311b9321188SToby Isaac       PetscEventPerfInfo *event_info;
312b9321188SToby Isaac 
313*dff009beSToby Isaac       PetscCall(PetscLogHandlerGetEventPerfInfo(nested->handler, 0, event_id, &event_info));
314b9321188SToby Isaac       perf[node] = *event_info;
315b9321188SToby Isaac     } else {
316b9321188SToby Isaac       PetscCall(PetscArrayzero(&perf[node], 1));
317b9321188SToby Isaac     }
318b9321188SToby Isaac   }
319b9321188SToby Isaac   *tree_p = tree;
320b9321188SToby Isaac   *perf_p = perf;
321b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
322b9321188SToby Isaac }
323b9321188SToby Isaac 
324b9321188SToby Isaac static PetscErrorCode PetscLogHandlerView_Nested(PetscLogHandler handler, PetscViewer viewer)
325b9321188SToby Isaac {
326b9321188SToby Isaac   PetscLogHandler_Nested nested = (PetscLogHandler_Nested)handler->data;
327b9321188SToby Isaac   PetscNestedEventNode  *nodes;
328b9321188SToby Isaac   PetscEventPerfInfo    *perf;
329b9321188SToby Isaac   PetscLogGlobalNames    global_events;
330b9321188SToby Isaac   PetscNestedEventTree   tree;
331b9321188SToby Isaac   PetscViewerFormat      format;
332b9321188SToby Isaac   MPI_Comm               comm = PetscObjectComm((PetscObject)viewer);
333b9321188SToby Isaac 
334b9321188SToby Isaac   PetscFunctionBegin;
335b9321188SToby Isaac   PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, nested->state->registry, &global_events));
336b9321188SToby Isaac   PetscCall(PetscLogNestedCreatePerfNodes(comm, nested, global_events, &nodes, &perf));
337b9321188SToby Isaac   tree.comm          = comm;
338b9321188SToby Isaac   tree.global_events = global_events;
339b9321188SToby Isaac   tree.perf          = perf;
340b9321188SToby Isaac   tree.nodes         = nodes;
341b9321188SToby Isaac   PetscCall(PetscViewerGetFormat(viewer, &format));
342b9321188SToby Isaac   if (format == PETSC_VIEWER_ASCII_XML) {
343b9321188SToby Isaac     PetscCall(PetscLogHandlerView_Nested_XML(nested, &tree, viewer));
344b9321188SToby Isaac   } else if (format == PETSC_VIEWER_ASCII_FLAMEGRAPH) {
345b9321188SToby Isaac     PetscCall(PetscLogHandlerView_Nested_Flamegraph(nested, &tree, viewer));
346b9321188SToby Isaac   } else SETERRQ(comm, PETSC_ERR_ARG_INCOMP, "No nested viewer for this format");
347b9321188SToby Isaac   PetscCall(PetscLogGlobalNamesDestroy(&global_events));
348b9321188SToby Isaac   PetscCall(PetscFree(tree.nodes));
349b9321188SToby Isaac   PetscCall(PetscFree(tree.perf));
350b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
351b9321188SToby Isaac }
352b9321188SToby Isaac 
353b9321188SToby Isaac /*MC
354294de794SToby Isaac   PETSCLOGHANDLERNESTED - PETSCLOGHANDLERNESTED = "nested" -  A `PetscLogHandler` that collects data for PETSc's
355b665b14eSToby Isaac   XML and flamegraph profiling log viewers.  A log handler of this type is created and started by
356b665b14eSToby Isaac   by `PetscLogNestedBegin()`.
357b9321188SToby Isaac 
358b9321188SToby Isaac   Level: developer
359b9321188SToby Isaac 
360b9321188SToby Isaac .seealso: [](ch_profiling), `PetscLogHandler`
361b9321188SToby Isaac M*/
362b9321188SToby Isaac 
363b9321188SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Nested(PetscLogHandler handler)
364b9321188SToby Isaac {
365b9321188SToby Isaac   PetscFunctionBegin;
366b9321188SToby Isaac   PetscCall(PetscLogHandlerContextCreate_Nested(PetscObjectComm((PetscObject)handler), (PetscLogHandler_Nested *)&handler->data));
367b9321188SToby Isaac   handler->ops->destroy       = PetscLogHandlerDestroy_Nested;
368b9321188SToby Isaac   handler->ops->stagepush     = PetscLogHandlerStagePush_Nested;
369b9321188SToby Isaac   handler->ops->stagepop      = PetscLogHandlerStagePop_Nested;
370b9321188SToby Isaac   handler->ops->eventbegin    = PetscLogHandlerEventBegin_Nested;
371b9321188SToby Isaac   handler->ops->eventend      = PetscLogHandlerEventEnd_Nested;
372b9321188SToby Isaac   handler->ops->eventsync     = PetscLogHandlerEventSync_Nested;
373b9321188SToby Isaac   handler->ops->objectcreate  = PetscLogHandlerObjectCreate_Nested;
374b9321188SToby Isaac   handler->ops->objectdestroy = PetscLogHandlerObjectDestroy_Nested;
375b9321188SToby Isaac   handler->ops->view          = PetscLogHandlerView_Nested;
376b9321188SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
377b9321188SToby Isaac }
378