Actual source code: logdefault.c

  1: #include <petscviewer.h>
  2: #include <petscdevice.h>
  3: #include <petsc/private/logimpl.h>
  4: #include <petsc/private/loghandlerimpl.h>
  5: #include <petsc/private/deviceimpl.h>
  6: #include <petscconfiginfo.h>
  7: #include <petscmachineinfo.h>
  8: #include "logdefault.h"

 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: }

 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: }

 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: }

 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: }

 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: }

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: }

126: static PetscErrorCode PetscEventPerfInfoPause(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event)
127: {
128:   PetscFunctionBegin;
129:   PetscCall(PetscEventPerfInfoToc_Internal(eventInfo, time, logMemory, event, PETSC_TRUE));
130:   PetscFunctionReturn(PETSC_SUCCESS);
131: }

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: }

159: PETSC_LOG_RESIZABLE_ARRAY(EventPerfArray, PetscEventPerfInfo, PetscLogEvent, PetscEventPerfInfoInit, NULL, NULL)

161: /* --- PetscClassPerf --- */

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;

170: static PetscErrorCode PetscClassPerfInit(PetscClassPerf *classInfo)
171: {
172:   PetscFunctionBegin;
173:   PetscCall(PetscMemzero(classInfo, sizeof(*classInfo)));
174:   PetscFunctionReturn(PETSC_SUCCESS);
175: }

177: PETSC_LOG_RESIZABLE_ARRAY(ClassPerfArray, PetscClassPerf, PetscLogClass, PetscClassPerfInit, NULL, NULL)

179: /* --- PetscStagePerf --- */

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;

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: }

198: static PetscErrorCode PetscStageInfoReset(PetscStagePerf *stageInfo)
199: {
200:   PetscFunctionBegin;
201:   PetscCall(PetscLogEventPerfArrayDestroy(&stageInfo->eventLog));
202:   PetscCall(PetscLogClassPerfArrayDestroy(&stageInfo->classLog));
203:   PetscFunctionReturn(PETSC_SUCCESS);
204: }

206: PETSC_LOG_RESIZABLE_ARRAY(StageInfoArray, PetscStagePerf, PetscLogStage, PetscStageInfoInit, PetscStageInfoReset, NULL)

208: /* --- Action --- */

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;

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:   PetscObjectId      id1, id2, id3; /* The ids of associated objects */
227: } Action;

229: PETSC_LOG_RESIZABLE_ARRAY(ActionArray, Action, PetscLogEvent, NULL, NULL, NULL)

231: /* --- Object --- */

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;

242: PETSC_LOG_RESIZABLE_ARRAY(ObjectArray, Object, PetscObject, NULL, NULL, NULL)

244: /* Map from (threadid,stage,event) to perfInfo data struct */
245: #include <petsc/private/hashmapijk.h>

247: PETSC_HASH_MAP(HMapEvent, PetscHashIJKKey, PetscEventPerfInfo *, PetscHashIJKKeyHash, PetscHashIJKKeyEqual, NULL)

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: };

264: /* --- PetscLogHandler_Default --- */

266: static PetscErrorCode PetscLogHandlerContextCreate_Default(PetscLogHandler_Default *def_p)
267: {
268:   PetscLogHandler_Default def;

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:   PetscCall(PetscSpinlockCreate(&def->lock));

278:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_include_actions", &def->petsc_logActions, NULL));
279:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_include_objects", &def->petsc_logObjects, NULL));
280:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_handler_default_use_threadsafe_events", &def->use_threadsafe, NULL));
281:   if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { PetscCall(PetscHMapEventCreate(&def->eventInfoMap_th)); }
282:   PetscFunctionReturn(PETSC_SUCCESS);
283: }

285: static PetscErrorCode PetscLogHandlerDestroy_Default(PetscLogHandler h)
286: {
287:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;

289:   PetscFunctionBegin;
290:   PetscCall(PetscLogStageInfoArrayDestroy(&def->stages));
291:   PetscCall(PetscLogActionArrayDestroy(&def->petsc_actions));
292:   PetscCall(PetscLogObjectArrayDestroy(&def->petsc_objects));
293:   PetscCall(PetscSpinlockDestroy(&def->lock));
294:   if (def->eventInfoMap_th) {
295:     PetscEventPerfInfo **array;
296:     PetscInt             n, off = 0;

298:     PetscCall(PetscHMapEventGetSize(def->eventInfoMap_th, &n));
299:     PetscCall(PetscMalloc1(n, &array));
300:     PetscCall(PetscHMapEventGetVals(def->eventInfoMap_th, &off, array));
301:     for (PetscInt i = 0; i < n; i++) PetscCall(PetscFree(array[i]));
302:     PetscCall(PetscFree(array));
303:     PetscCall(PetscHMapEventDestroy(&def->eventInfoMap_th));
304:   }
305:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerGetEventPerfInfo_C", NULL));
306:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerGetStagePerfInfo_C", NULL));
307:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerSetLogActions_C", NULL));
308:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerSetLogObjects_C", NULL));
309:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerLogObjectState_C", NULL));
310:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerGetNumObjects_C", NULL));
311:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventDeactivatePush_C", NULL));
312:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventDeactivatePop_C", NULL));
313:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventsPause_C", NULL));
314:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventsResume_C", NULL));
315:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerDump_C", NULL));
316:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerStageSetVisible_C", NULL));
317:   PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerStageGetVisible_C", NULL));
318:   PetscCall(PetscFree(def));
319:   PetscFunctionReturn(PETSC_SUCCESS);
320: }

322: static PetscErrorCode PetscLogHandlerDefaultGetStageInfo(PetscLogHandler handler, PetscLogStage stage, PetscStagePerf **stage_info_p)
323: {
324:   PetscStagePerf         *stage_info = NULL;
325:   PetscLogHandler_Default def        = (PetscLogHandler_Default)handler->data;

327:   PetscFunctionBegin;
328:   PetscCall(PetscLogStageInfoArrayResize(def->stages, stage + 1));
329:   PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
330:   *stage_info_p = stage_info;
331:   PetscFunctionReturn(PETSC_SUCCESS);
332: }

334: static PetscErrorCode PetscLogHandlerGetEventPerfInfo_Default(PetscLogHandler handler, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **event_info_p)
335: {
336:   PetscEventPerfInfo    *event_info = NULL;
337:   PetscStagePerf        *stage_info = NULL;
338:   PetscLogEventPerfArray event_log;

340:   PetscFunctionBegin;
341:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(handler->state, &stage));
342:   PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info));
343:   event_log = stage_info->eventLog;
344:   PetscCall(PetscLogEventPerfArrayResize(event_log, event + 1));
345:   PetscCall(PetscLogEventPerfArrayGetRef(event_log, event, &event_info));
346:   event_info->id = event;
347:   *event_info_p  = event_info;
348:   PetscFunctionReturn(PETSC_SUCCESS);
349: }

351: static PetscErrorCode PetscLogHandlerGetStagePerfInfo_Default(PetscLogHandler handler, PetscLogStage stage, PetscEventPerfInfo **stage_info_p)
352: {
353:   PetscStagePerf *stage_perf_info = NULL;

355:   PetscFunctionBegin;
356:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(handler->state, &stage));
357:   PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_perf_info));
358:   *stage_info_p = &stage_perf_info->perfInfo;
359:   PetscFunctionReturn(PETSC_SUCCESS);
360: }

362: static PetscErrorCode PetscLogHandlerDefaultGetClassPerf(PetscLogHandler handler, PetscLogStage stage, PetscLogClass clss, PetscClassPerf **class_info)
363: {
364:   PetscLogClassPerfArray class_log;
365:   PetscStagePerf        *stage_info;

367:   PetscFunctionBegin;
368:   PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info));
369:   class_log = stage_info->classLog;
370:   PetscCall(PetscLogClassPerfArrayResize(class_log, clss + 1));
371:   PetscCall(PetscLogClassPerfArrayGetRef(class_log, clss, class_info));
372:   PetscFunctionReturn(PETSC_SUCCESS);
373: }

375: static PetscErrorCode PetscLogHandlerObjectCreate_Default(PetscLogHandler h, PetscObject obj)
376: {
377:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
378:   PetscLogState           state;
379:   PetscLogStage           stage;
380:   PetscClassPerf         *classInfo;
381:   int                     oclass = 0;

383:   PetscFunctionBegin;
384:   PetscCall(PetscLogHandlerGetState(h, &state));
385:   PetscCall(PetscSpinlockLock(&def->lock));
386:   /* Record stage info */
387:   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
388:   PetscCall(PetscLogStateGetClassFromClassId(state, obj->classid, &oclass));
389:   PetscCall(PetscLogHandlerDefaultGetClassPerf(h, stage, oclass, &classInfo));
390:   classInfo->creations++;
391:   /* Record the creation action */
392:   if (def->petsc_logActions) {
393:     Action new_action;

395:     PetscCall(PetscTime(&new_action.time));
396:     new_action.time -= petsc_BaseTime;
397:     new_action.action  = PETSC_LOG_ACTION_CREATE;
398:     new_action.event   = -1;
399:     new_action.classid = obj->classid;
400:     new_action.id1     = obj->id;
401:     new_action.id2     = -1;
402:     new_action.id3     = -1;
403:     new_action.flops   = petsc_TotalFlops;
404:     PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
405:     PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
406:     PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
407:   }
408:   /* We don't just use obj->id to count all objects that are created
409:      because PetscLogHandlers are objects and PetscLogObjectDestroy() will not
410:      be called for them: the number of objects created and destroyed as counted
411:      here and below would have an imbalance */
412:   def->petsc_numObjectsCreated++;
413:   /* Record the object */
414:   if (def->petsc_logObjects) {
415:     Object   new_object;
416:     PetscInt objid;

418:     new_object.parent = -1;
419:     new_object.obj    = obj;
420:     new_object.mem    = 0;

422:     PetscCall(PetscMemzero(new_object.name, sizeof(new_object.name)));
423:     PetscCall(PetscMemzero(new_object.info, sizeof(new_object.info)));
424:     PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1");
425:     PetscCall(PetscIntCast(obj->id, &objid));
426:     PetscCall(PetscLogObjectArrayResize(def->petsc_objects, objid));
427:     PetscCall(PetscLogObjectArraySet(def->petsc_objects, objid - 1, new_object));
428:   }
429:   PetscCall(PetscSpinlockUnlock(&def->lock));
430:   PetscFunctionReturn(PETSC_SUCCESS);
431: }

