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