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, &registry->events));
 50:   PetscCall(PetscLogStageArrayCreate(8, &registry->stages));
 51:   PetscCall(PetscLogClassArrayCreate(128, &registry->classes));
 52:   PetscFunctionReturn(PETSC_SUCCESS);
 53: }

 55: PETSC_INTERN PetscErrorCode PetscLogRegistryDestroy(PetscLogRegistry registry)
 56: {
 57:   PetscFunctionBegin;
 58:   PetscCall(PetscLogEventArrayDestroy(&registry->events));
 59:   PetscCall(PetscLogClassArrayDestroy(&registry->classes));
 60:   PetscCall(PetscLogStageArrayDestroy(&registry->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: 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)
212: {
213:   PetscMPIInt size, rank;
214:   PetscInt    num_names_global          = 0;
215:   PetscInt    num_names_local_remaining = num_names_local;
216:   PetscBool  *local_name_seen;
217:   PetscInt   *global_index_to_local_index = NULL;
218:   PetscInt   *local_index_to_global_index = NULL;
219:   PetscInt    max_name_len                = 0;
220:   char       *str_buffer;
221:   char      **global_names = NULL;
222:   PetscMPIInt p;

224:   PetscFunctionBegin;
225:   PetscCallMPI(MPI_Comm_size(comm, &size));
226:   if (size == 1) {
227:     PetscCall(PetscMalloc1(num_names_local, &global_index_to_local_index));
228:     PetscCall(PetscMalloc1(num_names_local, &local_index_to_global_index));
229:     PetscCall(PetscMalloc1(num_names_local, &global_names));
230:     for (PetscInt i = 0; i < num_names_local; i++) {
231:       global_index_to_local_index[i] = i;
232:       local_index_to_global_index[i] = i;
233:       PetscCall(PetscStrallocpy(names[i], &global_names[i]));
234:     }
235:     *num_names_global_p            = num_names_local;
236:     *global_index_to_local_index_p = global_index_to_local_index;
237:     *local_index_to_global_index_p = local_index_to_global_index;
238:     *global_names_p                = (const char **)global_names;
239:     PetscFunctionReturn(PETSC_SUCCESS);
240:   }
241:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
242:   PetscCall(PetscCalloc1(num_names_local, &local_name_seen));
243:   PetscCall(PetscMalloc1(num_names_local, &local_index_to_global_index));

245:   for (PetscInt i = 0; i < num_names_local; i++) {
246:     size_t i_len;
247:     PetscCall(PetscStrlen(names[i], &i_len));
248:     max_name_len = PetscMax(max_name_len, (PetscInt)i_len);
249:   }
250:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &max_name_len, 1, MPIU_INT, MPI_MAX, comm));
251:   PetscCall(PetscCalloc1(max_name_len + 1, &str_buffer));

253:   p = 0;
254:   while (p < size) {
255:     PetscInt my_loc, next_loc;
256:     PetscInt num_to_add;

258:     my_loc = num_names_local_remaining > 0 ? rank : PETSC_MPI_INT_MAX;
259:     PetscCallMPI(MPIU_Allreduce(&my_loc, &next_loc, 1, MPIU_INT, MPI_MIN, comm));
260:     if (next_loc == PETSC_MPI_INT_MAX) break;
261:     PetscAssert(next_loc >= p, comm, PETSC_ERR_PLIB, "Failed invariant, expected increasing next process");
262:     p          = next_loc;
263:     num_to_add = (rank == p) ? num_names_local_remaining : -1;
264:     PetscCallMPI(MPI_Bcast(&num_to_add, 1, MPIU_INT, p, comm));
265:     {
266:       PetscInt  new_num_names_global = num_names_global + num_to_add;
267:       PetscInt *new_global_index_to_local_index;
268:       char    **new_global_names;

270:       PetscCall(PetscMalloc1(new_num_names_global, &new_global_index_to_local_index));
271:       PetscCall(PetscArraycpy(new_global_index_to_local_index, global_index_to_local_index, num_names_global));
272:       for (PetscInt i = num_names_global; i < new_num_names_global; i++) new_global_index_to_local_index[i] = -1;
273:       PetscCall(PetscFree(global_index_to_local_index));
274:       global_index_to_local_index = new_global_index_to_local_index;

276:       PetscCall(PetscCalloc1(new_num_names_global, &new_global_names));
277:       PetscCall(PetscArraycpy(new_global_names, global_names, num_names_global));
278:       PetscCall(PetscFree(global_names));
279:       global_names = new_global_names;
280:     }

282:     if (rank == p) {
283:       for (PetscInt s = 0; s < num_names_local; s++) {
284:         if (local_name_seen[s]) continue;
285:         local_name_seen[s] = PETSC_TRUE;
286:         PetscCall(PetscArrayzero(str_buffer, max_name_len + 1));
287:         PetscCall(PetscStrallocpy(names[s], &global_names[num_names_global]));
288:         PetscCall(PetscStrncpy(str_buffer, names[s], max_name_len + 1));
289:         PetscCallMPI(MPI_Bcast(str_buffer, max_name_len + 1, MPI_CHAR, p, comm));
290:         local_index_to_global_index[s]                  = num_names_global;
291:         global_index_to_local_index[num_names_global++] = s;
292:         num_names_local_remaining--;
293:       }
294:     } else {
295:       for (PetscInt i = 0; i < num_to_add; i++) {
296:         PetscInt s;
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) {
313:           global_index_to_local_index[num_names_global] = -1; // this name is not present on this process
314:         }
315:         num_names_global++;
316:       }
317:     }
318:   }

