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