Actual source code: stagelog.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  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:
 23:     Inline since called for EACH PetscEventLogBeginDefault() and PetscEventLogEndDefault()

 25: .keywords: log, stage
 26: .seealso: PetscStageLogCreate()
 27: @*/
 28: PetscErrorCode PetscLogGetStageLog(PetscStageLog *stageLog)
 29: {
 32:   if (!petsc_stageLog) {
 33:     fprintf(stderr, "PETSC ERROR: Logging has not been enabled.\nYou might have forgotten to call PetscInitialize().\n");
 34:     MPI_Abort(MPI_COMM_WORLD, PETSC_ERR_SUP);
 35:   }
 36:   *stageLog = petsc_stageLog;
 37:   return(0);
 38: }

 40: /*@C
 41:   PetscStageLogGetCurrent - This function returns the stage from the top of the stack.

 43:   Not Collective

 45:   Input Parameter:
 46: . stageLog - The PetscStageLog

 48:   Output Parameter:
 49: . stage    - The current stage

 51:   Notes:
 52:   If no stage is currently active, stage is set to -1.

 54:   Level: developer

 56:   Developer Notes:
 57:     Inline since called for EACH PetscEventLogBeginDefault() and PetscEventLogEndDefault()

 59: .keywords: log, stage
 60: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
 61: @*/
 62: PetscErrorCode  PetscStageLogGetCurrent(PetscStageLog stageLog, int *stage)
 63: {
 64:   PetscBool      empty;

 68:   PetscIntStackEmpty(stageLog->stack, &empty);
 69:   if (empty) {
 70:     *stage = -1;
 71:   } else {
 72:     PetscIntStackTop(stageLog->stack, stage);
 73:   }
 74: #ifdef PETSC_USE_DEBUG
 75:   if (*stage != stageLog->curStage) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Inconsistency in stage log: stage %d should be %d", *stage, stageLog->curStage);
 76: #endif
 77:   return(0);
 78: }

 80: /*@C
 81:   PetscStageLogGetEventPerfLog - This function returns the PetscEventPerfLog for the given stage.

 83:   Not Collective

 85:   Input Parameters:
 86: + stageLog - The PetscStageLog
 87: - stage    - The stage

 89:   Output Parameter:
 90: . eventLog - The PetscEventPerfLog

 92:   Level: developer

 94:   Developer Notes:
 95:     Inline since called for EACH PetscEventLogBeginDefault() and PetscEventLogEndDefault()

 97: .keywords: log, stage
 98: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
 99: @*/
100: PetscErrorCode  PetscStageLogGetEventPerfLog(PetscStageLog stageLog, int stage, PetscEventPerfLog *eventLog)
101: {
104:   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);
105:   *eventLog = stageLog->stageInfo[stage].eventLog;
106:   return(0);
107: }

109: /*@C
110:   PetscStageInfoDestroy - This destroys a PetscStageInfo object.

112:   Not collective

114:   Input Paramter:
115: . stageInfo - The PetscStageInfo

117:   Level: developer

119: .keywords: log, stage, destroy
120: .seealso: PetscStageLogCreate()
121: @*/
122: PetscErrorCode  PetscStageInfoDestroy(PetscStageInfo *stageInfo)
123: {

127:   PetscFree(stageInfo->name);
128:   PetscEventPerfLogDestroy(stageInfo->eventLog);
129:   PetscClassPerfLogDestroy(stageInfo->classLog);
130:   return(0);
131: }

133: /*@C
134:   PetscStageLogDestroy - This destroys a PetscStageLog object.

136:   Not collective

138:   Input Paramter:
139: . stageLog - The PetscStageLog

141:   Level: developer

143: .keywords: log, stage, destroy
144: .seealso: PetscStageLogCreate()
145: @*/
146: PetscErrorCode  PetscStageLogDestroy(PetscStageLog stageLog)
147: {
148:   int            stage;

152:   if (!stageLog) return(0);
153:   PetscIntStackDestroy(stageLog->stack);
154:   PetscEventRegLogDestroy(stageLog->eventLog);
155:   PetscClassRegLogDestroy(stageLog->classLog);
156:   for (stage = 0; stage < stageLog->numStages; stage++) {
157:     PetscStageInfoDestroy(&stageLog->stageInfo[stage]);
158:   }
159:   PetscFree(stageLog->stageInfo);
160:   PetscFree(stageLog);
161:   return(0);
162: }

