Actual source code: logdefault.c

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

 10: static PetscErrorCode PetscEventPerfInfoInit(PetscEventPerfInfo *eventInfo)
 11: {
 12:   PetscFunctionBegin;
 13:   PetscCall(PetscMemzero(eventInfo, sizeof(*eventInfo)));
 14:   eventInfo->visible   = PETSC_TRUE;
 15:   eventInfo->id        = -1;
 16:   eventInfo->dof[0]    = -1.0;
 17:   eventInfo->dof[1]    = -1.0;
 18:   eventInfo->dof[2]    = -1.0;
 19:   eventInfo->dof[3]    = -1.0;
 20:   eventInfo->dof[4]    = -1.0;
 21:   eventInfo->dof[5]    = -1.0;
 22:   eventInfo->dof[6]    = -1.0;
 23:   eventInfo->dof[7]    = -1.0;
 24:   eventInfo->errors[0] = -1.0;
 25:   eventInfo->errors[1] = -1.0;
 26:   eventInfo->errors[2] = -1.0;
 27:   eventInfo->errors[3] = -1.0;
 28:   eventInfo->errors[4] = -1.0;
 29:   eventInfo->errors[5] = -1.0;
 30:   eventInfo->errors[6] = -1.0;
 31:   eventInfo->errors[7] = -1.0;
 32:   PetscFunctionReturn(PETSC_SUCCESS);
 33: }

 35: static PetscErrorCode PetscEventPerfInfoTic_Internal(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event, PetscBool resume)
 36: {
 37:   PetscFunctionBegin;
 38:   if (resume) {
 39:     eventInfo->timeTmp -= time;
 40:     eventInfo->flopsTmp -= petsc_TotalFlops_th;
 41:   } else {
 42:     eventInfo->timeTmp  = -time;
 43:     eventInfo->flopsTmp = -petsc_TotalFlops_th;
 44:   }
 45:   eventInfo->numMessages -= petsc_irecv_ct_th + petsc_isend_ct_th + petsc_recv_ct_th + petsc_send_ct_th;
 46:   eventInfo->messageLength -= petsc_irecv_len_th + petsc_isend_len_th + petsc_recv_len_th + petsc_send_len_th;
 47:   eventInfo->numReductions -= petsc_allreduce_ct_th + petsc_gather_ct_th + petsc_scatter_ct_th;
 48: #if defined(PETSC_HAVE_DEVICE)
 49:   eventInfo->CpuToGpuCount -= petsc_ctog_ct_th;
 50:   eventInfo->GpuToCpuCount -= petsc_gtoc_ct_th;
 51:   eventInfo->CpuToGpuSize -= petsc_ctog_sz_th;
 52:   eventInfo->GpuToCpuSize -= petsc_gtoc_sz_th;
 53:   eventInfo->GpuFlops -= petsc_gflops_th;
 54:   eventInfo->GpuTime -= petsc_gtime;
 55: #endif
 56:   if (logMemory) {
 57:     PetscLogDouble usage;
 58:     PetscCall(PetscMemoryGetCurrentUsage(&usage));
 59:     eventInfo->memIncrease -= usage;
 60:     PetscCall(PetscMallocGetCurrentUsage(&usage));
 61:     eventInfo->mallocSpace -= usage;
 62:     PetscCall(PetscMallocGetMaximumUsage(&usage));
 63:     eventInfo->mallocIncrease -= usage;
 64:     PetscCall(PetscMallocPushMaximumUsage(event));
 65:   }
 66:   PetscFunctionReturn(PETSC_SUCCESS);
 67: }

 69: static PetscErrorCode PetscEventPerfInfoTic(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event)
 70: {
 71:   PetscFunctionBegin;
 72:   PetscCall(PetscEventPerfInfoTic_Internal(eventInfo, time, logMemory, event, PETSC_FALSE));
 73:   PetscFunctionReturn(PETSC_SUCCESS);
 74: }

 76: static PetscErrorCode PetscEventPerfInfoResume(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event)
 77: {
 78:   PetscFunctionBegin;
 79:   PetscCall(PetscEventPerfInfoTic_Internal(eventInfo, time, logMemory, event, PETSC_TRUE));
 80:   PetscFunctionReturn(PETSC_SUCCESS);
 81: }

 83: static PetscErrorCode PetscEventPerfInfoToc_Internal(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event, PetscBool pause)
 84: {
 85:   PetscFunctionBegin;
 86:   eventInfo->timeTmp += time;
 87:   eventInfo->flopsTmp += petsc_TotalFlops_th;
 88:   if (!pause) {
 89:     eventInfo->time += eventInfo->timeTmp;
 90:     eventInfo->time2 += eventInfo->timeTmp * eventInfo->timeTmp;
 91:     eventInfo->flops += eventInfo->flopsTmp;
 92:     eventInfo->flops2 += eventInfo->flopsTmp * eventInfo->flopsTmp;
 93:   }
 94:   eventInfo->numMessages += petsc_irecv_ct_th + petsc_isend_ct_th + petsc_recv_ct_th + petsc_send_ct_th;
 95:   eventInfo->messageLength += petsc_irecv_len_th + petsc_isend_len_th + petsc_recv_len + petsc_send_len_th;
 96:   eventInfo->numReductions += petsc_allreduce_ct_th + petsc_gather_ct_th + petsc_scatter_ct_th;
 97: #if defined(PETSC_HAVE_DEVICE)
 98:   eventInfo->CpuToGpuCount += petsc_ctog_ct_th;
 99:   eventInfo->GpuToCpuCount += petsc_gtoc_ct_th;
100:   eventInfo->CpuToGpuSize += petsc_ctog_sz_th;
101:   eventInfo->GpuToCpuSize += petsc_gtoc_sz_th;
102:   eventInfo->GpuFlops += petsc_gflops_th;
103:   eventInfo->GpuTime += petsc_gtime;
104: #endif
105:   if (logMemory) {
106:     PetscLogDouble usage, musage;
107:     PetscCall(PetscMemoryGetCurrentUsage(&usage)); /* the comments below match the column labels printed in PetscLogView_Default() */
108:     eventInfo->memIncrease += usage;               /* RMI */
109:     PetscCall(PetscMallocGetCurrentUsage(&usage));
110:     eventInfo->mallocSpace += usage; /* Malloc */
111:     PetscCall(PetscMallocPopMaximumUsage((int)event, &musage));
112:     eventInfo->mallocIncreaseEvent = PetscMax(musage - usage, eventInfo->mallocIncreaseEvent); /* EMalloc */
113:     PetscCall(PetscMallocGetMaximumUsage(&usage));
114:     eventInfo->mallocIncrease += usage; /* MMalloc */
115:   }
116:   PetscFunctionReturn(PETSC_SUCCESS);
117: }

119: static PetscErrorCode PetscEventPerfInfoToc(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event)
120: {
121:   PetscFunctionBegin;
122:   PetscCall(PetscEventPerfInfoToc_Internal(eventInfo, time, logMemory, event, PETSC_FALSE));
123:   PetscFunctionReturn(PETSC_SUCCESS);
124: }

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

133: static PetscErrorCode PetscEventPerfInfoAdd_Internal(const PetscEventPerfInfo *eventInfo, PetscEventPerfInfo *outInfo)
134: {
135:   PetscFunctionBegin;
136:   outInfo->count += eventInfo->count;
137:   outInfo->time += eventInfo->time;
138:   outInfo->time2 += eventInfo->time2;
139:   outInfo->flops += eventInfo->flops;
140:   outInfo->flops2 += eventInfo->flops2;
141:   outInfo->numMessages += eventInfo->numMessages;
142:   outInfo->messageLength += eventInfo->messageLength;
143:   outInfo->numReductions += eventInfo->numReductions;
144: #if defined(PETSC_HAVE_DEVICE)
145:   outInfo->CpuToGpuCount += eventInfo->CpuToGpuCount;
146:   outInfo->GpuToCpuCount += eventInfo->GpuToCpuCount;
147:   outInfo->CpuToGpuSize += eventInfo->CpuToGpuSize;
148:   outInfo->GpuToCpuSize += eventInfo->GpuToCpuSize;
149:   outInfo->GpuFlops += eventInfo->GpuFlops;
150:   outInfo->GpuTime += eventInfo->GpuTime;
151: #endif
152:   outInfo->memIncrease += eventInfo->memIncrease;
153:   outInfo->mallocSpace += eventInfo->mallocSpace;
154:   outInfo->mallocIncreaseEvent += eventInfo->mallocIncreaseEvent;
155:   outInfo->mallocIncrease += eventInfo->mallocIncrease;
156:   PetscFunctionReturn(PETSC_SUCCESS);
157: }

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

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

163: typedef struct {
164:   int            creations;    /* The number of objects of this class created */
165:   int            destructions; /* The number of objects of this class destroyed */
166:   PetscLogDouble mem;          /* The total memory allocated by objects of this class; this is completely wrong and should possibly be removed */
167:   PetscLogDouble descMem;      /* The total memory allocated by descendents of these objects; this is completely wrong and should possibly be removed */
168: } PetscClassPerf;

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

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

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

181: typedef struct _PetscStagePerf {
182:   PetscBool              used;     /* The stage was pushed on this processor */
183:   PetscEventPerfInfo     perfInfo; /* The stage performance information */
184:   PetscLogEventPerfArray eventLog; /* The event information for this stage */
185:   PetscLogClassPerfArray classLog; /* The class information for this stage */
186: } PetscStagePerf;

188: static PetscErrorCode PetscStageInfoInit(PetscStagePerf *stageInfo)
189: {
190:   PetscFunctionBegin;
191:   PetscCall(PetscMemzero(stageInfo, sizeof(*stageInfo)));
192:   PetscCall(PetscLogEventPerfArrayCreate(128, &(stageInfo->eventLog)));
193:   PetscCall(PetscLogClassPerfArrayCreate(128, &(stageInfo->classLog)));
194:   PetscCall(PetscEventPerfInfoInit(&stageInfo->perfInfo));
195:   PetscFunctionReturn(PETSC_SUCCESS);
196: }

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

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

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

