Actual source code: plog.c

petsc-3.4.5 2014-06-29
  2: /*
  3:       PETSc code to log object creation and destruction and PETSc events.

  5:       This provides the public API used by the rest of PETSc and by users.

  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.

 10: */
 11: #include <petsc-private/logimpl.h>        /*I    "petscsys.h"   I*/
 12: #include <petsctime.h>
 13: #include <petscviewer.h>
 14: #include <petscthreadcomm.h>

 16: PetscLogEvent PETSC_LARGEST_EVENT = PETSC_EVENT;

 18: #if defined(PETSC_CLANGUAGE_CXX)
 19: std::map<std::string,PETSc::LogEvent> PETSc::Log::event_registry;
 20: std::map<std::string,PETSc::LogStage> PETSc::Log::stage_registry;
 21: #endif

 23: #if defined(PETSC_USE_LOG)
 24: #include <petscmachineinfo.h>
 25: #include <petscconfiginfo.h>

 27: /* used in the MPI_XXX() count macros in petsclog.h */

 29: /* Action and object logging variables */
 30: Action    *petsc_actions            = NULL;
 31: Object    *petsc_objects            = NULL;
 32: PetscBool petsc_logActions          = PETSC_FALSE;
 33: PetscBool petsc_logObjects          = PETSC_FALSE;
 34: int       petsc_numActions          = 0, petsc_maxActions = 100;
 35: int       petsc_numObjects          = 0, petsc_maxObjects = 100;
 36: int       petsc_numObjectsDestroyed = 0;

 38: /* Global counters */
 39: PetscLogDouble petsc_BaseTime        = 0.0;
 40: PetscLogDouble petsc_TotalFlops      = 0.0;  /* The number of flops */
 41: PetscLogDouble petsc_tmp_flops       = 0.0;  /* The incremental number of flops */
 42: PetscLogDouble petsc_send_ct         = 0.0;  /* The number of sends */
 43: PetscLogDouble petsc_recv_ct         = 0.0;  /* The number of receives */
 44: PetscLogDouble petsc_send_len        = 0.0;  /* The total length of all sent messages */
 45: PetscLogDouble petsc_recv_len        = 0.0;  /* The total length of all received messages */
 46: PetscLogDouble petsc_isend_ct        = 0.0;  /* The number of immediate sends */
 47: PetscLogDouble petsc_irecv_ct        = 0.0;  /* The number of immediate receives */
 48: PetscLogDouble petsc_isend_len       = 0.0;  /* The total length of all immediate send messages */
 49: PetscLogDouble petsc_irecv_len       = 0.0;  /* The total length of all immediate receive messages */
 50: PetscLogDouble petsc_wait_ct         = 0.0;  /* The number of waits */
 51: PetscLogDouble petsc_wait_any_ct     = 0.0;  /* The number of anywaits */
 52: PetscLogDouble petsc_wait_all_ct     = 0.0;  /* The number of waitalls */
 53: PetscLogDouble petsc_sum_of_waits_ct = 0.0;  /* The total number of waits */
 54: PetscLogDouble petsc_allreduce_ct    = 0.0;  /* The number of reductions */
 55: PetscLogDouble petsc_gather_ct       = 0.0;  /* The number of gathers and gathervs */
 56: PetscLogDouble petsc_scatter_ct      = 0.0;  /* The number of scatters and scattervs */

 58: /* Logging functions */
 59: PetscErrorCode (*PetscLogPHC)(PetscObject) = NULL;
 60: PetscErrorCode (*PetscLogPHD)(PetscObject) = NULL;
 61: PetscErrorCode (*PetscLogPLB)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject) = NULL;
 62: PetscErrorCode (*PetscLogPLE)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject) = NULL;

 64: /* Tracing event logging variables */
 65: FILE             *petsc_tracefile            = NULL;
 66: int              petsc_tracelevel            = 0;
 67: const char       *petsc_traceblanks          = "                                                                                                    ";
 68: char             petsc_tracespace[128]       = " ";
 69: PetscLogDouble   petsc_tracetime             = 0.0;
 70: static PetscBool PetscLogBegin_PrivateCalled = PETSC_FALSE;

 72: /*---------------------------------------------- General Functions --------------------------------------------------*/
 75: /*@C
 76:   PetscLogDestroy - Destroys the object and event logging data and resets the global counters.

 78:   Not Collective

 80:   Notes:
 81:   This routine should not usually be used by programmers. Instead employ
 82:   PetscLogStagePush() and PetscLogStagePop().

 84:   Level: developer

 86: .keywords: log, destroy
 87: .seealso: PetscLogDump(), PetscLogAllBegin(), PetscLogView(), PetscLogStagePush(), PlogStagePop()
 88: @*/
 89: PetscErrorCode  PetscLogDestroy(void)
 90: {
 91:   PetscStageLog  stageLog;

 95:   PetscFree(petsc_actions);
 96:   PetscFree(petsc_objects);
 97:   PetscLogSet(NULL, NULL);

 99:   /* Resetting phase */
100:   PetscLogGetStageLog(&stageLog);
101:   PetscStageLogDestroy(stageLog);

103:   petsc_TotalFlops            = 0.0;
104:   petsc_numActions            = 0;
105:   petsc_numObjects            = 0;
106:   petsc_numObjectsDestroyed   = 0;
107:   petsc_maxActions            = 100;
108:   petsc_maxObjects            = 100;
109:   petsc_actions               = NULL;
110:   petsc_objects               = NULL;
111:   petsc_logActions            = PETSC_FALSE;
112:   petsc_logObjects            = PETSC_FALSE;
113:   petsc_BaseTime              = 0.0;
114:   petsc_TotalFlops            = 0.0;
115:   petsc_tmp_flops             = 0.0;
116:   petsc_send_ct               = 0.0;
117:   petsc_recv_ct               = 0.0;
118:   petsc_send_len              = 0.0;
119:   petsc_recv_len              = 0.0;
120:   petsc_isend_ct              = 0.0;
121:   petsc_irecv_ct              = 0.0;
122:   petsc_isend_len             = 0.0;
123:   petsc_irecv_len             = 0.0;
124:   petsc_wait_ct               = 0.0;
125:   petsc_wait_any_ct           = 0.0;
126:   petsc_wait_all_ct           = 0.0;
127:   petsc_sum_of_waits_ct       = 0.0;
128:   petsc_allreduce_ct          = 0.0;
129:   petsc_gather_ct             = 0.0;
130:   petsc_scatter_ct            = 0.0;
131:   PETSC_LARGEST_EVENT         = PETSC_EVENT;
132:   PetscLogPHC                 = NULL;
133:   PetscLogPHD                 = NULL;
134:   petsc_tracefile             = NULL;
135:   petsc_tracelevel            = 0;
136:   petsc_traceblanks           = "                                                                                                    ";
137:   petsc_tracespace[0]         = ' '; petsc_tracespace[1] = 0;
138:   petsc_tracetime             = 0.0;
139:   PETSC_LARGEST_CLASSID       = PETSC_SMALLEST_CLASSID;
140:   PETSC_OBJECT_CLASSID        = 0;
141:   petsc_stageLog              = 0;
142:   PetscLogBegin_PrivateCalled = PETSC_FALSE;
143:   return(0);
144: }

148: /*@C
149:   PetscLogSet - Sets the logging functions called at the beginning and ending of every event.

151:   Not Collective

153:   Input Parameters:
154: + b - The function called at beginning of event
155: - e - The function called at end of event

157:   Level: developer

159: .seealso: PetscLogDump(), PetscLogBegin(), PetscLogAllBegin(), PetscLogTraceBegin()
160: @*/
161: PetscErrorCode  PetscLogSet(PetscErrorCode (*b)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject),
162:                             PetscErrorCode (*e)(PetscLogEvent, int, PetscObject, PetscObject, PetscObject, PetscObject))
163: {
165:   PetscLogPLB = b;
166:   PetscLogPLE = e;
167:   return(0);
168: }

170: #if defined(PETSC_HAVE_CHUD)
171: #include <CHUD/CHUD.h>
172: #endif
173: #if defined(PETSC_HAVE_PAPI)
174: #include <papi.h>
175: int PAPIEventSet = PAPI_NULL;
176: #endif

