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(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, 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, 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:   /* Pause stages in reverse of the order they were pushed */
703:   for (PetscInt stage = num_stages - 1; stage >= 0; stage--) {
704:     PetscStagePerf *stage_info = NULL;
705:     PetscInt        num_events;

707:     PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
708:     PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL));
709:     /* Pause events in reverse of the order they were pushed */
710:     for (PetscInt event = num_events - 1; event >= 0; event--) {
711:       PetscEventPerfInfo *event_info = NULL;
712:       PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info));
713:       if (event_info->depth > 0) {
714:         event_info->depth *= -1;
715:         PetscCall(PetscEventPerfInfoPause(event_info, time, PetscLogMemory, event));
716:       }
717:     }
718:     if (stage > 0 && stage_info->perfInfo.depth > 0) {
719:       stage_info->perfInfo.depth *= -1;
720:       PetscCall(PetscEventPerfInfoPause(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2)));
721:     }
722:   }
723:   PetscFunctionReturn(PETSC_SUCCESS);
724: }

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

732:   PetscFunctionBegin;
733:   if (--def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
734:   PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL));
735:   PetscCall(PetscTime(&time));
736:   /* Unpause stages in the same order they were pushed */
737:   for (PetscInt stage = 0; stage < num_stages; stage++) {
738:     PetscStagePerf *stage_info = NULL;
739:     PetscInt        num_events;

741:     PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
742:     PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL));
743:     /* Unpause events in the same order they were pushed */
744:     for (PetscInt event = 0; event < num_events; event++) {
745:       PetscEventPerfInfo *event_info = NULL;
746:       PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info));
747:       if (event_info->depth < 0) {
748:         event_info->depth *= -1;
749:         PetscCall(PetscEventPerfInfoResume(event_info, time, PetscLogMemory, event));
750:       }
751:     }
752:     if (stage > 0 && stage_info->perfInfo.depth < 0) {
753:       stage_info->perfInfo.depth *= -1;
754:       PetscCall(PetscEventPerfInfoResume(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2)));
755:     }
756:   }
757:   PetscFunctionReturn(PETSC_SUCCESS);
758: }

760: static PetscErrorCode PetscLogHandlerStagePush_Default(PetscLogHandler h, PetscLogStage new_stage)
761: {
762:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
763:   PetscLogDouble          time;
764:   PetscLogState           state;
765:   PetscLogStage           current_stage;
766:   PetscStagePerf         *new_stage_info;

768:   PetscFunctionBegin;
769:   if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
770:   PetscCall(PetscLogHandlerGetState(h, &state));
771:   current_stage = state->current_stage;
772:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, new_stage, &new_stage_info));
773:   PetscCall(PetscTime(&time));

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

791: static PetscErrorCode PetscLogHandlerStagePop_Default(PetscLogHandler h, PetscLogStage old_stage)
792: {
793:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
794:   PetscLogStage           current_stage;
795:   PetscStagePerf         *old_stage_info;
796:   PetscLogState           state;
797:   PetscLogDouble          time;

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

817: static PetscErrorCode PetscLogHandlerStageSetVisible_Default(PetscLogHandler h, PetscLogStage stage, PetscBool is_visible)
818: {
819:   PetscStagePerf *stage_info;

821:   PetscFunctionBegin;
822:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
823:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info));
824:   stage_info->perfInfo.visible = is_visible;
825:   PetscFunctionReturn(PETSC_SUCCESS);
826: }

828: static PetscErrorCode PetscLogHandlerStageGetVisible_Default(PetscLogHandler h, PetscLogStage stage, PetscBool *is_visible)
829: {
830:   PetscStagePerf *stage_info;

832:   PetscFunctionBegin;
833:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
834:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info));
835:   *is_visible = stage_info->perfInfo.visible;
836:   PetscFunctionReturn(PETSC_SUCCESS);
837: }

839: static PetscErrorCode PetscLogHandlerSetLogActions_Default(PetscLogHandler handler, PetscBool flag)
840: {
841:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;

843:   PetscFunctionBegin;
844:   def->petsc_logActions = flag;
845:   PetscFunctionReturn(PETSC_SUCCESS);
846: }

