Actual source code: logdefault.c
1: #include <petscviewer.h>
2: #include <petscdevice.h>
3: #include <petsc/private/logimpl.h>
4: #include <petsc/private/loghandlerimpl.h>
5: #include <petsc/private/deviceimpl.h>
6: #include <petscconfiginfo.h>
7: #include <petscmachineinfo.h>
8: #include "logdefault.h"
10: static PetscErrorCode PetscEventPerfInfoInit(PetscEventPerfInfo *eventInfo)
11: {
12: PetscFunctionBegin;
13: PetscCall(PetscMemzero(eventInfo, sizeof(*eventInfo)));
14: eventInfo->visible = PETSC_TRUE;
15: eventInfo->id = -1;
16: eventInfo->dof[0] = -1.0;
17: eventInfo->dof[1] = -1.0;
18: eventInfo->dof[2] = -1.0;
19: eventInfo->dof[3] = -1.0;
20: eventInfo->dof[4] = -1.0;
21: eventInfo->dof[5] = -1.0;
22: eventInfo->dof[6] = -1.0;
23: eventInfo->dof[7] = -1.0;
24: eventInfo->errors[0] = -1.0;
25: eventInfo->errors[1] = -1.0;
26: eventInfo->errors[2] = -1.0;
27: eventInfo->errors[3] = -1.0;
28: eventInfo->errors[4] = -1.0;
29: eventInfo->errors[5] = -1.0;
30: eventInfo->errors[6] = -1.0;
31: eventInfo->errors[7] = -1.0;
32: PetscFunctionReturn(PETSC_SUCCESS);
33: }
35: static PetscErrorCode PetscEventPerfInfoTic_Internal(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event, PetscBool resume)
36: {
37: PetscFunctionBegin;
38: if (resume) {
39: eventInfo->timeTmp -= time;
40: eventInfo->flopsTmp -= petsc_TotalFlops_th;
41: } else {
42: eventInfo->timeTmp = -time;
43: eventInfo->flopsTmp = -petsc_TotalFlops_th;
44: }
45: eventInfo->numMessages -= petsc_irecv_ct_th + petsc_isend_ct_th + petsc_recv_ct_th + petsc_send_ct_th;
46: eventInfo->messageLength -= petsc_irecv_len_th + petsc_isend_len_th + petsc_recv_len_th + petsc_send_len_th;
47: eventInfo->numReductions -= petsc_allreduce_ct_th + petsc_gather_ct_th + petsc_scatter_ct_th;
48: #if defined(PETSC_HAVE_DEVICE)
49: eventInfo->CpuToGpuCount -= petsc_ctog_ct_th;
50: eventInfo->GpuToCpuCount -= petsc_gtoc_ct_th;
51: eventInfo->CpuToGpuSize -= petsc_ctog_sz_th;
52: eventInfo->GpuToCpuSize -= petsc_gtoc_sz_th;
53: eventInfo->GpuFlops -= petsc_gflops_th;
54: eventInfo->GpuTime -= petsc_gtime;
55: #endif
56: if (logMemory) {
57: PetscLogDouble usage;
58: PetscCall(PetscMemoryGetCurrentUsage(&usage));
59: eventInfo->memIncrease -= usage;
60: PetscCall(PetscMallocGetCurrentUsage(&usage));
61: eventInfo->mallocSpace -= usage;
62: PetscCall(PetscMallocGetMaximumUsage(&usage));
63: eventInfo->mallocIncrease -= usage;
64: PetscCall(PetscMallocPushMaximumUsage(event));
65: }
66: PetscFunctionReturn(PETSC_SUCCESS);
67: }
69: static PetscErrorCode PetscEventPerfInfoTic(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event)
70: {
71: PetscFunctionBegin;
72: PetscCall(PetscEventPerfInfoTic_Internal(eventInfo, time, logMemory, event, PETSC_FALSE));
73: PetscFunctionReturn(PETSC_SUCCESS);
74: }
76: static PetscErrorCode PetscEventPerfInfoResume(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event)
77: {
78: PetscFunctionBegin;
79: PetscCall(PetscEventPerfInfoTic_Internal(eventInfo, time, logMemory, event, PETSC_TRUE));
80: PetscFunctionReturn(PETSC_SUCCESS);
81: }
83: static PetscErrorCode PetscEventPerfInfoToc_Internal(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event, PetscBool pause)
84: {
85: PetscFunctionBegin;
86: eventInfo->timeTmp += time;
87: eventInfo->flopsTmp += petsc_TotalFlops_th;
88: if (!pause) {
89: eventInfo->time += eventInfo->timeTmp;
90: eventInfo->time2 += eventInfo->timeTmp * eventInfo->timeTmp;
91: eventInfo->flops += eventInfo->flopsTmp;
92: eventInfo->flops2 += eventInfo->flopsTmp * eventInfo->flopsTmp;
93: }
94: eventInfo->numMessages += petsc_irecv_ct_th + petsc_isend_ct_th + petsc_recv_ct_th + petsc_send_ct_th;
95: eventInfo->messageLength += petsc_irecv_len_th + petsc_isend_len_th + petsc_recv_len + petsc_send_len_th;
96: eventInfo->numReductions += petsc_allreduce_ct_th + petsc_gather_ct_th + petsc_scatter_ct_th;
97: #if defined(PETSC_HAVE_DEVICE)
98: eventInfo->CpuToGpuCount += petsc_ctog_ct_th;
99: eventInfo->GpuToCpuCount += petsc_gtoc_ct_th;
100: eventInfo->CpuToGpuSize += petsc_ctog_sz_th;
101: eventInfo->GpuToCpuSize += petsc_gtoc_sz_th;
102: eventInfo->GpuFlops += petsc_gflops_th;
103: eventInfo->GpuTime += petsc_gtime;
104: #endif
105: if (logMemory) {
106: PetscLogDouble usage, musage;
107: PetscCall(PetscMemoryGetCurrentUsage(&usage)); /* the comments below match the column labels printed in PetscLogView_Default() */
108: eventInfo->memIncrease += usage; /* RMI */
109: PetscCall(PetscMallocGetCurrentUsage(&usage));
110: eventInfo->mallocSpace += usage; /* Malloc */
111: PetscCall(PetscMallocPopMaximumUsage((int)event, &musage));
112: eventInfo->mallocIncreaseEvent = PetscMax(musage - usage, eventInfo->mallocIncreaseEvent); /* EMalloc */
113: PetscCall(PetscMallocGetMaximumUsage(&usage));
114: eventInfo->mallocIncrease += usage; /* MMalloc */
115: }
116: PetscFunctionReturn(PETSC_SUCCESS);
117: }
119: static PetscErrorCode PetscEventPerfInfoToc(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event)
120: {
121: PetscFunctionBegin;
122: PetscCall(PetscEventPerfInfoToc_Internal(eventInfo, time, logMemory, event, PETSC_FALSE));
123: PetscFunctionReturn(PETSC_SUCCESS);
124: }
126: static PetscErrorCode PetscEventPerfInfoPause(PetscEventPerfInfo *eventInfo, PetscLogDouble time, PetscBool logMemory, int event)
127: {
128: PetscFunctionBegin;
129: PetscCall(PetscEventPerfInfoToc_Internal(eventInfo, time, logMemory, event, PETSC_TRUE));
130: PetscFunctionReturn(PETSC_SUCCESS);
131: }
133: static PetscErrorCode PetscEventPerfInfoAdd_Internal(const PetscEventPerfInfo *eventInfo, PetscEventPerfInfo *outInfo)
134: {
135: PetscFunctionBegin;
136: outInfo->count += eventInfo->count;
137: outInfo->time += eventInfo->time;
138: outInfo->time2 += eventInfo->time2;
139: outInfo->flops += eventInfo->flops;
140: outInfo->flops2 += eventInfo->flops2;
141: outInfo->numMessages += eventInfo->numMessages;
142: outInfo->messageLength += eventInfo->messageLength;
143: outInfo->numReductions += eventInfo->numReductions;
144: #if defined(PETSC_HAVE_DEVICE)
145: outInfo->CpuToGpuCount += eventInfo->CpuToGpuCount;
146: outInfo->GpuToCpuCount += eventInfo->GpuToCpuCount;
147: outInfo->CpuToGpuSize += eventInfo->CpuToGpuSize;
148: outInfo->GpuToCpuSize += eventInfo->GpuToCpuSize;
149: outInfo->GpuFlops += eventInfo->GpuFlops;
150: outInfo->GpuTime += eventInfo->GpuTime;
151: #endif
152: outInfo->memIncrease += eventInfo->memIncrease;
153: outInfo->mallocSpace += eventInfo->mallocSpace;
154: outInfo->mallocIncreaseEvent += eventInfo->mallocIncreaseEvent;
155: outInfo->mallocIncrease += eventInfo->mallocIncrease;
156: PetscFunctionReturn(PETSC_SUCCESS);
157: }
159: PETSC_LOG_RESIZABLE_ARRAY(EventPerfArray, PetscEventPerfInfo, PetscLogEvent, PetscEventPerfInfoInit, NULL, NULL)
161: /* --- PetscClassPerf --- */
163: typedef struct {
164: int creations; /* The number of objects of this class created */
165: int destructions; /* The number of objects of this class destroyed */
166: PetscLogDouble mem; /* The total memory allocated by objects of this class; this is completely wrong and should possibly be removed */
167: PetscLogDouble descMem; /* The total memory allocated by descendents of these objects; this is completely wrong and should possibly be removed */
168: } PetscClassPerf;
170: static PetscErrorCode PetscClassPerfInit(PetscClassPerf *classInfo)
171: {
172: PetscFunctionBegin;
173: PetscCall(PetscMemzero(classInfo, sizeof(*classInfo)));
174: PetscFunctionReturn(PETSC_SUCCESS);
175: }
177: PETSC_LOG_RESIZABLE_ARRAY(ClassPerfArray, PetscClassPerf, PetscLogClass, PetscClassPerfInit, NULL, NULL)
179: /* --- PetscStagePerf --- */
181: typedef struct _PetscStagePerf {
182: PetscBool used; /* The stage was pushed on this processor */
183: PetscEventPerfInfo perfInfo; /* The stage performance information */
184: PetscLogEventPerfArray eventLog; /* The event information for this stage */
185: PetscLogClassPerfArray classLog; /* The class information for this stage */
186: } PetscStagePerf;
188: static PetscErrorCode PetscStageInfoInit(PetscStagePerf *stageInfo)
189: {
190: PetscFunctionBegin;
191: PetscCall(PetscMemzero(stageInfo, sizeof(*stageInfo)));
192: PetscCall(PetscLogEventPerfArrayCreate(128, &stageInfo->eventLog));
193: PetscCall(PetscLogClassPerfArrayCreate(128, &stageInfo->classLog));
194: PetscCall(PetscEventPerfInfoInit(&stageInfo->perfInfo));
195: PetscFunctionReturn(PETSC_SUCCESS);
196: }
198: static PetscErrorCode PetscStageInfoReset(PetscStagePerf *stageInfo)
199: {
200: PetscFunctionBegin;
201: PetscCall(PetscLogEventPerfArrayDestroy(&stageInfo->eventLog));
202: PetscCall(PetscLogClassPerfArrayDestroy(&stageInfo->classLog));
203: PetscFunctionReturn(PETSC_SUCCESS);
204: }
206: PETSC_LOG_RESIZABLE_ARRAY(StageInfoArray, PetscStagePerf, PetscLogStage, PetscStageInfoInit, PetscStageInfoReset, NULL)
208: /* --- Action --- */
210: /* The structure for action logging */
211: typedef enum {
212: PETSC_LOG_ACTION_CREATE,
213: PETSC_LOG_ACTION_DESTROY,
214: PETSC_LOG_ACTION_BEGIN,
215: PETSC_LOG_ACTION_END,
216: } PetscLogActionType;
218: typedef struct _Action {
219: PetscLogActionType action; /* The type of execution */
220: PetscLogEvent event; /* The event number */
221: PetscClassId classid; /* The event class id */
222: PetscLogDouble time; /* The time of occurrence */
223: PetscLogDouble flops; /* The cumulative flops */
224: PetscLogDouble mem; /* The current memory usage */
225: PetscLogDouble maxmem; /* The maximum memory usage */
226: int id1, id2, id3; /* The ids of associated objects */
227: } Action;
229: PETSC_LOG_RESIZABLE_ARRAY(ActionArray, Action, PetscLogEvent, NULL, NULL, NULL)
231: /* --- Object --- */
233: /* The structure for object logging */
234: typedef struct _Object {
235: PetscObject obj; /* The associated PetscObject */
236: int parent; /* The parent id */
237: PetscLogDouble mem; /* The memory associated with the object */
238: char name[64]; /* The object name */
239: char info[64]; /* The information string */
240: } Object;
242: PETSC_LOG_RESIZABLE_ARRAY(ObjectArray, Object, PetscObject, NULL, NULL, NULL)
244: /* Map from (threadid,stage,event) to perfInfo data struct */
245: #include <petsc/private/hashmapijk.h>
247: PETSC_HASH_MAP(HMapEvent, PetscHashIJKKey, PetscEventPerfInfo *, PetscHashIJKKeyHash, PetscHashIJKKeyEqual, NULL)
249: typedef struct _n_PetscLogHandler_Default *PetscLogHandler_Default;
250: struct _n_PetscLogHandler_Default {
251: PetscLogStageInfoArray stages;
252: PetscSpinlock lock;
253: PetscLogActionArray petsc_actions;
254: PetscLogObjectArray petsc_objects;
255: PetscBool petsc_logActions;
256: PetscBool petsc_logObjects;
257: int petsc_numObjectsCreated;
258: int petsc_numObjectsDestroyed;
259: PetscHMapEvent eventInfoMap_th;
260: int pause_depth;
261: PetscBool use_threadsafe;
262: };
264: /* --- PetscLogHandler_Default --- */
266: static PetscErrorCode PetscLogHandlerContextCreate_Default(PetscLogHandler_Default *def_p)
267: {
268: PetscLogHandler_Default def;
270: PetscFunctionBegin;
271: PetscCall(PetscNew(def_p));
272: def = *def_p;
273: PetscCall(PetscLogStageInfoArrayCreate(8, &def->stages));
274: PetscCall(PetscLogActionArrayCreate(64, &def->petsc_actions));
275: PetscCall(PetscLogObjectArrayCreate(64, &def->petsc_objects));
276: PetscCall(PetscSpinlockCreate(&def->lock));
278: PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_include_actions", &def->petsc_logActions, NULL));
279: PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_include_objects", &def->petsc_logObjects, NULL));
280: PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_handler_default_use_threadsafe_events", &def->use_threadsafe, NULL));
281: if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { PetscCall(PetscHMapEventCreate(&def->eventInfoMap_th)); }
282: PetscFunctionReturn(PETSC_SUCCESS);
283: }
285: static PetscErrorCode PetscLogHandlerDestroy_Default(PetscLogHandler h)
286: {
287: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
289: PetscFunctionBegin;
290: PetscCall(PetscLogStageInfoArrayDestroy(&def->stages));
291: PetscCall(PetscLogActionArrayDestroy(&def->petsc_actions));
292: PetscCall(PetscLogObjectArrayDestroy(&def->petsc_objects));
293: PetscCall(PetscSpinlockDestroy(&def->lock));
294: if (def->eventInfoMap_th) {
295: PetscEventPerfInfo **array;
296: PetscInt n, off = 0;
298: PetscCall(PetscHMapEventGetSize(def->eventInfoMap_th, &n));
299: PetscCall(PetscMalloc1(n, &array));
300: PetscCall(PetscHMapEventGetVals(def->eventInfoMap_th, &off, array));
301: for (PetscInt i = 0; i < n; i++) PetscCall(PetscFree(array[i]));
302: PetscCall(PetscFree(array));
303: PetscCall(PetscHMapEventDestroy(&def->eventInfoMap_th));
304: }
305: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerGetEventPerfInfo_C", NULL));
306: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerGetStagePerfInfo_C", NULL));
307: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerSetLogActions_C", NULL));
308: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerSetLogObjects_C", NULL));
309: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerLogObjectState_C", NULL));
310: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerGetNumObjects_C", NULL));
311: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventDeactivatePush_C", NULL));
312: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventDeactivatePop_C", NULL));
313: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventsPause_C", NULL));
314: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventsResume_C", NULL));
315: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerDump_C", NULL));
316: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerStageSetVisible_C", NULL));
317: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerStageGetVisible_C", NULL));
318: PetscCall(PetscFree(def));
319: PetscFunctionReturn(PETSC_SUCCESS);
320: }
322: static PetscErrorCode PetscLogHandlerDefaultGetStageInfo(PetscLogHandler handler, PetscLogStage stage, PetscStagePerf **stage_info_p)
323: {
324: PetscStagePerf *stage_info = NULL;
325: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
327: PetscFunctionBegin;
328: PetscCall(PetscLogStageInfoArrayResize(def->stages, stage + 1));
329: PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
330: *stage_info_p = stage_info;
331: PetscFunctionReturn(PETSC_SUCCESS);
332: }
334: static PetscErrorCode PetscLogHandlerGetEventPerfInfo_Default(PetscLogHandler handler, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **event_info_p)
335: {
336: PetscEventPerfInfo *event_info = NULL;
337: PetscStagePerf *stage_info = NULL;
338: PetscLogEventPerfArray event_log;
340: PetscFunctionBegin;
341: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(handler->state, &stage));
342: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info));
343: event_log = stage_info->eventLog;
344: PetscCall(PetscLogEventPerfArrayResize(event_log, event + 1));
345: PetscCall(PetscLogEventPerfArrayGetRef(event_log, event, &event_info));
346: event_info->id = event;
347: *event_info_p = event_info;
348: PetscFunctionReturn(PETSC_SUCCESS);
349: }
351: static PetscErrorCode PetscLogHandlerGetStagePerfInfo_Default(PetscLogHandler handler, PetscLogStage stage, PetscEventPerfInfo **stage_info_p)
352: {
353: PetscStagePerf *stage_perf_info = NULL;
355: PetscFunctionBegin;
356: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(handler->state, &stage));
357: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_perf_info));
358: *stage_info_p = &stage_perf_info->perfInfo;
359: PetscFunctionReturn(PETSC_SUCCESS);
360: }
362: static PetscErrorCode PetscLogHandlerDefaultGetClassPerf(PetscLogHandler handler, PetscLogStage stage, PetscLogClass clss, PetscClassPerf **class_info)
363: {
364: PetscLogClassPerfArray class_log;
365: PetscStagePerf *stage_info;
367: PetscFunctionBegin;
368: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info));
369: class_log = stage_info->classLog;
370: PetscCall(PetscLogClassPerfArrayResize(class_log, clss + 1));
371: PetscCall(PetscLogClassPerfArrayGetRef(class_log, clss, class_info));
372: PetscFunctionReturn(PETSC_SUCCESS);
373: }
375: static PetscErrorCode PetscLogHandlerObjectCreate_Default(PetscLogHandler h, PetscObject obj)
376: {
377: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
378: PetscLogState state;
379: PetscLogStage stage;
380: PetscClassPerf *classInfo;
381: int oclass = 0;
383: PetscFunctionBegin;
384: PetscCall(PetscLogHandlerGetState(h, &state));
385: PetscCall(PetscSpinlockLock(&def->lock));
386: /* Record stage info */
387: PetscCall(PetscLogStateGetCurrentStage(state, &stage));
388: PetscCall(PetscLogStateGetClassFromClassId(state, obj->classid, &oclass));
389: PetscCall(PetscLogHandlerDefaultGetClassPerf(h, stage, oclass, &classInfo));
390: classInfo->creations++;
391: /* Record the creation action */
392: if (def->petsc_logActions) {
393: Action new_action;
395: PetscCall(PetscTime(&new_action.time));
396: new_action.time -= petsc_BaseTime;
397: new_action.action = PETSC_LOG_ACTION_CREATE;
398: new_action.event = -1;
399: new_action.classid = obj->classid;
400: new_action.id1 = obj->id;
401: new_action.id2 = -1;
402: new_action.id3 = -1;
403: new_action.flops = petsc_TotalFlops;
404: PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
405: PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
406: PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
407: }
408: /* We don't just use obj->id to count all objects that are created
409: because PetscLogHandlers are objects and PetscLogObjectDestroy() will not
410: be called for them: the number of objects created and destroyed as counted
411: here and below would have an imbalance */
412: def->petsc_numObjectsCreated++;
413: /* Record the object */
414: if (def->petsc_logObjects) {
415: Object new_object;
417: new_object.parent = -1;
418: new_object.obj = obj;
419: new_object.mem = 0;
421: PetscCall(PetscMemzero(new_object.name, sizeof(new_object.name)));
422: PetscCall(PetscMemzero(new_object.info, sizeof(new_object.info)));
423: PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1");
424: PetscCall(PetscLogObjectArrayResize(def->petsc_objects, obj->id));
425: PetscCall(PetscLogObjectArraySet(def->petsc_objects, obj->id - 1, new_object));
426: }
427: PetscCall(PetscSpinlockUnlock(&def->lock));
428: PetscFunctionReturn(PETSC_SUCCESS);
429: }
431: static PetscErrorCode PetscLogHandlerObjectDestroy_Default(PetscLogHandler h, PetscObject obj)
432: {
433: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
434: PetscLogState state;
435: PetscLogStage stage;
436: PetscClassPerf *classInfo;
437: int oclass = 0;
439: PetscFunctionBegin;
440: PetscCall(PetscLogHandlerGetState(h, &state));
441: /* Record stage info */
442: PetscCall(PetscSpinlockLock(&def->lock));
443: PetscCall(PetscLogStateGetCurrentStage(state, &stage));
444: if (stage >= 0) {
445: /* stage < 0 can happen if the log summary is output before some things are destroyed */
446: PetscCall(PetscLogStateGetClassFromClassId(state, obj->classid, &oclass));
447: PetscCall(PetscLogHandlerDefaultGetClassPerf(h, stage, oclass, &classInfo));
448: classInfo->destructions++;
449: }
450: /* Cannot Credit all ancestors with your memory because they may have already been destroyed*/
451: def->petsc_numObjectsDestroyed++;
452: /* Dynamically enlarge logging structures */
453: /* Record the destruction action */
454: if (def->petsc_logActions) {
455: Action new_action;
457: PetscCall(PetscTime(&new_action.time));
458: new_action.time -= petsc_BaseTime;
459: new_action.event = -1;
460: new_action.action = PETSC_LOG_ACTION_DESTROY;
461: new_action.classid = obj->classid;
462: new_action.id1 = obj->id;
463: new_action.id2 = -1;
464: new_action.id3 = -1;
465: new_action.flops = petsc_TotalFlops;
466: PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
467: PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
468: PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
469: }
470: if (def->petsc_logObjects) {
471: Object *obj_entry = NULL;
473: PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1");
474: PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, obj->id - 1, &obj_entry));
475: if (obj->name) PetscCall(PetscStrncpy(obj_entry->name, obj->name, 64));
476: obj_entry->obj = NULL;
477: }
478: PetscCall(PetscSpinlockUnlock(&def->lock));
479: PetscFunctionReturn(PETSC_SUCCESS);
480: }
482: static PetscErrorCode PetscLogHandlerEventSync_Default(PetscLogHandler h, PetscLogEvent event, MPI_Comm comm)
483: {
484: PetscLogState state;
485: PetscLogEventInfo event_info;
486: PetscEventPerfInfo *event_perf_info;
487: int stage;
488: PetscLogDouble time = 0.0;
490: PetscFunctionBegin;
491: if (!PetscLogSyncOn || comm == MPI_COMM_NULL) PetscFunctionReturn(PETSC_SUCCESS);
492: PetscCall(PetscLogHandlerGetState(h, &state));
493: PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
494: if (!event_info.collective) PetscFunctionReturn(PETSC_SUCCESS);
495: PetscCall(PetscLogStateGetCurrentStage(state, &stage));
496: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
497: if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
499: PetscCall(PetscTimeSubtract(&time));
500: PetscCallMPI(MPI_Barrier(comm));
501: PetscCall(PetscTimeAdd(&time));
502: event_perf_info->syncTime += time;
503: PetscFunctionReturn(PETSC_SUCCESS);
504: }
506: static PetscErrorCode PetscLogGetStageEventPerfInfo_threaded(PetscLogHandler_Default def, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **eventInfo)
507: {
508: PetscEventPerfInfo *leventInfo = NULL;
509: PetscHashIJKKey key;
511: PetscFunctionBegin;
512: #if PetscDefined(HAVE_THREADSAFETY)
513: key.i = PetscLogGetTid();
514: #else
515: key.i = 0;
516: #endif
517: key.j = stage;
518: key.k = event;
519: PetscCall(PetscSpinlockLock(&def->lock));
520: PetscCall(PetscHMapEventGet(def->eventInfoMap_th, key, &leventInfo));
521: if (!leventInfo) {
522: PetscCall(PetscNew(&leventInfo));
523: leventInfo->id = event;
524: PetscCall(PetscHMapEventSet(def->eventInfoMap_th, key, leventInfo));
525: }
526: PetscCall(PetscSpinlockUnlock(&def->lock));
527: *eventInfo = leventInfo;
528: PetscFunctionReturn(PETSC_SUCCESS);
529: }
531: static PetscErrorCode PetscLogHandlerEventBegin_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
532: {
533: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
534: PetscEventPerfInfo *event_perf_info = NULL;
535: PetscLogEventInfo event_info;
536: PetscLogDouble time;
537: PetscLogState state;
538: PetscLogStage stage;
540: PetscFunctionBegin;
541: PetscCall(PetscLogHandlerGetState(h, &state));
542: if (PetscDefined(USE_DEBUG)) {
543: PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
548: if (event_info.collective && o1) {
549: PetscInt64 b1[2], b2[2];
551: b1[0] = -o1->cidx;
552: b1[1] = o1->cidx;
553: PetscCall(MPIU_Allreduce(b1, b2, 2, MPIU_INT64, MPI_MAX, PetscObjectComm(o1)));
554: PetscCheck(-b2[0] == b2[1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Collective event %s not called collectively %" PetscInt64_FMT " != %" PetscInt64_FMT, event_info.name, -b2[0], b2[1]);
555: }
556: }
557: /* Synchronization */
558: PetscCall(PetscLogHandlerEventSync_Default(h, event, PetscObjectComm(o1)));
559: PetscCall(PetscLogStateGetCurrentStage(state, &stage));
560: if (def->pause_depth > 0) stage = 0; // in pause-mode, all events run on the main stage
561: if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
562: PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info));
563: if (event_perf_info->depth == 0) { PetscCall(PetscEventPerfInfoInit(event_perf_info)); }
564: } else {
565: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
566: }
567: PetscCheck(event_perf_info->depth >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to begin a paused event, this is not allowed");
568: event_perf_info->depth++;
569: /* Check for double counting */
570: if (event_perf_info->depth > 1) PetscFunctionReturn(PETSC_SUCCESS);
571: PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
572: /* Log the performance info */
573: event_perf_info->count++;
574: PetscCall(PetscTime(&time));
575: PetscCall(PetscEventPerfInfoTic(event_perf_info, time, PetscLogMemory, (int)event));
576: if (def->petsc_logActions) {
577: PetscLogDouble curTime;
578: Action new_action;
580: PetscCall(PetscTime(&curTime));
581: new_action.time = curTime - petsc_BaseTime;
582: new_action.action = PETSC_LOG_ACTION_BEGIN;
583: new_action.event = event;
584: new_action.classid = event_info.classid;
585: new_action.id1 = o1 ? o1->id : -1;
586: new_action.id2 = o2 ? o2->id : -1;
587: new_action.id3 = o3 ? o3->id : -1;
588: new_action.flops = petsc_TotalFlops;
589: PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
590: PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
591: PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
592: }
593: PetscFunctionReturn(PETSC_SUCCESS);
594: }
596: static PetscErrorCode PetscLogHandlerEventEnd_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
597: {
598: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
599: PetscEventPerfInfo *event_perf_info = NULL;
600: PetscLogDouble time;
601: PetscLogState state;
602: int stage;
603: PetscLogEventInfo event_info;
605: PetscFunctionBegin;
606: PetscCall(PetscLogHandlerGetState(h, &state));
607: if (PetscDefined(USE_DEBUG)) {
608: PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
613: if (event_info.collective && o1) {
614: PetscInt64 b1[2], b2[2];
616: b1[0] = -o1->cidx;
617: b1[1] = o1->cidx;
618: PetscCall(MPIU_Allreduce(b1, b2, 2, MPIU_INT64, MPI_MAX, PetscObjectComm(o1)));
619: PetscCheck(-b2[0] == b2[1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Collective event %s not called collectively %" PetscInt64_FMT " != %" PetscInt64_FMT, event_info.name, -b2[0], b2[1]);
620: }
621: }
622: if (def->petsc_logActions) {
623: PetscLogDouble curTime;
624: Action new_action;
626: PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
627: PetscCall(PetscTime(&curTime));
628: new_action.time = curTime - petsc_BaseTime;
629: new_action.action = PETSC_LOG_ACTION_END;
630: new_action.event = event;
631: new_action.classid = event_info.classid;
632: new_action.id1 = o1 ? o1->id : -1;
633: new_action.id2 = o2 ? o2->id : -2;
634: new_action.id3 = o3 ? o3->id : -3;
635: new_action.flops = petsc_TotalFlops;
636: PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
637: PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
638: PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
639: }
640: PetscCall(PetscLogStateGetCurrentStage(state, &stage));
641: if (def->pause_depth > 0) stage = 0; // all events run on the main stage in pause-mode
642: if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
643: PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info));
644: } else {
645: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
646: }
647: PetscCheck(event_perf_info->depth > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to end paused event, not allowed");
648: event_perf_info->depth--;
649: /* Check for double counting */
650: if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
651: else PetscCheck(event_perf_info->depth == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Logging event had unbalanced begin/end pairs");
653: /* Log performance info */
654: PetscCall(PetscTime(&time));
655: PetscCall(PetscEventPerfInfoToc(event_perf_info, time, PetscLogMemory, (int)event));
656: if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
657: PetscEventPerfInfo *event_perf_info_global;
658: PetscCall(PetscSpinlockLock(&def->lock));
659: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info_global));
660: PetscCall(PetscEventPerfInfoAdd_Internal(event_perf_info, event_perf_info_global));
661: PetscCall(PetscSpinlockUnlock(&def->lock));
662: }
663: PetscFunctionReturn(PETSC_SUCCESS);
664: }
666: static PetscErrorCode PetscLogHandlerEventDeactivatePush_Default(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event)
667: {
668: PetscEventPerfInfo *event_perf_info;
670: PetscFunctionBegin;
671: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
672: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
673: event_perf_info->depth++;
674: PetscFunctionReturn(PETSC_SUCCESS);
675: }
677: static PetscErrorCode PetscLogHandlerEventDeactivatePop_Default(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event)
678: {
679: PetscEventPerfInfo *event_perf_info;
681: PetscFunctionBegin;
682: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
683: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
684: event_perf_info->depth--;
685: PetscFunctionReturn(PETSC_SUCCESS);
686: }
688: static PetscErrorCode PetscLogHandlerEventsPause_Default(PetscLogHandler h)
689: {
690: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
691: PetscLogDouble time;
692: PetscInt num_stages;
694: PetscFunctionBegin;
695: if (def->pause_depth++ > 0) PetscFunctionReturn(PETSC_SUCCESS);
696: PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL));
697: PetscCall(PetscTime(&time));
698: for (PetscInt stage = 0; stage < num_stages; stage++) {
699: PetscStagePerf *stage_info = NULL;
700: PetscInt num_events;
702: PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
703: PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL));
704: for (PetscInt event = 0; event < num_events; event++) {
705: PetscEventPerfInfo *event_info = NULL;
706: PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info));
707: if (event_info->depth > 0) {
708: event_info->depth *= -1;
709: PetscCall(PetscEventPerfInfoPause(event_info, time, PetscLogMemory, event));
710: }
711: }
712: if (stage > 0 && stage_info->perfInfo.depth > 0) {
713: stage_info->perfInfo.depth *= -1;
714: PetscCall(PetscEventPerfInfoPause(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2)));
715: }
716: }
717: PetscFunctionReturn(PETSC_SUCCESS);
718: }
720: static PetscErrorCode PetscLogHandlerEventsResume_Default(PetscLogHandler h)
721: {
722: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
723: PetscLogDouble time;
724: PetscInt num_stages;
726: PetscFunctionBegin;
727: if (--def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
728: PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL));
729: PetscCall(PetscTime(&time));
730: for (PetscInt stage = 0; stage < num_stages; stage++) {
731: PetscStagePerf *stage_info = NULL;
732: PetscInt num_events;
734: PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
735: PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL));
736: for (PetscInt event = 0; event < num_events; event++) {
737: PetscEventPerfInfo *event_info = NULL;
738: PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info));
739: if (event_info->depth < 0) {
740: event_info->depth *= -1;
741: PetscCall(PetscEventPerfInfoResume(event_info, time, PetscLogMemory, event));
742: }
743: }
744: if (stage > 0 && stage_info->perfInfo.depth < 0) {
745: stage_info->perfInfo.depth *= -1;
746: PetscCall(PetscEventPerfInfoResume(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2)));
747: }
748: }
749: PetscFunctionReturn(PETSC_SUCCESS);
750: }
752: static PetscErrorCode PetscLogHandlerStagePush_Default(PetscLogHandler h, PetscLogStage new_stage)
753: {
754: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
755: PetscLogDouble time;
756: PetscLogState state;
757: PetscLogStage current_stage;
758: PetscStagePerf *new_stage_info;
760: PetscFunctionBegin;
761: if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
762: PetscCall(PetscLogHandlerGetState(h, &state));
763: current_stage = state->current_stage;
764: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, new_stage, &new_stage_info));
765: PetscCall(PetscTime(&time));
767: /* Record flops/time of previous stage */
768: if (current_stage >= 0) {
769: if (PetscBTLookup(state->active, current_stage)) {
770: PetscStagePerf *current_stage_info;
771: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, ¤t_stage_info));
772: PetscCall(PetscEventPerfInfoToc(¤t_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2)));
773: }
774: }
775: new_stage_info->used = PETSC_TRUE;
776: new_stage_info->perfInfo.count++;
777: new_stage_info->perfInfo.depth++;
778: /* Subtract current quantities so that we obtain the difference when we pop */
779: if (PetscBTLookup(state->active, new_stage)) PetscCall(PetscEventPerfInfoTic(&new_stage_info->perfInfo, time, PetscLogMemory, (int)-(new_stage + 2)));
780: PetscFunctionReturn(PETSC_SUCCESS);
781: }
783: static PetscErrorCode PetscLogHandlerStagePop_Default(PetscLogHandler h, PetscLogStage old_stage)
784: {
785: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
786: PetscLogStage current_stage;
787: PetscStagePerf *old_stage_info;
788: PetscLogState state;
789: PetscLogDouble time;
791: PetscFunctionBegin;
792: if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
793: PetscCall(PetscLogHandlerGetState(h, &state));
794: current_stage = state->current_stage;
795: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, old_stage, &old_stage_info));
796: PetscCall(PetscTime(&time));
797: old_stage_info->perfInfo.depth--;
798: if (PetscBTLookup(state->active, old_stage)) { PetscCall(PetscEventPerfInfoToc(&old_stage_info->perfInfo, time, PetscLogMemory, (int)-(old_stage + 2))); }
799: if (current_stage >= 0) {
800: if (PetscBTLookup(state->active, current_stage)) {
801: PetscStagePerf *current_stage_info;
802: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, ¤t_stage_info));
803: PetscCall(PetscEventPerfInfoTic(¤t_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2)));
804: }
805: }
806: PetscFunctionReturn(PETSC_SUCCESS);
807: }
809: static PetscErrorCode PetscLogHandlerStageSetVisible_Default(PetscLogHandler h, PetscLogStage stage, PetscBool is_visible)
810: {
811: PetscStagePerf *stage_info;
813: PetscFunctionBegin;
814: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
815: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info));
816: stage_info->perfInfo.visible = is_visible;
817: PetscFunctionReturn(PETSC_SUCCESS);
818: }
820: static PetscErrorCode PetscLogHandlerStageGetVisible_Default(PetscLogHandler h, PetscLogStage stage, PetscBool *is_visible)
821: {
822: PetscStagePerf *stage_info;
824: PetscFunctionBegin;
825: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
826: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info));
827: *is_visible = stage_info->perfInfo.visible;
828: PetscFunctionReturn(PETSC_SUCCESS);
829: }
831: static PetscErrorCode PetscLogHandlerSetLogActions_Default(PetscLogHandler handler, PetscBool flag)
832: {
833: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
835: PetscFunctionBegin;
836: def->petsc_logActions = flag;
837: PetscFunctionReturn(PETSC_SUCCESS);
838: }
840: static PetscErrorCode PetscLogHandlerSetLogObjects_Default(PetscLogHandler handler, PetscBool flag)
841: {
842: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
844: PetscFunctionBegin;
845: def->petsc_logObjects = flag;
846: PetscFunctionReturn(PETSC_SUCCESS);
847: }
849: static PetscErrorCode PetscLogHandlerLogObjectState_Default(PetscLogHandler handler, PetscObject obj, const char format[], va_list Argp)
850: {
851: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
852: size_t fullLength;
854: PetscFunctionBegin;
855: if (def->petsc_logObjects) {
856: Object *obj_entry = NULL;
858: PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, obj->id - 1, &obj_entry));
859: PetscCall(PetscVSNPrintf(obj_entry->info, 64, format, &fullLength, Argp));
860: }
861: PetscFunctionReturn(PETSC_SUCCESS);
862: }
864: static PetscErrorCode PetscLogHandlerGetNumObjects_Default(PetscLogHandler handler, PetscInt *num_objects)
865: {
866: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
868: PetscFunctionBegin;
869: PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, num_objects, NULL));
870: PetscFunctionReturn(PETSC_SUCCESS);
871: }
873: static PetscErrorCode PetscLogHandlerDump_Default(PetscLogHandler handler, const char sname[])
874: {
875: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
876: FILE *fd;
877: char file[PETSC_MAX_PATH_LEN], fname[PETSC_MAX_PATH_LEN];
878: PetscLogDouble flops, _TotalTime;
879: PetscMPIInt rank;
880: int curStage;
881: PetscLogState state;
882: PetscInt num_events;
883: PetscLogEvent event;
885: PetscFunctionBegin;
886: /* Calculate the total elapsed time */
887: PetscCall(PetscTime(&_TotalTime));
888: _TotalTime -= petsc_BaseTime;
889: /* Open log file */
890: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)handler), &rank));
891: PetscCall(PetscSNPrintf(file, PETSC_STATIC_ARRAY_LENGTH(file), "%s.%d", sname && sname[0] ? sname : "Log", rank));
892: PetscCall(PetscFixFilename(file, fname));
893: PetscCall(PetscFOpen(PETSC_COMM_SELF, fname, "w", &fd));
894: PetscCheck(!(rank == 0) || !(!fd), PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file: %s", fname);
895: /* Output totals */
896: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Total Flop %14e %16.8e\n", petsc_TotalFlops, _TotalTime));
897: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Clock Resolution %g\n", 0.0));
898: /* Output actions */
899: if (def->petsc_logActions) {
900: PetscInt num_actions;
901: PetscCall(PetscLogActionArrayGetSize(def->petsc_actions, &num_actions, NULL));
902: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Actions accomplished %d\n", (int)num_actions));
903: for (int a = 0; a < num_actions; a++) {
904: Action *action;
906: PetscCall(PetscLogActionArrayGetRef(def->petsc_actions, a, &action));
907: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "%g %d %d %d %d %d %d %g %g %g\n", action->time, action->action, (int)action->event, (int)action->classid, action->id1, action->id2, action->id3, action->flops, action->mem, action->maxmem));
908: }
909: }
910: /* Output objects */
911: if (def->petsc_logObjects) {
912: PetscInt num_objects;
914: PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
915: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Objects created %d destroyed %d\n", def->petsc_numObjectsCreated, def->petsc_numObjectsDestroyed));
916: for (int o = 0; o < num_objects; o++) {
917: Object *object = NULL;
919: PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, o, &object));
920: if (object->parent != -1) continue; // object with this id wasn't logged, probably a PetscLogHandler
921: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Parent ID: %d Memory: %d\n", object->parent, (int)object->mem));
922: if (!object->name[0]) {
923: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Name\n"));
924: } else {
925: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Name: %s\n", object->name));
926: }
927: if (!object->info[0]) {
928: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Info\n"));
929: } else {
930: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Info: %s\n", object->info));
931: }
932: }
933: }
934: /* Output events */
935: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Event log:\n"));
936: PetscCall(PetscLogHandlerGetState(handler, &state));
937: PetscCall(PetscLogStateGetNumEvents(state, &num_events));
938: PetscCall(PetscLogStateGetCurrentStage(state, &curStage));
939: for (event = 0; event < num_events; event++) {
940: PetscEventPerfInfo *event_info;
942: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, curStage, event, &event_info));
943: if (event_info->time != 0.0) flops = event_info->flops / event_info->time;
944: else flops = 0.0;
945: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "%d %16d %16g %16g %16g\n", event, event_info->count, event_info->flops, event_info->time, flops));
946: }
947: PetscCall(PetscFClose(PETSC_COMM_SELF, fd));
948: PetscFunctionReturn(PETSC_SUCCESS);
949: }
951: /*
952: PetscLogView_Detailed - Each process prints the times for its own events
954: */
955: static PetscErrorCode PetscLogHandlerView_Default_Detailed(PetscLogHandler handler, PetscViewer viewer)
956: {
957: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
958: PetscLogDouble locTotalTime, numRed, maxMem;
959: PetscInt numStages, numEvents;
960: MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
961: PetscMPIInt rank, size;
962: PetscLogGlobalNames global_stages, global_events;
963: PetscLogState state;
964: PetscEventPerfInfo zero_info;
966: PetscFunctionBegin;
967: PetscCall(PetscLogHandlerGetState(handler, &state));
968: PetscCallMPI(MPI_Comm_size(comm, &size));
969: PetscCallMPI(MPI_Comm_rank(comm, &rank));
970: /* Must preserve reduction count before we go on */
971: numRed = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
972: /* Get the total elapsed time */
973: PetscCall(PetscTime(&locTotalTime));
974: locTotalTime -= petsc_BaseTime;
975: PetscCall(PetscViewerASCIIPrintf(viewer, "size = %d\n", size));
976: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalTimes = {}\n"));
977: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessages = {}\n"));
978: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessageLens = {}\n"));
979: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalReductions = {}\n"));
980: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalFlop = {}\n"));
981: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalObjects = {}\n"));
982: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMemory = {}\n"));
983: PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
984: PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
985: PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
986: PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
987: PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
988: PetscCall(PetscViewerASCIIPrintf(viewer, "Stages = {}\n"));
989: for (PetscInt stage = 0; stage < numStages; stage++) {
990: PetscInt stage_id;
991: const char *stage_name;
993: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
994: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
995: PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"] = {}\n", stage_name));
996: PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"summary\"] = {}\n", stage_name));
997: for (PetscInt event = 0; event < numEvents; event++) {
998: PetscEventPerfInfo *eventInfo = &zero_info;
999: PetscBool is_zero = PETSC_FALSE;
1000: PetscInt event_id;
1001: const char *event_name;
1003: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1004: PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1005: if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1006: is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE;
1007: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1008: if (!is_zero) { PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"%s\"] = {}\n", stage_name, event_name)); }
1009: }
1010: }
1011: PetscCall(PetscMallocGetMaximumUsage(&maxMem));
1012: PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1013: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalTimes[%d] = %g\n", rank, locTotalTime));
1014: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessages[%d] = %g\n", rank, (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct)));
1015: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessageLens[%d] = %g\n", rank, (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len)));
1016: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalReductions[%d] = %g\n", rank, numRed));
1017: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalFlop[%d] = %g\n", rank, petsc_TotalFlops));
1018: {
1019: PetscInt num_objects;
1021: PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
1022: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalObjects[%d] = %d\n", rank, (int)num_objects));
1023: }
1024: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMemory[%d] = %g\n", rank, maxMem));
1025: PetscCall(PetscViewerFlush(viewer));
1026: for (PetscInt stage = 0; stage < numStages; stage++) {
1027: PetscEventPerfInfo *stage_perf_info = &zero_info;
1028: PetscInt stage_id;
1029: const char *stage_name;
1031: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1032: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1033: if (stage_id >= 0) {
1034: PetscStagePerf *stage_info;
1035: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info));
1036: stage_perf_info = &stage_info->perfInfo;
1037: }
1038: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Stages[\"%s\"][\"summary\"][%d] = {\"time\" : %g, \"numMessages\" : %g, \"messageLength\" : %g, \"numReductions\" : %g, \"flop\" : %g}\n", stage_name, rank, stage_perf_info->time,
1039: stage_perf_info->numMessages, stage_perf_info->messageLength, stage_perf_info->numReductions, stage_perf_info->flops));
1040: for (PetscInt event = 0; event < numEvents; event++) {
1041: PetscEventPerfInfo *eventInfo = &zero_info;
1042: PetscBool is_zero = PETSC_FALSE;
1043: PetscInt event_id;
1044: const char *event_name;
1046: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1047: PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1048: if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1049: is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE;
1050: PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero));
1051: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1052: if (!is_zero) {
1053: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Stages[\"%s\"][\"%s\"][%d] = {\"count\" : %d, \"time\" : %g, \"syncTime\" : %g, \"numMessages\" : %g, \"messageLength\" : %g, \"numReductions\" : %g, \"flop\" : %g", stage_name, event_name, rank,
1054: eventInfo->count, eventInfo->time, eventInfo->syncTime, eventInfo->numMessages, eventInfo->messageLength, eventInfo->numReductions, eventInfo->flops));
1055: if (eventInfo->dof[0] >= 0.) {
1056: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"dof\" : ["));
1057: for (PetscInt d = 0; d < 8; ++d) {
1058: if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", "));
1059: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->dof[d]));
1060: }
1061: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]"));
1062: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"error\" : ["));
1063: for (PetscInt e = 0; e < 8; ++e) {
1064: if (e > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", "));
1065: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->errors[e]));
1066: }
1067: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]"));
1068: }
1069: }
1070: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "}\n"));
1071: }
1072: }
1073: PetscCall(PetscViewerFlush(viewer));
1074: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1075: PetscCall(PetscLogGlobalNamesDestroy(&global_events));
1076: PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1077: PetscFunctionReturn(PETSC_SUCCESS);
1078: }
1080: /*
1081: PetscLogView_CSV - Each process prints the times for its own events in Comma-Separated Value Format
1082: */
1083: static PetscErrorCode PetscLogHandlerView_Default_CSV(PetscLogHandler handler, PetscViewer viewer)
1084: {
1085: PetscLogDouble locTotalTime, maxMem;
1086: PetscInt numStages, numEvents, stage, event;
1087: MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
1088: PetscMPIInt rank, size;
1089: PetscLogGlobalNames global_stages, global_events;
1090: PetscLogState state;
1091: PetscEventPerfInfo zero_info;
1093: PetscFunctionBegin;
1094: PetscCall(PetscLogHandlerGetState(handler, &state));
1095: PetscCallMPI(MPI_Comm_size(comm, &size));
1096: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1097: /* Must preserve reduction count before we go on */
1098: /* Get the total elapsed time */
1099: PetscCall(PetscTime(&locTotalTime));
1100: locTotalTime -= petsc_BaseTime;
1101: PetscCall(PetscMallocGetMaximumUsage(&maxMem));
1102: PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
1103: PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
1104: PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
1105: PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
1106: PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
1107: PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1108: PetscCall(PetscViewerASCIIPrintf(viewer, "Stage Name,Event Name,Rank,Count,Time,Num Messages,Message Length,Num Reductions,FLOP,dof0,dof1,dof2,dof3,dof4,dof5,dof6,dof7,e0,e1,e2,e3,e4,e5,e6,e7,%d\n", size));
1109: PetscCall(PetscViewerFlush(viewer));
1110: for (stage = 0; stage < numStages; stage++) {
1111: PetscEventPerfInfo *stage_perf_info;
1112: PetscInt stage_id;
1113: const char *stage_name;
1115: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1116: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1117: stage_perf_info = &zero_info;
1118: if (stage_id >= 0) {
1119: PetscStagePerf *stage_info;
1120: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info));
1121: stage_perf_info = &stage_info->perfInfo;
1122: }
1123: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s,summary,%d,1,%g,%g,%g,%g,%g\n", stage_name, rank, stage_perf_info->time, stage_perf_info->numMessages, stage_perf_info->messageLength, stage_perf_info->numReductions, stage_perf_info->flops));
1124: for (event = 0; event < numEvents; event++) {
1125: PetscEventPerfInfo *eventInfo = &zero_info;
1126: PetscBool is_zero = PETSC_FALSE;
1127: PetscInt event_id;
1128: const char *event_name;
1130: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1131: PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1132: if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1133: PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero));
1134: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1135: if (!is_zero) {
1136: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s,%s,%d,%d,%g,%g,%g,%g,%g", stage_name, event_name, rank, eventInfo->count, eventInfo->time, eventInfo->numMessages, eventInfo->messageLength, eventInfo->numReductions, eventInfo->flops));
1137: if (eventInfo->dof[0] >= 0.) {
1138: for (PetscInt d = 0; d < 8; ++d) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->dof[d]));
1139: for (PetscInt e = 0; e < 8; ++e) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->errors[e]));
1140: }
1141: }
1142: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
1143: }
1144: }
1145: PetscCall(PetscViewerFlush(viewer));
1146: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1147: PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1148: PetscCall(PetscLogGlobalNamesDestroy(&global_events));
1149: PetscFunctionReturn(PETSC_SUCCESS);
1150: }
1152: static PetscErrorCode PetscLogViewWarnSync(PetscViewer viewer)
1153: {
1154: PetscFunctionBegin;
1155: if (!PetscLogSyncOn) PetscFunctionReturn(PETSC_SUCCESS);
1156: PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1157: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n"));
1158: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1159: PetscCall(PetscViewerASCIIPrintf(viewer, " # WARNING!!! #\n"));
1160: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1161: PetscCall(PetscViewerASCIIPrintf(viewer, " # This program was run with logging synchronization. #\n"));
1162: PetscCall(PetscViewerASCIIPrintf(viewer, " # This option provides more meaningful imbalance #\n"));
1163: PetscCall(PetscViewerASCIIPrintf(viewer, " # figures at the expense of slowing things down and #\n"));
1164: PetscCall(PetscViewerASCIIPrintf(viewer, " # providing a distorted view of the overall runtime. #\n"));
1165: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1166: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n\n\n"));
1167: PetscFunctionReturn(PETSC_SUCCESS);
1168: }
1170: static PetscErrorCode PetscLogViewWarnDebugging(PetscViewer viewer)
1171: {
1172: PetscFunctionBegin;
1173: if (PetscDefined(USE_DEBUG)) {
1174: PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1175: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n"));
1176: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1177: PetscCall(PetscViewerASCIIPrintf(viewer, " # WARNING!!! #\n"));
1178: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1179: PetscCall(PetscViewerASCIIPrintf(viewer, " # This code was compiled with a debugging option. #\n"));
1180: PetscCall(PetscViewerASCIIPrintf(viewer, " # To get timing results run ./configure #\n"));
1181: PetscCall(PetscViewerASCIIPrintf(viewer, " # using --with-debugging=no, the performance will #\n"));
1182: PetscCall(PetscViewerASCIIPrintf(viewer, " # be generally two or three times faster. #\n"));
1183: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1184: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n\n\n"));
1185: }
1186: PetscFunctionReturn(PETSC_SUCCESS);
1187: }
1189: static PetscErrorCode PetscLogViewWarnNoGpuAwareMpi(PetscViewer viewer)
1190: {
1191: #if defined(PETSC_HAVE_DEVICE)
1192: PetscMPIInt size;
1193: PetscBool deviceInitialized = PETSC_FALSE;
1195: PetscFunctionBegin;
1196: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
1197: for (int i = PETSC_DEVICE_HOST + 1; i < PETSC_DEVICE_MAX; ++i) {
1198: const PetscDeviceType dtype = PetscDeviceTypeCast(i);
1199: if (PetscDeviceInitialized(dtype)) { /* a non-host device was initialized */
1200: deviceInitialized = PETSC_TRUE;
1201: break;
1202: }
1203: }
1204: /* the last condition says petsc is configured with device but it is a pure CPU run, so don't print misleading warnings */
1205: if (use_gpu_aware_mpi || size == 1 || !deviceInitialized) PetscFunctionReturn(PETSC_SUCCESS);
1206: PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1207: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n"));
1208: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1209: PetscCall(PetscViewerASCIIPrintf(viewer, " # WARNING!!! #\n"));
1210: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1211: PetscCall(PetscViewerASCIIPrintf(viewer, " # This code was compiled with GPU support and you've #\n"));
1212: PetscCall(PetscViewerASCIIPrintf(viewer, " # created PETSc/GPU objects, but you intentionally #\n"));
1213: PetscCall(PetscViewerASCIIPrintf(viewer, " # used -use_gpu_aware_mpi 0, requiring PETSc to copy #\n"));
1214: PetscCall(PetscViewerASCIIPrintf(viewer, " # additional data between the GPU and CPU. To obtain #\n"));
1215: PetscCall(PetscViewerASCIIPrintf(viewer, " # meaningful timing results on multi-rank runs, use #\n"));
1216: PetscCall(PetscViewerASCIIPrintf(viewer, " # GPU-aware MPI instead. #\n"));
1217: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1218: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n\n\n"));
1219: PetscFunctionReturn(PETSC_SUCCESS);
1220: #else
1221: return PETSC_SUCCESS;
1222: #endif
1223: }
1225: static PetscErrorCode PetscLogViewWarnGpuTime(PetscViewer viewer)
1226: {
1227: #if defined(PETSC_HAVE_DEVICE) && !defined(PETSC_HAVE_KOKKOS_WITHOUT_GPU)
1229: PetscFunctionBegin;
1230: if (!PetscLogGpuTimeFlag || petsc_gflops == 0) PetscFunctionReturn(PETSC_SUCCESS);
1231: PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1232: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n"));
1233: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1234: PetscCall(PetscViewerASCIIPrintf(viewer, " # WARNING!!! #\n"));
1235: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1236: PetscCall(PetscViewerASCIIPrintf(viewer, " # This code was run with -log_view_gpu_time #\n"));
1237: PetscCall(PetscViewerASCIIPrintf(viewer, " # This provides accurate timing within the GPU kernels #\n"));
1238: PetscCall(PetscViewerASCIIPrintf(viewer, " # but can slow down the entire computation by a #\n"));
1239: PetscCall(PetscViewerASCIIPrintf(viewer, " # measurable amount. For fastest runs we recommend #\n"));
1240: PetscCall(PetscViewerASCIIPrintf(viewer, " # not using this option. #\n"));
1241: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1242: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n\n\n"));
1243: PetscFunctionReturn(PETSC_SUCCESS);
1244: #else
1245: return PETSC_SUCCESS;
1246: #endif
1247: }
1249: static PetscErrorCode PetscLogHandlerView_Default_Info(PetscLogHandler handler, PetscViewer viewer)
1250: {
1251: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
1252: char arch[128], hostname[128], username[128], pname[PETSC_MAX_PATH_LEN], date[128];
1253: PetscLogDouble locTotalTime, TotalTime, TotalFlops;
1254: PetscLogDouble numMessages, messageLength, avgMessLen, numReductions;
1255: PetscLogDouble stageTime, flops, flopr, mem, mess, messLen, red;
1256: PetscLogDouble fracTime, fracFlops, fracMessages, fracLength, fracReductions, fracMess, fracMessLen, fracRed;
1257: PetscLogDouble fracStageTime, fracStageFlops, fracStageMess, fracStageMessLen, fracStageRed;
1258: PetscLogDouble min, max, tot, ratio, avg, x, y;
1259: PetscLogDouble minf, maxf, totf, ratf, mint, maxt, tott, ratt, ratC, totm, totml, totr, mal, malmax, emalmax;
1260: #if defined(PETSC_HAVE_DEVICE)
1261: PetscLogEvent KSP_Solve, SNES_Solve, TS_Step, TAO_Solve; /* These need to be fixed to be some events registered with certain objects */
1262: PetscLogDouble cct, gct, csz, gsz, gmaxt, gflops, gflopr, fracgflops;
1263: #endif
1264: PetscMPIInt minC, maxC;
1265: PetscMPIInt size, rank;
1266: PetscBool *localStageUsed, *stageUsed;
1267: PetscBool *localStageVisible, *stageVisible;
1268: PetscInt numStages, numEvents;
1269: int stage, oclass;
1270: PetscLogEvent event;
1271: char version[256];
1272: MPI_Comm comm;
1273: #if defined(PETSC_HAVE_DEVICE) && !defined(PETSC_HAVE_KOKKOS_WITHOUT_GPU)
1274: PetscInt64 nas = 0x7FF0000000000002;
1275: #endif
1276: PetscLogGlobalNames global_stages, global_events;
1277: PetscEventPerfInfo zero_info;
1278: PetscLogState state;
1280: PetscFunctionBegin;
1281: PetscCall(PetscLogHandlerGetState(handler, &state));
1282: PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1283: PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
1284: PetscCallMPI(MPI_Comm_size(comm, &size));
1285: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1286: /* Get the total elapsed time */
1287: PetscCall(PetscTime(&locTotalTime));
1288: locTotalTime -= petsc_BaseTime;
1290: PetscCall(PetscViewerASCIIPrintf(viewer, "****************************************************************************************************************************************************************\n"));
1291: PetscCall(PetscViewerASCIIPrintf(viewer, "*** WIDEN YOUR WINDOW TO 160 CHARACTERS. Use 'enscript -r -fCourier9' to print this document ***\n"));
1292: PetscCall(PetscViewerASCIIPrintf(viewer, "****************************************************************************************************************************************************************\n"));
1293: PetscCall(PetscViewerASCIIPrintf(viewer, "\n------------------------------------------------------------------ PETSc Performance Summary: ------------------------------------------------------------------\n\n"));
1294: PetscCall(PetscLogViewWarnSync(viewer));
1295: PetscCall(PetscLogViewWarnDebugging(viewer));
1296: PetscCall(PetscLogViewWarnNoGpuAwareMpi(viewer));
1297: PetscCall(PetscLogViewWarnGpuTime(viewer));
1298: PetscCall(PetscGetArchType(arch, sizeof(arch)));
1299: PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
1300: PetscCall(PetscGetUserName(username, sizeof(username)));
1301: PetscCall(PetscGetProgramName(pname, sizeof(pname)));
1302: PetscCall(PetscGetDate(date, sizeof(date)));
1303: PetscCall(PetscGetVersion(version, sizeof(version)));
1305: #if defined(PETSC_HAVE_CUPM)
1306: const char *cupm = PetscDefined(HAVE_CUDA) ? "CUDA" : "HIP";
1307: if (PetscDeviceCUPMRuntimeArch)
1308: PetscCall(PetscViewerASCIIPrintf(viewer, "%s on a %s named %s with %d process%s and %s architecture %d, by %s on %s\n", pname, arch, hostname, size, size > 1 ? "es" : "", cupm, PetscDeviceCUPMRuntimeArch, username, date));
1309: else
1310: #endif
1311: PetscCall(PetscViewerASCIIPrintf(viewer, "%s on a %s named %s with %d process%s, by %s on %s\n", pname, arch, hostname, size, size > 1 ? "es" : "", username, date));
1313: #if defined(PETSC_HAVE_OPENMP)
1314: PetscCall(PetscViewerASCIIPrintf(viewer, "Using %" PetscInt_FMT " OpenMP threads\n", PetscNumOMPThreads));
1315: #endif
1316: PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s\n", version));
1318: /* Must preserve reduction count before we go on */
1319: red = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
1321: /* Calculate summary information */
1322: PetscCall(PetscViewerASCIIPrintf(viewer, "\n Max Max/Min Avg Total\n"));
1323: /* Time */
1324: PetscCall(MPIU_Allreduce(&locTotalTime, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1325: PetscCall(MPIU_Allreduce(&locTotalTime, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1326: PetscCall(MPIU_Allreduce(&locTotalTime, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1327: avg = tot / ((PetscLogDouble)size);
1328: if (min != 0.0) ratio = max / min;
1329: else ratio = 0.0;
1330: PetscCall(PetscViewerASCIIPrintf(viewer, "Time (sec): %5.3e %7.3f %5.3e\n", max, ratio, avg));
1331: TotalTime = tot;
1332: /* Objects */
1333: {
1334: PetscInt num_objects;
1336: PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
1337: avg = (PetscLogDouble)num_objects;
1338: }
1339: PetscCall(MPIU_Allreduce(&avg, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1340: PetscCall(MPIU_Allreduce(&avg, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1341: PetscCall(MPIU_Allreduce(&avg, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1342: avg = tot / ((PetscLogDouble)size);
1343: if (min != 0.0) ratio = max / min;
1344: else ratio = 0.0;
1345: PetscCall(PetscViewerASCIIPrintf(viewer, "Objects: %5.3e %7.3f %5.3e\n", max, ratio, avg));
1346: /* Flops */
1347: PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1348: PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1349: PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1350: avg = tot / ((PetscLogDouble)size);
1351: if (min != 0.0) ratio = max / min;
1352: else ratio = 0.0;
1353: PetscCall(PetscViewerASCIIPrintf(viewer, "Flops: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot));
1354: TotalFlops = tot;
1355: /* Flops/sec -- Must talk to Barry here */
1356: if (locTotalTime != 0.0) flops = petsc_TotalFlops / locTotalTime;
1357: else flops = 0.0;
1358: PetscCall(MPIU_Allreduce(&flops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1359: PetscCall(MPIU_Allreduce(&flops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1360: PetscCall(MPIU_Allreduce(&flops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1361: avg = tot / ((PetscLogDouble)size);
1362: if (min != 0.0) ratio = max / min;
1363: else ratio = 0.0;
1364: PetscCall(PetscViewerASCIIPrintf(viewer, "Flops/sec: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot));
1365: /* Memory */
1366: PetscCall(PetscMallocGetMaximumUsage(&mem));
1367: if (mem > 0.0) {
1368: PetscCall(MPIU_Allreduce(&mem, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1369: PetscCall(MPIU_Allreduce(&mem, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1370: PetscCall(MPIU_Allreduce(&mem, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1371: avg = tot / ((PetscLogDouble)size);
1372: if (min != 0.0) ratio = max / min;
1373: else ratio = 0.0;
1374: PetscCall(PetscViewerASCIIPrintf(viewer, "Memory (bytes): %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot));
1375: }
1376: /* Messages */
1377: mess = 0.5 * (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct);
1378: PetscCall(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1379: PetscCall(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1380: PetscCall(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1381: avg = tot / ((PetscLogDouble)size);
1382: if (min != 0.0) ratio = max / min;
1383: else ratio = 0.0;
1384: PetscCall(PetscViewerASCIIPrintf(viewer, "MPI Msg Count: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot));
1385: numMessages = tot;
1386: /* Message Lengths */
1387: mess = 0.5 * (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len);
1388: PetscCall(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1389: PetscCall(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1390: PetscCall(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1391: if (numMessages != 0) avg = tot / numMessages;
1392: else avg = 0.0;
1393: if (min != 0.0) ratio = max / min;
1394: else ratio = 0.0;
1395: PetscCall(PetscViewerASCIIPrintf(viewer, "MPI Msg Len (bytes): %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot));
1396: messageLength = tot;
1397: /* Reductions */
1398: PetscCall(MPIU_Allreduce(&red, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1399: PetscCall(MPIU_Allreduce(&red, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1400: PetscCall(MPIU_Allreduce(&red, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1401: if (min != 0.0) ratio = max / min;
1402: else ratio = 0.0;
1403: PetscCall(PetscViewerASCIIPrintf(viewer, "MPI Reductions: %5.3e %7.3f\n", max, ratio));
1404: numReductions = red; /* wrong because uses count from process zero */
1405: PetscCall(PetscViewerASCIIPrintf(viewer, "\nFlop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)\n"));
1406: PetscCall(PetscViewerASCIIPrintf(viewer, " e.g., VecAXPY() for real vectors of length N --> 2N flops\n"));
1407: PetscCall(PetscViewerASCIIPrintf(viewer, " and VecAXPY() for complex vectors of length N --> 8N flops\n"));
1409: PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
1410: PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
1411: PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
1412: PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
1413: PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
1414: PetscCall(PetscMalloc1(numStages, &localStageUsed));
1415: PetscCall(PetscMalloc1(numStages, &stageUsed));
1416: PetscCall(PetscMalloc1(numStages, &localStageVisible));
1417: PetscCall(PetscMalloc1(numStages, &stageVisible));
1418: if (numStages > 0) {
1419: for (stage = 0; stage < numStages; stage++) {
1420: PetscInt stage_id;
1422: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1423: if (stage_id >= 0) {
1424: PetscStagePerf *stage_info;
1426: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info));
1427: localStageUsed[stage] = stage_info->used;
1428: localStageVisible[stage] = stage_info->perfInfo.visible;
1429: } else {
1430: localStageUsed[stage] = PETSC_FALSE;
1431: localStageVisible[stage] = PETSC_TRUE;
1432: }
1433: }
1434: PetscCall(MPIU_Allreduce(localStageUsed, stageUsed, numStages, MPIU_BOOL, MPI_LOR, comm));
1435: PetscCall(MPIU_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm));
1436: for (stage = 0; stage < numStages; stage++) {
1437: if (stageUsed[stage] && stageVisible[stage]) {
1438: PetscCall(PetscViewerASCIIPrintf(viewer, "\nSummary of Stages: ----- Time ------ ----- Flop ------ --- Messages --- -- Message Lengths -- -- Reductions --\n"));
1439: PetscCall(PetscViewerASCIIPrintf(viewer, " Avg %%Total Avg %%Total Count %%Total Avg %%Total Count %%Total\n"));
1440: break;
1441: }
1442: }
1443: for (stage = 0; stage < numStages; stage++) {
1444: PetscInt stage_id;
1445: PetscEventPerfInfo *stage_info;
1446: const char *stage_name;
1448: if (!(stageUsed[stage] && stageVisible[stage])) continue;
1449: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1450: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1451: stage_info = &zero_info;
1452: if (localStageUsed[stage]) {
1453: PetscStagePerf *stage_perf_info;
1455: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_perf_info));
1456: stage_info = &stage_perf_info->perfInfo;
1457: }
1458: PetscCall(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1459: PetscCall(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1460: PetscCall(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1461: PetscCall(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1462: PetscCall(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1463: mess *= 0.5;
1464: messLen *= 0.5;
1465: red /= size;
1466: if (TotalTime != 0.0) fracTime = stageTime / TotalTime;
1467: else fracTime = 0.0;
1468: if (TotalFlops != 0.0) fracFlops = flops / TotalFlops;
1469: else fracFlops = 0.0;
1470: /* Talk to Barry if (stageTime != 0.0) flops = (size*flops)/stageTime; else flops = 0.0; */
1471: if (numMessages != 0.0) fracMessages = mess / numMessages;
1472: else fracMessages = 0.0;
1473: if (mess != 0.0) avgMessLen = messLen / mess;
1474: else avgMessLen = 0.0;
1475: if (messageLength != 0.0) fracLength = messLen / messageLength;
1476: else fracLength = 0.0;
1477: if (numReductions != 0.0) fracReductions = red / numReductions;
1478: else fracReductions = 0.0;
1479: PetscCall(PetscViewerASCIIPrintf(viewer, "%2d: %15s: %6.4e %5.1f%% %6.4e %5.1f%% %5.3e %5.1f%% %5.3e %5.1f%% %5.3e %5.1f%%\n", stage, stage_name, stageTime / size, 100.0 * fracTime, flops, 100.0 * fracFlops, mess, 100.0 * fracMessages, avgMessLen, 100.0 * fracLength, red, 100.0 * fracReductions));
1480: }
1481: }
1483: PetscCall(PetscViewerASCIIPrintf(viewer, "\n------------------------------------------------------------------------------------------------------------------------\n"));
1484: PetscCall(PetscViewerASCIIPrintf(viewer, "See the 'Profiling' chapter of the users' manual for details on interpreting output.\n"));
1485: PetscCall(PetscViewerASCIIPrintf(viewer, "Phase summary info:\n"));
1486: PetscCall(PetscViewerASCIIPrintf(viewer, " Count: number of times phase was executed\n"));
1487: PetscCall(PetscViewerASCIIPrintf(viewer, " Time and Flop: Max - maximum over all processors\n"));
1488: PetscCall(PetscViewerASCIIPrintf(viewer, " Ratio - ratio of maximum to minimum over all processors\n"));
1489: PetscCall(PetscViewerASCIIPrintf(viewer, " Mess: number of messages sent\n"));
1490: PetscCall(PetscViewerASCIIPrintf(viewer, " AvgLen: average message length (bytes)\n"));
1491: PetscCall(PetscViewerASCIIPrintf(viewer, " Reduct: number of global reductions\n"));
1492: PetscCall(PetscViewerASCIIPrintf(viewer, " Global: entire computation\n"));
1493: PetscCall(PetscViewerASCIIPrintf(viewer, " Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().\n"));
1494: PetscCall(PetscViewerASCIIPrintf(viewer, " %%T - percent time in this phase %%F - percent flop in this phase\n"));
1495: PetscCall(PetscViewerASCIIPrintf(viewer, " %%M - percent messages in this phase %%L - percent message lengths in this phase\n"));
1496: PetscCall(PetscViewerASCIIPrintf(viewer, " %%R - percent reductions in this phase\n"));
1497: PetscCall(PetscViewerASCIIPrintf(viewer, " Total Mflop/s: 10e-6 * (sum of flop over all processors)/(max time over all processors)\n"));
1498: if (PetscLogMemory) {
1499: PetscCall(PetscViewerASCIIPrintf(viewer, " Memory usage is summed over all MPI processes, it is given in mega-bytes\n"));
1500: PetscCall(PetscViewerASCIIPrintf(viewer, " Malloc Mbytes: Memory allocated and kept during event (sum over all calls to event). May be negative\n"));
1501: PetscCall(PetscViewerASCIIPrintf(viewer, " EMalloc Mbytes: extra memory allocated during event and then freed (maximum over all calls to events). Never negative\n"));
1502: PetscCall(PetscViewerASCIIPrintf(viewer, " MMalloc Mbytes: Increase in high water mark of allocated memory (sum over all calls to event). Never negative\n"));
1503: PetscCall(PetscViewerASCIIPrintf(viewer, " RMI Mbytes: Increase in resident memory (sum over all calls to event)\n"));
1504: }
1505: #if defined(PETSC_HAVE_DEVICE)
1506: PetscCall(PetscViewerASCIIPrintf(viewer, " GPU Mflop/s: 10e-6 * (sum of flop on GPU over all processors)/(max GPU time over all processors)\n"));
1507: PetscCall(PetscViewerASCIIPrintf(viewer, " CpuToGpu Count: total number of CPU to GPU copies per processor\n"));
1508: PetscCall(PetscViewerASCIIPrintf(viewer, " CpuToGpu Size (Mbytes): 10e-6 * (total size of CPU to GPU copies per processor)\n"));
1509: PetscCall(PetscViewerASCIIPrintf(viewer, " GpuToCpu Count: total number of GPU to CPU copies per processor\n"));
1510: PetscCall(PetscViewerASCIIPrintf(viewer, " GpuToCpu Size (Mbytes): 10e-6 * (total size of GPU to CPU copies per processor)\n"));
1511: PetscCall(PetscViewerASCIIPrintf(viewer, " GPU %%F: percent flops on GPU in this event\n"));
1512: #endif
1513: PetscCall(PetscViewerASCIIPrintf(viewer, "------------------------------------------------------------------------------------------------------------------------\n"));
1515: PetscCall(PetscLogViewWarnDebugging(viewer));
1517: /* Report events */
1518: PetscCall(PetscViewerASCIIPrintf(viewer, "Event Count Time (sec) Flop --- Global --- --- Stage ---- Total"));
1519: if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, " Malloc EMalloc MMalloc RMI"));
1520: #if defined(PETSC_HAVE_DEVICE)
1521: PetscCall(PetscViewerASCIIPrintf(viewer, " GPU - CpuToGpu - - GpuToCpu - GPU"));
1522: #endif
1523: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1524: PetscCall(PetscViewerASCIIPrintf(viewer, " Max Ratio Max Ratio Max Ratio Mess AvgLen Reduct %%T %%F %%M %%L %%R %%T %%F %%M %%L %%R Mflop/s"));
1525: if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, " Mbytes Mbytes Mbytes Mbytes"));
1526: #if defined(PETSC_HAVE_DEVICE)
1527: PetscCall(PetscViewerASCIIPrintf(viewer, " Mflop/s Count Size Count Size %%F"));
1528: #endif
1529: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1530: PetscCall(PetscViewerASCIIPrintf(viewer, "------------------------------------------------------------------------------------------------------------------------"));
1531: if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, "-----------------------------"));
1532: #if defined(PETSC_HAVE_DEVICE)
1533: PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------------"));
1534: #endif
1535: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1537: #if defined(PETSC_HAVE_DEVICE)
1538: /* this indirect way of accessing these values is needed when PETSc is build with multiple libraries since the symbols are not in libpetscsys */
1539: PetscCall(PetscLogStateGetEventFromName(state, "TaoSolve", &TAO_Solve));
1540: PetscCall(PetscLogStateGetEventFromName(state, "TSStep", &TS_Step));
1541: PetscCall(PetscLogStateGetEventFromName(state, "SNESSolve", &SNES_Solve));
1542: PetscCall(PetscLogStateGetEventFromName(state, "KSPSolve", &KSP_Solve));
1543: #endif
1545: for (stage = 0; stage < numStages; stage++) {
1546: PetscInt stage_id;
1547: PetscEventPerfInfo *stage_info;
1548: const char *stage_name;
1550: if (!(stageVisible[stage] && stageUsed[stage])) continue;
1551: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1552: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1553: PetscCall(PetscViewerASCIIPrintf(viewer, "\n--- Event Stage %d: %s\n\n", stage, stage_name));
1554: stage_info = &zero_info;
1555: if (localStageUsed[stage]) {
1556: PetscStagePerf *stage_perf_info;
1558: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_perf_info));
1559: stage_info = &stage_perf_info->perfInfo;
1560: }
1561: PetscCall(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1562: PetscCall(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1563: PetscCall(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1564: PetscCall(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1565: PetscCall(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1566: mess *= 0.5;
1567: messLen *= 0.5;
1568: red /= size;
1570: for (event = 0; event < numEvents; event++) {
1571: PetscInt event_id;
1572: PetscEventPerfInfo *event_info = &zero_info;
1573: PetscBool is_zero = PETSC_FALSE;
1574: const char *event_name;
1576: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1577: PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1578: if (event_id >= 0 && stage_id >= 0) { PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &event_info)); }
1579: PetscCall(PetscMemcmp(event_info, &zero_info, sizeof(zero_info), &is_zero));
1580: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1581: if (!is_zero) {
1582: flopr = event_info->flops;
1583: PetscCall(MPIU_Allreduce(&flopr, &minf, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1584: PetscCall(MPIU_Allreduce(&flopr, &maxf, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1585: PetscCall(MPIU_Allreduce(&event_info->flops, &totf, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1586: PetscCall(MPIU_Allreduce(&event_info->time, &mint, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1587: PetscCall(MPIU_Allreduce(&event_info->time, &maxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1588: PetscCall(MPIU_Allreduce(&event_info->time, &tott, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1589: PetscCall(MPIU_Allreduce(&event_info->numMessages, &totm, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1590: PetscCall(MPIU_Allreduce(&event_info->messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1591: PetscCall(MPIU_Allreduce(&event_info->numReductions, &totr, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1592: PetscCall(MPIU_Allreduce(&event_info->count, &minC, 1, MPI_INT, MPI_MIN, comm));
1593: PetscCall(MPIU_Allreduce(&event_info->count, &maxC, 1, MPI_INT, MPI_MAX, comm));
1594: if (PetscLogMemory) {
1595: PetscCall(MPIU_Allreduce(&event_info->memIncrease, &mem, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1596: PetscCall(MPIU_Allreduce(&event_info->mallocSpace, &mal, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1597: PetscCall(MPIU_Allreduce(&event_info->mallocIncrease, &malmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1598: PetscCall(MPIU_Allreduce(&event_info->mallocIncreaseEvent, &emalmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1599: }
1600: #if defined(PETSC_HAVE_DEVICE)
1601: PetscCall(MPIU_Allreduce(&event_info->CpuToGpuCount, &cct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1602: PetscCall(MPIU_Allreduce(&event_info->GpuToCpuCount, &gct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1603: PetscCall(MPIU_Allreduce(&event_info->CpuToGpuSize, &csz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1604: PetscCall(MPIU_Allreduce(&event_info->GpuToCpuSize, &gsz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1605: PetscCall(MPIU_Allreduce(&event_info->GpuFlops, &gflops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1606: PetscCall(MPIU_Allreduce(&event_info->GpuTime, &gmaxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1607: #endif
1608: if (mint < 0.0) {
1609: PetscCall(PetscViewerASCIIPrintf(viewer, "WARNING!!! Minimum time %g over all processors for %s is negative! This happens\n on some machines whose times cannot handle too rapid calls.!\n artificially changing minimum to zero.\n", mint, event_name));
1610: mint = 0;
1611: }
1612: PetscCheck(minf >= 0.0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Minimum flop %g over all processors for %s is negative! Not possible!", minf, event_name);
1613: #if defined(PETSC_HAVE_DEVICE) && !defined(PETSC_HAVE_KOKKOS_WITHOUT_GPU)
1614: /* Put NaN into the time for all events that may not be time accurately since they may happen asynchronously on the GPU */
1615: if (!PetscLogGpuTimeFlag && petsc_gflops > 0) {
1616: memcpy(&gmaxt, &nas, sizeof(PetscLogDouble));
1617: if (event_id != SNES_Solve && event_id != KSP_Solve && event_id != TS_Step && event_id != TAO_Solve) {
1618: memcpy(&mint, &nas, sizeof(PetscLogDouble));
1619: memcpy(&maxt, &nas, sizeof(PetscLogDouble));
1620: }
1621: }
1622: #endif
1623: totm *= 0.5;
1624: totml *= 0.5;
1625: totr /= size;
1627: if (maxC != 0) {
1628: if (minC != 0) ratC = ((PetscLogDouble)maxC) / minC;
1629: else ratC = 0.0;
1630: if (mint != 0.0) ratt = maxt / mint;
1631: else ratt = 0.0;
1632: if (minf != 0.0) ratf = maxf / minf;
1633: else ratf = 0.0;
1634: if (TotalTime != 0.0) fracTime = tott / TotalTime;
1635: else fracTime = 0.0;
1636: if (TotalFlops != 0.0) fracFlops = totf / TotalFlops;
1637: else fracFlops = 0.0;
1638: if (stageTime != 0.0) fracStageTime = tott / stageTime;
1639: else fracStageTime = 0.0;
1640: if (flops != 0.0) fracStageFlops = totf / flops;
1641: else fracStageFlops = 0.0;
1642: if (numMessages != 0.0) fracMess = totm / numMessages;
1643: else fracMess = 0.0;
1644: if (messageLength != 0.0) fracMessLen = totml / messageLength;
1645: else fracMessLen = 0.0;
1646: if (numReductions != 0.0) fracRed = totr / numReductions;
1647: else fracRed = 0.0;
1648: if (mess != 0.0) fracStageMess = totm / mess;
1649: else fracStageMess = 0.0;
1650: if (messLen != 0.0) fracStageMessLen = totml / messLen;
1651: else fracStageMessLen = 0.0;
1652: if (red != 0.0) fracStageRed = totr / red;
1653: else fracStageRed = 0.0;
1654: if (totm != 0.0) totml /= totm;
1655: else totml = 0.0;
1656: if (maxt != 0.0) flopr = totf / maxt;
1657: else flopr = 0.0;
1658: if (fracStageTime > 1.0 || fracStageFlops > 1.0 || fracStageMess > 1.0 || fracStageMessLen > 1.0 || fracStageRed > 1.0)
1659: PetscCall(PetscViewerASCIIPrintf(viewer, "%-16s %7d %3.1f %5.4e %3.1f %3.2e %3.1f %2.1e %2.1e %2.1e %2.0f %2.0f %2.0f %2.0f %2.0f Multiple stages %5.0f", event_name, maxC, ratC, maxt, ratt, maxf, ratf, totm, totml, totr, 100.0 * fracTime, 100.0 * fracFlops, 100.0 * fracMess, 100.0 * fracMessLen, 100.0 * fracRed, PetscAbs(flopr) / 1.0e6));
1660: else {
1661: if (PetscIsNanReal((PetscReal)maxt)) { // when maxt, ratt, flopr are NaN (i.e., run with GPUs but without -log_view_gpu_time), replace the confusing "nan" with "n/a"
1662: PetscCall(PetscViewerASCIIPrintf(viewer, "%-16s %7d %3.1f n/a n/a %3.2e %3.1f %2.1e %2.1e %2.1e %2.0f %2.0f %2.0f %2.0f %2.0f %3.0f %2.0f %2.0f %2.0f %2.0f n/a", event_name, maxC, ratC, maxf, ratf, totm, totml, totr, 100.0 * fracTime, 100.0 * fracFlops, 100.0 * fracMess, 100.0 * fracMessLen, 100.0 * fracRed, 100.0 * fracStageTime, 100.0 * fracStageFlops, 100.0 * fracStageMess, 100.0 * fracStageMessLen, 100.0 * fracStageRed));
1663: } else {
1664: PetscCall(PetscViewerASCIIPrintf(viewer, "%-16s %7d %3.1f %5.4e %3.1f %3.2e %3.1f %2.1e %2.1e %2.1e %2.0f %2.0f %2.0f %2.0f %2.0f %3.0f %2.0f %2.0f %2.0f %2.0f %5.0f", event_name, maxC, ratC, maxt, ratt, maxf, ratf, totm, totml, totr, 100.0 * fracTime, 100.0 * fracFlops, 100.0 * fracMess, 100.0 * fracMessLen, 100.0 * fracRed, 100.0 * fracStageTime, 100.0 * fracStageFlops, 100.0 * fracStageMess, 100.0 * fracStageMessLen, 100.0 * fracStageRed, PetscAbs(flopr) / 1.0e6));
1665: }
1666: }
1667: if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, " %5.0f %5.0f %5.0f %5.0f", mal / 1.0e6, emalmax / 1.0e6, malmax / 1.0e6, mem / 1.0e6));
1668: #if defined(PETSC_HAVE_DEVICE)
1669: if (totf != 0.0) fracgflops = gflops / totf;
1670: else fracgflops = 0.0;
1671: if (gmaxt != 0.0) gflopr = gflops / gmaxt;
1672: else gflopr = 0.0;
1673: if (PetscIsNanReal((PetscReal)gflopr)) {
1674: PetscCall(PetscViewerASCIIPrintf(viewer, " n/a %4.0f %3.2e %4.0f %3.2e % 2.0f", cct / size, csz / (1.0e6 * size), gct / size, gsz / (1.0e6 * size), 100.0 * fracgflops));
1675: } else {
1676: PetscCall(PetscViewerASCIIPrintf(viewer, " %5.0f %4.0f %3.2e %4.0f %3.2e % 2.0f", PetscAbs(gflopr) / 1.0e6, cct / size, csz / (1.0e6 * size), gct / size, gsz / (1.0e6 * size), 100.0 * fracgflops));
1677: }
1678: #endif
1679: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1680: }
1681: }
1682: }
1683: }
1685: /* Memory usage and object creation */
1686: PetscCall(PetscViewerASCIIPrintf(viewer, "------------------------------------------------------------------------------------------------------------------------"));
1687: if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, "-----------------------------"));
1688: #if defined(PETSC_HAVE_DEVICE)
1689: PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------------"));
1690: #endif
1691: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1692: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1694: /* Right now, only stages on the first processor are reported here, meaning only objects associated with
1695: the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then
1696: stats for stages local to processor sets.
1697: */
1698: /* We should figure out the longest object name here (now 20 characters) */
1699: PetscCall(PetscViewerASCIIPrintf(viewer, "Object Type Creations Destructions. Reports information only for process 0.\n"));
1700: for (stage = 0; stage < numStages; stage++) {
1701: const char *stage_name;
1703: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1704: PetscCall(PetscViewerASCIIPrintf(viewer, "\n--- Event Stage %d: %s\n\n", stage, stage_name));
1705: if (localStageUsed[stage]) {
1706: PetscInt num_classes;
1708: PetscCall(PetscLogStateGetNumClasses(state, &num_classes));
1709: for (oclass = 0; oclass < num_classes; oclass++) {
1710: PetscClassPerf *class_perf_info;
1712: PetscCall(PetscLogHandlerDefaultGetClassPerf(handler, stage, oclass, &class_perf_info));
1713: if ((class_perf_info->creations > 0) || (class_perf_info->destructions > 0)) {
1714: PetscLogClassInfo class_reg_info;
1715: PetscBool flg;
1717: PetscCall(PetscLogStateClassGetInfo(state, oclass, &class_reg_info));
1718: if (stage == 0 && oclass == num_classes - 1) {
1719: PetscCall(PetscStrcmp(class_reg_info.name, "Viewer", &flg));
1720: PetscCheck(flg && class_perf_info->creations == 1 && class_perf_info->destructions == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "The last PetscObject type of the main PetscLogStage should be PetscViewer with a single creation and no destruction");
1721: } else PetscCall(PetscViewerASCIIPrintf(viewer, "%20s %5d %5d\n", class_reg_info.name, class_perf_info->creations, class_perf_info->destructions));
1722: }
1723: }
1724: }
1725: }
1727: PetscCall(PetscFree(localStageUsed));
1728: PetscCall(PetscFree(stageUsed));
1729: PetscCall(PetscFree(localStageVisible));
1730: PetscCall(PetscFree(stageVisible));
1731: PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1732: PetscCall(PetscLogGlobalNamesDestroy(&global_events));
1734: /* Information unrelated to this particular run */
1735: PetscCall(PetscViewerASCIIPrintf(viewer, "========================================================================================================================\n"));
1736: PetscCall(PetscTime(&y));
1737: PetscCall(PetscTime(&x));
1738: PetscCall(PetscTime(&y));
1739: PetscCall(PetscTime(&y));
1740: PetscCall(PetscTime(&y));
1741: PetscCall(PetscTime(&y));
1742: PetscCall(PetscTime(&y));
1743: PetscCall(PetscTime(&y));
1744: PetscCall(PetscTime(&y));
1745: PetscCall(PetscTime(&y));
1746: PetscCall(PetscTime(&y));
1747: PetscCall(PetscTime(&y));
1748: PetscCall(PetscViewerASCIIPrintf(viewer, "Average time to get PetscTime(): %g\n", (y - x) / 10.0));
1749: /* MPI information */
1750: if (size > 1) {
1751: MPI_Status status;
1752: PetscMPIInt tag;
1753: MPI_Comm newcomm;
1755: PetscCallMPI(MPI_Barrier(comm));
1756: PetscCall(PetscTime(&x));
1757: PetscCallMPI(MPI_Barrier(comm));
1758: PetscCallMPI(MPI_Barrier(comm));
1759: PetscCallMPI(MPI_Barrier(comm));
1760: PetscCallMPI(MPI_Barrier(comm));
1761: PetscCallMPI(MPI_Barrier(comm));
1762: PetscCall(PetscTime(&y));
1763: PetscCall(PetscViewerASCIIPrintf(viewer, "Average time for MPI_Barrier(): %g\n", (y - x) / 5.0));
1764: PetscCall(PetscCommDuplicate(comm, &newcomm, &tag));
1765: PetscCallMPI(MPI_Barrier(comm));
1766: if (rank) {
1767: PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, rank - 1, tag, newcomm, &status));
1768: PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, (rank + 1) % size, tag, newcomm));
1769: } else {
1770: PetscCall(PetscTime(&x));
1771: PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, 1, tag, newcomm));
1772: PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, size - 1, tag, newcomm, &status));
1773: PetscCall(PetscTime(&y));
1774: PetscCall(PetscViewerASCIIPrintf(viewer, "Average time for zero size MPI_Send(): %g\n", (y - x) / size));
1775: }
1776: PetscCall(PetscCommDestroy(&newcomm));
1777: }
1778: PetscCall(PetscOptionsView(NULL, viewer));
1780: /* Machine and compile information */
1781: if (PetscDefined(USE_FORTRAN_KERNELS)) {
1782: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with FORTRAN kernels\n"));
1783: } else {
1784: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled without FORTRAN kernels\n"));
1785: }
1786: if (PetscDefined(USE_64BIT_INDICES)) {
1787: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 64-bit PetscInt\n"));
1788: } else if (PetscDefined(USE___FLOAT128)) {
1789: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 32-bit PetscInt\n"));
1790: }
1791: if (PetscDefined(USE_REAL_SINGLE)) {
1792: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with single precision PetscScalar and PetscReal\n"));
1793: } else if (PetscDefined(USE___FLOAT128)) {
1794: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 128 bit precision PetscScalar and PetscReal\n"));
1795: }
1796: if (PetscDefined(USE_REAL_MAT_SINGLE)) {
1797: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with single precision matrices\n"));
1798: } else {
1799: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with full precision matrices (default)\n"));
1800: }
1801: PetscCall(PetscViewerASCIIPrintf(viewer, "sizeof(short) %d sizeof(int) %d sizeof(long) %d sizeof(void*) %d sizeof(PetscScalar) %d sizeof(PetscInt) %d\n", (int)sizeof(short), (int)sizeof(int), (int)sizeof(long), (int)sizeof(void *), (int)sizeof(PetscScalar), (int)sizeof(PetscInt)));
1803: PetscCall(PetscViewerASCIIPrintf(viewer, "Configure options: %s", petscconfigureoptions));
1804: PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petscmachineinfo));
1805: PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsccompilerinfo));
1806: PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsccompilerflagsinfo));
1807: PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsclinkerinfo));
1809: /* Cleanup */
1810: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1811: PetscCall(PetscLogViewWarnNoGpuAwareMpi(viewer));
1812: PetscCall(PetscLogViewWarnDebugging(viewer));
1813: PetscCall(PetscFPTrapPop());
1814: PetscFunctionReturn(PETSC_SUCCESS);
1815: }
1817: static PetscErrorCode PetscLogHandlerView_Default(PetscLogHandler handler, PetscViewer viewer)
1818: {
1819: PetscViewerFormat format;
1821: PetscFunctionBegin;
1822: PetscCall(PetscViewerGetFormat(viewer, &format));
1823: if (format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_ASCII_INFO) {
1824: PetscCall(PetscLogHandlerView_Default_Info(handler, viewer));
1825: } else if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1826: PetscCall(PetscLogHandlerView_Default_Detailed(handler, viewer));
1827: } else if (format == PETSC_VIEWER_ASCII_CSV) {
1828: PetscCall(PetscLogHandlerView_Default_CSV(handler, viewer));
1829: }
1830: PetscFunctionReturn(PETSC_SUCCESS);
1831: }
1833: /*MC
1834: PETSCLOGHANDLERDEFAULT - PETSCLOGHANDLERDEFAULT = "default" - A `PetscLogHandler` that collects data for PETSc
1835: default profiling log viewers (`PetscLogView()` and `PetscLogDump()`). A log handler of this type is
1836: created and started (`PetscLogHandlerStart()`) by `PetscLogDefaultBegin()`.
1838: Options Database Keys:
1839: + -log_include_actions - include a growing list of actions (event beginnings and endings, object creations and destructions) in `PetscLogDump()` (`PetscLogActions()`).
1840: - -log_include_objects - include a growing list of object creations and destructions in `PetscLogDump()` (`PetscLogObjects()`).
1842: Level: developer
1844: .seealso: [](ch_profiling), `PetscLogHandler`
1845: M*/
1847: PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Default(PetscLogHandler handler)
1848: {
1849: PetscFunctionBegin;
1850: PetscCall(PetscLogHandlerContextCreate_Default((PetscLogHandler_Default *)&handler->data));
1851: handler->ops->destroy = PetscLogHandlerDestroy_Default;
1852: handler->ops->eventbegin = PetscLogHandlerEventBegin_Default;
1853: handler->ops->eventend = PetscLogHandlerEventEnd_Default;
1854: handler->ops->eventsync = PetscLogHandlerEventSync_Default;
1855: handler->ops->objectcreate = PetscLogHandlerObjectCreate_Default;
1856: handler->ops->objectdestroy = PetscLogHandlerObjectDestroy_Default;
1857: handler->ops->stagepush = PetscLogHandlerStagePush_Default;
1858: handler->ops->stagepop = PetscLogHandlerStagePop_Default;
1859: handler->ops->view = PetscLogHandlerView_Default;
1860: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetEventPerfInfo_C", PetscLogHandlerGetEventPerfInfo_Default));
1861: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetStagePerfInfo_C", PetscLogHandlerGetStagePerfInfo_Default));
1862: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerSetLogActions_C", PetscLogHandlerSetLogActions_Default));
1863: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerSetLogObjects_C", PetscLogHandlerSetLogObjects_Default));
1864: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerLogObjectState_C", PetscLogHandlerLogObjectState_Default));
1865: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetNumObjects_C", PetscLogHandlerGetNumObjects_Default));
1866: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventDeactivatePush_C", PetscLogHandlerEventDeactivatePush_Default));
1867: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventDeactivatePop_C", PetscLogHandlerEventDeactivatePop_Default));
1868: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventsPause_C", PetscLogHandlerEventsPause_Default));
1869: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventsResume_C", PetscLogHandlerEventsResume_Default));
1870: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerDump_C", PetscLogHandlerDump_Default));
1871: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerStageSetVisible_C", PetscLogHandlerStageSetVisible_Default));
1872: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerStageGetVisible_C", PetscLogHandlerStageGetVisible_Default));
1873: PetscFunctionReturn(PETSC_SUCCESS);
1874: }