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