xref: /petsc/src/sys/logging/plog.c (revision a16fd2c93c02146fccd68469496ac02ca99b9ebe)
1 
2 /*
3       PETSc code to log object creation and destruction and PETSc events.
4 
5       This provides the public API used by the rest of PETSc and by users.
6 
7       These routines use a private API that is not used elsewhere in PETSc and is not
8       accessible to users. The private API is defined in logimpl.h and the utils directory.
9 
10 */
11 #include <petsc/private/logimpl.h> /*I    "petscsys.h"   I*/
12 #include <petsctime.h>
13 #include <petscviewer.h>
14 
15 /* this is not consistently used and is difficult to keep correct if done manually, possibly it should be removed */
16 PetscErrorCode PetscLogObjectParent(PetscObject p, PetscObject c) {
17   if (!c || !p) return 0;
18   c->parent   = p;
19   c->parentid = p->id;
20   return 0;
21 }
22 
23 /*@C
24    PetscLogObjectMemory - Adds to an object a count of additional amount of memory that is used by the object.
25 
26    Not collective.
27 
28    Input Parameters:
29 +  obj  - the PETSc object
30 -  mem  - the amount of memory that is being added to the object
31 
32    Level: developer
33 
34    Developer Notes:
35    This is not used consistently. It is very difficult to manually track the memory usage per object so this should
36    likely be removed and replaced with an automated system.
37 
38 .seealso: `PetscFinalize()`, `PetscInitializeFortran()`, `PetscGetArgs()`, `PetscInitializeNoArguments()`
39 
40 @*/
41 PetscErrorCode PetscLogObjectMemory(PetscObject p, PetscLogDouble m) {
42   if (!p) return 0;
43   p->mem += m;
44   return 0;
45 }
46 
47 PetscLogEvent PETSC_LARGEST_EVENT = PETSC_EVENT;
48 
49 #if defined(PETSC_USE_LOG)
50 #include <petscmachineinfo.h>
51 #include <petscconfiginfo.h>
52 
53 /* used in the MPI_XXX() count macros in petsclog.h */
54 
55 /* Action and object logging variables */
56 Action   *petsc_actions    = NULL;
57 Object   *petsc_objects    = NULL;
58 PetscBool petsc_logActions = PETSC_FALSE;
59 PetscBool petsc_logObjects = PETSC_FALSE;
60 int       petsc_numActions = 0, petsc_maxActions = 100;
61 int       petsc_numObjects = 0, petsc_maxObjects = 100;
62 int       petsc_numObjectsDestroyed = 0;
63 
64 /* Global counters */
65 PetscLogDouble petsc_BaseTime        = 0.0;
66 PetscLogDouble petsc_TotalFlops      = 0.0; /* The number of flops */
67 PetscLogDouble petsc_tmp_flops       = 0.0; /* The incremental number of flops */
68 PetscLogDouble petsc_send_ct         = 0.0; /* The number of sends */
69 PetscLogDouble petsc_recv_ct         = 0.0; /* The number of receives */
70 PetscLogDouble petsc_send_len        = 0.0; /* The total length of all sent messages */
71 PetscLogDouble petsc_recv_len        = 0.0; /* The total length of all received messages */
72 PetscLogDouble petsc_isend_ct        = 0.0; /* The number of immediate sends */
73 PetscLogDouble petsc_irecv_ct        = 0.0; /* The number of immediate receives */
74 PetscLogDouble petsc_isend_len       = 0.0; /* The total length of all immediate send messages */
75 PetscLogDouble petsc_irecv_len       = 0.0; /* The total length of all immediate receive messages */
76 PetscLogDouble petsc_wait_ct         = 0.0; /* The number of waits */
77 PetscLogDouble petsc_wait_any_ct     = 0.0; /* The number of anywaits */
78 PetscLogDouble petsc_wait_all_ct     = 0.0; /* The number of waitalls */
79 PetscLogDouble petsc_sum_of_waits_ct = 0.0; /* The total number of waits */
80 PetscLogDouble petsc_allreduce_ct    = 0.0; /* The number of reductions */
81 PetscLogDouble petsc_gather_ct       = 0.0; /* The number of gathers and gathervs */
82 PetscLogDouble petsc_scatter_ct      = 0.0; /* The number of scatters and scattervs */
83 #if defined(PETSC_HAVE_DEVICE)
84 PetscLogDouble petsc_ctog_ct        = 0.0; /* The total number of CPU to GPU copies */
85 PetscLogDouble petsc_gtoc_ct        = 0.0; /* The total number of GPU to CPU copies */
86 PetscLogDouble petsc_ctog_sz        = 0.0; /* The total size of CPU to GPU copies */
87 PetscLogDouble petsc_gtoc_sz        = 0.0; /* The total size of GPU to CPU copies */
88 PetscLogDouble petsc_ctog_ct_scalar = 0.0; /* The total number of CPU to GPU copies */
89 PetscLogDouble petsc_gtoc_ct_scalar = 0.0; /* The total number of GPU to CPU copies */
90 PetscLogDouble petsc_ctog_sz_scalar = 0.0; /* The total size of CPU to GPU copies */
91 PetscLogDouble petsc_gtoc_sz_scalar = 0.0; /* The total size of GPU to CPU copies */
92 PetscLogDouble petsc_gflops         = 0.0; /* The flops done on a GPU */
93 PetscLogDouble petsc_gtime          = 0.0; /* The time spent on a GPU */
94 #endif
95 
96 /* Logging functions */
97 PetscErrorCode (*PetscLogPHC)(PetscObject)                                                            = NULL;
98 PetscErrorCode (*PetscLogPHD)(PetscObject)                                                            = NULL;
99 PetscErrorCode (*PetscLogPLB)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject) = NULL;
100 PetscErrorCode (*PetscLogPLE)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject) = NULL;
101 
102 /* Tracing event logging variables */
103 FILE            *petsc_tracefile          = NULL;
104 int              petsc_tracelevel         = 0;
105 const char      *petsc_traceblanks        = "                                                                                                    ";
106 char             petsc_tracespace[128]    = " ";
107 PetscLogDouble   petsc_tracetime          = 0.0;
108 static PetscBool PetscLogInitializeCalled = PETSC_FALSE;
109 
110 static PetscIntStack current_log_event_stack = NULL;
111 
112 PETSC_INTERN PetscErrorCode PetscLogInitialize(void) {
113   int       stage;
114   PetscBool opt;
115 
116   PetscFunctionBegin;
117   if (PetscLogInitializeCalled) PetscFunctionReturn(0);
118   PetscLogInitializeCalled = PETSC_TRUE;
119 
120   PetscCall(PetscIntStackCreate(&current_log_event_stack));
121   PetscCall(PetscOptionsHasName(NULL, NULL, "-log_exclude_actions", &opt));
122   if (opt) petsc_logActions = PETSC_FALSE;
123   PetscCall(PetscOptionsHasName(NULL, NULL, "-log_exclude_objects", &opt));
124   if (opt) petsc_logObjects = PETSC_FALSE;
125   if (petsc_logActions) PetscCall(PetscMalloc1(petsc_maxActions, &petsc_actions));
126   if (petsc_logObjects) PetscCall(PetscMalloc1(petsc_maxObjects, &petsc_objects));
127   PetscLogPHC = PetscLogObjCreateDefault;
128   PetscLogPHD = PetscLogObjDestroyDefault;
129   /* Setup default logging structures */
130   PetscCall(PetscStageLogCreate(&petsc_stageLog));
131   PetscCall(PetscStageLogRegister(petsc_stageLog, "Main Stage", &stage));
132 
133   /* All processors sync here for more consistent logging */
134   PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
135   PetscTime(&petsc_BaseTime);
136   PetscCall(PetscLogStagePush(stage));
137   PetscFunctionReturn(0);
138 }
139 
140 PETSC_INTERN PetscErrorCode PetscLogFinalize(void) {
141   PetscStageLog stageLog;
142 
143   PetscFunctionBegin;
144   PetscCall(PetscFree(petsc_actions));
145   PetscCall(PetscFree(petsc_objects));
146   PetscCall(PetscLogNestedEnd());
147   PetscCall(PetscLogSet(NULL, NULL));
148 
149   /* Resetting phase */
150   PetscCall(PetscLogGetStageLog(&stageLog));
151   PetscCall(PetscStageLogDestroy(stageLog));
152   PetscCall(PetscIntStackDestroy(current_log_event_stack));
153   current_log_event_stack = NULL;
154 
155   petsc_TotalFlops          = 0.0;
156   petsc_numActions          = 0;
157   petsc_numObjects          = 0;
158   petsc_numObjectsDestroyed = 0;
159   petsc_maxActions          = 100;
160   petsc_maxObjects          = 100;
161   petsc_actions             = NULL;
162   petsc_objects             = NULL;
163   petsc_logActions          = PETSC_FALSE;
164   petsc_logObjects          = PETSC_FALSE;
165   petsc_BaseTime            = 0.0;
166   petsc_TotalFlops          = 0.0;
167   petsc_tmp_flops           = 0.0;
168   petsc_send_ct             = 0.0;
169   petsc_recv_ct             = 0.0;
170   petsc_send_len            = 0.0;
171   petsc_recv_len            = 0.0;
172   petsc_isend_ct            = 0.0;
173   petsc_irecv_ct            = 0.0;
174   petsc_isend_len           = 0.0;
175   petsc_irecv_len           = 0.0;
176   petsc_wait_ct             = 0.0;
177   petsc_wait_any_ct         = 0.0;
178   petsc_wait_all_ct         = 0.0;
179   petsc_sum_of_waits_ct     = 0.0;
180   petsc_allreduce_ct        = 0.0;
181   petsc_gather_ct           = 0.0;
182   petsc_scatter_ct          = 0.0;
183 #if defined(PETSC_HAVE_DEVICE)
184   petsc_ctog_ct = 0.0;
185   petsc_gtoc_ct = 0.0;
186   petsc_ctog_sz = 0.0;
187   petsc_gtoc_sz = 0.0;
188   petsc_gflops  = 0.0;
189   petsc_gtime   = 0.0;
190 #endif
191   PETSC_LARGEST_EVENT      = PETSC_EVENT;
192   PetscLogPHC              = NULL;
193   PetscLogPHD              = NULL;
194   petsc_tracefile          = NULL;
195   petsc_tracelevel         = 0;
196   petsc_traceblanks        = "                                                                                                    ";
197   petsc_tracespace[0]      = ' ';
198   petsc_tracespace[1]      = 0;
199   petsc_tracetime          = 0.0;
200   PETSC_LARGEST_CLASSID    = PETSC_SMALLEST_CLASSID;
201   PETSC_OBJECT_CLASSID     = 0;
202   petsc_stageLog           = NULL;
203   PetscLogInitializeCalled = PETSC_FALSE;
204   PetscFunctionReturn(0);
205 }
206 
207 /*@C
208   PetscLogSet - Sets the logging functions called at the beginning and ending of every event.
209 
210   Not Collective
211 
212   Input Parameters:
213 + b - The function called at beginning of event
214 - e - The function called at end of event
215 
216   Level: developer
217 
218 .seealso: `PetscLogDump()`, `PetscLogDefaultBegin()`, `PetscLogAllBegin()`, `PetscLogTraceBegin()`
219 @*/
220 PetscErrorCode PetscLogSet(PetscErrorCode (*b)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject), PetscErrorCode (*e)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject)) {
221   PetscFunctionBegin;
222   PetscLogPLB = b;
223   PetscLogPLE = e;
224   PetscFunctionReturn(0);
225 }
226 
227 /*@C
228   PetscLogIsActive - Check if logging is currently in progress.
229 
230   Not Collective
231 
232   Output Parameter:
233 . isActive - PETSC_TRUE if logging is in progress, PETSC_FALSE otherwise
234 
235   Level: beginner
236 
237 .seealso: `PetscLogDefaultBegin()`, `PetscLogAllBegin()`, `PetscLogSet()`
238 @*/
239 PetscErrorCode PetscLogIsActive(PetscBool *isActive) {
240   PetscFunctionBegin;
241   *isActive = (PetscLogPLB && PetscLogPLE) ? PETSC_TRUE : PETSC_FALSE;
242   PetscFunctionReturn(0);
243 }
244 
245 /*@C
246   PetscLogDefaultBegin - Turns on logging of objects and events. This logs flop
247   rates and object creation and should not slow programs down too much.
248   This routine may be called more than once.
249 
250   Logically Collective over PETSC_COMM_WORLD
251 
252   Options Database Keys:
253 . -log_view [viewertype:filename:viewerformat] - Prints summary of flop and timing information to the
254                   screen (for code configured with --with-log=1 (which is the default))
255 
256   Usage:
257 .vb
258       PetscInitialize(...);
259       PetscLogDefaultBegin();
260        ... code ...
261       PetscLogView(viewer); or PetscLogDump();
262       PetscFinalize();
263 .ve
264 
265   Notes:
266   PetscLogView(viewer) or PetscLogDump() actually cause the printing of
267   the logging information.
268 
269   Level: advanced
270 
271 .seealso: `PetscLogDump()`, `PetscLogAllBegin()`, `PetscLogView()`, `PetscLogTraceBegin()`
272 @*/
273 PetscErrorCode PetscLogDefaultBegin(void) {
274   PetscFunctionBegin;
275   PetscCall(PetscLogSet(PetscLogEventBeginDefault, PetscLogEventEndDefault));
276   PetscFunctionReturn(0);
277 }
278 
279 /*@C
280   PetscLogAllBegin - Turns on extensive logging of objects and events. Logs
281   all events. This creates large log files and slows the program down.
282 
283   Logically Collective on PETSC_COMM_WORLD
284 
285   Options Database Keys:
286 . -log_all - Prints extensive log information
287 
288   Usage:
289 .vb
290      PetscInitialize(...);
291      PetscLogAllBegin();
292      ... code ...
293      PetscLogDump(filename);
294      PetscFinalize();
295 .ve
296 
297   Notes:
298   A related routine is PetscLogDefaultBegin() (with the options key -log), which is
299   intended for production runs since it logs only flop rates and object
300   creation (and shouldn't significantly slow the programs).
301 
302   Level: advanced
303 
304 .seealso: `PetscLogDump()`, `PetscLogDefaultBegin()`, `PetscLogTraceBegin()`
305 @*/
306 PetscErrorCode PetscLogAllBegin(void) {
307   PetscFunctionBegin;
308   PetscCall(PetscLogSet(PetscLogEventBeginComplete, PetscLogEventEndComplete));
309   PetscFunctionReturn(0);
310 }
311 
312 /*@C
313   PetscLogTraceBegin - Activates trace logging.  Every time a PETSc event
314   begins or ends, the event name is printed.
315 
316   Logically Collective on PETSC_COMM_WORLD
317 
318   Input Parameter:
319 . file - The file to print trace in (e.g. stdout)
320 
321   Options Database Key:
322 . -log_trace [filename] - Activates PetscLogTraceBegin()
323 
324   Notes:
325   PetscLogTraceBegin() prints the processor number, the execution time (sec),
326   then "Event begin:" or "Event end:" followed by the event name.
327 
328   PetscLogTraceBegin() allows tracing of all PETSc calls, which is useful
329   to determine where a program is hanging without running in the
330   debugger.  Can be used in conjunction with the -info option.
331 
332   Level: intermediate
333 
334 .seealso: `PetscLogDump()`, `PetscLogAllBegin()`, `PetscLogView()`, `PetscLogDefaultBegin()`
335 @*/
336 PetscErrorCode PetscLogTraceBegin(FILE *file) {
337   PetscFunctionBegin;
338   petsc_tracefile = file;
339 
340   PetscCall(PetscLogSet(PetscLogEventBeginTrace, PetscLogEventEndTrace));
341   PetscFunctionReturn(0);
342 }
343 
344 /*@
345   PetscLogActions - Determines whether actions are logged for the graphical viewer.
346 
347   Not Collective
348 
349   Input Parameter:
350 . flag - PETSC_TRUE if actions are to be logged
351 
352   Level: intermediate
353 
354   Note: Logging of actions continues to consume more memory as the program
355   runs. Long running programs should consider turning this feature off.
356 
357   Options Database Keys:
358 . -log_exclude_actions - Turns off actions logging
359 
360 .seealso: `PetscLogStagePush()`, `PetscLogStagePop()`
361 @*/
362 PetscErrorCode PetscLogActions(PetscBool flag) {
363   PetscFunctionBegin;
364   petsc_logActions = flag;
365   PetscFunctionReturn(0);
366 }
367 
368 /*@
369   PetscLogObjects - Determines whether objects are logged for the graphical viewer.
370 
371   Not Collective
372 
373   Input Parameter:
374 . flag - PETSC_TRUE if objects are to be logged
375 
376   Level: intermediate
377 
378   Note: Logging of objects continues to consume more memory as the program
379   runs. Long running programs should consider turning this feature off.
380 
381   Options Database Keys:
382 . -log_exclude_objects - Turns off objects logging
383 
384 .seealso: `PetscLogStagePush()`, `PetscLogStagePop()`
385 @*/
386 PetscErrorCode PetscLogObjects(PetscBool flag) {
387   PetscFunctionBegin;
388   petsc_logObjects = flag;
389   PetscFunctionReturn(0);
390 }
391 
392 /*------------------------------------------------ Stage Functions --------------------------------------------------*/
393 /*@C
394   PetscLogStageRegister - Attaches a character string name to a logging stage.
395 
396   Not Collective
397 
398   Input Parameter:
399 . sname - The name to associate with that stage
400 
401   Output Parameter:
402 . stage - The stage number
403 
404   Level: intermediate
405 
406 .seealso: `PetscLogStagePush()`, `PetscLogStagePop()`
407 @*/
408 PetscErrorCode PetscLogStageRegister(const char sname[], PetscLogStage *stage) {
409   PetscStageLog stageLog;
410   PetscLogEvent event;
411 
412   PetscFunctionBegin;
413   PetscCall(PetscLogGetStageLog(&stageLog));
414   PetscCall(PetscStageLogRegister(stageLog, sname, stage));
415   /* Copy events already changed in the main stage, this sucks */
416   PetscCall(PetscEventPerfLogEnsureSize(stageLog->stageInfo[*stage].eventLog, stageLog->eventLog->numEvents));
417   for (event = 0; event < stageLog->eventLog->numEvents; event++) PetscCall(PetscEventPerfInfoCopy(&stageLog->stageInfo[0].eventLog->eventInfo[event], &stageLog->stageInfo[*stage].eventLog->eventInfo[event]));
418   PetscCall(PetscClassPerfLogEnsureSize(stageLog->stageInfo[*stage].classLog, stageLog->classLog->numClasses));
419   PetscFunctionReturn(0);
420 }
421 
422 /*@C
423   PetscLogStagePush - This function pushes a stage on the stack.
424 
425   Not Collective
426 
427   Input Parameter:
428 . stage - The stage on which to log
429 
430   Usage:
431   If the option -log_sumary is used to run the program containing the
432   following code, then 2 sets of summary data will be printed during
433   PetscFinalize().
434 .vb
435       PetscInitialize(int *argc,char ***args,0,0);
436       [stage 0 of code]
437       PetscLogStagePush(1);
438       [stage 1 of code]
439       PetscLogStagePop();
440       PetscBarrier(...);
441       [more stage 0 of code]
442       PetscFinalize();
443 .ve
444 
445   Notes:
446   Use PetscLogStageRegister() to register a stage.
447 
448   Level: intermediate
449 
450 .seealso: `PetscLogStagePop()`, `PetscLogStageRegister()`, `PetscBarrier()`
451 @*/
452 PetscErrorCode PetscLogStagePush(PetscLogStage stage) {
453   PetscStageLog stageLog;
454 
455   PetscFunctionBegin;
456   PetscCall(PetscLogGetStageLog(&stageLog));
457   PetscCall(PetscStageLogPush(stageLog, stage));
458   PetscFunctionReturn(0);
459 }
460 
461 /*@C
462   PetscLogStagePop - This function pops a stage from the stack.
463 
464   Not Collective
465 
466   Usage:
467   If the option -log_sumary is used to run the program containing the
468   following code, then 2 sets of summary data will be printed during
469   PetscFinalize().
470 .vb
471       PetscInitialize(int *argc,char ***args,0,0);
472       [stage 0 of code]
473       PetscLogStagePush(1);
474       [stage 1 of code]
475       PetscLogStagePop();
476       PetscBarrier(...);
477       [more stage 0 of code]
478       PetscFinalize();
479 .ve
480 
481   Notes:
482   Use PetscLogStageRegister() to register a stage.
483 
484   Level: intermediate
485 
486 .seealso: `PetscLogStagePush()`, `PetscLogStageRegister()`, `PetscBarrier()`
487 @*/
488 PetscErrorCode PetscLogStagePop(void) {
489   PetscStageLog stageLog;
490 
491   PetscFunctionBegin;
492   PetscCall(PetscLogGetStageLog(&stageLog));
493   PetscCall(PetscStageLogPop(stageLog));
494   PetscFunctionReturn(0);
495 }
496 
497 /*@
498   PetscLogStageSetActive - Determines stage activity for PetscLogEventBegin() and PetscLogEventEnd().
499 
500   Not Collective
501 
502   Input Parameters:
503 + stage    - The stage
504 - isActive - The activity flag, PETSC_TRUE for logging, else PETSC_FALSE (defaults to PETSC_TRUE)
505 
506   Level: intermediate
507 
508 .seealso: `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`, `PetscPreLoadStage()`
509 @*/
510 PetscErrorCode PetscLogStageSetActive(PetscLogStage stage, PetscBool isActive) {
511   PetscStageLog stageLog;
512 
513   PetscFunctionBegin;
514   PetscCall(PetscLogGetStageLog(&stageLog));
515   PetscCall(PetscStageLogSetActive(stageLog, stage, isActive));
516   PetscFunctionReturn(0);
517 }
518 
519 /*@
520   PetscLogStageGetActive - Returns stage activity for PetscLogEventBegin() and PetscLogEventEnd().
521 
522   Not Collective
523 
524   Input Parameter:
525 . stage    - The stage
526 
527   Output Parameter:
528 . isActive - The activity flag, PETSC_TRUE for logging, else PETSC_FALSE (defaults to PETSC_TRUE)
529 
530   Level: intermediate
531 
532 .seealso: `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`, `PetscPreLoadStage()`
533 @*/
534 PetscErrorCode PetscLogStageGetActive(PetscLogStage stage, PetscBool *isActive) {
535   PetscStageLog stageLog;
536 
537   PetscFunctionBegin;
538   PetscCall(PetscLogGetStageLog(&stageLog));
539   PetscCall(PetscStageLogGetActive(stageLog, stage, isActive));
540   PetscFunctionReturn(0);
541 }
542 
543 /*@
544   PetscLogStageSetVisible - Determines stage visibility in PetscLogView()
545 
546   Not Collective
547 
548   Input Parameters:
549 + stage     - The stage
550 - isVisible - The visibility flag, PETSC_TRUE to print, else PETSC_FALSE (defaults to PETSC_TRUE)
551 
552   Level: intermediate
553 
554 .seealso: `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscLogView()`
555 @*/
556 PetscErrorCode PetscLogStageSetVisible(PetscLogStage stage, PetscBool isVisible) {
557   PetscStageLog stageLog;
558 
559   PetscFunctionBegin;
560   PetscCall(PetscLogGetStageLog(&stageLog));
561   PetscCall(PetscStageLogSetVisible(stageLog, stage, isVisible));
562   PetscFunctionReturn(0);
563 }
564 
565 /*@
566   PetscLogStageGetVisible - Returns stage visibility in PetscLogView()
567 
568   Not Collective
569 
570   Input Parameter:
571 . stage     - The stage
572 
573   Output Parameter:
574 . isVisible - The visibility flag, PETSC_TRUE to print, else PETSC_FALSE (defaults to PETSC_TRUE)
575 
576   Level: intermediate
577 
578 .seealso: `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscLogView()`
579 @*/
580 PetscErrorCode PetscLogStageGetVisible(PetscLogStage stage, PetscBool *isVisible) {
581   PetscStageLog stageLog;
582 
583   PetscFunctionBegin;
584   PetscCall(PetscLogGetStageLog(&stageLog));
585   PetscCall(PetscStageLogGetVisible(stageLog, stage, isVisible));
586   PetscFunctionReturn(0);
587 }
588 
589 /*@C
590   PetscLogStageGetId - Returns the stage id when given the stage name.
591 
592   Not Collective
593 
594   Input Parameter:
595 . name  - The stage name
596 
597   Output Parameter:
598 . stage - The stage, , or -1 if no stage with that name exists
599 
600   Level: intermediate
601 
602 .seealso: `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`, `PetscPreLoadStage()`
603 @*/
604 PetscErrorCode PetscLogStageGetId(const char name[], PetscLogStage *stage) {
605   PetscStageLog stageLog;
606 
607   PetscFunctionBegin;
608   PetscCall(PetscLogGetStageLog(&stageLog));
609   PetscCall(PetscStageLogGetStage(stageLog, name, stage));
610   PetscFunctionReturn(0);
611 }
612 
613 /*------------------------------------------------ Event Functions --------------------------------------------------*/
614 
615 /*@C
616   PetscLogEventRegister - Registers an event name for logging operations in an application code.
617 
618   Not Collective
619 
620   Input Parameters:
621 + name   - The name associated with the event
622 - classid - The classid associated to the class for this event, obtain either with
623            PetscClassIdRegister() or use a predefined one such as KSP_CLASSID, SNES_CLASSID, the predefined ones
624            are only available in C code
625 
626   Output Parameter:
627 . event - The event id for use with PetscLogEventBegin() and PetscLogEventEnd().
628 
629   Example of Usage:
630 .vb
631       PetscLogEvent USER_EVENT;
632       PetscClassId classid;
633       PetscLogDouble user_event_flops;
634       PetscClassIdRegister("class name",&classid);
635       PetscLogEventRegister("User event name",classid,&USER_EVENT);
636       PetscLogEventBegin(USER_EVENT,0,0,0,0);
637          [code segment to monitor]
638          PetscLogFlops(user_event_flops);
639       PetscLogEventEnd(USER_EVENT,0,0,0,0);
640 .ve
641 
642   Notes:
643   PETSc automatically logs library events if the code has been
644   configured with --with-log (which is the default) and
645   -log_view or -log_all is specified.  PetscLogEventRegister() is
646   intended for logging user events to supplement this PETSc
647   information.
648 
649   PETSc can gather data for use with the utilities Jumpshot
650   (part of the MPICH distribution).  If PETSc has been compiled
651   with flag -DPETSC_HAVE_MPE (MPE is an additional utility within
652   MPICH), the user can employ another command line option, -log_mpe,
653   to create a logfile, "mpe.log", which can be visualized
654   Jumpshot.
655 
656   The classid is associated with each event so that classes of events
657   can be disabled simultaneously, such as all matrix events. The user
658   can either use an existing classid, such as MAT_CLASSID, or create
659   their own as shown in the example.
660 
661   If an existing event with the same name exists, its event handle is
662   returned instead of creating a new event.
663 
664   Level: intermediate
665 
666 .seealso: `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscLogFlops()`,
667           `PetscLogEventActivate()`, `PetscLogEventDeactivate()`, `PetscClassIdRegister()`
668 @*/
669 PetscErrorCode PetscLogEventRegister(const char name[], PetscClassId classid, PetscLogEvent *event) {
670   PetscStageLog stageLog;
671   int           stage;
672 
673   PetscFunctionBegin;
674   *event = PETSC_DECIDE;
675   PetscCall(PetscLogGetStageLog(&stageLog));
676   PetscCall(PetscEventRegLogGetEvent(stageLog->eventLog, name, event));
677   if (*event > 0) PetscFunctionReturn(0);
678   PetscCall(PetscEventRegLogRegister(stageLog->eventLog, name, classid, event));
679   for (stage = 0; stage < stageLog->numStages; stage++) {
680     PetscCall(PetscEventPerfLogEnsureSize(stageLog->stageInfo[stage].eventLog, stageLog->eventLog->numEvents));
681     PetscCall(PetscClassPerfLogEnsureSize(stageLog->stageInfo[stage].classLog, stageLog->classLog->numClasses));
682   }
683   PetscFunctionReturn(0);
684 }
685 
686 /*@
687   PetscLogEventSetCollective - Indicates that a particular event is collective.
688 
689   Not Collective
690 
691   Input Parameters:
692 + event - The event id
693 - collective - Bolean flag indicating whether a particular event is collective
694 
695   Note:
696   New events returned from PetscLogEventRegister() are collective by default.
697 
698   Level: developer
699 
700 .seealso: `PetscLogEventRegister()`
701 @*/
702 PetscErrorCode PetscLogEventSetCollective(PetscLogEvent event, PetscBool collective) {
703   PetscStageLog    stageLog;
704   PetscEventRegLog eventRegLog;
705 
706   PetscFunctionBegin;
707   PetscCall(PetscLogGetStageLog(&stageLog));
708   PetscCall(PetscStageLogGetEventRegLog(stageLog, &eventRegLog));
709   PetscCheck(event >= 0 && event <= eventRegLog->numEvents, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid event id");
710   eventRegLog->eventInfo[event].collective = collective;
711   PetscFunctionReturn(0);
712 }
713 
714 /*@
715   PetscLogEventIncludeClass - Activates event logging for a PETSc object class in every stage.
716 
717   Not Collective
718 
719   Input Parameter:
720 . classid - The object class, for example MAT_CLASSID, SNES_CLASSID, etc.
721 
722   Level: developer
723 
724 .seealso: `PetscLogEventActivateClass()`, `PetscLogEventDeactivateClass()`, `PetscLogEventActivate()`, `PetscLogEventDeactivate()`
725 @*/
726 PetscErrorCode PetscLogEventIncludeClass(PetscClassId classid) {
727   PetscStageLog stageLog;
728   int           stage;
729 
730   PetscFunctionBegin;
731   PetscCall(PetscLogGetStageLog(&stageLog));
732   for (stage = 0; stage < stageLog->numStages; stage++) PetscCall(PetscEventPerfLogActivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid));
733   PetscFunctionReturn(0);
734 }
735 
736 /*@
737   PetscLogEventExcludeClass - Deactivates event logging for a PETSc object class in every stage.
738 
739   Not Collective
740 
741   Input Parameter:
742 . classid - The object class, for example MAT_CLASSID, SNES_CLASSID, etc.
743 
744   Level: developer
745 
746 .seealso: `PetscLogEventDeactivateClass()`, `PetscLogEventActivateClass()`, `PetscLogEventDeactivate()`, `PetscLogEventActivate()`
747 @*/
748 PetscErrorCode PetscLogEventExcludeClass(PetscClassId classid) {
749   PetscStageLog stageLog;
750   int           stage;
751 
752   PetscFunctionBegin;
753   PetscCall(PetscLogGetStageLog(&stageLog));
754   for (stage = 0; stage < stageLog->numStages; stage++) PetscCall(PetscEventPerfLogDeactivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid));
755   PetscFunctionReturn(0);
756 }
757 
758 /*@
759   PetscLogEventActivate - Indicates that a particular event should be logged.
760 
761   Not Collective
762 
763   Input Parameter:
764 . event - The event id
765 
766   Usage:
767 .vb
768       PetscLogEventDeactivate(VEC_SetValues);
769         [code where you do not want to log VecSetValues()]
770       PetscLogEventActivate(VEC_SetValues);
771         [code where you do want to log VecSetValues()]
772 .ve
773 
774   Note:
775   The event may be either a pre-defined PETSc event (found in include/petsclog.h)
776   or an event number obtained with PetscLogEventRegister().
777 
778   Level: advanced
779 
780 .seealso: `PlogEventDeactivate()`, `PlogEventDeactivatePush()`, `PetscLogEventDeactivatePop()`
781 @*/
782 PetscErrorCode PetscLogEventActivate(PetscLogEvent event) {
783   PetscStageLog stageLog;
784   int           stage;
785 
786   PetscFunctionBegin;
787   PetscCall(PetscLogGetStageLog(&stageLog));
788   PetscCall(PetscStageLogGetCurrent(stageLog, &stage));
789   PetscCall(PetscEventPerfLogActivate(stageLog->stageInfo[stage].eventLog, event));
790   PetscFunctionReturn(0);
791 }
792 
793 /*@
794   PetscLogEventDeactivate - Indicates that a particular event should not be logged.
795 
796   Not Collective
797 
798   Input Parameter:
799 . event - The event id
800 
801   Usage:
802 .vb
803       PetscLogEventDeactivate(VEC_SetValues);
804         [code where you do not want to log VecSetValues()]
805       PetscLogEventActivate(VEC_SetValues);
806         [code where you do want to log VecSetValues()]
807 .ve
808 
809   Note:
810   The event may be either a pre-defined PETSc event (found in
811   include/petsclog.h) or an event number obtained with PetscLogEventRegister()).
812 
813   Level: advanced
814 
815 .seealso: `PetscLogEventActivate()`, `PetscLogEventDeactivatePush()`, `PetscLogEventDeactivatePop()`
816 @*/
817 PetscErrorCode PetscLogEventDeactivate(PetscLogEvent event) {
818   PetscStageLog stageLog;
819   int           stage;
820 
821   PetscFunctionBegin;
822   PetscCall(PetscLogGetStageLog(&stageLog));
823   PetscCall(PetscStageLogGetCurrent(stageLog, &stage));
824   PetscCall(PetscEventPerfLogDeactivate(stageLog->stageInfo[stage].eventLog, event));
825   PetscFunctionReturn(0);
826 }
827 
828 /*@
829   PetscLogEventDeactivatePush - Indicates that a particular event should not be logged.
830 
831   Not Collective
832 
833   Input Parameter:
834 . event - The event id
835 
836   Usage:
837 .vb
838       PetscLogEventDeactivatePush(VEC_SetValues);
839         [code where you do not want to log VecSetValues()]
840       PetscLogEventDeactivatePop(VEC_SetValues);
841         [code where you do want to log VecSetValues()]
842 .ve
843 
844   Note:
845   The event may be either a pre-defined PETSc event (found in
846   include/petsclog.h) or an event number obtained with PetscLogEventRegister()).
847 
848   Level: advanced
849 
850 .seealso: `PetscLogEventActivate()`, `PetscLogEventDeactivatePop()`
851 @*/
852 PetscErrorCode PetscLogEventDeactivatePush(PetscLogEvent event) {
853   PetscStageLog stageLog;
854   int           stage;
855 
856   PetscFunctionBegin;
857   PetscCall(PetscLogGetStageLog(&stageLog));
858   PetscCall(PetscStageLogGetCurrent(stageLog, &stage));
859   PetscCall(PetscEventPerfLogDeactivatePush(stageLog->stageInfo[stage].eventLog, event));
860   PetscFunctionReturn(0);
861 }
862 
863 /*@
864   PetscLogEventDeactivatePop - Indicates that a particular event should be logged.
865 
866   Not Collective
867 
868   Input Parameter:
869 . event - The event id
870 
871   Usage:
872 .vb
873       PetscLogEventDeactivatePush(VEC_SetValues);
874         [code where you do not want to log VecSetValues()]
875       PetscLogEventDeactivatePop(VEC_SetValues);
876         [code where you do want to log VecSetValues()]
877 .ve
878 
879   Note:
880   The event may be either a pre-defined PETSc event (found in
881   include/petsclog.h) or an event number obtained with PetscLogEventRegister()).
882 
883   Level: advanced
884 
885 .seealso: `PetscLogEventActivate()`, `PetscLogEventDeactivatePush()`
886 @*/
887 PetscErrorCode PetscLogEventDeactivatePop(PetscLogEvent event) {
888   PetscStageLog stageLog;
889   int           stage;
890 
891   PetscFunctionBegin;
892   PetscCall(PetscLogGetStageLog(&stageLog));
893   PetscCall(PetscStageLogGetCurrent(stageLog, &stage));
894   PetscCall(PetscEventPerfLogDeactivatePop(stageLog->stageInfo[stage].eventLog, event));
895   PetscFunctionReturn(0);
896 }
897 
898 /*@
899   PetscLogEventSetActiveAll - Sets the event activity in every stage.
900 
901   Not Collective
902 
903   Input Parameters:
904 + event    - The event id
905 - isActive - The activity flag determining whether the event is logged
906 
907   Level: advanced
908 
909 .seealso: `PlogEventActivate()`, `PlogEventDeactivate()`
910 @*/
911 PetscErrorCode PetscLogEventSetActiveAll(PetscLogEvent event, PetscBool isActive) {
912   PetscStageLog stageLog;
913   int           stage;
914 
915   PetscFunctionBegin;
916   PetscCall(PetscLogGetStageLog(&stageLog));
917   for (stage = 0; stage < stageLog->numStages; stage++) {
918     if (isActive) {
919       PetscCall(PetscEventPerfLogActivate(stageLog->stageInfo[stage].eventLog, event));
920     } else {
921       PetscCall(PetscEventPerfLogDeactivate(stageLog->stageInfo[stage].eventLog, event));
922     }
923   }
924   PetscFunctionReturn(0);
925 }
926 
927 /*@
928   PetscLogEventActivateClass - Activates event logging for a PETSc object class.
929 
930   Not Collective
931 
932   Input Parameter:
933 . classid - The event class, for example MAT_CLASSID, SNES_CLASSID, etc.
934 
935   Level: developer
936 
937 .seealso: `PetscLogEventDeactivateClass()`, `PetscLogEventActivate()`, `PetscLogEventDeactivate()`
938 @*/
939 PetscErrorCode PetscLogEventActivateClass(PetscClassId classid) {
940   PetscStageLog stageLog;
941   int           stage;
942 
943   PetscFunctionBegin;
944   PetscCall(PetscLogGetStageLog(&stageLog));
945   PetscCall(PetscStageLogGetCurrent(stageLog, &stage));
946   PetscCall(PetscEventPerfLogActivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid));
947   PetscFunctionReturn(0);
948 }
949 
950 /*@
951   PetscLogEventDeactivateClass - Deactivates event logging for a PETSc object class.
952 
953   Not Collective
954 
955   Input Parameter:
956 . classid - The event class, for example MAT_CLASSID, SNES_CLASSID, etc.
957 
958   Level: developer
959 
960 .seealso: `PetscLogEventActivateClass()`, `PetscLogEventActivate()`, `PetscLogEventDeactivate()`
961 @*/
962 PetscErrorCode PetscLogEventDeactivateClass(PetscClassId classid) {
963   PetscStageLog stageLog;
964   int           stage;
965 
966   PetscFunctionBegin;
967   PetscCall(PetscLogGetStageLog(&stageLog));
968   PetscCall(PetscStageLogGetCurrent(stageLog, &stage));
969   PetscCall(PetscEventPerfLogDeactivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid));
970   PetscFunctionReturn(0);
971 }
972 
973 /*MC
974    PetscLogEventSync - Synchronizes the beginning of a user event.
975 
976    Synopsis:
977    #include <petsclog.h>
978    PetscErrorCode PetscLogEventSync(int e,MPI_Comm comm)
979 
980    Collective
981 
982    Input Parameters:
983 +  e - integer associated with the event obtained from PetscLogEventRegister()
984 -  comm - an MPI communicator
985 
986    Usage:
987 .vb
988      PetscLogEvent USER_EVENT;
989      PetscLogEventRegister("User event",0,&USER_EVENT);
990      PetscLogEventSync(USER_EVENT,PETSC_COMM_WORLD);
991      PetscLogEventBegin(USER_EVENT,0,0,0,0);
992         [code segment to monitor]
993      PetscLogEventEnd(USER_EVENT,0,0,0,0);
994 .ve
995 
996    Notes:
997    This routine should be called only if there is not a
998    PetscObject available to pass to PetscLogEventBegin().
999 
1000    Level: developer
1001 
1002 .seealso: `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`
1003 
1004 M*/
1005 
1006 /*MC
1007    PetscLogEventBegin - Logs the beginning of a user event.
1008 
1009    Synopsis:
1010    #include <petsclog.h>
1011    PetscErrorCode PetscLogEventBegin(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
1012 
1013    Not Collective
1014 
1015    Input Parameters:
1016 +  e - integer associated with the event obtained from PetscLogEventRegister()
1017 -  o1,o2,o3,o4 - objects associated with the event, or 0
1018 
1019    Fortran Synopsis:
1020    void PetscLogEventBegin(int e,PetscErrorCode ierr)
1021 
1022    Usage:
1023 .vb
1024      PetscLogEvent USER_EVENT;
1025      PetscLogDouble user_event_flops;
1026      PetscLogEventRegister("User event",0,&USER_EVENT);
1027      PetscLogEventBegin(USER_EVENT,0,0,0,0);
1028         [code segment to monitor]
1029         PetscLogFlops(user_event_flops);
1030      PetscLogEventEnd(USER_EVENT,0,0,0,0);
1031 .ve
1032 
1033    Notes:
1034    You need to register each integer event with the command
1035    PetscLogEventRegister().
1036 
1037    Developer Notes:
1038      PetscLogEventBegin() and PetscLogEventBegin() return error codes instead of explicitly handling the
1039      errors that occur in the macro directly because other packages that use this macros have used them in their
1040      own functions or methods that do not return error codes and it would be disruptive to change the current
1041      behavior.
1042 
1043    Level: intermediate
1044 
1045 .seealso: `PetscLogEventRegister()`, `PetscLogEventEnd()`, `PetscLogFlops()`
1046 
1047 M*/
1048 
1049 /*MC
1050    PetscLogEventEnd - Log the end of a user event.
1051 
1052    Synopsis:
1053    #include <petsclog.h>
1054    PetscErrorCode PetscLogEventEnd(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
1055 
1056    Not Collective
1057 
1058    Input Parameters:
1059 +  e - integer associated with the event obtained with PetscLogEventRegister()
1060 -  o1,o2,o3,o4 - objects associated with the event, or 0
1061 
1062    Fortran Synopsis:
1063    void PetscLogEventEnd(int e,PetscErrorCode ierr)
1064 
1065    Usage:
1066 .vb
1067      PetscLogEvent USER_EVENT;
1068      PetscLogDouble user_event_flops;
1069      PetscLogEventRegister("User event",0,&USER_EVENT,);
1070      PetscLogEventBegin(USER_EVENT,0,0,0,0);
1071         [code segment to monitor]
1072         PetscLogFlops(user_event_flops);
1073      PetscLogEventEnd(USER_EVENT,0,0,0,0);
1074 .ve
1075 
1076    Notes:
1077    You should also register each additional integer event with the command
1078    PetscLogEventRegister().
1079 
1080    Level: intermediate
1081 
1082 .seealso: `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogFlops()`
1083 
1084 M*/
1085 
1086 /*@C
1087   PetscLogEventGetId - Returns the event id when given the event name.
1088 
1089   Not Collective
1090 
1091   Input Parameter:
1092 . name  - The event name
1093 
1094   Output Parameter:
1095 . event - The event, or -1 if no event with that name exists
1096 
1097   Level: intermediate
1098 
1099 .seealso: `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscLogStageGetId()`
1100 @*/
1101 PetscErrorCode PetscLogEventGetId(const char name[], PetscLogEvent *event) {
1102   PetscStageLog stageLog;
1103 
1104   PetscFunctionBegin;
1105   PetscCall(PetscLogGetStageLog(&stageLog));
1106   PetscCall(PetscEventRegLogGetEvent(stageLog->eventLog, name, event));
1107   PetscFunctionReturn(0);
1108 }
1109 
1110 PetscErrorCode PetscLogPushCurrentEvent_Internal(PetscLogEvent event) {
1111   PetscFunctionBegin;
1112   PetscCall(PetscIntStackPush(current_log_event_stack, event));
1113   PetscFunctionReturn(0);
1114 }
1115 
1116 PetscErrorCode PetscLogPopCurrentEvent_Internal(void) {
1117   PetscFunctionBegin;
1118   PetscCall(PetscIntStackPop(current_log_event_stack, NULL));
1119   PetscFunctionReturn(0);
1120 }
1121 
1122 PetscErrorCode PetscLogGetCurrentEvent_Internal(PetscLogEvent *event) {
1123   PetscBool empty;
1124 
1125   PetscFunctionBegin;
1126   PetscValidIntPointer(event, 1);
1127   *event = PETSC_DECIDE;
1128   PetscCall(PetscIntStackEmpty(current_log_event_stack, &empty));
1129   if (!empty) PetscCall(PetscIntStackTop(current_log_event_stack, event));
1130   PetscFunctionReturn(0);
1131 }
1132 
1133 PetscErrorCode PetscLogEventPause_Internal(PetscLogEvent event) {
1134   PetscFunctionBegin;
1135   if (event != PETSC_DECIDE) PetscCall(PetscLogEventEnd(event, NULL, NULL, NULL, NULL));
1136   PetscFunctionReturn(0);
1137 }
1138 
1139 PetscErrorCode PetscLogEventResume_Internal(PetscLogEvent event) {
1140   PetscStageLog     stageLog;
1141   PetscEventPerfLog eventLog;
1142   int               stage;
1143 
1144   PetscFunctionBegin;
1145   if (event == PETSC_DECIDE) PetscFunctionReturn(0);
1146   PetscCall(PetscLogEventBegin(event, NULL, NULL, NULL, NULL));
1147   PetscCall(PetscLogGetStageLog(&stageLog));
1148   PetscCall(PetscStageLogGetCurrent(stageLog, &stage));
1149   PetscCall(PetscStageLogGetEventPerfLog(stageLog, stage, &eventLog));
1150   eventLog->eventInfo[event].count--;
1151   PetscFunctionReturn(0);
1152 }
1153 
1154 /*------------------------------------------------ Output Functions -------------------------------------------------*/
1155 /*@C
1156   PetscLogDump - Dumps logs of objects to a file. This file is intended to
1157   be read by bin/petscview. This program no longer exists.
1158 
1159   Collective on PETSC_COMM_WORLD
1160 
1161   Input Parameter:
1162 . name - an optional file name
1163 
1164   Usage:
1165 .vb
1166      PetscInitialize(...);
1167      PetscLogDefaultBegin(); or PetscLogAllBegin();
1168      ... code ...
1169      PetscLogDump(filename);
1170      PetscFinalize();
1171 .ve
1172 
1173   Notes:
1174   The default file name is
1175 $    Log.<rank>
1176   where <rank> is the processor number. If no name is specified,
1177   this file will be used.
1178 
1179   Level: advanced
1180 
1181 .seealso: `PetscLogDefaultBegin()`, `PetscLogAllBegin()`, `PetscLogView()`
1182 @*/
1183 PetscErrorCode PetscLogDump(const char sname[]) {
1184   PetscStageLog       stageLog;
1185   PetscEventPerfInfo *eventInfo;
1186   FILE               *fd;
1187   char                file[PETSC_MAX_PATH_LEN], fname[PETSC_MAX_PATH_LEN];
1188   PetscLogDouble      flops, _TotalTime;
1189   PetscMPIInt         rank;
1190   int                 action, object, curStage;
1191   PetscLogEvent       event;
1192 
1193   PetscFunctionBegin;
1194   /* Calculate the total elapsed time */
1195   PetscTime(&_TotalTime);
1196   _TotalTime -= petsc_BaseTime;
1197   /* Open log file */
1198   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1199   if (sname && sname[0]) sprintf(file, "%s.%d", sname, rank);
1200   else sprintf(file, "Log.%d", rank);
1201   PetscCall(PetscFixFilename(file, fname));
1202   PetscCall(PetscFOpen(PETSC_COMM_WORLD, fname, "w", &fd));
1203   PetscCheck(!(rank == 0) || !(!fd), PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file: %s", fname);
1204   /* Output totals */
1205   PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "Total Flop %14e %16.8e\n", petsc_TotalFlops, _TotalTime));
1206   PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "Clock Resolution %g\n", 0.0));
1207   /* Output actions */
1208   if (petsc_logActions) {
1209     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "Actions accomplished %d\n", petsc_numActions));
1210     for (action = 0; action < petsc_numActions; action++) {
1211       PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%g %d %d %d %d %d %d %g %g %g\n", petsc_actions[action].time, petsc_actions[action].action, (int)petsc_actions[action].event, (int)petsc_actions[action].classid, petsc_actions[action].id1,
1212                              petsc_actions[action].id2, petsc_actions[action].id3, petsc_actions[action].flops, petsc_actions[action].mem, petsc_actions[action].maxmem));
1213     }
1214   }
1215   /* Output objects */
1216   if (petsc_logObjects) {
1217     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "Objects created %d destroyed %d\n", petsc_numObjects, petsc_numObjectsDestroyed));
1218     for (object = 0; object < petsc_numObjects; object++) {
1219       PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "Parent ID: %d Memory: %d\n", petsc_objects[object].parent, (int)petsc_objects[object].mem));
1220       if (!petsc_objects[object].name[0]) {
1221         PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "No Name\n"));
1222       } else {
1223         PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "Name: %s\n", petsc_objects[object].name));
1224       }
1225       if (petsc_objects[object].info[0] != 0) {
1226         PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "No Info\n"));
1227       } else {
1228         PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "Info: %s\n", petsc_objects[object].info));
1229       }
1230     }
1231   }
1232   /* Output events */
1233   PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "Event log:\n"));
1234   PetscCall(PetscLogGetStageLog(&stageLog));
1235   PetscCall(PetscIntStackTop(stageLog->stack, &curStage));
1236   eventInfo = stageLog->stageInfo[curStage].eventLog->eventInfo;
1237   for (event = 0; event < stageLog->stageInfo[curStage].eventLog->numEvents; event++) {
1238     if (eventInfo[event].time != 0.0) flops = eventInfo[event].flops / eventInfo[event].time;
1239     else flops = 0.0;
1240     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fd, "%d %16d %16g %16g %16g\n", event, eventInfo[event].count, eventInfo[event].flops, eventInfo[event].time, flops));
1241   }
1242   PetscCall(PetscFClose(PETSC_COMM_WORLD, fd));
1243   PetscFunctionReturn(0);
1244 }
1245 
1246 /*
1247   PetscLogView_Detailed - Each process prints the times for its own events
1248 
1249 */
1250 PetscErrorCode PetscLogView_Detailed(PetscViewer viewer) {
1251   PetscStageLog       stageLog;
1252   PetscEventPerfInfo *eventInfo = NULL, *stageInfo = NULL;
1253   PetscLogDouble      locTotalTime, numRed, maxMem;
1254   int                 numStages, numEvents, stage, event;
1255   MPI_Comm            comm = PetscObjectComm((PetscObject)viewer);
1256   PetscMPIInt         rank, size;
1257 
1258   PetscFunctionBegin;
1259   PetscCallMPI(MPI_Comm_size(comm, &size));
1260   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1261   /* Must preserve reduction count before we go on */
1262   numRed = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
1263   /* Get the total elapsed time */
1264   PetscTime(&locTotalTime);
1265   locTotalTime -= petsc_BaseTime;
1266   PetscCall(PetscViewerASCIIPrintf(viewer, "size = %d\n", size));
1267   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalTimes = {}\n"));
1268   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessages = {}\n"));
1269   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessageLens = {}\n"));
1270   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalReductions = {}\n"));
1271   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalFlop = {}\n"));
1272   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalObjects = {}\n"));
1273   PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMemory = {}\n"));
1274   PetscCall(PetscLogGetStageLog(&stageLog));
1275   PetscCallMPI(MPI_Allreduce(&stageLog->numStages, &numStages, 1, MPI_INT, MPI_MAX, comm));
1276   PetscCall(PetscViewerASCIIPrintf(viewer, "Stages = {}\n"));
1277   for (stage = 0; stage < numStages; stage++) {
1278     PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"] = {}\n", stageLog->stageInfo[stage].name));
1279     PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"summary\"] = {}\n", stageLog->stageInfo[stage].name));
1280     PetscCallMPI(MPI_Allreduce(&stageLog->stageInfo[stage].eventLog->numEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm));
1281     for (event = 0; event < numEvents; event++) PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"%s\"] = {}\n", stageLog->stageInfo[stage].name, stageLog->eventLog->eventInfo[event].name));
1282   }
1283   PetscCall(PetscMallocGetMaximumUsage(&maxMem));
1284   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1285   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalTimes[%d] = %g\n", rank, locTotalTime));
1286   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessages[%d] = %g\n", rank, (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct)));
1287   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessageLens[%d] = %g\n", rank, (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len)));
1288   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalReductions[%d] = %g\n", rank, numRed));
1289   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalFlop[%d] = %g\n", rank, petsc_TotalFlops));
1290   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalObjects[%d] = %d\n", rank, petsc_numObjects));
1291   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMemory[%d] = %g\n", rank, maxMem));
1292   PetscCall(PetscViewerFlush(viewer));
1293   for (stage = 0; stage < numStages; stage++) {
1294     stageInfo = &stageLog->stageInfo[stage].perfInfo;
1295     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Stages[\"%s\"][\"summary\"][%d] = {\"time\" : %g, \"numMessages\" : %g, \"messageLength\" : %g, \"numReductions\" : %g, \"flop\" : %g}\n", stageLog->stageInfo[stage].name, rank, stageInfo->time,
1296                                                  stageInfo->numMessages, stageInfo->messageLength, stageInfo->numReductions, stageInfo->flops));
1297     PetscCallMPI(MPI_Allreduce(&stageLog->stageInfo[stage].eventLog->numEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm));
1298     for (event = 0; event < numEvents; event++) {
1299       eventInfo = &stageLog->stageInfo[stage].eventLog->eventInfo[event];
1300       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Stages[\"%s\"][\"%s\"][%d] = {\"count\" : %d, \"time\" : %g, \"syncTime\" : %g, \"numMessages\" : %g, \"messageLength\" : %g, \"numReductions\" : %g, \"flop\" : %g",
1301                                                    stageLog->stageInfo[stage].name, stageLog->eventLog->eventInfo[event].name, rank, eventInfo->count, eventInfo->time, eventInfo->syncTime, eventInfo->numMessages, eventInfo->messageLength, eventInfo->numReductions,
1302                                                    eventInfo->flops));
1303       if (eventInfo->dof[0] >= 0.) {
1304         PetscInt d, e;
1305 
1306         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"dof\" : ["));
1307         for (d = 0; d < 8; ++d) {
1308           if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", "));
1309           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->dof[d]));
1310         }
1311         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]"));
1312         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"error\" : ["));
1313         for (e = 0; e < 8; ++e) {
1314           if (e > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", "));
1315           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->errors[e]));
1316         }
1317         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]"));
1318       }
1319       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "}\n"));
1320     }
1321   }
1322   PetscCall(PetscViewerFlush(viewer));
1323   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1324   PetscFunctionReturn(0);
1325 }
1326 
1327 /*
1328   PetscLogView_CSV - Each process prints the times for its own events in Comma-Separated Value Format
1329 */
1330 PetscErrorCode PetscLogView_CSV(PetscViewer viewer) {
1331   PetscStageLog       stageLog;
1332   PetscEventPerfInfo *eventInfo = NULL;
1333   PetscLogDouble      locTotalTime, maxMem;
1334   int                 numStages, numEvents, stage, event;
1335   MPI_Comm            comm = PetscObjectComm((PetscObject)viewer);
1336   PetscMPIInt         rank, size;
1337 
1338   PetscFunctionBegin;
1339   PetscCallMPI(MPI_Comm_size(comm, &size));
1340   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1341   /* Must preserve reduction count before we go on */
1342   /* Get the total elapsed time */
1343   PetscTime(&locTotalTime);
1344   locTotalTime -= petsc_BaseTime;
1345   PetscCall(PetscLogGetStageLog(&stageLog));
1346   PetscCallMPI(MPI_Allreduce(&stageLog->numStages, &numStages, 1, MPI_INT, MPI_MAX, comm));
1347   PetscCall(PetscMallocGetMaximumUsage(&maxMem));
1348   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1349   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));
1350   PetscCall(PetscViewerFlush(viewer));
1351   for (stage = 0; stage < numStages; stage++) {
1352     PetscEventPerfInfo *stageInfo = &stageLog->stageInfo[stage].perfInfo;
1353 
1354     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s,summary,%d,1,%g,%g,%g,%g,%g\n", stageLog->stageInfo[stage].name, rank, stageInfo->time, stageInfo->numMessages, stageInfo->messageLength, stageInfo->numReductions, stageInfo->flops));
1355     PetscCallMPI(MPI_Allreduce(&stageLog->stageInfo[stage].eventLog->numEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm));
1356     for (event = 0; event < numEvents; event++) {
1357       eventInfo = &stageLog->stageInfo[stage].eventLog->eventInfo[event];
1358       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s,%s,%d,%d,%g,%g,%g,%g,%g", stageLog->stageInfo[stage].name, stageLog->eventLog->eventInfo[event].name, rank, eventInfo->count, eventInfo->time, eventInfo->numMessages, eventInfo->messageLength,
1359                                                    eventInfo->numReductions, eventInfo->flops));
1360       if (eventInfo->dof[0] >= 0.) {
1361         PetscInt d, e;
1362 
1363         for (d = 0; d < 8; ++d) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->dof[d]));
1364         for (e = 0; e < 8; ++e) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->errors[e]));
1365       }
1366       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
1367     }
1368   }
1369   PetscCall(PetscViewerFlush(viewer));
1370   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1371   PetscFunctionReturn(0);
1372 }
1373 
1374 static PetscErrorCode PetscLogViewWarnSync(MPI_Comm comm, FILE *fd) {
1375   PetscFunctionBegin;
1376   if (!PetscLogSyncOn) PetscFunctionReturn(0);
1377   PetscCall(PetscFPrintf(comm, fd, "\n\n"));
1378   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n"));
1379   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1380   PetscCall(PetscFPrintf(comm, fd, "      #                       WARNING!!!                       #\n"));
1381   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1382   PetscCall(PetscFPrintf(comm, fd, "      #   This program was run with logging synchronization.   #\n"));
1383   PetscCall(PetscFPrintf(comm, fd, "      #   This option provides more meaningful imbalance       #\n"));
1384   PetscCall(PetscFPrintf(comm, fd, "      #   figures at the expense of slowing things down and    #\n"));
1385   PetscCall(PetscFPrintf(comm, fd, "      #   providing a distorted view of the overall runtime.   #\n"));
1386   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1387   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n\n\n"));
1388   PetscFunctionReturn(0);
1389 }
1390 
1391 static PetscErrorCode PetscLogViewWarnDebugging(MPI_Comm comm, FILE *fd) {
1392   PetscFunctionBegin;
1393   if (PetscDefined(USE_DEBUG)) {
1394     PetscCall(PetscFPrintf(comm, fd, "\n\n"));
1395     PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n"));
1396     PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1397     PetscCall(PetscFPrintf(comm, fd, "      #                       WARNING!!!                       #\n"));
1398     PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1399     PetscCall(PetscFPrintf(comm, fd, "      #   This code was compiled with a debugging option.      #\n"));
1400     PetscCall(PetscFPrintf(comm, fd, "      #   To get timing results run ./configure                #\n"));
1401     PetscCall(PetscFPrintf(comm, fd, "      #   using --with-debugging=no, the performance will      #\n"));
1402     PetscCall(PetscFPrintf(comm, fd, "      #   be generally two or three times faster.              #\n"));
1403     PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1404     PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n\n\n"));
1405   }
1406   PetscFunctionReturn(0);
1407 }
1408 
1409 static PetscErrorCode PetscLogViewWarnNoGpuAwareMpi(MPI_Comm comm, FILE *fd) {
1410 #if defined(PETSC_HAVE_DEVICE)
1411   PetscMPIInt size;
1412 
1413   PetscFunctionBegin;
1414   PetscCallMPI(MPI_Comm_size(comm, &size));
1415   if (use_gpu_aware_mpi || size == 1) PetscFunctionReturn(0);
1416   PetscCall(PetscFPrintf(comm, fd, "\n\n"));
1417   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n"));
1418   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1419   PetscCall(PetscFPrintf(comm, fd, "      #                       WARNING!!!                       #\n"));
1420   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1421   PetscCall(PetscFPrintf(comm, fd, "      #   This code was compiled with GPU support and you've   #\n"));
1422   PetscCall(PetscFPrintf(comm, fd, "      #   created PETSc/GPU objects, but you intentionally     #\n"));
1423   PetscCall(PetscFPrintf(comm, fd, "      #   used -use_gpu_aware_mpi 0, requiring PETSc to copy   #\n"));
1424   PetscCall(PetscFPrintf(comm, fd, "      #   additional data between the GPU and CPU. To obtain   #\n"));
1425   PetscCall(PetscFPrintf(comm, fd, "      #   meaningful timing results on multi-rank runs, use    #\n"));
1426   PetscCall(PetscFPrintf(comm, fd, "      #   GPU-aware MPI instead.                               #\n"));
1427   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1428   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n\n\n"));
1429   PetscFunctionReturn(0);
1430 #else
1431   return 0;
1432 #endif
1433 }
1434 
1435 static PetscErrorCode PetscLogViewWarnGpuTime(MPI_Comm comm, FILE *fd) {
1436 #if defined(PETSC_HAVE_DEVICE)
1437 
1438   PetscFunctionBegin;
1439   if (!PetscLogGpuTimeFlag || petsc_gflops == 0) PetscFunctionReturn(0);
1440   PetscCall(PetscFPrintf(comm, fd, "\n\n"));
1441   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n"));
1442   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1443   PetscCall(PetscFPrintf(comm, fd, "      #                       WARNING!!!                       #\n"));
1444   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1445   PetscCall(PetscFPrintf(comm, fd, "      #   This code was run with -log_view_gpu_time            #\n"));
1446   PetscCall(PetscFPrintf(comm, fd, "      #   This provides accurate timing within the GPU kernels #\n"));
1447   PetscCall(PetscFPrintf(comm, fd, "      #   but can slow down the entire computation by a        #\n"));
1448   PetscCall(PetscFPrintf(comm, fd, "      #   measurable amount. For fastest runs we recommend     #\n"));
1449   PetscCall(PetscFPrintf(comm, fd, "      #   not using this option.                               #\n"));
1450   PetscCall(PetscFPrintf(comm, fd, "      #                                                        #\n"));
1451   PetscCall(PetscFPrintf(comm, fd, "      ##########################################################\n\n\n"));
1452   PetscFunctionReturn(0);
1453 #else
1454   return 0;
1455 #endif
1456 }
1457 
1458 PetscErrorCode PetscLogView_Default(PetscViewer viewer) {
1459   FILE               *fd;
1460   PetscLogDouble      zero = 0.0;
1461   PetscStageLog       stageLog;
1462   PetscStageInfo     *stageInfo = NULL;
1463   PetscEventPerfInfo *eventInfo = NULL;
1464   PetscClassPerfInfo *classInfo;
1465   char                arch[128], hostname[128], username[128], pname[PETSC_MAX_PATH_LEN], date[128];
1466   const char         *name;
1467   PetscLogDouble      locTotalTime, TotalTime, TotalFlops;
1468   PetscLogDouble      numMessages, messageLength, avgMessLen, numReductions;
1469   PetscLogDouble      stageTime, flops, flopr, mem, mess, messLen, red;
1470   PetscLogDouble      fracTime, fracFlops, fracMessages, fracLength, fracReductions, fracMess, fracMessLen, fracRed;
1471   PetscLogDouble      fracStageTime, fracStageFlops, fracStageMess, fracStageMessLen, fracStageRed;
1472   PetscLogDouble      min, max, tot, ratio, avg, x, y;
1473   PetscLogDouble      minf, maxf, totf, ratf, mint, maxt, tott, ratt, ratC, totm, totml, totr, mal, malmax, emalmax;
1474 #if defined(PETSC_HAVE_DEVICE)
1475   PetscLogEvent  KSP_Solve, SNES_Solve, TS_Step, TAO_Solve; /* These need to be fixed to be some events registered with certain objects */
1476   PetscLogDouble cct, gct, csz, gsz, gmaxt, gflops, gflopr, fracgflops;
1477 #endif
1478   PetscMPIInt    minC, maxC;
1479   PetscMPIInt    size, rank;
1480   PetscBool     *localStageUsed, *stageUsed;
1481   PetscBool     *localStageVisible, *stageVisible;
1482   int            numStages, localNumEvents, numEvents;
1483   int            stage, oclass;
1484   PetscLogEvent  event;
1485   PetscErrorCode ierr = 0;
1486   char           version[256];
1487   MPI_Comm       comm;
1488 #if defined(PETSC_HAVE_DEVICE)
1489   PetscLogEvent eventid;
1490   PetscInt64    nas = 0x7FF0000000000002;
1491 #endif
1492 
1493   PetscFunctionBegin;
1494   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1495   PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
1496   PetscCall(PetscViewerASCIIGetPointer(viewer, &fd));
1497   PetscCallMPI(MPI_Comm_size(comm, &size));
1498   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1499   /* Get the total elapsed time */
1500   PetscTime(&locTotalTime);
1501   locTotalTime -= petsc_BaseTime;
1502 
1503   PetscCall(PetscFPrintf(comm, fd, "****************************************************************************************************************************************************************\n"));
1504   PetscCall(PetscFPrintf(comm, fd, "***                                WIDEN YOUR WINDOW TO 160 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document                                 ***\n"));
1505   PetscCall(PetscFPrintf(comm, fd, "****************************************************************************************************************************************************************\n"));
1506   PetscCall(PetscFPrintf(comm, fd, "\n------------------------------------------------------------------ PETSc Performance Summary: ------------------------------------------------------------------\n\n"));
1507   PetscCall(PetscLogViewWarnSync(comm, fd));
1508   PetscCall(PetscLogViewWarnDebugging(comm, fd));
1509   PetscCall(PetscLogViewWarnNoGpuAwareMpi(comm, fd));
1510   PetscCall(PetscLogViewWarnGpuTime(comm, fd));
1511   PetscCall(PetscGetArchType(arch, sizeof(arch)));
1512   PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
1513   PetscCall(PetscGetUserName(username, sizeof(username)));
1514   PetscCall(PetscGetProgramName(pname, sizeof(pname)));
1515   PetscCall(PetscGetDate(date, sizeof(date)));
1516   PetscCall(PetscGetVersion(version, sizeof(version)));
1517   if (size == 1) {
1518     PetscCall(PetscFPrintf(comm, fd, "%s on a %s named %s with %d processor, by %s %s\n", pname, arch, hostname, size, username, date));
1519   } else {
1520     PetscCall(PetscFPrintf(comm, fd, "%s on a %s named %s with %d processors, by %s %s\n", pname, arch, hostname, size, username, date));
1521   }
1522 #if defined(PETSC_HAVE_OPENMP)
1523   PetscCall(PetscFPrintf(comm, fd, "Using %" PetscInt_FMT " OpenMP threads\n", PetscNumOMPThreads));
1524 #endif
1525   PetscCall(PetscFPrintf(comm, fd, "Using %s\n", version));
1526 
1527   /* Must preserve reduction count before we go on */
1528   red = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
1529 
1530   /* Calculate summary information */
1531   PetscCall(PetscFPrintf(comm, fd, "\n                         Max       Max/Min     Avg       Total\n"));
1532   /*   Time */
1533   PetscCallMPI(MPI_Allreduce(&locTotalTime, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1534   PetscCallMPI(MPI_Allreduce(&locTotalTime, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1535   PetscCallMPI(MPI_Allreduce(&locTotalTime, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1536   avg = tot / ((PetscLogDouble)size);
1537   if (min != 0.0) ratio = max / min;
1538   else ratio = 0.0;
1539   PetscCall(PetscFPrintf(comm, fd, "Time (sec):           %5.3e   %7.3f   %5.3e\n", max, ratio, avg));
1540   TotalTime = tot;
1541   /*   Objects */
1542   avg       = (PetscLogDouble)petsc_numObjects;
1543   PetscCallMPI(MPI_Allreduce(&avg, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1544   PetscCallMPI(MPI_Allreduce(&avg, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1545   PetscCallMPI(MPI_Allreduce(&avg, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1546   avg = tot / ((PetscLogDouble)size);
1547   if (min != 0.0) ratio = max / min;
1548   else ratio = 0.0;
1549   PetscCall(PetscFPrintf(comm, fd, "Objects:              %5.3e   %7.3f   %5.3e\n", max, ratio, avg));
1550   /*   Flops */
1551   PetscCallMPI(MPI_Allreduce(&petsc_TotalFlops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1552   PetscCallMPI(MPI_Allreduce(&petsc_TotalFlops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1553   PetscCallMPI(MPI_Allreduce(&petsc_TotalFlops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1554   avg = tot / ((PetscLogDouble)size);
1555   if (min != 0.0) ratio = max / min;
1556   else ratio = 0.0;
1557   PetscCall(PetscFPrintf(comm, fd, "Flops:                %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1558   TotalFlops = tot;
1559   /*   Flops/sec -- Must talk to Barry here */
1560   if (locTotalTime != 0.0) flops = petsc_TotalFlops / locTotalTime;
1561   else flops = 0.0;
1562   PetscCallMPI(MPI_Allreduce(&flops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1563   PetscCallMPI(MPI_Allreduce(&flops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1564   PetscCallMPI(MPI_Allreduce(&flops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1565   avg = tot / ((PetscLogDouble)size);
1566   if (min != 0.0) ratio = max / min;
1567   else ratio = 0.0;
1568   PetscCall(PetscFPrintf(comm, fd, "Flops/sec:            %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1569   /*   Memory */
1570   PetscCall(PetscMallocGetMaximumUsage(&mem));
1571   if (mem > 0.0) {
1572     PetscCallMPI(MPI_Allreduce(&mem, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1573     PetscCallMPI(MPI_Allreduce(&mem, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1574     PetscCallMPI(MPI_Allreduce(&mem, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1575     avg = tot / ((PetscLogDouble)size);
1576     if (min != 0.0) ratio = max / min;
1577     else ratio = 0.0;
1578     PetscCall(PetscFPrintf(comm, fd, "Memory (bytes):       %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1579   }
1580   /*   Messages */
1581   mess = 0.5 * (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct);
1582   PetscCallMPI(MPI_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1583   PetscCallMPI(MPI_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1584   PetscCallMPI(MPI_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1585   avg = tot / ((PetscLogDouble)size);
1586   if (min != 0.0) ratio = max / min;
1587   else ratio = 0.0;
1588   PetscCall(PetscFPrintf(comm, fd, "MPI Msg Count:        %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1589   numMessages = tot;
1590   /*   Message Lengths */
1591   mess        = 0.5 * (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len);
1592   PetscCallMPI(MPI_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1593   PetscCallMPI(MPI_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1594   PetscCallMPI(MPI_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1595   if (numMessages != 0) avg = tot / numMessages;
1596   else avg = 0.0;
1597   if (min != 0.0) ratio = max / min;
1598   else ratio = 0.0;
1599   PetscCall(PetscFPrintf(comm, fd, "MPI Msg Len (bytes):  %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot));
1600   messageLength = tot;
1601   /*   Reductions */
1602   PetscCallMPI(MPI_Allreduce(&red, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1603   PetscCallMPI(MPI_Allreduce(&red, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1604   PetscCallMPI(MPI_Allreduce(&red, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1605   if (min != 0.0) ratio = max / min;
1606   else ratio = 0.0;
1607   PetscCall(PetscFPrintf(comm, fd, "MPI Reductions:       %5.3e   %7.3f\n", max, ratio));
1608   numReductions = red; /* wrong because uses count from process zero */
1609   PetscCall(PetscFPrintf(comm, fd, "\nFlop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)\n"));
1610   PetscCall(PetscFPrintf(comm, fd, "                            e.g., VecAXPY() for real vectors of length N --> 2N flops\n"));
1611   PetscCall(PetscFPrintf(comm, fd, "                            and VecAXPY() for complex vectors of length N --> 8N flops\n"));
1612 
1613   /* Get total number of stages --
1614        Currently, a single processor can register more stages than another, but stages must all be registered in order.
1615        We can removed this requirement if necessary by having a global stage numbering and indirection on the stage ID.
1616        This seems best accomplished by assoicating a communicator with each stage.
1617   */
1618   PetscCall(PetscLogGetStageLog(&stageLog));
1619   PetscCallMPI(MPI_Allreduce(&stageLog->numStages, &numStages, 1, MPI_INT, MPI_MAX, comm));
1620   PetscCall(PetscMalloc1(numStages, &localStageUsed));
1621   PetscCall(PetscMalloc1(numStages, &stageUsed));
1622   PetscCall(PetscMalloc1(numStages, &localStageVisible));
1623   PetscCall(PetscMalloc1(numStages, &stageVisible));
1624   if (numStages > 0) {
1625     stageInfo = stageLog->stageInfo;
1626     for (stage = 0; stage < numStages; stage++) {
1627       if (stage < stageLog->numStages) {
1628         localStageUsed[stage]    = stageInfo[stage].used;
1629         localStageVisible[stage] = stageInfo[stage].perfInfo.visible;
1630       } else {
1631         localStageUsed[stage]    = PETSC_FALSE;
1632         localStageVisible[stage] = PETSC_TRUE;
1633       }
1634     }
1635     PetscCallMPI(MPI_Allreduce(localStageUsed, stageUsed, numStages, MPIU_BOOL, MPI_LOR, comm));
1636     PetscCallMPI(MPI_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm));
1637     for (stage = 0; stage < numStages; stage++) {
1638       if (stageUsed[stage]) {
1639         PetscCall(PetscFPrintf(comm, fd, "\nSummary of Stages:   ----- Time ------  ----- Flop ------  --- Messages ---  -- Message Lengths --  -- Reductions --\n"));
1640         PetscCall(PetscFPrintf(comm, fd, "                        Avg     %%Total     Avg     %%Total    Count   %%Total     Avg         %%Total    Count   %%Total\n"));
1641         break;
1642       }
1643     }
1644     for (stage = 0; stage < numStages; stage++) {
1645       if (!stageUsed[stage]) continue;
1646       /* CANNOT use MPI_Allreduce() since it might fail the line number check */
1647       if (localStageUsed[stage]) {
1648         PetscCallMPI(MPI_Allreduce(&stageInfo[stage].perfInfo.time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1649         PetscCallMPI(MPI_Allreduce(&stageInfo[stage].perfInfo.flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1650         PetscCallMPI(MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1651         PetscCallMPI(MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1652         PetscCallMPI(MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1653         name = stageInfo[stage].name;
1654       } else {
1655         PetscCallMPI(MPI_Allreduce(&zero, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1656         PetscCallMPI(MPI_Allreduce(&zero, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1657         PetscCallMPI(MPI_Allreduce(&zero, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1658         PetscCallMPI(MPI_Allreduce(&zero, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1659         PetscCallMPI(MPI_Allreduce(&zero, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1660         name = "";
1661       }
1662       mess *= 0.5;
1663       messLen *= 0.5;
1664       red /= size;
1665       if (TotalTime != 0.0) fracTime = stageTime / TotalTime;
1666       else fracTime = 0.0;
1667       if (TotalFlops != 0.0) fracFlops = flops / TotalFlops;
1668       else fracFlops = 0.0;
1669       /* Talk to Barry if (stageTime     != 0.0) flops          = (size*flops)/stageTime; else flops          = 0.0; */
1670       if (numMessages != 0.0) fracMessages = mess / numMessages;
1671       else fracMessages = 0.0;
1672       if (mess != 0.0) avgMessLen = messLen / mess;
1673       else avgMessLen = 0.0;
1674       if (messageLength != 0.0) fracLength = messLen / messageLength;
1675       else fracLength = 0.0;
1676       if (numReductions != 0.0) fracReductions = red / numReductions;
1677       else fracReductions = 0.0;
1678       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, name, stageTime / size, 100.0 * fracTime, flops, 100.0 * fracFlops, mess, 100.0 * fracMessages, avgMessLen, 100.0 * fracLength, red, 100.0 * fracReductions));
1679     }
1680   }
1681 
1682   PetscCall(PetscFPrintf(comm, fd, "\n------------------------------------------------------------------------------------------------------------------------\n"));
1683   PetscCall(PetscFPrintf(comm, fd, "See the 'Profiling' chapter of the users' manual for details on interpreting output.\n"));
1684   PetscCall(PetscFPrintf(comm, fd, "Phase summary info:\n"));
1685   PetscCall(PetscFPrintf(comm, fd, "   Count: number of times phase was executed\n"));
1686   PetscCall(PetscFPrintf(comm, fd, "   Time and Flop: Max - maximum over all processors\n"));
1687   PetscCall(PetscFPrintf(comm, fd, "                  Ratio - ratio of maximum to minimum over all processors\n"));
1688   PetscCall(PetscFPrintf(comm, fd, "   Mess: number of messages sent\n"));
1689   PetscCall(PetscFPrintf(comm, fd, "   AvgLen: average message length (bytes)\n"));
1690   PetscCall(PetscFPrintf(comm, fd, "   Reduct: number of global reductions\n"));
1691   PetscCall(PetscFPrintf(comm, fd, "   Global: entire computation\n"));
1692   PetscCall(PetscFPrintf(comm, fd, "   Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().\n"));
1693   PetscCall(PetscFPrintf(comm, fd, "      %%T - percent time in this phase         %%F - percent flop in this phase\n"));
1694   PetscCall(PetscFPrintf(comm, fd, "      %%M - percent messages in this phase     %%L - percent message lengths in this phase\n"));
1695   PetscCall(PetscFPrintf(comm, fd, "      %%R - percent reductions in this phase\n"));
1696   PetscCall(PetscFPrintf(comm, fd, "   Total Mflop/s: 10e-6 * (sum of flop over all processors)/(max time over all processors)\n"));
1697   if (PetscLogMemory) {
1698     PetscCall(PetscFPrintf(comm, fd, "   Malloc Mbytes: Memory allocated and kept during event (sum over all calls to event). May be negative\n"));
1699     PetscCall(PetscFPrintf(comm, fd, "   EMalloc Mbytes: extra memory allocated during event and then freed (maximum over all calls to events). Never negative\n"));
1700     PetscCall(PetscFPrintf(comm, fd, "   MMalloc Mbytes: Increase in high water mark of allocated memory (sum over all calls to event). Never negative\n"));
1701     PetscCall(PetscFPrintf(comm, fd, "   RMI Mbytes: Increase in resident memory (sum over all calls to event)\n"));
1702   }
1703 #if defined(PETSC_HAVE_DEVICE)
1704   PetscCall(PetscFPrintf(comm, fd, "   GPU Mflop/s: 10e-6 * (sum of flop on GPU over all processors)/(max GPU time over all processors)\n"));
1705   PetscCall(PetscFPrintf(comm, fd, "   CpuToGpu Count: total number of CPU to GPU copies per processor\n"));
1706   PetscCall(PetscFPrintf(comm, fd, "   CpuToGpu Size (Mbytes): 10e-6 * (total size of CPU to GPU copies per processor)\n"));
1707   PetscCall(PetscFPrintf(comm, fd, "   GpuToCpu Count: total number of GPU to CPU copies per processor\n"));
1708   PetscCall(PetscFPrintf(comm, fd, "   GpuToCpu Size (Mbytes): 10e-6 * (total size of GPU to CPU copies per processor)\n"));
1709   PetscCall(PetscFPrintf(comm, fd, "   GPU %%F: percent flops on GPU in this event\n"));
1710 #endif
1711   PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------\n"));
1712 
1713   PetscCall(PetscLogViewWarnDebugging(comm, fd));
1714 
1715   /* Report events */
1716   PetscCall(PetscFPrintf(comm, fd, "Event                Count      Time (sec)     Flop                              --- Global ---  --- Stage ----  Total"));
1717   if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, "  Malloc EMalloc MMalloc RMI"));
1718 #if defined(PETSC_HAVE_DEVICE)
1719   PetscCall(PetscFPrintf(comm, fd, "   GPU    - CpuToGpu -   - GpuToCpu - GPU"));
1720 #endif
1721   PetscCall(PetscFPrintf(comm, fd, "\n"));
1722   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"));
1723   if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, " Mbytes Mbytes Mbytes Mbytes"));
1724 #if defined(PETSC_HAVE_DEVICE)
1725   PetscCall(PetscFPrintf(comm, fd, " Mflop/s Count   Size   Count   Size  %%F"));
1726 #endif
1727   PetscCall(PetscFPrintf(comm, fd, "\n"));
1728   PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------"));
1729   if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, "-----------------------------"));
1730 #if defined(PETSC_HAVE_DEVICE)
1731   PetscCall(PetscFPrintf(comm, fd, "---------------------------------------"));
1732 #endif
1733   PetscCall(PetscFPrintf(comm, fd, "\n"));
1734 
1735 #if defined(PETSC_HAVE_DEVICE)
1736   /* this indirect way of accessing these values is needed when PETSc is build with multiple libraries since the symbols are not in libpetscsys */
1737   PetscCall(PetscEventRegLogGetEvent(stageLog->eventLog, "TAOSolve", &TAO_Solve));
1738   PetscCall(PetscEventRegLogGetEvent(stageLog->eventLog, "TSStep", &TS_Step));
1739   PetscCall(PetscEventRegLogGetEvent(stageLog->eventLog, "SNESSolve", &SNES_Solve));
1740   PetscCall(PetscEventRegLogGetEvent(stageLog->eventLog, "KSPSolve", &KSP_Solve));
1741 #endif
1742 
1743   /* Problem: The stage name will not show up unless the stage executed on proc 1 */
1744   for (stage = 0; stage < numStages; stage++) {
1745     if (!stageVisible[stage]) continue;
1746     /* CANNOT use MPI_Allreduce() since it might fail the line number check */
1747     if (localStageUsed[stage]) {
1748       PetscCall(PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stageInfo[stage].name));
1749       PetscCallMPI(MPI_Allreduce(&stageInfo[stage].perfInfo.time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1750       PetscCallMPI(MPI_Allreduce(&stageInfo[stage].perfInfo.flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1751       PetscCallMPI(MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1752       PetscCallMPI(MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1753       PetscCallMPI(MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1754     } else {
1755       PetscCall(PetscFPrintf(comm, fd, "\n--- Event Stage %d: Unknown\n\n", stage));
1756       PetscCallMPI(MPI_Allreduce(&zero, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1757       PetscCallMPI(MPI_Allreduce(&zero, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1758       PetscCallMPI(MPI_Allreduce(&zero, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1759       PetscCallMPI(MPI_Allreduce(&zero, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1760       PetscCallMPI(MPI_Allreduce(&zero, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1761     }
1762     mess *= 0.5;
1763     messLen *= 0.5;
1764     red /= size;
1765 
1766     /* Get total number of events in this stage --
1767        Currently, a single processor can register more events than another, but events must all be registered in order,
1768        just like stages. We can removed this requirement if necessary by having a global event numbering and indirection
1769        on the event ID. This seems best accomplished by associating a communicator with each stage.
1770 
1771        Problem: If the event did not happen on proc 1, its name will not be available.
1772        Problem: Event visibility is not implemented
1773     */
1774     if (localStageUsed[stage]) {
1775       eventInfo      = stageLog->stageInfo[stage].eventLog->eventInfo;
1776       localNumEvents = stageLog->stageInfo[stage].eventLog->numEvents;
1777     } else localNumEvents = 0;
1778     PetscCallMPI(MPI_Allreduce(&localNumEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm));
1779     for (event = 0; event < numEvents; event++) {
1780       /* CANNOT use MPI_Allreduce() since it might fail the line number check */
1781       if (localStageUsed[stage] && (event < stageLog->stageInfo[stage].eventLog->numEvents) && (eventInfo[event].depth == 0)) {
1782         if ((eventInfo[event].count > 0) && (eventInfo[event].time > 0.0)) flopr = eventInfo[event].flops;
1783         else flopr = 0.0;
1784         PetscCallMPI(MPI_Allreduce(&flopr, &minf, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1785         PetscCallMPI(MPI_Allreduce(&flopr, &maxf, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1786         PetscCallMPI(MPI_Allreduce(&eventInfo[event].flops, &totf, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1787         PetscCallMPI(MPI_Allreduce(&eventInfo[event].time, &mint, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1788         PetscCallMPI(MPI_Allreduce(&eventInfo[event].time, &maxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1789         PetscCallMPI(MPI_Allreduce(&eventInfo[event].time, &tott, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1790         PetscCallMPI(MPI_Allreduce(&eventInfo[event].numMessages, &totm, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1791         PetscCallMPI(MPI_Allreduce(&eventInfo[event].messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1792         PetscCallMPI(MPI_Allreduce(&eventInfo[event].numReductions, &totr, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1793         PetscCallMPI(MPI_Allreduce(&eventInfo[event].count, &minC, 1, MPI_INT, MPI_MIN, comm));
1794         PetscCallMPI(MPI_Allreduce(&eventInfo[event].count, &maxC, 1, MPI_INT, MPI_MAX, comm));
1795         if (PetscLogMemory) {
1796           PetscCallMPI(MPI_Allreduce(&eventInfo[event].memIncrease, &mem, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1797           PetscCallMPI(MPI_Allreduce(&eventInfo[event].mallocSpace, &mal, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1798           PetscCallMPI(MPI_Allreduce(&eventInfo[event].mallocIncrease, &malmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1799           PetscCallMPI(MPI_Allreduce(&eventInfo[event].mallocIncreaseEvent, &emalmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1800         }
1801 #if defined(PETSC_HAVE_DEVICE)
1802         PetscCallMPI(MPI_Allreduce(&eventInfo[event].CpuToGpuCount, &cct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1803         PetscCallMPI(MPI_Allreduce(&eventInfo[event].GpuToCpuCount, &gct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1804         PetscCallMPI(MPI_Allreduce(&eventInfo[event].CpuToGpuSize, &csz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1805         PetscCallMPI(MPI_Allreduce(&eventInfo[event].GpuToCpuSize, &gsz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1806         PetscCallMPI(MPI_Allreduce(&eventInfo[event].GpuFlops, &gflops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1807         PetscCallMPI(MPI_Allreduce(&eventInfo[event].GpuTime, &gmaxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1808 #endif
1809         name = stageLog->eventLog->eventInfo[event].name;
1810       } else {
1811         flopr = 0.0;
1812         PetscCallMPI(MPI_Allreduce(&flopr, &minf, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1813         PetscCallMPI(MPI_Allreduce(&flopr, &maxf, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1814         PetscCallMPI(MPI_Allreduce(&zero, &totf, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1815         PetscCallMPI(MPI_Allreduce(&zero, &mint, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1816         PetscCallMPI(MPI_Allreduce(&zero, &maxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1817         PetscCallMPI(MPI_Allreduce(&zero, &tott, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1818         PetscCallMPI(MPI_Allreduce(&zero, &totm, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1819         PetscCallMPI(MPI_Allreduce(&zero, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1820         PetscCallMPI(MPI_Allreduce(&zero, &totr, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1821         PetscCallMPI(MPI_Allreduce(&ierr, &minC, 1, MPI_INT, MPI_MIN, comm));
1822         PetscCallMPI(MPI_Allreduce(&ierr, &maxC, 1, MPI_INT, MPI_MAX, comm));
1823         if (PetscLogMemory) {
1824           PetscCallMPI(MPI_Allreduce(&zero, &mem, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1825           PetscCallMPI(MPI_Allreduce(&zero, &mal, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1826           PetscCallMPI(MPI_Allreduce(&zero, &malmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1827           PetscCallMPI(MPI_Allreduce(&zero, &emalmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1828         }
1829 #if defined(PETSC_HAVE_DEVICE)
1830         PetscCallMPI(MPI_Allreduce(&zero, &cct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1831         PetscCallMPI(MPI_Allreduce(&zero, &gct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1832         PetscCallMPI(MPI_Allreduce(&zero, &csz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1833         PetscCallMPI(MPI_Allreduce(&zero, &gsz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1834         PetscCallMPI(MPI_Allreduce(&zero, &gflops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1835         PetscCallMPI(MPI_Allreduce(&zero, &gmaxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1836 #endif
1837         name = "";
1838       }
1839       if (mint < 0.0) {
1840         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, name));
1841         mint = 0;
1842       }
1843       PetscCheck(minf >= 0.0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Minimum flop %g over all processors for %s is negative! Not possible!", minf, name);
1844 /* Put NaN into the time for all events that may not be time accurately since they may happen asynchronously on the GPU */
1845 #if defined(PETSC_HAVE_DEVICE)
1846       if (!PetscLogGpuTimeFlag && petsc_gflops > 0) {
1847         memcpy(&gmaxt, &nas, sizeof(PetscLogDouble));
1848         PetscCall(PetscEventRegLogGetEvent(stageLog->eventLog, name, &eventid));
1849         if (eventid != SNES_Solve && eventid != KSP_Solve && eventid != TS_Step && eventid != TAO_Solve) {
1850           memcpy(&mint, &nas, sizeof(PetscLogDouble));
1851           memcpy(&maxt, &nas, sizeof(PetscLogDouble));
1852         }
1853       }
1854 #endif
1855       totm *= 0.5;
1856       totml *= 0.5;
1857       totr /= size;
1858 
1859       if (maxC != 0) {
1860         if (minC != 0) ratC = ((PetscLogDouble)maxC) / minC;
1861         else ratC = 0.0;
1862         if (mint != 0.0) ratt = maxt / mint;
1863         else ratt = 0.0;
1864         if (minf != 0.0) ratf = maxf / minf;
1865         else ratf = 0.0;
1866         if (TotalTime != 0.0) fracTime = tott / TotalTime;
1867         else fracTime = 0.0;
1868         if (TotalFlops != 0.0) fracFlops = totf / TotalFlops;
1869         else fracFlops = 0.0;
1870         if (stageTime != 0.0) fracStageTime = tott / stageTime;
1871         else fracStageTime = 0.0;
1872         if (flops != 0.0) fracStageFlops = totf / flops;
1873         else fracStageFlops = 0.0;
1874         if (numMessages != 0.0) fracMess = totm / numMessages;
1875         else fracMess = 0.0;
1876         if (messageLength != 0.0) fracMessLen = totml / messageLength;
1877         else fracMessLen = 0.0;
1878         if (numReductions != 0.0) fracRed = totr / numReductions;
1879         else fracRed = 0.0;
1880         if (mess != 0.0) fracStageMess = totm / mess;
1881         else fracStageMess = 0.0;
1882         if (messLen != 0.0) fracStageMessLen = totml / messLen;
1883         else fracStageMessLen = 0.0;
1884         if (red != 0.0) fracStageRed = totr / red;
1885         else fracStageRed = 0.0;
1886         if (totm != 0.0) totml /= totm;
1887         else totml = 0.0;
1888         if (maxt != 0.0) flopr = totf / maxt;
1889         else flopr = 0.0;
1890         if (fracStageTime > 1.00) PetscCall(PetscFPrintf(comm, fd, "Warning -- total time of event greater than time of entire stage -- something is wrong with the timer\n"));
1891         PetscCall(PetscFPrintf(comm, fd, "%-16s %7d%4.1f %5.4e%4.1f %3.2e%4.1f %2.1e %2.1e %2.1e%3.0f%3.0f%3.0f%3.0f%3.0f %3.0f%3.0f%3.0f%3.0f%3.0f %5.0f", 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));
1892         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));
1893 #if defined(PETSC_HAVE_DEVICE)
1894         if (totf != 0.0) fracgflops = gflops / totf;
1895         else fracgflops = 0.0;
1896         if (gmaxt != 0.0) gflopr = gflops / gmaxt;
1897         else gflopr = 0.0;
1898         PetscCall(PetscFPrintf(comm, fd, "   %5.0f   %4.0f %3.2e %4.0f %3.2e% 3.0f", PetscAbs(gflopr) / 1.0e6, cct / size, csz / (1.0e6 * size), gct / size, gsz / (1.0e6 * size), 100.0 * fracgflops));
1899 #endif
1900         PetscCall(PetscFPrintf(comm, fd, "\n"));
1901       }
1902     }
1903   }
1904 
1905   /* Memory usage and object creation */
1906   PetscCall(PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------"));
1907   if (PetscLogMemory) PetscCall(PetscFPrintf(comm, fd, "-----------------------------"));
1908 #if defined(PETSC_HAVE_DEVICE)
1909   PetscCall(PetscFPrintf(comm, fd, "---------------------------------------"));
1910 #endif
1911   PetscCall(PetscFPrintf(comm, fd, "\n"));
1912   PetscCall(PetscFPrintf(comm, fd, "\n"));
1913 
1914   /* Right now, only stages on the first processor are reported here, meaning only objects associated with
1915      the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then
1916      stats for stages local to processor sets.
1917   */
1918   /* We should figure out the longest object name here (now 20 characters) */
1919   PetscCall(PetscFPrintf(comm, fd, "Object Type          Creations   Destructions. Reports information only for process 0.\n"));
1920   for (stage = 0; stage < numStages; stage++) {
1921     if (localStageUsed[stage]) {
1922       classInfo = stageLog->stageInfo[stage].classLog->classInfo;
1923       PetscCall(PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stageInfo[stage].name));
1924       for (oclass = 0; oclass < stageLog->stageInfo[stage].classLog->numClasses; oclass++) {
1925         if ((classInfo[oclass].creations > 0) || (classInfo[oclass].destructions > 0)) {
1926           PetscCall(PetscFPrintf(comm, fd, "%20s %5d          %5d\n", stageLog->classLog->classInfo[oclass].name, classInfo[oclass].creations, classInfo[oclass].destructions));
1927         }
1928       }
1929     } else {
1930       if (!localStageVisible[stage]) continue;
1931       PetscCall(PetscFPrintf(comm, fd, "\n--- Event Stage %d: Unknown\n\n", stage));
1932     }
1933   }
1934 
1935   PetscCall(PetscFree(localStageUsed));
1936   PetscCall(PetscFree(stageUsed));
1937   PetscCall(PetscFree(localStageVisible));
1938   PetscCall(PetscFree(stageVisible));
1939 
1940   /* Information unrelated to this particular run */
1941   PetscCall(PetscFPrintf(comm, fd, "========================================================================================================================\n"));
1942   PetscTime(&y);
1943   PetscTime(&x);
1944   PetscTime(&y);
1945   PetscTime(&y);
1946   PetscTime(&y);
1947   PetscTime(&y);
1948   PetscTime(&y);
1949   PetscTime(&y);
1950   PetscTime(&y);
1951   PetscTime(&y);
1952   PetscTime(&y);
1953   PetscTime(&y);
1954   PetscCall(PetscFPrintf(comm, fd, "Average time to get PetscTime(): %g\n", (y - x) / 10.0));
1955   /* MPI information */
1956   if (size > 1) {
1957     MPI_Status  status;
1958     PetscMPIInt tag;
1959     MPI_Comm    newcomm;
1960 
1961     PetscCallMPI(MPI_Barrier(comm));
1962     PetscTime(&x);
1963     PetscCallMPI(MPI_Barrier(comm));
1964     PetscCallMPI(MPI_Barrier(comm));
1965     PetscCallMPI(MPI_Barrier(comm));
1966     PetscCallMPI(MPI_Barrier(comm));
1967     PetscCallMPI(MPI_Barrier(comm));
1968     PetscTime(&y);
1969     PetscCall(PetscFPrintf(comm, fd, "Average time for MPI_Barrier(): %g\n", (y - x) / 5.0));
1970     PetscCall(PetscCommDuplicate(comm, &newcomm, &tag));
1971     PetscCallMPI(MPI_Barrier(comm));
1972     if (rank) {
1973       PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, rank - 1, tag, newcomm, &status));
1974       PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, (rank + 1) % size, tag, newcomm));
1975     } else {
1976       PetscTime(&x);
1977       PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, 1, tag, newcomm));
1978       PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, size - 1, tag, newcomm, &status));
1979       PetscTime(&y);
1980       PetscCall(PetscFPrintf(comm, fd, "Average time for zero size MPI_Send(): %g\n", (y - x) / size));
1981     }
1982     PetscCall(PetscCommDestroy(&newcomm));
1983   }
1984   PetscCall(PetscOptionsView(NULL, viewer));
1985 
1986   /* Machine and compile information */
1987 #if defined(PETSC_USE_FORTRAN_KERNELS)
1988   PetscCall(PetscFPrintf(comm, fd, "Compiled with FORTRAN kernels\n"));
1989 #else
1990   PetscCall(PetscFPrintf(comm, fd, "Compiled without FORTRAN kernels\n"));
1991 #endif
1992 #if defined(PETSC_USE_64BIT_INDICES)
1993   PetscCall(PetscFPrintf(comm, fd, "Compiled with 64 bit PetscInt\n"));
1994 #elif defined(PETSC_USE___FLOAT128)
1995   PetscCall(PetscFPrintf(comm, fd, "Compiled with 32 bit PetscInt\n"));
1996 #endif
1997 #if defined(PETSC_USE_REAL_SINGLE)
1998   PetscCall(PetscFPrintf(comm, fd, "Compiled with single precision PetscScalar and PetscReal\n"));
1999 #elif defined(PETSC_USE___FLOAT128)
2000   PetscCall(PetscFPrintf(comm, fd, "Compiled with 128 bit precision PetscScalar and PetscReal\n"));
2001 #endif
2002 #if defined(PETSC_USE_REAL_MAT_SINGLE)
2003   PetscCall(PetscFPrintf(comm, fd, "Compiled with single precision matrices\n"));
2004 #else
2005   PetscCall(PetscFPrintf(comm, fd, "Compiled with full precision matrices (default)\n"));
2006 #endif
2007   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)));
2008 
2009   PetscCall(PetscFPrintf(comm, fd, "Configure options: %s", petscconfigureoptions));
2010   PetscCall(PetscFPrintf(comm, fd, "%s", petscmachineinfo));
2011   PetscCall(PetscFPrintf(comm, fd, "%s", petsccompilerinfo));
2012   PetscCall(PetscFPrintf(comm, fd, "%s", petsccompilerflagsinfo));
2013   PetscCall(PetscFPrintf(comm, fd, "%s", petsclinkerinfo));
2014 
2015   /* Cleanup */
2016   PetscCall(PetscFPrintf(comm, fd, "\n"));
2017   PetscCall(PetscLogViewWarnNoGpuAwareMpi(comm, fd));
2018   PetscCall(PetscLogViewWarnDebugging(comm, fd));
2019   PetscCall(PetscFPTrapPop());
2020   PetscFunctionReturn(0);
2021 }
2022 
2023 /*@C
2024   PetscLogView - Prints a summary of the logging.
2025 
2026   Collective over MPI_Comm
2027 
2028   Input Parameter:
2029 .  viewer - an ASCII viewer
2030 
2031   Options Database Keys:
2032 +  -log_view [:filename] - Prints summary of log information
2033 .  -log_view :filename.py:ascii_info_detail - Saves logging information from each process as a Python file
2034 .  -log_view :filename.xml:ascii_xml - Saves a summary of the logging information in a nested format (see below for how to view it)
2035 .  -log_view :filename.txt:ascii_flamegraph - Saves logging information in a format suitable for visualising as a Flame Graph (see below for how to view it)
2036 .  -log_view_memory - Also display memory usage in each event
2037 .  -log_view_gpu_time - Also display time in each event for GPU kernels (Note this may slow the computation)
2038 .  -log_all - Saves a file Log.rank for each MPI process with details of each step of the computation
2039 -  -log_trace [filename] - Displays a trace of what each process is doing
2040 
2041   Notes:
2042   It is possible to control the logging programatically but we recommend using the options database approach whenever possible
2043   By default the summary is printed to stdout.
2044 
2045   Before calling this routine you must have called either PetscLogDefaultBegin() or PetscLogNestedBegin()
2046 
2047   If PETSc is configured with --with-logging=0 then this functionality is not available
2048 
2049   To view the nested XML format filename.xml first copy  ${PETSC_DIR}/share/petsc/xml/performance_xml2html.xsl to the current
2050   directory then open filename.xml with your browser. Specific notes for certain browsers
2051 $    Firefox and Internet explorer - simply open the file
2052 $    Google Chrome - you must start up Chrome with the option --allow-file-access-from-files
2053 $    Safari - see https://ccm.net/faq/36342-safari-how-to-enable-local-file-access
2054   or one can use the package http://xmlsoft.org/XSLT/xsltproc2.html to translate the xml file to html and then open it with
2055   your browser.
2056   Alternatively, use the script ${PETSC_DIR}/lib/petsc/bin/petsc-performance-view to automatically open a new browser
2057   window and render the XML log file contents.
2058 
2059   The nested XML format was kindly donated by Koos Huijssen and Christiaan M. Klaij  MARITIME  RESEARCH  INSTITUTE  NETHERLANDS
2060 
2061   The Flame Graph output can be visualised using either the original Flame Graph script (https://github.com/brendangregg/FlameGraph)
2062   or using speedscope (https://www.speedscope.app).
2063   Old XML profiles may be converted into this format using the script ${PETSC_DIR}/lib/petsc/bin/xml2flamegraph.py.
2064 
2065   Level: beginner
2066 
2067 .seealso: `PetscLogDefaultBegin()`, `PetscLogDump()`
2068 @*/
2069 PetscErrorCode PetscLogView(PetscViewer viewer) {
2070   PetscBool         isascii;
2071   PetscViewerFormat format;
2072   int               stage, lastStage;
2073   PetscStageLog     stageLog;
2074 
2075   PetscFunctionBegin;
2076   PetscCheck(PetscLogPLB, PETSC_COMM_SELF, PETSC_ERR_SUP, "Must use -log_view or PetscLogDefaultBegin() before calling this routine");
2077   /* Pop off any stages the user forgot to remove */
2078   lastStage = 0;
2079   PetscCall(PetscLogGetStageLog(&stageLog));
2080   PetscCall(PetscStageLogGetCurrent(stageLog, &stage));
2081   while (stage >= 0) {
2082     lastStage = stage;
2083     PetscCall(PetscStageLogPop(stageLog));
2084     PetscCall(PetscStageLogGetCurrent(stageLog, &stage));
2085   }
2086   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2087   PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Currently can only view logging to ASCII");
2088   PetscCall(PetscViewerGetFormat(viewer, &format));
2089   if (format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_ASCII_INFO) {
2090     PetscCall(PetscLogView_Default(viewer));
2091   } else if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
2092     PetscCall(PetscLogView_Detailed(viewer));
2093   } else if (format == PETSC_VIEWER_ASCII_CSV) {
2094     PetscCall(PetscLogView_CSV(viewer));
2095   } else if (format == PETSC_VIEWER_ASCII_XML) {
2096     PetscCall(PetscLogView_Nested(viewer));
2097   } else if (format == PETSC_VIEWER_ASCII_FLAMEGRAPH) {
2098     PetscCall(PetscLogView_Flamegraph(viewer));
2099   }
2100   PetscCall(PetscStageLogPush(stageLog, lastStage));
2101   PetscFunctionReturn(0);
2102 }
2103 
2104 /*@C
2105   PetscLogViewFromOptions - Processes command line options to determine if/how a PetscLog is to be viewed.
2106 
2107   Collective on PETSC_COMM_WORLD
2108 
2109   Not normally called by user
2110 
2111   Level: intermediate
2112 
2113 @*/
2114 PetscErrorCode PetscLogViewFromOptions(void) {
2115   PetscViewer       viewer;
2116   PetscBool         flg;
2117   PetscViewerFormat format;
2118 
2119   PetscFunctionBegin;
2120   PetscCall(PetscOptionsGetViewer(PETSC_COMM_WORLD, NULL, NULL, "-log_view", &viewer, &format, &flg));
2121   if (flg) {
2122     PetscCall(PetscViewerPushFormat(viewer, format));
2123     PetscCall(PetscLogView(viewer));
2124     PetscCall(PetscViewerPopFormat(viewer));
2125     PetscCall(PetscViewerDestroy(&viewer));
2126   }
2127   PetscFunctionReturn(0);
2128 }
2129 
2130 /*----------------------------------------------- Counter Functions -------------------------------------------------*/
2131 /*@C
2132    PetscGetFlops - Returns the number of flops used on this processor
2133    since the program began.
2134 
2135    Not Collective
2136 
2137    Output Parameter:
2138    flops - number of floating point operations
2139 
2140    Notes:
2141    A global counter logs all PETSc flop counts.  The user can use
2142    PetscLogFlops() to increment this counter to include flops for the
2143    application code.
2144 
2145    Level: intermediate
2146 
2147 .seealso: `PetscTime()`, `PetscLogFlops()`
2148 @*/
2149 PetscErrorCode PetscGetFlops(PetscLogDouble *flops) {
2150   PetscFunctionBegin;
2151   *flops = petsc_TotalFlops;
2152   PetscFunctionReturn(0);
2153 }
2154 
2155 PetscErrorCode PetscLogObjectState(PetscObject obj, const char format[], ...) {
2156   size_t  fullLength;
2157   va_list Argp;
2158 
2159   PetscFunctionBegin;
2160   if (!petsc_logObjects) PetscFunctionReturn(0);
2161   va_start(Argp, format);
2162   PetscCall(PetscVSNPrintf(petsc_objects[obj->id].info, 64, format, &fullLength, Argp));
2163   va_end(Argp);
2164   PetscFunctionReturn(0);
2165 }
2166 
2167 /*MC
2168    PetscLogFlops - Adds floating point operations to the global counter.
2169 
2170    Synopsis:
2171    #include <petsclog.h>
2172    PetscErrorCode PetscLogFlops(PetscLogDouble f)
2173 
2174    Not Collective
2175 
2176    Input Parameter:
2177 .  f - flop counter
2178 
2179    Usage:
2180 .vb
2181      PetscLogEvent USER_EVENT;
2182      PetscLogEventRegister("User event",0,&USER_EVENT);
2183      PetscLogEventBegin(USER_EVENT,0,0,0,0);
2184         [code segment to monitor]
2185         PetscLogFlops(user_flops)
2186      PetscLogEventEnd(USER_EVENT,0,0,0,0);
2187 .ve
2188 
2189    Notes:
2190    A global counter logs all PETSc flop counts.  The user can use
2191    PetscLogFlops() to increment this counter to include flops for the
2192    application code.
2193 
2194    Level: intermediate
2195 
2196 .seealso: `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscGetFlops()`
2197 
2198 M*/
2199 
2200 /*MC
2201    PetscPreLoadBegin - Begin a segment of code that may be preloaded (run twice)
2202     to get accurate timings
2203 
2204    Synopsis:
2205    #include <petsclog.h>
2206    void PetscPreLoadBegin(PetscBool  flag,char *name);
2207 
2208    Not Collective
2209 
2210    Input Parameters:
2211 +   flag - PETSC_TRUE to run twice, PETSC_FALSE to run once, may be overridden
2212            with command line option -preload true or -preload false
2213 -   name - name of first stage (lines of code timed separately with -log_view) to
2214            be preloaded
2215 
2216    Usage:
2217 .vb
2218      PetscPreLoadBegin(PETSC_TRUE,"first stage);
2219        lines of code
2220        PetscPreLoadStage("second stage");
2221        lines of code
2222      PetscPreLoadEnd();
2223 .ve
2224 
2225    Notes:
2226     Only works in C/C++, not Fortran
2227 
2228      Flags available within the macro.
2229 +    PetscPreLoadingUsed - true if we are or have done preloading
2230 .    PetscPreLoadingOn - true if it is CURRENTLY doing preload
2231 .    PetscPreLoadIt - 0 for the first computation (with preloading turned off it is only 0) 1 for the second
2232 -    PetscPreLoadMax - number of times it will do the computation, only one when preloading is turned on
2233      The first two variables are available throughout the program, the second two only between the PetscPreLoadBegin()
2234      and PetscPreLoadEnd()
2235 
2236    Level: intermediate
2237 
2238 .seealso: `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadEnd()`, `PetscPreLoadStage()`
2239 
2240 M*/
2241 
2242 /*MC
2243    PetscPreLoadEnd - End a segment of code that may be preloaded (run twice)
2244     to get accurate timings
2245 
2246    Synopsis:
2247    #include <petsclog.h>
2248    void PetscPreLoadEnd(void);
2249 
2250    Not Collective
2251 
2252    Usage:
2253 .vb
2254      PetscPreLoadBegin(PETSC_TRUE,"first stage);
2255        lines of code
2256        PetscPreLoadStage("second stage");
2257        lines of code
2258      PetscPreLoadEnd();
2259 .ve
2260 
2261    Notes:
2262     only works in C/C++ not fortran
2263 
2264    Level: intermediate
2265 
2266 .seealso: `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadBegin()`, `PetscPreLoadStage()`
2267 
2268 M*/
2269 
2270 /*MC
2271    PetscPreLoadStage - Start a new segment of code to be timed separately.
2272     to get accurate timings
2273 
2274    Synopsis:
2275    #include <petsclog.h>
2276    void PetscPreLoadStage(char *name);
2277 
2278    Not Collective
2279 
2280    Usage:
2281 .vb
2282      PetscPreLoadBegin(PETSC_TRUE,"first stage);
2283        lines of code
2284        PetscPreLoadStage("second stage");
2285        lines of code
2286      PetscPreLoadEnd();
2287 .ve
2288 
2289    Notes:
2290     only works in C/C++ not fortran
2291 
2292    Level: intermediate
2293 
2294 .seealso: `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`
2295 
2296 M*/
2297 
2298 #if PetscDefined(HAVE_DEVICE)
2299 #include <petsc/private/deviceimpl.h>
2300 
2301 PetscBool PetscLogGpuTimeFlag = PETSC_FALSE;
2302 
2303 /*
2304    This cannot be called by users between PetscInitialize() and PetscFinalize() at any random location in the code
2305    because it will result in timing results that cannot be interpreted.
2306 */
2307 static PetscErrorCode PetscLogGpuTime_Off(void) {
2308   PetscLogGpuTimeFlag = PETSC_FALSE;
2309   return 0;
2310 }
2311 
2312 /*@C
2313      PetscLogGpuTime - turn on the logging of GPU time for GPU kernels
2314 
2315   Options Database:
2316 .   -log_view_gpu_time - provide the GPU times in the -log_view output
2317 
2318   Notes:
2319     Because the logging of GPU time requires blocking the CPU execution for each kernel, turning on the timing of the
2320     GPU kernels can slow down the entire computation and should only be used when studying the performance
2321     of operations on GPU such as vector operations and matrix-vector operations.
2322 
2323     This routine should only be called once near the beginning of the program. Once it is started it cannot be turned off.
2324 
2325    Level: advanced
2326 
2327 .seealso: `PetscLogView()`, `PetscLogGpuFlops()`, `PetscLogGpuTimeEnd()`, `PetscLogGpuTimeBegin()`
2328 @*/
2329 PetscErrorCode PetscLogGpuTime(void) {
2330   if (!PetscLogGpuTimeFlag) PetscCall(PetscRegisterFinalize(PetscLogGpuTime_Off));
2331   PetscLogGpuTimeFlag = PETSC_TRUE;
2332   return 0;
2333 }
2334 
2335 /*@C
2336   PetscLogGpuTimeBegin - Start timer for device
2337 
2338   Notes:
2339     When CUDA or HIP is enabled, the timer is run on the GPU, it is a separate logging of time devoted to GPU computations (excluding kernel launch times).
2340     When CUDA or HIP is not available, the timer is run on the CPU, it is a separate logging of time devoted to GPU computations (including kernel launch times).
2341     There is no need to call WaitForCUDA() or WaitForHIP() between PetscLogGpuTimeBegin and PetscLogGpuTimeEnd
2342     This timer should NOT include times for data transfers between the GPU and CPU, nor setup actions such as allocating space.
2343     The regular logging captures the time for data transfers and any CPU activites during the event
2344     It is used to compute the flop rate on the GPU as it is actively engaged in running a kernel.
2345 
2346   Developer Notes:
2347     The GPU event timer captures the execution time of all the kernels launched in the default stream by the CPU between PetscLogGpuTimeBegin() and PetsLogGpuTimeEnd().
2348     PetscLogGpuTimeBegin() and PetsLogGpuTimeEnd() insert the begin and end events into the default stream (stream 0). The device will record a time stamp for the event when it reaches that event in the stream. The function xxxEventSynchronize() is called in PetsLogGpuTimeEnd() to block CPU execution, but not continued GPU excution, until the timer event is recorded.
2349 
2350   Level: intermediate
2351 
2352 .seealso: `PetscLogView()`, `PetscLogGpuFlops()`, `PetscLogGpuTimeEnd()`, `PetscLogGpuTime()`
2353 @*/
2354 PetscErrorCode PetscLogGpuTimeBegin(void) {
2355   PetscFunctionBegin;
2356   if (!PetscLogPLB || !PetscLogGpuTimeFlag) PetscFunctionReturn(0);
2357   if (PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP)) {
2358     PetscDeviceContext dctx;
2359 
2360     PetscCall(PetscDeviceContextGetCurrentContext(&dctx));
2361     PetscCall(PetscDeviceContextBeginTimer_Internal(dctx));
2362   } else {
2363     PetscCall(PetscTimeSubtract(&petsc_gtime));
2364   }
2365   PetscFunctionReturn(0);
2366 }
2367 
2368 /*@C
2369   PetscLogGpuTimeEnd - Stop timer for device
2370 
2371   Level: intermediate
2372 
2373 .seealso: `PetscLogView()`, `PetscLogGpuFlops()`, `PetscLogGpuTimeBegin()`
2374 @*/
2375 PetscErrorCode PetscLogGpuTimeEnd(void) {
2376   PetscFunctionBegin;
2377   if (!PetscLogPLE || !PetscLogGpuTimeFlag) PetscFunctionReturn(0);
2378   if (PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP)) {
2379     PetscDeviceContext dctx;
2380     PetscLogDouble     elapsed;
2381 
2382     PetscCall(PetscDeviceContextGetCurrentContext(&dctx));
2383     PetscCall(PetscDeviceContextEndTimer_Internal(dctx, &elapsed));
2384     petsc_gtime += (elapsed / 1000.0);
2385   } else {
2386     PetscCall(PetscTimeAdd(&petsc_gtime));
2387   }
2388   PetscFunctionReturn(0);
2389 }
2390 #endif /* end of PETSC_HAVE_DEVICE */
2391 
2392 #else /* end of -DPETSC_USE_LOG section */
2393 
2394 PetscErrorCode PetscLogObjectState(PetscObject obj, const char format[], ...) {
2395   PetscFunctionBegin;
2396   PetscFunctionReturn(0);
2397 }
2398 
2399 #endif /* PETSC_USE_LOG*/
2400 
2401 PetscClassId PETSC_LARGEST_CLASSID = PETSC_SMALLEST_CLASSID;
2402 PetscClassId PETSC_OBJECT_CLASSID  = 0;
2403 
2404 /*@C
2405   PetscClassIdRegister - Registers a new class name for objects and logging operations in an application code.
2406 
2407   Not Collective
2408 
2409   Input Parameter:
2410 . name   - The class name
2411 
2412   Output Parameter:
2413 . oclass - The class id or classid
2414 
2415   Level: developer
2416 
2417 @*/
2418 PetscErrorCode PetscClassIdRegister(const char name[], PetscClassId *oclass) {
2419 #if defined(PETSC_USE_LOG)
2420   PetscStageLog stageLog;
2421   PetscInt      stage;
2422 #endif
2423 
2424   PetscFunctionBegin;
2425   *oclass = ++PETSC_LARGEST_CLASSID;
2426 #if defined(PETSC_USE_LOG)
2427   PetscCall(PetscLogGetStageLog(&stageLog));
2428   PetscCall(PetscClassRegLogRegister(stageLog->classLog, name, *oclass));
2429   for (stage = 0; stage < stageLog->numStages; stage++) PetscCall(PetscClassPerfLogEnsureSize(stageLog->stageInfo[stage].classLog, stageLog->classLog->numClasses));
2430 #endif
2431   PetscFunctionReturn(0);
2432 }
2433 
2434 #if defined(PETSC_USE_LOG) && defined(PETSC_HAVE_MPE)
2435 #include <mpe.h>
2436 
2437 PetscBool PetscBeganMPE = PETSC_FALSE;
2438 
2439 PETSC_INTERN PetscErrorCode PetscLogEventBeginMPE(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject);
2440 PETSC_INTERN PetscErrorCode PetscLogEventEndMPE(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject);
2441 
2442 /*@C
2443    PetscLogMPEBegin - Turns on MPE logging of events. This creates large log files
2444    and slows the program down.
2445 
2446    Collective over PETSC_COMM_WORLD
2447 
2448    Options Database Keys:
2449 . -log_mpe - Prints extensive log information
2450 
2451    Notes:
2452    A related routine is PetscLogDefaultBegin() (with the options key -log_view), which is
2453    intended for production runs since it logs only flop rates and object
2454    creation (and should not significantly slow the programs).
2455 
2456    Level: advanced
2457 
2458 .seealso: `PetscLogDump()`, `PetscLogDefaultBegin()`, `PetscLogAllBegin()`, `PetscLogEventActivate()`,
2459           `PetscLogEventDeactivate()`
2460 @*/
2461 PetscErrorCode PetscLogMPEBegin(void) {
2462   PetscFunctionBegin;
2463   /* Do MPE initialization */
2464   if (!MPE_Initialized_logging()) { /* This function exists in mpich 1.1.2 and higher */
2465     PetscCall(PetscInfo(0, "Initializing MPE.\n"));
2466     PetscCall(MPE_Init_log());
2467 
2468     PetscBeganMPE = PETSC_TRUE;
2469   } else {
2470     PetscCall(PetscInfo(0, "MPE already initialized. Not attempting to reinitialize.\n"));
2471   }
2472   PetscCall(PetscLogSet(PetscLogEventBeginMPE, PetscLogEventEndMPE));
2473   PetscFunctionReturn(0);
2474 }
2475 
2476 /*@C
2477    PetscLogMPEDump - Dumps the MPE logging info to file for later use with Jumpshot.
2478 
2479    Collective over PETSC_COMM_WORLD
2480 
2481    Level: advanced
2482 
2483 .seealso: `PetscLogDump()`, `PetscLogAllBegin()`, `PetscLogMPEBegin()`
2484 @*/
2485 PetscErrorCode PetscLogMPEDump(const char sname[]) {
2486   char name[PETSC_MAX_PATH_LEN];
2487 
2488   PetscFunctionBegin;
2489   if (PetscBeganMPE) {
2490     PetscCall(PetscInfo(0, "Finalizing MPE.\n"));
2491     if (sname) {
2492       PetscCall(PetscStrcpy(name, sname));
2493     } else {
2494       PetscCall(PetscGetProgramName(name, sizeof(name)));
2495     }
2496     PetscCall(MPE_Finish_log(name));
2497   } else {
2498     PetscCall(PetscInfo(0, "Not finalizing MPE (not started by PETSc).\n"));
2499   }
2500   PetscFunctionReturn(0);
2501 }
2502 
2503 #define PETSC_RGB_COLORS_MAX 39
2504 static const char *PetscLogMPERGBColors[PETSC_RGB_COLORS_MAX] = {"OliveDrab:      ", "BlueViolet:     ", "CadetBlue:      ", "CornflowerBlue: ", "DarkGoldenrod:  ", "DarkGreen:      ", "DarkKhaki:      ", "DarkOliveGreen: ",
2505                                                                  "DarkOrange:     ", "DarkOrchid:     ", "DarkSeaGreen:   ", "DarkSlateGray:  ", "DarkTurquoise:  ", "DeepPink:       ", "DarkKhaki:      ", "DimGray:        ",
2506                                                                  "DodgerBlue:     ", "GreenYellow:    ", "HotPink:        ", "IndianRed:      ", "LavenderBlush:  ", "LawnGreen:      ", "LemonChiffon:   ", "LightCoral:     ",
2507                                                                  "LightCyan:      ", "LightPink:      ", "LightSalmon:    ", "LightSlateGray: ", "LightYellow:    ", "LimeGreen:      ", "MediumPurple:   ", "MediumSeaGreen: ",
2508                                                                  "MediumSlateBlue:", "MidnightBlue:   ", "MintCream:      ", "MistyRose:      ", "NavajoWhite:    ", "NavyBlue:       ", "OliveDrab:      "};
2509 
2510 /*@C
2511   PetscLogMPEGetRGBColor - This routine returns a rgb color useable with PetscLogEventRegister()
2512 
2513   Not collective. Maybe it should be?
2514 
2515   Output Parameter:
2516 . str - character string representing the color
2517 
2518   Level: developer
2519 
2520 .seealso: `PetscLogEventRegister`
2521 @*/
2522 PetscErrorCode PetscLogMPEGetRGBColor(const char *str[]) {
2523   static int idx = 0;
2524 
2525   PetscFunctionBegin;
2526   *str = PetscLogMPERGBColors[idx];
2527   idx  = (idx + 1) % PETSC_RGB_COLORS_MAX;
2528   PetscFunctionReturn(0);
2529 }
2530 
2531 #endif /* PETSC_USE_LOG && PETSC_HAVE_MPE */
2532