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