433: static PetscErrorCode PetscLogHandlerObjectDestroy_Default(PetscLogHandler h, PetscObject obj)
434: {
435:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
436:   PetscLogState           state;
437:   PetscLogStage           stage;
438:   PetscClassPerf         *classInfo;
439:   int                     oclass = 0;

441:   PetscFunctionBegin;
442:   PetscCall(PetscLogHandlerGetState(h, &state));
443:   /* Record stage info */
444:   PetscCall(PetscSpinlockLock(&def->lock));
445:   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
446:   if (stage >= 0) {
447:     /* stage < 0 can happen if the log summary is output before some things are destroyed */
448:     PetscCall(PetscLogStateGetClassFromClassId(state, obj->classid, &oclass));
449:     PetscCall(PetscLogHandlerDefaultGetClassPerf(h, stage, oclass, &classInfo));
450:     classInfo->destructions++;
451:   }
452:   /* Cannot Credit all ancestors with your memory because they may have already been destroyed*/
453:   def->petsc_numObjectsDestroyed++;
454:   /* Dynamically enlarge logging structures */
455:   /* Record the destruction action */
456:   if (def->petsc_logActions) {
457:     Action new_action;

459:     PetscCall(PetscTime(&new_action.time));
460:     new_action.time -= petsc_BaseTime;
461:     new_action.event   = -1;
462:     new_action.action  = PETSC_LOG_ACTION_DESTROY;
463:     new_action.classid = obj->classid;
464:     new_action.id1     = obj->id;
465:     new_action.id2     = -1;
466:     new_action.id3     = -1;
467:     new_action.flops   = petsc_TotalFlops;
468:     PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
469:     PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
470:     PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
471:   }
472:   if (def->petsc_logObjects) {
473:     Object  *obj_entry = NULL;
474:     PetscInt objid;

476:     PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1");
477:     PetscCall(PetscIntCast(obj->id - 1, &objid));
478:     PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, objid, &obj_entry));
479:     if (obj->name) PetscCall(PetscStrncpy(obj_entry->name, obj->name, 64));
480:     obj_entry->obj = NULL;
481:   }
482:   PetscCall(PetscSpinlockUnlock(&def->lock));
483:   PetscFunctionReturn(PETSC_SUCCESS);
484: }

486: static PetscErrorCode PetscLogHandlerEventSync_Default(PetscLogHandler h, PetscLogEvent event, MPI_Comm comm)
487: {
488:   PetscLogState       state;
489:   PetscLogEventInfo   event_info;
490:   PetscEventPerfInfo *event_perf_info;
491:   int                 stage;
492:   PetscLogDouble      time = 0.0;

494:   PetscFunctionBegin;
495:   if (!PetscLogSyncOn || comm == MPI_COMM_NULL) PetscFunctionReturn(PETSC_SUCCESS);
496:   PetscCall(PetscLogHandlerGetState(h, &state));
497:   PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
498:   if (!event_info.collective) PetscFunctionReturn(PETSC_SUCCESS);
499:   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
500:   PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
501:   if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS);

503:   PetscCall(PetscTimeSubtract(&time));
504:   PetscCallMPI(MPI_Barrier(comm));
505:   PetscCall(PetscTimeAdd(&time));
506:   event_perf_info->syncTime += time;
507:   PetscFunctionReturn(PETSC_SUCCESS);
508: }

510: static PetscErrorCode PetscLogGetStageEventPerfInfo_threaded(PetscLogHandler_Default def, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **eventInfo)
511: {
512:   PetscEventPerfInfo *leventInfo = NULL;
513:   PetscHashIJKKey     key;

515:   PetscFunctionBegin;
516: #if PetscDefined(HAVE_THREADSAFETY)
517:   key.i = PetscLogGetTid();
518: #else
519:   key.i = 0;
520: #endif
521:   key.j = stage;
522:   key.k = event;
523:   PetscCall(PetscSpinlockLock(&def->lock));
524:   PetscCall(PetscHMapEventGet(def->eventInfoMap_th, key, &leventInfo));
525:   if (!leventInfo) {
526:     PetscCall(PetscNew(&leventInfo));
527:     leventInfo->id = event;
528:     PetscCall(PetscHMapEventSet(def->eventInfoMap_th, key, leventInfo));
529:   }
530:   PetscCall(PetscSpinlockUnlock(&def->lock));
531:   *eventInfo = leventInfo;
532:   PetscFunctionReturn(PETSC_SUCCESS);
533: }

535: static PetscErrorCode PetscLogHandlerEventBegin_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
536: {
537:   PetscLogHandler_Default def             = (PetscLogHandler_Default)h->data;
538:   PetscEventPerfInfo     *event_perf_info = NULL;
539:   PetscLogEventInfo       event_info;
540:   PetscLogDouble          time;
541:   PetscLogState           state;
542:   PetscLogStage           stage;

544:   PetscFunctionBegin;
545:   PetscCall(PetscLogHandlerGetState(h, &state));
546:   if (PetscDefined(USE_DEBUG)) {
547:     PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
552:     if (event_info.collective && o1) {
553:       PetscInt64 b1[2], b2[2];

555:       b1[0] = -o1->cidx;
556:       b1[1] = o1->cidx;
557:       PetscCallMPI(MPIU_Allreduce(b1, b2, 2, MPIU_INT64, MPI_MAX, PetscObjectComm(o1)));
558:       PetscCheck(-b2[0] == b2[1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Collective event %s not called collectively %" PetscInt64_FMT " != %" PetscInt64_FMT, event_info.name, -b2[0], b2[1]);
559:     }
560:   }
561:   /* Synchronization */
562:   PetscCall(PetscLogHandlerEventSync_Default(h, event, PetscObjectComm(o1)));
563:   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
564:   if (def->pause_depth > 0) stage = 0; // in pause-mode, all events run on the main stage
565:   if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
566:     PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info));
567:     if (event_perf_info->depth == 0) { PetscCall(PetscEventPerfInfoInit(event_perf_info)); }
568:   } else {
569:     PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
570:   }
571:   PetscCheck(event_perf_info->depth >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to begin a paused event, this is not allowed");
572:   event_perf_info->depth++;
573:   /* Check for double counting */
574:   if (event_perf_info->depth > 1) PetscFunctionReturn(PETSC_SUCCESS);
575:   PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
576:   /* Log the performance info */
577:   event_perf_info->count++;
578:   PetscCall(PetscTime(&time));
579:   PetscCall(PetscEventPerfInfoTic(event_perf_info, time, PetscLogMemory, (int)event));
580:   if (def->petsc_logActions) {
581:     PetscLogDouble curTime;
582:     Action         new_action;

584:     PetscCall(PetscTime(&curTime));
585:     new_action.time    = curTime - petsc_BaseTime;
586:     new_action.action  = PETSC_LOG_ACTION_BEGIN;
587:     new_action.event   = event;
588:     new_action.classid = event_info.classid;
589:     new_action.id1     = o1 ? o1->id : -1;
590:     new_action.id2     = o2 ? o2->id : -1;
591:     new_action.id3     = o3 ? o3->id : -1;
592:     new_action.flops   = petsc_TotalFlops;
593:     PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
594:     PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
595:     PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
596:   }
597:   PetscFunctionReturn(PETSC_SUCCESS);
598: }

600: static PetscErrorCode PetscLogHandlerEventEnd_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
601: {
602:   PetscLogHandler_Default def             = (PetscLogHandler_Default)h->data;
603:   PetscEventPerfInfo     *event_perf_info = NULL;
604:   PetscLogDouble          time;
605:   PetscLogState           state;
606:   int                     stage;
607:   PetscLogEventInfo       event_info;

609:   PetscFunctionBegin;
610:   PetscCall(PetscLogHandlerGetState(h, &state));
611:   if (PetscDefined(USE_DEBUG)) {
612:     PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
617:     if (event_info.collective && o1) {
618:       PetscInt64 b1[2], b2[2];

620:       b1[0] = -o1->cidx;
621:       b1[1] = o1->cidx;
622:       PetscCallMPI(MPIU_Allreduce(b1, b2, 2, MPIU_INT64, MPI_MAX, PetscObjectComm(o1)));
623:       PetscCheck(-b2[0] == b2[1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Collective event %s not called collectively %" PetscInt64_FMT " != %" PetscInt64_FMT, event_info.name, -b2[0], b2[1]);
624:     }
625:   }
626:   if (def->petsc_logActions) {
627:     PetscLogDouble curTime;
628:     Action         new_action;

630:     PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
631:     PetscCall(PetscTime(&curTime));
632:     new_action.time    = curTime - petsc_BaseTime;
633:     new_action.action  = PETSC_LOG_ACTION_END;
634:     new_action.event   = event;
635:     new_action.classid = event_info.classid;
636:     new_action.id1     = o1 ? o1->id : -1;
637:     new_action.id2     = o2 ? o2->id : -2;
638:     new_action.id3     = o3 ? o3->id : -3;
639:     new_action.flops   = petsc_TotalFlops;
640:     PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
641:     PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
642:     PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
643:   }
644:   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
645:   if (def->pause_depth > 0) stage = 0; // all events run on the main stage in pause-mode
646:   if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
647:     PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info));
648:   } else {
649:     PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
650:   }
651:   PetscCheck(event_perf_info->depth > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to end paused event, not allowed");
652:   event_perf_info->depth--;
653:   /* Check for double counting */
654:   if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
655:   else PetscCheck(event_perf_info->depth == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Logging event had unbalanced begin/end pairs");

657:   /* Log performance info */
658:   PetscCall(PetscTime(&time));
659:   PetscCall(PetscEventPerfInfoToc(event_perf_info, time, PetscLogMemory, (int)event));
660:   if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
661:     PetscEventPerfInfo *event_perf_info_global;
662:     PetscCall(PetscSpinlockLock(&def->lock));
663:     PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info_global));
664:     PetscCall(PetscEventPerfInfoAdd_Internal(event_perf_info, event_perf_info_global));
665:     PetscCall(PetscSpinlockUnlock(&def->lock));
666:   }
667:   PetscFunctionReturn(PETSC_SUCCESS);
668: }

670: static PetscErrorCode PetscLogHandlerEventDeactivatePush_Default(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event)
671: {
672:   PetscEventPerfInfo *event_perf_info;

674:   PetscFunctionBegin;
675:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
676:   PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
677:   event_perf_info->depth++;
678:   PetscFunctionReturn(PETSC_SUCCESS);
679: }

681: static PetscErrorCode PetscLogHandlerEventDeactivatePop_Default(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event)
682: {
683:   PetscEventPerfInfo *event_perf_info;

685:   PetscFunctionBegin;
686:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
687:   PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
688:   event_perf_info->depth--;
689:   PetscFunctionReturn(PETSC_SUCCESS);
690: }