210: /* The structure for action logging */
211: typedef enum {
212:   PETSC_LOG_ACTION_CREATE,
213:   PETSC_LOG_ACTION_DESTROY,
214:   PETSC_LOG_ACTION_BEGIN,
215:   PETSC_LOG_ACTION_END,
216: } PetscLogActionType;

218: typedef struct _Action {
219:   PetscLogActionType action;        /* The type of execution */
220:   PetscLogEvent      event;         /* The event number */
221:   PetscClassId       classid;       /* The event class id */
222:   PetscLogDouble     time;          /* The time of occurrence */
223:   PetscLogDouble     flops;         /* The cumulative flops */
224:   PetscLogDouble     mem;           /* The current memory usage */
225:   PetscLogDouble     maxmem;        /* The maximum memory usage */
226:   int                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)));

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

415:     new_object.parent = -1;
416:     new_object.obj    = obj;
417:     new_object.mem    = 0;

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

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

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

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

471:     PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1");
472:     PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, obj->id - 1, &obj_entry));
473:     if (obj->name) PetscCall(PetscStrncpy(obj_entry->name, obj->name, 64));
474:     obj_entry->obj = NULL;
475:   }
476:   PetscCall(PetscSpinlockUnlock(&def->lock));
477:   PetscFunctionReturn(PETSC_SUCCESS);
478: }

480: static PetscErrorCode PetscLogHandlerEventSync_Default(PetscLogHandler h, PetscLogEvent event, MPI_Comm comm)
481: {
482:   PetscLogState       state;
483:   PetscLogEventInfo   event_info;
484:   PetscEventPerfInfo *event_perf_info;
485:   int                 stage;
486:   PetscLogDouble      time = 0.0;

488:   PetscFunctionBegin;
489:   if (!(PetscLogSyncOn) || comm == MPI_COMM_NULL) PetscFunctionReturn(PETSC_SUCCESS);
490:   PetscCall(PetscLogHandlerGetState(h, &state));
491:   PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
492:   if (!event_info.collective) PetscFunctionReturn(PETSC_SUCCESS);
493:   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
494:   PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
495:   if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS);

497:   PetscCall(PetscTimeSubtract(&time));
498:   PetscCallMPI(MPI_Barrier(comm));
499:   PetscCall(PetscTimeAdd(&time));
500:   event_perf_info->syncTime += time;
501:   PetscFunctionReturn(PETSC_SUCCESS);
502: }

504: static PetscErrorCode PetscLogGetStageEventPerfInfo_threaded(PetscLogHandler_Default def, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **eventInfo)
505: {
506:   PetscEventPerfInfo *leventInfo = NULL;
507:   PetscHashIJKKey     key;

509:   PetscFunctionBegin;

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

530: static PetscErrorCode PetscLogHandlerEventBegin_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
531: {
532:   PetscLogHandler_Default def             = (PetscLogHandler_Default)h->data;
533:   PetscEventPerfInfo     *event_perf_info = NULL;
534:   PetscLogEventInfo       event_info;
535:   PetscLogDouble          time;
536:   PetscLogState           state;
537:   PetscLogStage           stage;

539:   PetscFunctionBegin;
540:   PetscCall(PetscLogHandlerGetState(h, &state));
541:   /* Synchronization */
542:   PetscCall(PetscLogHandlerEventSync_Default(h, event, PetscObjectComm(o1)));
543:   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
544:   if (def->pause_depth > 0) stage = 0; // in pause-mode, all events run on the main stage
545:   if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
546:     PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info));
547:     if (event_perf_info->depth == 0) { PetscCall(PetscEventPerfInfoInit(event_perf_info)); }
548:   } else {
549:     PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
550:   }
551:   PetscCheck(event_perf_info->depth >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to begin a paused event, this is not allowed");
552:   event_perf_info->depth++;
553:   /* Check for double counting */
554:   if (event_perf_info->depth > 1) PetscFunctionReturn(PETSC_SUCCESS);
555:   PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
556:   /* Log the performance info */
557:   event_perf_info->count++;
558:   PetscCall(PetscTime(&time));
559:   PetscCall(PetscEventPerfInfoTic(event_perf_info, time, PetscLogMemory, (int)event));
560:   if (def->petsc_logActions) {
561:     PetscLogDouble curTime;
562:     Action         new_action;

564:     PetscCall(PetscTime(&curTime));
565:     new_action.time    = curTime - petsc_BaseTime;
566:     new_action.action  = PETSC_LOG_ACTION_BEGIN;
567:     new_action.event   = event;
568:     new_action.classid = event_info.classid;
569:     new_action.id1     = o1 ? o1->id : -1;
570:     new_action.id2     = o2 ? o2->id : -1;
571:     new_action.id3     = o3 ? o3->id : -1;
572:     new_action.flops   = petsc_TotalFlops;
573:     PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
574:     PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
575:     PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
576:   }
577:   PetscFunctionReturn(PETSC_SUCCESS);
578: }

580: static PetscErrorCode PetscLogHandlerEventEnd_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
581: {
582:   PetscLogHandler_Default def             = (PetscLogHandler_Default)h->data;
583:   PetscEventPerfInfo     *event_perf_info = NULL;
584:   PetscLogDouble          time;
585:   PetscLogState           state;
586:   int                     stage;

588:   PetscFunctionBegin;
589:   PetscCall(PetscLogHandlerGetState(h, &state));
590:   if (def->petsc_logActions) {
591:     PetscLogEventInfo event_info;
592:     PetscLogDouble    curTime;
593:     Action            new_action;

595:     PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
596:     PetscCall(PetscTime(&curTime));
597:     new_action.time    = curTime - petsc_BaseTime;
598:     new_action.action  = PETSC_LOG_ACTION_END;
599:     new_action.event   = event;
600:     new_action.classid = event_info.classid;
601:     new_action.id1     = o1 ? o1->id : -1;
602:     new_action.id2     = o2 ? o2->id : -2;
603:     new_action.id3     = o3 ? o3->id : -3;
604:     new_action.flops   = petsc_TotalFlops;
605:     PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
606:     PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
607:     PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
608:   }
609:   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
610:   if (def->pause_depth > 0) stage = 0; // all events run on the main stage in pause-mode
611:   if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
612:     PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info));
613:   } else {
614:     PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
615:   }
616:   PetscCheck(event_perf_info->depth > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to end paused event, not allowed");
617:   event_perf_info->depth--;
618:   /* Check for double counting */
619:   if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
620:   else PetscCheck(event_perf_info->depth == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Logging event had unbalanced begin/end pairs");

622:   /* Log performance info */
623:   PetscCall(PetscTime(&time));
624:   PetscCall(PetscEventPerfInfoToc(event_perf_info, time, PetscLogMemory, (int)event));
625:   if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
626:     PetscEventPerfInfo *event_perf_info_global;
627:     PetscCall(PetscSpinlockLock(&def->lock));
628:     PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info_global));
629:     PetscCall(PetscEventPerfInfoAdd_Internal(event_perf_info, event_perf_info_global));
630:     PetscCall(PetscSpinlockUnlock(&def->lock));
631:   }
632:   PetscFunctionReturn(PETSC_SUCCESS);
633: }

635: static PetscErrorCode PetscLogHandlerEventDeactivatePush_Default(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event)
636: {
637:   PetscEventPerfInfo *event_perf_info;

639:   PetscFunctionBegin;
640:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
641:   PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
642:   event_perf_info->depth++;
643:   PetscFunctionReturn(PETSC_SUCCESS);
644: }

646: static PetscErrorCode PetscLogHandlerEventDeactivatePop_Default(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event)
647: {
648:   PetscEventPerfInfo *event_perf_info;

650:   PetscFunctionBegin;
651:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
652:   PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
653:   event_perf_info->depth--;
654:   PetscFunctionReturn(PETSC_SUCCESS);
655: }

657: static PetscErrorCode PetscLogHandlerEventsPause_Default(PetscLogHandler h)
658: {
659:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
660:   PetscLogDouble          time;
661:   PetscInt                num_stages;

663:   PetscFunctionBegin;
664:   if (def->pause_depth++ > 0) PetscFunctionReturn(PETSC_SUCCESS);
665:   PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL));
666:   PetscCall(PetscTime(&time));
667:   for (PetscInt stage = 0; stage < num_stages; stage++) {
668:     PetscStagePerf *stage_info = NULL;
669:     PetscInt        num_events;

671:     PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
672:     PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL));
673:     for (PetscInt event = 0; event < num_events; event++) {
674:       PetscEventPerfInfo *event_info = NULL;
675:       PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info));
676:       if (event_info->depth > 0) {
677:         event_info->depth *= -1;
678:         PetscCall(PetscEventPerfInfoPause(event_info, time, PetscLogMemory, event));
679:       }
680:     }
681:     if (stage > 0 && stage_info->perfInfo.depth > 0) {
682:       stage_info->perfInfo.depth *= -1;
683:       PetscCall(PetscEventPerfInfoPause(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2)));
684:     }
685:   }
686:   PetscFunctionReturn(PETSC_SUCCESS);
687: }

689: static PetscErrorCode PetscLogHandlerEventsResume_Default(PetscLogHandler h)
690: {
691:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
692:   PetscLogDouble          time;
693:   PetscInt                num_stages;

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

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

721: static PetscErrorCode PetscLogHandlerStagePush_Default(PetscLogHandler h, PetscLogStage new_stage)
722: {
723:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
724:   PetscLogDouble          time;
725:   PetscLogState           state;
726:   PetscLogStage           current_stage;
727:   PetscStagePerf         *new_stage_info;

729:   PetscFunctionBegin;
730:   if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
731:   PetscCall(PetscLogHandlerGetState(h, &state));
732:   current_stage = state->current_stage;
733:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, new_stage, &new_stage_info));
734:   PetscCall(PetscTime(&time));

