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