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