736:   /* Record flops/time of previous stage */
737:   if (current_stage >= 0) {
738:     if (PetscBTLookup(state->active, current_stage)) {
739:       PetscStagePerf *current_stage_info;
740:       PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, &current_stage_info));
741:       PetscCall(PetscEventPerfInfoToc(&current_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2)));
742:     }
743:   }
744:   new_stage_info->used = PETSC_TRUE;
745:   new_stage_info->perfInfo.count++;
746:   new_stage_info->perfInfo.depth++;
747:   /* Subtract current quantities so that we obtain the difference when we pop */
748:   if (PetscBTLookup(state->active, new_stage)) PetscCall(PetscEventPerfInfoTic(&new_stage_info->perfInfo, time, PetscLogMemory, (int)-(new_stage + 2)));
749:   PetscFunctionReturn(PETSC_SUCCESS);
750: }

752: static PetscErrorCode PetscLogHandlerStagePop_Default(PetscLogHandler h, PetscLogStage old_stage)
753: {
754:   PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
755:   PetscLogStage           current_stage;
756:   PetscStagePerf         *old_stage_info;
757:   PetscLogState           state;
758:   PetscLogDouble          time;

760:   PetscFunctionBegin;
761:   if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
762:   PetscCall(PetscLogHandlerGetState(h, &state));
763:   current_stage = state->current_stage;
764:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, old_stage, &old_stage_info));
765:   PetscCall(PetscTime(&time));
766:   old_stage_info->perfInfo.depth--;
767:   if (PetscBTLookup(state->active, old_stage)) { PetscCall(PetscEventPerfInfoToc(&old_stage_info->perfInfo, time, PetscLogMemory, (int)-(old_stage + 2))); }
768:   if (current_stage >= 0) {
769:     if (PetscBTLookup(state->active, current_stage)) {
770:       PetscStagePerf *current_stage_info;
771:       PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, &current_stage_info));
772:       PetscCall(PetscEventPerfInfoTic(&current_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2)));
773:     }
774:   }
775:   PetscFunctionReturn(PETSC_SUCCESS);
776: }

778: static PetscErrorCode PetscLogHandlerStageSetVisible_Default(PetscLogHandler h, PetscLogStage stage, PetscBool is_visible)
779: {
780:   PetscStagePerf *stage_info;

782:   PetscFunctionBegin;
783:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
784:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info));
785:   stage_info->perfInfo.visible = is_visible;
786:   PetscFunctionReturn(PETSC_SUCCESS);
787: }

789: static PetscErrorCode PetscLogHandlerStageGetVisible_Default(PetscLogHandler h, PetscLogStage stage, PetscBool *is_visible)
790: {
791:   PetscStagePerf *stage_info;

793:   PetscFunctionBegin;
794:   if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
795:   PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info));
796:   *is_visible = stage_info->perfInfo.visible;
797:   PetscFunctionReturn(PETSC_SUCCESS);
798: }

800: static PetscErrorCode PetscLogHandlerSetLogActions_Default(PetscLogHandler handler, PetscBool flag)
801: {
802:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;

804:   PetscFunctionBegin;
805:   def->petsc_logActions = flag;
806:   PetscFunctionReturn(PETSC_SUCCESS);
807: }

809: static PetscErrorCode PetscLogHandlerSetLogObjects_Default(PetscLogHandler handler, PetscBool flag)
810: {
811:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;

813:   PetscFunctionBegin;
814:   def->petsc_logObjects = flag;
815:   PetscFunctionReturn(PETSC_SUCCESS);
816: }

818: static PetscErrorCode PetscLogHandlerLogObjectState_Default(PetscLogHandler handler, PetscObject obj, const char format[], va_list Argp)
819: {
820:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
821:   size_t                  fullLength;

823:   PetscFunctionBegin;
824:   if (def->petsc_logObjects) {
825:     Object *obj_entry = NULL;

827:     PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, obj->id - 1, &obj_entry));
828:     PetscCall(PetscVSNPrintf(obj_entry->info, 64, format, &fullLength, Argp));
829:   }
830:   PetscFunctionReturn(PETSC_SUCCESS);
831: }

833: static PetscErrorCode PetscLogHandlerGetNumObjects_Default(PetscLogHandler handler, PetscInt *num_objects)
834: {
835:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;

837:   PetscFunctionBegin;
838:   PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, num_objects, NULL));
839:   PetscFunctionReturn(PETSC_SUCCESS);
840: }

842: static PetscErrorCode PetscLogHandlerDump_Default(PetscLogHandler handler, const char sname[])
843: {
844:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
845:   FILE                   *fd;
846:   char                    file[PETSC_MAX_PATH_LEN], fname[PETSC_MAX_PATH_LEN];
847:   PetscLogDouble          flops, _TotalTime;
848:   PetscMPIInt             rank;
849:   int                     curStage;
850:   PetscLogState           state;
851:   PetscInt                num_events;
852:   PetscLogEvent           event;

854:   PetscFunctionBegin;
855:   /* Calculate the total elapsed time */
856:   PetscCall(PetscTime(&_TotalTime));
857:   _TotalTime -= petsc_BaseTime;
858:   /* Open log file */
859:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)handler), &rank));
860:   PetscCall(PetscSNPrintf(file, PETSC_STATIC_ARRAY_LENGTH(file), "%s.%d", sname && sname[0] ? sname : "Log", rank));
861:   PetscCall(PetscFixFilename(file, fname));
862:   PetscCall(PetscFOpen(PETSC_COMM_SELF, fname, "w", &fd));
863:   PetscCheck(!(rank == 0) || !(!fd), PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file: %s", fname);
864:   /* Output totals */
865:   PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Total Flop %14e %16.8e\n", petsc_TotalFlops, _TotalTime));
866:   PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Clock Resolution %g\n", 0.0));
867:   /* Output actions */
868:   if (def->petsc_logActions) {
869:     PetscInt num_actions;
870:     PetscCall(PetscLogActionArrayGetSize(def->petsc_actions, &num_actions, NULL));
871:     PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Actions accomplished %d\n", (int)num_actions));
872:     for (int a = 0; a < num_actions; a++) {
873:       Action *action;

875:       PetscCall(PetscLogActionArrayGetRef(def->petsc_actions, a, &action));
876:       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "%g %d %d %d %d %d %d %g %g %g\n", action->time, action->action, (int)action->event, (int)action->classid, action->id1, action->id2, action->id3, action->flops, action->mem, action->maxmem));
877:     }
878:   }
879:   /* Output objects */
880:   if (def->petsc_logObjects) {
881:     PetscInt num_objects;

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

888:       PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, o, &object));
889:       if (object->parent != -1) continue; // object with this id wasn't logged, probably a PetscLogHandler
890:       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Parent ID: %d Memory: %d\n", object->parent, (int)object->mem));
891:       if (!object->name[0]) {
892:         PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Name\n"));
893:       } else {
894:         PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Name: %s\n", object->name));
895:       }
896:       if (!object->info[0]) {
897:         PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Info\n"));
898:       } else {
899:         PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Info: %s\n", object->info));
900:       }
901:     }
902:   }
903:   /* Output events */
904:   PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Event log:\n"));
905:   PetscCall(PetscLogHandlerGetState(handler, &state));
906:   PetscCall(PetscLogStateGetNumEvents(state, &num_events));
907:   PetscCall(PetscLogStateGetCurrentStage(state, &curStage));
908:   for (event = 0; event < num_events; event++) {
909:     PetscEventPerfInfo *event_info;

911:     PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, curStage, event, &event_info));
912:     if (event_info->time != 0.0) flops = event_info->flops / event_info->time;
913:     else flops = 0.0;
914:     PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "%d %16d %16g %16g %16g\n", event, event_info->count, event_info->flops, event_info->time, flops));
915:   }
916:   PetscCall(PetscFClose(PETSC_COMM_SELF, fd));
917:   PetscFunctionReturn(PETSC_SUCCESS);
918: }

