xref: /petsc/src/sys/logging/state/logregistry.c (revision 7f031e8bf1f008cdd443cdad0cd45837cb20997c)
16873511fSToby Isaac #include <petsc/private/logimpl.h> /*I "petsclog.h" I*/
26873511fSToby Isaac 
36873511fSToby Isaac #define PETSC_LOG_RESIZABLE_ARRAY_HAS_NAME(Container, Entry, Key, Equal) \
46873511fSToby Isaac   static inline PETSC_UNUSED PetscErrorCode PetscLog##Container##Destructor(Entry *entry) \
56873511fSToby Isaac   { \
66873511fSToby Isaac     PetscFunctionBegin; \
76873511fSToby Isaac     PetscCall(PetscFree(entry->name)); \
86873511fSToby Isaac     PetscFunctionReturn(PETSC_SUCCESS); \
96873511fSToby Isaac   } \
106873511fSToby Isaac   PETSC_LOG_RESIZABLE_ARRAY(Container, Entry, Key, NULL, PetscLog##Container##Destructor, Equal)
116873511fSToby Isaac 
126873511fSToby Isaac #define PETSC_LOG_RESIZABLE_ARRAY_KEY_BY_NAME(Container, Entry) \
136873511fSToby Isaac   static inline PETSC_UNUSED PetscErrorCode PetscLog##Container##Equal(Entry *entry, const char *name, PetscBool *is_equal) \
146873511fSToby Isaac   { \
156873511fSToby Isaac     PetscFunctionBegin; \
166873511fSToby Isaac     PetscCall(PetscStrcmp(entry->name, name, is_equal)); \
176873511fSToby Isaac     PetscFunctionReturn(PETSC_SUCCESS); \
186873511fSToby Isaac   } \
196873511fSToby Isaac   PETSC_LOG_RESIZABLE_ARRAY_HAS_NAME(Container, Entry, const char *, PetscLog##Container##Equal)
206873511fSToby Isaac 
PetscLogClassArrayEqual(PetscLogClassInfo * class_info,PetscLogClassInfo * key,PetscBool * is_equal)216873511fSToby Isaac static PetscErrorCode PetscLogClassArrayEqual(PetscLogClassInfo *class_info, PetscLogClassInfo *key, PetscBool *is_equal)
226873511fSToby Isaac {
236873511fSToby Isaac   PetscFunctionBegin;
246873511fSToby Isaac   if (key->name) {
256873511fSToby Isaac     PetscCall(PetscStrcmp(class_info->name, key->name, is_equal));
266873511fSToby Isaac   } else {
276873511fSToby Isaac     *is_equal = (class_info->classid == key->classid) ? PETSC_TRUE : PETSC_FALSE;
286873511fSToby Isaac   }
296873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
306873511fSToby Isaac }
316873511fSToby Isaac 
326873511fSToby Isaac PETSC_LOG_RESIZABLE_ARRAY_KEY_BY_NAME(EventArray, PetscLogEventInfo)
336873511fSToby Isaac PETSC_LOG_RESIZABLE_ARRAY_KEY_BY_NAME(StageArray, PetscLogStageInfo)
346873511fSToby Isaac PETSC_LOG_RESIZABLE_ARRAY_HAS_NAME(ClassArray, PetscLogClassInfo, PetscLogClassInfo *, PetscLogClassArrayEqual)
356873511fSToby Isaac 
366873511fSToby Isaac struct _n_PetscLogRegistry {
376873511fSToby Isaac   PetscLogEventArray events;
386873511fSToby Isaac   PetscLogClassArray classes;
396873511fSToby Isaac   PetscLogStageArray stages;
406873511fSToby Isaac };
416873511fSToby Isaac 
PetscLogRegistryCreate(PetscLogRegistry * registry_p)426873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryCreate(PetscLogRegistry *registry_p)
436873511fSToby Isaac {
446873511fSToby Isaac   PetscLogRegistry registry;
456873511fSToby Isaac 
466873511fSToby Isaac   PetscFunctionBegin;
476873511fSToby Isaac   PetscCall(PetscNew(registry_p));
486873511fSToby Isaac   registry = *registry_p;
496873511fSToby Isaac   PetscCall(PetscLogEventArrayCreate(128, &registry->events));
506873511fSToby Isaac   PetscCall(PetscLogStageArrayCreate(8, &registry->stages));
516873511fSToby Isaac   PetscCall(PetscLogClassArrayCreate(128, &registry->classes));
526873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
536873511fSToby Isaac }
546873511fSToby Isaac 
PetscLogRegistryDestroy(PetscLogRegistry registry)556873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryDestroy(PetscLogRegistry registry)
566873511fSToby Isaac {
576873511fSToby Isaac   PetscFunctionBegin;
586873511fSToby Isaac   PetscCall(PetscLogEventArrayDestroy(&registry->events));
596873511fSToby Isaac   PetscCall(PetscLogClassArrayDestroy(&registry->classes));
606873511fSToby Isaac   PetscCall(PetscLogStageArrayDestroy(&registry->stages));
616873511fSToby Isaac   PetscCall(PetscFree(registry));
626873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
636873511fSToby Isaac }
646873511fSToby Isaac 
PetscLogRegistryGetNumEvents(PetscLogRegistry registry,PetscInt * num_events,PetscInt * max_events)656873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryGetNumEvents(PetscLogRegistry registry, PetscInt *num_events, PetscInt *max_events)
666873511fSToby Isaac {
676873511fSToby Isaac   PetscFunctionBegin;
686873511fSToby Isaac   PetscCall(PetscLogEventArrayGetSize(registry->events, num_events, max_events));
696873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
706873511fSToby Isaac }
716873511fSToby Isaac 
PetscLogRegistryGetNumStages(PetscLogRegistry registry,PetscInt * num_stages,PetscInt * max_stages)726873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryGetNumStages(PetscLogRegistry registry, PetscInt *num_stages, PetscInt *max_stages)
736873511fSToby Isaac {
746873511fSToby Isaac   PetscFunctionBegin;
756873511fSToby Isaac   PetscCall(PetscLogStageArrayGetSize(registry->stages, num_stages, max_stages));
766873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
776873511fSToby Isaac }
786873511fSToby Isaac 
PetscLogRegistryGetNumClasses(PetscLogRegistry registry,PetscInt * num_classes,PetscInt * max_classes)796873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryGetNumClasses(PetscLogRegistry registry, PetscInt *num_classes, PetscInt *max_classes)
806873511fSToby Isaac {
816873511fSToby Isaac   PetscFunctionBegin;
826873511fSToby Isaac   PetscCall(PetscLogClassArrayGetSize(registry->classes, num_classes, max_classes));
836873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
846873511fSToby Isaac }
856873511fSToby Isaac 
PetscLogRegistryStageRegister(PetscLogRegistry registry,const char sname[],int * stage)866873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryStageRegister(PetscLogRegistry registry, const char sname[], int *stage)
876873511fSToby Isaac {
886873511fSToby Isaac   int               idx;
896873511fSToby Isaac   PetscLogStageInfo stage_info;
906873511fSToby Isaac 
916873511fSToby Isaac   PetscFunctionBegin;
926873511fSToby Isaac   PetscCall(PetscLogStageArrayFind(registry->stages, sname, &idx));
936873511fSToby Isaac   PetscCheck(idx == -1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "An event named %s is already registered", sname);
946873511fSToby Isaac   *stage = registry->stages->num_entries;
956873511fSToby Isaac   PetscCall(PetscStrallocpy(sname, &stage_info.name));
966873511fSToby Isaac   PetscCall(PetscLogStageArrayPush(registry->stages, stage_info));
976873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
986873511fSToby Isaac }
996873511fSToby Isaac 
PetscLogRegistryEventRegister(PetscLogRegistry registry,const char name[],PetscClassId classid,PetscLogEvent * event)1006873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryEventRegister(PetscLogRegistry registry, const char name[], PetscClassId classid, PetscLogEvent *event)
1016873511fSToby Isaac {
1026873511fSToby Isaac   PetscLogEventInfo new_info;
1036873511fSToby Isaac 
1046873511fSToby Isaac   PetscFunctionBegin;
1056873511fSToby Isaac   PetscCall(PetscLogRegistryGetEventFromName(registry, name, event));
1066873511fSToby Isaac   if (*event >= 0) PetscFunctionReturn(PETSC_SUCCESS);
1076873511fSToby Isaac   *event              = registry->events->num_entries;
1086873511fSToby Isaac   new_info.classid    = classid;
1096873511fSToby Isaac   new_info.collective = PETSC_TRUE;
1106873511fSToby Isaac   PetscCall(PetscStrallocpy(name, &new_info.name));
1116873511fSToby Isaac   PetscCall(PetscLogEventArrayPush(registry->events, new_info));
1126873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1136873511fSToby Isaac }
1146873511fSToby Isaac 
PetscLogRegistryClassRegister(PetscLogRegistry registry,const char name[],PetscClassId classid,PetscLogClass * clss)1156873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryClassRegister(PetscLogRegistry registry, const char name[], PetscClassId classid, PetscLogClass *clss)
1166873511fSToby Isaac {
1176873511fSToby Isaac   PetscLogClassInfo new_info;
1186873511fSToby Isaac 
1196873511fSToby Isaac   PetscFunctionBegin;
1206873511fSToby Isaac   PetscCall(PetscLogRegistryGetClassFromClassId(registry, classid, clss));
1216873511fSToby Isaac   if (*clss >= 0) PetscFunctionReturn(PETSC_SUCCESS);
1226873511fSToby Isaac   *clss            = registry->classes->num_entries;
1236873511fSToby Isaac   new_info.classid = classid;
1246873511fSToby Isaac   PetscCall(PetscStrallocpy(name, &new_info.name));
1256873511fSToby Isaac   PetscCall(PetscLogClassArrayPush(registry->classes, new_info));
1266873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1276873511fSToby Isaac }
1286873511fSToby Isaac 
PetscLogRegistryGetEventFromName(PetscLogRegistry registry,const char name[],PetscLogEvent * event)1296873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryGetEventFromName(PetscLogRegistry registry, const char name[], PetscLogEvent *event)
1306873511fSToby Isaac {
1316873511fSToby Isaac   PetscFunctionBegin;
1326873511fSToby Isaac   PetscCall(PetscLogEventArrayFind(registry->events, name, event));
1336873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1346873511fSToby Isaac }
1356873511fSToby Isaac 
PetscLogRegistryGetStageFromName(PetscLogRegistry registry,const char name[],PetscLogStage * stage)1366873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryGetStageFromName(PetscLogRegistry registry, const char name[], PetscLogStage *stage)
1376873511fSToby Isaac {
1386873511fSToby Isaac   PetscFunctionBegin;
1396873511fSToby Isaac   PetscCall(PetscLogStageArrayFind(registry->stages, name, stage));
1406873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1416873511fSToby Isaac }
1426873511fSToby Isaac 
PetscLogRegistryGetClassFromClassId(PetscLogRegistry registry,PetscClassId classid,PetscLogStage * clss)1436873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryGetClassFromClassId(PetscLogRegistry registry, PetscClassId classid, PetscLogStage *clss)
1446873511fSToby Isaac {
1456873511fSToby Isaac   PetscLogClassInfo key;
1466873511fSToby Isaac 
1476873511fSToby Isaac   PetscFunctionBegin;
1486873511fSToby Isaac   key.name    = NULL;
1496873511fSToby Isaac   key.classid = classid;
1506873511fSToby Isaac   PetscCall(PetscLogClassArrayFind(registry->classes, &key, clss));
1516873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1526873511fSToby Isaac }
1536873511fSToby Isaac 
PetscLogRegistryGetClassFromName(PetscLogRegistry registry,const char name[],PetscLogStage * clss)1546873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryGetClassFromName(PetscLogRegistry registry, const char name[], PetscLogStage *clss)
1556873511fSToby Isaac {
1566873511fSToby Isaac   PetscLogClassInfo key;
1576873511fSToby Isaac 
1586873511fSToby Isaac   PetscFunctionBegin;
1596873511fSToby Isaac   key.name    = (char *)name;
1606873511fSToby Isaac   key.classid = -1;
1616873511fSToby Isaac   PetscCall(PetscLogClassArrayFind(registry->classes, &key, clss));
1626873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1636873511fSToby Isaac }
1646873511fSToby Isaac 
PetscLogRegistryEventGetInfo(PetscLogRegistry registry,PetscLogEvent event,PetscLogEventInfo * event_info)1656873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryEventGetInfo(PetscLogRegistry registry, PetscLogEvent event, PetscLogEventInfo *event_info)
1666873511fSToby Isaac {
1676873511fSToby Isaac   PetscFunctionBegin;
1686873511fSToby Isaac   PetscCall(PetscLogEventArrayGet(registry->events, event, event_info));
1696873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1706873511fSToby Isaac }
1716873511fSToby Isaac 
PetscLogRegistryStageGetInfo(PetscLogRegistry registry,PetscLogStage stage,PetscLogStageInfo * stage_info)1726873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryStageGetInfo(PetscLogRegistry registry, PetscLogStage stage, PetscLogStageInfo *stage_info)
1736873511fSToby Isaac {
1746873511fSToby Isaac   PetscFunctionBegin;
1756873511fSToby Isaac   PetscCall(PetscLogStageArrayGet(registry->stages, stage, stage_info));
1766873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1776873511fSToby Isaac }
1786873511fSToby Isaac 
PetscLogRegistryClassGetInfo(PetscLogRegistry registry,PetscLogClass clss,PetscLogClassInfo * class_info)1796873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryClassGetInfo(PetscLogRegistry registry, PetscLogClass clss, PetscLogClassInfo *class_info)
1806873511fSToby Isaac {
1816873511fSToby Isaac   PetscFunctionBegin;
1826873511fSToby Isaac   PetscCall(PetscLogClassArrayGet(registry->classes, clss, class_info));
1836873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1846873511fSToby Isaac }
1856873511fSToby Isaac 
PetscLogRegistryEventSetCollective(PetscLogRegistry registry,PetscLogEvent event,PetscBool collective)1866873511fSToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryEventSetCollective(PetscLogRegistry registry, PetscLogEvent event, PetscBool collective)
1876873511fSToby Isaac {
1886873511fSToby Isaac   PetscLogEventInfo *event_info;
1896873511fSToby Isaac 
1906873511fSToby Isaac   PetscFunctionBegin;
1916873511fSToby Isaac   PetscCall(PetscLogEventArrayGetRef(registry->events, event, &event_info));
1926873511fSToby Isaac   event_info->collective = collective;
1936873511fSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1946873511fSToby Isaac }
195176456b8SToby Isaac 
196176456b8SToby Isaac /* Given a list of strings on each process, create a global numbering.  Order
197176456b8SToby Isaac  them by their order on the first process, then the remaining by their order
198176456b8SToby Isaac  on the second process, etc.  The expectation is that most processes have the
199176456b8SToby Isaac  same names in the same order so it shouldn't take too many rounds to figure
200176456b8SToby Isaac  out */
201176456b8SToby Isaac 
202176456b8SToby Isaac struct _n_PetscLogGlobalNames {
203176456b8SToby Isaac   MPI_Comm     comm;
204176456b8SToby Isaac   PetscInt     count_global;
205176456b8SToby Isaac   PetscInt     count_local;
206176456b8SToby Isaac   const char **names;
207176456b8SToby Isaac   PetscInt    *global_to_local;
208176456b8SToby Isaac   PetscInt    *local_to_global;
209176456b8SToby Isaac };
210176456b8SToby Isaac 
211*31a765c4SPierre Jolivet PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wconversion")
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)212176456b8SToby Isaac 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)
213176456b8SToby Isaac {
214176456b8SToby Isaac   PetscMPIInt size, rank;
215176456b8SToby Isaac   PetscInt    num_names_global          = 0;
216176456b8SToby Isaac   PetscInt    num_names_local_remaining = num_names_local;
217176456b8SToby Isaac   PetscBool  *local_name_seen;
218176456b8SToby Isaac   PetscInt   *global_index_to_local_index = NULL;
219176456b8SToby Isaac   PetscInt   *local_index_to_global_index = NULL;
220176456b8SToby Isaac   PetscInt    max_name_len                = 0;
221176456b8SToby Isaac   char       *str_buffer;
222176456b8SToby Isaac   char      **global_names = NULL;
223176456b8SToby Isaac   PetscMPIInt p;
224176456b8SToby Isaac 
225176456b8SToby Isaac   PetscFunctionBegin;
226176456b8SToby Isaac   PetscCallMPI(MPI_Comm_size(comm, &size));
227176456b8SToby Isaac   if (size == 1) {
228176456b8SToby Isaac     PetscCall(PetscMalloc1(num_names_local, &global_index_to_local_index));
229176456b8SToby Isaac     PetscCall(PetscMalloc1(num_names_local, &local_index_to_global_index));
230176456b8SToby Isaac     PetscCall(PetscMalloc1(num_names_local, &global_names));
231176456b8SToby Isaac     for (PetscInt i = 0; i < num_names_local; i++) {
232176456b8SToby Isaac       global_index_to_local_index[i] = i;
233176456b8SToby Isaac       local_index_to_global_index[i] = i;
234176456b8SToby Isaac       PetscCall(PetscStrallocpy(names[i], &global_names[i]));
235176456b8SToby Isaac     }
236176456b8SToby Isaac     *num_names_global_p            = num_names_local;
237176456b8SToby Isaac     *global_index_to_local_index_p = global_index_to_local_index;
238176456b8SToby Isaac     *local_index_to_global_index_p = local_index_to_global_index;
239176456b8SToby Isaac     *global_names_p                = (const char **)global_names;
240176456b8SToby Isaac     PetscFunctionReturn(PETSC_SUCCESS);
241176456b8SToby Isaac   }
242176456b8SToby Isaac   PetscCallMPI(MPI_Comm_rank(comm, &rank));
243176456b8SToby Isaac   PetscCall(PetscCalloc1(num_names_local, &local_name_seen));
244176456b8SToby Isaac   PetscCall(PetscMalloc1(num_names_local, &local_index_to_global_index));
245176456b8SToby Isaac 
246176456b8SToby Isaac   for (PetscInt i = 0; i < num_names_local; i++) {
247176456b8SToby Isaac     size_t i_len;
248176456b8SToby Isaac     PetscCall(PetscStrlen(names[i], &i_len));
249176456b8SToby Isaac     max_name_len = PetscMax(max_name_len, (PetscInt)i_len);
250176456b8SToby Isaac   }
251176456b8SToby Isaac   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &max_name_len, 1, MPIU_INT, MPI_MAX, comm));
252176456b8SToby Isaac   PetscCall(PetscCalloc1(max_name_len + 1, &str_buffer));
253176456b8SToby Isaac 
254176456b8SToby Isaac   p = 0;
255176456b8SToby Isaac   while (p < size) {
256176456b8SToby Isaac     PetscInt my_loc, next_loc;
257176456b8SToby Isaac     PetscInt num_to_add;
258176456b8SToby Isaac 
259176456b8SToby Isaac     my_loc = num_names_local_remaining > 0 ? rank : PETSC_MPI_INT_MAX;
260176456b8SToby Isaac     PetscCallMPI(MPIU_Allreduce(&my_loc, &next_loc, 1, MPIU_INT, MPI_MIN, comm));
261176456b8SToby Isaac     if (next_loc == PETSC_MPI_INT_MAX) break;
262176456b8SToby Isaac     PetscAssert(next_loc >= p, comm, PETSC_ERR_PLIB, "Failed invariant, expected increasing next process");
263*31a765c4SPierre Jolivet     PetscCall(PetscCIntCast(next_loc, &p));
264176456b8SToby Isaac     num_to_add = (rank == p) ? num_names_local_remaining : -1;
265176456b8SToby Isaac     PetscCallMPI(MPI_Bcast(&num_to_add, 1, MPIU_INT, p, comm));
266176456b8SToby Isaac     {
267176456b8SToby Isaac       PetscInt  new_num_names_global = num_names_global + num_to_add;
268176456b8SToby Isaac       PetscInt *new_global_index_to_local_index;
269176456b8SToby Isaac       char    **new_global_names;
270176456b8SToby Isaac 
271176456b8SToby Isaac       PetscCall(PetscMalloc1(new_num_names_global, &new_global_index_to_local_index));
272176456b8SToby Isaac       PetscCall(PetscArraycpy(new_global_index_to_local_index, global_index_to_local_index, num_names_global));
273176456b8SToby Isaac       for (PetscInt i = num_names_global; i < new_num_names_global; i++) new_global_index_to_local_index[i] = -1;
274176456b8SToby Isaac       PetscCall(PetscFree(global_index_to_local_index));
275176456b8SToby Isaac       global_index_to_local_index = new_global_index_to_local_index;
276176456b8SToby Isaac 
277176456b8SToby Isaac       PetscCall(PetscCalloc1(new_num_names_global, &new_global_names));
278176456b8SToby Isaac       PetscCall(PetscArraycpy(new_global_names, global_names, num_names_global));
279176456b8SToby Isaac       PetscCall(PetscFree(global_names));
280176456b8SToby Isaac       global_names = new_global_names;
281176456b8SToby Isaac     }
282176456b8SToby Isaac 
283176456b8SToby Isaac     if (rank == p) {
284176456b8SToby Isaac       for (PetscInt s = 0; s < num_names_local; s++) {
285176456b8SToby Isaac         if (local_name_seen[s]) continue;
286176456b8SToby Isaac         local_name_seen[s] = PETSC_TRUE;
287176456b8SToby Isaac         PetscCall(PetscArrayzero(str_buffer, max_name_len + 1));
288176456b8SToby Isaac         PetscCall(PetscStrallocpy(names[s], &global_names[num_names_global]));
289176456b8SToby Isaac         PetscCall(PetscStrncpy(str_buffer, names[s], max_name_len + 1));
290176456b8SToby Isaac         PetscCallMPI(MPI_Bcast(str_buffer, max_name_len + 1, MPI_CHAR, p, comm));
291176456b8SToby Isaac         local_index_to_global_index[s]                  = num_names_global;
292176456b8SToby Isaac         global_index_to_local_index[num_names_global++] = s;
293176456b8SToby Isaac         num_names_local_remaining--;
294176456b8SToby Isaac       }
295176456b8SToby Isaac     } else {
296*31a765c4SPierre Jolivet       for (PetscInt i = 0, s; i < num_to_add; i++) {
297176456b8SToby Isaac         PetscCallMPI(MPI_Bcast(str_buffer, max_name_len + 1, MPI_CHAR, p, comm));
298176456b8SToby Isaac         PetscCall(PetscStrallocpy(str_buffer, &global_names[num_names_global]));
299176456b8SToby Isaac         for (s = 0; s < num_names_local; s++) {
300176456b8SToby Isaac           PetscBool same;
301176456b8SToby Isaac 
302176456b8SToby Isaac           if (local_name_seen[s]) continue;
303176456b8SToby Isaac           PetscCall(PetscStrncmp(names[s], str_buffer, max_name_len + 1, &same));
304176456b8SToby Isaac           if (same) {
305176456b8SToby Isaac             local_name_seen[s]                            = PETSC_TRUE;
306176456b8SToby Isaac             global_index_to_local_index[num_names_global] = s;
307176456b8SToby Isaac             local_index_to_global_index[s]                = num_names_global;
308176456b8SToby Isaac             num_names_local_remaining--;
309176456b8SToby Isaac             break;
310176456b8SToby Isaac           }
311176456b8SToby Isaac         }
3123a7d0413SPierre Jolivet         if (s == num_names_local) global_index_to_local_index[num_names_global] = -1; // this name is not present on this process
313176456b8SToby Isaac         num_names_global++;
314176456b8SToby Isaac       }
315176456b8SToby Isaac     }
316176456b8SToby Isaac   }
317176456b8SToby Isaac 
318176456b8SToby Isaac   PetscCall(PetscFree(str_buffer));
319176456b8SToby Isaac   PetscCall(PetscFree(local_name_seen));
320176456b8SToby Isaac   *num_names_global_p            = num_names_global;
321176456b8SToby Isaac   *global_index_to_local_index_p = global_index_to_local_index;
322176456b8SToby Isaac   *local_index_to_global_index_p = local_index_to_global_index;
323176456b8SToby Isaac   *global_names_p                = (const char **)global_names;
324176456b8SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
325176456b8SToby Isaac }
PETSC_PRAGMA_DIAGNOSTIC_IGNORED_END()326*31a765c4SPierre Jolivet PETSC_PRAGMA_DIAGNOSTIC_IGNORED_END()
327176456b8SToby Isaac 
328176456b8SToby Isaac PETSC_INTERN PetscErrorCode PetscLogGlobalNamesCreate(MPI_Comm comm, PetscInt num_names_local, const char **local_names, PetscLogGlobalNames *global_names_p)
329176456b8SToby Isaac {
330176456b8SToby Isaac   PetscLogGlobalNames global_names;
331176456b8SToby Isaac 
332176456b8SToby Isaac   PetscFunctionBegin;
333176456b8SToby Isaac   PetscCall(PetscNew(&global_names));
334176456b8SToby Isaac   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));
335176456b8SToby Isaac   global_names->count_local = num_names_local;
336176456b8SToby Isaac   *global_names_p           = global_names;
337176456b8SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
338176456b8SToby Isaac }
339176456b8SToby Isaac 
PetscLogGlobalNamesDestroy(PetscLogGlobalNames * global_names_p)340176456b8SToby Isaac PETSC_INTERN PetscErrorCode PetscLogGlobalNamesDestroy(PetscLogGlobalNames *global_names_p)
341176456b8SToby Isaac {
342176456b8SToby Isaac   PetscLogGlobalNames global_names;
343176456b8SToby Isaac 
344176456b8SToby Isaac   PetscFunctionBegin;
345176456b8SToby Isaac   global_names    = *global_names_p;
346176456b8SToby Isaac   *global_names_p = NULL;
347176456b8SToby Isaac   PetscCall(PetscFree(global_names->global_to_local));
348176456b8SToby Isaac   PetscCall(PetscFree(global_names->local_to_global));
3493a7d0413SPierre Jolivet   for (PetscInt i = 0; i < global_names->count_global; i++) PetscCall(PetscFree(global_names->names[i]));
350176456b8SToby Isaac   PetscCall(PetscFree(global_names->names));
351176456b8SToby Isaac   PetscCall(PetscFree(global_names));
352176456b8SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
353176456b8SToby Isaac }
354176456b8SToby Isaac 
PetscLogGlobalNamesGlobalGetName(PetscLogGlobalNames global_names,PetscInt idx,const char ** name)355176456b8SToby Isaac PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGlobalGetName(PetscLogGlobalNames global_names, PetscInt idx, const char **name)
356176456b8SToby Isaac {
357176456b8SToby Isaac   PetscFunctionBegin;
358835f2295SStefano Zampini   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);
359176456b8SToby Isaac   *name = global_names->names[idx];
360176456b8SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
361176456b8SToby Isaac }
362176456b8SToby Isaac 
PetscLogGlobalNamesGlobalGetLocal(PetscLogGlobalNames global_names,PetscInt idx,PetscInt * local_idx)363176456b8SToby Isaac PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGlobalGetLocal(PetscLogGlobalNames global_names, PetscInt idx, PetscInt *local_idx)
364176456b8SToby Isaac {
365176456b8SToby Isaac   PetscFunctionBegin;
366835f2295SStefano Zampini   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);
367176456b8SToby Isaac   *local_idx = global_names->global_to_local[idx];
368176456b8SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
369176456b8SToby Isaac }
370176456b8SToby Isaac 
PetscLogGlobalNamesLocalGetGlobal(PetscLogGlobalNames global_names,PetscInt local_idx,PetscInt * idx)371176456b8SToby Isaac PETSC_INTERN PetscErrorCode PetscLogGlobalNamesLocalGetGlobal(PetscLogGlobalNames global_names, PetscInt local_idx, PetscInt *idx)
372176456b8SToby Isaac {
373176456b8SToby Isaac   PetscFunctionBegin;
374835f2295SStefano Zampini   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);
375176456b8SToby Isaac   *idx = global_names->local_to_global[local_idx];
376176456b8SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
377176456b8SToby Isaac }
378176456b8SToby Isaac 
PetscLogGlobalNamesGetSize(PetscLogGlobalNames global_names,PetscInt * local_size,PetscInt * global_size)379176456b8SToby Isaac PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGetSize(PetscLogGlobalNames global_names, PetscInt *local_size, PetscInt *global_size)
380176456b8SToby Isaac {
381176456b8SToby Isaac   PetscFunctionBegin;
382176456b8SToby Isaac   if (local_size) *local_size = global_names->count_local;
383176456b8SToby Isaac   if (global_size) *global_size = global_names->count_global;
384176456b8SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
385176456b8SToby Isaac }
386176456b8SToby Isaac 
PetscLogRegistryCreateGlobalStageNames(MPI_Comm comm,PetscLogRegistry registry,PetscLogGlobalNames * global_names_p)387176456b8SToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryCreateGlobalStageNames(MPI_Comm comm, PetscLogRegistry registry, PetscLogGlobalNames *global_names_p)
388176456b8SToby Isaac {
389176456b8SToby Isaac   PetscInt     num_stages_local;
390176456b8SToby Isaac   const char **names;
391176456b8SToby Isaac 
392176456b8SToby Isaac   PetscFunctionBegin;
393176456b8SToby Isaac   PetscCall(PetscLogStageArrayGetSize(registry->stages, &num_stages_local, NULL));
394176456b8SToby Isaac   PetscCall(PetscMalloc1(num_stages_local, &names));
395*31a765c4SPierre Jolivet   for (PetscLogEvent i = 0; i < num_stages_local; i++) {
396176456b8SToby Isaac     PetscLogStageInfo stage_info = {NULL};
397176456b8SToby Isaac     PetscCall(PetscLogRegistryStageGetInfo(registry, i, &stage_info));
398176456b8SToby Isaac     names[i] = stage_info.name;
399176456b8SToby Isaac   }
400176456b8SToby Isaac   PetscCall(PetscLogGlobalNamesCreate(comm, num_stages_local, names, global_names_p));
401176456b8SToby Isaac   PetscCall(PetscFree(names));
402176456b8SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
403176456b8SToby Isaac }
404176456b8SToby Isaac 
PetscLogRegistryCreateGlobalEventNames(MPI_Comm comm,PetscLogRegistry registry,PetscLogGlobalNames * global_names_p)405176456b8SToby Isaac PETSC_INTERN PetscErrorCode PetscLogRegistryCreateGlobalEventNames(MPI_Comm comm, PetscLogRegistry registry, PetscLogGlobalNames *global_names_p)
406176456b8SToby Isaac {
407176456b8SToby Isaac   PetscInt     num_events_local;
408176456b8SToby Isaac   const char **names;
409176456b8SToby Isaac 
410176456b8SToby Isaac   PetscFunctionBegin;
411176456b8SToby Isaac   PetscCall(PetscLogEventArrayGetSize(registry->events, &num_events_local, NULL));
412176456b8SToby Isaac   PetscCall(PetscMalloc1(num_events_local, &names));
413*31a765c4SPierre Jolivet   for (PetscLogEvent i = 0; i < num_events_local; i++) {
414176456b8SToby Isaac     PetscLogEventInfo event_info = {NULL, 0, PETSC_FALSE};
415176456b8SToby Isaac 
416176456b8SToby Isaac     PetscCall(PetscLogRegistryEventGetInfo(registry, i, &event_info));
417176456b8SToby Isaac     names[i] = event_info.name;
418176456b8SToby Isaac   }
419176456b8SToby Isaac   PetscCall(PetscLogGlobalNamesCreate(comm, num_events_local, names, global_names_p));
420176456b8SToby Isaac   PetscCall(PetscFree(names));
421176456b8SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
422176456b8SToby Isaac }
423