Actual source code: ex48.c
1: static char help[] = "Tests HDF5 attribute I/O.\n\n";
3: #include <petscviewerhdf5.h>
4: #include <petscvec.h>
6: static PetscInt n = 5; /* testing vector size */
7: static PetscBool verbose = PETSC_FALSE;
8: #define SLEN 128
10: /* sequence of unique absolute paths */
11: #define nap 9
12: static const char *apaths[nap] = {
13: /* 0 */
14: "/", "/g1", "/g1/g2", "/g1/nonExistingGroup1", "/g1/g3",
15: /* 5 */
16: "/g1/g3/g4", "/g1/nonExistingGroup2", "/g1/nonExistingGroup2/g5", "/g1/g6/g7"};
18: #define np 21
19: /* sequence of paths (absolute or relative); "<" encodes Pop */
20: static const char *paths[np] = {
21: /* 0 */
22: "/",
23: "/g1",
24: "/g1/g2",
25: "/g1/nonExistingGroup1",
26: "<",
27: /* 5 */
28: ".", /* /g1/g2 */
29: "<",
30: "<",
31: "g3", /* /g1/g3 */
32: "g4", /* /g1/g3/g4 */
33: /* 10 */
34: "<",
35: "<",
36: ".", /* /g1 */
37: "<",
38: "nonExistingGroup2", /* /g1/nonExistingG2 */
39: /* 15 */
40: "g5", /* /g1/nonExistingG2/g5 */
41: "<",
42: "<",
43: "g6/g7", /* /g1/g6/g7 */
44: "<",
45: /* 20 */
46: "<",
47: };
48: /* corresponding expected absolute paths - positions in abspath */
49: static const PetscInt paths2apaths[np] = {
50: /* 0 */
51: 0,
52: 1,
53: 2,
54: 3,
55: 2,
56: /* 5 */
57: 2,
58: 2,
59: 1,
60: 4,
61: 5,
62: /* 10 */
63: 4,
64: 1,
65: 1,
66: 1,
67: 6,
68: /* 15 */
69: 7,
70: 6,
71: 1,
72: 8,
73: 1,
74: /* 20 */
75: 0,
76: };
78: #define ns 4
79: /* for "" attribute will be stored to group, otherwise to given dataset */
80: static const char *datasets[ns] = {"", "x", "nonExistingVec", "y"};
82: /* beware this yields PETSC_FALSE for "" but group "" is interpreted as "/" */
83: static inline PetscErrorCode shouldExist(const char name[], PetscBool emptyExists, PetscBool *has)
84: {
85: size_t len = 0;
87: PetscFunctionBegin;
88: PetscCall(PetscStrlen(name, &len));
89: *has = emptyExists;
90: if (len) {
91: char *loc = NULL;
92: PetscCall(PetscStrstr(name, "nonExisting", &loc));
93: *has = PetscNot(loc);
94: }
95: PetscFunctionReturn(PETSC_SUCCESS);
96: }
98: static inline PetscErrorCode isPop(const char path[], PetscBool *has)
99: {
100: PetscFunctionBegin;
101: PetscCall(PetscStrcmp(path, "<", has));
102: PetscFunctionReturn(PETSC_SUCCESS);
103: }
105: static inline PetscErrorCode isDot(const char path[], PetscBool *has)
106: {
107: PetscFunctionBegin;
108: PetscCall(PetscStrcmp(path, ".", has));
109: PetscFunctionReturn(PETSC_SUCCESS);
110: }
112: static inline PetscErrorCode isRoot(const char path[], PetscBool *flg)
113: {
114: size_t len;
116: PetscFunctionBegin;
117: PetscCall(PetscStrlen(path, &len));
118: *flg = PetscNot(len);
119: if (!*flg) PetscCall(PetscStrcmp(path, "/", flg));
120: PetscFunctionReturn(PETSC_SUCCESS);
121: }
123: static inline PetscErrorCode compare(PetscDataType dt, void *ptr0, void *ptr1, PetscBool *flg)
124: {
125: PetscFunctionBegin;
126: switch (dt) {
127: case PETSC_INT:
128: *flg = (PetscBool)(*(PetscInt *)ptr0 == *(PetscInt *)ptr1);
129: if (verbose) {
130: if (*flg) {
131: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT, *(PetscInt *)ptr0));
132: } else {
133: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " != %" PetscInt_FMT "\n", *(PetscInt *)ptr0, *(PetscInt *)ptr1));
134: }
135: }
136: break;
137: case PETSC_REAL:
138: *flg = (PetscBool)(*(PetscReal *)ptr0 == *(PetscReal *)ptr1);
139: if (verbose) {
140: if (*flg) {
141: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%f", *(PetscReal *)ptr0));
142: } else {
143: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%f != %f\n", *(PetscReal *)ptr0, *(PetscReal *)ptr1));
144: }
145: }
146: break;
147: case PETSC_BOOL:
148: *flg = (PetscBool)(*(PetscBool *)ptr0 == *(PetscBool *)ptr1);
149: if (verbose) {
150: if (*flg) {
151: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%s", PetscBools[*(PetscBool *)ptr0]));
152: } else {
153: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%s != %s\n", PetscBools[*(PetscBool *)ptr0], PetscBools[*(PetscBool *)ptr1]));
154: }
155: }
156: break;
157: case PETSC_STRING:
158: PetscCall(PetscStrcmp((const char *)ptr0, (const char *)ptr1, flg));
159: if (verbose) {
160: if (*flg) {
161: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%s", (char *)ptr0));
162: } else {
163: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%s != %s\n", (char *)ptr0, (char *)ptr1));
164: }
165: }
166: break;
167: default:
168: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscDataType %s not handled here", PetscDataTypes[dt]);
169: }
170: PetscFunctionReturn(PETSC_SUCCESS);
171: }
173: static inline PetscErrorCode alterString(const char oldstr[], char str[])
174: {
175: size_t i, n;
177: PetscFunctionBegin;
178: PetscCall(PetscStrlen(oldstr, &n));
179: PetscCall(PetscStrncpy(str, oldstr, n + 1));
180: for (i = 0; i < n; i++) {
181: if (('A' <= str[i] && str[i] < 'Z') || ('a' <= str[i] && str[i] < 'z')) {
182: str[i]++;
183: break;
184: }
185: }
186: PetscFunctionReturn(PETSC_SUCCESS);
187: }
189: /* if name given, check dataset with this name exists under current group, otherwise just check current group exists */
190: /* flg: 0 doesn't exist, 1 group, 2 dataset */
191: static PetscErrorCode hasGroupOrDataset(PetscViewer viewer, const char path[], int *flg)
192: {
193: PetscBool has;
195: PetscFunctionBegin;
196: *flg = 0;
197: PetscCall(PetscViewerHDF5HasGroup(viewer, path, &has));
198: if (has) *flg = 1;
199: else {
200: PetscCall(PetscViewerHDF5HasDataset(viewer, path, &has));
201: if (has) *flg = 2;
202: }
203: PetscFunctionReturn(PETSC_SUCCESS);
204: }
206: #define nt 5 /* number of datatypes */
207: typedef struct _n_Capsule *Capsule;
208: struct _n_Capsule {
209: char names[nt][SLEN];
210: PetscDataType types[nt];
211: char typeNames[nt][SLEN];
212: size_t sizes[nt];
213: void *vals[nt];
214: PetscInt id, ntypes;
215: };
217: static PetscErrorCode CapsuleCreate(Capsule old, Capsule *newcapsule)
218: {
219: Capsule c;
220: PetscBool bool0 = PETSC_TRUE;
221: PetscInt int0 = -1;
222: PetscReal real0 = -1.1;
223: char str0[] = "Test String";
224: char nestr0[] = "NONEXISTING STRING"; /* this attribute shall be skipped for writing */
225: void *vals[nt] = {&bool0, &int0, &real0, str0, nestr0};
226: size_t sizes[nt] = {sizeof(bool0), sizeof(int0), sizeof(real0), sizeof(str0), sizeof(str0)};
227: PetscDataType types[nt] = {PETSC_BOOL, PETSC_INT, PETSC_REAL, PETSC_STRING, PETSC_STRING};
228: const char *tNames[nt] = {"bool", "int", "real", "str", "nonExisting"};
229: PetscInt t;
231: PetscFunctionBegin;
232: PetscCall(PetscNew(&c));
233: c->id = 0;
234: c->ntypes = nt;
235: if (old) {
236: /* alter values */
237: t = 0;
238: bool0 = PetscNot(*((PetscBool *)old->vals[t]));
239: t++;
240: int0 = *((PetscInt *)old->vals[t]) * -2;
241: t++;
242: real0 = *((PetscReal *)old->vals[t]) * -2.0;
243: t++;
244: PetscCall(alterString((const char *)old->vals[t], str0));
245: t++;
246: c->id = old->id + 1;
247: }
248: for (t = 0; t < nt; t++) {
249: c->sizes[t] = sizes[t];
250: c->types[t] = types[t];
251: PetscCall(PetscStrncpy(c->typeNames[t], tNames[t], sizeof(c->typeNames[t])));
252: PetscCall(PetscSNPrintf(c->names[t], SLEN, "attr_%" PetscInt_FMT "_%s", c->id, tNames[t]));
253: PetscCall(PetscMalloc(sizes[t], &c->vals[t]));
254: PetscCall(PetscMemcpy(c->vals[t], vals[t], sizes[t]));
255: }
256: *newcapsule = c;
257: PetscFunctionReturn(PETSC_SUCCESS);
258: }
259: #undef nt
261: static PetscErrorCode CapsuleWriteAttributes(Capsule c, PetscViewer v, const char parent[])
262: {
263: PetscInt t;
264: PetscBool flg = PETSC_FALSE;
266: PetscFunctionBegin;
267: for (t = 0; t < c->ntypes; t++) {
268: PetscCall(shouldExist(c->names[t], PETSC_FALSE, &flg));
269: if (!flg) continue;
270: PetscCall(PetscViewerHDF5WriteAttribute(v, parent, c->names[t], c->types[t], c->vals[t]));
271: }
272: PetscFunctionReturn(PETSC_SUCCESS);
273: }
275: static PetscErrorCode CapsuleReadAndCompareAttributes(Capsule c, PetscViewer v, const char parent[])
276: {
277: char *group;
278: int gd = 0;
279: PetscInt t;
280: PetscBool flg = PETSC_FALSE, hasAttr = PETSC_FALSE;
281: MPI_Comm comm;
283: PetscFunctionBegin;
284: PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
285: PetscCall(PetscViewerHDF5GetGroup(v, NULL, &group));
286: PetscCall(hasGroupOrDataset(v, parent, &gd));
287: /* check correct existence of attributes */
288: for (t = 0; t < c->ntypes; t++) {
289: const char *attribute = c->names[t];
290: PetscCall(shouldExist(attribute, PETSC_FALSE, &flg));
291: PetscCall(PetscViewerHDF5HasAttribute(v, parent, attribute, &hasAttr));
292: if (verbose) {
293: PetscCall(PetscPrintf(comm, " %-24s = ", attribute));
294: if (!hasAttr) PetscCall(PetscPrintf(comm, "---"));
295: }
296: PetscCheck(gd || !hasAttr, comm, PETSC_ERR_PLIB, "Attribute %s/%s/%s exists while its parent %s/%s doesn't exist", group, parent, attribute, group, parent);
297: PetscCheck(flg == hasAttr, comm, PETSC_ERR_PLIB, "Attribute %s/%s should exist? %s Exists? %s", parent, attribute, PetscBools[flg], PetscBools[hasAttr]);
299: /* check loaded attributes are the same as original */
300: if (hasAttr) {
301: char buffer[SLEN];
302: char *str;
303: void *ptr0;
304: /* check the stored data is the same as original */
305: //TODO datatype should better be output arg, not input
306: //TODO string attributes should probably have a separate function since the handling is different;
307: //TODO or maybe it should just accept string buffer rather than pointer to string
308: if (c->types[t] == PETSC_STRING) {
309: PetscCall(PetscViewerHDF5ReadAttribute(v, parent, attribute, c->types[t], NULL, &str));
310: ptr0 = str;
311: } else {
312: PetscCall(PetscViewerHDF5ReadAttribute(v, parent, attribute, c->types[t], NULL, &buffer));
313: ptr0 = &buffer;
314: }
315: PetscCall(compare(c->types[t], ptr0, c->vals[t], &flg));
316: PetscCheck(flg, comm, PETSC_ERR_PLIB, "Value of attribute %s/%s/%s is not equal to the original value", group, parent, attribute);
317: if (verbose) PetscCall(PetscPrintf(comm, " (=)"));
318: if (c->types[t] == PETSC_STRING) PetscCall(PetscFree(str));
319: }
320: if (verbose && gd) PetscCall(PetscPrintf(comm, "\n"));
321: }
322: PetscCall(PetscFree(group));
323: PetscFunctionReturn(PETSC_SUCCESS);
324: }
326: static PetscErrorCode CapsuleDestroy(Capsule *c)
327: {
328: PetscInt t;
330: PetscFunctionBegin;
331: if (!*c) PetscFunctionReturn(PETSC_SUCCESS);
332: for (t = 0; t < (*c)->ntypes; t++) PetscCall(PetscFree((*c)->vals[t]));
333: PetscCall(PetscFree(*c));
334: PetscFunctionReturn(PETSC_SUCCESS);
335: }
337: static PetscErrorCode testGroupsDatasets(PetscViewer viewer)
338: {
339: char buf[PETSC_MAX_PATH_LEN];
340: Vec vecs[nap][ns];
341: PetscInt p, s;
342: PetscBool flg = PETSC_FALSE, flg1 = PETSC_FALSE, flg2 = PETSC_FALSE;
343: PetscRandom rand;
344: const char *filename;
345: MPI_Comm comm;
347: PetscFunctionBegin;
348: PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
349: PetscCall(PetscViewerFileGetName(viewer, &filename));
350: if (verbose) PetscCall(PetscPrintf(comm, "# TEST testGroupsDatasets\n"));
351: /* store random vectors */
352: PetscCall(PetscRandomCreate(comm, &rand));
353: PetscCall(PetscRandomSetInterval(rand, 0.0, 10.0));
354: PetscCall(PetscRandomSetFromOptions(rand));
355: PetscCall(PetscMemzero(vecs, nap * ns * sizeof(Vec)));
357: /* test dataset writing */
358: if (verbose) PetscCall(PetscPrintf(comm, "## WRITE PHASE\n"));
359: for (p = 0; p < np; p++) {
360: PetscCall(isPop(paths[p], &flg));
361: PetscCall(isDot(paths[p], &flg1));
362: PetscCall(shouldExist(apaths[paths2apaths[p]], PETSC_FALSE, &flg2));
363: if (flg) {
364: PetscCall(PetscViewerHDF5PopGroup(viewer));
365: } else {
366: PetscCall(PetscViewerHDF5PushGroup(viewer, paths[p]));
367: }
368: if (verbose) PetscCall(PetscPrintf(comm, "%-32s => %4s => %-32s should exist? %s\n", paths[p], flg ? "pop" : "push", apaths[paths2apaths[p]], PetscBools[flg2]));
369: if (flg || flg1 || !flg2) continue;
371: for (s = 0; s < ns; s++) {
372: Vec v;
374: PetscCall(shouldExist(datasets[s], PETSC_FALSE, &flg));
375: if (!flg) continue;
377: PetscCall(VecCreate(comm, &v));
378: PetscCall(PetscObjectSetName((PetscObject)v, datasets[s]));
379: PetscCall(VecSetSizes(v, n, PETSC_DECIDE));
380: PetscCall(VecSetFromOptions(v));
381: PetscCall(VecSetRandom(v, rand));
382: if (verbose) {
383: PetscReal min, max;
384: PetscCall(VecMin(v, NULL, &min));
385: PetscCall(VecMax(v, NULL, &max));
386: PetscCall(PetscPrintf(comm, " Create dataset %s/%s, keep in memory in vecs[%" PetscInt_FMT "][%" PetscInt_FMT "], min %.3e max %.3e\n", apaths[paths2apaths[p]], datasets[s], paths2apaths[p], s, min, max));
387: }
389: PetscCall(VecView(v, viewer));
390: vecs[paths2apaths[p]][s] = v;
391: }
392: }
393: PetscCall(PetscViewerFlush(viewer));
394: PetscCall(PetscRandomDestroy(&rand));
396: if (verbose) PetscCall(PetscPrintf(comm, "\n## READ PHASE\n"));
397: /* check correct existence of groups in file */
398: for (p = 0; p < np; p++) {
399: char *group;
400: const char *expected = apaths[paths2apaths[p]];
402: /* check Push/Pop is correct */
403: PetscCall(isPop(paths[p], &flg));
404: if (flg) {
405: PetscCall(PetscViewerHDF5PopGroup(viewer));
406: } else {
407: PetscCall(PetscViewerHDF5PushGroup(viewer, paths[p]));
408: }
409: PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
410: PetscCall(PetscViewerHDF5HasGroup(viewer, NULL, &flg1));
411: if (verbose) PetscCall(PetscPrintf(comm, "%-32s => %4s => %-32s exists? %s\n", paths[p], flg ? "pop" : "push", group, PetscBools[flg1]));
412: PetscCall(PetscStrcmp(group, expected, &flg2));
413: PetscCheck(flg2, comm, PETSC_ERR_PLIB, "Current group %s not equal to expected %s", group, expected);
414: PetscCall(shouldExist(group, PETSC_TRUE, &flg2));
415: PetscCheck(flg1 == flg2, comm, PETSC_ERR_PLIB, "Group %s should exist? %s Exists in %s? %s", group, PetscBools[flg2], filename, PetscBools[flg1]);
416: PetscCall(PetscFree(group));
417: }
419: /* check existence of datasets; compare loaded vectors with original ones */
420: for (p = 0; p < np; p++) {
421: char *group;
423: /* check Push/Pop is correct */
424: PetscCall(isPop(paths[p], &flg));
425: if (flg) {
426: PetscCall(PetscViewerHDF5PopGroup(viewer));
427: } else {
428: PetscCall(PetscViewerHDF5PushGroup(viewer, paths[p]));
429: }
430: PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
431: PetscCall(PetscViewerHDF5HasGroup(viewer, NULL, &flg));
432: if (verbose) PetscCall(PetscPrintf(comm, "Has %s group? %s\n", group, PetscBools[flg]));
433: for (s = 0; s < ns; s++) {
434: const char *name = datasets[s];
435: char *fullname = buf;
437: /* check correct existence of datasets in file */
438: PetscCall(PetscSNPrintf(fullname, sizeof(buf), "%s/%s", group, name));
439: PetscCall(shouldExist(name, PETSC_FALSE, &flg1));
440: flg1 = (PetscBool)(flg && flg1); /* both group and dataset need to exist */
441: PetscCall(PetscViewerHDF5HasDataset(viewer, name, &flg2));
442: if (verbose) PetscCall(PetscPrintf(comm, " %s dataset? %s", fullname, PetscBools[flg2]));
443: PetscCheck(flg2 == flg1, comm, PETSC_ERR_PLIB, "Dataset %s should exist? %s Exists in %s? %s", fullname, PetscBools[flg1], filename, PetscBools[flg2]);
445: if (flg2) {
446: Vec v;
447: /* check loaded Vec is the same as original */
448: PetscCall(VecCreate(comm, &v));
449: PetscCall(PetscObjectSetName((PetscObject)v, name));
450: PetscCall(VecLoad(v, viewer));
451: PetscCall(VecEqual(v, vecs[paths2apaths[p]][s], &flg1));
452: PetscCheck(flg1, comm, PETSC_ERR_PLIB, "Dataset %s in %s is not equal to the original Vec", fullname, filename);
453: if (verbose) PetscCall(PetscPrintf(comm, " (=)"));
454: PetscCall(VecDestroy(&v));
455: }
456: if (verbose) PetscCall(PetscPrintf(comm, "\n"));
457: }
458: PetscCall(PetscFree(group));
459: }
460: PetscCall(PetscViewerFlush(viewer));
461: for (p = 0; p < nap; p++)
462: for (s = 0; s < ns; s++) PetscCall(VecDestroy(&vecs[p][s]));
463: if (verbose) PetscCall(PetscPrintf(comm, "# END testGroupsDatasets\n\n"));
464: PetscFunctionReturn(PETSC_SUCCESS);
465: }
467: static inline PetscErrorCode formPath(PetscBool relativize, const char path[], const char dataset[], char buf[], size_t bufsize)
468: {
469: PetscBool isroot = PETSC_FALSE;
471: PetscFunctionBegin;
472: PetscCall(isRoot(path, &isroot));
473: if (relativize) {
474: if (isroot) {
475: PetscCall(PetscStrncpy(buf, dataset, bufsize));
476: } else {
477: /* skip initial '/' in paths[p] if prefix given */
478: PetscCall(PetscSNPrintf(buf, bufsize, "%s/%s", path + 1, dataset));
479: }
480: } else {
481: PetscCall(PetscSNPrintf(buf, bufsize, "%s/%s", isroot ? "" : path, dataset));
482: }
483: PetscFunctionReturn(PETSC_SUCCESS);
484: }
486: /* test attribute writing, existence checking and reading, use absolute paths */
487: static PetscErrorCode testAttributesAbsolutePath(PetscViewer viewer, const char prefix[])
488: {
489: char buf[PETSC_MAX_PATH_LEN];
490: Capsule capsules[nap][ns], c = NULL, old = NULL;
491: PetscInt p, s;
492: PetscBool flg = PETSC_FALSE, flg1 = PETSC_FALSE;
493: MPI_Comm comm;
495: PetscFunctionBegin;
496: PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
497: if (verbose) {
498: if (prefix) {
499: PetscCall(PetscPrintf(comm, "# TEST testAttributesAbsolutePath, prefix=\"%s\"\n", prefix));
500: } else {
501: PetscCall(PetscPrintf(comm, "# TEST testAttributesAbsolutePath\n"));
502: }
503: PetscCall(PetscPrintf(comm, "## WRITE PHASE\n"));
504: }
505: PetscCall(PetscMemzero(capsules, nap * ns * sizeof(Capsule)));
507: /* test attribute writing */
508: if (prefix) PetscCall(PetscViewerHDF5PushGroup(viewer, prefix));
509: for (p = 0; p < np; p++)
510: for (s = 0; s < ns; s++) {
511: /* we test only absolute paths here */
512: PetscCall(PetscViewerHDF5PathIsRelative(paths[p], PETSC_FALSE, &flg));
513: if (flg) continue;
514: {
515: char *group;
517: PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
518: PetscCall(PetscStrcmp(group, prefix, &flg));
519: PetscCheck(flg, comm, PETSC_ERR_PLIB, "prefix %s not equal to pushed group %s", prefix, group);
520: PetscCall(PetscFree(group));
521: }
522: PetscCall(formPath((PetscBool) !!prefix, paths[p], datasets[s], buf, sizeof(buf)));
523: PetscCall(shouldExist(buf, PETSC_TRUE, &flg));
524: if (!flg) continue;
526: if (verbose) {
527: if (prefix) {
528: PetscCall(PetscPrintf(comm, "Write attributes to %s/%s\n", prefix, buf));
529: } else {
530: PetscCall(PetscPrintf(comm, "Write attributes to %s\n", buf));
531: }
532: }
534: PetscCall(CapsuleCreate(old, &c));
535: PetscCall(CapsuleWriteAttributes(c, viewer, buf));
536: PetscCheck(!capsules[paths2apaths[p]][s], comm, PETSC_ERR_PLIB, "capsules[%" PetscInt_FMT "][%" PetscInt_FMT "] gets overwritten for %s", paths2apaths[p], s, buf);
537: capsules[paths2apaths[p]][s] = c;
538: old = c;
539: }
540: if (prefix) PetscCall(PetscViewerHDF5PopGroup(viewer));
541: PetscCall(PetscViewerFlush(viewer));
543: if (verbose) PetscCall(PetscPrintf(comm, "\n## READ PHASE\n"));
544: if (prefix) PetscCall(PetscViewerHDF5PushGroup(viewer, prefix));
545: for (p = 0; p < np; p++)
546: for (s = 0; s < ns; s++) {
547: /* we test only absolute paths here */
548: PetscCall(PetscViewerHDF5PathIsRelative(paths[p], PETSC_FALSE, &flg));
549: if (flg) continue;
551: /* check existence of given group/dataset */
552: PetscCall(formPath((PetscBool) !!prefix, paths[p], datasets[s], buf, sizeof(buf)));
553: PetscCall(shouldExist(buf, PETSC_TRUE, &flg));
554: if (verbose) {
555: if (prefix) {
556: PetscCall(PetscPrintf(comm, "Has %s/%s? %s\n", prefix, buf, PetscBools[flg]));
557: } else {
558: PetscCall(PetscPrintf(comm, "Has %s? %s\n", buf, PetscBools[flg]));
559: }
560: }
562: /* check attribute capsule has been created for given path */
563: c = capsules[paths2apaths[p]][s];
564: flg1 = (PetscBool) !!c;
565: PetscCheck(flg == flg1, comm, PETSC_ERR_PLIB, "Capsule should exist for %s? %s Exists? %s", buf, PetscBools[flg], PetscBools[flg1]);
566: if (!flg) continue;
568: /* check correct existence and fidelity of attributes in file */
569: PetscCall(CapsuleReadAndCompareAttributes(c, viewer, buf));
570: }
571: if (prefix) PetscCall(PetscViewerHDF5PopGroup(viewer));
572: PetscCall(PetscViewerFlush(viewer));
573: for (p = 0; p < nap; p++)
574: for (s = 0; s < ns; s++) PetscCall(CapsuleDestroy(&capsules[p][s]));
575: if (verbose) PetscCall(PetscPrintf(comm, "# END testAttributesAbsolutePath\n\n"));
576: PetscFunctionReturn(PETSC_SUCCESS);
577: }
579: /* test attribute writing, existence checking and reading, use group push/pop */
580: static PetscErrorCode testAttributesPushedPath(PetscViewer viewer)
581: {
582: Capsule capsules[nap][ns], c = NULL, old = NULL;
583: PetscInt p, s;
584: int gd;
585: PetscBool flg = PETSC_FALSE, flg1 = PETSC_FALSE;
586: MPI_Comm comm;
588: PetscFunctionBegin;
589: PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
590: if (verbose) {
591: PetscCall(PetscPrintf(comm, "# TEST testAttributesPushedPath\n"));
592: PetscCall(PetscPrintf(comm, "## WRITE PHASE\n"));
593: }
594: PetscCall(PetscMemzero(capsules, nap * ns * sizeof(Capsule)));
596: /* test attribute writing */
597: for (p = 0; p < np; p++) {
598: PetscCall(isPop(paths[p], &flg));
599: PetscCall(isDot(paths[p], &flg1));
600: if (flg) {
601: PetscCall(PetscViewerHDF5PopGroup(viewer));
602: } else {
603: PetscCall(PetscViewerHDF5PushGroup(viewer, paths[p]));
604: }
605: /* < and . have been already visited => skip */
606: if (flg || flg1) continue;
608: /* assume here that groups and datasets are already in the file */
609: for (s = 0; s < ns; s++) {
610: PetscCall(hasGroupOrDataset(viewer, datasets[s], &gd));
611: if (!gd) continue;
612: if (verbose) PetscCall(PetscPrintf(comm, "Write attributes to %s/%s\n", apaths[paths2apaths[p]], datasets[s]));
613: PetscCall(CapsuleCreate(old, &c));
614: PetscCall(CapsuleWriteAttributes(c, viewer, datasets[s]));
615: PetscCheck(!capsules[paths2apaths[p]][s], comm, PETSC_ERR_PLIB, "capsules[%" PetscInt_FMT "][%" PetscInt_FMT "] gets overwritten for %s/%s", paths2apaths[p], s, paths[p], datasets[s]);
616: capsules[paths2apaths[p]][s] = c;
617: old = c;
618: }
619: }
620: PetscCall(PetscViewerFlush(viewer));
622: if (verbose) PetscCall(PetscPrintf(comm, "\n## READ PHASE\n"));
623: for (p = 0; p < np; p++) {
624: char *group;
626: PetscCall(isPop(paths[p], &flg1));
627: if (flg1) {
628: PetscCall(PetscViewerHDF5PopGroup(viewer));
629: } else {
630: PetscCall(PetscViewerHDF5PushGroup(viewer, paths[p]));
631: }
632: PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
633: for (s = 0; s < ns; s++) {
634: PetscCall(hasGroupOrDataset(viewer, datasets[s], &gd));
635: if (verbose) PetscCall(PetscPrintf(comm, "%s/%s %s\n", group, datasets[s], gd ? (gd == 1 ? "is group" : "is dataset") : "does not exist"));
637: /* check attribute capsule has been created for given path */
638: c = capsules[paths2apaths[p]][s];
639: flg = (PetscBool) !!gd;
640: flg1 = (PetscBool) !!c;
641: PetscCheck(flg == flg1, comm, PETSC_ERR_PLIB, "Capsule should exist for %s/%s? %s Exists? %s", group, datasets[s], PetscBools[flg], PetscBools[flg1]);
642: if (!flg) continue;
644: /* check correct existence of attributes in file */
645: PetscCall(CapsuleReadAndCompareAttributes(c, viewer, datasets[s]));
646: }
647: PetscCall(PetscFree(group));
648: }
649: PetscCall(PetscViewerFlush(viewer));
650: for (p = 0; p < nap; p++)
651: for (s = 0; s < ns; s++) PetscCall(CapsuleDestroy(&capsules[p][s]));
652: if (verbose) PetscCall(PetscPrintf(comm, "# END testAttributesPushedPath\n\n"));
653: PetscFunctionReturn(PETSC_SUCCESS);
654: }
656: /* test attribute writing, existence checking and reading, use group push/pop */
657: static PetscErrorCode testObjectAttributes(PetscViewer viewer)
658: {
659: Capsule capsules[nap][ns], c = NULL, old = NULL;
660: PetscInt p, s;
661: PetscBool flg = PETSC_FALSE, flg1 = PETSC_FALSE;
662: MPI_Comm comm;
664: PetscFunctionBegin;
665: PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
666: if (verbose) {
667: PetscCall(PetscPrintf(comm, "# TEST testObjectAttributes\n"));
668: PetscCall(PetscPrintf(comm, "## WRITE PHASE\n"));
669: }
670: PetscCall(PetscMemzero(capsules, nap * ns * sizeof(Capsule)));
672: /* test attribute writing */
673: for (p = 0; p < np; p++) {
674: PetscCall(isPop(paths[p], &flg));
675: PetscCall(isDot(paths[p], &flg1));
676: if (flg) {
677: PetscCall(PetscViewerHDF5PopGroup(viewer));
678: } else {
679: PetscCall(PetscViewerHDF5PushGroup(viewer, paths[p]));
680: }
681: /* < and . have been already visited => skip */
682: if (flg || flg1) continue;
684: /* assume here that groups and datasets are already in the file */
685: for (s = 0; s < ns; s++) {
686: Vec v;
687: size_t len;
688: const char *name = datasets[s];
690: PetscCall(PetscStrlen(name, &len));
691: if (!len) continue;
692: PetscCall(VecCreate(comm, &v));
693: PetscCall(PetscObjectSetName((PetscObject)v, name));
694: PetscCall(PetscViewerHDF5HasObject(viewer, (PetscObject)v, &flg));
695: if (flg) {
696: if (verbose) PetscCall(PetscPrintf(comm, "Write attributes to %s/%s\n", apaths[paths2apaths[p]], name));
697: PetscCall(CapsuleCreate(old, &c));
698: PetscCall(CapsuleWriteAttributes(c, viewer, name));
699: PetscCheck(!capsules[paths2apaths[p]][s], comm, PETSC_ERR_PLIB, "capsules[%" PetscInt_FMT "][%" PetscInt_FMT "] gets overwritten for %s/%s", paths2apaths[p], s, paths[p], name);
700: capsules[paths2apaths[p]][s] = c;
701: old = c;
702: }
703: PetscCall(VecDestroy(&v));
704: }
705: }
706: PetscCall(PetscViewerFlush(viewer));
708: if (verbose) PetscCall(PetscPrintf(comm, "\n## READ PHASE\n"));
709: for (p = 0; p < np; p++) {
710: char *group;
712: PetscCall(isPop(paths[p], &flg));
713: if (flg) {
714: PetscCall(PetscViewerHDF5PopGroup(viewer));
715: } else {
716: PetscCall(PetscViewerHDF5PushGroup(viewer, paths[p]));
717: }
718: PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
719: for (s = 0; s < ns; s++) {
720: Vec v;
721: size_t len;
722: const char *name = datasets[s];
724: PetscCall(PetscStrlen(name, &len));
725: if (!len) continue;
726: PetscCall(VecCreate(comm, &v));
727: PetscCall(PetscObjectSetName((PetscObject)v, name));
728: PetscCall(PetscViewerHDF5HasObject(viewer, (PetscObject)v, &flg));
729: if (verbose) PetscCall(PetscPrintf(comm, "Is %s/%s dataset? %s\n", group, name, PetscBools[flg]));
731: /* check attribute capsule has been created for given path */
732: c = capsules[paths2apaths[p]][s];
733: flg1 = (PetscBool) !!c;
734: PetscCheck(flg == flg1, comm, PETSC_ERR_PLIB, "Capsule should exist for %s/%s? %s Exists? %s", group, name, PetscBools[flg], PetscBools[flg1]);
736: /* check correct existence of attributes in file */
737: if (flg) PetscCall(CapsuleReadAndCompareAttributes(c, viewer, name));
738: PetscCall(VecDestroy(&v));
739: }
740: PetscCall(PetscFree(group));
741: }
742: PetscCall(PetscViewerFlush(viewer));
743: for (p = 0; p < nap; p++)
744: for (s = 0; s < ns; s++) PetscCall(CapsuleDestroy(&capsules[p][s]));
745: if (verbose) PetscCall(PetscPrintf(comm, "# END testObjectAttributes\n\n"));
746: PetscFunctionReturn(PETSC_SUCCESS);
747: }
749: static PetscErrorCode testAttributesDefaultValue(PetscViewer viewer)
750: {
751: #define nv 4
752: PetscBool bools[nv];
753: PetscInt ints[nv];
754: PetscReal reals[nv];
755: char *strings[nv];
756: PetscBool flg;
757: PetscInt i;
758: MPI_Comm comm;
760: PetscFunctionBegin;
761: PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
762: if (verbose) PetscCall(PetscPrintf(comm, "# TEST testAttributesDefaultValue\n"));
764: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_0_bool", PETSC_BOOL, NULL, &bools[0]));
765: bools[1] = PetscNot(bools[0]);
766: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_0_bool", PETSC_BOOL, &bools[1], &bools[2]));
767: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_nonExisting_bool", PETSC_BOOL, &bools[1], &bools[3]));
768: PetscCheck(bools[2] == bools[0], comm, PETSC_ERR_PLIB, "%s = bools[2] != bools[0] = %s", PetscBools[bools[2]], PetscBools[bools[0]]);
769: PetscCheck(bools[3] == bools[1], comm, PETSC_ERR_PLIB, "%s = bools[3] != bools[1] = %s", PetscBools[bools[3]], PetscBools[bools[1]]);
771: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_0_int", PETSC_INT, NULL, &ints[0]));
772: ints[1] = ints[0] * -333;
773: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_0_int", PETSC_INT, &ints[1], &ints[2]));
774: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_nonExisting_int", PETSC_INT, &ints[1], &ints[3]));
775: PetscCheck(ints[2] == ints[0], comm, PETSC_ERR_PLIB, "%" PetscInt_FMT " = ints[2] != ints[0] = %" PetscInt_FMT, ints[2], ints[0]);
776: PetscCheck(ints[3] == ints[1], comm, PETSC_ERR_PLIB, "%" PetscInt_FMT " = ints[3] != ints[1] = %" PetscInt_FMT, ints[3], ints[1]);
777: if (verbose) PetscCall(PetscIntView(nv, ints, PETSC_VIEWER_STDOUT_WORLD));
779: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_0_real", PETSC_REAL, NULL, &reals[0]));
780: reals[1] = reals[0] * -11.1;
781: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_0_real", PETSC_REAL, &reals[1], &reals[2]));
782: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_nonExisting_real", PETSC_REAL, &reals[1], &reals[3]));
783: PetscCheck(reals[2] == reals[0], comm, PETSC_ERR_PLIB, "%f = reals[2] != reals[0] = %f", reals[2], reals[0]);
784: PetscCheck(reals[3] == reals[1], comm, PETSC_ERR_PLIB, "%f = reals[3] != reals[1] = %f", reals[3], reals[1]);
785: if (verbose) PetscCall(PetscRealView(nv, reals, PETSC_VIEWER_STDOUT_WORLD));
787: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_0_str", PETSC_STRING, NULL, &strings[0]));
788: PetscCall(PetscStrallocpy(strings[0], &strings[1]));
789: PetscCall(alterString(strings[0], strings[1]));
790: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_0_str", PETSC_STRING, &strings[1], &strings[2]));
791: PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", "attr_nonExisting_str", PETSC_STRING, &strings[1], &strings[3]));
792: PetscCall(PetscStrcmp(strings[2], strings[0], &flg));
793: PetscCheck(flg, comm, PETSC_ERR_PLIB, "%s = strings[2] != strings[0] = %s", strings[2], strings[0]);
794: PetscCall(PetscStrcmp(strings[3], strings[1], &flg));
795: PetscCheck(flg, comm, PETSC_ERR_PLIB, "%s = strings[3] != strings[1] = %s", strings[3], strings[1]);
796: for (i = 0; i < nv; i++) PetscCall(PetscFree(strings[i]));
798: PetscCall(PetscViewerFlush(viewer));
799: if (verbose) PetscCall(PetscPrintf(comm, "# END testAttributesDefaultValue\n"));
800: #undef nv
801: PetscFunctionReturn(PETSC_SUCCESS);
802: }
804: int main(int argc, char **argv)
805: {
806: static char filename[PETSC_MAX_PATH_LEN] = "ex48.h5";
807: PetscMPIInt rank;
808: MPI_Comm comm;
809: PetscViewer viewer;
811: PetscFunctionBegin;
812: PetscFunctionBeginUser;
813: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
814: comm = PETSC_COMM_WORLD;
815: PetscCallMPI(MPI_Comm_rank(comm, &rank));
816: PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
817: PetscCall(PetscOptionsGetBool(NULL, NULL, "-verbose", &verbose, NULL));
818: PetscCall(PetscOptionsGetString(NULL, NULL, "-filename", filename, sizeof(filename), NULL));
819: if (verbose) PetscCall(PetscPrintf(comm, "np ns " PetscStringize(np) " " PetscStringize(ns) "\n"));
821: PetscCall(PetscViewerHDF5Open(comm, filename, FILE_MODE_WRITE, &viewer));
822: PetscCall(testGroupsDatasets(viewer));
823: PetscCall(testAttributesAbsolutePath(viewer, "/"));
824: PetscCall(testAttributesAbsolutePath(viewer, "/prefix"));
825: PetscCall(PetscViewerDestroy(&viewer));
827: /* test reopening in update mode */
828: PetscCall(PetscViewerHDF5Open(comm, filename, FILE_MODE_UPDATE, &viewer));
829: PetscCall(testAttributesPushedPath(viewer));
830: PetscCall(testObjectAttributes(viewer));
831: PetscCall(testAttributesDefaultValue(viewer));
832: PetscCall(PetscViewerDestroy(&viewer));
833: PetscCall(PetscFinalize());
834: return 0;
835: }
837: /*TEST
839: build:
840: requires: hdf5
842: test:
843: nsize: {{1 4}}
845: TEST*/