Actual source code: eventlog.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2: /*
  3:      This defines part of the private API for logging performance information. It is intended to be used only by the
  4:    PETSc PetscLog...() interface and not elsewhere, nor by users. Hence the prototypes for these functions are NOT
  5:    in the public PETSc include files.

  7: */
  8:  #include <petsc/private/logimpl.h>

 10: PetscBool PetscLogSyncOn = PETSC_FALSE;

 12: /*----------------------------------------------- Creation Functions -------------------------------------------------*/
 13: /* Note: these functions do not have prototypes in a public directory, so they are considered "internal" and not exported. */

 15: /*@C
 16:   PetscEventRegLogCreate - This creates a PetscEventRegLog object.

 18:   Not collective

 20:   Input Parameter:
 21: . eventLog - The PetscEventRegLog

 23:   Level: developer

 25: .keywords: log, event, create
 26: .seealso: PetscEventRegLogDestroy(), PetscStageLogCreate()
 27: @*/
 28: PetscErrorCode PetscEventRegLogCreate(PetscEventRegLog *eventLog)
 29: {
 30:   PetscEventRegLog l;
 31:   PetscErrorCode   ierr;

 34:   PetscNew(&l);
 35:   l->numEvents = 0;
 36:   l->maxEvents = 100;
 37:   PetscMalloc1(l->maxEvents,&l->eventInfo);
 38:   *eventLog    = l;
 39:   return(0);
 40: }

 42: /*@C
 43:   PetscEventRegLogDestroy - This destroys a PetscEventRegLog object.

 45:   Not collective

 47:   Input Paramter:
 48: . eventLog - The PetscEventRegLog

 50:   Level: developer

 52: .keywords: log, event, destroy
 53: .seealso: PetscEventRegLogCreate()
 54: @*/
 55: PetscErrorCode PetscEventRegLogDestroy(PetscEventRegLog eventLog)
 56: {
 57:   int            e;

 61:   for (e = 0; e < eventLog->numEvents; e++) {
 62:     PetscFree(eventLog->eventInfo[e].name);
 63:   }
 64:   PetscFree(eventLog->eventInfo);
 65:   PetscFree(eventLog);
 66:   return(0);
 67: }

 69: /*@C
 70:   PetscEventPerfLogCreate - This creates a PetscEventPerfLog object.

 72:   Not collective

 74:   Input Parameter:
 75: . eventLog - The PetscEventPerfLog

 77:   Level: developer

 79: .keywords: log, event, create
 80: .seealso: PetscEventPerfLogDestroy(), PetscStageLogCreate()
 81: @*/
 82: PetscErrorCode PetscEventPerfLogCreate(PetscEventPerfLog *eventLog)
 83: {
 84:   PetscEventPerfLog l;
 85:   PetscErrorCode    ierr;

 88:   PetscNew(&l);
 89:   l->numEvents = 0;
 90:   l->maxEvents = 100;
 91:   PetscMalloc1(l->maxEvents,&l->eventInfo);
 92:   *eventLog    = l;
 93:   return(0);
 94: }

 96: /*@C
 97:   PetscEventPerfLogDestroy - This destroys a PetscEventPerfLog object.

 99:   Not collective

101:   Input Paramter:
102: . eventLog - The PetscEventPerfLog

104:   Level: developer

106: .keywords: log, event, destroy
107: .seealso: PetscEventPerfLogCreate()
108: @*/
109: PetscErrorCode PetscEventPerfLogDestroy(PetscEventPerfLog eventLog)
110: {

114:   PetscFree(eventLog->eventInfo);
115:   PetscFree(eventLog);
116:   return(0);
117: }

