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 Notes:
 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:   Notes:
 50:   If no stage is currently active, stage is set to -1.

 52:   Level: developer

 54:   Developer Notes:
 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 Notes:
 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++) {
141:     PetscStageInfoDestroy(&stageLog->stageInfo[stage]);
142:   }
143:   PetscFree(stageLog->stageInfo);
144:   PetscFree(stageLog);
145:   return 0;
146: }

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

151:   Not Collective

153:   Input Parameters:
154: + stageLog - The PetscStageLog
155: - sname    - the name to associate with that stage

157:   Output Parameter:
158: . stage    - The stage index

160:   Level: developer

162: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscStageLogCreate()
163: @*/
164: PetscErrorCode  PetscStageLogRegister(PetscStageLog stageLog, const char sname[], int *stage)
165: {
166:   PetscStageInfo *stageInfo;
167:   int            s;

171:   /* Check stage already registered */
172:   for (s = 0; s < stageLog->numStages; ++s) {
173:     PetscBool same;

175:     PetscStrcmp(stageLog->stageInfo[s].name, sname, &same);
177:   }
178:   /* Create new stage */
179:   s = stageLog->numStages++;
180:   if (stageLog->numStages > stageLog->maxStages) {
181:     PetscMalloc1(stageLog->maxStages*2, &stageInfo);
182:     PetscArraycpy(stageInfo, stageLog->stageInfo, stageLog->maxStages);
183:     PetscFree(stageLog->stageInfo);
184:     stageLog->stageInfo  = stageInfo;
185:     stageLog->maxStages *= 2;
186:   }
187:   /* Setup new stage info */
188:   stageInfo = &stageLog->stageInfo[s];
189:   PetscMemzero(stageInfo,sizeof(PetscStageInfo));
190:   PetscStrallocpy(sname,&stageInfo->name);
191:   stageInfo->used             = PETSC_FALSE;
192:   stageInfo->perfInfo.active  = PETSC_TRUE;
193:   stageInfo->perfInfo.visible = PETSC_TRUE;
194:   PetscEventPerfLogCreate(&stageInfo->eventLog);
195:   PetscClassPerfLogCreate(&stageInfo->classLog);
196:   *stage = s;
197:   return 0;
198: }

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

203:   Not Collective

205:   Input Parameters:
206: + stageLog   - The PetscStageLog
207: - stage - The stage to log

209:   Database Options:
210: . -log_view - Activates logging

212:   Usage:
213:   If the option -log_sumary is used to run the program containing the
214:   following code, then 2 sets of summary data will be printed during
215:   PetscFinalize().
216: .vb
217:       PetscInitialize(int *argc,char ***args,0,0);
218:       [stage 0 of code]
219:       PetscStageLogPush(stageLog,1);
220:       [stage 1 of code]
221:       PetscStageLogPop(stageLog);
222:       PetscBarrier(...);
223:       [more stage 0 of code]
224:       PetscFinalize();
225: .ve

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

231:   Level: developer

