Actual source code: eventlog.c


  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;
 11: PetscBool PetscLogMemory = PETSC_FALSE;
 12: #if defined(PETSC_HAVE_DEVICE)
 13: PetscBool PetscLogGpuTraffic = PETSC_FALSE;
 14: #endif

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

 19: /*@C
 20:   PetscEventRegLogCreate - This creates a PetscEventRegLog object.

 22:   Not collective

 24:   Input Parameter:
 25: . eventLog - The PetscEventRegLog

 27:   Level: developer

 29: .seealso: PetscEventRegLogDestroy(), PetscStageLogCreate()
 30: @*/
 31: PetscErrorCode PetscEventRegLogCreate(PetscEventRegLog *eventLog)
 32: {
 33:   PetscEventRegLog l;
 34:   PetscErrorCode   ierr;

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

 45: /*@C
 46:   PetscEventRegLogDestroy - This destroys a PetscEventRegLog object.

 48:   Not collective

 50:   Input Parameter:
 51: . eventLog - The PetscEventRegLog

 53:   Level: developer

 55: .seealso: PetscEventRegLogCreate()
 56: @*/
 57: PetscErrorCode PetscEventRegLogDestroy(PetscEventRegLog eventLog)
 58: {
 59:   int            e;

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

 71: /*@C
 72:   PetscEventPerfLogCreate - This creates a PetscEventPerfLog object.

 74:   Not collective

 76:   Input Parameter:
 77: . eventLog - The PetscEventPerfLog

 79:   Level: developer

 81: .seealso: PetscEventPerfLogDestroy(), PetscStageLogCreate()
 82: @*/
 83: PetscErrorCode PetscEventPerfLogCreate(PetscEventPerfLog *eventLog)
 84: {
 85:   PetscEventPerfLog l;
 86:   PetscErrorCode    ierr;

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

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

100:   Not collective

102:   Input Parameter:
103: . eventLog - The PetscEventPerfLog

105:   Level: developer

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 Parameter:
126: . eventInfo - The PetscEventPerfInfo

128:   Level: developer

130: .seealso: PetscEventPerfLogCreate()
131: @*/
132: PetscErrorCode PetscEventPerfInfoClear(PetscEventPerfInfo *eventInfo)
133: {
135:   eventInfo->id            = -1;
136:   eventInfo->active        = PETSC_TRUE;
137:   eventInfo->visible       = PETSC_TRUE;
138:   eventInfo->depth         = 0;
139:   eventInfo->count         = 0;
140:   eventInfo->flops         = 0.0;
141:   eventInfo->flops2        = 0.0;
142:   eventInfo->flopsTmp      = 0.0;
143:   eventInfo->time          = 0.0;
144:   eventInfo->time2         = 0.0;
145:   eventInfo->timeTmp       = 0.0;
146:   eventInfo->syncTime      = 0.0;
147:   eventInfo->dof[0]        = -1.0;
148:   eventInfo->dof[1]        = -1.0;
149:   eventInfo->dof[2]        = -1.0;
150:   eventInfo->dof[3]        = -1.0;
151:   eventInfo->dof[4]        = -1.0;
152:   eventInfo->dof[5]        = -1.0;
153:   eventInfo->dof[6]        = -1.0;
154:   eventInfo->dof[7]        = -1.0;
155:   eventInfo->errors[0]     = -1.0;
156:   eventInfo->errors[1]     = -1.0;
157:   eventInfo->errors[2]     = -1.0;
158:   eventInfo->errors[3]     = -1.0;
159:   eventInfo->errors[4]     = -1.0;
160:   eventInfo->errors[5]     = -1.0;
161:   eventInfo->errors[6]     = -1.0;
162:   eventInfo->errors[7]     = -1.0;
163:   eventInfo->numMessages   = 0.0;
164:   eventInfo->messageLength = 0.0;
165:   eventInfo->numReductions = 0.0;
166:   #if defined(PETSC_HAVE_DEVICE)
167:   eventInfo->CpuToGpuCount = 0.0;
168:   eventInfo->GpuToCpuCount = 0.0;
169:   eventInfo->CpuToGpuSize  = 0.0;
170:   eventInfo->GpuToCpuSize  = 0.0;
171:   eventInfo->GpuFlops      = 0.0;
172:   eventInfo->GpuTime       = 0.0;
173:   #endif
174:   return(0);
175: }

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

180:   Not collective

182:   Input Parameter:
183: . eventInfo - The input PetscEventPerfInfo

185:   Output Parameter:
186: . outInfo   - The output PetscEventPerfInfo

188:   Level: developer

190: .seealso: PetscEventPerfInfoClear()
191: @*/
192: PetscErrorCode PetscEventPerfInfoCopy(PetscEventPerfInfo *eventInfo,PetscEventPerfInfo *outInfo)
193: {
195:   outInfo->id      = eventInfo->id;
196:   outInfo->active  = eventInfo->active;
197:   outInfo->visible = eventInfo->visible;
198:   return(0);
199: }

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

204:   Not collective

206:   Input Parameters:
207: + eventLog - The PetscEventPerfLog
208: - size     - The size

210:   Level: developer

212: .seealso: PetscEventPerfLogCreate()
213: @*/
214: PetscErrorCode PetscEventPerfLogEnsureSize(PetscEventPerfLog eventLog,int size)
215: {
216:   PetscEventPerfInfo *eventInfo;
217:   PetscErrorCode     ierr;

220:   while (size > eventLog->maxEvents) {
221:     PetscCalloc1(eventLog->maxEvents*2,&eventInfo);
222:     PetscArraycpy(eventInfo,eventLog->eventInfo,eventLog->maxEvents);
223:     PetscFree(eventLog->eventInfo);
224:     eventLog->eventInfo  = eventInfo;
225:     eventLog->maxEvents *= 2;
226:   }
227:   while (eventLog->numEvents < size) {
228:     PetscEventPerfInfoClear(&eventLog->eventInfo[eventLog->numEvents++]);
229:   }
230:   return(0);
231: }

233: #if defined(PETSC_HAVE_MPE)
234: #include <mpe.h>
235: PETSC_INTERN PetscErrorCode PetscLogMPEGetRGBColor(const char*[]);
236: PetscErrorCode PetscLogEventBeginMPE(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
237: {
238:   PetscErrorCode    ierr;

241:   MPE_Log_event(petsc_stageLog->eventLog->eventInfo[event].mpe_id_begin,0,NULL);
242:   return(0);
243: }

245: PetscErrorCode PetscLogEventEndMPE(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
246: {
247:   PetscErrorCode    ierr;

250:   MPE_Log_event(petsc_stageLog->eventLog->eventInfo[event].mpe_id_end,0,NULL);
251:   return(0);
252: }
253: #endif

255: /*--------------------------------------------- Registration Functions ----------------------------------------------*/
256: /*@C
257:   PetscEventRegLogRegister - Registers an event for logging operations in an application code.

259:   Not Collective

261:   Input Parameters:
262: + eventLog - The PetscEventLog
263: . ename    - The name associated with the event
264: - classid   - The classid associated to the class for this event

266:   Output Parameter:
267: . event    - The event

269:   Example of Usage:
270: .vb
271:       int USER_EVENT;
272:       PetscLogDouble user_event_flops;
273:       PetscLogEventRegister("User event name",0,&USER_EVENT);
274:       PetscLogEventBegin(USER_EVENT,0,0,0,0);
275:          [code segment to monitor]
276:          PetscLogFlops(user_event_flops);
277:       PetscLogEventEnd(USER_EVENT,0,0,0,0);
278: .ve

280:   Notes:

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

289:   Level: developer

291: .seealso: PetscLogEventBegin(), PetscLogEventEnd(), PetscLogFlops(),
292:           PetscEventLogActivate(), PetscEventLogDeactivate()
293: @*/
294: PetscErrorCode PetscEventRegLogRegister(PetscEventRegLog eventLog,const char ename[],PetscClassId classid,PetscLogEvent *event)
295: {
296:   PetscEventRegInfo *eventInfo;
297:   char              *str;
298:   int               e;
299:   PetscErrorCode    ierr;

304:   /* Should check classid I think */
305:   e = eventLog->numEvents++;
306:   if (eventLog->numEvents > eventLog->maxEvents) {
307:     PetscCalloc1(eventLog->maxEvents*2,&eventInfo);
308:     PetscArraycpy(eventInfo,eventLog->eventInfo,eventLog->maxEvents);
309:     PetscFree(eventLog->eventInfo);
310:     eventLog->eventInfo  = eventInfo;
311:     eventLog->maxEvents *= 2;
312:   }
313:   PetscStrallocpy(ename,&str);

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

324:     beginID = MPE_Log_get_event_number();
325:     endID   = MPE_Log_get_event_number();

327:     eventLog->eventInfo[e].mpe_id_begin = beginID;
328:     eventLog->eventInfo[e].mpe_id_end   = endID;

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

341: /*---------------------------------------------- Activation Functions -----------------------------------------------*/
342: /*@C
343:   PetscEventPerfLogActivate - Indicates that a particular event should be logged.

345:   Not Collective

347:   Input Parameters:
348: + eventLog - The PetscEventPerfLog
349: - event    - The event

351:    Usage:
352: .vb
353:       PetscEventPerfLogDeactivate(log, VEC_SetValues);
354:         [code where you do not want to log VecSetValues()]
355:       PetscEventPerfLogActivate(log, VEC_SetValues);
356:         [code where you do want to log VecSetValues()]
357: .ve

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

363:   Level: developer

365: .seealso: PetscEventPerfLogDeactivate(), PetscEventPerfLogDeactivatePop(), PetscEventPerfLogDeactivatePush()
366: @*/
367: PetscErrorCode PetscEventPerfLogActivate(PetscEventPerfLog eventLog,PetscLogEvent event)
368: {
370:   eventLog->eventInfo[event].active = PETSC_TRUE;
371:   return(0);
372: }

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

377:   Not Collective

379:   Input Parameters:
380: + eventLog - The PetscEventPerfLog
381: - event    - The event

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

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

395:   Level: developer

397: .seealso: PetscEventPerfLogActivate(), PetscEventPerfLogDeactivatePop(), PetscEventPerfLogDeactivatePush()
398: @*/
399: PetscErrorCode PetscEventPerfLogDeactivate(PetscEventPerfLog eventLog,PetscLogEvent event)
400: {
402:   eventLog->eventInfo[event].active = PETSC_FALSE;
403:   return(0);
404: }

406: /*@C
407:   PetscEventPerfLogDeactivatePush - Indicates that a particular event should not be logged.

409:   Not Collective

411:   Input Parameters:
412: + eventLog - The PetscEventPerfLog
413: - event    - The event

415:    Usage:
416: .vb
417:       PetscEventPerfLogDeactivatePush(log, VEC_SetValues);
418:         [code where you do not want to log VecSetValues()]
419:       PetscEventPerfLogDeactivatePop(log, VEC_SetValues);
420:         [code where you do want to log VecSetValues()]
421: .ve

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

427:   Level: developer

429: .seealso: PetscEventPerfLogDeactivate(), PetscEventPerfLogActivate(), PetscEventPerfLogDeactivatePop()
430: @*/
431: PetscErrorCode PetscEventPerfLogDeactivatePush(PetscEventPerfLog eventLog,PetscLogEvent event)
432: {
434:   eventLog->eventInfo[event].depth++;
435:   return(0);
436: }

438: /*@C
439:   PetscEventPerfLogDeactivatePop - Indicates that a particular event should  be logged.

441:   Not Collective

443:   Input Parameters:
444: + eventLog - The PetscEventPerfLog
445: - event    - The event

447:    Usage:
448: .vb
449:       PetscEventPerfLogDeactivatePush(log, VEC_SetValues);
450:         [code where you do not want to log VecSetValues()]
451:       PetscEventPerfLogDeactivatePop(log, VEC_SetValues);
452:         [code where you do want to log VecSetValues()]
453: .ve

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

459:   Level: developer

461: .seealso: PetscEventPerfLogDeactivate(), PetscEventPerfLogActivate(), PetscEventPerfLogDeactivatePush()
462: @*/
463: PetscErrorCode PetscEventPerfLogDeactivatePop(PetscEventPerfLog eventLog,PetscLogEvent event)
464: {
466:   eventLog->eventInfo[event].depth--;
467:   return(0);
468: }

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

473:   Not Collective

475:   Input Parameters:
476: + eventLog    - The PetscEventPerfLog
477: . eventRegLog - The PetscEventRegLog
478: - classid      - The class id, for example MAT_CLASSID, SNES_CLASSID,

480:   Level: developer

482: .seealso: PetscEventPerfLogDeactivateClass(), PetscEventPerfLogActivate(), PetscEventPerfLogDeactivate()
483: @*/
484: PetscErrorCode PetscEventPerfLogActivateClass(PetscEventPerfLog eventLog,PetscEventRegLog eventRegLog,PetscClassId classid)
485: {
486:   int e;

489:   for (e = 0; e < eventLog->numEvents; e++) {
490:     int c = eventRegLog->eventInfo[e].classid;
491:     if (c == classid) eventLog->eventInfo[e].active = PETSC_TRUE;
492:   }
493:   return(0);
494: }

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

499:   Not Collective

501:   Input Parameters:
502: + eventLog    - The PetscEventPerfLog
503: . eventRegLog - The PetscEventRegLog
504: - classid - The class id, for example MAT_CLASSID, SNES_CLASSID,

506:   Level: developer

508: .seealso: PetscEventPerfLogDeactivateClass(), PetscEventPerfLogDeactivate(), PetscEventPerfLogActivate()
509: @*/
510: PetscErrorCode PetscEventPerfLogDeactivateClass(PetscEventPerfLog eventLog,PetscEventRegLog eventRegLog,PetscClassId classid)
511: {
512:   int e;

515:   for (e = 0; e < eventLog->numEvents; e++) {
516:     int c = eventRegLog->eventInfo[e].classid;
517:     if (c == classid) eventLog->eventInfo[e].active = PETSC_FALSE;
518:   }
519:   return(0);
520: }

522: /*------------------------------------------------ Query Functions --------------------------------------------------*/
523: /*@C
524:   PetscEventRegLogGetEvent - This function returns the event id given the event name.

526:   Not Collective

528:   Input Parameters:
529: + eventLog - The PetscEventRegLog
530: - name     - The stage name

532:   Output Parameter:
533: . event    - The event id, or -1 if not found

535:   Level: developer

537: .seealso: PetscEventRegLogRegister()
538: @*/
539: PetscErrorCode  PetscEventRegLogGetEvent(PetscEventRegLog eventLog,const char name[],PetscLogEvent *event)
540: {
541:   PetscBool      match;
542:   int            e;

548:   *event = -1;
549:   for (e = 0; e < eventLog->numEvents; e++) {
550:     PetscStrcasecmp(eventLog->eventInfo[e].name,name,&match);
551:     if (match) {
552:       *event = e;
553:       break;
554:     }
555:   }
556:   return(0);
557: }

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

562:   Not Collective

564:   Input Parameters:
565: + eventLog  - The PetscEventPerfLog
566: . event     - The event to log
567: - isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)

569:   Database Options:
570: . -log_view - Activates log summary

572:   Level: developer

574: .seealso: PetscEventPerfLogGetVisible(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
575: @*/
576: PetscErrorCode PetscEventPerfLogSetVisible(PetscEventPerfLog eventLog,PetscLogEvent event,PetscBool isVisible)
577: {
579:   eventLog->eventInfo[event].visible = isVisible;
580:   return(0);
581: }

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

586:   Not Collective

588:   Input Parameters:
589: + eventLog  - The PetscEventPerfLog
590: - event     - The event id to log

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

595:   Database Options:
596: . -log_view - Activates log summary

598:   Level: developer

600: .seealso: PetscEventPerfLogSetVisible(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
601: @*/
602: PetscErrorCode PetscEventPerfLogGetVisible(PetscEventPerfLog eventLog,PetscLogEvent event,PetscBool  *isVisible)
603: {
606:   *isVisible = eventLog->eventInfo[event].visible;
607:   return(0);
608: }

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

613:   Input Parameters:
614: + stage - The stage number or PETSC_DETERMINE for the current stage
615: - event - The event number

617:   Output Parameters:
618: . info - This structure is filled with the performance information

620:   Level: Intermediate

622: .seealso: PetscLogEventGetFlops()
623: @*/
624: PetscErrorCode PetscLogEventGetPerfInfo(int stage,PetscLogEvent event,PetscEventPerfInfo *info)
625: {
626:   PetscStageLog     stageLog;
627:   PetscEventPerfLog eventLog = NULL;
628:   PetscErrorCode    ierr;

632:   if (!PetscLogPLB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must use -log_view or PetscLogDefaultBegin() before calling this routine");
633:   PetscLogGetStageLog(&stageLog);
634:   if (stage < 0) {PetscStageLogGetCurrent(stageLog,&stage);}
635:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
636:   *info = eventLog->eventInfo[event];
637:   return(0);
638: }

640: PetscErrorCode PetscLogEventGetFlops(PetscLogEvent event,PetscLogDouble *flops)
641: {
642:   PetscStageLog     stageLog;
643:   PetscEventPerfLog eventLog = NULL;
644:   int               stage;
645:   PetscErrorCode    ierr;

648:   if (!PetscLogPLB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must use -log_view or PetscLogDefaultBegin() before calling this routine");
649:   PetscLogGetStageLog(&stageLog);
650:   PetscStageLogGetCurrent(stageLog,&stage);
651:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
652:   *flops = eventLog->eventInfo[event].flops;
653:   return(0);
654: }

656: PetscErrorCode PetscLogEventZeroFlops(PetscLogEvent event)
657: {
658:   PetscStageLog     stageLog;
659:   PetscEventPerfLog eventLog = NULL;
660:   int               stage;
661:   PetscErrorCode    ierr;

664:   PetscLogGetStageLog(&stageLog);
665:   PetscStageLogGetCurrent(stageLog,&stage);
666:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);

668:   eventLog->eventInfo[event].flops    = 0.0;
669:   eventLog->eventInfo[event].flops2   = 0.0;
670:   eventLog->eventInfo[event].flopsTmp = 0.0;
671:   return(0);
672: }

674: PetscErrorCode PetscLogEventSynchronize(PetscLogEvent event,MPI_Comm comm)
675: {
676:   PetscStageLog     stageLog;
677:   PetscEventRegLog  eventRegLog;
678:   PetscEventPerfLog eventLog = NULL;
679:   int               stage;
680:   PetscLogDouble    time = 0.0;
681:   PetscErrorCode    ierr;

684:   if (!PetscLogSyncOn || comm == MPI_COMM_NULL) return(0);
685:   PetscLogGetStageLog(&stageLog);
686:   PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
687:   if (!eventRegLog->eventInfo[event].collective) return(0);
688:   PetscStageLogGetCurrent(stageLog,&stage);
689:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
690:   if (eventLog->eventInfo[event].depth > 0) return(0);

692:   PetscTimeSubtract(&time);
693:   MPI_Barrier(comm);
694:   PetscTimeAdd(&time);
695:   eventLog->eventInfo[event].syncTime += time;
696:   return(0);
697: }

699: PetscErrorCode PetscLogEventBeginDefault(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
700: {
701:   PetscStageLog     stageLog;
702:   PetscEventPerfLog eventLog = NULL;
703:   int               stage;
704:   PetscErrorCode    ierr;

707:   PetscLogGetStageLog(&stageLog);
708:   PetscStageLogGetCurrent(stageLog,&stage);
709:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
710:   /* Synchronization */
711:   PetscLogEventSynchronize(event,PetscObjectComm(o1));
712:   /* Check for double counting */
713:   eventLog->eventInfo[event].depth++;
714:   if (eventLog->eventInfo[event].depth > 1) return(0);
715:   /* Log the performance info */
716:   eventLog->eventInfo[event].count++;
717:   eventLog->eventInfo[event].timeTmp = 0.0;
718:   PetscTimeSubtract(&eventLog->eventInfo[event].timeTmp);
719:   eventLog->eventInfo[event].flopsTmp       = -petsc_TotalFlops;
720:   eventLog->eventInfo[event].numMessages   -= petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
721:   eventLog->eventInfo[event].messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
722:   eventLog->eventInfo[event].numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
723:   if (PetscLogMemory) {
724:     PetscLogDouble usage;
725:     PetscMemoryGetCurrentUsage(&usage);
726:     eventLog->eventInfo[event].memIncrease -= usage;
727:     PetscMallocGetCurrentUsage(&usage);
728:     eventLog->eventInfo[event].mallocSpace -= usage;
729:     PetscMallocGetMaximumUsage(&usage);
730:     eventLog->eventInfo[event].mallocIncrease -= usage;
731:     PetscMallocPushMaximumUsage((int)event);
732:   }
733:   #if defined(PETSC_HAVE_DEVICE)
734:   eventLog->eventInfo[event].CpuToGpuCount -= petsc_ctog_ct;
735:   eventLog->eventInfo[event].GpuToCpuCount -= petsc_gtoc_ct;
736:   eventLog->eventInfo[event].CpuToGpuSize  -= petsc_ctog_sz;
737:   eventLog->eventInfo[event].GpuToCpuSize  -= petsc_gtoc_sz;
738:   eventLog->eventInfo[event].GpuFlops      -= petsc_gflops;
739:   eventLog->eventInfo[event].GpuTime       -= petsc_gtime;
740:   #endif
741:   return(0);
742: }

744: PetscErrorCode PetscLogEventEndDefault(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
745: {
746:   PetscStageLog     stageLog;
747:   PetscEventPerfLog eventLog = NULL;
748:   int               stage;
749:   PetscErrorCode    ierr;

752:   PetscLogGetStageLog(&stageLog);
753:   PetscStageLogGetCurrent(stageLog,&stage);
754:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
755:   /* Check for double counting */
756:   eventLog->eventInfo[event].depth--;
757:   if (eventLog->eventInfo[event].depth > 0) return(0);
758:   else if (eventLog->eventInfo[event].depth < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Logging event had unbalanced begin/end pairs");
759:   /* Log performance info */
760:   PetscTimeAdd(&eventLog->eventInfo[event].timeTmp);
761:   eventLog->eventInfo[event].time          += eventLog->eventInfo[event].timeTmp;
762:   eventLog->eventInfo[event].time2         += eventLog->eventInfo[event].timeTmp*eventLog->eventInfo[event].timeTmp;
763:   eventLog->eventInfo[event].flopsTmp      += petsc_TotalFlops;
764:   eventLog->eventInfo[event].flops         += eventLog->eventInfo[event].flopsTmp;
765:   eventLog->eventInfo[event].flops2        += eventLog->eventInfo[event].flopsTmp*eventLog->eventInfo[event].flopsTmp;
766:   eventLog->eventInfo[event].numMessages   += petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
767:   eventLog->eventInfo[event].messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
768:   eventLog->eventInfo[event].numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
769:   if (PetscLogMemory) {
770:     PetscLogDouble usage,musage;
771:     PetscMemoryGetCurrentUsage(&usage);
772:     eventLog->eventInfo[event].memIncrease += usage;
773:     PetscMallocGetCurrentUsage(&usage);
774:     eventLog->eventInfo[event].mallocSpace += usage;
775:     PetscMallocPopMaximumUsage((int)event,&musage);
776:     eventLog->eventInfo[event].mallocIncreaseEvent = PetscMax(musage-usage,eventLog->eventInfo[event].mallocIncreaseEvent);
777:     PetscMallocGetMaximumUsage(&usage);
778:     eventLog->eventInfo[event].mallocIncrease += usage;
779:   }
780:   #if defined(PETSC_HAVE_DEVICE)
781:   eventLog->eventInfo[event].CpuToGpuCount += petsc_ctog_ct;
782:   eventLog->eventInfo[event].GpuToCpuCount += petsc_gtoc_ct;
783:   eventLog->eventInfo[event].CpuToGpuSize  += petsc_ctog_sz;
784:   eventLog->eventInfo[event].GpuToCpuSize  += petsc_gtoc_sz;
785:   eventLog->eventInfo[event].GpuFlops      += petsc_gflops;
786:   eventLog->eventInfo[event].GpuTime       += petsc_gtime;
787:   #endif
788:   return(0);
789: }

791: PetscErrorCode PetscLogEventBeginComplete(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
792: {
793:   PetscStageLog     stageLog;
794:   PetscEventRegLog  eventRegLog;
795:   PetscEventPerfLog eventPerfLog = NULL;
796:   Action            *tmpAction;
797:   PetscLogDouble    start,end;
798:   PetscLogDouble    curTime;
799:   int               stage;
800:   PetscErrorCode    ierr;

803:   /* Dynamically enlarge logging structures */
804:   if (petsc_numActions >= petsc_maxActions) {
805:     PetscTime(&start);
806:     PetscCalloc1(petsc_maxActions*2,&tmpAction);
807:     PetscArraycpy(tmpAction,petsc_actions,petsc_maxActions);
808:     PetscFree(petsc_actions);

810:     petsc_actions     = tmpAction;
811:     petsc_maxActions *= 2;
812:     PetscTime(&end);
813:     petsc_BaseTime += (end - start);
814:   }
815:   /* Record the event */
816:   PetscLogGetStageLog(&stageLog);
817:   PetscStageLogGetCurrent(stageLog,&stage);
818:   PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
819:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
820:   PetscTime(&curTime);
821:   if (petsc_logActions) {
822:     petsc_actions[petsc_numActions].time    = curTime - petsc_BaseTime;
823:     petsc_actions[petsc_numActions].action  = ACTIONBEGIN;
824:     petsc_actions[petsc_numActions].event   = event;
825:     petsc_actions[petsc_numActions].classid = eventRegLog->eventInfo[event].classid;
826:     if (o1) petsc_actions[petsc_numActions].id1 = o1->id;
827:     else petsc_actions[petsc_numActions].id1 = -1;
828:     if (o2) petsc_actions[petsc_numActions].id2 = o2->id;
829:     else petsc_actions[petsc_numActions].id2 = -1;
830:     if (o3) petsc_actions[petsc_numActions].id3 = o3->id;
831:     else petsc_actions[petsc_numActions].id3 = -1;
832:     petsc_actions[petsc_numActions].flops = petsc_TotalFlops;

834:     PetscMallocGetCurrentUsage(&petsc_actions[petsc_numActions].mem);
835:     PetscMallocGetMaximumUsage(&petsc_actions[petsc_numActions].maxmem);
836:     petsc_numActions++;
837:   }
838:   /* Check for double counting */
839:   eventPerfLog->eventInfo[event].depth++;
840:   if (eventPerfLog->eventInfo[event].depth > 1) return(0);
841:   /* Log the performance info */
842:   eventPerfLog->eventInfo[event].count++;
843:   eventPerfLog->eventInfo[event].time          -= curTime;
844:   eventPerfLog->eventInfo[event].flops         -= petsc_TotalFlops;
845:   eventPerfLog->eventInfo[event].numMessages   -= petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
846:   eventPerfLog->eventInfo[event].messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
847:   eventPerfLog->eventInfo[event].numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
848:   return(0);
849: }

851: PetscErrorCode PetscLogEventEndComplete(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
852: {
853:   PetscStageLog     stageLog;
854:   PetscEventRegLog  eventRegLog;
855:   PetscEventPerfLog eventPerfLog = NULL;
856:   Action            *tmpAction;
857:   PetscLogDouble    start,end;
858:   PetscLogDouble    curTime;
859:   int               stage;
860:   PetscErrorCode    ierr;

863:   /* Dynamically enlarge logging structures */
864:   if (petsc_numActions >= petsc_maxActions) {
865:     PetscTime(&start);
866:     PetscCalloc1(petsc_maxActions*2,&tmpAction);
867:     PetscArraycpy(tmpAction,petsc_actions,petsc_maxActions);
868:     PetscFree(petsc_actions);

870:     petsc_actions     = tmpAction;
871:     petsc_maxActions *= 2;
872:     PetscTime(&end);
873:     petsc_BaseTime += (end - start);
874:   }
875:   /* Record the event */
876:   PetscLogGetStageLog(&stageLog);
877:   PetscStageLogGetCurrent(stageLog,&stage);
878:   PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
879:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
880:   PetscTime(&curTime);
881:   if (petsc_logActions) {
882:     petsc_actions[petsc_numActions].time    = curTime - petsc_BaseTime;
883:     petsc_actions[petsc_numActions].action  = ACTIONEND;
884:     petsc_actions[petsc_numActions].event   = event;
885:     petsc_actions[petsc_numActions].classid = eventRegLog->eventInfo[event].classid;
886:     if (o1) petsc_actions[petsc_numActions].id1 = o1->id;
887:     else petsc_actions[petsc_numActions].id1 = -1;
888:     if (o2) petsc_actions[petsc_numActions].id2 = o2->id;
889:     else petsc_actions[petsc_numActions].id2 = -1;
890:     if (o3) petsc_actions[petsc_numActions].id3 = o3->id;
891:     else petsc_actions[petsc_numActions].id3 = -1;
892:     petsc_actions[petsc_numActions].flops = petsc_TotalFlops;

894:     PetscMallocGetCurrentUsage(&petsc_actions[petsc_numActions].mem);
895:     PetscMallocGetMaximumUsage(&petsc_actions[petsc_numActions].maxmem);
896:     petsc_numActions++;
897:   }
898:   /* Check for double counting */
899:   eventPerfLog->eventInfo[event].depth--;
900:   if (eventPerfLog->eventInfo[event].depth > 0) return(0);
901:   else if (eventPerfLog->eventInfo[event].depth < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Logging event had unbalanced begin/end pairs");
902:   /* Log the performance info */
903:   eventPerfLog->eventInfo[event].count++;
904:   eventPerfLog->eventInfo[event].time          += curTime;
905:   eventPerfLog->eventInfo[event].flops         += petsc_TotalFlops;
906:   eventPerfLog->eventInfo[event].numMessages   += petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
907:   eventPerfLog->eventInfo[event].messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
908:   eventPerfLog->eventInfo[event].numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
909:   return(0);
910: }

912: PetscErrorCode PetscLogEventBeginTrace(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
913: {
914:   PetscStageLog     stageLog;
915:   PetscEventRegLog  eventRegLog;
916:   PetscEventPerfLog eventPerfLog = NULL;
917:   PetscLogDouble    cur_time;
918:   PetscMPIInt       rank;
919:   int               stage,err;
920:   PetscErrorCode    ierr;

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

925:   petsc_tracelevel++;
926:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
927:   PetscLogGetStageLog(&stageLog);
928:   PetscStageLogGetCurrent(stageLog,&stage);
929:   PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
930:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
931:   /* Check for double counting */
932:   eventPerfLog->eventInfo[event].depth++;
933:   if (eventPerfLog->eventInfo[event].depth > 1) return(0);
934:   /* Log performance info */
935:   PetscTime(&cur_time);
936:   PetscFPrintf(PETSC_COMM_SELF,petsc_tracefile,"%s[%d] %g Event begin: %s\n",petsc_tracespace,rank,cur_time-petsc_tracetime,eventRegLog->eventInfo[event].name);
937:   PetscStrncpy(petsc_tracespace,petsc_traceblanks,2*petsc_tracelevel);

939:   petsc_tracespace[2*petsc_tracelevel] = 0;

941:   err = fflush(petsc_tracefile);
942:   if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
943:   return(0);
944: }

946: PetscErrorCode PetscLogEventEndTrace(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
947: {
948:   PetscStageLog     stageLog;
949:   PetscEventRegLog  eventRegLog;
950:   PetscEventPerfLog eventPerfLog = NULL;
951:   PetscLogDouble    cur_time;
952:   int               stage,err;
953:   PetscMPIInt       rank;
954:   PetscErrorCode    ierr;

957:   petsc_tracelevel--;
958:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
959:   PetscLogGetStageLog(&stageLog);
960:   PetscStageLogGetCurrent(stageLog,&stage);
961:   PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
962:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
963:   /* Check for double counting */
964:   eventPerfLog->eventInfo[event].depth--;
965:   if (eventPerfLog->eventInfo[event].depth > 0) return(0);
966:   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");

968:   /* Log performance info */
969:   if (petsc_tracelevel) {
970:     PetscStrncpy(petsc_tracespace,petsc_traceblanks,2*petsc_tracelevel);
971:   }
972:   petsc_tracespace[2*petsc_tracelevel] = 0;
973:   PetscTime(&cur_time);
974:   PetscFPrintf(PETSC_COMM_SELF,petsc_tracefile,"%s[%d] %g Event end: %s\n",petsc_tracespace,rank,cur_time-petsc_tracetime,eventRegLog->eventInfo[event].name);
975:   err  = fflush(petsc_tracefile);
976:   if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
977:   return(0);
978: }

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

983:   Not Collective

985:   Input Parameters:
986: + event - The event id to log
987: . n     - The dof index, in [0, 8)
988: - dof   - The number of dofs

990:   Database Options:
991: . -log_view - Activates log summary

993:   Note: This is to enable logging of convergence

995:   Level: developer

997: .seealso: PetscLogEventSetError(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
998: @*/
999: PetscErrorCode PetscLogEventSetDof(PetscLogEvent event, PetscInt n, PetscLogDouble dof)
1000: {
1001:   PetscStageLog     stageLog;
1002:   PetscEventPerfLog eventLog = NULL;
1003:   int               stage;
1004:   PetscErrorCode    ierr;

1007:   if ((n < 0) || (n > 7)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Error index %D is not in [0, 8)", n);
1008:   PetscLogGetStageLog(&stageLog);
1009:   PetscStageLogGetCurrent(stageLog,&stage);
1010:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
1011:   eventLog->eventInfo[event].dof[n] = dof;
1012:   return(0);
1013: }

1015: /*@C
1016:   PetscLogEventSetError - Set the nth error associated with this event

1018:   Not Collective

1020:   Input Parameters:
1021: + event - The event id to log
1022: . n     - The error index, in [0, 8)
1023: - error - The error

1025:   Database Options:
1026: . -log_view - Activates log summary

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

1031:   Level: developer

1033: .seealso: PetscLogEventSetDof(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
1034: @*/
1035: PetscErrorCode PetscLogEventSetError(PetscLogEvent event, PetscInt n, PetscLogDouble error)
1036: {
1037:   PetscStageLog     stageLog;
1038:   PetscEventPerfLog eventLog = NULL;
1039:   int               stage;
1040:   PetscErrorCode    ierr;

1043:   if ((n < 0) || (n > 7)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Error index %D is not in [0, 8)", n);
1044:   PetscLogGetStageLog(&stageLog);
1045:   PetscStageLogGetCurrent(stageLog,&stage);
1046:   PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
1047:   eventLog->eventInfo[event].errors[n] = error;
1048:   return(0);
1049: }