Actual source code: logmpe.c

  1: #include <petsc/private/logimpl.h>
  2: #include <petsc/private/loghandlerimpl.h>
  3: #include <mpe.h>

  5: typedef struct _n_PetscEventMPE {
  6:   int start;
  7:   int final;
  8: } PetscEventMPE;

 10: PETSC_LOG_RESIZABLE_ARRAY(MPEArray, PetscEventMPE, void *, NULL, NULL, NULL);

 12: typedef struct _n_PetscLogHandler_MPE *PetscLogHandler_MPE;

 14: struct _n_PetscLogHandler_MPE {
 15:   PetscLogMPEArray events;
 16: };

 18: static PetscErrorCode PetscLogHandlerContextCreate_MPE(PetscLogHandler_MPE *mpe_p)
 19: {
 20:   PetscLogHandler_MPE mpe;

 22:   PetscFunctionBegin;
 23:   PetscCall(PetscNew(mpe_p));
 24:   mpe = *mpe_p;
 25:   PetscCall(PetscLogMPEArrayCreate(128, &mpe->events));
 26:   PetscFunctionReturn(PETSC_SUCCESS);
 27: }

 29: static PetscErrorCode PetscLogHandlerDestroy_MPE(PetscLogHandler h)
 30: {
 31:   PetscLogHandler_MPE mpe = (PetscLogHandler_MPE)h->data;

 33:   PetscFunctionBegin;
 34:   PetscCall(PetscLogMPEArrayDestroy(&mpe->events));
 35:   PetscCall(PetscFree(mpe));
 36:   PetscFunctionReturn(PETSC_SUCCESS);
 37: }

 39: #define PETSC_RGB_COLORS_MAX 39
 40: static const char *PetscLogMPERGBColors[PETSC_RGB_COLORS_MAX] = {"OliveDrab:      ", "BlueViolet:     ", "CadetBlue:      ", "CornflowerBlue: ", "DarkGoldenrod:  ", "DarkGreen:      ", "DarkKhaki:      ", "DarkOliveGreen: ",
 41:                                                                  "DarkOrange:     ", "DarkOrchid:     ", "DarkSeaGreen:   ", "DarkSlateGray:  ", "DarkTurquoise:  ", "DeepPink:       ", "DarkKhaki:      ", "DimGray:        ",
 42:                                                                  "DodgerBlue:     ", "GreenYellow:    ", "HotPink:        ", "IndianRed:      ", "LavenderBlush:  ", "LawnGreen:      ", "LemonChiffon:   ", "LightCoral:     ",
 43:                                                                  "LightCyan:      ", "LightPink:      ", "LightSalmon:    ", "LightSlateGray: ", "LightYellow:    ", "LimeGreen:      ", "MediumPurple:   ", "MediumSeaGreen: ",
 44:                                                                  "MediumSlateBlue:", "MidnightBlue:   ", "MintCream:      ", "MistyRose:      ", "NavajoWhite:    ", "NavyBlue:       ", "OliveDrab:      "};

 46: static PetscErrorCode PetscLogMPEGetRGBColor_Internal(const char *str[])
 47: {
 48:   static int idx = 0;

 50:   PetscFunctionBegin;
 51:   *str = PetscLogMPERGBColors[idx];
 52:   idx  = (idx + 1) % PETSC_RGB_COLORS_MAX;
 53:   PetscFunctionReturn(PETSC_SUCCESS);
 54: }

 56: static PetscErrorCode PetscLogHandlerMPECreateEvent(const char name[], PetscLogMPEArray array)
 57: {
 58:   PetscEventMPE mpe_event;
 59:   PetscMPIInt   rank;

 61:   PetscFunctionBegin;
 62:   MPE_Log_get_state_eventIDs(&mpe_event.start, &mpe_event.final);
 63:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
 64:   if (rank == 0) {
 65:     const char *color;

 67:     PetscCall(PetscLogMPEGetRGBColor_Internal(&color));
 68:     MPE_Describe_state(mpe_event.start, mpe_event.final, name, (char *)color);
 69:   }
 70:   PetscCall(PetscLogMPEArrayPush(array, mpe_event));
 71:   PetscFunctionReturn(PETSC_SUCCESS);
 72: }

 74: static PetscErrorCode PetscLogHandlerMPEUpdate(PetscLogHandler h)
 75: {
 76:   PetscLogHandler_MPE mpe = (PetscLogHandler_MPE)h->data;
 77:   PetscLogState       state;
 78:   PetscInt            num_events, num_events_old;

 80:   PetscFunctionBegin;
 81:   PetscCall(PetscLogHandlerGetState(h, &state));
 82:   PetscCall(PetscLogStateGetNumEvents(state, &num_events));
 83:   PetscCall(PetscLogMPEArrayGetSize(mpe->events, &num_events_old, NULL));
 84:   for (PetscInt i = num_events_old; i < num_events; i++) {
 85:     PetscLogEventInfo event_info;

 87:     PetscCall(PetscLogStateEventGetInfo(state, (PetscLogEvent)i, &event_info));
 88:     PetscCall(PetscLogHandlerMPECreateEvent(event_info.name, mpe->events));
 89:   }
 90:   PetscFunctionReturn(PETSC_SUCCESS);
 91: }

 93: static PetscErrorCode PetscLogHandlerEventBegin_MPE(PetscLogHandler handler, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
 94: {
 95:   PetscLogHandler_MPE mpe = (PetscLogHandler_MPE)handler->data;
 96:   PetscEventMPE       mpe_event;

 98:   PetscFunctionBegin;
 99:   PetscCall(PetscLogHandlerMPEUpdate(handler));
100:   PetscCall(PetscLogMPEArrayGet(mpe->events, event, &mpe_event));
101:   PetscCall(MPE_Log_event(mpe_event.start, 0, NULL));
102:   PetscFunctionReturn(PETSC_SUCCESS);
103: }

105: static PetscErrorCode PetscLogHandlerEventEnd_MPE(PetscLogHandler handler, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
106: {
107:   PetscLogHandler_MPE mpe = (PetscLogHandler_MPE)handler->data;
108:   PetscEventMPE       mpe_event;

110:   PetscFunctionBegin;
111:   PetscCall(PetscLogHandlerMPEUpdate(handler));
112:   PetscCall(PetscLogMPEArrayGet(mpe->events, event, &mpe_event));
113:   PetscCall(MPE_Log_event(mpe_event.final, 0, NULL));
114:   PetscFunctionReturn(PETSC_SUCCESS);
115: }

117: /*MC
118:   PETSCLOGHANDLERMPE - PETSCLOGHANDLERMPE = "mpe" -  A
119:   `PetscLogHandler` that collects data for MPE, the MPI Parallel Environment for
120:   performance visualization.  A log handler of this type is created and started
121:   by `PetscLogMPEBegin()`.

123:   Level: developer

125: .seealso: [](ch_profiling), `PetscLogHandler`
126: M*/

128: PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_MPE(PetscLogHandler handler)
129: {
130:   PetscFunctionBegin;
131:   PetscCall(PetscLogHandlerContextCreate_MPE((PetscLogHandler_MPE *)&handler->data));
132:   handler->ops->destroy    = PetscLogHandlerDestroy_MPE;
133:   handler->ops->eventbegin = PetscLogHandlerEventBegin_MPE;
134:   handler->ops->eventend   = PetscLogHandlerEventEnd_MPE;
135:   PetscFunctionReturn(PETSC_SUCCESS);
136: }