233: .seealso: PetscStageLogPop(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
234: @*/
235: PetscErrorCode  PetscStageLogPush(PetscStageLog stageLog, int stage)
236: {
237:   int            curStage = 0;
238:   PetscBool      empty;


242:   /* Record flops/time of previous stage */
243:   PetscIntStackEmpty(stageLog->stack, &empty);
244:   if (!empty) {
245:     PetscIntStackTop(stageLog->stack, &curStage);
246:     if (stageLog->stageInfo[curStage].perfInfo.active) {
247:       PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time);
248:       stageLog->stageInfo[curStage].perfInfo.flops         += petsc_TotalFlops;
249:       stageLog->stageInfo[curStage].perfInfo.numMessages   += petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
250:       stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
251:       stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
252:     }
253:   }
254:   /* Activate the stage */
255:   PetscIntStackPush(stageLog->stack, stage);

257:   stageLog->stageInfo[stage].used = PETSC_TRUE;
258:   stageLog->stageInfo[stage].perfInfo.count++;
259:   stageLog->curStage = stage;
260:   /* Subtract current quantities so that we obtain the difference when we pop */
261:   if (stageLog->stageInfo[stage].perfInfo.active) {
262:     PetscTimeSubtract(&stageLog->stageInfo[stage].perfInfo.time);
263:     stageLog->stageInfo[stage].perfInfo.flops         -= petsc_TotalFlops;
264:     stageLog->stageInfo[stage].perfInfo.numMessages   -= petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
265:     stageLog->stageInfo[stage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
266:     stageLog->stageInfo[stage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
267:   }
268:   return 0;
269: }

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

274:   Not Collective

276:   Input Parameter:
277: . stageLog - The PetscStageLog

279:   Usage:
280:   If the option -log_sumary is used to run the program containing the
281:   following code, then 2 sets of summary data will be printed during
282:   PetscFinalize().
283: .vb
284:       PetscInitialize(int *argc,char ***args,0,0);
285:       [stage 0 of code]
286:       PetscStageLogPush(stageLog,1);
287:       [stage 1 of code]
288:       PetscStageLogPop(stageLog);
289:       PetscBarrier(...);
290:       [more stage 0 of code]
291:       PetscFinalize();
292: .ve

294:   Notes:
295:   Use PetscStageLogRegister() to register a stage.

297:   Level: developer

299: .seealso: PetscStageLogPush(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
300: @*/
301: PetscErrorCode  PetscStageLogPop(PetscStageLog stageLog)
302: {
303:   int            curStage;
304:   PetscBool      empty;

306:   /* Record flops/time of current stage */
307:   PetscIntStackPop(stageLog->stack, &curStage);
308:   if (stageLog->stageInfo[curStage].perfInfo.active) {
309:     PetscTimeAdd(&stageLog->stageInfo[curStage].perfInfo.time);
310:     stageLog->stageInfo[curStage].perfInfo.flops         += petsc_TotalFlops;
311:     stageLog->stageInfo[curStage].perfInfo.numMessages   += petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
312:     stageLog->stageInfo[curStage].perfInfo.messageLength += petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
313:     stageLog->stageInfo[curStage].perfInfo.numReductions += petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
314:   }
315:   PetscIntStackEmpty(stageLog->stack, &empty);
316:   if (!empty) {
317:     /* Subtract current quantities so that we obtain the difference when we pop */
318:     PetscIntStackTop(stageLog->stack, &curStage);
319:     if (stageLog->stageInfo[curStage].perfInfo.active) {
320:       PetscTimeSubtract(&stageLog->stageInfo[curStage].perfInfo.time);
321:       stageLog->stageInfo[curStage].perfInfo.flops         -= petsc_TotalFlops;
322:       stageLog->stageInfo[curStage].perfInfo.numMessages   -= petsc_irecv_ct  + petsc_isend_ct  + petsc_recv_ct  + petsc_send_ct;
323:       stageLog->stageInfo[curStage].perfInfo.messageLength -= petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len;
324:       stageLog->stageInfo[curStage].perfInfo.numReductions -= petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
325:     }
326:     stageLog->curStage = curStage;
327:   } else stageLog->curStage = -1;
328:   return 0;
329: }

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

334:   Not Collective

336:   Input Parameters:
337: . stageLog - The PetscStageLog

339:   Output Parameter:
340: . classLog - The PetscClassRegLog

342:   Level: developer

344: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
345: @*/
346: PetscErrorCode  PetscStageLogGetClassRegLog(PetscStageLog stageLog, PetscClassRegLog *classLog)
347: {
349:   *classLog = stageLog->classLog;
350:   return 0;
351: }

353: /*@C
354:   PetscStageLogGetEventRegLog - This function returns the PetscEventRegLog.

356:   Not Collective

358:   Input Parameters:
359: . stageLog - The PetscStageLog

361:   Output Parameter:
362: . eventLog - The PetscEventRegLog

364:   Level: developer

366: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
367: @*/
368: PetscErrorCode  PetscStageLogGetEventRegLog(PetscStageLog stageLog, PetscEventRegLog *eventLog)
369: {
371:   *eventLog = stageLog->eventLog;
372:   return 0;
373: }

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

378:   Not Collective

380:   Input Parameters:
381: + stageLog - The PetscStageLog
382: - stage    - The stage

384:   Output Parameter:
385: . classLog - The PetscClassPerfLog

387:   Level: developer

389: .seealso: PetscStageLogPush(), PetscStageLogPop(), PetscLogGetStageLog()
390: @*/
391: PetscErrorCode  PetscStageLogGetClassPerfLog(PetscStageLog stageLog, int stage, PetscClassPerfLog *classLog)
392: {
395:   *classLog = stageLog->stageInfo[stage].classLog;
396:   return 0;
397: }

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

402:   Not Collective

404:   Input Parameters:
405: + stageLog - The PetscStageLog
406: . stage    - The stage to log
407: - isActive - The activity flag, PETSC_TRUE for logging, otherwise PETSC_FALSE (default is PETSC_TRUE)

409:   Level: developer

411: .seealso: PetscStageLogGetActive(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
412: @*/
413: PetscErrorCode  PetscStageLogSetActive(PetscStageLog stageLog, int stage, PetscBool isActive)
414: {
416:   stageLog->stageInfo[stage].perfInfo.active = isActive;
417:   return 0;
418: }

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

423:   Not Collective

425:   Input Parameters:
426: + stageLog - The PetscStageLog
427: - stage    - The stage to log

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

432:   Level: developer

434: .seealso: PetscStageLogSetActive(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
435: @*/
436: PetscErrorCode  PetscStageLogGetActive(PetscStageLog stageLog, int stage, PetscBool  *isActive)
437: {
440:   *isActive = stageLog->stageInfo[stage].perfInfo.active;
441:   return 0;
442: }

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

447:   Not Collective

449:   Input Parameters:
450: + stageLog  - The PetscStageLog
451: . stage     - The stage to log
452: - isVisible - The visibility flag, PETSC_TRUE for printing, otherwise PETSC_FALSE (default is PETSC_TRUE)

454:   Database Options:
455: . -log_view - Activates log summary

457:   Level: developer

459: .seealso: PetscStageLogGetVisible(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
460: @*/
461: PetscErrorCode  PetscStageLogSetVisible(PetscStageLog stageLog, int stage, PetscBool isVisible)
462: {
464:   stageLog->stageInfo[stage].perfInfo.visible = isVisible;
465:   return 0;
466: }

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

471:   Not Collective

473:   Input Parameters:
474: + stageLog  - The PetscStageLog
475: - stage     - The stage to log

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

480:   Database Options:
481: . -log_view - Activates log summary

483:   Level: developer

485: .seealso: PetscStageLogSetVisible(), PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
486: @*/
487: PetscErrorCode  PetscStageLogGetVisible(PetscStageLog stageLog, int stage, PetscBool  *isVisible)
488: {
491:   *isVisible = stageLog->stageInfo[stage].perfInfo.visible;
492:   return 0;
493: }

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

498:   Not Collective

500:   Input Parameters:
501: + stageLog - The PetscStageLog
502: - name     - The stage name

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

507:   Level: developer

509: .seealso: PetscStageLogGetCurrent(), PetscStageLogRegister(), PetscLogGetStageLog()
510: @*/
511: PetscErrorCode  PetscStageLogGetStage(PetscStageLog stageLog, const char name[], PetscLogStage *stage)
512: {
513:   PetscBool      match;
514:   int            s;

518:   *stage = -1;
519:   for (s = 0; s < stageLog->numStages; s++) {
520:     PetscStrcasecmp(stageLog->stageInfo[s].name, name, &match);
521:     if (match) {
522:       *stage = s;
523:       break;
524:     }
525:   }
526:   return 0;
527: }

529: /*@C
530:   PetscStageLogCreate - This creates a PetscStageLog object.

532:   Not collective

534:   Input Parameter:
535: . stageLog - The PetscStageLog

537:   Level: developer

539: .seealso: PetscStageLogCreate()
540: @*/
541: PetscErrorCode  PetscStageLogCreate(PetscStageLog *stageLog)
542: {
543:   PetscStageLog  l;

545:   PetscNew(&l);

547:   l->numStages = 0;
548:   l->maxStages = 10;
549:   l->curStage  = -1;

551:   PetscIntStackCreate(&l->stack);
552:   PetscMalloc1(l->maxStages, &l->stageInfo);
553:   PetscEventRegLogCreate(&l->eventLog);
554:   PetscClassRegLogCreate(&l->classLog);

556:   *stageLog = l;
557:   return 0;
558: }