920: /*
921:   PetscLogView_Detailed - Each process prints the times for its own events

923: */
924: static PetscErrorCode PetscLogHandlerView_Default_Detailed(PetscLogHandler handler, PetscViewer viewer)
925: {
926:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
927:   PetscLogDouble          locTotalTime, numRed, maxMem;
928:   PetscInt                numStages, numEvents;
929:   MPI_Comm                comm = PetscObjectComm((PetscObject)viewer);
930:   PetscMPIInt             rank, size;
931:   PetscLogGlobalNames     global_stages, global_events;
932:   PetscLogState           state;
933:   PetscEventPerfInfo      zero_info;

935:   PetscFunctionBegin;
936:   PetscCall(PetscLogHandlerGetState(handler, &state));
937:   PetscCallMPI(MPI_Comm_size(comm, &size));
938:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
939:   /* Must preserve reduction count before we go on */
940:   numRed = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
941:   /* Get the total elapsed time */
942:   PetscCall(PetscTime(&locTotalTime));
943:   locTotalTime -= petsc_BaseTime;
944:   PetscCall(PetscViewerASCIIPrintf(viewer, "size = %d\n", size));
945:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalTimes = {}\n"));
946:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessages = {}\n"));
947:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessageLens = {}\n"));
948:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalReductions = {}\n"));
949:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalFlop = {}\n"));
950:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalObjects = {}\n"));
951:   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMemory = {}\n"));
952:   PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
953:   PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
954:   PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
955:   PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
956:   PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
957:   PetscCall(PetscViewerASCIIPrintf(viewer, "Stages = {}\n"));
958:   for (PetscInt stage = 0; stage < numStages; stage++) {
959:     PetscInt    stage_id;
960:     const char *stage_name;

962:     PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
963:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
964:     PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"] = {}\n", stage_name));
965:     PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"summary\"] = {}\n", stage_name));
966:     for (PetscInt event = 0; event < numEvents; event++) {
967:       PetscEventPerfInfo *eventInfo = &zero_info;
968:       PetscBool           is_zero   = PETSC_FALSE;
969:       PetscInt            event_id;
970:       const char         *event_name;

972:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
973:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
974:       if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
975:       is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE;
976:       PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
977:       if (!is_zero) { PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"%s\"] = {}\n", stage_name, event_name)); }
978:     }
979:   }
980:   PetscCall(PetscMallocGetMaximumUsage(&maxMem));
981:   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
982:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalTimes[%d] = %g\n", rank, locTotalTime));
983:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessages[%d] = %g\n", rank, (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct)));
984:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessageLens[%d] = %g\n", rank, (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len)));
985:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalReductions[%d] = %g\n", rank, numRed));
986:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalFlop[%d] = %g\n", rank, petsc_TotalFlops));
987:   {
988:     PetscInt num_objects;

990:     PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
991:     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalObjects[%d] = %d\n", rank, (int)num_objects));
992:   }
993:   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMemory[%d] = %g\n", rank, maxMem));
994:   PetscCall(PetscViewerFlush(viewer));
995:   for (PetscInt stage = 0; stage < numStages; stage++) {
996:     PetscEventPerfInfo *stage_perf_info = &zero_info;
997:     PetscInt            stage_id;
998:     const char         *stage_name;

1000:     PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1001:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1002:     if (stage_id >= 0) {
1003:       PetscStagePerf *stage_info;
1004:       PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info));
1005:       stage_perf_info = &stage_info->perfInfo;
1006:     }
1007:     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,
1008:                                                  stage_perf_info->numMessages, stage_perf_info->messageLength, stage_perf_info->numReductions, stage_perf_info->flops));
1009:     for (PetscInt event = 0; event < numEvents; event++) {
1010:       PetscEventPerfInfo *eventInfo = &zero_info;
1011:       PetscBool           is_zero   = PETSC_FALSE;
1012:       PetscInt            event_id;
1013:       const char         *event_name;

1015:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1016:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1017:       if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1018:       is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE;
1019:       PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero));
1020:       PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1021:       if (!is_zero) {
1022:         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,
1023:                                                      eventInfo->count, eventInfo->time, eventInfo->syncTime, eventInfo->numMessages, eventInfo->messageLength, eventInfo->numReductions, eventInfo->flops));
1024:         if (eventInfo->dof[0] >= 0.) {
1025:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"dof\" : ["));
1026:           for (PetscInt d = 0; d < 8; ++d) {
1027:             if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", "));
1028:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->dof[d]));
1029:           }
1030:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]"));
1031:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"error\" : ["));
1032:           for (PetscInt e = 0; e < 8; ++e) {
1033:             if (e > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", "));
1034:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->errors[e]));
1035:           }
1036:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]"));
1037:         }
1038:       }
1039:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "}\n"));
1040:     }
1041:   }
1042:   PetscCall(PetscViewerFlush(viewer));
1043:   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1044:   PetscCall(PetscLogGlobalNamesDestroy(&global_events));
1045:   PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1046:   PetscFunctionReturn(PETSC_SUCCESS);
1047: }

1049: /*
1050:   PetscLogView_CSV - Each process prints the times for its own events in Comma-Separated Value Format
1051: */
1052: static PetscErrorCode PetscLogHandlerView_Default_CSV(PetscLogHandler handler, PetscViewer viewer)
1053: {
1054:   PetscLogDouble      locTotalTime, maxMem;
1055:   PetscInt            numStages, numEvents, stage, event;
1056:   MPI_Comm            comm = PetscObjectComm((PetscObject)viewer);
1057:   PetscMPIInt         rank, size;
1058:   PetscLogGlobalNames global_stages, global_events;
1059:   PetscLogState       state;
1060:   PetscEventPerfInfo  zero_info;

1062:   PetscFunctionBegin;
1063:   PetscCall(PetscLogHandlerGetState(handler, &state));
1064:   PetscCallMPI(MPI_Comm_size(comm, &size));
1065:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1066:   /* Must preserve reduction count before we go on */
1067:   /* Get the total elapsed time */
1068:   PetscCall(PetscTime(&locTotalTime));
1069:   locTotalTime -= petsc_BaseTime;
1070:   PetscCall(PetscMallocGetMaximumUsage(&maxMem));
1071:   PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
1072:   PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
1073:   PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
1074:   PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
1075:   PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
1076:   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1077:   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));
1078:   PetscCall(PetscViewerFlush(viewer));
1079:   for (stage = 0; stage < numStages; stage++) {
1080:     PetscEventPerfInfo *stage_perf_info;
1081:     PetscInt            stage_id;
1082:     const char         *stage_name;

1084:     PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1085:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1086:     stage_perf_info = &zero_info;
1087:     if (stage_id >= 0) {
1088:       PetscStagePerf *stage_info;
1089:       PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info));
1090:       stage_perf_info = &stage_info->perfInfo;
1091:     }
1092:     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));
1093:     for (event = 0; event < numEvents; event++) {
1094:       PetscEventPerfInfo *eventInfo = &zero_info;
1095:       PetscBool           is_zero   = PETSC_FALSE;
1096:       PetscInt            event_id;
1097:       const char         *event_name;

1099:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1100:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1101:       if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1102:       PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero));
1103:       PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1104:       if (!is_zero) {
1105:         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));
1106:         if (eventInfo->dof[0] >= 0.) {
1107:           for (PetscInt d = 0; d < 8; ++d) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->dof[d]));
1108:           for (PetscInt e = 0; e < 8; ++e) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->errors[e]));
1109:         }
1110:       }
1111:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
1112:     }
1113:   }
1114:   PetscCall(PetscViewerFlush(viewer));
1115:   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1116:   PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1117:   PetscCall(PetscLogGlobalNamesDestroy(&global_events));
1118:   PetscFunctionReturn(PETSC_SUCCESS);
1119: }

1121: static PetscErrorCode PetscLogViewWarnSync(MPI_Comm comm, FILE *fd)
1122: {
1123:   PetscFunctionBegin;
1124:   if (!PetscLogSyncOn) PetscFunctionReturn(PETSC_SUCCESS);
1125:   PetscCall(PetscFPrintf(comm, fd, "\n\n"));
1126:   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n"));
1127:   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1128:   PetscCall(PetscFPrintf(comm, fd, "      #                       WARNING!!!                       #\n"));
1129:   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1130:   PetscCall(PetscFPrintf(comm, fd, "      #   This program was run with logging synchronization.   #\n"));
1131:   PetscCall(PetscFPrintf(comm, fd, "      #   This option provides more meaningful imbalance       #\n"));
1132:   PetscCall(PetscFPrintf(comm, fd, "      #   figures at the expense of slowing things down and    #\n"));
1133:   PetscCall(PetscFPrintf(comm, fd, "      #   providing a distorted view of the overall runtime.   #\n"));
1134:   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1135:   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n\n\n"));
1136:   PetscFunctionReturn(PETSC_SUCCESS);
1137: }

1139: static PetscErrorCode PetscLogViewWarnDebugging(MPI_Comm comm, FILE *fd)
1140: {
1141:   PetscFunctionBegin;
1142:   if (PetscDefined(USE_DEBUG)) {
1143:     PetscCall(PetscFPrintf(comm, fd, "\n\n"));
1144:     PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n"));
1145:     PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1146:     PetscCall(PetscFPrintf(comm, fd, "      #                       WARNING!!!                       #\n"));
1147:     PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1148:     PetscCall(PetscFPrintf(comm, fd, "      #   This code was compiled with a debugging option.      #\n"));
1149:     PetscCall(PetscFPrintf(comm, fd, "      #   To get timing results run ./configure                #\n"));
1150:     PetscCall(PetscFPrintf(comm, fd, "      #   using --with-debugging=no, the performance will      #\n"));
1151:     PetscCall(PetscFPrintf(comm, fd, "      #   be generally two or three times faster.              #\n"));
1152:     PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1153:     PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n\n\n"));
1154:   }
1155:   PetscFunctionReturn(PETSC_SUCCESS);
1156: }

1158: static PetscErrorCode PetscLogViewWarnNoGpuAwareMpi(MPI_Comm comm, FILE *fd)
1159: {
1160: #if defined(PETSC_HAVE_DEVICE)
1161:   PetscMPIInt size;
1162:   PetscBool   deviceInitialized = PETSC_FALSE;

1164:   PetscFunctionBegin;
1165:   PetscCallMPI(MPI_Comm_size(comm, &size));
1166:   for (int i = PETSC_DEVICE_HOST + 1; i < PETSC_DEVICE_MAX; ++i) {
1167:     const PetscDeviceType dtype = PetscDeviceTypeCast(i);
1168:     if (PetscDeviceInitialized(dtype)) { /* a non-host device was initialized */
1169:       deviceInitialized = PETSC_TRUE;
1170:       break;
1171:     }
1172:   }
1173:   /* the last condition says petsc is configured with device but it is a pure CPU run, so don't print misleading warnings */
1174:   if (use_gpu_aware_mpi || size == 1 || !deviceInitialized) PetscFunctionReturn(PETSC_SUCCESS);
1175:   PetscCall(PetscFPrintf(comm, fd, "\n\n"));
1176:   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n"));
1177:   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1178:   PetscCall(PetscFPrintf(comm, fd, "      #                       WARNING!!!                       #\n"));
1179:   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1180:   PetscCall(PetscFPrintf(comm, fd, "      #   This code was compiled with GPU support and you've   #\n"));
1181:   PetscCall(PetscFPrintf(comm, fd, "      #   created PETSc/GPU objects, but you intentionally     #\n"));
1182:   PetscCall(PetscFPrintf(comm, fd, "      #   used -use_gpu_aware_mpi 0, requiring PETSc to copy   #\n"));
1183:   PetscCall(PetscFPrintf(comm, fd, "      #   additional data between the GPU and CPU. To obtain   #\n"));
1184:   PetscCall(PetscFPrintf(comm, fd, "      #   meaningful timing results on multi-rank runs, use    #\n"));
1185:   PetscCall(PetscFPrintf(comm, fd, "      #   GPU-aware MPI instead.                               #\n"));
1186:   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1187:   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n\n\n"));
1188:   PetscFunctionReturn(PETSC_SUCCESS);
1189: #else
1190:   return PETSC_SUCCESS;
1191: #endif
1192: }

