1 #include <petscviewer.h> 2 #include <petscdevice.h> 3 #include <petsc/private/logimpl.h> /*I "petscsys.h" I*/ 4 #include <petsc/private/loghandlerimpl.h> 5 #include <petsc/private/deviceimpl.h> 6 #include <petscconfiginfo.h> 7 #include <petscmachineinfo.h> 8 #include "logdefault.h" 9 10 static PetscErrorCode PetscEventPerfInfoInit(PetscEventPerfInfo *eventInfo) 11 { 12 PetscFunctionBegin; 13 PetscCall(PetscMemzero(eventInfo, sizeof(*eventInfo))); 14 eventInfo->visible = PETSC_TRUE; 15 eventInfo->id = -1; 16 eventInfo->dof[0] = -1.0; 17 eventInfo->dof[1] = -1.0; 18 eventInfo->dof[2] = -1.0; 19 eventInfo->dof[3] = -1.0; 20 eventInfo->dof[4] = -1.0; 21 eventInfo->dof[5] = -1.0; 22 eventInfo->dof[6] = -1.0; 23 eventInfo->dof[7] = -1.0; 24 eventInfo->errors[0] = -1.0; 25 eventInfo->errors[1] = -1.0; 26 eventInfo->errors[2] = -1.0; 27 eventInfo->errors[3] = -1.0; 28 eventInfo->errors[4] = -1.0; 29 eventInfo->errors[5] = -1.0; 30 eventInfo->errors[6] = -1.0; 31 eventInfo->errors[7] = -1.0; 32 PetscFunctionReturn(PETSC_SUCCESS); 33 } 34 35 static PetscErrorCode PetscEventPerfInfoTic_Internal(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event, PetscBool resume) 36 { 37 PetscFunctionBegin; 38 if (resume) { 39 eventInfo->timeTmp -= time; 40 eventInfo->flopsTmp -= petsc_TotalFlops_th; 41 } else { 42 eventInfo->timeTmp = -time; 43 eventInfo->flopsTmp = -petsc_TotalFlops_th; 44 } 45 eventInfo->numMessages -= petsc_irecv_ct_th + petsc_isend_ct_th + petsc_recv_ct_th + petsc_send_ct_th; 46 eventInfo->messageLength -= petsc_irecv_len_th + petsc_isend_len_th + petsc_recv_len_th + petsc_send_len_th; 47 eventInfo->numReductions -= petsc_allreduce_ct_th + petsc_gather_ct_th + petsc_scatter_ct_th; 48 #if defined(PETSC_HAVE_DEVICE) 49 eventInfo->CpuToGpuCount -= petsc_ctog_ct_th; 50 eventInfo->GpuToCpuCount -= petsc_gtoc_ct_th; 51 eventInfo->CpuToGpuSize -= petsc_ctog_sz_th; 52 eventInfo->GpuToCpuSize -= petsc_gtoc_sz_th; 53 eventInfo->GpuFlops -= petsc_gflops_th; 54 eventInfo->GpuTime -= petsc_gtime; 55 #endif 56 if (logMemory) { 57 PetscLogDouble usage; 58 PetscCall(PetscMemoryGetCurrentUsage(&usage)); 59 eventInfo->memIncrease -= usage; 60 PetscCall(PetscMallocGetCurrentUsage(&usage)); 61 eventInfo->mallocSpace -= usage; 62 PetscCall(PetscMallocGetMaximumUsage(&usage)); 63 eventInfo->mallocIncrease -= usage; 64 PetscCall(PetscMallocPushMaximumUsage(event)); 65 } 66 PetscFunctionReturn(PETSC_SUCCESS); 67 } 68 69 static PetscErrorCode PetscEventPerfInfoTic(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event) 70 { 71 PetscFunctionBegin; 72 PetscCall(PetscEventPerfInfoTic_Internal(eventInfo, time, logMemory, event, PETSC_FALSE)); 73 PetscFunctionReturn(PETSC_SUCCESS); 74 } 75 76 static PetscErrorCode PetscEventPerfInfoResume(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event) 77 { 78 PetscFunctionBegin; 79 PetscCall(PetscEventPerfInfoTic_Internal(eventInfo, time, logMemory, event, PETSC_TRUE)); 80 PetscFunctionReturn(PETSC_SUCCESS); 81 } 82 83 static PetscErrorCode PetscEventPerfInfoToc_Internal(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event, PetscBool pause) 84 { 85 PetscFunctionBegin; 86 eventInfo->timeTmp += time; 87 eventInfo->flopsTmp += petsc_TotalFlops_th; 88 if (!pause) { 89 eventInfo->time += eventInfo->timeTmp; 90 eventInfo->time2 += eventInfo->timeTmp * eventInfo->timeTmp; 91 eventInfo->flops += eventInfo->flopsTmp; 92 eventInfo->flops2 += eventInfo->flopsTmp * eventInfo->flopsTmp; 93 } 94 eventInfo->numMessages += petsc_irecv_ct_th + petsc_isend_ct_th + petsc_recv_ct_th + petsc_send_ct_th; 95 eventInfo->messageLength += petsc_irecv_len_th + petsc_isend_len_th + petsc_recv_len + petsc_send_len_th; 96 eventInfo->numReductions += petsc_allreduce_ct_th + petsc_gather_ct_th + petsc_scatter_ct_th; 97 #if defined(PETSC_HAVE_DEVICE) 98 eventInfo->CpuToGpuCount += petsc_ctog_ct_th; 99 eventInfo->GpuToCpuCount += petsc_gtoc_ct_th; 100 eventInfo->CpuToGpuSize += petsc_ctog_sz_th; 101 eventInfo->GpuToCpuSize += petsc_gtoc_sz_th; 102 eventInfo->GpuFlops += petsc_gflops_th; 103 eventInfo->GpuTime += petsc_gtime; 104 #endif 105 if (logMemory) { 106 PetscLogDouble usage, musage; 107 PetscCall(PetscMemoryGetCurrentUsage(&usage)); /* the comments below match the column labels printed in PetscLogView_Default() */ 108 eventInfo->memIncrease += usage; /* RMI */ 109 PetscCall(PetscMallocGetCurrentUsage(&usage)); 110 eventInfo->mallocSpace += usage; /* Malloc */ 111 PetscCall(PetscMallocPopMaximumUsage((int)event, &musage)); 112 eventInfo->mallocIncreaseEvent = PetscMax(musage - usage, eventInfo->mallocIncreaseEvent); /* EMalloc */ 113 PetscCall(PetscMallocGetMaximumUsage(&usage)); 114 eventInfo->mallocIncrease += usage; /* MMalloc */ 115 } 116 PetscFunctionReturn(PETSC_SUCCESS); 117 } 118 119 static PetscErrorCode PetscEventPerfInfoToc(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event) 120 { 121 PetscFunctionBegin; 122 PetscCall(PetscEventPerfInfoToc_Internal(eventInfo, time, logMemory, event, PETSC_FALSE)); 123 PetscFunctionReturn(PETSC_SUCCESS); 124 } 125 126 static PetscErrorCode PetscEventPerfInfoPause(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event) 127 { 128 PetscFunctionBegin; 129 PetscCall(PetscEventPerfInfoToc_Internal(eventInfo, time, logMemory, event, PETSC_FALSE)); 130 PetscFunctionReturn(PETSC_SUCCESS); 131 } 132 133 static PetscErrorCode PetscEventPerfInfoAdd_Internal(const PetscEventPerfInfo *eventInfo, PetscEventPerfInfo *outInfo) 134 { 135 PetscFunctionBegin; 136 outInfo->count += eventInfo->count; 137 outInfo->time += eventInfo->time; 138 outInfo->time2 += eventInfo->time2; 139 outInfo->flops += eventInfo->flops; 140 outInfo->flops2 += eventInfo->flops2; 141 outInfo->numMessages += eventInfo->numMessages; 142 outInfo->messageLength += eventInfo->messageLength; 143 outInfo->numReductions += eventInfo->numReductions; 144 #if defined(PETSC_HAVE_DEVICE) 145 outInfo->CpuToGpuCount += eventInfo->CpuToGpuCount; 146 outInfo->GpuToCpuCount += eventInfo->GpuToCpuCount; 147 outInfo->CpuToGpuSize += eventInfo->CpuToGpuSize; 148 outInfo->GpuToCpuSize += eventInfo->GpuToCpuSize; 149 outInfo->GpuFlops += eventInfo->GpuFlops; 150 outInfo->GpuTime += eventInfo->GpuTime; 151 #endif 152 outInfo->memIncrease += eventInfo->memIncrease; 153 outInfo->mallocSpace += eventInfo->mallocSpace; 154 outInfo->mallocIncreaseEvent += eventInfo->mallocIncreaseEvent; 155 outInfo->mallocIncrease += eventInfo->mallocIncrease; 156 PetscFunctionReturn(PETSC_SUCCESS); 157 } 158 159 PETSC_LOG_RESIZABLE_ARRAY(EventPerfArray, PetscEventPerfInfo, PetscLogEvent, PetscEventPerfInfoInit, NULL, NULL) 160 161 /* --- PetscClassPerf --- */ 162 163 typedef struct { 164 int creations; /* The number of objects of this class created */ 165 int destructions; /* The number of objects of this class destroyed */ 166 PetscLogDouble mem; /* The total memory allocated by objects of this class; this is completely wrong and should possibly be removed */ 167 PetscLogDouble descMem; /* The total memory allocated by descendents of these objects; this is completely wrong and should possibly be removed */ 168 } PetscClassPerf; 169 170 static PetscErrorCode PetscClassPerfInit(PetscClassPerf *classInfo) 171 { 172 PetscFunctionBegin; 173 PetscCall(PetscMemzero(classInfo, sizeof(*classInfo))); 174 PetscFunctionReturn(PETSC_SUCCESS); 175 } 176 177 PETSC_LOG_RESIZABLE_ARRAY(ClassPerfArray, PetscClassPerf, PetscLogClass, PetscClassPerfInit, NULL, NULL) 178 179 /* --- PetscStagePerf --- */ 180 181 typedef struct _PetscStagePerf { 182 PetscBool used; /* The stage was pushed on this processor */ 183 PetscEventPerfInfo perfInfo; /* The stage performance information */ 184 PetscLogEventPerfArray eventLog; /* The event information for this stage */ 185 PetscLogClassPerfArray classLog; /* The class information for this stage */ 186 } PetscStagePerf; 187 188 static PetscErrorCode PetscStageInfoInit(PetscStagePerf *stageInfo) 189 { 190 PetscFunctionBegin; 191 PetscCall(PetscMemzero(stageInfo, sizeof(*stageInfo))); 192 PetscCall(PetscLogEventPerfArrayCreate(128, &(stageInfo->eventLog))); 193 PetscCall(PetscLogClassPerfArrayCreate(128, &(stageInfo->classLog))); 194 PetscCall(PetscEventPerfInfoInit(&stageInfo->perfInfo)); 195 PetscFunctionReturn(PETSC_SUCCESS); 196 } 197 198 static PetscErrorCode PetscStageInfoReset(PetscStagePerf *stageInfo) 199 { 200 PetscFunctionBegin; 201 PetscCall(PetscLogEventPerfArrayDestroy(&(stageInfo->eventLog))); 202 PetscCall(PetscLogClassPerfArrayDestroy(&(stageInfo->classLog))); 203 PetscFunctionReturn(PETSC_SUCCESS); 204 } 205 206 PETSC_LOG_RESIZABLE_ARRAY(StageInfoArray, PetscStagePerf, PetscLogStage, PetscStageInfoInit, PetscStageInfoReset, NULL) 207 208 /* --- Action --- */ 209 210 /* The structure for action logging */ 211 typedef enum { 212 PETSC_LOG_ACTION_CREATE, 213 PETSC_LOG_ACTION_DESTROY, 214 PETSC_LOG_ACTION_BEGIN, 215 PETSC_LOG_ACTION_END, 216 } PetscLogActionType; 217 218 typedef struct _Action { 219 PetscLogActionType action; /* The type of execution */ 220 PetscLogEvent event; /* The event number */ 221 PetscClassId classid; /* The event class id */ 222 PetscLogDouble time; /* The time of occurrence */ 223 PetscLogDouble flops; /* The cumulative flops */ 224 PetscLogDouble mem; /* The current memory usage */ 225 PetscLogDouble maxmem; /* The maximum memory usage */ 226 int id1, id2, id3; /* The ids of associated objects */ 227 } Action; 228 229 PETSC_LOG_RESIZABLE_ARRAY(ActionArray, Action, PetscLogEvent, NULL, NULL, NULL) 230 231 /* --- Object --- */ 232 233 /* The structure for object logging */ 234 typedef struct _Object { 235 PetscObject obj; /* The associated PetscObject */ 236 int parent; /* The parent id */ 237 PetscLogDouble mem; /* The memory associated with the object */ 238 char name[64]; /* The object name */ 239 char info[64]; /* The information string */ 240 } Object; 241 242 PETSC_LOG_RESIZABLE_ARRAY(ObjectArray, Object, PetscObject, NULL, NULL, NULL) 243 244 /* Map from (threadid,stage,event) to perfInfo data struct */ 245 #include <petsc/private/hashmapijk.h> 246 247 PETSC_HASH_MAP(HMapEvent, PetscHashIJKKey, PetscEventPerfInfo *, PetscHashIJKKeyHash, PetscHashIJKKeyEqual, NULL) 248 249 typedef struct _n_PetscLogHandler_Default *PetscLogHandler_Default; 250 struct _n_PetscLogHandler_Default { 251 PetscLogStageInfoArray stages; 252 PetscSpinlock lock; 253 PetscLogActionArray petsc_actions; 254 PetscLogObjectArray petsc_objects; 255 PetscBool petsc_logActions; 256 PetscBool petsc_logObjects; 257 int petsc_numObjectsCreated; 258 int petsc_numObjectsDestroyed; 259 PetscHMapEvent eventInfoMap_th; 260 int pause_depth; 261 PetscBool use_threadsafe; 262 }; 263 264 /* --- PetscLogHandler_Default --- */ 265 266 static PetscErrorCode PetscLogHandlerContextCreate_Default(PetscLogHandler_Default *def_p) 267 { 268 PetscLogHandler_Default def; 269 270 PetscFunctionBegin; 271 PetscCall(PetscNew(def_p)); 272 def = *def_p; 273 PetscCall(PetscLogStageInfoArrayCreate(8, &(def->stages))); 274 PetscCall(PetscLogActionArrayCreate(64, &(def->petsc_actions))); 275 PetscCall(PetscLogObjectArrayCreate(64, &(def->petsc_objects))); 276 277 PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_include_actions", &def->petsc_logActions, NULL)); 278 PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_include_objects", &def->petsc_logObjects, NULL)); 279 PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_handler_default_use_threadsafe_events", &def->use_threadsafe, NULL)); 280 if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { PetscCall(PetscHMapEventCreate(&def->eventInfoMap_th)); } 281 PetscFunctionReturn(PETSC_SUCCESS); 282 } 283 284 static PetscErrorCode PetscLogHandlerDestroy_Default(PetscLogHandler h) 285 { 286 PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 287 288 PetscFunctionBegin; 289 PetscCall(PetscLogStageInfoArrayDestroy(&def->stages)); 290 PetscCall(PetscLogActionArrayDestroy(&def->petsc_actions)); 291 PetscCall(PetscLogObjectArrayDestroy(&def->petsc_objects)); 292 if (def->eventInfoMap_th) { 293 PetscEventPerfInfo **array; 294 PetscInt n, off = 0; 295 296 PetscCall(PetscHMapEventGetSize(def->eventInfoMap_th, &n)); 297 PetscCall(PetscMalloc1(n, &array)); 298 PetscCall(PetscHMapEventGetVals(def->eventInfoMap_th, &off, array)); 299 for (PetscInt i = 0; i < n; i++) PetscCall(PetscFree(array[i])); 300 PetscCall(PetscFree(array)); 301 PetscCall(PetscHMapEventDestroy(&def->eventInfoMap_th)); 302 } 303 PetscCall(PetscFree(def)); 304 PetscFunctionReturn(PETSC_SUCCESS); 305 } 306 307 static PetscErrorCode PetscLogHandlerDefaultGetStageInfo(PetscLogHandler handler, PetscLogStage stage, PetscStagePerf **stage_info_p) 308 { 309 PetscStagePerf *stage_info = NULL; 310 PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 311 312 PetscFunctionBegin; 313 PetscCall(PetscLogStageInfoArrayResize(def->stages, stage + 1)); 314 PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info)); 315 *stage_info_p = stage_info; 316 PetscFunctionReturn(PETSC_SUCCESS); 317 } 318 319 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultGetEventPerfInfo(PetscLogHandler handler, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **event_info_p) 320 { 321 PetscEventPerfInfo *event_info = NULL; 322 PetscStagePerf *stage_info = NULL; 323 PetscLogEventPerfArray event_log; 324 325 PetscFunctionBegin; 326 if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(handler->state, &stage)); 327 PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info)); 328 event_log = stage_info->eventLog; 329 PetscCall(PetscLogEventPerfArrayResize(event_log, event + 1)); 330 PetscCall(PetscLogEventPerfArrayGetRef(event_log, event, &event_info)); 331 event_info->id = event; 332 *event_info_p = event_info; 333 PetscFunctionReturn(PETSC_SUCCESS); 334 } 335 336 static PetscErrorCode PetscLogHandlerDefaultGetClassPerf(PetscLogHandler handler, PetscLogStage stage, PetscLogClass clss, PetscClassPerf **class_info) 337 { 338 PetscLogClassPerfArray class_log; 339 PetscStagePerf *stage_info; 340 341 PetscFunctionBegin; 342 PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info)); 343 class_log = stage_info->classLog; 344 PetscCall(PetscLogClassPerfArrayResize(class_log, clss + 1)); 345 PetscCall(PetscLogClassPerfArrayGetRef(class_log, clss, class_info)); 346 PetscFunctionReturn(PETSC_SUCCESS); 347 } 348 349 static PetscErrorCode PetscLogHandlerObjectCreate_Default(PetscLogHandler h, PetscObject obj) 350 { 351 PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 352 PetscLogState state; 353 PetscLogStage stage; 354 PetscClassPerf *classInfo; 355 int oclass = 0; 356 357 PetscFunctionBegin; 358 PetscCall(PetscLogHandlerGetState(h, &state)); 359 PetscCall(PetscSpinlockLock(&def->lock)); 360 /* Record stage info */ 361 PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 362 PetscCall(PetscLogStateGetClassFromClassId(state, obj->classid, &oclass)); 363 PetscCall(PetscLogHandlerDefaultGetClassPerf(h, stage, oclass, &classInfo)); 364 classInfo->creations++; 365 /* Record the creation action */ 366 if (def->petsc_logActions) { 367 Action new_action; 368 369 PetscCall(PetscTime(&new_action.time)); 370 new_action.time -= petsc_BaseTime; 371 new_action.action = PETSC_LOG_ACTION_CREATE; 372 new_action.event = -1; 373 new_action.classid = obj->classid; 374 new_action.id1 = obj->id; 375 new_action.id2 = -1; 376 new_action.id3 = -1; 377 new_action.flops = petsc_TotalFlops; 378 PetscCall(PetscMallocGetCurrentUsage(&new_action.mem)); 379 PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem)); 380 PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action)); 381 } 382 /* We don't just use obj->id to count all objects that are created 383 because PetscLogHandlers are objects and PetscLogObjectDestroy() will not 384 be called for them: the number of objects created and destroyed as counted 385 here and below would have an imbalance */ 386 def->petsc_numObjectsCreated++; 387 /* Record the object */ 388 if (def->petsc_logObjects) { 389 Object new_object; 390 391 new_object.parent = -1; 392 new_object.obj = obj; 393 new_object.mem = 0; 394 395 PetscCall(PetscMemzero(new_object.name, sizeof(new_object.name))); 396 PetscCall(PetscMemzero(new_object.info, sizeof(new_object.info))); 397 PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1"); 398 PetscCall(PetscLogObjectArrayResize(def->petsc_objects, obj->id)); 399 PetscCall(PetscLogObjectArraySet(def->petsc_objects, obj->id - 1, new_object)); 400 } 401 PetscCall(PetscSpinlockUnlock(&def->lock)); 402 PetscFunctionReturn(PETSC_SUCCESS); 403 } 404 405 static PetscErrorCode PetscLogHandlerObjectDestroy_Default(PetscLogHandler h, PetscObject obj) 406 { 407 PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 408 PetscLogState state; 409 PetscLogStage stage; 410 PetscClassPerf *classInfo; 411 int oclass = 0; 412 413 PetscFunctionBegin; 414 PetscCall(PetscLogHandlerGetState(h, &state)); 415 /* Record stage info */ 416 PetscCall(PetscSpinlockLock(&def->lock)); 417 PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 418 if (stage >= 0) { 419 /* stage < 0 can happen if the log summary is output before some things are destroyed */ 420 PetscCall(PetscLogStateGetClassFromClassId(state, obj->classid, &oclass)); 421 PetscCall(PetscLogHandlerDefaultGetClassPerf(h, stage, oclass, &classInfo)); 422 classInfo->destructions++; 423 } 424 /* Cannot Credit all ancestors with your memory because they may have already been destroyed*/ 425 def->petsc_numObjectsDestroyed++; 426 /* Dynamically enlarge logging structures */ 427 /* Record the destruction action */ 428 if (def->petsc_logActions) { 429 Action new_action; 430 431 PetscCall(PetscTime(&new_action.time)); 432 new_action.time -= petsc_BaseTime; 433 new_action.event = -1; 434 new_action.action = PETSC_LOG_ACTION_DESTROY; 435 new_action.classid = obj->classid; 436 new_action.id1 = obj->id; 437 new_action.id2 = -1; 438 new_action.id3 = -1; 439 new_action.flops = petsc_TotalFlops; 440 PetscCall(PetscMallocGetCurrentUsage(&new_action.mem)); 441 PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem)); 442 PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action)); 443 } 444 if (def->petsc_logObjects) { 445 Object *obj_entry = NULL; 446 447 PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1"); 448 PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, obj->id - 1, &obj_entry)); 449 if (obj->name) PetscCall(PetscStrncpy(obj_entry->name, obj->name, 64)); 450 obj_entry->obj = NULL; 451 } 452 PetscCall(PetscSpinlockUnlock(&def->lock)); 453 PetscFunctionReturn(PETSC_SUCCESS); 454 } 455 456 static PetscErrorCode PetscLogHandlerEventSync_Default(PetscLogHandler h, PetscLogEvent event, MPI_Comm comm) 457 { 458 PetscLogState state; 459 PetscLogEventInfo event_info; 460 PetscEventPerfInfo *event_perf_info; 461 int stage; 462 PetscLogDouble time = 0.0; 463 464 PetscFunctionBegin; 465 if (!(PetscLogSyncOn) || comm == MPI_COMM_NULL) PetscFunctionReturn(PETSC_SUCCESS); 466 PetscCall(PetscLogHandlerGetState(h, &state)); 467 PetscCall(PetscLogStateEventGetInfo(state, event, &event_info)); 468 if (!event_info.collective) PetscFunctionReturn(PETSC_SUCCESS); 469 PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 470 PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info)); 471 if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS); 472 473 PetscCall(PetscTimeSubtract(&time)); 474 PetscCallMPI(MPI_Barrier(comm)); 475 PetscCall(PetscTimeAdd(&time)); 476 event_perf_info->syncTime += time; 477 PetscFunctionReturn(PETSC_SUCCESS); 478 } 479 480 static PetscErrorCode PetscLogGetStageEventPerfInfo_threaded(PetscLogHandler_Default def, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **eventInfo) 481 { 482 PetscEventPerfInfo *leventInfo = NULL; 483 PetscHashIJKKey key; 484 485 PetscFunctionBegin; 486 487 #if PetscDefined(HAVE_THREADSAFETY) 488 key.i = PetscLogGetTid(); 489 #else 490 key.i = 0; 491 #endif 492 key.j = stage; 493 key.k = event; 494 PetscCall(PetscSpinlockLock(&def->lock)); 495 PetscCall(PetscHMapEventGet(def->eventInfoMap_th, key, &leventInfo)); 496 if (!leventInfo) { 497 PetscCall(PetscNew(&leventInfo)); 498 leventInfo->id = event; 499 PetscCall(PetscHMapEventSet(def->eventInfoMap_th, key, leventInfo)); 500 } 501 PetscCall(PetscSpinlockUnlock(&def->lock)); 502 *eventInfo = leventInfo; 503 PetscFunctionReturn(PETSC_SUCCESS); 504 } 505 506 #if defined(PETSC_HAVE_CUDA) 507 #include <nvToolsExt.h> 508 #endif 509 static PetscErrorCode PetscLogHandlerEventBegin_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4) 510 { 511 PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 512 PetscEventPerfInfo *event_perf_info = NULL; 513 PetscLogEventInfo event_info; 514 PetscLogDouble time; 515 PetscLogState state; 516 PetscLogStage stage; 517 518 PetscFunctionBegin; 519 PetscCall(PetscLogHandlerGetState(h, &state)); 520 /* Synchronization */ 521 PetscCall(PetscLogHandlerEventSync_Default(h, event, PetscObjectComm(o1))); 522 PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 523 if (def->pause_depth > 0) stage = 0; // in pause-mode, all events run on the main stage 524 if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { 525 PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info)); 526 if (event_perf_info->depth == 0) { PetscCall(PetscEventPerfInfoInit(event_perf_info)); } 527 } else { 528 PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info)); 529 } 530 PetscCheck(event_perf_info->depth >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to begin a paused event, this is not allowed"); 531 event_perf_info->depth++; 532 /* Check for double counting */ 533 if (event_perf_info->depth > 1) PetscFunctionReturn(PETSC_SUCCESS); 534 #if defined(PETSC_HAVE_CUDA) 535 if (PetscDeviceInitialized(PETSC_DEVICE_CUDA)) { 536 PetscLogEventInfo event_reg_info; 537 538 PetscCall(PetscLogStateEventGetInfo(state, event, &event_reg_info)); 539 nvtxRangePushA(event_reg_info.name); 540 } 541 #endif 542 PetscCall(PetscLogStateEventGetInfo(state, event, &event_info)); 543 /* Log the performance info */ 544 event_perf_info->count++; 545 PetscCall(PetscTime(&time)); 546 PetscCall(PetscEventPerfInfoTic(event_perf_info, time, PetscLogMemory, (int)event)); 547 if (def->petsc_logActions) { 548 PetscLogDouble curTime; 549 Action new_action; 550 551 PetscCall(PetscTime(&curTime)); 552 new_action.time = curTime - petsc_BaseTime; 553 new_action.action = PETSC_LOG_ACTION_BEGIN; 554 new_action.event = event; 555 new_action.classid = event_info.classid; 556 new_action.id1 = o1 ? o1->id : -1; 557 new_action.id2 = o2 ? o2->id : -1; 558 new_action.id3 = o3 ? o3->id : -1; 559 new_action.flops = petsc_TotalFlops; 560 PetscCall(PetscMallocGetCurrentUsage(&new_action.mem)); 561 PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem)); 562 PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action)); 563 } 564 PetscFunctionReturn(PETSC_SUCCESS); 565 } 566 567 static PetscErrorCode PetscLogHandlerEventEnd_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4) 568 { 569 PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 570 PetscEventPerfInfo *event_perf_info = NULL; 571 PetscLogDouble time; 572 PetscLogState state; 573 int stage; 574 575 PetscFunctionBegin; 576 PetscCall(PetscLogHandlerGetState(h, &state)); 577 if (def->petsc_logActions) { 578 PetscLogEventInfo event_info; 579 PetscLogDouble curTime; 580 Action new_action; 581 582 PetscCall(PetscLogStateEventGetInfo(state, event, &event_info)); 583 PetscCall(PetscTime(&curTime)); 584 new_action.time = curTime - petsc_BaseTime; 585 new_action.action = PETSC_LOG_ACTION_END; 586 new_action.event = event; 587 new_action.classid = event_info.classid; 588 new_action.id1 = o1 ? o1->id : -1; 589 new_action.id2 = o2 ? o2->id : -2; 590 new_action.id3 = o3 ? o3->id : -3; 591 new_action.flops = petsc_TotalFlops; 592 PetscCall(PetscMallocGetCurrentUsage(&new_action.mem)); 593 PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem)); 594 PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action)); 595 } 596 PetscCall(PetscLogStateGetCurrentStage(state, &stage)); 597 if (def->pause_depth > 0) stage = 0; // all events run on the main stage in pause-mode 598 if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { 599 PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info)); 600 } else { 601 PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info)); 602 } 603 PetscCheck(event_perf_info->depth > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to end paused event, not allowed"); 604 event_perf_info->depth--; 605 /* Check for double counting */ 606 if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS); 607 else PetscCheck(event_perf_info->depth == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Logging event had unbalanced begin/end pairs"); 608 609 /* Log performance info */ 610 PetscCall(PetscTime(&time)); 611 PetscCall(PetscEventPerfInfoToc(event_perf_info, time, PetscLogMemory, (int)event)); 612 if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { 613 PetscEventPerfInfo *event_perf_info_global; 614 PetscCall(PetscSpinlockLock(&def->lock)); 615 PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info_global)); 616 PetscCall(PetscEventPerfInfoAdd_Internal(event_perf_info, event_perf_info_global)); 617 PetscCall(PetscSpinlockUnlock(&def->lock)); 618 } 619 #if defined(PETSC_HAVE_CUDA) 620 if (PetscDeviceInitialized(PETSC_DEVICE_CUDA)) nvtxRangePop(); 621 #endif 622 PetscFunctionReturn(PETSC_SUCCESS); 623 } 624 625 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultDeactivatePush(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event) 626 { 627 PetscEventPerfInfo *event_perf_info; 628 629 PetscFunctionBegin; 630 if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage)); 631 PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info)); 632 event_perf_info->depth++; 633 PetscFunctionReturn(PETSC_SUCCESS); 634 } 635 636 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultDeactivatePop(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event) 637 { 638 PetscEventPerfInfo *event_perf_info; 639 640 PetscFunctionBegin; 641 if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage)); 642 PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(h, stage, event, &event_perf_info)); 643 event_perf_info->depth--; 644 PetscFunctionReturn(PETSC_SUCCESS); 645 } 646 647 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultEventsPause(PetscLogHandler h) 648 { 649 PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 650 PetscLogDouble time; 651 PetscInt num_stages; 652 653 PetscFunctionBegin; 654 if (def->pause_depth++ > 0) PetscFunctionReturn(PETSC_SUCCESS); 655 PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL)); 656 PetscCall(PetscTime(&time)); 657 for (PetscInt stage = 0; stage < num_stages; stage++) { 658 PetscStagePerf *stage_info = NULL; 659 PetscInt num_events; 660 661 PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info)); 662 PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL)); 663 for (PetscInt event = 0; event < num_events; event++) { 664 PetscEventPerfInfo *event_info = NULL; 665 PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info)); 666 if (event_info->depth > 0) { 667 event_info->depth *= -1; 668 PetscCall(PetscEventPerfInfoPause(event_info, time, PetscLogMemory, event)); 669 } 670 } 671 if (stage > 0 && stage_info->perfInfo.depth > 0) { 672 stage_info->perfInfo.depth *= -1; 673 PetscCall(PetscEventPerfInfoPause(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2))); 674 } 675 } 676 PetscFunctionReturn(PETSC_SUCCESS); 677 } 678 679 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultEventsResume(PetscLogHandler h) 680 { 681 PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 682 PetscLogDouble time; 683 PetscInt num_stages; 684 685 PetscFunctionBegin; 686 if (--def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS); 687 PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL)); 688 PetscCall(PetscTime(&time)); 689 for (PetscInt stage = 0; stage < num_stages; stage++) { 690 PetscStagePerf *stage_info = NULL; 691 PetscInt num_events; 692 693 PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info)); 694 PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL)); 695 for (PetscInt event = 0; event < num_events; event++) { 696 PetscEventPerfInfo *event_info = NULL; 697 PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info)); 698 if (event_info->depth < 0) { 699 event_info->depth *= -1; 700 PetscCall(PetscEventPerfInfoResume(event_info, time, PetscLogMemory, event)); 701 } 702 } 703 if (stage > 0 && stage_info->perfInfo.depth < 0) { 704 stage_info->perfInfo.depth *= -1; 705 PetscCall(PetscEventPerfInfoResume(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2))); 706 } 707 } 708 PetscFunctionReturn(PETSC_SUCCESS); 709 } 710 711 static PetscErrorCode PetscLogHandlerStagePush_Default(PetscLogHandler h, PetscLogStage new_stage) 712 { 713 PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 714 PetscLogDouble time; 715 PetscLogState state; 716 PetscLogStage current_stage; 717 PetscStagePerf *new_stage_info; 718 719 PetscFunctionBegin; 720 if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS); 721 PetscCall(PetscLogHandlerGetState(h, &state)); 722 current_stage = state->current_stage; 723 PetscCall(PetscLogHandlerDefaultGetStageInfo(h, new_stage, &new_stage_info)); 724 PetscCall(PetscTime(&time)); 725 726 /* Record flops/time of previous stage */ 727 if (current_stage >= 0) { 728 if (PetscBTLookup(state->active, current_stage)) { 729 PetscStagePerf *current_stage_info; 730 PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, ¤t_stage_info)); 731 PetscCall(PetscEventPerfInfoToc(¤t_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2))); 732 } 733 } 734 new_stage_info->used = PETSC_TRUE; 735 new_stage_info->perfInfo.count++; 736 new_stage_info->perfInfo.depth++; 737 /* Subtract current quantities so that we obtain the difference when we pop */ 738 if (PetscBTLookup(state->active, new_stage)) PetscCall(PetscEventPerfInfoTic(&new_stage_info->perfInfo, time, PetscLogMemory, (int)-(new_stage + 2))); 739 PetscFunctionReturn(PETSC_SUCCESS); 740 } 741 742 static PetscErrorCode PetscLogHandlerStagePop_Default(PetscLogHandler h, PetscLogStage old_stage) 743 { 744 PetscLogHandler_Default def = (PetscLogHandler_Default)h->data; 745 PetscLogStage current_stage; 746 PetscStagePerf *old_stage_info; 747 PetscLogState state; 748 PetscLogDouble time; 749 750 PetscFunctionBegin; 751 if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS); 752 PetscCall(PetscLogHandlerGetState(h, &state)); 753 current_stage = state->current_stage; 754 PetscCall(PetscLogHandlerDefaultGetStageInfo(h, old_stage, &old_stage_info)); 755 PetscCall(PetscTime(&time)); 756 old_stage_info->perfInfo.depth--; 757 if (PetscBTLookup(state->active, old_stage)) { PetscCall(PetscEventPerfInfoToc(&old_stage_info->perfInfo, time, PetscLogMemory, (int)-(old_stage + 2))); } 758 if (current_stage >= 0) { 759 if (PetscBTLookup(state->active, current_stage)) { 760 PetscStagePerf *current_stage_info; 761 PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, ¤t_stage_info)); 762 PetscCall(PetscEventPerfInfoTic(¤t_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2))); 763 } 764 } 765 PetscFunctionReturn(PETSC_SUCCESS); 766 } 767 768 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultStageSetVisible(PetscLogHandler h, PetscLogStage stage, PetscBool is_visible) 769 { 770 PetscStagePerf *stage_info; 771 772 PetscFunctionBegin; 773 if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage)); 774 PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info)); 775 stage_info->perfInfo.visible = is_visible; 776 PetscFunctionReturn(PETSC_SUCCESS); 777 } 778 779 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultStageGetVisible(PetscLogHandler h, PetscLogStage stage, PetscBool *is_visible) 780 { 781 PetscStagePerf *stage_info; 782 783 PetscFunctionBegin; 784 if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage)); 785 PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info)); 786 *is_visible = stage_info->perfInfo.visible; 787 PetscFunctionReturn(PETSC_SUCCESS); 788 } 789 790 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultSetLogActions(PetscLogHandler handler, PetscBool flag) 791 { 792 PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 793 794 PetscFunctionBegin; 795 def->petsc_logActions = flag; 796 PetscFunctionReturn(PETSC_SUCCESS); 797 } 798 799 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultSetLogObjects(PetscLogHandler handler, PetscBool flag) 800 { 801 PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 802 803 PetscFunctionBegin; 804 def->petsc_logObjects = flag; 805 PetscFunctionReturn(PETSC_SUCCESS); 806 } 807 808 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultLogObjectState(PetscLogHandler handler, PetscObject obj, const char format[], va_list Argp) 809 { 810 PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 811 size_t fullLength; 812 813 PetscFunctionBegin; 814 if (def->petsc_logObjects) { 815 Object *obj_entry = NULL; 816 817 PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, obj->id - 1, &obj_entry)); 818 PetscCall(PetscVSNPrintf(obj_entry->info, 64, format, &fullLength, Argp)); 819 } 820 PetscFunctionReturn(PETSC_SUCCESS); 821 } 822 823 PETSC_INTERN PetscErrorCode PetscLogHandlerDefaultGetNumObjects(PetscLogHandler handler, PetscInt *num_objects) 824 { 825 PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 826 827 PetscFunctionBegin; 828 PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, num_objects, NULL)); 829 PetscFunctionReturn(PETSC_SUCCESS); 830 } 831 832 PETSC_INTERN PetscErrorCode PetscLogHandlerDump_Default(PetscLogHandler handler, const char sname[]) 833 { 834 PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 835 FILE *fd; 836 char file[PETSC_MAX_PATH_LEN], fname[PETSC_MAX_PATH_LEN]; 837 PetscLogDouble flops, _TotalTime; 838 PetscMPIInt rank; 839 int curStage; 840 PetscLogState state; 841 PetscInt num_events; 842 PetscLogEvent event; 843 844 PetscFunctionBegin; 845 /* Calculate the total elapsed time */ 846 PetscCall(PetscTime(&_TotalTime)); 847 _TotalTime -= petsc_BaseTime; 848 /* Open log file */ 849 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)handler), &rank)); 850 PetscCall(PetscSNPrintf(file, PETSC_STATIC_ARRAY_LENGTH(file), "%s.%d", sname && sname[0] ? sname : "Log", rank)); 851 PetscCall(PetscFixFilename(file, fname)); 852 PetscCall(PetscFOpen(PETSC_COMM_SELF, fname, "w", &fd)); 853 PetscCheck(!(rank == 0) || !(!fd), PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file: %s", fname); 854 /* Output totals */ 855 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Total Flop %14e %16.8e\n", petsc_TotalFlops, _TotalTime)); 856 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Clock Resolution %g\n", 0.0)); 857 /* Output actions */ 858 if (def->petsc_logActions) { 859 PetscInt num_actions; 860 PetscCall(PetscLogActionArrayGetSize(def->petsc_actions, &num_actions, NULL)); 861 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Actions accomplished %d\n", (int)num_actions)); 862 for (int a = 0; a < num_actions; a++) { 863 Action *action; 864 865 PetscCall(PetscLogActionArrayGetRef(def->petsc_actions, a, &action)); 866 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)); 867 } 868 } 869 /* Output objects */ 870 if (def->petsc_logObjects) { 871 PetscInt num_objects; 872 873 PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL)); 874 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Objects created %d destroyed %d\n", def->petsc_numObjectsCreated, def->petsc_numObjectsDestroyed)); 875 for (int o = 0; o < num_objects; o++) { 876 Object *object = NULL; 877 878 PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, o, &object)); 879 if (object->parent != -1) continue; // object with this id wasn't logged, probably a PetscLogHandler 880 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Parent ID: %d Memory: %d\n", object->parent, (int)object->mem)); 881 if (!object->name[0]) { 882 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Name\n")); 883 } else { 884 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Name: %s\n", object->name)); 885 } 886 if (!object->info[0]) { 887 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Info\n")); 888 } else { 889 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Info: %s\n", object->info)); 890 } 891 } 892 } 893 /* Output events */ 894 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Event log:\n")); 895 PetscCall(PetscLogHandlerGetState(handler, &state)); 896 PetscCall(PetscLogStateGetNumEvents(state, &num_events)); 897 PetscCall(PetscLogStateGetCurrentStage(state, &curStage)); 898 for (event = 0; event < num_events; event++) { 899 PetscEventPerfInfo *event_info; 900 901 PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(handler, curStage, event, &event_info)); 902 if (event_info->time != 0.0) flops = event_info->flops / event_info->time; 903 else flops = 0.0; 904 PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "%d %16d %16g %16g %16g\n", event, event_info->count, event_info->flops, event_info->time, flops)); 905 } 906 PetscCall(PetscFClose(PETSC_COMM_SELF, fd)); 907 PetscFunctionReturn(PETSC_SUCCESS); 908 } 909 910 /* 911 PetscLogView_Detailed - Each process prints the times for its own events 912 913 */ 914 static PetscErrorCode PetscLogHandlerView_Default_Detailed(PetscLogHandler handler, PetscViewer viewer) 915 { 916 PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 917 PetscLogDouble locTotalTime, numRed, maxMem; 918 PetscInt numStages, numEvents; 919 MPI_Comm comm = PetscObjectComm((PetscObject)viewer); 920 PetscMPIInt rank, size; 921 PetscLogGlobalNames global_stages, global_events; 922 PetscLogState state; 923 PetscEventPerfInfo zero_info; 924 925 PetscFunctionBegin; 926 PetscCall(PetscLogHandlerGetState(handler, &state)); 927 PetscCallMPI(MPI_Comm_size(comm, &size)); 928 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 929 /* Must preserve reduction count before we go on */ 930 numRed = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct; 931 /* Get the total elapsed time */ 932 PetscCall(PetscTime(&locTotalTime)); 933 locTotalTime -= petsc_BaseTime; 934 PetscCall(PetscViewerASCIIPrintf(viewer, "size = %d\n", size)); 935 PetscCall(PetscViewerASCIIPrintf(viewer, "LocalTimes = {}\n")); 936 PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessages = {}\n")); 937 PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessageLens = {}\n")); 938 PetscCall(PetscViewerASCIIPrintf(viewer, "LocalReductions = {}\n")); 939 PetscCall(PetscViewerASCIIPrintf(viewer, "LocalFlop = {}\n")); 940 PetscCall(PetscViewerASCIIPrintf(viewer, "LocalObjects = {}\n")); 941 PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMemory = {}\n")); 942 PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages)); 943 PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events)); 944 PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages)); 945 PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents)); 946 PetscCall(PetscMemzero(&zero_info, sizeof(zero_info))); 947 PetscCall(PetscViewerASCIIPrintf(viewer, "Stages = {}\n")); 948 for (PetscInt stage = 0; stage < numStages; stage++) { 949 PetscInt stage_id; 950 const char *stage_name; 951 952 PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 953 PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 954 PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"] = {}\n", stage_name)); 955 PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"summary\"] = {}\n", stage_name)); 956 for (PetscInt event = 0; event < numEvents; event++) { 957 PetscEventPerfInfo *eventInfo = &zero_info; 958 PetscBool is_zero = PETSC_FALSE; 959 PetscInt event_id; 960 const char *event_name; 961 962 PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id)); 963 PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name)); 964 if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(handler, stage_id, event_id, &eventInfo)); 965 is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE; 966 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm)); 967 if (!is_zero) { PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"%s\"] = {}\n", stage_name, event_name)); } 968 } 969 } 970 PetscCall(PetscMallocGetMaximumUsage(&maxMem)); 971 PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 972 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalTimes[%d] = %g\n", rank, locTotalTime)); 973 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessages[%d] = %g\n", rank, (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct))); 974 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessageLens[%d] = %g\n", rank, (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len))); 975 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalReductions[%d] = %g\n", rank, numRed)); 976 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalFlop[%d] = %g\n", rank, petsc_TotalFlops)); 977 { 978 PetscInt num_objects; 979 980 PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL)); 981 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalObjects[%d] = %d\n", rank, (int)num_objects)); 982 } 983 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMemory[%d] = %g\n", rank, maxMem)); 984 PetscCall(PetscViewerFlush(viewer)); 985 for (PetscInt stage = 0; stage < numStages; stage++) { 986 PetscEventPerfInfo *stage_perf_info = &zero_info; 987 PetscInt stage_id; 988 const char *stage_name; 989 990 PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 991 PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 992 if (stage_id >= 0) { 993 PetscStagePerf *stage_info; 994 PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info)); 995 stage_perf_info = &stage_info->perfInfo; 996 } 997 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, 998 stage_perf_info->numMessages, stage_perf_info->messageLength, stage_perf_info->numReductions, stage_perf_info->flops)); 999 for (PetscInt event = 0; event < numEvents; event++) { 1000 PetscEventPerfInfo *eventInfo = &zero_info; 1001 PetscBool is_zero = PETSC_FALSE; 1002 PetscInt event_id; 1003 const char *event_name; 1004 1005 PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id)); 1006 PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name)); 1007 if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(handler, stage_id, event_id, &eventInfo)); 1008 is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE; 1009 PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero)); 1010 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm)); 1011 if (!is_zero) { 1012 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, 1013 eventInfo->count, eventInfo->time, eventInfo->syncTime, eventInfo->numMessages, eventInfo->messageLength, eventInfo->numReductions, eventInfo->flops)); 1014 if (eventInfo->dof[0] >= 0.) { 1015 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"dof\" : [")); 1016 for (PetscInt d = 0; d < 8; ++d) { 1017 if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", ")); 1018 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->dof[d])); 1019 } 1020 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]")); 1021 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"error\" : [")); 1022 for (PetscInt e = 0; e < 8; ++e) { 1023 if (e > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", ")); 1024 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->errors[e])); 1025 } 1026 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]")); 1027 } 1028 } 1029 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "}\n")); 1030 } 1031 } 1032 PetscCall(PetscViewerFlush(viewer)); 1033 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1034 PetscCall(PetscLogGlobalNamesDestroy(&global_events)); 1035 PetscCall(PetscLogGlobalNamesDestroy(&global_stages)); 1036 PetscFunctionReturn(PETSC_SUCCESS); 1037 } 1038 1039 /* 1040 PetscLogView_CSV - Each process prints the times for its own events in Comma-Separated Value Format 1041 */ 1042 static PetscErrorCode PetscLogHandlerView_Default_CSV(PetscLogHandler handler, PetscViewer viewer) 1043 { 1044 PetscLogDouble locTotalTime, maxMem; 1045 PetscInt numStages, numEvents, stage, event; 1046 MPI_Comm comm = PetscObjectComm((PetscObject)viewer); 1047 PetscMPIInt rank, size; 1048 PetscLogGlobalNames global_stages, global_events; 1049 PetscLogState state; 1050 PetscEventPerfInfo zero_info; 1051 1052 PetscFunctionBegin; 1053 PetscCall(PetscLogHandlerGetState(handler, &state)); 1054 PetscCallMPI(MPI_Comm_size(comm, &size)); 1055 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1056 /* Must preserve reduction count before we go on */ 1057 /* Get the total elapsed time */ 1058 PetscCall(PetscTime(&locTotalTime)); 1059 locTotalTime -= petsc_BaseTime; 1060 PetscCall(PetscMallocGetMaximumUsage(&maxMem)); 1061 PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages)); 1062 PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events)); 1063 PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages)); 1064 PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents)); 1065 PetscCall(PetscMemzero(&zero_info, sizeof(zero_info))); 1066 PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1067 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)); 1068 PetscCall(PetscViewerFlush(viewer)); 1069 for (stage = 0; stage < numStages; stage++) { 1070 PetscEventPerfInfo *stage_perf_info; 1071 PetscInt stage_id; 1072 const char *stage_name; 1073 1074 PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 1075 PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 1076 stage_perf_info = &zero_info; 1077 if (stage_id >= 0) { 1078 PetscStagePerf *stage_info; 1079 PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info)); 1080 stage_perf_info = &stage_info->perfInfo; 1081 } 1082 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)); 1083 for (event = 0; event < numEvents; event++) { 1084 PetscEventPerfInfo *eventInfo = &zero_info; 1085 PetscBool is_zero = PETSC_FALSE; 1086 PetscInt event_id; 1087 const char *event_name; 1088 1089 PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id)); 1090 PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name)); 1091 if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(handler, stage_id, event_id, &eventInfo)); 1092 PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero)); 1093 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm)); 1094 if (!is_zero) { 1095 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)); 1096 if (eventInfo->dof[0] >= 0.) { 1097 for (PetscInt d = 0; d < 8; ++d) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->dof[d])); 1098 for (PetscInt e = 0; e < 8; ++e) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->errors[e])); 1099 } 1100 } 1101 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 1102 } 1103 } 1104 PetscCall(PetscViewerFlush(viewer)); 1105 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1106 PetscCall(PetscLogGlobalNamesDestroy(&global_stages)); 1107 PetscCall(PetscLogGlobalNamesDestroy(&global_events)); 1108 PetscFunctionReturn(PETSC_SUCCESS); 1109 } 1110 1111 static PetscErrorCode PetscLogViewWarnSync(MPI_Comm comm, FILE *fd) 1112 { 1113 PetscFunctionBegin; 1114 if (!PetscLogSyncOn) PetscFunctionReturn(PETSC_SUCCESS); 1115 PetscCall(PetscFPrintf(comm, fd, "\n\n")); 1116 PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n")); 1117 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1118 PetscCall(PetscFPrintf(comm, fd, " # WARNING!!! #\n")); 1119 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1120 PetscCall(PetscFPrintf(comm, fd, " # This program was run with logging synchronization. #\n")); 1121 PetscCall(PetscFPrintf(comm, fd, " # This option provides more meaningful imbalance #\n")); 1122 PetscCall(PetscFPrintf(comm, fd, " # figures at the expense of slowing things down and #\n")); 1123 PetscCall(PetscFPrintf(comm, fd, " # providing a distorted view of the overall runtime. #\n")); 1124 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1125 PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n\n\n")); 1126 PetscFunctionReturn(PETSC_SUCCESS); 1127 } 1128 1129 static PetscErrorCode PetscLogViewWarnDebugging(MPI_Comm comm, FILE *fd) 1130 { 1131 PetscFunctionBegin; 1132 if (PetscDefined(USE_DEBUG)) { 1133 PetscCall(PetscFPrintf(comm, fd, "\n\n")); 1134 PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n")); 1135 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1136 PetscCall(PetscFPrintf(comm, fd, " # WARNING!!! #\n")); 1137 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1138 PetscCall(PetscFPrintf(comm, fd, " # This code was compiled with a debugging option. #\n")); 1139 PetscCall(PetscFPrintf(comm, fd, " # To get timing results run ./configure #\n")); 1140 PetscCall(PetscFPrintf(comm, fd, " # using --with-debugging=no, the performance will #\n")); 1141 PetscCall(PetscFPrintf(comm, fd, " # be generally two or three times faster. #\n")); 1142 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1143 PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n\n\n")); 1144 } 1145 PetscFunctionReturn(PETSC_SUCCESS); 1146 } 1147 1148 static PetscErrorCode PetscLogViewWarnNoGpuAwareMpi(MPI_Comm comm, FILE *fd) 1149 { 1150 #if defined(PETSC_HAVE_DEVICE) 1151 PetscMPIInt size; 1152 PetscBool deviceInitialized = PETSC_FALSE; 1153 1154 PetscFunctionBegin; 1155 PetscCallMPI(MPI_Comm_size(comm, &size)); 1156 for (int i = PETSC_DEVICE_HOST + 1; i < PETSC_DEVICE_MAX; ++i) { 1157 const PetscDeviceType dtype = PetscDeviceTypeCast(i); 1158 if (PetscDeviceInitialized(dtype)) { /* a non-host device was initialized */ 1159 deviceInitialized = PETSC_TRUE; 1160 break; 1161 } 1162 } 1163 /* the last condition says petsc is configured with device but it is a pure CPU run, so don't print misleading warnings */ 1164 if (use_gpu_aware_mpi || size == 1 || !deviceInitialized) PetscFunctionReturn(PETSC_SUCCESS); 1165 PetscCall(PetscFPrintf(comm, fd, "\n\n")); 1166 PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n")); 1167 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1168 PetscCall(PetscFPrintf(comm, fd, " # WARNING!!! #\n")); 1169 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1170 PetscCall(PetscFPrintf(comm, fd, " # This code was compiled with GPU support and you've #\n")); 1171 PetscCall(PetscFPrintf(comm, fd, " # created PETSc/GPU objects, but you intentionally #\n")); 1172 PetscCall(PetscFPrintf(comm, fd, " # used -use_gpu_aware_mpi 0, requiring PETSc to copy #\n")); 1173 PetscCall(PetscFPrintf(comm, fd, " # additional data between the GPU and CPU. To obtain #\n")); 1174 PetscCall(PetscFPrintf(comm, fd, " # meaningful timing results on multi-rank runs, use #\n")); 1175 PetscCall(PetscFPrintf(comm, fd, " # GPU-aware MPI instead. #\n")); 1176 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1177 PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n\n\n")); 1178 PetscFunctionReturn(PETSC_SUCCESS); 1179 #else 1180 return PETSC_SUCCESS; 1181 #endif 1182 } 1183 1184 static PetscErrorCode PetscLogViewWarnGpuTime(MPI_Comm comm, FILE *fd) 1185 { 1186 #if defined(PETSC_HAVE_DEVICE) 1187 1188 PetscFunctionBegin; 1189 if (!PetscLogGpuTimeFlag || petsc_gflops == 0) PetscFunctionReturn(PETSC_SUCCESS); 1190 PetscCall(PetscFPrintf(comm, fd, "\n\n")); 1191 PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n")); 1192 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1193 PetscCall(PetscFPrintf(comm, fd, " # WARNING!!! #\n")); 1194 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1195 PetscCall(PetscFPrintf(comm, fd, " # This code was run with -log_view_gpu_time #\n")); 1196 PetscCall(PetscFPrintf(comm, fd, " # This provides accurate timing within the GPU kernels #\n")); 1197 PetscCall(PetscFPrintf(comm, fd, " # but can slow down the entire computation by a #\n")); 1198 PetscCall(PetscFPrintf(comm, fd, " # measurable amount. For fastest runs we recommend #\n")); 1199 PetscCall(PetscFPrintf(comm, fd, " # not using this option. #\n")); 1200 PetscCall(PetscFPrintf(comm, fd, " # #\n")); 1201 PetscCall(PetscFPrintf(comm, fd, " ##########################################################\n\n\n")); 1202 PetscFunctionReturn(PETSC_SUCCESS); 1203 #else 1204 return PETSC_SUCCESS; 1205 #endif 1206 } 1207 1208 static PetscErrorCode PetscLogHandlerView_Default_Info(PetscLogHandler handler, PetscViewer viewer) 1209 { 1210 FILE *fd; 1211 PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data; 1212 char arch[128], hostname[128], username[128], pname[PETSC_MAX_PATH_LEN], date[128]; 1213 PetscLogDouble locTotalTime, TotalTime, TotalFlops; 1214 PetscLogDouble numMessages, messageLength, avgMessLen, numReductions; 1215 PetscLogDouble stageTime, flops, flopr, mem, mess, messLen, red; 1216 PetscLogDouble fracTime, fracFlops, fracMessages, fracLength, fracReductions, fracMess, fracMessLen, fracRed; 1217 PetscLogDouble fracStageTime, fracStageFlops, fracStageMess, fracStageMessLen, fracStageRed; 1218 PetscLogDouble min, max, tot, ratio, avg, x, y; 1219 PetscLogDouble minf, maxf, totf, ratf, mint, maxt, tott, ratt, ratC, totm, totml, totr, mal, malmax, emalmax; 1220 #if defined(PETSC_HAVE_DEVICE) 1221 PetscLogEvent KSP_Solve, SNES_Solve, TS_Step, TAO_Solve; /* These need to be fixed to be some events registered with certain objects */ 1222 PetscLogDouble cct, gct, csz, gsz, gmaxt, gflops, gflopr, fracgflops; 1223 #endif 1224 PetscMPIInt minC, maxC; 1225 PetscMPIInt size, rank; 1226 PetscBool *localStageUsed, *stageUsed; 1227 PetscBool *localStageVisible, *stageVisible; 1228 PetscInt numStages, numEvents; 1229 int stage, oclass; 1230 PetscLogEvent event; 1231 char version[256]; 1232 MPI_Comm comm; 1233 #if defined(PETSC_HAVE_DEVICE) 1234 PetscInt64 nas = 0x7FF0000000000002; 1235 #endif 1236 PetscLogGlobalNames global_stages, global_events; 1237 PetscEventPerfInfo zero_info; 1238 PetscLogState state; 1239 1240 PetscFunctionBegin; 1241 PetscCall(PetscLogHandlerGetState(handler, &state)); 1242 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 1243 PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm)); 1244 PetscCall(PetscViewerASCIIGetPointer(viewer, &fd)); 1245 PetscCallMPI(MPI_Comm_size(comm, &size)); 1246 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1247 /* Get the total elapsed time */ 1248 PetscCall(PetscTime(&locTotalTime)); 1249 locTotalTime -= petsc_BaseTime; 1250 1251 PetscCall(PetscFPrintf(comm, fd, "****************************************************************************************************************************************************************\n")); 1252 PetscCall(PetscFPrintf(comm, fd, "*** WIDEN YOUR WINDOW TO 160 CHARACTERS. Use 'enscript -r -fCourier9' to print this document ***\n")); 1253 PetscCall(PetscFPrintf(comm, fd, "****************************************************************************************************************************************************************\n")); 1254 PetscCall(PetscFPrintf(comm, fd, "\n------------------------------------------------------------------ PETSc Performance Summary: ------------------------------------------------------------------\n\n")); 1255 PetscCall(PetscLogViewWarnSync(comm, fd)); 1256 PetscCall(PetscLogViewWarnDebugging(comm, fd)); 1257 PetscCall(PetscLogViewWarnNoGpuAwareMpi(comm, fd)); 1258 PetscCall(PetscLogViewWarnGpuTime(comm, fd)); 1259 PetscCall(PetscGetArchType(arch, sizeof(arch))); 1260 PetscCall(PetscGetHostName(hostname, sizeof(hostname))); 1261 PetscCall(PetscGetUserName(username, sizeof(username))); 1262 PetscCall(PetscGetProgramName(pname, sizeof(pname))); 1263 PetscCall(PetscGetDate(date, sizeof(date))); 1264 PetscCall(PetscGetVersion(version, sizeof(version))); 1265 if (size == 1) { 1266 PetscCall(PetscFPrintf(comm, fd, "%s on a %s named %s with %d processor, by %s %s\n", pname, arch, hostname, size, username, date)); 1267 } else { 1268 PetscCall(PetscFPrintf(comm, fd, "%s on a %s named %s with %d processors, by %s %s\n", pname, arch, hostname, size, username, date)); 1269 } 1270 #if defined(PETSC_HAVE_OPENMP) 1271 PetscCall(PetscFPrintf(comm, fd, "Using %" PetscInt_FMT " OpenMP threads\n", PetscNumOMPThreads)); 1272 #endif 1273 PetscCall(PetscFPrintf(comm, fd, "Using %s\n", version)); 1274 1275 /* Must preserve reduction count before we go on */ 1276 red = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct; 1277 1278 /* Calculate summary information */ 1279 PetscCall(PetscFPrintf(comm, fd, "\n Max Max/Min Avg Total\n")); 1280 /* Time */ 1281 PetscCall(MPIU_Allreduce(&locTotalTime, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 1282 PetscCall(MPIU_Allreduce(&locTotalTime, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1283 PetscCall(MPIU_Allreduce(&locTotalTime, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1284 avg = tot / ((PetscLogDouble)size); 1285 if (min != 0.0) ratio = max / min; 1286 else ratio = 0.0; 1287 PetscCall(PetscFPrintf(comm, fd, "Time (sec): %5.3e %7.3f %5.3e\n", max, ratio, avg)); 1288 TotalTime = tot; 1289 /* Objects */ 1290 { 1291 PetscInt num_objects; 1292 1293 PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL)); 1294 avg = (PetscLogDouble)num_objects; 1295 } 1296 PetscCall(MPIU_Allreduce(&avg, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 1297 PetscCall(MPIU_Allreduce(&avg, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1298 PetscCall(MPIU_Allreduce(&avg, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1299 avg = tot / ((PetscLogDouble)size); 1300 if (min != 0.0) ratio = max / min; 1301 else ratio = 0.0; 1302 PetscCall(PetscFPrintf(comm, fd, "Objects: %5.3e %7.3f %5.3e\n", max, ratio, avg)); 1303 /* Flops */ 1304 PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 1305 PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1306 PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1307 avg = tot / ((PetscLogDouble)size); 1308 if (min != 0.0) ratio = max / min; 1309 else ratio = 0.0; 1310 PetscCall(PetscFPrintf(comm, fd, "Flops: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot)); 1311 TotalFlops = tot; 1312 /* Flops/sec -- Must talk to Barry here */ 1313 if (locTotalTime != 0.0) flops = petsc_TotalFlops / locTotalTime; 1314 else flops = 0.0; 1315 PetscCall(MPIU_Allreduce(&flops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 1316 PetscCall(MPIU_Allreduce(&flops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1317 PetscCall(MPIU_Allreduce(&flops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1318 avg = tot / ((PetscLogDouble)size); 1319 if (min != 0.0) ratio = max / min; 1320 else ratio = 0.0; 1321 PetscCall(PetscFPrintf(comm, fd, "Flops/sec: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot)); 1322 /* Memory */ 1323 PetscCall(PetscMallocGetMaximumUsage(&mem)); 1324 if (mem > 0.0) { 1325 PetscCall(MPIU_Allreduce(&mem, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 1326 PetscCall(MPIU_Allreduce(&mem, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1327 PetscCall(MPIU_Allreduce(&mem, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1328 avg = tot / ((PetscLogDouble)size); 1329 if (min != 0.0) ratio = max / min; 1330 else ratio = 0.0; 1331 PetscCall(PetscFPrintf(comm, fd, "Memory (bytes): %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot)); 1332 } 1333 /* Messages */ 1334 mess = 0.5 * (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct); 1335 PetscCall(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 1336 PetscCall(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1337 PetscCall(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1338 avg = tot / ((PetscLogDouble)size); 1339 if (min != 0.0) ratio = max / min; 1340 else ratio = 0.0; 1341 PetscCall(PetscFPrintf(comm, fd, "MPI Msg Count: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot)); 1342 numMessages = tot; 1343 /* Message Lengths */ 1344 mess = 0.5 * (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len); 1345 PetscCall(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 1346 PetscCall(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1347 PetscCall(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1348 if (numMessages != 0) avg = tot / numMessages; 1349 else avg = 0.0; 1350 if (min != 0.0) ratio = max / min; 1351 else ratio = 0.0; 1352 PetscCall(PetscFPrintf(comm, fd, "MPI Msg Len (bytes): %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot)); 1353 messageLength = tot; 1354 /* Reductions */ 1355 PetscCall(MPIU_Allreduce(&red, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 1356 PetscCall(MPIU_Allreduce(&red, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1357 PetscCall(MPIU_Allreduce(&red, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1358 if (min != 0.0) ratio = max / min; 1359 else ratio = 0.0; 1360 PetscCall(PetscFPrintf(comm, fd, "MPI Reductions: %5.3e %7.3f\n", max, ratio)); 1361 numReductions = red; /* wrong because uses count from process zero */ 1362 PetscCall(PetscFPrintf(comm, fd, "\nFlop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)\n")); 1363 PetscCall(PetscFPrintf(comm, fd, " e.g., VecAXPY() for real vectors of length N --> 2N flops\n")); 1364 PetscCall(PetscFPrintf(comm, fd, " and VecAXPY() for complex vectors of length N --> 8N flops\n")); 1365 1366 PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages)); 1367 PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events)); 1368 PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages)); 1369 PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents)); 1370 PetscCall(PetscMemzero(&zero_info, sizeof(zero_info))); 1371 PetscCall(PetscMalloc1(numStages, &localStageUsed)); 1372 PetscCall(PetscMalloc1(numStages, &stageUsed)); 1373 PetscCall(PetscMalloc1(numStages, &localStageVisible)); 1374 PetscCall(PetscMalloc1(numStages, &stageVisible)); 1375 if (numStages > 0) { 1376 for (stage = 0; stage < numStages; stage++) { 1377 PetscInt stage_id; 1378 1379 PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 1380 if (stage_id >= 0) { 1381 PetscStagePerf *stage_info; 1382 1383 PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info)); 1384 localStageUsed[stage] = stage_info->used; 1385 localStageVisible[stage] = stage_info->perfInfo.visible; 1386 } else { 1387 localStageUsed[stage] = PETSC_FALSE; 1388 localStageVisible[stage] = PETSC_TRUE; 1389 } 1390 } 1391 PetscCall(MPIU_Allreduce(localStageUsed, stageUsed, numStages, MPIU_BOOL, MPI_LOR, comm)); 1392 PetscCall(MPIU_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm)); 1393 for (stage = 0; stage < numStages; stage++) { 1394 if (stageUsed[stage] && stageVisible[stage]) { 1395 PetscCall(PetscFPrintf(comm, fd, "\nSummary of Stages: ----- Time ------ ----- Flop ------ --- Messages --- -- Message Lengths -- -- Reductions --\n")); 1396 PetscCall(PetscFPrintf(comm, fd, " Avg %%Total Avg %%Total Count %%Total Avg %%Total Count %%Total\n")); 1397 break; 1398 } 1399 } 1400 for (stage = 0; stage < numStages; stage++) { 1401 PetscInt stage_id; 1402 PetscEventPerfInfo *stage_info; 1403 const char *stage_name; 1404 1405 if (!(stageUsed[stage] && stageVisible[stage])) continue; 1406 PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 1407 PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 1408 stage_info = &zero_info; 1409 if (localStageUsed[stage]) { 1410 PetscStagePerf *stage_perf_info; 1411 1412 PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_perf_info)); 1413 stage_info = &stage_perf_info->perfInfo; 1414 } 1415 PetscCall(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1416 PetscCall(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1417 PetscCall(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1418 PetscCall(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1419 PetscCall(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1420 mess *= 0.5; 1421 messLen *= 0.5; 1422 red /= size; 1423 if (TotalTime != 0.0) fracTime = stageTime / TotalTime; 1424 else fracTime = 0.0; 1425 if (TotalFlops != 0.0) fracFlops = flops / TotalFlops; 1426 else fracFlops = 0.0; 1427 /* Talk to Barry if (stageTime != 0.0) flops = (size*flops)/stageTime; else flops = 0.0; */ 1428 if (numMessages != 0.0) fracMessages = mess / numMessages; 1429 else fracMessages = 0.0; 1430 if (mess != 0.0) avgMessLen = messLen / mess; 1431 else avgMessLen = 0.0; 1432 if (messageLength != 0.0) fracLength = messLen / messageLength; 1433 else fracLength = 0.0; 1434 if (numReductions != 0.0) fracReductions = red / numReductions; 1435 else fracReductions = 0.0; 1436 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)); 1437 } 1438 } 1439 1440 PetscCall(PetscFPrintf(comm, fd, "\n------------------------------------------------------------------------------------------------------------------------\n")); 1441 PetscCall(PetscFPrintf(comm, fd, "See the 'Profiling' chapter of the users' manual for details on interpreting output.\n")); 1442 PetscCall(PetscFPrintf(comm, fd, "Phase summary info:\n")); 1443 PetscCall(PetscFPrintf(comm, fd, " Count: number of times phase was executed\n")); 1444 PetscCall(PetscFPrintf(comm, fd, " Time and Flop: Max - maximum over all processors\n")); 1445 PetscCall(PetscFPrintf(comm, fd, " Ratio - ratio of maximum to minimum over all processors\n")); 1446 PetscCall(PetscFPrintf(comm, fd, " Mess: number of messages sent\n")); 1447 PetscCall(PetscFPrintf(comm, fd, " AvgLen: average message length (bytes)\n")); 1448 PetscCall(PetscFPrintf(comm, fd, " Reduct: number of global reductions\n")); 1449 PetscCall(PetscFPrintf(comm, fd, " Global: entire computation\n")); 1450 PetscCall(PetscFPrintf(comm, fd, " Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().\n")); 1451 PetscCall(PetscFPrintf(comm, fd, " %%T - percent time in this phase %%F - percent flop in this phase\n")); 1452 PetscCall(PetscFPrintf(comm, fd, " %%M - percent messages in this phase %%L - percent message lengths in this phase\n")); 1453 PetscCall(PetscFPrintf(comm, fd, " %%R - percent reductions in this phase\n")); 1454 PetscCall(PetscFPrintf(comm, fd, " Total Mflop/s: 10e-6 * (sum of flop over all processors)/(max time over all processors)\n")); 1455 if (PetscLogMemory) { 1456 PetscCall(PetscFPrintf(comm, fd, " Memory usage is summed over all MPI processes, it is given in mega-bytes\n")); 1457 PetscCall(PetscFPrintf(comm, fd, " Malloc Mbytes: Memory allocated and kept during event (sum over all calls to event). May be negative\n")); 1458 PetscCall(PetscFPrintf(comm, fd, " EMalloc Mbytes: extra memory allocated during event and then freed (maximum over all calls to events). Never negative\n")); 1459 PetscCall(PetscFPrintf(comm, fd, " MMalloc Mbytes: Increase in high water mark of allocated memory (sum over all calls to event). Never negative\n")); 1460 PetscCall(PetscFPrintf(comm, fd, " RMI Mbytes: Increase in resident memory (sum over all calls to event)\n")); 1461 } 1462 #if defined(PETSC_HAVE_DEVICE) 1463 PetscCall(PetscFPrintf(comm, fd, " GPU Mflop/s: 10e-6 * (sum of flop on GPU over all processors)/(max GPU time over all processors)\n")); 1464 PetscCall(PetscFPrintf(comm, fd, " CpuToGpu Count: total number of CPU to GPU copies per processor\n")); 1465 PetscCall(PetscFPrintf(comm, fd, " CpuToGpu Size (Mbytes): 10e-6 * (total size of CPU to GPU copies per processor)\n")); 1466 PetscCall(PetscFPrintf(comm, fd, " GpuToCpu Count: total number of GPU to CPU copies per processor\n")); 1467 PetscCall(PetscFPrintf(comm, fd, " GpuToCpu Size (Mbytes): 10e-6 * (total size of GPU to CPU copies per processor)\n")); 1468 PetscCall(PetscFPrintf(comm, fd, " GPU %%F: percent flops on GPU in this event\n")); 1469 #endif 1470 PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------\n")); 1471 1472 PetscCall(PetscLogViewWarnDebugging(comm, fd)); 1473 1474 /* Report events */ 1475 PetscCall(PetscFPrintf(comm, fd, "Event Count Time (sec) Flop --- Global --- --- Stage ---- Total")); 1476 if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, " Malloc EMalloc MMalloc RMI")); 1477 #if defined(PETSC_HAVE_DEVICE) 1478 PetscCall(PetscFPrintf(comm, fd, " GPU - CpuToGpu - - GpuToCpu - GPU")); 1479 #endif 1480 PetscCall(PetscFPrintf(comm, fd, "\n")); 1481 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")); 1482 if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, " Mbytes Mbytes Mbytes Mbytes")); 1483 #if defined(PETSC_HAVE_DEVICE) 1484 PetscCall(PetscFPrintf(comm, fd, " Mflop/s Count Size Count Size %%F")); 1485 #endif 1486 PetscCall(PetscFPrintf(comm, fd, "\n")); 1487 PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------")); 1488 if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, "-----------------------------")); 1489 #if defined(PETSC_HAVE_DEVICE) 1490 PetscCall(PetscFPrintf(comm, fd, "---------------------------------------")); 1491 #endif 1492 PetscCall(PetscFPrintf(comm, fd, "\n")); 1493 1494 #if defined(PETSC_HAVE_DEVICE) 1495 /* this indirect way of accessing these values is needed when PETSc is build with multiple libraries since the symbols are not in libpetscsys */ 1496 PetscCall(PetscLogStateGetEventFromName(state, "TAOSolve", &TAO_Solve)); 1497 PetscCall(PetscLogStateGetEventFromName(state, "TSStep", &TS_Step)); 1498 PetscCall(PetscLogStateGetEventFromName(state, "SNESSolve", &SNES_Solve)); 1499 PetscCall(PetscLogStateGetEventFromName(state, "KSPSolve", &KSP_Solve)); 1500 #endif 1501 1502 for (stage = 0; stage < numStages; stage++) { 1503 PetscInt stage_id; 1504 PetscEventPerfInfo *stage_info; 1505 const char *stage_name; 1506 1507 if (!(stageVisible[stage] && stageUsed[stage])) continue; 1508 PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id)); 1509 PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 1510 PetscCall(PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stage_name)); 1511 stage_info = &zero_info; 1512 if (localStageUsed[stage]) { 1513 PetscStagePerf *stage_perf_info; 1514 1515 PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_perf_info)); 1516 stage_info = &stage_perf_info->perfInfo; 1517 } 1518 PetscCall(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1519 PetscCall(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1520 PetscCall(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1521 PetscCall(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1522 PetscCall(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1523 mess *= 0.5; 1524 messLen *= 0.5; 1525 red /= size; 1526 1527 for (event = 0; event < numEvents; event++) { 1528 PetscInt event_id; 1529 PetscEventPerfInfo *event_info = &zero_info; 1530 PetscBool is_zero = PETSC_FALSE; 1531 const char *event_name; 1532 1533 PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id)); 1534 PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name)); 1535 if (event_id >= 0 && stage_id >= 0) { PetscCall(PetscLogHandlerDefaultGetEventPerfInfo(handler, stage_id, event_id, &event_info)); } 1536 PetscCall(PetscMemcmp(event_info, &zero_info, sizeof(zero_info), &is_zero)); 1537 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm)); 1538 if (!is_zero) { 1539 flopr = event_info->flops; 1540 PetscCall(MPIU_Allreduce(&flopr, &minf, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 1541 PetscCall(MPIU_Allreduce(&flopr, &maxf, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1542 PetscCall(MPIU_Allreduce(&event_info->flops, &totf, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1543 PetscCall(MPIU_Allreduce(&event_info->time, &mint, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm)); 1544 PetscCall(MPIU_Allreduce(&event_info->time, &maxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1545 PetscCall(MPIU_Allreduce(&event_info->time, &tott, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1546 PetscCall(MPIU_Allreduce(&event_info->numMessages, &totm, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1547 PetscCall(MPIU_Allreduce(&event_info->messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1548 PetscCall(MPIU_Allreduce(&event_info->numReductions, &totr, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1549 PetscCall(MPIU_Allreduce(&event_info->count, &minC, 1, MPI_INT, MPI_MIN, comm)); 1550 PetscCall(MPIU_Allreduce(&event_info->count, &maxC, 1, MPI_INT, MPI_MAX, comm)); 1551 if (PetscLogMemory) { 1552 PetscCall(MPIU_Allreduce(&event_info->memIncrease, &mem, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1553 PetscCall(MPIU_Allreduce(&event_info->mallocSpace, &mal, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1554 PetscCall(MPIU_Allreduce(&event_info->mallocIncrease, &malmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1555 PetscCall(MPIU_Allreduce(&event_info->mallocIncreaseEvent, &emalmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1556 } 1557 #if defined(PETSC_HAVE_DEVICE) 1558 PetscCall(MPIU_Allreduce(&event_info->CpuToGpuCount, &cct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1559 PetscCall(MPIU_Allreduce(&event_info->GpuToCpuCount, &gct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1560 PetscCall(MPIU_Allreduce(&event_info->CpuToGpuSize, &csz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1561 PetscCall(MPIU_Allreduce(&event_info->GpuToCpuSize, &gsz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1562 PetscCall(MPIU_Allreduce(&event_info->GpuFlops, &gflops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm)); 1563 PetscCall(MPIU_Allreduce(&event_info->GpuTime, &gmaxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm)); 1564 #endif 1565 if (mint < 0.0) { 1566 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)); 1567 mint = 0; 1568 } 1569 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); 1570 #if defined(PETSC_HAVE_DEVICE) 1571 /* Put NaN into the time for all events that may not be time accurately since they may happen asynchronously on the GPU */ 1572 if (!PetscLogGpuTimeFlag && petsc_gflops > 0) { 1573 memcpy(&gmaxt, &nas, sizeof(PetscLogDouble)); 1574 if (event_id != SNES_Solve && event_id != KSP_Solve && event_id != TS_Step && event_id != TAO_Solve) { 1575 memcpy(&mint, &nas, sizeof(PetscLogDouble)); 1576 memcpy(&maxt, &nas, sizeof(PetscLogDouble)); 1577 } 1578 } 1579 #endif 1580 totm *= 0.5; 1581 totml *= 0.5; 1582 totr /= size; 1583 1584 if (maxC != 0) { 1585 if (minC != 0) ratC = ((PetscLogDouble)maxC) / minC; 1586 else ratC = 0.0; 1587 if (mint != 0.0) ratt = maxt / mint; 1588 else ratt = 0.0; 1589 if (minf != 0.0) ratf = maxf / minf; 1590 else ratf = 0.0; 1591 if (TotalTime != 0.0) fracTime = tott / TotalTime; 1592 else fracTime = 0.0; 1593 if (TotalFlops != 0.0) fracFlops = totf / TotalFlops; 1594 else fracFlops = 0.0; 1595 if (stageTime != 0.0) fracStageTime = tott / stageTime; 1596 else fracStageTime = 0.0; 1597 if (flops != 0.0) fracStageFlops = totf / flops; 1598 else fracStageFlops = 0.0; 1599 if (numMessages != 0.0) fracMess = totm / numMessages; 1600 else fracMess = 0.0; 1601 if (messageLength != 0.0) fracMessLen = totml / messageLength; 1602 else fracMessLen = 0.0; 1603 if (numReductions != 0.0) fracRed = totr / numReductions; 1604 else fracRed = 0.0; 1605 if (mess != 0.0) fracStageMess = totm / mess; 1606 else fracStageMess = 0.0; 1607 if (messLen != 0.0) fracStageMessLen = totml / messLen; 1608 else fracStageMessLen = 0.0; 1609 if (red != 0.0) fracStageRed = totr / red; 1610 else fracStageRed = 0.0; 1611 if (totm != 0.0) totml /= totm; 1612 else totml = 0.0; 1613 if (maxt != 0.0) flopr = totf / maxt; 1614 else flopr = 0.0; 1615 if (fracStageTime > 1.0 || fracStageFlops > 1.0 || fracStageMess > 1.0 || fracStageMessLen > 1.0 || fracStageRed > 1.0) 1616 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)); 1617 else 1618 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)); 1619 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)); 1620 #if defined(PETSC_HAVE_DEVICE) 1621 if (totf != 0.0) fracgflops = gflops / totf; 1622 else fracgflops = 0.0; 1623 if (gmaxt != 0.0) gflopr = gflops / gmaxt; 1624 else gflopr = 0.0; 1625 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)); 1626 #endif 1627 PetscCall(PetscFPrintf(comm, fd, "\n")); 1628 } 1629 } 1630 } 1631 } 1632 1633 /* Memory usage and object creation */ 1634 PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------")); 1635 if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, "-----------------------------")); 1636 #if defined(PETSC_HAVE_DEVICE) 1637 PetscCall(PetscFPrintf(comm, fd, "---------------------------------------")); 1638 #endif 1639 PetscCall(PetscFPrintf(comm, fd, "\n")); 1640 PetscCall(PetscFPrintf(comm, fd, "\n")); 1641 1642 /* Right now, only stages on the first processor are reported here, meaning only objects associated with 1643 the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then 1644 stats for stages local to processor sets. 1645 */ 1646 /* We should figure out the longest object name here (now 20 characters) */ 1647 PetscCall(PetscFPrintf(comm, fd, "Object Type Creations Destructions. Reports information only for process 0.\n")); 1648 for (stage = 0; stage < numStages; stage++) { 1649 const char *stage_name; 1650 1651 PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name)); 1652 PetscCall(PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stage_name)); 1653 if (localStageUsed[stage]) { 1654 PetscInt num_classes; 1655 1656 PetscCall(PetscLogStateGetNumClasses(state, &num_classes)); 1657 for (oclass = 0; oclass < num_classes; oclass++) { 1658 PetscClassPerf *class_perf_info; 1659 1660 PetscCall(PetscLogHandlerDefaultGetClassPerf(handler, stage, oclass, &class_perf_info)); 1661 if ((class_perf_info->creations > 0) || (class_perf_info->destructions > 0)) { 1662 PetscLogClassInfo class_reg_info; 1663 1664 PetscCall(PetscLogStateClassGetInfo(state, oclass, &class_reg_info)); 1665 PetscCall(PetscFPrintf(comm, fd, "%20s %5d %5d\n", class_reg_info.name, class_perf_info->creations, class_perf_info->destructions)); 1666 } 1667 } 1668 } 1669 } 1670 1671 PetscCall(PetscFree(localStageUsed)); 1672 PetscCall(PetscFree(stageUsed)); 1673 PetscCall(PetscFree(localStageVisible)); 1674 PetscCall(PetscFree(stageVisible)); 1675 PetscCall(PetscLogGlobalNamesDestroy(&global_stages)); 1676 PetscCall(PetscLogGlobalNamesDestroy(&global_events)); 1677 1678 /* Information unrelated to this particular run */ 1679 PetscCall(PetscFPrintf(comm, fd, "========================================================================================================================\n")); 1680 PetscCall(PetscTime(&y)); 1681 PetscCall(PetscTime(&x)); 1682 PetscCall(PetscTime(&y)); 1683 PetscCall(PetscTime(&y)); 1684 PetscCall(PetscTime(&y)); 1685 PetscCall(PetscTime(&y)); 1686 PetscCall(PetscTime(&y)); 1687 PetscCall(PetscTime(&y)); 1688 PetscCall(PetscTime(&y)); 1689 PetscCall(PetscTime(&y)); 1690 PetscCall(PetscTime(&y)); 1691 PetscCall(PetscTime(&y)); 1692 PetscCall(PetscFPrintf(comm, fd, "Average time to get PetscTime(): %g\n", (y - x) / 10.0)); 1693 /* MPI information */ 1694 if (size > 1) { 1695 MPI_Status status; 1696 PetscMPIInt tag; 1697 MPI_Comm newcomm; 1698 1699 PetscCallMPI(MPI_Barrier(comm)); 1700 PetscCall(PetscTime(&x)); 1701 PetscCallMPI(MPI_Barrier(comm)); 1702 PetscCallMPI(MPI_Barrier(comm)); 1703 PetscCallMPI(MPI_Barrier(comm)); 1704 PetscCallMPI(MPI_Barrier(comm)); 1705 PetscCallMPI(MPI_Barrier(comm)); 1706 PetscCall(PetscTime(&y)); 1707 PetscCall(PetscFPrintf(comm, fd, "Average time for MPI_Barrier(): %g\n", (y - x) / 5.0)); 1708 PetscCall(PetscCommDuplicate(comm, &newcomm, &tag)); 1709 PetscCallMPI(MPI_Barrier(comm)); 1710 if (rank) { 1711 PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, rank - 1, tag, newcomm, &status)); 1712 PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, (rank + 1) % size, tag, newcomm)); 1713 } else { 1714 PetscCall(PetscTime(&x)); 1715 PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, 1, tag, newcomm)); 1716 PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, size - 1, tag, newcomm, &status)); 1717 PetscCall(PetscTime(&y)); 1718 PetscCall(PetscFPrintf(comm, fd, "Average time for zero size MPI_Send(): %g\n", (y - x) / size)); 1719 } 1720 PetscCall(PetscCommDestroy(&newcomm)); 1721 } 1722 PetscCall(PetscOptionsView(NULL, viewer)); 1723 1724 /* Machine and compile information */ 1725 if (PetscDefined(USE_FORTRAN_KERNELS)) { 1726 PetscCall(PetscFPrintf(comm, fd, "Compiled with FORTRAN kernels\n")); 1727 } else { 1728 PetscCall(PetscFPrintf(comm, fd, "Compiled without FORTRAN kernels\n")); 1729 } 1730 if (PetscDefined(USE_64BIT_INDICES)) { 1731 PetscCall(PetscFPrintf(comm, fd, "Compiled with 64-bit PetscInt\n")); 1732 } else if (PetscDefined(USE___FLOAT128)) { 1733 PetscCall(PetscFPrintf(comm, fd, "Compiled with 32-bit PetscInt\n")); 1734 } 1735 if (PetscDefined(USE_REAL_SINGLE)) { 1736 PetscCall(PetscFPrintf(comm, fd, "Compiled with single precision PetscScalar and PetscReal\n")); 1737 } else if (PetscDefined(USE___FLOAT128)) { 1738 PetscCall(PetscFPrintf(comm, fd, "Compiled with 128 bit precision PetscScalar and PetscReal\n")); 1739 } 1740 if (PetscDefined(USE_REAL_MAT_SINGLE)) { 1741 PetscCall(PetscFPrintf(comm, fd, "Compiled with single precision matrices\n")); 1742 } else { 1743 PetscCall(PetscFPrintf(comm, fd, "Compiled with full precision matrices (default)\n")); 1744 } 1745 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))); 1746 1747 PetscCall(PetscFPrintf(comm, fd, "Configure options: %s", petscconfigureoptions)); 1748 PetscCall(PetscFPrintf(comm, fd, "%s", petscmachineinfo)); 1749 PetscCall(PetscFPrintf(comm, fd, "%s", petsccompilerinfo)); 1750 PetscCall(PetscFPrintf(comm, fd, "%s", petsccompilerflagsinfo)); 1751 PetscCall(PetscFPrintf(comm, fd, "%s", petsclinkerinfo)); 1752 1753 /* Cleanup */ 1754 PetscCall(PetscFPrintf(comm, fd, "\n")); 1755 PetscCall(PetscLogViewWarnNoGpuAwareMpi(comm, fd)); 1756 PetscCall(PetscLogViewWarnDebugging(comm, fd)); 1757 PetscCall(PetscFPTrapPop()); 1758 PetscFunctionReturn(PETSC_SUCCESS); 1759 } 1760 1761 static PetscErrorCode PetscLogHandlerView_Default(PetscLogHandler handler, PetscViewer viewer) 1762 { 1763 PetscViewerFormat format; 1764 PetscFunctionBegin; 1765 PetscCall(PetscViewerGetFormat(viewer, &format)); 1766 if (format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_ASCII_INFO) { 1767 PetscCall(PetscLogHandlerView_Default_Info(handler, viewer)); 1768 } else if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1769 PetscCall(PetscLogHandlerView_Default_Detailed(handler, viewer)); 1770 } else if (format == PETSC_VIEWER_ASCII_CSV) { 1771 PetscCall(PetscLogHandlerView_Default_CSV(handler, viewer)); 1772 } 1773 PetscFunctionReturn(PETSC_SUCCESS); 1774 } 1775 1776 /*MC 1777 PETSCLOGHANDLERDEFAULT - PETSCLOGHANDLERDEFAULT = "default" - A `PetscLogHandler` that collects data for PETSc 1778 default profiling log viewers (`PetscLogView()` and `PetscLogDump()`). A log handler of this type is 1779 created and started (`PetscLogHandlerStart()`) by `PetscLogDefaultBegin()`. 1780 1781 Options Database Keys: 1782 + -log_include_actions - include a growing list of actions (event beginnings and endings, object creations and destructions) in `PetscLogDump()` (`PetscLogActions()`). 1783 - -log_include_objects - include a growing list of object creations and destructions in `PetscLogDump()` (`PetscLogObjects()`). 1784 1785 Level: developer 1786 1787 .seealso: [](ch_profiling), `PetscLogHandler` 1788 M*/ 1789 1790 PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Default(PetscLogHandler handler) 1791 { 1792 PetscFunctionBegin; 1793 PetscCall(PetscLogHandlerContextCreate_Default((PetscLogHandler_Default *)&handler->data)); 1794 handler->ops->destroy = PetscLogHandlerDestroy_Default; 1795 handler->ops->eventbegin = PetscLogHandlerEventBegin_Default; 1796 handler->ops->eventend = PetscLogHandlerEventEnd_Default; 1797 handler->ops->eventsync = PetscLogHandlerEventSync_Default; 1798 handler->ops->objectcreate = PetscLogHandlerObjectCreate_Default; 1799 handler->ops->objectdestroy = PetscLogHandlerObjectDestroy_Default; 1800 handler->ops->stagepush = PetscLogHandlerStagePush_Default; 1801 handler->ops->stagepop = PetscLogHandlerStagePop_Default; 1802 handler->ops->view = PetscLogHandlerView_Default; 1803 PetscFunctionReturn(PETSC_SUCCESS); 1804 } 1805