164: /*@C
165:   PetscStageLogRegister - Registers a stage name for logging operations in an application code.

167:   Not Collective

169:   Input Parameter:
170: + stageLog - The PetscStageLog
171: - sname    - the name to associate with that stage

173:   Output Parameter:
174: . stage    - The stage index

176:   Level: developer

178: .keywords: log, stage, register
179: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscStageLogCreate()
180: @*/
181: PetscErrorCode  PetscStageLogRegister(PetscStageLog stageLog, const char sname[], int *stage)
182: {
183:   PetscStageInfo *stageInfo;
184:   int            s;

190:   /* Check stage already registered */
191:   for (s = 0; s < stageLog->numStages; ++s) {
192:     PetscBool same;

194:     PetscStrcmp(stageLog->stageInfo[s].name, sname, &same);
195:     if (same) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Duplicate stage name given: %s", sname);
196:   }
197:   /* Create new stage */
198:   s = stageLog->numStages++;
199:   if (stageLog->numStages > stageLog->maxStages) {
200:     PetscMalloc1(stageLog->maxStages*2, &stageInfo);
201:     PetscMemcpy(stageInfo, stageLog->stageInfo, stageLog->maxStages * sizeof(PetscStageInfo));
202:     PetscFree(stageLog->stageInfo);
203:     stageLog->stageInfo  = stageInfo;
204:     stageLog->maxStages *= 2;
205:   }
206:   /* Setup new stage info */
207:   stageInfo = &stageLog->stageInfo[s];
208:   PetscMemzero(stageInfo,sizeof(PetscStageInfo));
209:   PetscStrallocpy(sname,&stageInfo->name);
210:   stageInfo->used             = PETSC_FALSE;
211:   stageInfo->perfInfo.active  = PETSC_TRUE;
212:   stageInfo->perfInfo.visible = PETSC_TRUE;
213:   PetscEventPerfLogCreate(&stageInfo->eventLog);
214:   PetscClassPerfLogCreate(&stageInfo->classLog);
215:   *stage = s;
216:   return(0);
217: }

219: /*@C
220:   PetscStageLogPush - This function pushes a stage on the stack.

222:   Not Collective

224:   Input Parameters:
225: + stageLog   - The PetscStageLog
226: - stage - The stage to log

228:   Database Options:
229: . -log_view - Activates logging

231:   Usage:
232:   If the option -log_sumary is used to run the program containing the
233:   following code, then 2 sets of summary data will be printed during
234:   PetscFinalize().
235: .vb
236:       PetscInitialize(int *argc,char ***args,0,0);
237:       [stage 0 of code]
238:       PetscStageLogPush(stageLog,1);
239:       [stage 1 of code]
240:       PetscStageLogPop(stageLog);
241:       PetscBarrier(...);
242:       [more stage 0 of code]
243:       PetscFinalize();
244: .ve

246:   Notes:
247:   Use PetscLogStageRegister() to register a stage. All previous stages are
248:   accumulating time and flops, but events will only be logged in this stage.

250:   Level: developer

252: .keywords: log, push, stage
253: .seealso: PetscStageLogPop(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
254: @*/
255: PetscErrorCode  PetscStageLogPush(PetscStageLog stageLog, int stage)
256: {
257:   int            curStage = 0;
258:   PetscBool      empty;

262:   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);

264:   /* Record flops/time of previous stage */
265:   PetscIntStackEmpty(stageLog->stack, &empty);
266:   if (!empty) {
267:     PetscIntStackTop(stageLog->stack, &curStage);
268:     if (stageLog->stageInfo[curStage].perfInfo.active) {
269:       PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time);
270:       stageLog->stageInfo[curStage].perfInfo.flops         += petsc_TotalFlops;
271:       stageLog->stageInfo[curStage].perfInfo.numMessages   += petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
272:       stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
273:       stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
274:     }
275:   }
276:   /* Activate the stage */
277:   PetscIntStackPush(stageLog->stack, stage);