848: static PetscErrorCode PetscLogHandlerSetLogObjects_Default(PetscLogHandler handler, PetscBool flag)
849: {
850:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;

852:   PetscFunctionBegin;
853:   def->petsc_logObjects = flag;
854:   PetscFunctionReturn(PETSC_SUCCESS);
855: }

857: static PetscErrorCode PetscLogHandlerLogObjectState_Default(PetscLogHandler handler, PetscObject obj, const char format[], va_list Argp)
858: {
859:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
860:   size_t                  fullLength;

862:   PetscFunctionBegin;
863:   if (def->petsc_logObjects) {
864:     Object  *obj_entry = NULL;
865:     PetscInt objid;

867:     PetscCall(PetscIntCast(obj->id - 1, &objid));
868:     PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, objid, &obj_entry));
869:     PetscCall(PetscVSNPrintf(obj_entry->info, 64, format, &fullLength, Argp));
870:   }
871:   PetscFunctionReturn(PETSC_SUCCESS);
872: }

874: static PetscErrorCode PetscLogHandlerGetNumObjects_Default(PetscLogHandler handler, PetscInt *num_objects)
875: {
876:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;

878:   PetscFunctionBegin;
879:   PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, num_objects, NULL));
880:   PetscFunctionReturn(PETSC_SUCCESS);
881: }

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

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

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

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

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

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

962: /*
963:   PetscLogView_Detailed - Each process prints the times for its own events

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

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

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

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

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

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

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

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

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

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

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

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

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

1200: static PetscErrorCode PetscLogViewWarnNoGpuAwareMpi(PetscViewer viewer)
1201: {
1202: #if defined(PETSC_HAVE_DEVICE)
1203:   PetscMPIInt size;
1204:   PetscBool   deviceInitialized = PETSC_FALSE;

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

1236: static PetscErrorCode PetscLogViewWarnGpuTime(PetscViewer viewer)
1237: {
1238: #if defined(PETSC_HAVE_DEVICE) && !defined(PETSC_HAVE_KOKKOS_WITHOUT_GPU)

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

1260: PETSC_INTERN int    PetscGlobalArgc;
1261: PETSC_INTERN char **PetscGlobalArgs;

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

1294:   PetscFunctionBegin;
1295:   PetscCall(PetscLogHandlerGetState(handler, &state));
1296:   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1297:   PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
1298:   PetscCallMPI(MPI_Comm_size(comm, &size));
1299:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1300:   /* Get the total elapsed time */
1301:   PetscCall(PetscTime(&locTotalTime));
1302:   locTotalTime -= petsc_BaseTime;

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

1319: #if defined(PETSC_HAVE_CUPM)
1320:   const char *cupm = PetscDefined(HAVE_CUDA) ? "CUDA" : "HIP";
1321:   if (PetscDeviceCUPMRuntimeArch)
1322:     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));
1323:   else
1324: #endif
1325:     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));

1327: #if defined(PETSC_HAVE_OPENMP)
1328:   PetscCall(PetscViewerASCIIPrintf(viewer, "Using %" PetscInt_FMT " OpenMP threads\n", PetscNumOMPThreads));
1329: #endif
1330:   PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s\n", version));

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

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

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

