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));
277: PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_include_actions", &def->petsc_logActions, NULL));
278: PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_include_objects", &def->petsc_logObjects, NULL));
279: PetscCall(PetscOptionsGetBool(NULL, NULL, "-log_handler_default_use_threadsafe_events", &def->use_threadsafe, NULL));
280: if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) { PetscCall(PetscHMapEventCreate(&def->eventInfoMap_th)); }
281: PetscFunctionReturn(PETSC_SUCCESS);
282: }
284: static PetscErrorCode PetscLogHandlerDestroy_Default(PetscLogHandler h)
285: {
286: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
288: PetscFunctionBegin;
289: PetscCall(PetscLogStageInfoArrayDestroy(&def->stages));
290: PetscCall(PetscLogActionArrayDestroy(&def->petsc_actions));
291: PetscCall(PetscLogObjectArrayDestroy(&def->petsc_objects));
292: if (def->eventInfoMap_th) {
293: PetscEventPerfInfo **array;
294: PetscInt n, off = 0;
296: PetscCall(PetscHMapEventGetSize(def->eventInfoMap_th, &n));
297: PetscCall(PetscMalloc1(n, &array));
298: PetscCall(PetscHMapEventGetVals(def->eventInfoMap_th, &off, array));
299: for (PetscInt i = 0; i < n; i++) PetscCall(PetscFree(array[i]));
300: PetscCall(PetscFree(array));
301: PetscCall(PetscHMapEventDestroy(&def->eventInfoMap_th));
302: }
303: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerGetEventPerfInfo_C", NULL));
304: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerGetStagePerfInfo_C", NULL));
305: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerSetLogActions_C", NULL));
306: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerSetLogObjects_C", NULL));
307: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerLogObjectState_C", NULL));
308: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerGetNumObjects_C", NULL));
309: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventDeactivatePush_C", NULL));
310: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventDeactivatePop_C", NULL));
311: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventsPause_C", NULL));
312: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerEventsResume_C", NULL));
313: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerDump_C", NULL));
314: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerStageSetVisible_C", NULL));
315: PetscCall(PetscObjectComposeFunction((PetscObject)h, "PetscLogHandlerStageGetVisible_C", NULL));
316: PetscCall(PetscFree(def));
317: PetscFunctionReturn(PETSC_SUCCESS);
318: }
320: static PetscErrorCode PetscLogHandlerDefaultGetStageInfo(PetscLogHandler handler, PetscLogStage stage, PetscStagePerf **stage_info_p)
321: {
322: PetscStagePerf *stage_info = NULL;
323: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
325: PetscFunctionBegin;
326: PetscCall(PetscLogStageInfoArrayResize(def->stages, stage + 1));
327: PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
328: *stage_info_p = stage_info;
329: PetscFunctionReturn(PETSC_SUCCESS);
330: }
332: static PetscErrorCode PetscLogHandlerGetEventPerfInfo_Default(PetscLogHandler handler, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **event_info_p)
333: {
334: PetscEventPerfInfo *event_info = NULL;
335: PetscStagePerf *stage_info = NULL;
336: PetscLogEventPerfArray event_log;
338: PetscFunctionBegin;
339: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(handler->state, &stage));
340: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info));
341: event_log = stage_info->eventLog;
342: PetscCall(PetscLogEventPerfArrayResize(event_log, event + 1));
343: PetscCall(PetscLogEventPerfArrayGetRef(event_log, event, &event_info));
344: event_info->id = event;
345: *event_info_p = event_info;
346: PetscFunctionReturn(PETSC_SUCCESS);
347: }
349: static PetscErrorCode PetscLogHandlerGetStagePerfInfo_Default(PetscLogHandler handler, PetscLogStage stage, PetscEventPerfInfo **stage_info_p)
350: {
351: PetscStagePerf *stage_perf_info = NULL;
353: PetscFunctionBegin;
354: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(handler->state, &stage));
355: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_perf_info));
356: *stage_info_p = &stage_perf_info->perfInfo;
357: PetscFunctionReturn(PETSC_SUCCESS);
358: }
360: static PetscErrorCode PetscLogHandlerDefaultGetClassPerf(PetscLogHandler handler, PetscLogStage stage, PetscLogClass clss, PetscClassPerf **class_info)
361: {
362: PetscLogClassPerfArray class_log;
363: PetscStagePerf *stage_info;
365: PetscFunctionBegin;
366: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info));
367: class_log = stage_info->classLog;
368: PetscCall(PetscLogClassPerfArrayResize(class_log, clss + 1));
369: PetscCall(PetscLogClassPerfArrayGetRef(class_log, clss, class_info));
370: PetscFunctionReturn(PETSC_SUCCESS);
371: }
373: static PetscErrorCode PetscLogHandlerObjectCreate_Default(PetscLogHandler h, PetscObject obj)
374: {
375: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
376: PetscLogState state;
377: PetscLogStage stage;
378: PetscClassPerf *classInfo;
379: int oclass = 0;
381: PetscFunctionBegin;
382: PetscCall(PetscLogHandlerGetState(h, &state));
383: PetscCall(PetscSpinlockLock(&def->lock));
384: /* Record stage info */
385: PetscCall(PetscLogStateGetCurrentStage(state, &stage));
386: PetscCall(PetscLogStateGetClassFromClassId(state, obj->classid, &oclass));
387: PetscCall(PetscLogHandlerDefaultGetClassPerf(h, stage, oclass, &classInfo));
388: classInfo->creations++;
389: /* Record the creation action */
390: if (def->petsc_logActions) {
391: Action new_action;
393: PetscCall(PetscTime(&new_action.time));
394: new_action.time -= petsc_BaseTime;
395: new_action.action = PETSC_LOG_ACTION_CREATE;
396: new_action.event = -1;
397: new_action.classid = obj->classid;
398: new_action.id1 = obj->id;
399: new_action.id2 = -1;
400: new_action.id3 = -1;
401: new_action.flops = petsc_TotalFlops;
402: PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
403: PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
404: PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
405: }
406: /* We don't just use obj->id to count all objects that are created
407: because PetscLogHandlers are objects and PetscLogObjectDestroy() will not
408: be called for them: the number of objects created and destroyed as counted
409: here and below would have an imbalance */
410: def->petsc_numObjectsCreated++;
411: /* Record the object */
412: if (def->petsc_logObjects) {
413: Object new_object;
415: new_object.parent = -1;
416: new_object.obj = obj;
417: new_object.mem = 0;
419: PetscCall(PetscMemzero(new_object.name, sizeof(new_object.name)));
420: PetscCall(PetscMemzero(new_object.info, sizeof(new_object.info)));
421: PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1");
422: PetscCall(PetscLogObjectArrayResize(def->petsc_objects, obj->id));
423: PetscCall(PetscLogObjectArraySet(def->petsc_objects, obj->id - 1, new_object));
424: }
425: PetscCall(PetscSpinlockUnlock(&def->lock));
426: PetscFunctionReturn(PETSC_SUCCESS);
427: }
429: static PetscErrorCode PetscLogHandlerObjectDestroy_Default(PetscLogHandler h, PetscObject obj)
430: {
431: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
432: PetscLogState state;
433: PetscLogStage stage;
434: PetscClassPerf *classInfo;
435: int oclass = 0;
437: PetscFunctionBegin;
438: PetscCall(PetscLogHandlerGetState(h, &state));
439: /* Record stage info */
440: PetscCall(PetscSpinlockLock(&def->lock));
441: PetscCall(PetscLogStateGetCurrentStage(state, &stage));
442: if (stage >= 0) {
443: /* stage < 0 can happen if the log summary is output before some things are destroyed */
444: PetscCall(PetscLogStateGetClassFromClassId(state, obj->classid, &oclass));
445: PetscCall(PetscLogHandlerDefaultGetClassPerf(h, stage, oclass, &classInfo));
446: classInfo->destructions++;
447: }
448: /* Cannot Credit all ancestors with your memory because they may have already been destroyed*/
449: def->petsc_numObjectsDestroyed++;
450: /* Dynamically enlarge logging structures */
451: /* Record the destruction action */
452: if (def->petsc_logActions) {
453: Action new_action;
455: PetscCall(PetscTime(&new_action.time));
456: new_action.time -= petsc_BaseTime;
457: new_action.event = -1;
458: new_action.action = PETSC_LOG_ACTION_DESTROY;
459: new_action.classid = obj->classid;
460: new_action.id1 = obj->id;
461: new_action.id2 = -1;
462: new_action.id3 = -1;
463: new_action.flops = petsc_TotalFlops;
464: PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
465: PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
466: PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
467: }
468: if (def->petsc_logObjects) {
469: Object *obj_entry = NULL;
471: PetscAssert(obj->id >= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Object ids from PetscObjectNewId_Internal() start at 1");
472: PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, obj->id - 1, &obj_entry));
473: if (obj->name) PetscCall(PetscStrncpy(obj_entry->name, obj->name, 64));
474: obj_entry->obj = NULL;
475: }
476: PetscCall(PetscSpinlockUnlock(&def->lock));
477: PetscFunctionReturn(PETSC_SUCCESS);
478: }
480: static PetscErrorCode PetscLogHandlerEventSync_Default(PetscLogHandler h, PetscLogEvent event, MPI_Comm comm)
481: {
482: PetscLogState state;
483: PetscLogEventInfo event_info;
484: PetscEventPerfInfo *event_perf_info;
485: int stage;
486: PetscLogDouble time = 0.0;
488: PetscFunctionBegin;
489: if (!PetscLogSyncOn || comm == MPI_COMM_NULL) PetscFunctionReturn(PETSC_SUCCESS);
490: PetscCall(PetscLogHandlerGetState(h, &state));
491: PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
492: if (!event_info.collective) PetscFunctionReturn(PETSC_SUCCESS);
493: PetscCall(PetscLogStateGetCurrentStage(state, &stage));
494: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
495: if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
497: PetscCall(PetscTimeSubtract(&time));
498: PetscCallMPI(MPI_Barrier(comm));
499: PetscCall(PetscTimeAdd(&time));
500: event_perf_info->syncTime += time;
501: PetscFunctionReturn(PETSC_SUCCESS);
502: }
504: static PetscErrorCode PetscLogGetStageEventPerfInfo_threaded(PetscLogHandler_Default def, PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo **eventInfo)
505: {
506: PetscEventPerfInfo *leventInfo = NULL;
507: PetscHashIJKKey key;
509: PetscFunctionBegin;
510: #if PetscDefined(HAVE_THREADSAFETY)
511: key.i = PetscLogGetTid();
512: #else
513: key.i = 0;
514: #endif
515: key.j = stage;
516: key.k = event;
517: PetscCall(PetscSpinlockLock(&def->lock));
518: PetscCall(PetscHMapEventGet(def->eventInfoMap_th, key, &leventInfo));
519: if (!leventInfo) {
520: PetscCall(PetscNew(&leventInfo));
521: leventInfo->id = event;
522: PetscCall(PetscHMapEventSet(def->eventInfoMap_th, key, leventInfo));
523: }
524: PetscCall(PetscSpinlockUnlock(&def->lock));
525: *eventInfo = leventInfo;
526: PetscFunctionReturn(PETSC_SUCCESS);
527: }
529: static PetscErrorCode PetscLogHandlerEventBegin_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
530: {
531: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
532: PetscEventPerfInfo *event_perf_info = NULL;
533: PetscLogEventInfo event_info;
534: PetscLogDouble time;
535: PetscLogState state;
536: PetscLogStage stage;
538: PetscFunctionBegin;
539: PetscCall(PetscLogHandlerGetState(h, &state));
540: if (PetscDefined(USE_DEBUG)) {
541: PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
546: if (event_info.collective && o1) {
547: PetscInt64 b1[2], b2[2];
549: b1[0] = -o1->cidx;
550: b1[1] = o1->cidx;
551: PetscCall(MPIU_Allreduce(b1, b2, 2, MPIU_INT64, MPI_MAX, PetscObjectComm(o1)));
552: PetscCheck(-b2[0] == b2[1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Collective event %s not called collectively %" PetscInt64_FMT " != %" PetscInt64_FMT, event_info.name, b1[1], b2[1]);
553: }
554: }
555: /* Synchronization */
556: PetscCall(PetscLogHandlerEventSync_Default(h, event, PetscObjectComm(o1)));
557: PetscCall(PetscLogStateGetCurrentStage(state, &stage));
558: if (def->pause_depth > 0) stage = 0; // in pause-mode, all events run on the main stage
559: if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
560: PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info));
561: if (event_perf_info->depth == 0) { PetscCall(PetscEventPerfInfoInit(event_perf_info)); }
562: } else {
563: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
564: }
565: PetscCheck(event_perf_info->depth >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to begin a paused event, this is not allowed");
566: event_perf_info->depth++;
567: /* Check for double counting */
568: if (event_perf_info->depth > 1) PetscFunctionReturn(PETSC_SUCCESS);
569: PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
570: /* Log the performance info */
571: event_perf_info->count++;
572: PetscCall(PetscTime(&time));
573: PetscCall(PetscEventPerfInfoTic(event_perf_info, time, PetscLogMemory, (int)event));
574: if (def->petsc_logActions) {
575: PetscLogDouble curTime;
576: Action new_action;
578: PetscCall(PetscTime(&curTime));
579: new_action.time = curTime - petsc_BaseTime;
580: new_action.action = PETSC_LOG_ACTION_BEGIN;
581: new_action.event = event;
582: new_action.classid = event_info.classid;
583: new_action.id1 = o1 ? o1->id : -1;
584: new_action.id2 = o2 ? o2->id : -1;
585: new_action.id3 = o3 ? o3->id : -1;
586: new_action.flops = petsc_TotalFlops;
587: PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
588: PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
589: PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
590: }
591: PetscFunctionReturn(PETSC_SUCCESS);
592: }
594: static PetscErrorCode PetscLogHandlerEventEnd_Default(PetscLogHandler h, PetscLogEvent event, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
595: {
596: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
597: PetscEventPerfInfo *event_perf_info = NULL;
598: PetscLogDouble time;
599: PetscLogState state;
600: int stage;
601: PetscLogEventInfo event_info;
603: PetscFunctionBegin;
604: PetscCall(PetscLogHandlerGetState(h, &state));
605: if (PetscDefined(USE_DEBUG)) {
606: PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
611: if (event_info.collective && o1) {
612: PetscInt64 b1[2], b2[2];
614: b1[0] = -o1->cidx;
615: b1[1] = o1->cidx;
616: PetscCall(MPIU_Allreduce(b1, b2, 2, MPIU_INT64, MPI_MAX, PetscObjectComm(o1)));
617: PetscCheck(-b2[0] == b2[1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Collective event %s not called collectively %" PetscInt64_FMT " != %" PetscInt64_FMT, event_info.name, b1[1], b2[1]);
618: }
619: }
620: if (def->petsc_logActions) {
621: PetscLogDouble curTime;
622: Action new_action;
624: PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
625: PetscCall(PetscTime(&curTime));
626: new_action.time = curTime - petsc_BaseTime;
627: new_action.action = PETSC_LOG_ACTION_END;
628: new_action.event = event;
629: new_action.classid = event_info.classid;
630: new_action.id1 = o1 ? o1->id : -1;
631: new_action.id2 = o2 ? o2->id : -2;
632: new_action.id3 = o3 ? o3->id : -3;
633: new_action.flops = petsc_TotalFlops;
634: PetscCall(PetscMallocGetCurrentUsage(&new_action.mem));
635: PetscCall(PetscMallocGetMaximumUsage(&new_action.maxmem));
636: PetscCall(PetscLogActionArrayPush(def->petsc_actions, new_action));
637: }
638: PetscCall(PetscLogStateGetCurrentStage(state, &stage));
639: if (def->pause_depth > 0) stage = 0; // all events run on the main stage in pause-mode
640: if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
641: PetscCall(PetscLogGetStageEventPerfInfo_threaded(def, stage, event, &event_perf_info));
642: } else {
643: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
644: }
645: PetscCheck(event_perf_info->depth > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Trying to end paused event, not allowed");
646: event_perf_info->depth--;
647: /* Check for double counting */
648: if (event_perf_info->depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
649: else PetscCheck(event_perf_info->depth == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Logging event had unbalanced begin/end pairs");
651: /* Log performance info */
652: PetscCall(PetscTime(&time));
653: PetscCall(PetscEventPerfInfoToc(event_perf_info, time, PetscLogMemory, (int)event));
654: if (PetscDefined(HAVE_THREADSAFETY) || def->use_threadsafe) {
655: PetscEventPerfInfo *event_perf_info_global;
656: PetscCall(PetscSpinlockLock(&def->lock));
657: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info_global));
658: PetscCall(PetscEventPerfInfoAdd_Internal(event_perf_info, event_perf_info_global));
659: PetscCall(PetscSpinlockUnlock(&def->lock));
660: }
661: PetscFunctionReturn(PETSC_SUCCESS);
662: }
664: static PetscErrorCode PetscLogHandlerEventDeactivatePush_Default(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event)
665: {
666: PetscEventPerfInfo *event_perf_info;
668: PetscFunctionBegin;
669: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
670: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
671: event_perf_info->depth++;
672: PetscFunctionReturn(PETSC_SUCCESS);
673: }
675: static PetscErrorCode PetscLogHandlerEventDeactivatePop_Default(PetscLogHandler h, PetscLogStage stage, PetscLogEvent event)
676: {
677: PetscEventPerfInfo *event_perf_info;
679: PetscFunctionBegin;
680: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
681: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(h, stage, event, &event_perf_info));
682: event_perf_info->depth--;
683: PetscFunctionReturn(PETSC_SUCCESS);
684: }
686: static PetscErrorCode PetscLogHandlerEventsPause_Default(PetscLogHandler h)
687: {
688: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
689: PetscLogDouble time;
690: PetscInt num_stages;
692: PetscFunctionBegin;
693: if (def->pause_depth++ > 0) PetscFunctionReturn(PETSC_SUCCESS);
694: PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL));
695: PetscCall(PetscTime(&time));
696: for (PetscInt stage = 0; stage < num_stages; stage++) {
697: PetscStagePerf *stage_info = NULL;
698: PetscInt num_events;
700: PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
701: PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL));
702: for (PetscInt event = 0; event < num_events; event++) {
703: PetscEventPerfInfo *event_info = NULL;
704: PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info));
705: if (event_info->depth > 0) {
706: event_info->depth *= -1;
707: PetscCall(PetscEventPerfInfoPause(event_info, time, PetscLogMemory, event));
708: }
709: }
710: if (stage > 0 && stage_info->perfInfo.depth > 0) {
711: stage_info->perfInfo.depth *= -1;
712: PetscCall(PetscEventPerfInfoPause(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2)));
713: }
714: }
715: PetscFunctionReturn(PETSC_SUCCESS);
716: }
718: static PetscErrorCode PetscLogHandlerEventsResume_Default(PetscLogHandler h)
719: {
720: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
721: PetscLogDouble time;
722: PetscInt num_stages;
724: PetscFunctionBegin;
725: if (--def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
726: PetscCall(PetscLogStageInfoArrayGetSize(def->stages, &num_stages, NULL));
727: PetscCall(PetscTime(&time));
728: for (PetscInt stage = 0; stage < num_stages; stage++) {
729: PetscStagePerf *stage_info = NULL;
730: PetscInt num_events;
732: PetscCall(PetscLogStageInfoArrayGetRef(def->stages, stage, &stage_info));
733: PetscCall(PetscLogEventPerfArrayGetSize(stage_info->eventLog, &num_events, NULL));
734: for (PetscInt event = 0; event < num_events; event++) {
735: PetscEventPerfInfo *event_info = NULL;
736: PetscCall(PetscLogEventPerfArrayGetRef(stage_info->eventLog, event, &event_info));
737: if (event_info->depth < 0) {
738: event_info->depth *= -1;
739: PetscCall(PetscEventPerfInfoResume(event_info, time, PetscLogMemory, event));
740: }
741: }
742: if (stage > 0 && stage_info->perfInfo.depth < 0) {
743: stage_info->perfInfo.depth *= -1;
744: PetscCall(PetscEventPerfInfoResume(&stage_info->perfInfo, time, PetscLogMemory, -(stage + 2)));
745: }
746: }
747: PetscFunctionReturn(PETSC_SUCCESS);
748: }
750: static PetscErrorCode PetscLogHandlerStagePush_Default(PetscLogHandler h, PetscLogStage new_stage)
751: {
752: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
753: PetscLogDouble time;
754: PetscLogState state;
755: PetscLogStage current_stage;
756: PetscStagePerf *new_stage_info;
758: PetscFunctionBegin;
759: if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
760: PetscCall(PetscLogHandlerGetState(h, &state));
761: current_stage = state->current_stage;
762: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, new_stage, &new_stage_info));
763: PetscCall(PetscTime(&time));
765: /* Record flops/time of previous stage */
766: if (current_stage >= 0) {
767: if (PetscBTLookup(state->active, current_stage)) {
768: PetscStagePerf *current_stage_info;
769: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, ¤t_stage_info));
770: PetscCall(PetscEventPerfInfoToc(¤t_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2)));
771: }
772: }
773: new_stage_info->used = PETSC_TRUE;
774: new_stage_info->perfInfo.count++;
775: new_stage_info->perfInfo.depth++;
776: /* Subtract current quantities so that we obtain the difference when we pop */
777: if (PetscBTLookup(state->active, new_stage)) PetscCall(PetscEventPerfInfoTic(&new_stage_info->perfInfo, time, PetscLogMemory, (int)-(new_stage + 2)));
778: PetscFunctionReturn(PETSC_SUCCESS);
779: }
781: static PetscErrorCode PetscLogHandlerStagePop_Default(PetscLogHandler h, PetscLogStage old_stage)
782: {
783: PetscLogHandler_Default def = (PetscLogHandler_Default)h->data;
784: PetscLogStage current_stage;
785: PetscStagePerf *old_stage_info;
786: PetscLogState state;
787: PetscLogDouble time;
789: PetscFunctionBegin;
790: if (def->pause_depth > 0) PetscFunctionReturn(PETSC_SUCCESS);
791: PetscCall(PetscLogHandlerGetState(h, &state));
792: current_stage = state->current_stage;
793: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, old_stage, &old_stage_info));
794: PetscCall(PetscTime(&time));
795: old_stage_info->perfInfo.depth--;
796: if (PetscBTLookup(state->active, old_stage)) { PetscCall(PetscEventPerfInfoToc(&old_stage_info->perfInfo, time, PetscLogMemory, (int)-(old_stage + 2))); }
797: if (current_stage >= 0) {
798: if (PetscBTLookup(state->active, current_stage)) {
799: PetscStagePerf *current_stage_info;
800: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, current_stage, ¤t_stage_info));
801: PetscCall(PetscEventPerfInfoTic(¤t_stage_info->perfInfo, time, PetscLogMemory, (int)-(current_stage + 2)));
802: }
803: }
804: PetscFunctionReturn(PETSC_SUCCESS);
805: }
807: static PetscErrorCode PetscLogHandlerStageSetVisible_Default(PetscLogHandler h, PetscLogStage stage, PetscBool is_visible)
808: {
809: PetscStagePerf *stage_info;
811: PetscFunctionBegin;
812: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
813: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info));
814: stage_info->perfInfo.visible = is_visible;
815: PetscFunctionReturn(PETSC_SUCCESS);
816: }
818: static PetscErrorCode PetscLogHandlerStageGetVisible_Default(PetscLogHandler h, PetscLogStage stage, PetscBool *is_visible)
819: {
820: PetscStagePerf *stage_info;
822: PetscFunctionBegin;
823: if (stage < 0) PetscCall(PetscLogStateGetCurrentStage(h->state, &stage));
824: PetscCall(PetscLogHandlerDefaultGetStageInfo(h, stage, &stage_info));
825: *is_visible = stage_info->perfInfo.visible;
826: PetscFunctionReturn(PETSC_SUCCESS);
827: }
829: static PetscErrorCode PetscLogHandlerSetLogActions_Default(PetscLogHandler handler, PetscBool flag)
830: {
831: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
833: PetscFunctionBegin;
834: def->petsc_logActions = flag;
835: PetscFunctionReturn(PETSC_SUCCESS);
836: }
838: static PetscErrorCode PetscLogHandlerSetLogObjects_Default(PetscLogHandler handler, PetscBool flag)
839: {
840: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
842: PetscFunctionBegin;
843: def->petsc_logObjects = flag;
844: PetscFunctionReturn(PETSC_SUCCESS);
845: }
847: static PetscErrorCode PetscLogHandlerLogObjectState_Default(PetscLogHandler handler, PetscObject obj, const char format[], va_list Argp)
848: {
849: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
850: size_t fullLength;
852: PetscFunctionBegin;
853: if (def->petsc_logObjects) {
854: Object *obj_entry = NULL;
856: PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, obj->id - 1, &obj_entry));
857: PetscCall(PetscVSNPrintf(obj_entry->info, 64, format, &fullLength, Argp));
858: }
859: PetscFunctionReturn(PETSC_SUCCESS);
860: }
862: static PetscErrorCode PetscLogHandlerGetNumObjects_Default(PetscLogHandler handler, PetscInt *num_objects)
863: {
864: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
866: PetscFunctionBegin;
867: PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, num_objects, NULL));
868: PetscFunctionReturn(PETSC_SUCCESS);
869: }
871: static PetscErrorCode PetscLogHandlerDump_Default(PetscLogHandler handler, const char sname[])
872: {
873: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
874: FILE *fd;
875: char file[PETSC_MAX_PATH_LEN], fname[PETSC_MAX_PATH_LEN];
876: PetscLogDouble flops, _TotalTime;
877: PetscMPIInt rank;
878: int curStage;
879: PetscLogState state;
880: PetscInt num_events;
881: PetscLogEvent event;
883: PetscFunctionBegin;
884: /* Calculate the total elapsed time */
885: PetscCall(PetscTime(&_TotalTime));
886: _TotalTime -= petsc_BaseTime;
887: /* Open log file */
888: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)handler), &rank));
889: PetscCall(PetscSNPrintf(file, PETSC_STATIC_ARRAY_LENGTH(file), "%s.%d", sname && sname[0] ? sname : "Log", rank));
890: PetscCall(PetscFixFilename(file, fname));
891: PetscCall(PetscFOpen(PETSC_COMM_SELF, fname, "w", &fd));
892: PetscCheck(!(rank == 0) || !(!fd), PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file: %s", fname);
893: /* Output totals */
894: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Total Flop %14e %16.8e\n", petsc_TotalFlops, _TotalTime));
895: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Clock Resolution %g\n", 0.0));
896: /* Output actions */
897: if (def->petsc_logActions) {
898: PetscInt num_actions;
899: PetscCall(PetscLogActionArrayGetSize(def->petsc_actions, &num_actions, NULL));
900: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Actions accomplished %d\n", (int)num_actions));
901: for (int a = 0; a < num_actions; a++) {
902: Action *action;
904: PetscCall(PetscLogActionArrayGetRef(def->petsc_actions, a, &action));
905: 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));
906: }
907: }
908: /* Output objects */
909: if (def->petsc_logObjects) {
910: PetscInt num_objects;
912: PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
913: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Objects created %d destroyed %d\n", def->petsc_numObjectsCreated, def->petsc_numObjectsDestroyed));
914: for (int o = 0; o < num_objects; o++) {
915: Object *object = NULL;
917: PetscCall(PetscLogObjectArrayGetRef(def->petsc_objects, o, &object));
918: if (object->parent != -1) continue; // object with this id wasn't logged, probably a PetscLogHandler
919: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Parent ID: %d Memory: %d\n", object->parent, (int)object->mem));
920: if (!object->name[0]) {
921: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Name\n"));
922: } else {
923: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Name: %s\n", object->name));
924: }
925: if (!object->info[0]) {
926: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "No Info\n"));
927: } else {
928: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Info: %s\n", object->info));
929: }
930: }
931: }
932: /* Output events */
933: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "Event log:\n"));
934: PetscCall(PetscLogHandlerGetState(handler, &state));
935: PetscCall(PetscLogStateGetNumEvents(state, &num_events));
936: PetscCall(PetscLogStateGetCurrentStage(state, &curStage));
937: for (event = 0; event < num_events; event++) {
938: PetscEventPerfInfo *event_info;
940: PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, curStage, event, &event_info));
941: if (event_info->time != 0.0) flops = event_info->flops / event_info->time;
942: else flops = 0.0;
943: PetscCall(PetscFPrintf(PETSC_COMM_SELF, fd, "%d %16d %16g %16g %16g\n", event, event_info->count, event_info->flops, event_info->time, flops));
944: }
945: PetscCall(PetscFClose(PETSC_COMM_SELF, fd));
946: PetscFunctionReturn(PETSC_SUCCESS);
947: }
949: /*
950: PetscLogView_Detailed - Each process prints the times for its own events
952: */
953: static PetscErrorCode PetscLogHandlerView_Default_Detailed(PetscLogHandler handler, PetscViewer viewer)
954: {
955: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
956: PetscLogDouble locTotalTime, numRed, maxMem;
957: PetscInt numStages, numEvents;
958: MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
959: PetscMPIInt rank, size;
960: PetscLogGlobalNames global_stages, global_events;
961: PetscLogState state;
962: PetscEventPerfInfo zero_info;
964: PetscFunctionBegin;
965: PetscCall(PetscLogHandlerGetState(handler, &state));
966: PetscCallMPI(MPI_Comm_size(comm, &size));
967: PetscCallMPI(MPI_Comm_rank(comm, &rank));
968: /* Must preserve reduction count before we go on */
969: numRed = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
970: /* Get the total elapsed time */
971: PetscCall(PetscTime(&locTotalTime));
972: locTotalTime -= petsc_BaseTime;
973: PetscCall(PetscViewerASCIIPrintf(viewer, "size = %d\n", size));
974: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalTimes = {}\n"));
975: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessages = {}\n"));
976: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMessageLens = {}\n"));
977: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalReductions = {}\n"));
978: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalFlop = {}\n"));
979: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalObjects = {}\n"));
980: PetscCall(PetscViewerASCIIPrintf(viewer, "LocalMemory = {}\n"));
981: PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
982: PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
983: PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
984: PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
985: PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
986: PetscCall(PetscViewerASCIIPrintf(viewer, "Stages = {}\n"));
987: for (PetscInt stage = 0; stage < numStages; stage++) {
988: PetscInt stage_id;
989: const char *stage_name;
991: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
992: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
993: PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"] = {}\n", stage_name));
994: PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"summary\"] = {}\n", stage_name));
995: for (PetscInt event = 0; event < numEvents; event++) {
996: PetscEventPerfInfo *eventInfo = &zero_info;
997: PetscBool is_zero = PETSC_FALSE;
998: PetscInt event_id;
999: const char *event_name;
1001: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1002: PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1003: if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1004: is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE;
1005: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1006: if (!is_zero) { PetscCall(PetscViewerASCIIPrintf(viewer, "Stages[\"%s\"][\"%s\"] = {}\n", stage_name, event_name)); }
1007: }
1008: }
1009: PetscCall(PetscMallocGetMaximumUsage(&maxMem));
1010: PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1011: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalTimes[%d] = %g\n", rank, locTotalTime));
1012: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessages[%d] = %g\n", rank, (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct)));
1013: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMessageLens[%d] = %g\n", rank, (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len)));
1014: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalReductions[%d] = %g\n", rank, numRed));
1015: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalFlop[%d] = %g\n", rank, petsc_TotalFlops));
1016: {
1017: PetscInt num_objects;
1019: PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
1020: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalObjects[%d] = %d\n", rank, (int)num_objects));
1021: }
1022: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "LocalMemory[%d] = %g\n", rank, maxMem));
1023: PetscCall(PetscViewerFlush(viewer));
1024: for (PetscInt stage = 0; stage < numStages; stage++) {
1025: PetscEventPerfInfo *stage_perf_info = &zero_info;
1026: PetscInt stage_id;
1027: const char *stage_name;
1029: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1030: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1031: if (stage_id >= 0) {
1032: PetscStagePerf *stage_info;
1033: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info));
1034: stage_perf_info = &stage_info->perfInfo;
1035: }
1036: 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,
1037: stage_perf_info->numMessages, stage_perf_info->messageLength, stage_perf_info->numReductions, stage_perf_info->flops));
1038: for (PetscInt event = 0; event < numEvents; event++) {
1039: PetscEventPerfInfo *eventInfo = &zero_info;
1040: PetscBool is_zero = PETSC_FALSE;
1041: PetscInt event_id;
1042: const char *event_name;
1044: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1045: PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1046: if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1047: is_zero = eventInfo->count == 0 ? PETSC_TRUE : PETSC_FALSE;
1048: PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero));
1049: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1050: if (!is_zero) {
1051: 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,
1052: eventInfo->count, eventInfo->time, eventInfo->syncTime, eventInfo->numMessages, eventInfo->messageLength, eventInfo->numReductions, eventInfo->flops));
1053: if (eventInfo->dof[0] >= 0.) {
1054: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"dof\" : ["));
1055: for (PetscInt d = 0; d < 8; ++d) {
1056: if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", "));
1057: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->dof[d]));
1058: }
1059: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]"));
1060: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", \"error\" : ["));
1061: for (PetscInt e = 0; e < 8; ++e) {
1062: if (e > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ", "));
1063: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", eventInfo->errors[e]));
1064: }
1065: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "]"));
1066: }
1067: }
1068: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "}\n"));
1069: }
1070: }
1071: PetscCall(PetscViewerFlush(viewer));
1072: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1073: PetscCall(PetscLogGlobalNamesDestroy(&global_events));
1074: PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1075: PetscFunctionReturn(PETSC_SUCCESS);
1076: }
1078: /*
1079: PetscLogView_CSV - Each process prints the times for its own events in Comma-Separated Value Format
1080: */
1081: static PetscErrorCode PetscLogHandlerView_Default_CSV(PetscLogHandler handler, PetscViewer viewer)
1082: {
1083: PetscLogDouble locTotalTime, maxMem;
1084: PetscInt numStages, numEvents, stage, event;
1085: MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
1086: PetscMPIInt rank, size;
1087: PetscLogGlobalNames global_stages, global_events;
1088: PetscLogState state;
1089: PetscEventPerfInfo zero_info;
1091: PetscFunctionBegin;
1092: PetscCall(PetscLogHandlerGetState(handler, &state));
1093: PetscCallMPI(MPI_Comm_size(comm, &size));
1094: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1095: /* Must preserve reduction count before we go on */
1096: /* Get the total elapsed time */
1097: PetscCall(PetscTime(&locTotalTime));
1098: locTotalTime -= petsc_BaseTime;
1099: PetscCall(PetscMallocGetMaximumUsage(&maxMem));
1100: PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
1101: PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
1102: PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
1103: PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
1104: PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
1105: PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1106: 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));
1107: PetscCall(PetscViewerFlush(viewer));
1108: for (stage = 0; stage < numStages; stage++) {
1109: PetscEventPerfInfo *stage_perf_info;
1110: PetscInt stage_id;
1111: const char *stage_name;
1113: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1114: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1115: stage_perf_info = &zero_info;
1116: if (stage_id >= 0) {
1117: PetscStagePerf *stage_info;
1118: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_info));
1119: stage_perf_info = &stage_info->perfInfo;
1120: }
1121: 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));
1122: for (event = 0; event < numEvents; event++) {
1123: PetscEventPerfInfo *eventInfo = &zero_info;
1124: PetscBool is_zero = PETSC_FALSE;
1125: PetscInt event_id;
1126: const char *event_name;
1128: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1129: PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1130: if (event_id >= 0 && stage_id >= 0) PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &eventInfo));
1131: PetscCall(PetscMemcmp(eventInfo, &zero_info, sizeof(zero_info), &is_zero));
1132: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1133: if (!is_zero) {
1134: 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));
1135: if (eventInfo->dof[0] >= 0.) {
1136: for (PetscInt d = 0; d < 8; ++d) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->dof[d]));
1137: for (PetscInt e = 0; e < 8; ++e) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ",%g", eventInfo->errors[e]));
1138: }
1139: }
1140: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
1141: }
1142: }
1143: PetscCall(PetscViewerFlush(viewer));
1144: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1145: PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1146: PetscCall(PetscLogGlobalNamesDestroy(&global_events));
1147: PetscFunctionReturn(PETSC_SUCCESS);
1148: }
1150: static PetscErrorCode PetscLogViewWarnSync(PetscViewer viewer)
1151: {
1152: PetscFunctionBegin;
1153: if (!PetscLogSyncOn) PetscFunctionReturn(PETSC_SUCCESS);
1154: PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1155: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n"));
1156: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1157: PetscCall(PetscViewerASCIIPrintf(viewer, " # WARNING!!! #\n"));
1158: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1159: PetscCall(PetscViewerASCIIPrintf(viewer, " # This program was run with logging synchronization. #\n"));
1160: PetscCall(PetscViewerASCIIPrintf(viewer, " # This option provides more meaningful imbalance #\n"));
1161: PetscCall(PetscViewerASCIIPrintf(viewer, " # figures at the expense of slowing things down and #\n"));
1162: PetscCall(PetscViewerASCIIPrintf(viewer, " # providing a distorted view of the overall runtime. #\n"));
1163: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1164: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n\n\n"));
1165: PetscFunctionReturn(PETSC_SUCCESS);
1166: }
1168: static PetscErrorCode PetscLogViewWarnDebugging(PetscViewer viewer)
1169: {
1170: PetscFunctionBegin;
1171: if (PetscDefined(USE_DEBUG)) {
1172: PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1173: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n"));
1174: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1175: PetscCall(PetscViewerASCIIPrintf(viewer, " # WARNING!!! #\n"));
1176: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1177: PetscCall(PetscViewerASCIIPrintf(viewer, " # This code was compiled with a debugging option. #\n"));
1178: PetscCall(PetscViewerASCIIPrintf(viewer, " # To get timing results run ./configure #\n"));
1179: PetscCall(PetscViewerASCIIPrintf(viewer, " # using --with-debugging=no, the performance will #\n"));
1180: PetscCall(PetscViewerASCIIPrintf(viewer, " # be generally two or three times faster. #\n"));
1181: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1182: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n\n\n"));
1183: }
1184: PetscFunctionReturn(PETSC_SUCCESS);
1185: }
1187: static PetscErrorCode PetscLogViewWarnNoGpuAwareMpi(PetscViewer viewer)
1188: {
1189: #if defined(PETSC_HAVE_DEVICE)
1190: PetscMPIInt size;
1191: PetscBool deviceInitialized = PETSC_FALSE;
1193: PetscFunctionBegin;
1194: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
1195: for (int i = PETSC_DEVICE_HOST + 1; i < PETSC_DEVICE_MAX; ++i) {
1196: const PetscDeviceType dtype = PetscDeviceTypeCast(i);
1197: if (PetscDeviceInitialized(dtype)) { /* a non-host device was initialized */
1198: deviceInitialized = PETSC_TRUE;
1199: break;
1200: }
1201: }
1202: /* the last condition says petsc is configured with device but it is a pure CPU run, so don't print misleading warnings */
1203: if (use_gpu_aware_mpi || size == 1 || !deviceInitialized) PetscFunctionReturn(PETSC_SUCCESS);
1204: PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1205: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n"));
1206: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1207: PetscCall(PetscViewerASCIIPrintf(viewer, " # WARNING!!! #\n"));
1208: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1209: PetscCall(PetscViewerASCIIPrintf(viewer, " # This code was compiled with GPU support and you've #\n"));
1210: PetscCall(PetscViewerASCIIPrintf(viewer, " # created PETSc/GPU objects, but you intentionally #\n"));
1211: PetscCall(PetscViewerASCIIPrintf(viewer, " # used -use_gpu_aware_mpi 0, requiring PETSc to copy #\n"));
1212: PetscCall(PetscViewerASCIIPrintf(viewer, " # additional data between the GPU and CPU. To obtain #\n"));
1213: PetscCall(PetscViewerASCIIPrintf(viewer, " # meaningful timing results on multi-rank runs, use #\n"));
1214: PetscCall(PetscViewerASCIIPrintf(viewer, " # GPU-aware MPI instead. #\n"));
1215: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1216: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n\n\n"));
1217: PetscFunctionReturn(PETSC_SUCCESS);
1218: #else
1219: return PETSC_SUCCESS;
1220: #endif
1221: }
1223: static PetscErrorCode PetscLogViewWarnGpuTime(PetscViewer viewer)
1224: {
1225: #if defined(PETSC_HAVE_DEVICE)
1227: PetscFunctionBegin;
1228: if (!PetscLogGpuTimeFlag || petsc_gflops == 0) PetscFunctionReturn(PETSC_SUCCESS);
1229: PetscCall(PetscViewerASCIIPrintf(viewer, "\n\n"));
1230: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n"));
1231: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1232: PetscCall(PetscViewerASCIIPrintf(viewer, " # WARNING!!! #\n"));
1233: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1234: PetscCall(PetscViewerASCIIPrintf(viewer, " # This code was run with -log_view_gpu_time #\n"));
1235: PetscCall(PetscViewerASCIIPrintf(viewer, " # This provides accurate timing within the GPU kernels #\n"));
1236: PetscCall(PetscViewerASCIIPrintf(viewer, " # but can slow down the entire computation by a #\n"));
1237: PetscCall(PetscViewerASCIIPrintf(viewer, " # measurable amount. For fastest runs we recommend #\n"));
1238: PetscCall(PetscViewerASCIIPrintf(viewer, " # not using this option. #\n"));
1239: PetscCall(PetscViewerASCIIPrintf(viewer, " # #\n"));
1240: PetscCall(PetscViewerASCIIPrintf(viewer, " ##########################################################\n\n\n"));
1241: PetscFunctionReturn(PETSC_SUCCESS);
1242: #else
1243: return PETSC_SUCCESS;
1244: #endif
1245: }
1247: static PetscErrorCode PetscLogHandlerView_Default_Info(PetscLogHandler handler, PetscViewer viewer)
1248: {
1249: PetscLogHandler_Default def = (PetscLogHandler_Default)handler->data;
1250: char arch[128], hostname[128], username[128], pname[PETSC_MAX_PATH_LEN], date[128];
1251: PetscLogDouble locTotalTime, TotalTime, TotalFlops;
1252: PetscLogDouble numMessages, messageLength, avgMessLen, numReductions;
1253: PetscLogDouble stageTime, flops, flopr, mem, mess, messLen, red;
1254: PetscLogDouble fracTime, fracFlops, fracMessages, fracLength, fracReductions, fracMess, fracMessLen, fracRed;
1255: PetscLogDouble fracStageTime, fracStageFlops, fracStageMess, fracStageMessLen, fracStageRed;
1256: PetscLogDouble min, max, tot, ratio, avg, x, y;
1257: PetscLogDouble minf, maxf, totf, ratf, mint, maxt, tott, ratt, ratC, totm, totml, totr, mal, malmax, emalmax;
1258: #if defined(PETSC_HAVE_DEVICE)
1259: PetscLogEvent KSP_Solve, SNES_Solve, TS_Step, TAO_Solve; /* These need to be fixed to be some events registered with certain objects */
1260: PetscLogDouble cct, gct, csz, gsz, gmaxt, gflops, gflopr, fracgflops;
1261: #endif
1262: PetscMPIInt minC, maxC;
1263: PetscMPIInt size, rank;
1264: PetscBool *localStageUsed, *stageUsed;
1265: PetscBool *localStageVisible, *stageVisible;
1266: PetscInt numStages, numEvents;
1267: int stage, oclass;
1268: PetscLogEvent event;
1269: char version[256];
1270: MPI_Comm comm;
1271: #if defined(PETSC_HAVE_DEVICE)
1272: PetscInt64 nas = 0x7FF0000000000002;
1273: #endif
1274: PetscLogGlobalNames global_stages, global_events;
1275: PetscEventPerfInfo zero_info;
1276: PetscLogState state;
1278: PetscFunctionBegin;
1279: PetscCall(PetscLogHandlerGetState(handler, &state));
1280: PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1281: PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
1282: PetscCallMPI(MPI_Comm_size(comm, &size));
1283: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1284: /* Get the total elapsed time */
1285: PetscCall(PetscTime(&locTotalTime));
1286: locTotalTime -= petsc_BaseTime;
1288: PetscCall(PetscViewerASCIIPrintf(viewer, "****************************************************************************************************************************************************************\n"));
1289: PetscCall(PetscViewerASCIIPrintf(viewer, "*** WIDEN YOUR WINDOW TO 160 CHARACTERS. Use 'enscript -r -fCourier9' to print this document ***\n"));
1290: PetscCall(PetscViewerASCIIPrintf(viewer, "****************************************************************************************************************************************************************\n"));
1291: PetscCall(PetscViewerASCIIPrintf(viewer, "\n------------------------------------------------------------------ PETSc Performance Summary: ------------------------------------------------------------------\n\n"));
1292: PetscCall(PetscLogViewWarnSync(viewer));
1293: PetscCall(PetscLogViewWarnDebugging(viewer));
1294: PetscCall(PetscLogViewWarnNoGpuAwareMpi(viewer));
1295: PetscCall(PetscLogViewWarnGpuTime(viewer));
1296: PetscCall(PetscGetArchType(arch, sizeof(arch)));
1297: PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
1298: PetscCall(PetscGetUserName(username, sizeof(username)));
1299: PetscCall(PetscGetProgramName(pname, sizeof(pname)));
1300: PetscCall(PetscGetDate(date, sizeof(date)));
1301: PetscCall(PetscGetVersion(version, sizeof(version)));
1303: #if defined(PETSC_HAVE_CUPM)
1304: const char *cupm = PetscDefined(HAVE_CUDA) ? "CUDA" : "HIP";
1305: if (PetscDeviceCUPMRuntimeArch) PetscCall(PetscViewerASCIIPrintf(viewer, "%s on a %s named %s with %d processor(s) and %s architecture %d, by %s on %s\n", pname, arch, hostname, size, cupm, PetscDeviceCUPMRuntimeArch, username, date));
1306: else
1307: #endif
1308: PetscCall(PetscViewerASCIIPrintf(viewer, "%s on a %s named %s with %d processor(s), by %s on %s\n", pname, arch, hostname, size, username, date));
1310: #if defined(PETSC_HAVE_OPENMP)
1311: PetscCall(PetscViewerASCIIPrintf(viewer, "Using %" PetscInt_FMT " OpenMP threads\n", PetscNumOMPThreads));
1312: #endif
1313: PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s\n", version));
1315: /* Must preserve reduction count before we go on */
1316: red = petsc_allreduce_ct + petsc_gather_ct + petsc_scatter_ct;
1318: /* Calculate summary information */
1319: PetscCall(PetscViewerASCIIPrintf(viewer, "\n Max Max/Min Avg Total\n"));
1320: /* Time */
1321: PetscCall(MPIU_Allreduce(&locTotalTime, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1322: PetscCall(MPIU_Allreduce(&locTotalTime, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1323: PetscCall(MPIU_Allreduce(&locTotalTime, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1324: avg = tot / ((PetscLogDouble)size);
1325: if (min != 0.0) ratio = max / min;
1326: else ratio = 0.0;
1327: PetscCall(PetscViewerASCIIPrintf(viewer, "Time (sec): %5.3e %7.3f %5.3e\n", max, ratio, avg));
1328: TotalTime = tot;
1329: /* Objects */
1330: {
1331: PetscInt num_objects;
1333: PetscCall(PetscLogObjectArrayGetSize(def->petsc_objects, &num_objects, NULL));
1334: avg = (PetscLogDouble)num_objects;
1335: }
1336: PetscCall(MPIU_Allreduce(&avg, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1337: PetscCall(MPIU_Allreduce(&avg, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1338: PetscCall(MPIU_Allreduce(&avg, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1339: avg = tot / ((PetscLogDouble)size);
1340: if (min != 0.0) ratio = max / min;
1341: else ratio = 0.0;
1342: PetscCall(PetscViewerASCIIPrintf(viewer, "Objects: %5.3e %7.3f %5.3e\n", max, ratio, avg));
1343: /* Flops */
1344: PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1345: PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1346: PetscCall(MPIU_Allreduce(&petsc_TotalFlops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1347: avg = tot / ((PetscLogDouble)size);
1348: if (min != 0.0) ratio = max / min;
1349: else ratio = 0.0;
1350: PetscCall(PetscViewerASCIIPrintf(viewer, "Flops: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot));
1351: TotalFlops = tot;
1352: /* Flops/sec -- Must talk to Barry here */
1353: if (locTotalTime != 0.0) flops = petsc_TotalFlops / locTotalTime;
1354: else flops = 0.0;
1355: PetscCall(MPIU_Allreduce(&flops, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1356: PetscCall(MPIU_Allreduce(&flops, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1357: PetscCall(MPIU_Allreduce(&flops, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1358: avg = tot / ((PetscLogDouble)size);
1359: if (min != 0.0) ratio = max / min;
1360: else ratio = 0.0;
1361: PetscCall(PetscViewerASCIIPrintf(viewer, "Flops/sec: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot));
1362: /* Memory */
1363: PetscCall(PetscMallocGetMaximumUsage(&mem));
1364: if (mem > 0.0) {
1365: PetscCall(MPIU_Allreduce(&mem, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1366: PetscCall(MPIU_Allreduce(&mem, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1367: PetscCall(MPIU_Allreduce(&mem, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1368: avg = tot / ((PetscLogDouble)size);
1369: if (min != 0.0) ratio = max / min;
1370: else ratio = 0.0;
1371: PetscCall(PetscViewerASCIIPrintf(viewer, "Memory (bytes): %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot));
1372: }
1373: /* Messages */
1374: mess = 0.5 * (petsc_irecv_ct + petsc_isend_ct + petsc_recv_ct + petsc_send_ct);
1375: PetscCall(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1376: PetscCall(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1377: PetscCall(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1378: avg = tot / ((PetscLogDouble)size);
1379: if (min != 0.0) ratio = max / min;
1380: else ratio = 0.0;
1381: PetscCall(PetscViewerASCIIPrintf(viewer, "MPI Msg Count: %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot));
1382: numMessages = tot;
1383: /* Message Lengths */
1384: mess = 0.5 * (petsc_irecv_len + petsc_isend_len + petsc_recv_len + petsc_send_len);
1385: PetscCall(MPIU_Allreduce(&mess, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1386: PetscCall(MPIU_Allreduce(&mess, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1387: PetscCall(MPIU_Allreduce(&mess, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1388: if (numMessages != 0) avg = tot / numMessages;
1389: else avg = 0.0;
1390: if (min != 0.0) ratio = max / min;
1391: else ratio = 0.0;
1392: PetscCall(PetscViewerASCIIPrintf(viewer, "MPI Msg Len (bytes): %5.3e %7.3f %5.3e %5.3e\n", max, ratio, avg, tot));
1393: messageLength = tot;
1394: /* Reductions */
1395: PetscCall(MPIU_Allreduce(&red, &min, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1396: PetscCall(MPIU_Allreduce(&red, &max, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1397: PetscCall(MPIU_Allreduce(&red, &tot, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1398: if (min != 0.0) ratio = max / min;
1399: else ratio = 0.0;
1400: PetscCall(PetscViewerASCIIPrintf(viewer, "MPI Reductions: %5.3e %7.3f\n", max, ratio));
1401: numReductions = red; /* wrong because uses count from process zero */
1402: PetscCall(PetscViewerASCIIPrintf(viewer, "\nFlop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)\n"));
1403: PetscCall(PetscViewerASCIIPrintf(viewer, " e.g., VecAXPY() for real vectors of length N --> 2N flops\n"));
1404: PetscCall(PetscViewerASCIIPrintf(viewer, " and VecAXPY() for complex vectors of length N --> 8N flops\n"));
1406: PetscCall(PetscLogRegistryCreateGlobalStageNames(comm, state->registry, &global_stages));
1407: PetscCall(PetscLogRegistryCreateGlobalEventNames(comm, state->registry, &global_events));
1408: PetscCall(PetscLogGlobalNamesGetSize(global_stages, NULL, &numStages));
1409: PetscCall(PetscLogGlobalNamesGetSize(global_events, NULL, &numEvents));
1410: PetscCall(PetscMemzero(&zero_info, sizeof(zero_info)));
1411: PetscCall(PetscMalloc1(numStages, &localStageUsed));
1412: PetscCall(PetscMalloc1(numStages, &stageUsed));
1413: PetscCall(PetscMalloc1(numStages, &localStageVisible));
1414: PetscCall(PetscMalloc1(numStages, &stageVisible));
1415: if (numStages > 0) {
1416: for (stage = 0; stage < numStages; stage++) {
1417: PetscInt stage_id;
1419: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1420: if (stage_id >= 0) {
1421: PetscStagePerf *stage_info;
1423: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_info));
1424: localStageUsed[stage] = stage_info->used;
1425: localStageVisible[stage] = stage_info->perfInfo.visible;
1426: } else {
1427: localStageUsed[stage] = PETSC_FALSE;
1428: localStageVisible[stage] = PETSC_TRUE;
1429: }
1430: }
1431: PetscCall(MPIU_Allreduce(localStageUsed, stageUsed, numStages, MPIU_BOOL, MPI_LOR, comm));
1432: PetscCall(MPIU_Allreduce(localStageVisible, stageVisible, numStages, MPIU_BOOL, MPI_LAND, comm));
1433: for (stage = 0; stage < numStages; stage++) {
1434: if (stageUsed[stage] && stageVisible[stage]) {
1435: PetscCall(PetscViewerASCIIPrintf(viewer, "\nSummary of Stages: ----- Time ------ ----- Flop ------ --- Messages --- -- Message Lengths -- -- Reductions --\n"));
1436: PetscCall(PetscViewerASCIIPrintf(viewer, " Avg %%Total Avg %%Total Count %%Total Avg %%Total Count %%Total\n"));
1437: break;
1438: }
1439: }
1440: for (stage = 0; stage < numStages; stage++) {
1441: PetscInt stage_id;
1442: PetscEventPerfInfo *stage_info;
1443: const char *stage_name;
1445: if (!(stageUsed[stage] && stageVisible[stage])) continue;
1446: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1447: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1448: stage_info = &zero_info;
1449: if (localStageUsed[stage]) {
1450: PetscStagePerf *stage_perf_info;
1452: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage, &stage_perf_info));
1453: stage_info = &stage_perf_info->perfInfo;
1454: }
1455: PetscCall(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1456: PetscCall(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1457: PetscCall(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1458: PetscCall(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1459: PetscCall(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1460: mess *= 0.5;
1461: messLen *= 0.5;
1462: red /= size;
1463: if (TotalTime != 0.0) fracTime = stageTime / TotalTime;
1464: else fracTime = 0.0;
1465: if (TotalFlops != 0.0) fracFlops = flops / TotalFlops;
1466: else fracFlops = 0.0;
1467: /* Talk to Barry if (stageTime != 0.0) flops = (size*flops)/stageTime; else flops = 0.0; */
1468: if (numMessages != 0.0) fracMessages = mess / numMessages;
1469: else fracMessages = 0.0;
1470: if (mess != 0.0) avgMessLen = messLen / mess;
1471: else avgMessLen = 0.0;
1472: if (messageLength != 0.0) fracLength = messLen / messageLength;
1473: else fracLength = 0.0;
1474: if (numReductions != 0.0) fracReductions = red / numReductions;
1475: else fracReductions = 0.0;
1476: 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));
1477: }
1478: }
1480: PetscCall(PetscViewerASCIIPrintf(viewer, "\n------------------------------------------------------------------------------------------------------------------------\n"));
1481: PetscCall(PetscViewerASCIIPrintf(viewer, "See the 'Profiling' chapter of the users' manual for details on interpreting output.\n"));
1482: PetscCall(PetscViewerASCIIPrintf(viewer, "Phase summary info:\n"));
1483: PetscCall(PetscViewerASCIIPrintf(viewer, " Count: number of times phase was executed\n"));
1484: PetscCall(PetscViewerASCIIPrintf(viewer, " Time and Flop: Max - maximum over all processors\n"));
1485: PetscCall(PetscViewerASCIIPrintf(viewer, " Ratio - ratio of maximum to minimum over all processors\n"));
1486: PetscCall(PetscViewerASCIIPrintf(viewer, " Mess: number of messages sent\n"));
1487: PetscCall(PetscViewerASCIIPrintf(viewer, " AvgLen: average message length (bytes)\n"));
1488: PetscCall(PetscViewerASCIIPrintf(viewer, " Reduct: number of global reductions\n"));
1489: PetscCall(PetscViewerASCIIPrintf(viewer, " Global: entire computation\n"));
1490: PetscCall(PetscViewerASCIIPrintf(viewer, " Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().\n"));
1491: PetscCall(PetscViewerASCIIPrintf(viewer, " %%T - percent time in this phase %%F - percent flop in this phase\n"));
1492: PetscCall(PetscViewerASCIIPrintf(viewer, " %%M - percent messages in this phase %%L - percent message lengths in this phase\n"));
1493: PetscCall(PetscViewerASCIIPrintf(viewer, " %%R - percent reductions in this phase\n"));
1494: PetscCall(PetscViewerASCIIPrintf(viewer, " Total Mflop/s: 10e-6 * (sum of flop over all processors)/(max time over all processors)\n"));
1495: if (PetscLogMemory) {
1496: PetscCall(PetscViewerASCIIPrintf(viewer, " Memory usage is summed over all MPI processes, it is given in mega-bytes\n"));
1497: PetscCall(PetscViewerASCIIPrintf(viewer, " Malloc Mbytes: Memory allocated and kept during event (sum over all calls to event). May be negative\n"));
1498: PetscCall(PetscViewerASCIIPrintf(viewer, " EMalloc Mbytes: extra memory allocated during event and then freed (maximum over all calls to events). Never negative\n"));
1499: PetscCall(PetscViewerASCIIPrintf(viewer, " MMalloc Mbytes: Increase in high water mark of allocated memory (sum over all calls to event). Never negative\n"));
1500: PetscCall(PetscViewerASCIIPrintf(viewer, " RMI Mbytes: Increase in resident memory (sum over all calls to event)\n"));
1501: }
1502: #if defined(PETSC_HAVE_DEVICE)
1503: PetscCall(PetscViewerASCIIPrintf(viewer, " GPU Mflop/s: 10e-6 * (sum of flop on GPU over all processors)/(max GPU time over all processors)\n"));
1504: PetscCall(PetscViewerASCIIPrintf(viewer, " CpuToGpu Count: total number of CPU to GPU copies per processor\n"));
1505: PetscCall(PetscViewerASCIIPrintf(viewer, " CpuToGpu Size (Mbytes): 10e-6 * (total size of CPU to GPU copies per processor)\n"));
1506: PetscCall(PetscViewerASCIIPrintf(viewer, " GpuToCpu Count: total number of GPU to CPU copies per processor\n"));
1507: PetscCall(PetscViewerASCIIPrintf(viewer, " GpuToCpu Size (Mbytes): 10e-6 * (total size of GPU to CPU copies per processor)\n"));
1508: PetscCall(PetscViewerASCIIPrintf(viewer, " GPU %%F: percent flops on GPU in this event\n"));
1509: #endif
1510: PetscCall(PetscViewerASCIIPrintf(viewer, "------------------------------------------------------------------------------------------------------------------------\n"));
1512: PetscCall(PetscLogViewWarnDebugging(viewer));
1514: /* Report events */
1515: PetscCall(PetscViewerASCIIPrintf(viewer, "Event Count Time (sec) Flop --- Global --- --- Stage ---- Total"));
1516: if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, " Malloc EMalloc MMalloc RMI"));
1517: #if defined(PETSC_HAVE_DEVICE)
1518: PetscCall(PetscViewerASCIIPrintf(viewer, " GPU - CpuToGpu - - GpuToCpu - GPU"));
1519: #endif
1520: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1521: PetscCall(PetscViewerASCIIPrintf(viewer, " Max Ratio Max Ratio Max Ratio Mess AvgLen Reduct %%T %%F %%M %%L %%R %%T %%F %%M %%L %%R Mflop/s"));
1522: if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, " Mbytes Mbytes Mbytes Mbytes"));
1523: #if defined(PETSC_HAVE_DEVICE)
1524: PetscCall(PetscViewerASCIIPrintf(viewer, " Mflop/s Count Size Count Size %%F"));
1525: #endif
1526: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1527: PetscCall(PetscViewerASCIIPrintf(viewer, "------------------------------------------------------------------------------------------------------------------------"));
1528: if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, "-----------------------------"));
1529: #if defined(PETSC_HAVE_DEVICE)
1530: PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------------"));
1531: #endif
1532: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1534: #if defined(PETSC_HAVE_DEVICE)
1535: /* this indirect way of accessing these values is needed when PETSc is build with multiple libraries since the symbols are not in libpetscsys */
1536: PetscCall(PetscLogStateGetEventFromName(state, "TaoSolve", &TAO_Solve));
1537: PetscCall(PetscLogStateGetEventFromName(state, "TSStep", &TS_Step));
1538: PetscCall(PetscLogStateGetEventFromName(state, "SNESSolve", &SNES_Solve));
1539: PetscCall(PetscLogStateGetEventFromName(state, "KSPSolve", &KSP_Solve));
1540: #endif
1542: for (stage = 0; stage < numStages; stage++) {
1543: PetscInt stage_id;
1544: PetscEventPerfInfo *stage_info;
1545: const char *stage_name;
1547: if (!(stageVisible[stage] && stageUsed[stage])) continue;
1548: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_stages, stage, &stage_id));
1549: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1550: PetscCall(PetscViewerASCIIPrintf(viewer, "\n--- Event Stage %d: %s\n\n", stage, stage_name));
1551: stage_info = &zero_info;
1552: if (localStageUsed[stage]) {
1553: PetscStagePerf *stage_perf_info;
1555: PetscCall(PetscLogHandlerDefaultGetStageInfo(handler, stage_id, &stage_perf_info));
1556: stage_info = &stage_perf_info->perfInfo;
1557: }
1558: PetscCall(MPIU_Allreduce(&stage_info->time, &stageTime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1559: PetscCall(MPIU_Allreduce(&stage_info->flops, &flops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1560: PetscCall(MPIU_Allreduce(&stage_info->numMessages, &mess, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1561: PetscCall(MPIU_Allreduce(&stage_info->messageLength, &messLen, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1562: PetscCall(MPIU_Allreduce(&stage_info->numReductions, &red, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1563: mess *= 0.5;
1564: messLen *= 0.5;
1565: red /= size;
1567: for (event = 0; event < numEvents; event++) {
1568: PetscInt event_id;
1569: PetscEventPerfInfo *event_info = &zero_info;
1570: PetscBool is_zero = PETSC_FALSE;
1571: const char *event_name;
1573: PetscCall(PetscLogGlobalNamesGlobalGetLocal(global_events, event, &event_id));
1574: PetscCall(PetscLogGlobalNamesGlobalGetName(global_events, event, &event_name));
1575: if (event_id >= 0 && stage_id >= 0) { PetscCall(PetscLogHandlerGetEventPerfInfo_Default(handler, stage_id, event_id, &event_info)); }
1576: PetscCall(PetscMemcmp(event_info, &zero_info, sizeof(zero_info), &is_zero));
1577: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &is_zero, 1, MPIU_BOOL, MPI_LAND, comm));
1578: if (!is_zero) {
1579: flopr = event_info->flops;
1580: PetscCall(MPIU_Allreduce(&flopr, &minf, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1581: PetscCall(MPIU_Allreduce(&flopr, &maxf, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1582: PetscCall(MPIU_Allreduce(&event_info->flops, &totf, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1583: PetscCall(MPIU_Allreduce(&event_info->time, &mint, 1, MPIU_PETSCLOGDOUBLE, MPI_MIN, comm));
1584: PetscCall(MPIU_Allreduce(&event_info->time, &maxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1585: PetscCall(MPIU_Allreduce(&event_info->time, &tott, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1586: PetscCall(MPIU_Allreduce(&event_info->numMessages, &totm, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1587: PetscCall(MPIU_Allreduce(&event_info->messageLength, &totml, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1588: PetscCall(MPIU_Allreduce(&event_info->numReductions, &totr, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1589: PetscCall(MPIU_Allreduce(&event_info->count, &minC, 1, MPI_INT, MPI_MIN, comm));
1590: PetscCall(MPIU_Allreduce(&event_info->count, &maxC, 1, MPI_INT, MPI_MAX, comm));
1591: if (PetscLogMemory) {
1592: PetscCall(MPIU_Allreduce(&event_info->memIncrease, &mem, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1593: PetscCall(MPIU_Allreduce(&event_info->mallocSpace, &mal, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1594: PetscCall(MPIU_Allreduce(&event_info->mallocIncrease, &malmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1595: PetscCall(MPIU_Allreduce(&event_info->mallocIncreaseEvent, &emalmax, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1596: }
1597: #if defined(PETSC_HAVE_DEVICE)
1598: PetscCall(MPIU_Allreduce(&event_info->CpuToGpuCount, &cct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1599: PetscCall(MPIU_Allreduce(&event_info->GpuToCpuCount, &gct, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1600: PetscCall(MPIU_Allreduce(&event_info->CpuToGpuSize, &csz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1601: PetscCall(MPIU_Allreduce(&event_info->GpuToCpuSize, &gsz, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1602: PetscCall(MPIU_Allreduce(&event_info->GpuFlops, &gflops, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, comm));
1603: PetscCall(MPIU_Allreduce(&event_info->GpuTime, &gmaxt, 1, MPIU_PETSCLOGDOUBLE, MPI_MAX, comm));
1604: #endif
1605: if (mint < 0.0) {
1606: 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));
1607: mint = 0;
1608: }
1609: 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);
1610: #if defined(PETSC_HAVE_DEVICE)
1611: /* Put NaN into the time for all events that may not be time accurately since they may happen asynchronously on the GPU */
1612: if (!PetscLogGpuTimeFlag && petsc_gflops > 0) {
1613: memcpy(&gmaxt, &nas, sizeof(PetscLogDouble));
1614: if (event_id != SNES_Solve && event_id != KSP_Solve && event_id != TS_Step && event_id != TAO_Solve) {
1615: memcpy(&mint, &nas, sizeof(PetscLogDouble));
1616: memcpy(&maxt, &nas, sizeof(PetscLogDouble));
1617: }
1618: }
1619: #endif
1620: totm *= 0.5;
1621: totml *= 0.5;
1622: totr /= size;
1624: if (maxC != 0) {
1625: if (minC != 0) ratC = ((PetscLogDouble)maxC) / minC;
1626: else ratC = 0.0;
1627: if (mint != 0.0) ratt = maxt / mint;
1628: else ratt = 0.0;
1629: if (minf != 0.0) ratf = maxf / minf;
1630: else ratf = 0.0;
1631: if (TotalTime != 0.0) fracTime = tott / TotalTime;
1632: else fracTime = 0.0;
1633: if (TotalFlops != 0.0) fracFlops = totf / TotalFlops;
1634: else fracFlops = 0.0;
1635: if (stageTime != 0.0) fracStageTime = tott / stageTime;
1636: else fracStageTime = 0.0;
1637: if (flops != 0.0) fracStageFlops = totf / flops;
1638: else fracStageFlops = 0.0;
1639: if (numMessages != 0.0) fracMess = totm / numMessages;
1640: else fracMess = 0.0;
1641: if (messageLength != 0.0) fracMessLen = totml / messageLength;
1642: else fracMessLen = 0.0;
1643: if (numReductions != 0.0) fracRed = totr / numReductions;
1644: else fracRed = 0.0;
1645: if (mess != 0.0) fracStageMess = totm / mess;
1646: else fracStageMess = 0.0;
1647: if (messLen != 0.0) fracStageMessLen = totml / messLen;
1648: else fracStageMessLen = 0.0;
1649: if (red != 0.0) fracStageRed = totr / red;
1650: else fracStageRed = 0.0;
1651: if (totm != 0.0) totml /= totm;
1652: else totml = 0.0;
1653: if (maxt != 0.0) flopr = totf / maxt;
1654: else flopr = 0.0;
1655: if (fracStageTime > 1.0 || fracStageFlops > 1.0 || fracStageMess > 1.0 || fracStageMessLen > 1.0 || fracStageRed > 1.0)
1656: 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));
1657: else {
1658: 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"
1659: 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));
1660: } else {
1661: 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));
1662: }
1663: }
1664: 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));
1665: #if defined(PETSC_HAVE_DEVICE)
1666: if (totf != 0.0) fracgflops = gflops / totf;
1667: else fracgflops = 0.0;
1668: if (gmaxt != 0.0) gflopr = gflops / gmaxt;
1669: else gflopr = 0.0;
1670: if (PetscIsNanReal((PetscReal)gflopr)) {
1671: 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));
1672: } else {
1673: 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));
1674: }
1675: #endif
1676: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1677: }
1678: }
1679: }
1680: }
1682: /* Memory usage and object creation */
1683: PetscCall(PetscViewerASCIIPrintf(viewer, "------------------------------------------------------------------------------------------------------------------------"));
1684: if (PetscLogMemory) PetscCall(PetscViewerASCIIPrintf(viewer, "-----------------------------"));
1685: #if defined(PETSC_HAVE_DEVICE)
1686: PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------------"));
1687: #endif
1688: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1689: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1691: /* Right now, only stages on the first processor are reported here, meaning only objects associated with
1692: the global communicator, or MPI_COMM_SELF for proc 1. We really should report global stats and then
1693: stats for stages local to processor sets.
1694: */
1695: /* We should figure out the longest object name here (now 20 characters) */
1696: PetscCall(PetscViewerASCIIPrintf(viewer, "Object Type Creations Destructions. Reports information only for process 0.\n"));
1697: for (stage = 0; stage < numStages; stage++) {
1698: const char *stage_name;
1700: PetscCall(PetscLogGlobalNamesGlobalGetName(global_stages, stage, &stage_name));
1701: PetscCall(PetscViewerASCIIPrintf(viewer, "\n--- Event Stage %d: %s\n\n", stage, stage_name));
1702: if (localStageUsed[stage]) {
1703: PetscInt num_classes;
1705: PetscCall(PetscLogStateGetNumClasses(state, &num_classes));
1706: for (oclass = 0; oclass < num_classes; oclass++) {
1707: PetscClassPerf *class_perf_info;
1709: PetscCall(PetscLogHandlerDefaultGetClassPerf(handler, stage, oclass, &class_perf_info));
1710: if ((class_perf_info->creations > 0) || (class_perf_info->destructions > 0)) {
1711: PetscLogClassInfo class_reg_info;
1712: PetscBool flg;
1714: PetscCall(PetscLogStateClassGetInfo(state, oclass, &class_reg_info));
1715: if (stage == 0 && oclass == num_classes - 1) {
1716: PetscCall(PetscStrcmp(class_reg_info.name, "Viewer", &flg));
1717: 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");
1718: } else PetscCall(PetscViewerASCIIPrintf(viewer, "%20s %5d %5d\n", class_reg_info.name, class_perf_info->creations, class_perf_info->destructions));
1719: }
1720: }
1721: }
1722: }
1724: PetscCall(PetscFree(localStageUsed));
1725: PetscCall(PetscFree(stageUsed));
1726: PetscCall(PetscFree(localStageVisible));
1727: PetscCall(PetscFree(stageVisible));
1728: PetscCall(PetscLogGlobalNamesDestroy(&global_stages));
1729: PetscCall(PetscLogGlobalNamesDestroy(&global_events));
1731: /* Information unrelated to this particular run */
1732: PetscCall(PetscViewerASCIIPrintf(viewer, "========================================================================================================================\n"));
1733: PetscCall(PetscTime(&y));
1734: PetscCall(PetscTime(&x));
1735: PetscCall(PetscTime(&y));
1736: PetscCall(PetscTime(&y));
1737: PetscCall(PetscTime(&y));
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(PetscViewerASCIIPrintf(viewer, "Average time to get PetscTime(): %g\n", (y - x) / 10.0));
1746: /* MPI information */
1747: if (size > 1) {
1748: MPI_Status status;
1749: PetscMPIInt tag;
1750: MPI_Comm newcomm;
1752: PetscCallMPI(MPI_Barrier(comm));
1753: PetscCall(PetscTime(&x));
1754: PetscCallMPI(MPI_Barrier(comm));
1755: PetscCallMPI(MPI_Barrier(comm));
1756: PetscCallMPI(MPI_Barrier(comm));
1757: PetscCallMPI(MPI_Barrier(comm));
1758: PetscCallMPI(MPI_Barrier(comm));
1759: PetscCall(PetscTime(&y));
1760: PetscCall(PetscViewerASCIIPrintf(viewer, "Average time for MPI_Barrier(): %g\n", (y - x) / 5.0));
1761: PetscCall(PetscCommDuplicate(comm, &newcomm, &tag));
1762: PetscCallMPI(MPI_Barrier(comm));
1763: if (rank) {
1764: PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, rank - 1, tag, newcomm, &status));
1765: PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, (rank + 1) % size, tag, newcomm));
1766: } else {
1767: PetscCall(PetscTime(&x));
1768: PetscCallMPI(MPI_Send(NULL, 0, MPI_INT, 1, tag, newcomm));
1769: PetscCallMPI(MPI_Recv(NULL, 0, MPI_INT, size - 1, tag, newcomm, &status));
1770: PetscCall(PetscTime(&y));
1771: PetscCall(PetscViewerASCIIPrintf(viewer, "Average time for zero size MPI_Send(): %g\n", (y - x) / size));
1772: }
1773: PetscCall(PetscCommDestroy(&newcomm));
1774: }
1775: PetscCall(PetscOptionsView(NULL, viewer));
1777: /* Machine and compile information */
1778: if (PetscDefined(USE_FORTRAN_KERNELS)) {
1779: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with FORTRAN kernels\n"));
1780: } else {
1781: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled without FORTRAN kernels\n"));
1782: }
1783: if (PetscDefined(USE_64BIT_INDICES)) {
1784: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 64-bit PetscInt\n"));
1785: } else if (PetscDefined(USE___FLOAT128)) {
1786: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 32-bit PetscInt\n"));
1787: }
1788: if (PetscDefined(USE_REAL_SINGLE)) {
1789: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with single precision PetscScalar and PetscReal\n"));
1790: } else if (PetscDefined(USE___FLOAT128)) {
1791: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with 128 bit precision PetscScalar and PetscReal\n"));
1792: }
1793: if (PetscDefined(USE_REAL_MAT_SINGLE)) {
1794: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with single precision matrices\n"));
1795: } else {
1796: PetscCall(PetscViewerASCIIPrintf(viewer, "Compiled with full precision matrices (default)\n"));
1797: }
1798: 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)));
1800: PetscCall(PetscViewerASCIIPrintf(viewer, "Configure options: %s", petscconfigureoptions));
1801: PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petscmachineinfo));
1802: PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsccompilerinfo));
1803: PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsccompilerflagsinfo));
1804: PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsclinkerinfo));
1806: /* Cleanup */
1807: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1808: PetscCall(PetscLogViewWarnNoGpuAwareMpi(viewer));
1809: PetscCall(PetscLogViewWarnDebugging(viewer));
1810: PetscCall(PetscFPTrapPop());
1811: PetscFunctionReturn(PETSC_SUCCESS);
1812: }
1814: static PetscErrorCode PetscLogHandlerView_Default(PetscLogHandler handler, PetscViewer viewer)
1815: {
1816: PetscViewerFormat format;
1818: PetscFunctionBegin;
1819: PetscCall(PetscViewerGetFormat(viewer, &format));
1820: if (format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_ASCII_INFO) {
1821: PetscCall(PetscLogHandlerView_Default_Info(handler, viewer));
1822: } else if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1823: PetscCall(PetscLogHandlerView_Default_Detailed(handler, viewer));
1824: } else if (format == PETSC_VIEWER_ASCII_CSV) {
1825: PetscCall(PetscLogHandlerView_Default_CSV(handler, viewer));
1826: }
1827: PetscFunctionReturn(PETSC_SUCCESS);
1828: }
1830: /*MC
1831: PETSCLOGHANDLERDEFAULT - PETSCLOGHANDLERDEFAULT = "default" - A `PetscLogHandler` that collects data for PETSc
1832: default profiling log viewers (`PetscLogView()` and `PetscLogDump()`). A log handler of this type is
1833: created and started (`PetscLogHandlerStart()`) by `PetscLogDefaultBegin()`.
1835: Options Database Keys:
1836: + -log_include_actions - include a growing list of actions (event beginnings and endings, object creations and destructions) in `PetscLogDump()` (`PetscLogActions()`).
1837: - -log_include_objects - include a growing list of object creations and destructions in `PetscLogDump()` (`PetscLogObjects()`).
1839: Level: developer
1841: .seealso: [](ch_profiling), `PetscLogHandler`
1842: M*/
1844: PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Default(PetscLogHandler handler)
1845: {
1846: PetscFunctionBegin;
1847: PetscCall(PetscLogHandlerContextCreate_Default((PetscLogHandler_Default *)&handler->data));
1848: handler->ops->destroy = PetscLogHandlerDestroy_Default;
1849: handler->ops->eventbegin = PetscLogHandlerEventBegin_Default;
1850: handler->ops->eventend = PetscLogHandlerEventEnd_Default;
1851: handler->ops->eventsync = PetscLogHandlerEventSync_Default;
1852: handler->ops->objectcreate = PetscLogHandlerObjectCreate_Default;
1853: handler->ops->objectdestroy = PetscLogHandlerObjectDestroy_Default;
1854: handler->ops->stagepush = PetscLogHandlerStagePush_Default;
1855: handler->ops->stagepop = PetscLogHandlerStagePop_Default;
1856: handler->ops->view = PetscLogHandlerView_Default;
1857: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetEventPerfInfo_C", PetscLogHandlerGetEventPerfInfo_Default));
1858: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetStagePerfInfo_C", PetscLogHandlerGetStagePerfInfo_Default));
1859: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerSetLogActions_C", PetscLogHandlerSetLogActions_Default));
1860: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerSetLogObjects_C", PetscLogHandlerSetLogObjects_Default));
1861: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerLogObjectState_C", PetscLogHandlerLogObjectState_Default));
1862: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerGetNumObjects_C", PetscLogHandlerGetNumObjects_Default));
1863: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventDeactivatePush_C", PetscLogHandlerEventDeactivatePush_Default));
1864: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventDeactivatePop_C", PetscLogHandlerEventDeactivatePop_Default));
1865: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventsPause_C", PetscLogHandlerEventsPause_Default));
1866: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerEventsResume_C", PetscLogHandlerEventsResume_Default));
1867: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerDump_C", PetscLogHandlerDump_Default));
1868: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerStageSetVisible_C", PetscLogHandlerStageSetVisible_Default));
1869: PetscCall(PetscObjectComposeFunction((PetscObject)handler, "PetscLogHandlerStageGetVisible_C", PetscLogHandlerStageGetVisible_Default));
1870: PetscFunctionReturn(PETSC_SUCCESS);
1871: }