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