692: static PetscErrorCode PetscLogHandlerEventsPause_Default(PetscLogHandler h)
693: {
694:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
695:   PetscLogDouble          time;
696:   PetscInt                num_stages;

698:   PetscFunctionBegin;
699:   if (def->pause_depth++ > 0) PetscFunctionReturn(PETSC_SUCCESS);
700:   PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL));
701:   PetscCall(PetscTime(&time));
702:   for (PetscInt stage = 0; stage < num_stages; stage++) {
703:     PetscStagePerf *stage_info = NULL;
704:     PetscInt        num_events;

706:     PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
707:     PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL));
708:     for (PetscInt event = 0; event < num_events; event++) {
709:       PetscEventPerfInfo *event_info = NULL;
710:       PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info));
711:       if (event_info->depth > 0) {
712:         event_info->depth *= -1;
713:         PetscCall(PetscEventPerfInfoPause(event_info, time, PetscLogMemory, event));
714:       }
715:     }
716:     if (stage > 0 && stage_info->perfInfo.depth > 0) {
717:       stage_info->perfInfo.depth *= -1;
718:       PetscCall(PetscEventPerfInfoPause(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2)));
719:     }
720:   }
721:   PetscFunctionReturn(PETSC_SUCCESS);
722: }

724: static PetscErrorCode PetscLogHandlerEventsResume_Default(PetscLogHandler h)
725: {
726:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
727:   PetscLogDouble          time;
728:   PetscInt                num_stages;

730:   PetscFunctionBegin;
731:   if (--def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
732:   PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL));
733:   PetscCall(PetscTime(&time));
734:   for (PetscInt stage = 0; stage < num_stages; stage++) {
735:     PetscStagePerf *stage_info = NULL;
736:     PetscInt        num_events;

738:     PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
739:     PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL));
740:     for (PetscInt event = 0; event < num_events; event++) {
741:       PetscEventPerfInfo *event_info = NULL;
742:       PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info));
743:       if (event_info->depth < 0) {
744:         event_info->depth *= -1;
745:         PetscCall(PetscEventPerfInfoResume(event_info, time, PetscLogMemory, event));
746:       }
747:     }
748:     if (stage > 0 && stage_info->perfInfo.depth < 0) {
749:       stage_info->perfInfo.depth *= -1;
750:       PetscCall(PetscEventPerfInfoResume(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2)));
751:     }
752:   }
753:   PetscFunctionReturn(PETSC_SUCCESS);
754: }

756: static PetscErrorCode PetscLogHandlerStagePush_Default(PetscLogHandler h, PetscLogStage new_stage)
757: {
758:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
759:   PetscLogDouble          time;
760:   PetscLogState           state;
761:   PetscLogStage           current_stage;
762:   PetscStagePerf         *new_stage_info;

764:   PetscFunctionBegin;
765:   if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
766:   PetscCall(PetscLogHandlerGetState(h, &state));
767:   current_stage = state->current_stage;
768:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, new_stage, &new_stage_info));
769:   PetscCall(PetscTime(&time));

771:   /* Record flops/time of previous stage */
772:   if (current_stage >= 0) {
773:     if (PetscBTLookup(state->active, current_stage)) {
774:       PetscStagePerf *current_stage_info;
775:       PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, &current_stage_info));
776:       PetscCall(PetscEventPerfInfoToc(&current_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2)));
777:     }
778:   }
779:   new_stage_info->used = PETSC_TRUE;
780:   new_stage_info->perfInfo.count++;
781:   new_stage_info->perfInfo.depth++;
782:   /* Subtract current quantities so that we obtain the difference when we pop */
783:   if (PetscBTLookup(state->active, new_stage)) PetscCall(PetscEventPerfInfoTic(&new_stage_info->perfInfo, time, PetscLogMemory, (int)-(new_stage + 2)));
784:   PetscFunctionReturn(PETSC_SUCCESS);
785: }

787: static PetscErrorCode PetscLogHandlerStagePop_Default(PetscLogHandler h, PetscLogStage old_stage)
788: {
789:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
790:   PetscLogStage           current_stage;
791:   PetscStagePerf         *old_stage_info;
792:   PetscLogState           state;
793:   PetscLogDouble          time;

795:   PetscFunctionBegin;
796:   if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
797:   PetscCall(PetscLogHandlerGetState(h, &state));
798:   current_stage = state->current_stage;
799:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, old_stage, &old_stage_info));
800:   PetscCall(PetscTime(&time));
801:   old_stage_info->perfInfo.depth--;
802:   if (PetscBTLookup(state->active, old_stage)) { PetscCall(PetscEventPerfInfoToc(&old_stage_info->perfInfo, time, PetscLogMemory, (int)-(old_stage + 2))); }
803:   if (current_stage >= 0) {
804:     if (PetscBTLookup(state->active, current_stage)) {
805:       PetscStagePerf *current_stage_info;
806:       PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, &current_stage_info));
807:       PetscCall(PetscEventPerfInfoTic(&current_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2)));
808:     }
809:   }
810:   PetscFunctionReturn(PETSC_SUCCESS);
811: }

813: static PetscErrorCode PetscLogHandlerStageSetVisible_Default(PetscLogHandler h, PetscLogStage stage, PetscBool is_visible)
814: {
815:   PetscStagePerf *stage_info;

817:   PetscFunctionBegin;
818:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
819:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info));
820:   stage_info->perfInfo.visible = is_visible;
821:   PetscFunctionReturn(PETSC_SUCCESS);
822: }

824: static PetscErrorCode PetscLogHandlerStageGetVisible_Default(PetscLogHandler h, PetscLogStage stage, PetscBool *is_visible)
825: {
826:   PetscStagePerf *stage_info;

828:   PetscFunctionBegin;
829:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
830:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info));
831:   *is_visible = stage_info->perfInfo.visible;
832:   PetscFunctionReturn(PETSC_SUCCESS);
833: }

835: static PetscErrorCode PetscLogHandlerSetLogActions_Default(PetscLogHandler handler, PetscBool flag)
836: {
837:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;

839:   PetscFunctionBegin;
840:   def->petsc_logActions = flag;
841:   PetscFunctionReturn(PETSC_SUCCESS);
842: }

844: static PetscErrorCode PetscLogHandlerSetLogObjects_Default(PetscLogHandler handler, PetscBool flag)
845: {
846:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;

848:   PetscFunctionBegin;
849:   def->petsc_logObjects = flag;
850:   PetscFunctionReturn(PETSC_SUCCESS);
851: }

853: static PetscErrorCode PetscLogHandlerLogObjectState_Default(PetscLogHandler handler, PetscObject obj, const char format[], va_list Argp)
854: {
855:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
856:   size_t                  fullLength;

858:   PetscFunctionBegin;
859:   if (def->petsc_logObjects) {
860:     Object  *obj_entry = NULL;
861:     PetscInt objid;

863:     PetscCall(PetscIntCast(obj->id - 1, &objid));
864:     PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, objid, &obj_entry));
865:     PetscCall(PetscVSNPrintf(obj_entry->info, 64, format, &fullLength, Argp));
866:   }
867:   PetscFunctionReturn(PETSC_SUCCESS);
868: }

870: static PetscErrorCode PetscLogHandlerGetNumObjects_Default(PetscLogHandler handler, PetscInt *num_objects)
871: {
872:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;

874:   PetscFunctionBegin;
875:   PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, num_objects, NULL));
876:   PetscFunctionReturn(PETSC_SUCCESS);
877: }

879: static PetscErrorCode PetscLogHandlerDump_Default(PetscLogHandler handler, const char sname[])
880: {
881:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
882:   FILE                   *fd;
883:   char                    file[PETSC_MAX_PATH_LEN], fname[PETSC_MAX_PATH_LEN];
884:   PetscLogDouble          flops, _TotalTime;
885:   PetscMPIInt             rank;
886:   int                     curStage;
887:   PetscLogState           state;
888:   PetscInt                num_events;
889:   PetscLogEvent           event;

891:   PetscFunctionBegin;
892:   /* Calculate the total elapsed time */
893:   PetscCall(PetscTime(&_TotalTime));
894:   _TotalTime -= petsc_BaseTime;
895:   /* Open log file */
896:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)handler), &rank));
897:   PetscCall(PetscSNPrintf(file, PETSC_STATIC_ARRAY_LENGTH(file), "%s.%d", sname && sname[0] ? sname : "Log", rank));
898:   PetscCall(PetscFixFilename(file, fname));
899:   PetscCall(PetscFOpen(PETSC_COMM_SELF, fname, "w", &fd));
900:   PetscCheck(!(rank == 0) || !(!fd), PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file: %s", fname);
901:   /* Output totals */
902:   PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Total Flop %14e %16.8e\n", petsc_TotalFlops, _TotalTime));
903:   PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Clock Resolution %g\n", 0.0));
904:   /* Output actions */
905:   if (def->petsc_logActions) {
906:     PetscInt num_actions;
907:     PetscCall(PetscLogActionArrayGetSize(def->petsc_actions, &num_actions, NULL));
908:     PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Actions accomplished %d\n", (int)num_actions));
909:     for (int a = 0; a < num_actions; a++) {
910:       Action *action;

912:       PetscCall(PetscLogActionArrayGetRef(def->petsc_actions, a, &action));
913:       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "%g %d %d %d  %" PetscInt64_FMT " %" PetscInt64_FMT " %" PetscInt64_FMT " %g %g %g\n", action->time, action->action, (int)action->event, (int)action->classid, action->id1, action->id2, action->id3,
914:                              action->flops, action->mem, action->maxmem));
915:     }
916:   }
917:   /* Output objects */
918:   if (def->petsc_logObjects) {
919:     PetscInt num_objects;

921:     PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
922:     PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Objects created %d destroyed %d\n", def->petsc_numObjectsCreated, def->petsc_numObjectsDestroyed));
923:     for (int o = 0; o < num_objects; o++) {
924:       Object *object = NULL;

926:       PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, o, &object));
927:       if (object->parent != -1) continue; // object with this id wasn't logged, probably a PetscLogHandler
928:       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Parent ID: %d Memory: %d\n", object->parent, (int)object->mem));
929:       if (!object->name[0]) {
930:         PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Name\n"));
931:       } else {
932:         PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Name: %s\n", object->name));
933:       }
934:       if (!object->info[0]) {
935:         PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Info\n"));
936:       } else {
937:         PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Info: %s\n", object->info));
938:       }
939:     }
940:   }
941:   /* Output events */
942:   PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Event log:\n"));
943:   PetscCall(PetscLogHandlerGetState(handler, &state));
944:   PetscCall(PetscLogStateGetNumEvents(state, &num_events));
945:   PetscCall(PetscLogStateGetCurrentStage(state, &curStage));
946:   for (event = 0; event < num_events; event++) {
947:     PetscEventPerfInfo *event_info;

949:     PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, curStage, event, &event_info));
950:     if (event_info->time != 0.0) flops = event_info->flops / event_info->time;
951:     else flops = 0.0;
952:     PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "%d %16d %16g %16g %16g\n", event, event_info->count, event_info->flops, event_info->time, flops));
953:   }
954:   PetscCall(PetscFClose(PETSC_COMM_SELF, fd));
955:   PetscFunctionReturn(PETSC_SUCCESS);
956: }

