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