1*6467efc9SToby Isaac #include <petsc/private/logimpl.h> /*I "petscsys.h" I*/ 2*6467efc9SToby Isaac #include <petsc/private/loghandlerimpl.h> 3*6467efc9SToby Isaac 4*6467efc9SToby Isaac typedef struct _n_PetscLogHandler_Trace *PetscLogHandler_Trace; 5*6467efc9SToby Isaac struct _n_PetscLogHandler_Trace { 6*6467efc9SToby Isaac FILE *petsc_tracefile; 7*6467efc9SToby Isaac size_t petsc_tracelevel; 8*6467efc9SToby Isaac char petsc_tracespace[128]; 9*6467efc9SToby Isaac PetscLogDouble petsc_tracetime; 10*6467efc9SToby Isaac }; 11*6467efc9SToby Isaac 12*6467efc9SToby Isaac static PetscErrorCode PetscLogHandlerEventBegin_Trace(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4) 13*6467efc9SToby Isaac { 14*6467efc9SToby Isaac PetscLogHandler_Trace tr = (PetscLogHandler_Trace)h->data; 15*6467efc9SToby Isaac PetscLogEventInfo event_info; 16*6467efc9SToby Isaac PetscLogDouble cur_time; 17*6467efc9SToby Isaac PetscMPIInt rank; 18*6467efc9SToby Isaac PetscLogState state; 19*6467efc9SToby Isaac PetscLogStage stage; 20*6467efc9SToby Isaac 21*6467efc9SToby Isaac PetscFunctionBegin; 22*6467efc9SToby Isaac if (!tr->petsc_tracetime) PetscCall(PetscTime(&tr->petsc_tracetime)); 23*6467efc9SToby Isaac tr->petsc_tracelevel++; 24*6467efc9SToby Isaac PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)h), &rank)); 25*6467efc9SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state)); 26*6467efc9SToby Isaac PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 27*6467efc9SToby Isaac /* Log performance info */ 28*6467efc9SToby Isaac PetscCall(PetscTime(&cur_time)); 29*6467efc9SToby Isaac PetscCall(PetscLogStateEventGetInfo(state, event, &event_info)); 30*6467efc9SToby Isaac 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*6467efc9SToby Isaac for (size_t i = 0; i < PetscMin(sizeof(tr->petsc_tracespace), 2 * tr->petsc_tracelevel); i++) tr->petsc_tracespace[i] = ' '; 32*6467efc9SToby Isaac tr->petsc_tracespace[PetscMin(127, 2 * tr->petsc_tracelevel)] = '\0'; 33*6467efc9SToby Isaac PetscCall(PetscFFlush(tr->petsc_tracefile)); 34*6467efc9SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 35*6467efc9SToby Isaac } 36*6467efc9SToby Isaac 37*6467efc9SToby Isaac static PetscErrorCode PetscLogHandlerEventEnd_Trace(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4) 38*6467efc9SToby Isaac { 39*6467efc9SToby Isaac PetscLogHandler_Trace tr = (PetscLogHandler_Trace)h->data; 40*6467efc9SToby Isaac PetscLogEventInfo event_info; 41*6467efc9SToby Isaac PetscLogDouble cur_time; 42*6467efc9SToby Isaac PetscLogState state; 43*6467efc9SToby Isaac PetscLogStage stage; 44*6467efc9SToby Isaac PetscMPIInt rank; 45*6467efc9SToby Isaac 46*6467efc9SToby Isaac PetscFunctionBegin; 47*6467efc9SToby Isaac tr->petsc_tracelevel--; 48*6467efc9SToby Isaac PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)h), &rank)); 49*6467efc9SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state)); 50*6467efc9SToby Isaac PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 51*6467efc9SToby Isaac /* Log performance info */ 52*6467efc9SToby Isaac for (size_t i = 0; i < PetscMin(sizeof(tr->petsc_tracespace), 2 * tr->petsc_tracelevel); i++) tr->petsc_tracespace[i] = ' '; 53*6467efc9SToby Isaac tr->petsc_tracespace[PetscMin(127, 2 * tr->petsc_tracelevel)] = '\0'; 54*6467efc9SToby Isaac PetscCall(PetscTime(&cur_time)); 55*6467efc9SToby Isaac PetscCall(PetscLogStateEventGetInfo(state, event, &event_info)); 56*6467efc9SToby Isaac 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*6467efc9SToby Isaac PetscCall(PetscFFlush(tr->petsc_tracefile)); 58*6467efc9SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 59*6467efc9SToby Isaac } 60*6467efc9SToby Isaac 61*6467efc9SToby Isaac static PetscErrorCode PetscLogHandlerDestroy_Trace(PetscLogHandler h) 62*6467efc9SToby Isaac { 63*6467efc9SToby Isaac PetscFunctionBegin; 64*6467efc9SToby Isaac PetscCall(PetscFree(h->data)); 65*6467efc9SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 66*6467efc9SToby Isaac } 67*6467efc9SToby Isaac 68*6467efc9SToby Isaac /*MC 69*6467efc9SToby Isaac PETSC_LOG_HANDLER_TRACE - PETSC_LOG_HANDLER_TRACE = "trace" - A 70*6467efc9SToby Isaac `PetscLogHandler` that collects data for PETSc's tracing log viewer. 71*6467efc9SToby Isaac 72*6467efc9SToby Isaac Level: developer 73*6467efc9SToby Isaac 74*6467efc9SToby Isaac .seealso: [](ch_profiling), `PetscLogHandler`, `PetscLogHandlerCreateTrace()` 75*6467efc9SToby Isaac M*/ 76*6467efc9SToby Isaac 77*6467efc9SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Trace(PetscLogHandler handler) 78*6467efc9SToby Isaac { 79*6467efc9SToby Isaac PetscLogHandler_Trace tr; 80*6467efc9SToby Isaac 81*6467efc9SToby Isaac PetscFunctionBegin; 82*6467efc9SToby Isaac PetscCall(PetscNew(&tr)); 83*6467efc9SToby Isaac handler->data = (void *)tr; 84*6467efc9SToby Isaac handler->ops->eventbegin = PetscLogHandlerEventBegin_Trace; 85*6467efc9SToby Isaac handler->ops->eventend = PetscLogHandlerEventEnd_Trace; 86*6467efc9SToby Isaac handler->ops->destroy = PetscLogHandlerDestroy_Trace; 87*6467efc9SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 88*6467efc9SToby Isaac } 89*6467efc9SToby Isaac 90*6467efc9SToby Isaac /*@C 91*6467efc9SToby Isaac PetscLogHandlerCreateTrace - Create a logger that traces events and stages to a given file descriptor 92*6467efc9SToby Isaac 93*6467efc9SToby Isaac Collective 94*6467efc9SToby Isaac 95*6467efc9SToby Isaac Input Parameters: 96*6467efc9SToby Isaac + comm - an MPI communicator 97*6467efc9SToby Isaac - file - a file descriptor 98*6467efc9SToby Isaac 99*6467efc9SToby Isaac Output Parameters: 100*6467efc9SToby Isaac . handler - a `PetscLogHandler of type `PETSC_LOG_HANDLER_TRACE` 101*6467efc9SToby Isaac 102*6467efc9SToby Isaac Level: developer 103*6467efc9SToby Isaac 104*6467efc9SToby Isaac Notes: 105*6467efc9SToby Isaac Most users can just use `PetscLogTraceBegin()` to create and immediately start (`PetscLogHandlerStart()`) a tracing log handler 106*6467efc9SToby Isaac 107*6467efc9SToby Isaac .seealso: [](ch_profiling), `PetscLogHandler`, `PetscLogHandlerTraceBegin()` 108*6467efc9SToby Isaac @*/ 109*6467efc9SToby Isaac PetscErrorCode PetscLogHandlerCreateTrace(MPI_Comm comm, FILE *file, PetscLogHandler *handler) 110*6467efc9SToby Isaac { 111*6467efc9SToby Isaac PetscLogHandler h; 112*6467efc9SToby Isaac PetscLogHandler_Trace tr; 113*6467efc9SToby Isaac 114*6467efc9SToby Isaac PetscFunctionBegin; 115*6467efc9SToby Isaac PetscCall(PetscLogHandlerCreate(comm, handler)); 116*6467efc9SToby Isaac h = *handler; 117*6467efc9SToby Isaac PetscCall(PetscLogHandlerSetType(h, PETSC_LOG_HANDLER_TRACE)); 118*6467efc9SToby Isaac tr = (PetscLogHandler_Trace)(h->data); 119*6467efc9SToby Isaac tr->petsc_tracefile = file; 120*6467efc9SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 121*6467efc9SToby Isaac } 122