1194: static PetscErrorCode PetscLogViewWarnGpuTime(MPI_Comm comm, FILE *fd)
1195: {
1196: #if defined(PETSC_HAVE_DEVICE)

1198:   PetscFunctionBegin;
1199:   if (!PetscLogGpuTimeFlag || petsc_gflops == 0) PetscFunctionReturn(PETSC_SUCCESS);
1200:   PetscCall(PetscFPrintf(comm, fd, "\n\n"));
1201:   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n"));
1202:   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1203:   PetscCall(PetscFPrintf(comm, fd, "      #                       WARNING!!!                       #\n"));
1204:   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1205:   PetscCall(PetscFPrintf(comm, fd, "      #   This code was run with -log_view_gpu_time            #\n"));
1206:   PetscCall(PetscFPrintf(comm, fd, "      #   This provides accurate timing within the GPU kernels #\n"));
1207:   PetscCall(PetscFPrintf(comm, fd, "      #   but can slow down the entire computation by a        #\n"));
1208:   PetscCall(PetscFPrintf(comm, fd, "      #   measurable amount. For fastest runs we recommend     #\n"));
1209:   PetscCall(PetscFPrintf(comm, fd, "      #   not using this option.                               #\n"));
1210:   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1211:   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n\n\n"));
1212:   PetscFunctionReturn(PETSC_SUCCESS);
1213: #else
1214:   return PETSC_SUCCESS;
1215: #endif
1216: }

