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