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