178f1b9b4SToby Isaac #include <petscviewer.h> 278f1b9b4SToby Isaac #include <petscdevice.h> 378f1b9b4SToby Isaac #include <petsc/private/logimpl.h> /*I "petscsys.h" I*/ 478f1b9b4SToby Isaac #include <petsc/private/loghandlerimpl.h> 578f1b9b4SToby Isaac #include <petsc/private/deviceimpl.h> 678f1b9b4SToby Isaac #include <petscconfiginfo.h> 778f1b9b4SToby Isaac #include <petscmachineinfo.h> 878f1b9b4SToby Isaac #include "logdefault.h" 978f1b9b4SToby Isaac 1078f1b9b4SToby Isaac static PetscErrorCode PetscEventPerfInfoInit(PetscEventPerfInfo *eventInfo) 1178f1b9b4SToby Isaac { 1278f1b9b4SToby Isaac PetscFunctionBegin; 1378f1b9b4SToby Isaac PetscCall(PetscMemzero(eventInfo, sizeof(*eventInfo))); 1478f1b9b4SToby Isaac eventInfo->visible = PETSC_TRUE; 1578f1b9b4SToby Isaac eventInfo->id = -1; 1678f1b9b4SToby Isaac eventInfo->dof[0] = -1.0; 1778f1b9b4SToby Isaac eventInfo->dof[1] = -1.0; 1878f1b9b4SToby Isaac eventInfo->dof[2] = -1.0; 1978f1b9b4SToby Isaac eventInfo->dof[3] = -1.0; 2078f1b9b4SToby Isaac eventInfo->dof[4] = -1.0; 2178f1b9b4SToby Isaac eventInfo->dof[5] = -1.0; 2278f1b9b4SToby Isaac eventInfo->dof[6] = -1.0; 2378f1b9b4SToby Isaac eventInfo->dof[7] = -1.0; 2478f1b9b4SToby Isaac eventInfo->errors[0] = -1.0; 2578f1b9b4SToby Isaac eventInfo->errors[1] = -1.0; 2678f1b9b4SToby Isaac eventInfo->errors[2] = -1.0; 2778f1b9b4SToby Isaac eventInfo->errors[3] = -1.0; 2878f1b9b4SToby Isaac eventInfo->errors[4] = -1.0; 2978f1b9b4SToby Isaac eventInfo->errors[5] = -1.0; 3078f1b9b4SToby Isaac eventInfo->errors[6] = -1.0; 3178f1b9b4SToby Isaac eventInfo->errors[7] = -1.0; 3278f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 3378f1b9b4SToby Isaac } 3478f1b9b4SToby Isaac 3578f1b9b4SToby Isaac static PetscErrorCode PetscEventPerfInfoTic_Internal(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event, PetscBool resume) 3678f1b9b4SToby Isaac { 3778f1b9b4SToby Isaac PetscFunctionBegin; 3878f1b9b4SToby Isaac if (resume) { 3978f1b9b4SToby Isaac eventInfo->timeTmp -= time; 4078f1b9b4SToby Isaac eventInfo->flopsTmp -= petsc_TotalFlops_th; 4178f1b9b4SToby Isaac } else { 4278f1b9b4SToby Isaac eventInfo->timeTmp = -time; 4378f1b9b4SToby Isaac eventInfo->flopsTmp = -petsc_TotalFlops_th; 4478f1b9b4SToby Isaac } 4578f1b9b4SToby Isaac eventInfo->numMessages -= petsc_irecv_ct_th + petsc_isend_ct_th + petsc_recv_ct_th + petsc_send_ct_th; 4678f1b9b4SToby Isaac eventInfo->messageLength -= petsc_irecv_len_th + petsc_isend_len_th + petsc_recv_len_th + petsc_send_len_th; 4778f1b9b4SToby Isaac eventInfo->numReductions -= petsc_allreduce_ct_th + petsc_gather_ct_th + petsc_scatter_ct_th; 4878f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 4978f1b9b4SToby Isaac eventInfo->CpuToGpuCount -= petsc_ctog_ct_th; 5078f1b9b4SToby Isaac eventInfo->GpuToCpuCount -= petsc_gtoc_ct_th; 5178f1b9b4SToby Isaac eventInfo->CpuToGpuSize -= petsc_ctog_sz_th; 5278f1b9b4SToby Isaac eventInfo->GpuToCpuSize -= petsc_gtoc_sz_th; 5378f1b9b4SToby Isaac eventInfo->GpuFlops -= petsc_gflops_th; 5478f1b9b4SToby Isaac eventInfo->GpuTime -= petsc_gtime; 5578f1b9b4SToby Isaac #endif 5678f1b9b4SToby Isaac if (logMemory) { 5778f1b9b4SToby Isaac PetscLogDouble usage; 5878f1b9b4SToby Isaac PetscCall(PetscMemoryGetCurrentUsage(&usage)); 5978f1b9b4SToby Isaac eventInfo->memIncrease -= usage; 6078f1b9b4SToby Isaac PetscCall(PetscMallocGetCurrentUsage(&usage)); 6178f1b9b4SToby Isaac eventInfo->mallocSpace -= usage; 6278f1b9b4SToby Isaac PetscCall(PetscMallocGetMaximumUsage(&usage)); 6378f1b9b4SToby Isaac eventInfo->mallocIncrease -= usage; 6478f1b9b4SToby Isaac PetscCall(PetscMallocPushMaximumUsage(event)); 6578f1b9b4SToby Isaac } 6678f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 6778f1b9b4SToby Isaac } 6878f1b9b4SToby Isaac 6978f1b9b4SToby Isaac static PetscErrorCode PetscEventPerfInfoTic(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event) 7078f1b9b4SToby Isaac { 7178f1b9b4SToby Isaac PetscFunctionBegin; 7278f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoTic_Internal(eventInfo, time, logMemory, event, PETSC_FALSE)); 7378f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 7478f1b9b4SToby Isaac } 7578f1b9b4SToby Isaac 7678f1b9b4SToby Isaac static PetscErrorCode PetscEventPerfInfoResume(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event) 7778f1b9b4SToby Isaac { 7878f1b9b4SToby Isaac PetscFunctionBegin; 7978f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoTic_Internal(eventInfo, time, logMemory, event, PETSC_TRUE)); 8078f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 8178f1b9b4SToby Isaac } 8278f1b9b4SToby Isaac 8378f1b9b4SToby Isaac static PetscErrorCode PetscEventPerfInfoToc_Internal(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event, PetscBool pause) 8478f1b9b4SToby Isaac { 8578f1b9b4SToby Isaac PetscFunctionBegin; 8678f1b9b4SToby Isaac eventInfo->timeTmp += time; 8778f1b9b4SToby Isaac eventInfo->flopsTmp += petsc_TotalFlops_th; 8878f1b9b4SToby Isaac if (!pause) { 8978f1b9b4SToby Isaac eventInfo->time += eventInfo->timeTmp; 9078f1b9b4SToby Isaac eventInfo->time2 += eventInfo->timeTmp * eventInfo->timeTmp; 9178f1b9b4SToby Isaac eventInfo->flops += eventInfo->flopsTmp; 9278f1b9b4SToby Isaac eventInfo->flops2 += eventInfo->flopsTmp * eventInfo->flopsTmp; 9378f1b9b4SToby Isaac } 9478f1b9b4SToby Isaac eventInfo->numMessages += petsc_irecv_ct_th + petsc_isend_ct_th + petsc_recv_ct_th + petsc_send_ct_th; 9578f1b9b4SToby Isaac eventInfo->messageLength += petsc_irecv_len_th + petsc_isend_len_th + petsc_recv_len + petsc_send_len_th; 9678f1b9b4SToby Isaac eventInfo->numReductions += petsc_allreduce_ct_th + petsc_gather_ct_th + petsc_scatter_ct_th; 9778f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 9878f1b9b4SToby Isaac eventInfo->CpuToGpuCount += petsc_ctog_ct_th; 9978f1b9b4SToby Isaac eventInfo->GpuToCpuCount += petsc_gtoc_ct_th; 10078f1b9b4SToby Isaac eventInfo->CpuToGpuSize += petsc_ctog_sz_th; 10178f1b9b4SToby Isaac eventInfo->GpuToCpuSize += petsc_gtoc_sz_th; 10278f1b9b4SToby Isaac eventInfo->GpuFlops += petsc_gflops_th; 10378f1b9b4SToby Isaac eventInfo->GpuTime += petsc_gtime; 10478f1b9b4SToby Isaac #endif 10578f1b9b4SToby Isaac if (logMemory) { 10678f1b9b4SToby Isaac PetscLogDouble usage, musage; 10778f1b9b4SToby Isaac PetscCall(PetscMemoryGetCurrentUsage(&usage)); /* the comments below match the column labels printed in PetscLogView_Default() */ 10878f1b9b4SToby Isaac eventInfo->memIncrease += usage; /* RMI */ 10978f1b9b4SToby Isaac PetscCall(PetscMallocGetCurrentUsage(&usage)); 11078f1b9b4SToby Isaac eventInfo->mallocSpace += usage; /* Malloc */ 11178f1b9b4SToby Isaac PetscCall(PetscMallocPopMaximumUsage((int)event, &musage)); 11278f1b9b4SToby Isaac eventInfo->mallocIncreaseEvent = PetscMax(musage - usage, eventInfo->mallocIncreaseEvent); /* EMalloc */ 11378f1b9b4SToby Isaac PetscCall(PetscMallocGetMaximumUsage(&usage)); 11478f1b9b4SToby Isaac eventInfo->mallocIncrease += usage; /* MMalloc */ 11578f1b9b4SToby Isaac } 11678f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 11778f1b9b4SToby Isaac } 11878f1b9b4SToby Isaac 11978f1b9b4SToby Isaac static PetscErrorCode PetscEventPerfInfoToc(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event) 12078f1b9b4SToby Isaac { 12178f1b9b4SToby Isaac PetscFunctionBegin; 12278f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoToc_Internal(eventInfo, time, logMemory, event, PETSC_FALSE)); 12378f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 12478f1b9b4SToby Isaac } 12578f1b9b4SToby Isaac 12678f1b9b4SToby Isaac static PetscErrorCode PetscEventPerfInfoPause(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event) 12778f1b9b4SToby Isaac { 12878f1b9b4SToby Isaac PetscFunctionBegin; 12978f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoToc_Internal(eventInfo, time, logMemory, event, PETSC_FALSE)); 13078f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 13178f1b9b4SToby Isaac } 13278f1b9b4SToby Isaac 13378f1b9b4SToby Isaac static PetscErrorCode PetscEventPerfInfoAdd_Internal(const PetscEventPerfInfo *eventInfo, PetscEventPerfInfo *outInfo) 13478f1b9b4SToby Isaac { 13578f1b9b4SToby Isaac PetscFunctionBegin; 13678f1b9b4SToby Isaac outInfo->count += eventInfo->count; 13778f1b9b4SToby Isaac outInfo->time += eventInfo->time; 13878f1b9b4SToby Isaac outInfo->time2 += eventInfo->time2; 13978f1b9b4SToby Isaac outInfo->flops += eventInfo->flops; 14078f1b9b4SToby Isaac outInfo->flops2 += eventInfo->flops2; 14178f1b9b4SToby Isaac outInfo->numMessages += eventInfo->numMessages; 14278f1b9b4SToby Isaac outInfo->messageLength += eventInfo->messageLength; 14378f1b9b4SToby Isaac outInfo->numReductions += eventInfo->numReductions; 14478f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 14578f1b9b4SToby Isaac outInfo->CpuToGpuCount += eventInfo->CpuToGpuCount; 14678f1b9b4SToby Isaac outInfo->GpuToCpuCount += eventInfo->GpuToCpuCount; 14778f1b9b4SToby Isaac outInfo->CpuToGpuSize += eventInfo->CpuToGpuSize; 14878f1b9b4SToby Isaac outInfo->GpuToCpuSize += eventInfo->GpuToCpuSize; 14978f1b9b4SToby Isaac outInfo->GpuFlops += eventInfo->GpuFlops; 15078f1b9b4SToby Isaac outInfo->GpuTime += eventInfo->GpuTime; 15178f1b9b4SToby Isaac #endif 15278f1b9b4SToby Isaac outInfo->memIncrease += eventInfo->memIncrease; 15378f1b9b4SToby Isaac outInfo->mallocSpace += eventInfo->mallocSpace; 15478f1b9b4SToby Isaac outInfo->mallocIncreaseEvent += eventInfo->mallocIncreaseEvent; 15578f1b9b4SToby Isaac outInfo->mallocIncrease += eventInfo->mallocIncrease; 15678f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 15778f1b9b4SToby Isaac } 15878f1b9b4SToby Isaac 15978f1b9b4SToby Isaac PETSC_LOG_RESIZABLE_ARRAY(EventPerfArray, PetscEventPerfInfo, PetscLogEvent, PetscEventPerfInfoInit, NULL, NULL) 16078f1b9b4SToby Isaac 16178f1b9b4SToby Isaac /* --- PetscClassPerf --- */ 16278f1b9b4SToby Isaac 16378f1b9b4SToby Isaac typedef struct { 16478f1b9b4SToby Isaac int creations; /* The number of objects of this class created */ 16578f1b9b4SToby Isaac int destructions; /* The number of objects of this class destroyed */ 16678f1b9b4SToby Isaac PetscLogDouble mem; /* The total memory allocated by objects of this class; this is completely wrong and should possibly be removed */ 16778f1b9b4SToby Isaac PetscLogDouble descMem; /* The total memory allocated by descendents of these objects; this is completely wrong and should possibly be removed */ 16878f1b9b4SToby Isaac } PetscClassPerf; 16978f1b9b4SToby Isaac 17078f1b9b4SToby Isaac static PetscErrorCode PetscClassPerfInit(PetscClassPerf *classInfo) 17178f1b9b4SToby Isaac { 17278f1b9b4SToby Isaac PetscFunctionBegin; 17378f1b9b4SToby Isaac PetscCall(PetscMemzero(classInfo, sizeof(*classInfo))); 17478f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 17578f1b9b4SToby Isaac } 17678f1b9b4SToby Isaac 17778f1b9b4SToby Isaac PETSC_LOG_RESIZABLE_ARRAY(ClassPerfArray, PetscClassPerf, PetscLogClass, PetscClassPerfInit, NULL, NULL) 17878f1b9b4SToby Isaac 17978f1b9b4SToby Isaac /* --- PetscStagePerf --- */ 18078f1b9b4SToby Isaac 18178f1b9b4SToby Isaac typedef struct _PetscStagePerf { 18278f1b9b4SToby Isaac PetscBool used; /* The stage was pushed on this processor */ 18378f1b9b4SToby Isaac PetscEventPerfInfo perfInfo; /* The stage performance information */ 18478f1b9b4SToby Isaac PetscLogEventPerfArray eventLog; /* The event information for this stage */ 18578f1b9b4SToby Isaac PetscLogClassPerfArray classLog; /* The class information for this stage */ 18678f1b9b4SToby Isaac } PetscStagePerf; 18778f1b9b4SToby Isaac 18878f1b9b4SToby Isaac static PetscErrorCode PetscStageInfoInit(PetscStagePerf *stageInfo) 18978f1b9b4SToby Isaac { 19078f1b9b4SToby Isaac PetscFunctionBegin; 19178f1b9b4SToby Isaac PetscCall(PetscMemzero(stageInfo, sizeof(*stageInfo))); 19278f1b9b4SToby Isaac PetscCall(PetscLogEventPerfArrayCreate(128, &(stageInfo->eventLog))); 19378f1b9b4SToby Isaac PetscCall(PetscLogClassPerfArrayCreate(128, &(stageInfo->classLog))); 19478f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoInit(&stageInfo->perfInfo)); 19578f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 19678f1b9b4SToby Isaac } 19778f1b9b4SToby Isaac 19878f1b9b4SToby Isaac static PetscErrorCode PetscStageInfoReset(PetscStagePerf *stageInfo) 19978f1b9b4SToby Isaac { 20078f1b9b4SToby Isaac PetscFunctionBegin; 20178f1b9b4SToby Isaac PetscCall(PetscLogEventPerfArrayDestroy(&(stageInfo->eventLog))); 20278f1b9b4SToby Isaac PetscCall(PetscLogClassPerfArrayDestroy(&(stageInfo->classLog))); 20378f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 20478f1b9b4SToby Isaac } 20578f1b9b4SToby Isaac 20678f1b9b4SToby Isaac PETSC_LOG_RESIZABLE_ARRAY(StageInfoArray, PetscStagePerf, PetscLogStage, PetscStageInfoInit, PetscStageInfoReset, NULL) 20778f1b9b4SToby Isaac 20878f1b9b4SToby Isaac /* --- Action --- */ 20978f1b9b4SToby Isaac 21078f1b9b4SToby Isaac /* The structure for action logging */ 21178f1b9b4SToby Isaac typedef enum { 21278f1b9b4SToby Isaac PETSC_LOG_ACTION_CREATE, 21378f1b9b4SToby Isaac PETSC_LOG_ACTION_DESTROY, 21478f1b9b4SToby Isaac PETSC_LOG_ACTION_BEGIN, 21578f1b9b4SToby Isaac PETSC_LOG_ACTION_END, 21678f1b9b4SToby Isaac } PetscLogActionType; 21778f1b9b4SToby Isaac 218b665b14eSToby Isaac typedef struct _Action { 219b665b14eSToby Isaac PetscLogActionType action; /* The type of execution */ 220b665b14eSToby Isaac PetscLogEvent event; /* The event number */ 221b665b14eSToby Isaac PetscClassId classid; /* The event class id */ 222b665b14eSToby Isaac PetscLogDouble time; /* The time of occurrence */ 223b665b14eSToby Isaac PetscLogDouble flops; /* The cumulative flops */ 224b665b14eSToby Isaac PetscLogDouble mem; /* The current memory usage */ 225b665b14eSToby Isaac PetscLogDouble maxmem; /* The maximum memory usage */ 226b665b14eSToby Isaac int id1, id2, id3; /* The ids of associated objects */ 227b665b14eSToby Isaac } Action; 228b665b14eSToby Isaac 22978f1b9b4SToby Isaac PETSC_LOG_RESIZABLE_ARRAY(ActionArray, Action, PetscLogEvent, NULL, NULL, NULL) 23078f1b9b4SToby Isaac 23178f1b9b4SToby Isaac /* --- Object --- */ 23278f1b9b4SToby Isaac 233b665b14eSToby Isaac /* The structure for object logging */ 234b665b14eSToby Isaac typedef struct _Object { 235b665b14eSToby Isaac PetscObject obj; /* The associated PetscObject */ 236b665b14eSToby Isaac int parent; /* The parent id */ 237b665b14eSToby Isaac PetscLogDouble mem; /* The memory associated with the object */ 238b665b14eSToby Isaac char name[64]; /* The object name */ 239b665b14eSToby Isaac char info[64]; /* The information string */ 240b665b14eSToby Isaac } Object; 241b665b14eSToby Isaac 24278f1b9b4SToby Isaac PETSC_LOG_RESIZABLE_ARRAY(ObjectArray, Object, PetscObject, NULL, NULL, NULL) 24378f1b9b4SToby Isaac 24478f1b9b4SToby Isaac /* Map from (threadid,stage,event) to perfInfo data struct */ 24578f1b9b4SToby Isaac #include <petsc/private/hashmapijk.h> 24678f1b9b4SToby Isaac 24778f1b9b4SToby Isaac PETSC_HASH_MAP(HMapEvent, PetscHashIJKKey, PetscEventPerfInfo *, PetscHashIJKKeyHash, PetscHashIJKKeyEqual, NULL) 24878f1b9b4SToby Isaac 24978f1b9b4SToby Isaac typedef struct _n_PetscLogHandler_Default *PetscLogHandler_Default; 25078f1b9b4SToby Isaac struct _n_PetscLogHandler_Default { 25178f1b9b4SToby Isaac PetscLogStageInfoArray stages; 25278f1b9b4SToby Isaac PetscSpinlock lock; 25378f1b9b4SToby Isaac PetscLogActionArray petsc_actions; 25478f1b9b4SToby Isaac PetscLogObjectArray petsc_objects; 25578f1b9b4SToby Isaac PetscBool petsc_logActions; 25678f1b9b4SToby Isaac PetscBool petsc_logObjects; 25778f1b9b4SToby Isaac int petsc_numObjectsCreated; 25878f1b9b4SToby Isaac int petsc_numObjectsDestroyed; 25978f1b9b4SToby Isaac PetscHMapEvent eventInfoMap_th; 26078f1b9b4SToby Isaac int pause_depth; 261461318ebSToby Isaac PetscBool use_threadsafe; 26278f1b9b4SToby Isaac }; 26378f1b9b4SToby Isaac 26478f1b9b4SToby Isaac /* --- PetscLogHandler_Default --- */ 26578f1b9b4SToby Isaac 26678f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerContextCreate_Default(PetscLogHandler_Default *def_p) 26778f1b9b4SToby Isaac { 26878f1b9b4SToby Isaac PetscLogHandler_Default def; 26978f1b9b4SToby Isaac 27078f1b9b4SToby Isaac PetscFunctionBegin; 27178f1b9b4SToby Isaac PetscCall(PetscNew(def_p)); 27278f1b9b4SToby Isaac def = *def_p; 27378f1b9b4SToby Isaac PetscCall(PetscLogStageInfoArrayCreate(8, &(def->stages))); 27478f1b9b4SToby Isaac PetscCall(PetscLogActionArrayCreate(64, &(def->petsc_actions))); 27578f1b9b4SToby Isaac PetscCall(PetscLogObjectArrayCreate(64, &(def->petsc_objects))); 27678f1b9b4SToby Isaac 27778f1b9b4SToby Isaac PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_include_actions", &def->petsc_logActions, NULL)); 27878f1b9b4SToby Isaac PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_include_objects", &def->petsc_logObjects, NULL)); 279461318ebSToby Isaac PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_handler_default_use_threadsafe_events", &def->use_threadsafe, NULL)); 280461318ebSToby Isaac if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { PetscCall(PetscHMapEventCreate(&def->eventInfoMap_th)); } 28178f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 28278f1b9b4SToby Isaac } 28378f1b9b4SToby Isaac 28478f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerDestroy_Default(PetscLogHandler h) 28578f1b9b4SToby Isaac { 28678f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 28778f1b9b4SToby Isaac 28878f1b9b4SToby Isaac PetscFunctionBegin; 28978f1b9b4SToby Isaac PetscCall(PetscLogStageInfoArrayDestroy(&def->stages)); 29078f1b9b4SToby Isaac PetscCall(PetscLogActionArrayDestroy(&def->petsc_actions)); 29178f1b9b4SToby Isaac PetscCall(PetscLogObjectArrayDestroy(&def->petsc_objects)); 29278f1b9b4SToby Isaac if (def->eventInfoMap_th) { 29378f1b9b4SToby Isaac PetscEventPerfInfo **array; 29478f1b9b4SToby Isaac PetscInt n, off = 0; 29578f1b9b4SToby Isaac 29678f1b9b4SToby Isaac PetscCall(PetscHMapEventGetSize(def->eventInfoMap_th, &n)); 29778f1b9b4SToby Isaac PetscCall(PetscMalloc1(n, &array)); 29878f1b9b4SToby Isaac PetscCall(PetscHMapEventGetVals(def->eventInfoMap_th, &off, array)); 29978f1b9b4SToby Isaac for (PetscInt i = 0; i < n; i++) PetscCall(PetscFree(array[i])); 30078f1b9b4SToby Isaac PetscCall(PetscFree(array)); 30178f1b9b4SToby Isaac PetscCall(PetscHMapEventDestroy(&def->eventInfoMap_th)); 30278f1b9b4SToby Isaac } 30378f1b9b4SToby Isaac PetscCall(PetscFree(def)); 30478f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 30578f1b9b4SToby Isaac } 30678f1b9b4SToby Isaac 30778f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerDefaultGetStageInfo(PetscLogHandler handler, PetscLogStage stage, PetscStagePerf **stage_info_p) 30878f1b9b4SToby Isaac { 30978f1b9b4SToby Isaac PetscStagePerf *stage_info = NULL; 31078f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 31178f1b9b4SToby Isaac 31278f1b9b4SToby Isaac PetscFunctionBegin; 31378f1b9b4SToby Isaac PetscCall(PetscLogStageInfoArrayResize(def->stages, stage + 1)); 31478f1b9b4SToby Isaac PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info)); 31578f1b9b4SToby Isaac *stage_info_p = stage_info; 31678f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 31778f1b9b4SToby Isaac } 31878f1b9b4SToby Isaac 31978f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultGetEventPerfInfo(PetscLogHandler handler, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **event_info_p) 32078f1b9b4SToby Isaac { 32178f1b9b4SToby Isaac PetscEventPerfInfo *event_info = NULL; 32278f1b9b4SToby Isaac PetscStagePerf *stage_info = NULL; 32378f1b9b4SToby Isaac PetscLogEventPerfArray event_log; 32478f1b9b4SToby Isaac 32578f1b9b4SToby Isaac PetscFunctionBegin; 32678f1b9b4SToby Isaac if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(handler->state, &stage)); 32778f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info)); 32878f1b9b4SToby Isaac event_log = stage_info->eventLog; 32978f1b9b4SToby Isaac PetscCall(PetscLogEventPerfArrayResize(event_log, event + 1)); 33078f1b9b4SToby Isaac PetscCall(PetscLogEventPerfArrayGetRef(event_log, event, &event_info)); 33178f1b9b4SToby Isaac event_info->id = event; 33278f1b9b4SToby Isaac *event_info_p = event_info; 33378f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 33478f1b9b4SToby Isaac } 33578f1b9b4SToby Isaac 33678f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerDefaultGetClassPerf(PetscLogHandler handler, PetscLogStage stage, PetscLogClass clss, PetscClassPerf **class_info) 33778f1b9b4SToby Isaac { 33878f1b9b4SToby Isaac PetscLogClassPerfArray class_log; 33978f1b9b4SToby Isaac PetscStagePerf *stage_info; 34078f1b9b4SToby Isaac 34178f1b9b4SToby Isaac PetscFunctionBegin; 34278f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info)); 34378f1b9b4SToby Isaac class_log = stage_info->classLog; 34478f1b9b4SToby Isaac PetscCall(PetscLogClassPerfArrayResize(class_log, clss + 1)); 34578f1b9b4SToby Isaac PetscCall(PetscLogClassPerfArrayGetRef(class_log, clss, class_info)); 34678f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 34778f1b9b4SToby Isaac } 34878f1b9b4SToby Isaac 34978f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerObjectCreate_Default(PetscLogHandler h, PetscObject obj) 35078f1b9b4SToby Isaac { 35178f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 35278f1b9b4SToby Isaac PetscLogState state; 35378f1b9b4SToby Isaac PetscLogStage stage; 35478f1b9b4SToby Isaac PetscClassPerf *classInfo; 35578f1b9b4SToby Isaac int oclass = 0; 35678f1b9b4SToby Isaac 35778f1b9b4SToby Isaac PetscFunctionBegin; 35878f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state)); 35978f1b9b4SToby Isaac PetscCall(PetscSpinlockLock(&def->lock)); 36078f1b9b4SToby Isaac /* Record stage info */ 36178f1b9b4SToby Isaac PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 36278f1b9b4SToby Isaac PetscCall(PetscLogStateGetClassFromClassId(state, obj->classid, &oclass)); 36378f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetClassPerf(h, stage, oclass, &classInfo)); 36478f1b9b4SToby Isaac classInfo->creations++; 36578f1b9b4SToby Isaac /* Record the creation action */ 36678f1b9b4SToby Isaac if (def->petsc_logActions) { 36778f1b9b4SToby Isaac Action new_action; 36878f1b9b4SToby Isaac 36978f1b9b4SToby Isaac PetscCall(PetscTime(&new_action.time)); 37078f1b9b4SToby Isaac new_action.time -= petsc_BaseTime; 37178f1b9b4SToby Isaac new_action.action = PETSC_LOG_ACTION_CREATE; 37278f1b9b4SToby Isaac new_action.event = -1; 37378f1b9b4SToby Isaac new_action.classid = obj->classid; 37478f1b9b4SToby Isaac new_action.id1 = obj->id; 37578f1b9b4SToby Isaac new_action.id2 = -1; 37678f1b9b4SToby Isaac new_action.id3 = -1; 37778f1b9b4SToby Isaac new_action.flops = petsc_TotalFlops; 37878f1b9b4SToby Isaac PetscCall(PetscMallocGetCurrentUsage(&new_action.mem)); 37978f1b9b4SToby Isaac PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem)); 38078f1b9b4SToby Isaac PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action)); 38178f1b9b4SToby Isaac } 38278f1b9b4SToby Isaac /* We don't just use obj->id to count all objects that are created 38378f1b9b4SToby Isaac because PetscLogHandlers are objects and PetscLogObjectDestroy() will not 38478f1b9b4SToby Isaac be called for them: the number of objects created and destroyed as counted 38578f1b9b4SToby Isaac here and below would have an imbalance */ 38678f1b9b4SToby Isaac def->petsc_numObjectsCreated++; 38778f1b9b4SToby Isaac /* Record the object */ 38878f1b9b4SToby Isaac if (def->petsc_logObjects) { 38978f1b9b4SToby Isaac Object new_object; 39078f1b9b4SToby Isaac 39178f1b9b4SToby Isaac new_object.parent = -1; 39278f1b9b4SToby Isaac new_object.obj = obj; 39378f1b9b4SToby Isaac new_object.mem = 0; 39478f1b9b4SToby Isaac 39578f1b9b4SToby Isaac PetscCall(PetscMemzero(new_object.name, sizeof(new_object.name))); 39678f1b9b4SToby Isaac PetscCall(PetscMemzero(new_object.info, sizeof(new_object.info))); 39778f1b9b4SToby Isaac PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1"); 39878f1b9b4SToby Isaac PetscCall(PetscLogObjectArrayResize(def->petsc_objects, obj->id)); 39978f1b9b4SToby Isaac PetscCall(PetscLogObjectArraySet(def->petsc_objects, obj->id - 1, new_object)); 40078f1b9b4SToby Isaac } 40178f1b9b4SToby Isaac PetscCall(PetscSpinlockUnlock(&def->lock)); 40278f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 40378f1b9b4SToby Isaac } 40478f1b9b4SToby Isaac 40578f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerObjectDestroy_Default(PetscLogHandler h, PetscObject obj) 40678f1b9b4SToby Isaac { 40778f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 40878f1b9b4SToby Isaac PetscLogState state; 40978f1b9b4SToby Isaac PetscLogStage stage; 41078f1b9b4SToby Isaac PetscClassPerf *classInfo; 41178f1b9b4SToby Isaac int oclass = 0; 41278f1b9b4SToby Isaac 41378f1b9b4SToby Isaac PetscFunctionBegin; 41478f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state)); 41578f1b9b4SToby Isaac /* Record stage info */ 41678f1b9b4SToby Isaac PetscCall(PetscSpinlockLock(&def->lock)); 41778f1b9b4SToby Isaac PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 41878f1b9b4SToby Isaac if (stage >= 0) { 41978f1b9b4SToby Isaac /* stage < 0 can happen if the log summary is output before some things are destroyed */ 42078f1b9b4SToby Isaac PetscCall(PetscLogStateGetClassFromClassId(state, obj->classid, &oclass)); 42178f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetClassPerf(h, stage, oclass, &classInfo)); 42278f1b9b4SToby Isaac classInfo->destructions++; 42378f1b9b4SToby Isaac } 42478f1b9b4SToby Isaac /* Cannot Credit all ancestors with your memory because they may have already been destroyed*/ 42578f1b9b4SToby Isaac def->petsc_numObjectsDestroyed++; 42678f1b9b4SToby Isaac /* Dynamically enlarge logging structures */ 42778f1b9b4SToby Isaac /* Record the destruction action */ 42878f1b9b4SToby Isaac if (def->petsc_logActions) { 42978f1b9b4SToby Isaac Action new_action; 43078f1b9b4SToby Isaac 43178f1b9b4SToby Isaac PetscCall(PetscTime(&new_action.time)); 43278f1b9b4SToby Isaac new_action.time -= petsc_BaseTime; 43378f1b9b4SToby Isaac new_action.event = -1; 43478f1b9b4SToby Isaac new_action.action = PETSC_LOG_ACTION_DESTROY; 43578f1b9b4SToby Isaac new_action.classid = obj->classid; 43678f1b9b4SToby Isaac new_action.id1 = obj->id; 43778f1b9b4SToby Isaac new_action.id2 = -1; 43878f1b9b4SToby Isaac new_action.id3 = -1; 43978f1b9b4SToby Isaac new_action.flops = petsc_TotalFlops; 44078f1b9b4SToby Isaac PetscCall(PetscMallocGetCurrentUsage(&new_action.mem)); 44178f1b9b4SToby Isaac PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem)); 44278f1b9b4SToby Isaac PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action)); 44378f1b9b4SToby Isaac } 44478f1b9b4SToby Isaac if (def->petsc_logObjects) { 44578f1b9b4SToby Isaac Object *obj_entry = NULL; 44678f1b9b4SToby Isaac 44778f1b9b4SToby Isaac PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1"); 44878f1b9b4SToby Isaac PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, obj->id - 1, &obj_entry)); 44978f1b9b4SToby Isaac if (obj->name) PetscCall(PetscStrncpy(obj_entry->name, obj->name, 64)); 45078f1b9b4SToby Isaac obj_entry->obj = NULL; 45178f1b9b4SToby Isaac } 45278f1b9b4SToby Isaac PetscCall(PetscSpinlockUnlock(&def->lock)); 45378f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 45478f1b9b4SToby Isaac } 45578f1b9b4SToby Isaac 45678f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerEventSync_Default(PetscLogHandler h, PetscLogEvent event, MPI_Comm comm) 45778f1b9b4SToby Isaac { 45878f1b9b4SToby Isaac PetscLogState state; 45978f1b9b4SToby Isaac PetscLogEventInfo event_info; 46078f1b9b4SToby Isaac PetscEventPerfInfo *event_perf_info; 46178f1b9b4SToby Isaac int stage; 46278f1b9b4SToby Isaac PetscLogDouble time = 0.0; 46378f1b9b4SToby Isaac 46478f1b9b4SToby Isaac PetscFunctionBegin; 46578f1b9b4SToby Isaac if (!(PetscLogSyncOn) || comm == MPI_COMM_NULL) PetscFunctionReturn(PETSC_SUCCESS); 46678f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state)); 46778f1b9b4SToby Isaac PetscCall(PetscLogStateEventGetInfo(state, event, &event_info)); 46878f1b9b4SToby Isaac if (!event_info.collective) PetscFunctionReturn(PETSC_SUCCESS); 46978f1b9b4SToby Isaac PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 47078f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info)); 47178f1b9b4SToby Isaac if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS); 47278f1b9b4SToby Isaac 47378f1b9b4SToby Isaac PetscCall(PetscTimeSubtract(&time)); 47478f1b9b4SToby Isaac PetscCallMPI(MPI_Barrier(comm)); 47578f1b9b4SToby Isaac PetscCall(PetscTimeAdd(&time)); 47678f1b9b4SToby Isaac event_perf_info->syncTime += time; 47778f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 47878f1b9b4SToby Isaac } 47978f1b9b4SToby Isaac 48078f1b9b4SToby Isaac static PetscErrorCode PetscLogGetStageEventPerfInfo_threaded(PetscLogHandler_Default def, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **eventInfo) 48178f1b9b4SToby Isaac { 48278f1b9b4SToby Isaac PetscEventPerfInfo *leventInfo = NULL; 48378f1b9b4SToby Isaac PetscHashIJKKey key; 48478f1b9b4SToby Isaac 48578f1b9b4SToby Isaac PetscFunctionBegin; 48678f1b9b4SToby Isaac 48778f1b9b4SToby Isaac #if PetscDefined(HAVE_THREADSAFETY) 48878f1b9b4SToby Isaac key.i = PetscLogGetTid(); 48978f1b9b4SToby Isaac #else 49078f1b9b4SToby Isaac key.i = 0; 49178f1b9b4SToby Isaac #endif 49278f1b9b4SToby Isaac key.j = stage; 49378f1b9b4SToby Isaac key.k = event; 49478f1b9b4SToby Isaac PetscCall(PetscSpinlockLock(&def->lock)); 49578f1b9b4SToby Isaac PetscCall(PetscHMapEventGet(def->eventInfoMap_th, key, &leventInfo)); 49678f1b9b4SToby Isaac if (!leventInfo) { 49778f1b9b4SToby Isaac PetscCall(PetscNew(&leventInfo)); 49878f1b9b4SToby Isaac leventInfo->id = event; 49978f1b9b4SToby Isaac PetscCall(PetscHMapEventSet(def->eventInfoMap_th, key, leventInfo)); 50078f1b9b4SToby Isaac } 50178f1b9b4SToby Isaac PetscCall(PetscSpinlockUnlock(&def->lock)); 50278f1b9b4SToby Isaac *eventInfo = leventInfo; 50378f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 50478f1b9b4SToby Isaac } 50578f1b9b4SToby Isaac 50678f1b9b4SToby Isaac #if defined(PETSC_HAVE_CUDA) 50778f1b9b4SToby Isaac #include <nvToolsExt.h> 50878f1b9b4SToby Isaac #endif 50978f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerEventBegin_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4) 51078f1b9b4SToby Isaac { 51178f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 51278f1b9b4SToby Isaac PetscEventPerfInfo *event_perf_info = NULL; 51378f1b9b4SToby Isaac PetscLogEventInfo event_info; 51478f1b9b4SToby Isaac PetscLogDouble time; 51578f1b9b4SToby Isaac PetscLogState state; 51678f1b9b4SToby Isaac PetscLogStage stage; 51778f1b9b4SToby Isaac 51878f1b9b4SToby Isaac PetscFunctionBegin; 51978f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state)); 52078f1b9b4SToby Isaac /* Synchronization */ 52178f1b9b4SToby Isaac PetscCall(PetscLogHandlerEventSync_Default(h, event, PetscObjectComm(o1))); 52278f1b9b4SToby Isaac PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 52378f1b9b4SToby Isaac if (def->pause_depth > 0) stage = 0; // in pause-mode, all events run on the main stage 524461318ebSToby Isaac if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { 52578f1b9b4SToby Isaac PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info)); 52678f1b9b4SToby Isaac if (event_perf_info->depth == 0) { PetscCall(PetscEventPerfInfoInit(event_perf_info)); } 52778f1b9b4SToby Isaac } else { 52878f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info)); 52978f1b9b4SToby Isaac } 53078f1b9b4SToby Isaac PetscCheck(event_perf_info->depth >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to begin a paused event, this is not allowed"); 53178f1b9b4SToby Isaac event_perf_info->depth++; 53278f1b9b4SToby Isaac /* Check for double counting */ 53378f1b9b4SToby Isaac if (event_perf_info->depth > 1) PetscFunctionReturn(PETSC_SUCCESS); 53478f1b9b4SToby Isaac #if defined(PETSC_HAVE_CUDA) 53578f1b9b4SToby Isaac if (PetscDeviceInitialized(PETSC_DEVICE_CUDA)) { 53678f1b9b4SToby Isaac PetscLogEventInfo event_reg_info; 53778f1b9b4SToby Isaac 53878f1b9b4SToby Isaac PetscCall(PetscLogStateEventGetInfo(state, event, &event_reg_info)); 53978f1b9b4SToby Isaac nvtxRangePushA(event_reg_info.name); 54078f1b9b4SToby Isaac } 54178f1b9b4SToby Isaac #endif 54278f1b9b4SToby Isaac PetscCall(PetscLogStateEventGetInfo(state, event, &event_info)); 54378f1b9b4SToby Isaac /* Log the performance info */ 54478f1b9b4SToby Isaac event_perf_info->count++; 54578f1b9b4SToby Isaac PetscCall(PetscTime(&time)); 54678f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoTic(event_perf_info, time, PetscLogMemory, (int)event)); 54778f1b9b4SToby Isaac if (def->petsc_logActions) { 54878f1b9b4SToby Isaac PetscLogDouble curTime; 54978f1b9b4SToby Isaac Action new_action; 55078f1b9b4SToby Isaac 55178f1b9b4SToby Isaac PetscCall(PetscTime(&curTime)); 55278f1b9b4SToby Isaac new_action.time = curTime - petsc_BaseTime; 55378f1b9b4SToby Isaac new_action.action = PETSC_LOG_ACTION_BEGIN; 55478f1b9b4SToby Isaac new_action.event = event; 55578f1b9b4SToby Isaac new_action.classid = event_info.classid; 55678f1b9b4SToby Isaac new_action.id1 = o1 ? o1->id : -1; 55778f1b9b4SToby Isaac new_action.id2 = o2 ? o2->id : -1; 55878f1b9b4SToby Isaac new_action.id3 = o3 ? o3->id : -1; 55978f1b9b4SToby Isaac new_action.flops = petsc_TotalFlops; 56078f1b9b4SToby Isaac PetscCall(PetscMallocGetCurrentUsage(&new_action.mem)); 56178f1b9b4SToby Isaac PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem)); 56278f1b9b4SToby Isaac PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action)); 56378f1b9b4SToby Isaac } 56478f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 56578f1b9b4SToby Isaac } 56678f1b9b4SToby Isaac 56778f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerEventEnd_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4) 56878f1b9b4SToby Isaac { 56978f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 57078f1b9b4SToby Isaac PetscEventPerfInfo *event_perf_info = NULL; 57178f1b9b4SToby Isaac PetscLogDouble time; 57278f1b9b4SToby Isaac PetscLogState state; 57378f1b9b4SToby Isaac int stage; 57478f1b9b4SToby Isaac 57578f1b9b4SToby Isaac PetscFunctionBegin; 57678f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state)); 57778f1b9b4SToby Isaac if (def->petsc_logActions) { 57878f1b9b4SToby Isaac PetscLogEventInfo event_info; 57978f1b9b4SToby Isaac PetscLogDouble curTime; 58078f1b9b4SToby Isaac Action new_action; 58178f1b9b4SToby Isaac 58278f1b9b4SToby Isaac PetscCall(PetscLogStateEventGetInfo(state, event, &event_info)); 58378f1b9b4SToby Isaac PetscCall(PetscTime(&curTime)); 58478f1b9b4SToby Isaac new_action.time = curTime - petsc_BaseTime; 58578f1b9b4SToby Isaac new_action.action = PETSC_LOG_ACTION_END; 58678f1b9b4SToby Isaac new_action.event = event; 58778f1b9b4SToby Isaac new_action.classid = event_info.classid; 58878f1b9b4SToby Isaac new_action.id1 = o1 ? o1->id : -1; 58978f1b9b4SToby Isaac new_action.id2 = o2 ? o2->id : -2; 59078f1b9b4SToby Isaac new_action.id3 = o3 ? o3->id : -3; 59178f1b9b4SToby Isaac new_action.flops = petsc_TotalFlops; 59278f1b9b4SToby Isaac PetscCall(PetscMallocGetCurrentUsage(&new_action.mem)); 59378f1b9b4SToby Isaac PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem)); 59478f1b9b4SToby Isaac PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action)); 59578f1b9b4SToby Isaac } 59678f1b9b4SToby Isaac PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 59778f1b9b4SToby Isaac if (def->pause_depth > 0) stage = 0; // all events run on the main stage in pause-mode 598461318ebSToby Isaac if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { 59978f1b9b4SToby Isaac PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info)); 60078f1b9b4SToby Isaac } else { 60178f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info)); 60278f1b9b4SToby Isaac } 60378f1b9b4SToby Isaac PetscCheck(event_perf_info->depth > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to end paused event, not allowed"); 60478f1b9b4SToby Isaac event_perf_info->depth--; 60578f1b9b4SToby Isaac /* Check for double counting */ 60678f1b9b4SToby Isaac if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS); 60778f1b9b4SToby Isaac else PetscCheck(event_perf_info->depth == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Logging event had unbalanced begin/end pairs"); 60878f1b9b4SToby Isaac 60978f1b9b4SToby Isaac /* Log performance info */ 61078f1b9b4SToby Isaac PetscCall(PetscTime(&time)); 61178f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoToc(event_perf_info, time, PetscLogMemory, (int)event)); 612461318ebSToby Isaac if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { 61378f1b9b4SToby Isaac PetscEventPerfInfo *event_perf_info_global; 61478f1b9b4SToby Isaac PetscCall(PetscSpinlockLock(&def->lock)); 61578f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info_global)); 61678f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoAdd_Internal(event_perf_info, event_perf_info_global)); 61778f1b9b4SToby Isaac PetscCall(PetscSpinlockUnlock(&def->lock)); 61878f1b9b4SToby Isaac } 61978f1b9b4SToby Isaac #if defined(PETSC_HAVE_CUDA) 62078f1b9b4SToby Isaac if (PetscDeviceInitialized(PETSC_DEVICE_CUDA)) nvtxRangePop(); 62178f1b9b4SToby Isaac #endif 62278f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 62378f1b9b4SToby Isaac } 62478f1b9b4SToby Isaac 62578f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultDeactivatePush(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event) 62678f1b9b4SToby Isaac { 62778f1b9b4SToby Isaac PetscEventPerfInfo *event_perf_info; 62878f1b9b4SToby Isaac 62978f1b9b4SToby Isaac PetscFunctionBegin; 63078f1b9b4SToby Isaac if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage)); 63178f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info)); 63278f1b9b4SToby Isaac event_perf_info->depth++; 63378f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 63478f1b9b4SToby Isaac } 63578f1b9b4SToby Isaac 63678f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultDeactivatePop(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event) 63778f1b9b4SToby Isaac { 63878f1b9b4SToby Isaac PetscEventPerfInfo *event_perf_info; 63978f1b9b4SToby Isaac 64078f1b9b4SToby Isaac PetscFunctionBegin; 64178f1b9b4SToby Isaac if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage)); 64278f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info)); 64378f1b9b4SToby Isaac event_perf_info->depth--; 64478f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 64578f1b9b4SToby Isaac } 64678f1b9b4SToby Isaac 64778f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultEventsPause(PetscLogHandler h) 64878f1b9b4SToby Isaac { 64978f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 65078f1b9b4SToby Isaac PetscLogDouble time; 65178f1b9b4SToby Isaac PetscInt num_stages; 65278f1b9b4SToby Isaac 65378f1b9b4SToby Isaac PetscFunctionBegin; 65478f1b9b4SToby Isaac if (def->pause_depth++ > 0) PetscFunctionReturn(PETSC_SUCCESS); 65578f1b9b4SToby Isaac PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL)); 65678f1b9b4SToby Isaac PetscCall(PetscTime(&time)); 65778f1b9b4SToby Isaac for (PetscInt stage = 0; stage < num_stages; stage++) { 65878f1b9b4SToby Isaac PetscStagePerf *stage_info = NULL; 65978f1b9b4SToby Isaac PetscInt num_events; 66078f1b9b4SToby Isaac 66178f1b9b4SToby Isaac PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info)); 66278f1b9b4SToby Isaac PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL)); 66378f1b9b4SToby Isaac for (PetscInt event = 0; event < num_events; event++) { 66478f1b9b4SToby Isaac PetscEventPerfInfo *event_info = NULL; 66578f1b9b4SToby Isaac PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info)); 66678f1b9b4SToby Isaac if (event_info->depth > 0) { 66778f1b9b4SToby Isaac event_info->depth *= -1; 66878f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoPause(event_info, time, PetscLogMemory, event)); 66978f1b9b4SToby Isaac } 67078f1b9b4SToby Isaac } 67178f1b9b4SToby Isaac if (stage > 0 && stage_info->perfInfo.depth > 0) { 67278f1b9b4SToby Isaac stage_info->perfInfo.depth *= -1; 67378f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoPause(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2))); 67478f1b9b4SToby Isaac } 67578f1b9b4SToby Isaac } 67678f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 67778f1b9b4SToby Isaac } 67878f1b9b4SToby Isaac 67978f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultEventsResume(PetscLogHandler h) 68078f1b9b4SToby Isaac { 68178f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 68278f1b9b4SToby Isaac PetscLogDouble time; 68378f1b9b4SToby Isaac PetscInt num_stages; 68478f1b9b4SToby Isaac 68578f1b9b4SToby Isaac PetscFunctionBegin; 68678f1b9b4SToby Isaac if (--def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS); 68778f1b9b4SToby Isaac PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL)); 68878f1b9b4SToby Isaac PetscCall(PetscTime(&time)); 68978f1b9b4SToby Isaac for (PetscInt stage = 0; stage < num_stages; stage++) { 69078f1b9b4SToby Isaac PetscStagePerf *stage_info = NULL; 69178f1b9b4SToby Isaac PetscInt num_events; 69278f1b9b4SToby Isaac 69378f1b9b4SToby Isaac PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info)); 69478f1b9b4SToby Isaac PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL)); 69578f1b9b4SToby Isaac for (PetscInt event = 0; event < num_events; event++) { 69678f1b9b4SToby Isaac PetscEventPerfInfo *event_info = NULL; 69778f1b9b4SToby Isaac PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info)); 69878f1b9b4SToby Isaac if (event_info->depth < 0) { 69978f1b9b4SToby Isaac event_info->depth *= -1; 70078f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoResume(event_info, time, PetscLogMemory, event)); 70178f1b9b4SToby Isaac } 70278f1b9b4SToby Isaac } 70378f1b9b4SToby Isaac if (stage > 0 && stage_info->perfInfo.depth < 0) { 70478f1b9b4SToby Isaac stage_info->perfInfo.depth *= -1; 70578f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoResume(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2))); 70678f1b9b4SToby Isaac } 70778f1b9b4SToby Isaac } 70878f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 70978f1b9b4SToby Isaac } 71078f1b9b4SToby Isaac 71178f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerStagePush_Default(PetscLogHandler h, PetscLogStage new_stage) 71278f1b9b4SToby Isaac { 71378f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 71478f1b9b4SToby Isaac PetscLogDouble time; 71578f1b9b4SToby Isaac PetscLogState state; 71678f1b9b4SToby Isaac PetscLogStage current_stage; 71778f1b9b4SToby Isaac PetscStagePerf *new_stage_info; 71878f1b9b4SToby Isaac 71978f1b9b4SToby Isaac PetscFunctionBegin; 72078f1b9b4SToby Isaac if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS); 72178f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state)); 72278f1b9b4SToby Isaac current_stage = state->current_stage; 72378f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(h, new_stage, &new_stage_info)); 72478f1b9b4SToby Isaac PetscCall(PetscTime(&time)); 72578f1b9b4SToby Isaac 72678f1b9b4SToby Isaac /* Record flops/time of previous stage */ 72778f1b9b4SToby Isaac if (current_stage >= 0) { 72878f1b9b4SToby Isaac if (PetscBTLookup(state->active, current_stage)) { 72978f1b9b4SToby Isaac PetscStagePerf *current_stage_info; 73078f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, ¤t_stage_info)); 73178f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoToc(¤t_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2))); 73278f1b9b4SToby Isaac } 73378f1b9b4SToby Isaac } 73478f1b9b4SToby Isaac new_stage_info->used = PETSC_TRUE; 73578f1b9b4SToby Isaac new_stage_info->perfInfo.count++; 73678f1b9b4SToby Isaac new_stage_info->perfInfo.depth++; 73778f1b9b4SToby Isaac /* Subtract current quantities so that we obtain the difference when we pop */ 73878f1b9b4SToby Isaac if (PetscBTLookup(state->active, new_stage)) PetscCall(PetscEventPerfInfoTic(&new_stage_info->perfInfo, time, PetscLogMemory, (int)-(new_stage + 2))); 73978f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 74078f1b9b4SToby Isaac } 74178f1b9b4SToby Isaac 74278f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerStagePop_Default(PetscLogHandler h, PetscLogStage old_stage) 74378f1b9b4SToby Isaac { 74478f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 74578f1b9b4SToby Isaac PetscLogStage current_stage; 74678f1b9b4SToby Isaac PetscStagePerf *old_stage_info; 74778f1b9b4SToby Isaac PetscLogState state; 74878f1b9b4SToby Isaac PetscLogDouble time; 74978f1b9b4SToby Isaac 75078f1b9b4SToby Isaac PetscFunctionBegin; 75178f1b9b4SToby Isaac if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS); 75278f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(h, &state)); 75378f1b9b4SToby Isaac current_stage = state->current_stage; 75478f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(h, old_stage, &old_stage_info)); 75578f1b9b4SToby Isaac PetscCall(PetscTime(&time)); 75678f1b9b4SToby Isaac old_stage_info->perfInfo.depth--; 75778f1b9b4SToby Isaac if (PetscBTLookup(state->active, old_stage)) { PetscCall(PetscEventPerfInfoToc(&old_stage_info->perfInfo, time, PetscLogMemory, (int)-(old_stage + 2))); } 75878f1b9b4SToby Isaac if (current_stage >= 0) { 75978f1b9b4SToby Isaac if (PetscBTLookup(state->active, current_stage)) { 76078f1b9b4SToby Isaac PetscStagePerf *current_stage_info; 76178f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, ¤t_stage_info)); 76278f1b9b4SToby Isaac PetscCall(PetscEventPerfInfoTic(¤t_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2))); 76378f1b9b4SToby Isaac } 76478f1b9b4SToby Isaac } 76578f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 76678f1b9b4SToby Isaac } 76778f1b9b4SToby Isaac 76878f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultStageSetVisible(PetscLogHandler h, PetscLogStage stage, PetscBool is_visible) 76978f1b9b4SToby Isaac { 77078f1b9b4SToby Isaac PetscStagePerf *stage_info; 77178f1b9b4SToby Isaac 77278f1b9b4SToby Isaac PetscFunctionBegin; 77378f1b9b4SToby Isaac if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage)); 77478f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info)); 77578f1b9b4SToby Isaac stage_info->perfInfo.visible = is_visible; 77678f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 77778f1b9b4SToby Isaac } 77878f1b9b4SToby Isaac 77978f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultStageGetVisible(PetscLogHandler h, PetscLogStage stage, PetscBool *is_visible) 78078f1b9b4SToby Isaac { 78178f1b9b4SToby Isaac PetscStagePerf *stage_info; 78278f1b9b4SToby Isaac 78378f1b9b4SToby Isaac PetscFunctionBegin; 78478f1b9b4SToby Isaac if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage)); 78578f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info)); 78678f1b9b4SToby Isaac *is_visible = stage_info->perfInfo.visible; 78778f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 78878f1b9b4SToby Isaac } 78978f1b9b4SToby Isaac 79078f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultSetLogActions(PetscLogHandler handler, PetscBool flag) 79178f1b9b4SToby Isaac { 79278f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 79378f1b9b4SToby Isaac 79478f1b9b4SToby Isaac PetscFunctionBegin; 79578f1b9b4SToby Isaac def->petsc_logActions = flag; 79678f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 79778f1b9b4SToby Isaac } 79878f1b9b4SToby Isaac 79978f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultSetLogObjects(PetscLogHandler handler, PetscBool flag) 80078f1b9b4SToby Isaac { 80178f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 80278f1b9b4SToby Isaac 80378f1b9b4SToby Isaac PetscFunctionBegin; 80478f1b9b4SToby Isaac def->petsc_logObjects = flag; 80578f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 80678f1b9b4SToby Isaac } 80778f1b9b4SToby Isaac 80878f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultLogObjectState(PetscLogHandler handler, PetscObject obj, const char format[], va_list Argp) 80978f1b9b4SToby Isaac { 81078f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 81178f1b9b4SToby Isaac size_t fullLength; 81278f1b9b4SToby Isaac 81378f1b9b4SToby Isaac PetscFunctionBegin; 81478f1b9b4SToby Isaac if (def->petsc_logObjects) { 81578f1b9b4SToby Isaac Object *obj_entry = NULL; 81678f1b9b4SToby Isaac 81778f1b9b4SToby Isaac PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, obj->id - 1, &obj_entry)); 81878f1b9b4SToby Isaac PetscCall(PetscVSNPrintf(obj_entry->info, 64, format, &fullLength, Argp)); 81978f1b9b4SToby Isaac } 82078f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 82178f1b9b4SToby Isaac } 82278f1b9b4SToby Isaac 82378f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultGetNumObjects(PetscLogHandler handler, PetscInt *num_objects) 82478f1b9b4SToby Isaac { 82578f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 82678f1b9b4SToby Isaac 82778f1b9b4SToby Isaac PetscFunctionBegin; 82878f1b9b4SToby Isaac PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, num_objects, NULL)); 82978f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 83078f1b9b4SToby Isaac } 83178f1b9b4SToby Isaac 83278f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerDump_Default(PetscLogHandler handler, const char sname[]) 83378f1b9b4SToby Isaac { 83478f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 83578f1b9b4SToby Isaac FILE *fd; 83678f1b9b4SToby Isaac char file[PETSC_MAX_PATH_LEN], fname[PETSC_MAX_PATH_LEN]; 83778f1b9b4SToby Isaac PetscLogDouble flops, _TotalTime; 83878f1b9b4SToby Isaac PetscMPIInt rank; 83978f1b9b4SToby Isaac int curStage; 84078f1b9b4SToby Isaac PetscLogState state; 84178f1b9b4SToby Isaac PetscInt num_events; 84278f1b9b4SToby Isaac PetscLogEvent event; 84378f1b9b4SToby Isaac 84478f1b9b4SToby Isaac PetscFunctionBegin; 84578f1b9b4SToby Isaac /* Calculate the total elapsed time */ 84678f1b9b4SToby Isaac PetscCall(PetscTime(&_TotalTime)); 84778f1b9b4SToby Isaac _TotalTime -= petsc_BaseTime; 84878f1b9b4SToby Isaac /* Open log file */ 84978f1b9b4SToby Isaac PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)handler), &rank)); 85078f1b9b4SToby Isaac PetscCall(PetscSNPrintf(file, PETSC_STATIC_ARRAY_LENGTH(file), "%s.%d", sname && sname[0] ? sname : "Log", rank)); 85178f1b9b4SToby Isaac PetscCall(PetscFixFilename(file, fname)); 85278f1b9b4SToby Isaac PetscCall(PetscFOpen(PETSC_COMM_SELF, fname, "w", &fd)); 85378f1b9b4SToby Isaac PetscCheck(!(rank == 0) || !(!fd), PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file: %s", fname); 85478f1b9b4SToby Isaac /* Output totals */ 85578f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Total Flop %14e %16.8e\n", petsc_TotalFlops, _TotalTime)); 85678f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Clock Resolution %g\n", 0.0)); 85778f1b9b4SToby Isaac /* Output actions */ 85878f1b9b4SToby Isaac if (def->petsc_logActions) { 85978f1b9b4SToby Isaac PetscInt num_actions; 86078f1b9b4SToby Isaac PetscCall(PetscLogActionArrayGetSize(def->petsc_actions, &num_actions, NULL)); 86178f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Actions accomplished %d\n", (int)num_actions)); 86278f1b9b4SToby Isaac for (int a = 0; a < num_actions; a++) { 86378f1b9b4SToby Isaac Action *action; 86478f1b9b4SToby Isaac 86578f1b9b4SToby Isaac PetscCall(PetscLogActionArrayGetRef(def->petsc_actions, a, &action)); 86678f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "%g %d %d %d %d %d %d %g %g %g\n", action->time, action->action, (int)action->event, (int)action->classid, action->id1, action->id2, action->id3, action->flops, action->mem, action->maxmem)); 86778f1b9b4SToby Isaac } 86878f1b9b4SToby Isaac } 86978f1b9b4SToby Isaac /* Output objects */ 87078f1b9b4SToby Isaac if (def->petsc_logObjects) { 87178f1b9b4SToby Isaac PetscInt num_objects; 87278f1b9b4SToby Isaac 87378f1b9b4SToby Isaac PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL)); 87478f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Objects created %d destroyed %d\n", def->petsc_numObjectsCreated, def->petsc_numObjectsDestroyed)); 87578f1b9b4SToby Isaac for (int o = 0; o < num_objects; o++) { 87678f1b9b4SToby Isaac Object *object = NULL; 87778f1b9b4SToby Isaac 87878f1b9b4SToby Isaac PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, o, &object)); 87978f1b9b4SToby Isaac if (object->parent != -1) continue; // object with this id wasn't logged, probably a PetscLogHandler 88078f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Parent ID: %d Memory: %d\n", object->parent, (int)object->mem)); 88178f1b9b4SToby Isaac if (!object->name[0]) { 88278f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Name\n")); 88378f1b9b4SToby Isaac } else { 88478f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Name: %s\n", object->name)); 88578f1b9b4SToby Isaac } 88678f1b9b4SToby Isaac if (!object->info[0]) { 88778f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Info\n")); 88878f1b9b4SToby Isaac } else { 88978f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Info: %s\n", object->info)); 89078f1b9b4SToby Isaac } 89178f1b9b4SToby Isaac } 89278f1b9b4SToby Isaac } 89378f1b9b4SToby Isaac /* Output events */ 89478f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Event log:\n")); 89578f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(handler, &state)); 89678f1b9b4SToby Isaac PetscCall(PetscLogStateGetNumEvents(state, &num_events)); 89778f1b9b4SToby Isaac PetscCall(PetscLogStateGetCurrentStage(state, &curStage)); 89878f1b9b4SToby Isaac for (event = 0; event < num_events; event++) { 89978f1b9b4SToby Isaac PetscEventPerfInfo *event_info; 90078f1b9b4SToby Isaac 90178f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(handler, curStage, event, &event_info)); 90278f1b9b4SToby Isaac if (event_info->time != 0.0) flops = event_info->flops / event_info->time; 90378f1b9b4SToby Isaac else flops = 0.0; 90478f1b9b4SToby Isaac PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "%d %16d %16g %16g %16g\n", event, event_info->count, event_info->flops, event_info->time, flops)); 90578f1b9b4SToby Isaac } 90678f1b9b4SToby Isaac PetscCall(PetscFClose(PETSC_COMM_SELF, fd)); 90778f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 90878f1b9b4SToby Isaac } 90978f1b9b4SToby Isaac 91078f1b9b4SToby Isaac /* 91178f1b9b4SToby Isaac PetscLogView_Detailed - Each process prints the times for its own events 91278f1b9b4SToby Isaac 91378f1b9b4SToby Isaac */ 91478f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerView_Default_Detailed(PetscLogHandler handler, PetscViewer viewer) 91578f1b9b4SToby Isaac { 91678f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 91778f1b9b4SToby Isaac PetscLogDouble locTotalTime, numRed, maxMem; 91878f1b9b4SToby Isaac PetscInt numStages, numEvents; 91978f1b9b4SToby Isaac MPI_Comm comm = PetscObjectComm((PetscObject)viewer); 92078f1b9b4SToby Isaac PetscMPIInt rank, size; 92178f1b9b4SToby Isaac PetscLogGlobalNames global_stages, global_events; 92278f1b9b4SToby Isaac PetscLogState state; 92378f1b9b4SToby Isaac PetscEventPerfInfo zero_info; 92478f1b9b4SToby Isaac 92578f1b9b4SToby Isaac PetscFunctionBegin; 92678f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(handler, &state)); 92778f1b9b4SToby Isaac PetscCallMPI(MPI_Comm_size(comm, &size)); 92878f1b9b4SToby Isaac PetscCallMPI(MPI_Comm_rank(comm, &rank)); 92978f1b9b4SToby Isaac /* Must preserve reduction count before we go on */ 93078f1b9b4SToby Isaac numRed = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct; 93178f1b9b4SToby Isaac /* Get the total elapsed time */ 93278f1b9b4SToby Isaac PetscCall(PetscTime(&locTotalTime)); 93378f1b9b4SToby Isaac locTotalTime -= petsc_BaseTime; 93478f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "size = %d\n", size)); 93578f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "LocalTimes = {}\n")); 93678f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessages = {}\n")); 93778f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessageLens = {}\n")); 93878f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "LocalReductions = {}\n")); 93978f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "LocalFlop = {}\n")); 94078f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "LocalObjects = {}\n")); 94178f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMemory = {}\n")); 94278f1b9b4SToby Isaac PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages)); 94378f1b9b4SToby Isaac PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events)); 94478f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages)); 94578f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents)); 94678f1b9b4SToby Isaac PetscCall(PetscMemzero(&zero_info, sizeof(zero_info))); 94778f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "Stages = {}\n")); 94878f1b9b4SToby Isaac for (PetscInt stage = 0; stage < numStages; stage++) { 94978f1b9b4SToby Isaac PetscInt stage_id; 95078f1b9b4SToby Isaac const char *stage_name; 95178f1b9b4SToby Isaac 95278f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 95378f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 95478f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"] = {}\n", stage_name)); 95578f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"summary\"] = {}\n", stage_name)); 95678f1b9b4SToby Isaac for (PetscInt event = 0; event < numEvents; event++) { 95778f1b9b4SToby Isaac PetscEventPerfInfo *eventInfo = &zero_info; 95878f1b9b4SToby Isaac PetscBool is_zero = PETSC_FALSE; 95978f1b9b4SToby Isaac PetscInt event_id; 96078f1b9b4SToby Isaac const char *event_name; 96178f1b9b4SToby Isaac 96278f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id)); 96378f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name)); 96478f1b9b4SToby Isaac if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(handler, stage_id, event_id, &eventInfo)); 96578f1b9b4SToby Isaac is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE; 96678f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm)); 96778f1b9b4SToby Isaac if (!is_zero) { PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"%s\"] = {}\n", stage_name, event_name)); } 96878f1b9b4SToby Isaac } 96978f1b9b4SToby Isaac } 97078f1b9b4SToby Isaac PetscCall(PetscMallocGetMaximumUsage(&maxMem)); 97178f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 97278f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalTimes[%d] = %g\n", rank, locTotalTime)); 97378f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessages[%d] = %g\n", rank, (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct))); 97478f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessageLens[%d] = %g\n", rank, (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len))); 97578f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalReductions[%d] = %g\n", rank, numRed)); 97678f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalFlop[%d] = %g\n", rank, petsc_TotalFlops)); 97778f1b9b4SToby Isaac { 97878f1b9b4SToby Isaac PetscInt num_objects; 97978f1b9b4SToby Isaac 98078f1b9b4SToby Isaac PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL)); 98178f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalObjects[%d] = %d\n", rank, (int)num_objects)); 98278f1b9b4SToby Isaac } 98378f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMemory[%d] = %g\n", rank, maxMem)); 98478f1b9b4SToby Isaac PetscCall(PetscViewerFlush(viewer)); 98578f1b9b4SToby Isaac for (PetscInt stage = 0; stage < numStages; stage++) { 98678f1b9b4SToby Isaac PetscEventPerfInfo *stage_perf_info = &zero_info; 98778f1b9b4SToby Isaac PetscInt stage_id; 98878f1b9b4SToby Isaac const char *stage_name; 98978f1b9b4SToby Isaac 99078f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 99178f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 99278f1b9b4SToby Isaac if (stage_id >= 0) { 99378f1b9b4SToby Isaac PetscStagePerf *stage_info; 99478f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info)); 99578f1b9b4SToby Isaac stage_perf_info = &stage_info->perfInfo; 99678f1b9b4SToby Isaac } 99778f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Stages[\"%s\"][\"summary\"][%d] = {\"time\" : %g, \"numMessages\" : %g, \"messageLength\" : %g, \"numReductions\" : %g, \"flop\" : %g}\n", stage_name, rank, stage_perf_info->time, 99878f1b9b4SToby Isaac stage_perf_info->numMessages, stage_perf_info->messageLength, stage_perf_info->numReductions, stage_perf_info->flops)); 99978f1b9b4SToby Isaac for (PetscInt event = 0; event < numEvents; event++) { 100078f1b9b4SToby Isaac PetscEventPerfInfo *eventInfo = &zero_info; 100178f1b9b4SToby Isaac PetscBool is_zero = PETSC_FALSE; 100278f1b9b4SToby Isaac PetscInt event_id; 100378f1b9b4SToby Isaac const char *event_name; 100478f1b9b4SToby Isaac 100578f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id)); 100678f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name)); 100778f1b9b4SToby Isaac if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(handler, stage_id, event_id, &eventInfo)); 100878f1b9b4SToby Isaac is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE; 100978f1b9b4SToby Isaac PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero)); 101078f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm)); 101178f1b9b4SToby Isaac if (!is_zero) { 101278f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Stages[\"%s\"][\"%s\"][%d] = {\"count\" : %d, \"time\" : %g, \"syncTime\" : %g, \"numMessages\" : %g, \"messageLength\" : %g, \"numReductions\" : %g, \"flop\" : %g", stage_name, event_name, rank, 101378f1b9b4SToby Isaac eventInfo->count, eventInfo->time, eventInfo->syncTime, eventInfo->numMessages, eventInfo->messageLength, eventInfo->numReductions, eventInfo->flops)); 101478f1b9b4SToby Isaac if (eventInfo->dof[0] >= 0.) { 101578f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"dof\" : [")); 101678f1b9b4SToby Isaac for (PetscInt d = 0; d < 8; ++d) { 101778f1b9b4SToby Isaac if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", ")); 101878f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->dof[d])); 101978f1b9b4SToby Isaac } 102078f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]")); 102178f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"error\" : [")); 102278f1b9b4SToby Isaac for (PetscInt e = 0; e < 8; ++e) { 102378f1b9b4SToby Isaac if (e > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", ")); 102478f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->errors[e])); 102578f1b9b4SToby Isaac } 102678f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]")); 102778f1b9b4SToby Isaac } 102878f1b9b4SToby Isaac } 102978f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "}\n")); 103078f1b9b4SToby Isaac } 103178f1b9b4SToby Isaac } 103278f1b9b4SToby Isaac PetscCall(PetscViewerFlush(viewer)); 103378f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 103478f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesDestroy(&global_events)); 103578f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesDestroy(&global_stages)); 103678f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 103778f1b9b4SToby Isaac } 103878f1b9b4SToby Isaac 103978f1b9b4SToby Isaac /* 104078f1b9b4SToby Isaac PetscLogView_CSV - Each process prints the times for its own events in Comma-Separated Value Format 104178f1b9b4SToby Isaac */ 104278f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerView_Default_CSV(PetscLogHandler handler, PetscViewer viewer) 104378f1b9b4SToby Isaac { 104478f1b9b4SToby Isaac PetscLogDouble locTotalTime, maxMem; 104578f1b9b4SToby Isaac PetscInt numStages, numEvents, stage, event; 104678f1b9b4SToby Isaac MPI_Comm comm = PetscObjectComm((PetscObject)viewer); 104778f1b9b4SToby Isaac PetscMPIInt rank, size; 104878f1b9b4SToby Isaac PetscLogGlobalNames global_stages, global_events; 104978f1b9b4SToby Isaac PetscLogState state; 105078f1b9b4SToby Isaac PetscEventPerfInfo zero_info; 105178f1b9b4SToby Isaac 105278f1b9b4SToby Isaac PetscFunctionBegin; 105378f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(handler, &state)); 105478f1b9b4SToby Isaac PetscCallMPI(MPI_Comm_size(comm, &size)); 105578f1b9b4SToby Isaac PetscCallMPI(MPI_Comm_rank(comm, &rank)); 105678f1b9b4SToby Isaac /* Must preserve reduction count before we go on */ 105778f1b9b4SToby Isaac /* Get the total elapsed time */ 105878f1b9b4SToby Isaac PetscCall(PetscTime(&locTotalTime)); 105978f1b9b4SToby Isaac locTotalTime -= petsc_BaseTime; 106078f1b9b4SToby Isaac PetscCall(PetscMallocGetMaximumUsage(&maxMem)); 106178f1b9b4SToby Isaac PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages)); 106278f1b9b4SToby Isaac PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events)); 106378f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages)); 106478f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents)); 106578f1b9b4SToby Isaac PetscCall(PetscMemzero(&zero_info, sizeof(zero_info))); 106678f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 106778f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "Stage Name,Event Name,Rank,Count,Time,Num Messages,Message Length,Num Reductions,FLOP,dof0,dof1,dof2,dof3,dof4,dof5,dof6,dof7,e0,e1,e2,e3,e4,e5,e6,e7,%d\n", size)); 106878f1b9b4SToby Isaac PetscCall(PetscViewerFlush(viewer)); 106978f1b9b4SToby Isaac for (stage = 0; stage < numStages; stage++) { 107078f1b9b4SToby Isaac PetscEventPerfInfo *stage_perf_info; 107178f1b9b4SToby Isaac PetscInt stage_id; 107278f1b9b4SToby Isaac const char *stage_name; 107378f1b9b4SToby Isaac 107478f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 107578f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 107678f1b9b4SToby Isaac stage_perf_info = &zero_info; 107778f1b9b4SToby Isaac if (stage_id >= 0) { 107878f1b9b4SToby Isaac PetscStagePerf *stage_info; 107978f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info)); 108078f1b9b4SToby Isaac stage_perf_info = &stage_info->perfInfo; 108178f1b9b4SToby Isaac } 108278f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s,summary,%d,1,%g,%g,%g,%g,%g\n", stage_name, rank, stage_perf_info->time, stage_perf_info->numMessages, stage_perf_info->messageLength, stage_perf_info->numReductions, stage_perf_info->flops)); 108378f1b9b4SToby Isaac for (event = 0; event < numEvents; event++) { 108478f1b9b4SToby Isaac PetscEventPerfInfo *eventInfo = &zero_info; 108578f1b9b4SToby Isaac PetscBool is_zero = PETSC_FALSE; 108678f1b9b4SToby Isaac PetscInt event_id; 108778f1b9b4SToby Isaac const char *event_name; 108878f1b9b4SToby Isaac 108978f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id)); 109078f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name)); 109178f1b9b4SToby Isaac if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(handler, stage_id, event_id, &eventInfo)); 109278f1b9b4SToby Isaac PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero)); 109378f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm)); 109478f1b9b4SToby Isaac if (!is_zero) { 109578f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s,%s,%d,%d,%g,%g,%g,%g,%g", stage_name, event_name, rank, eventInfo->count, eventInfo->time, eventInfo->numMessages, eventInfo->messageLength, eventInfo->numReductions, eventInfo->flops)); 109678f1b9b4SToby Isaac if (eventInfo->dof[0] >= 0.) { 109778f1b9b4SToby Isaac for (PetscInt d = 0; d < 8; ++d) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->dof[d])); 109878f1b9b4SToby Isaac for (PetscInt e = 0; e < 8; ++e) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->errors[e])); 109978f1b9b4SToby Isaac } 110078f1b9b4SToby Isaac } 110178f1b9b4SToby Isaac PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 110278f1b9b4SToby Isaac } 110378f1b9b4SToby Isaac } 110478f1b9b4SToby Isaac PetscCall(PetscViewerFlush(viewer)); 110578f1b9b4SToby Isaac PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 110678f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesDestroy(&global_stages)); 110778f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesDestroy(&global_events)); 110878f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 110978f1b9b4SToby Isaac } 111078f1b9b4SToby Isaac 111178f1b9b4SToby Isaac static PetscErrorCode PetscLogViewWarnSync(MPI_Comm comm, FILE *fd) 111278f1b9b4SToby Isaac { 111378f1b9b4SToby Isaac PetscFunctionBegin; 111478f1b9b4SToby Isaac if (!PetscLogSyncOn) PetscFunctionReturn(PETSC_SUCCESS); 111578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n\n")); 111678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n")); 111778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 111878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # WARNING!!! #\n")); 111978f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 112078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # This program was run with logging synchronization. #\n")); 112178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # This option provides more meaningful imbalance #\n")); 112278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # figures at the expense of slowing things down and #\n")); 112378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # providing a distorted view of the overall runtime. #\n")); 112478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 112578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n\n\n")); 112678f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 112778f1b9b4SToby Isaac } 112878f1b9b4SToby Isaac 112978f1b9b4SToby Isaac static PetscErrorCode PetscLogViewWarnDebugging(MPI_Comm comm, FILE *fd) 113078f1b9b4SToby Isaac { 113178f1b9b4SToby Isaac PetscFunctionBegin; 113278f1b9b4SToby Isaac if (PetscDefined(USE_DEBUG)) { 113378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n\n")); 113478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n")); 113578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 113678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # WARNING!!! #\n")); 113778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 113878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # This code was compiled with a debugging option. #\n")); 113978f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # To get timing results run ./configure #\n")); 114078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # using --with-debugging=no, the performance will #\n")); 114178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # be generally two or three times faster. #\n")); 114278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 114378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n\n\n")); 114478f1b9b4SToby Isaac } 114578f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 114678f1b9b4SToby Isaac } 114778f1b9b4SToby Isaac 114878f1b9b4SToby Isaac static PetscErrorCode PetscLogViewWarnNoGpuAwareMpi(MPI_Comm comm, FILE *fd) 114978f1b9b4SToby Isaac { 115078f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 115178f1b9b4SToby Isaac PetscMPIInt size; 115278f1b9b4SToby Isaac PetscBool deviceInitialized = PETSC_FALSE; 115378f1b9b4SToby Isaac 115478f1b9b4SToby Isaac PetscFunctionBegin; 115578f1b9b4SToby Isaac PetscCallMPI(MPI_Comm_size(comm, &size)); 115678f1b9b4SToby Isaac for (int i = PETSC_DEVICE_HOST + 1; i < PETSC_DEVICE_MAX; ++i) { 115778f1b9b4SToby Isaac const PetscDeviceType dtype = PetscDeviceTypeCast(i); 115878f1b9b4SToby Isaac if (PetscDeviceInitialized(dtype)) { /* a non-host device was initialized */ 115978f1b9b4SToby Isaac deviceInitialized = PETSC_TRUE; 116078f1b9b4SToby Isaac break; 116178f1b9b4SToby Isaac } 116278f1b9b4SToby Isaac } 116378f1b9b4SToby Isaac /* the last condition says petsc is configured with device but it is a pure CPU run, so don't print misleading warnings */ 116478f1b9b4SToby Isaac if (use_gpu_aware_mpi || size == 1 || !deviceInitialized) PetscFunctionReturn(PETSC_SUCCESS); 116578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n\n")); 116678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n")); 116778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 116878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # WARNING!!! #\n")); 116978f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 117078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # This code was compiled with GPU support and you've #\n")); 117178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # created PETSc/GPU objects, but you intentionally #\n")); 117278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # used -use_gpu_aware_mpi 0, requiring PETSc to copy #\n")); 117378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # additional data between the GPU and CPU. To obtain #\n")); 117478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # meaningful timing results on multi-rank runs, use #\n")); 117578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # GPU-aware MPI instead. #\n")); 117678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 117778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n\n\n")); 117878f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 117978f1b9b4SToby Isaac #else 118078f1b9b4SToby Isaac return PETSC_SUCCESS; 118178f1b9b4SToby Isaac #endif 118278f1b9b4SToby Isaac } 118378f1b9b4SToby Isaac 118478f1b9b4SToby Isaac static PetscErrorCode PetscLogViewWarnGpuTime(MPI_Comm comm, FILE *fd) 118578f1b9b4SToby Isaac { 118678f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 118778f1b9b4SToby Isaac 118878f1b9b4SToby Isaac PetscFunctionBegin; 118978f1b9b4SToby Isaac if (!PetscLogGpuTimeFlag || petsc_gflops == 0) PetscFunctionReturn(PETSC_SUCCESS); 119078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n\n")); 119178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n")); 119278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 119378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # WARNING!!! #\n")); 119478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 119578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # This code was run with -log_view_gpu_time #\n")); 119678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # This provides accurate timing within the GPU kernels #\n")); 119778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # but can slow down the entire computation by a #\n")); 119878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # measurable amount. For fastest runs we recommend #\n")); 119978f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # not using this option. #\n")); 120078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " # #\n")); 120178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n\n\n")); 120278f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 120378f1b9b4SToby Isaac #else 120478f1b9b4SToby Isaac return PETSC_SUCCESS; 120578f1b9b4SToby Isaac #endif 120678f1b9b4SToby Isaac } 120778f1b9b4SToby Isaac 120878f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerView_Default_Info(PetscLogHandler handler, PetscViewer viewer) 120978f1b9b4SToby Isaac { 121078f1b9b4SToby Isaac FILE *fd; 121178f1b9b4SToby Isaac PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 121278f1b9b4SToby Isaac char arch[128], hostname[128], username[128], pname[PETSC_MAX_PATH_LEN], date[128]; 121378f1b9b4SToby Isaac PetscLogDouble locTotalTime, TotalTime, TotalFlops; 121478f1b9b4SToby Isaac PetscLogDouble numMessages, messageLength, avgMessLen, numReductions; 121578f1b9b4SToby Isaac PetscLogDouble stageTime, flops, flopr, mem, mess, messLen, red; 121678f1b9b4SToby Isaac PetscLogDouble fracTime, fracFlops, fracMessages, fracLength, fracReductions, fracMess, fracMessLen, fracRed; 121778f1b9b4SToby Isaac PetscLogDouble fracStageTime, fracStageFlops, fracStageMess, fracStageMessLen, fracStageRed; 121878f1b9b4SToby Isaac PetscLogDouble min, max, tot, ratio, avg, x, y; 121978f1b9b4SToby Isaac PetscLogDouble minf, maxf, totf, ratf, mint, maxt, tott, ratt, ratC, totm, totml, totr, mal, malmax, emalmax; 122078f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 122178f1b9b4SToby Isaac PetscLogEvent KSP_Solve, SNES_Solve, TS_Step, TAO_Solve; /* These need to be fixed to be some events registered with certain objects */ 122278f1b9b4SToby Isaac PetscLogDouble cct, gct, csz, gsz, gmaxt, gflops, gflopr, fracgflops; 122378f1b9b4SToby Isaac #endif 122478f1b9b4SToby Isaac PetscMPIInt minC, maxC; 122578f1b9b4SToby Isaac PetscMPIInt size, rank; 122678f1b9b4SToby Isaac PetscBool *localStageUsed, *stageUsed; 122778f1b9b4SToby Isaac PetscBool *localStageVisible, *stageVisible; 122878f1b9b4SToby Isaac PetscInt numStages, numEvents; 122978f1b9b4SToby Isaac int stage, oclass; 123078f1b9b4SToby Isaac PetscLogEvent event; 123178f1b9b4SToby Isaac char version[256]; 123278f1b9b4SToby Isaac MPI_Comm comm; 123378f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 123478f1b9b4SToby Isaac PetscInt64 nas = 0x7FF0000000000002; 123578f1b9b4SToby Isaac #endif 123678f1b9b4SToby Isaac PetscLogGlobalNames global_stages, global_events; 123778f1b9b4SToby Isaac PetscEventPerfInfo zero_info; 123878f1b9b4SToby Isaac PetscLogState state; 123978f1b9b4SToby Isaac 124078f1b9b4SToby Isaac PetscFunctionBegin; 124178f1b9b4SToby Isaac PetscCall(PetscLogHandlerGetState(handler, &state)); 124278f1b9b4SToby Isaac PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 124378f1b9b4SToby Isaac PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm)); 124478f1b9b4SToby Isaac PetscCall(PetscViewerASCIIGetPointer(viewer, &fd)); 124578f1b9b4SToby Isaac PetscCallMPI(MPI_Comm_size(comm, &size)); 124678f1b9b4SToby Isaac PetscCallMPI(MPI_Comm_rank(comm, &rank)); 124778f1b9b4SToby Isaac /* Get the total elapsed time */ 124878f1b9b4SToby Isaac PetscCall(PetscTime(&locTotalTime)); 124978f1b9b4SToby Isaac locTotalTime -= petsc_BaseTime; 125078f1b9b4SToby Isaac 125178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "****************************************************************************************************************************************************************\n")); 125278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "*** WIDEN YOUR WINDOW TO 160 CHARACTERS. Use 'enscript -r -fCourier9' to print this document ***\n")); 125378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "****************************************************************************************************************************************************************\n")); 125478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n------------------------------------------------------------------ PETSc Performance Summary: ------------------------------------------------------------------\n\n")); 125578f1b9b4SToby Isaac PetscCall(PetscLogViewWarnSync(comm, fd)); 125678f1b9b4SToby Isaac PetscCall(PetscLogViewWarnDebugging(comm, fd)); 125778f1b9b4SToby Isaac PetscCall(PetscLogViewWarnNoGpuAwareMpi(comm, fd)); 125878f1b9b4SToby Isaac PetscCall(PetscLogViewWarnGpuTime(comm, fd)); 125978f1b9b4SToby Isaac PetscCall(PetscGetArchType(arch, sizeof(arch))); 126078f1b9b4SToby Isaac PetscCall(PetscGetHostName(hostname, sizeof(hostname))); 126178f1b9b4SToby Isaac PetscCall(PetscGetUserName(username, sizeof(username))); 126278f1b9b4SToby Isaac PetscCall(PetscGetProgramName(pname, sizeof(pname))); 126378f1b9b4SToby Isaac PetscCall(PetscGetDate(date, sizeof(date))); 126478f1b9b4SToby Isaac PetscCall(PetscGetVersion(version, sizeof(version))); 126578f1b9b4SToby Isaac if (size == 1) { 126678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "%s on a %s named %s with %d processor, by %s %s\n", pname, arch, hostname, size, username, date)); 126778f1b9b4SToby Isaac } else { 126878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "%s on a %s named %s with %d processors, by %s %s\n", pname, arch, hostname, size, username, date)); 126978f1b9b4SToby Isaac } 127078f1b9b4SToby Isaac #if defined(PETSC_HAVE_OPENMP) 127178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Using %" PetscInt_FMT " OpenMP threads\n", PetscNumOMPThreads)); 127278f1b9b4SToby Isaac #endif 127378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Using %s\n", version)); 127478f1b9b4SToby Isaac 127578f1b9b4SToby Isaac /* Must preserve reduction count before we go on */ 127678f1b9b4SToby Isaac red = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct; 127778f1b9b4SToby Isaac 127878f1b9b4SToby Isaac /* Calculate summary information */ 127978f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n Max Max/Min Avg Total\n")); 128078f1b9b4SToby Isaac /* Time */ 128178f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&locTotalTime, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 128278f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&locTotalTime, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 128378f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&locTotalTime, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 128478f1b9b4SToby Isaac avg = tot / ((PetscLogDouble)size); 128578f1b9b4SToby Isaac if (min != 0.0) ratio = max / min; 128678f1b9b4SToby Isaac else ratio = 0.0; 128778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Time (sec): %5.3e %7.3f %5.3e\n", max, ratio, avg)); 128878f1b9b4SToby Isaac TotalTime = tot; 128978f1b9b4SToby Isaac /* Objects */ 129078f1b9b4SToby Isaac { 129178f1b9b4SToby Isaac PetscInt num_objects; 129278f1b9b4SToby Isaac 129378f1b9b4SToby Isaac PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL)); 129478f1b9b4SToby Isaac avg = (PetscLogDouble)num_objects; 129578f1b9b4SToby Isaac } 129678f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&avg, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 129778f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&avg, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 129878f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&avg, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 129978f1b9b4SToby Isaac avg = tot / ((PetscLogDouble)size); 130078f1b9b4SToby Isaac if (min != 0.0) ratio = max / min; 130178f1b9b4SToby Isaac else ratio = 0.0; 130278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Objects: %5.3e %7.3f %5.3e\n", max, ratio, avg)); 130378f1b9b4SToby Isaac /* Flops */ 130478f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 130578f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 130678f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 130778f1b9b4SToby Isaac avg = tot / ((PetscLogDouble)size); 130878f1b9b4SToby Isaac if (min != 0.0) ratio = max / min; 130978f1b9b4SToby Isaac else ratio = 0.0; 131078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Flops: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot)); 131178f1b9b4SToby Isaac TotalFlops = tot; 131278f1b9b4SToby Isaac /* Flops/sec -- Must talk to Barry here */ 131378f1b9b4SToby Isaac if (locTotalTime != 0.0) flops = petsc_TotalFlops / locTotalTime; 131478f1b9b4SToby Isaac else flops = 0.0; 131578f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&flops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 131678f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&flops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 131778f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&flops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 131878f1b9b4SToby Isaac avg = tot / ((PetscLogDouble)size); 131978f1b9b4SToby Isaac if (min != 0.0) ratio = max / min; 132078f1b9b4SToby Isaac else ratio = 0.0; 132178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Flops/sec: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot)); 132278f1b9b4SToby Isaac /* Memory */ 132378f1b9b4SToby Isaac PetscCall(PetscMallocGetMaximumUsage(&mem)); 132478f1b9b4SToby Isaac if (mem > 0.0) { 132578f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&mem, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 132678f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&mem, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 132778f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&mem, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 132878f1b9b4SToby Isaac avg = tot / ((PetscLogDouble)size); 132978f1b9b4SToby Isaac if (min != 0.0) ratio = max / min; 133078f1b9b4SToby Isaac else ratio = 0.0; 133178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Memory (bytes): %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot)); 133278f1b9b4SToby Isaac } 133378f1b9b4SToby Isaac /* Messages */ 133478f1b9b4SToby Isaac mess = 0.5 * (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct); 133578f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 133678f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 133778f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 133878f1b9b4SToby Isaac avg = tot / ((PetscLogDouble)size); 133978f1b9b4SToby Isaac if (min != 0.0) ratio = max / min; 134078f1b9b4SToby Isaac else ratio = 0.0; 134178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "MPI Msg Count: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot)); 134278f1b9b4SToby Isaac numMessages = tot; 134378f1b9b4SToby Isaac /* Message Lengths */ 134478f1b9b4SToby Isaac mess = 0.5 * (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len); 134578f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 134678f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 134778f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 134878f1b9b4SToby Isaac if (numMessages != 0) avg = tot / numMessages; 134978f1b9b4SToby Isaac else avg = 0.0; 135078f1b9b4SToby Isaac if (min != 0.0) ratio = max / min; 135178f1b9b4SToby Isaac else ratio = 0.0; 135278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "MPI Msg Len (bytes): %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot)); 135378f1b9b4SToby Isaac messageLength = tot; 135478f1b9b4SToby Isaac /* Reductions */ 135578f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&red, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 135678f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&red, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 135778f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&red, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 135878f1b9b4SToby Isaac if (min != 0.0) ratio = max / min; 135978f1b9b4SToby Isaac else ratio = 0.0; 136078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "MPI Reductions: %5.3e %7.3f\n", max, ratio)); 136178f1b9b4SToby Isaac numReductions = red; /* wrong because uses count from process zero */ 136278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\nFlop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)\n")); 136378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " e.g., VecAXPY() for real vectors of length N --> 2N flops\n")); 136478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " and VecAXPY() for complex vectors of length N --> 8N flops\n")); 136578f1b9b4SToby Isaac 136678f1b9b4SToby Isaac PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages)); 136778f1b9b4SToby Isaac PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events)); 136878f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages)); 136978f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents)); 137078f1b9b4SToby Isaac PetscCall(PetscMemzero(&zero_info, sizeof(zero_info))); 137178f1b9b4SToby Isaac PetscCall(PetscMalloc1(numStages, &localStageUsed)); 137278f1b9b4SToby Isaac PetscCall(PetscMalloc1(numStages, &stageUsed)); 137378f1b9b4SToby Isaac PetscCall(PetscMalloc1(numStages, &localStageVisible)); 137478f1b9b4SToby Isaac PetscCall(PetscMalloc1(numStages, &stageVisible)); 137578f1b9b4SToby Isaac if (numStages > 0) { 137678f1b9b4SToby Isaac for (stage = 0; stage < numStages; stage++) { 137778f1b9b4SToby Isaac PetscInt stage_id; 137878f1b9b4SToby Isaac 137978f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 138078f1b9b4SToby Isaac if (stage_id >= 0) { 138178f1b9b4SToby Isaac PetscStagePerf *stage_info; 138278f1b9b4SToby Isaac 138378f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info)); 138478f1b9b4SToby Isaac localStageUsed[stage] = stage_info->used; 138578f1b9b4SToby Isaac localStageVisible[stage] = stage_info->perfInfo.visible; 138678f1b9b4SToby Isaac } else { 138778f1b9b4SToby Isaac localStageUsed[stage] = PETSC_FALSE; 138878f1b9b4SToby Isaac localStageVisible[stage] = PETSC_TRUE; 138978f1b9b4SToby Isaac } 139078f1b9b4SToby Isaac } 139178f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(localStageUsed, stageUsed, numStages, MPIU_BOOL, MPI_LOR, comm)); 139278f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm)); 139378f1b9b4SToby Isaac for (stage = 0; stage < numStages; stage++) { 139478f1b9b4SToby Isaac if (stageUsed[stage] && stageVisible[stage]) { 139578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\nSummary of Stages: ----- Time ------ ----- Flop ------ --- Messages --- -- Message Lengths -- -- Reductions --\n")); 139678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Avg %%Total Avg %%Total Count %%Total Avg %%Total Count %%Total\n")); 139778f1b9b4SToby Isaac break; 139878f1b9b4SToby Isaac } 139978f1b9b4SToby Isaac } 140078f1b9b4SToby Isaac for (stage = 0; stage < numStages; stage++) { 140178f1b9b4SToby Isaac PetscInt stage_id; 140278f1b9b4SToby Isaac PetscEventPerfInfo *stage_info; 140378f1b9b4SToby Isaac const char *stage_name; 140478f1b9b4SToby Isaac 140578f1b9b4SToby Isaac if (!(stageUsed[stage] && stageVisible[stage])) continue; 140678f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 140778f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 140878f1b9b4SToby Isaac stage_info = &zero_info; 140978f1b9b4SToby Isaac if (localStageUsed[stage]) { 141078f1b9b4SToby Isaac PetscStagePerf *stage_perf_info; 141178f1b9b4SToby Isaac 141278f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_perf_info)); 141378f1b9b4SToby Isaac stage_info = &stage_perf_info->perfInfo; 141478f1b9b4SToby Isaac } 141578f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 141678f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 141778f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 141878f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 141978f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 142078f1b9b4SToby Isaac mess *= 0.5; 142178f1b9b4SToby Isaac messLen *= 0.5; 142278f1b9b4SToby Isaac red /= size; 142378f1b9b4SToby Isaac if (TotalTime != 0.0) fracTime = stageTime / TotalTime; 142478f1b9b4SToby Isaac else fracTime = 0.0; 142578f1b9b4SToby Isaac if (TotalFlops != 0.0) fracFlops = flops / TotalFlops; 142678f1b9b4SToby Isaac else fracFlops = 0.0; 142778f1b9b4SToby Isaac /* Talk to Barry if (stageTime != 0.0) flops = (size*flops)/stageTime; else flops = 0.0; */ 142878f1b9b4SToby Isaac if (numMessages != 0.0) fracMessages = mess / numMessages; 142978f1b9b4SToby Isaac else fracMessages = 0.0; 143078f1b9b4SToby Isaac if (mess != 0.0) avgMessLen = messLen / mess; 143178f1b9b4SToby Isaac else avgMessLen = 0.0; 143278f1b9b4SToby Isaac if (messageLength != 0.0) fracLength = messLen / messageLength; 143378f1b9b4SToby Isaac else fracLength = 0.0; 143478f1b9b4SToby Isaac if (numReductions != 0.0) fracReductions = red / numReductions; 143578f1b9b4SToby Isaac else fracReductions = 0.0; 143678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "%2d: %15s: %6.4e %5.1f%% %6.4e %5.1f%% %5.3e %5.1f%% %5.3e %5.1f%% %5.3e %5.1f%%\n", stage, stage_name, stageTime / size, 100.0 * fracTime, flops, 100.0 * fracFlops, mess, 100.0 * fracMessages, avgMessLen, 100.0 * fracLength, red, 100.0 * fracReductions)); 143778f1b9b4SToby Isaac } 143878f1b9b4SToby Isaac } 143978f1b9b4SToby Isaac 144078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n------------------------------------------------------------------------------------------------------------------------\n")); 144178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "See the 'Profiling' chapter of the users' manual for details on interpreting output.\n")); 144278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Phase summary info:\n")); 144378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Count: number of times phase was executed\n")); 144478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Time and Flop: Max - maximum over all processors\n")); 144578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Ratio - ratio of maximum to minimum over all processors\n")); 144678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Mess: number of messages sent\n")); 144778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " AvgLen: average message length (bytes)\n")); 144878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Reduct: number of global reductions\n")); 144978f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Global: entire computation\n")); 145078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().\n")); 145178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " %%T - percent time in this phase %%F - percent flop in this phase\n")); 145278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " %%M - percent messages in this phase %%L - percent message lengths in this phase\n")); 145378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " %%R - percent reductions in this phase\n")); 145478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Total Mflop/s: 10e-6 * (sum of flop over all processors)/(max time over all processors)\n")); 145578f1b9b4SToby Isaac if (PetscLogMemory) { 145678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Memory usage is summed over all MPI processes, it is given in mega-bytes\n")); 145778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Malloc Mbytes: Memory allocated and kept during event (sum over all calls to event). May be negative\n")); 145878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " EMalloc Mbytes: extra memory allocated during event and then freed (maximum over all calls to events). Never negative\n")); 145978f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " MMalloc Mbytes: Increase in high water mark of allocated memory (sum over all calls to event). Never negative\n")); 146078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " RMI Mbytes: Increase in resident memory (sum over all calls to event)\n")); 146178f1b9b4SToby Isaac } 146278f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 146378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " GPU Mflop/s: 10e-6 * (sum of flop on GPU over all processors)/(max GPU time over all processors)\n")); 146478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " CpuToGpu Count: total number of CPU to GPU copies per processor\n")); 146578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " CpuToGpu Size (Mbytes): 10e-6 * (total size of CPU to GPU copies per processor)\n")); 146678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " GpuToCpu Count: total number of GPU to CPU copies per processor\n")); 146778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " GpuToCpu Size (Mbytes): 10e-6 * (total size of GPU to CPU copies per processor)\n")); 146878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " GPU %%F: percent flops on GPU in this event\n")); 146978f1b9b4SToby Isaac #endif 147078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------\n")); 147178f1b9b4SToby Isaac 147278f1b9b4SToby Isaac PetscCall(PetscLogViewWarnDebugging(comm, fd)); 147378f1b9b4SToby Isaac 147478f1b9b4SToby Isaac /* Report events */ 147578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Event Count Time (sec) Flop --- Global --- --- Stage ---- Total")); 147678f1b9b4SToby Isaac if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, " Malloc EMalloc MMalloc RMI")); 147778f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 147878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " GPU - CpuToGpu - - GpuToCpu - GPU")); 147978f1b9b4SToby Isaac #endif 148078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n")); 148178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Max Ratio Max Ratio Max Ratio Mess AvgLen Reduct %%T %%F %%M %%L %%R %%T %%F %%M %%L %%R Mflop/s")); 148278f1b9b4SToby Isaac if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, " Mbytes Mbytes Mbytes Mbytes")); 148378f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 148478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " Mflop/s Count Size Count Size %%F")); 148578f1b9b4SToby Isaac #endif 148678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n")); 148778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------")); 148878f1b9b4SToby Isaac if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, "-----------------------------")); 148978f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 149078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "---------------------------------------")); 149178f1b9b4SToby Isaac #endif 149278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n")); 149378f1b9b4SToby Isaac 149478f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 149578f1b9b4SToby Isaac /* this indirect way of accessing these values is needed when PETSc is build with multiple libraries since the symbols are not in libpetscsys */ 149678f1b9b4SToby Isaac PetscCall(PetscLogStateGetEventFromName(state, "TAOSolve", &TAO_Solve)); 149778f1b9b4SToby Isaac PetscCall(PetscLogStateGetEventFromName(state, "TSStep", &TS_Step)); 149878f1b9b4SToby Isaac PetscCall(PetscLogStateGetEventFromName(state, "SNESSolve", &SNES_Solve)); 149978f1b9b4SToby Isaac PetscCall(PetscLogStateGetEventFromName(state, "KSPSolve", &KSP_Solve)); 150078f1b9b4SToby Isaac #endif 150178f1b9b4SToby Isaac 150278f1b9b4SToby Isaac for (stage = 0; stage < numStages; stage++) { 150378f1b9b4SToby Isaac PetscInt stage_id; 150478f1b9b4SToby Isaac PetscEventPerfInfo *stage_info; 150578f1b9b4SToby Isaac const char *stage_name; 150678f1b9b4SToby Isaac 150778f1b9b4SToby Isaac if (!(stageVisible[stage] && stageUsed[stage])) continue; 150878f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 150978f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 151078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stage_name)); 151178f1b9b4SToby Isaac stage_info = &zero_info; 151278f1b9b4SToby Isaac if (localStageUsed[stage]) { 151378f1b9b4SToby Isaac PetscStagePerf *stage_perf_info; 151478f1b9b4SToby Isaac 151578f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_perf_info)); 151678f1b9b4SToby Isaac stage_info = &stage_perf_info->perfInfo; 151778f1b9b4SToby Isaac } 151878f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 151978f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 152078f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 152178f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 152278f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 152378f1b9b4SToby Isaac mess *= 0.5; 152478f1b9b4SToby Isaac messLen *= 0.5; 152578f1b9b4SToby Isaac red /= size; 152678f1b9b4SToby Isaac 152778f1b9b4SToby Isaac for (event = 0; event < numEvents; event++) { 152878f1b9b4SToby Isaac PetscInt event_id; 152978f1b9b4SToby Isaac PetscEventPerfInfo *event_info = &zero_info; 153078f1b9b4SToby Isaac PetscBool is_zero = PETSC_FALSE; 153178f1b9b4SToby Isaac const char *event_name; 153278f1b9b4SToby Isaac 153378f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id)); 153478f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name)); 153578f1b9b4SToby Isaac if (event_id >= 0 && stage_id >= 0) { PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(handler, stage_id, event_id, &event_info)); } 153678f1b9b4SToby Isaac PetscCall(PetscMemcmp(event_info, &zero_info, sizeof(zero_info), &is_zero)); 153778f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm)); 153878f1b9b4SToby Isaac if (!is_zero) { 153978f1b9b4SToby Isaac flopr = event_info->flops; 154078f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&flopr, &minf, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 154178f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&flopr, &maxf, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 154278f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->flops, &totf, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 154378f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->time, &mint, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 154478f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->time, &maxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 154578f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->time, &tott, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 154678f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->numMessages, &totm, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 154778f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 154878f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->numReductions, &totr, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 154978f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->count, &minC, 1, MPI_INT, MPI_MIN, comm)); 155078f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->count, &maxC, 1, MPI_INT, MPI_MAX, comm)); 155178f1b9b4SToby Isaac if (PetscLogMemory) { 155278f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->memIncrease, &mem, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 155378f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->mallocSpace, &mal, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 155478f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->mallocIncrease, &malmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 155578f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->mallocIncreaseEvent, &emalmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 155678f1b9b4SToby Isaac } 155778f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 155878f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->CpuToGpuCount, &cct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 155978f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->GpuToCpuCount, &gct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 156078f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->CpuToGpuSize, &csz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 156178f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->GpuToCpuSize, &gsz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 156278f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->GpuFlops, &gflops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 156378f1b9b4SToby Isaac PetscCall(MPIU_Allreduce(&event_info->GpuTime, &gmaxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 156478f1b9b4SToby Isaac #endif 156578f1b9b4SToby Isaac if (mint < 0.0) { 156678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "WARNING!!! Minimum time %g over all processors for %s is negative! This happens\n on some machines whose times cannot handle too rapid calls.!\n artificially changing minimum to zero.\n", mint, event_name)); 156778f1b9b4SToby Isaac mint = 0; 156878f1b9b4SToby Isaac } 156978f1b9b4SToby Isaac PetscCheck(minf >= 0.0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Minimum flop %g over all processors for %s is negative! Not possible!", minf, event_name); 157078f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 157178f1b9b4SToby Isaac /* Put NaN into the time for all events that may not be time accurately since they may happen asynchronously on the GPU */ 157278f1b9b4SToby Isaac if (!PetscLogGpuTimeFlag && petsc_gflops > 0) { 157378f1b9b4SToby Isaac memcpy(&gmaxt, &nas, sizeof(PetscLogDouble)); 157478f1b9b4SToby Isaac if (event_id != SNES_Solve && event_id != KSP_Solve && event_id != TS_Step && event_id != TAO_Solve) { 157578f1b9b4SToby Isaac memcpy(&mint, &nas, sizeof(PetscLogDouble)); 157678f1b9b4SToby Isaac memcpy(&maxt, &nas, sizeof(PetscLogDouble)); 157778f1b9b4SToby Isaac } 157878f1b9b4SToby Isaac } 157978f1b9b4SToby Isaac #endif 158078f1b9b4SToby Isaac totm *= 0.5; 158178f1b9b4SToby Isaac totml *= 0.5; 158278f1b9b4SToby Isaac totr /= size; 158378f1b9b4SToby Isaac 158478f1b9b4SToby Isaac if (maxC != 0) { 158578f1b9b4SToby Isaac if (minC != 0) ratC = ((PetscLogDouble)maxC) / minC; 158678f1b9b4SToby Isaac else ratC = 0.0; 158778f1b9b4SToby Isaac if (mint != 0.0) ratt = maxt / mint; 158878f1b9b4SToby Isaac else ratt = 0.0; 158978f1b9b4SToby Isaac if (minf != 0.0) ratf = maxf / minf; 159078f1b9b4SToby Isaac else ratf = 0.0; 159178f1b9b4SToby Isaac if (TotalTime != 0.0) fracTime = tott / TotalTime; 159278f1b9b4SToby Isaac else fracTime = 0.0; 159378f1b9b4SToby Isaac if (TotalFlops != 0.0) fracFlops = totf / TotalFlops; 159478f1b9b4SToby Isaac else fracFlops = 0.0; 159578f1b9b4SToby Isaac if (stageTime != 0.0) fracStageTime = tott / stageTime; 159678f1b9b4SToby Isaac else fracStageTime = 0.0; 159778f1b9b4SToby Isaac if (flops != 0.0) fracStageFlops = totf / flops; 159878f1b9b4SToby Isaac else fracStageFlops = 0.0; 159978f1b9b4SToby Isaac if (numMessages != 0.0) fracMess = totm / numMessages; 160078f1b9b4SToby Isaac else fracMess = 0.0; 160178f1b9b4SToby Isaac if (messageLength != 0.0) fracMessLen = totml / messageLength; 160278f1b9b4SToby Isaac else fracMessLen = 0.0; 160378f1b9b4SToby Isaac if (numReductions != 0.0) fracRed = totr / numReductions; 160478f1b9b4SToby Isaac else fracRed = 0.0; 160578f1b9b4SToby Isaac if (mess != 0.0) fracStageMess = totm / mess; 160678f1b9b4SToby Isaac else fracStageMess = 0.0; 160778f1b9b4SToby Isaac if (messLen != 0.0) fracStageMessLen = totml / messLen; 160878f1b9b4SToby Isaac else fracStageMessLen = 0.0; 160978f1b9b4SToby Isaac if (red != 0.0) fracStageRed = totr / red; 161078f1b9b4SToby Isaac else fracStageRed = 0.0; 161178f1b9b4SToby Isaac if (totm != 0.0) totml /= totm; 161278f1b9b4SToby Isaac else totml = 0.0; 161378f1b9b4SToby Isaac if (maxt != 0.0) flopr = totf / maxt; 161478f1b9b4SToby Isaac else flopr = 0.0; 161578f1b9b4SToby Isaac if (fracStageTime > 1.0 || fracStageFlops > 1.0 || fracStageMess > 1.0 || fracStageMessLen > 1.0 || fracStageRed > 1.0) 161678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "%-16s %7d %3.1f %5.4e %3.1f %3.2e %3.1f %2.1e %2.1e %2.1e %2.0f %2.0f %2.0f %2.0f %2.0f Multiple stages %5.0f", event_name, maxC, ratC, maxt, ratt, maxf, ratf, totm, totml, totr, 100.0 * fracTime, 100.0 * fracFlops, 100.0 * fracMess, 100.0 * fracMessLen, 100.0 * fracRed, PetscAbs(flopr) / 1.0e6)); 161778f1b9b4SToby Isaac else 161878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "%-16s %7d %3.1f %5.4e %3.1f %3.2e %3.1f %2.1e %2.1e %2.1e %2.0f %2.0f %2.0f %2.0f %2.0f %3.0f %2.0f %2.0f %2.0f %2.0f %5.0f", event_name, maxC, ratC, maxt, ratt, maxf, ratf, totm, totml, totr, 100.0 * fracTime, 100.0 * fracFlops, 100.0 * fracMess, 100.0 * fracMessLen, 100.0 * fracRed, 100.0 * fracStageTime, 100.0 * fracStageFlops, 100.0 * fracStageMess, 100.0 * fracStageMessLen, 100.0 * fracStageRed, PetscAbs(flopr) / 1.0e6)); 161978f1b9b4SToby Isaac if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, " %5.0f %5.0f %5.0f %5.0f", mal / 1.0e6, emalmax / 1.0e6, malmax / 1.0e6, mem / 1.0e6)); 162078f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 162178f1b9b4SToby Isaac if (totf != 0.0) fracgflops = gflops / totf; 162278f1b9b4SToby Isaac else fracgflops = 0.0; 162378f1b9b4SToby Isaac if (gmaxt != 0.0) gflopr = gflops / gmaxt; 162478f1b9b4SToby Isaac else gflopr = 0.0; 162578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, " %5.0f %4.0f %3.2e %4.0f %3.2e % 2.0f", PetscAbs(gflopr) / 1.0e6, cct / size, csz / (1.0e6 * size), gct / size, gsz / (1.0e6 * size), 100.0 * fracgflops)); 162678f1b9b4SToby Isaac #endif 162778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n")); 162878f1b9b4SToby Isaac } 162978f1b9b4SToby Isaac } 163078f1b9b4SToby Isaac } 163178f1b9b4SToby Isaac } 163278f1b9b4SToby Isaac 163378f1b9b4SToby Isaac /* Memory usage and object creation */ 163478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------")); 163578f1b9b4SToby Isaac if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, "-----------------------------")); 163678f1b9b4SToby Isaac #if defined(PETSC_HAVE_DEVICE) 163778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "---------------------------------------")); 163878f1b9b4SToby Isaac #endif 163978f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n")); 164078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n")); 164178f1b9b4SToby Isaac 164278f1b9b4SToby Isaac /* Right now, only stages on the first processor are reported here, meaning only objects associated with 164378f1b9b4SToby Isaac the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then 164478f1b9b4SToby Isaac stats for stages local to processor sets. 164578f1b9b4SToby Isaac */ 164678f1b9b4SToby Isaac /* We should figure out the longest object name here (now 20 characters) */ 164778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Object Type Creations Destructions. Reports information only for process 0.\n")); 164878f1b9b4SToby Isaac for (stage = 0; stage < numStages; stage++) { 164978f1b9b4SToby Isaac const char *stage_name; 165078f1b9b4SToby Isaac 165178f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 165278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stage_name)); 165378f1b9b4SToby Isaac if (localStageUsed[stage]) { 165478f1b9b4SToby Isaac PetscInt num_classes; 165578f1b9b4SToby Isaac 165678f1b9b4SToby Isaac PetscCall(PetscLogStateGetNumClasses(state, &num_classes)); 165778f1b9b4SToby Isaac for (oclass = 0; oclass < num_classes; oclass++) { 165878f1b9b4SToby Isaac PetscClassPerf *class_perf_info; 165978f1b9b4SToby Isaac 166078f1b9b4SToby Isaac PetscCall(PetscLogHandlerDefaultGetClassPerf(handler, stage, oclass, &class_perf_info)); 166178f1b9b4SToby Isaac if ((class_perf_info->creations > 0) || (class_perf_info->destructions > 0)) { 166278f1b9b4SToby Isaac PetscLogClassInfo class_reg_info; 166378f1b9b4SToby Isaac 166478f1b9b4SToby Isaac PetscCall(PetscLogStateClassGetInfo(state, oclass, &class_reg_info)); 166578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "%20s %5d %5d\n", class_reg_info.name, class_perf_info->creations, class_perf_info->destructions)); 166678f1b9b4SToby Isaac } 166778f1b9b4SToby Isaac } 166878f1b9b4SToby Isaac } 166978f1b9b4SToby Isaac } 167078f1b9b4SToby Isaac 167178f1b9b4SToby Isaac PetscCall(PetscFree(localStageUsed)); 167278f1b9b4SToby Isaac PetscCall(PetscFree(stageUsed)); 167378f1b9b4SToby Isaac PetscCall(PetscFree(localStageVisible)); 167478f1b9b4SToby Isaac PetscCall(PetscFree(stageVisible)); 167578f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesDestroy(&global_stages)); 167678f1b9b4SToby Isaac PetscCall(PetscLogGlobalNamesDestroy(&global_events)); 167778f1b9b4SToby Isaac 167878f1b9b4SToby Isaac /* Information unrelated to this particular run */ 167978f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "========================================================================================================================\n")); 168078f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 168178f1b9b4SToby Isaac PetscCall(PetscTime(&x)); 168278f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 168378f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 168478f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 168578f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 168678f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 168778f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 168878f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 168978f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 169078f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 169178f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 169278f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Average time to get PetscTime(): %g\n", (y - x) / 10.0)); 169378f1b9b4SToby Isaac /* MPI information */ 169478f1b9b4SToby Isaac if (size > 1) { 169578f1b9b4SToby Isaac MPI_Status status; 169678f1b9b4SToby Isaac PetscMPIInt tag; 169778f1b9b4SToby Isaac MPI_Comm newcomm; 169878f1b9b4SToby Isaac 169978f1b9b4SToby Isaac PetscCallMPI(MPI_Barrier(comm)); 170078f1b9b4SToby Isaac PetscCall(PetscTime(&x)); 170178f1b9b4SToby Isaac PetscCallMPI(MPI_Barrier(comm)); 170278f1b9b4SToby Isaac PetscCallMPI(MPI_Barrier(comm)); 170378f1b9b4SToby Isaac PetscCallMPI(MPI_Barrier(comm)); 170478f1b9b4SToby Isaac PetscCallMPI(MPI_Barrier(comm)); 170578f1b9b4SToby Isaac PetscCallMPI(MPI_Barrier(comm)); 170678f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 170778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Average time for MPI_Barrier(): %g\n", (y - x) / 5.0)); 170878f1b9b4SToby Isaac PetscCall(PetscCommDuplicate(comm, &newcomm, &tag)); 170978f1b9b4SToby Isaac PetscCallMPI(MPI_Barrier(comm)); 171078f1b9b4SToby Isaac if (rank) { 171178f1b9b4SToby Isaac PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, rank - 1, tag, newcomm, &status)); 171278f1b9b4SToby Isaac PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, (rank + 1) % size, tag, newcomm)); 171378f1b9b4SToby Isaac } else { 171478f1b9b4SToby Isaac PetscCall(PetscTime(&x)); 171578f1b9b4SToby Isaac PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, 1, tag, newcomm)); 171678f1b9b4SToby Isaac PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, size - 1, tag, newcomm, &status)); 171778f1b9b4SToby Isaac PetscCall(PetscTime(&y)); 171878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Average time for zero size MPI_Send(): %g\n", (y - x) / size)); 171978f1b9b4SToby Isaac } 172078f1b9b4SToby Isaac PetscCall(PetscCommDestroy(&newcomm)); 172178f1b9b4SToby Isaac } 172278f1b9b4SToby Isaac PetscCall(PetscOptionsView(NULL, viewer)); 172378f1b9b4SToby Isaac 172478f1b9b4SToby Isaac /* Machine and compile information */ 172578f1b9b4SToby Isaac if (PetscDefined(USE_FORTRAN_KERNELS)) { 172678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Compiled with FORTRAN kernels\n")); 172778f1b9b4SToby Isaac } else { 172878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Compiled without FORTRAN kernels\n")); 172978f1b9b4SToby Isaac } 173078f1b9b4SToby Isaac if (PetscDefined(USE_64BIT_INDICES)) { 173178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Compiled with 64-bit PetscInt\n")); 173278f1b9b4SToby Isaac } else if (PetscDefined(USE___FLOAT128)) { 173378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Compiled with 32-bit PetscInt\n")); 173478f1b9b4SToby Isaac } 173578f1b9b4SToby Isaac if (PetscDefined(USE_REAL_SINGLE)) { 173678f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Compiled with single precision PetscScalar and PetscReal\n")); 173778f1b9b4SToby Isaac } else if (PetscDefined(USE___FLOAT128)) { 173878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Compiled with 128 bit precision PetscScalar and PetscReal\n")); 173978f1b9b4SToby Isaac } 174078f1b9b4SToby Isaac if (PetscDefined(USE_REAL_MAT_SINGLE)) { 174178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Compiled with single precision matrices\n")); 174278f1b9b4SToby Isaac } else { 174378f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Compiled with full precision matrices (default)\n")); 174478f1b9b4SToby Isaac } 174578f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "sizeof(short) %d sizeof(int) %d sizeof(long) %d sizeof(void*) %d sizeof(PetscScalar) %d sizeof(PetscInt) %d\n", (int)sizeof(short), (int)sizeof(int), (int)sizeof(long), (int)sizeof(void *), (int)sizeof(PetscScalar), (int)sizeof(PetscInt))); 174678f1b9b4SToby Isaac 174778f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "Configure options: %s", petscconfigureoptions)); 174878f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "%s", petscmachineinfo)); 174978f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "%s", petsccompilerinfo)); 175078f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "%s", petsccompilerflagsinfo)); 175178f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "%s", petsclinkerinfo)); 175278f1b9b4SToby Isaac 175378f1b9b4SToby Isaac /* Cleanup */ 175478f1b9b4SToby Isaac PetscCall(PetscFPrintf(comm, fd, "\n")); 175578f1b9b4SToby Isaac PetscCall(PetscLogViewWarnNoGpuAwareMpi(comm, fd)); 175678f1b9b4SToby Isaac PetscCall(PetscLogViewWarnDebugging(comm, fd)); 175778f1b9b4SToby Isaac PetscCall(PetscFPTrapPop()); 175878f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 175978f1b9b4SToby Isaac } 176078f1b9b4SToby Isaac 176178f1b9b4SToby Isaac static PetscErrorCode PetscLogHandlerView_Default(PetscLogHandler handler, PetscViewer viewer) 176278f1b9b4SToby Isaac { 176378f1b9b4SToby Isaac PetscViewerFormat format; 176478f1b9b4SToby Isaac PetscFunctionBegin; 176578f1b9b4SToby Isaac PetscCall(PetscViewerGetFormat(viewer, &format)); 176678f1b9b4SToby Isaac if (format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_ASCII_INFO) { 176778f1b9b4SToby Isaac PetscCall(PetscLogHandlerView_Default_Info(handler, viewer)); 176878f1b9b4SToby Isaac } else if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 176978f1b9b4SToby Isaac PetscCall(PetscLogHandlerView_Default_Detailed(handler, viewer)); 177078f1b9b4SToby Isaac } else if (format == PETSC_VIEWER_ASCII_CSV) { 177178f1b9b4SToby Isaac PetscCall(PetscLogHandlerView_Default_CSV(handler, viewer)); 177278f1b9b4SToby Isaac } 177378f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 177478f1b9b4SToby Isaac } 177578f1b9b4SToby Isaac 177678f1b9b4SToby Isaac /*MC 1777*294de794SToby Isaac PETSCLOGHANDLERDEFAULT - PETSCLOGHANDLERDEFAULT = "default" - A `PetscLogHandler` that collects data for PETSc 1778b665b14eSToby Isaac default profiling log viewers (`PetscLogView()` and `PetscLogDump()`). A log handler of this type is 1779b665b14eSToby Isaac created and started (`PetscLogHandlerStart()`) by `PetscLogDefaultBegin()`. 178078f1b9b4SToby Isaac 178178f1b9b4SToby Isaac Options Database Keys: 178278f1b9b4SToby Isaac + -log_include_actions - include a growing list of actions (event beginnings and endings, object creations and destructions) in `PetscLogDump()` (`PetscLogActions()`). 178378f1b9b4SToby Isaac - -log_include_objects - include a growing list of object creations and destructions in `PetscLogDump()` (`PetscLogObjects()`). 178478f1b9b4SToby Isaac 178578f1b9b4SToby Isaac Level: developer 178678f1b9b4SToby Isaac 178778f1b9b4SToby Isaac .seealso: [](ch_profiling), `PetscLogHandler` 178878f1b9b4SToby Isaac M*/ 178978f1b9b4SToby Isaac 179078f1b9b4SToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Default(PetscLogHandler handler) 179178f1b9b4SToby Isaac { 179278f1b9b4SToby Isaac PetscFunctionBegin; 179378f1b9b4SToby Isaac PetscCall(PetscLogHandlerContextCreate_Default((PetscLogHandler_Default *)&handler->data)); 179478f1b9b4SToby Isaac handler->ops->destroy = PetscLogHandlerDestroy_Default; 179578f1b9b4SToby Isaac handler->ops->eventbegin = PetscLogHandlerEventBegin_Default; 179678f1b9b4SToby Isaac handler->ops->eventend = PetscLogHandlerEventEnd_Default; 179778f1b9b4SToby Isaac handler->ops->eventsync = PetscLogHandlerEventSync_Default; 179878f1b9b4SToby Isaac handler->ops->objectcreate = PetscLogHandlerObjectCreate_Default; 179978f1b9b4SToby Isaac handler->ops->objectdestroy = PetscLogHandlerObjectDestroy_Default; 180078f1b9b4SToby Isaac handler->ops->stagepush = PetscLogHandlerStagePush_Default; 180178f1b9b4SToby Isaac handler->ops->stagepop = PetscLogHandlerStagePop_Default; 180278f1b9b4SToby Isaac handler->ops->view = PetscLogHandlerView_Default; 180378f1b9b4SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 180478f1b9b4SToby Isaac } 1805