Actual source code: stagelog.c
petsc-3.5.4 2015-05-23
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> /*I "petscsys.h" I*/
10: PetscStageLog petsc_stageLog = 0;
14: /*@C
15: PetscLogGetStageLog - This function returns the default stage logging object.
17: Not collective
19: Output Parameter:
20: . stageLog - The default PetscStageLog
22: Level: developer
24: Developer Notes: Inline since called for EACH PetscEventLogBeginDefault() and PetscEventLogEndDefault()
26: .keywords: log, stage
27: .seealso: PetscStageLogCreate()
28: @*/
29: PetscErrorCode PetscLogGetStageLog(PetscStageLog *stageLog)
30: {
33: if (!petsc_stageLog) {
34: fprintf(stderr, "PETSC ERROR: Logging has not been enabled.\nYou might have forgotten to call PetscInitialize().\n");
35: MPI_Abort(MPI_COMM_WORLD, PETSC_ERR_SUP);
36: }
37: *stageLog = petsc_stageLog;
38: return(0);
39: }
43: /*@C
44: PetscStageLogGetCurrent - This function returns the stage from the top of the stack.
46: Not Collective
48: Input Parameter:
49: . stageLog - The PetscStageLog
51: Output Parameter:
52: . stage - The current stage
54: Notes:
55: If no stage is currently active, stage is set to -1.
57: Level: developer
59: Developer Notes: Inline since called for EACH PetscEventLogBeginDefault() and PetscEventLogEndDefault()
61: .keywords: log, stage
62: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
63: @*/
64: PetscErrorCode PetscStageLogGetCurrent(PetscStageLog stageLog, int *stage)
65: {
66: PetscBool empty;
70: PetscIntStackEmpty(stageLog->stack, &empty);
71: if (empty) {
72: *stage = -1;
73: } else {
74: PetscIntStackTop(stageLog->stack, stage);
75: }
76: #ifdef PETSC_USE_DEBUG
77: if (*stage != stageLog->curStage) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Inconsistency in stage log: stage %d should be %d", *stage, stageLog->curStage);
78: #endif
79: return(0);
80: }
84: /*@C
85: PetscStageLogGetEventPerfLog - This function returns the PetscEventPerfLog for the given stage.
87: Not Collective
89: Input Parameters:
90: + stageLog - The PetscStageLog
91: - stage - The stage
93: Output Parameter:
94: . eventLog - The PetscEventPerfLog
96: Level: developer
98: Developer Notes: Inline since called for EACH PetscEventLogBeginDefault() and PetscEventLogEndDefault()
100: .keywords: log, stage
101: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
102: @*/
103: PetscErrorCode PetscStageLogGetEventPerfLog(PetscStageLog stageLog, int stage, PetscEventPerfLog *eventLog)
104: {
107: 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);
108: *eventLog = stageLog->stageInfo[stage].eventLog;
109: return(0);
110: }
114: /*@C
115: PetscStageInfoDestroy - This destroys a PetscStageInfo object.
117: Not collective
119: Input Paramter:
120: . stageInfo - The PetscStageInfo
122: Level: developer
124: .keywords: log, stage, destroy
125: .seealso: PetscStageLogCreate()
126: @*/
127: PetscErrorCode PetscStageInfoDestroy(PetscStageInfo *stageInfo)
128: {
132: PetscFree(stageInfo->name);
133: EventPerfLogDestroy(stageInfo->eventLog);
134: ClassPerfLogDestroy(stageInfo->classLog);
135: return(0);
136: }
140: /*@C
141: PetscStageLogDestroy - This destroys a PetscStageLog object.
143: Not collective
145: Input Paramter:
146: . stageLog - The PetscStageLog
148: Level: developer
150: .keywords: log, stage, destroy
151: .seealso: PetscStageLogCreate()
152: @*/
153: PetscErrorCode PetscStageLogDestroy(PetscStageLog stageLog)
154: {
155: int stage;
159: if (!stageLog) return(0);
160: PetscIntStackDestroy(stageLog->stack);
161: EventRegLogDestroy(stageLog->eventLog);
162: PetscClassRegLogDestroy(stageLog->classLog);
163: for (stage = 0; stage < stageLog->numStages; stage++) {
164: PetscStageInfoDestroy(&stageLog->stageInfo[stage]);
165: }
166: PetscFree(stageLog->stageInfo);
167: PetscFree(stageLog);
168: return(0);
169: }
173: /*@C
174: PetscStageLogRegister - Registers a stage name for logging operations in an application code.
176: Not Collective
178: Input Parameter:
179: + stageLog - The PetscStageLog
180: - sname - the name to associate with that stage
182: Output Parameter:
183: . stage - The stage index
185: Level: developer
187: .keywords: log, stage, register
188: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscStageLogCreate()
189: @*/
190: PetscErrorCode PetscStageLogRegister(PetscStageLog stageLog, const char sname[], int *stage)
191: {
192: PetscStageInfo *stageInfo;
193: char *str;
194: int s, st;
200: for (st = 0; st < stageLog->numStages; ++st) {
201: PetscBool same;
203: PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);
204: if (same) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Duplicate stage name given: %s", sname);
205: }
206: s = stageLog->numStages++;
207: if (stageLog->numStages > stageLog->maxStages) {
208: PetscMalloc1(stageLog->maxStages*2, &stageInfo);
209: PetscMemcpy(stageInfo, stageLog->stageInfo, stageLog->maxStages * sizeof(PetscStageInfo));
210: PetscFree(stageLog->stageInfo);
212: stageLog->stageInfo = stageInfo;
213: stageLog->maxStages *= 2;
214: }
215: /* Setup stage */
216: PetscStrallocpy(sname, &str);
218: stageLog->stageInfo[s].name = str;
219: stageLog->stageInfo[s].used = PETSC_FALSE;
220: stageLog->stageInfo[s].perfInfo.active = PETSC_TRUE;
221: stageLog->stageInfo[s].perfInfo.visible = PETSC_TRUE;
222: stageLog->stageInfo[s].perfInfo.count = 0;
223: stageLog->stageInfo[s].perfInfo.flops = 0.0;
224: stageLog->stageInfo[s].perfInfo.time = 0.0;
225: stageLog->stageInfo[s].perfInfo.numMessages = 0.0;
226: stageLog->stageInfo[s].perfInfo.messageLength = 0.0;
227: stageLog->stageInfo[s].perfInfo.numReductions = 0.0;
229: EventPerfLogCreate(&stageLog->stageInfo[s].eventLog);
230: ClassPerfLogCreate(&stageLog->stageInfo[s].classLog);
231: *stage = s;
232: return(0);
233: }
237: /*@C
238: PetscStageLogPush - This function pushes a stage on the stack.
240: Not Collective
242: Input Parameters:
243: + stageLog - The PetscStageLog
244: - stage - The stage to log
246: Database Options:
247: . -log_summary - Activates logging
249: Usage:
250: If the option -log_sumary is used to run the program containing the
251: following code, then 2 sets of summary data will be printed during
252: PetscFinalize().
253: .vb
254: PetscInitialize(int *argc,char ***args,0,0);
255: [stage 0 of code]
256: PetscStageLogPush(stageLog,1);
257: [stage 1 of code]
258: PetscStageLogPop(stageLog);
259: PetscBarrier(...);
260: [more stage 0 of code]
261: PetscFinalize();
262: .ve
264: Notes:
265: Use PetscLogStageRegister() to register a stage. All previous stages are
266: accumulating time and flops, but events will only be logged in this stage.
268: Level: developer
270: .keywords: log, push, stage
271: .seealso: PetscStageLogPop(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
272: @*/
273: PetscErrorCode PetscStageLogPush(PetscStageLog stageLog, int stage)
274: {
275: int curStage = 0;
276: PetscBool empty;
280: 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);
282: /* Record flops/time of previous stage */
283: PetscIntStackEmpty(stageLog->stack, &empty);
284: if (!empty) {
285: PetscIntStackTop(stageLog->stack, &curStage);
286: if (stageLog->stageInfo[curStage].perfInfo.active) {
287: PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time);
288: stageLog->stageInfo[curStage].perfInfo.flops += petsc_TotalFlops;
289: stageLog->stageInfo[curStage].perfInfo.numMessages += petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
290: stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
291: stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
292: }
293: }
294: /* Activate the stage */
295: PetscIntStackPush(stageLog->stack, stage);
297: stageLog->stageInfo[stage].used = PETSC_TRUE;
298: stageLog->stageInfo[stage].perfInfo.count++;
299: stageLog->curStage = stage;
300: /* Subtract current quantities so that we obtain the difference when we pop */
301: if (stageLog->stageInfo[stage].perfInfo.active) {
302: PetscTimeSubtract(&stageLog->stageInfo[stage].perfInfo.time);
303: stageLog->stageInfo[stage].perfInfo.flops -= petsc_TotalFlops;
304: stageLog->stageInfo[stage].perfInfo.numMessages -= petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
305: stageLog->stageInfo[stage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
306: stageLog->stageInfo[stage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
307: }
308: return(0);
309: }
313: /*@C
314: PetscStageLogPop - This function pops a stage from the stack.
316: Not Collective
318: Input Parameter:
319: . stageLog - The PetscStageLog
321: Usage:
322: If the option -log_sumary is used to run the program containing the
323: following code, then 2 sets of summary data will be printed during
324: PetscFinalize().
325: .vb
326: PetscInitialize(int *argc,char ***args,0,0);
327: [stage 0 of code]
328: PetscStageLogPush(stageLog,1);
329: [stage 1 of code]
330: PetscStageLogPop(stageLog);
331: PetscBarrier(...);
332: [more stage 0 of code]
333: PetscFinalize();
334: .ve
336: Notes:
337: Use PetscStageLogRegister() to register a stage.
339: Level: developer
341: .keywords: log, pop, stage
342: .seealso: PetscStageLogPush(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
343: @*/
344: PetscErrorCode PetscStageLogPop(PetscStageLog stageLog)
345: {
346: int curStage;
347: PetscBool empty;
351: /* Record flops/time of current stage */
352: PetscIntStackPop(stageLog->stack, &curStage);
353: if (stageLog->stageInfo[curStage].perfInfo.active) {
354: PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time);
355: stageLog->stageInfo[curStage].perfInfo.flops += petsc_TotalFlops;
356: stageLog->stageInfo[curStage].perfInfo.numMessages += petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
357: stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
358: stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
359: }
360: PetscIntStackEmpty(stageLog->stack, &empty);
361: if (!empty) {
362: /* Subtract current quantities so that we obtain the difference when we pop */
363: PetscIntStackTop(stageLog->stack, &curStage);
364: if (stageLog->stageInfo[curStage].perfInfo.active) {
365: PetscTimeSubtract(&stageLog->stageInfo[curStage].perfInfo.time);
366: stageLog->stageInfo[curStage].perfInfo.flops -= petsc_TotalFlops;
367: stageLog->stageInfo[curStage].perfInfo.numMessages -= petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct;
368: stageLog->stageInfo[curStage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
369: stageLog->stageInfo[curStage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
370: }
371: stageLog->curStage = curStage;
372: } else stageLog->curStage = -1;
373: return(0);
374: }
379: /*@C
380: PetscStageLogGetClassRegLog - This function returns the PetscClassRegLog for the given stage.
382: Not Collective
384: Input Parameters:
385: . stageLog - The PetscStageLog
387: Output Parameter:
388: . classLog - The PetscClassRegLog
390: Level: developer
392: .keywords: log, stage
393: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
394: @*/
395: PetscErrorCode PetscStageLogGetClassRegLog(PetscStageLog stageLog, PetscClassRegLog *classLog)
396: {
399: *classLog = stageLog->classLog;
400: return(0);
401: }
405: /*@C
406: PetscStageLogGetEventRegLog - This function returns the PetscEventRegLog.
408: Not Collective
410: Input Parameters:
411: . stageLog - The PetscStageLog
413: Output Parameter:
414: . eventLog - The PetscEventRegLog
416: Level: developer
418: .keywords: log, stage
419: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
420: @*/
421: PetscErrorCode PetscStageLogGetEventRegLog(PetscStageLog stageLog, PetscEventRegLog *eventLog)
422: {
425: *eventLog = stageLog->eventLog;
426: return(0);
427: }
431: /*@C
432: PetscStageLogGetClassPerfLog - This function returns the ClassPerfLog for the given stage.
434: Not Collective
436: Input Parameters:
437: + stageLog - The PetscStageLog
438: - stage - The stage
440: Output Parameter:
441: . classLog - The ClassPerfLog
443: Level: developer
445: .keywords: log, stage
446: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
447: @*/
448: PetscErrorCode PetscStageLogGetClassPerfLog(PetscStageLog stageLog, int stage, PetscClassPerfLog *classLog)
449: {
452: 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);
453: *classLog = stageLog->stageInfo[stage].classLog;
454: return(0);
455: }
460: /*@C
461: PetscStageLogSetActive - This function determines whether events will be logged during this state.
463: Not Collective
465: Input Parameters:
466: + stageLog - The PetscStageLog
467: . stage - The stage to log
468: - isActive - The activity flag, PETSC_TRUE for logging, otherwise PETSC_FALSE (default is PETSC_TRUE)
470: Level: developer
472: .keywords: log, active, stage
473: .seealso: PetscStageLogGetActive(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
474: @*/
475: PetscErrorCode PetscStageLogSetActive(PetscStageLog stageLog, int stage, PetscBool isActive)
476: {
478: 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);
479: stageLog->stageInfo[stage].perfInfo.active = isActive;
480: return(0);
481: }
485: /*@C
486: PetscStageLogGetActive - This function returns whether events will be logged suring this stage.
488: Not Collective
490: Input Parameters:
491: + stageLog - The PetscStageLog
492: - stage - The stage to log
494: Output Parameter:
495: . isActive - The activity flag, PETSC_TRUE for logging, otherwise PETSC_FALSE (default is PETSC_TRUE)
497: Level: developer
499: .keywords: log, visible, stage
500: .seealso: PetscStageLogSetActive(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
501: @*/
502: PetscErrorCode PetscStageLogGetActive(PetscStageLog stageLog, int stage, PetscBool *isActive)
503: {
505: 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: *isActive = stageLog->stageInfo[stage].perfInfo.active;
508: return(0);
509: }
513: /*@C
514: PetscStageLogSetVisible - This function determines whether a stage is printed during PetscLogView()
516: Not Collective
518: Input Parameters:
519: + stageLog - The PetscStageLog
520: . stage - The stage to log
521: - isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)
523: Database Options:
524: . -log_summary - Activates log summary
526: Level: developer
528: .keywords: log, visible, stage
529: .seealso: PetscStageLogGetVisible(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
530: @*/
531: PetscErrorCode PetscStageLogSetVisible(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);
535: stageLog->stageInfo[stage].perfInfo.visible = isVisible;
536: return(0);
537: }
541: /*@C
542: PetscStageLogGetVisible - This function returns whether a stage is printed during PetscLogView()
544: Not Collective
546: Input Parameters:
547: + stageLog - The PetscStageLog
548: - stage - The stage to log
550: Output Parameter:
551: . isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)
553: Database Options:
554: . -log_summary - Activates log summary
556: Level: developer
558: .keywords: log, visible, stage
559: .seealso: PetscStageLogSetVisible(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
560: @*/
561: PetscErrorCode PetscStageLogGetVisible(PetscStageLog stageLog, int stage, PetscBool *isVisible)
562: {
564: 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);
566: *isVisible = stageLog->stageInfo[stage].perfInfo.visible;
567: return(0);
568: }
572: /*@C
573: PetscStageLogGetStage - This function returns the stage id given the stage name.
575: Not Collective
577: Input Parameters:
578: + stageLog - The PetscStageLog
579: - name - The stage name
581: Output Parameter:
582: . stage - The stage id
584: Level: developer
586: .keywords: log, stage
587: .seealso: PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
588: @*/
589: PetscErrorCode PetscStageLogGetStage(PetscStageLog stageLog, const char name[], int *stage)
590: {
591: PetscBool match;
592: int s;
598: *stage = -1;
599: for (s = 0; s < stageLog->numStages; s++) {
600: PetscStrcasecmp(stageLog->stageInfo[s].name, name, &match);
601: if (match) break;
602: }
603: if (s == stageLog->numStages) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "No stage named %s", name);
604: *stage = s;
605: return(0);
606: }
610: /*@C
611: PetscStageLogCreate - This creates a PetscStageLog object.
613: Not collective
615: Input Parameter:
616: . stageLog - The PetscStageLog
618: Level: developer
620: .keywords: log, stage, create
621: .seealso: PetscStageLogCreate()
622: @*/
623: PetscErrorCode PetscStageLogCreate(PetscStageLog *stageLog)
624: {
625: PetscStageLog l;
629: PetscNew(&l);
631: l->numStages = 0;
632: l->maxStages = 10;
633: l->curStage = -1;
635: PetscIntStackCreate(&l->stack);
636: PetscMalloc1(l->maxStages, &l->stageInfo);
637: EventRegLogCreate(&l->eventLog);
638: PetscClassRegLogCreate(&l->classLog);
640: *stageLog = l;
641: return(0);
642: }