1423:   PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
1424:   PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
1425:   PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
1426:   PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
1427:   PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
1428:   PetscCall(PetscMalloc1(numStages, &localStageUsed));
1429:   PetscCall(PetscMalloc1(numStages, &stageUsed));
1430:   PetscCall(PetscMalloc1(numStages, &localStageVisible));
1431:   PetscCall(PetscMalloc1(numStages, &stageVisible));
1432:   if (numStages > 0) {
1433:     for (stage = 0; stage < numStages; stage++) {
1434:       PetscInt stage_id;

1436:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1437:       if (stage_id >= 0) {
1438:         PetscStagePerf *stage_info;

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

1462:       if (!(stageUsed[stage] && stageVisible[stage])) continue;
1463:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1464:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1465:       stage_info = &zero_info;
1466:       if (localStageUsed[stage]) {
1467:         PetscStagePerf *stage_perf_info;

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

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

1529:   PetscCall(PetscLogViewWarnDebugging(viewer));

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

1551: #if defined(PETSC_HAVE_DEVICE)
1552:   /* this indirect way of accessing these values is needed when PETSc is build with multiple libraries since the symbols are not in libpetscsys */
1553:   PetscCall(PetscLogStateGetEventFromName(state, "TaoSolve", &TAO_Solve));
1554:   PetscCall(PetscLogStateGetEventFromName(state, "TSStep", &TS_Step));
1555:   PetscCall(PetscLogStateGetEventFromName(state, "SNESSolve", &SNES_Solve));
1556:   PetscCall(PetscLogStateGetEventFromName(state, "KSPSolve", &KSP_Solve));
1557: #endif

1559:   for (stage = 0; stage < numStages; stage++) {
1560:     PetscInt            stage_id;
1561:     PetscEventPerfInfo *stage_info;
1562:     const char         *stage_name;

1564:     if (!(stageVisible[stage] && stageUsed[stage])) continue;
1565:     PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1566:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1567:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n--- Event Stage %d: %s\n\n", stage, stage_name));
1568:     stage_info = &zero_info;
1569:     if (localStageUsed[stage]) {
1570:       PetscStagePerf *stage_perf_info;

1572:       PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_perf_info));
1573:       stage_info = &stage_perf_info->perfInfo;
1574:     }
1575:     PetscCallMPI(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1576:     PetscCallMPI(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1577:     PetscCallMPI(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1578:     PetscCallMPI(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1579:     PetscCallMPI(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1580:     mess *= 0.5;
1581:     messLen *= 0.5;
1582:     red /= size;

1584:     for (event = 0; event < numEvents; event++) {
1585:       PetscInt            event_id;
1586:       PetscEventPerfInfo *event_info = &zero_info;
1587:       PetscBool           is_zero    = PETSC_FALSE;
1588:       const char         *event_name;

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

1641:         if (maxC != 0) {
1642:           if (minC != 0) ratC = ((PetscLogDouble)maxC) / minC;
1643:           else ratC = 0.0;
1644:           if (mint != 0.0) ratt = maxt / mint;
1645:           else ratt = 0.0;
1646:           if (minf != 0.0) ratf = maxf / minf;
1647:           else ratf = 0.0;
1648:           if (TotalTime != 0.0) fracTime = tott / TotalTime;
1649:           else fracTime = 0.0;
1650:           if (TotalFlops != 0.0) fracFlops = totf / TotalFlops;
1651:           else fracFlops = 0.0;
1652:           if (stageTime != 0.0) fracStageTime = tott / stageTime;
1653:           else fracStageTime = 0.0;
1654:           if (flops != 0.0) fracStageFlops = totf / flops;
1655:           else fracStageFlops = 0.0;
1656:           if (numMessages != 0.0) fracMess = totm / numMessages;
1657:           else fracMess = 0.0;
1658:           if (messageLength != 0.0) fracMessLen = totml / messageLength;
1659:           else fracMessLen = 0.0;
1660:           if (numReductions != 0.0) fracRed = totr / numReductions;
1661:           else fracRed = 0.0;
1662:           if (mess != 0.0) fracStageMess = totm / mess;
1663:           else fracStageMess = 0.0;
1664:           if (messLen != 0.0) fracStageMessLen = totml / messLen;
1665:           else fracStageMessLen = 0.0;
1666:           if (red != 0.0) fracStageRed = totr / red;
1667:           else fracStageRed = 0.0;
1668:           if (totm != 0.0) totml /= totm;
1669:           else totml = 0.0;
1670:           if (maxt != 0.0) flopr = totf / maxt;
1671:           else flopr = 0.0;
1672:           if (fracStageTime > 1.0 || fracStageFlops > 1.0 || fracStageMess > 1.0 || fracStageMessLen > 1.0 || fracStageRed > 1.0) {
1673:             if (PetscIsNanReal(maxt))
1674:               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));
1675:             else
1676:               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));
1677:           } else {
1678:             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"
1679:               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));
1680:             } else {
1681:               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));
1682:             }
1683:           }
1684:           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));
1685: #if defined(PETSC_HAVE_DEVICE)
1686:           if (totf != 0.0) fracgflops = gflops / totf;
1687:           else fracgflops = 0.0;
1688:           if (gmaxt != 0.0) gflopr = gflops / gmaxt;
1689:           else gflopr = 0.0;
1690:           if (PetscIsNanReal(gflopr)) {
1691:             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));
1692:           } else {
1693:             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));
1694:           }
1695: #endif
1696:           PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1697:         }
1698:       }
1699:     }
1700:   }

