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