279:   stageLog->stageInfo[stage].used = PETSC_TRUE;
280:   stageLog->stageInfo[stage].perfInfo.count++;
281:   stageLog->curStage = stage;
282:   /* Subtract current quantities so that we obtain the difference when we pop */
283:   if (stageLog->stageInfo[stage].perfInfo.active) {
284:     PetscTimeSubtract(&stageLog->stageInfo[stage].perfInfo.time);
285:     stageLog->stageInfo[stage].perfInfo.flops         -= petsc_TotalFlops;
286:     stageLog->stageInfo[stage].perfInfo.numMessages   -= petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
287:     stageLog->stageInfo[stage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
288:     stageLog->stageInfo[stage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
289:   }
290:   return(0);
291: }

293: /*@C
294:   PetscStageLogPop - This function pops a stage from the stack.

296:   Not Collective

298:   Input Parameter:
299: . stageLog - The PetscStageLog

301:   Usage:
302:   If the option -log_sumary is used to run the program containing the
303:   following code, then 2 sets of summary data will be printed during
304:   PetscFinalize().
305: .vb
306:       PetscInitialize(int *argc,char ***args,0,0);
307:       [stage 0 of code]
308:       PetscStageLogPush(stageLog,1);
309:       [stage 1 of code]
310:       PetscStageLogPop(stageLog);
311:       PetscBarrier(...);
312:       [more stage 0 of code]
313:       PetscFinalize();
314: .ve

316:   Notes:
317:   Use PetscStageLogRegister() to register a stage.

319:   Level: developer

321: .keywords: log, pop, stage
322: .seealso: PetscStageLogPush(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
323: @*/
324: PetscErrorCode  PetscStageLogPop(PetscStageLog stageLog)
325: {
326:   int            curStage;
327:   PetscBool      empty;

331:   /* Record flops/time of current stage */
332:   PetscIntStackPop(stageLog->stack, &curStage);
333:   if (stageLog->stageInfo[curStage].perfInfo.active) {
334:     PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time);
335:     stageLog->stageInfo[curStage].perfInfo.flops         += petsc_TotalFlops;
336:     stageLog->stageInfo[curStage].perfInfo.numMessages   += petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
337:     stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
338:     stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
339:   }
340:   PetscIntStackEmpty(stageLog->stack, &empty);
341:   if (!empty) {
342:     /* Subtract current quantities so that we obtain the difference when we pop */
343:     PetscIntStackTop(stageLog->stack, &curStage);
344:     if (stageLog->stageInfo[curStage].perfInfo.active) {
345:       PetscTimeSubtract(&stageLog->stageInfo[curStage].perfInfo.time);
346:       stageLog->stageInfo[curStage].perfInfo.flops         -= petsc_TotalFlops;
347:       stageLog->stageInfo[curStage].perfInfo.numMessages   -= petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
348:       stageLog->stageInfo[curStage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
349:       stageLog->stageInfo[curStage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
350:     }
351:     stageLog->curStage = curStage;
352:   } else stageLog->curStage = -1;
353:   return(0);
354: }


357: /*@C
358:   PetscStageLogGetClassRegLog - This function returns the PetscClassRegLog for the given stage.

360:   Not Collective

362:   Input Parameters:
363: . stageLog - The PetscStageLog

365:   Output Parameter:
366: . classLog - The PetscClassRegLog

368:   Level: developer

370: .keywords: log, stage
371: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
372: @*/
373: PetscErrorCode  PetscStageLogGetClassRegLog(PetscStageLog stageLog, PetscClassRegLog *classLog)
374: {
377:   *classLog = stageLog->classLog;
378:   return(0);
379: }

381: /*@C
382:   PetscStageLogGetEventRegLog - This function returns the PetscEventRegLog.

384:   Not Collective

386:   Input Parameters:
387: . stageLog - The PetscStageLog

389:   Output Parameter:
390: . eventLog - The PetscEventRegLog

392:   Level: developer

394: .keywords: log, stage
395: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
396: @*/
397: PetscErrorCode  PetscStageLogGetEventRegLog(PetscStageLog stageLog, PetscEventRegLog *eventLog)
398: {
401:   *eventLog = stageLog->eventLog;
402:   return(0);
403: }

405: /*@C
406:   PetscStageLogGetClassPerfLog - This function returns the PetscClassPerfLog for the given stage.

408:   Not Collective

410:   Input Parameters:
411: + stageLog - The PetscStageLog
412: - stage    - The stage

414:   Output Parameter:
415: . classLog - The PetscClassPerfLog

417:   Level: developer

419: .keywords: log, stage
420: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
421: @*/
422: PetscErrorCode  PetscStageLogGetClassPerfLog(PetscStageLog stageLog, int stage, PetscClassPerfLog *classLog)
423: {
426:   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);
427:   *classLog = stageLog->stageInfo[stage].classLog;
428:   return(0);
429: }


432: /*@C
433:   PetscStageLogSetActive - This function determines whether events will be logged during this state.

435:   Not Collective

437:   Input Parameters:
438: + stageLog - The PetscStageLog
439: . stage    - The stage to log
440: - isActive - The activity flag, PETSC_TRUE for logging, otherwise PETSC_FALSE (default is PETSC_TRUE)

442:   Level: developer

444: .keywords: log, active, stage
445: .seealso: PetscStageLogGetActive(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
446: @*/
447: PetscErrorCode  PetscStageLogSetActive(PetscStageLog stageLog, int stage, PetscBool isActive)
448: {
450:   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);
451:   stageLog->stageInfo[stage].perfInfo.active = isActive;
452:   return(0);
453: }

