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