119: /*------------------------------------------------ General Functions -------------------------------------------------*/
120: /*@C
121:   PetscEventPerfInfoClear - This clears a PetscEventPerfInfo object.

123:   Not collective

125:   Input Paramter:
126: . eventInfo - The PetscEventPerfInfo

128:   Level: developer

130: .keywords: log, event, destroy
131: .seealso: PetscEventPerfLogCreate()
132: @*/
133: PetscErrorCode PetscEventPerfInfoClear(PetscEventPerfInfo *eventInfo)
134: {
136:   eventInfo->id            = -1;
137:   eventInfo->active        = PETSC_TRUE;
138:   eventInfo->visible       = PETSC_TRUE;
139:   eventInfo->depth         = 0;
140:   eventInfo->count         = 0;
141:   eventInfo->flops         = 0.0;
142:   eventInfo->flops2        = 0.0;
143:   eventInfo->flopsTmp      = 0.0;
144:   eventInfo->time          = 0.0;
145:   eventInfo->time2         = 0.0;
146:   eventInfo->timeTmp       = 0.0;
147:   eventInfo->syncTime      = 0.0;
148:   eventInfo->dof[0]        = -1.0;
149:   eventInfo->dof[1]        = -1.0;
150:   eventInfo->dof[2]        = -1.0;
151:   eventInfo->dof[3]        = -1.0;
152:   eventInfo->dof[4]        = -1.0;
153:   eventInfo->dof[5]        = -1.0;
154:   eventInfo->dof[6]        = -1.0;
155:   eventInfo->dof[7]        = -1.0;
156:   eventInfo->errors[0]     = -1.0;
157:   eventInfo->errors[1]     = -1.0;
158:   eventInfo->errors[2]     = -1.0;
159:   eventInfo->errors[3]     = -1.0;
160:   eventInfo->errors[4]     = -1.0;
161:   eventInfo->errors[5]     = -1.0;
162:   eventInfo->errors[6]     = -1.0;
163:   eventInfo->errors[7]     = -1.0;
164:   eventInfo->numMessages   = 0.0;
165:   eventInfo->messageLength = 0.0;
166:   eventInfo->numReductions = 0.0;
167:   return(0);
168: }

170: /*@C
171:   PetscEventPerfInfoCopy - Copy the activity and visibility data in eventInfo to outInfo

173:   Not collective

175:   Input Paramter:
176: . eventInfo - The input PetscEventPerfInfo

178:   Output Paramter:
179: . outInfo   - The output PetscEventPerfInfo

181:   Level: developer

183: .keywords: log, event, copy
184: .seealso: PetscEventPerfInfoClear()
185: @*/
186: PetscErrorCode PetscEventPerfInfoCopy(PetscEventPerfInfo *eventInfo,PetscEventPerfInfo *outInfo)
187: {
189:   outInfo->id      = eventInfo->id;
190:   outInfo->active  = eventInfo->active;
191:   outInfo->visible = eventInfo->visible;
192:   return(0);
193: }

195: /*@C
196:   PetscEventPerfLogEnsureSize - This ensures that a PetscEventPerfLog is at least of a certain size.

198:   Not collective

200:   Input Paramters:
201: + eventLog - The PetscEventPerfLog
202: - size     - The size

204:   Level: developer

206: .keywords: log, event, size, ensure
207: .seealso: PetscEventPerfLogCreate()
208: @*/
209: PetscErrorCode PetscEventPerfLogEnsureSize(PetscEventPerfLog eventLog,int size)
210: {
211:   PetscEventPerfInfo *eventInfo;
212:   PetscErrorCode     ierr;

215:   while (size > eventLog->maxEvents) {
216:     PetscMalloc1(eventLog->maxEvents*2,&eventInfo);
217:     PetscMemcpy(eventInfo,eventLog->eventInfo,eventLog->maxEvents * sizeof(PetscEventPerfInfo));
218:     PetscFree(eventLog->eventInfo);
219:     eventLog->eventInfo  = eventInfo;
220:     eventLog->maxEvents *= 2;
221:   }
222:   while (eventLog->numEvents < size) {
223:     PetscEventPerfInfoClear(&eventLog->eventInfo[eventLog->numEvents++]);
224:   }
225:   return(0);
226: }

