xref: /petsc/src/sys/logging/plog.c (revision 2c9ec6dfe874b911fa49ef7e759d29a8430d6aff)
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   If an existing event with the same name exists, its event handle is
733   returned instead of creating a new event.
734 
735   Level: intermediate
736 
737 .keywords: log, event, register
738 .seealso: PetscLogEventBegin(), PetscLogEventEnd(), PetscLogFlops(),
739           PetscLogEventMPEActivate(), PetscLogEventMPEDeactivate(),
740           PetscLogEventActivate(), PetscLogEventDeactivate(), PetscClassIdRegister()
741 @*/
742 PetscErrorCode  PetscLogEventRegister(const char name[],PetscClassId classid,PetscLogEvent *event)
743 {
744   PetscStageLog  stageLog;
745   int            stage;
746   PetscErrorCode ierr;
747 
748   PetscFunctionBegin;
749   *event = PETSC_DECIDE;
750   ierr   = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
751   ierr   = EventRegLogGetEvent(stageLog->eventLog, name, event);CHKERRQ(ierr);
752   if (*event > 0) PetscFunctionReturn(0);
753   ierr   = EventRegLogRegister(stageLog->eventLog, name, classid, event);CHKERRQ(ierr);
754   for (stage = 0; stage < stageLog->numStages; stage++) {
755     ierr = EventPerfLogEnsureSize(stageLog->stageInfo[stage].eventLog, stageLog->eventLog->numEvents);CHKERRQ(ierr);
756     ierr = ClassPerfLogEnsureSize(stageLog->stageInfo[stage].classLog, stageLog->classLog->numClasses);CHKERRQ(ierr);
757   }
758   PetscFunctionReturn(0);
759 }
760 
761 #undef __FUNCT__
762 #define __FUNCT__ "PetscLogEventActivate"
763 /*@
764   PetscLogEventActivate - Indicates that a particular event should be logged.
765 
766   Not Collective
767 
768   Input Parameter:
769 . event - The event id
770 
771   Usage:
772 .vb
773       PetscLogEventDeactivate(VEC_SetValues);
774         [code where you do not want to log VecSetValues()]
775       PetscLogEventActivate(VEC_SetValues);
776         [code where you do want to log VecSetValues()]
777 .ve
778 
779   Note:
780   The event may be either a pre-defined PETSc event (found in include/petsclog.h)
781   or an event number obtained with PetscLogEventRegister().
782 
783   Level: advanced
784 
785 .keywords: log, event, activate
786 .seealso: PetscLogEventMPEDeactivate(),PetscLogEventMPEActivate(),PlogEventDeactivate()
787 @*/
788 PetscErrorCode  PetscLogEventActivate(PetscLogEvent event)
789 {
790   PetscStageLog  stageLog;
791   int            stage;
792   PetscErrorCode ierr;
793 
794   PetscFunctionBegin;
795   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
796   ierr = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
797   ierr = EventPerfLogActivate(stageLog->stageInfo[stage].eventLog, event);CHKERRQ(ierr);
798   PetscFunctionReturn(0);
799 }
800 
801 #undef __FUNCT__
802 #define __FUNCT__ "PetscLogEventDeactivate"
803 /*@
804   PetscLogEventDeactivate - Indicates that a particular event should not be logged.
805 
806   Not Collective
807 
808   Input Parameter:
809 . event - The event id
810 
811   Usage:
812 .vb
813       PetscLogEventDeactivate(VEC_SetValues);
814         [code where you do not want to log VecSetValues()]
815       PetscLogEventActivate(VEC_SetValues);
816         [code where you do want to log VecSetValues()]
817 .ve
818 
819   Note:
820   The event may be either a pre-defined PETSc event (found in
821   include/petsclog.h) or an event number obtained with PetscLogEventRegister()).
822 
823   Level: advanced
824 
825 .keywords: log, event, deactivate
826 .seealso: PetscLogEventMPEDeactivate(),PetscLogEventMPEActivate(),PlogEventActivate()
827 @*/
828 PetscErrorCode  PetscLogEventDeactivate(PetscLogEvent event)
829 {
830   PetscStageLog  stageLog;
831   int            stage;
832   PetscErrorCode ierr;
833 
834   PetscFunctionBegin;
835   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
836   ierr = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
837   ierr = EventPerfLogDeactivate(stageLog->stageInfo[stage].eventLog, event);CHKERRQ(ierr);
838   PetscFunctionReturn(0);
839 }
840 
841 #undef __FUNCT__
842 #define __FUNCT__ "PetscLogEventSetActiveAll"
843 /*@
844   PetscLogEventSetActiveAll - Sets the event activity in every stage.
845 
846   Not Collective
847 
848   Input Parameters:
849 + event    - The event id
850 - isActive - The activity flag determining whether the event is logged
851 
852   Level: advanced
853 
854 .keywords: log, event, activate
855 .seealso: PetscLogEventMPEDeactivate(),PetscLogEventMPEActivate(),PlogEventActivate(),PlogEventDeactivate()
856 @*/
857 PetscErrorCode  PetscLogEventSetActiveAll(PetscLogEvent event, PetscBool isActive)
858 {
859   PetscStageLog  stageLog;
860   int            stage;
861   PetscErrorCode ierr;
862 
863   PetscFunctionBegin;
864   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
865   for (stage = 0; stage < stageLog->numStages; stage++) {
866     if (isActive) {
867       ierr = EventPerfLogActivate(stageLog->stageInfo[stage].eventLog, event);CHKERRQ(ierr);
868     } else {
869       ierr = EventPerfLogDeactivate(stageLog->stageInfo[stage].eventLog, event);CHKERRQ(ierr);
870     }
871   }
872   PetscFunctionReturn(0);
873 }
874 
875 #undef __FUNCT__
876 #define __FUNCT__ "PetscLogEventActivateClass"
877 /*@
878   PetscLogEventActivateClass - Activates event logging for a PETSc object class.
879 
880   Not Collective
881 
882   Input Parameter:
883 . classid - The event class, for example MAT_CLASSID, SNES_CLASSID, etc.
884 
885   Level: developer
886 
887 .keywords: log, event, activate, class
888 .seealso: PetscInfoActivate(),PetscInfo(),PetscInfoAllow(),PetscLogEventDeactivateClass(), PetscLogEventActivate(),PetscLogEventDeactivate()
889 @*/
890 PetscErrorCode  PetscLogEventActivateClass(PetscClassId classid)
891 {
892   PetscStageLog  stageLog;
893   int            stage;
894   PetscErrorCode ierr;
895 
896   PetscFunctionBegin;
897   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
898   ierr = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
899   ierr = EventPerfLogActivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid);CHKERRQ(ierr);
900   PetscFunctionReturn(0);
901 }
902 
903 #undef __FUNCT__
904 #define __FUNCT__ "PetscLogEventDeactivateClass"
905 /*@
906   PetscLogEventDeactivateClass - Deactivates event logging for a PETSc object class.
907 
908   Not Collective
909 
910   Input Parameter:
911 . classid - The event class, for example MAT_CLASSID, SNES_CLASSID, etc.
912 
913   Level: developer
914 
915 .keywords: log, event, deactivate, class
916 .seealso: PetscInfoActivate(),PetscInfo(),PetscInfoAllow(),PetscLogEventActivateClass(), PetscLogEventActivate(),PetscLogEventDeactivate()
917 @*/
918 PetscErrorCode  PetscLogEventDeactivateClass(PetscClassId classid)
919 {
920   PetscStageLog  stageLog;
921   int            stage;
922   PetscErrorCode ierr;
923 
924   PetscFunctionBegin;
925   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
926   ierr = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
927   ierr = EventPerfLogDeactivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid);CHKERRQ(ierr);
928   PetscFunctionReturn(0);
929 }
930 
931 /*MC
932    PetscLogEventBegin - Logs the beginning of a user event.
933 
934    Synopsis:
935    #include <petsclog.h>
936    PetscErrorCode PetscLogEventBegin(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
937 
938    Not Collective
939 
940    Input Parameters:
941 +  e - integer associated with the event obtained from PetscLogEventRegister()
942 -  o1,o2,o3,o4 - objects associated with the event, or 0
943 
944 
945    Fortran Synopsis:
946    void PetscLogEventBegin(int e,PetscErrorCode ierr)
947 
948    Usage:
949 .vb
950      PetscLogEvent USER_EVENT;
951      PetscLogDouble user_event_flops;
952      PetscLogEventRegister("User event",0,&USER_EVENT);
953      PetscLogEventBegin(USER_EVENT,0,0,0,0);
954         [code segment to monitor]
955         PetscLogFlops(user_event_flops);
956      PetscLogEventEnd(USER_EVENT,0,0,0,0);
957 .ve
958 
959    Notes:
960    You need to register each integer event with the command
961    PetscLogEventRegister().  The source code must be compiled with
962    -DPETSC_USE_LOG, which is the default.
963 
964    PETSc automatically logs library events if the code has been
965    compiled with -DPETSC_USE_LOG, and -log, -log_summary, or -log_all are
966    specified.  PetscLogEventBegin() is intended for logging user events
967    to supplement this PETSc information.
968 
969    Level: intermediate
970 
971 .seealso: PetscLogEventRegister(), PetscLogEventEnd(), PetscLogFlops()
972 
973 .keywords: log, event, begin
974 M*/
975 
976 /*MC
977    PetscLogEventEnd - Log the end of a user event.
978 
979    Synopsis:
980    #include <petsclog.h>
981    PetscErrorCode PetscLogEventEnd(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
982 
983    Not Collective
984 
985    Input Parameters:
986 +  e - integer associated with the event obtained with PetscLogEventRegister()
987 -  o1,o2,o3,o4 - objects associated with the event, or 0
988 
989 
990    Fortran Synopsis:
991    void PetscLogEventEnd(int e,PetscErrorCode ierr)
992 
993    Usage:
994 .vb
995      PetscLogEvent USER_EVENT;
996      PetscLogDouble user_event_flops;
997      PetscLogEventRegister("User event",0,&USER_EVENT,);
998      PetscLogEventBegin(USER_EVENT,0,0,0,0);
999         [code segment to monitor]
1000         PetscLogFlops(user_event_flops);
1001      PetscLogEventEnd(USER_EVENT,0,0,0,0);
1002 .ve
1003 
1004    Notes:
1005    You should also register each additional integer event with the command
1006    PetscLogEventRegister(). Source code must be compiled with
1007    -DPETSC_USE_LOG, which is the default.
1008 
1009    PETSc automatically logs library events if the code has been
1010    compiled with -DPETSC_USE_LOG, and -log, -log_summary, or -log_all are
1011    specified.  PetscLogEventEnd() is intended for logging user events
1012    to supplement this PETSc information.
1013 
1014    Level: intermediate
1015 
1016 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogFlops()
1017 
1018 .keywords: log, event, end
1019 M*/
1020 
1021 /*MC
1022    PetscLogEventBarrierBegin - Logs the time in a barrier before an event.
1023 
1024    Synopsis:
1025    #include <petsclog.h>
1026    PetscErrorCode PetscLogEventBarrierBegin(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4,MPI_Comm comm)
1027 
1028    Not Collective
1029 
1030    Input Parameters:
1031 .  e - integer associated with the event obtained from PetscLogEventRegister()
1032 .  o1,o2,o3,o4 - objects associated with the event, or 0
1033 .  comm - communicator the barrier takes place over
1034 
1035 
1036    Usage:
1037 .vb
1038      PetscLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
1039        MPI_Allreduce()
1040      PetscLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm);
1041 .ve
1042 
1043    Notes:
1044    This is for logging the amount of time spent in a barrier for an event
1045    that requires synchronization.
1046 
1047    Additional Notes:
1048    Synchronization events always come in pairs; for example, VEC_NormBarrier and
1049    VEC_NormComm = VEC_NormBarrier + 1
1050 
1051    Level: advanced
1052 
1053 .seealso: PetscLogEventRegister(), PetscLogEventEnd(), PetscLogFlops(), PetscLogEventBegin(),
1054           PetscLogEventBarrierEnd()
1055 
1056 .keywords: log, event, begin, barrier
1057 M*/
1058 
1059 /*MC
1060    PetscLogEventBarrierEnd - Logs the time in a barrier before an event.
1061 
1062    Synopsis:
1063    #include <petsclog.h>
1064    PetscErrorCode PetscLogEventBarrierEnd(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4,MPI_Comm comm)
1065 
1066    Logically Collective on MPI_Comm
1067 
1068    Input Parameters:
1069 .  e - integer associated with the event obtained from PetscLogEventRegister()
1070 .  o1,o2,o3,o4 - objects associated with the event, or 0
1071 .  comm - communicator the barrier takes place over
1072 
1073 
1074     Usage:
1075 .vb
1076      PetscLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
1077        MPI_Allreduce()
1078      PetscLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm);
1079 .ve
1080 
1081    Notes:
1082    This is for logging the amount of time spent in a barrier for an event
1083    that requires synchronization.
1084 
1085    Additional Notes:
1086    Synchronization events always come in pairs; for example, VEC_NormBarrier and
1087    VEC_NormComm = VEC_NormBarrier + 1
1088 
1089    Level: advanced
1090 
1091 .seealso: PetscLogEventRegister(), PetscLogEventEnd(), PetscLogFlops(), PetscLogEventBegin(),
1092           PetscLogEventBarrierBegin()
1093 
1094 .keywords: log, event, begin, barrier
1095 M*/
1096 
1097 #undef __FUNCT__
1098 #define __FUNCT__ "PetscLogEventGetId"
1099 /*@C
1100   PetscLogEventGetId - Returns the event id when given the event name.
1101 
1102   Not Collective
1103 
1104   Input Parameter:
1105 . name  - The event name
1106 
1107   Output Parameter:
1108 . event - The event, or -1 if no event with that name exists
1109 
1110   Level: intermediate
1111 
1112 .seealso: PetscLogEventBegin(), PetscLogEventEnd(), PetscLogStageGetId()
1113 @*/
1114 PetscErrorCode  PetscLogEventGetId(const char name[], PetscLogEvent *event)
1115 {
1116   PetscStageLog  stageLog;
1117   PetscErrorCode ierr;
1118 
1119   PetscFunctionBegin;
1120   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1121   ierr = EventRegLogGetEvent(stageLog->eventLog, name, event);CHKERRQ(ierr);
1122   PetscFunctionReturn(0);
1123 }
1124 
1125 
1126 /*------------------------------------------------ Output Functions -------------------------------------------------*/
1127 #undef __FUNCT__
1128 #define __FUNCT__ "PetscLogDump"
1129 /*@C
1130   PetscLogDump - Dumps logs of objects to a file. This file is intended to
1131   be read by bin/petscview. This program no longer exists.
1132 
1133   Collective on PETSC_COMM_WORLD
1134 
1135   Input Parameter:
1136 . name - an optional file name
1137 
1138   Options Database Keys:
1139 + -log     - Prints basic log information (for code compiled with PETSC_USE_LOG)
1140 - -log_all - Prints extensive log information (for code compiled with PETSC_USE_LOG)
1141 
1142   Usage:
1143 .vb
1144      PetscInitialize(...);
1145      PetscLogBegin(); or PetscLogAllBegin();
1146      ... code ...
1147      PetscLogDump(filename);
1148      PetscFinalize();
1149 .ve
1150 
1151   Notes:
1152   The default file name is
1153 $    Log.<rank>
1154   where <rank> is the processor number. If no name is specified,
1155   this file will be used.
1156 
1157   Level: advanced
1158 
1159 .keywords: log, dump
1160 .seealso: PetscLogBegin(), PetscLogAllBegin(), PetscLogView()
1161 @*/
1162 PetscErrorCode  PetscLogDump(const char sname[])
1163 {
1164   PetscStageLog      stageLog;
1165   PetscEventPerfInfo *eventInfo;
1166   FILE               *fd;
1167   char               file[PETSC_MAX_PATH_LEN], fname[PETSC_MAX_PATH_LEN];
1168   PetscLogDouble     flops, _TotalTime;
1169   PetscMPIInt        rank;
1170   int                action, object, curStage;
1171   PetscLogEvent      event;
1172   PetscErrorCode     ierr;
1173 
1174   PetscFunctionBegin;
1175   /* Calculate the total elapsed time */
1176   PetscTime(&_TotalTime);
1177   _TotalTime -= petsc_BaseTime;
1178   /* Open log file */
1179   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
1180   if (sname) sprintf(file, "%s.%d", sname, rank);
1181   else sprintf(file, "Log.%d", rank);
1182   ierr = PetscFixFilename(file, fname);CHKERRQ(ierr);
1183   ierr = PetscFOpen(PETSC_COMM_WORLD, fname, "w", &fd);CHKERRQ(ierr);
1184   if ((!rank) && (!fd)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN, "Cannot open file: %s", fname);
1185   /* Output totals */
1186   ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Total Flops %14e %16.8e\n", petsc_TotalFlops, _TotalTime);CHKERRQ(ierr);
1187   ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Clock Resolution %g\n", 0.0);CHKERRQ(ierr);
1188   /* Output actions */
1189   if (petsc_logActions) {
1190     ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Actions accomplished %d\n", petsc_numActions);CHKERRQ(ierr);
1191     for (action = 0; action < petsc_numActions; action++) {
1192       ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "%g %d %d %d %d %d %d %g %g %g\n",
1193                           petsc_actions[action].time, petsc_actions[action].action, (int)petsc_actions[action].event, (int)petsc_actions[action].classid, petsc_actions[action].id1,
1194                           petsc_actions[action].id2, petsc_actions[action].id3, petsc_actions[action].flops, petsc_actions[action].mem, petsc_actions[action].maxmem);CHKERRQ(ierr);
1195     }
1196   }
1197   /* Output objects */
1198   if (petsc_logObjects) {
1199     ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Objects created %d destroyed %d\n", petsc_numObjects, petsc_numObjectsDestroyed);CHKERRQ(ierr);
1200     for (object = 0; object < petsc_numObjects; object++) {
1201       ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Parent ID: %d Memory: %d\n", petsc_objects[object].parent, (int) petsc_objects[object].mem);CHKERRQ(ierr);
1202       if (!petsc_objects[object].name[0]) {
1203         ierr = PetscFPrintf(PETSC_COMM_WORLD, fd,"No Name\n");CHKERRQ(ierr);
1204       } else {
1205         ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Name: %s\n", petsc_objects[object].name);CHKERRQ(ierr);
1206       }
1207       if (petsc_objects[object].info[0] != 0) {
1208         ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "No Info\n");CHKERRQ(ierr);
1209       } else {
1210         ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Info: %s\n", petsc_objects[object].info);CHKERRQ(ierr);
1211       }
1212     }
1213   }
1214   /* Output events */
1215   ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Event log:\n");CHKERRQ(ierr);
1216   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1217   ierr = PetscIntStackTop(stageLog->stack, &curStage);CHKERRQ(ierr);
1218   eventInfo = stageLog->stageInfo[curStage].eventLog->eventInfo;
1219   for (event = 0; event < stageLog->stageInfo[curStage].eventLog->numEvents; event++) {
1220     if (eventInfo[event].time != 0.0) flops = eventInfo[event].flops/eventInfo[event].time;
1221     else flops = 0.0;
1222     ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "%d %16d %16g %16g %16g\n", event, eventInfo[event].count,
1223                         eventInfo[event].flops, eventInfo[event].time, flops);CHKERRQ(ierr);
1224   }
1225   ierr = PetscFClose(PETSC_COMM_WORLD, fd);CHKERRQ(ierr);
1226   PetscFunctionReturn(0);
1227 }
1228 
1229 #undef __FUNCT__
1230 #define __FUNCT__ "PetscLogView_Detailed"
1231 /*
1232   PetscLogView_Detailed - Each process prints the times for its own events
1233 
1234 */
1235 PetscErrorCode  PetscLogView_Detailed(PetscViewer viewer)
1236 {
1237   MPI_Comm           comm       = PetscObjectComm((PetscObject) viewer);
1238   PetscEventPerfInfo *eventInfo = NULL;
1239   PetscLogDouble     locTotalTime, numRed, maxMem;
1240   PetscStageLog      stageLog;
1241   int                numStages,numEvents,stage,event;
1242   PetscMPIInt        rank,size;
1243   PetscErrorCode     ierr;
1244 
1245   PetscFunctionBegin;
1246   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1247   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1248   /* Must preserve reduction count before we go on */
1249   numRed = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
1250   /* Get the total elapsed time */
1251   PetscTime(&locTotalTime);  locTotalTime -= petsc_BaseTime;
1252   ierr = PetscViewerASCIIPrintf(viewer,"numProcs   = %d\n",size);CHKERRQ(ierr);
1253   ierr = PetscViewerASCIIPrintf(viewer,"LocalTimes = {}\n");CHKERRQ(ierr);
1254   ierr = PetscViewerASCIIPrintf(viewer,"LocalFlops = {}\n");CHKERRQ(ierr);
1255   ierr = PetscViewerASCIIPrintf(viewer,"LocalMessageLens = {}\n");CHKERRQ(ierr);
1256   ierr = PetscViewerASCIIPrintf(viewer,"LocalMessages = {}\n");CHKERRQ(ierr);
1257   ierr = PetscViewerASCIIPrintf(viewer,"LocalReductions = {}\n");CHKERRQ(ierr);
1258   ierr = PetscViewerASCIIPrintf(viewer,"LocalObjects = {}\n");CHKERRQ(ierr);
1259   ierr = PetscViewerASCIIPrintf(viewer,"LocalMemory = {}\n");CHKERRQ(ierr);
1260   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1261   ierr = MPI_Allreduce(&stageLog->numStages, &numStages, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1262   ierr = PetscViewerASCIIPrintf(viewer,"Stages = {}\n");CHKERRQ(ierr);
1263   for (stage=0; stage<numStages; stage++) {
1264     ierr = PetscViewerASCIIPrintf(viewer,"Stages[\"%s\"] = {}\n",stageLog->stageInfo[stage].name);CHKERRQ(ierr);
1265     ierr = MPI_Allreduce(&stageLog->stageInfo[stage].eventLog->numEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1266     for (event = 0; event < numEvents; event++) {
1267       ierr = PetscViewerASCIIPrintf(viewer,"Stages[\"%s\"][\"%s\"] = {}\n",stageLog->stageInfo[stage].name,stageLog->eventLog->eventInfo[event].name);CHKERRQ(ierr);
1268     }
1269   }
1270   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
1271   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalTimes[%d] = %g\n",rank,locTotalTime);CHKERRQ(ierr);
1272   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalFlops[%d] = %g\n",rank,petsc_TotalFlops);CHKERRQ(ierr);
1273   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalMessageLens[%d] = %g\n",rank,(petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len));CHKERRQ(ierr);
1274   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalMessages[%d] = %g\n",rank,(petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct));CHKERRQ(ierr);
1275   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalReductions[%d] = %g\n",rank,numRed);CHKERRQ(ierr);
1276   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalObjects[%d] = %g\n",rank,petsc_numObjects);CHKERRQ(ierr);
1277   ierr = PetscMallocGetMaximumUsage(&maxMem);CHKERRQ(ierr);
1278   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalMemory[%d] = %g\n",rank,maxMem);CHKERRQ(ierr);
1279   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1280   for (stage=0; stage<numStages; stage++) {
1281     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Stages[\"%s\"][\"summary\"][%d] = {\"time\" : %g, \"numMessages\" : %g, \"messageLength\" : %g, \"numReductions\" : %g, \"flops\" : %g}\n",
1282                                               stageLog->stageInfo[stage].name,rank,
1283                                               stageLog->stageInfo[stage].perfInfo.time,stageLog->stageInfo[stage].perfInfo.numMessages,stageLog->stageInfo[stage].perfInfo.messageLength,
1284                                               stageLog->stageInfo[stage].perfInfo.numReductions,stageLog->stageInfo[stage].perfInfo.flops);CHKERRQ(ierr);
1285     ierr = MPI_Allreduce(&stageLog->stageInfo[stage].eventLog->numEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1286     for (event = 0; event < numEvents; event++) {
1287       eventInfo = stageLog->stageInfo[stage].eventLog->eventInfo;
1288       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,
1289                                                 eventInfo[event].count, eventInfo[event].time,eventInfo[event].numMessages, eventInfo[event].messageLength,
1290                                                 eventInfo[event].numReductions,eventInfo[event].flops);CHKERRQ(ierr);
1291     }
1292   }
1293   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1294   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
1295   PetscFunctionReturn(0);
1296 }
1297 
1298 #undef __FUNCT__
1299 #define __FUNCT__ "PetscLogView_Default"
1300 PetscErrorCode  PetscLogView_Default(PetscViewer viewer)
1301 {
1302   FILE               *fd;
1303   PetscLogDouble     zero       = 0.0;
1304   PetscStageLog      stageLog;
1305   PetscStageInfo     *stageInfo = NULL;
1306   PetscEventPerfInfo *eventInfo = NULL;
1307   PetscClassPerfInfo *classInfo;
1308   char               arch[128],hostname[128],username[128],pname[PETSC_MAX_PATH_LEN],date[128];
1309   const char         *name;
1310   PetscLogDouble     locTotalTime, TotalTime, TotalFlops;
1311   PetscLogDouble     numMessages, messageLength, avgMessLen, numReductions;
1312   PetscLogDouble     stageTime, flops, flopr, mem, mess, messLen, red;
1313   PetscLogDouble     fracTime, fracFlops, fracMessages, fracLength, fracReductions, fracMess, fracMessLen, fracRed;
1314   PetscLogDouble     fracStageTime, fracStageFlops, fracStageMess, fracStageMessLen, fracStageRed;
1315   PetscLogDouble     min, max, tot, ratio, avg, x, y;
1316   PetscLogDouble     minf, maxf, totf, ratf, mint, maxt, tott, ratt, ratCt, totm, totml, totr;
1317   PetscMPIInt        minCt, maxCt;
1318   PetscMPIInt        size, rank;
1319   PetscBool          *localStageUsed,    *stageUsed;
1320   PetscBool          *localStageVisible, *stageVisible;
1321   int                numStages, localNumEvents, numEvents;
1322   int                stage, oclass;
1323   PetscLogEvent      event;
1324   PetscErrorCode     ierr;
1325   char               version[256];
1326   MPI_Comm           comm;
1327   PetscInt           nthreads;
1328 
1329   PetscFunctionBegin;
1330   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
1331   ierr = PetscViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr);
1332   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1333   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1334   /* Get the total elapsed time */
1335   PetscTime(&locTotalTime);  locTotalTime -= petsc_BaseTime;
1336 
1337   ierr = PetscFPrintf(comm, fd, "************************************************************************************************************************\n");CHKERRQ(ierr);
1338   ierr = PetscFPrintf(comm, fd, "***             WIDEN YOUR WINDOW TO 120 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document            ***\n");CHKERRQ(ierr);
1339   ierr = PetscFPrintf(comm, fd, "************************************************************************************************************************\n");CHKERRQ(ierr);
1340   ierr = PetscFPrintf(comm, fd, "\n---------------------------------------------- PETSc Performance Summary: ----------------------------------------------\n\n");CHKERRQ(ierr);
1341   ierr = PetscGetArchType(arch,sizeof(arch));CHKERRQ(ierr);
1342   ierr = PetscGetHostName(hostname,sizeof(hostname));CHKERRQ(ierr);
1343   ierr = PetscGetUserName(username,sizeof(username));CHKERRQ(ierr);
1344   ierr = PetscGetProgramName(pname,sizeof(pname));CHKERRQ(ierr);
1345   ierr = PetscGetDate(date,sizeof(date));CHKERRQ(ierr);
1346   ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
1347   if (size == 1) {
1348     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);
1349   } else {
1350     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);
1351   }
1352   ierr = PetscThreadCommGetNThreads(PETSC_COMM_WORLD,&nthreads);CHKERRQ(ierr);
1353   if (nthreads > 1) {
1354     ierr = PetscFPrintf(comm,fd,"With %d threads per MPI_Comm\n", (int)nthreads);CHKERRQ(ierr);
1355   }
1356 
1357   ierr = PetscFPrintf(comm, fd, "Using %s\n", version);CHKERRQ(ierr);
1358 
1359   /* Must preserve reduction count before we go on */
1360   red = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
1361 
1362   /* Calculate summary information */
1363   ierr = PetscFPrintf(comm, fd, "\n                         Max       Max/Min        Avg      Total \n");CHKERRQ(ierr);
1364   /*   Time */
1365   ierr = MPI_Allreduce(&locTotalTime, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1366   ierr = MPI_Allreduce(&locTotalTime, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1367   ierr = MPI_Allreduce(&locTotalTime, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1368   avg  = (tot)/((PetscLogDouble) size);
1369   if (min != 0.0) ratio = max/min;
1370   else ratio = 0.0;
1371   ierr = PetscFPrintf(comm, fd, "Time (sec):           %5.3e   %10.5f   %5.3e\n", max, ratio, avg);CHKERRQ(ierr);
1372   TotalTime = tot;
1373   /*   Objects */
1374   avg  = (PetscLogDouble) petsc_numObjects;
1375   ierr = MPI_Allreduce(&avg,          &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1376   ierr = MPI_Allreduce(&avg,          &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1377   ierr = MPI_Allreduce(&avg,          &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1378   avg  = (tot)/((PetscLogDouble) size);
1379   if (min != 0.0) ratio = max/min;
1380   else ratio = 0.0;
1381   ierr = PetscFPrintf(comm, fd, "Objects:              %5.3e   %10.5f   %5.3e\n", max, ratio, avg);CHKERRQ(ierr);
1382   /*   Flops */
1383   ierr = MPI_Allreduce(&petsc_TotalFlops,  &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1384   ierr = MPI_Allreduce(&petsc_TotalFlops,  &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1385   ierr = MPI_Allreduce(&petsc_TotalFlops,  &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1386   avg  = (tot)/((PetscLogDouble) size);
1387   if (min != 0.0) ratio = max/min;
1388   else ratio = 0.0;
1389   ierr = PetscFPrintf(comm, fd, "Flops:                %5.3e   %10.5f   %5.3e  %5.3e\n", max, ratio, avg, tot);CHKERRQ(ierr);
1390   TotalFlops = tot;
1391   /*   Flops/sec -- Must talk to Barry here */
1392   if (locTotalTime != 0.0) flops = petsc_TotalFlops/locTotalTime;
1393   else flops = 0.0;
1394   ierr = MPI_Allreduce(&flops,        &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1395   ierr = MPI_Allreduce(&flops,        &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1396   ierr = MPI_Allreduce(&flops,        &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1397   avg  = (tot)/((PetscLogDouble) size);
1398   if (min != 0.0) ratio = max/min;
1399   else ratio = 0.0;
1400   ierr = PetscFPrintf(comm, fd, "Flops/sec:            %5.3e   %10.5f   %5.3e  %5.3e\n", max, ratio, avg, tot);CHKERRQ(ierr);
1401   /*   Memory */
1402   ierr = PetscMallocGetMaximumUsage(&mem);CHKERRQ(ierr);
1403   if (mem > 0.0) {
1404     ierr = MPI_Allreduce(&mem,          &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1405     ierr = MPI_Allreduce(&mem,          &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1406     ierr = MPI_Allreduce(&mem,          &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1407     avg  = (tot)/((PetscLogDouble) size);
1408     if (min != 0.0) ratio = max/min;
1409     else ratio = 0.0;
1410     ierr = PetscFPrintf(comm, fd, "Memory:               %5.3e   %10.5f              %5.3e\n", max, ratio, tot);CHKERRQ(ierr);
1411   }
1412   /*   Messages */
1413   mess = 0.5*(petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct);
1414   ierr = MPI_Allreduce(&mess,         &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1415   ierr = MPI_Allreduce(&mess,         &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1416   ierr = MPI_Allreduce(&mess,         &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1417   avg  = (tot)/((PetscLogDouble) size);
1418   if (min != 0.0) ratio = max/min;
1419   else ratio = 0.0;
1420   ierr = PetscFPrintf(comm, fd, "MPI Messages:         %5.3e   %10.5f   %5.3e  %5.3e\n", max, ratio, avg, tot);CHKERRQ(ierr);
1421   numMessages = tot;
1422   /*   Message Lengths */
1423   mess = 0.5*(petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len);
1424   ierr = MPI_Allreduce(&mess,         &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1425   ierr = MPI_Allreduce(&mess,         &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1426   ierr = MPI_Allreduce(&mess,         &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1427   if (numMessages != 0) avg = (tot)/(numMessages);
1428   else avg = 0.0;
1429   if (min != 0.0) ratio = max/min;
1430   else ratio = 0.0;
1431   ierr = PetscFPrintf(comm, fd, "MPI Message Lengths:  %5.3e   %10.5f   %5.3e  %5.3e\n", max, ratio, avg, tot);CHKERRQ(ierr);
1432   messageLength = tot;
1433   /*   Reductions */
1434   ierr = MPI_Allreduce(&red,          &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1435   ierr = MPI_Allreduce(&red,          &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1436   ierr = MPI_Allreduce(&red,          &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1437   if (min != 0.0) ratio = max/min;
1438   else ratio = 0.0;
1439   ierr = PetscFPrintf(comm, fd, "MPI Reductions:       %5.3e   %10.5f\n", max, ratio);CHKERRQ(ierr);
1440   numReductions = red; /* wrong because uses count from process zero */
1441   ierr = PetscFPrintf(comm, fd, "\nFlop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)\n");CHKERRQ(ierr);
1442   ierr = PetscFPrintf(comm, fd, "                            e.g., VecAXPY() for real vectors of length N --> 2N flops\n");CHKERRQ(ierr);
1443   ierr = PetscFPrintf(comm, fd, "                            and VecAXPY() for complex vectors of length N --> 8N flops\n");CHKERRQ(ierr);
1444 
1445   /* Get total number of stages --
1446        Currently, a single processor can register more stages than another, but stages must all be registered in order.
1447        We can removed this requirement if necessary by having a global stage numbering and indirection on the stage ID.
1448        This seems best accomplished by assoicating a communicator with each stage.
1449   */
1450   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1451   ierr = MPI_Allreduce(&stageLog->numStages, &numStages, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1452   ierr = PetscMalloc1(numStages, &localStageUsed);CHKERRQ(ierr);
1453   ierr = PetscMalloc1(numStages, &stageUsed);CHKERRQ(ierr);
1454   ierr = PetscMalloc1(numStages, &localStageVisible);CHKERRQ(ierr);
1455   ierr = PetscMalloc1(numStages, &stageVisible);CHKERRQ(ierr);
1456   if (numStages > 0) {
1457     stageInfo = stageLog->stageInfo;
1458     for (stage = 0; stage < numStages; stage++) {
1459       if (stage < stageLog->numStages) {
1460         localStageUsed[stage]    = stageInfo[stage].used;
1461         localStageVisible[stage] = stageInfo[stage].perfInfo.visible;
1462       } else {
1463         localStageUsed[stage]    = PETSC_FALSE;
1464         localStageVisible[stage] = PETSC_TRUE;
1465       }
1466     }
1467     ierr = MPI_Allreduce(localStageUsed,    stageUsed,    numStages, MPIU_BOOL, MPI_LOR,  comm);CHKERRQ(ierr);
1468     ierr = MPI_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr);
1469     for (stage = 0; stage < numStages; stage++) {
1470       if (stageUsed[stage]) {
1471         ierr = PetscFPrintf(comm, fd, "\nSummary of Stages:   ----- Time ------  ----- Flops -----  --- Messages ---  -- Message Lengths --  -- Reductions --\n");CHKERRQ(ierr);
1472         ierr = PetscFPrintf(comm, fd, "                        Avg     %%Total     Avg     %%Total   counts   %%Total     Avg         %%Total   counts   %%Total \n");CHKERRQ(ierr);
1473         break;
1474       }
1475     }
1476     for (stage = 0; stage < numStages; stage++) {
1477       if (!stageUsed[stage]) continue;
1478       if (localStageUsed[stage]) {
1479         ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.time,          &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1480         ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.flops,         &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1481         ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages,   &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1482         ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1483         ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1484         name = stageInfo[stage].name;
1485       } else {
1486         ierr = MPI_Allreduce(&zero,                           &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1487         ierr = MPI_Allreduce(&zero,                           &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1488         ierr = MPI_Allreduce(&zero,                           &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1489         ierr = MPI_Allreduce(&zero,                           &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1490         ierr = MPI_Allreduce(&zero,                           &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1491         name = "";
1492       }
1493       mess *= 0.5; messLen *= 0.5; red /= size;
1494       if (TotalTime     != 0.0) fracTime       = stageTime/TotalTime;    else fracTime       = 0.0;
1495       if (TotalFlops    != 0.0) fracFlops      = flops/TotalFlops;       else fracFlops      = 0.0;
1496       /* Talk to Barry if (stageTime     != 0.0) flops          = (size*flops)/stageTime; else flops          = 0.0; */
1497       if (numMessages   != 0.0) fracMessages   = mess/numMessages;       else fracMessages   = 0.0;
1498       if (numMessages   != 0.0) avgMessLen     = messLen/numMessages;    else avgMessLen     = 0.0;
1499       if (messageLength != 0.0) fracLength     = messLen/messageLength;  else fracLength     = 0.0;
1500       if (numReductions != 0.0) fracReductions = red/numReductions;      else fracReductions = 0.0;
1501       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",
1502                           stage, name, stageTime/size, 100.0*fracTime, flops, 100.0*fracFlops,
1503                           mess, 100.0*fracMessages, avgMessLen, 100.0*fracLength, red, 100.0*fracReductions);CHKERRQ(ierr);
1504     }
1505   }
1506 
1507   ierr = PetscFPrintf(comm, fd,"\n------------------------------------------------------------------------------------------------------------------------\n");CHKERRQ(ierr);
1508   ierr = PetscFPrintf(comm, fd, "See the 'Profiling' chapter of the users' manual for details on interpreting output.\n");CHKERRQ(ierr);
1509   ierr = PetscFPrintf(comm, fd, "Phase summary info:\n");CHKERRQ(ierr);
1510   ierr = PetscFPrintf(comm, fd, "   Count: number of times phase was executed\n");CHKERRQ(ierr);
1511   ierr = PetscFPrintf(comm, fd, "   Time and Flops: Max - maximum over all processors\n");CHKERRQ(ierr);
1512   ierr = PetscFPrintf(comm, fd, "                   Ratio - ratio of maximum to minimum over all processors\n");CHKERRQ(ierr);
1513   ierr = PetscFPrintf(comm, fd, "   Mess: number of messages sent\n");CHKERRQ(ierr);
1514   ierr = PetscFPrintf(comm, fd, "   Avg. len: average message length (bytes)\n");CHKERRQ(ierr);
1515   ierr = PetscFPrintf(comm, fd, "   Reduct: number of global reductions\n");CHKERRQ(ierr);
1516   ierr = PetscFPrintf(comm, fd, "   Global: entire computation\n");CHKERRQ(ierr);
1517   ierr = PetscFPrintf(comm, fd, "   Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().\n");CHKERRQ(ierr);
1518   ierr = PetscFPrintf(comm, fd, "      %%T - percent time in this phase         %%F - percent flops in this phase\n");CHKERRQ(ierr);
1519   ierr = PetscFPrintf(comm, fd, "      %%M - percent messages in this phase     %%L - percent message lengths in this phase\n");CHKERRQ(ierr);
1520   ierr = PetscFPrintf(comm, fd, "      %%R - percent reductions in this phase\n");CHKERRQ(ierr);
1521   ierr = PetscFPrintf(comm, fd, "   Total Mflop/s: 10e-6 * (sum of flops over all processors)/(max time over all processors)\n");CHKERRQ(ierr);
1522   ierr = PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------\n");CHKERRQ(ierr);
1523 
1524 #if defined(PETSC_USE_DEBUG)
1525   ierr = PetscFPrintf(comm, fd, "\n\n");CHKERRQ(ierr);
1526   ierr = PetscFPrintf(comm, fd, "      ##########################################################\n");CHKERRQ(ierr);
1527   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1528   ierr = PetscFPrintf(comm, fd, "      #                          WARNING!!!                    #\n");CHKERRQ(ierr);
1529   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1530   ierr = PetscFPrintf(comm, fd, "      #   This code was compiled with a debugging option,      #\n");CHKERRQ(ierr);
1531   ierr = PetscFPrintf(comm, fd, "      #   To get timing results run ./configure                #\n");CHKERRQ(ierr);
1532   ierr = PetscFPrintf(comm, fd, "      #   using --with-debugging=no, the performance will      #\n");CHKERRQ(ierr);
1533   ierr = PetscFPrintf(comm, fd, "      #   be generally two or three times faster.              #\n");CHKERRQ(ierr);
1534   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1535   ierr = PetscFPrintf(comm, fd, "      ##########################################################\n\n\n");CHKERRQ(ierr);
1536 #endif
1537 #if defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_FORTRAN_KERNELS)
1538   ierr = PetscFPrintf(comm, fd, "\n\n");CHKERRQ(ierr);
1539   ierr = PetscFPrintf(comm, fd, "      ##########################################################\n");CHKERRQ(ierr);
1540   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1541   ierr = PetscFPrintf(comm, fd, "      #                          WARNING!!!                    #\n");CHKERRQ(ierr);
1542   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1543   ierr = PetscFPrintf(comm, fd, "      #   The code for various complex numbers numerical       #\n");CHKERRQ(ierr);
1544   ierr = PetscFPrintf(comm, fd, "      #   kernels uses C++, which generally is not well        #\n");CHKERRQ(ierr);
1545   ierr = PetscFPrintf(comm, fd, "      #   optimized.  For performance that is about 4-5 times  #\n");CHKERRQ(ierr);
1546   ierr = PetscFPrintf(comm, fd, "      #   faster, specify --with-fortran-kernels=1             #\n");CHKERRQ(ierr);
1547   ierr = PetscFPrintf(comm, fd, "      #   when running ./configure.py.                         #\n");CHKERRQ(ierr);
1548   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1549   ierr = PetscFPrintf(comm, fd, "      ##########################################################\n\n\n");CHKERRQ(ierr);
1550 #endif
1551 
1552   /* Report events */
1553   ierr = PetscFPrintf(comm, fd,"Event                Count      Time (sec)     Flops                             --- Global ---  --- Stage ---   Total\n");CHKERRQ(ierr);
1554   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);
1555   ierr = PetscFPrintf(comm,fd,"------------------------------------------------------------------------------------------------------------------------\n");CHKERRQ(ierr);
1556 
1557   /* Problem: The stage name will not show up unless the stage executed on proc 1 */
1558   for (stage = 0; stage < numStages; stage++) {
1559     if (!stageVisible[stage]) continue;
1560     if (localStageUsed[stage]) {
1561       ierr = PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stageInfo[stage].name);CHKERRQ(ierr);
1562       ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.time,          &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1563       ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.flops,         &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1564       ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages,   &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1565       ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1566       ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1567     } else {
1568       ierr = PetscFPrintf(comm, fd, "\n--- Event Stage %d: Unknown\n\n", stage);CHKERRQ(ierr);
1569       ierr = MPI_Allreduce(&zero,                           &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1570       ierr = MPI_Allreduce(&zero,                           &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1571       ierr = MPI_Allreduce(&zero,                           &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1572       ierr = MPI_Allreduce(&zero,                           &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1573       ierr = MPI_Allreduce(&zero,                           &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1574     }
1575     mess *= 0.5; messLen *= 0.5; red /= size;
1576 
1577     /* Get total number of events in this stage --
1578        Currently, a single processor can register more events than another, but events must all be registered in order,
1579        just like stages. We can removed this requirement if necessary by having a global event numbering and indirection
1580        on the event ID. This seems best accomplished by assoicating a communicator with each stage.
1581 
1582        Problem: If the event did not happen on proc 1, its name will not be available.
1583        Problem: Event visibility is not implemented
1584     */
1585     if (localStageUsed[stage]) {
1586       eventInfo      = stageLog->stageInfo[stage].eventLog->eventInfo;
1587       localNumEvents = stageLog->stageInfo[stage].eventLog->numEvents;
1588     } else localNumEvents = 0;
1589     ierr = MPI_Allreduce(&localNumEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1590     for (event = 0; event < numEvents; event++) {
1591       if (localStageUsed[stage] && (event < stageLog->stageInfo[stage].eventLog->numEvents) && (eventInfo[event].depth == 0)) {
1592         if ((eventInfo[event].count > 0) && (eventInfo[event].time > 0.0)) flopr = eventInfo[event].flops;
1593         else flopr = 0.0;
1594 
1595         ierr = MPI_Allreduce(&flopr,                          &minf,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1596         ierr = MPI_Allreduce(&flopr,                          &maxf,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1597         ierr = MPI_Allreduce(&eventInfo[event].flops,         &totf,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1598         ierr = MPI_Allreduce(&eventInfo[event].time,          &mint,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1599         ierr = MPI_Allreduce(&eventInfo[event].time,          &maxt,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1600         ierr = MPI_Allreduce(&eventInfo[event].time,          &tott,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1601         ierr = MPI_Allreduce(&eventInfo[event].numMessages,   &totm,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1602         ierr = MPI_Allreduce(&eventInfo[event].messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1603         ierr = MPI_Allreduce(&eventInfo[event].numReductions, &totr,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1604         ierr = MPI_Allreduce(&eventInfo[event].count,         &minCt, 1, MPI_INT,             MPI_MIN, comm);CHKERRQ(ierr);
1605         ierr = MPI_Allreduce(&eventInfo[event].count,         &maxCt, 1, MPI_INT,             MPI_MAX, comm);CHKERRQ(ierr);
1606         name = stageLog->eventLog->eventInfo[event].name;
1607       } else {
1608         flopr = 0.0;
1609         ierr  = MPI_Allreduce(&flopr,                          &minf,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1610         ierr  = MPI_Allreduce(&flopr,                          &maxf,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1611         ierr  = MPI_Allreduce(&zero,                           &totf,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1612         ierr  = MPI_Allreduce(&zero,                           &mint,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1613         ierr  = MPI_Allreduce(&zero,                           &maxt,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1614         ierr  = MPI_Allreduce(&zero,                           &tott,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1615         ierr  = MPI_Allreduce(&zero,                           &totm,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1616         ierr  = MPI_Allreduce(&zero,                           &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1617         ierr  = MPI_Allreduce(&zero,                           &totr,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1618         ierr  = MPI_Allreduce(&ierr,                           &minCt, 1, MPI_INT,             MPI_MIN, comm);CHKERRQ(ierr);
1619         ierr  = MPI_Allreduce(&ierr,                           &maxCt, 1, MPI_INT,             MPI_MAX, comm);CHKERRQ(ierr);
1620         name  = "";
1621       }
1622       if (mint < 0.0) {
1623         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);
1624         mint = 0;
1625       }
1626       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);
1627       totm *= 0.5; totml *= 0.5; totr /= size;
1628 
1629       if (maxCt != 0) {
1630         if (minCt         != 0)   ratCt            = ((PetscLogDouble) maxCt)/minCt; else ratCt            = 0.0;
1631         if (mint          != 0.0) ratt             = maxt/mint;                  else ratt             = 0.0;
1632         if (minf          != 0.0) ratf             = maxf/minf;                  else ratf             = 0.0;
1633         if (TotalTime     != 0.0) fracTime         = tott/TotalTime;             else fracTime         = 0.0;
1634         if (TotalFlops    != 0.0) fracFlops        = totf/TotalFlops;            else fracFlops        = 0.0;
1635         if (stageTime     != 0.0) fracStageTime    = tott/stageTime;             else fracStageTime    = 0.0;
1636         if (flops         != 0.0) fracStageFlops   = totf/flops;                 else fracStageFlops   = 0.0;
1637         if (numMessages   != 0.0) fracMess         = totm/numMessages;           else fracMess         = 0.0;
1638         if (messageLength != 0.0) fracMessLen      = totml/messageLength;        else fracMessLen      = 0.0;
1639         if (numReductions != 0.0) fracRed          = totr/numReductions;         else fracRed          = 0.0;
1640         if (mess          != 0.0) fracStageMess    = totm/mess;                  else fracStageMess    = 0.0;
1641         if (messLen       != 0.0) fracStageMessLen = totml/messLen;              else fracStageMessLen = 0.0;
1642         if (red           != 0.0) fracStageRed     = totr/red;                   else fracStageRed     = 0.0;
1643         if (totm          != 0.0) totml           /= totm;                       else totml            = 0.0;
1644         if (maxt          != 0.0) flopr            = totf/maxt;                  else flopr            = 0.0;
1645         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);
1646         ierr = PetscFPrintf(comm, fd,
1647           "%-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",
1648                             name, maxCt, ratCt, maxt, ratt, maxf, ratf, totm, totml, totr,
1649                             100.0*fracTime, 100.0*fracFlops, 100.0*fracMess, 100.0*fracMessLen, 100.0*fracRed,
1650                             100.0*fracStageTime, 100.0*fracStageFlops, 100.0*fracStageMess, 100.0*fracStageMessLen, 100.0*fracStageRed,
1651                             PetscAbsReal(flopr/1.0e6));CHKERRQ(ierr);
1652       }
1653     }
1654   }
1655 
1656   /* Memory usage and object creation */
1657   ierr = PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------\n");CHKERRQ(ierr);
1658   ierr = PetscFPrintf(comm, fd, "\n");CHKERRQ(ierr);
1659   ierr = PetscFPrintf(comm, fd, "Memory usage is given in bytes:\n\n");CHKERRQ(ierr);
1660 
1661   /* Right now, only stages on the first processor are reported here, meaning only objects associated with
1662      the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then
1663      stats for stages local to processor sets.
1664   */
1665   /* We should figure out the longest object name here (now 20 characters) */
1666   ierr = PetscFPrintf(comm, fd, "Object Type          Creations   Destructions     Memory  Descendants' Mem.\n");CHKERRQ(ierr);
1667   ierr = PetscFPrintf(comm, fd, "Reports information only for process 0.\n");CHKERRQ(ierr);
1668   for (stage = 0; stage < numStages; stage++) {
1669     if (localStageUsed[stage]) {
1670       classInfo = stageLog->stageInfo[stage].classLog->classInfo;
1671       ierr = PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stageInfo[stage].name);CHKERRQ(ierr);
1672       for (oclass = 0; oclass < stageLog->stageInfo[stage].classLog->numClasses; oclass++) {
1673         if ((classInfo[oclass].creations > 0) || (classInfo[oclass].destructions > 0)) {
1674           ierr = PetscFPrintf(comm, fd, "%20s %5d          %5d  %11.0f     %g\n", stageLog->classLog->classInfo[oclass].name,
1675                               classInfo[oclass].creations, classInfo[oclass].destructions, classInfo[oclass].mem,
1676                               classInfo[oclass].descMem);CHKERRQ(ierr);
1677         }
1678       }
1679     } else {
1680       ierr = PetscFPrintf(comm, fd, "\n--- Event Stage %d: Unknown\n\n", stage);CHKERRQ(ierr);
1681     }
1682   }
1683 
1684   ierr = PetscFree(localStageUsed);CHKERRQ(ierr);
1685   ierr = PetscFree(stageUsed);CHKERRQ(ierr);
1686   ierr = PetscFree(localStageVisible);CHKERRQ(ierr);
1687   ierr = PetscFree(stageVisible);CHKERRQ(ierr);
1688 
1689   /* Information unrelated to this particular run */
1690   ierr = PetscFPrintf(comm, fd, "========================================================================================================================\n");CHKERRQ(ierr);
1691   PetscTime(&y);
1692   PetscTime(&x);
1693   PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y);
1694   PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y);
1695   ierr = PetscFPrintf(comm,fd,"Average time to get PetscTime(): %g\n", (y-x)/10.0);CHKERRQ(ierr);
1696   /* MPI information */
1697   if (size > 1) {
1698     MPI_Status  status;
1699     PetscMPIInt tag;
1700     MPI_Comm    newcomm;
1701 
1702     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1703     PetscTime(&x);
1704     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1705     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1706     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1707     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1708     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1709     PetscTime(&y);
1710     ierr = PetscFPrintf(comm, fd, "Average time for MPI_Barrier(): %g\n", (y-x)/5.0);CHKERRQ(ierr);
1711     ierr = PetscCommDuplicate(comm,&newcomm, &tag);CHKERRQ(ierr);
1712     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1713     if (rank) {
1714       ierr = MPI_Recv(0, 0, MPI_INT, rank-1,            tag, newcomm, &status);CHKERRQ(ierr);
1715       ierr = MPI_Send(0, 0, MPI_INT, (rank+1)%size, tag, newcomm);CHKERRQ(ierr);
1716     } else {
1717       PetscTime(&x);
1718       ierr = MPI_Send(0, 0, MPI_INT, 1,          tag, newcomm);CHKERRQ(ierr);
1719       ierr = MPI_Recv(0, 0, MPI_INT, size-1, tag, newcomm, &status);CHKERRQ(ierr);
1720       PetscTime(&y);
1721       ierr = PetscFPrintf(comm,fd,"Average time for zero size MPI_Send(): %g\n", (y-x)/size);CHKERRQ(ierr);
1722     }
1723     ierr = PetscCommDestroy(&newcomm);CHKERRQ(ierr);
1724   }
1725   ierr = PetscOptionsView(viewer);CHKERRQ(ierr);
1726 
1727   /* Machine and compile information */
1728 #if defined(PETSC_USE_FORTRAN_KERNELS)
1729   ierr = PetscFPrintf(comm, fd, "Compiled with FORTRAN kernels\n");CHKERRQ(ierr);
1730 #else
1731   ierr = PetscFPrintf(comm, fd, "Compiled without FORTRAN kernels\n");CHKERRQ(ierr);
1732 #endif
1733 #if defined(PETSC_USE_REAL_SINGLE)
1734   ierr = PetscFPrintf(comm, fd, "Compiled with single precision PetscScalar and PetscReal\n");CHKERRQ(ierr);
1735 #elif defined(PETSC_USE_LONGDOUBLE)
1736   ierr = PetscFPrintf(comm, fd, "Compiled with long double precision PetscScalar and PetscReal\n");CHKERRQ(ierr);
1737 #endif
1738 
1739 #if defined(PETSC_USE_REAL_MAT_SINGLE)
1740   ierr = PetscFPrintf(comm, fd, "Compiled with single precision matrices\n");CHKERRQ(ierr);
1741 #else
1742   ierr = PetscFPrintf(comm, fd, "Compiled with full precision matrices (default)\n");CHKERRQ(ierr);
1743 #endif
1744   ierr = PetscFPrintf(comm, fd, "sizeof(short) %d sizeof(int) %d sizeof(long) %d sizeof(void*) %d sizeof(PetscScalar) %d sizeof(PetscInt) %d\n",
1745                       (int) sizeof(short), (int) sizeof(int), (int) sizeof(long), (int) sizeof(void*),(int) sizeof(PetscScalar),(int) sizeof(PetscInt));CHKERRQ(ierr);
1746 
1747   ierr = PetscFPrintf(comm, fd, "Configure options: %s",petscconfigureoptions);CHKERRQ(ierr);
1748   ierr = PetscFPrintf(comm, fd, "%s", petscmachineinfo);CHKERRQ(ierr);
1749   ierr = PetscFPrintf(comm, fd, "%s", petsccompilerinfo);CHKERRQ(ierr);
1750   ierr = PetscFPrintf(comm, fd, "%s", petsccompilerflagsinfo);CHKERRQ(ierr);
1751   ierr = PetscFPrintf(comm, fd, "%s", petsclinkerinfo);CHKERRQ(ierr);
1752 
1753   /* Cleanup */
1754   ierr = PetscFPrintf(comm, fd, "\n");CHKERRQ(ierr);
1755   PetscFunctionReturn(0);
1756 }
1757 
1758 #undef __FUNCT__
1759 #define __FUNCT__ "PetscLogView"
1760 /*@C
1761   PetscLogView - Prints a summary of the logging.
1762 
1763   Collective over MPI_Comm
1764 
1765   Input Parameter:
1766 .  viewer - an ASCII viewer
1767 
1768   Options Database Keys:
1769 . -log_view [viewertype[:filename[:format]]] - Prints summary of log information (for code compiled with PETSC_USE_LOG)
1770 
1771   Usage:
1772 .vb
1773      PetscInitialize(...);
1774      PetscLogBegin();
1775      ... code ...
1776      PetscLogView(PetscViewer);
1777      PetscFinalize(...);
1778 .ve
1779 
1780   Notes:
1781   By default the summary is printed to stdout.
1782 
1783   Level: beginner
1784 
1785 .keywords: log, dump, print
1786 .seealso: PetscLogBegin(), PetscLogDump()
1787 @*/
1788 PetscErrorCode  PetscLogView(PetscViewer viewer)
1789 {
1790   PetscErrorCode    ierr;
1791   PetscBool         isascii;
1792   PetscViewerFormat format;
1793   int               stage, lastStage;
1794   PetscStageLog     stageLog;
1795 
1796   PetscFunctionBegin;
1797   if (!PetscLogBegin_PrivateCalled) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "No call to PetscLogBegin() before PetscLogView()");
1798   /* Pop off any stages the user forgot to remove */
1799   lastStage = 0;
1800   ierr      = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1801   ierr      = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
1802   while (stage >= 0) {
1803     lastStage = stage;
1804     ierr      = PetscStageLogPop(stageLog);CHKERRQ(ierr);
1805     ierr      = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
1806   }
1807   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1808   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Currently can only view logging to ASCII");
1809   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1810   if (format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_ASCII_INFO) {
1811     ierr = PetscLogView_Default(viewer);CHKERRQ(ierr);
1812   } else if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1813     ierr = PetscLogView_Detailed(viewer);CHKERRQ(ierr);
1814   }
1815   ierr = PetscStageLogPush(stageLog, lastStage);CHKERRQ(ierr);
1816   PetscFunctionReturn(0);
1817 }
1818 
1819 #undef __FUNCT__
1820 #define __FUNCT__ "PetscLogViewFromOptions"
1821 /*@C
1822   PetscLogViewFromOptions - Processes command line options to determine if/how a PetscLog is to be viewed.
1823 
1824   Collective on PETSC_COMM_WORLD
1825 
1826   Not normally called by user
1827 
1828   Level: intermediate
1829 
1830 @*/
1831 PetscErrorCode PetscLogViewFromOptions(void)
1832 {
1833   PetscErrorCode    ierr;
1834   PetscViewer       viewer;
1835   PetscBool         flg;
1836   PetscViewerFormat format;
1837 
1838   PetscFunctionBegin;
1839   ierr   = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,"-log_view",&viewer,&format,&flg);CHKERRQ(ierr);
1840   if (flg) {
1841     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1842     ierr = PetscLogView(viewer);CHKERRQ(ierr);
1843     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1844     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1845   }
1846   PetscFunctionReturn(0);
1847 }
1848 
1849 
1850 
1851 /*----------------------------------------------- Counter Functions -------------------------------------------------*/
1852 #undef __FUNCT__
1853 #define __FUNCT__ "PetscGetFlops"
1854 /*@C
1855    PetscGetFlops - Returns the number of flops used on this processor
1856    since the program began.
1857 
1858    Not Collective
1859 
1860    Output Parameter:
1861    flops - number of floating point operations
1862 
1863    Notes:
1864    A global counter logs all PETSc flop counts.  The user can use
1865    PetscLogFlops() to increment this counter to include flops for the
1866    application code.
1867 
1868    PETSc automatically logs library events if the code has been
1869    compiled with -DPETSC_USE_LOG (which is the default), and -log,
1870    -log_summary, or -log_all are specified.  PetscLogFlops() is
1871    intended for logging user flops to supplement this PETSc
1872    information.
1873 
1874    Level: intermediate
1875 
1876 .keywords: log, flops, floating point operations
1877 
1878 .seealso: PetscTime(), PetscLogFlops()
1879 @*/
1880 PetscErrorCode  PetscGetFlops(PetscLogDouble *flops)
1881 {
1882   PetscFunctionBegin;
1883   *flops = petsc_TotalFlops;
1884   PetscFunctionReturn(0);
1885 }
1886 
1887 #undef __FUNCT__
1888 #define __FUNCT__ "PetscLogObjectState"
1889 PetscErrorCode  PetscLogObjectState(PetscObject obj, const char format[], ...)
1890 {
1891   PetscErrorCode ierr;
1892   size_t         fullLength;
1893   va_list        Argp;
1894 
1895   PetscFunctionBegin;
1896   if (!petsc_logObjects) PetscFunctionReturn(0);
1897   va_start(Argp, format);
1898   ierr = PetscVSNPrintf(petsc_objects[obj->id].info, 64,format,&fullLength, Argp);CHKERRQ(ierr);
1899   va_end(Argp);
1900   PetscFunctionReturn(0);
1901 }
1902 
1903 
1904 /*MC
1905    PetscLogFlops - Adds floating point operations to the global counter.
1906 
1907    Synopsis:
1908    #include <petsclog.h>
1909    PetscErrorCode PetscLogFlops(PetscLogDouble f)
1910 
1911    Not Collective
1912 
1913    Input Parameter:
1914 .  f - flop counter
1915 
1916 
1917    Usage:
1918 .vb
1919      PetscLogEvent USER_EVENT;
1920      PetscLogEventRegister("User event",0,&USER_EVENT);
1921      PetscLogEventBegin(USER_EVENT,0,0,0,0);
1922         [code segment to monitor]
1923         PetscLogFlops(user_flops)
1924      PetscLogEventEnd(USER_EVENT,0,0,0,0);
1925 .ve
1926 
1927    Notes:
1928    A global counter logs all PETSc flop counts.  The user can use
1929    PetscLogFlops() to increment this counter to include flops for the
1930    application code.
1931 
1932    PETSc automatically logs library events if the code has been
1933    compiled with -DPETSC_USE_LOG (which is the default), and -log,
1934    -log_summary, or -log_all are specified.  PetscLogFlops() is
1935    intended for logging user flops to supplement this PETSc
1936    information.
1937 
1938    Level: intermediate
1939 
1940 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscGetFlops()
1941 
1942 .keywords: log, flops, floating point operations
1943 M*/
1944 
1945 /*MC
1946    PetscPreLoadBegin - Begin a segment of code that may be preloaded (run twice)
1947     to get accurate timings
1948 
1949    Synopsis:
1950    #include <petsclog.h>
1951    void PetscPreLoadBegin(PetscBool  flag,char *name);
1952 
1953    Not Collective
1954 
1955    Input Parameter:
1956 +   flag - PETSC_TRUE to run twice, PETSC_FALSE to run once, may be overridden
1957            with command line option -preload true or -preload false
1958 -   name - name of first stage (lines of code timed separately with -log_summary) to
1959            be preloaded
1960 
1961    Usage:
1962 .vb
1963      PetscPreLoadBegin(PETSC_TRUE,"first stage);
1964        lines of code
1965        PetscPreLoadStage("second stage");
1966        lines of code
1967      PetscPreLoadEnd();
1968 .ve
1969 
1970    Notes: Only works in C/C++, not Fortran
1971 
1972      Flags available within the macro.
1973 +    PetscPreLoadingUsed - true if we are or have done preloading
1974 .    PetscPreLoadingOn - true if it is CURRENTLY doing preload
1975 .    PetscPreLoadIt - 0 for the first computation (with preloading turned off it is only 0) 1 for the second
1976 -    PetscPreLoadMax - number of times it will do the computation, only one when preloading is turned on
1977      The first two variables are available throughout the program, the second two only between the PetscPreLoadBegin()
1978      and PetscPreLoadEnd()
1979 
1980    Level: intermediate
1981 
1982 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadEnd(), PetscPreLoadStage()
1983 
1984    Concepts: preloading
1985    Concepts: timing^accurate
1986    Concepts: paging^eliminating effects of
1987 
1988 
1989 M*/
1990 
1991 /*MC
1992    PetscPreLoadEnd - End a segment of code that may be preloaded (run twice)
1993     to get accurate timings
1994 
1995    Synopsis:
1996    #include <petsclog.h>
1997    void PetscPreLoadEnd(void);
1998 
1999    Not Collective
2000 
2001    Usage:
2002 .vb
2003      PetscPreLoadBegin(PETSC_TRUE,"first stage);
2004        lines of code
2005        PetscPreLoadStage("second stage");
2006        lines of code
2007      PetscPreLoadEnd();
2008 .ve
2009 
2010    Notes: only works in C/C++ not fortran
2011 
2012    Level: intermediate
2013 
2014 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadBegin(), PetscPreLoadStage()
2015 
2016 M*/
2017 
2018 /*MC
2019    PetscPreLoadStage - Start a new segment of code to be timed separately.
2020     to get accurate timings
2021 
2022    Synopsis:
2023    #include <petsclog.h>
2024    void PetscPreLoadStage(char *name);
2025 
2026    Not Collective
2027 
2028    Usage:
2029 .vb
2030      PetscPreLoadBegin(PETSC_TRUE,"first stage);
2031        lines of code
2032        PetscPreLoadStage("second stage");
2033        lines of code
2034      PetscPreLoadEnd();
2035 .ve
2036 
2037    Notes: only works in C/C++ not fortran
2038 
2039    Level: intermediate
2040 
2041 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadBegin(), PetscPreLoadEnd()
2042 
2043 M*/
2044 
2045 
2046 #else /* end of -DPETSC_USE_LOG section */
2047 
2048 #undef __FUNCT__
2049 #define __FUNCT__ "PetscLogObjectState"
2050 PetscErrorCode  PetscLogObjectState(PetscObject obj, const char format[], ...)
2051 {
2052   PetscFunctionBegin;
2053   PetscFunctionReturn(0);
2054 }
2055 
2056 #endif /* PETSC_USE_LOG*/
2057 
2058 
2059 PetscClassId PETSC_LARGEST_CLASSID = PETSC_SMALLEST_CLASSID;
2060 PetscClassId PETSC_OBJECT_CLASSID  = 0;
2061 
2062 #undef __FUNCT__
2063 #define __FUNCT__ "PetscClassIdRegister"
2064 /*@C
2065   PetscClassIdRegister - Registers a new class name for objects and logging operations in an application code.
2066 
2067   Not Collective
2068 
2069   Input Parameter:
2070 . name   - The class name
2071 
2072   Output Parameter:
2073 . oclass - The class id or classid
2074 
2075   Level: developer
2076 
2077 .keywords: log, class, register
2078 
2079 @*/
2080 PetscErrorCode  PetscClassIdRegister(const char name[],PetscClassId *oclass)
2081 {
2082 #if defined(PETSC_USE_LOG)
2083   PetscStageLog  stageLog;
2084   PetscInt       stage;
2085   PetscErrorCode ierr;
2086 #endif
2087 
2088   PetscFunctionBegin;
2089   *oclass = ++PETSC_LARGEST_CLASSID;
2090 #if defined(PETSC_USE_LOG)
2091   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
2092   ierr = PetscClassRegLogRegister(stageLog->classLog, name, *oclass);CHKERRQ(ierr);
2093   for (stage = 0; stage < stageLog->numStages; stage++) {
2094     ierr = ClassPerfLogEnsureSize(stageLog->stageInfo[stage].classLog, stageLog->classLog->numClasses);CHKERRQ(ierr);
2095   }
2096 #endif
2097   PetscFunctionReturn(0);
2098 }
2099 
2100 #if defined(PETSC_USE_LOG) && defined(PETSC_HAVE_MPE)
2101 #include <mpe.h>
2102 
2103 PetscBool PetscBeganMPE = PETSC_FALSE;
2104 
2105 PETSC_INTERN PetscErrorCode PetscLogEventBeginMPE(PetscLogEvent,int,PetscObject,PetscObject,PetscObject,PetscObject);
2106 PETSC_INTERN PetscErrorCode PetscLogEventEndMPE(PetscLogEvent,int,PetscObject,PetscObject,PetscObject,PetscObject);
2107 
2108 #undef __FUNCT__
2109 #define __FUNCT__ "PetscLogMPEBegin"
2110 /*@C
2111    PetscLogMPEBegin - Turns on MPE logging of events. This creates large log files
2112    and slows the program down.
2113 
2114    Collective over PETSC_COMM_WORLD
2115 
2116    Options Database Keys:
2117 . -log_mpe - Prints extensive log information (for code compiled with PETSC_USE_LOG)
2118 
2119    Notes:
2120    A related routine is PetscLogBegin() (with the options key -log_summary), which is
2121    intended for production runs since it logs only flop rates and object
2122    creation (and should not significantly slow the programs).
2123 
2124    Level: advanced
2125 
2126    Concepts: logging^MPE
2127    Concepts: logging^message passing
2128 
2129 .seealso: PetscLogDump(), PetscLogBegin(), PetscLogAllBegin(), PetscLogEventActivate(),
2130           PetscLogEventDeactivate()
2131 @*/
2132 PetscErrorCode  PetscLogMPEBegin(void)
2133 {
2134   PetscErrorCode ierr;
2135 
2136   PetscFunctionBegin;
2137   /* Do MPE initialization */
2138   if (!MPE_Initialized_logging()) { /* This function exists in mpich 1.1.2 and higher */
2139     ierr = PetscInfo(0,"Initializing MPE.\n");CHKERRQ(ierr);
2140     ierr = MPE_Init_log();CHKERRQ(ierr);
2141 
2142     PetscBeganMPE = PETSC_TRUE;
2143   } else {
2144     ierr = PetscInfo(0,"MPE already initialized. Not attempting to reinitialize.\n");CHKERRQ(ierr);
2145   }
2146   ierr = PetscLogSet(PetscLogEventBeginMPE, PetscLogEventEndMPE);CHKERRQ(ierr);
2147   PetscFunctionReturn(0);
2148 }
2149 
2150 #undef __FUNCT__
2151 #define __FUNCT__ "PetscLogMPEDump"
2152 /*@C
2153    PetscLogMPEDump - Dumps the MPE logging info to file for later use with Jumpshot.
2154 
2155    Collective over PETSC_COMM_WORLD
2156 
2157    Level: advanced
2158 
2159 .seealso: PetscLogDump(), PetscLogAllBegin(), PetscLogMPEBegin()
2160 @*/
2161 PetscErrorCode  PetscLogMPEDump(const char sname[])
2162 {
2163   char           name[PETSC_MAX_PATH_LEN];
2164   PetscErrorCode ierr;
2165 
2166   PetscFunctionBegin;
2167   if (PetscBeganMPE) {
2168     ierr = PetscInfo(0,"Finalizing MPE.\n");CHKERRQ(ierr);
2169     if (sname) {
2170       ierr = PetscStrcpy(name,sname);CHKERRQ(ierr);
2171     } else {
2172       ierr = PetscGetProgramName(name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
2173     }
2174     ierr = MPE_Finish_log(name);CHKERRQ(ierr);
2175   } else {
2176     ierr = PetscInfo(0,"Not finalizing MPE (not started by PETSc).\n");CHKERRQ(ierr);
2177   }
2178   PetscFunctionReturn(0);
2179 }
2180 
2181 #define PETSC_RGB_COLORS_MAX 39
2182 static const char *PetscLogMPERGBColors[PETSC_RGB_COLORS_MAX] = {
2183   "OliveDrab:      ",
2184   "BlueViolet:     ",
2185   "CadetBlue:      ",
2186   "CornflowerBlue: ",
2187   "DarkGoldenrod:  ",
2188   "DarkGreen:      ",
2189   "DarkKhaki:      ",
2190   "DarkOliveGreen: ",
2191   "DarkOrange:     ",
2192   "DarkOrchid:     ",
2193   "DarkSeaGreen:   ",
2194   "DarkSlateGray:  ",
2195   "DarkTurquoise:  ",
2196   "DeepPink:       ",
2197   "DarkKhaki:      ",
2198   "DimGray:        ",
2199   "DodgerBlue:     ",
2200   "GreenYellow:    ",
2201   "HotPink:        ",
2202   "IndianRed:      ",
2203   "LavenderBlush:  ",
2204   "LawnGreen:      ",
2205   "LemonChiffon:   ",
2206   "LightCoral:     ",
2207   "LightCyan:      ",
2208   "LightPink:      ",
2209   "LightSalmon:    ",
2210   "LightSlateGray: ",
2211   "LightYellow:    ",
2212   "LimeGreen:      ",
2213   "MediumPurple:   ",
2214   "MediumSeaGreen: ",
2215   "MediumSlateBlue:",
2216   "MidnightBlue:   ",
2217   "MintCream:      ",
2218   "MistyRose:      ",
2219   "NavajoWhite:    ",
2220   "NavyBlue:       ",
2221   "OliveDrab:      "
2222 };
2223 
2224 #undef __FUNCT__
2225 #define __FUNCT__ "PetscLogMPEGetRGBColor"
2226 /*@C
2227   PetscLogMPEGetRGBColor - This routine returns a rgb color useable with PetscLogEventRegister()
2228 
2229   Not collective. Maybe it should be?
2230 
2231   Output Parameter
2232 . str - character string representing the color
2233 
2234   Level: developer
2235 
2236 .keywords: log, mpe , color
2237 .seealso: PetscLogEventRegister
2238 @*/
2239 PetscErrorCode  PetscLogMPEGetRGBColor(const char *str[])
2240 {
2241   static int idx = 0;
2242 
2243   PetscFunctionBegin;
2244   *str = PetscLogMPERGBColors[idx];
2245   idx  = (idx + 1)% PETSC_RGB_COLORS_MAX;
2246   PetscFunctionReturn(0);
2247 }
2248 
2249 #endif /* PETSC_USE_LOG && PETSC_HAVE_MPE */
2250