1702:   /* Memory usage and object creation */
1703:   PetscCall(PetscViewerASCIIPrintf(viewer, "------------------------------------------------------------------------------------------------------------------------"));
1704:   if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, "-----------------------------"));
1705: #if defined(PETSC_HAVE_DEVICE)
1706:   PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------------"));
1707: #endif
1708:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1709:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));

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

1720:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1721:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n--- Event Stage %d: %s\n\n", stage, stage_name));
1722:     if (localStageUsed[stage]) {
1723:       PetscInt num_classes;

1725:       PetscCall(PetscLogStateGetNumClasses(state, &num_classes));
1726:       for (oclass = 0; oclass < num_classes; oclass++) {
1727:         PetscClassPerf *class_perf_info;

1729:         PetscCall(PetscLogHandlerDefaultGetClassPerf(handler, stage, oclass, &class_perf_info));
1730:         if (class_perf_info->creations > 0 || class_perf_info->destructions > 0) {
1731:           PetscLogClassInfo class_reg_info;
1732:           PetscBool         flg = PETSC_FALSE;

1734:           PetscCall(PetscLogStateClassGetInfo(state, oclass, &class_reg_info));
1735:           if (stage == 0 && oclass == num_classes - 1) {
1736:             if (PetscGlobalArgc == 0 && PetscGlobalArgs == NULL) {
1737:               PetscCall(PetscStrcmp(class_reg_info.name, "Viewer", &flg));
1738:               if (flg && class_perf_info->creations == PetscLogNumViewersCreated && class_perf_info->destructions == PetscLogNumViewersDestroyed) continue;
1739:             }
1740:           }
1741:           PetscCall(PetscViewerASCIIPrintf(viewer, "%20s %5d          %5d\n", class_reg_info.name, class_perf_info->creations, class_perf_info->destructions));
1742:         }
1743:       }
1744:     }
1745:   }

1747:   PetscCall(PetscFree(localStageUsed));
1748:   PetscCall(PetscFree(stageUsed));
1749:   PetscCall(PetscFree(localStageVisible));
1750:   PetscCall(PetscFree(stageVisible));
1751:   PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1752:   PetscCall(PetscLogGlobalNamesDestroy(&global_events));

