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