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