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