Actual source code: stagelog.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: PetscStageLog petsc_stageLog = NULL;
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 Note:
23: Inline since called for EACH `PetscEventLogBeginDefault()` and `PetscEventLogEndDefault()`
25: .seealso: `PetscStageLogCreate()`
26: @*/
27: PetscErrorCode PetscLogGetStageLog(PetscStageLog *stageLog)
28: {
30: if (!petsc_stageLog) {
31: fprintf(stderr, "PETSC ERROR: Logging has not been enabled.\nYou might have forgotten to call PetscInitialize().\n");
32: PETSCABORT(MPI_COMM_WORLD, PETSC_ERR_SUP);
33: }
34: *stageLog = petsc_stageLog;
35: return 0;
36: }
38: /*@C
39: PetscStageLogGetCurrent - This function returns the stage from the top of the stack.
41: Not Collective
43: Input Parameter:
44: . stageLog - The `PetscStageLog`
46: Output Parameter:
47: . stage - The current stage
49: Note:
50: If no stage is currently active, stage is set to -1.
52: Level: developer
54: Developer Note:
55: Inline since called for EACH `PetscEventLogBeginDefault()` and `PetscEventLogEndDefault()`
57: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscLogGetStageLog()`
58: @*/
59: PetscErrorCode PetscStageLogGetCurrent(PetscStageLog stageLog, int *stage)
60: {
61: PetscBool empty;
63: PetscIntStackEmpty(stageLog->stack, &empty);
64: if (empty) {
65: *stage = -1;
66: } else {
67: PetscIntStackTop(stageLog->stack, stage);
68: }
70: return 0;
71: }
73: /*@C
74: PetscStageLogGetEventPerfLog - This function returns the `PetscEventPerfLog` for the given stage.
76: Not Collective
78: Input Parameters:
79: + stageLog - The `PetscStageLog`
80: - stage - The stage
82: Output Parameter:
83: . eventLog - The `PetscEventPerfLog`
85: Level: developer
87: Developer Note:
88: Inline since called for EACH `PetscEventLogBeginDefault()` and `PetscEventLogEndDefault()`
90: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscLogGetStageLog()`
91: @*/
92: PetscErrorCode PetscStageLogGetEventPerfLog(PetscStageLog stageLog, int stage, PetscEventPerfLog *eventLog)
93: {
96: *eventLog = stageLog->stageInfo[stage].eventLog;
97: return 0;
98: }
100: /*@C
101: PetscStageInfoDestroy - This destroys a `PetscStageInfo` object.
103: Not collective
105: Input Parameter:
106: . stageInfo - The `PetscStageInfo`
108: Level: developer
110: .seealso: `PetscStageLogCreate()`
111: @*/
112: PetscErrorCode PetscStageInfoDestroy(PetscStageInfo *stageInfo)
113: {
114: PetscFree(stageInfo->name);
115: PetscEventPerfLogDestroy(stageInfo->eventLog);
116: PetscClassPerfLogDestroy(stageInfo->classLog);
117: return 0;
118: }
120: /*@C
121: PetscStageLogDestroy - This destroys a `PetscStageLog` object.
123: Not collective
125: Input Parameter:
126: . stageLog - The `PetscStageLog`
128: Level: developer
130: .seealso: `PetscStageLogCreate()`
131: @*/
132: PetscErrorCode PetscStageLogDestroy(PetscStageLog stageLog)
133: {
134: int stage;
136: if (!stageLog) return 0;
137: PetscIntStackDestroy(stageLog->stack);
138: PetscEventRegLogDestroy(stageLog->eventLog);
139: PetscClassRegLogDestroy(stageLog->classLog);
140: for (stage = 0; stage < stageLog->numStages; stage++) PetscStageInfoDestroy(&stageLog->stageInfo[stage]);
141: PetscFree(stageLog->stageInfo);
142: PetscFree(stageLog);
143: return 0;
144: }
146: /*@C
147: PetscStageLogRegister - Registers a stage name for logging operations in an application code.
149: Not Collective
151: Input Parameters:
152: + stageLog - The `PetscStageLog`
153: - sname - the name to associate with that stage
155: Output Parameter:
156: . stage - The stage index
158: Level: developer
160: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscStageLogCreate()`
161: @*/
162: PetscErrorCode PetscStageLogRegister(PetscStageLog stageLog, const char sname[], int *stage)
163: {
164: PetscStageInfo *stageInfo;
165: int s;
169: /* Check stage already registered */
170: for (s = 0; s < stageLog->numStages; ++s) {
171: PetscBool same;
173: PetscStrcmp(stageLog->stageInfo[s].name, sname, &same);
175: }
176: /* Create new stage */
177: s = stageLog->numStages++;
178: if (stageLog->numStages > stageLog->maxStages) {
179: PetscMalloc1(stageLog->maxStages * 2, &stageInfo);
180: PetscArraycpy(stageInfo, stageLog->stageInfo, stageLog->maxStages);
181: PetscFree(stageLog->stageInfo);
182: stageLog->stageInfo = stageInfo;
183: stageLog->maxStages *= 2;
184: }
185: /* Setup new stage info */
186: stageInfo = &stageLog->stageInfo[s];
187: PetscMemzero(stageInfo, sizeof(PetscStageInfo));
188: PetscStrallocpy(sname, &stageInfo->name);
189: stageInfo->used = PETSC_FALSE;
190: stageInfo->perfInfo.active = PETSC_TRUE;
191: stageInfo->perfInfo.visible = PETSC_TRUE;
192: PetscEventPerfLogCreate(&stageInfo->eventLog);
193: PetscClassPerfLogCreate(&stageInfo->classLog);
194: *stage = s;
195: return 0;
196: }
198: /*@C
199: PetscStageLogPush - This function pushes a stage on the stack.
201: Not Collective
203: Input Parameters:
204: + stageLog - The `PetscStageLog`
205: - stage - The stage to log
207: Database Options:
208: . -log_view - Activates logging
210: Usage:
211: If the option -log_view is used to run the program containing the
212: following code, then 2 sets of summary data will be printed during
213: `PetscFinalize()`.
214: .vb
215: PetscInitialize(int *argc,char ***args,0,0);
216: [stage 0 of code]
217: PetscStageLogPush(stageLog,1);
218: [stage 1 of code]
219: PetscStageLogPop(stageLog);
220: PetscBarrier(...);
221: [more stage 0 of code]
222: PetscFinalize();
223: .ve
225: Note;
226: Use `PetscLogStageRegister()` to register a stage. All previous stages are
227: accumulating time and flops, but events will only be logged in this stage.
229: Level: developer
231: .seealso: `PetscStageLogPop()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
232: @*/
233: PetscErrorCode PetscStageLogPush(PetscStageLog stageLog, int stage)
234: {
235: int curStage = 0;
236: PetscBool empty;
240: /* Record flops/time of previous stage */
241: PetscIntStackEmpty(stageLog->stack, &empty);
242: if (!empty) {
243: PetscIntStackTop(stageLog->stack, &curStage);
244: if (stageLog->stageInfo[curStage].perfInfo.active) {
245: PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time);
246: stageLog->stageInfo[curStage].perfInfo.flops += petsc_TotalFlops;
247: stageLog->stageInfo[curStage].perfInfo.numMessages += petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
248: stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
249: stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
250: }
251: }
252: /* Activate the stage */
253: PetscIntStackPush(stageLog->stack, stage);
255: stageLog->stageInfo[stage].used = PETSC_TRUE;
256: stageLog->stageInfo[stage].perfInfo.count++;
257: stageLog->curStage = stage;
258: /* Subtract current quantities so that we obtain the difference when we pop */
259: if (stageLog->stageInfo[stage].perfInfo.active) {
260: PetscTimeSubtract(&stageLog->stageInfo[stage].perfInfo.time);
261: stageLog->stageInfo[stage].perfInfo.flops -= petsc_TotalFlops;
262: stageLog->stageInfo[stage].perfInfo.numMessages -= petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
263: stageLog->stageInfo[stage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
264: stageLog->stageInfo[stage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
265: }
266: return 0;
267: }
269: /*@C
270: PetscStageLogPop - This function pops a stage from the stack.
272: Not Collective
274: Input Parameter:
275: . stageLog - The `PetscStageLog`
277: Usage:
278: If the option -log_view is used to run the program containing the
279: following code, then 2 sets of summary data will be printed during
280: PetscFinalize().
281: .vb
282: PetscInitialize(int *argc,char ***args,0,0);
283: [stage 0 of code]
284: PetscStageLogPush(stageLog,1);
285: [stage 1 of code]
286: PetscStageLogPop(stageLog);
287: PetscBarrier(...);
288: [more stage 0 of code]
289: PetscFinalize();
290: .ve
292: Note:
293: Use `PetscStageLogRegister()` to register a stage.
295: Level: developer
297: .seealso: `PetscStageLogPush()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
298: @*/
299: PetscErrorCode PetscStageLogPop(PetscStageLog stageLog)
300: {
301: int curStage;
302: PetscBool empty;
304: /* Record flops/time of current stage */
305: PetscIntStackPop(stageLog->stack, &curStage);
306: if (stageLog->stageInfo[curStage].perfInfo.active) {
307: PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time);
308: stageLog->stageInfo[curStage].perfInfo.flops += petsc_TotalFlops;
309: stageLog->stageInfo[curStage].perfInfo.numMessages += petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
310: stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
311: stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
312: }
313: PetscIntStackEmpty(stageLog->stack, &empty);
314: if (!empty) {
315: /* Subtract current quantities so that we obtain the difference when we pop */
316: PetscIntStackTop(stageLog->stack, &curStage);
317: if (stageLog->stageInfo[curStage].perfInfo.active) {
318: PetscTimeSubtract(&stageLog->stageInfo[curStage].perfInfo.time);
319: stageLog->stageInfo[curStage].perfInfo.flops -= petsc_TotalFlops;
320: stageLog->stageInfo[curStage].perfInfo.numMessages -= petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
321: stageLog->stageInfo[curStage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
322: stageLog->stageInfo[curStage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
323: }
324: stageLog->curStage = curStage;
325: } else stageLog->curStage = -1;
326: return 0;
327: }
329: /*@C
330: PetscStageLogGetClassRegLog - This function returns the PetscClassRegLog for the given stage.
332: Not Collective
334: Input Parameters:
335: . stageLog - The `PetscStageLog`
337: Output Parameter:
338: . classLog - The `PetscClassRegLog`
340: Level: developer
342: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscLogGetStageLog()`
343: @*/
344: PetscErrorCode PetscStageLogGetClassRegLog(PetscStageLog stageLog, PetscClassRegLog *classLog)
345: {
347: *classLog = stageLog->classLog;
348: return 0;
349: }
351: /*@C
352: PetscStageLogGetEventRegLog - This function returns the `PetscEventRegLog`.
354: Not Collective
356: Input Parameters:
357: . stageLog - The `PetscStageLog`
359: Output Parameter:
360: . eventLog - The `PetscEventRegLog`
362: Level: developer
364: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscLogGetStageLog()`
365: @*/
366: PetscErrorCode PetscStageLogGetEventRegLog(PetscStageLog stageLog, PetscEventRegLog *eventLog)
367: {
369: *eventLog = stageLog->eventLog;
370: return 0;
371: }
373: /*@C
374: PetscStageLogGetClassPerfLog - This function returns the `PetscClassPerfLog` for the given stage.
376: Not Collective
378: Input Parameters:
379: + stageLog - The `PetscStageLog`
380: - stage - The stage
382: Output Parameter:
383: . classLog - The `PetscClassPerfLog`
385: Level: developer
387: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscLogGetStageLog()`
388: @*/
389: PetscErrorCode PetscStageLogGetClassPerfLog(PetscStageLog stageLog, int stage, PetscClassPerfLog *classLog)
390: {
393: *classLog = stageLog->stageInfo[stage].classLog;
394: return 0;
395: }
397: /*@C
398: PetscStageLogSetActive - This function determines whether events will be logged during this state.
400: Not Collective
402: Input Parameters:
403: + stageLog - The `PetscStageLog`
404: . stage - The stage to log
405: - isActive - The activity flag, `PETSC_TRUE` for logging, otherwise `PETSC_FALSE` (default is `PETSC_TRUE`)
407: Level: developer
409: .seealso: `PetscStageLogGetActive()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
410: @*/
411: PetscErrorCode PetscStageLogSetActive(PetscStageLog stageLog, int stage, PetscBool isActive)
412: {
414: stageLog->stageInfo[stage].perfInfo.active = isActive;
415: return 0;
416: }
418: /*@C
419: PetscStageLogGetActive - This function returns whether events will be logged suring this stage.
421: Not Collective
423: Input Parameters:
424: + stageLog - The `PetscStageLog`
425: - stage - The stage to log
427: Output Parameter:
428: . isActive - The activity flag, `PETSC_TRUE` for logging, otherwise `PETSC_FALSE` (default is `PETSC_TRUE`)
430: Level: developer
432: .seealso: `PetscStageLogSetActive()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
433: @*/
434: PetscErrorCode PetscStageLogGetActive(PetscStageLog stageLog, int stage, PetscBool *isActive)
435: {
438: *isActive = stageLog->stageInfo[stage].perfInfo.active;
439: return 0;
440: }
442: /*@C
443: PetscStageLogSetVisible - This function determines whether a stage is printed during `PetscLogView()`
445: Not Collective
447: Input Parameters:
448: + stageLog - The `PetscStageLog`
449: . stage - The stage to log
450: - isVisible - The visibility flag, `PETSC_TRUE` for printing, otherwise `PETSC_FALSE` (default is `PETSC_TRUE`)
452: Database Options:
453: . -log_view - Activates log summary
455: Level: developer
457: .seealso: `PetscStageLogGetVisible()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
458: @*/
459: PetscErrorCode PetscStageLogSetVisible(PetscStageLog stageLog, int stage, PetscBool isVisible)
460: {
462: stageLog->stageInfo[stage].perfInfo.visible = isVisible;
463: return 0;
464: }
466: /*@C
467: PetscStageLogGetVisible - This function returns whether a stage is printed during `PetscLogView()`
469: Not Collective
471: Input Parameters:
472: + stageLog - The `PetscStageLog`
473: - stage - The stage to log
475: Output Parameter:
476: . isVisible - The visibility flag, `PETSC_TRUE` for printing, otherwise `PETSC_FALSE` (default is `PETSC_TRUE`)
478: Database Options:
479: . -log_view - Activates log summary
481: Level: developer
483: .seealso: `PetscStageLogSetVisible()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
484: @*/
485: PetscErrorCode PetscStageLogGetVisible(PetscStageLog stageLog, int stage, PetscBool *isVisible)
486: {
489: *isVisible = stageLog->stageInfo[stage].perfInfo.visible;
490: return 0;
491: }
493: /*@C
494: PetscStageLogGetStage - This function returns the stage id given the stage name.
496: Not Collective
498: Input Parameters:
499: + stageLog - The `PetscStageLog`
500: - name - The stage name
502: Output Parameter:
503: . stage - The stage id, or -1 if it does not exist
505: Level: developer
507: .seealso: `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
508: @*/
509: PetscErrorCode PetscStageLogGetStage(PetscStageLog stageLog, const char name[], PetscLogStage *stage)
510: {
511: PetscBool match;
512: int s;
516: *stage = -1;
517: for (s = 0; s < stageLog->numStages; s++) {
518: PetscStrcasecmp(stageLog->stageInfo[s].name, name, &match);
519: if (match) {
520: *stage = s;
521: break;
522: }
523: }
524: return 0;
525: }
527: /*@C
528: PetscStageLogCreate - This creates a `PetscStageLog` object.
530: Not collective
532: Output Parameter:
533: . stageLog - The `PetscStageLog`
535: Level: developer
537: .seealso: `PetscStageLogCreate()`
538: @*/
539: PetscErrorCode PetscStageLogCreate(PetscStageLog *stageLog)
540: {
541: PetscStageLog l;
543: PetscNew(&l);
545: l->numStages = 0;
546: l->maxStages = 10;
547: l->curStage = -1;
549: PetscIntStackCreate(&l->stack);
550: PetscMalloc1(l->maxStages, &l->stageInfo);
551: PetscEventRegLogCreate(&l->eventLog);
552: PetscClassRegLogCreate(&l->classLog);
554: *stageLog = l;
555: return 0;
556: }