xref: /petsc/src/sys/logging/state/logregistry.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
1 
2 #include <petsc/private/logimpl.h> /*I "petsclog.h" I*/
3 
4 #define PETSC_LOG_RESIZABLE_ARRAY_HAS_NAME(Container, Entry, Key, Equal) \
5   static inline PETSC_UNUSED PetscErrorCode PetscLog##Container##Destructor(Entry *entry) \
6   { \
7     PetscFunctionBegin; \
8     PetscCall(PetscFree(entry->name)); \
9     PetscFunctionReturn(PETSC_SUCCESS); \
10   } \
11   PETSC_LOG_RESIZABLE_ARRAY(Container, Entry, Key, NULL, PetscLog##Container##Destructor, Equal)
12 
13 #define PETSC_LOG_RESIZABLE_ARRAY_KEY_BY_NAME(Container, Entry) \
14   static inline PETSC_UNUSED PetscErrorCode PetscLog##Container##Equal(Entry *entry, const char *name, PetscBool *is_equal) \
15   { \
16     PetscFunctionBegin; \
17     PetscCall(PetscStrcmp(entry->name, name, is_equal)); \
18     PetscFunctionReturn(PETSC_SUCCESS); \
19   } \
20   PETSC_LOG_RESIZABLE_ARRAY_HAS_NAME(Container, Entry, const char *, PetscLog##Container##Equal)
21 
22 static PetscErrorCode PetscLogClassArrayEqual(PetscLogClassInfo *class_info, PetscLogClassInfo *key, PetscBool *is_equal)
23 {
24   PetscFunctionBegin;
25   if (key->name) {
26     PetscCall(PetscStrcmp(class_info->name, key->name, is_equal));
27   } else {
28     *is_equal = (class_info->classid == key->classid) ? PETSC_TRUE : PETSC_FALSE;
29   }
30   PetscFunctionReturn(PETSC_SUCCESS);
31 }
32 
33 PETSC_LOG_RESIZABLE_ARRAY_KEY_BY_NAME(EventArray, PetscLogEventInfo)
34 PETSC_LOG_RESIZABLE_ARRAY_KEY_BY_NAME(StageArray, PetscLogStageInfo)
35 PETSC_LOG_RESIZABLE_ARRAY_HAS_NAME(ClassArray, PetscLogClassInfo, PetscLogClassInfo *, PetscLogClassArrayEqual)
36 
37 struct _n_PetscLogRegistry {
38   PetscLogEventArray events;
39   PetscLogClassArray classes;
40   PetscLogStageArray stages;
41 };
42 
43 PETSC_INTERN PetscErrorCode PetscLogRegistryCreate(PetscLogRegistry *registry_p)
44 {
45   PetscLogRegistry registry;
46 
47   PetscFunctionBegin;
48   PetscCall(PetscNew(registry_p));
49   registry = *registry_p;
50   PetscCall(PetscLogEventArrayCreate(128, &registry->events));
51   PetscCall(PetscLogStageArrayCreate(8, &registry->stages));
52   PetscCall(PetscLogClassArrayCreate(128, &registry->classes));
53   PetscFunctionReturn(PETSC_SUCCESS);
54 }
55 
56 PETSC_INTERN PetscErrorCode PetscLogRegistryDestroy(PetscLogRegistry registry)
57 {
58   PetscFunctionBegin;
59   PetscCall(PetscLogEventArrayDestroy(&registry->events));
60   PetscCall(PetscLogClassArrayDestroy(&registry->classes));
61   PetscCall(PetscLogStageArrayDestroy(&registry->stages));
62   PetscCall(PetscFree(registry));
63   PetscFunctionReturn(PETSC_SUCCESS);
64 }
65 
66 PETSC_INTERN PetscErrorCode PetscLogRegistryGetNumEvents(PetscLogRegistry registry, PetscInt *num_events, PetscInt *max_events)
67 {
68   PetscFunctionBegin;
69   PetscCall(PetscLogEventArrayGetSize(registry->events, num_events, max_events));
70   PetscFunctionReturn(PETSC_SUCCESS);
71 }
72 
73 PETSC_INTERN PetscErrorCode PetscLogRegistryGetNumStages(PetscLogRegistry registry, PetscInt *num_stages, PetscInt *max_stages)
74 {
75   PetscFunctionBegin;
76   PetscCall(PetscLogStageArrayGetSize(registry->stages, num_stages, max_stages));
77   PetscFunctionReturn(PETSC_SUCCESS);
78 }
79 
80 PETSC_INTERN PetscErrorCode PetscLogRegistryGetNumClasses(PetscLogRegistry registry, PetscInt *num_classes, PetscInt *max_classes)
81 {
82   PetscFunctionBegin;
83   PetscCall(PetscLogClassArrayGetSize(registry->classes, num_classes, max_classes));
84   PetscFunctionReturn(PETSC_SUCCESS);
85 }
86 
87 PETSC_INTERN PetscErrorCode PetscLogRegistryStageRegister(PetscLogRegistry registry, const char sname[], int *stage)
88 {
89   int               idx;
90   PetscLogStageInfo stage_info;
91 
92   PetscFunctionBegin;
93   PetscCall(PetscLogStageArrayFind(registry->stages, sname, &idx));
94   PetscCheck(idx == -1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "An event named %s is already registered", sname);
95   *stage = registry->stages->num_entries;
96   PetscCall(PetscStrallocpy(sname, &stage_info.name));
97   PetscCall(PetscLogStageArrayPush(registry->stages, stage_info));
98   PetscFunctionReturn(PETSC_SUCCESS);
99 }
100 
101 PETSC_INTERN PetscErrorCode PetscLogRegistryEventRegister(PetscLogRegistry registry, const char name[], PetscClassId classid, PetscLogEvent *event)
102 {
103   PetscLogEventInfo new_info;
104 
105   PetscFunctionBegin;
106   PetscCall(PetscLogRegistryGetEventFromName(registry, name, event));
107   if (*event >= 0) PetscFunctionReturn(PETSC_SUCCESS);
108   *event              = registry->events->num_entries;
109   new_info.classid    = classid;
110   new_info.collective = PETSC_TRUE;
111   PetscCall(PetscStrallocpy(name, &new_info.name));
112   PetscCall(PetscLogEventArrayPush(registry->events, new_info));
113   PetscFunctionReturn(PETSC_SUCCESS);
114 }
115 
116 PETSC_INTERN PetscErrorCode PetscLogRegistryClassRegister(PetscLogRegistry registry, const char name[], PetscClassId classid, PetscLogClass *clss)
117 {
118   PetscLogClassInfo new_info;
119 
120   PetscFunctionBegin;
121   PetscCall(PetscLogRegistryGetClassFromClassId(registry, classid, clss));
122   if (*clss >= 0) PetscFunctionReturn(PETSC_SUCCESS);
123   *clss            = registry->classes->num_entries;
124   new_info.classid = classid;
125   PetscCall(PetscStrallocpy(name, &new_info.name));
126   PetscCall(PetscLogClassArrayPush(registry->classes, new_info));
127   PetscFunctionReturn(PETSC_SUCCESS);
128 }
129 
130 PETSC_INTERN PetscErrorCode PetscLogRegistryGetEventFromName(PetscLogRegistry registry, const char name[], PetscLogEvent *event)
131 {
132   PetscFunctionBegin;
133   PetscCall(PetscLogEventArrayFind(registry->events, name, event));
134   PetscFunctionReturn(PETSC_SUCCESS);
135 }
136 
137 PETSC_INTERN PetscErrorCode PetscLogRegistryGetStageFromName(PetscLogRegistry registry, const char name[], PetscLogStage *stage)
138 {
139   PetscFunctionBegin;
140   PetscCall(PetscLogStageArrayFind(registry->stages, name, stage));
141   PetscFunctionReturn(PETSC_SUCCESS);
142 }
143 
144 PETSC_INTERN PetscErrorCode PetscLogRegistryGetClassFromClassId(PetscLogRegistry registry, PetscClassId classid, PetscLogStage *clss)
145 {
146   PetscLogClassInfo key;
147 
148   PetscFunctionBegin;
149   key.name    = NULL;
150   key.classid = classid;
151   PetscCall(PetscLogClassArrayFind(registry->classes, &key, clss));
152   PetscFunctionReturn(PETSC_SUCCESS);
153 }
154 
155 PETSC_INTERN PetscErrorCode PetscLogRegistryGetClassFromName(PetscLogRegistry registry, const char name[], PetscLogStage *clss)
156 {
157   PetscLogClassInfo key;
158 
159   PetscFunctionBegin;
160   key.name    = (char *)name;
161   key.classid = -1;
162   PetscCall(PetscLogClassArrayFind(registry->classes, &key, clss));
163   PetscFunctionReturn(PETSC_SUCCESS);
164 }
165 
166 PETSC_INTERN PetscErrorCode PetscLogRegistryEventGetInfo(PetscLogRegistry registry, PetscLogEvent event, PetscLogEventInfo *event_info)
167 {
168   PetscFunctionBegin;
169   PetscCall(PetscLogEventArrayGet(registry->events, event, event_info));
170   PetscFunctionReturn(PETSC_SUCCESS);
171 }
172 
173 PETSC_INTERN PetscErrorCode PetscLogRegistryStageGetInfo(PetscLogRegistry registry, PetscLogStage stage, PetscLogStageInfo *stage_info)
174 {
175   PetscFunctionBegin;
176   PetscCall(PetscLogStageArrayGet(registry->stages, stage, stage_info));
177   PetscFunctionReturn(PETSC_SUCCESS);
178 }
179 
180 PETSC_INTERN PetscErrorCode PetscLogRegistryClassGetInfo(PetscLogRegistry registry, PetscLogClass clss, PetscLogClassInfo *class_info)
181 {
182   PetscFunctionBegin;
183   PetscCall(PetscLogClassArrayGet(registry->classes, clss, class_info));
184   PetscFunctionReturn(PETSC_SUCCESS);
185 }
186 
187 PETSC_INTERN PetscErrorCode PetscLogRegistryEventSetCollective(PetscLogRegistry registry, PetscLogEvent event, PetscBool collective)
188 {
189   PetscLogEventInfo *event_info;
190 
191   PetscFunctionBegin;
192   PetscCall(PetscLogEventArrayGetRef(registry->events, event, &event_info));
193   event_info->collective = collective;
194   PetscFunctionReturn(PETSC_SUCCESS);
195 }
196 
197 /* Given a list of strings on each process, create a global numbering.  Order
198  them by their order on the first process, then the remaining by their order
199  on the second process, etc.  The expectation is that most processes have the
200  same names in the same order so it shouldn't take too many rounds to figure
201  out */
202 
203 struct _n_PetscLogGlobalNames {
204   MPI_Comm     comm;
205   PetscInt     count_global;
206   PetscInt     count_local;
207   const char **names;
208   PetscInt    *global_to_local;
209   PetscInt    *local_to_global;
210 };
211 
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;
224 
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));
245 
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));
253 
254   p = 0;
255   while (p < size) {
256     PetscInt my_loc, next_loc;
257     PetscInt num_to_add;
258 
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     p          = next_loc;
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;
270 
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;
276 
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     }
282 
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; i < num_to_add; i++) {
297         PetscInt s;
298         PetscCallMPI(MPI_Bcast(str_buffer, max_name_len + 1, MPI_CHAR, p, comm));
299         PetscCall(PetscStrallocpy(str_buffer, &global_names[num_names_global]));
300         for (s = 0; s < num_names_local; s++) {
301           PetscBool same;
302 
303           if (local_name_seen[s]) continue;
304           PetscCall(PetscStrncmp(names[s], str_buffer, max_name_len + 1, &same));
305           if (same) {
306             local_name_seen[s]                            = PETSC_TRUE;
307             global_index_to_local_index[num_names_global] = s;
308             local_index_to_global_index[s]                = num_names_global;
309             num_names_local_remaining--;
310             break;
311           }
312         }
313         if (s == num_names_local) {
314           global_index_to_local_index[num_names_global] = -1; // this name is not present on this process
315         }
316         num_names_global++;
317       }
318     }
319   }
320 
321   PetscCall(PetscFree(str_buffer));
322   PetscCall(PetscFree(local_name_seen));
323   *num_names_global_p            = num_names_global;
324   *global_index_to_local_index_p = global_index_to_local_index;
325   *local_index_to_global_index_p = local_index_to_global_index;
326   *global_names_p                = (const char **)global_names;
327 
328   PetscFunctionReturn(PETSC_SUCCESS);
329 }
330 
331 PETSC_INTERN PetscErrorCode PetscLogGlobalNamesCreate(MPI_Comm comm, PetscInt num_names_local, const char **local_names, PetscLogGlobalNames *global_names_p)
332 {
333   PetscLogGlobalNames global_names;
334 
335   PetscFunctionBegin;
336   PetscCall(PetscNew(&global_names));
337   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));
338   global_names->count_local = num_names_local;
339   *global_names_p           = global_names;
340   PetscFunctionReturn(PETSC_SUCCESS);
341 }
342 
343 PETSC_INTERN PetscErrorCode PetscLogGlobalNamesDestroy(PetscLogGlobalNames *global_names_p)
344 {
345   PetscLogGlobalNames global_names;
346 
347   PetscFunctionBegin;
348   global_names    = *global_names_p;
349   *global_names_p = NULL;
350   PetscCall(PetscFree(global_names->global_to_local));
351   PetscCall(PetscFree(global_names->local_to_global));
352   for (PetscInt i = 0; i < global_names->count_global; i++) { PetscCall(PetscFree(global_names->names[i])); }
353   PetscCall(PetscFree(global_names->names));
354   PetscCall(PetscFree(global_names));
355   PetscFunctionReturn(PETSC_SUCCESS);
356 }
357 
358 PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGlobalGetName(PetscLogGlobalNames global_names, PetscInt idx, const char **name)
359 {
360   PetscFunctionBegin;
361   PetscCheck(idx >= 0 && idx < global_names->count_global, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Idx %d not in range [0,%d)", (int)idx, (int)global_names->count_global);
362   *name = global_names->names[idx];
363   PetscFunctionReturn(PETSC_SUCCESS);
364 }
365 
366 PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGlobalGetLocal(PetscLogGlobalNames global_names, PetscInt idx, PetscInt *local_idx)
367 {
368   PetscFunctionBegin;
369   PetscCheck(idx >= 0 && idx < global_names->count_global, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Idx %d not in range [0,%d)", (int)idx, (int)global_names->count_global);
370   *local_idx = global_names->global_to_local[idx];
371   PetscFunctionReturn(PETSC_SUCCESS);
372 }
373 
374 PETSC_INTERN PetscErrorCode PetscLogGlobalNamesLocalGetGlobal(PetscLogGlobalNames global_names, PetscInt local_idx, PetscInt *idx)
375 {
376   PetscFunctionBegin;
377   PetscCheck(local_idx >= 0 && local_idx < global_names->count_local, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Idx %d not in range [0,%d)", (int)local_idx, (int)global_names->count_local);
378   *idx = global_names->local_to_global[local_idx];
379   PetscFunctionReturn(PETSC_SUCCESS);
380 }
381 
382 PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGetSize(PetscLogGlobalNames global_names, PetscInt *local_size, PetscInt *global_size)
383 {
384   PetscFunctionBegin;
385   if (local_size) *local_size = global_names->count_local;
386   if (global_size) *global_size = global_names->count_global;
387   PetscFunctionReturn(PETSC_SUCCESS);
388 }
389 
390 PETSC_INTERN PetscErrorCode PetscLogRegistryCreateGlobalStageNames(MPI_Comm comm, PetscLogRegistry registry, PetscLogGlobalNames *global_names_p)
391 {
392   PetscInt     num_stages_local;
393   const char **names;
394 
395   PetscFunctionBegin;
396   PetscCall(PetscLogStageArrayGetSize(registry->stages, &num_stages_local, NULL));
397   PetscCall(PetscMalloc1(num_stages_local, &names));
398   for (PetscInt i = 0; i < num_stages_local; i++) {
399     PetscLogStageInfo stage_info = {NULL};
400     PetscCall(PetscLogRegistryStageGetInfo(registry, i, &stage_info));
401     names[i] = stage_info.name;
402   }
403   PetscCall(PetscLogGlobalNamesCreate(comm, num_stages_local, names, global_names_p));
404   PetscCall(PetscFree(names));
405   PetscFunctionReturn(PETSC_SUCCESS);
406 }
407 
408 PETSC_INTERN PetscErrorCode PetscLogRegistryCreateGlobalEventNames(MPI_Comm comm, PetscLogRegistry registry, PetscLogGlobalNames *global_names_p)
409 {
410   PetscInt     num_events_local;
411   const char **names;
412 
413   PetscFunctionBegin;
414   PetscCall(PetscLogEventArrayGetSize(registry->events, &num_events_local, NULL));
415   PetscCall(PetscMalloc1(num_events_local, &names));
416   for (PetscInt i = 0; i < num_events_local; i++) {
417     PetscLogEventInfo event_info = {NULL, 0, PETSC_FALSE};
418 
419     PetscCall(PetscLogRegistryEventGetInfo(registry, i, &event_info));
420     names[i] = event_info.name;
421   }
422   PetscCall(PetscLogGlobalNamesCreate(comm, num_events_local, names, global_names_p));
423   PetscCall(PetscFree(names));
424   PetscFunctionReturn(PETSC_SUCCESS);
425 }
426