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