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