320:   PetscCall(PetscFree(str_buffer));
321:   PetscCall(PetscFree(local_name_seen));
322:   *num_names_global_p            = num_names_global;
323:   *global_index_to_local_index_p = global_index_to_local_index;
324:   *local_index_to_global_index_p = local_index_to_global_index;
325:   *global_names_p                = (const char **)global_names;
326:   PetscFunctionReturn(PETSC_SUCCESS);
327: }

329: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesCreate(MPI_Comm comm, PetscInt num_names_local, const char **local_names, PetscLogGlobalNames *global_names_p)
330: {
331:   PetscLogGlobalNames global_names;

333:   PetscFunctionBegin;
334:   PetscCall(PetscNew(&global_names));
335:   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));
336:   global_names->count_local = num_names_local;
337:   *global_names_p           = global_names;
338:   PetscFunctionReturn(PETSC_SUCCESS);
339: }

341: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesDestroy(PetscLogGlobalNames *global_names_p)
342: {
343:   PetscLogGlobalNames global_names;

345:   PetscFunctionBegin;
346:   global_names    = *global_names_p;
347:   *global_names_p = NULL;
348:   PetscCall(PetscFree(global_names->global_to_local));
349:   PetscCall(PetscFree(global_names->local_to_global));
350:   for (PetscInt i = 0; i < global_names->count_global; i++) { PetscCall(PetscFree(global_names->names[i])); }
351:   PetscCall(PetscFree(global_names->names));
352:   PetscCall(PetscFree(global_names));
353:   PetscFunctionReturn(PETSC_SUCCESS);
354: }

356: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGlobalGetName(PetscLogGlobalNames global_names, PetscInt idx, const char **name)
357: {
358:   PetscFunctionBegin;
359:   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);
360:   *name = global_names->names[idx];
361:   PetscFunctionReturn(PETSC_SUCCESS);
362: }

364: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGlobalGetLocal(PetscLogGlobalNames global_names, PetscInt idx, PetscInt *local_idx)
365: {
366:   PetscFunctionBegin;
367:   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);
368:   *local_idx = global_names->global_to_local[idx];
369:   PetscFunctionReturn(PETSC_SUCCESS);
370: }

372: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesLocalGetGlobal(PetscLogGlobalNames global_names, PetscInt local_idx, PetscInt *idx)
373: {
374:   PetscFunctionBegin;
375:   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);
376:   *idx = global_names->local_to_global[local_idx];
377:   PetscFunctionReturn(PETSC_SUCCESS);
378: }

380: PETSC_INTERN PetscErrorCode PetscLogGlobalNamesGetSize(PetscLogGlobalNames global_names, PetscInt *local_size, PetscInt *global_size)
381: {
382:   PetscFunctionBegin;
383:   if (local_size) *local_size = global_names->count_local;
384:   if (global_size) *global_size = global_names->count_global;
385:   PetscFunctionReturn(PETSC_SUCCESS);
386: }

388: PETSC_INTERN PetscErrorCode PetscLogRegistryCreateGlobalStageNames(MPI_Comm comm, PetscLogRegistry registry, PetscLogGlobalNames *global_names_p)
389: {
390:   PetscInt     num_stages_local;
391:   const char **names;

393:   PetscFunctionBegin;
394:   PetscCall(PetscLogStageArrayGetSize(registry->stages, &num_stages_local, NULL));
395:   PetscCall(PetscMalloc1(num_stages_local, &names));
396:   for (PetscInt i = 0; i < num_stages_local; i++) {
397:     PetscLogStageInfo stage_info = {NULL};
398:     PetscCall(PetscLogRegistryStageGetInfo(registry, i, &stage_info));
399:     names[i] = stage_info.name;
400:   }
401:   PetscCall(PetscLogGlobalNamesCreate(comm, num_stages_local, names, global_names_p));
402:   PetscCall(PetscFree(names));
403:   PetscFunctionReturn(PETSC_SUCCESS);
404: }

406: PETSC_INTERN PetscErrorCode PetscLogRegistryCreateGlobalEventNames(MPI_Comm comm, PetscLogRegistry registry, PetscLogGlobalNames *global_names_p)
407: {
408:   PetscInt     num_events_local;
409:   const char **names;

411:   PetscFunctionBegin;
412:   PetscCall(PetscLogEventArrayGetSize(registry->events, &num_events_local, NULL));
413:   PetscCall(PetscMalloc1(num_events_local, &names));
414:   for (PetscInt i = 0; i < num_events_local; i++) {
415:     PetscLogEventInfo event_info = {NULL, 0, PETSC_FALSE};

417:     PetscCall(PetscLogRegistryEventGetInfo(registry, i, &event_info));
418:     names[i] = event_info.name;
419:   }
420:   PetscCall(PetscLogGlobalNamesCreate(comm, num_events_local, names, global_names_p));
421:   PetscCall(PetscFree(names));
422:   PetscFunctionReturn(PETSC_SUCCESS);
423: }