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