958: /*
959:   PetscLogView_Detailed - Each process prints the times for its own events

961: */
962: static PetscErrorCode PetscLogHandlerView_Default_Detailed(PetscLogHandler handler, PetscViewer viewer)
963: {
964:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
965:   PetscLogDouble          locTotalTime, numRed, maxMem;
966:   PetscInt                numStages, numEvents;
967:   MPI_Comm                comm = PetscObjectComm((PetscObject)viewer);
968:   PetscMPIInt             rank, size;
969:   PetscLogGlobalNames     global_stages, global_events;
970:   PetscLogState           state;
971:   PetscEventPerfInfo      zero_info;

973:   PetscFunctionBegin;
974:   PetscCall(PetscLogHandlerGetState(handler, &state));
975:   PetscCallMPI(MPI_Comm_size(comm, &size));
976:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
977:   /* Must preserve reduction count before we go on */
978:   numRed = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
979:   /* Get the total elapsed time */
980:   PetscCall(PetscTime(&locTotalTime));
981:   locTotalTime -= petsc_BaseTime;
982:   PetscCall(PetscViewerASCIIPrintf(viewer, "size = %d\n", size));
983:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalTimes = {}\n"));
984:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessages = {}\n"));
985:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessageLens = {}\n"));
986:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalReductions = {}\n"));
987:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalFlop = {}\n"));
988:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalObjects = {}\n"));
989:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMemory = {}\n"));
990:   PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
991:   PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
992:   PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
993:   PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
994:   PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
995:   PetscCall(PetscViewerASCIIPrintf(viewer, "Stages = {}\n"));
996:   for (PetscInt stage = 0; stage < numStages; stage++) {
997:     PetscInt    stage_id;
998:     const char *stage_name;

1000:     PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1001:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1002:     PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"] = {}\n", stage_name));
1003:     PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"summary\"] = {}\n", stage_name));
1004:     for (PetscInt event = 0; event < numEvents; event++) {
1005:       PetscEventPerfInfo *eventInfo = &zero_info;
1006:       PetscBool           is_zero   = PETSC_FALSE;
1007:       PetscInt            event_id;
1008:       const char         *event_name;

1010:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1011:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1012:       if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1013:       is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE;
1014:       PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1015:       if (!is_zero) { PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"%s\"] = {}\n", stage_name, event_name)); }
1016:     }
1017:   }
1018:   PetscCall(PetscMallocGetMaximumUsage(&maxMem));
1019:   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1020:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalTimes[%d] = %g\n", rank, locTotalTime));
1021:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessages[%d] = %g\n", rank, petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct));
1022:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessageLens[%d] = %g\n", rank, petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len));
1023:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalReductions[%d] = %g\n", rank, numRed));
1024:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalFlop[%d] = %g\n", rank, petsc_TotalFlops));
1025:   {
1026:     PetscInt num_objects;

1028:     PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
1029:     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalObjects[%d] = %d\n", rank, (int)num_objects));
1030:   }
1031:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMemory[%d] = %g\n", rank, maxMem));
1032:   PetscCall(PetscViewerFlush(viewer));
1033:   for (PetscInt stage = 0; stage < numStages; stage++) {
1034:     PetscEventPerfInfo *stage_perf_info = &zero_info;
1035:     PetscInt            stage_id;
1036:     const char         *stage_name;

1038:     PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1039:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1040:     if (stage_id >= 0) {
1041:       PetscStagePerf *stage_info;
1042:       PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info));
1043:       stage_perf_info = &stage_info->perfInfo;
1044:     }
1045:     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,
1046:                                                  stage_perf_info->numMessages, stage_perf_info->messageLength, stage_perf_info->numReductions, stage_perf_info->flops));
1047:     for (PetscInt event = 0; event < numEvents; event++) {
1048:       PetscEventPerfInfo *eventInfo = &zero_info;
1049:       PetscBool           is_zero   = PETSC_FALSE;
1050:       PetscInt            event_id;
1051:       const char         *event_name;

1053:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1054:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1055:       if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1056:       is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE;
1057:       PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero));
1058:       PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1059:       if (!is_zero) {
1060:         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,
1061:                                                      eventInfo->count, eventInfo->time, eventInfo->syncTime, eventInfo->numMessages, eventInfo->messageLength, eventInfo->numReductions, eventInfo->flops));
1062:         if (eventInfo->dof[0] >= 0.) {
1063:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"dof\" : ["));
1064:           for (PetscInt d = 0; d < 8; ++d) {
1065:             if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", "));
1066:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->dof[d]));
1067:           }
1068:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]"));
1069:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"error\" : ["));
1070:           for (PetscInt e = 0; e < 8; ++e) {
1071:             if (e > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", "));
1072:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->errors[e]));
1073:           }
1074:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]"));
1075:         }
1076:       }
1077:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "}\n"));
1078:     }
1079:   }
1080:   PetscCall(PetscViewerFlush(viewer));
1081:   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1082:   PetscCall(PetscLogGlobalNamesDestroy(&global_events));
1083:   PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1084:   PetscFunctionReturn(PETSC_SUCCESS);
1085: }

1087: /*
1088:   PetscLogView_CSV - Each process prints the times for its own events in Comma-Separated Value Format
1089: */
1090: static PetscErrorCode PetscLogHandlerView_Default_CSV(PetscLogHandler handler, PetscViewer viewer)
1091: {
1092:   PetscLogDouble      locTotalTime, maxMem;
1093:   PetscInt            numStages, numEvents, stage, event;
1094:   MPI_Comm            comm = PetscObjectComm((PetscObject)viewer);
1095:   PetscMPIInt         rank, size;
1096:   PetscLogGlobalNames global_stages, global_events;
1097:   PetscLogState       state;
1098:   PetscEventPerfInfo  zero_info;

1100:   PetscFunctionBegin;
1101:   PetscCall(PetscLogHandlerGetState(handler, &state));
1102:   PetscCallMPI(MPI_Comm_size(comm, &size));
1103:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1104:   /* Must preserve reduction count before we go on */
1105:   /* Get the total elapsed time */
1106:   PetscCall(PetscTime(&locTotalTime));
1107:   locTotalTime -= petsc_BaseTime;
1108:   PetscCall(PetscMallocGetMaximumUsage(&maxMem));
1109:   PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
1110:   PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
1111:   PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
1112:   PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
1113:   PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
1114:   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1115:   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));
1116:   PetscCall(PetscViewerFlush(viewer));
1117:   for (stage = 0; stage < numStages; stage++) {
1118:     PetscEventPerfInfo *stage_perf_info;
1119:     PetscInt            stage_id;
1120:     const char         *stage_name;

1122:     PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1123:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1124:     stage_perf_info = &zero_info;
1125:     if (stage_id >= 0) {
1126:       PetscStagePerf *stage_info;
1127:       PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info));
1128:       stage_perf_info = &stage_info->perfInfo;
1129:     }
1130:     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));
1131:     for (event = 0; event < numEvents; event++) {
1132:       PetscEventPerfInfo *eventInfo = &zero_info;
1133:       PetscBool           is_zero   = PETSC_FALSE;
1134:       PetscInt            event_id;
1135:       const char         *event_name;

1137:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1138:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1139:       if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1140:       PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero));
1141:       PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1142:       if (!is_zero) {
1143:         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));
1144:         if (eventInfo->dof[0] >= 0.) {
1145:           for (PetscInt d = 0; d < 8; ++d) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->dof[d]));
1146:           for (PetscInt e = 0; e < 8; ++e) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->errors[e]));
1147:         }
1148:       }
1149:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
1150:     }
1151:   }
1152:   PetscCall(PetscViewerFlush(viewer));
1153:   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1154:   PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1155:   PetscCall(PetscLogGlobalNamesDestroy(&global_events));
1156:   PetscFunctionReturn(PETSC_SUCCESS);
1157: }

1159: static PetscErrorCode PetscLogViewWarnSync(PetscViewer viewer)
1160: {
1161:   PetscFunctionBegin;
1162:   if (!PetscLogSyncOn) PetscFunctionReturn(PETSC_SUCCESS);
1163:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1164:   PetscCall(PetscViewerASCIIPrintf(viewer, "      ##########################################################\n"));
1165:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1166:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                       WARNING!!!                       #\n"));
1167:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1168:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   This program was run with logging synchronization.   #\n"));
1169:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   This option provides more meaningful imbalance       #\n"));
1170:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   figures at the expense of slowing things down and    #\n"));
1171:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   providing a distorted view of the overall runtime.   #\n"));
1172:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1173:   PetscCall(PetscViewerASCIIPrintf(viewer, "      ##########################################################\n\n\n"));
1174:   PetscFunctionReturn(PETSC_SUCCESS);
1175: }

1177: static PetscErrorCode PetscLogViewWarnDebugging(PetscViewer viewer)
1178: {
1179:   PetscFunctionBegin;
1180:   if (PetscDefined(USE_DEBUG)) {
1181:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1182:     PetscCall(PetscViewerASCIIPrintf(viewer, "      ##########################################################\n"));
1183:     PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1184:     PetscCall(PetscViewerASCIIPrintf(viewer, "      #                       WARNING!!!                       #\n"));
1185:     PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1186:     PetscCall(PetscViewerASCIIPrintf(viewer, "      #   This code was compiled with a debugging option.      #\n"));
1187:     PetscCall(PetscViewerASCIIPrintf(viewer, "      #   To get timing results run ./configure                #\n"));
1188:     PetscCall(PetscViewerASCIIPrintf(viewer, "      #   using --with-debugging=no, the performance will      #\n"));
1189:     PetscCall(PetscViewerASCIIPrintf(viewer, "      #   be generally two or three times faster.              #\n"));
1190:     PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1191:     PetscCall(PetscViewerASCIIPrintf(viewer, "      ##########################################################\n\n\n"));
1192:   }
1193:   PetscFunctionReturn(PETSC_SUCCESS);
1194: }