1218: static PetscErrorCode PetscLogHandlerView_Default_Info(PetscLogHandler handler, PetscViewer viewer)
1219: {
1220:   FILE                   *fd;
1221:   PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
1222:   char                    arch[128], hostname[128], username[128], pname[PETSC_MAX_PATH_LEN], date[128];
1223:   PetscLogDouble          locTotalTime, TotalTime, TotalFlops;
1224:   PetscLogDouble          numMessages, messageLength, avgMessLen, numReductions;
1225:   PetscLogDouble          stageTime, flops, flopr, mem, mess, messLen, red;
1226:   PetscLogDouble          fracTime, fracFlops, fracMessages, fracLength, fracReductions, fracMess, fracMessLen, fracRed;
1227:   PetscLogDouble          fracStageTime, fracStageFlops, fracStageMess, fracStageMessLen, fracStageRed;
1228:   PetscLogDouble          min, max, tot, ratio, avg, x, y;
1229:   PetscLogDouble          minf, maxf, totf, ratf, mint, maxt, tott, ratt, ratC, totm, totml, totr, mal, malmax, emalmax;
1230: #if defined(PETSC_HAVE_DEVICE)
1231:   PetscLogEvent  KSP_Solve, SNES_Solve, TS_Step, TAO_Solve; /* These need to be fixed to be some events registered with certain objects */
1232:   PetscLogDouble cct, gct, csz, gsz, gmaxt, gflops, gflopr, fracgflops;
1233: #endif
1234:   PetscMPIInt   minC, maxC;
1235:   PetscMPIInt   size, rank;
1236:   PetscBool    *localStageUsed, *stageUsed;
1237:   PetscBool    *localStageVisible, *stageVisible;
1238:   PetscInt      numStages, numEvents;
1239:   int           stage, oclass;
1240:   PetscLogEvent event;
1241:   char          version[256];
1242:   MPI_Comm      comm;
1243: #if defined(PETSC_HAVE_DEVICE)
1244:   PetscInt64 nas = 0x7FF0000000000002;
1245: #endif
1246:   PetscLogGlobalNames global_stages, global_events;
1247:   PetscEventPerfInfo  zero_info;
1248:   PetscLogState       state;

1250:   PetscFunctionBegin;
1251:   PetscCall(PetscLogHandlerGetState(handler, &state));
1252:   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1253:   PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
1254:   PetscCall(PetscViewerASCIIGetPointer(viewer, &fd));
1255:   PetscCallMPI(MPI_Comm_size(comm, &size));
1256:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1257:   /* Get the total elapsed time */
1258:   PetscCall(PetscTime(&locTotalTime));
1259:   locTotalTime -= petsc_BaseTime;

1261:   PetscCall(PetscFPrintf(comm, fd, "****************************************************************************************************************************************************************\n"));
1262:   PetscCall(PetscFPrintf(comm, fd, "***                                WIDEN YOUR WINDOW TO 160 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document                                 ***\n"));
1263:   PetscCall(PetscFPrintf(comm, fd, "****************************************************************************************************************************************************************\n"));
1264:   PetscCall(PetscFPrintf(comm, fd, "\n------------------------------------------------------------------ PETSc Performance Summary: ------------------------------------------------------------------\n\n"));
1265:   PetscCall(PetscLogViewWarnSync(comm, fd));
1266:   PetscCall(PetscLogViewWarnDebugging(comm, fd));
1267:   PetscCall(PetscLogViewWarnNoGpuAwareMpi(comm, fd));
1268:   PetscCall(PetscLogViewWarnGpuTime(comm, fd));
1269:   PetscCall(PetscGetArchType(arch, sizeof(arch)));
1270:   PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
1271:   PetscCall(PetscGetUserName(username, sizeof(username)));
1272:   PetscCall(PetscGetProgramName(pname, sizeof(pname)));
1273:   PetscCall(PetscGetDate(date, sizeof(date)));
1274:   PetscCall(PetscGetVersion(version, sizeof(version)));
1275:   if (size == 1) {
1276:     PetscCall(PetscFPrintf(comm, fd, "%s on a %s named %s with %d processor, by %s %s\n", pname, arch, hostname, size, username, date));
1277:   } else {
1278:     PetscCall(PetscFPrintf(comm, fd, "%s on a %s named %s with %d processors, by %s %s\n", pname, arch, hostname, size, username, date));
1279:   }
1280: #if defined(PETSC_HAVE_OPENMP)
1281:   PetscCall(PetscFPrintf(comm, fd, "Using %" PetscInt_FMT " OpenMP threads\n", PetscNumOMPThreads));
1282: #endif
1283:   PetscCall(PetscFPrintf(comm, fd, "Using %s\n", version));

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

1288:   /* Calculate summary information */
1289:   PetscCall(PetscFPrintf(comm, fd, "\n                         Max       Max/Min     Avg       Total\n"));
1290:   /*   Time */
1291:   PetscCall(MPIU_Allreduce(&locTotalTime, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1292:   PetscCall(MPIU_Allreduce(&locTotalTime, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1293:   PetscCall(MPIU_Allreduce(&locTotalTime, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1294:   avg = tot / ((PetscLogDouble)size);
1295:   if (min != 0.0) ratio = max / min;
1296:   else ratio = 0.0;
1297:   PetscCall(PetscFPrintf(comm, fd, "Time (sec):           %5.3e   %7.3f   %5.3e\n", max, ratio, avg));
1298:   TotalTime = tot;
1299:   /*   Objects */
1300:   {
1301:     PetscInt num_objects;

1303:     PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
1304:     avg = (PetscLogDouble)num_objects;
1305:   }
1306:   PetscCall(MPIU_Allreduce(&avg, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1307:   PetscCall(MPIU_Allreduce(&avg, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1308:   PetscCall(MPIU_Allreduce(&avg, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1309:   avg = tot / ((PetscLogDouble)size);
1310:   if (min != 0.0) ratio = max / min;
1311:   else ratio = 0.0;
1312:   PetscCall(PetscFPrintf(comm, fd, "Objects:              %5.3e   %7.3f   %5.3e\n", max, ratio, avg));
1313:   /*   Flops */
1314:   PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1315:   PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1316:   PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1317:   avg = tot / ((PetscLogDouble)size);
1318:   if (min != 0.0) ratio = max / min;
1319:   else ratio = 0.0;
1320:   PetscCall(PetscFPrintf(comm, fd, "Flops:                %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1321:   TotalFlops = tot;
1322:   /*   Flops/sec -- Must talk to Barry here */
1323:   if (locTotalTime != 0.0) flops = petsc_TotalFlops / locTotalTime;
1324:   else flops = 0.0;
1325:   PetscCall(MPIU_Allreduce(&flops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1326:   PetscCall(MPIU_Allreduce(&flops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1327:   PetscCall(MPIU_Allreduce(&flops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1328:   avg = tot / ((PetscLogDouble)size);
1329:   if (min != 0.0) ratio = max / min;
1330:   else ratio = 0.0;
1331:   PetscCall(PetscFPrintf(comm, fd, "Flops/sec:            %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1332:   /*   Memory */
1333:   PetscCall(PetscMallocGetMaximumUsage(&mem));
1334:   if (mem > 0.0) {
1335:     PetscCall(MPIU_Allreduce(&mem, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1336:     PetscCall(MPIU_Allreduce(&mem, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1337:     PetscCall(MPIU_Allreduce(&mem, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1338:     avg = tot / ((PetscLogDouble)size);
1339:     if (min != 0.0) ratio = max / min;
1340:     else ratio = 0.0;
1341:     PetscCall(PetscFPrintf(comm, fd, "Memory (bytes):       %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1342:   }
1343:   /*   Messages */
1344:   mess = 0.5 * (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct);
1345:   PetscCall(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1346:   PetscCall(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1347:   PetscCall(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1348:   avg = tot / ((PetscLogDouble)size);
1349:   if (min != 0.0) ratio = max / min;
1350:   else ratio = 0.0;
1351:   PetscCall(PetscFPrintf(comm, fd, "MPI Msg Count:        %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1352:   numMessages = tot;
1353:   /*   Message Lengths */
1354:   mess = 0.5 * (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len);
1355:   PetscCall(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1356:   PetscCall(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1357:   PetscCall(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1358:   if (numMessages != 0) avg = tot / numMessages;
1359:   else avg = 0.0;
1360:   if (min != 0.0) ratio = max / min;
1361:   else ratio = 0.0;
1362:   PetscCall(PetscFPrintf(comm, fd, "MPI Msg Len (bytes):  %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1363:   messageLength = tot;
1364:   /*   Reductions */
1365:   PetscCall(MPIU_Allreduce(&red, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1366:   PetscCall(MPIU_Allreduce(&red, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1367:   PetscCall(MPIU_Allreduce(&red, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1368:   if (min != 0.0) ratio = max / min;
1369:   else ratio = 0.0;
1370:   PetscCall(PetscFPrintf(comm, fd, "MPI Reductions:       %5.3e   %7.3f\n", max, ratio));
1371:   numReductions = red; /* wrong because uses count from process zero */
1372:   PetscCall(PetscFPrintf(comm, fd, "\nFlop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)\n"));
1373:   PetscCall(PetscFPrintf(comm, fd, "                            e.g., VecAXPY() for real vectors of length N --> 2N flops\n"));
1374:   PetscCall(PetscFPrintf(comm, fd, "                            and VecAXPY() for complex vectors of length N --> 8N flops\n"));

1376:   PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
1377:   PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
1378:   PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
1379:   PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
1380:   PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
1381:   PetscCall(PetscMalloc1(numStages, &localStageUsed));
1382:   PetscCall(PetscMalloc1(numStages, &stageUsed));
1383:   PetscCall(PetscMalloc1(numStages, &localStageVisible));
1384:   PetscCall(PetscMalloc1(numStages, &stageVisible));
1385:   if (numStages > 0) {
1386:     for (stage = 0; stage < numStages; stage++) {
1387:       PetscInt stage_id;

1389:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1390:       if (stage_id >= 0) {
1391:         PetscStagePerf *stage_info;

1393:         PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info));
1394:         localStageUsed[stage]    = stage_info->used;
1395:         localStageVisible[stage] = stage_info->perfInfo.visible;
1396:       } else {
1397:         localStageUsed[stage]    = PETSC_FALSE;
1398:         localStageVisible[stage] = PETSC_TRUE;
1399:       }
1400:     }
1401:     PetscCall(MPIU_Allreduce(localStageUsed, stageUsed, numStages, MPIU_BOOL, MPI_LOR, comm));
1402:     PetscCall(MPIU_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm));
1403:     for (stage = 0; stage < numStages; stage++) {
1404:       if (stageUsed[stage] && stageVisible[stage]) {
1405:         PetscCall(PetscFPrintf(comm, fd, "\nSummary of Stages:   ----- Time ------  ----- Flop ------  --- Messages ---  -- Message Lengths --  -- Reductions --\n"));
1406:         PetscCall(PetscFPrintf(comm, fd, "                        Avg     %%Total     Avg     %%Total    Count   %%Total     Avg         %%Total    Count   %%Total\n"));
1407:         break;
1408:       }
1409:     }
1410:     for (stage = 0; stage < numStages; stage++) {
1411:       PetscInt            stage_id;
1412:       PetscEventPerfInfo *stage_info;
1413:       const char         *stage_name;

1415:       if (!(stageUsed[stage] && stageVisible[stage])) continue;
1416:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1417:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1418:       stage_info = &zero_info;
1419:       if (localStageUsed[stage]) {
1420:         PetscStagePerf *stage_perf_info;

1422:         PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_perf_info));
1423:         stage_info = &stage_perf_info->perfInfo;
1424:       }
1425:       PetscCall(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1426:       PetscCall(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1427:       PetscCall(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1428:       PetscCall(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1429:       PetscCall(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1430:       mess *= 0.5;
1431:       messLen *= 0.5;
1432:       red /= size;
1433:       if (TotalTime != 0.0) fracTime = stageTime / TotalTime;
1434:       else fracTime = 0.0;
1435:       if (TotalFlops != 0.0) fracFlops = flops / TotalFlops;
1436:       else fracFlops = 0.0;
1437:       /* Talk to Barry if (stageTime     != 0.0) flops          = (size*flops)/stageTime; else flops          = 0.0; */
1438:       if (numMessages != 0.0) fracMessages = mess / numMessages;
1439:       else fracMessages = 0.0;
1440:       if (mess != 0.0) avgMessLen = messLen / mess;
1441:       else avgMessLen = 0.0;
1442:       if (messageLength != 0.0) fracLength = messLen / messageLength;
1443:       else fracLength = 0.0;
1444:       if (numReductions != 0.0) fracReductions = red / numReductions;
1445:       else fracReductions = 0.0;
1446:       PetscCall(PetscFPrintf(comm, fd, "%2d: %15s: %6.4e %5.1f%%  %6.4e %5.1f%%  %5.3e %5.1f%%  %5.3e      %5.1f%%  %5.3e %5.1f%%\n", stage, stage_name, stageTime / size, 100.0 * fracTime, flops, 100.0 * fracFlops, mess, 100.0 * fracMessages, avgMessLen, 100.0 * fracLength, red, 100.0 * fracReductions));
1447:     }
1448:   }

1450:   PetscCall(PetscFPrintf(comm, fd, "\n------------------------------------------------------------------------------------------------------------------------\n"));
1451:   PetscCall(PetscFPrintf(comm, fd, "See the 'Profiling' chapter of the users' manual for details on interpreting output.\n"));
1452:   PetscCall(PetscFPrintf(comm, fd, "Phase summary info:\n"));
1453:   PetscCall(PetscFPrintf(comm, fd, "   Count: number of times phase was executed\n"));
1454:   PetscCall(PetscFPrintf(comm, fd, "   Time and Flop: Max - maximum over all processors\n"));
1455:   PetscCall(PetscFPrintf(comm, fd, "                  Ratio - ratio of maximum to minimum over all processors\n"));
1456:   PetscCall(PetscFPrintf(comm, fd, "   Mess: number of messages sent\n"));
1457:   PetscCall(PetscFPrintf(comm, fd, "   AvgLen: average message length (bytes)\n"));
1458:   PetscCall(PetscFPrintf(comm, fd, "   Reduct: number of global reductions\n"));
1459:   PetscCall(PetscFPrintf(comm, fd, "   Global: entire computation\n"));
1460:   PetscCall(PetscFPrintf(comm, fd, "   Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().\n"));
1461:   PetscCall(PetscFPrintf(comm, fd, "      %%T - percent time in this phase         %%F - percent flop in this phase\n"));
1462:   PetscCall(PetscFPrintf(comm, fd, "      %%M - percent messages in this phase     %%L - percent message lengths in this phase\n"));
1463:   PetscCall(PetscFPrintf(comm, fd, "      %%R - percent reductions in this phase\n"));
1464:   PetscCall(PetscFPrintf(comm, fd, "   Total Mflop/s: 10e-6 * (sum of flop over all processors)/(max time over all processors)\n"));
1465:   if (PetscLogMemory) {
1466:     PetscCall(PetscFPrintf(comm, fd, "   Memory usage is summed over all MPI processes, it is given in mega-bytes\n"));
1467:     PetscCall(PetscFPrintf(comm, fd, "   Malloc Mbytes: Memory allocated and kept during event (sum over all calls to event). May be negative\n"));
1468:     PetscCall(PetscFPrintf(comm, fd, "   EMalloc Mbytes: extra memory allocated during event and then freed (maximum over all calls to events). Never negative\n"));
1469:     PetscCall(PetscFPrintf(comm, fd, "   MMalloc Mbytes: Increase in high water mark of allocated memory (sum over all calls to event). Never negative\n"));
1470:     PetscCall(PetscFPrintf(comm, fd, "   RMI Mbytes: Increase in resident memory (sum over all calls to event)\n"));
1471:   }
1472: #if defined(PETSC_HAVE_DEVICE)
1473:   PetscCall(PetscFPrintf(comm, fd, "   GPU Mflop/s: 10e-6 * (sum of flop on GPU over all processors)/(max GPU time over all processors)\n"));
1474:   PetscCall(PetscFPrintf(comm, fd, "   CpuToGpu Count: total number of CPU to GPU copies per processor\n"));
1475:   PetscCall(PetscFPrintf(comm, fd, "   CpuToGpu Size (Mbytes): 10e-6 * (total size of CPU to GPU copies per processor)\n"));
1476:   PetscCall(PetscFPrintf(comm, fd, "   GpuToCpu Count: total number of GPU to CPU copies per processor\n"));
1477:   PetscCall(PetscFPrintf(comm, fd, "   GpuToCpu Size (Mbytes): 10e-6 * (total size of GPU to CPU copies per processor)\n"));
1478:   PetscCall(PetscFPrintf(comm, fd, "   GPU %%F: percent flops on GPU in this event\n"));
1479: #endif
1480:   PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------\n"));

1482:   PetscCall(PetscLogViewWarnDebugging(comm, fd));

1484:   /* Report events */
1485:   PetscCall(PetscFPrintf(comm, fd, "Event                Count      Time (sec)     Flop                              --- Global ---  --- Stage ----  Total"));
1486:   if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, "  Malloc EMalloc MMalloc RMI"));
1487: #if defined(PETSC_HAVE_DEVICE)
1488:   PetscCall(PetscFPrintf(comm, fd, "   GPU    - CpuToGpu -   - GpuToCpu - GPU"));
1489: #endif
1490:   PetscCall(PetscFPrintf(comm, fd, "\n"));
1491:   PetscCall(PetscFPrintf(comm, fd, "                   Max Ratio  Max     Ratio   Max  Ratio  Mess   AvgLen  Reduct  %%T %%F %%M %%L %%R  %%T %%F %%M %%L %%R Mflop/s"));
1492:   if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, " Mbytes Mbytes Mbytes Mbytes"));
1493: #if defined(PETSC_HAVE_DEVICE)
1494:   PetscCall(PetscFPrintf(comm, fd, " Mflop/s Count   Size   Count   Size  %%F"));
1495: #endif
1496:   PetscCall(PetscFPrintf(comm, fd, "\n"));
1497:   PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------"));
1498:   if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, "-----------------------------"));
1499: #if defined(PETSC_HAVE_DEVICE)
1500:   PetscCall(PetscFPrintf(comm, fd, "---------------------------------------"));
1501: #endif
1502:   PetscCall(PetscFPrintf(comm, fd, "\n"));

1504: #if defined(PETSC_HAVE_DEVICE)
1505:   /* this indirect way of accessing these values is needed when PETSc is build with multiple libraries since the symbols are not in libpetscsys */
1506:   PetscCall(PetscLogStateGetEventFromName(state, "TaoSolve", &TAO_Solve));
1507:   PetscCall(PetscLogStateGetEventFromName(state, "TSStep", &TS_Step));
1508:   PetscCall(PetscLogStateGetEventFromName(state, "SNESSolve", &SNES_Solve));
1509:   PetscCall(PetscLogStateGetEventFromName(state, "KSPSolve", &KSP_Solve));
1510: #endif

1512:   for (stage = 0; stage < numStages; stage++) {
1513:     PetscInt            stage_id;
1514:     PetscEventPerfInfo *stage_info;
1515:     const char         *stage_name;

1517:     if (!(stageVisible[stage] && stageUsed[stage])) continue;
1518:     PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1519:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1520:     PetscCall(PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stage_name));
1521:     stage_info = &zero_info;
1522:     if (localStageUsed[stage]) {
1523:       PetscStagePerf *stage_perf_info;

1525:       PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_perf_info));
1526:       stage_info = &stage_perf_info->perfInfo;
1527:     }
1528:     PetscCall(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1529:     PetscCall(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1530:     PetscCall(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1531:     PetscCall(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1532:     PetscCall(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1533:     mess *= 0.5;
1534:     messLen *= 0.5;
1535:     red /= size;

1537:     for (event = 0; event < numEvents; event++) {
1538:       PetscInt            event_id;
1539:       PetscEventPerfInfo *event_info = &zero_info;
1540:       PetscBool           is_zero    = PETSC_FALSE;
1541:       const char         *event_name;

1543:       PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1544:       PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1545:       if (event_id >= 0 && stage_id >= 0) { PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &event_info)); }
1546:       PetscCall(PetscMemcmp(event_info, &zero_info, sizeof(zero_info), &is_zero));
1547:       PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1548:       if (!is_zero) {
1549:         flopr = event_info->flops;
1550:         PetscCall(MPIU_Allreduce(&flopr, &minf, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1551:         PetscCall(MPIU_Allreduce(&flopr, &maxf, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1552:         PetscCall(MPIU_Allreduce(&event_info->flops, &totf, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1553:         PetscCall(MPIU_Allreduce(&event_info->time, &mint, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1554:         PetscCall(MPIU_Allreduce(&event_info->time, &maxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1555:         PetscCall(MPIU_Allreduce(&event_info->time, &tott, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1556:         PetscCall(MPIU_Allreduce(&event_info->numMessages, &totm, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1557:         PetscCall(MPIU_Allreduce(&event_info->messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1558:         PetscCall(MPIU_Allreduce(&event_info->numReductions, &totr, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1559:         PetscCall(MPIU_Allreduce(&event_info->count, &minC, 1, MPI_INT, MPI_MIN, comm));
1560:         PetscCall(MPIU_Allreduce(&event_info->count, &maxC, 1, MPI_INT, MPI_MAX, comm));
1561:         if (PetscLogMemory) {
1562:           PetscCall(MPIU_Allreduce(&event_info->memIncrease, &mem, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1563:           PetscCall(MPIU_Allreduce(&event_info->mallocSpace, &mal, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1564:           PetscCall(MPIU_Allreduce(&event_info->mallocIncrease, &malmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1565:           PetscCall(MPIU_Allreduce(&event_info->mallocIncreaseEvent, &emalmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1566:         }
1567: #if defined(PETSC_HAVE_DEVICE)
1568:         PetscCall(MPIU_Allreduce(&event_info->CpuToGpuCount, &cct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1569:         PetscCall(MPIU_Allreduce(&event_info->GpuToCpuCount, &gct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1570:         PetscCall(MPIU_Allreduce(&event_info->CpuToGpuSize, &csz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1571:         PetscCall(MPIU_Allreduce(&event_info->GpuToCpuSize, &gsz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1572:         PetscCall(MPIU_Allreduce(&event_info->GpuFlops, &gflops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1573:         PetscCall(MPIU_Allreduce(&event_info->GpuTime, &gmaxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1574: #endif
1575:         if (mint < 0.0) {
1576:           PetscCall(PetscFPrintf(comm, fd, "WARNING!!! Minimum time %g over all processors for %s is negative! This happens\n on some machines whose times cannot handle too rapid calls.!\n artificially changing minimum to zero.\n", mint, event_name));
1577:           mint = 0;
1578:         }
1579:         PetscCheck(minf >= 0.0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Minimum flop %g over all processors for %s is negative! Not possible!", minf, event_name);
1580: #if defined(PETSC_HAVE_DEVICE)
1581:         /* Put NaN into the time for all events that may not be time accurately since they may happen asynchronously on the GPU */
1582:         if (!PetscLogGpuTimeFlag && petsc_gflops > 0) {
1583:           memcpy(&gmaxt, &nas, sizeof(PetscLogDouble));
1584:           if (event_id != SNES_Solve && event_id != KSP_Solve && event_id != TS_Step && event_id != TAO_Solve) {
1585:             memcpy(&mint, &nas, sizeof(PetscLogDouble));
1586:             memcpy(&maxt, &nas, sizeof(PetscLogDouble));
1587:           }
1588:         }
1589: #endif
1590:         totm *= 0.5;
1591:         totml *= 0.5;
1592:         totr /= size;

1594:         if (maxC != 0) {
1595:           if (minC != 0) ratC = ((PetscLogDouble)maxC) / minC;
1596:           else ratC = 0.0;
1597:           if (mint != 0.0) ratt = maxt / mint;
1598:           else ratt = 0.0;
1599:           if (minf != 0.0) ratf = maxf / minf;
1600:           else ratf = 0.0;
1601:           if (TotalTime != 0.0) fracTime = tott / TotalTime;
1602:           else fracTime = 0.0;
1603:           if (TotalFlops != 0.0) fracFlops = totf / TotalFlops;
1604:           else fracFlops = 0.0;
1605:           if (stageTime != 0.0) fracStageTime = tott / stageTime;
1606:           else fracStageTime = 0.0;
1607:           if (flops != 0.0) fracStageFlops = totf / flops;
1608:           else fracStageFlops = 0.0;
1609:           if (numMessages != 0.0) fracMess = totm / numMessages;
1610:           else fracMess = 0.0;
1611:           if (messageLength != 0.0) fracMessLen = totml / messageLength;
1612:           else fracMessLen = 0.0;
1613:           if (numReductions != 0.0) fracRed = totr / numReductions;
1614:           else fracRed = 0.0;
1615:           if (mess != 0.0) fracStageMess = totm / mess;
1616:           else fracStageMess = 0.0;
1617:           if (messLen != 0.0) fracStageMessLen = totml / messLen;
1618:           else fracStageMessLen = 0.0;
1619:           if (red != 0.0) fracStageRed = totr / red;
1620:           else fracStageRed = 0.0;
1621:           if (totm != 0.0) totml /= totm;
1622:           else totml = 0.0;
1623:           if (maxt != 0.0) flopr = totf / maxt;
1624:           else flopr = 0.0;
1625:           if (fracStageTime > 1.0 || fracStageFlops > 1.0 || fracStageMess > 1.0 || fracStageMessLen > 1.0 || fracStageRed > 1.0)
1626:             PetscCall(PetscFPrintf(comm, fd, "%-16s %7d %3.1f %5.4e %3.1f %3.2e %3.1f %2.1e %2.1e %2.1e %2.0f %2.0f %2.0f %2.0f %2.0f Multiple stages %5.0f", event_name, maxC, ratC, maxt, ratt, maxf, ratf, totm, totml, totr, 100.0 * fracTime, 100.0 * fracFlops, 100.0 * fracMess, 100.0 * fracMessLen, 100.0 * fracRed, PetscAbs(flopr) / 1.0e6));
1627:           else
1628:             PetscCall(PetscFPrintf(comm, fd, "%-16s %7d %3.1f %5.4e %3.1f %3.2e %3.1f %2.1e %2.1e %2.1e %2.0f %2.0f %2.0f %2.0f %2.0f %3.0f %2.0f %2.0f %2.0f %2.0f %5.0f", event_name, maxC, ratC, maxt, ratt, maxf, ratf, totm, totml, totr, 100.0 * fracTime, 100.0 * fracFlops, 100.0 * fracMess, 100.0 * fracMessLen, 100.0 * fracRed, 100.0 * fracStageTime, 100.0 * fracStageFlops, 100.0 * fracStageMess, 100.0 * fracStageMessLen, 100.0 * fracStageRed, PetscAbs(flopr) / 1.0e6));
1629:           if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, " %5.0f   %5.0f   %5.0f   %5.0f", mal / 1.0e6, emalmax / 1.0e6, malmax / 1.0e6, mem / 1.0e6));
1630: #if defined(PETSC_HAVE_DEVICE)
1631:           if (totf != 0.0) fracgflops = gflops / totf;
1632:           else fracgflops = 0.0;
1633:           if (gmaxt != 0.0) gflopr = gflops / gmaxt;
1634:           else gflopr = 0.0;
1635:           PetscCall(PetscFPrintf(comm, fd, "   %5.0f   %4.0f %3.2e %4.0f %3.2e % 2.0f", PetscAbs(gflopr) / 1.0e6, cct / size, csz / (1.0e6 * size), gct / size, gsz / (1.0e6 * size), 100.0 * fracgflops));
1636: #endif
1637:           PetscCall(PetscFPrintf(comm, fd, "\n"));
1638:         }
1639:       }
1640:     }
1641:   }

1643:   /* Memory usage and object creation */
1644:   PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------"));
1645:   if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, "-----------------------------"));
1646: #if defined(PETSC_HAVE_DEVICE)
1647:   PetscCall(PetscFPrintf(comm, fd, "---------------------------------------"));
1648: #endif
1649:   PetscCall(PetscFPrintf(comm, fd, "\n"));
1650:   PetscCall(PetscFPrintf(comm, fd, "\n"));

1652:   /* Right now, only stages on the first processor are reported here, meaning only objects associated with
1653:      the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then
1654:      stats for stages local to processor sets.
1655:   */
1656:   /* We should figure out the longest object name here (now 20 characters) */
1657:   PetscCall(PetscFPrintf(comm, fd, "Object Type          Creations   Destructions. Reports information only for process 0.\n"));
1658:   for (stage = 0; stage < numStages; stage++) {
1659:     const char *stage_name;

1661:     PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1662:     PetscCall(PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stage_name));
1663:     if (localStageUsed[stage]) {
1664:       PetscInt num_classes;

1666:       PetscCall(PetscLogStateGetNumClasses(state, &num_classes));
1667:       for (oclass = 0; oclass < num_classes; oclass++) {
1668:         PetscClassPerf *class_perf_info;

1670:         PetscCall(PetscLogHandlerDefaultGetClassPerf(handler, stage, oclass, &class_perf_info));
1671:         if ((class_perf_info->creations > 0) || (class_perf_info->destructions > 0)) {
1672:           PetscLogClassInfo class_reg_info;
1673:           PetscBool         flg;

1675:           PetscCall(PetscLogStateClassGetInfo(state, oclass, &class_reg_info));
1676:           if (stage == 0 && oclass == num_classes - 1) {
1677:             PetscCall(PetscStrcmp(class_reg_info.name, "Viewer", &flg));
1678:             PetscCheck(flg && class_perf_info->creations == 1 && class_perf_info->destructions == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "The last PetscObject type of the main PetscLogStage should be PetscViewer with a single creation and no destruction");
1679:           } else PetscCall(PetscFPrintf(comm, fd, "%20s %5d          %5d\n", class_reg_info.name, class_perf_info->creations, class_perf_info->destructions));
1680:         }
1681:       }
1682:     }
1683:   }

1685:   PetscCall(PetscFree(localStageUsed));
1686:   PetscCall(PetscFree(stageUsed));
1687:   PetscCall(PetscFree(localStageVisible));
1688:   PetscCall(PetscFree(stageVisible));
1689:   PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1690:   PetscCall(PetscLogGlobalNamesDestroy(&global_events));

1692:   /* Information unrelated to this particular run */
1693:   PetscCall(PetscFPrintf(comm, fd, "========================================================================================================================\n"));
1694:   PetscCall(PetscTime(&y));
1695:   PetscCall(PetscTime(&x));
1696:   PetscCall(PetscTime(&y));
1697:   PetscCall(PetscTime(&y));
1698:   PetscCall(PetscTime(&y));
1699:   PetscCall(PetscTime(&y));
1700:   PetscCall(PetscTime(&y));
1701:   PetscCall(PetscTime(&y));
1702:   PetscCall(PetscTime(&y));
1703:   PetscCall(PetscTime(&y));
1704:   PetscCall(PetscTime(&y));
1705:   PetscCall(PetscTime(&y));
1706:   PetscCall(PetscFPrintf(comm, fd, "Average time to get PetscTime(): %g\n", (y - x) / 10.0));
1707:   /* MPI information */
1708:   if (size > 1) {
1709:     MPI_Status  status;
1710:     PetscMPIInt tag;
1711:     MPI_Comm    newcomm;

1713:     PetscCallMPI(MPI_Barrier(comm));
1714:     PetscCall(PetscTime(&x));
1715:     PetscCallMPI(MPI_Barrier(comm));
1716:     PetscCallMPI(MPI_Barrier(comm));
1717:     PetscCallMPI(MPI_Barrier(comm));
1718:     PetscCallMPI(MPI_Barrier(comm));
1719:     PetscCallMPI(MPI_Barrier(comm));
1720:     PetscCall(PetscTime(&y));
1721:     PetscCall(PetscFPrintf(comm, fd, "Average time for MPI_Barrier(): %g\n", (y - x) / 5.0));
1722:     PetscCall(PetscCommDuplicate(comm, &newcomm, &tag));
1723:     PetscCallMPI(MPI_Barrier(comm));
1724:     if (rank) {
1725:       PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, rank - 1, tag, newcomm, &status));
1726:       PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, (rank + 1) % size, tag, newcomm));
1727:     } else {
1728:       PetscCall(PetscTime(&x));
1729:       PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, 1, tag, newcomm));
1730:       PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, size - 1, tag, newcomm, &status));
1731:       PetscCall(PetscTime(&y));
1732:       PetscCall(PetscFPrintf(comm, fd, "Average time for zero size MPI_Send(): %g\n", (y - x) / size));
1733:     }
1734:     PetscCall(PetscCommDestroy(&newcomm));
1735:   }
1736:   PetscCall(PetscOptionsView(NULL, viewer));

1738:   /* Machine and compile information */
1739:   if (PetscDefined(USE_FORTRAN_KERNELS)) {
1740:     PetscCall(PetscFPrintf(comm, fd, "Compiled with FORTRAN kernels\n"));
1741:   } else {
1742:     PetscCall(PetscFPrintf(comm, fd, "Compiled without FORTRAN kernels\n"));
1743:   }
1744:   if (PetscDefined(USE_64BIT_INDICES)) {
1745:     PetscCall(PetscFPrintf(comm, fd, "Compiled with 64-bit PetscInt\n"));
1746:   } else if (PetscDefined(USE___FLOAT128)) {
1747:     PetscCall(PetscFPrintf(comm, fd, "Compiled with 32-bit PetscInt\n"));
1748:   }
1749:   if (PetscDefined(USE_REAL_SINGLE)) {
1750:     PetscCall(PetscFPrintf(comm, fd, "Compiled with single precision PetscScalar and PetscReal\n"));
1751:   } else if (PetscDefined(USE___FLOAT128)) {
1752:     PetscCall(PetscFPrintf(comm, fd, "Compiled with 128 bit precision PetscScalar and PetscReal\n"));
1753:   }
1754:   if (PetscDefined(USE_REAL_MAT_SINGLE)) {
1755:     PetscCall(PetscFPrintf(comm, fd, "Compiled with single precision matrices\n"));
1756:   } else {
1757:     PetscCall(PetscFPrintf(comm, fd, "Compiled with full precision matrices (default)\n"));
1758:   }
1759:   PetscCall(PetscFPrintf(comm, fd, "sizeof(short) %d sizeof(int) %d sizeof(long) %d sizeof(void*) %d sizeof(PetscScalar) %d sizeof(PetscInt) %d\n", (int)sizeof(short), (int)sizeof(int), (int)sizeof(long), (int)sizeof(void *), (int)sizeof(PetscScalar), (int)sizeof(PetscInt)));

1761:   PetscCall(PetscFPrintf(comm, fd, "Configure options: %s", petscconfigureoptions));
1762:   PetscCall(PetscFPrintf(comm, fd, "%s", petscmachineinfo));
1763:   PetscCall(PetscFPrintf(comm, fd, "%s", petsccompilerinfo));
1764:   PetscCall(PetscFPrintf(comm, fd, "%s", petsccompilerflagsinfo));
1765:   PetscCall(PetscFPrintf(comm, fd, "%s", petsclinkerinfo));

1767:   /* Cleanup */
1768:   PetscCall(PetscFPrintf(comm, fd, "\n"));
1769:   PetscCall(PetscLogViewWarnNoGpuAwareMpi(comm, fd));
1770:   PetscCall(PetscLogViewWarnDebugging(comm, fd));
1771:   PetscCall(PetscFPTrapPop());
1772:   PetscFunctionReturn(PETSC_SUCCESS);
1773: }

1775: static PetscErrorCode PetscLogHandlerView_Default(PetscLogHandler handler, PetscViewer viewer)
1776: {
1777:   PetscViewerFormat format;
1778:   PetscFunctionBegin;
1779:   PetscCall(PetscViewerGetFormat(viewer, &format));
1780:   if (format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_ASCII_INFO) {
1781:     PetscCall(PetscLogHandlerView_Default_Info(handler, viewer));
1782:   } else if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1783:     PetscCall(PetscLogHandlerView_Default_Detailed(handler, viewer));
1784:   } else if (format == PETSC_VIEWER_ASCII_CSV) {
1785:     PetscCall(PetscLogHandlerView_Default_CSV(handler, viewer));
1786:   }
1787:   PetscFunctionReturn(PETSC_SUCCESS);
1788: }

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

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

1799:   Level: developer

1801: .seealso: [](ch_profiling), `PetscLogHandler`
1802: M*/

1804: PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Default(PetscLogHandler handler)
1805: {
1806:   PetscFunctionBegin;
1807:   PetscCall(PetscLogHandlerContextCreate_Default((PetscLogHandler_Default *)&handler->data));
1808:   handler->ops->destroy       = PetscLogHandlerDestroy_Default;
1809:   handler->ops->eventbegin    = PetscLogHandlerEventBegin_Default;
1810:   handler->ops->eventend      = PetscLogHandlerEventEnd_Default;
1811:   handler->ops->eventsync     = PetscLogHandlerEventSync_Default;
1812:   handler->ops->objectcreate  = PetscLogHandlerObjectCreate_Default;
1813:   handler->ops->objectdestroy = PetscLogHandlerObjectDestroy_Default;
1814:   handler->ops->stagepush     = PetscLogHandlerStagePush_Default;
1815:   handler->ops->stagepop      = PetscLogHandlerStagePop_Default;
1816:   handler->ops->view          = PetscLogHandlerView_Default;
1817:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetEventPerfInfo_C", PetscLogHandlerGetEventPerfInfo_Default));
1818:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetStagePerfInfo_C", PetscLogHandlerGetStagePerfInfo_Default));
1819:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerSetLogActions_C", PetscLogHandlerSetLogActions_Default));
1820:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerSetLogObjects_C", PetscLogHandlerSetLogObjects_Default));
1821:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerLogObjectState_C", PetscLogHandlerLogObjectState_Default));
1822:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetNumObjects_C", PetscLogHandlerGetNumObjects_Default));
1823:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventDeactivatePush_C", PetscLogHandlerEventDeactivatePush_Default));
1824:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventDeactivatePop_C", PetscLogHandlerEventDeactivatePop_Default));
1825:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventsPause_C", PetscLogHandlerEventsPause_Default));
1826:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventsResume_C", PetscLogHandlerEventsResume_Default));
1827:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerDump_C", PetscLogHandlerDump_Default));
1828:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerStageSetVisible_C", PetscLogHandlerStageSetVisible_Default));
1829:   PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerStageGetVisible_C", PetscLogHandlerStageGetVisible_Default));
1830:   PetscFunctionReturn(PETSC_SUCCESS);
1831: }