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: {
29: PetscFunctionBegin;
31: if (!petsc_stageLog) {
32: fprintf(stderr, "PETSC ERROR: Logging has not been enabled.\nYou might have forgotten to call PetscInitialize().\n");
33: PETSCABORT(MPI_COMM_WORLD, PETSC_ERR_SUP);
34: }
35: *stageLog = petsc_stageLog;
36: PetscFunctionReturn(PETSC_SUCCESS);
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: Note:
51: If no stage is currently active, stage is set to -1.
53: Level: developer
55: Developer Note:
56: Inline since called for EACH `PetscEventLogBeginDefault()` and `PetscEventLogEndDefault()`
58: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscLogGetStageLog()`
59: @*/
60: PetscErrorCode PetscStageLogGetCurrent(PetscStageLog stageLog, int *stage)
61: {
62: PetscBool empty;
64: PetscFunctionBegin;
65: PetscCall(PetscIntStackEmpty(stageLog->stack, &empty));
66: if (empty) {
67: *stage = -1;
68: } else {
69: PetscCall(PetscIntStackTop(stageLog->stack, stage));
70: }
71: PetscCheck(*stage == stageLog->curStage, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in stage log: stage %d should be %d", *stage, stageLog->curStage);
72: PetscFunctionReturn(PETSC_SUCCESS);
73: }
75: /*@C
76: PetscStageLogGetEventPerfLog - This function returns the `PetscEventPerfLog` for the given stage.
78: Not Collective
80: Input Parameters:
81: + stageLog - The `PetscStageLog`
82: - stage - The stage
84: Output Parameter:
85: . eventLog - The `PetscEventPerfLog`
87: Level: developer
89: Developer Note:
90: Inline since called for EACH `PetscEventLogBeginDefault()` and `PetscEventLogEndDefault()`
92: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscLogGetStageLog()`
93: @*/
94: PetscErrorCode PetscStageLogGetEventPerfLog(PetscStageLog stageLog, int stage, PetscEventPerfLog *eventLog)
95: {
96: PetscFunctionBegin;
98: PetscCheck(!(stage < 0) && !(stage >= stageLog->numStages), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
99: *eventLog = stageLog->stageInfo[stage].eventLog;
100: PetscFunctionReturn(PETSC_SUCCESS);
101: }
103: /*@C
104: PetscStageInfoDestroy - This destroys a `PetscStageInfo` object.
106: Not collective
108: Input Parameter:
109: . stageInfo - The `PetscStageInfo`
111: Level: developer
113: .seealso: `PetscStageLogCreate()`
114: @*/
115: PetscErrorCode PetscStageInfoDestroy(PetscStageInfo *stageInfo)
116: {
117: PetscFunctionBegin;
118: PetscCall(PetscFree(stageInfo->name));
119: PetscCall(PetscEventPerfLogDestroy(stageInfo->eventLog));
120: PetscCall(PetscClassPerfLogDestroy(stageInfo->classLog));
121: PetscFunctionReturn(PETSC_SUCCESS);
122: }
124: /*@C
125: PetscStageLogDestroy - This destroys a `PetscStageLog` object.
127: Not collective
129: Input Parameter:
130: . stageLog - The `PetscStageLog`
132: Level: developer
134: .seealso: `PetscStageLogCreate()`
135: @*/
136: PetscErrorCode PetscStageLogDestroy(PetscStageLog stageLog)
137: {
138: int stage;
140: PetscFunctionBegin;
141: if (!stageLog) PetscFunctionReturn(PETSC_SUCCESS);
142: PetscCall(PetscIntStackDestroy(stageLog->stack));
143: PetscCall(PetscEventRegLogDestroy(stageLog->eventLog));
144: PetscCall(PetscClassRegLogDestroy(stageLog->classLog));
145: for (stage = 0; stage < stageLog->numStages; stage++) PetscCall(PetscStageInfoDestroy(&stageLog->stageInfo[stage]));
146: PetscCall(PetscFree(stageLog->stageInfo));
147: PetscCall(PetscFree(stageLog));
148: PetscFunctionReturn(PETSC_SUCCESS);
149: }
151: /*@C
152: PetscStageLogRegister - Registers a stage name for logging operations in an application code.
154: Not Collective
156: Input Parameters:
157: + stageLog - The `PetscStageLog`
158: - sname - the name to associate with that stage
160: Output Parameter:
161: . stage - The stage index
163: Level: developer
165: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscStageLogCreate()`
166: @*/
167: PetscErrorCode PetscStageLogRegister(PetscStageLog stageLog, const char sname[], int *stage)
168: {
169: PetscStageInfo *stageInfo;
170: int s;
172: PetscFunctionBegin;
175: /* Check stage already registered */
176: for (s = 0; s < stageLog->numStages; ++s) {
177: PetscBool same;
179: PetscCall(PetscStrcmp(stageLog->stageInfo[s].name, sname, &same));
180: PetscCheck(!same, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Duplicate stage name given: %s", sname);
181: }
182: /* Create new stage */
183: s = stageLog->numStages++;
184: if (stageLog->numStages > stageLog->maxStages) {
185: PetscCall(PetscMalloc1(stageLog->maxStages * 2, &stageInfo));
186: PetscCall(PetscArraycpy(stageInfo, stageLog->stageInfo, stageLog->maxStages));
187: PetscCall(PetscFree(stageLog->stageInfo));
188: stageLog->stageInfo = stageInfo;
189: stageLog->maxStages *= 2;
190: }
191: /* Setup new stage info */
192: stageInfo = &stageLog->stageInfo[s];
193: PetscCall(PetscMemzero(stageInfo, sizeof(PetscStageInfo)));
194: PetscCall(PetscStrallocpy(sname, &stageInfo->name));
195: stageInfo->used = PETSC_FALSE;
196: stageInfo->perfInfo.active = PETSC_TRUE;
197: stageInfo->perfInfo.visible = PETSC_TRUE;
198: PetscCall(PetscEventPerfLogCreate(&stageInfo->eventLog));
199: PetscCall(PetscClassPerfLogCreate(&stageInfo->classLog));
200: *stage = s;
201: PetscFunctionReturn(PETSC_SUCCESS);
202: }
204: /*@C
205: PetscStageLogPush - This function pushes a stage on the stack.
207: Not Collective
209: Input Parameters:
210: + stageLog - The `PetscStageLog`
211: - stage - The stage to log
213: Options Database Key:
214: . -log_view - Activates logging
216: Usage:
217: If the option -log_view is used to run the program containing the
218: following code, then 2 sets of summary data will be printed during
219: `PetscFinalize()`.
220: .vb
221: PetscInitialize(int *argc,char ***args,0,0);
222: [stage 0 of code]
223: PetscStageLogPush(stageLog,1);
224: [stage 1 of code]
225: PetscStageLogPop(stageLog);
226: PetscBarrier(...);
227: [more stage 0 of code]
228: PetscFinalize();
229: .ve
231: Note;
232: Use `PetscLogStageRegister()` to register a stage. All previous stages are
233: accumulating time and flops, but events will only be logged in this stage.
235: Level: developer
237: .seealso: `PetscStageLogPop()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
238: @*/
239: PetscErrorCode PetscStageLogPush(PetscStageLog stageLog, int stage)
240: {
241: int curStage = 0;
242: PetscBool empty;
244: PetscFunctionBegin;
245: PetscCheck(!(stage < 0) && !(stage >= stageLog->numStages), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
247: /* Record flops/time of previous stage */
248: PetscCall(PetscIntStackEmpty(stageLog->stack, &empty));
249: if (!empty) {
250: PetscCall(PetscIntStackTop(stageLog->stack, &curStage));
251: if (stageLog->stageInfo[curStage].perfInfo.active) {
252: PetscCall(PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time));
253: stageLog->stageInfo[curStage].perfInfo.flops += petsc_TotalFlops;
254: stageLog->stageInfo[curStage].perfInfo.numMessages += petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
255: stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
256: stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
257: }
258: }
259: /* Activate the stage */
260: PetscCall(PetscIntStackPush(stageLog->stack, stage));
262: stageLog->stageInfo[stage].used = PETSC_TRUE;
263: stageLog->stageInfo[stage].perfInfo.count++;
264: stageLog->curStage = stage;
265: /* Subtract current quantities so that we obtain the difference when we pop */
266: if (stageLog->stageInfo[stage].perfInfo.active) {
267: PetscCall(PetscTimeSubtract(&stageLog->stageInfo[stage].perfInfo.time));
268: stageLog->stageInfo[stage].perfInfo.flops -= petsc_TotalFlops;
269: stageLog->stageInfo[stage].perfInfo.numMessages -= petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
270: stageLog->stageInfo[stage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
271: stageLog->stageInfo[stage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
272: }
273: PetscFunctionReturn(PETSC_SUCCESS);
274: }
276: /*@C
277: PetscStageLogPop - This function pops a stage from the stack.
279: Not Collective
281: Input Parameter:
282: . stageLog - The `PetscStageLog`
284: Usage:
285: If the option -log_view is used to run the program containing the
286: following code, then 2 sets of summary data will be printed during
287: PetscFinalize().
288: .vb
289: PetscInitialize(int *argc,char ***args,0,0);
290: [stage 0 of code]
291: PetscStageLogPush(stageLog,1);
292: [stage 1 of code]
293: PetscStageLogPop(stageLog);
294: PetscBarrier(...);
295: [more stage 0 of code]
296: PetscFinalize();
297: .ve
299: Note:
300: Use `PetscStageLogRegister()` to register a stage.
302: Level: developer
304: .seealso: `PetscStageLogPush()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
305: @*/
306: PetscErrorCode PetscStageLogPop(PetscStageLog stageLog)
307: {
308: int curStage;
309: PetscBool empty;
311: PetscFunctionBegin;
312: /* Record flops/time of current stage */
313: PetscCall(PetscIntStackPop(stageLog->stack, &curStage));
314: if (stageLog->stageInfo[curStage].perfInfo.active) {
315: PetscCall(PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time));
316: stageLog->stageInfo[curStage].perfInfo.flops += petsc_TotalFlops;
317: stageLog->stageInfo[curStage].perfInfo.numMessages += petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
318: stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
319: stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
320: }
321: PetscCall(PetscIntStackEmpty(stageLog->stack, &empty));
322: if (!empty) {
323: /* Subtract current quantities so that we obtain the difference when we pop */
324: PetscCall(PetscIntStackTop(stageLog->stack, &curStage));
325: if (stageLog->stageInfo[curStage].perfInfo.active) {
326: PetscCall(PetscTimeSubtract(&stageLog->stageInfo[curStage].perfInfo.time));
327: stageLog->stageInfo[curStage].perfInfo.flops -= petsc_TotalFlops;
328: stageLog->stageInfo[curStage].perfInfo.numMessages -= petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
329: stageLog->stageInfo[curStage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
330: stageLog->stageInfo[curStage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
331: }
332: stageLog->curStage = curStage;
333: } else stageLog->curStage = -1;
334: PetscFunctionReturn(PETSC_SUCCESS);
335: }
337: /*@C
338: PetscStageLogGetClassRegLog - This function returns the PetscClassRegLog for the given stage.
340: Not Collective
342: Input Parameter:
343: . stageLog - The `PetscStageLog`
345: Output Parameter:
346: . classLog - The `PetscClassRegLog`
348: Level: developer
350: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscLogGetStageLog()`
351: @*/
352: PetscErrorCode PetscStageLogGetClassRegLog(PetscStageLog stageLog, PetscClassRegLog *classLog)
353: {
354: PetscFunctionBegin;
356: *classLog = stageLog->classLog;
357: PetscFunctionReturn(PETSC_SUCCESS);
358: }
360: /*@C
361: PetscStageLogGetEventRegLog - This function returns the `PetscEventRegLog`.
363: Not Collective
365: Input Parameter:
366: . stageLog - The `PetscStageLog`
368: Output Parameter:
369: . eventLog - The `PetscEventRegLog`
371: Level: developer
373: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscLogGetStageLog()`
374: @*/
375: PetscErrorCode PetscStageLogGetEventRegLog(PetscStageLog stageLog, PetscEventRegLog *eventLog)
376: {
377: PetscFunctionBegin;
379: *eventLog = stageLog->eventLog;
380: PetscFunctionReturn(PETSC_SUCCESS);
381: }
383: /*@C
384: PetscStageLogGetClassPerfLog - This function returns the `PetscClassPerfLog` for the given stage.
386: Not Collective
388: Input Parameters:
389: + stageLog - The `PetscStageLog`
390: - stage - The stage
392: Output Parameter:
393: . classLog - The `PetscClassPerfLog`
395: Level: developer
397: .seealso: `PetscStageLogPush()`, `PetscStageLogPop()`, `PetscLogGetStageLog()`
398: @*/
399: PetscErrorCode PetscStageLogGetClassPerfLog(PetscStageLog stageLog, int stage, PetscClassPerfLog *classLog)
400: {
401: PetscFunctionBegin;
403: PetscCheck(!(stage < 0) && !(stage >= stageLog->numStages), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
404: *classLog = stageLog->stageInfo[stage].classLog;
405: PetscFunctionReturn(PETSC_SUCCESS);
406: }
408: /*@C
409: PetscStageLogSetActive - This function determines whether events will be logged during this state.
411: Not Collective
413: Input Parameters:
414: + stageLog - The `PetscStageLog`
415: . stage - The stage to log
416: - isActive - The activity flag, `PETSC_TRUE` for logging, otherwise `PETSC_FALSE` (default is `PETSC_TRUE`)
418: Level: developer
420: .seealso: `PetscStageLogGetActive()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
421: @*/
422: PetscErrorCode PetscStageLogSetActive(PetscStageLog stageLog, int stage, PetscBool isActive)
423: {
424: PetscFunctionBegin;
425: PetscCheck(!(stage < 0) && !(stage >= stageLog->numStages), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
426: stageLog->stageInfo[stage].perfInfo.active = isActive;
427: PetscFunctionReturn(PETSC_SUCCESS);
428: }
430: /*@C
431: PetscStageLogGetActive - This function returns whether events will be logged suring this stage.
433: Not Collective
435: Input Parameters:
436: + stageLog - The `PetscStageLog`
437: - stage - The stage to log
439: Output Parameter:
440: . isActive - The activity flag, `PETSC_TRUE` for logging, otherwise `PETSC_FALSE` (default is `PETSC_TRUE`)
442: Level: developer
444: .seealso: `PetscStageLogSetActive()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
445: @*/
446: PetscErrorCode PetscStageLogGetActive(PetscStageLog stageLog, int stage, PetscBool *isActive)
447: {
448: PetscFunctionBegin;
449: PetscCheck(!(stage < 0) && !(stage >= stageLog->numStages), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
451: *isActive = stageLog->stageInfo[stage].perfInfo.active;
452: PetscFunctionReturn(PETSC_SUCCESS);
453: }
455: /*@C
456: PetscStageLogSetVisible - This function determines whether a stage is printed during `PetscLogView()`
458: Not Collective
460: Input Parameters:
461: + stageLog - The `PetscStageLog`
462: . stage - The stage to log
463: - isVisible - The visibility flag, `PETSC_TRUE` for printing, otherwise `PETSC_FALSE` (default is `PETSC_TRUE`)
465: Options Database Key:
466: . -log_view - Activates log summary
468: Level: developer
470: .seealso: `PetscStageLogGetVisible()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
471: @*/
472: PetscErrorCode PetscStageLogSetVisible(PetscStageLog stageLog, int stage, PetscBool isVisible)
473: {
474: PetscFunctionBegin;
475: PetscCheck(!(stage < 0) && !(stage >= stageLog->numStages), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
476: stageLog->stageInfo[stage].perfInfo.visible = isVisible;
477: PetscFunctionReturn(PETSC_SUCCESS);
478: }
480: /*@C
481: PetscStageLogGetVisible - This function returns whether a stage is printed during `PetscLogView()`
483: Not Collective
485: Input Parameters:
486: + stageLog - The `PetscStageLog`
487: - stage - The stage to log
489: Output Parameter:
490: . isVisible - The visibility flag, `PETSC_TRUE` for printing, otherwise `PETSC_FALSE` (default is `PETSC_TRUE`)
492: Options Database Key:
493: . -log_view - Activates log summary
495: Level: developer
497: .seealso: `PetscStageLogSetVisible()`, `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
498: @*/
499: PetscErrorCode PetscStageLogGetVisible(PetscStageLog stageLog, int stage, PetscBool *isVisible)
500: {
501: PetscFunctionBegin;
502: PetscCheck(!(stage < 0) && !(stage >= stageLog->numStages), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid stage %d should be in [0,%d)", stage, stageLog->numStages);
504: *isVisible = stageLog->stageInfo[stage].perfInfo.visible;
505: PetscFunctionReturn(PETSC_SUCCESS);
506: }
508: /*@C
509: PetscStageLogGetStage - This function returns the stage id given the stage name.
511: Not Collective
513: Input Parameters:
514: + stageLog - The `PetscStageLog`
515: - name - The stage name
517: Output Parameter:
518: . stage - The stage id, or -1 if it does not exist
520: Level: developer
522: .seealso: `PetscStageLogGetCurrent()`, `PetscStageLogRegister()`, `PetscLogGetStageLog()`
523: @*/
524: PetscErrorCode PetscStageLogGetStage(PetscStageLog stageLog, const char name[], PetscLogStage *stage)
525: {
526: PetscBool match;
527: int s;
529: PetscFunctionBegin;
532: *stage = -1;
533: for (s = 0; s < stageLog->numStages; s++) {
534: PetscCall(PetscStrcasecmp(stageLog->stageInfo[s].name, name, &match));
535: if (match) {
536: *stage = s;
537: break;
538: }
539: }
540: PetscFunctionReturn(PETSC_SUCCESS);
541: }
543: /*@C
544: PetscStageLogCreate - This creates a `PetscStageLog` object.
546: Not collective
548: Output Parameter:
549: . stageLog - The `PetscStageLog`
551: Level: developer
553: .seealso: `PetscStageLogCreate()`
554: @*/
555: PetscErrorCode PetscStageLogCreate(PetscStageLog *stageLog)
556: {
557: PetscStageLog l;
559: PetscFunctionBegin;
560: PetscCall(PetscNew(&l));
562: l->numStages = 0;
563: l->maxStages = 10;
564: l->curStage = -1;
566: PetscCall(PetscIntStackCreate(&l->stack));
567: PetscCall(PetscMalloc1(l->maxStages, &l->stageInfo));
568: PetscCall(PetscEventRegLogCreate(&l->eventLog));
569: PetscCall(PetscClassRegLogCreate(&l->classLog));
571: *stageLog = l;
572: PetscFunctionReturn(PETSC_SUCCESS);
573: }