178: /*------------------------------------------- Initialization Functions ----------------------------------------------*/
181: PetscErrorCode  PetscLogBegin_Private(void)
182: {
183:   int            stage;
184:   PetscBool      opt;

188:   if (PetscLogBegin_PrivateCalled) return(0);
189:   PetscLogBegin_PrivateCalled = PETSC_TRUE;

191:   PetscOptionsHasName(NULL, "-log_exclude_actions", &opt);
192:   if (opt) petsc_logActions = PETSC_FALSE;
193:   PetscOptionsHasName(NULL, "-log_exclude_objects", &opt);
194:   if (opt) petsc_logObjects = PETSC_FALSE;
195:   if (petsc_logActions) {
196:     PetscMalloc(petsc_maxActions * sizeof(Action), &petsc_actions);
197:   }
198:   if (petsc_logObjects) {
199:     PetscMalloc(petsc_maxObjects * sizeof(Object), &petsc_objects);
200:   }
201:   PetscLogPHC = PetscLogObjCreateDefault;
202:   PetscLogPHD = PetscLogObjDestroyDefault;
203:   /* Setup default logging structures */
204:   PetscStageLogCreate(&petsc_stageLog);
205:   PetscStageLogRegister(petsc_stageLog, "Main Stage", &stage);
206: #if defined(PETSC_HAVE_CHUD)
207:   chudInitialize();
208:   chudAcquireSamplingFacility(CHUD_BLOCKING);
209:   chudSetSamplingDevice(chudCPU1Dev);
210:   chudSetStartDelay(0,chudNanoSeconds);
211:   chudClearPMCMode(chudCPU1Dev,chudUnused);
212:   chudClearPMCs();
213:   /* chudSetPMCMuxPosition(chudCPU1Dev,0,0); */
214:   printf("%s\n",chudGetEventName(chudCPU1Dev,PMC_1,193));
215:   printf("%s\n",chudGetEventDescription(chudCPU1Dev,PMC_1,193));
216:   printf("%s\n",chudGetEventNotes(chudCPU1Dev,PMC_1,193));
217:   chudSetPMCEvent(chudCPU1Dev,PMC_1,193);
218:   chudSetPMCMode(chudCPU1Dev,PMC_1,chudCounter);
219:   chudSetPrivilegeFilter(chudCPU1Dev,PMC_1,chudCountUserEvents);
220:   chudSetPMCEventMask(chudCPU1Dev,PMC_1,0xFE);
221:   if (!chudIsEventValid(chudCPU1Dev,PMC_1,193)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Event is not valid %d",193);
222:   chudStartPMCs();
223: #endif
224: #if defined(PETSC_HAVE_PAPI)
225:   PAPI_library_init(PAPI_VER_CURRENT);
226:   if (ierr != PAPI_VER_CURRENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot initialize PAPI");
227:   PAPI_query_event(PAPI_FP_INS);
228:   PAPI_create_eventset(&PAPIEventSet);
229:   PAPI_add_event(PAPIEventSet,PAPI_FP_INS);
230:   PAPI_start(PAPIEventSet);
231: #endif

233:   /* All processors sync here for more consistent logging */
234:   MPI_Barrier(PETSC_COMM_WORLD);
235:   PetscTime(&petsc_BaseTime);
236:   PetscLogStagePush(stage);
237:   return(0);
238: }

242: /*@C
243:   PetscLogBegin - Turns on logging of objects and events. This logs flop
244:   rates and object creation and should not slow programs down too much.
245:   This routine may be called more than once.

247:   Logically Collective over PETSC_COMM_WORLD

249:   Options Database Keys:
250: + -log_summary - Prints summary of flop and timing information to the
251:                   screen (for code compiled with PETSC_USE_LOG)
252: - -log - Prints detailed log information (for code compiled with PETSC_USE_LOG)

254:   Usage:
255: .vb
256:       PetscInitialize(...);
257:       PetscLogBegin();
258:        ... code ...
259:       PetscLogView(viewer); or PetscLogDump();
260:       PetscFinalize();
261: .ve

263:   Notes:
264:   PetscLogView(viewer) or PetscLogDump() actually cause the printing of
265:   the logging information.

267:   Level: advanced

269: .keywords: log, begin
270: .seealso: PetscLogDump(), PetscLogAllBegin(), PetscLogView(), PetscLogTraceBegin()
271: @*/
272: PetscErrorCode  PetscLogBegin(void)
273: {

277:   PetscLogSet(PetscLogEventBeginDefault, PetscLogEventEndDefault);
278:   PetscLogBegin_Private();
279:   return(0);
280: }

284: /*@C
285:   PetscLogAllBegin - Turns on extensive logging of objects and events. Logs
286:   all events. This creates large log files and slows the program down.

288:   Logically Collective on PETSC_COMM_WORLD

290:   Options Database Keys:
291: . -log_all - Prints extensive log information (for code compiled with PETSC_USE_LOG)

293:   Usage:
294: .vb
295:      PetscInitialize(...);
296:      PetscLogAllBegin();
297:      ... code ...
298:      PetscLogDump(filename);
299:      PetscFinalize();
300: .ve

302:   Notes:
303:   A related routine is PetscLogBegin() (with the options key -log), which is
304:   intended for production runs since it logs only flop rates and object
305:   creation (and shouldn't significantly slow the programs).

307:   Level: advanced

309: .keywords: log, all, begin
310: .seealso: PetscLogDump(), PetscLogBegin(), PetscLogTraceBegin()
311: @*/
312: PetscErrorCode  PetscLogAllBegin(void)
313: {

317:   PetscLogSet(PetscLogEventBeginComplete, PetscLogEventEndComplete);
318:   PetscLogBegin_Private();
319:   return(0);
320: }

324: /*@
325:   PetscLogTraceBegin - Activates trace logging.  Every time a PETSc event
326:   begins or ends, the event name is printed.

328:   Logically Collective on PETSC_COMM_WORLD

330:   Input Parameter:
331: . file - The file to print trace in (e.g. stdout)

333:   Options Database Key:
334: . -log_trace [filename] - Activates PetscLogTraceBegin()

336:   Notes:
337:   PetscLogTraceBegin() prints the processor number, the execution time (sec),
338:   then "Event begin:" or "Event end:" followed by the event name.

340:   PetscLogTraceBegin() allows tracing of all PETSc calls, which is useful
341:   to determine where a program is hanging without running in the
342:   debugger.  Can be used in conjunction with the -info option.

344:   Level: intermediate

346: .seealso: PetscLogDump(), PetscLogAllBegin(), PetscLogView(), PetscLogBegin()
347: @*/
348: PetscErrorCode  PetscLogTraceBegin(FILE *file)
349: {

353:   petsc_tracefile = file;

355:   PetscLogSet(PetscLogEventBeginTrace, PetscLogEventEndTrace);
356:   PetscLogBegin_Private();
357:   return(0);
358: }

362: /*@
363:   PetscLogActions - Determines whether actions are logged for the graphical viewer.

365:   Not Collective

367:   Input Parameter:
368: . flag - PETSC_TRUE if actions are to be logged

370:   Level: intermediate

372:   Note: Logging of actions continues to consume more memory as the program
373:   runs. Long running programs should consider turning this feature off.

375:   Options Database Keys:
376: . -log_exclude_actions - Turns off actions logging

378: .keywords: log, stage, register
379: .seealso: PetscLogStagePush(), PetscLogStagePop()
380: @*/
381: PetscErrorCode  PetscLogActions(PetscBool flag)
382: {
384:   petsc_logActions = flag;
385:   return(0);
386: }

390: /*@
391:   PetscLogObjects - Determines whether objects are logged for the graphical viewer.

393:   Not Collective

395:   Input Parameter:
396: . flag - PETSC_TRUE if objects are to be logged

398:   Level: intermediate

400:   Note: Logging of objects continues to consume more memory as the program
401:   runs. Long running programs should consider turning this feature off.

403:   Options Database Keys:
404: . -log_exclude_objects - Turns off objects logging

406: .keywords: log, stage, register
407: .seealso: PetscLogStagePush(), PetscLogStagePop()
408: @*/
409: PetscErrorCode  PetscLogObjects(PetscBool flag)
410: {
412:   petsc_logObjects = flag;
413:   return(0);
414: }

416: /*------------------------------------------------ Stage Functions --------------------------------------------------*/
419: /*@C
420:   PetscLogStageRegister - Attaches a charactor string name to a logging stage.

422:   Not Collective

424:   Input Parameter:
425: . sname - The name to associate with that stage

427:   Output Parameter:
428: . stage - The stage number

430:   Level: intermediate

432: .keywords: log, stage, register
433: .seealso: PetscLogStagePush(), PetscLogStagePop()
434: @*/
435: PetscErrorCode  PetscLogStageRegister(const char sname[],PetscLogStage *stage)
436: {
437:   PetscStageLog  stageLog;
438:   PetscLogEvent  event;

442:   PetscLogGetStageLog(&stageLog);
443:   PetscStageLogRegister(stageLog, sname, stage);
444:   /* Copy events already changed in the main stage, this sucks */
445:   EventPerfLogEnsureSize(stageLog->stageInfo[*stage].eventLog, stageLog->eventLog->numEvents);
446:   for (event = 0; event < stageLog->eventLog->numEvents; event++) {
447:     EventPerfInfoCopy(&stageLog->stageInfo[0].eventLog->eventInfo[event],&stageLog->stageInfo[*stage].eventLog->eventInfo[event]);
448:   }
449:   ClassPerfLogEnsureSize(stageLog->stageInfo[*stage].classLog, stageLog->classLog->numClasses);
450:   return(0);
451: }

455: /*@C
456:   PetscLogStagePush - This function pushes a stage on the stack.

458:   Not Collective

460:   Input Parameter:
461: . stage - The stage on which to log

463:   Usage:
464:   If the option -log_sumary is used to run the program containing the
465:   following code, then 2 sets of summary data will be printed during
466:   PetscFinalize().
467: .vb
468:       PetscInitialize(int *argc,char ***args,0,0);
469:       [stage 0 of code]
470:       PetscLogStagePush(1);
471:       [stage 1 of code]
472:       PetscLogStagePop();
473:       PetscBarrier(...);
474:       [more stage 0 of code]
475:       PetscFinalize();
476: .ve

478:   Notes:
479:   Use PetscLogStageRegister() to register a stage.

481:   Level: intermediate

483: .keywords: log, push, stage
484: .seealso: PetscLogStagePop(), PetscLogStageRegister(), PetscBarrier()
485: @*/
486: PetscErrorCode  PetscLogStagePush(PetscLogStage stage)
487: {
488:   PetscStageLog  stageLog;

492:   PetscLogGetStageLog(&stageLog);
493:   PetscStageLogPush(stageLog, stage);
494:   return(0);
495: }

499: /*@C
500:   PetscLogStagePop - This function pops a stage from the stack.

502:   Not Collective

504:   Usage:
505:   If the option -log_sumary is used to run the program containing the
506:   following code, then 2 sets of summary data will be printed during
507:   PetscFinalize().
508: .vb
509:       PetscInitialize(int *argc,char ***args,0,0);
510:       [stage 0 of code]
511:       PetscLogStagePush(1);
512:       [stage 1 of code]
513:       PetscLogStagePop();
514:       PetscBarrier(...);
515:       [more stage 0 of code]
516:       PetscFinalize();
517: .ve

519:   Notes:
520:   Use PetscLogStageRegister() to register a stage.

522:   Level: intermediate

524: .keywords: log, pop, stage
525: .seealso: PetscLogStagePush(), PetscLogStageRegister(), PetscBarrier()
526: @*/
527: PetscErrorCode  PetscLogStagePop(void)
528: {
529:   PetscStageLog  stageLog;

533:   PetscLogGetStageLog(&stageLog);
534:   PetscStageLogPop(stageLog);
535:   return(0);
536: }

540: /*@
541:   PetscLogStageSetActive - Determines stage activity for PetscLogEventBegin() and PetscLogEventEnd().

543:   Not Collective

545:   Input Parameters:
546: + stage    - The stage
547: - isActive - The activity flag, PETSC_TRUE for logging, else PETSC_FALSE (defaults to PETSC_TRUE)

549:   Level: intermediate

551: .seealso: PetscLogStagePush(), PetscLogStagePop(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadBegin(), PetscPreLoadEnd(), PetscPreLoadStage()
552: @*/
553: PetscErrorCode  PetscLogStageSetActive(PetscLogStage stage, PetscBool isActive)
554: {
555:   PetscStageLog  stageLog;

559:   PetscLogGetStageLog(&stageLog);
560:   PetscStageLogSetActive(stageLog, stage, isActive);
561:   return(0);
562: }

566: /*@
567:   PetscLogStageGetActive - Returns stage activity for PetscLogEventBegin() and PetscLogEventEnd().

569:   Not Collective

571:   Input Parameter:
572: . stage    - The stage

574:   Output Parameter:
575: . isActive - The activity flag, PETSC_TRUE for logging, else PETSC_FALSE (defaults to PETSC_TRUE)

577:   Level: intermediate

579: .seealso: PetscLogStagePush(), PetscLogStagePop(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadBegin(), PetscPreLoadEnd(), PetscPreLoadStage()
580: @*/
581: PetscErrorCode  PetscLogStageGetActive(PetscLogStage stage, PetscBool  *isActive)
582: {
583:   PetscStageLog  stageLog;

587:   PetscLogGetStageLog(&stageLog);
588:   PetscStageLogGetActive(stageLog, stage, isActive);
589:   return(0);
590: }

594: /*@
595:   PetscLogStageSetVisible - Determines stage visibility in PetscLogView()

597:   Not Collective

599:   Input Parameters:
600: + stage     - The stage
601: - isVisible - The visibility flag, PETSC_TRUE to print, else PETSC_FALSE (defaults to PETSC_TRUE)

603:   Level: intermediate

605: .seealso: PetscLogStagePush(), PetscLogStagePop(), PetscLogView()
606: @*/
607: PetscErrorCode  PetscLogStageSetVisible(PetscLogStage stage, PetscBool isVisible)
608: {
609:   PetscStageLog  stageLog;

613:   PetscLogGetStageLog(&stageLog);
614:   PetscStageLogSetVisible(stageLog, stage, isVisible);
615:   return(0);
616: }

620: /*@
621:   PetscLogStageGetVisible - Returns stage visibility in PetscLogView()

623:   Not Collective

625:   Input Parameter:
626: . stage     - The stage

628:   Output Parameter:
629: . isVisible - The visibility flag, PETSC_TRUE to print, else PETSC_FALSE (defaults to PETSC_TRUE)

631:   Level: intermediate

633: .seealso: PetscLogStagePush(), PetscLogStagePop(), PetscLogView()
634: @*/
635: PetscErrorCode  PetscLogStageGetVisible(PetscLogStage stage, PetscBool  *isVisible)
636: {
637:   PetscStageLog  stageLog;

641:   PetscLogGetStageLog(&stageLog);
642:   PetscStageLogGetVisible(stageLog, stage, isVisible);
643:   return(0);
644: }

648: /*@C
649:   PetscLogStageGetId - Returns the stage id when given the stage name.

651:   Not Collective

653:   Input Parameter:
654: . name  - The stage name

656:   Output Parameter:
657: . stage - The stage

659:   Level: intermediate

661: .seealso: PetscLogStagePush(), PetscLogStagePop(), PetscPreLoadBegin(), PetscPreLoadEnd(), PetscPreLoadStage()
662: @*/
663: PetscErrorCode  PetscLogStageGetId(const char name[], PetscLogStage *stage)
664: {
665:   PetscStageLog  stageLog;

669:   PetscLogGetStageLog(&stageLog);
670:   PetscStageLogGetStage(stageLog, name, stage);
671:   return(0);
672: }

674: /*------------------------------------------------ Event Functions --------------------------------------------------*/
677: /*@C
678:   PetscLogEventRegister - Registers an event name for logging operations in an application code.

680:   Not Collective

682:   Input Parameter:
683: + name   - The name associated with the event
684: - classid - The classid associated to the class for this event, obtain either with
685:            PetscClassIdRegister() or use a predefined one such as KSP_CLASSID, SNES_CLASSID, the predefined ones
686:            are only available in C code

688:   Output Parameter:
689: . event - The event id for use with PetscLogEventBegin() and PetscLogEventEnd().

691:   Example of Usage:
692: .vb
693:       PetscLogEvent USER_EVENT;
694:       PetscClassId classid;
695:       PetscLogDouble user_event_flops;
696:       PetscClassIdRegister("class name",&classid);
697:       PetscLogEventRegister("User event name",classid,&USER_EVENT);
698:       PetscLogEventBegin(USER_EVENT,0,0,0,0);
699:          [code segment to monitor]
700:          PetscLogFlops(user_event_flops);
701:       PetscLogEventEnd(USER_EVENT,0,0,0,0);
702: .ve

704:   Notes:
705:   PETSc automatically logs library events if the code has been
706:   compiled with -DPETSC_USE_LOG (which is the default) and -log,
707:   -log_summary, or -log_all are specified.  PetscLogEventRegister() is
708:   intended for logging user events to supplement this PETSc
709:   information.

711:   PETSc can gather data for use with the utilities Jumpshot
712:   (part of the MPICH distribution).  If PETSc has been compiled
713:   with flag -DPETSC_HAVE_MPE (MPE is an additional utility within
714:   MPICH), the user can employ another command line option, -log_mpe,
715:   to create a logfile, "mpe.log", which can be visualized
716:   Jumpshot.

718:   The classid is associated with each event so that classes of events
719:   can be disabled simultaneously, such as all matrix events. The user
720:   can either use an existing classid, such as MAT_CLASSID, or create
721:   their own as shown in the example.

723:   Level: intermediate

725: .keywords: log, event, register
726: .seealso: PetscLogEventBegin(), PetscLogEventEnd(), PetscLogFlops(),
727:           PetscLogEventMPEActivate(), PetscLogEventMPEDeactivate(),
728:           PetscLogEventActivate(), PetscLogEventDeactivate(), PetscClassIdRegister()
729: @*/
730: PetscErrorCode  PetscLogEventRegister(const char name[],PetscClassId classid,PetscLogEvent *event)
731: {
732:   PetscStageLog  stageLog;
733:   int            stage;

737:   *event = PETSC_DECIDE;
738:   PetscLogGetStageLog(&stageLog);
739:   EventRegLogRegister(stageLog->eventLog, name, classid, event);
740:   for (stage = 0; stage < stageLog->numStages; stage++) {
741:     EventPerfLogEnsureSize(stageLog->stageInfo[stage].eventLog, stageLog->eventLog->numEvents);
742:     ClassPerfLogEnsureSize(stageLog->stageInfo[stage].classLog, stageLog->classLog->numClasses);
743:   }
744:   return(0);
745: }

749: /*@
750:   PetscLogEventActivate - Indicates that a particular event should be logged.

752:   Not Collective

754:   Input Parameter:
755: . event - The event id

757:   Usage:
758: .vb
759:       PetscLogEventDeactivate(VEC_SetValues);
760:         [code where you do not want to log VecSetValues()]
761:       PetscLogEventActivate(VEC_SetValues);
762:         [code where you do want to log VecSetValues()]
763: .ve

765:   Note:
766:   The event may be either a pre-defined PETSc event (found in include/petsclog.h)
767:   or an event number obtained with PetscLogEventRegister().

769:   Level: advanced

771: .keywords: log, event, activate
772: .seealso: PetscLogEventMPEDeactivate(),PetscLogEventMPEActivate(),PlogEventDeactivate()
773: @*/
774: PetscErrorCode  PetscLogEventActivate(PetscLogEvent event)
775: {
776:   PetscStageLog  stageLog;
777:   int            stage;

781:   PetscLogGetStageLog(&stageLog);
782:   PetscStageLogGetCurrent(stageLog, &stage);
783:   EventPerfLogActivate(stageLog->stageInfo[stage].eventLog, event);
784:   return(0);
785: }

789: /*@
790:   PetscLogEventDeactivate - Indicates that a particular event should not be logged.

792:   Not Collective

794:   Input Parameter:
795: . event - The event id

797:   Usage:
798: .vb
799:       PetscLogEventDeactivate(VEC_SetValues);
800:         [code where you do not want to log VecSetValues()]
801:       PetscLogEventActivate(VEC_SetValues);
802:         [code where you do want to log VecSetValues()]
803: .ve

805:   Note:
806:   The event may be either a pre-defined PETSc event (found in
807:   include/petsclog.h) or an event number obtained with PetscLogEventRegister()).

809:   Level: advanced

811: .keywords: log, event, deactivate
812: .seealso: PetscLogEventMPEDeactivate(),PetscLogEventMPEActivate(),PlogEventActivate()
813: @*/
814: PetscErrorCode  PetscLogEventDeactivate(PetscLogEvent event)
815: {
816:   PetscStageLog  stageLog;
817:   int            stage;

821:   PetscLogGetStageLog(&stageLog);
822:   PetscStageLogGetCurrent(stageLog, &stage);
823:   EventPerfLogDeactivate(stageLog->stageInfo[stage].eventLog, event);
824:   return(0);
825: }

829: /*@
830:   PetscLogEventSetActiveAll - Sets the event activity in every stage.

832:   Not Collective

834:   Input Parameters:
835: + event    - The event id
836: - isActive - The activity flag determining whether the event is logged

838:   Level: advanced

840: .keywords: log, event, activate
841: .seealso: PetscLogEventMPEDeactivate(),PetscLogEventMPEActivate(),PlogEventActivate(),PlogEventDeactivate()
842: @*/
843: PetscErrorCode  PetscLogEventSetActiveAll(PetscLogEvent event, PetscBool isActive)
844: {
845:   PetscStageLog  stageLog;
846:   int            stage;

850:   PetscLogGetStageLog(&stageLog);
851:   for (stage = 0; stage < stageLog->numStages; stage++) {
852:     if (isActive) {
853:       EventPerfLogActivate(stageLog->stageInfo[stage].eventLog, event);
854:     } else {
855:       EventPerfLogDeactivate(stageLog->stageInfo[stage].eventLog, event);
856:     }
857:   }
858:   return(0);
859: }

863: /*@
864:   PetscLogEventActivateClass - Activates event logging for a PETSc object class.

866:   Not Collective

868:   Input Parameter:
869: . classid - The event class, for example MAT_CLASSID, SNES_CLASSID, etc.

871:   Level: developer

873: .keywords: log, event, activate, class
874: .seealso: PetscInfoActivate(),PetscInfo(),PetscInfoAllow(),PetscLogEventDeactivateClass(), PetscLogEventActivate(),PetscLogEventDeactivate()
875: @*/
876: PetscErrorCode  PetscLogEventActivateClass(PetscClassId classid)
877: {
878:   PetscStageLog  stageLog;
879:   int            stage;

883:   PetscLogGetStageLog(&stageLog);
884:   PetscStageLogGetCurrent(stageLog, &stage);
885:   EventPerfLogActivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid);
886:   return(0);
887: }

891: /*@
892:   PetscLogEventDeactivateClass - Deactivates event logging for a PETSc object class.

894:   Not Collective

896:   Input Parameter:
897: . classid - The event class, for example MAT_CLASSID, SNES_CLASSID, etc.

899:   Level: developer

901: .keywords: log, event, deactivate, class
902: .seealso: PetscInfoActivate(),PetscInfo(),PetscInfoAllow(),PetscLogEventActivateClass(), PetscLogEventActivate(),PetscLogEventDeactivate()
903: @*/
904: PetscErrorCode  PetscLogEventDeactivateClass(PetscClassId classid)
905: {
906:   PetscStageLog  stageLog;
907:   int            stage;

911:   PetscLogGetStageLog(&stageLog);
912:   PetscStageLogGetCurrent(stageLog, &stage);
913:   EventPerfLogDeactivateClass(stageLog->stageInfo[stage].eventLog, stageLog->eventLog, classid);
914:   return(0);
915: }

917: /*MC
918:    PetscLogEventBegin - Logs the beginning of a user event.

920:    Synopsis:
921:    #include "petsclog.h"
922:    PetscErrorCode PetscLogEventBegin(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)

924:    Not Collective

926:    Input Parameters:
927: +  e - integer associated with the event obtained from PetscLogEventRegister()
928: -  o1,o2,o3,o4 - objects associated with the event, or 0


931:    Fortran Synopsis:
932:    void PetscLogEventBegin(int e,PetscErrorCode ierr)

934:    Usage:
935: .vb
936:      PetscLogEvent USER_EVENT;
937:      PetscLogDouble user_event_flops;
938:      PetscLogEventRegister("User event",0,&USER_EVENT);
939:      PetscLogEventBegin(USER_EVENT,0,0,0,0);
940:         [code segment to monitor]
941:         PetscLogFlops(user_event_flops);
942:      PetscLogEventEnd(USER_EVENT,0,0,0,0);
943: .ve

945:    Notes:
946:    You need to register each integer event with the command
947:    PetscLogEventRegister().  The source code must be compiled with
948:    -DPETSC_USE_LOG, which is the default.

950:    PETSc automatically logs library events if the code has been
951:    compiled with -DPETSC_USE_LOG, and -log, -log_summary, or -log_all are
952:    specified.  PetscLogEventBegin() is intended for logging user events
953:    to supplement this PETSc information.

955:    Level: intermediate

957: .seealso: PetscLogEventRegister(), PetscLogEventEnd(), PetscLogFlops()

959: .keywords: log, event, begin
960: M*/

962: /*MC
963:    PetscLogEventEnd - Log the end of a user event.

965:    Synopsis:
966:    #include "petsclog.h"
967:    PetscErrorCode PetscLogEventEnd(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)

969:    Not Collective

971:    Input Parameters:
972: +  e - integer associated with the event obtained with PetscLogEventRegister()
973: -  o1,o2,o3,o4 - objects associated with the event, or 0


976:    Fortran Synopsis:
977:    void PetscLogEventEnd(int e,PetscErrorCode ierr)

979:    Usage:
980: .vb
981:      PetscLogEvent USER_EVENT;
982:      PetscLogDouble user_event_flops;
983:      PetscLogEventRegister("User event",0,&USER_EVENT,);
984:      PetscLogEventBegin(USER_EVENT,0,0,0,0);
985:         [code segment to monitor]
986:         PetscLogFlops(user_event_flops);
987:      PetscLogEventEnd(USER_EVENT,0,0,0,0);
988: .ve

990:    Notes:
991:    You should also register each additional integer event with the command
992:    PetscLogEventRegister(). Source code must be compiled with
993:    -DPETSC_USE_LOG, which is the default.

995:    PETSc automatically logs library events if the code has been
996:    compiled with -DPETSC_USE_LOG, and -log, -log_summary, or -log_all are
997:    specified.  PetscLogEventEnd() is intended for logging user events
998:    to supplement this PETSc information.

1000:    Level: intermediate

1002: .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogFlops()

1004: .keywords: log, event, end
1005: M*/

1007: /*MC
1008:    PetscLogEventBarrierBegin - Logs the time in a barrier before an event.

1010:    Synopsis:
1011:    #include "petsclog.h"
1012:    PetscErrorCode PetscLogEventBarrierBegin(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4,MPI_Comm comm)

1014:    Not Collective

1016:    Input Parameters:
1017: .  e - integer associated with the event obtained from PetscLogEventRegister()
1018: .  o1,o2,o3,o4 - objects associated with the event, or 0
1019: .  comm - communicator the barrier takes place over


1022:    Usage:
1023: .vb
1024:      PetscLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
1025:        MPI_Allreduce()
1026:      PetscLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm);
1027: .ve

1029:    Notes:
1030:    This is for logging the amount of time spent in a barrier for an event
1031:    that requires synchronization.

1033:    Additional Notes:
1034:    Synchronization events always come in pairs; for example, VEC_NormBarrier and
1035:    VEC_NormComm = VEC_NormBarrier + 1

1037:    Level: advanced

1039: .seealso: PetscLogEventRegister(), PetscLogEventEnd(), PetscLogFlops(), PetscLogEventBegin(),
1040:           PetscLogEventBarrierEnd()

1042: .keywords: log, event, begin, barrier
1043: M*/

1045: /*MC
1046:    PetscLogEventBarrierEnd - Logs the time in a barrier before an event.

1048:    Synopsis:
1049:    #include "petsclog.h"
1050:    PetscErrorCode PetscLogEventBarrierEnd(int e,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4,MPI_Comm comm)

1052:    Logically Collective on MPI_Comm

1054:    Input Parameters:
1055: .  e - integer associated with the event obtained from PetscLogEventRegister()
1056: .  o1,o2,o3,o4 - objects associated with the event, or 0
1057: .  comm - communicator the barrier takes place over


1060:     Usage:
1061: .vb
1062:      PetscLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
1063:        MPI_Allreduce()
1064:      PetscLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm);
1065: .ve

1067:    Notes:
1068:    This is for logging the amount of time spent in a barrier for an event
1069:    that requires synchronization.

1071:    Additional Notes:
1072:    Synchronization events always come in pairs; for example, VEC_NormBarrier and
1073:    VEC_NormComm = VEC_NormBarrier + 1

1075:    Level: advanced

1077: .seealso: PetscLogEventRegister(), PetscLogEventEnd(), PetscLogFlops(), PetscLogEventBegin(),
1078:           PetscLogEventBarrierBegin()

1080: .keywords: log, event, begin, barrier
1081: M*/

1085: /*@C
1086:   PetscLogEventGetId - Returns the event id when given the event name.

1088:   Not Collective

1090:   Input Parameter:
1091: . name  - The event name

1093:   Output Parameter:
1094: . event - The event

1096:   Level: intermediate

1098: .seealso: PetscLogEventBegin(), PetscLogEventEnd(), PetscLogStageGetId()
1099: @*/
1100: PetscErrorCode  PetscLogEventGetId(const char name[], PetscLogEvent *event)
1101: {
1102:   PetscStageLog  stageLog;

1106:   PetscLogGetStageLog(&stageLog);
1107:   EventRegLogGetEvent(stageLog->eventLog, name, event);
1108:   return(0);
1109: }


1112: /*------------------------------------------------ Output Functions -------------------------------------------------*/
1115: /*@C
1116:   PetscLogDump - Dumps logs of objects to a file. This file is intended to
1117:   be read by bin/petscview. This program no longer exists.

1119:   Collective on PETSC_COMM_WORLD

1121:   Input Parameter:
1122: . name - an optional file name

1124:   Options Database Keys:
1125: + -log     - Prints basic log information (for code compiled with PETSC_USE_LOG)
1126: - -log_all - Prints extensive log information (for code compiled with PETSC_USE_LOG)

1128:   Usage:
1129: .vb
1130:      PetscInitialize(...);
1131:      PetscLogBegin(); or PetscLogAllBegin();
1132:      ... code ...
1133:      PetscLogDump(filename);
1134:      PetscFinalize();
1135: .ve

1137:   Notes:
1138:   The default file name is
1139: $    Log.<rank>
1140:   where <rank> is the processor number. If no name is specified,
1141:   this file will be used.

1143:   Level: advanced

1145: .keywords: log, dump
1146: .seealso: PetscLogBegin(), PetscLogAllBegin(), PetscLogView()
1147: @*/
1148: PetscErrorCode  PetscLogDump(const char sname[])
1149: {
1150:   PetscStageLog      stageLog;
1151:   PetscEventPerfInfo *eventInfo;
1152:   FILE               *fd;
1153:   char               file[PETSC_MAX_PATH_LEN], fname[PETSC_MAX_PATH_LEN];
1154:   PetscLogDouble     flops, _TotalTime;
1155:   PetscMPIInt        rank;
1156:   int                action, object, curStage;
1157:   PetscLogEvent      event;
1158:   PetscErrorCode     ierr;

1161:   /* Calculate the total elapsed time */
1162:   PetscTime(&_TotalTime);
1163:   _TotalTime -= petsc_BaseTime;
1164:   /* Open log file */
1165:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
1166:   if (sname) sprintf(file, "%s.%d", sname, rank);
1167:   else sprintf(file, "Log.%d", rank);
1168:   PetscFixFilename(file, fname);
1169:   PetscFOpen(PETSC_COMM_WORLD, fname, "w", &fd);
1170:   if ((!rank) && (!fd)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN, "Cannot open file: %s", fname);
1171:   /* Output totals */
1172:   PetscFPrintf(PETSC_COMM_WORLD, fd, "Total Flops %14e %16.8e\n", petsc_TotalFlops, _TotalTime);
1173:   PetscFPrintf(PETSC_COMM_WORLD, fd, "Clock Resolution %g\n", 0.0);
1174:   /* Output actions */
1175:   if (petsc_logActions) {
1176:     PetscFPrintf(PETSC_COMM_WORLD, fd, "Actions accomplished %d\n", petsc_numActions);
1177:     for (action = 0; action < petsc_numActions; action++) {
1178:       PetscFPrintf(PETSC_COMM_WORLD, fd, "%g %d %d %d %d %d %d %g %g %g\n",
1179:                           petsc_actions[action].time, petsc_actions[action].action, (int)petsc_actions[action].event, (int)petsc_actions[action].classid, petsc_actions[action].id1,
1180:                           petsc_actions[action].id2, petsc_actions[action].id3, petsc_actions[action].flops, petsc_actions[action].mem, petsc_actions[action].maxmem);
1181:     }
1182:   }
1183:   /* Output objects */
1184:   if (petsc_logObjects) {
1185:     PetscFPrintf(PETSC_COMM_WORLD, fd, "Objects created %d destroyed %d\n", petsc_numObjects, petsc_numObjectsDestroyed);
1186:     for (object = 0; object < petsc_numObjects; object++) {
1187:       PetscFPrintf(PETSC_COMM_WORLD, fd, "Parent ID: %d Memory: %d\n", petsc_objects[object].parent, (int) petsc_objects[object].mem);
1188:       if (!petsc_objects[object].name[0]) {
1189:         PetscFPrintf(PETSC_COMM_WORLD, fd,"No Name\n");
1190:       } else {
1191:         PetscFPrintf(PETSC_COMM_WORLD, fd, "Name: %s\n", petsc_objects[object].name);
1192:       }
1193:       if (petsc_objects[object].info[0] != 0) {
1194:         PetscFPrintf(PETSC_COMM_WORLD, fd, "No Info\n");
1195:       } else {
1196:         PetscFPrintf(PETSC_COMM_WORLD, fd, "Info: %s\n", petsc_objects[object].info);
1197:       }
1198:     }
1199:   }
1200:   /* Output events */
1201:   PetscFPrintf(PETSC_COMM_WORLD, fd, "Event log:\n");
1202:   PetscLogGetStageLog(&stageLog);
1203:   PetscIntStackTop(stageLog->stack, &curStage);
1204:   eventInfo = stageLog->stageInfo[curStage].eventLog->eventInfo;
1205:   for (event = 0; event < stageLog->stageInfo[curStage].eventLog->numEvents; event++) {
1206:     if (eventInfo[event].time != 0.0) flops = eventInfo[event].flops/eventInfo[event].time;
1207:     else flops = 0.0;
1208:     PetscFPrintf(PETSC_COMM_WORLD, fd, "%d %16d %16g %16g %16g\n", event, eventInfo[event].count,
1209:                         eventInfo[event].flops, eventInfo[event].time, flops);
1210:   }
1211:   PetscFClose(PETSC_COMM_WORLD, fd);
1212:   return(0);
1213: }

1217: /*@C
1218:   PetscLogViewer - Prints a summary of the logging.

1220:   Collective over MPI_Comm

1222:   Input Parameter:
1223: .  viewer - an ASCII viewer

1225:   Options Database Keys:
1226: . -log_summary - Prints summary of log information (for code compiled with PETSC_USE_LOG)

1228:   Usage:
1229: .vb
1230:      PetscInitialize(...);
1231:      PetscLogBegin();
1232:      ... code ...
1233:      PetscLogView(PetscViewer);
1234:      PetscFinalize(...);
1235: .ve

1237:   Notes:
1238:   By default the summary is printed to stdout.

1240:   Level: beginner

1242: .keywords: log, dump, print
1243: .seealso: PetscLogBegin(), PetscLogDump()
1244: @*/
1245: PetscErrorCode  PetscLogView(PetscViewer viewer)
1246: {
1247:   FILE               *fd;
1248:   PetscLogDouble     zero       = 0.0;
1249:   PetscStageLog      stageLog;
1250:   PetscStageInfo     *stageInfo = NULL;
1251:   PetscEventPerfInfo *eventInfo = NULL;
1252:   PetscClassPerfInfo *classInfo;
1253:   char               arch[128],hostname[128],username[128],pname[PETSC_MAX_PATH_LEN],date[128];
1254:   const char         *name;
1255:   PetscLogDouble     locTotalTime, TotalTime, TotalFlops;
1256:   PetscLogDouble     numMessages, messageLength, avgMessLen, numReductions;
1257:   PetscLogDouble     stageTime, flops, flopr, mem, mess, messLen, red;
1258:   PetscLogDouble     fracTime, fracFlops, fracMessages, fracLength, fracReductions, fracMess, fracMessLen, fracRed;
1259:   PetscLogDouble     fracStageTime, fracStageFlops, fracStageMess, fracStageMessLen, fracStageRed;
1260:   PetscLogDouble     min, max, tot, ratio, avg, x, y;
1261:   PetscLogDouble     minf, maxf, totf, ratf, mint, maxt, tott, ratt, ratCt, totm, totml, totr;
1262:   PetscMPIInt        minCt, maxCt;
1263:   PetscMPIInt        size, rank;
1264:   PetscBool          *localStageUsed,    *stageUsed;
1265:   PetscBool          *localStageVisible, *stageVisible;
1266:   int                numStages, localNumEvents, numEvents;
1267:   int                stage, lastStage, oclass;
1268:   PetscLogEvent      event;
1269:   PetscErrorCode     ierr;
1270:   char               version[256];
1271:   MPI_Comm           comm;
1272:   PetscInt           nthreads;

1275:   PetscObjectGetComm((PetscObject)viewer,&comm);
1276:   if (!PetscLogBegin_PrivateCalled) SETERRQ(comm, PETSC_ERR_ORDER, "No call to PetscLogBegin() before PetscLogView()");
1277:   PetscViewerASCIIGetPointer(viewer,&fd);
1278:   MPI_Comm_size(comm, &size);
1279:   MPI_Comm_rank(comm, &rank);
1280:   /* Pop off any stages the user forgot to remove */
1281:   lastStage = 0;
1282:   PetscLogGetStageLog(&stageLog);
1283:   PetscStageLogGetCurrent(stageLog, &stage);
1284:   while (stage >= 0) {
1285:     lastStage = stage;
1286:     PetscStageLogPop(stageLog);
1287:     PetscStageLogGetCurrent(stageLog, &stage);
1288:   }
1289:   /* Get the total elapsed time */
1290:   PetscTime(&locTotalTime);  locTotalTime -= petsc_BaseTime;

1292:   PetscFPrintf(comm, fd, "************************************************************************************************************************\n");
1293:   PetscFPrintf(comm, fd, "***             WIDEN YOUR WINDOW TO 120 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document            ***\n");
1294:   PetscFPrintf(comm, fd, "************************************************************************************************************************\n");
1295:   PetscFPrintf(comm, fd, "\n---------------------------------------------- PETSc Performance Summary: ----------------------------------------------\n\n");
1296:   PetscGetArchType(arch,sizeof(arch));
1297:   PetscGetHostName(hostname,sizeof(hostname));
1298:   PetscGetUserName(username,sizeof(username));
1299:   PetscGetProgramName(pname,sizeof(pname));
1300:   PetscGetDate(date,sizeof(date));
1301:   PetscGetVersion(version,sizeof(version));
1302:   if (size == 1) {
1303:     PetscFPrintf(comm,fd,"%s on a %s named %s with %d processor, by %s %s\n", pname, arch, hostname, size, username, date);
1304:   } else {
1305:     PetscFPrintf(comm,fd,"%s on a %s named %s with %d processors, by %s %s\n", pname, arch, hostname, size, username, date);
1306:   }
1307:   PetscThreadCommGetNThreads(PETSC_COMM_WORLD,&nthreads);
1308:   if (nthreads > 1) {
1309:     PetscFPrintf(comm,fd,"With %d threads per MPI_Comm\n", (int)nthreads);
1310:   }

1312:   PetscFPrintf(comm, fd, "Using %s\n", version);

1314:   /* Must preserve reduction count before we go on */
1315:   red = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;

1317:   /* Calculate summary information */
1318:   PetscFPrintf(comm, fd, "\n                         Max       Max/Min        Avg      Total \n");
1319:   /*   Time */
1320:   MPI_Allreduce(&locTotalTime, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1321:   MPI_Allreduce(&locTotalTime, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1322:   MPI_Allreduce(&locTotalTime, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1323:   avg  = (tot)/((PetscLogDouble) size);
1324:   if (min != 0.0) ratio = max/min;
1325:   else ratio = 0.0;
1326:   PetscFPrintf(comm, fd, "Time (sec):           %5.3e   %10.5f   %5.3e\n", max, ratio, avg);
1327:   TotalTime = tot;
1328:   /*   Objects */
1329:   avg  = (PetscLogDouble) petsc_numObjects;
1330:   MPI_Allreduce(&avg,          &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1331:   MPI_Allreduce(&avg,          &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1332:   MPI_Allreduce(&avg,          &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1333:   avg  = (tot)/((PetscLogDouble) size);
1334:   if (min != 0.0) ratio = max/min;
1335:   else ratio = 0.0;
1336:   PetscFPrintf(comm, fd, "Objects:              %5.3e   %10.5f   %5.3e\n", max, ratio, avg);
1337:   /*   Flops */
1338:   MPI_Allreduce(&petsc_TotalFlops,  &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1339:   MPI_Allreduce(&petsc_TotalFlops,  &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1340:   MPI_Allreduce(&petsc_TotalFlops,  &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1341:   avg  = (tot)/((PetscLogDouble) size);
1342:   if (min != 0.0) ratio = max/min;
1343:   else ratio = 0.0;
1344:   PetscFPrintf(comm, fd, "Flops:                %5.3e   %10.5f   %5.3e  %5.3e\n", max, ratio, avg, tot);
1345:   TotalFlops = tot;
1346:   /*   Flops/sec -- Must talk to Barry here */
1347:   if (locTotalTime != 0.0) flops = petsc_TotalFlops/locTotalTime;
1348:   else flops = 0.0;
1349:   MPI_Allreduce(&flops,        &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1350:   MPI_Allreduce(&flops,        &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1351:   MPI_Allreduce(&flops,        &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1352:   avg  = (tot)/((PetscLogDouble) size);
1353:   if (min != 0.0) ratio = max/min;
1354:   else ratio = 0.0;
1355:   PetscFPrintf(comm, fd, "Flops/sec:            %5.3e   %10.5f   %5.3e  %5.3e\n", max, ratio, avg, tot);
1356:   /*   Memory */
1357:   PetscMallocGetMaximumUsage(&mem);
1358:   if (mem > 0.0) {
1359:     MPI_Allreduce(&mem,          &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1360:     MPI_Allreduce(&mem,          &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1361:     MPI_Allreduce(&mem,          &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1362:     avg  = (tot)/((PetscLogDouble) size);
1363:     if (min != 0.0) ratio = max/min;
1364:     else ratio = 0.0;
1365:     PetscFPrintf(comm, fd, "Memory:               %5.3e   %10.5f              %5.3e\n", max, ratio, tot);
1366:   }
1367:   /*   Messages */
1368:   mess = 0.5*(petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct);
1369:   MPI_Allreduce(&mess,         &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1370:   MPI_Allreduce(&mess,         &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1371:   MPI_Allreduce(&mess,         &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1372:   avg  = (tot)/((PetscLogDouble) size);
1373:   if (min != 0.0) ratio = max/min;
1374:   else ratio = 0.0;
1375:   PetscFPrintf(comm, fd, "MPI Messages:         %5.3e   %10.5f   %5.3e  %5.3e\n", max, ratio, avg, tot);
1376:   numMessages = tot;
1377:   /*   Message Lengths */
1378:   mess = 0.5*(petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len);
1379:   MPI_Allreduce(&mess,         &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1380:   MPI_Allreduce(&mess,         &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1381:   MPI_Allreduce(&mess,         &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1382:   if (numMessages != 0) avg = (tot)/(numMessages);
1383:   else avg = 0.0;
1384:   if (min != 0.0) ratio = max/min;
1385:   else ratio = 0.0;
1386:   PetscFPrintf(comm, fd, "MPI Message Lengths:  %5.3e   %10.5f   %5.3e  %5.3e\n", max, ratio, avg, tot);
1387:   messageLength = tot;
1388:   /*   Reductions */
1389:   MPI_Allreduce(&red,          &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1390:   MPI_Allreduce(&red,          &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1391:   MPI_Allreduce(&red,          &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1392:   if (min != 0.0) ratio = max/min;
1393:   else ratio = 0.0;
1394:   PetscFPrintf(comm, fd, "MPI Reductions:       %5.3e   %10.5f\n", max, ratio);
1395:   numReductions = red; /* wrong because uses count from process zero */
1396:   PetscFPrintf(comm, fd, "\nFlop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)\n");
1397:   PetscFPrintf(comm, fd, "                            e.g., VecAXPY() for real vectors of length N --> 2N flops\n");
1398:   PetscFPrintf(comm, fd, "                            and VecAXPY() for complex vectors of length N --> 8N flops\n");

1400:   /* Get total number of stages --
1401:        Currently, a single processor can register more stages than another, but stages must all be registered in order.
1402:        We can removed this requirement if necessary by having a global stage numbering and indirection on the stage ID.
1403:        This seems best accomplished by assoicating a communicator with each stage.
1404:   */
1405:   MPI_Allreduce(&stageLog->numStages, &numStages, 1, MPI_INT, MPI_MAX, comm);
1406:   PetscMalloc(numStages * sizeof(PetscBool), &localStageUsed);
1407:   PetscMalloc(numStages * sizeof(PetscBool), &stageUsed);
1408:   PetscMalloc(numStages * sizeof(PetscBool), &localStageVisible);
1409:   PetscMalloc(numStages * sizeof(PetscBool), &stageVisible);
1410:   if (numStages > 0) {
1411:     stageInfo = stageLog->stageInfo;
1412:     for (stage = 0; stage < numStages; stage++) {
1413:       if (stage < stageLog->numStages) {
1414:         localStageUsed[stage]    = stageInfo[stage].used;
1415:         localStageVisible[stage] = stageInfo[stage].perfInfo.visible;
1416:       } else {
1417:         localStageUsed[stage]    = PETSC_FALSE;
1418:         localStageVisible[stage] = PETSC_TRUE;
1419:       }
1420:     }
1421:     MPI_Allreduce(localStageUsed,    stageUsed,    numStages, MPIU_BOOL, MPI_LOR,  comm);
1422:     MPI_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm);
1423:     for (stage = 0; stage < numStages; stage++) {
1424:       if (stageUsed[stage]) {
1425:         PetscFPrintf(comm, fd, "\nSummary of Stages:   ----- Time ------  ----- Flops -----  --- Messages ---  -- Message Lengths --  -- Reductions --\n");
1426:         PetscFPrintf(comm, fd, "                        Avg     %%Total     Avg     %%Total   counts   %%Total     Avg         %%Total   counts   %%Total \n");
1427:         break;
1428:       }
1429:     }
1430:     for (stage = 0; stage < numStages; stage++) {
1431:       if (!stageUsed[stage]) continue;
1432:       if (localStageUsed[stage]) {
1433:         MPI_Allreduce(&stageInfo[stage].perfInfo.time,          &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1434:         MPI_Allreduce(&stageInfo[stage].perfInfo.flops,         &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1435:         MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages,   &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1436:         MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1437:         MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1438:         name = stageInfo[stage].name;
1439:       } else {
1440:         MPI_Allreduce(&zero,                           &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1441:         MPI_Allreduce(&zero,                           &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1442:         MPI_Allreduce(&zero,                           &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1443:         MPI_Allreduce(&zero,                           &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1444:         MPI_Allreduce(&zero,                           &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1445:         name = "";
1446:       }
1447:       mess *= 0.5; messLen *= 0.5; red /= size;
1448:       if (TotalTime     != 0.0) fracTime       = stageTime/TotalTime;    else fracTime       = 0.0;
1449:       if (TotalFlops    != 0.0) fracFlops      = flops/TotalFlops;       else fracFlops      = 0.0;
1450:       /* Talk to Barry if (stageTime     != 0.0) flops          = (size*flops)/stageTime; else flops          = 0.0; */
1451:       if (numMessages   != 0.0) fracMessages   = mess/numMessages;       else fracMessages   = 0.0;
1452:       if (numMessages   != 0.0) avgMessLen     = messLen/numMessages;    else avgMessLen     = 0.0;
1453:       if (messageLength != 0.0) fracLength     = messLen/messageLength;  else fracLength     = 0.0;
1454:       if (numReductions != 0.0) fracReductions = red/numReductions;      else fracReductions = 0.0;
1455:       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",
1456:                           stage, name, stageTime/size, 100.0*fracTime, flops, 100.0*fracFlops,
1457:                           mess, 100.0*fracMessages, avgMessLen, 100.0*fracLength, red, 100.0*fracReductions);
1458:     }
1459:   }

1461:   PetscFPrintf(comm, fd,"\n------------------------------------------------------------------------------------------------------------------------\n");
1462:   PetscFPrintf(comm, fd, "See the 'Profiling' chapter of the users' manual for details on interpreting output.\n");
1463:   PetscFPrintf(comm, fd, "Phase summary info:\n");
1464:   PetscFPrintf(comm, fd, "   Count: number of times phase was executed\n");
1465:   PetscFPrintf(comm, fd, "   Time and Flops: Max - maximum over all processors\n");
1466:   PetscFPrintf(comm, fd, "                   Ratio - ratio of maximum to minimum over all processors\n");
1467:   PetscFPrintf(comm, fd, "   Mess: number of messages sent\n");
1468:   PetscFPrintf(comm, fd, "   Avg. len: average message length (bytes)\n");
1469:   PetscFPrintf(comm, fd, "   Reduct: number of global reductions\n");
1470:   PetscFPrintf(comm, fd, "   Global: entire computation\n");
1471:   PetscFPrintf(comm, fd, "   Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().\n");
1472:   PetscFPrintf(comm, fd, "      %%T - percent time in this phase         %%F - percent flops in this phase\n");
1473:   PetscFPrintf(comm, fd, "      %%M - percent messages in this phase     %%L - percent message lengths in this phase\n");
1474:   PetscFPrintf(comm, fd, "      %%R - percent reductions in this phase\n");
1475:   PetscFPrintf(comm, fd, "   Total Mflop/s: 10e-6 * (sum of flops over all processors)/(max time over all processors)\n");
1476:   PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------\n");

1478: #if defined(PETSC_USE_DEBUG)
1479:   PetscFPrintf(comm, fd, "\n\n");
1480:   PetscFPrintf(comm, fd, "      ##########################################################\n");
1481:   PetscFPrintf(comm, fd, "      #                                                        #\n");
1482:   PetscFPrintf(comm, fd, "      #                          WARNING!!!                    #\n");
1483:   PetscFPrintf(comm, fd, "      #                                                        #\n");
1484:   PetscFPrintf(comm, fd, "      #   This code was compiled with a debugging option,      #\n");
1485:   PetscFPrintf(comm, fd, "      #   To get timing results run ./configure                #\n");
1486:   PetscFPrintf(comm, fd, "      #   using --with-debugging=no, the performance will      #\n");
1487:   PetscFPrintf(comm, fd, "      #   be generally two or three times faster.              #\n");
1488:   PetscFPrintf(comm, fd, "      #                                                        #\n");
1489:   PetscFPrintf(comm, fd, "      ##########################################################\n\n\n");
1490: #endif
1491: #if defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_FORTRAN_KERNELS)
1492:   PetscFPrintf(comm, fd, "\n\n");
1493:   PetscFPrintf(comm, fd, "      ##########################################################\n");
1494:   PetscFPrintf(comm, fd, "      #                                                        #\n");
1495:   PetscFPrintf(comm, fd, "      #                          WARNING!!!                    #\n");
1496:   PetscFPrintf(comm, fd, "      #                                                        #\n");
1497:   PetscFPrintf(comm, fd, "      #   The code for various complex numbers numerical       #\n");
1498:   PetscFPrintf(comm, fd, "      #   kernels uses C++, which generally is not well        #\n");
1499:   PetscFPrintf(comm, fd, "      #   optimized.  For performance that is about 4-5 times  #\n");
1500:   PetscFPrintf(comm, fd, "      #   faster, specify --with-fortran-kernels=1             #\n");
1501:   PetscFPrintf(comm, fd, "      #   when running ./configure.py.                         #\n");
1502:   PetscFPrintf(comm, fd, "      #                                                        #\n");
1503:   PetscFPrintf(comm, fd, "      ##########################################################\n\n\n");
1504: #endif

1506:   /* Report events */
1507:   PetscFPrintf(comm, fd,"Event                Count      Time (sec)     Flops                             --- Global ---  --- Stage ---   Total\n");
1508:   PetscFPrintf(comm, fd,"                   Max Ratio  Max     Ratio   Max  Ratio  Mess   Avg len Reduct  %%T %%F %%M %%L %%R  %%T %%F %%M %%L %%R Mflop/s\n");
1509:   PetscFPrintf(comm,fd,"------------------------------------------------------------------------------------------------------------------------\n");

1511:   /* Problem: The stage name will not show up unless the stage executed on proc 1 */
1512:   for (stage = 0; stage < numStages; stage++) {
1513:     if (!stageVisible[stage]) continue;
1514:     if (localStageUsed[stage]) {
1515:       PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stageInfo[stage].name);
1516:       MPI_Allreduce(&stageInfo[stage].perfInfo.time,          &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1517:       MPI_Allreduce(&stageInfo[stage].perfInfo.flops,         &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1518:       MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages,   &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1519:       MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1520:       MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1521:     } else {
1522:       PetscFPrintf(comm, fd, "\n--- Event Stage %d: Unknown\n\n", stage);
1523:       MPI_Allreduce(&zero,                           &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1524:       MPI_Allreduce(&zero,                           &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1525:       MPI_Allreduce(&zero,                           &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1526:       MPI_Allreduce(&zero,                           &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1527:       MPI_Allreduce(&zero,                           &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1528:     }
1529:     mess *= 0.5; messLen *= 0.5; red /= size;

1531:     /* Get total number of events in this stage --
1532:        Currently, a single processor can register more events than another, but events must all be registered in order,
1533:        just like stages. We can removed this requirement if necessary by having a global event numbering and indirection
1534:        on the event ID. This seems best accomplished by assoicating a communicator with each stage.

1536:        Problem: If the event did not happen on proc 1, its name will not be available.
1537:        Problem: Event visibility is not implemented
1538:     */
1539:     if (localStageUsed[stage]) {
1540:       eventInfo      = stageLog->stageInfo[stage].eventLog->eventInfo;
1541:       localNumEvents = stageLog->stageInfo[stage].eventLog->numEvents;
1542:     } else localNumEvents = 0;
1543:     MPI_Allreduce(&localNumEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm);
1544:     for (event = 0; event < numEvents; event++) {
1545:       if (localStageUsed[stage] && (event < stageLog->stageInfo[stage].eventLog->numEvents) && (eventInfo[event].depth == 0)) {
1546:         if ((eventInfo[event].count > 0) && (eventInfo[event].time > 0.0)) flopr = eventInfo[event].flops;
1547:         else flopr = 0.0;

1549:         MPI_Allreduce(&flopr,                          &minf,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1550:         MPI_Allreduce(&flopr,                          &maxf,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1551:         MPI_Allreduce(&eventInfo[event].flops,         &totf,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1552:         MPI_Allreduce(&eventInfo[event].time,          &mint,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1553:         MPI_Allreduce(&eventInfo[event].time,          &maxt,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1554:         MPI_Allreduce(&eventInfo[event].time,          &tott,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1555:         MPI_Allreduce(&eventInfo[event].numMessages,   &totm,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1556:         MPI_Allreduce(&eventInfo[event].messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1557:         MPI_Allreduce(&eventInfo[event].numReductions, &totr,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1558:         MPI_Allreduce(&eventInfo[event].count,         &minCt, 1, MPI_INT,             MPI_MIN, comm);
1559:         MPI_Allreduce(&eventInfo[event].count,         &maxCt, 1, MPI_INT,             MPI_MAX, comm);
1560:         name = stageLog->eventLog->eventInfo[event].name;
1561:       } else {
1562:         flopr = 0.0;
1563:         MPI_Allreduce(&flopr,                          &minf,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1564:         MPI_Allreduce(&flopr,                          &maxf,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1565:         MPI_Allreduce(&zero,                           &totf,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1566:         MPI_Allreduce(&zero,                           &mint,  1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm);
1567:         MPI_Allreduce(&zero,                           &maxt,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm);
1568:         MPI_Allreduce(&zero,                           &tott,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1569:         MPI_Allreduce(&zero,                           &totm,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1570:         MPI_Allreduce(&zero,                           &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1571:         MPI_Allreduce(&zero,                           &totr,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
1572:         MPI_Allreduce(&ierr,                           &minCt, 1, MPI_INT,             MPI_MIN, comm);
1573:         MPI_Allreduce(&ierr,                           &maxCt, 1, MPI_INT,             MPI_MAX, comm);
1574:         name  = "";
1575:       }
1576:       if (mint < 0.0) {
1577:         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);
1578:         mint = 0;
1579:       }
1580:       if (minf < 0.0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Minimum flops %g over all processors for %s is negative! Not possible!",minf,name);
1581:       totm *= 0.5; totml *= 0.5; totr /= size;

1583:       if (maxCt != 0) {
1584:         if (minCt         != 0)   ratCt            = ((PetscLogDouble) maxCt)/minCt; else ratCt            = 0.0;
1585:         if (mint          != 0.0) ratt             = maxt/mint;                  else ratt             = 0.0;
1586:         if (minf          != 0.0) ratf             = maxf/minf;                  else ratf             = 0.0;
1587:         if (TotalTime     != 0.0) fracTime         = tott/TotalTime;             else fracTime         = 0.0;
1588:         if (TotalFlops    != 0.0) fracFlops        = totf/TotalFlops;            else fracFlops        = 0.0;
1589:         if (stageTime     != 0.0) fracStageTime    = tott/stageTime;             else fracStageTime    = 0.0;
1590:         if (flops         != 0.0) fracStageFlops   = totf/flops;                 else fracStageFlops   = 0.0;
1591:         if (numMessages   != 0.0) fracMess         = totm/numMessages;           else fracMess         = 0.0;
1592:         if (messageLength != 0.0) fracMessLen      = totml/messageLength;        else fracMessLen      = 0.0;
1593:         if (numReductions != 0.0) fracRed          = totr/numReductions;         else fracRed          = 0.0;
1594:         if (mess          != 0.0) fracStageMess    = totm/mess;                  else fracStageMess    = 0.0;
1595:         if (messLen       != 0.0) fracStageMessLen = totml/messLen;              else fracStageMessLen = 0.0;
1596:         if (red           != 0.0) fracStageRed     = totr/red;                   else fracStageRed     = 0.0;
1597:         if (totm          != 0.0) totml           /= totm;                       else totml            = 0.0;
1598:         if (maxt          != 0.0) flopr            = totf/maxt;                  else flopr            = 0.0;
1599:         if (fracStageTime > 1.00)  PetscFPrintf(comm, fd,"Warning -- total time of even greater than time of entire stage -- something is wrong with the timer\n");
1600:         PetscFPrintf(comm, fd,
1601:           "%-16s %7d%4.1f %5.4e%4.1f %3.2e%4.1f %2.1e %2.1e %2.1e%3.0f%3.0f%3.0f%3.0f%3.0f %3.0f%3.0f%3.0f%3.0f%3.0f %5.0f\n",
1602:                             name, maxCt, ratCt, maxt, ratt, maxf, ratf, totm, totml, totr,
1603:                             100.0*fracTime, 100.0*fracFlops, 100.0*fracMess, 100.0*fracMessLen, 100.0*fracRed,
1604:                             100.0*fracStageTime, 100.0*fracStageFlops, 100.0*fracStageMess, 100.0*fracStageMessLen, 100.0*fracStageRed,
1605:                             PetscAbsReal(flopr/1.0e6));
1606:       }
1607:     }
1608:   }

1610:   /* Memory usage and object creation */
1611:   PetscFPrintf(comm, fd, "------------------------------------------------------------------------------------------------------------------------\n");
1612:   PetscFPrintf(comm, fd, "\n");
1613:   PetscFPrintf(comm, fd, "Memory usage is given in bytes:\n\n");

1615:   /* Right now, only stages on the first processor are reported here, meaning only objects associated with
1616:      the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then
1617:      stats for stages local to processor sets.
1618:   */
1619:   /* We should figure out the longest object name here (now 20 characters) */
1620:   PetscFPrintf(comm, fd, "Object Type          Creations   Destructions     Memory  Descendants' Mem.\n");
1621:   PetscFPrintf(comm, fd, "Reports information only for process 0.\n");
1622:   for (stage = 0; stage < numStages; stage++) {
1623:     if (localStageUsed[stage]) {
1624:       classInfo = stageLog->stageInfo[stage].classLog->classInfo;
1625:       PetscFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stageInfo[stage].name);
1626:       for (oclass = 0; oclass < stageLog->stageInfo[stage].classLog->numClasses; oclass++) {
1627:         if ((classInfo[oclass].creations > 0) || (classInfo[oclass].destructions > 0)) {
1628:           PetscFPrintf(comm, fd, "%20s %5d          %5d  %11.0f     %g\n", stageLog->classLog->classInfo[oclass].name,
1629:                               classInfo[oclass].creations, classInfo[oclass].destructions, classInfo[oclass].mem,
1630:                               classInfo[oclass].descMem);
1631:         }
1632:       }
1633:     } else {
1634:       PetscFPrintf(comm, fd, "\n--- Event Stage %d: Unknown\n\n", stage);
1635:     }
1636:   }

1638:   PetscFree(localStageUsed);
1639:   PetscFree(stageUsed);
1640:   PetscFree(localStageVisible);
1641:   PetscFree(stageVisible);

1643:   /* Information unrelated to this particular run */
1644:   PetscFPrintf(comm, fd, "========================================================================================================================\n");
1645:   PetscTime(&y);
1646:   PetscTime(&x);
1647:   PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y);
1648:   PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y);
1649:   PetscFPrintf(comm,fd,"Average time to get PetscTime(): %g\n", (y-x)/10.0);
1650:   /* MPI information */
1651:   if (size > 1) {
1652:     MPI_Status  status;
1653:     PetscMPIInt tag;
1654:     MPI_Comm    newcomm;

1656:     MPI_Barrier(comm);
1657:     PetscTime(&x);
1658:     MPI_Barrier(comm);
1659:     MPI_Barrier(comm);
1660:     MPI_Barrier(comm);
1661:     MPI_Barrier(comm);
1662:     MPI_Barrier(comm);
1663:     PetscTime(&y);
1664:     PetscFPrintf(comm, fd, "Average time for MPI_Barrier(): %g\n", (y-x)/5.0);
1665:     PetscCommDuplicate(comm,&newcomm, &tag);
1666:     MPI_Barrier(comm);
1667:     if (rank) {
1668:       MPI_Recv(0, 0, MPI_INT, rank-1,            tag, newcomm, &status);
1669:       MPI_Send(0, 0, MPI_INT, (rank+1)%size, tag, newcomm);
1670:     } else {
1671:       PetscTime(&x);
1672:       MPI_Send(0, 0, MPI_INT, 1,          tag, newcomm);
1673:       MPI_Recv(0, 0, MPI_INT, size-1, tag, newcomm, &status);
1674:       PetscTime(&y);
1675:       PetscFPrintf(comm,fd,"Average time for zero size MPI_Send(): %g\n", (y-x)/size);
1676:     }
1677:     PetscCommDestroy(&newcomm);
1678:   }
1679:   PetscOptionsView(viewer);

1681:   /* Machine and compile information */
1682: #if defined(PETSC_USE_FORTRAN_KERNELS)
1683:   PetscFPrintf(comm, fd, "Compiled with FORTRAN kernels\n");
1684: #else
1685:   PetscFPrintf(comm, fd, "Compiled without FORTRAN kernels\n");
1686: #endif
1687: #if defined(PETSC_USE_REAL_SINGLE)
1688:   PetscFPrintf(comm, fd, "Compiled with single precision PetscScalar and PetscReal\n");
1689: #elif defined(PETSC_USE_LONGDOUBLE)
1690:   PetscFPrintf(comm, fd, "Compiled with long double precision PetscScalar and PetscReal\n");
1691: #endif

1693: #if defined(PETSC_USE_REAL_MAT_SINGLE)
1694:   PetscFPrintf(comm, fd, "Compiled with single precision matrices\n");
1695: #else
1696:   PetscFPrintf(comm, fd, "Compiled with full precision matrices (default)\n");
1697: #endif
1698:   PetscFPrintf(comm, fd, "sizeof(short) %d sizeof(int) %d sizeof(long) %d sizeof(void*) %d sizeof(PetscScalar) %d sizeof(PetscInt) %d\n",
1699:                       (int) sizeof(short), (int) sizeof(int), (int) sizeof(long), (int) sizeof(void*),(int) sizeof(PetscScalar),(int) sizeof(PetscInt));

1701:   PetscFPrintf(comm, fd, "Configure run at: %s\n",petscconfigureruntime);
1702:   PetscFPrintf(comm, fd, "Configure options: %s",petscconfigureoptions);
1703:   PetscFPrintf(comm, fd, "%s", petscmachineinfo);
1704:   PetscFPrintf(comm, fd, "%s", petsccompilerinfo);
1705:   PetscFPrintf(comm, fd, "%s", petsccompilerflagsinfo);
1706:   PetscFPrintf(comm, fd, "%s", petsclinkerinfo);

1708:   /* Cleanup */
1709:   PetscFPrintf(comm, fd, "\n");
1710:   PetscStageLogPush(stageLog, lastStage);
1711:   return(0);
1712: }

1716: /*@C
1717:   PetscLogPrintDetailed - Each process prints the times for its own events

1719:   Collective over MPI_Comm

1721:   Input Parameter:
1722: + comm - The MPI communicator (only one processor prints output)
1723: - file - [Optional] The output file name

1725:   Options Database Keys:
1726: . -log_summary_detailed - Prints summary of log information (for code compiled with PETSC_USE_LOG)

1728:   Usage:
1729: .vb
1730:      PetscInitialize(...);
1731:      PetscLogBegin();
1732:      ... code ...
1733:      PetscLogPrintDetailed(MPI_Comm,filename);
1734:      PetscFinalize(...);
1735: .ve

1737:   Notes:
1738:   By default the summary is printed to stdout.

1740:   Level: beginner

1742: .keywords: log, dump, print
1743: .seealso: PetscLogBegin(), PetscLogDump(), PetscLogView()
1744: @*/
1745: PetscErrorCode  PetscLogPrintDetailed(MPI_Comm comm, const char filename[])
1746: {
1747:   FILE               *fd = PETSC_STDOUT;
1748:   PetscStageLog      stageLog;
1749:   PetscStageInfo     *stageInfo = NULL;
1750:   PetscEventPerfInfo *eventInfo = NULL;
1751:   const char         *name      = NULL;
1752:   PetscLogDouble     TotalTime;
1753:   PetscLogDouble     stageTime, flops, flopr, mess, messLen, red;
1754:   PetscLogDouble     maxf, totf, maxt, tott, totm, totml, totr = 0.0;
1755:   PetscMPIInt        maxCt;
1756:   PetscBool          *stageUsed;
1757:   PetscBool          *stageVisible;
1758:   int                numStages, numEvents;
1759:   int                stage;
1760:   PetscLogEvent      event;
1761:   PetscErrorCode     ierr;

1764:   if (!PetscLogBegin_PrivateCalled) SETERRQ(comm, PETSC_ERR_ORDER, "No call to PetscLogBegin() before PetscLogPrintDetailed()");
1765:   /* Pop off any stages the user forgot to remove */
1766:   PetscLogGetStageLog(&stageLog);
1767:   PetscStageLogGetCurrent(stageLog, &stage);
1768:   while (stage >= 0) {
1769:     PetscStageLogPop(stageLog);
1770:     PetscStageLogGetCurrent(stageLog, &stage);
1771:   }
1772:   /* Get the total elapsed time */
1773:   PetscTime(&TotalTime);  TotalTime -= petsc_BaseTime;
1774:   /* Open the summary file */
1775:   if (filename) {
1776:     PetscFOpen(comm, filename, "w", &fd);
1777:   }

1779:   PetscFPrintf(comm, fd, "************************************************************************************************************************\n");
1780:   PetscFPrintf(comm, fd, "***             WIDEN YOUR WINDOW TO 120 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document            ***\n");
1781:   PetscFPrintf(comm, fd, "************************************************************************************************************************\n");


1784:   numStages = stageLog->numStages;
1785:   PetscMalloc(numStages * sizeof(PetscBool), &stageUsed);
1786:   PetscMalloc(numStages * sizeof(PetscBool), &stageVisible);
1787:   if (numStages > 0) {
1788:     stageInfo = stageLog->stageInfo;
1789:     for (stage = 0; stage < numStages; stage++) {
1790:       if (stage < stageLog->numStages) {
1791:         stageUsed[stage]    = stageInfo[stage].used;
1792:         stageVisible[stage] = stageInfo[stage].perfInfo.visible;
1793:       } else {
1794:         stageUsed[stage]    = PETSC_FALSE;
1795:         stageVisible[stage] = PETSC_TRUE;
1796:       }
1797:     }
1798:   }

1800:   /* Report events */
1801:   PetscFPrintf(comm, fd,"Event                Count      Time (sec)     Flops/sec                          \n");
1802:   PetscFPrintf(comm, fd,"                                                            Mess   Avg len Reduct \n");
1803:   PetscFPrintf(comm,fd,"-----------------------------------------------------------------------------------\n");
1804:   /* Problem: The stage name will not show up unless the stage executed on proc 1 */
1805:   for (stage = 0; stage < numStages; stage++) {
1806:     if (!stageVisible[stage]) continue;
1807:     if (stageUsed[stage]) {
1808:       PetscSynchronizedFPrintf(comm, fd, "\n--- Event Stage %d: %s\n\n", stage, stageInfo[stage].name);
1809:       MPI_Allreduce(&stageInfo[stage].perfInfo.time,          &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_SELF);
1810:       MPI_Allreduce(&stageInfo[stage].perfInfo.flops,         &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_SELF);
1811:       MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages,   &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_SELF);
1812:       MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_SELF);
1813:       MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_SELF);
1814:     }
1815:     mess *= 0.5; messLen *= 0.5;

1817:     /* Get total number of events in this stage --
1818:     */
1819:     if (stageUsed[stage]) {
1820:       eventInfo = stageLog->stageInfo[stage].eventLog->eventInfo;
1821:       numEvents = stageLog->stageInfo[stage].eventLog->numEvents;
1822:     } else numEvents = 0;

1824:     for (event = 0; event < numEvents; event++) {
1825:       if (stageUsed[stage] && (event < stageLog->stageInfo[stage].eventLog->numEvents)) {
1826:         if ((eventInfo[event].count > 0) && (eventInfo[event].time > 0.0)) flopr = eventInfo[event].flops/eventInfo[event].time;
1827:         else flopr = 0.0;

1829:         MPI_Allreduce(&flopr,                          &maxf,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, PETSC_COMM_SELF);
1830:         MPI_Allreduce(&eventInfo[event].flops,         &totf,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_SELF);
1831:         MPI_Allreduce(&eventInfo[event].time,          &maxt,  1, MPIU_PETSCLOGDOUBLE, MPI_MAX, PETSC_COMM_SELF);
1832:         MPI_Allreduce(&eventInfo[event].time,          &tott,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_SELF);
1833:         MPI_Allreduce(&eventInfo[event].numMessages,   &totm,  1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_SELF);
1834:         MPI_Allreduce(&eventInfo[event].messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_SELF);
1835:         totr  = eventInfo[event].numReductions;
1836:         MPI_Allreduce(&eventInfo[event].count,         &maxCt, 1, MPI_INT,             MPI_MAX, PETSC_COMM_SELF);
1837:         name  = stageLog->eventLog->eventInfo[event].name;
1838:         totm *= 0.5; totml *= 0.5;
1839:       }

1841:       if (maxCt != 0) {
1842:         if (totm != 0.0) totml /= totm;
1843:         else totml = 0.0;
1844:         PetscSynchronizedFPrintf(comm, fd,"%-16s %7d      %5.4e      %3.2e      %2.1e %2.1e %2.1e\n",name, maxCt,  maxt,  maxf, totm, totml, totr);
1845:       }
1846:     }
1847:   }
1848:   PetscSynchronizedFlush(comm);

1850:   PetscFree(stageUsed);
1851:   PetscFree(stageVisible);

1853:   PetscFClose(comm, fd);
1854:   return(0);
1855: }

1857: /*----------------------------------------------- Counter Functions -------------------------------------------------*/
1860: /*@C
1861:    PetscGetFlops - Returns the number of flops used on this processor
1862:    since the program began.

1864:    Not Collective

1866:    Output Parameter:
1867:    flops - number of floating point operations

1869:    Notes:
1870:    A global counter logs all PETSc flop counts.  The user can use
1871:    PetscLogFlops() to increment this counter to include flops for the
1872:    application code.

1874:    PETSc automatically logs library events if the code has been
1875:    compiled with -DPETSC_USE_LOG (which is the default), and -log,
1876:    -log_summary, or -log_all are specified.  PetscLogFlops() is
1877:    intended for logging user flops to supplement this PETSc
1878:    information.

1880:    Level: intermediate

1882: .keywords: log, flops, floating point operations

1884: .seealso: PetscTime(), PetscLogFlops()
1885: @*/
1886: PetscErrorCode  PetscGetFlops(PetscLogDouble *flops)
1887: {
1889:   *flops = petsc_TotalFlops;
1890:   return(0);
1891: }

1895: PetscErrorCode  PetscLogObjectState(PetscObject obj, const char format[], ...)
1896: {
1898:   size_t         fullLength;
1899:   va_list        Argp;

1902:   if (!petsc_logObjects) return(0);
1903:   va_start(Argp, format);
1904:   PetscVSNPrintf(petsc_objects[obj->id].info, 64,format,&fullLength, Argp);
1905:   va_end(Argp);
1906:   return(0);
1907: }


1910: /*MC
1911:    PetscLogFlops - Adds floating point operations to the global counter.

1913:    Synopsis:
1914:    #include "petsclog.h"
1915:    PetscErrorCode PetscLogFlops(PetscLogDouble f)

1917:    Not Collective

1919:    Input Parameter:
1920: .  f - flop counter


1923:    Usage:
1924: .vb
1925:      PetscLogEvent USER_EVENT;
1926:      PetscLogEventRegister("User event",0,&USER_EVENT);
1927:      PetscLogEventBegin(USER_EVENT,0,0,0,0);
1928:         [code segment to monitor]
1929:         PetscLogFlops(user_flops)
1930:      PetscLogEventEnd(USER_EVENT,0,0,0,0);
1931: .ve

1933:    Notes:
1934:    A global counter logs all PETSc flop counts.  The user can use
1935:    PetscLogFlops() to increment this counter to include flops for the
1936:    application code.

1938:    PETSc automatically logs library events if the code has been
1939:    compiled with -DPETSC_USE_LOG (which is the default), and -log,
1940:    -log_summary, or -log_all are specified.  PetscLogFlops() is
1941:    intended for logging user flops to supplement this PETSc
1942:    information.

1944:    Level: intermediate

1946: .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscGetFlops()

1948: .keywords: log, flops, floating point operations
1949: M*/

1951: /*MC
1952:    PetscPreLoadBegin - Begin a segment of code that may be preloaded (run twice)
1953:     to get accurate timings

1955:    Synopsis:
1956:    #include "petsclog.h"
1957:    void PetscPreLoadBegin(PetscBool  flag,char *name);

1959:    Not Collective

1961:    Input Parameter:
1962: +   flag - PETSC_TRUE to run twice, PETSC_FALSE to run once, may be overridden
1963:            with command line option -preload true or -preload false
1964: -   name - name of first stage (lines of code timed separately with -log_summary) to
1965:            be preloaded

1967:    Usage:
1968: .vb
1969:      PetscPreLoadBegin(PETSC_TRUE,"first stage);
1970:        lines of code
1971:        PetscPreLoadStage("second stage");
1972:        lines of code
1973:      PetscPreLoadEnd();
1974: .ve

1976:    Notes: Only works in C/C++, not Fortran

1978:      Flags available within the macro.
1979: +    PetscPreLoadingUsed - true if we are or have done preloading
1980: .    PetscPreLoadingOn - true if it is CURRENTLY doing preload
1981: .    PetscPreLoadIt - 0 for the first computation (with preloading turned off it is only 0) 1 for the second
1982: -    PetscPreLoadMax - number of times it will do the computation, only one when preloading is turned on
1983:      The first two variables are available throughout the program, the second two only between the PetscPreLoadBegin()
1984:      and PetscPreLoadEnd()

1986:    Level: intermediate

1988: .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadEnd(), PetscPreLoadStage()

1990:    Concepts: preloading
1991:    Concepts: timing^accurate
1992:    Concepts: paging^eliminating effects of


1995: M*/

1997: /*MC
1998:    PetscPreLoadEnd - End a segment of code that may be preloaded (run twice)
1999:     to get accurate timings

2001:    Synopsis:
2002:    #include "petsclog.h"
2003:    void PetscPreLoadEnd(void);

2005:    Not Collective

2007:    Usage:
2008: .vb
2009:      PetscPreLoadBegin(PETSC_TRUE,"first stage);
2010:        lines of code
2011:        PetscPreLoadStage("second stage");
2012:        lines of code
2013:      PetscPreLoadEnd();
2014: .ve

2016:    Notes: only works in C/C++ not fortran

2018:    Level: intermediate

2020: .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadBegin(), PetscPreLoadStage()

2022: M*/

2024: /*MC
2025:    PetscPreLoadStage - Start a new segment of code to be timed separately.
2026:     to get accurate timings

2028:    Synopsis:
2029:    #include "petsclog.h"
2030:    void PetscPreLoadStage(char *name);

2032:    Not Collective

2034:    Usage:
2035: .vb
2036:      PetscPreLoadBegin(PETSC_TRUE,"first stage);
2037:        lines of code
2038:        PetscPreLoadStage("second stage");
2039:        lines of code
2040:      PetscPreLoadEnd();
2041: .ve

2043:    Notes: only works in C/C++ not fortran

2045:    Level: intermediate

2047: .seealso: PetscLogEventRegister(), PetscLogEventBegin(), PetscLogEventEnd(), PetscPreLoadBegin(), PetscPreLoadEnd()

2049: M*/

2053: /*@
2054:    PetscLogViewPython - Saves logging information in a Python format.

2056:    Collective on PetscViewer

2058:    Input Paramter:
2059: .   viewer - viewer to save Python data

2061:   Level: intermediate

2063: @*/
2064: PetscErrorCode  PetscLogViewPython(PetscViewer viewer)
2065: {
2066:   FILE               *fd;
2067:   PetscLogDouble     zero                    = 0.0;
2068:   PetscStageLog      stageLog;
2069:   PetscStageInfo     *stageInfo              = NULL;
2070:   PetscEventPerfInfo *eventInfo              = NULL;
2071:   const char         *name;
2072:   char               stageName[2048];
2073:   char               eventName[2048];
2074:   PetscLogDouble     locTotalTime, TotalTime = 0, TotalFlops = 0;
2075:   PetscLogDouble     numMessages             = 0, messageLength = 0, avgMessLen, numReductions = 0;
2076:   PetscLogDouble     stageTime, flops, mem, mess, messLen, red;
2077:   PetscLogDouble     fracTime, fracFlops, fracMessages, fracLength;
2078:   PetscLogDouble     fracReductions;
2079:   PetscLogDouble     tot,avg,x,y,*mydata;
2080:   PetscMPIInt        maxCt;
2081:   PetscMPIInt        size, rank, *mycount;
2082:   PetscBool          *localStageUsed,    *stageUsed;
2083:   PetscBool          *localStageVisible, *stageVisible;
2084:   int                numStages, localNumEvents, numEvents;
2085:   int                stage, lastStage;
2086:   PetscLogEvent      event;
2087:   PetscErrorCode     ierr;
2088:   PetscInt           i;
2089:   MPI_Comm           comm;

2092:   PetscObjectGetComm((PetscObject)viewer,&comm);
2093:   if (!PetscLogBegin_PrivateCalled) SETERRQ(comm, PETSC_ERR_ORDER, "No call to PetscLogBegin() before PetscLogViewPython()");
2094:   MPI_Comm_size(comm, &size);
2095:   MPI_Comm_rank(comm, &rank);
2096:   PetscMalloc(size*sizeof(PetscLogDouble), &mydata);
2097:   PetscMalloc(size*sizeof(PetscMPIInt), &mycount);
2098:   PetscViewerASCIIGetPointer(viewer,&fd);

2100:   /* Pop off any stages the user forgot to remove */
2101:   lastStage = 0;
2102:   PetscLogGetStageLog(&stageLog);
2103:   PetscStageLogGetCurrent(stageLog, &stage);
2104:   while (stage >= 0) {
2105:     lastStage = stage;
2106:     PetscStageLogPop(stageLog);
2107:     PetscStageLogGetCurrent(stageLog, &stage);
2108:   }
2109:   /* Get the total elapsed time */
2110:   PetscTime(&locTotalTime);  locTotalTime -= petsc_BaseTime;

2112:   PetscFPrintf(comm, fd, "\n#------ PETSc Performance Summary ----------\n\n");
2113:   PetscFPrintf(comm, fd, "Nproc = %d\n",size);

2115:   /* Must preserve reduction count before we go on */
2116:   red = (petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct)/((PetscLogDouble) size);

2118:   /* Calculate summary information */

2120:   /*   Time */
2121:   MPI_Gather(&locTotalTime,1,MPIU_PETSCLOGDOUBLE,mydata,1,MPIU_PETSCLOGDOUBLE,0,comm);
2122:   if (!rank) {
2123:     PetscFPrintf(comm, fd, "Time = [ ");
2124:     tot  = 0.0;
2125:     for (i=0; i<size; i++) {
2126:       tot += mydata[i];
2127:       PetscFPrintf(comm, fd, "  %5.3e,",mydata[i]);
2128:     }
2129:     PetscFPrintf(comm, fd, "]\n");
2130:     avg       = (tot)/((PetscLogDouble) size);
2131:     TotalTime = tot;
2132:   }

2134:   /*   Objects */
2135:   avg  = (PetscLogDouble) petsc_numObjects;
2136:   MPI_Gather(&avg,1,MPIU_PETSCLOGDOUBLE,mydata,1,MPIU_PETSCLOGDOUBLE,0,comm);
2137:   if (!rank) {
2138:     PetscFPrintf(comm, fd, "Objects = [ ");
2139:     for (i=0; i<size; i++) {
2140:       PetscFPrintf(comm, fd, "  %5.3e,",mydata[i]);
2141:     }
2142:     PetscFPrintf(comm, fd, "]\n");
2143:   }

2145:   /*   Flops */
2146:   MPI_Gather(&petsc_TotalFlops,1,MPIU_PETSCLOGDOUBLE,mydata,1,MPIU_PETSCLOGDOUBLE,0,comm);
2147:   if (!rank) {
2148:     PetscFPrintf(comm, fd, "Flops = [ ");
2149:     tot  = 0.0;
2150:     for (i=0; i<size; i++) {
2151:       tot += mydata[i];
2152:       PetscFPrintf(comm, fd, "  %5.3e,",mydata[i]);
2153:     }
2154:     PetscFPrintf(comm, fd, "]\n");
2155:     TotalFlops = tot;
2156:   }

2158:   /*   Memory */
2159:   PetscMallocGetMaximumUsage(&mem);
2160:   MPI_Gather(&mem,1,MPIU_PETSCLOGDOUBLE,mydata,1,MPIU_PETSCLOGDOUBLE,0,comm);
2161:   if (!rank) {
2162:     PetscFPrintf(comm, fd, "Memory = [ ");
2163:     for (i=0; i<size; i++) {
2164:       PetscFPrintf(comm, fd, "  %5.3e,",mydata[i]);
2165:     }
2166:     PetscFPrintf(comm, fd, "]\n");
2167:   }

2169:   /*   Messages */
2170:   mess = 0.5*(petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct);
2171:   MPI_Gather(&mess,1,MPIU_PETSCLOGDOUBLE,mydata,1,MPIU_PETSCLOGDOUBLE,0,comm);
2172:   if (!rank) {
2173:     PetscFPrintf(comm, fd, "MPIMessages = [ ");
2174:     tot  = 0.0;
2175:     for (i=0; i<size; i++) {
2176:       tot += mydata[i];
2177:       PetscFPrintf(comm, fd, "  %5.3e,",mydata[i]);
2178:     }
2179:     PetscFPrintf(comm, fd, "]\n");
2180:     numMessages = tot;
2181:   }

2183:   /*   Message Lengths */
2184:   mess = 0.5*(petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len);
2185:   MPI_Gather(&mess,1,MPIU_PETSCLOGDOUBLE,mydata,1,MPIU_PETSCLOGDOUBLE,0,comm);
2186:   if (!rank) {
2187:     PetscFPrintf(comm, fd, "MPIMessageLengths = [ ");
2188:     tot  = 0.0;
2189:     for (i=0; i<size; i++) {
2190:       tot += mydata[i];
2191:       PetscFPrintf(comm, fd, "  %5.3e,",mydata[i]);
2192:     }
2193:     PetscFPrintf(comm, fd, "]\n");
2194:     messageLength = tot;
2195:   }

2197:   /*   Reductions */
2198:   MPI_Gather(&red,1,MPIU_PETSCLOGDOUBLE,mydata,1,MPIU_PETSCLOGDOUBLE,0,comm);
2199:   if (!rank) {
2200:     PetscFPrintf(comm, fd, "MPIReductions = [ ");
2201:     tot  = 0.0;
2202:     for (i=0; i<size; i++) {
2203:       tot += mydata[i];
2204:       PetscFPrintf(comm, fd, "  %5.3e,",mydata[i]);
2205:     }
2206:     PetscFPrintf(comm, fd, "]\n");
2207:     numReductions = tot;
2208:   }

2210:   /* Get total number of stages --
2211:        Currently, a single processor can register more stages than another, but stages must all be registered in order.
2212:        We can removed this requirement if necessary by having a global stage numbering and indirection on the stage ID.
2213:        This seems best accomplished by assoicating a communicator with each stage.
2214:   */
2215:   MPI_Allreduce(&stageLog->numStages, &numStages, 1, MPI_INT, MPI_MAX, comm);
2216:   PetscMalloc(numStages * sizeof(PetscBool), &localStageUsed);
2217:   PetscMalloc(numStages * sizeof(PetscBool), &stageUsed);
2218:   PetscMalloc(numStages * sizeof(PetscBool), &localStageVisible);
2219:   PetscMalloc(numStages * sizeof(PetscBool), &stageVisible);
2220:   if (numStages > 0) {
2221:     stageInfo = stageLog->stageInfo;
2222:     for (stage = 0; stage < numStages; stage++) {
2223:       if (stage < stageLog->numStages) {
2224:         localStageUsed[stage]    = stageInfo[stage].used;
2225:         localStageVisible[stage] = stageInfo[stage].perfInfo.visible;
2226:       } else {
2227:         localStageUsed[stage]    = PETSC_FALSE;
2228:         localStageVisible[stage] = PETSC_TRUE;
2229:       }
2230:     }
2231:     MPI_Allreduce(localStageUsed,    stageUsed,    numStages, MPIU_BOOL, MPI_LOR,  comm);
2232:     MPI_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm);
2233:     for (stage = 0; stage < numStages; stage++) {
2234:       if (stageUsed[stage]) {
2235:         PetscFPrintf(comm, fd, "\n#Summary of Stages:   ----- Time ------  ----- Flops -----  --- Messages ---  -- Message Lengths --  -- Reductions --\n");
2236:         PetscFPrintf(comm, fd, "#                       Avg     %%Total     Avg     %%Total   counts   %%Total     Avg         %%Total   counts   %%Total \n");
2237:         break;
2238:       }
2239:     }
2240:     for (stage = 0; stage < numStages; stage++) {
2241:       if (!stageUsed[stage]) continue;
2242:       if (localStageUsed[stage]) {
2243:         MPI_Allreduce(&stageInfo[stage].perfInfo.time,          &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2244:         MPI_Allreduce(&stageInfo[stage].perfInfo.flops,         &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2245:         MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages,   &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2246:         MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2247:         MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2248:         name = stageInfo[stage].name;
2249:       } else {
2250:         MPI_Allreduce(&zero,                           &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2251:         MPI_Allreduce(&zero,                           &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2252:         MPI_Allreduce(&zero,                           &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2253:         MPI_Allreduce(&zero,                           &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2254:         MPI_Allreduce(&zero,                           &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2255:         name = "";
2256:       }
2257:       mess *= 0.5; messLen *= 0.5; red /= size;
2258:       if (TotalTime     != 0.0) fracTime       = stageTime/TotalTime;    else fracTime       = 0.0;
2259:       if (TotalFlops    != 0.0) fracFlops      = flops/TotalFlops;       else fracFlops      = 0.0;
2260:       /* Talk to Barry if (stageTime     != 0.0) flops          = (size*flops)/stageTime; else flops          = 0.0; */
2261:       if (numMessages   != 0.0) fracMessages   = mess/numMessages;       else fracMessages   = 0.0;
2262:       if (numMessages   != 0.0) avgMessLen     = messLen/numMessages;    else avgMessLen     = 0.0;
2263:       if (messageLength != 0.0) fracLength     = messLen/messageLength;  else fracLength     = 0.0;
2264:       if (numReductions != 0.0) fracReductions = red/numReductions;      else fracReductions = 0.0;
2265:       PetscFPrintf(comm, fd, "# ");
2266:       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",
2267:                           stage, name, stageTime/size, 100.0*fracTime, flops, 100.0*fracFlops,
2268:                           mess, 100.0*fracMessages, avgMessLen, 100.0*fracLength, red, 100.0*fracReductions);
2269:     }
2270:   }

2272:   /* Report events */
2273:   PetscFPrintf(comm,fd,"\n# Event\n");
2274:   PetscFPrintf(comm,fd,"# ------------------------------------------------------\n");
2275:   PetscFPrintf(comm,fd,"class Stage(object):\n");
2276:   PetscFPrintf(comm,fd,"    def __init__(self, name, time, flops, numMessages, messageLength, numReductions):\n");
2277:   PetscFPrintf(comm,fd,"        # The time and flops represent totals across processes, whereas reductions are only counted once\n");
2278:   PetscFPrintf(comm,fd,"        self.name          = name\n");
2279:   PetscFPrintf(comm,fd,"        self.time          = time\n");
2280:   PetscFPrintf(comm,fd,"        self.flops         = flops\n");
2281:   PetscFPrintf(comm,fd,"        self.numMessages   = numMessages\n");
2282:   PetscFPrintf(comm,fd,"        self.messageLength = messageLength\n");
2283:   PetscFPrintf(comm,fd,"        self.numReductions = numReductions\n");
2284:   PetscFPrintf(comm,fd,"        self.event         = {}\n");
2285:   PetscFPrintf(comm,fd, "class Dummy(object):\n");
2286:   PetscFPrintf(comm,fd, "    pass\n");
2287:   /* Problem: The stage name will not show up unless the stage executed on proc 1 */
2288:   for (stage = 0; stage < numStages; stage++) {
2289:     if (!stageVisible[stage]) continue;
2290:     if (localStageUsed[stage]) {
2291:       MPI_Allreduce(&stageInfo[stage].perfInfo.time,          &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2292:       MPI_Allreduce(&stageInfo[stage].perfInfo.flops,         &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2293:       MPI_Allreduce(&stageInfo[stage].perfInfo.numMessages,   &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2294:       MPI_Allreduce(&stageInfo[stage].perfInfo.messageLength, &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2295:       MPI_Allreduce(&stageInfo[stage].perfInfo.numReductions, &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2296:     } else {
2297:       PetscFPrintf(comm, fd, "\n--- Event Stage %d: Unknown\n\n", stage);
2298:       MPI_Allreduce(&zero,                           &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2299:       MPI_Allreduce(&zero,                           &flops,     1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2300:       MPI_Allreduce(&zero,                           &mess,      1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2301:       MPI_Allreduce(&zero,                           &messLen,   1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2302:       MPI_Allreduce(&zero,                           &red,       1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm);
2303:     }
2304:     mess *= 0.5; messLen *= 0.5; red /= size;

2306:     /* Get total number of events in this stage --
2307:        Currently, a single processor can register more events than another, but events must all be registered in order,
2308:        just like stages. We can removed this requirement if necessary by having a global event numbering and indirection
2309:        on the event ID. This seems best accomplished by assoicating a communicator with each stage.

2311:        Problem: If the event did not happen on proc 1, its name will not be available.
2312:        Problem: Event visibility is not implemented
2313:     */

2315:     {
2316:       size_t len, c;

2318:       PetscStrcpy(stageName, stageInfo[stage].name);
2319:       PetscStrlen(stageName, &len);
2320:       for (c = 0; c < len; ++c) {
2321:         if (stageName[c] == ' ') stageName[c] = '_';
2322:       }
2323:     }
2324:     if (!rank) {
2325:       PetscFPrintf(comm, fd, "%s = Stage('%s', %g, %g, %g, %g, %g)\n", stageName, stageName, stageTime, flops, mess, messLen, red);
2326:     }

2328:     if (localStageUsed[stage]) {
2329:       eventInfo      = stageLog->stageInfo[stage].eventLog->eventInfo;
2330:       localNumEvents = stageLog->stageInfo[stage].eventLog->numEvents;
2331:     } else localNumEvents = 0;

2333:     MPI_Allreduce(&localNumEvents, &numEvents, 1, MPI_INT, MPI_MAX, comm);
2334:     for (event = 0; event < numEvents; event++) {
2335:       PetscBool      hasEvent = PETSC_TRUE;
2336:       PetscMPIInt    tmpI;
2337:       PetscLogDouble tmpR;

2339:       if (localStageUsed[stage] && (event < stageLog->stageInfo[stage].eventLog->numEvents) && (eventInfo[event].depth == 0)) {
2340:         size_t len, c;

2342:         MPI_Allreduce(&eventInfo[event].count, &maxCt, 1, MPI_INT, MPI_MAX, comm);
2343:         PetscStrcpy(eventName, stageLog->eventLog->eventInfo[event].name);
2344:         PetscStrlen(eventName, &len);
2345:         for (c = 0; c < len; ++c) {
2346:           if (eventName[c] == ' ') eventName[c] = '_';
2347:         }
2348:       } else {
2349:         MPI_Allreduce(&ierr, &maxCt, 1, MPI_INT, MPI_MAX, comm);
2350:         eventName[0] = 0;
2351:         hasEvent     = PETSC_FALSE;
2352:       }

2354:       if (maxCt != 0) {
2355:         PetscFPrintf(comm, fd,"#\n");
2356:         if (!rank) {
2357:           PetscFPrintf(comm, fd, "%s = Dummy()\n",eventName);
2358:           PetscFPrintf(comm, fd, "%s.event['%s'] = %s\n",stageName,eventName,eventName);
2359:         }
2360:         /* Count */
2361:         if (hasEvent) tmpI = eventInfo[event].count;
2362:         else          tmpI = 0;
2363:         MPI_Gather(&tmpI,1, MPI_INT, mycount, 1, MPI_INT, 0, comm);
2364:         PetscFPrintf(comm, fd, "%s.Count = [ ", eventName);
2365:         if (!rank) {
2366:           for (i=0; i<size; i++) {
2367:             PetscFPrintf(comm, fd, "  %7d,",mycount[i]);
2368:           }
2369:           PetscFPrintf(comm, fd, "]\n");
2370:         }
2371:         /* Time */
2372:         if (hasEvent) tmpR = eventInfo[event].time;
2373:         else          tmpR = 0.0;
2374:         MPI_Gather(&tmpR, 1, MPIU_PETSCLOGDOUBLE, mydata, 1, MPIU_PETSCLOGDOUBLE, 0, comm);
2375:         if (!rank) {
2376:           PetscFPrintf(comm, fd, "%s.Time  = [ ", eventName);
2377:           for (i=0; i<size; i++) {
2378:             PetscFPrintf(comm, fd, "  %5.3e,",mydata[i]);
2379:           }
2380:           PetscFPrintf(comm, fd, "]\n");
2381:         }
2382:         if (hasEvent) tmpR = eventInfo[event].time2;
2383:         else          tmpR = 0.0;
2384:         MPI_Gather(&tmpR, 1, MPIU_PETSCLOGDOUBLE, mydata, 1, MPIU_PETSCLOGDOUBLE, 0, comm);
2385:         if (!rank) {
2386:           PetscFPrintf(comm, fd, "%s.Time2 = [ ", eventName);
2387:           for (i=0; i<size; i++) {
2388:             PetscFPrintf(comm, fd, "  %5.3e,", mydata[i]);
2389:           }
2390:           PetscFPrintf(comm, fd, "]\n");
2391:         }
2392:         /* Flops */
2393:         if (hasEvent) tmpR = eventInfo[event].flops;
2394:         else          tmpR = 0.0;
2395:         MPI_Gather(&tmpR, 1, MPIU_PETSCLOGDOUBLE, mydata, 1, MPIU_PETSCLOGDOUBLE, 0, comm);
2396:         if (!rank) {
2397:           PetscFPrintf(comm, fd, "%s.Flops = [ ", eventName);
2398:           for (i=0; i<size; i++) {
2399:             PetscFPrintf(comm, fd, "  %5.3e,",mydata[i]);
2400:           }
2401:           PetscFPrintf(comm, fd, "]\n");
2402:         }
2403:         if (hasEvent) tmpR = eventInfo[event].flops2;
2404:         else          tmpR = 0.0;
2405:         MPI_Gather(&tmpR, 1, MPIU_PETSCLOGDOUBLE, mydata, 1, MPIU_PETSCLOGDOUBLE, 0, comm);
2406:         if (!rank) {
2407:           PetscFPrintf(comm, fd, "%s.Flops2 = [ ", eventName);
2408:           for (i=0; i<size; i++) {
2409:             PetscFPrintf(comm, fd, "  %5.3e,", mydata[i]);
2410:           }
2411:           PetscFPrintf(comm, fd, "]\n");
2412:         }
2413:         /* Num Messages */
2414:         if (hasEvent) tmpR = eventInfo[event].numMessages;
2415:         else          tmpR = 0.0;
2416:         MPI_Gather(&tmpR, 1, MPIU_PETSCLOGDOUBLE, mydata, 1, MPIU_PETSCLOGDOUBLE, 0, comm);
2417:         PetscFPrintf(comm, fd, "%s.NumMessages = [ ", eventName);
2418:         if (!rank) {
2419:           for (i=0; i<size; i++) {
2420:             PetscFPrintf(comm, fd, "  %7.1e,",mydata[i]);
2421:           }
2422:           PetscFPrintf(comm, fd, "]\n");
2423:         }
2424:         /* Message Length */
2425:         if (hasEvent) tmpR = eventInfo[event].messageLength;
2426:         else          tmpR = 0.0;
2427:         MPI_Gather(&tmpR, 1, MPIU_PETSCLOGDOUBLE, mydata, 1, MPIU_PETSCLOGDOUBLE, 0, comm);
2428:         if (!rank) {
2429:           PetscFPrintf(comm, fd, "%s.MessageLength = [ ", eventName);
2430:           for (i=0; i<size; i++) {
2431:             PetscFPrintf(comm, fd, "  %5.3e,",mydata[i]);
2432:           }
2433:           PetscFPrintf(comm, fd, "]\n");
2434:         }
2435:         /* Num Reductions */
2436:         if (hasEvent) tmpR = eventInfo[event].numReductions;
2437:         else          tmpR = 0.0;
2438:         MPI_Gather(&tmpR, 1, MPIU_PETSCLOGDOUBLE, mydata, 1, MPIU_PETSCLOGDOUBLE, 0, comm);
2439:         PetscFPrintf(comm, fd, "%s.NumReductions = [ ", eventName);
2440:         if (!rank) {
2441:           for (i=0; i<size; i++) {
2442:             PetscFPrintf(comm, fd, "  %7.1e,",mydata[i]);
2443:           }
2444:           PetscFPrintf(comm, fd, "]\n");
2445:         }
2446:       }
2447:     }
2448:   }

2450:   /* Right now, only stages on the first processor are reported here, meaning only objects associated with
2451:      the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then
2452:      stats for stages local to processor sets.
2453:   */
2454:   for (stage = 0; stage < numStages; stage++) {
2455:     if (!localStageUsed[stage]) {
2456:       PetscFPrintf(comm, fd, "\n--- Event Stage %d: Unknown\n\n", stage);
2457:     }
2458:   }

2460:   PetscFree(localStageUsed);
2461:   PetscFree(stageUsed);
2462:   PetscFree(localStageVisible);
2463:   PetscFree(stageVisible);
2464:   PetscFree(mydata);
2465:   PetscFree(mycount);

2467:   /* Information unrelated to this particular run */
2468:   PetscFPrintf(comm, fd, "# ========================================================================================================================\n");
2469:   PetscTime(&y);
2470:   PetscTime(&x);
2471:   PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y);
2472:   PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y); PetscTime(&y);
2473:   PetscFPrintf(comm,fd,"AveragetimetogetPetscTime = %g\n", (y-x)/10.0);
2474:   /* MPI information */
2475:   if (size > 1) {
2476:     MPI_Status  status;
2477:     PetscMPIInt tag;
2478:     MPI_Comm    newcomm;

2480:     MPI_Barrier(comm);
2481:     PetscTime(&x);
2482:     MPI_Barrier(comm);
2483:     MPI_Barrier(comm);
2484:     MPI_Barrier(comm);
2485:     MPI_Barrier(comm);
2486:     MPI_Barrier(comm);
2487:     PetscTime(&y);
2488:     PetscFPrintf(comm, fd, "AveragetimeforMPI_Barrier = %g\n", (y-x)/5.0);
2489:     PetscCommDuplicate(comm,&newcomm, &tag);
2490:     MPI_Barrier(comm);
2491:     if (rank) {
2492:       MPI_Recv(0, 0, MPI_INT, rank-1,            tag, newcomm, &status);
2493:       MPI_Send(0, 0, MPI_INT, (rank+1)%size, tag, newcomm);
2494:     } else {
2495:       PetscTime(&x);
2496:       MPI_Send(0, 0, MPI_INT, 1,          tag, newcomm);
2497:       MPI_Recv(0, 0, MPI_INT, size-1, tag, newcomm, &status);
2498:       PetscTime(&y);
2499:       PetscFPrintf(comm,fd,"AveragetimforzerosizeMPI_Send = %g\n", (y-x)/size);
2500:     }
2501:     PetscCommDestroy(&newcomm);
2502:   }

2504:   /* Cleanup */
2505:   PetscFPrintf(comm, fd, "\n");
2506:   PetscStageLogPush(stageLog, lastStage);
2507:   return(0);
2508: }

2510: #else /* end of -DPETSC_USE_LOG section */

2514: PetscErrorCode  PetscLogObjectState(PetscObject obj, const char format[], ...)
2515: {
2517:   return(0);
2518: }

2520: #endif /* PETSC_USE_LOG*/


2523: PetscClassId PETSC_LARGEST_CLASSID = PETSC_SMALLEST_CLASSID;
2524: PetscClassId PETSC_OBJECT_CLASSID  = 0;

2528: /*@C
2529:   PetscClassIdRegister - Registers a new class name for objects and logging operations in an application code.

2531:   Not Collective

2533:   Input Parameter:
2534: . name   - The class name

2536:   Output Parameter:
2537: . oclass - The class id or classid

2539:   Level: developer

2541: .keywords: log, class, register

2543: @*/
2544: PetscErrorCode  PetscClassIdRegister(const char name[],PetscClassId *oclass)
2545: {
2546: #if defined(PETSC_USE_LOG)
2547:   PetscStageLog  stageLog;
2548:   PetscInt       stage;
2550: #endif

2553:   *oclass = ++PETSC_LARGEST_CLASSID;
2554: #if defined(PETSC_USE_LOG)
2555:   PetscLogGetStageLog(&stageLog);
2556:   PetscClassRegLogRegister(stageLog->classLog, name, *oclass);
2557:   for (stage = 0; stage < stageLog->numStages; stage++) {
2558:     ClassPerfLogEnsureSize(stageLog->stageInfo[stage].classLog, stageLog->classLog->numClasses);
2559:   }
2560: #endif
2561:   return(0);
2562: }

2564: #if defined(PETSC_USE_LOG) && defined(PETSC_HAVE_MPE)
2565: #include <mpe.h>

2567: PetscBool PetscBeganMPE = PETSC_FALSE;

2569: PETSC_INTERN PetscErrorCode PetscLogEventBeginMPE(PetscLogEvent,int,PetscObject,PetscObject,PetscObject,PetscObject);
2570: PETSC_INTERN PetscErrorCode PetscLogEventEndMPE(PetscLogEvent,int,PetscObject,PetscObject,PetscObject,PetscObject);

2574: /*@C
2575:    PetscLogMPEBegin - Turns on MPE logging of events. This creates large log files
2576:    and slows the program down.

2578:    Collective over PETSC_COMM_WORLD

2580:    Options Database Keys:
2581: . -log_mpe - Prints extensive log information (for code compiled with PETSC_USE_LOG)

2583:    Notes:
2584:    A related routine is PetscLogBegin() (with the options key -log_summary), which is
2585:    intended for production runs since it logs only flop rates and object
2586:    creation (and should not significantly slow the programs).

2588:    Level: advanced

2590:    Concepts: logging^MPE
2591:    Concepts: logging^message passing

2593: .seealso: PetscLogDump(), PetscLogBegin(), PetscLogAllBegin(), PetscLogEventActivate(),
2594:           PetscLogEventDeactivate()
2595: @*/
2596: PetscErrorCode  PetscLogMPEBegin(void)
2597: {

2601:   /* Do MPE initialization */
2602:   if (!MPE_Initialized_logging()) { /* This function exists in mpich 1.1.2 and higher */
2603:     PetscInfo(0,"Initializing MPE.\n");
2604:     MPE_Init_log();

2606:     PetscBeganMPE = PETSC_TRUE;
2607:   } else {
2608:     PetscInfo(0,"MPE already initialized. Not attempting to reinitialize.\n");
2609:   }
2610:   PetscLogSet(PetscLogEventBeginMPE, PetscLogEventEndMPE);
2611:   return(0);
2612: }

2616: /*@C
2617:    PetscLogMPEDump - Dumps the MPE logging info to file for later use with Jumpshot.

2619:    Collective over PETSC_COMM_WORLD

2621:    Level: advanced

2623: .seealso: PetscLogDump(), PetscLogAllBegin(), PetscLogMPEBegin()
2624: @*/
2625: PetscErrorCode  PetscLogMPEDump(const char sname[])
2626: {
2627:   char           name[PETSC_MAX_PATH_LEN];

2631:   if (PetscBeganMPE) {
2632:     PetscInfo(0,"Finalizing MPE.\n");
2633:     if (sname) {
2634:       PetscStrcpy(name,sname);
2635:     } else {
2636:       PetscGetProgramName(name,PETSC_MAX_PATH_LEN);
2637:     }
2638:     MPE_Finish_log(name);
2639:   } else {
2640:     PetscInfo(0,"Not finalizing MPE (not started by PETSc).\n");
2641:   }
2642:   return(0);
2643: }

2645: #define PETSC_RGB_COLORS_MAX 39
2646: static const char *PetscLogMPERGBColors[PETSC_RGB_COLORS_MAX] = {
2647:   "OliveDrab:      ",
2648:   "BlueViolet:     ",
2649:   "CadetBlue:      ",
2650:   "CornflowerBlue: ",
2651:   "DarkGoldenrod:  ",
2652:   "DarkGreen:      ",
2653:   "DarkKhaki:      ",
2654:   "DarkOliveGreen: ",
2655:   "DarkOrange:     ",
2656:   "DarkOrchid:     ",
2657:   "DarkSeaGreen:   ",
2658:   "DarkSlateGray:  ",
2659:   "DarkTurquoise:  ",
2660:   "DeepPink:       ",
2661:   "DarkKhaki:      ",
2662:   "DimGray:        ",
2663:   "DodgerBlue:     ",
2664:   "GreenYellow:    ",
2665:   "HotPink:        ",
2666:   "IndianRed:      ",
2667:   "LavenderBlush:  ",
2668:   "LawnGreen:      ",
2669:   "LemonChiffon:   ",
2670:   "LightCoral:     ",
2671:   "LightCyan:      ",
2672:   "LightPink:      ",
2673:   "LightSalmon:    ",
2674:   "LightSlateGray: ",
2675:   "LightYellow:    ",
2676:   "LimeGreen:      ",
2677:   "MediumPurple:   ",
2678:   "MediumSeaGreen: ",
2679:   "MediumSlateBlue:",
2680:   "MidnightBlue:   ",
2681:   "MintCream:      ",
2682:   "MistyRose:      ",
2683:   "NavajoWhite:    ",
2684:   "NavyBlue:       ",
2685:   "OliveDrab:      "
2686: };

2690: /*@C
2691:   PetscLogMPEGetRGBColor - This routine returns a rgb color useable with PetscLogEventRegister()

2693:   Not collective. Maybe it should be?

2695:   Output Parameter
2696: . str - character string representing the color

2698:   Level: developer

2700: .keywords: log, mpe , color
2701: .seealso: PetscLogEventRegister
2702: @*/
2703: PetscErrorCode  PetscLogMPEGetRGBColor(const char *str[])
2704: {
2705:   static int idx = 0;

2708:   *str = PetscLogMPERGBColors[idx];
2709:   idx  = (idx + 1)% PETSC_RGB_COLORS_MAX;
2710:   return(0);
2711: }

2713: #endif /* PETSC_USE_LOG && PETSC_HAVE_MPE */