xref: /petsc/src/sys/logging/handler/impls/mpe/logmpe.c (revision 856bee69f0e0908e75ff837867b1777dfb1ced96)
1*856bee69SToby Isaac #include <petsc/private/logimpl.h> /*I "petscsys.h" I*/
2*856bee69SToby Isaac #include <petsc/private/loghandlerimpl.h>
3*856bee69SToby Isaac #include <mpe.h>
4*856bee69SToby Isaac 
5*856bee69SToby Isaac typedef struct _n_PetscEventMPE {
6*856bee69SToby Isaac   int start;
7*856bee69SToby Isaac   int final;
8*856bee69SToby Isaac } PetscEventMPE;
9*856bee69SToby Isaac 
10*856bee69SToby Isaac PETSC_LOG_RESIZABLE_ARRAY(MPEArray, PetscEventMPE, void *, NULL, NULL, NULL);
11*856bee69SToby Isaac 
12*856bee69SToby Isaac typedef struct _n_PetscLogHandler_MPE *PetscLogHandler_MPE;
13*856bee69SToby Isaac 
14*856bee69SToby Isaac struct _n_PetscLogHandler_MPE {
15*856bee69SToby Isaac   PetscLogMPEArray events;
16*856bee69SToby Isaac };
17*856bee69SToby Isaac 
18*856bee69SToby Isaac static PetscErrorCode PetscLogHandlerContextCreate_MPE(PetscLogHandler_MPE *mpe_p)
19*856bee69SToby Isaac {
20*856bee69SToby Isaac   PetscLogHandler_MPE mpe;
21*856bee69SToby Isaac 
22*856bee69SToby Isaac   PetscFunctionBegin;
23*856bee69SToby Isaac   PetscCall(PetscNew(mpe_p));
24*856bee69SToby Isaac   mpe = *mpe_p;
25*856bee69SToby Isaac   PetscCall(PetscLogMPEArrayCreate(128, &mpe->events));
26*856bee69SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
27*856bee69SToby Isaac }
28*856bee69SToby Isaac 
29*856bee69SToby Isaac static PetscErrorCode PetscLogHandlerDestroy_MPE(PetscLogHandler h)
30*856bee69SToby Isaac {
31*856bee69SToby Isaac   PetscLogHandler_MPE mpe = (PetscLogHandler_MPE)h->data;
32*856bee69SToby Isaac 
33*856bee69SToby Isaac   PetscFunctionBegin;
34*856bee69SToby Isaac   PetscCall(PetscLogMPEArrayDestroy(&mpe->events));
35*856bee69SToby Isaac   PetscCall(PetscFree(mpe));
36*856bee69SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
37*856bee69SToby Isaac }
38*856bee69SToby Isaac 
39*856bee69SToby Isaac #define PETSC_RGB_COLORS_MAX 39
40*856bee69SToby Isaac static const char *PetscLogMPERGBColors[PETSC_RGB_COLORS_MAX] = {"OliveDrab:      ", "BlueViolet:     ", "CadetBlue:      ", "CornflowerBlue: ", "DarkGoldenrod:  ", "DarkGreen:      ", "DarkKhaki:      ", "DarkOliveGreen: ",
41*856bee69SToby Isaac                                                                  "DarkOrange:     ", "DarkOrchid:     ", "DarkSeaGreen:   ", "DarkSlateGray:  ", "DarkTurquoise:  ", "DeepPink:       ", "DarkKhaki:      ", "DimGray:        ",
42*856bee69SToby Isaac                                                                  "DodgerBlue:     ", "GreenYellow:    ", "HotPink:        ", "IndianRed:      ", "LavenderBlush:  ", "LawnGreen:      ", "LemonChiffon:   ", "LightCoral:     ",
43*856bee69SToby Isaac                                                                  "LightCyan:      ", "LightPink:      ", "LightSalmon:    ", "LightSlateGray: ", "LightYellow:    ", "LimeGreen:      ", "MediumPurple:   ", "MediumSeaGreen: ",
44*856bee69SToby Isaac                                                                  "MediumSlateBlue:", "MidnightBlue:   ", "MintCream:      ", "MistyRose:      ", "NavajoWhite:    ", "NavyBlue:       ", "OliveDrab:      "};
45*856bee69SToby Isaac 
46*856bee69SToby Isaac static PetscErrorCode PetscLogMPEGetRGBColor_Internal(const char *str[])
47*856bee69SToby Isaac {
48*856bee69SToby Isaac   static int idx = 0;
49*856bee69SToby Isaac 
50*856bee69SToby Isaac   PetscFunctionBegin;
51*856bee69SToby Isaac   *str = PetscLogMPERGBColors[idx];
52*856bee69SToby Isaac   idx  = (idx + 1) % PETSC_RGB_COLORS_MAX;
53*856bee69SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
54*856bee69SToby Isaac }
55*856bee69SToby Isaac 
56*856bee69SToby Isaac static PetscErrorCode PetscLogHandlerMPECreateEvent(const char name[], PetscLogMPEArray array)
57*856bee69SToby Isaac {
58*856bee69SToby Isaac   PetscEventMPE mpe_event;
59*856bee69SToby Isaac   PetscMPIInt   rank;
60*856bee69SToby Isaac 
61*856bee69SToby Isaac   PetscFunctionBegin;
62*856bee69SToby Isaac   MPE_Log_get_state_eventIDs(&mpe_event.start, &mpe_event.final);
63*856bee69SToby Isaac   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
64*856bee69SToby Isaac   if (rank == 0) {
65*856bee69SToby Isaac     const char *color;
66*856bee69SToby Isaac 
67*856bee69SToby Isaac     PetscCall(PetscLogMPEGetRGBColor_Internal(&color));
68*856bee69SToby Isaac     MPE_Describe_state(mpe_event.start, mpe_event.final, name, (char *)color);
69*856bee69SToby Isaac   }
70*856bee69SToby Isaac   PetscCall(PetscLogMPEArrayPush(array, mpe_event));
71*856bee69SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
72*856bee69SToby Isaac }
73*856bee69SToby Isaac 
74*856bee69SToby Isaac static PetscErrorCode PetscLogHandlerMPEUpdate(PetscLogHandler h)
75*856bee69SToby Isaac {
76*856bee69SToby Isaac   PetscLogHandler_MPE mpe = (PetscLogHandler_MPE)h->data;
77*856bee69SToby Isaac   PetscLogState       state;
78*856bee69SToby Isaac   PetscInt            num_events, num_events_old;
79*856bee69SToby Isaac 
80*856bee69SToby Isaac   PetscFunctionBegin;
81*856bee69SToby Isaac   PetscCall(PetscLogHandlerGetState(h, &state));
82*856bee69SToby Isaac   PetscCall(PetscLogStateGetNumEvents(state, &num_events));
83*856bee69SToby Isaac   PetscCall(PetscLogMPEArrayGetSize(mpe->events, &num_events_old, NULL));
84*856bee69SToby Isaac   for (PetscInt i = num_events_old; i < num_events; i++) {
85*856bee69SToby Isaac     PetscLogEventInfo event_info;
86*856bee69SToby Isaac 
87*856bee69SToby Isaac     PetscCall(PetscLogStateEventGetInfo(state, (PetscLogEvent)i, &event_info));
88*856bee69SToby Isaac     PetscCall(PetscLogHandlerMPECreateEvent(event_info.name, mpe->events));
89*856bee69SToby Isaac   }
90*856bee69SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
91*856bee69SToby Isaac }
92*856bee69SToby Isaac 
93*856bee69SToby Isaac static PetscErrorCode PetscLogHandlerEventBegin_MPE(PetscLogHandler handler, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
94*856bee69SToby Isaac {
95*856bee69SToby Isaac   PetscLogHandler_MPE mpe = (PetscLogHandler_MPE)handler->data;
96*856bee69SToby Isaac   PetscEventMPE       mpe_event;
97*856bee69SToby Isaac 
98*856bee69SToby Isaac   PetscFunctionBegin;
99*856bee69SToby Isaac   PetscCall(PetscLogHandlerMPEUpdate(handler));
100*856bee69SToby Isaac   PetscCall(PetscLogMPEArrayGet(mpe->events, event, &mpe_event));
101*856bee69SToby Isaac   PetscCall(MPE_Log_event(mpe_event.start, 0, NULL));
102*856bee69SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
103*856bee69SToby Isaac }
104*856bee69SToby Isaac 
105*856bee69SToby Isaac static PetscErrorCode PetscLogHandlerEventEnd_MPE(PetscLogHandler handler, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
106*856bee69SToby Isaac {
107*856bee69SToby Isaac   PetscLogHandler_MPE mpe = (PetscLogHandler_MPE)handler->data;
108*856bee69SToby Isaac   PetscEventMPE       mpe_event;
109*856bee69SToby Isaac 
110*856bee69SToby Isaac   PetscFunctionBegin;
111*856bee69SToby Isaac   PetscCall(PetscLogHandlerMPEUpdate(handler));
112*856bee69SToby Isaac   PetscCall(PetscLogMPEArrayGet(mpe->events, event, &mpe_event));
113*856bee69SToby Isaac   PetscCall(MPE_Log_event(mpe_event.final, 0, NULL));
114*856bee69SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
115*856bee69SToby Isaac }
116*856bee69SToby Isaac 
117*856bee69SToby Isaac /*MC
118*856bee69SToby Isaac   PETSC_LOG_HANDLER_MPE - PETSC_LOG_HANDLER_MPE = "mpe" -  A
119*856bee69SToby Isaac   `PetscLogHandler` that collects data for MPE, the MPI Parallel Enviornment for
120*856bee69SToby Isaac   performance visualization.
121*856bee69SToby Isaac 
122*856bee69SToby Isaac   Level: developer
123*856bee69SToby Isaac 
124*856bee69SToby Isaac .seealso: [](ch_profiling), `PetscLogHandler`
125*856bee69SToby Isaac M*/
126*856bee69SToby Isaac 
127*856bee69SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_MPE(PetscLogHandler handler)
128*856bee69SToby Isaac {
129*856bee69SToby Isaac   PetscFunctionBegin;
130*856bee69SToby Isaac   PetscCall(PetscLogHandlerContextCreate_MPE((PetscLogHandler_MPE *)&handler->data));
131*856bee69SToby Isaac   handler->ops->destroy    = PetscLogHandlerDestroy_MPE;
132*856bee69SToby Isaac   handler->ops->eventbegin = PetscLogHandlerEventBegin_MPE;
133*856bee69SToby Isaac   handler->ops->eventend   = PetscLogHandlerEventEnd_MPE;
134*856bee69SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
135*856bee69SToby Isaac }
136