1754:   /* Information unrelated to this particular run */
1755:   PetscCall(PetscViewerASCIIPrintf(viewer, "========================================================================================================================\n"));
1756:   PetscCall(PetscTime(&y));
1757:   PetscCall(PetscTime(&x));
1758:   PetscCall(PetscTime(&y));
1759:   PetscCall(PetscTime(&y));
1760:   PetscCall(PetscTime(&y));
1761:   PetscCall(PetscTime(&y));
1762:   PetscCall(PetscTime(&y));
1763:   PetscCall(PetscTime(&y));
1764:   PetscCall(PetscTime(&y));
1765:   PetscCall(PetscTime(&y));
1766:   PetscCall(PetscTime(&y));
1767:   PetscCall(PetscTime(&y));
1768:   PetscCall(PetscViewerASCIIPrintf(viewer, "Average time to get PetscTime(): %g\n", (y - x) / 10.0));
1769:   /* MPI information */
1770:   if (size > 1) {
1771:     MPI_Status  status;
1772:     PetscMPIInt tag;
1773:     MPI_Comm    newcomm;

1775:     PetscCallMPI(MPI_Barrier(comm));
1776:     PetscCall(PetscTime(&x));
1777:     PetscCallMPI(MPI_Barrier(comm));
1778:     PetscCallMPI(MPI_Barrier(comm));
1779:     PetscCallMPI(MPI_Barrier(comm));
1780:     PetscCallMPI(MPI_Barrier(comm));
1781:     PetscCallMPI(MPI_Barrier(comm));
1782:     PetscCall(PetscTime(&y));
1783:     PetscCall(PetscViewerASCIIPrintf(viewer, "Average time for MPI_Barrier(): %g\n", (y - x) / 5.0));
1784:     PetscCall(PetscCommDuplicate(comm, &newcomm, &tag));
1785:     PetscCallMPI(MPI_Barrier(comm));
1786:     if (rank) {
1787:       PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, rank - 1, tag, newcomm, &status));
1788:       PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, (rank + 1) % size, tag, newcomm));
1789:     } else {
1790:       PetscCall(PetscTime(&x));
1791:       PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, 1, tag, newcomm));
1792:       PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, size - 1, tag, newcomm, &status));
1793:       PetscCall(PetscTime(&y));
1794:       PetscCall(PetscViewerASCIIPrintf(viewer, "Average time for zero size MPI_Send(): %g\n", (y - x) / size));
1795:     }
1796:     PetscCall(PetscCommDestroy(&newcomm));
1797:   }
1798:   PetscCall(PetscOptionsView(NULL, viewer));

1800:   /* Machine and compile information */
1801:   if (PetscDefined(USE_FORTRAN_KERNELS)) {
1802:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with FORTRAN kernels\n"));
1803:   } else {
1804:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled without FORTRAN kernels\n"));
1805:   }
1806:   if (PetscDefined(USE_64BIT_INDICES)) {
1807:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 64-bit PetscInt\n"));
1808:   } else if (PetscDefined(USE___FLOAT128)) {
1809:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 32-bit PetscInt\n"));
1810:   }
1811:   if (PetscDefined(USE_REAL_SINGLE)) {
1812:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with single precision PetscScalar and PetscReal\n"));
1813:   } else if (PetscDefined(USE___FLOAT128)) {
1814:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 128 bit precision PetscScalar and PetscReal\n"));
1815:   }
1816:   if (PetscDefined(USE_REAL_MAT_SINGLE)) {
1817:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with single precision matrices\n"));
1818:   } else {
1819:     PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with full precision matrices (default)\n"));
1820:   }
1821:   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)));

1823:   PetscCall(PetscViewerASCIIPrintf(viewer, "Configure options: %s", petscconfigureoptions));
1824:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petscmachineinfo));
1825:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsccompilerinfo));
1826:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsccompilerflagsinfo));
1827:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsclinkerinfo));

1829:   /* Cleanup */
1830:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1831:   PetscCall(PetscLogViewWarnNoGpuAwareMpi(viewer));
1832:   PetscCall(PetscLogViewWarnDebugging(viewer));
1833:   PetscCall(PetscFPTrapPop());
1834:   PetscFunctionReturn(PETSC_SUCCESS);
1835: }

1837: static PetscErrorCode PetscLogHandlerView_Default(PetscLogHandler handler, PetscViewer viewer)
1838: {
1839:   PetscViewerFormat format;

1841:   PetscFunctionBegin;
1842:   PetscCall(PetscViewerGetFormat(viewer, &format));
1843:   if (format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_ASCII_INFO) {
1844:     PetscCall(PetscLogHandlerView_Default_Info(handler, viewer));
1845:   } else if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1846:     PetscCall(PetscLogHandlerView_Default_Detailed(handler, viewer));
1847:   } else if (format == PETSC_VIEWER_ASCII_CSV) {
1848:     PetscCall(PetscLogHandlerView_Default_CSV(handler, viewer));
1849:   }
1850:   PetscFunctionReturn(PETSC_SUCCESS);
1851: }

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

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

1862:   Level: developer

1864: .seealso: [](ch_profiling), `PetscLogHandler`
1865: M*/

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