1196: static PetscErrorCode PetscLogViewWarnNoGpuAwareMpi(PetscViewer viewer)
1197: {
1198: #if defined(PETSC_HAVE_DEVICE)
1199:   PetscMPIInt size;
1200:   PetscBool   deviceInitialized = PETSC_FALSE;

1202:   PetscFunctionBegin;
1203:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
1204:   for (int i = PETSC_DEVICE_HOST + 1; i < PETSC_DEVICE_MAX; ++i) {
1205:     const PetscDeviceType dtype = PetscDeviceTypeCast(i);
1206:     if (PetscDeviceInitialized(dtype)) { /* a non-host device was initialized */
1207:       deviceInitialized = PETSC_TRUE;
1208:       break;
1209:     }
1210:   }
1211:   /* the last condition says petsc is configured with device but it is a pure CPU run, so don't print misleading warnings */
1212:   if (use_gpu_aware_mpi || size == 1 || !deviceInitialized) PetscFunctionReturn(PETSC_SUCCESS);
1213:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1214:   PetscCall(PetscViewerASCIIPrintf(viewer, "      ##########################################################\n"));
1215:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1216:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                       WARNING!!!                       #\n"));
1217:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1218:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   This code was compiled with GPU support and you've   #\n"));
1219:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   created PETSc/GPU objects, but you intentionally     #\n"));
1220:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   used -use_gpu_aware_mpi 0, requiring PETSc to copy   #\n"));
1221:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   additional data between the GPU and CPU. To obtain   #\n"));
1222:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   meaningful timing results on multi-rank runs, use    #\n"));
1223:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   GPU-aware MPI instead.                               #\n"));
1224:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1225:   PetscCall(PetscViewerASCIIPrintf(viewer, "      ##########################################################\n\n\n"));
1226:   PetscFunctionReturn(PETSC_SUCCESS);
1227: #else
1228:   return PETSC_SUCCESS;
1229: #endif
1230: }

1232: static PetscErrorCode PetscLogViewWarnGpuTime(PetscViewer viewer)
1233: {
1234: #if defined(PETSC_HAVE_DEVICE) && !defined(PETSC_HAVE_KOKKOS_WITHOUT_GPU)

1236:   PetscFunctionBegin;
1237:   if (!PetscLogGpuTimeFlag || petsc_gflops == 0) PetscFunctionReturn(PETSC_SUCCESS);
1238:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1239:   PetscCall(PetscViewerASCIIPrintf(viewer, "      ##########################################################\n"));
1240:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1241:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                       WARNING!!!                       #\n"));
1242:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1243:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   This code was run with -log_view_gpu_time            #\n"));
1244:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   This provides accurate timing within the GPU kernels #\n"));
1245:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   but can slow down the entire computation by a        #\n"));
1246:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   measurable amount. For fastest runs we recommend     #\n"));
1247:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #   not using this option.                               #\n"));
1248:   PetscCall(PetscViewerASCIIPrintf(viewer, "      #                                                        #\n"));
1249:   PetscCall(PetscViewerASCIIPrintf(viewer, "      ##########################################################\n\n\n"));
1250:   PetscFunctionReturn(PETSC_SUCCESS);
1251: #else
1252:   return PETSC_SUCCESS;
1253: #endif
1254: }