228: #if defined(PETSC_HAVE_MPE)
229: #include <mpe.h>
230: PETSC_INTERN PetscErrorCode PetscLogMPEGetRGBColor(const char*[]);
231: PetscErrorCode PetscLogEventBeginMPE(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
232: {
233:   PetscErrorCode    ierr;

236:   MPE_Log_event(petsc_stageLog->eventLog->eventInfo[event].mpe_id_begin,0,NULL);
237:   return(0);
238: }

240: PetscErrorCode PetscLogEventEndMPE(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
241: {
242:   PetscErrorCode    ierr;

245:   MPE_Log_event(petsc_stageLog->eventLog->eventInfo[event].mpe_id_end,0,NULL);
246:   return(0);
247: }
248: #endif

250: /*--------------------------------------------- Registration Functions ----------------------------------------------*/
251: /*@C
252:   PetscEventRegLogRegister - Registers an event for logging operations in an application code.

254:   Not Collective

256:   Input Parameters:
257: + eventLog - The PetscEventLog
258: . ename    - The name associated with the event
259: - classid   - The classid associated to the class for this event

261:   Output Parameter:
262: . event    - The event

264:   Example of Usage:
265: .vb
266:       int USER_EVENT;
267:       PetscLogDouble user_event_flops;
268:       PetscLogEventRegister("User event name",0,&USER_EVENT);
269:       PetscLogEventBegin(USER_EVENT,0,0,0,0);
270:          [code segment to monitor]
271:          PetscLogFlops(user_event_flops);
272:       PetscLogEventEnd(USER_EVENT,0,0,0,0);
273: .ve

275:   Notes:

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

284:   Level: developer

286: .keywords: log, event, register
287: .seealso: PetscLogEventBegin(), PetscLogEventEnd(), PetscLogFlops(), 
288:           PetscEventLogActivate(), PetscEventLogDeactivate()
289: @*/
290: PetscErrorCode PetscEventRegLogRegister(PetscEventRegLog eventLog,const char ename[],PetscClassId classid,PetscLogEvent *event)
291: {
292:   PetscEventRegInfo *eventInfo;
293:   char              *str;
294:   int               e;
295:   PetscErrorCode    ierr;

300:   /* Should check classid I think */
301:   e = eventLog->numEvents++;
302:   if (eventLog->numEvents > eventLog->maxEvents) {
303:     PetscMalloc1(eventLog->maxEvents*2,&eventInfo);
304:     PetscMemcpy(eventInfo,eventLog->eventInfo,eventLog->maxEvents * sizeof(PetscEventRegInfo));
305:     PetscFree(eventLog->eventInfo);
306:     eventLog->eventInfo  = eventInfo;
307:     eventLog->maxEvents *= 2;
308:   }
309:   PetscStrallocpy(ename,&str);

311:   eventLog->eventInfo[e].name       = str;
312:   eventLog->eventInfo[e].classid    = classid;
313:   eventLog->eventInfo[e].collective = PETSC_TRUE;
314: #if defined(PETSC_HAVE_MPE)
315:   if (PetscLogPLB == PetscLogEventBeginMPE) {
316:     const char  *color;
317:     PetscMPIInt rank;
318:     int         beginID,endID;

320:     beginID = MPE_Log_get_event_number();
321:     endID   = MPE_Log_get_event_number();

323:     eventLog->eventInfo[e].mpe_id_begin = beginID;
324:     eventLog->eventInfo[e].mpe_id_end   = endID;

326:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
327:     if (!rank) {
328:       PetscLogMPEGetRGBColor(&color);
329:       MPE_Describe_state(beginID,endID,str,(char*)color);
330:     }
331:   }
332: #endif
333:   *event = e;
334:   return(0);
335: }

337: /*---------------------------------------------- Activation Functions -----------------------------------------------*/
338: /*@C
339:   PetscEventPerfLogActivate - Indicates that a particular event should be logged.

341:   Not Collective

343:   Input Parameters:
344: + eventLog - The PetscEventPerfLog
345: - event    - The event

347:    Usage:
348: .vb
349:       PetscEventPerfLogDeactivate(log, VEC_SetValues);
350:         [code where you do not want to log VecSetValues()]
351:       PetscEventPerfLogActivate(log, VEC_SetValues);
352:         [code where you do want to log VecSetValues()]
353: .ve

355:   Note:
356:   The event may be either a pre-defined PETSc event (found in
357:   include/petsclog.h) or an event number obtained with PetscEventRegLogRegister().

359:   Level: developer

361: .keywords: log, event, activate
362: .seealso: PetscEventPerfLogDeactivate()
363: @*/
364: PetscErrorCode PetscEventPerfLogActivate(PetscEventPerfLog eventLog,PetscLogEvent event)
365: {
367:   eventLog->eventInfo[event].active = PETSC_TRUE;
368:   return(0);
369: }

371: /*@C
372:   PetscEventPerfLogDeactivate - Indicates that a particular event should not be logged.

374:   Not Collective

376:   Input Parameters:
377: + eventLog - The PetscEventPerfLog
378: - event    - The event

380:    Usage:
381: .vb
382:       PetscEventPerfLogDeactivate(log, VEC_SetValues);
383:         [code where you do not want to log VecSetValues()]
384:       PetscEventPerfLogActivate(log, VEC_SetValues);
385:         [code where you do want to log VecSetValues()]
386: .ve

388:   Note:
389:   The event may be either a pre-defined PETSc event (found in
390:   include/petsclog.h) or an event number obtained with PetscEventRegLogRegister().

392:   Level: developer

394: .keywords: log, event, activate
395: .seealso: PetscEventPerfLogActivate()
396: @*/
397: PetscErrorCode PetscEventPerfLogDeactivate(PetscEventPerfLog eventLog,PetscLogEvent event)
398: {
400:   eventLog->eventInfo[event].active = PETSC_FALSE;
401:   return(0);
402: }

404: /*@C
405:   PetscEventPerfLogActivateClass - Activates event logging for a PETSc object class.

407:   Not Collective

409:   Input Parameters:
410: + eventLog    - The PetscEventPerfLog
411: . eventRegLog - The PetscEventRegLog
412: - classid      - The class id, for example MAT_CLASSID, SNES_CLASSID,

414:   Level: developer

416: .seealso: PetscEventPerfLogDeactivateClass(), PetscEventPerfLogActivate(), PetscEventPerfLogDeactivate()
417: @*/
418: PetscErrorCode PetscEventPerfLogActivateClass(PetscEventPerfLog eventLog,PetscEventRegLog eventRegLog,PetscClassId classid)
419: {
420:   int e;

423:   for (e = 0; e < eventLog->numEvents; e++) {
424:     int c = eventRegLog->eventInfo[e].classid;
425:     if (c == classid) eventLog->eventInfo[e].active = PETSC_TRUE;
426:   }
427:   return(0);
428: }

430: /*@C
431:   PetscEventPerfLogDeactivateClass - Deactivates event logging for a PETSc object class.

433:   Not Collective

435:   Input Parameters:
436: + eventLog    - The PetscEventPerfLog
437: . eventRegLog - The PetscEventRegLog
438: - classid - The class id, for example MAT_CLASSID, SNES_CLASSID,

440:   Level: developer

442: .seealso: PetscEventPerfLogDeactivateClass(), PetscEventPerfLogDeactivate(), PetscEventPerfLogActivate()
443: @*/
444: PetscErrorCode PetscEventPerfLogDeactivateClass(PetscEventPerfLog eventLog,PetscEventRegLog eventRegLog,PetscClassId classid)
445: {
446:   int e;

449:   for (e = 0; e < eventLog->numEvents; e++) {
450:     int c = eventRegLog->eventInfo[e].classid;
451:     if (c == classid) eventLog->eventInfo[e].active = PETSC_FALSE;
452:   }
453:   return(0);
454: }

456: /*------------------------------------------------ Query Functions --------------------------------------------------*/
457: /*@C
458:   PetscEventRegLogGetEvent - This function returns the event id given the event name.

460:   Not Collective

462:   Input Parameters:
463: + eventLog - The PetscEventRegLog
464: - name     - The stage name

466:   Output Parameter:
467: . event    - The event id, or -1 if not found

469:   Level: developer

471: .keywords: log, stage
472: .seealso: PetscEventRegLogRegister()
473: @*/
474: PetscErrorCode  PetscEventRegLogGetEvent(PetscEventRegLog eventLog,const char name[],PetscLogEvent *event)
475: {
476:   PetscBool      match;
477:   int            e;

483:   *event = -1;
484:   for (e = 0; e < eventLog->numEvents; e++) {
485:     PetscStrcasecmp(eventLog->eventInfo[e].name,name,&match);
486:     if (match) {
487:       *event = e;
488:       break;
489:     }
490:   }
491:   return(0);
492: }

494: /*@C
495:   PetscEventPerfLogSetVisible - This function determines whether an event is printed during PetscLogView()

497:   Not Collective

499:   Input Parameters:
500: + eventLog  - The PetscEventPerfLog
501: . event     - The event to log
502: - isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)

504:   Database Options:
505: . -log_view - Activates log summary

507:   Level: developer

509: .keywords: log, visible, event
510: .seealso: PetscEventPerfLogGetVisible(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
511: @*/
512: PetscErrorCode PetscEventPerfLogSetVisible(PetscEventPerfLog eventLog,PetscLogEvent event,PetscBool isVisible)
513: {
515:   eventLog->eventInfo[event].visible = isVisible;
516:   return(0);
517: }

519: /*@C
520:   PetscEventPerfLogGetVisible - This function returns whether an event is printed during PetscLogView()

522:   Not Collective

524:   Input Parameters:
525: + eventLog  - The PetscEventPerfLog
526: - event     - The event id to log

528:   Output Parameter:
529: . isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)

531:   Database Options:
532: . -log_view - Activates log summary

534:   Level: developer

536: .keywords: log, visible, event
537: .seealso: PetscEventPerfLogSetVisible(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
538: @*/
539: PetscErrorCode PetscEventPerfLogGetVisible(PetscEventPerfLog eventLog,PetscLogEvent event,PetscBool  *isVisible)
540: {
543:   *isVisible = eventLog->eventInfo[event].visible;
544:   return(0);
545: }

547: /*@C
548:   PetscLogEventGetPerfInfo - Return the performance information about the given event in the given stage

550:   Input Parameters:
551: + stage - The stage number or PETSC_DETERMINE for the current stage
552: - event - The event number

554:   Output Parameters:
555: . info - This structure is filled with the performance information

557:   Level: Intermediate

559: .seealso: PetscLogEventGetFlops()
560: @*/
561: PetscErrorCode PetscLogEventGetPerfInfo(int stage,PetscLogEvent event,PetscEventPerfInfo *info)
562: {
563:   PetscStageLog     stageLog;
564:   PetscEventPerfLog eventLog = NULL;
565:   PetscErrorCode    ierr;

569:   if (!PetscLogPLB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must use -log_view or PetscLogDefaultBegin() before calling this routine");
570:   PetscLogGetStageLog(&stageLog);
571:   if (stage < 0) {PetscStageLogGetCurrent(stageLog,&stage);}
572:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
573:   *info = eventLog->eventInfo[event];
574:   return(0);
575: }

577: PetscErrorCode PetscLogEventGetFlops(PetscLogEvent event,PetscLogDouble *flops)
578: {
579:   PetscStageLog     stageLog;
580:   PetscEventPerfLog eventLog = NULL;
581:   int               stage;
582:   PetscErrorCode    ierr;

585:   if (!PetscLogPLB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must use -log_view or PetscLogDefaultBegin() before calling this routine");
586:   PetscLogGetStageLog(&stageLog);
587:   PetscStageLogGetCurrent(stageLog,&stage);
588:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
589:   *flops = eventLog->eventInfo[event].flops;
590:   return(0);
591: }

593: PetscErrorCode PetscLogEventZeroFlops(PetscLogEvent event)
594: {
595:   PetscStageLog     stageLog;
596:   PetscEventPerfLog eventLog = NULL;
597:   int               stage;
598:   PetscErrorCode    ierr;

601:   PetscLogGetStageLog(&stageLog);
602:   PetscStageLogGetCurrent(stageLog,&stage);
603:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);

605:   eventLog->eventInfo[event].flops    = 0.0;
606:   eventLog->eventInfo[event].flops2   = 0.0;
607:   eventLog->eventInfo[event].flopsTmp = 0.0;
608:   return(0);
609: }

611: PetscErrorCode PetscLogEventSynchronize(PetscLogEvent event,MPI_Comm comm)
612: {
613:   PetscStageLog     stageLog;
614:   PetscEventRegLog  eventRegLog;
615:   PetscEventPerfLog eventLog = NULL;
616:   int               stage;
617:   PetscLogDouble    time = 0.0;
618:   PetscErrorCode    ierr;

621:   if (!PetscLogSyncOn || comm == MPI_COMM_NULL) return(0);
622:   PetscLogGetStageLog(&stageLog);
623:   PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
624:   if (!eventRegLog->eventInfo[event].collective) return(0);
625:   PetscStageLogGetCurrent(stageLog,&stage);
626:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
627:   if (eventLog->eventInfo[event].depth > 0) return(0);

629:   PetscTimeSubtract(&time);
630:   MPI_Barrier(comm);
631:   PetscTimeAdd(&time);
632:   eventLog->eventInfo[event].syncTime += time;
633:   return(0);
634: }

636: PetscErrorCode PetscLogEventBeginDefault(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
637: {
638:   PetscStageLog     stageLog;
639:   PetscEventPerfLog eventLog = NULL;
640:   int               stage;
641:   PetscErrorCode    ierr;

644:   PetscLogGetStageLog(&stageLog);
645:   PetscStageLogGetCurrent(stageLog,&stage);
646:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
647:   /* Synchronization */
648:   PetscLogEventSynchronize(event,PetscObjectComm(o1));
649:   /* Check for double counting */
650:   eventLog->eventInfo[event].depth++;
651:   if (eventLog->eventInfo[event].depth > 1) return(0);
652:   /* Log the performance info */
653:   eventLog->eventInfo[event].count++;
654:   eventLog->eventInfo[event].timeTmp = 0.0;
655:   PetscTimeSubtract(&eventLog->eventInfo[event].timeTmp);
656:   eventLog->eventInfo[event].flopsTmp       = 0.0;
657:   eventLog->eventInfo[event].flopsTmp      -= petsc_TotalFlops;
658:   eventLog->eventInfo[event].numMessages   -= petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
659:   eventLog->eventInfo[event].messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
660:   eventLog->eventInfo[event].numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
661:   return(0);
662: }

664: PetscErrorCode PetscLogEventEndDefault(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
665: {
666:   PetscStageLog     stageLog;
667:   PetscEventPerfLog eventLog = NULL;
668:   int               stage;
669:   PetscErrorCode    ierr;

672:   PetscLogGetStageLog(&stageLog);
673:   PetscStageLogGetCurrent(stageLog,&stage);
674:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
675:   /* Check for double counting */
676:   eventLog->eventInfo[event].depth--;
677:   if (eventLog->eventInfo[event].depth > 0) return(0);
678:   else if (eventLog->eventInfo[event].depth < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Logging event had unbalanced begin/end pairs");
679:   /* Log performance info */
680:   PetscTimeAdd(&eventLog->eventInfo[event].timeTmp);
681:   eventLog->eventInfo[event].time          += eventLog->eventInfo[event].timeTmp;
682:   eventLog->eventInfo[event].time2         += eventLog->eventInfo[event].timeTmp*eventLog->eventInfo[event].timeTmp;
683:   eventLog->eventInfo[event].flopsTmp      += petsc_TotalFlops;
684:   eventLog->eventInfo[event].flops         += eventLog->eventInfo[event].flopsTmp;
685:   eventLog->eventInfo[event].flops2        += eventLog->eventInfo[event].flopsTmp*eventLog->eventInfo[event].flopsTmp;
686:   eventLog->eventInfo[event].numMessages   += petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
687:   eventLog->eventInfo[event].messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
688:   eventLog->eventInfo[event].numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
689:   return(0);
690: }

692: PetscErrorCode PetscLogEventBeginComplete(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
693: {
694:   PetscStageLog     stageLog;
695:   PetscEventRegLog  eventRegLog;
696:   PetscEventPerfLog eventPerfLog = NULL;
697:   Action            *tmpAction;
698:   PetscLogDouble    start,end;
699:   PetscLogDouble    curTime;
700:   int               stage;
701:   PetscErrorCode    ierr;

704:   /* Dynamically enlarge logging structures */
705:   if (petsc_numActions >= petsc_maxActions) {
706:     PetscTime(&start);
707:     PetscMalloc1(petsc_maxActions*2,&tmpAction);
708:     PetscMemcpy(tmpAction,petsc_actions,petsc_maxActions * sizeof(Action));
709:     PetscFree(petsc_actions);

711:     petsc_actions     = tmpAction;
712:     petsc_maxActions *= 2;
713:     PetscTime(&end);
714:     petsc_BaseTime += (end - start);
715:   }
716:   /* Record the event */
717:   PetscLogGetStageLog(&stageLog);
718:   PetscStageLogGetCurrent(stageLog,&stage);
719:   PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
720:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
721:   PetscTime(&curTime);
722:   if (petsc_logActions) {
723:     petsc_actions[petsc_numActions].time    = curTime - petsc_BaseTime;
724:     petsc_actions[petsc_numActions].action  = ACTIONBEGIN;
725:     petsc_actions[petsc_numActions].event   = event;
726:     petsc_actions[petsc_numActions].classid = eventRegLog->eventInfo[event].classid;
727:     if (o1) petsc_actions[petsc_numActions].id1 = o1->id;
728:     else petsc_actions[petsc_numActions].id1 = -1;
729:     if (o2) petsc_actions[petsc_numActions].id2 = o2->id;
730:     else petsc_actions[petsc_numActions].id2 = -1;
731:     if (o3) petsc_actions[petsc_numActions].id3 = o3->id;
732:     else petsc_actions[petsc_numActions].id3 = -1;
733:     petsc_actions[petsc_numActions].flops = petsc_TotalFlops;

735:     PetscMallocGetCurrentUsage(&petsc_actions[petsc_numActions].mem);
736:     PetscMallocGetMaximumUsage(&petsc_actions[petsc_numActions].maxmem);
737:     petsc_numActions++;
738:   }
739:   /* Check for double counting */
740:   eventPerfLog->eventInfo[event].depth++;
741:   if (eventPerfLog->eventInfo[event].depth > 1) return(0);
742:   /* Log the performance info */
743:   eventPerfLog->eventInfo[event].count++;
744:   eventPerfLog->eventInfo[event].time          -= curTime;
745:   eventPerfLog->eventInfo[event].flops         -= petsc_TotalFlops;
746:   eventPerfLog->eventInfo[event].numMessages   -= petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
747:   eventPerfLog->eventInfo[event].messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
748:   eventPerfLog->eventInfo[event].numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
749:   return(0);
750: }

752: PetscErrorCode PetscLogEventEndComplete(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
753: {
754:   PetscStageLog     stageLog;
755:   PetscEventRegLog  eventRegLog;
756:   PetscEventPerfLog eventPerfLog = NULL;
757:   Action            *tmpAction;
758:   PetscLogDouble    start,end;
759:   PetscLogDouble    curTime;
760:   int               stage;
761:   PetscErrorCode    ierr;

764:   /* Dynamically enlarge logging structures */
765:   if (petsc_numActions >= petsc_maxActions) {
766:     PetscTime(&start);
767:     PetscMalloc1(petsc_maxActions*2,&tmpAction);
768:     PetscMemcpy(tmpAction,petsc_actions,petsc_maxActions * sizeof(Action));
769:     PetscFree(petsc_actions);

771:     petsc_actions     = tmpAction;
772:     petsc_maxActions *= 2;
773:     PetscTime(&end);
774:     petsc_BaseTime += (end - start);
775:   }
776:   /* Record the event */
777:   PetscLogGetStageLog(&stageLog);
778:   PetscStageLogGetCurrent(stageLog,&stage);
779:   PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
780:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
781:   PetscTime(&curTime);
782:   if (petsc_logActions) {
783:     petsc_actions[petsc_numActions].time    = curTime - petsc_BaseTime;
784:     petsc_actions[petsc_numActions].action  = ACTIONEND;
785:     petsc_actions[petsc_numActions].event   = event;
786:     petsc_actions[petsc_numActions].classid = eventRegLog->eventInfo[event].classid;
787:     if (o1) petsc_actions[petsc_numActions].id1 = o1->id;
788:     else petsc_actions[petsc_numActions].id1 = -1;
789:     if (o2) petsc_actions[petsc_numActions].id2 = o2->id;
790:     else petsc_actions[petsc_numActions].id2 = -1;
791:     if (o3) petsc_actions[petsc_numActions].id3 = o3->id;
792:     else petsc_actions[petsc_numActions].id3 = -1;
793:     petsc_actions[petsc_numActions].flops = petsc_TotalFlops;

795:     PetscMallocGetCurrentUsage(&petsc_actions[petsc_numActions].mem);
796:     PetscMallocGetMaximumUsage(&petsc_actions[petsc_numActions].maxmem);
797:     petsc_numActions++;
798:   }
799:   /* Check for double counting */
800:   eventPerfLog->eventInfo[event].depth--;
801:   if (eventPerfLog->eventInfo[event].depth > 0) return(0);
802:   else if (eventPerfLog->eventInfo[event].depth < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Logging event had unbalanced begin/end pairs");
803:   /* Log the performance info */
804:   eventPerfLog->eventInfo[event].count++;
805:   eventPerfLog->eventInfo[event].time          += curTime;
806:   eventPerfLog->eventInfo[event].flops         += petsc_TotalFlops;
807:   eventPerfLog->eventInfo[event].numMessages   += petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
808:   eventPerfLog->eventInfo[event].messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
809:   eventPerfLog->eventInfo[event].numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
810:   return(0);
811: }

813: PetscErrorCode PetscLogEventBeginTrace(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
814: {
815:   PetscStageLog     stageLog;
816:   PetscEventRegLog  eventRegLog;
817:   PetscEventPerfLog eventPerfLog = NULL;
818:   PetscLogDouble    cur_time;
819:   PetscMPIInt       rank;
820:   int               stage,err;
821:   PetscErrorCode    ierr;

824:   if (!petsc_tracetime) PetscTime(&petsc_tracetime);

826:   petsc_tracelevel++;
827:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
828:   PetscLogGetStageLog(&stageLog);
829:   PetscStageLogGetCurrent(stageLog,&stage);
830:   PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
831:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
832:   /* Check for double counting */
833:   eventPerfLog->eventInfo[event].depth++;
834:   if (eventPerfLog->eventInfo[event].depth > 1) return(0);
835:   /* Log performance info */
836:   PetscTime(&cur_time);
837:   PetscFPrintf(PETSC_COMM_SELF,petsc_tracefile,"%s[%d] %g Event begin: %s\n",petsc_tracespace,rank,cur_time-petsc_tracetime,eventRegLog->eventInfo[event].name);
838:   PetscStrncpy(petsc_tracespace,petsc_traceblanks,2*petsc_tracelevel);

840:   petsc_tracespace[2*petsc_tracelevel] = 0;

842:   err = fflush(petsc_tracefile);
843:   if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
844:   return(0);
845: }

847: PetscErrorCode PetscLogEventEndTrace(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
848: {
849:   PetscStageLog     stageLog;
850:   PetscEventRegLog  eventRegLog;
851:   PetscEventPerfLog eventPerfLog = NULL;
852:   PetscLogDouble    cur_time;
853:   int               stage,err;
854:   PetscMPIInt       rank;
855:   PetscErrorCode    ierr;

858:   petsc_tracelevel--;
859:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
860:   PetscLogGetStageLog(&stageLog);
861:   PetscStageLogGetCurrent(stageLog,&stage);
862:   PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
863:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
864:   /* Check for double counting */
865:   eventPerfLog->eventInfo[event].depth--;
866:   if (eventPerfLog->eventInfo[event].depth > 0) return(0);
867:   else if (eventPerfLog->eventInfo[event].depth < 0 || petsc_tracelevel < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Logging event had unbalanced begin/end pairs");

869:   /* Log performance info */
870:   if (petsc_tracelevel) {
871:     PetscStrncpy(petsc_tracespace,petsc_traceblanks,2*petsc_tracelevel);
872:   }
873:   petsc_tracespace[2*petsc_tracelevel] = 0;
874:   PetscTime(&cur_time);
875:   PetscFPrintf(PETSC_COMM_SELF,petsc_tracefile,"%s[%d] %g Event end: %s\n",petsc_tracespace,rank,cur_time-petsc_tracetime,eventRegLog->eventInfo[event].name);
876:   err  = fflush(petsc_tracefile);
877:   if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
878:   return(0);
879: }

881: /*@C
882:   PetscLogEventSetDof - Set the nth number of degrees of freedom associated with this event

884:   Not Collective

886:   Input Parameters:
887: + event - The event id to log
888: . n     - The dof index, in [0, 8)
889: - dof   - The number of dofs

891:   Database Options:
892: . -log_view - Activates log summary

894:   Note: This is to enable logging of convergence

896:   Level: developer

898: .keywords: log, visible, event
899: .seealso: PetscLogEventSetError(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
900: @*/
901: PetscErrorCode PetscLogEventSetDof(PetscLogEvent event, PetscInt n, PetscLogDouble dof)
902: {
903:   PetscStageLog     stageLog;
904:   PetscEventPerfLog eventLog = NULL;
905:   int               stage;
906:   PetscErrorCode    ierr;

909:   if ((n < 0) || (n > 7)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Error index %D is not in [0, 8)", n);
910:   PetscLogGetStageLog(&stageLog);
911:   PetscStageLogGetCurrent(stageLog,&stage);
912:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
913:   eventLog->eventInfo[event].dof[n] = dof;
914:   return(0);
915: }

917: /*@C
918:   PetscLogEventSetError - Set the nth error associated with this event

920:   Not Collective

922:   Input Parameters:
923: + event - The event id to log
924: . n     - The error index, in [0, 8)
925: . error - The error

927:   Database Options:
928: . -log_view - Activates log summary

930:   Note: This is to enable logging of convergence, and enable users to interpret the errors as they wish. For example,
931:   as different norms, or as errors for different fields

933:   Level: developer

935: .keywords: log, visible, event
936: .seealso: PetscLogEventSetDof(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
937: @*/
938: PetscErrorCode PetscLogEventSetError(PetscLogEvent event, PetscInt n, PetscLogDouble error)
939: {
940:   PetscStageLog     stageLog;
941:   PetscEventPerfLog eventLog = NULL;
942:   int               stage;
943:   PetscErrorCode    ierr;

946:   if ((n < 0) || (n > 7)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Error index %D is not in [0, 8)", n);
947:   PetscLogGetStageLog(&stageLog);
948:   PetscStageLogGetCurrent(stageLog,&stage);
949:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
950:   eventLog->eventInfo[event].errors[n] = error;
951:   return(0);
952: }