Actual source code: logtrace.c

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

  4: typedef struct _n_PetscLogHandler_Trace *PetscLogHandler_Trace;
  5: struct _n_PetscLogHandler_Trace {
  6:   FILE          *petsc_tracefile;
  7:   size_t         petsc_tracelevel;
  8:   char           petsc_tracespace[128];
  9:   PetscLogDouble petsc_tracetime;
 10: };

 12: static PetscErrorCode PetscLogHandlerEventBegin_Trace(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
 13: {
 14:   PetscLogHandler_Trace tr = (PetscLogHandler_Trace)h->data;
 15:   PetscLogEventInfo     event_info;
 16:   PetscLogDouble        cur_time;
 17:   PetscMPIInt           rank;
 18:   PetscLogState         state;
 19:   PetscLogStage         stage;

 21:   PetscFunctionBegin;
 22:   if (!tr->petsc_tracetime) PetscCall(PetscTime(&tr->petsc_tracetime));
 23:   tr->petsc_tracelevel++;
 24:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)h), &rank));
 25:   PetscCall(PetscLogHandlerGetState(h, &state));
 26:   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
 27:   /* Log performance info */
 28:   PetscCall(PetscTime(&cur_time));
 29:   PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
 30:   PetscCall(PetscFPrintf(PETSC_COMM_SELF, tr->petsc_tracefile, "%s[%d] %g Event begin: %s\n", tr->petsc_tracespace, rank, cur_time - tr->petsc_tracetime, event_info.name));
 31:   for (size_t i = 0; i < PetscMin(sizeof(tr->petsc_tracespace), 2 * tr->petsc_tracelevel); i++) tr->petsc_tracespace[i] = ' ';
 32:   tr->petsc_tracespace[PetscMin(127, 2 * tr->petsc_tracelevel)] = '\0';
 33:   PetscCall(PetscFFlush(tr->petsc_tracefile));
 34:   PetscFunctionReturn(PETSC_SUCCESS);
 35: }

 37: static PetscErrorCode PetscLogHandlerEventEnd_Trace(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
 38: {
 39:   PetscLogHandler_Trace tr = (PetscLogHandler_Trace)h->data;
 40:   PetscLogEventInfo     event_info;
 41:   PetscLogDouble        cur_time;
 42:   PetscLogState         state;
 43:   PetscLogStage         stage;
 44:   PetscMPIInt           rank;

 46:   PetscFunctionBegin;
 47:   tr->petsc_tracelevel--;
 48:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)h), &rank));
 49:   PetscCall(PetscLogHandlerGetState(h, &state));
 50:   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
 51:   /* Log performance info */
 52:   for (size_t i = 0; i < PetscMin(sizeof(tr->petsc_tracespace), 2 * tr->petsc_tracelevel); i++) tr->petsc_tracespace[i] = ' ';
 53:   tr->petsc_tracespace[PetscMin(127, 2 * tr->petsc_tracelevel)] = '\0';
 54:   PetscCall(PetscTime(&cur_time));
 55:   PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
 56:   PetscCall(PetscFPrintf(PETSC_COMM_SELF, tr->petsc_tracefile, "%s[%d] %g Event end: %s\n", tr->petsc_tracespace, rank, cur_time - tr->petsc_tracetime, event_info.name));
 57:   PetscCall(PetscFFlush(tr->petsc_tracefile));
 58:   PetscFunctionReturn(PETSC_SUCCESS);
 59: }

 61: static PetscErrorCode PetscLogHandlerDestroy_Trace(PetscLogHandler h)
 62: {
 63:   PetscFunctionBegin;
 64:   PetscCall(PetscFree(h->data));
 65:   PetscFunctionReturn(PETSC_SUCCESS);
 66: }

 68: /*MC
 69:   PETSCLOGHANDLERTRACE - PETSCLOGHANDLERTRACE = "trace" -  A
 70:   `PetscLogHandler` that collects data for PETSc's tracing log viewer.
 71:   A log handler of this type is created and started by `PetscLogTraceBegin()`.

 73:   Level: developer

 75: .seealso: [](ch_profiling), `PetscLogHandler`, `PetscLogHandlerCreateTrace()`
 76: M*/

 78: PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Trace(PetscLogHandler handler)
 79: {
 80:   PetscLogHandler_Trace tr;

 82:   PetscFunctionBegin;
 83:   PetscCall(PetscNew(&tr));
 84:   handler->data            = (void *)tr;
 85:   handler->ops->eventbegin = PetscLogHandlerEventBegin_Trace;
 86:   handler->ops->eventend   = PetscLogHandlerEventEnd_Trace;
 87:   handler->ops->destroy    = PetscLogHandlerDestroy_Trace;
 88:   PetscFunctionReturn(PETSC_SUCCESS);
 89: }

 91: /*@C
 92:   PetscLogHandlerCreateTrace - Create a logger that traces events and stages to a given file descriptor

 94:   Collective

 96:   Input Parameters:
 97: + comm - an MPI communicator
 98: - file - a file descriptor

100:   Output Parameters:
101: . handler - a `PetscLogHandler of type `PETSCLOGHANDLERTRACE`

103:   Level: developer

105:   Notes:
106:   Most users can just use `PetscLogTraceBegin()` to create and immediately start (`PetscLogHandlerStart()`) a tracing log handler

108: .seealso: [](ch_profiling), `PetscLogHandler`, `PetscLogHandlerTraceBegin()`
109: @*/
110: PetscErrorCode PetscLogHandlerCreateTrace(MPI_Comm comm, FILE *file, PetscLogHandler *handler)
111: {
112:   PetscLogHandler       h;
113:   PetscLogHandler_Trace tr;

115:   PetscFunctionBegin;
116:   PetscCall(PetscLogHandlerCreate(comm, handler));
117:   h = *handler;
118:   PetscCall(PetscLogHandlerSetType(h, PETSCLOGHANDLERTRACE));
119:   tr                  = (PetscLogHandler_Trace)h->data;
120:   tr->petsc_tracefile = file;
121:   PetscFunctionReturn(PETSC_SUCCESS);
122: }