1256: static PetscErrorCode PetscLogHandlerView_Default_Info(PetscLogHandler handler, PetscViewer viewer)
1257: {
1258:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
1259:   char                    arch[128], hostname[128], username[128], pname[PETSC_MAX_PATH_LEN], date[128];
1260:   PetscLogDouble          locTotalTime, TotalTime, TotalFlops;
1261:   PetscLogDouble          numMessages, messageLength, avgMessLen, numReductions;
1262:   PetscLogDouble          stageTime, flops, flopr, mem, mess, messLen, red;
1263:   PetscLogDouble          fracTime, fracFlops, fracMessages, fracLength, fracReductions, fracMess, fracMessLen, fracRed;
1264:   PetscLogDouble          fracStageTime, fracStageFlops, fracStageMess, fracStageMessLen, fracStageRed;
1265:   PetscLogDouble          min, max, tot, ratio, avg, x, y;
1266:   PetscLogDouble          minf, maxf, totf, ratf, mint, maxt, tott, ratt, ratC, totm, totml, totr, mal, malmax, emalmax;
1267: #if defined(PETSC_HAVE_DEVICE)
1268:   PetscLogEvent  KSP_Solve, SNES_Solve, TS_Step, TAO_Solve; /* These need to be fixed to be some events registered with certain objects */
1269:   PetscLogDouble cct, gct, csz, gsz, gmaxt, gflops, gflopr, fracgflops;
1270: #endif
1271:   PetscMPIInt   minC, maxC;
1272:   PetscMPIInt   size, rank;
1273:   PetscBool    *localStageUsed, *stageUsed;
1274:   PetscBool    *localStageVisible, *stageVisible;
1275:   PetscInt      numStages, numEvents;
1276:   int           stage, oclass;
1277:   PetscLogEvent event;
1278:   char          version[256];
1279:   MPI_Comm      comm;
1280: #if defined(PETSC_HAVE_DEVICE) && !defined(PETSC_HAVE_KOKKOS_WITHOUT_GPU)
1281:   PetscInt64 nas = 0x7FF0000000000002;
1282: #endif
1283:   PetscLogGlobalNames global_stages, global_events;
1284:   PetscEventPerfInfo  zero_info;
1285:   PetscLogState       state;

1287:   PetscFunctionBegin;
1288:   PetscCall(PetscLogHandlerGetState(handler, &state));
1289:   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1290:   PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
1291:   PetscCallMPI(MPI_Comm_size(comm, &size));
1292:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1293:   /* Get the total elapsed time */
1294:   PetscCall(PetscTime(&locTotalTime));
1295:   locTotalTime -= petsc_BaseTime;

1297:   PetscCall(PetscViewerASCIIPrintf(viewer, "****************************************************************************************************************************************************************\n"));
1298:   PetscCall(PetscViewerASCIIPrintf(viewer, "***                                WIDEN YOUR WINDOW TO 160 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document                                 ***\n"));
1299:   PetscCall(PetscViewerASCIIPrintf(viewer, "****************************************************************************************************************************************************************\n"));
1300:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n------------------------------------------------------------------ PETSc Performance Summary: ------------------------------------------------------------------\n\n"));
1301:   PetscCall(PetscLogViewWarnSync(viewer));
1302:   PetscCall(PetscLogViewWarnDebugging(viewer));
1303:   PetscCall(PetscLogViewWarnNoGpuAwareMpi(viewer));
1304:   PetscCall(PetscLogViewWarnGpuTime(viewer));
1305:   PetscCall(PetscGetArchType(arch, sizeof(arch)));
1306:   PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
1307:   PetscCall(PetscGetUserName(username, sizeof(username)));
1308:   PetscCall(PetscGetProgramName(pname, sizeof(pname)));
1309:   PetscCall(PetscGetDate(date, sizeof(date)));
1310:   PetscCall(PetscGetVersion(version, sizeof(version)));

1312: #if defined(PETSC_HAVE_CUPM)
1313:   const char *cupm = PetscDefined(HAVE_CUDA) ? "CUDA" : "HIP";
1314:   if (PetscDeviceCUPMRuntimeArch)
1315:     PetscCall(PetscViewerASCIIPrintf(viewer, "%s on a %s named %s with %d process%s and %s architecture %d, by %s on %s\n", pname, arch, hostname, size, size > 1 ? "es" : "", cupm, PetscDeviceCUPMRuntimeArch, username, date));
1316:   else
1317: #endif
1318:     PetscCall(PetscViewerASCIIPrintf(viewer, "%s on a %s named %s with %d process%s, by %s on %s\n", pname, arch, hostname, size, size > 1 ? "es" : "", username, date));

1320: #if defined(PETSC_HAVE_OPENMP)
1321:   PetscCall(PetscViewerASCIIPrintf(viewer, "Using %" PetscInt_FMT " OpenMP threads\n", PetscNumOMPThreads));
1322: #endif
1323:   PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s\n", version));

1325:   /* Must preserve reduction count before we go on */
1326:   red = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;

1328:   /* Calculate summary information */
1329:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n                         Max       Max/Min     Avg       Total\n"));
1330:   /*   Time */
1331:   PetscCallMPI(MPIU_Allreduce(&locTotalTime, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1332:   PetscCallMPI(MPIU_Allreduce(&locTotalTime, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1333:   PetscCallMPI(MPIU_Allreduce(&locTotalTime, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1334:   avg = tot / ((PetscLogDouble)size);
1335:   if (min != 0.0) ratio = max / min;
1336:   else ratio = 0.0;
1337:   PetscCall(PetscViewerASCIIPrintf(viewer, "Time (sec):           %5.3e   %7.3f   %5.3e\n", max, ratio, avg));
1338:   TotalTime = tot;
1339:   /*   Objects */
1340:   {
1341:     PetscInt num_objects;

1343:     PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
1344:     avg = (PetscLogDouble)num_objects;
1345:   }
1346:   PetscCallMPI(MPIU_Allreduce(&avg, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1347:   PetscCallMPI(MPIU_Allreduce(&avg, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1348:   PetscCallMPI(MPIU_Allreduce(&avg, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1349:   avg = tot / ((PetscLogDouble)size);
1350:   if (min != 0.0) ratio = max / min;
1351:   else ratio = 0.0;
1352:   PetscCall(PetscViewerASCIIPrintf(viewer, "Objects:              %5.3e   %7.3f   %5.3e\n", max, ratio, avg));
1353:   /*   Flops */
1354:   PetscCallMPI(MPIU_Allreduce(&petsc_TotalFlops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1355:   PetscCallMPI(MPIU_Allreduce(&petsc_TotalFlops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1356:   PetscCallMPI(MPIU_Allreduce(&petsc_TotalFlops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1357:   avg = tot / ((PetscLogDouble)size);
1358:   if (min != 0.0) ratio = max / min;
1359:   else ratio = 0.0;
1360:   PetscCall(PetscViewerASCIIPrintf(viewer, "Flops:                %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1361:   TotalFlops = tot;
1362:   /*   Flops/sec -- Must talk to Barry here */
1363:   if (locTotalTime != 0.0) flops = petsc_TotalFlops / locTotalTime;
1364:   else flops = 0.0;
1365:   PetscCallMPI(MPIU_Allreduce(&flops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1366:   PetscCallMPI(MPIU_Allreduce(&flops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1367:   PetscCallMPI(MPIU_Allreduce(&flops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1368:   avg = tot / ((PetscLogDouble)size);
1369:   if (min != 0.0) ratio = max / min;
1370:   else ratio = 0.0;
1371:   PetscCall(PetscViewerASCIIPrintf(viewer, "Flops/sec:            %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1372:   /*   Memory */
1373:   PetscCall(PetscMallocGetMaximumUsage(&mem));
1374:   if (mem > 0.0) {
1375:     PetscCallMPI(MPIU_Allreduce(&mem, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1376:     PetscCallMPI(MPIU_Allreduce(&mem, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1377:     PetscCallMPI(MPIU_Allreduce(&mem, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1378:     avg = tot / ((PetscLogDouble)size);
1379:     if (min != 0.0) ratio = max / min;
1380:     else ratio = 0.0;
1381:     PetscCall(PetscViewerASCIIPrintf(viewer, "Memory (bytes):       %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1382:   }
1383:   /*   Messages */
1384:   mess = 0.5 * (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct);
1385:   PetscCallMPI(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1386:   PetscCallMPI(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1387:   PetscCallMPI(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1388:   avg = tot / ((PetscLogDouble)size);
1389:   if (min != 0.0) ratio = max / min;
1390:   else ratio = 0.0;
1391:   PetscCall(PetscViewerASCIIPrintf(viewer, "MPI Msg Count:        %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1392:   numMessages = tot;
1393:   /*   Message Lengths */
1394:   mess = 0.5 * (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len);
1395:   PetscCallMPI(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1396:   PetscCallMPI(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1397:   PetscCallMPI(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1398:   if (numMessages != 0) avg = tot / numMessages;
1399:   else avg = 0.0;
1400:   if (min != 0.0) ratio = max / min;
1401:   else ratio = 0.0;
1402:   PetscCall(PetscViewerASCIIPrintf(viewer, "MPI Msg Len (bytes):  %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1403:   messageLength = tot;
1404:   /*   Reductions */
1405:   PetscCallMPI(MPIU_Allreduce(&red, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1406:   PetscCallMPI(MPIU_Allreduce(&red, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1407:   PetscCallMPI(MPIU_Allreduce(&red, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1408:   if (min != 0.0) ratio = max / min;
1409:   else ratio = 0.0;
1410:   PetscCall(PetscViewerASCIIPrintf(viewer, "MPI Reductions:       %5.3e   %7.3f\n", max, ratio));
1411:   numReductions = red; /* wrong because uses count from process zero */
1412:   PetscCall(PetscViewerASCIIPrintf(viewer, "\nFlop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)\n"));
1413:   PetscCall(PetscViewerASCIIPrintf(viewer, "                            e.g., VecAXPY() for real vectors of length N --> 2N flops\n"));
1414:   PetscCall(PetscViewerASCIIPrintf(viewer, "                            and VecAXPY() for complex vectors of length N --> 8N flops\n"));

1416:   PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
1417:   PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
1418:   PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
1419:   PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
1420:   PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
1421:   PetscCall(PetscMalloc1(numStages, &localStageUsed));
1422:   PetscCall(PetscMalloc1(numStages, &stageUsed));
1423:   PetscCall(PetscMalloc1(numStages, &localStageVisible));
1424:   PetscCall(PetscMalloc1(numStages, &stageVisible));
1425:   if (numStages > 0) {
1426:     for (stage = 0; stage < numStages; stage++) {
1427:       PetscInt stage_id;

1429:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1430:       if (stage_id >= 0) {
1431:         PetscStagePerf *stage_info;

1433:         PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info));
1434:         localStageUsed[stage]    = stage_info->used;
1435:         localStageVisible[stage] = stage_info->perfInfo.visible;
1436:       } else {
1437:         localStageUsed[stage]    = PETSC_FALSE;
1438:         localStageVisible[stage] = PETSC_TRUE;
1439:       }
1440:     }
1441:     PetscCallMPI(MPIU_Allreduce(localStageUsed, stageUsed, numStages, MPIU_BOOL, MPI_LOR, comm));
1442:     PetscCallMPI(MPIU_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm));
1443:     for (stage = 0; stage < numStages; stage++) {
1444:       if (stageUsed[stage] && stageVisible[stage]) {
1445:         PetscCall(PetscViewerASCIIPrintf(viewer, "\nSummary of Stages:   ----- Time ------  ----- Flop ------  --- Messages ---  -- Message Lengths --  -- Reductions --\n"));
1446:         PetscCall(PetscViewerASCIIPrintf(viewer, "                        Avg     %%Total     Avg     %%Total    Count   %%Total     Avg         %%Total    Count   %%Total\n"));
1447:         break;
1448:       }
1449:     }
1450:     for (stage = 0; stage < numStages; stage++) {
1451:       PetscInt            stage_id;
1452:       PetscEventPerfInfo *stage_info;
1453:       const char         *stage_name;

1455:       if (!(stageUsed[stage] && stageVisible[stage])) continue;
1456:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1457:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1458:       stage_info = &zero_info;
1459:       if (localStageUsed[stage]) {
1460:         PetscStagePerf *stage_perf_info;

1462:         PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_perf_info));
1463:         stage_info = &stage_perf_info->perfInfo;
1464:       }
1465:       PetscCallMPI(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1466:       PetscCallMPI(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1467:       PetscCallMPI(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1468:       PetscCallMPI(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1469:       PetscCallMPI(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1470:       mess *= 0.5;
1471:       messLen *= 0.5;
1472:       red /= size;
1473:       if (TotalTime != 0.0) fracTime = stageTime / TotalTime;
1474:       else fracTime = 0.0;
1475:       if (TotalFlops != 0.0) fracFlops = flops / TotalFlops;
1476:       else fracFlops = 0.0;
1477:       /* Talk to Barry if (stageTime     != 0.0) flops          = (size*flops)/stageTime; else flops          = 0.0; */
1478:       if (numMessages != 0.0) fracMessages = mess / numMessages;
1479:       else fracMessages = 0.0;
1480:       if (mess != 0.0) avgMessLen = messLen / mess;
1481:       else avgMessLen = 0.0;
1482:       if (messageLength != 0.0) fracLength = messLen / messageLength;
1483:       else fracLength = 0.0;
1484:       if (numReductions != 0.0) fracReductions = red / numReductions;
1485:       else fracReductions = 0.0;
1486:       PetscCall(PetscViewerASCIIPrintf(viewer, "%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));
1487:     }
1488:   }

1490:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n------------------------------------------------------------------------------------------------------------------------\n"));
1491:   PetscCall(PetscViewerASCIIPrintf(viewer, "See the 'Profiling' chapter of the users' manual for details on interpreting output.\n"));
1492:   PetscCall(PetscViewerASCIIPrintf(viewer, "Phase summary info:\n"));
1493:   PetscCall(PetscViewerASCIIPrintf(viewer, "   Count: number of times phase was executed\n"));
1494:   PetscCall(PetscViewerASCIIPrintf(viewer, "   Time and Flop: Max - maximum over all processors\n"));
1495:   PetscCall(PetscViewerASCIIPrintf(viewer, "                  Ratio - ratio of maximum to minimum over all processors\n"));
1496:   PetscCall(PetscViewerASCIIPrintf(viewer, "   Mess: number of messages sent\n"));
1497:   PetscCall(PetscViewerASCIIPrintf(viewer, "   AvgLen: average message length (bytes)\n"));
1498:   PetscCall(PetscViewerASCIIPrintf(viewer, "   Reduct: number of global reductions\n"));
1499:   PetscCall(PetscViewerASCIIPrintf(viewer, "   Global: entire computation\n"));
1500:   PetscCall(PetscViewerASCIIPrintf(viewer, "   Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().\n"));
1501:   PetscCall(PetscViewerASCIIPrintf(viewer, "      %%T - percent time in this phase         %%F - percent flop in this phase\n"));
1502:   PetscCall(PetscViewerASCIIPrintf(viewer, "      %%M - percent messages in this phase     %%L - percent message lengths in this phase\n"));
1503:   PetscCall(PetscViewerASCIIPrintf(viewer, "      %%R - percent reductions in this phase\n"));
1504:   PetscCall(PetscViewerASCIIPrintf(viewer, "   Total Mflop/s: 10e-6 * (sum of flop over all processors)/(max time over all processors)\n"));
1505:   if (PetscLogMemory) {
1506:     PetscCall(PetscViewerASCIIPrintf(viewer, "   Memory usage is summed over all MPI processes, it is given in mega-bytes\n"));
1507:     PetscCall(PetscViewerASCIIPrintf(viewer, "   Malloc Mbytes: Memory allocated and kept during event (sum over all calls to event). May be negative\n"));
1508:     PetscCall(PetscViewerASCIIPrintf(viewer, "   EMalloc Mbytes: extra memory allocated during event and then freed (maximum over all calls to events). Never negative\n"));
1509:     PetscCall(PetscViewerASCIIPrintf(viewer, "   MMalloc Mbytes: Increase in high water mark of allocated memory (sum over all calls to event). Never negative\n"));
1510:     PetscCall(PetscViewerASCIIPrintf(viewer, "   RMI Mbytes: Increase in resident memory (sum over all calls to event)\n"));
1511:   }
1512: #if defined(PETSC_HAVE_DEVICE)
1513:   PetscCall(PetscViewerASCIIPrintf(viewer, "   GPU Mflop/s: 10e-6 * (sum of flop on GPU over all processors)/(max GPU time over all processors)\n"));
1514:   PetscCall(PetscViewerASCIIPrintf(viewer, "   CpuToGpu Count: total number of CPU to GPU copies per processor\n"));
1515:   PetscCall(PetscViewerASCIIPrintf(viewer, "   CpuToGpu Size (Mbytes): 10e-6 * (total size of CPU to GPU copies per processor)\n"));
1516:   PetscCall(PetscViewerASCIIPrintf(viewer, "   GpuToCpu Count: total number of GPU to CPU copies per processor\n"));
1517:   PetscCall(PetscViewerASCIIPrintf(viewer, "   GpuToCpu Size (Mbytes): 10e-6 * (total size of GPU to CPU copies per processor)\n"));
1518:   PetscCall(PetscViewerASCIIPrintf(viewer, "   GPU %%F: percent flops on GPU in this event\n"));
1519: #endif
1520:   PetscCall(PetscViewerASCIIPrintf(viewer, "------------------------------------------------------------------------------------------------------------------------\n"));

1522:   PetscCall(PetscLogViewWarnDebugging(viewer));

1524:   /* Report events */
1525:   PetscCall(PetscViewerASCIIPrintf(viewer, "Event                Count      Time (sec)     Flop                              --- Global ---  --- Stage ----  Total"));
1526:   if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, "  Malloc EMalloc MMalloc RMI"));
1527: #if defined(PETSC_HAVE_DEVICE)
1528:   PetscCall(PetscViewerASCIIPrintf(viewer, "   GPU    - CpuToGpu -   - GpuToCpu - GPU"));
1529: #endif
1530:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1531:   PetscCall(PetscViewerASCIIPrintf(viewer, "                   Max Ratio  Max     Ratio   Max  Ratio  Mess   AvgLen  Reduct  %%T %%F %%M %%L %%R  %%T %%F %%M %%L %%R Mflop/s"));
1532:   if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, " Mbytes Mbytes Mbytes Mbytes"));
1533: #if defined(PETSC_HAVE_DEVICE)
1534:   PetscCall(PetscViewerASCIIPrintf(viewer, " Mflop/s Count   Size   Count   Size  %%F"));
1535: #endif
1536:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1537:   PetscCall(PetscViewerASCIIPrintf(viewer, "------------------------------------------------------------------------------------------------------------------------"));
1538:   if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, "-----------------------------"));
1539: #if defined(PETSC_HAVE_DEVICE)
1540:   PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------------"));
1541: #endif
1542:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));

1544: #if defined(PETSC_HAVE_DEVICE)
1545:   /* this indirect way of accessing these values is needed when PETSc is build with multiple libraries since the symbols are not in libpetscsys */
1546:   PetscCall(PetscLogStateGetEventFromName(state, "TaoSolve", &TAO_Solve));
1547:   PetscCall(PetscLogStateGetEventFromName(state, "TSStep", &TS_Step));
1548:   PetscCall(PetscLogStateGetEventFromName(state, "SNESSolve", &SNES_Solve));
1549:   PetscCall(PetscLogStateGetEventFromName(state, "KSPSolve", &KSP_Solve));
1550: #endif

1552:   for (stage = 0; stage < numStages; stage++) {
1553:     PetscInt            stage_id;
1554:     PetscEventPerfInfo *stage_info;
1555:     const char         *stage_name;

1557:     if (!(stageVisible[stage] && stageUsed[stage])) continue;
1558:     PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1559:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1560:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n--- Event Stage %d: %s\n\n", stage, stage_name));
1561:     stage_info = &zero_info;
1562:     if (localStageUsed[stage]) {
1563:       PetscStagePerf *stage_perf_info;

1565:       PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_perf_info));
1566:       stage_info = &stage_perf_info->perfInfo;
1567:     }
1568:     PetscCallMPI(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1569:     PetscCallMPI(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1570:     PetscCallMPI(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1571:     PetscCallMPI(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1572:     PetscCallMPI(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1573:     mess *= 0.5;
1574:     messLen *= 0.5;
1575:     red /= size;

1577:     for (event = 0; event < numEvents; event++) {
1578:       PetscInt            event_id;
1579:       PetscEventPerfInfo *event_info = &zero_info;
1580:       PetscBool           is_zero    = PETSC_FALSE;
1581:       const char         *event_name;

1583:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1584:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1585:       if (event_id >= 0 && stage_id >= 0) { PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &event_info)); }
1586:       PetscCall(PetscMemcmp(event_info, &zero_info, sizeof(zero_info), &is_zero));
1587:       PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1588:       if (!is_zero) {
1589:         flopr = event_info->flops;
1590:         PetscCallMPI(MPIU_Allreduce(&flopr, &minf, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1591:         PetscCallMPI(MPIU_Allreduce(&flopr, &maxf, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1592:         PetscCallMPI(MPIU_Allreduce(&event_info->flops, &totf, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1593:         PetscCallMPI(MPIU_Allreduce(&event_info->time, &mint, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1594:         PetscCallMPI(MPIU_Allreduce(&event_info->time, &maxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1595:         PetscCallMPI(MPIU_Allreduce(&event_info->time, &tott, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1596:         PetscCallMPI(MPIU_Allreduce(&event_info->numMessages, &totm, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1597:         PetscCallMPI(MPIU_Allreduce(&event_info->messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1598:         PetscCallMPI(MPIU_Allreduce(&event_info->numReductions, &totr, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1599:         PetscCallMPI(MPIU_Allreduce(&event_info->count, &minC, 1, MPI_INT, MPI_MIN, comm));
1600:         PetscCallMPI(MPIU_Allreduce(&event_info->count, &maxC, 1, MPI_INT, MPI_MAX, comm));
1601:         if (PetscLogMemory) {
1602:           PetscCallMPI(MPIU_Allreduce(&event_info->memIncrease, &mem, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1603:           PetscCallMPI(MPIU_Allreduce(&event_info->mallocSpace, &mal, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1604:           PetscCallMPI(MPIU_Allreduce(&event_info->mallocIncrease, &malmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1605:           PetscCallMPI(MPIU_Allreduce(&event_info->mallocIncreaseEvent, &emalmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1606:         }
1607: #if defined(PETSC_HAVE_DEVICE)
1608:         PetscCallMPI(MPIU_Allreduce(&event_info->CpuToGpuCount, &cct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1609:         PetscCallMPI(MPIU_Allreduce(&event_info->GpuToCpuCount, &gct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1610:         PetscCallMPI(MPIU_Allreduce(&event_info->CpuToGpuSize, &csz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1611:         PetscCallMPI(MPIU_Allreduce(&event_info->GpuToCpuSize, &gsz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1612:         PetscCallMPI(MPIU_Allreduce(&event_info->GpuFlops, &gflops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1613:         PetscCallMPI(MPIU_Allreduce(&event_info->GpuTime, &gmaxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1614: #endif
1615:         if (mint < 0.0) {
1616:           PetscCall(PetscViewerASCIIPrintf(viewer, "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));
1617:           mint = 0;
1618:         }
1619:         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);
1620: #if defined(PETSC_HAVE_DEVICE) && !defined(PETSC_HAVE_KOKKOS_WITHOUT_GPU)
1621:         /* Put NaN into the time for all events that may not be time accurately since they may happen asynchronously on the GPU */
1622:         if (!PetscLogGpuTimeFlag && petsc_gflops > 0) {
1623:           memcpy(&gmaxt, &nas, sizeof(PetscLogDouble));
1624:           if (event_id != SNES_Solve && event_id != KSP_Solve && event_id != TS_Step && event_id != TAO_Solve) {
1625:             memcpy(&mint, &nas, sizeof(PetscLogDouble));
1626:             memcpy(&maxt, &nas, sizeof(PetscLogDouble));
1627:           }
1628:         }
1629: #endif
1630:         totm *= 0.5;
1631:         totml *= 0.5;
1632:         totr /= size;

1634:         if (maxC != 0) {
1635:           if (minC != 0) ratC = ((PetscLogDouble)maxC) / minC;
1636:           else ratC = 0.0;
1637:           if (mint != 0.0) ratt = maxt / mint;
1638:           else ratt = 0.0;
1639:           if (minf != 0.0) ratf = maxf / minf;
1640:           else ratf = 0.0;
1641:           if (TotalTime != 0.0) fracTime = tott / TotalTime;
1642:           else fracTime = 0.0;
1643:           if (TotalFlops != 0.0) fracFlops = totf / TotalFlops;
1644:           else fracFlops = 0.0;
1645:           if (stageTime != 0.0) fracStageTime = tott / stageTime;
1646:           else fracStageTime = 0.0;
1647:           if (flops != 0.0) fracStageFlops = totf / flops;
1648:           else fracStageFlops = 0.0;
1649:           if (numMessages != 0.0) fracMess = totm / numMessages;
1650:           else fracMess = 0.0;
1651:           if (messageLength != 0.0) fracMessLen = totml / messageLength;
1652:           else fracMessLen = 0.0;
1653:           if (numReductions != 0.0) fracRed = totr / numReductions;
1654:           else fracRed = 0.0;
1655:           if (mess != 0.0) fracStageMess = totm / mess;
1656:           else fracStageMess = 0.0;
1657:           if (messLen != 0.0) fracStageMessLen = totml / messLen;
1658:           else fracStageMessLen = 0.0;
1659:           if (red != 0.0) fracStageRed = totr / red;
1660:           else fracStageRed = 0.0;
1661:           if (totm != 0.0) totml /= totm;
1662:           else totml = 0.0;
1663:           if (maxt != 0.0) flopr = totf / maxt;
1664:           else flopr = 0.0;
1665:           if (fracStageTime > 1.0 || fracStageFlops > 1.0 || fracStageMess > 1.0 || fracStageMessLen > 1.0 || fracStageRed > 1.0) {
1666:             if (PetscIsNanReal(maxt))
1667:               PetscCall(PetscViewerASCIIPrintf(viewer, "%-16s %7d %3.1f  n/a     n/a   %3.2e %3.1f %2.1e %2.1e %2.1e %2.0f %2.0f %2.0f %2.0f %2.0f Multiple stages n/a", event_name, maxC, ratC, maxf, ratf, totm, totml, totr, 100.0 * fracTime, 100.0 * fracFlops, 100.0 * fracMess, 100.0 * fracMessLen, 100.0 * fracRed));
1668:             else
1669:               PetscCall(PetscViewerASCIIPrintf(viewer, "%-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));
1670:           } else {
1671:             if (PetscIsNanReal(maxt)) { // when maxt, ratt, flopr are NaN (i.e., run with GPUs but without -log_view_gpu_time), replace the confusing "nan" with "n/a"
1672:               PetscCall(PetscViewerASCIIPrintf(viewer, "%-16s %7d %3.1f  n/a     n/a   %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  n/a", event_name, maxC, ratC, 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));
1673:             } else {
1674:               PetscCall(PetscViewerASCIIPrintf(viewer, "%-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));
1675:             }
1676:           }
1677:           if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, " %5.0f   %5.0f   %5.0f   %5.0f", mal / 1.0e6, emalmax / 1.0e6, malmax / 1.0e6, mem / 1.0e6));
1678: #if defined(PETSC_HAVE_DEVICE)
1679:           if (totf != 0.0) fracgflops = gflops / totf;
1680:           else fracgflops = 0.0;
1681:           if (gmaxt != 0.0) gflopr = gflops / gmaxt;
1682:           else gflopr = 0.0;
1683:           if (PetscIsNanReal(gflopr)) {
1684:             PetscCall(PetscViewerASCIIPrintf(viewer, "    n/a    %4.0f %3.2e %4.0f %3.2e % 2.0f", cct / size, csz / (1.0e6 * size), gct / size, gsz / (1.0e6 * size), 100.0 * fracgflops));
1685:           } else {
1686:             PetscCall(PetscViewerASCIIPrintf(viewer, "   %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));
1687:           }
1688: #endif
1689:           PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1690:         }
1691:       }
1692:     }
1693:   }

1695:   /* Memory usage and object creation */
1696:   PetscCall(PetscViewerASCIIPrintf(viewer, "------------------------------------------------------------------------------------------------------------------------"));
1697:   if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, "-----------------------------"));
1698: #if defined(PETSC_HAVE_DEVICE)
1699:   PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------------"));
1700: #endif
1701:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1702:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));

1704:   /* Right now, only stages on the first processor are reported here, meaning only objects associated with
1705:      the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then
1706:      stats for stages local to processor sets.
1707:   */
1708:   /* We should figure out the longest object name here (now 20 characters) */
1709:   PetscCall(PetscViewerASCIIPrintf(viewer, "Object Type          Creations   Destructions. Reports information only for process 0.\n"));
1710:   for (stage = 0; stage < numStages; stage++) {
1711:     const char *stage_name;

1713:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1714:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n--- Event Stage %d: %s\n\n", stage, stage_name));
1715:     if (localStageUsed[stage]) {
1716:       PetscInt num_classes;

1718:       PetscCall(PetscLogStateGetNumClasses(state, &num_classes));
1719:       for (oclass = 0; oclass < num_classes; oclass++) {
1720:         PetscClassPerf *class_perf_info;

1722:         PetscCall(PetscLogHandlerDefaultGetClassPerf(handler, stage, oclass, &class_perf_info));
1723:         if ((class_perf_info->creations > 0) || (class_perf_info->destructions > 0)) {
1724:           PetscLogClassInfo class_reg_info;
1725:           PetscBool         flg;

1727:           PetscCall(PetscLogStateClassGetInfo(state, oclass, &class_reg_info));
1728:           if (stage == 0 && oclass == num_classes - 1) {
1729:             PetscCall(PetscStrcmp(class_reg_info.name, "Viewer", &flg));
1730:             PetscCheck(flg && class_perf_info->creations == 1 && class_perf_info->destructions == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "The last PetscObject type of the main PetscLogStage should be PetscViewer with a single creation and no destruction");
1731:           } else PetscCall(PetscViewerASCIIPrintf(viewer, "%20s %5d          %5d\n", class_reg_info.name, class_perf_info->creations, class_perf_info->destructions));
1732:         }
1733:       }
1734:     }
1735:   }

1737:   PetscCall(PetscFree(localStageUsed));
1738:   PetscCall(PetscFree(stageUsed));
1739:   PetscCall(PetscFree(localStageVisible));
1740:   PetscCall(PetscFree(stageVisible));
1741:   PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1742:   PetscCall(PetscLogGlobalNamesDestroy(&global_events));

1744:   /* Information unrelated to this particular run */
1745:   PetscCall(PetscViewerASCIIPrintf(viewer, "========================================================================================================================\n"));
1746:   PetscCall(PetscTime(&y));
1747:   PetscCall(PetscTime(&x));
1748:   PetscCall(PetscTime(&y));
1749:   PetscCall(PetscTime(&y));
1750:   PetscCall(PetscTime(&y));
1751:   PetscCall(PetscTime(&y));
1752:   PetscCall(PetscTime(&y));
1753:   PetscCall(PetscTime(&y));
1754:   PetscCall(PetscTime(&y));
1755:   PetscCall(PetscTime(&y));
1756:   PetscCall(PetscTime(&y));
1757:   PetscCall(PetscTime(&y));
1758:   PetscCall(PetscViewerASCIIPrintf(viewer, "Average time to get PetscTime(): %g\n", (y - x) / 10.0));
1759:   /* MPI information */
1760:   if (size > 1) {
1761:     MPI_Status  status;
1762:     PetscMPIInt tag;
1763:     MPI_Comm    newcomm;

1765:     PetscCallMPI(MPI_Barrier(comm));
1766:     PetscCall(PetscTime(&x));
1767:     PetscCallMPI(MPI_Barrier(comm));
1768:     PetscCallMPI(MPI_Barrier(comm));
1769:     PetscCallMPI(MPI_Barrier(comm));
1770:     PetscCallMPI(MPI_Barrier(comm));
1771:     PetscCallMPI(MPI_Barrier(comm));
1772:     PetscCall(PetscTime(&y));
1773:     PetscCall(PetscViewerASCIIPrintf(viewer, "Average time for MPI_Barrier(): %g\n", (y - x) / 5.0));
1774:     PetscCall(PetscCommDuplicate(comm, &newcomm, &tag));
1775:     PetscCallMPI(MPI_Barrier(comm));
1776:     if (rank) {
1777:       PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, rank - 1, tag, newcomm, &status));
1778:       PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, (rank + 1) % size, tag, newcomm));
1779:     } else {
1780:       PetscCall(PetscTime(&x));
1781:       PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, 1, tag, newcomm));
1782:       PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, size - 1, tag, newcomm, &status));
1783:       PetscCall(PetscTime(&y));
1784:       PetscCall(PetscViewerASCIIPrintf(viewer, "Average time for zero size MPI_Send(): %g\n", (y - x) / size));
1785:     }
1786:     PetscCall(PetscCommDestroy(&newcomm));
1787:   }
1788:   PetscCall(PetscOptionsView(NULL, viewer));

1790:   /* Machine and compile information */
1791:   if (PetscDefined(USE_FORTRAN_KERNELS)) {
1792:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with FORTRAN kernels\n"));
1793:   } else {
1794:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled without FORTRAN kernels\n"));
1795:   }
1796:   if (PetscDefined(USE_64BIT_INDICES)) {
1797:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 64-bit PetscInt\n"));
1798:   } else if (PetscDefined(USE___FLOAT128)) {
1799:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 32-bit PetscInt\n"));
1800:   }
1801:   if (PetscDefined(USE_REAL_SINGLE)) {
1802:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with single precision PetscScalar and PetscReal\n"));
1803:   } else if (PetscDefined(USE___FLOAT128)) {
1804:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 128 bit precision PetscScalar and PetscReal\n"));
1805:   }
1806:   if (PetscDefined(USE_REAL_MAT_SINGLE)) {
1807:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with single precision matrices\n"));
1808:   } else {
1809:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with full precision matrices (default)\n"));
1810:   }
1811:   PetscCall(PetscViewerASCIIPrintf(viewer, "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)));

1813:   PetscCall(PetscViewerASCIIPrintf(viewer, "Configure options: %s", petscconfigureoptions));
1814:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petscmachineinfo));
1815:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsccompilerinfo));
1816:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsccompilerflagsinfo));
1817:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsclinkerinfo));

1819:   /* Cleanup */
1820:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1821:   PetscCall(PetscLogViewWarnNoGpuAwareMpi(viewer));
1822:   PetscCall(PetscLogViewWarnDebugging(viewer));
1823:   PetscCall(PetscFPTrapPop());
1824:   PetscFunctionReturn(PETSC_SUCCESS);
1825: }

1827: static PetscErrorCode PetscLogHandlerView_Default(PetscLogHandler handler, PetscViewer viewer)
1828: {
1829:   PetscViewerFormat format;

1831:   PetscFunctionBegin;
1832:   PetscCall(PetscViewerGetFormat(viewer, &format));
1833:   if (format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_ASCII_INFO) {
1834:     PetscCall(PetscLogHandlerView_Default_Info(handler, viewer));
1835:   } else if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1836:     PetscCall(PetscLogHandlerView_Default_Detailed(handler, viewer));
1837:   } else if (format == PETSC_VIEWER_ASCII_CSV) {
1838:     PetscCall(PetscLogHandlerView_Default_CSV(handler, viewer));
1839:   }
1840:   PetscFunctionReturn(PETSC_SUCCESS);
1841: }

1843: /*MC
1844:   PETSCLOGHANDLERDEFAULT - PETSCLOGHANDLERDEFAULT = "default" -  A `PetscLogHandler` that collects data for PETSc
1845:   default profiling log viewers (`PetscLogView()` and `PetscLogDump()`).  A log handler of this type is
1846:   created and started (`PetscLogHandlerStart()`) by `PetscLogDefaultBegin()`.

1848:   Options Database Keys:
1849: + -log_include_actions - include a growing list of actions (event beginnings and endings, object creations and destructions) in `PetscLogDump()` (`PetscLogActions()`).
1850: - -log_include_objects - include a growing list of object creations and destructions in `PetscLogDump()` (`PetscLogObjects()`).

1852:   Level: developer

1854: .seealso: [](ch_profiling), `PetscLogHandler`
1855: M*/

1857: PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Default(PetscLogHandler handler)
1858: {
1859:   PetscFunctionBegin;
1860:   PetscCall(PetscLogHandlerContextCreate_Default((PetscLogHandler_Default *)&handler->data));
1861:   handler->ops->destroy       = PetscLogHandlerDestroy_Default;
1862:   handler->ops->eventbegin    = PetscLogHandlerEventBegin_Default;
1863:   handler->ops->eventend      = PetscLogHandlerEventEnd_Default;
1864:   handler->ops->eventsync     = PetscLogHandlerEventSync_Default;
1865:   handler->ops->objectcreate  = PetscLogHandlerObjectCreate_Default;
1866:   handler->ops->objectdestroy = PetscLogHandlerObjectDestroy_Default;
1867:   handler->ops->stagepush     = PetscLogHandlerStagePush_Default;
1868:   handler->ops->stagepop      = PetscLogHandlerStagePop_Default;
1869:   handler->ops->view          = PetscLogHandlerView_Default;
1870:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetEventPerfInfo_C", PetscLogHandlerGetEventPerfInfo_Default));
1871:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetStagePerfInfo_C", PetscLogHandlerGetStagePerfInfo_Default));
1872:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerSetLogActions_C", PetscLogHandlerSetLogActions_Default));
1873:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerSetLogObjects_C", PetscLogHandlerSetLogObjects_Default));
1874:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerLogObjectState_C", PetscLogHandlerLogObjectState_Default));
1875:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetNumObjects_C", PetscLogHandlerGetNumObjects_Default));
1876:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventDeactivatePush_C", PetscLogHandlerEventDeactivatePush_Default));
1877:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventDeactivatePop_C", PetscLogHandlerEventDeactivatePop_Default));
1878:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventsPause_C", PetscLogHandlerEventsPause_Default));
1879:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventsResume_C", PetscLogHandlerEventsResume_Default));
1880:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerDump_C", PetscLogHandlerDump_Default));
1881:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerStageSetVisible_C", PetscLogHandlerStageSetVisible_Default));
1882:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerStageGetVisible_C", PetscLogHandlerStageGetVisible_Default));
1883:   PetscFunctionReturn(PETSC_SUCCESS);
1884: }