xref: /petsc/src/sys/logging/plog.c (revision 5a856986583887c326abe5dfd149e8184a29cd80) !
1 
2 /*
3       PETSc code to log object creation and destruction and PETSc events.
4 
5       This provides the public API used by the rest of PETSc and by users.
6 
7       These routines use a private API that is not used elsewhere in PETSc and is not
8       accessible to users. The private API is defined in logimpl.h and the utils directory.
9 
10 */
11 #include <petsc/private/logimpl.h>        /*I    "petscsys.h"   I*/
12 #include <petsctime.h>
13 #include <petscviewer.h>
14 
15 PetscErrorCode PetscLogObjectParent(PetscObject p,PetscObject c)
16 {
17   if (!c || !p) return 0;
18   c->parent   = p;
19   c->parentid = p->id;
20   return 0;
21 }
22 
23 /*@C
24    PetscLogObjectMemory - Adds to an object a count of additional amount of memory that is used by the object.
25 
26    Not collective.
27 
28    Input Parameters:
29 +  obj  - the PETSc object
30 -  mem  - the amount of memory that is being added to the object
31 
32    Level: developer
33 
34    Developer Notes:
35     Currently we do not always do a good job of associating all memory allocations with an object.
36 
37 .seealso: PetscFinalize(), PetscInitializeFortran(), PetscGetArgs(), PetscInitializeNoArguments()
38 
39 @*/
40 PetscErrorCode PetscLogObjectMemory(PetscObject p,PetscLogDouble m)
41 {
42   if (!p) return 0;
43   p->mem += m;
44   return 0;
45 }
46 
47 PetscLogEvent PETSC_LARGEST_EVENT = PETSC_EVENT;
48 
49 #if defined(PETSC_USE_LOG)
50 #include <petscmachineinfo.h>
51 #include <petscconfiginfo.h>
52 
53 /* used in the MPI_XXX() count macros in petsclog.h */
54 
55 /* Action and object logging variables */
56 Action    *petsc_actions            = NULL;
57 Object    *petsc_objects            = NULL;
58 PetscBool petsc_logActions          = PETSC_FALSE;
59 PetscBool petsc_logObjects          = PETSC_FALSE;
60 int       petsc_numActions          = 0, petsc_maxActions = 100;
61 int       petsc_numObjects          = 0, petsc_maxObjects = 100;
62 int       petsc_numObjectsDestroyed = 0;
63 
64 /* Global counters */
65 PetscLogDouble petsc_BaseTime        = 0.0;
66 PetscLogDouble petsc_TotalFlops      = 0.0;  /* The number of flops */
67 PetscLogDouble petsc_tmp_flops       = 0.0;  /* The incremental number of flops */
68 PetscLogDouble petsc_send_ct         = 0.0;  /* The number of sends */
69 PetscLogDouble petsc_recv_ct         = 0.0;  /* The number of receives */
70 PetscLogDouble petsc_send_len        = 0.0;  /* The total length of all sent messages */
71 PetscLogDouble petsc_recv_len        = 0.0;  /* The total length of all received messages */
72 PetscLogDouble petsc_isend_ct        = 0.0;  /* The number of immediate sends */
73 PetscLogDouble petsc_irecv_ct        = 0.0;  /* The number of immediate receives */
74 PetscLogDouble petsc_isend_len       = 0.0;  /* The total length of all immediate send messages */
75 PetscLogDouble petsc_irecv_len       = 0.0;  /* The total length of all immediate receive messages */
76 PetscLogDouble petsc_wait_ct         = 0.0;  /* The number of waits */
77 PetscLogDouble petsc_wait_any_ct     = 0.0;  /* The number of anywaits */
78 PetscLogDouble petsc_wait_all_ct     = 0.0;  /* The number of waitalls */
79 PetscLogDouble petsc_sum_of_waits_ct = 0.0;  /* The total number of waits */
80 PetscLogDouble petsc_allreduce_ct    = 0.0;  /* The number of reductions */
81 PetscLogDouble petsc_gather_ct       = 0.0;  /* The number of gathers and gathervs */
82 PetscLogDouble petsc_scatter_ct      = 0.0;  /* The number of scatters and scattervs */
83 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
84 PetscLogDouble petsc_ctog_ct         = 0.0;  /* The total number of CPU to GPU copies */
85 PetscLogDouble petsc_gtoc_ct         = 0.0;  /* The total number of GPU to CPU copies */
86 PetscLogDouble petsc_ctog_sz         = 0.0;  /* The total size of CPU to GPU copies */
87 PetscLogDouble petsc_gtoc_sz         = 0.0;  /* The total size of GPU to CPU copies */
88 #endif
89 
90 /* Logging functions */
91 PetscErrorCode (*PetscLogPHC)(PetscObject) = NULL;
92 PetscErrorCode (*PetscLogPHD)(PetscObject) = NULL;
93 PetscErrorCode (*PetscLogPLB)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject) = NULL;
94 PetscErrorCode (*PetscLogPLE)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject) = NULL;
95 
96 /* Tracing event logging variables */
97 FILE             *petsc_tracefile            = NULL;
98 int              petsc_tracelevel            = 0;
99 const char       *petsc_traceblanks          = "                                                                                                    ";
100 char             petsc_tracespace[128]       = " ";
101 PetscLogDouble   petsc_tracetime             = 0.0;
102 static PetscBool PetscLogInitializeCalled = PETSC_FALSE;
103 
104 PETSC_INTERN PetscErrorCode PetscLogInitialize(void)
105 {
106   int            stage;
107   PetscBool      opt;
108   PetscErrorCode ierr;
109 
110   PetscFunctionBegin;
111   if (PetscLogInitializeCalled) PetscFunctionReturn(0);
112   PetscLogInitializeCalled = PETSC_TRUE;
113 
114   ierr = PetscOptionsHasName(NULL,NULL, "-log_exclude_actions", &opt);CHKERRQ(ierr);
115   if (opt) petsc_logActions = PETSC_FALSE;
116   ierr = PetscOptionsHasName(NULL,NULL, "-log_exclude_objects", &opt);CHKERRQ(ierr);
117   if (opt) petsc_logObjects = PETSC_FALSE;
118   if (petsc_logActions) {
119     ierr = PetscMalloc1(petsc_maxActions, &petsc_actions);CHKERRQ(ierr);
120   }
121   if (petsc_logObjects) {
122     ierr = PetscMalloc1(petsc_maxObjects, &petsc_objects);CHKERRQ(ierr);
123   }
124   PetscLogPHC = PetscLogObjCreateDefault;
125   PetscLogPHD = PetscLogObjDestroyDefault;
126   /* Setup default logging structures */
127   ierr = PetscStageLogCreate(&petsc_stageLog);CHKERRQ(ierr);
128   ierr = PetscStageLogRegister(petsc_stageLog, "Main Stage", &stage);CHKERRQ(ierr);
129 
130   /* All processors sync here for more consistent logging */
131   ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr);
132   PetscTime(&petsc_BaseTime);
133   ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
134   PetscFunctionReturn(0);
135 }
136 
137 PETSC_INTERN PetscErrorCode PetscLogFinalize(void)
138 {
139   PetscStageLog  stageLog;
140   PetscErrorCode ierr;
141 
142   PetscFunctionBegin;
143   ierr = PetscFree(petsc_actions);CHKERRQ(ierr);
144   ierr = PetscFree(petsc_objects);CHKERRQ(ierr);
145   ierr = PetscLogNestedEnd();CHKERRQ(ierr);
146   ierr = PetscLogSet(NULL, NULL);CHKERRQ(ierr);
147 
148   /* Resetting phase */
149   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
150   ierr = PetscStageLogDestroy(stageLog);CHKERRQ(ierr);
151 
152   petsc_TotalFlops            = 0.0;
153   petsc_numActions            = 0;
154   petsc_numObjects            = 0;
155   petsc_numObjectsDestroyed   = 0;
156   petsc_maxActions            = 100;
157   petsc_maxObjects            = 100;
158   petsc_actions               = NULL;
159   petsc_objects               = NULL;
160   petsc_logActions            = PETSC_FALSE;
161   petsc_logObjects            = PETSC_FALSE;
162   petsc_BaseTime              = 0.0;
163   petsc_TotalFlops            = 0.0;
164   petsc_tmp_flops             = 0.0;
165   petsc_send_ct               = 0.0;
166   petsc_recv_ct               = 0.0;
167   petsc_send_len              = 0.0;
168   petsc_recv_len              = 0.0;
169   petsc_isend_ct              = 0.0;
170   petsc_irecv_ct              = 0.0;
171   petsc_isend_len             = 0.0;
172   petsc_irecv_len             = 0.0;
173   petsc_wait_ct               = 0.0;
174   petsc_wait_any_ct           = 0.0;
175   petsc_wait_all_ct           = 0.0;
176   petsc_sum_of_waits_ct       = 0.0;
177   petsc_allreduce_ct          = 0.0;
178   petsc_gather_ct             = 0.0;
179   petsc_scatter_ct            = 0.0;
180   #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
181   petsc_ctog_ct               = 0.0;
182   petsc_gtoc_ct               = 0.0;
183   petsc_ctog_sz               = 0.0;
184   petsc_gtoc_sz               = 0.0;
185   #endif
186   PETSC_LARGEST_EVENT         = PETSC_EVENT;
187   PetscLogPHC                 = NULL;
188   PetscLogPHD                 = NULL;
189   petsc_tracefile             = NULL;
190   petsc_tracelevel            = 0;
191   petsc_traceblanks           = "                                                                                                    ";
192   petsc_tracespace[0]         = ' '; petsc_tracespace[1] = 0;
193   petsc_tracetime             = 0.0;
194   PETSC_LARGEST_CLASSID       = PETSC_SMALLEST_CLASSID;
195   PETSC_OBJECT_CLASSID        = 0;
196   petsc_stageLog              = 0;
197   PetscLogInitializeCalled    = PETSC_FALSE;
198   PetscFunctionReturn(0);
199 }
200 
201 /*@C
202   PetscLogSet - Sets the logging functions called at the beginning and ending of every event.
203 
204   Not Collective
205 
206   Input Parameters:
207 + b - The function called at beginning of event
208 - e - The function called at end of event
209 
210   Level: developer
211 
212 .seealso: PetscLogDump(), PetscLogDefaultBegin(), PetscLogAllBegin(), PetscLogTraceBegin()
213 @*/
214 PetscErrorCode  PetscLogSet(PetscErrorCode (*b)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject),
215                             PetscErrorCode (*e)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject))
216 {
217   PetscFunctionBegin;
218   PetscLogPLB = b;
219   PetscLogPLE = e;
220   PetscFunctionReturn(0);
221 }
222 
223 /*@C
224   PetscLogDefaultBegin - Turns on logging of objects and events. This logs flop
225   rates and object creation and should not slow programs down too much.
226   This routine may be called more than once.
227 
228   Logically Collective over PETSC_COMM_WORLD
229 
230   Options Database Keys:
231 . -log_view [viewertype:filename:viewerformat] - Prints summary of flop and timing information to the
232                   screen (for code configured with --with-log=1 (which is the default))
233 
234   Usage:
235 .vb
236       PetscInitialize(...);
237       PetscLogDefaultBegin();
238        ... code ...
239       PetscLogView(viewer); or PetscLogDump();
240       PetscFinalize();
241 .ve
242 
243   Notes:
244   PetscLogView(viewer) or PetscLogDump() actually cause the printing of
245   the logging information.
246 
247   Level: advanced
248 
249 .seealso: PetscLogDump(), PetscLogAllBegin(), PetscLogView(), PetscLogTraceBegin()
250 @*/
251 PetscErrorCode  PetscLogDefaultBegin(void)
252 {
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256   ierr = PetscLogSet(PetscLogEventBeginDefault, PetscLogEventEndDefault);CHKERRQ(ierr);
257   PetscFunctionReturn(0);
258 }
259 
260 /*@C
261   PetscLogAllBegin - Turns on extensive logging of objects and events. Logs
262   all events. This creates large log files and slows the program down.
263 
264   Logically Collective on PETSC_COMM_WORLD
265 
266   Options Database Keys:
267 . -log_all - Prints extensive log information
268 
269   Usage:
270 .vb
271      PetscInitialize(...);
272      PetscLogAllBegin();
273      ... code ...
274      PetscLogDump(filename);
275      PetscFinalize();
276 .ve
277 
278   Notes:
279   A related routine is PetscLogDefaultBegin() (with the options key -log), which is
280   intended for production runs since it logs only flop rates and object
281   creation (and shouldn't significantly slow the programs).
282 
283   Level: advanced
284 
285 .seealso: PetscLogDump(), PetscLogDefaultBegin(), PetscLogTraceBegin()
286 @*/
287 PetscErrorCode  PetscLogAllBegin(void)
288 {
289   PetscErrorCode ierr;
290 
291   PetscFunctionBegin;
292   ierr = PetscLogSet(PetscLogEventBeginComplete, PetscLogEventEndComplete);CHKERRQ(ierr);
293   PetscFunctionReturn(0);
294 }
295 
296 /*@
297   PetscLogTraceBegin - Activates trace logging.  Every time a PETSc event
298   begins or ends, the event name is printed.
299 
300   Logically Collective on PETSC_COMM_WORLD
301 
302   Input Parameter:
303 . file - The file to print trace in (e.g. stdout)
304 
305   Options Database Key:
306 . -log_trace [filename] - Activates PetscLogTraceBegin()
307 
308   Notes:
309   PetscLogTraceBegin() prints the processor number, the execution time (sec),
310   then "Event begin:" or "Event end:" followed by the event name.
311 
312   PetscLogTraceBegin() allows tracing of all PETSc calls, which is useful
313   to determine where a program is hanging without running in the
314   debugger.  Can be used in conjunction with the -info option.
315 
316   Level: intermediate
317 
318 .seealso: PetscLogDump(), PetscLogAllBegin(), PetscLogView(), PetscLogDefaultBegin()
319 @*/
320 PetscErrorCode  PetscLogTraceBegin(FILE *file)
321 {
322   PetscErrorCode ierr;
323 
324   PetscFunctionBegin;
325   petsc_tracefile = file;
326 
327   ierr = PetscLogSet(PetscLogEventBeginTrace, PetscLogEventEndTrace);CHKERRQ(ierr);
328   PetscFunctionReturn(0);
329 }
330 
331 /*@
332   PetscLogActions - Determines whether actions are logged for the graphical viewer.
333 
334   Not Collective
335 
336   Input Parameter:
337 . flag - PETSC_TRUE if actions are to be logged
338 
339   Level: intermediate
340 
341   Note: Logging of actions continues to consume more memory as the program
342   runs. Long running programs should consider turning this feature off.
343 
344   Options Database Keys:
345 . -log_exclude_actions - Turns off actions logging
346 
347 .seealso: PetscLogStagePush(), PetscLogStagePop()
348 @*/
349 PetscErrorCode  PetscLogActions(PetscBool flag)
350 {
351   PetscFunctionBegin;
352   petsc_logActions = flag;
353   PetscFunctionReturn(0);
354 }
355 
356 /*@
357   PetscLogObjects - Determines whether objects are logged for the graphical viewer.
358 
359   Not Collective
360 
361   Input Parameter:
362 . flag - PETSC_TRUE if objects are to be logged
363 
364   Level: intermediate
365 
366   Note: Logging of objects continues to consume more memory as the program
367   runs. Long running programs should consider turning this feature off.
368 
369   Options Database Keys:
370 . -log_exclude_objects - Turns off objects logging
371 
372 .seealso: PetscLogStagePush(), PetscLogStagePop()
373 @*/
374 PetscErrorCode  PetscLogObjects(PetscBool flag)
375 {
376   PetscFunctionBegin;
377   petsc_logObjects = flag;
378   PetscFunctionReturn(0);
379 }
380 
381 /*------------------------------------------------ Stage Functions --------------------------------------------------*/
382 /*@C
383   PetscLogStageRegister - Attaches a character string name to a logging stage.
384 
385   Not Collective
386 
387   Input Parameter:
388 . sname - The name to associate with that stage
389 
390   Output Parameter:
391 . stage - The stage number
392 
393   Level: intermediate
394 
395 .seealso: PetscLogStagePush(), PetscLogStagePop()
396 @*/
397 PetscErrorCode  PetscLogStageRegister(const char sname[],PetscLogStage *stage)
398 {
399   PetscStageLog  stageLog;
400   PetscLogEvent  event;
401   PetscErrorCode ierr;
402 
403   PetscFunctionBegin;
404   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
405   ierr = PetscStageLogRegister(stageLog, sname, stage);CHKERRQ(ierr);
406   /* Copy events already changed in the main stage, this sucks */
407   ierr = PetscEventPerfLogEnsureSize(stageLog->stageInfo[*stage].eventLog, stageLog->eventLog->numEvents);CHKERRQ(ierr);
408   for (event = 0; event < stageLog->eventLog->numEvents; event++) {
409     ierr = PetscEventPerfInfoCopy(&stageLog->stageInfo[0].eventLog->eventInfo[event],&stageLog->stageInfo[*stage].eventLog->eventInfo[event]);CHKERRQ(ierr);
410   }
411   ierr = PetscClassPerfLogEnsureSize(stageLog->stageInfo[*stage].classLog, stageLog->classLog->numClasses);CHKERRQ(ierr);
412   PetscFunctionReturn(0);
413 }
414 
415 /*@C
416   PetscLogStagePush - This function pushes a stage on the stack.
417 
418   Not Collective
419 
420   Input Parameter:
421 . stage - The stage on which to log
422 
423   Usage:
424   If the option -log_sumary is used to run the program containing the
425   following code, then 2 sets of summary data will be printed during
426   PetscFinalize().
427 .vb
428       PetscInitialize(int *argc,char ***args,0,0);
429       [stage 0 of code]
430       PetscLogStagePush(1);
431       [stage 1 of code]
432       PetscLogStagePop();
433       PetscBarrier(...);
434       [more stage 0 of code]
435       PetscFinalize();
436 .ve
437 
438   Notes:
439   Use PetscLogStageRegister() to register a stage.
440 
441   Level: intermediate
442 
443 .seealso: PetscLogStagePop(), PetscLogStageRegister(), PetscBarrier()
444 @*/
445 PetscErrorCode  PetscLogStagePush(PetscLogStage stage)
446 {
447   PetscStageLog  stageLog;
448   PetscErrorCode ierr;
449 
450   PetscFunctionBegin;
451   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
452   ierr = PetscStageLogPush(stageLog, stage);CHKERRQ(ierr);
453   PetscFunctionReturn(0);
454 }
455 
456 /*@C
457   PetscLogStagePop - This function pops a stage from the stack.
458 
459   Not Collective
460 
461   Usage:
462   If the option -log_sumary is used to run the program containing the
463   following code, then 2 sets of summary data will be printed during
464   PetscFinalize().
465 .vb
466       PetscInitialize(int *argc,char ***args,0,0);
467       [stage 0 of code]
468       PetscLogStagePush(1);
469       [stage 1 of code]
470       PetscLogStagePop();
471       PetscBarrier(...);
472       [more stage 0 of code]
473       PetscFinalize();
474 .ve
475 
476   Notes:
477   Use PetscLogStageRegister() to register a stage.
478 
479   Level: intermediate
480 
481 .seealso: PetscLogStagePush(), PetscLogStageRegister(), PetscBarrier()
482 @*/
483 PetscErrorCode  PetscLogStagePop(void)
484 {
485   PetscStageLog  stageLog;
486   PetscErrorCode ierr;
487 
488   PetscFunctionBegin;
489   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
490   ierr = PetscStageLogPop(stageLog);CHKERRQ(ierr);
491   PetscFunctionReturn(0);
492 }
493 
494 /*@
495   PetscLogStageSetActive - Determines stage activity for PetscLogEventBegin() and PetscLogEventEnd().
496 
497   Not Collective
498 
499   Input Parameters:
500 + stage    - The stage
501 - isActive - The activity flag, PETSC_TRUE for logging, else PETSC_FALSE (defaults to PETSC_TRUE)
502 
503   Level: intermediate
504 
505 .seealso: PetscLogStagePush(), PetscLogStagePop(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadBegin(), PetscPreLoadEnd(), PetscPreLoadStage()
506 @*/
507 PetscErrorCode  PetscLogStageSetActive(PetscLogStage stage, PetscBool isActive)
508 {
509   PetscStageLog  stageLog;
510   PetscErrorCode ierr;
511 
512   PetscFunctionBegin;
513   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
514   ierr = PetscStageLogSetActive(stageLog, stage, isActive);CHKERRQ(ierr);
515   PetscFunctionReturn(0);
516 }
517 
518 /*@
519   PetscLogStageGetActive - Returns stage activity for PetscLogEventBegin() and PetscLogEventEnd().
520 
521   Not Collective
522 
523   Input Parameter:
524 . stage    - The stage
525 
526   Output Parameter:
527 . isActive - The activity flag, PETSC_TRUE for logging, else PETSC_FALSE (defaults to PETSC_TRUE)
528 
529   Level: intermediate
530 
531 .seealso: PetscLogStagePush(), PetscLogStagePop(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadBegin(), PetscPreLoadEnd(), PetscPreLoadStage()
532 @*/
533 PetscErrorCode  PetscLogStageGetActive(PetscLogStage stage, PetscBool  *isActive)
534 {
535   PetscStageLog  stageLog;
536   PetscErrorCode ierr;
537 
538   PetscFunctionBegin;
539   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
540   ierr = PetscStageLogGetActive(stageLog, stage, isActive);CHKERRQ(ierr);
541   PetscFunctionReturn(0);
542 }
543 
544 /*@
545   PetscLogStageSetVisible - Determines stage visibility in PetscLogView()
546 
547   Not Collective
548 
549   Input Parameters:
550 + stage     - The stage
551 - isVisible - The visibility flag, PETSC_TRUE to print, else PETSC_FALSE (defaults to PETSC_TRUE)
552 
553   Level: intermediate
554 
555 .seealso: PetscLogStagePush(), PetscLogStagePop(), PetscLogView()
556 @*/
557 PetscErrorCode  PetscLogStageSetVisible(PetscLogStage stage, PetscBool isVisible)
558 {
559   PetscStageLog  stageLog;
560   PetscErrorCode ierr;
561 
562   PetscFunctionBegin;
563   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
564   ierr = PetscStageLogSetVisible(stageLog, stage, isVisible);CHKERRQ(ierr);
565   PetscFunctionReturn(0);
566 }
567 
568 /*@
569   PetscLogStageGetVisible - Returns stage visibility in PetscLogView()
570 
571   Not Collective
572 
573   Input Parameter:
574 . stage     - The stage
575 
576   Output Parameter:
577 . isVisible - The visibility flag, PETSC_TRUE to print, else PETSC_FALSE (defaults to PETSC_TRUE)
578 
579   Level: intermediate
580 
581 .seealso: PetscLogStagePush(), PetscLogStagePop(), PetscLogView()
582 @*/
583 PetscErrorCode  PetscLogStageGetVisible(PetscLogStage stage, PetscBool  *isVisible)
584 {
585   PetscStageLog  stageLog;
586   PetscErrorCode ierr;
587 
588   PetscFunctionBegin;
589   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
590   ierr = PetscStageLogGetVisible(stageLog, stage, isVisible);CHKERRQ(ierr);
591   PetscFunctionReturn(0);
592 }
593 
594 /*@C
595   PetscLogStageGetId - Returns the stage id when given the stage name.
596 
597   Not Collective
598 
599   Input Parameter:
600 . name  - The stage name
601 
602   Output Parameter:
603 . stage - The stage, , or -1 if no stage with that name exists
604 
605   Level: intermediate
606 
607 .seealso: PetscLogStagePush(), PetscLogStagePop(), PetscPreLoadBegin(), PetscPreLoadEnd(), PetscPreLoadStage()
608 @*/
609 PetscErrorCode  PetscLogStageGetId(const char name[], PetscLogStage *stage)
610 {
611   PetscStageLog  stageLog;
612   PetscErrorCode ierr;
613 
614   PetscFunctionBegin;
615   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
616   ierr = PetscStageLogGetStage(stageLog, name, stage);CHKERRQ(ierr);
617   PetscFunctionReturn(0);
618 }
619 
620 /*------------------------------------------------ Event Functions --------------------------------------------------*/
621 /*@C
622   PetscLogEventRegister - Registers an event name for logging operations in an application code.
623 
624   Not Collective
625 
626   Input Parameter:
627 + name   - The name associated with the event
628 - classid - The classid associated to the class for this event, obtain either with
629            PetscClassIdRegister() or use a predefined one such as KSP_CLASSID, SNES_CLASSID, the predefined ones
630            are only available in C code
631 
632   Output Parameter:
633 . event - The event id for use with PetscLogEventBegin() and PetscLogEventEnd().
634 
635   Example of Usage:
636 .vb
637       PetscLogEvent USER_EVENT;
638       PetscClassId classid;
639       PetscLogDouble user_event_flops;
640       PetscClassIdRegister("class name",&classid);
641       PetscLogEventRegister("User event name",classid,&USER_EVENT);
642       PetscLogEventBegin(USER_EVENT,0,0,0,0);
643          [code segment to monitor]
644          PetscLogFlops(user_event_flops);
645       PetscLogEventEnd(USER_EVENT,0,0,0,0);
646 .ve
647 
648   Notes:
649   PETSc automatically logs library events if the code has been
650   configured with --with-log (which is the default) and
651   -log_view or -log_all is specified.  PetscLogEventRegister() is
652   intended for logging user events to supplement this PETSc
653   information.
654 
655   PETSc can gather data for use with the utilities Jumpshot
656   (part of the MPICH distribution).  If PETSc has been compiled
657   with flag -DPETSC_HAVE_MPE (MPE is an additional utility within
658   MPICH), the user can employ another command line option, -log_mpe,
659   to create a logfile, "mpe.log", which can be visualized
660   Jumpshot.
661 
662   The classid is associated with each event so that classes of events
663   can be disabled simultaneously, such as all matrix events. The user
664   can either use an existing classid, such as MAT_CLASSID, or create
665   their own as shown in the example.
666 
667   If an existing event with the same name exists, its event handle is
668   returned instead of creating a new event.
669 
670   Level: intermediate
671 
672 .seealso: PetscLogEventBegin(), PetscLogEventEnd(), PetscLogFlops(),
673           PetscLogEventActivate(), PetscLogEventDeactivate(), PetscClassIdRegister()
674 @*/
675 PetscErrorCode  PetscLogEventRegister(const char name[],PetscClassId classid,PetscLogEvent *event)
676 {
677   PetscStageLog  stageLog;
678   int            stage;
679   PetscErrorCode ierr;
680 
681   PetscFunctionBegin;
682   *event = PETSC_DECIDE;
683   ierr   = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
684   ierr   = PetscEventRegLogGetEvent(stageLog->eventLog, name, event);CHKERRQ(ierr);
685   if (*event > 0) PetscFunctionReturn(0);
686   ierr   = PetscEventRegLogRegister(stageLog->eventLog, name, classid, event);CHKERRQ(ierr);
687   for (stage = 0; stage < stageLog->numStages; stage++) {
688     ierr = PetscEventPerfLogEnsureSize(stageLog->stageInfo[stage].eventLog, stageLog->eventLog->numEvents);CHKERRQ(ierr);
689     ierr = PetscClassPerfLogEnsureSize(stageLog->stageInfo[stage].classLog, stageLog->classLog->numClasses);CHKERRQ(ierr);
690   }
691   PetscFunctionReturn(0);
692 }
693 
694 /*@
695   PetscLogEventSetCollective - Indicates that a particular event is collective.
696 
697   Not Collective
698 
699   Input Parameter:
700 + event - The event id
701 - collective - Bolean flag indicating whether a particular event is collective
702 
703   Note:
704   New events returned from PetscLogEventRegister() are collective by default.
705 
706   Level: developer
707 
708 .seealso: PetscLogEventRegister()
709 @*/
710 PetscErrorCode PetscLogEventSetCollective(PetscLogEvent event,PetscBool collective)
711 {
712   PetscStageLog    stageLog;
713   PetscEventRegLog eventRegLog;
714   PetscErrorCode   ierr;
715 
716   PetscFunctionBegin;
717   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
718   ierr = PetscStageLogGetEventRegLog(stageLog,&eventRegLog);CHKERRQ(ierr);
719   if (event < 0 || event > eventRegLog->numEvents) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid event id");
720   eventRegLog->eventInfo[event].collective = collective;
721   PetscFunctionReturn(0);
722 }
723 
724 /*@
725   PetscLogEventIncludeClass - Activates event logging for a PETSc object class in every stage.
726 
727   Not Collective
728 
729   Input Parameter:
730 . classid - The object class, for example MAT_CLASSID, SNES_CLASSID, etc.
731 
732   Level: developer
733 
734 .seealso: PetscLogEventActivateClass(),PetscLogEventDeactivateClass(),PetscLogEventActivate(),PetscLogEventDeactivate()
735 @*/
736 PetscErrorCode  PetscLogEventIncludeClass(PetscClassId classid)
737 {
738   PetscStageLog  stageLog;
739   int            stage;
740   PetscErrorCode ierr;
741 
742   PetscFunctionBegin;
743   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
744   for (stage = 0; stage < stageLog->numStages; stage++) {
745     ierr = PetscEventPerfLogActivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid);CHKERRQ(ierr);
746   }
747   PetscFunctionReturn(0);
748 }
749 
750 /*@
751   PetscLogEventExcludeClass - Deactivates event logging for a PETSc object class in every stage.
752 
753   Not Collective
754 
755   Input Parameter:
756 . classid - The object class, for example MAT_CLASSID, SNES_CLASSID, etc.
757 
758   Level: developer
759 
760 .seealso: PetscLogEventDeactivateClass(),PetscLogEventActivateClass(),PetscLogEventDeactivate(),PetscLogEventActivate()
761 @*/
762 PetscErrorCode  PetscLogEventExcludeClass(PetscClassId classid)
763 {
764   PetscStageLog  stageLog;
765   int            stage;
766   PetscErrorCode ierr;
767 
768   PetscFunctionBegin;
769   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
770   for (stage = 0; stage < stageLog->numStages; stage++) {
771     ierr = PetscEventPerfLogDeactivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid);CHKERRQ(ierr);
772   }
773   PetscFunctionReturn(0);
774 }
775 
776 /*@
777   PetscLogEventActivate - Indicates that a particular event should be logged.
778 
779   Not Collective
780 
781   Input Parameter:
782 . event - The event id
783 
784   Usage:
785 .vb
786       PetscLogEventDeactivate(VEC_SetValues);
787         [code where you do not want to log VecSetValues()]
788       PetscLogEventActivate(VEC_SetValues);
789         [code where you do want to log VecSetValues()]
790 .ve
791 
792   Note:
793   The event may be either a pre-defined PETSc event (found in include/petsclog.h)
794   or an event number obtained with PetscLogEventRegister().
795 
796   Level: advanced
797 
798 .seealso: PlogEventDeactivate()
799 @*/
800 PetscErrorCode  PetscLogEventActivate(PetscLogEvent event)
801 {
802   PetscStageLog  stageLog;
803   int            stage;
804   PetscErrorCode ierr;
805 
806   PetscFunctionBegin;
807   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
808   ierr = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
809   ierr = PetscEventPerfLogActivate(stageLog->stageInfo[stage].eventLog, event);CHKERRQ(ierr);
810   PetscFunctionReturn(0);
811 }
812 
813 /*@
814   PetscLogEventDeactivate - Indicates that a particular event should not be logged.
815 
816   Not Collective
817 
818   Input Parameter:
819 . event - The event id
820 
821   Usage:
822 .vb
823       PetscLogEventDeactivate(VEC_SetValues);
824         [code where you do not want to log VecSetValues()]
825       PetscLogEventActivate(VEC_SetValues);
826         [code where you do want to log VecSetValues()]
827 .ve
828 
829   Note:
830   The event may be either a pre-defined PETSc event (found in
831   include/petsclog.h) or an event number obtained with PetscLogEventRegister()).
832 
833   Level: advanced
834 
835 .seealso: PlogEventActivate()
836 @*/
837 PetscErrorCode  PetscLogEventDeactivate(PetscLogEvent event)
838 {
839   PetscStageLog  stageLog;
840   int            stage;
841   PetscErrorCode ierr;
842 
843   PetscFunctionBegin;
844   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
845   ierr = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
846   ierr = PetscEventPerfLogDeactivate(stageLog->stageInfo[stage].eventLog, event);CHKERRQ(ierr);
847   PetscFunctionReturn(0);
848 }
849 
850 /*@
851   PetscLogEventSetActiveAll - Sets the event activity in every stage.
852 
853   Not Collective
854 
855   Input Parameters:
856 + event    - The event id
857 - isActive - The activity flag determining whether the event is logged
858 
859   Level: advanced
860 
861 .seealso: PlogEventActivate(),PlogEventDeactivate()
862 @*/
863 PetscErrorCode  PetscLogEventSetActiveAll(PetscLogEvent event, PetscBool isActive)
864 {
865   PetscStageLog  stageLog;
866   int            stage;
867   PetscErrorCode ierr;
868 
869   PetscFunctionBegin;
870   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
871   for (stage = 0; stage < stageLog->numStages; stage++) {
872     if (isActive) {
873       ierr = PetscEventPerfLogActivate(stageLog->stageInfo[stage].eventLog, event);CHKERRQ(ierr);
874     } else {
875       ierr = PetscEventPerfLogDeactivate(stageLog->stageInfo[stage].eventLog, event);CHKERRQ(ierr);
876     }
877   }
878   PetscFunctionReturn(0);
879 }
880 
881 /*@
882   PetscLogEventActivateClass - Activates event logging for a PETSc object class.
883 
884   Not Collective
885 
886   Input Parameter:
887 . classid - The event class, for example MAT_CLASSID, SNES_CLASSID, etc.
888 
889   Level: developer
890 
891 .seealso: PetscLogEventDeactivateClass(),PetscLogEventActivate(),PetscLogEventDeactivate()
892 @*/
893 PetscErrorCode  PetscLogEventActivateClass(PetscClassId classid)
894 {
895   PetscStageLog  stageLog;
896   int            stage;
897   PetscErrorCode ierr;
898 
899   PetscFunctionBegin;
900   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
901   ierr = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
902   ierr = PetscEventPerfLogActivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid);CHKERRQ(ierr);
903   PetscFunctionReturn(0);
904 }
905 
906 /*@
907   PetscLogEventDeactivateClass - Deactivates event logging for a PETSc object class.
908 
909   Not Collective
910 
911   Input Parameter:
912 . classid - The event class, for example MAT_CLASSID, SNES_CLASSID, etc.
913 
914   Level: developer
915 
916 .seealso: PetscLogEventActivateClass(),PetscLogEventActivate(),PetscLogEventDeactivate()
917 @*/
918 PetscErrorCode  PetscLogEventDeactivateClass(PetscClassId classid)
919 {
920   PetscStageLog  stageLog;
921   int            stage;
922   PetscErrorCode ierr;
923 
924   PetscFunctionBegin;
925   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
926   ierr = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
927   ierr = PetscEventPerfLogDeactivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid);CHKERRQ(ierr);
928   PetscFunctionReturn(0);
929 }
930 
931 /*MC
932    PetscLogEventSync - Synchronizes the beginning of a user event.
933 
934    Synopsis:
935    #include <petsclog.h>
936    PetscErrorCode PetscLogEventSync(int e,MPI_Comm comm)
937 
938    Collective
939 
940    Input Parameters:
941 +  e - integer associated with the event obtained from PetscLogEventRegister()
942 -  comm - an MPI communicator
943 
944    Usage:
945 .vb
946      PetscLogEvent USER_EVENT;
947      PetscLogEventRegister("User event",0,&USER_EVENT);
948      PetscLogEventSync(USER_EVENT,PETSC_COMM_WORLD);
949      PetscLogEventBegin(USER_EVENT,0,0,0,0);
950         [code segment to monitor]
951      PetscLogEventEnd(USER_EVENT,0,0,0,0);
952 .ve
953 
954    Notes:
955    This routine should be called only if there is not a
956    PetscObject available to pass to PetscLogEventBegin().
957 
958    Level: developer
959 
960 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd()
961 
962 M*/
963 
964 /*MC
965    PetscLogEventBegin - Logs the beginning of a user event.
966 
967    Synopsis:
968    #include <petsclog.h>
969    PetscErrorCode PetscLogEventBegin(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
970 
971    Not Collective
972 
973    Input Parameters:
974 +  e - integer associated with the event obtained from PetscLogEventRegister()
975 -  o1,o2,o3,o4 - objects associated with the event, or 0
976 
977 
978    Fortran Synopsis:
979    void PetscLogEventBegin(int e,PetscErrorCode ierr)
980 
981    Usage:
982 .vb
983      PetscLogEvent USER_EVENT;
984      PetscLogDouble user_event_flops;
985      PetscLogEventRegister("User event",0,&USER_EVENT);
986      PetscLogEventBegin(USER_EVENT,0,0,0,0);
987         [code segment to monitor]
988         PetscLogFlops(user_event_flops);
989      PetscLogEventEnd(USER_EVENT,0,0,0,0);
990 .ve
991 
992    Notes:
993    You need to register each integer event with the command
994    PetscLogEventRegister().
995 
996    Level: intermediate
997 
998 .seealso: PetscLogEventRegister(), PetscLogEventEnd(), PetscLogFlops()
999 
1000 M*/
1001 
1002 /*MC
1003    PetscLogEventEnd - Log the end of a user event.
1004 
1005    Synopsis:
1006    #include <petsclog.h>
1007    PetscErrorCode PetscLogEventEnd(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
1008 
1009    Not Collective
1010 
1011    Input Parameters:
1012 +  e - integer associated with the event obtained with PetscLogEventRegister()
1013 -  o1,o2,o3,o4 - objects associated with the event, or 0
1014 
1015 
1016    Fortran Synopsis:
1017    void PetscLogEventEnd(int e,PetscErrorCode ierr)
1018 
1019    Usage:
1020 .vb
1021      PetscLogEvent USER_EVENT;
1022      PetscLogDouble user_event_flops;
1023      PetscLogEventRegister("User event",0,&USER_EVENT,);
1024      PetscLogEventBegin(USER_EVENT,0,0,0,0);
1025         [code segment to monitor]
1026         PetscLogFlops(user_event_flops);
1027      PetscLogEventEnd(USER_EVENT,0,0,0,0);
1028 .ve
1029 
1030    Notes:
1031    You should also register each additional integer event with the command
1032    PetscLogEventRegister().
1033 
1034    Level: intermediate
1035 
1036 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogFlops()
1037 
1038 M*/
1039 
1040 /*@C
1041   PetscLogEventGetId - Returns the event id when given the event name.
1042 
1043   Not Collective
1044 
1045   Input Parameter:
1046 . name  - The event name
1047 
1048   Output Parameter:
1049 . event - The event, or -1 if no event with that name exists
1050 
1051   Level: intermediate
1052 
1053 .seealso: PetscLogEventBegin(), PetscLogEventEnd(), PetscLogStageGetId()
1054 @*/
1055 PetscErrorCode  PetscLogEventGetId(const char name[], PetscLogEvent *event)
1056 {
1057   PetscStageLog  stageLog;
1058   PetscErrorCode ierr;
1059 
1060   PetscFunctionBegin;
1061   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1062   ierr = PetscEventRegLogGetEvent(stageLog->eventLog, name, event);CHKERRQ(ierr);
1063   PetscFunctionReturn(0);
1064 }
1065 
1066 
1067 /*------------------------------------------------ Output Functions -------------------------------------------------*/
1068 /*@C
1069   PetscLogDump - Dumps logs of objects to a file. This file is intended to
1070   be read by bin/petscview. This program no longer exists.
1071 
1072   Collective on PETSC_COMM_WORLD
1073 
1074   Input Parameter:
1075 . name - an optional file name
1076 
1077   Usage:
1078 .vb
1079      PetscInitialize(...);
1080      PetscLogDefaultBegin(); or PetscLogAllBegin();
1081      ... code ...
1082      PetscLogDump(filename);
1083      PetscFinalize();
1084 .ve
1085 
1086   Notes:
1087   The default file name is
1088 $    Log.<rank>
1089   where <rank> is the processor number. If no name is specified,
1090   this file will be used.
1091 
1092   Level: advanced
1093 
1094 .seealso: PetscLogDefaultBegin(), PetscLogAllBegin(), PetscLogView()
1095 @*/
1096 PetscErrorCode  PetscLogDump(const char sname[])
1097 {
1098   PetscStageLog      stageLog;
1099   PetscEventPerfInfo *eventInfo;
1100   FILE               *fd;
1101   char               file[PETSC_MAX_PATH_LEN], fname[PETSC_MAX_PATH_LEN];
1102   PetscLogDouble     flops, _TotalTime;
1103   PetscMPIInt        rank;
1104   int                action, object, curStage;
1105   PetscLogEvent      event;
1106   PetscErrorCode     ierr;
1107 
1108   PetscFunctionBegin;
1109   /* Calculate the total elapsed time */
1110   PetscTime(&_TotalTime);
1111   _TotalTime -= petsc_BaseTime;
1112   /* Open log file */
1113   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
1114   if (sname && sname[0]) sprintf(file, "%s.%d", sname, rank);
1115   else sprintf(file, "Log.%d", rank);
1116   ierr = PetscFixFilename(file, fname);CHKERRQ(ierr);
1117   ierr = PetscFOpen(PETSC_COMM_WORLD, fname, "w", &fd);CHKERRQ(ierr);
1118   if ((!rank) && (!fd)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN, "Cannot open file: %s", fname);
1119   /* Output totals */
1120   ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Total Flop %14e %16.8e\n", petsc_TotalFlops, _TotalTime);CHKERRQ(ierr);
1121   ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Clock Resolution %g\n", 0.0);CHKERRQ(ierr);
1122   /* Output actions */
1123   if (petsc_logActions) {
1124     ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Actions accomplished %d\n", petsc_numActions);CHKERRQ(ierr);
1125     for (action = 0; action < petsc_numActions; action++) {
1126       ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "%g %d %d %d %d %d %d %g %g %g\n",
1127                           petsc_actions[action].time, petsc_actions[action].action, (int)petsc_actions[action].event, (int)petsc_actions[action].classid, petsc_actions[action].id1,
1128                           petsc_actions[action].id2, petsc_actions[action].id3, petsc_actions[action].flops, petsc_actions[action].mem, petsc_actions[action].maxmem);CHKERRQ(ierr);
1129     }
1130   }
1131   /* Output objects */
1132   if (petsc_logObjects) {
1133     ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Objects created %d destroyed %d\n", petsc_numObjects, petsc_numObjectsDestroyed);CHKERRQ(ierr);
1134     for (object = 0; object < petsc_numObjects; object++) {
1135       ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Parent ID: %d Memory: %d\n", petsc_objects[object].parent, (int) petsc_objects[object].mem);CHKERRQ(ierr);
1136       if (!petsc_objects[object].name[0]) {
1137         ierr = PetscFPrintf(PETSC_COMM_WORLD, fd,"No Name\n");CHKERRQ(ierr);
1138       } else {
1139         ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Name: %s\n", petsc_objects[object].name);CHKERRQ(ierr);
1140       }
1141       if (petsc_objects[object].info[0] != 0) {
1142         ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "No Info\n");CHKERRQ(ierr);
1143       } else {
1144         ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Info: %s\n", petsc_objects[object].info);CHKERRQ(ierr);
1145       }
1146     }
1147   }
1148   /* Output events */
1149   ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "Event log:\n");CHKERRQ(ierr);
1150   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1151   ierr = PetscIntStackTop(stageLog->stack, &curStage);CHKERRQ(ierr);
1152   eventInfo = stageLog->stageInfo[curStage].eventLog->eventInfo;
1153   for (event = 0; event < stageLog->stageInfo[curStage].eventLog->numEvents; event++) {
1154     if (eventInfo[event].time != 0.0) flops = eventInfo[event].flops/eventInfo[event].time;
1155     else flops = 0.0;
1156     ierr = PetscFPrintf(PETSC_COMM_WORLD, fd, "%d %16d %16g %16g %16g\n", event, eventInfo[event].count,
1157                         eventInfo[event].flops, eventInfo[event].time, flops);CHKERRQ(ierr);
1158   }
1159   ierr = PetscFClose(PETSC_COMM_WORLD, fd);CHKERRQ(ierr);
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 /*
1164   PetscLogView_Detailed - Each process prints the times for its own events
1165 
1166 */
1167 PetscErrorCode  PetscLogView_Detailed(PetscViewer viewer)
1168 {
1169   PetscStageLog      stageLog;
1170   PetscEventPerfInfo *eventInfo = NULL, *stageInfo = NULL;
1171   PetscLogDouble     locTotalTime, numRed, maxMem;
1172   int                numStages,numEvents,stage,event;
1173   MPI_Comm           comm = PetscObjectComm((PetscObject) viewer);
1174   PetscMPIInt        rank,size;
1175   PetscErrorCode     ierr;
1176 
1177   PetscFunctionBegin;
1178   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1179   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1180   /* Must preserve reduction count before we go on */
1181   numRed = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
1182   /* Get the total elapsed time */
1183   PetscTime(&locTotalTime);  locTotalTime -= petsc_BaseTime;
1184   ierr = PetscViewerASCIIPrintf(viewer,"size = %d\n",size);CHKERRQ(ierr);
1185   ierr = PetscViewerASCIIPrintf(viewer,"LocalTimes = {}\n");CHKERRQ(ierr);
1186   ierr = PetscViewerASCIIPrintf(viewer,"LocalMessages = {}\n");CHKERRQ(ierr);
1187   ierr = PetscViewerASCIIPrintf(viewer,"LocalMessageLens = {}\n");CHKERRQ(ierr);
1188   ierr = PetscViewerASCIIPrintf(viewer,"LocalReductions = {}\n");CHKERRQ(ierr);
1189   ierr = PetscViewerASCIIPrintf(viewer,"LocalFlop = {}\n");CHKERRQ(ierr);
1190   ierr = PetscViewerASCIIPrintf(viewer,"LocalObjects = {}\n");CHKERRQ(ierr);
1191   ierr = PetscViewerASCIIPrintf(viewer,"LocalMemory = {}\n");CHKERRQ(ierr);
1192   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1193   ierr = MPIU_Allreduce(&stageLog->numStages, &numStages, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1194   ierr = PetscViewerASCIIPrintf(viewer,"Stages = {}\n");CHKERRQ(ierr);
1195   for (stage=0; stage<numStages; stage++) {
1196     ierr = PetscViewerASCIIPrintf(viewer,"Stages[\"%s\"] = {}\n",stageLog->stageInfo[stage].name);CHKERRQ(ierr);
1197     ierr = PetscViewerASCIIPrintf(viewer,"Stages[\"%s\"][\"summary\"] = {}\n",stageLog->stageInfo[stage].name);CHKERRQ(ierr);
1198     ierr = MPIU_Allreduce(&stageLog->stageInfo[stage].eventLog->numEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1199     for (event = 0; event < numEvents; event++) {
1200       ierr = PetscViewerASCIIPrintf(viewer,"Stages[\"%s\"][\"%s\"] = {}\n",stageLog->stageInfo[stage].name,stageLog->eventLog->eventInfo[event].name);CHKERRQ(ierr);
1201     }
1202   }
1203   ierr = PetscMallocGetMaximumUsage(&maxMem);CHKERRQ(ierr);
1204   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
1205   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalTimes[%d] = %g\n",rank,locTotalTime);CHKERRQ(ierr);
1206   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalMessages[%d] = %g\n",rank,(petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct));CHKERRQ(ierr);
1207   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalMessageLens[%d] = %g\n",rank,(petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len));CHKERRQ(ierr);
1208   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalReductions[%d] = %g\n",rank,numRed);CHKERRQ(ierr);
1209   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalFlop[%d] = %g\n",rank,petsc_TotalFlops);CHKERRQ(ierr);
1210   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalObjects[%d] = %d\n",rank,petsc_numObjects);CHKERRQ(ierr);
1211   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"LocalMemory[%d] = %g\n",rank,maxMem);CHKERRQ(ierr);
1212   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1213   for (stage=0; stage<numStages; stage++) {
1214     stageInfo = &stageLog->stageInfo[stage].perfInfo;
1215     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Stages[\"%s\"][\"summary\"][%d] = {\"time\" : %g, \"numMessages\" : %g, \"messageLength\" : %g, \"numReductions\" : %g, \"flop\" : %g}\n",
1216                                               stageLog->stageInfo[stage].name,rank,
1217                                               stageInfo->time,stageInfo->numMessages,stageInfo->messageLength,stageInfo->numReductions,stageInfo->flops);CHKERRQ(ierr);
1218     ierr = MPIU_Allreduce(&stageLog->stageInfo[stage].eventLog->numEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1219     for (event = 0; event < numEvents; event++) {
1220       eventInfo = &stageLog->stageInfo[stage].eventLog->eventInfo[event];
1221       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Stages[\"%s\"][\"%s\"][%d] = {\"count\" : %D, \"time\" : %g, \"syncTime\" : %g, \"numMessages\" : %g, \"messageLength\" : %g, \"numReductions\" : %g, \"flop\" : %g",
1222                                                 stageLog->stageInfo[stage].name,stageLog->eventLog->eventInfo[event].name,rank,
1223                                                 eventInfo->count,eventInfo->time,eventInfo->syncTime,eventInfo->numMessages,eventInfo->messageLength,eventInfo->numReductions,eventInfo->flops);CHKERRQ(ierr);
1224       if (eventInfo->dof[0] >= 0.) {
1225         PetscInt d, e;
1226 
1227         ierr = PetscViewerASCIISynchronizedPrintf(viewer, ", \"dof\" : [");CHKERRQ(ierr);
1228         for (d = 0; d < 8; ++d) {
1229           if (d > 0) {ierr = PetscViewerASCIISynchronizedPrintf(viewer, ", ");CHKERRQ(ierr);}
1230           ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->dof[d]);CHKERRQ(ierr);
1231         }
1232         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "]");CHKERRQ(ierr);
1233         ierr = PetscViewerASCIISynchronizedPrintf(viewer, ", \"error\" : [");CHKERRQ(ierr);
1234         for (e = 0; e < 8; ++e) {
1235           if (e > 0) {ierr = PetscViewerASCIISynchronizedPrintf(viewer, ", ");CHKERRQ(ierr);}
1236           ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->errors[e]);CHKERRQ(ierr);
1237         }
1238         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "]");CHKERRQ(ierr);
1239       }
1240       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"}\n");CHKERRQ(ierr);
1241     }
1242   }
1243   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1244   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
1245   PetscFunctionReturn(0);
1246 }
1247 
1248 /*
1249   PetscLogView_CSV - Each process prints the times for its own events in Comma-Separated Value Format
1250 */
1251 PetscErrorCode  PetscLogView_CSV(PetscViewer viewer)
1252 {
1253   PetscStageLog      stageLog;
1254   PetscEventPerfInfo *eventInfo = NULL;
1255   PetscLogDouble     locTotalTime, maxMem;
1256   int                numStages,numEvents,stage,event;
1257   MPI_Comm           comm = PetscObjectComm((PetscObject) viewer);
1258   PetscMPIInt        rank,size;
1259   PetscErrorCode     ierr;
1260 
1261   PetscFunctionBegin;
1262   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1263   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1264   /* Must preserve reduction count before we go on */
1265   /* Get the total elapsed time */
1266   PetscTime(&locTotalTime);  locTotalTime -= petsc_BaseTime;
1267   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1268   ierr = MPIU_Allreduce(&stageLog->numStages, &numStages, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1269   ierr = PetscMallocGetMaximumUsage(&maxMem);CHKERRQ(ierr);
1270   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
1271   ierr = PetscViewerASCIIPrintf(viewer,"Stage Name,Event Name,Rank,Time,Num Messages,Message Length,Num Reductions,FLOP,dof0,dof1,dof2,dof3,dof4,dof5,dof6,dof7,e0,e1,e2,e3,e4,e5,e6,e7,%d\n", size);
1272   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1273   for (stage=0; stage<numStages; stage++) {
1274     PetscEventPerfInfo *stageInfo = &stageLog->stageInfo[stage].perfInfo;
1275 
1276     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%s,summary,%d,%g,%g,%g,%g,%g\n",
1277                                               stageLog->stageInfo[stage].name,rank,stageInfo->time,stageInfo->numMessages,stageInfo->messageLength,stageInfo->numReductions,stageInfo->flops);CHKERRQ(ierr);
1278     ierr = MPIU_Allreduce(&stageLog->stageInfo[stage].eventLog->numEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1279     for (event = 0; event < numEvents; event++) {
1280       eventInfo = &stageLog->stageInfo[stage].eventLog->eventInfo[event];
1281       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%s,%s,%d,%g,%g,%g,%g,%g",stageLog->stageInfo[stage].name,
1282                                                 stageLog->eventLog->eventInfo[event].name,rank,eventInfo->time,eventInfo->numMessages,
1283                                                 eventInfo->messageLength,eventInfo->numReductions,eventInfo->flops);CHKERRQ(ierr);
1284       if (eventInfo->dof[0] >= 0.) {
1285         PetscInt d, e;
1286 
1287         for (d = 0; d < 8; ++d) {
1288           ierr = PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->dof[d]);CHKERRQ(ierr);
1289         }
1290         for (e = 0; e < 8; ++e) {
1291           ierr = PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->errors[e]);CHKERRQ(ierr);
1292         }
1293       }
1294       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1295     }
1296   }
1297   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1298   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
1299   PetscFunctionReturn(0);
1300 }
1301 
1302 static PetscErrorCode PetscLogViewWarnSync(MPI_Comm comm,FILE *fd)
1303 {
1304   PetscErrorCode ierr;
1305   PetscFunctionBegin;
1306   if (!PetscLogSyncOn) PetscFunctionReturn(0);
1307   ierr = PetscFPrintf(comm, fd, "\n\n");CHKERRQ(ierr);
1308   ierr = PetscFPrintf(comm, fd, "      ##########################################################\n");CHKERRQ(ierr);
1309   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1310   ierr = PetscFPrintf(comm, fd, "      #                       WARNING!!!                       #\n");CHKERRQ(ierr);
1311   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1312   ierr = PetscFPrintf(comm, fd, "      #   This program was run with logging synchronization.   #\n");CHKERRQ(ierr);
1313   ierr = PetscFPrintf(comm, fd, "      #   This option provides more meaningful imbalance       #\n");CHKERRQ(ierr);
1314   ierr = PetscFPrintf(comm, fd, "      #   figures at the expense of slowing things down and    #\n");CHKERRQ(ierr);
1315   ierr = PetscFPrintf(comm, fd, "      #   providing a distorted view of the overall runtime.   #\n");CHKERRQ(ierr);
1316   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1317   ierr = PetscFPrintf(comm, fd, "      ##########################################################\n\n\n");CHKERRQ(ierr);
1318   PetscFunctionReturn(0);
1319   PetscFunctionReturn(0);
1320 }
1321 
1322 static PetscErrorCode PetscLogViewWarnDebugging(MPI_Comm comm,FILE *fd)
1323 {
1324 #if defined(PETSC_USE_DEBUG)
1325   PetscErrorCode ierr;
1326 
1327   PetscFunctionBegin;
1328   ierr = PetscFPrintf(comm, fd, "\n\n");CHKERRQ(ierr);
1329   ierr = PetscFPrintf(comm, fd, "      ##########################################################\n");CHKERRQ(ierr);
1330   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1331   ierr = PetscFPrintf(comm, fd, "      #                       WARNING!!!                       #\n");CHKERRQ(ierr);
1332   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1333   ierr = PetscFPrintf(comm, fd, "      #   This code was compiled with a debugging option.      #\n");CHKERRQ(ierr);
1334   ierr = PetscFPrintf(comm, fd, "      #   To get timing results run ./configure                #\n");CHKERRQ(ierr);
1335   ierr = PetscFPrintf(comm, fd, "      #   using --with-debugging=no, the performance will      #\n");CHKERRQ(ierr);
1336   ierr = PetscFPrintf(comm, fd, "      #   be generally two or three times faster.              #\n");CHKERRQ(ierr);
1337   ierr = PetscFPrintf(comm, fd, "      #                                                        #\n");CHKERRQ(ierr);
1338   ierr = PetscFPrintf(comm, fd, "      ##########################################################\n\n\n");CHKERRQ(ierr);
1339   PetscFunctionReturn(0);
1340 #else
1341   return 0;
1342 #endif
1343 }
1344 
1345 #if defined(PETSC_HAVE_OPENMP)
1346 extern PetscInt PetscNumOMPThreads;
1347 #endif
1348 
1349 PetscErrorCode  PetscLogView_Default(PetscViewer viewer)
1350 {
1351   FILE               *fd;
1352   PetscLogDouble     zero       = 0.0;
1353   PetscStageLog      stageLog;
1354   PetscStageInfo     *stageInfo = NULL;
1355   PetscEventPerfInfo *eventInfo = NULL;
1356   PetscClassPerfInfo *classInfo;
1357   char               arch[128],hostname[128],username[128],pname[PETSC_MAX_PATH_LEN],date[128];
1358   const char         *name;
1359   PetscLogDouble     locTotalTime, TotalTime, TotalFlops;
1360   PetscLogDouble     numMessages, messageLength, avgMessLen, numReductions;
1361   PetscLogDouble     stageTime, flops, flopr, mem, mess, messLen, red;
1362   PetscLogDouble     fracTime, fracFlops, fracMessages, fracLength, fracReductions, fracMess, fracMessLen, fracRed;
1363   PetscLogDouble     fracStageTime, fracStageFlops, fracStageMess, fracStageMessLen, fracStageRed;
1364   PetscLogDouble     min, max, tot, ratio, avg, x, y;
1365   PetscLogDouble     minf, maxf, totf, ratf, mint, maxt, tott, ratt, ratC, totm, totml, totr, mal, malmax, emalmax;
1366   #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1367   PetscLogDouble     cct, gct, csz, gsz;
1368   #endif
1369   PetscMPIInt        minC, maxC;
1370   PetscMPIInt        size, rank;
1371   PetscBool          *localStageUsed,    *stageUsed;
1372   PetscBool          *localStageVisible, *stageVisible;
1373   int                numStages, localNumEvents, numEvents;
1374   int                stage, oclass;
1375   PetscLogEvent      event;
1376   PetscErrorCode     ierr;
1377   char               version[256];
1378   MPI_Comm           comm;
1379 
1380   PetscFunctionBegin;
1381   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
1382   ierr = PetscViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr);
1383   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1384   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1385   /* Get the total elapsed time */
1386   PetscTime(&locTotalTime);  locTotalTime -= petsc_BaseTime;
1387 
1388   ierr = PetscFPrintf(comm, fd, "************************************************************************************************************************\n");CHKERRQ(ierr);
1389   ierr = PetscFPrintf(comm, fd, "***             WIDEN YOUR WINDOW TO 120 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document            ***\n");CHKERRQ(ierr);
1390   ierr = PetscFPrintf(comm, fd, "************************************************************************************************************************\n");CHKERRQ(ierr);
1391   ierr = PetscFPrintf(comm, fd, "\n---------------------------------------------- PETSc Performance Summary: ----------------------------------------------\n\n");CHKERRQ(ierr);
1392   ierr = PetscLogViewWarnSync(comm,fd);CHKERRQ(ierr);
1393   ierr = PetscLogViewWarnDebugging(comm,fd);CHKERRQ(ierr);
1394   ierr = PetscGetArchType(arch,sizeof(arch));CHKERRQ(ierr);
1395   ierr = PetscGetHostName(hostname,sizeof(hostname));CHKERRQ(ierr);
1396   ierr = PetscGetUserName(username,sizeof(username));CHKERRQ(ierr);
1397   ierr = PetscGetProgramName(pname,sizeof(pname));CHKERRQ(ierr);
1398   ierr = PetscGetDate(date,sizeof(date));CHKERRQ(ierr);
1399   ierr = PetscGetVersion(version,sizeof(version));CHKERRQ(ierr);
1400   if (size == 1) {
1401     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);
1402   } else {
1403     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);
1404   }
1405 #if defined(PETSC_HAVE_OPENMP)
1406   ierr = PetscFPrintf(comm,fd,"Using %D OpenMP threads\n", PetscNumOMPThreads);CHKERRQ(ierr);
1407 #endif
1408   ierr = PetscFPrintf(comm, fd, "Using %s\n", version);CHKERRQ(ierr);
1409 
1410   /* Must preserve reduction count before we go on */
1411   red = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
1412 
1413   /* Calculate summary information */
1414   ierr = PetscFPrintf(comm, fd, "\n                         Max       Max/Min     Avg       Total \n");CHKERRQ(ierr);
1415   /*   Time */
1416   ierr = MPIU_Allreduce(&locTotalTime, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1417   ierr = MPIU_Allreduce(&locTotalTime, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1418   ierr = MPIU_Allreduce(&locTotalTime, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1419   avg  = tot/((PetscLogDouble) size);
1420   if (min != 0.0) ratio = max/min; else ratio = 0.0;
1421   ierr = PetscFPrintf(comm, fd, "Time (sec):           %5.3e   %7.3f   %5.3e\n", max, ratio, avg);CHKERRQ(ierr);
1422   TotalTime = tot;
1423   /*   Objects */
1424   avg  = (PetscLogDouble) petsc_numObjects;
1425   ierr = MPIU_Allreduce(&avg,          &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1426   ierr = MPIU_Allreduce(&avg,          &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1427   ierr = MPIU_Allreduce(&avg,          &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1428   avg  = tot/((PetscLogDouble) size);
1429   if (min != 0.0) ratio = max/min; else ratio = 0.0;
1430   ierr = PetscFPrintf(comm, fd, "Objects:              %5.3e   %7.3f   %5.3e\n", max, ratio, avg);CHKERRQ(ierr);
1431   /*   Flops */
1432   ierr = MPIU_Allreduce(&petsc_TotalFlops,  &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1433   ierr = MPIU_Allreduce(&petsc_TotalFlops,  &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1434   ierr = MPIU_Allreduce(&petsc_TotalFlops,  &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1435   avg  = tot/((PetscLogDouble) size);
1436   if (min != 0.0) ratio = max/min; else ratio = 0.0;
1437   ierr = PetscFPrintf(comm, fd, "Flop:                 %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot);CHKERRQ(ierr);
1438   TotalFlops = tot;
1439   /*   Flops/sec -- Must talk to Barry here */
1440   if (locTotalTime != 0.0) flops = petsc_TotalFlops/locTotalTime; else flops = 0.0;
1441   ierr = MPIU_Allreduce(&flops,        &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1442   ierr = MPIU_Allreduce(&flops,        &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1443   ierr = MPIU_Allreduce(&flops,        &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1444   avg  = tot/((PetscLogDouble) size);
1445   if (min != 0.0) ratio = max/min; else ratio = 0.0;
1446   ierr = PetscFPrintf(comm, fd, "Flop/sec:             %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot);CHKERRQ(ierr);
1447   /*   Memory */
1448   ierr = PetscMallocGetMaximumUsage(&mem);CHKERRQ(ierr);
1449   if (mem > 0.0) {
1450     ierr = MPIU_Allreduce(&mem,          &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1451     ierr = MPIU_Allreduce(&mem,          &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1452     ierr = MPIU_Allreduce(&mem,          &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1453     avg  = tot/((PetscLogDouble) size);
1454     if (min != 0.0) ratio = max/min; else ratio = 0.0;
1455     ierr = PetscFPrintf(comm, fd, "Memory:               %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot);CHKERRQ(ierr);
1456   }
1457   /*   Messages */
1458   mess = 0.5*(petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct);
1459   ierr = MPIU_Allreduce(&mess,         &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1460   ierr = MPIU_Allreduce(&mess,         &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1461   ierr = MPIU_Allreduce(&mess,         &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1462   avg  = tot/((PetscLogDouble) size);
1463   if (min != 0.0) ratio = max/min; else ratio = 0.0;
1464   ierr = PetscFPrintf(comm, fd, "MPI Messages:         %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot);CHKERRQ(ierr);
1465   numMessages = tot;
1466   /*   Message Lengths */
1467   mess = 0.5*(petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len);
1468   ierr = MPIU_Allreduce(&mess,         &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1469   ierr = MPIU_Allreduce(&mess,         &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1470   ierr = MPIU_Allreduce(&mess,         &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1471   if (numMessages != 0) avg = tot/numMessages; else avg = 0.0;
1472   if (min != 0.0) ratio = max/min; else ratio = 0.0;
1473   ierr = PetscFPrintf(comm, fd, "MPI Message Lengths:  %5.3e   %7.3f   %5.3e  %5.3e\n", max, ratio, avg, tot);CHKERRQ(ierr);
1474   messageLength = tot;
1475   /*   Reductions */
1476   ierr = MPIU_Allreduce(&red,          &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1477   ierr = MPIU_Allreduce(&red,          &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1478   ierr = MPIU_Allreduce(&red,          &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1479   if (min != 0.0) ratio = max/min; else ratio = 0.0;
1480   ierr = PetscFPrintf(comm, fd, "MPI Reductions:       %5.3e   %7.3f\n", max, ratio);CHKERRQ(ierr);
1481   numReductions = red; /* wrong because uses count from process zero */
1482   ierr = PetscFPrintf(comm, fd, "\nFlop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)\n");CHKERRQ(ierr);
1483   ierr = PetscFPrintf(comm, fd, "                            e.g., VecAXPY() for real vectors of length N --> 2N flop\n");CHKERRQ(ierr);
1484   ierr = PetscFPrintf(comm, fd, "                            and VecAXPY() for complex vectors of length N --> 8N flop\n");CHKERRQ(ierr);
1485 
1486   /* Get total number of stages --
1487        Currently, a single processor can register more stages than another, but stages must all be registered in order.
1488        We can removed this requirement if necessary by having a global stage numbering and indirection on the stage ID.
1489        This seems best accomplished by assoicating a communicator with each stage.
1490   */
1491   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1492   ierr = MPIU_Allreduce(&stageLog->numStages, &numStages, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1493   ierr = PetscMalloc1(numStages, &localStageUsed);CHKERRQ(ierr);
1494   ierr = PetscMalloc1(numStages, &stageUsed);CHKERRQ(ierr);
1495   ierr = PetscMalloc1(numStages, &localStageVisible);CHKERRQ(ierr);
1496   ierr = PetscMalloc1(numStages, &stageVisible);CHKERRQ(ierr);
1497   if (numStages > 0) {
1498     stageInfo = stageLog->stageInfo;
1499     for (stage = 0; stage < numStages; stage++) {
1500       if (stage < stageLog->numStages) {
1501         localStageUsed[stage]    = stageInfo[stage].used;
1502         localStageVisible[stage] = stageInfo[stage].perfInfo.visible;
1503       } else {
1504         localStageUsed[stage]    = PETSC_FALSE;
1505         localStageVisible[stage] = PETSC_TRUE;
1506       }
1507     }
1508     ierr = MPIU_Allreduce(localStageUsed,    stageUsed,    numStages, MPIU_BOOL, MPI_LOR,  comm);CHKERRQ(ierr);
1509     ierr = MPIU_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr);
1510     for (stage = 0; stage < numStages; stage++) {
1511       if (stageUsed[stage]) {
1512         ierr = PetscFPrintf(comm, fd, "\nSummary of Stages:   ----- Time ------  ----- Flop ------  --- Messages ---  -- Message Lengths --  -- Reductions --\n");CHKERRQ(ierr);
1513         ierr = PetscFPrintf(comm, fd, "                        Avg     %%Total     Avg     %%Total    Count   %%Total     Avg         %%Total    Count   %%Total \n");CHKERRQ(ierr);
1514         break;
1515       }
1516     }
1517     for (stage = 0; stage < numStages; stage++) {
1518       if (!stageUsed[stage]) continue;
1519       /* CANNOT use MPIU_Allreduce() since it might fail the line number check */
1520       if (localStageUsed[stage]) {
1521         ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.time,          &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1522         ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.flops,         &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1523         ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages,   &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1524         ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1525         ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1526         name = stageInfo[stage].name;
1527       } else {
1528         ierr = MPI_Allreduce(&zero,                           &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1529         ierr = MPI_Allreduce(&zero,                           &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1530         ierr = MPI_Allreduce(&zero,                           &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1531         ierr = MPI_Allreduce(&zero,                           &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1532         ierr = MPI_Allreduce(&zero,                           &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1533         name = "";
1534       }
1535       mess *= 0.5; messLen *= 0.5; red /= size;
1536       if (TotalTime     != 0.0) fracTime       = stageTime/TotalTime;    else fracTime       = 0.0;
1537       if (TotalFlops    != 0.0) fracFlops      = flops/TotalFlops;       else fracFlops      = 0.0;
1538       /* Talk to Barry if (stageTime     != 0.0) flops          = (size*flops)/stageTime; else flops          = 0.0; */
1539       if (numMessages   != 0.0) fracMessages   = mess/numMessages;       else fracMessages   = 0.0;
1540       if (mess          != 0.0) avgMessLen     = messLen/mess;           else avgMessLen     = 0.0;
1541       if (messageLength != 0.0) fracLength     = messLen/messageLength;  else fracLength     = 0.0;
1542       if (numReductions != 0.0) fracReductions = red/numReductions;      else fracReductions = 0.0;
1543       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",
1544                           stage, name, stageTime/size, 100.0*fracTime, flops, 100.0*fracFlops,
1545                           mess, 100.0*fracMessages, avgMessLen, 100.0*fracLength, red, 100.0*fracReductions);CHKERRQ(ierr);
1546     }
1547   }
1548 
1549   ierr = PetscFPrintf(comm, fd,"\n------------------------------------------------------------------------------------------------------------------------\n");CHKERRQ(ierr);
1550   ierr = PetscFPrintf(comm, fd, "See the 'Profiling' chapter of the users' manual for details on interpreting output.\n");CHKERRQ(ierr);
1551   ierr = PetscFPrintf(comm, fd, "Phase summary info:\n");CHKERRQ(ierr);
1552   ierr = PetscFPrintf(comm, fd, "   Count: number of times phase was executed\n");CHKERRQ(ierr);
1553   ierr = PetscFPrintf(comm, fd, "   Time and Flop: Max - maximum over all processors\n");CHKERRQ(ierr);
1554   ierr = PetscFPrintf(comm, fd, "                  Ratio - ratio of maximum to minimum over all processors\n");CHKERRQ(ierr);
1555   ierr = PetscFPrintf(comm, fd, "   Mess: number of messages sent\n");CHKERRQ(ierr);
1556   ierr = PetscFPrintf(comm, fd, "   AvgLen: average message length (bytes)\n");CHKERRQ(ierr);
1557   ierr = PetscFPrintf(comm, fd, "   Reduct: number of global reductions\n");CHKERRQ(ierr);
1558   ierr = PetscFPrintf(comm, fd, "   Global: entire computation\n");CHKERRQ(ierr);
1559   ierr = PetscFPrintf(comm, fd, "   Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().\n");CHKERRQ(ierr);
1560   ierr = PetscFPrintf(comm, fd, "      %%T - percent time in this phase         %%F - percent flop in this phase\n");CHKERRQ(ierr);
1561   ierr = PetscFPrintf(comm, fd, "      %%M - percent messages in this phase     %%L - percent message lengths in this phase\n");CHKERRQ(ierr);
1562   ierr = PetscFPrintf(comm, fd, "      %%R - percent reductions in this phase\n");CHKERRQ(ierr);
1563   ierr = PetscFPrintf(comm, fd, "   Total Mflop/s: 10e-6 * (sum of flop over all processors)/(max time over all processors)\n");CHKERRQ(ierr);
1564   if (PetscLogMemory) {
1565     ierr = PetscFPrintf(comm, fd, "   Malloc Mbytes: Memory allocated and kept during event (sum over all calls to event)\n");CHKERRQ(ierr);
1566     ierr = PetscFPrintf(comm, fd, "   EMalloc Mbytes: extra memory allocated during event and then freed (maximum over all calls to events)\n");CHKERRQ(ierr);
1567     ierr = PetscFPrintf(comm, fd, "   MMalloc Mbytes: Increase in high water mark of allocated memory (sum over all calls to event)\n");CHKERRQ(ierr);
1568     ierr = PetscFPrintf(comm, fd, "   RMI Mbytes: Increase in resident memory (sum over all calls to event)\n");CHKERRQ(ierr);
1569   }
1570   #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1571   ierr = PetscFPrintf(comm, fd, "   CpuToGpu Count: total number of CPU to GPU copies over all processors\n");CHKERRQ(ierr);
1572   ierr = PetscFPrintf(comm, fd, "   CpuToGpu Size (Mbytes): 10e-6 * (total size of CPU to GPU copies)\n");CHKERRQ(ierr);
1573   ierr = PetscFPrintf(comm, fd, "   GpuToCpu Count: total number of GPU to CPU copies over all processors\n");CHKERRQ(ierr);
1574   ierr = PetscFPrintf(comm, fd, "   GpuToCpu Size (Mbytes): 10e-6 * (total size of GPU to CPU copies)\n");CHKERRQ(ierr);
1575   #endif
1576   ierr = PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------\n");CHKERRQ(ierr);
1577 
1578   ierr = PetscLogViewWarnDebugging(comm,fd);CHKERRQ(ierr);
1579 
1580   /* Report events */
1581   ierr = PetscFPrintf(comm, fd,"Event                Count      Time (sec)     Flop                              --- Global ---  --- Stage ----  Total");CHKERRQ(ierr);
1582   if (PetscLogMemory) {
1583     ierr = PetscFPrintf(comm, fd,"  Malloc EMalloc MMalloc RMI");CHKERRQ(ierr);
1584   }
1585   #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1586   ierr = PetscFPrintf(comm, fd,"   - CpuToGpu -   - GpuToCpu -");CHKERRQ(ierr);
1587   #endif
1588   ierr = PetscFPrintf(comm, fd,"\n");CHKERRQ(ierr);
1589   ierr = PetscFPrintf(comm, fd,"                   Max Ratio  Max     Ratio   Max  Ratio  Mess   AvgLen  Reduct  %%T %%F %%M %%L %%R  %%T %%F %%M %%L %%R Mflop/s");CHKERRQ(ierr);
1590   if (PetscLogMemory) {
1591     ierr = PetscFPrintf(comm, fd," Mbytes Mbytes Mbytes Mbytes");CHKERRQ(ierr);
1592   }
1593   #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1594   ierr = PetscFPrintf(comm, fd,"  Count   Size   Count   Size");CHKERRQ(ierr);
1595   #endif
1596   ierr = PetscFPrintf(comm, fd,"\n");CHKERRQ(ierr);
1597   ierr = PetscFPrintf(comm, fd,"------------------------------------------------------------------------------------------------------------------------\n");CHKERRQ(ierr);
1598 
1599   /* Problem: The stage name will not show up unless the stage executed on proc 1 */
1600   for (stage = 0; stage < numStages; stage++) {
1601     if (!stageVisible[stage]) continue;
1602     /* CANNOT use MPIU_Allreduce() since it might fail the line number check */
1603     if (localStageUsed[stage]) {
1604       ierr = PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stageInfo[stage].name);CHKERRQ(ierr);
1605       ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.time,          &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1606       ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.flops,         &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1607       ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages,   &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1608       ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1609       ierr = MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1610     } else {
1611       ierr = PetscFPrintf(comm, fd, "\n--- Event Stage %d: Unknown\n\n", stage);CHKERRQ(ierr);
1612       ierr = MPI_Allreduce(&zero,                           &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1613       ierr = MPI_Allreduce(&zero,                           &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1614       ierr = MPI_Allreduce(&zero,                           &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1615       ierr = MPI_Allreduce(&zero,                           &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1616       ierr = MPI_Allreduce(&zero,                           &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1617     }
1618     mess *= 0.5; messLen *= 0.5; red /= size;
1619 
1620     /* Get total number of events in this stage --
1621        Currently, a single processor can register more events than another, but events must all be registered in order,
1622        just like stages. We can removed this requirement if necessary by having a global event numbering and indirection
1623        on the event ID. This seems best accomplished by associating a communicator with each stage.
1624 
1625        Problem: If the event did not happen on proc 1, its name will not be available.
1626        Problem: Event visibility is not implemented
1627     */
1628     if (localStageUsed[stage]) {
1629       eventInfo      = stageLog->stageInfo[stage].eventLog->eventInfo;
1630       localNumEvents = stageLog->stageInfo[stage].eventLog->numEvents;
1631     } else localNumEvents = 0;
1632     ierr = MPIU_Allreduce(&localNumEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
1633     for (event = 0; event < numEvents; event++) {
1634       /* CANNOT use MPIU_Allreduce() since it might fail the line number check */
1635       if (localStageUsed[stage] && (event < stageLog->stageInfo[stage].eventLog->numEvents) && (eventInfo[event].depth == 0)) {
1636         if ((eventInfo[event].count > 0) && (eventInfo[event].time > 0.0)) flopr = eventInfo[event].flops; else flopr = 0.0;
1637         ierr = MPI_Allreduce(&flopr,                          &minf,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1638         ierr = MPI_Allreduce(&flopr,                          &maxf,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1639         ierr = MPI_Allreduce(&eventInfo[event].flops,         &totf,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1640         ierr = MPI_Allreduce(&eventInfo[event].time,          &mint,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1641         ierr = MPI_Allreduce(&eventInfo[event].time,          &maxt,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1642         ierr = MPI_Allreduce(&eventInfo[event].time,          &tott,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1643         ierr = MPI_Allreduce(&eventInfo[event].numMessages,   &totm,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1644         ierr = MPI_Allreduce(&eventInfo[event].messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1645         ierr = MPI_Allreduce(&eventInfo[event].numReductions, &totr,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1646         ierr = MPI_Allreduce(&eventInfo[event].count,         &minC,  1, MPI_INT,             MPI_MIN, comm);CHKERRQ(ierr);
1647         ierr = MPI_Allreduce(&eventInfo[event].count,         &maxC,  1, MPI_INT,             MPI_MAX, comm);CHKERRQ(ierr);
1648         if (PetscLogMemory) {
1649           ierr  = MPI_Allreduce(&eventInfo[event].memIncrease,    &mem,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1650           ierr  = MPI_Allreduce(&eventInfo[event].mallocSpace,    &mal,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1651           ierr  = MPI_Allreduce(&eventInfo[event].mallocIncrease, &malmax,1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1652           ierr  = MPI_Allreduce(&eventInfo[event].mallocIncreaseEvent, &emalmax,1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1653         }
1654         #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1655         ierr  = MPI_Allreduce(&eventInfo[event].CpuToGpuCount,    &cct,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1656         ierr  = MPI_Allreduce(&eventInfo[event].GpuToCpuCount,    &gct,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1657         ierr  = MPI_Allreduce(&eventInfo[event].CpuToGpuSize,     &csz,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1658         ierr  = MPI_Allreduce(&eventInfo[event].GpuToCpuSize,     &gsz,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1659         #endif
1660         name = stageLog->eventLog->eventInfo[event].name;
1661       } else {
1662         flopr = 0.0;
1663         ierr  = MPI_Allreduce(&flopr,                         &minf,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1664         ierr  = MPI_Allreduce(&flopr,                         &maxf,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1665         ierr  = MPI_Allreduce(&zero,                          &totf,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1666         ierr  = MPI_Allreduce(&zero,                          &mint,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);CHKERRQ(ierr);
1667         ierr  = MPI_Allreduce(&zero,                          &maxt,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);CHKERRQ(ierr);
1668         ierr  = MPI_Allreduce(&zero,                          &tott,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1669         ierr  = MPI_Allreduce(&zero,                          &totm,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1670         ierr  = MPI_Allreduce(&zero,                          &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1671         ierr  = MPI_Allreduce(&zero,                          &totr,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1672         ierr  = MPI_Allreduce(&ierr,                          &minC,  1, MPI_INT,             MPI_MIN, comm);CHKERRQ(ierr);
1673         ierr  = MPI_Allreduce(&ierr,                          &maxC,  1, MPI_INT,             MPI_MAX, comm);CHKERRQ(ierr);
1674         if (PetscLogMemory) {
1675           ierr  = MPI_Allreduce(&zero,                        &mem,    1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1676           ierr  = MPI_Allreduce(&zero,                        &mal,    1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1677           ierr  = MPI_Allreduce(&zero,                        &malmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1678           ierr  = MPI_Allreduce(&zero,                        &emalmax,1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1679         }
1680         #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1681         ierr  = MPI_Allreduce(&zero,                          &cct,    1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1682         ierr  = MPI_Allreduce(&zero,                          &gct,    1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1683         ierr  = MPI_Allreduce(&zero,                          &csz,    1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1684         ierr  = MPI_Allreduce(&zero,                          &gsz,    1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);CHKERRQ(ierr);
1685         #endif
1686         name  = "";
1687       }
1688       if (mint < 0.0) {
1689         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);
1690         mint = 0;
1691       }
1692       if (minf < 0.0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Minimum flop %g over all processors for %s is negative! Not possible!",minf,name);
1693       totm *= 0.5; totml *= 0.5; totr /= size;
1694 
1695       if (maxC != 0) {
1696         if (minC          != 0)   ratC             = ((PetscLogDouble)maxC)/minC;else ratC             = 0.0;
1697         if (mint          != 0.0) ratt             = maxt/mint;                  else ratt             = 0.0;
1698         if (minf          != 0.0) ratf             = maxf/minf;                  else ratf             = 0.0;
1699         if (TotalTime     != 0.0) fracTime         = tott/TotalTime;             else fracTime         = 0.0;
1700         if (TotalFlops    != 0.0) fracFlops        = totf/TotalFlops;            else fracFlops        = 0.0;
1701         if (stageTime     != 0.0) fracStageTime    = tott/stageTime;             else fracStageTime    = 0.0;
1702         if (flops         != 0.0) fracStageFlops   = totf/flops;                 else fracStageFlops   = 0.0;
1703         if (numMessages   != 0.0) fracMess         = totm/numMessages;           else fracMess         = 0.0;
1704         if (messageLength != 0.0) fracMessLen      = totml/messageLength;        else fracMessLen      = 0.0;
1705         if (numReductions != 0.0) fracRed          = totr/numReductions;         else fracRed          = 0.0;
1706         if (mess          != 0.0) fracStageMess    = totm/mess;                  else fracStageMess    = 0.0;
1707         if (messLen       != 0.0) fracStageMessLen = totml/messLen;              else fracStageMessLen = 0.0;
1708         if (red           != 0.0) fracStageRed     = totr/red;                   else fracStageRed     = 0.0;
1709         if (totm          != 0.0) totml           /= totm;                       else totml            = 0.0;
1710         if (maxt          != 0.0) flopr            = totf/maxt;                  else flopr            = 0.0;
1711         if (fracStageTime > 1.00)  ierr = PetscFPrintf(comm, fd,"Warning -- total time of event greater than time of entire stage -- something is wrong with the timer\n");CHKERRQ(ierr);
1712         ierr = PetscFPrintf(comm, fd,
1713                             "%-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",
1714                             name, maxC, ratC, maxt, ratt, maxf, ratf, totm, totml, totr,
1715                             100.0*fracTime, 100.0*fracFlops, 100.0*fracMess, 100.0*fracMessLen, 100.0*fracRed,
1716                             100.0*fracStageTime, 100.0*fracStageFlops, 100.0*fracStageMess, 100.0*fracStageMessLen, 100.0*fracStageRed,
1717                             PetscAbs(flopr)/1.0e6);CHKERRQ(ierr);
1718         if (PetscLogMemory) {
1719           ierr = PetscFPrintf(comm, fd," %5.0f   %5.0f   %5.0f   %5.0f",mal/1.0e6,emalmax/1.0e6,malmax/1.0e6,mem/1.0e6);CHKERRQ(ierr);
1720         }
1721         #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1722         ierr = PetscFPrintf(comm, fd,"    %4.0f %3.2e %4.0f %3.2e", cct, csz/1.0e6, gct, gsz/1.0e6);CHKERRQ(ierr);
1723         #endif
1724         ierr = PetscFPrintf(comm, fd,"\n");CHKERRQ(ierr);
1725       }
1726     }
1727   }
1728 
1729   /* Memory usage and object creation */
1730   ierr = PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------\n");CHKERRQ(ierr);
1731   ierr = PetscFPrintf(comm, fd, "\n");CHKERRQ(ierr);
1732   ierr = PetscFPrintf(comm, fd, "Memory usage is given in bytes:\n\n");CHKERRQ(ierr);
1733 
1734   /* Right now, only stages on the first processor are reported here, meaning only objects associated with
1735      the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then
1736      stats for stages local to processor sets.
1737   */
1738   /* We should figure out the longest object name here (now 20 characters) */
1739   ierr = PetscFPrintf(comm, fd, "Object Type          Creations   Destructions     Memory  Descendants' Mem.\n");CHKERRQ(ierr);
1740   ierr = PetscFPrintf(comm, fd, "Reports information only for process 0.\n");CHKERRQ(ierr);
1741   for (stage = 0; stage < numStages; stage++) {
1742     if (localStageUsed[stage]) {
1743       classInfo = stageLog->stageInfo[stage].classLog->classInfo;
1744       ierr = PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stageInfo[stage].name);CHKERRQ(ierr);
1745       for (oclass = 0; oclass < stageLog->stageInfo[stage].classLog->numClasses; oclass++) {
1746         if ((classInfo[oclass].creations > 0) || (classInfo[oclass].destructions > 0)) {
1747           ierr = PetscFPrintf(comm, fd, "%20s %5d          %5d  %11.0f     %g\n", stageLog->classLog->classInfo[oclass].name,
1748                               classInfo[oclass].creations, classInfo[oclass].destructions, classInfo[oclass].mem,
1749                               classInfo[oclass].descMem);CHKERRQ(ierr);
1750         }
1751       }
1752     } else {
1753       if (!localStageVisible[stage]) continue;
1754       ierr = PetscFPrintf(comm, fd, "\n--- Event Stage %d: Unknown\n\n", stage);CHKERRQ(ierr);
1755     }
1756   }
1757 
1758   ierr = PetscFree(localStageUsed);CHKERRQ(ierr);
1759   ierr = PetscFree(stageUsed);CHKERRQ(ierr);
1760   ierr = PetscFree(localStageVisible);CHKERRQ(ierr);
1761   ierr = PetscFree(stageVisible);CHKERRQ(ierr);
1762 
1763   /* Information unrelated to this particular run */
1764   ierr = PetscFPrintf(comm, fd, "========================================================================================================================\n");CHKERRQ(ierr);
1765   PetscTime(&y);
1766   PetscTime(&x);
1767   PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y);
1768   PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y);
1769   ierr = PetscFPrintf(comm,fd,"Average time to get PetscTime(): %g\n", (y-x)/10.0);CHKERRQ(ierr);
1770   /* MPI information */
1771   if (size > 1) {
1772     MPI_Status  status;
1773     PetscMPIInt tag;
1774     MPI_Comm    newcomm;
1775 
1776     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1777     PetscTime(&x);
1778     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1779     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1780     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1781     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1782     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1783     PetscTime(&y);
1784     ierr = PetscFPrintf(comm, fd, "Average time for MPI_Barrier(): %g\n", (y-x)/5.0);CHKERRQ(ierr);
1785     ierr = PetscCommDuplicate(comm,&newcomm, &tag);CHKERRQ(ierr);
1786     ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1787     if (rank) {
1788       ierr = MPI_Recv(NULL, 0, MPI_INT, rank-1,            tag, newcomm, &status);CHKERRQ(ierr);
1789       ierr = MPI_Send(NULL, 0, MPI_INT, (rank+1)%size, tag, newcomm);CHKERRQ(ierr);
1790     } else {
1791       PetscTime(&x);
1792       ierr = MPI_Send(NULL, 0, MPI_INT, 1,          tag, newcomm);CHKERRQ(ierr);
1793       ierr = MPI_Recv(NULL, 0, MPI_INT, size-1, tag, newcomm, &status);CHKERRQ(ierr);
1794       PetscTime(&y);
1795       ierr = PetscFPrintf(comm,fd,"Average time for zero size MPI_Send(): %g\n", (y-x)/size);CHKERRQ(ierr);
1796     }
1797     ierr = PetscCommDestroy(&newcomm);CHKERRQ(ierr);
1798   }
1799   ierr = PetscOptionsView(NULL,viewer);CHKERRQ(ierr);
1800 
1801   /* Machine and compile information */
1802 #if defined(PETSC_USE_FORTRAN_KERNELS)
1803   ierr = PetscFPrintf(comm, fd, "Compiled with FORTRAN kernels\n");CHKERRQ(ierr);
1804 #else
1805   ierr = PetscFPrintf(comm, fd, "Compiled without FORTRAN kernels\n");CHKERRQ(ierr);
1806 #endif
1807 #if defined(PETSC_USE_64BIT_INDICES)
1808   ierr = PetscFPrintf(comm, fd, "Compiled with 64 bit PetscInt\n");CHKERRQ(ierr);
1809 #elif defined(PETSC_USE___FLOAT128)
1810   ierr = PetscFPrintf(comm, fd, "Compiled with 32 bit PetscInt\n");CHKERRQ(ierr);
1811 #endif
1812 #if defined(PETSC_USE_REAL_SINGLE)
1813   ierr = PetscFPrintf(comm, fd, "Compiled with single precision PetscScalar and PetscReal\n");CHKERRQ(ierr);
1814 #elif defined(PETSC_USE___FLOAT128)
1815   ierr = PetscFPrintf(comm, fd, "Compiled with 128 bit precision PetscScalar and PetscReal\n");CHKERRQ(ierr);
1816 #endif
1817 #if defined(PETSC_USE_REAL_MAT_SINGLE)
1818   ierr = PetscFPrintf(comm, fd, "Compiled with single precision matrices\n");CHKERRQ(ierr);
1819 #else
1820   ierr = PetscFPrintf(comm, fd, "Compiled with full precision matrices (default)\n");CHKERRQ(ierr);
1821 #endif
1822   ierr = PetscFPrintf(comm, fd, "sizeof(short) %d sizeof(int) %d sizeof(long) %d sizeof(void*) %d sizeof(PetscScalar) %d sizeof(PetscInt) %d\n",
1823                       (int) sizeof(short), (int) sizeof(int), (int) sizeof(long), (int) sizeof(void*),(int) sizeof(PetscScalar),(int) sizeof(PetscInt));CHKERRQ(ierr);
1824 
1825   ierr = PetscFPrintf(comm, fd, "Configure options: %s",petscconfigureoptions);CHKERRQ(ierr);
1826   ierr = PetscFPrintf(comm, fd, "%s", petscmachineinfo);CHKERRQ(ierr);
1827   ierr = PetscFPrintf(comm, fd, "%s", petsccompilerinfo);CHKERRQ(ierr);
1828   ierr = PetscFPrintf(comm, fd, "%s", petsccompilerflagsinfo);CHKERRQ(ierr);
1829   ierr = PetscFPrintf(comm, fd, "%s", petsclinkerinfo);CHKERRQ(ierr);
1830 
1831   /* Cleanup */
1832   ierr = PetscFPrintf(comm, fd, "\n");CHKERRQ(ierr);
1833   ierr = PetscLogViewWarnDebugging(comm,fd);CHKERRQ(ierr);
1834   PetscFunctionReturn(0);
1835 }
1836 
1837 /*@C
1838   PetscLogView - Prints a summary of the logging.
1839 
1840   Collective over MPI_Comm
1841 
1842   Input Parameter:
1843 .  viewer - an ASCII viewer
1844 
1845   Options Database Keys:
1846 +  -log_view [:filename] - Prints summary of log information
1847 .  -log_view :filename.py:ascii_info_detail - Saves logging information from each process as a Python file
1848 .  -log_view :filename.xml:ascii_xml - Saves a summary of the logging information in a nested format (see below for how to view it)
1849 .  -log_all - Saves a file Log.rank for each MPI process with details of each step of the computation
1850 -  -log_trace [filename] - Displays a trace of what each process is doing
1851 
1852   Notes:
1853   It is possible to control the logging programatically but we recommend using the options database approach whenever possible
1854   By default the summary is printed to stdout.
1855 
1856   Before calling this routine you must have called either PetscLogDefaultBegin() or PetscLogNestedBegin()
1857 
1858   If PETSc is configured with --with-logging=0 then this functionality is not available
1859 
1860   To view the nested XML format filename.xml first copy  ${PETSC_DIR}/share/petsc/xml/performance_xml2html.xsl to the current
1861   directory then open filename.xml with your browser. Specific notes for certain browsers
1862 $    Firefox and Internet explorer - simply open the file
1863 $    Google Chrome - you must start up Chrome with the option --allow-file-access-from-files
1864 $    Safari - see https://ccm.net/faq/36342-safari-how-to-enable-local-file-access
1865   or one can use the package http://xmlsoft.org/XSLT/xsltproc2.html to translate the xml file to html and then open it with
1866   your browser.
1867   Alternatively, use the script ${PETSC_DIR}/lib/petsc/bin/petsc-performance-view to automatically open a new browser
1868   window and render the XML log file contents.
1869 
1870   The nested XML format was kindly donated by Koos Huijssen and Christiaan M. Klaij  MARITIME  RESEARCH  INSTITUTE  NETHERLANDS
1871 
1872   Level: beginner
1873 
1874 .seealso: PetscLogDefaultBegin(), PetscLogDump()
1875 @*/
1876 PetscErrorCode  PetscLogView(PetscViewer viewer)
1877 {
1878   PetscErrorCode    ierr;
1879   PetscBool         isascii;
1880   PetscViewerFormat format;
1881   int               stage, lastStage;
1882   PetscStageLog     stageLog;
1883 
1884   PetscFunctionBegin;
1885   if (!PetscLogPLB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must use -log_view or PetscLogDefaultBegin() before calling this routine");
1886   /* Pop off any stages the user forgot to remove */
1887   lastStage = 0;
1888   ierr      = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
1889   ierr      = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
1890   while (stage >= 0) {
1891     lastStage = stage;
1892     ierr      = PetscStageLogPop(stageLog);CHKERRQ(ierr);
1893     ierr      = PetscStageLogGetCurrent(stageLog, &stage);CHKERRQ(ierr);
1894   }
1895   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1896   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Currently can only view logging to ASCII");
1897   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1898   if (format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_ASCII_INFO) {
1899     ierr = PetscLogView_Default(viewer);CHKERRQ(ierr);
1900   } else if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1901     ierr = PetscLogView_Detailed(viewer);CHKERRQ(ierr);
1902   } else if (format == PETSC_VIEWER_ASCII_CSV) {
1903     ierr = PetscLogView_CSV(viewer);CHKERRQ(ierr);
1904   } else if (format == PETSC_VIEWER_ASCII_XML) {
1905     ierr = PetscLogView_Nested(viewer);CHKERRQ(ierr);
1906   }
1907   ierr = PetscStageLogPush(stageLog, lastStage);CHKERRQ(ierr);
1908   PetscFunctionReturn(0);
1909 }
1910 
1911 /*@C
1912   PetscLogViewFromOptions - Processes command line options to determine if/how a PetscLog is to be viewed.
1913 
1914   Collective on PETSC_COMM_WORLD
1915 
1916   Not normally called by user
1917 
1918   Level: intermediate
1919 
1920 @*/
1921 PetscErrorCode PetscLogViewFromOptions(void)
1922 {
1923   PetscErrorCode    ierr;
1924   PetscViewer       viewer;
1925   PetscBool         flg;
1926   PetscViewerFormat format;
1927 
1928   PetscFunctionBegin;
1929   ierr   = PetscOptionsGetViewer(PETSC_COMM_WORLD,NULL,NULL,"-log_view",&viewer,&format,&flg);CHKERRQ(ierr);
1930   if (flg) {
1931     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1932     ierr = PetscLogView(viewer);CHKERRQ(ierr);
1933     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1934     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1935   }
1936   PetscFunctionReturn(0);
1937 }
1938 
1939 
1940 
1941 /*----------------------------------------------- Counter Functions -------------------------------------------------*/
1942 /*@C
1943    PetscGetFlops - Returns the number of flops used on this processor
1944    since the program began.
1945 
1946    Not Collective
1947 
1948    Output Parameter:
1949    flops - number of floating point operations
1950 
1951    Notes:
1952    A global counter logs all PETSc flop counts.  The user can use
1953    PetscLogFlops() to increment this counter to include flops for the
1954    application code.
1955 
1956    Level: intermediate
1957 
1958 .seealso: PetscTime(), PetscLogFlops()
1959 @*/
1960 PetscErrorCode  PetscGetFlops(PetscLogDouble *flops)
1961 {
1962   PetscFunctionBegin;
1963   *flops = petsc_TotalFlops;
1964   PetscFunctionReturn(0);
1965 }
1966 
1967 PetscErrorCode  PetscLogObjectState(PetscObject obj, const char format[], ...)
1968 {
1969   PetscErrorCode ierr;
1970   size_t         fullLength;
1971   va_list        Argp;
1972 
1973   PetscFunctionBegin;
1974   if (!petsc_logObjects) PetscFunctionReturn(0);
1975   va_start(Argp, format);
1976   ierr = PetscVSNPrintf(petsc_objects[obj->id].info, 64,format,&fullLength, Argp);CHKERRQ(ierr);
1977   va_end(Argp);
1978   PetscFunctionReturn(0);
1979 }
1980 
1981 
1982 /*MC
1983    PetscLogFlops - Adds floating point operations to the global counter.
1984 
1985    Synopsis:
1986    #include <petsclog.h>
1987    PetscErrorCode PetscLogFlops(PetscLogDouble f)
1988 
1989    Not Collective
1990 
1991    Input Parameter:
1992 .  f - flop counter
1993 
1994 
1995    Usage:
1996 .vb
1997      PetscLogEvent USER_EVENT;
1998      PetscLogEventRegister("User event",0,&USER_EVENT);
1999      PetscLogEventBegin(USER_EVENT,0,0,0,0);
2000         [code segment to monitor]
2001         PetscLogFlops(user_flops)
2002      PetscLogEventEnd(USER_EVENT,0,0,0,0);
2003 .ve
2004 
2005    Notes:
2006    A global counter logs all PETSc flop counts.  The user can use
2007    PetscLogFlops() to increment this counter to include flops for the
2008    application code.
2009 
2010    Level: intermediate
2011 
2012 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscGetFlops()
2013 
2014 M*/
2015 
2016 /*MC
2017    PetscPreLoadBegin - Begin a segment of code that may be preloaded (run twice)
2018     to get accurate timings
2019 
2020    Synopsis:
2021    #include <petsclog.h>
2022    void PetscPreLoadBegin(PetscBool  flag,char *name);
2023 
2024    Not Collective
2025 
2026    Input Parameter:
2027 +   flag - PETSC_TRUE to run twice, PETSC_FALSE to run once, may be overridden
2028            with command line option -preload true or -preload false
2029 -   name - name of first stage (lines of code timed separately with -log_view) to
2030            be preloaded
2031 
2032    Usage:
2033 .vb
2034      PetscPreLoadBegin(PETSC_TRUE,"first stage);
2035        lines of code
2036        PetscPreLoadStage("second stage");
2037        lines of code
2038      PetscPreLoadEnd();
2039 .ve
2040 
2041    Notes:
2042     Only works in C/C++, not Fortran
2043 
2044      Flags available within the macro.
2045 +    PetscPreLoadingUsed - true if we are or have done preloading
2046 .    PetscPreLoadingOn - true if it is CURRENTLY doing preload
2047 .    PetscPreLoadIt - 0 for the first computation (with preloading turned off it is only 0) 1 for the second
2048 -    PetscPreLoadMax - number of times it will do the computation, only one when preloading is turned on
2049      The first two variables are available throughout the program, the second two only between the PetscPreLoadBegin()
2050      and PetscPreLoadEnd()
2051 
2052    Level: intermediate
2053 
2054 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadEnd(), PetscPreLoadStage()
2055 
2056 
2057 M*/
2058 
2059 /*MC
2060    PetscPreLoadEnd - End a segment of code that may be preloaded (run twice)
2061     to get accurate timings
2062 
2063    Synopsis:
2064    #include <petsclog.h>
2065    void PetscPreLoadEnd(void);
2066 
2067    Not Collective
2068 
2069    Usage:
2070 .vb
2071      PetscPreLoadBegin(PETSC_TRUE,"first stage);
2072        lines of code
2073        PetscPreLoadStage("second stage");
2074        lines of code
2075      PetscPreLoadEnd();
2076 .ve
2077 
2078    Notes:
2079     only works in C/C++ not fortran
2080 
2081    Level: intermediate
2082 
2083 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadBegin(), PetscPreLoadStage()
2084 
2085 M*/
2086 
2087 /*MC
2088    PetscPreLoadStage - Start a new segment of code to be timed separately.
2089     to get accurate timings
2090 
2091    Synopsis:
2092    #include <petsclog.h>
2093    void PetscPreLoadStage(char *name);
2094 
2095    Not Collective
2096 
2097    Usage:
2098 .vb
2099      PetscPreLoadBegin(PETSC_TRUE,"first stage);
2100        lines of code
2101        PetscPreLoadStage("second stage");
2102        lines of code
2103      PetscPreLoadEnd();
2104 .ve
2105 
2106    Notes:
2107     only works in C/C++ not fortran
2108 
2109    Level: intermediate
2110 
2111 .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadBegin(), PetscPreLoadEnd()
2112 
2113 M*/
2114 
2115 
2116 #else /* end of -DPETSC_USE_LOG section */
2117 
2118 PetscErrorCode  PetscLogObjectState(PetscObject obj, const char format[], ...)
2119 {
2120   PetscFunctionBegin;
2121   PetscFunctionReturn(0);
2122 }
2123 
2124 #endif /* PETSC_USE_LOG*/
2125 
2126 
2127 PetscClassId PETSC_LARGEST_CLASSID = PETSC_SMALLEST_CLASSID;
2128 PetscClassId PETSC_OBJECT_CLASSID  = 0;
2129 
2130 /*@C
2131   PetscClassIdRegister - Registers a new class name for objects and logging operations in an application code.
2132 
2133   Not Collective
2134 
2135   Input Parameter:
2136 . name   - The class name
2137 
2138   Output Parameter:
2139 . oclass - The class id or classid
2140 
2141   Level: developer
2142 
2143 @*/
2144 PetscErrorCode  PetscClassIdRegister(const char name[],PetscClassId *oclass)
2145 {
2146 #if defined(PETSC_USE_LOG)
2147   PetscStageLog  stageLog;
2148   PetscInt       stage;
2149   PetscErrorCode ierr;
2150 #endif
2151 
2152   PetscFunctionBegin;
2153   *oclass = ++PETSC_LARGEST_CLASSID;
2154 #if defined(PETSC_USE_LOG)
2155   ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
2156   ierr = PetscClassRegLogRegister(stageLog->classLog, name, *oclass);CHKERRQ(ierr);
2157   for (stage = 0; stage < stageLog->numStages; stage++) {
2158     ierr = PetscClassPerfLogEnsureSize(stageLog->stageInfo[stage].classLog, stageLog->classLog->numClasses);CHKERRQ(ierr);
2159   }
2160 #endif
2161   PetscFunctionReturn(0);
2162 }
2163 
2164 #if defined(PETSC_USE_LOG) && defined(PETSC_HAVE_MPE)
2165 #include <mpe.h>
2166 
2167 PetscBool PetscBeganMPE = PETSC_FALSE;
2168 
2169 PETSC_INTERN PetscErrorCode PetscLogEventBeginMPE(PetscLogEvent,int,PetscObject,PetscObject,PetscObject,PetscObject);
2170 PETSC_INTERN PetscErrorCode PetscLogEventEndMPE(PetscLogEvent,int,PetscObject,PetscObject,PetscObject,PetscObject);
2171 
2172 /*@C
2173    PetscLogMPEBegin - Turns on MPE logging of events. This creates large log files
2174    and slows the program down.
2175 
2176    Collective over PETSC_COMM_WORLD
2177 
2178    Options Database Keys:
2179 . -log_mpe - Prints extensive log information
2180 
2181    Notes:
2182    A related routine is PetscLogDefaultBegin() (with the options key -log_view), which is
2183    intended for production runs since it logs only flop rates and object
2184    creation (and should not significantly slow the programs).
2185 
2186    Level: advanced
2187 
2188 
2189 .seealso: PetscLogDump(), PetscLogDefaultBegin(), PetscLogAllBegin(), PetscLogEventActivate(),
2190           PetscLogEventDeactivate()
2191 @*/
2192 PetscErrorCode  PetscLogMPEBegin(void)
2193 {
2194   PetscErrorCode ierr;
2195 
2196   PetscFunctionBegin;
2197   /* Do MPE initialization */
2198   if (!MPE_Initialized_logging()) { /* This function exists in mpich 1.1.2 and higher */
2199     ierr = PetscInfo(0,"Initializing MPE.\n");CHKERRQ(ierr);
2200     ierr = MPE_Init_log();CHKERRQ(ierr);
2201 
2202     PetscBeganMPE = PETSC_TRUE;
2203   } else {
2204     ierr = PetscInfo(0,"MPE already initialized. Not attempting to reinitialize.\n");CHKERRQ(ierr);
2205   }
2206   ierr = PetscLogSet(PetscLogEventBeginMPE, PetscLogEventEndMPE);CHKERRQ(ierr);
2207   PetscFunctionReturn(0);
2208 }
2209 
2210 /*@C
2211    PetscLogMPEDump - Dumps the MPE logging info to file for later use with Jumpshot.
2212 
2213    Collective over PETSC_COMM_WORLD
2214 
2215    Level: advanced
2216 
2217 .seealso: PetscLogDump(), PetscLogAllBegin(), PetscLogMPEBegin()
2218 @*/
2219 PetscErrorCode  PetscLogMPEDump(const char sname[])
2220 {
2221   char           name[PETSC_MAX_PATH_LEN];
2222   PetscErrorCode ierr;
2223 
2224   PetscFunctionBegin;
2225   if (PetscBeganMPE) {
2226     ierr = PetscInfo(0,"Finalizing MPE.\n");CHKERRQ(ierr);
2227     if (sname) {
2228       ierr = PetscStrcpy(name,sname);CHKERRQ(ierr);
2229     } else {
2230       ierr = PetscGetProgramName(name,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
2231     }
2232     ierr = MPE_Finish_log(name);CHKERRQ(ierr);
2233   } else {
2234     ierr = PetscInfo(0,"Not finalizing MPE (not started by PETSc).\n");CHKERRQ(ierr);
2235   }
2236   PetscFunctionReturn(0);
2237 }
2238 
2239 #define PETSC_RGB_COLORS_MAX 39
2240 static const char *PetscLogMPERGBColors[PETSC_RGB_COLORS_MAX] = {
2241   "OliveDrab:      ",
2242   "BlueViolet:     ",
2243   "CadetBlue:      ",
2244   "CornflowerBlue: ",
2245   "DarkGoldenrod:  ",
2246   "DarkGreen:      ",
2247   "DarkKhaki:      ",
2248   "DarkOliveGreen: ",
2249   "DarkOrange:     ",
2250   "DarkOrchid:     ",
2251   "DarkSeaGreen:   ",
2252   "DarkSlateGray:  ",
2253   "DarkTurquoise:  ",
2254   "DeepPink:       ",
2255   "DarkKhaki:      ",
2256   "DimGray:        ",
2257   "DodgerBlue:     ",
2258   "GreenYellow:    ",
2259   "HotPink:        ",
2260   "IndianRed:      ",
2261   "LavenderBlush:  ",
2262   "LawnGreen:      ",
2263   "LemonChiffon:   ",
2264   "LightCoral:     ",
2265   "LightCyan:      ",
2266   "LightPink:      ",
2267   "LightSalmon:    ",
2268   "LightSlateGray: ",
2269   "LightYellow:    ",
2270   "LimeGreen:      ",
2271   "MediumPurple:   ",
2272   "MediumSeaGreen: ",
2273   "MediumSlateBlue:",
2274   "MidnightBlue:   ",
2275   "MintCream:      ",
2276   "MistyRose:      ",
2277   "NavajoWhite:    ",
2278   "NavyBlue:       ",
2279   "OliveDrab:      "
2280 };
2281 
2282 /*@C
2283   PetscLogMPEGetRGBColor - This routine returns a rgb color useable with PetscLogEventRegister()
2284 
2285   Not collective. Maybe it should be?
2286 
2287   Output Parameter
2288 . str - character string representing the color
2289 
2290   Level: developer
2291 
2292 .seealso: PetscLogEventRegister
2293 @*/
2294 PetscErrorCode  PetscLogMPEGetRGBColor(const char *str[])
2295 {
2296   static int idx = 0;
2297 
2298   PetscFunctionBegin;
2299   *str = PetscLogMPERGBColors[idx];
2300   idx  = (idx + 1)% PETSC_RGB_COLORS_MAX;
2301   PetscFunctionReturn(0);
2302 }
2303 
2304 #endif /* PETSC_USE_LOG && PETSC_HAVE_MPE */
2305