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, ®istry->events)); 51 PetscCall(PetscLogStageArrayCreate(8, ®istry->stages)); 52 PetscCall(PetscLogClassArrayCreate(128, ®istry->classes)); 53 PetscFunctionReturn(PETSC_SUCCESS); 54 } 55 56 PETSC_INTERN PetscErrorCode PetscLogRegistryDestroy(PetscLogRegistry registry) 57 { 58 PetscFunctionBegin; 59 PetscCall(PetscLogEventArrayDestroy(®istry->events)); 60 PetscCall(PetscLogClassArrayDestroy(®istry->classes)); 61 PetscCall(PetscLogStageArrayDestroy(®istry->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