Actual source code: stagelog.c
petsc-3.8.4 2018-03-24
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: PetscStageLog petsc_stageLog = 0;
12: /*@C
13: PetscLogGetStageLog - This function returns the default stage logging object.
15: Not collective
17: Output Parameter:
18: . stageLog - The default PetscStageLog
20: Level: developer
22: Developer Notes: Inline since called for EACH PetscEventLogBeginDefault() and PetscEventLogEndDefault()
24: .keywords: log, stage
25: .seealso: PetscStageLogCreate()
26: @*/
27: PetscErrorCode PetscLogGetStageLog(PetscStageLog *stageLog)
28: {
31: if (!petsc_stageLog) {
32: fprintf(stderr, "PETSC ERROR: Logging has not been enabled.\nYou might have forgotten to call PetscInitialize().\n");
33: MPI_Abort(MPI_COMM_WORLD, PETSC_ERR_SUP);
34: }
35: *stageLog = petsc_stageLog;
36: return(0);
37: }
39: /*@C
40: PetscStageLogGetCurrent - This function returns the stage from the top of the stack.
42: Not Collective
44: Input Parameter:
45: . stageLog - The PetscStageLog
47: Output Parameter:
48: . stage - The current stage
50: Notes:
51: If no stage is currently active, stage is set to -1.
53: Level: developer
55: Developer Notes: Inline since called for EACH PetscEventLogBeginDefault() and PetscEventLogEndDefault()
57: .keywords: log, stage
58: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
59: @*/
60: PetscErrorCode PetscStageLogGetCurrent(PetscStageLog stageLog, int *stage)
61: {
62: PetscBool empty;
66: PetscIntStackEmpty(stageLog->stack, &empty);
67: if (empty) {
68: *stage = -1;
69: } else {
70: PetscIntStackTop(stageLog->stack, stage);
71: }
72: #ifdef PETSC_USE_DEBUG
73: if (*stage != stageLog->curStage) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Inconsistency in stage log: stage %d should be %d", *stage, stageLog->curStage);
74: #endif
75: return(0);
76: }
78: /*@C
79: PetscStageLogGetEventPerfLog - This function returns the PetscEventPerfLog for the given stage.
81: Not Collective
83: Input Parameters:
84: + stageLog - The PetscStageLog
85: - stage - The stage
87: Output Parameter:
88: . eventLog - The PetscEventPerfLog
90: Level: developer
92: Developer Notes: Inline since called for EACH PetscEventLogBeginDefault() and PetscEventLogEndDefault()
94: .keywords: log, stage
95: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
96: @*/
97: PetscErrorCode PetscStageLogGetEventPerfLog(PetscStageLog stageLog, int stage, PetscEventPerfLog *eventLog)
98: {
101: if ((stage < 0) || (stage >= stageLog->numStages)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
102: *eventLog = stageLog->stageInfo[stage].eventLog;
103: return(0);
104: }
106: /*@C
107: PetscStageInfoDestroy - This destroys a PetscStageInfo object.
109: Not collective
111: Input Paramter:
112: . stageInfo - The PetscStageInfo
114: Level: developer
116: .keywords: log, stage, destroy
117: .seealso: PetscStageLogCreate()
118: @*/
119: PetscErrorCode PetscStageInfoDestroy(PetscStageInfo *stageInfo)
120: {
124: PetscFree(stageInfo->name);
125: PetscEventPerfLogDestroy(stageInfo->eventLog);
126: PetscClassPerfLogDestroy(stageInfo->classLog);
127: return(0);
128: }
130: /*@C
131: PetscStageLogDestroy - This destroys a PetscStageLog object.
133: Not collective
135: Input Paramter:
136: . stageLog - The PetscStageLog
138: Level: developer
140: .keywords: log, stage, destroy
141: .seealso: PetscStageLogCreate()
142: @*/
143: PetscErrorCode PetscStageLogDestroy(PetscStageLog stageLog)
144: {
145: int stage;
149: if (!stageLog) return(0);
150: PetscIntStackDestroy(stageLog->stack);
151: PetscEventRegLogDestroy(stageLog->eventLog);
152: PetscClassRegLogDestroy(stageLog->classLog);
153: for (stage = 0; stage < stageLog->numStages; stage++) {
154: PetscStageInfoDestroy(&stageLog->stageInfo[stage]);
155: }
156: PetscFree(stageLog->stageInfo);
157: PetscFree(stageLog);
158: return(0);
159: }
161: /*@C
162: PetscStageLogRegister - Registers a stage name for logging operations in an application code.
164: Not Collective
166: Input Parameter:
167: + stageLog - The PetscStageLog
168: - sname - the name to associate with that stage
170: Output Parameter:
171: . stage - The stage index
173: Level: developer
175: .keywords: log, stage, register
176: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscStageLogCreate()
177: @*/
178: PetscErrorCode PetscStageLogRegister(PetscStageLog stageLog, const char sname[], int *stage)
179: {
180: PetscStageInfo *stageInfo;
181: char *str;
182: int s, st;
188: for (st = 0; st < stageLog->numStages; ++st) {
189: PetscBool same;
191: PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);
192: if (same) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Duplicate stage name given: %s", sname);
193: }
194: s = stageLog->numStages++;
195: if (stageLog->numStages > stageLog->maxStages) {
196: PetscMalloc1(stageLog->maxStages*2, &stageInfo);
197: PetscMemcpy(stageInfo, stageLog->stageInfo, stageLog->maxStages * sizeof(PetscStageInfo));
198: PetscFree(stageLog->stageInfo);
200: stageLog->stageInfo = stageInfo;
201: stageLog->maxStages *= 2;
202: }
203: /* Setup stage */
204: PetscStrallocpy(sname, &str);
206: stageLog->stageInfo[s].name = str;
207: stageLog->stageInfo[s].used = PETSC_FALSE;
208: stageLog->stageInfo[s].perfInfo.active = PETSC_TRUE;
209: stageLog->stageInfo[s].perfInfo.visible = PETSC_TRUE;
210: stageLog->stageInfo[s].perfInfo.count = 0;
211: stageLog->stageInfo[s].perfInfo.flops = 0.0;
212: stageLog->stageInfo[s].perfInfo.time = 0.0;
213: stageLog->stageInfo[s].perfInfo.numMessages = 0.0;
214: stageLog->stageInfo[s].perfInfo.messageLength = 0.0;
215: stageLog->stageInfo[s].perfInfo.numReductions = 0.0;
217: PetscEventPerfLogCreate(&stageLog->stageInfo[s].eventLog);
218: PetscClassPerfLogCreate(&stageLog->stageInfo[s].classLog);
219: *stage = s;
220: return(0);
221: }
223: /*@C
224: PetscStageLogPush - This function pushes a stage on the stack.
226: Not Collective
228: Input Parameters:
229: + stageLog - The PetscStageLog
230: - stage - The stage to log
232: Database Options:
233: . -log_view - Activates logging
235: Usage:
236: If the option -log_sumary is used to run the program containing the
237: following code, then 2 sets of summary data will be printed during
238: PetscFinalize().
239: .vb
240: PetscInitialize(int *argc,char ***args,0,0);
241: [stage 0 of code]
242: PetscStageLogPush(stageLog,1);
243: [stage 1 of code]
244: PetscStageLogPop(stageLog);
245: PetscBarrier(...);
246: [more stage 0 of code]
247: PetscFinalize();
248: .ve
250: Notes:
251: Use PetscLogStageRegister() to register a stage. All previous stages are
252: accumulating time and flops, but events will only be logged in this stage.
254: Level: developer
256: .keywords: log, push, stage
257: .seealso: PetscStageLogPop(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
258: @*/
259: PetscErrorCode PetscStageLogPush(PetscStageLog stageLog, int stage)
260: {
261: int curStage = 0;
262: PetscBool empty;
266: if ((stage < 0) || (stage >= stageLog->numStages)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
268: /* Record flops/time of previous stage */
269: PetscIntStackEmpty(stageLog->stack, &empty);
270: if (!empty) {
271: PetscIntStackTop(stageLog->stack, &curStage);
272: if (stageLog->stageInfo[curStage].perfInfo.active) {
273: PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time);
274: stageLog->stageInfo[curStage].perfInfo.flops += petsc_TotalFlops;
275: stageLog->stageInfo[curStage].perfInfo.numMessages += petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
276: stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
277: stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
278: }
279: }
280: /* Activate the stage */
281: PetscIntStackPush(stageLog->stack, stage);
283: stageLog->stageInfo[stage].used = PETSC_TRUE;
284: stageLog->stageInfo[stage].perfInfo.count++;
285: stageLog->curStage = stage;
286: /* Subtract current quantities so that we obtain the difference when we pop */
287: if (stageLog->stageInfo[stage].perfInfo.active) {
288: PetscTimeSubtract(&stageLog->stageInfo[stage].perfInfo.time);
289: stageLog->stageInfo[stage].perfInfo.flops -= petsc_TotalFlops;
290: stageLog->stageInfo[stage].perfInfo.numMessages -= petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
291: stageLog->stageInfo[stage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
292: stageLog->stageInfo[stage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
293: }
294: return(0);
295: }
297: /*@C
298: PetscStageLogPop - This function pops a stage from the stack.
300: Not Collective
302: Input Parameter:
303: . stageLog - The PetscStageLog
305: Usage:
306: If the option -log_sumary is used to run the program containing the
307: following code, then 2 sets of summary data will be printed during
308: PetscFinalize().
309: .vb
310: PetscInitialize(int *argc,char ***args,0,0);
311: [stage 0 of code]
312: PetscStageLogPush(stageLog,1);
313: [stage 1 of code]
314: PetscStageLogPop(stageLog);
315: PetscBarrier(...);
316: [more stage 0 of code]
317: PetscFinalize();
318: .ve
320: Notes:
321: Use PetscStageLogRegister() to register a stage.
323: Level: developer
325: .keywords: log, pop, stage
326: .seealso: PetscStageLogPush(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
327: @*/
328: PetscErrorCode PetscStageLogPop(PetscStageLog stageLog)
329: {
330: int curStage;
331: PetscBool empty;
335: /* Record flops/time of current stage */
336: PetscIntStackPop(stageLog->stack, &curStage);
337: if (stageLog->stageInfo[curStage].perfInfo.active) {
338: PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time);
339: stageLog->stageInfo[curStage].perfInfo.flops += petsc_TotalFlops;
340: stageLog->stageInfo[curStage].perfInfo.numMessages += petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
341: stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
342: stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
343: }
344: PetscIntStackEmpty(stageLog->stack, &empty);
345: if (!empty) {
346: /* Subtract current quantities so that we obtain the difference when we pop */
347: PetscIntStackTop(stageLog->stack, &curStage);
348: if (stageLog->stageInfo[curStage].perfInfo.active) {
349: PetscTimeSubtract(&stageLog->stageInfo[curStage].perfInfo.time);
350: stageLog->stageInfo[curStage].perfInfo.flops -= petsc_TotalFlops;
351: stageLog->stageInfo[curStage].perfInfo.numMessages -= petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
352: stageLog->stageInfo[curStage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
353: stageLog->stageInfo[curStage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
354: }
355: stageLog->curStage = curStage;
356: } else stageLog->curStage = -1;
357: return(0);
358: }
361: /*@C
362: PetscStageLogGetClassRegLog - This function returns the PetscClassRegLog for the given stage.
364: Not Collective
366: Input Parameters:
367: . stageLog - The PetscStageLog
369: Output Parameter:
370: . classLog - The PetscClassRegLog
372: Level: developer
374: .keywords: log, stage
375: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
376: @*/
377: PetscErrorCode PetscStageLogGetClassRegLog(PetscStageLog stageLog, PetscClassRegLog *classLog)
378: {
381: *classLog = stageLog->classLog;
382: return(0);
383: }
385: /*@C
386: PetscStageLogGetEventRegLog - This function returns the PetscEventRegLog.
388: Not Collective
390: Input Parameters:
391: . stageLog - The PetscStageLog
393: Output Parameter:
394: . eventLog - The PetscEventRegLog
396: Level: developer
398: .keywords: log, stage
399: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
400: @*/
401: PetscErrorCode PetscStageLogGetEventRegLog(PetscStageLog stageLog, PetscEventRegLog *eventLog)
402: {
405: *eventLog = stageLog->eventLog;
406: return(0);
407: }
409: /*@C
410: PetscStageLogGetClassPerfLog - This function returns the PetscClassPerfLog for the given stage.
412: Not Collective
414: Input Parameters:
415: + stageLog - The PetscStageLog
416: - stage - The stage
418: Output Parameter:
419: . classLog - The PetscClassPerfLog
421: Level: developer
423: .keywords: log, stage
424: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
425: @*/
426: PetscErrorCode PetscStageLogGetClassPerfLog(PetscStageLog stageLog, int stage, PetscClassPerfLog *classLog)
427: {
430: if ((stage < 0) || (stage >= stageLog->numStages)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
431: *classLog = stageLog->stageInfo[stage].classLog;
432: return(0);
433: }
436: /*@C
437: PetscStageLogSetActive - This function determines whether events will be logged during this state.
439: Not Collective
441: Input Parameters:
442: + stageLog - The PetscStageLog
443: . stage - The stage to log
444: - isActive - The activity flag, PETSC_TRUE for logging, otherwise PETSC_FALSE (default is PETSC_TRUE)
446: Level: developer
448: .keywords: log, active, stage
449: .seealso: PetscStageLogGetActive(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
450: @*/
451: PetscErrorCode PetscStageLogSetActive(PetscStageLog stageLog, int stage, PetscBool isActive)
452: {
454: if ((stage < 0) || (stage >= stageLog->numStages)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
455: stageLog->stageInfo[stage].perfInfo.active = isActive;
456: return(0);
457: }
459: /*@C
460: PetscStageLogGetActive - This function returns whether events will be logged suring this stage.
462: Not Collective
464: Input Parameters:
465: + stageLog - The PetscStageLog
466: - stage - The stage to log
468: Output Parameter:
469: . isActive - The activity flag, PETSC_TRUE for logging, otherwise PETSC_FALSE (default is PETSC_TRUE)
471: Level: developer
473: .keywords: log, visible, stage
474: .seealso: PetscStageLogSetActive(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
475: @*/
476: PetscErrorCode PetscStageLogGetActive(PetscStageLog stageLog, int stage, PetscBool *isActive)
477: {
479: if ((stage < 0) || (stage >= stageLog->numStages)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
481: *isActive = stageLog->stageInfo[stage].perfInfo.active;
482: return(0);
483: }
485: /*@C
486: PetscStageLogSetVisible - This function determines whether a stage is printed during PetscLogView()
488: Not Collective
490: Input Parameters:
491: + stageLog - The PetscStageLog
492: . stage - The stage to log
493: - isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)
495: Database Options:
496: . -log_view - Activates log summary
498: Level: developer
500: .keywords: log, visible, stage
501: .seealso: PetscStageLogGetVisible(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
502: @*/
503: PetscErrorCode PetscStageLogSetVisible(PetscStageLog stageLog, int stage, PetscBool isVisible)
504: {
506: if ((stage < 0) || (stage >= stageLog->numStages)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
507: stageLog->stageInfo[stage].perfInfo.visible = isVisible;
508: return(0);
509: }
511: /*@C
512: PetscStageLogGetVisible - This function returns whether a stage is printed during PetscLogView()
514: Not Collective
516: Input Parameters:
517: + stageLog - The PetscStageLog
518: - stage - The stage to log
520: Output Parameter:
521: . isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)
523: Database Options:
524: . -log_view - Activates log summary
526: Level: developer
528: .keywords: log, visible, stage
529: .seealso: PetscStageLogSetVisible(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
530: @*/
531: PetscErrorCode PetscStageLogGetVisible(PetscStageLog stageLog, int stage, PetscBool *isVisible)
532: {
534: if ((stage < 0) || (stage >= stageLog->numStages)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
536: *isVisible = stageLog->stageInfo[stage].perfInfo.visible;
537: return(0);
538: }
540: /*@C
541: PetscStageLogGetStage - This function returns the stage id given the stage name.
543: Not Collective
545: Input Parameters:
546: + stageLog - The PetscStageLog
547: - name - The stage name
549: Output Parameter:
550: . stage - The stage id, or -1 if it does not exist
552: Level: developer
554: .keywords: log, stage
555: .seealso: PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
556: @*/
557: PetscErrorCode PetscStageLogGetStage(PetscStageLog stageLog, const char name[], PetscLogStage *stage)
558: {
559: PetscBool match;
560: int s;
566: *stage = -1;
567: for (s = 0; s < stageLog->numStages; s++) {
568: PetscStrcasecmp(stageLog->stageInfo[s].name, name, &match);
569: if (match) {
570: *stage = s;
571: break;
572: }
573: }
574: return(0);
575: }
577: /*@C
578: PetscStageLogCreate - This creates a PetscStageLog object.
580: Not collective
582: Input Parameter:
583: . stageLog - The PetscStageLog
585: Level: developer
587: .keywords: log, stage, create
588: .seealso: PetscStageLogCreate()
589: @*/
590: PetscErrorCode PetscStageLogCreate(PetscStageLog *stageLog)
591: {
592: PetscStageLog l;
596: PetscNew(&l);
598: l->numStages = 0;
599: l->maxStages = 10;
600: l->curStage = -1;
602: PetscIntStackCreate(&l->stack);
603: PetscMalloc1(l->maxStages, &l->stageInfo);
604: PetscEventRegLogCreate(&l->eventLog);
605: PetscClassRegLogCreate(&l->classLog);
607: *stageLog = l;
608: return(0);
609: }