1 #include <petsc/private/logimpl.h> /*I "petscsys.h" I*/ 2 #include <petsc/private/loghandlerimpl.h> 3 #include <../src/sys/perfstubs/timer.h> 4 5 typedef struct _n_PetscEventPS { 6 void *timer; 7 int depth; 8 } PetscEventPS; 9 10 PETSC_LOG_RESIZABLE_ARRAY(PSArray, PetscEventPS, void *, NULL, NULL, NULL) 11 12 typedef struct _n_PetscLogHandler_Perfstubs *PetscLogHandler_Perfstubs; 13 14 struct _n_PetscLogHandler_Perfstubs { 15 PetscLogPSArray events; 16 PetscLogPSArray stages; 17 PetscBool started_perfstubs; 18 }; 19 20 static PetscErrorCode PetscLogHandlerContextCreate_Perfstubs(PetscLogHandler_Perfstubs *ps_p) 21 { 22 PetscLogHandler_Perfstubs ps; 23 24 PetscFunctionBegin; 25 PetscCall(PetscNew(ps_p)); 26 ps = *ps_p; 27 PetscCall(PetscLogPSArrayCreate(128, &ps->events)); 28 PetscCall(PetscLogPSArrayCreate(8, &ps->stages)); 29 PetscFunctionReturn(PETSC_SUCCESS); 30 } 31 32 static PetscErrorCode PetscLogHandlerDestroy_Perfstubs(PetscLogHandler h) 33 { 34 PetscInt num_events, num_stages; 35 PetscLogHandler_Perfstubs ps = (PetscLogHandler_Perfstubs)h->data; 36 37 PetscFunctionBegin; 38 PetscCall(PetscLogPSArrayGetSize(ps->events, &num_events, NULL)); 39 for (PetscInt i = 0; i < num_events; i++) { 40 PetscEventPS event = {NULL, 0}; 41 42 PetscCall(PetscLogPSArrayGet(ps->events, i, &event)); 43 PetscStackCallExternalVoid("ps_timer_destroy_", ps_timer_destroy_(event.timer)); 44 } 45 PetscCall(PetscLogPSArrayDestroy(&ps->events)); 46 47 PetscCall(PetscLogPSArrayGetSize(ps->stages, &num_stages, NULL)); 48 for (PetscInt i = 0; i < num_stages; i++) { 49 PetscEventPS stage = {NULL, 0}; 50 51 PetscCall(PetscLogPSArrayGet(ps->stages, i, &stage)); 52 PetscStackCallExternalVoid("ps_timer_destroy_", ps_timer_destroy_(stage.timer)); 53 } 54 PetscCall(PetscLogPSArrayDestroy(&ps->stages)); 55 56 if (ps->started_perfstubs) PetscStackCallExternalVoid("ps_finalize_", ps_finalize_()); 57 PetscCall(PetscFree(ps)); 58 PetscFunctionReturn(PETSC_SUCCESS); 59 } 60 61 static PetscErrorCode PetscLogHandlerPSUpdateEvents(PetscLogHandler h) 62 { 63 PetscLogHandler_Perfstubs ps = (PetscLogHandler_Perfstubs)h->data; 64 PetscLogState state; 65 PetscInt num_events, num_events_old; 66 67 PetscFunctionBegin; 68 PetscCall(PetscLogHandlerGetState(h, &state)); 69 PetscCall(PetscLogStateGetNumEvents(state, &num_events)); 70 PetscCall(PetscLogPSArrayGetSize(ps->events, &num_events_old, NULL)); 71 for (PetscInt i = num_events_old; i < num_events; i++) { 72 PetscLogEventInfo event_info = {NULL, -1, PETSC_FALSE}; 73 PetscEventPS ps_event = {NULL, 0}; 74 PetscLogEvent ei; 75 76 PetscCall(PetscMPIIntCast(i, &ei)); 77 PetscCall(PetscLogStateEventGetInfo(state, ei, &event_info)); 78 PetscStackCallExternalVoid("ps_timer_create_", ps_event.timer = ps_timer_create_(event_info.name)); 79 ps_event.depth = 0; 80 PetscCall(PetscLogPSArrayPush(ps->events, ps_event)); 81 } 82 PetscFunctionReturn(PETSC_SUCCESS); 83 } 84 85 static PetscErrorCode PetscLogHandlerPSUpdateStages(PetscLogHandler h) 86 { 87 PetscLogHandler_Perfstubs ps = (PetscLogHandler_Perfstubs)h->data; 88 PetscLogState state; 89 PetscInt num_stages, num_stages_old; 90 91 PetscFunctionBegin; 92 PetscCall(PetscLogHandlerGetState(h, &state)); 93 PetscCall(PetscLogStateGetNumStages(state, &num_stages)); 94 PetscCall(PetscLogPSArrayGetSize(ps->stages, &num_stages_old, NULL)); 95 for (PetscInt i = num_stages_old; i < num_stages; i++) { 96 PetscLogStageInfo stage_info = {NULL}; 97 PetscEventPS ps_stage = {NULL, 0}; 98 PetscLogEvent si; 99 100 PetscCall(PetscMPIIntCast(i, &si)); 101 PetscCall(PetscLogStateStageGetInfo(state, si, &stage_info)); 102 PetscStackCallExternalVoid("ps_timer_create_", ps_stage.timer = ps_timer_create_(stage_info.name)); 103 ps_stage.depth = 0; 104 PetscCall(PetscLogPSArrayPush(ps->stages, ps_stage)); 105 } 106 PetscFunctionReturn(PETSC_SUCCESS); 107 } 108 109 static PetscErrorCode PetscLogHandlerEventBegin_Perfstubs(PetscLogHandler handler, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4) 110 { 111 PetscLogHandler_Perfstubs ps = (PetscLogHandler_Perfstubs)handler->data; 112 PetscEventPS ps_event = {NULL, 0}; 113 114 PetscFunctionBegin; 115 if (event >= ps->events->num_entries) PetscCall(PetscLogHandlerPSUpdateEvents(handler)); 116 PetscCall(PetscLogPSArrayGet(ps->events, event, &ps_event)); 117 ps_event.depth++; 118 PetscCall(PetscLogPSArraySet(ps->events, event, ps_event)); 119 if (ps_event.depth == 1 && ps_event.timer != NULL) PetscStackCallExternalVoid("ps_timer_start_", ps_timer_start_(ps_event.timer)); 120 PetscFunctionReturn(PETSC_SUCCESS); 121 } 122 123 static PetscErrorCode PetscLogHandlerEventEnd_Perfstubs(PetscLogHandler handler, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4) 124 { 125 PetscLogHandler_Perfstubs ps = (PetscLogHandler_Perfstubs)handler->data; 126 PetscEventPS ps_event = {NULL, 0}; 127 128 PetscFunctionBegin; 129 if (event >= ps->events->num_entries) PetscCall(PetscLogHandlerPSUpdateEvents(handler)); 130 PetscCall(PetscLogPSArrayGet(ps->events, event, &ps_event)); 131 ps_event.depth--; 132 PetscCall(PetscLogPSArraySet(ps->events, event, ps_event)); 133 if (ps_event.depth == 0 && ps_event.timer != NULL) PetscStackCallExternalVoid("ps_timer_stop_", ps_timer_stop_(ps_event.timer)); 134 PetscFunctionReturn(PETSC_SUCCESS); 135 } 136 137 static PetscErrorCode PetscLogHandlerStagePush_Perfstubs(PetscLogHandler handler, PetscLogStage stage) 138 { 139 PetscLogHandler_Perfstubs ps = (PetscLogHandler_Perfstubs)handler->data; 140 PetscEventPS ps_event = {NULL, 0}; 141 142 PetscFunctionBegin; 143 if (stage >= ps->stages->num_entries) PetscCall(PetscLogHandlerPSUpdateStages(handler)); 144 PetscCall(PetscLogPSArrayGet(ps->stages, stage, &ps_event)); 145 if (ps_event.timer != NULL) PetscStackCallExternalVoid("ps_timer_start_", ps_timer_start_(ps_event.timer)); 146 PetscFunctionReturn(PETSC_SUCCESS); 147 } 148 149 static PetscErrorCode PetscLogHandlerStagePop_Perfstubs(PetscLogHandler handler, PetscLogStage stage) 150 { 151 PetscLogHandler_Perfstubs ps = (PetscLogHandler_Perfstubs)handler->data; 152 PetscEventPS ps_event = {NULL, 0}; 153 154 PetscFunctionBegin; 155 if (stage >= ps->stages->num_entries) PetscCall(PetscLogHandlerPSUpdateStages(handler)); 156 PetscCall(PetscLogPSArrayGet(ps->stages, stage, &ps_event)); 157 if (ps_event.timer != NULL) PetscStackCallExternalVoid("ps_timer_stop_", ps_timer_stop_(ps_event.timer)); 158 PetscFunctionReturn(PETSC_SUCCESS); 159 } 160 161 /*MC 162 PETSCLOGHANDLERPERFSTUBS - PETSCLOGHANDLERPERFSTUBS = "perfstubs" - A 163 `PetscLogHandler` that collects data for the PerfStubs/TAU instrumentation 164 library. A log handler of this type is created and started by 165 `PetscLogPerfstubsBegin()`. 166 167 Level: developer 168 169 .seealso: [](ch_profiling), `PetscLogHandler` 170 M*/ 171 172 PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Perfstubs(PetscLogHandler handler) 173 { 174 PetscBool started_perfstubs; 175 PetscLogHandler_Perfstubs lps; 176 177 PetscFunctionBegin; 178 if (perfstubs_initialized == PERFSTUBS_UNKNOWN) { 179 PetscStackCallExternalVoid("ps_initialize_", ps_initialize_()); 180 started_perfstubs = PETSC_TRUE; 181 } else { 182 started_perfstubs = PETSC_FALSE; 183 } 184 PetscCheck(perfstubs_initialized == PERFSTUBS_SUCCESS, PetscObjectComm((PetscObject)handler), PETSC_ERR_LIB, "perfstubs could not be initialized"); 185 PetscCall(PetscLogHandlerContextCreate_Perfstubs(&lps)); 186 lps->started_perfstubs = started_perfstubs; 187 handler->data = (void *)lps; 188 handler->ops->destroy = PetscLogHandlerDestroy_Perfstubs; 189 handler->ops->eventbegin = PetscLogHandlerEventBegin_Perfstubs; 190 handler->ops->eventend = PetscLogHandlerEventEnd_Perfstubs; 191 handler->ops->stagepush = PetscLogHandlerStagePush_Perfstubs; 192 handler->ops->stagepop = PetscLogHandlerStagePop_Perfstubs; 193 PetscFunctionReturn(PETSC_SUCCESS); 194 } 195