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