Actual source code: logregistry.c
1: #include <petsc/private/logimpl.h>
3: #define PETSC_LOG_RESIZABLE_ARRAY_HAS_NAME(Container, Entry, Key, Equal) \
4: static inline PETSC_UNUSED PetscErrorCode PetscLog##Container##Destructor(Entry *entry) \
5: { \
6: PetscFunctionBegin; \
7: PetscCall(PetscFree(entry->name)); \
8: PetscFunctionReturn(PETSC_SUCCESS); \
9: } \
10: PETSC_LOG_RESIZABLE_ARRAY(Container, Entry, Key, NULL, PetscLog##Container##Destructor, Equal)
12: #define PETSC_LOG_RESIZABLE_ARRAY_KEY_BY_NAME(Container, Entry) \
13: static inline PETSC_UNUSED PetscErrorCode PetscLog##Container##Equal(Entry *entry, const char *name, PetscBool *is_equal) \
14: { \
15: PetscFunctionBegin; \
16: PetscCall(PetscStrcmp(entry->name, name, is_equal)); \
17: PetscFunctionReturn(PETSC_SUCCESS); \
18: } \
19: PETSC_LOG_RESIZABLE_ARRAY_HAS_NAME(Container, Entry, const char *, PetscLog##Container##Equal)
21: static PetscErrorCode PetscLogClassArrayEqual(PetscLogClassInfo *class_info, PetscLogClassInfo *key, PetscBool *is_equal)
22: {
23: PetscFunctionBegin;
24: if (key->name) {
25: PetscCall(PetscStrcmp(class_info->name, key->name, is_equal));
26: } else {
27: *is_equal = (class_info->classid == key->classid) ? PETSC_TRUE : PETSC_FALSE;
28: }
29: PetscFunctionReturn(PETSC_SUCCESS);
30: }
32: PETSC_LOG_RESIZABLE_ARRAY_KEY_BY_NAME(EventArray, PetscLogEventInfo)
33: PETSC_LOG_RESIZABLE_ARRAY_KEY_BY_NAME(StageArray, PetscLogStageInfo)
34: PETSC_LOG_RESIZABLE_ARRAY_HAS_NAME(ClassArray, PetscLogClassInfo, PetscLogClassInfo *, PetscLogClassArrayEqual)
36: struct _n_PetscLogRegistry {
37: PetscLogEventArray events;
38: PetscLogClassArray classes;
39: PetscLogStageArray stages;
40: };
42: PETSC_INTERN PetscErrorCode PetscLogRegistryCreate(PetscLogRegistry *registry_p)
43: {
44: PetscLogRegistry registry;
46: PetscFunctionBegin;
47: PetscCall(PetscNew(registry_p));
48: registry = *registry_p;
49: PetscCall(PetscLogEventArrayCreate(128, ®istry->events));
50: PetscCall(PetscLogStageArrayCreate(8, ®istry->stages));
51: PetscCall(PetscLogClassArrayCreate(128, ®istry->classes));
52: PetscFunctionReturn(PETSC_SUCCESS);
53: }
55: PETSC_INTERN PetscErrorCode PetscLogRegistryDestroy(PetscLogRegistry registry)
56: {
57: PetscFunctionBegin;
58: PetscCall(PetscLogEventArrayDestroy(®istry->events));
59: PetscCall(PetscLogClassArrayDestroy(®istry->classes));
60: PetscCall(PetscLogStageArrayDestroy(®istry->stages));
61: PetscCall(PetscFree(registry));
62: PetscFunctionReturn(PETSC_SUCCESS);
63: }
65: PETSC_INTERN PetscErrorCode PetscLogRegistryGetNumEvents(PetscLogRegistry registry, PetscInt *num_events, PetscInt *max_events)
66: {
67: PetscFunctionBegin;
68: PetscCall(PetscLogEventArrayGetSize(registry->events, num_events, max_events));
69: PetscFunctionReturn(PETSC_SUCCESS);
70: }
72: PETSC_INTERN PetscErrorCode PetscLogRegistryGetNumStages(PetscLogRegistry registry, PetscInt *num_stages, PetscInt *max_stages)
73: {
74: PetscFunctionBegin;
75: PetscCall(PetscLogStageArrayGetSize(registry->stages, num_stages, max_stages));
76: PetscFunctionReturn(PETSC_SUCCESS);
77: }
79: PETSC_INTERN PetscErrorCode PetscLogRegistryGetNumClasses(PetscLogRegistry registry, PetscInt *num_classes, PetscInt *max_classes)
80: {
81: PetscFunctionBegin;
82: PetscCall(PetscLogClassArrayGetSize(registry->classes, num_classes, max_classes));
83: PetscFunctionReturn(PETSC_SUCCESS);
84: }
86: PETSC_INTERN PetscErrorCode PetscLogRegistryStageRegister(PetscLogRegistry registry, const char sname[], int *stage)
87: {
88: int idx;
89: PetscLogStageInfo stage_info;
91: PetscFunctionBegin;
92: PetscCall(PetscLogStageArrayFind(registry->stages, sname, &idx));
93: PetscCheck(idx == -1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "An event named %s is already registered", sname);
94: *stage = registry->stages->num_entries;
95: PetscCall(PetscStrallocpy(sname, &stage_info.name));
96: PetscCall(PetscLogStageArrayPush(registry->stages, stage_info));
97: PetscFunctionReturn(PETSC_SUCCESS);
98: }
100: PETSC_INTERN PetscErrorCode PetscLogRegistryEventRegister(PetscLogRegistry registry, const char name[], PetscClassId classid, PetscLogEvent *event)
101: {
102: PetscLogEventInfo new_info;
104: PetscFunctionBegin;
105: PetscCall(PetscLogRegistryGetEventFromName(registry, name, event));
106: if (*event >= 0) PetscFunctionReturn(PETSC_SUCCESS);
107: *event = registry->events->num_entries;
108: new_info.classid = classid;
109: new_info.collective = PETSC_TRUE;
110: PetscCall(PetscStrallocpy(name, &new_info.name));
111: PetscCall(PetscLogEventArrayPush(registry->events, new_info));
112: PetscFunctionReturn(PETSC_SUCCESS);
113: }
115: PETSC_INTERN PetscErrorCode PetscLogRegistryClassRegister(PetscLogRegistry registry, const char name[], PetscClassId classid, PetscLogClass *clss)
116: {
117: PetscLogClassInfo new_info;
119: PetscFunctionBegin;
120: PetscCall(PetscLogRegistryGetClassFromClassId(registry, classid, clss));
121: if (*clss >= 0) PetscFunctionReturn(PETSC_SUCCESS);
122: *clss = registry->classes->num_entries;
123: new_info.classid = classid;
124: PetscCall(PetscStrallocpy(name, &new_info.name));
125: PetscCall(PetscLogClassArrayPush(registry->classes, new_info));
126: PetscFunctionReturn(PETSC_SUCCESS);
127: }
129: PETSC_INTERN PetscErrorCode PetscLogRegistryGetEventFromName(PetscLogRegistry registry, const char name[], PetscLogEvent *event)
130: {
131: PetscFunctionBegin;
132: PetscCall(PetscLogEventArrayFind(registry->events, name, event));
133: PetscFunctionReturn(PETSC_SUCCESS);
134: }
136: PETSC_INTERN PetscErrorCode PetscLogRegistryGetStageFromName(PetscLogRegistry registry, const char name[], PetscLogStage *stage)
137: {
138: PetscFunctionBegin;
139: PetscCall(PetscLogStageArrayFind(registry->stages, name, stage));
140: PetscFunctionReturn(PETSC_SUCCESS);
141: }
143: PETSC_INTERN PetscErrorCode PetscLogRegistryGetClassFromClassId(PetscLogRegistry registry, PetscClassId classid, PetscLogStage *clss)
144: {
145: PetscLogClassInfo key;
147: PetscFunctionBegin;
148: key.name = NULL;
149: key.classid = classid;
150: PetscCall(PetscLogClassArrayFind(registry->classes, &key, clss));
151: PetscFunctionReturn(PETSC_SUCCESS);
152: }
154: PETSC_INTERN PetscErrorCode PetscLogRegistryGetClassFromName(PetscLogRegistry registry, const char name[], PetscLogStage *clss)
155: {
156: PetscLogClassInfo key;
158: PetscFunctionBegin;
159: key.name = (char *)name;
160: key.classid = -1;
161: PetscCall(PetscLogClassArrayFind(registry->classes, &key, clss));
162: PetscFunctionReturn(PETSC_SUCCESS);
163: }
165: PETSC_INTERN PetscErrorCode PetscLogRegistryEventGetInfo(PetscLogRegistry registry, PetscLogEvent event, PetscLogEventInfo *event_info)
166: {
167: PetscFunctionBegin;
168: PetscCall(PetscLogEventArrayGet(registry->events, event, event_info));
169: PetscFunctionReturn(PETSC_SUCCESS);
170: }
172: PETSC_INTERN PetscErrorCode PetscLogRegistryStageGetInfo(PetscLogRegistry registry, PetscLogStage stage, PetscLogStageInfo *stage_info)
173: {
174: PetscFunctionBegin;
175: PetscCall(PetscLogStageArrayGet(registry->stages, stage, stage_info));
176: PetscFunctionReturn(PETSC_SUCCESS);
177: }
179: PETSC_INTERN PetscErrorCode PetscLogRegistryClassGetInfo(PetscLogRegistry registry, PetscLogClass clss, PetscLogClassInfo *class_info)
180: {
181: PetscFunctionBegin;
182: PetscCall(PetscLogClassArrayGet(registry->classes, clss, class_info));
183: PetscFunctionReturn(PETSC_SUCCESS);
184: }
186: PETSC_INTERN PetscErrorCode PetscLogRegistryEventSetCollective(PetscLogRegistry registry, PetscLogEvent event, PetscBool collective)
187: {
188: PetscLogEventInfo *event_info;
190: PetscFunctionBegin;
191: PetscCall(PetscLogEventArrayGetRef(registry->events, event, &event_info));
192: event_info->collective = collective;
193: PetscFunctionReturn(PETSC_SUCCESS);
194: }
196: /* Given a list of strings on each process, create a global numbering. Order
197: them by their order on the first process, then the remaining by their order
198: on the second process, etc. The expectation is that most processes have the
199: same names in the same order so it shouldn't take too many rounds to figure
200: out */
202: struct _n_PetscLogGlobalNames {
203: MPI_Comm comm;
204: PetscInt count_global;
205: PetscInt count_local;
206: const char **names;
207: PetscInt *global_to_local;
208: PetscInt *local_to_global;
209: };
211: PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wconversion")
212: static PetscErrorCode PetscLogGlobalNamesCreate_Internal(MPI_Comm comm, PetscInt num_names_local, const char **names, PetscInt *num_names_global_p, PetscInt **global_index_to_local_index_p, PetscInt **local_index_to_global_index_p, const char ***global_names_p)
213: {
214: PetscMPIInt size, rank;
215: PetscInt num_names_global = 0;
216: PetscInt num_names_local_remaining = num_names_local;
217: PetscBool *local_name_seen;
218: PetscInt *global_index_to_local_index = NULL;
219: PetscInt *local_index_to_global_index = NULL;
220: PetscInt max_name_len = 0;
221: char *str_buffer;
222: char **global_names = NULL;
223: PetscMPIInt p;
225: PetscFunctionBegin;
226: PetscCallMPI(MPI_Comm_size(comm, &size));
227: if (size == 1) {
228: PetscCall(PetscMalloc1(num_names_local, &global_index_to_local_index));
229: PetscCall(PetscMalloc1(num_names_local, &local_index_to_global_index));
230: PetscCall(PetscMalloc1(num_names_local, &global_names));
231: for (PetscInt i = 0; i < num_names_local; i++) {
232: global_index_to_local_index[i] = i;
233: local_index_to_global_index[i] = i;
234: PetscCall(PetscStrallocpy(names[i], &global_names[i]));
235: }
236: *num_names_global_p = num_names_local;
237: *global_index_to_local_index_p = global_index_to_local_index;
238: *local_index_to_global_index_p = local_index_to_global_index;
239: *global_names_p = (const char **)global_names;
240: PetscFunctionReturn(PETSC_SUCCESS);
241: }
242: PetscCallMPI(MPI_Comm_rank(comm, &rank));
243: PetscCall(PetscCalloc1(num_names_local, &local_name_seen));
244: PetscCall(PetscMalloc1(num_names_local, &local_index_to_global_index));
246: for (PetscInt i = 0; i < num_names_local; i++) {
247: size_t i_len;
248: PetscCall(PetscStrlen(names[i], &i_len));
249: max_name_len = PetscMax(max_name_len, (PetscInt)i_len);
250: }
251: PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &max_name_len, 1, MPIU_INT, MPI_MAX, comm));
252: PetscCall(PetscCalloc1(max_name_len + 1, &str_buffer));
254: p = 0;
255: while (p < size) {
256: PetscInt my_loc, next_loc;
257: PetscInt num_to_add;
259: my_loc = num_names_local_remaining > 0 ? rank : PETSC_MPI_INT_MAX;
260: PetscCallMPI(MPIU_Allreduce(&my_loc, &next_loc, 1, MPIU_INT, MPI_MIN, comm));
261: if (next_loc == PETSC_MPI_INT_MAX) break;
262: PetscAssert(next_loc >= p, comm, PETSC_ERR_PLIB, "Failed invariant, expected increasing next process");
263: PetscCall(PetscCIntCast(next_loc, &p));
264: num_to_add = (rank == p) ? num_names_local_remaining : -1;
265: PetscCallMPI(MPI_Bcast(&num_to_add, 1, MPIU_INT, p, comm));
266: {
267: PetscInt new_num_names_global = num_names_global + num_to_add;
268: PetscInt *new_global_index_to_local_index;
269: char **new_global_names;
271: PetscCall(PetscMalloc1(new_num_names_global, &new_global_index_to_local_index));
272: PetscCall(PetscArraycpy(new_global_index_to_local_index, global_index_to_local_index, num_names_global));
273: for (PetscInt i = num_names_global; i < new_num_names_global; i++) new_global_index_to_local_index[i] = -1;
274: PetscCall(PetscFree(global_index_to_local_index));
275: global_index_to_local_index = new_global_index_to_local_index;
277: PetscCall(PetscCalloc1(new_num_names_global, &new_global_names));
278: PetscCall(PetscArraycpy(new_global_names, global_names, num_names_global));
279: PetscCall(PetscFree(global_names));
280: global_names = new_global_names;
281: }
283: if (rank == p) {
284: for (PetscInt s = 0; s < num_names_local; s++) {
285: if (local_name_seen[s]) continue;
286: local_name_seen[s] = PETSC_TRUE;
287: PetscCall(PetscArrayzero(str_buffer, max_name_len + 1));
288: PetscCall(PetscStrallocpy(names[s], &global_names[num_names_global]));
289: PetscCall(PetscStrncpy(str_buffer, names[s], max_name_len + 1));
290: PetscCallMPI(MPI_Bcast(str_buffer, max_name_len + 1, MPI_CHAR, p, comm));
291: local_index_to_global_index[s] = num_names_global;
292: global_index_to_local_index[num_names_global++] = s;
293: num_names_local_remaining--;
294: }
295: } else {
296: for (PetscInt i = 0, s; i < num_to_add; i++) {
297: PetscCallMPI(MPI_Bcast(str_buffer, max_name_len + 1, MPI_CHAR, p, comm));
298: PetscCall(PetscStrallocpy(str_buffer, &global_names[num_names_global]));
299: for (s = 0; s < num_names_local; s++) {
300: PetscBool same;
302: if (local_name_seen[s]) continue;
303: PetscCall(PetscStrncmp(names[s], str_buffer, max_name_len + 1, &same));
304: if (same) {
305: local_name_seen[s] = PETSC_TRUE;
306: global_index_to_local_index[num_names_global] = s;
307: local_index_to_global_index[s] = num_names_global;
308: num_names_local_remaining--;
309: break;
310: }
311: }
312: if (s == num_names_local) global_index_to_local_index[num_names_global] = -1; // this name is not present on this process
313: num_names_global++;
314: }
315: }
316: }
318: PetscCall(PetscFree(str_buffer));
319: PetscCall(PetscFree(local_name_seen));
320: *num_names_global_p = num_names_global;
321: *global_index_to_local_index_p = global_index_to_local_index;
322: *local_index_to_global_index_p = local_index_to_global_index;
323: *global_names_p = (const char **)global_names;
324: PetscFunctionReturn(PETSC_SUCCESS);
325: }
326: PETSC_PRAGMA_DIAGNOSTIC_IGNORED_END()
328: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesCreate(MPI_Comm comm, PetscInt num_names_local, const char **local_names, PetscLogGlobalNames *global_names_p)
329: {
330: PetscLogGlobalNames global_names;
332: PetscFunctionBegin;
333: PetscCall(PetscNew(&global_names));
334: PetscCall(PetscLogGlobalNamesCreate_Internal(comm, num_names_local, local_names, &global_names->count_global, &global_names->global_to_local, &global_names->local_to_global, &global_names->names));
335: global_names->count_local = num_names_local;
336: *global_names_p = global_names;
337: PetscFunctionReturn(PETSC_SUCCESS);
338: }
340: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesDestroy(PetscLogGlobalNames *global_names_p)
341: {
342: PetscLogGlobalNames global_names;
344: PetscFunctionBegin;
345: global_names = *global_names_p;
346: *global_names_p = NULL;
347: PetscCall(PetscFree(global_names->global_to_local));
348: PetscCall(PetscFree(global_names->local_to_global));
349: for (PetscInt i = 0; i < global_names->count_global; i++) PetscCall(PetscFree(global_names->names[i]));
350: PetscCall(PetscFree(global_names->names));
351: PetscCall(PetscFree(global_names));
352: PetscFunctionReturn(PETSC_SUCCESS);
353: }
355: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGlobalGetName(PetscLogGlobalNames global_names, PetscInt idx, const char **name)
356: {
357: PetscFunctionBegin;
358: PetscCheck(idx >= 0 && idx < global_names->count_global, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Idx %" PetscInt_FMT " not in range [0,%" PetscInt_FMT ")", idx, global_names->count_global);
359: *name = global_names->names[idx];
360: PetscFunctionReturn(PETSC_SUCCESS);
361: }
363: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGlobalGetLocal(PetscLogGlobalNames global_names, PetscInt idx, PetscInt *local_idx)
364: {
365: PetscFunctionBegin;
366: PetscCheck(idx >= 0 && idx < global_names->count_global, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Idx %" PetscInt_FMT " not in range [0,%" PetscInt_FMT ")", idx, global_names->count_global);
367: *local_idx = global_names->global_to_local[idx];
368: PetscFunctionReturn(PETSC_SUCCESS);
369: }
371: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesLocalGetGlobal(PetscLogGlobalNames global_names, PetscInt local_idx, PetscInt *idx)
372: {
373: PetscFunctionBegin;
374: PetscCheck(local_idx >= 0 && local_idx < global_names->count_local, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Idx %" PetscInt_FMT " not in range [0,%" PetscInt_FMT ")", local_idx, global_names->count_local);
375: *idx = global_names->local_to_global[local_idx];
376: PetscFunctionReturn(PETSC_SUCCESS);
377: }
379: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGetSize(PetscLogGlobalNames global_names, PetscInt *local_size, PetscInt *global_size)
380: {
381: PetscFunctionBegin;
382: if (local_size) *local_size = global_names->count_local;
383: if (global_size) *global_size = global_names->count_global;
384: PetscFunctionReturn(PETSC_SUCCESS);
385: }
387: PETSC_INTERN PetscErrorCode PetscLogRegistryCreateGlobalStageNames(MPI_Comm comm, PetscLogRegistry registry, PetscLogGlobalNames *global_names_p)
388: {
389: PetscInt num_stages_local;
390: const char **names;
392: PetscFunctionBegin;
393: PetscCall(PetscLogStageArrayGetSize(registry->stages, &num_stages_local, NULL));
394: PetscCall(PetscMalloc1(num_stages_local, &names));
395: for (PetscLogEvent i = 0; i < num_stages_local; i++) {
396: PetscLogStageInfo stage_info = {NULL};
397: PetscCall(PetscLogRegistryStageGetInfo(registry, i, &stage_info));
398: names[i] = stage_info.name;
399: }
400: PetscCall(PetscLogGlobalNamesCreate(comm, num_stages_local, names, global_names_p));
401: PetscCall(PetscFree(names));
402: PetscFunctionReturn(PETSC_SUCCESS);
403: }
405: PETSC_INTERN PetscErrorCode PetscLogRegistryCreateGlobalEventNames(MPI_Comm comm, PetscLogRegistry registry, PetscLogGlobalNames *global_names_p)
406: {
407: PetscInt num_events_local;
408: const char **names;
410: PetscFunctionBegin;
411: PetscCall(PetscLogEventArrayGetSize(registry->events, &num_events_local, NULL));
412: PetscCall(PetscMalloc1(num_events_local, &names));
413: for (PetscLogEvent i = 0; i < num_events_local; i++) {
414: PetscLogEventInfo event_info = {NULL, 0, PETSC_FALSE};
416: PetscCall(PetscLogRegistryEventGetInfo(registry, i, &event_info));
417: names[i] = event_info.name;
418: }
419: PetscCall(PetscLogGlobalNamesCreate(comm, num_events_local, names, global_names_p));
420: PetscCall(PetscFree(names));
421: PetscFunctionReturn(PETSC_SUCCESS);
422: }