Actual source code: eventlog.c
petsc-3.12.5 2020-03-29
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_VIENNACL) || defined(PETSC_HAVE_CUDA)
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 Paramter:
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 Paramter:
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 Paramter:
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_VIENNACL) || defined(PETSC_HAVE_CUDA)
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 Paramter:
183: . eventInfo - The input PetscEventPerfInfo
185: Output Paramter:
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 Paramters:
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()
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()
398: @*/
399: PetscErrorCode PetscEventPerfLogDeactivate(PetscEventPerfLog eventLog,PetscLogEvent event)
400: {
402: eventLog->eventInfo[event].active = PETSC_FALSE;
403: return(0);
404: }
406: /*@C
407: PetscEventPerfLogActivateClass - Activates event logging for a PETSc object class.
409: Not Collective
411: Input Parameters:
412: + eventLog - The PetscEventPerfLog
413: . eventRegLog - The PetscEventRegLog
414: - classid - The class id, for example MAT_CLASSID, SNES_CLASSID,
416: Level: developer
418: .seealso: PetscEventPerfLogDeactivateClass(), PetscEventPerfLogActivate(), PetscEventPerfLogDeactivate()
419: @*/
420: PetscErrorCode PetscEventPerfLogActivateClass(PetscEventPerfLog eventLog,PetscEventRegLog eventRegLog,PetscClassId classid)
421: {
422: int e;
425: for (e = 0; e < eventLog->numEvents; e++) {
426: int c = eventRegLog->eventInfo[e].classid;
427: if (c == classid) eventLog->eventInfo[e].active = PETSC_TRUE;
428: }
429: return(0);
430: }
432: /*@C
433: PetscEventPerfLogDeactivateClass - Deactivates event logging for a PETSc object class.
435: Not Collective
437: Input Parameters:
438: + eventLog - The PetscEventPerfLog
439: . eventRegLog - The PetscEventRegLog
440: - classid - The class id, for example MAT_CLASSID, SNES_CLASSID,
442: Level: developer
444: .seealso: PetscEventPerfLogDeactivateClass(), PetscEventPerfLogDeactivate(), PetscEventPerfLogActivate()
445: @*/
446: PetscErrorCode PetscEventPerfLogDeactivateClass(PetscEventPerfLog eventLog,PetscEventRegLog eventRegLog,PetscClassId classid)
447: {
448: int e;
451: for (e = 0; e < eventLog->numEvents; e++) {
452: int c = eventRegLog->eventInfo[e].classid;
453: if (c == classid) eventLog->eventInfo[e].active = PETSC_FALSE;
454: }
455: return(0);
456: }
458: /*------------------------------------------------ Query Functions --------------------------------------------------*/
459: /*@C
460: PetscEventRegLogGetEvent - This function returns the event id given the event name.
462: Not Collective
464: Input Parameters:
465: + eventLog - The PetscEventRegLog
466: - name - The stage name
468: Output Parameter:
469: . event - The event id, or -1 if not found
471: Level: developer
473: .seealso: PetscEventRegLogRegister()
474: @*/
475: PetscErrorCode PetscEventRegLogGetEvent(PetscEventRegLog eventLog,const char name[],PetscLogEvent *event)
476: {
477: PetscBool match;
478: int e;
484: *event = -1;
485: for (e = 0; e < eventLog->numEvents; e++) {
486: PetscStrcasecmp(eventLog->eventInfo[e].name,name,&match);
487: if (match) {
488: *event = e;
489: break;
490: }
491: }
492: return(0);
493: }
495: /*@C
496: PetscEventPerfLogSetVisible - This function determines whether an event is printed during PetscLogView()
498: Not Collective
500: Input Parameters:
501: + eventLog - The PetscEventPerfLog
502: . event - The event to log
503: - isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)
505: Database Options:
506: . -log_view - Activates log summary
508: Level: developer
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: .seealso: PetscEventPerfLogSetVisible(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
537: @*/
538: PetscErrorCode PetscEventPerfLogGetVisible(PetscEventPerfLog eventLog,PetscLogEvent event,PetscBool *isVisible)
539: {
542: *isVisible = eventLog->eventInfo[event].visible;
543: return(0);
544: }
546: /*@C
547: PetscLogEventGetPerfInfo - Return the performance information about the given event in the given stage
549: Input Parameters:
550: + stage - The stage number or PETSC_DETERMINE for the current stage
551: - event - The event number
553: Output Parameters:
554: . info - This structure is filled with the performance information
556: Level: Intermediate
558: .seealso: PetscLogEventGetFlops()
559: @*/
560: PetscErrorCode PetscLogEventGetPerfInfo(int stage,PetscLogEvent event,PetscEventPerfInfo *info)
561: {
562: PetscStageLog stageLog;
563: PetscEventPerfLog eventLog = NULL;
564: PetscErrorCode ierr;
568: if (!PetscLogPLB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must use -log_view or PetscLogDefaultBegin() before calling this routine");
569: PetscLogGetStageLog(&stageLog);
570: if (stage < 0) {PetscStageLogGetCurrent(stageLog,&stage);}
571: PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
572: *info = eventLog->eventInfo[event];
573: return(0);
574: }
576: PetscErrorCode PetscLogEventGetFlops(PetscLogEvent event,PetscLogDouble *flops)
577: {
578: PetscStageLog stageLog;
579: PetscEventPerfLog eventLog = NULL;
580: int stage;
581: PetscErrorCode ierr;
584: if (!PetscLogPLB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must use -log_view or PetscLogDefaultBegin() before calling this routine");
585: PetscLogGetStageLog(&stageLog);
586: PetscStageLogGetCurrent(stageLog,&stage);
587: PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
588: *flops = eventLog->eventInfo[event].flops;
589: return(0);
590: }
592: PetscErrorCode PetscLogEventZeroFlops(PetscLogEvent event)
593: {
594: PetscStageLog stageLog;
595: PetscEventPerfLog eventLog = NULL;
596: int stage;
597: PetscErrorCode ierr;
600: PetscLogGetStageLog(&stageLog);
601: PetscStageLogGetCurrent(stageLog,&stage);
602: PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
604: eventLog->eventInfo[event].flops = 0.0;
605: eventLog->eventInfo[event].flops2 = 0.0;
606: eventLog->eventInfo[event].flopsTmp = 0.0;
607: return(0);
608: }
610: PetscErrorCode PetscLogEventSynchronize(PetscLogEvent event,MPI_Comm comm)
611: {
612: PetscStageLog stageLog;
613: PetscEventRegLog eventRegLog;
614: PetscEventPerfLog eventLog = NULL;
615: int stage;
616: PetscLogDouble time = 0.0;
617: PetscErrorCode ierr;
620: if (!PetscLogSyncOn || comm == MPI_COMM_NULL) return(0);
621: PetscLogGetStageLog(&stageLog);
622: PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
623: if (!eventRegLog->eventInfo[event].collective) return(0);
624: PetscStageLogGetCurrent(stageLog,&stage);
625: PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
626: if (eventLog->eventInfo[event].depth > 0) return(0);
628: PetscTimeSubtract(&time);
629: MPI_Barrier(comm);
630: PetscTimeAdd(&time);
631: eventLog->eventInfo[event].syncTime += time;
632: return(0);
633: }
635: PetscErrorCode PetscLogEventBeginDefault(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
636: {
637: PetscStageLog stageLog;
638: PetscEventPerfLog eventLog = NULL;
639: int stage;
640: PetscErrorCode ierr;
643: PetscLogGetStageLog(&stageLog);
644: PetscStageLogGetCurrent(stageLog,&stage);
645: PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
646: /* Synchronization */
647: PetscLogEventSynchronize(event,PetscObjectComm(o1));
648: /* Check for double counting */
649: eventLog->eventInfo[event].depth++;
650: if (eventLog->eventInfo[event].depth > 1) return(0);
651: /* Log the performance info */
652: eventLog->eventInfo[event].count++;
653: eventLog->eventInfo[event].timeTmp = 0.0;
654: PetscTimeSubtract(&eventLog->eventInfo[event].timeTmp);
655: eventLog->eventInfo[event].flopsTmp = 0.0;
656: eventLog->eventInfo[event].flopsTmp -= petsc_TotalFlops;
657: eventLog->eventInfo[event].numMessages -= petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
658: eventLog->eventInfo[event].messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
659: eventLog->eventInfo[event].numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
660: if (PetscLogMemory) {
661: PetscLogDouble usage;
662: PetscMemoryGetCurrentUsage(&usage);
663: eventLog->eventInfo[event].memIncrease -= usage;
664: PetscMallocGetCurrentUsage(&usage);
665: eventLog->eventInfo[event].mallocSpace -= usage;
666: PetscMallocGetMaximumUsage(&usage);
667: eventLog->eventInfo[event].mallocIncrease -= usage;
668: PetscMallocPushMaximumUsage((int)event);
669: }
670: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
671: eventLog->eventInfo[event].CpuToGpuCount -= petsc_ctog_ct;
672: eventLog->eventInfo[event].GpuToCpuCount -= petsc_gtoc_ct;
673: eventLog->eventInfo[event].CpuToGpuSize -= petsc_ctog_sz;
674: eventLog->eventInfo[event].GpuToCpuSize -= petsc_gtoc_sz;
675: eventLog->eventInfo[event].GpuFlops -= petsc_gflops;
676: eventLog->eventInfo[event].GpuTime -= petsc_gtime;
677: #endif
678: return(0);
679: }
681: PetscErrorCode PetscLogEventEndDefault(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
682: {
683: PetscStageLog stageLog;
684: PetscEventPerfLog eventLog = NULL;
685: int stage;
686: PetscErrorCode ierr;
689: PetscLogGetStageLog(&stageLog);
690: PetscStageLogGetCurrent(stageLog,&stage);
691: PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
692: /* Check for double counting */
693: eventLog->eventInfo[event].depth--;
694: if (eventLog->eventInfo[event].depth > 0) return(0);
695: else if (eventLog->eventInfo[event].depth < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Logging event had unbalanced begin/end pairs");
696: /* Log performance info */
697: PetscTimeAdd(&eventLog->eventInfo[event].timeTmp);
698: eventLog->eventInfo[event].time += eventLog->eventInfo[event].timeTmp;
699: eventLog->eventInfo[event].time2 += eventLog->eventInfo[event].timeTmp*eventLog->eventInfo[event].timeTmp;
700: eventLog->eventInfo[event].flopsTmp += petsc_TotalFlops;
701: eventLog->eventInfo[event].flops += eventLog->eventInfo[event].flopsTmp;
702: eventLog->eventInfo[event].flops2 += eventLog->eventInfo[event].flopsTmp*eventLog->eventInfo[event].flopsTmp;
703: eventLog->eventInfo[event].numMessages += petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
704: eventLog->eventInfo[event].messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
705: eventLog->eventInfo[event].numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
706: if (PetscLogMemory) {
707: PetscLogDouble usage,musage;
708: PetscMemoryGetCurrentUsage(&usage);
709: eventLog->eventInfo[event].memIncrease += usage;
710: PetscMallocGetCurrentUsage(&usage);
711: eventLog->eventInfo[event].mallocSpace += usage;
712: PetscMallocPopMaximumUsage((int)event,&musage);
713: eventLog->eventInfo[event].mallocIncreaseEvent = PetscMax(musage-usage,eventLog->eventInfo[event].mallocIncreaseEvent);
714: PetscMallocGetMaximumUsage(&usage);
715: eventLog->eventInfo[event].mallocIncrease += usage;
716: }
717: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
718: eventLog->eventInfo[event].CpuToGpuCount += petsc_ctog_ct;
719: eventLog->eventInfo[event].GpuToCpuCount += petsc_gtoc_ct;
720: eventLog->eventInfo[event].CpuToGpuSize += petsc_ctog_sz;
721: eventLog->eventInfo[event].GpuToCpuSize += petsc_gtoc_sz;
722: eventLog->eventInfo[event].GpuFlops += petsc_gflops;
723: eventLog->eventInfo[event].GpuTime += petsc_gtime;
724: #endif
725: return(0);
726: }
728: PetscErrorCode PetscLogEventBeginComplete(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
729: {
730: PetscStageLog stageLog;
731: PetscEventRegLog eventRegLog;
732: PetscEventPerfLog eventPerfLog = NULL;
733: Action *tmpAction;
734: PetscLogDouble start,end;
735: PetscLogDouble curTime;
736: int stage;
737: PetscErrorCode ierr;
740: /* Dynamically enlarge logging structures */
741: if (petsc_numActions >= petsc_maxActions) {
742: PetscTime(&start);
743: PetscCalloc1(petsc_maxActions*2,&tmpAction);
744: PetscArraycpy(tmpAction,petsc_actions,petsc_maxActions);
745: PetscFree(petsc_actions);
747: petsc_actions = tmpAction;
748: petsc_maxActions *= 2;
749: PetscTime(&end);
750: petsc_BaseTime += (end - start);
751: }
752: /* Record the event */
753: PetscLogGetStageLog(&stageLog);
754: PetscStageLogGetCurrent(stageLog,&stage);
755: PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
756: PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
757: PetscTime(&curTime);
758: if (petsc_logActions) {
759: petsc_actions[petsc_numActions].time = curTime - petsc_BaseTime;
760: petsc_actions[petsc_numActions].action = ACTIONBEGIN;
761: petsc_actions[petsc_numActions].event = event;
762: petsc_actions[petsc_numActions].classid = eventRegLog->eventInfo[event].classid;
763: if (o1) petsc_actions[petsc_numActions].id1 = o1->id;
764: else petsc_actions[petsc_numActions].id1 = -1;
765: if (o2) petsc_actions[petsc_numActions].id2 = o2->id;
766: else petsc_actions[petsc_numActions].id2 = -1;
767: if (o3) petsc_actions[petsc_numActions].id3 = o3->id;
768: else petsc_actions[petsc_numActions].id3 = -1;
769: petsc_actions[petsc_numActions].flops = petsc_TotalFlops;
771: PetscMallocGetCurrentUsage(&petsc_actions[petsc_numActions].mem);
772: PetscMallocGetMaximumUsage(&petsc_actions[petsc_numActions].maxmem);
773: petsc_numActions++;
774: }
775: /* Check for double counting */
776: eventPerfLog->eventInfo[event].depth++;
777: if (eventPerfLog->eventInfo[event].depth > 1) return(0);
778: /* Log the performance info */
779: eventPerfLog->eventInfo[event].count++;
780: eventPerfLog->eventInfo[event].time -= curTime;
781: eventPerfLog->eventInfo[event].flops -= petsc_TotalFlops;
782: eventPerfLog->eventInfo[event].numMessages -= petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
783: eventPerfLog->eventInfo[event].messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
784: eventPerfLog->eventInfo[event].numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
785: return(0);
786: }
788: PetscErrorCode PetscLogEventEndComplete(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
789: {
790: PetscStageLog stageLog;
791: PetscEventRegLog eventRegLog;
792: PetscEventPerfLog eventPerfLog = NULL;
793: Action *tmpAction;
794: PetscLogDouble start,end;
795: PetscLogDouble curTime;
796: int stage;
797: PetscErrorCode ierr;
800: /* Dynamically enlarge logging structures */
801: if (petsc_numActions >= petsc_maxActions) {
802: PetscTime(&start);
803: PetscCalloc1(petsc_maxActions*2,&tmpAction);
804: PetscArraycpy(tmpAction,petsc_actions,petsc_maxActions);
805: PetscFree(petsc_actions);
807: petsc_actions = tmpAction;
808: petsc_maxActions *= 2;
809: PetscTime(&end);
810: petsc_BaseTime += (end - start);
811: }
812: /* Record the event */
813: PetscLogGetStageLog(&stageLog);
814: PetscStageLogGetCurrent(stageLog,&stage);
815: PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
816: PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
817: PetscTime(&curTime);
818: if (petsc_logActions) {
819: petsc_actions[petsc_numActions].time = curTime - petsc_BaseTime;
820: petsc_actions[petsc_numActions].action = ACTIONEND;
821: petsc_actions[petsc_numActions].event = event;
822: petsc_actions[petsc_numActions].classid = eventRegLog->eventInfo[event].classid;
823: if (o1) petsc_actions[petsc_numActions].id1 = o1->id;
824: else petsc_actions[petsc_numActions].id1 = -1;
825: if (o2) petsc_actions[petsc_numActions].id2 = o2->id;
826: else petsc_actions[petsc_numActions].id2 = -1;
827: if (o3) petsc_actions[petsc_numActions].id3 = o3->id;
828: else petsc_actions[petsc_numActions].id3 = -1;
829: petsc_actions[petsc_numActions].flops = petsc_TotalFlops;
831: PetscMallocGetCurrentUsage(&petsc_actions[petsc_numActions].mem);
832: PetscMallocGetMaximumUsage(&petsc_actions[petsc_numActions].maxmem);
833: petsc_numActions++;
834: }
835: /* Check for double counting */
836: eventPerfLog->eventInfo[event].depth--;
837: if (eventPerfLog->eventInfo[event].depth > 0) return(0);
838: else if (eventPerfLog->eventInfo[event].depth < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Logging event had unbalanced begin/end pairs");
839: /* Log the performance info */
840: eventPerfLog->eventInfo[event].count++;
841: eventPerfLog->eventInfo[event].time += curTime;
842: eventPerfLog->eventInfo[event].flops += petsc_TotalFlops;
843: eventPerfLog->eventInfo[event].numMessages += petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
844: eventPerfLog->eventInfo[event].messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
845: eventPerfLog->eventInfo[event].numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
846: return(0);
847: }
849: PetscErrorCode PetscLogEventBeginTrace(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
850: {
851: PetscStageLog stageLog;
852: PetscEventRegLog eventRegLog;
853: PetscEventPerfLog eventPerfLog = NULL;
854: PetscLogDouble cur_time;
855: PetscMPIInt rank;
856: int stage,err;
857: PetscErrorCode ierr;
860: if (!petsc_tracetime) PetscTime(&petsc_tracetime);
862: petsc_tracelevel++;
863: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
864: PetscLogGetStageLog(&stageLog);
865: PetscStageLogGetCurrent(stageLog,&stage);
866: PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
867: PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
868: /* Check for double counting */
869: eventPerfLog->eventInfo[event].depth++;
870: if (eventPerfLog->eventInfo[event].depth > 1) return(0);
871: /* Log performance info */
872: PetscTime(&cur_time);
873: PetscFPrintf(PETSC_COMM_SELF,petsc_tracefile,"%s[%d] %g Event begin: %s\n",petsc_tracespace,rank,cur_time-petsc_tracetime,eventRegLog->eventInfo[event].name);
874: PetscStrncpy(petsc_tracespace,petsc_traceblanks,2*petsc_tracelevel);
876: petsc_tracespace[2*petsc_tracelevel] = 0;
878: err = fflush(petsc_tracefile);
879: if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
880: return(0);
881: }
883: PetscErrorCode PetscLogEventEndTrace(PetscLogEvent event,int t,PetscObject o1,PetscObject o2,PetscObject o3,PetscObject o4)
884: {
885: PetscStageLog stageLog;
886: PetscEventRegLog eventRegLog;
887: PetscEventPerfLog eventPerfLog = NULL;
888: PetscLogDouble cur_time;
889: int stage,err;
890: PetscMPIInt rank;
891: PetscErrorCode ierr;
894: petsc_tracelevel--;
895: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
896: PetscLogGetStageLog(&stageLog);
897: PetscStageLogGetCurrent(stageLog,&stage);
898: PetscStageLogGetEventRegLog(stageLog,&eventRegLog);
899: PetscStageLogGetEventPerfLog(stageLog,stage,&eventPerfLog);
900: /* Check for double counting */
901: eventPerfLog->eventInfo[event].depth--;
902: if (eventPerfLog->eventInfo[event].depth > 0) return(0);
903: 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");
905: /* Log performance info */
906: if (petsc_tracelevel) {
907: PetscStrncpy(petsc_tracespace,petsc_traceblanks,2*petsc_tracelevel);
908: }
909: petsc_tracespace[2*petsc_tracelevel] = 0;
910: PetscTime(&cur_time);
911: PetscFPrintf(PETSC_COMM_SELF,petsc_tracefile,"%s[%d] %g Event end: %s\n",petsc_tracespace,rank,cur_time-petsc_tracetime,eventRegLog->eventInfo[event].name);
912: err = fflush(petsc_tracefile);
913: if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
914: return(0);
915: }
917: /*@C
918: PetscLogEventSetDof - Set the nth number of degrees of freedom associated with this event
920: Not Collective
922: Input Parameters:
923: + event - The event id to log
924: . n - The dof index, in [0, 8)
925: - dof - The number of dofs
927: Database Options:
928: . -log_view - Activates log summary
930: Note: This is to enable logging of convergence
932: Level: developer
934: .seealso: PetscLogEventSetError(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
935: @*/
936: PetscErrorCode PetscLogEventSetDof(PetscLogEvent event, PetscInt n, PetscLogDouble dof)
937: {
938: PetscStageLog stageLog;
939: PetscEventPerfLog eventLog = NULL;
940: int stage;
941: PetscErrorCode ierr;
944: if ((n < 0) || (n > 7)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Error index %D is not in [0, 8)", n);
945: PetscLogGetStageLog(&stageLog);
946: PetscStageLogGetCurrent(stageLog,&stage);
947: PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
948: eventLog->eventInfo[event].dof[n] = dof;
949: return(0);
950: }
952: /*@C
953: PetscLogEventSetError - Set the nth error associated with this event
955: Not Collective
957: Input Parameters:
958: + event - The event id to log
959: . n - The error index, in [0, 8)
960: - error - The error
962: Database Options:
963: . -log_view - Activates log summary
965: Note: This is to enable logging of convergence, and enable users to interpret the errors as they wish. For example,
966: as different norms, or as errors for different fields
968: Level: developer
970: .seealso: PetscLogEventSetDof(), PetscEventRegLogRegister(), PetscStageLogGetEventLog()
971: @*/
972: PetscErrorCode PetscLogEventSetError(PetscLogEvent event, PetscInt n, PetscLogDouble error)
973: {
974: PetscStageLog stageLog;
975: PetscEventPerfLog eventLog = NULL;
976: int stage;
977: PetscErrorCode ierr;
980: if ((n < 0) || (n > 7)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Error index %D is not in [0, 8)", n);
981: PetscLogGetStageLog(&stageLog);
982: PetscStageLogGetCurrent(stageLog,&stage);
983: PetscStageLogGetEventPerfLog(stageLog,stage,&eventLog);
984: eventLog->eventInfo[event].errors[n] = error;
985: return(0);
986: }