455: /*@C
456:   PetscStageLogGetActive - This function returns whether events will be logged suring this stage.

458:   Not Collective

460:   Input Parameters:
461: + stageLog - The PetscStageLog
462: - stage    - The stage to log

464:   Output Parameter:
465: . isActive - The activity flag, PETSC_TRUE for logging, otherwise PETSC_FALSE (default is PETSC_TRUE)

467:   Level: developer

469: .keywords: log, visible, stage
470: .seealso: PetscStageLogSetActive(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
471: @*/
472: PetscErrorCode  PetscStageLogGetActive(PetscStageLog stageLog, int stage, PetscBool  *isActive)
473: {
475:   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);
477:   *isActive = stageLog->stageInfo[stage].perfInfo.active;
478:   return(0);
479: }

481: /*@C
482:   PetscStageLogSetVisible - This function determines whether a stage is printed during PetscLogView()

484:   Not Collective

486:   Input Parameters:
487: + stageLog  - The PetscStageLog
488: . stage     - The stage to log
489: - isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)

491:   Database Options:
492: . -log_view - Activates log summary

494:   Level: developer

496: .keywords: log, visible, stage
497: .seealso: PetscStageLogGetVisible(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
498: @*/
499: PetscErrorCode  PetscStageLogSetVisible(PetscStageLog stageLog, int stage, PetscBool isVisible)
500: {
502:   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);
503:   stageLog->stageInfo[stage].perfInfo.visible = isVisible;
504:   return(0);
505: }

507: /*@C
508:   PetscStageLogGetVisible - This function returns whether a stage is printed during PetscLogView()

510:   Not Collective

512:   Input Parameters:
513: + stageLog  - The PetscStageLog
514: - stage     - The stage to log

516:   Output Parameter:
517: . isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)

519:   Database Options:
520: . -log_view - Activates log summary

522:   Level: developer

524: .keywords: log, visible, stage
525: .seealso: PetscStageLogSetVisible(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
526: @*/
527: PetscErrorCode  PetscStageLogGetVisible(PetscStageLog stageLog, int stage, PetscBool  *isVisible)
528: {
530:   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);
532:   *isVisible = stageLog->stageInfo[stage].perfInfo.visible;
533:   return(0);
534: }

536: /*@C
537:   PetscStageLogGetStage - This function returns the stage id given the stage name.

539:   Not Collective

541:   Input Parameters:
542: + stageLog - The PetscStageLog
543: - name     - The stage name

545:   Output Parameter:
546: . stage    - The stage id, or -1 if it does not exist

548:   Level: developer

550: .keywords: log, stage
551: .seealso: PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
552: @*/
553: PetscErrorCode  PetscStageLogGetStage(PetscStageLog stageLog, const char name[], PetscLogStage *stage)
554: {
555:   PetscBool      match;
556:   int            s;

562:   *stage = -1;
563:   for (s = 0; s < stageLog->numStages; s++) {
564:     PetscStrcasecmp(stageLog->stageInfo[s].name, name, &match);
565:     if (match) {
566:       *stage = s;
567:       break;
568:     }
569:   }
570:   return(0);
571: }

573: /*@C
574:   PetscStageLogCreate - This creates a PetscStageLog object.

576:   Not collective

578:   Input Parameter:
579: . stageLog - The PetscStageLog

581:   Level: developer

583: .keywords: log, stage, create
584: .seealso: PetscStageLogCreate()
585: @*/
586: PetscErrorCode  PetscStageLogCreate(PetscStageLog *stageLog)
587: {
588:   PetscStageLog  l;

592:   PetscNew(&l);

594:   l->numStages = 0;
595:   l->maxStages = 10;
596:   l->curStage  = -1;

598:   PetscIntStackCreate(&l->stack);
599:   PetscMalloc1(l->maxStages, &l->stageInfo);
600:   PetscEventRegLogCreate(&l->eventLog);
601:   PetscClassRegLogCreate(&l->classLog);

603:   *stageLog = l;
604:   return(0);
605: }