summaryrefslogtreecommitdiffstats
diff options
context:
space:
mode:
authorJerome Soumagne <jsoumagne@hdfgroup.org>2016-02-29 00:47:48 (GMT)
committerJerome Soumagne <jsoumagne@hdfgroup.org>2016-11-29 23:42:32 (GMT)
commita52c9c172d44cdb2f58c2000e655b1766d249485 (patch)
treec2c07fc612eb01719127d4df8e96c6a341b22568
parentb8d0d9ea97d4bf92020b89ddfb467fbad6b22f58 (diff)
downloadhdf5-a52c9c172d44cdb2f58c2000e655b1766d249485.zip
hdf5-a52c9c172d44cdb2f58c2000e655b1766d249485.tar.gz
hdf5-a52c9c172d44cdb2f58c2000e655b1766d249485.tar.bz2
Update H5Q after reference changes
Clean up query test
-rw-r--r--src/H5Q.c60
-rw-r--r--src/H5Qprivate.h4
-rw-r--r--test/query.c84
3 files changed, 73 insertions, 75 deletions
diff --git a/src/H5Q.c b/src/H5Q.c
index 369222f..c36b221 100644
--- a/src/H5Q.c
+++ b/src/H5Q.c
@@ -197,7 +197,7 @@ typedef enum H5Q_match_type_t { /* The different kinds of native types we can ma
typedef struct H5Q_ref_entry_t H5Q_ref_entry_t;
struct H5Q_ref_entry_t {
- href_t *ref;
+ href_t ref;
H5Q_QUEUE_ENTRY(H5Q_ref_entry_t) entry;
};
@@ -216,7 +216,7 @@ typedef struct {
} H5Q_apply_arg_t;
typedef struct {
- H5G_loc_t *loc;
+ const char *filename;
const char *loc_name;
H5G_loc_t *obj_loc;
H5Q_apply_arg_t *apply_args;
@@ -249,7 +249,7 @@ static herr_t H5Q__apply_object_link(H5G_loc_t *loc, const char *name,
static herr_t H5Q__apply_object_data(H5G_loc_t *loc, const char *name,
const H5O_info_t *oinfo, void *udata);
static herr_t H5Q__apply_object_attr(H5G_loc_t *loc, const char *name,
- const H5O_info_t *oinfo, void *udata);
+ void *udata);
static herr_t H5Q__apply_object_attr_iterate(H5A_t *attr, void *udata);
static herr_t H5Q__apply_object_attr_name(H5A_t *attr, void *udata);
static herr_t H5Q__apply_object_attr_value(H5A_t *attr, void *udata);
@@ -1836,7 +1836,6 @@ hid_t
H5Qapply(hid_t loc_id, hid_t query_id, unsigned *result, hid_t vcpl_id)
{
H5Q_t *query = NULL;
- H5G_loc_t loc;
H5G_t *ret_grp;
hid_t ret_value;
@@ -1844,8 +1843,6 @@ H5Qapply(hid_t loc_id, hid_t query_id, unsigned *result, hid_t vcpl_id)
H5TRACE4("i", "ii*Iui", loc_id, query_id, result, vcpl_id);
/* Check args and get the query objects */
- if (FAIL == H5G_loc(loc_id, &loc))
- HGOTO_ERROR(H5E_ARGS, H5E_BADTYPE, FAIL, "not a location");
if (NULL == (query = (H5Q_t *) H5I_object_verify(query_id, H5I_QUERY)))
HGOTO_ERROR(H5E_ARGS, H5E_BADTYPE, FAIL, "not a query ID");
if (!result)
@@ -1860,7 +1857,7 @@ H5Qapply(hid_t loc_id, hid_t query_id, unsigned *result, hid_t vcpl_id)
HGOTO_ERROR(H5E_ARGS, H5E_BADTYPE, FAIL, "not index access parms");
/* Apply query */
- if (NULL == (ret_grp = H5Q_apply(&loc, loc_id, query, result, vcpl_id)))
+ if (NULL == (ret_grp = H5Q_apply(loc_id, query, result, vcpl_id)))
HGOTO_ERROR(H5E_QUERY, H5E_CANTCOMPARE, FAIL, "unable to apply query");
if (FAIL == (ret_value = H5I_register(H5I_GROUP, ret_grp, TRUE)))
@@ -1880,8 +1877,8 @@ done:
*-------------------------------------------------------------------------
*/
H5G_t *
-H5Q_apply(const H5G_loc_t *loc, hid_t loc_id, const H5Q_t *query, unsigned *result,
- hid_t H5_ATTR_UNUSED vcpl_id)
+H5Q_apply(hid_t loc_id, const H5Q_t *query, unsigned *result,
+ hid_t H5_ATTR_UNUSED vcpl_id)
{
H5Q_apply_arg_t args;
H5Q_view_t view = H5Q_VIEW_INITIALIZER(view); /* Resulting view */
@@ -1895,7 +1892,6 @@ H5Q_apply(const H5G_loc_t *loc, hid_t loc_id, const H5Q_t *query, unsigned *resu
FUNC_ENTER_NOAPI_NOINIT
- HDassert(loc);
HDassert(query);
HDassert(result);
@@ -2086,7 +2082,7 @@ H5Q__apply_object(hid_t oid, const char *name, const H5O_info_t *oinfo,
break;
case H5Q_TYPE_ATTR_NAME:
case H5Q_TYPE_ATTR_VALUE:
- if (FAIL == H5Q__apply_object_attr(&loc, name, oinfo, udata))
+ if (FAIL == H5Q__apply_object_attr(&loc, name, udata))
HGOTO_ERROR(H5E_QUERY, H5E_CANTCOMPARE, FAIL, "can't apply data query to object");
break;
case H5Q_TYPE_MISC:
@@ -2112,7 +2108,7 @@ static herr_t
H5Q__apply_object_link(H5G_loc_t *loc, const char *name, const H5O_info_t *oinfo,
void *udata)
{
- href_t *ref;
+ href_t ref;
hbool_t result;
H5Q_apply_arg_t *args = (H5Q_apply_arg_t *) udata;
const char *link_name = NULL;
@@ -2142,7 +2138,7 @@ H5Q__apply_object_link(H5G_loc_t *loc, const char *name, const H5O_info_t *oinfo
H5Q_LOG_DEBUG("Match link name: %s (%s)\n", link_name, name);
/* Keep object reference */
- if (NULL == (ref = H5R_create_ext_object(loc, name, H5AC_dxpl_id)))
+ if (NULL == (ref = H5R_create_ext_object(H5F_OPEN_NAME(loc->oloc->file), name)))
HGOTO_ERROR(H5E_DATASET, H5E_CANTCREATE, FAIL, "can't create object reference");
if (FAIL == H5Q__view_append(args->view, H5R_EXT_OBJECT, ref))
HGOTO_ERROR(H5E_QUERY, H5E_CANTAPPEND, FAIL, "can't append object reference to view");
@@ -2164,7 +2160,7 @@ static herr_t
H5Q__apply_object_data(H5G_loc_t *loc, const char *name, const H5O_info_t *oinfo,
void *udata)
{
- href_t *ref;
+ href_t ref;
H5Q_apply_arg_t *args = (H5Q_apply_arg_t *) udata;
hid_t obj_id = FAIL;
H5S_t *dataspace = NULL;
@@ -2203,7 +2199,7 @@ H5Q__apply_object_data(H5G_loc_t *loc, const char *name, const H5O_info_t *oinfo
*(args->result) = H5Q_REF_REG;
/* Keep dataset region reference */
- if (NULL == (ref = H5R_create_ext_region(loc, name, H5AC_dxpl_id, dataspace)))
+ if (NULL == (ref = H5R_create_ext_region(H5F_OPEN_NAME(loc->oloc->file), name, dataspace)))
HGOTO_ERROR(H5E_DATASET, H5E_CANTGET, FAIL, "can't get buffer size for region reference");
if (FAIL == H5Q__view_append(args->view, H5R_EXT_REGION, ref))
HGOTO_ERROR(H5E_QUERY, H5E_CANTAPPEND, FAIL, "can't append region reference to view");
@@ -2225,8 +2221,7 @@ done:
*-------------------------------------------------------------------------
*/
static herr_t
-H5Q__apply_object_attr(H5G_loc_t *loc, const char *name, const H5O_info_t *oinfo,
- void *udata)
+H5Q__apply_object_attr(H5G_loc_t *loc, const char *name, void *udata)
{
H5Q_apply_attr_arg_t attr_args;
H5A_attr_iter_op_t attr_op; /* Attribute operator */
@@ -2239,11 +2234,10 @@ H5Q__apply_object_attr(H5G_loc_t *loc, const char *name, const H5O_info_t *oinfo
HDassert(loc);
HDassert(name);
- HDassert(oinfo);
HDassert(args);
/* Build attribute operator info */
- attr_args.loc = loc;
+ attr_args.filename = H5F_OPEN_NAME(loc->oloc->file);
attr_args.loc_name = name;
attr_args.apply_args = args;
attr_op.op_type = H5A_ATTR_OP_LIB;
@@ -2329,7 +2323,7 @@ H5Q__apply_object_attr_name(H5A_t *attr, void *udata)
H5Q_apply_attr_arg_t *args = (H5Q_apply_attr_arg_t *) udata;
char *name = NULL;
size_t name_len;
- href_t *ref;
+ href_t ref;
hbool_t result = FALSE;
herr_t ret_value = SUCCEED; /* Return value */
@@ -2356,7 +2350,7 @@ H5Q__apply_object_attr_name(H5A_t *attr, void *udata)
H5Q_LOG_DEBUG("Match attribute name: %s\n", (const char *) name);
/* Keep attribute reference */
- if (NULL == (ref = H5R_create_ext_attr(args->loc, args->loc_name, H5AC_dxpl_id, (const char *) name)))
+ if (NULL == (ref = H5R_create_ext_attr(args->filename, args->loc_name, (const char *) name)))
HGOTO_ERROR(H5E_DATASET, H5E_CANTGET, FAIL, "can't get buffer size for attribute reference");
if (FAIL == H5Q__view_append(args->apply_args->view, H5R_EXT_ATTR, ref))
HGOTO_ERROR(H5E_QUERY, H5E_CANTAPPEND, FAIL, "can't append attribute reference to view");
@@ -2390,7 +2384,7 @@ H5Q__apply_object_attr_value(H5A_t *iter_attr, void *udata)
size_t nelmts, elmt_size;
H5S_sel_iter_op_t iter_op; /* Operator for iteration */
H5Q_apply_attr_elem_arg_t iter_args;
- href_t *ref;
+ href_t ref;
hbool_t result = FALSE;
herr_t ret_value = SUCCEED; /* Return value */
@@ -2451,7 +2445,7 @@ H5Q__apply_object_attr_value(H5A_t *iter_attr, void *udata)
H5Q_LOG_DEBUG("Match value of attribute name: %s\n", (const char *) name);
/* Keep attribute reference */
- if (NULL == (ref = H5R_create_ext_attr(args->loc, args->loc_name, H5AC_dxpl_id, (const char *) name)))
+ if (NULL == (ref = H5R_create_ext_attr(args->filename, args->loc_name, (const char *) name)))
HGOTO_ERROR(H5E_DATASET, H5E_CANTGET, FAIL, "can't get buffer size for attribute reference");
if (FAIL == H5Q__view_append(args->apply_args->view, H5R_EXT_ATTR, ref))
HGOTO_ERROR(H5E_QUERY, H5E_CANTAPPEND, FAIL, "can't append attribute reference to view");
@@ -2684,7 +2678,7 @@ H5Q__view_write(H5G_t *grp, H5Q_view_t *view)
if (FAIL == H5S_select_hyperslab(space, H5S_SELECT_SET, &start, NULL, &count, NULL))
HGOTO_ERROR(H5E_DATASPACE, H5E_CANTINIT, FAIL, "unable to set hyperslab selection")
if (FAIL == H5D_write(dset, FALSE, ref_types[i], mem_space, space,
- H5P_DATASET_XFER_DEFAULT, ref_entry->ref))
+ H5P_DATASET_XFER_DEFAULT, &ref_entry->ref))
HGOTO_ERROR(H5E_DATASET, H5E_WRITEERROR, FAIL, "unable to write dataset");
if (FAIL == H5S_close(mem_space))
HGOTO_ERROR(H5E_DATASPACE, H5E_CANTCLOSEOBJ, FAIL, "unable to close dataspace");
@@ -2758,22 +2752,16 @@ H5Qapply_multi(size_t loc_count, hid_t *loc_ids, hid_t query_id,
unsigned *result, hid_t vcpl_id)
{
H5Q_t *query = NULL;
- H5G_loc_t *locs = NULL;
H5G_t *ret_grp;
hid_t ret_value;
- int i;
FUNC_ENTER_API(FAIL)
/* Check args and get the query objects */
if (!loc_count)
HGOTO_ERROR(H5E_ARGS, H5E_BADVALUE, FAIL, "loc_count cannot be NULL");
- if (NULL == (locs = (H5G_loc_t *) H5MM_malloc(loc_count * sizeof(H5G_loc_t))))
- HGOTO_ERROR(H5E_QUERY, H5E_CANTALLOC, FAIL, "can't allocate locs buffer");
- for (i = 0; i < loc_count; i++) {
- if (FAIL == H5G_loc(loc_ids[i], &locs[i]))
- HGOTO_ERROR(H5E_ARGS, H5E_BADTYPE, FAIL, "not a location");
- }
+ if (!loc_ids)
+ HGOTO_ERROR(H5E_ARGS, H5E_BADVALUE, FAIL, "loc_ids cannot be NULL");
if (NULL == (query = (H5Q_t *) H5I_object_verify(query_id, H5I_QUERY)))
HGOTO_ERROR(H5E_ARGS, H5E_BADTYPE, FAIL, "not a query ID");
if (!result)
@@ -2788,14 +2776,13 @@ H5Qapply_multi(size_t loc_count, hid_t *loc_ids, hid_t query_id,
HGOTO_ERROR(H5E_ARGS, H5E_BADTYPE, FAIL, "not index access parms");
/* Apply query */
- if (NULL == (ret_grp = H5Q_apply_multi(loc_count, locs, loc_ids, query, result, vcpl_id)))
+ if (NULL == (ret_grp = H5Q_apply_multi(loc_count, loc_ids, query, result, vcpl_id)))
HGOTO_ERROR(H5E_QUERY, H5E_CANTCOMPARE, FAIL, "unable to apply query");
if (FAIL == (ret_value = H5I_register(H5I_GROUP, ret_grp, TRUE)))
HGOTO_ERROR(H5E_ATOM, H5E_CANTREGISTER, FAIL, "unable to register group");
done:
- H5MM_free(locs);
FUNC_LEAVE_API(ret_value)
} /* end H5Qapply_multi() */
@@ -2809,8 +2796,8 @@ done:
*-------------------------------------------------------------------------
*/
H5G_t *
-H5Q_apply_multi(size_t loc_count, const H5G_loc_t *locs, hid_t *loc_ids,
- const H5Q_t *query, unsigned *result, hid_t H5_ATTR_UNUSED vcpl_id)
+H5Q_apply_multi(size_t loc_count, hid_t *loc_ids, const H5Q_t *query,
+ unsigned *result, hid_t H5_ATTR_UNUSED vcpl_id)
{
H5Q_view_t view = H5Q_VIEW_INITIALIZER(view); /* Resulting view */
H5Q_ref_head_t *refs[H5Q_VIEW_REF_NTYPES] = { &view.reg_refs, &view.obj_refs, &view.attr_refs };
@@ -2827,7 +2814,6 @@ H5Q_apply_multi(size_t loc_count, const H5G_loc_t *locs, hid_t *loc_ids,
FUNC_ENTER_NOAPI_NOINIT
HDassert(loc_count);
- HDassert(locs);
HDassert(loc_ids);
HDassert(query);
HDassert(result);
diff --git a/src/H5Qprivate.h b/src/H5Qprivate.h
index 2f68866..f6a90d3 100644
--- a/src/H5Qprivate.h
+++ b/src/H5Qprivate.h
@@ -107,7 +107,7 @@ H5_DLL H5Q_t *H5Q_decode(const unsigned char **buf);
/* Apply query */
H5_DLL herr_t H5Q_apply_atom(const H5Q_t *query, hbool_t *result, ...);
-H5_DLL H5G_t *H5Q_apply(const H5G_loc_t *loc, hid_t loc_id, const H5Q_t *query, unsigned *result, hid_t vcpl_id);
-H5_DLL H5G_t *H5Q_apply_multi(size_t loc_count, const H5G_loc_t *locs, hid_t *loc_ids, const H5Q_t *query, unsigned *result, hid_t vcpl_id);
+H5_DLL H5G_t *H5Q_apply(hid_t loc_id, const H5Q_t *query, unsigned *result, hid_t vcpl_id);
+H5_DLL H5G_t *H5Q_apply_multi(size_t loc_count, hid_t *loc_ids, const H5Q_t *query, unsigned *result, hid_t vcpl_id);
#endif /* _H5Qprivate_H */
diff --git a/test/query.c b/test/query.c
index 345f562..8bc4ddb 100644
--- a/test/query.c
+++ b/test/query.c
@@ -433,13 +433,14 @@ error:
/* Read region */
static herr_t
-test_query_read_selection(hid_t file, hid_t fapl, hid_t view, H5R_type_t rtype)
+test_query_read_selection(size_t file_count, const char *filenames[],
+ hid_t *files, hid_t view, H5R_type_t rtype)
{
hid_t refs = H5I_BADID, ref_type = H5I_BADID, ref_space = H5I_BADID;
size_t n_refs, ref_size, ref_buf_size;
- void *ref_buf= NULL;
+ void *ref_buf = NULL;
href_t *ref_ptr = NULL;
- const char *ref_path;
+ const char *ref_path = NULL;
hid_t obj = H5I_BADID, type = H5I_BADID, space = H5I_BADID, mem_space = H5I_BADID;
size_t n_elem, elem_size, buf_size;
float *buf = NULL;
@@ -459,7 +460,7 @@ test_query_read_selection(hid_t file, hid_t fapl, hid_t view, H5R_type_t rtype)
if (0 == (n_refs = (size_t) H5Sget_select_npoints(ref_space))) FAIL_STACK_ERROR;
printf("Found %zu reference(s)\n", n_refs);
if (0 == (ref_size = H5Tget_size(ref_type))) FAIL_STACK_ERROR;
- printf("Reference type size: %d\n", ref_size);
+ printf("Reference type size: %zu\n", ref_size);
/* Allocate buffer to hold data */
ref_buf_size = n_refs * ref_size;
@@ -474,21 +475,28 @@ test_query_read_selection(hid_t file, hid_t fapl, hid_t view, H5R_type_t rtype)
char filename[MAX_NAME];
hid_t loc = H5I_BADID;
- if (file == H5I_BADID) {
- if (H5Rget_filename(ref_ptr, filename, MAX_NAME) < 0) FAIL_STACK_ERROR;
- printf("Found reference from file: %s\n", filename);
- if ((loc = H5Fopen(filename, H5F_ACC_RDONLY, fapl)) < 0) FAIL_STACK_ERROR;
+ if (H5Rget_file_name(ref_ptr[i], filename, MAX_NAME) < 0) FAIL_STACK_ERROR;
+ printf("Found reference from file: %s\n", filename);
+ if (file_count > 1) {
+ unsigned int j;
+ for (j = 0; j < file_count; j++) {
+ if (0 == HDstrcmp(filename, filenames[j])) {
+ loc = files[j];
+ break;
+ }
+ }
} else {
- loc = file;
+ if (0 != HDstrcmp(filename, filenames[0])) FAIL_STACK_ERROR;
+ loc = files[0];
}
- if ((obj = H5Rdereference(loc, H5P_DEFAULT, ref_ptr)) < 0) FAIL_STACK_ERROR;
- if (H5Iget_name(obj, obj_path, MAX_NAME) < 0) FAIL_STACK_ERROR;
+ if (H5Rget_obj_name(loc, ref_ptr[i], obj_path, MAX_NAME) < 0) FAIL_STACK_ERROR;
printf("Found reference from object: %s\n", obj_path);
+ if ((obj = H5Rget_object(loc, H5P_DEFAULT, ref_ptr[i])) < 0) FAIL_STACK_ERROR;
if (rtype == H5R_REGION) {
- int j;
+ unsigned int j;
- if ((space = H5Rget_region(loc, ref_ptr)) < 0) FAIL_STACK_ERROR;
+ if ((space = H5Rget_region2(loc, ref_ptr[i])) < 0) FAIL_STACK_ERROR;
if ((type = H5Dget_type(obj)) < 0) FAIL_STACK_ERROR;
if (0 == (n_elem = (size_t) H5Sget_select_npoints(space))) FAIL_STACK_ERROR;
if (0 == (elem_size = H5Tget_size(type))) FAIL_STACK_ERROR;
@@ -518,19 +526,14 @@ test_query_read_selection(hid_t file, hid_t fapl, hid_t view, H5R_type_t rtype)
if (rtype == H5R_ATTR) {
char attr_name[MAX_NAME];
- if (H5Rget_name(obj, ref_ptr, attr_name, MAX_NAME) < 0) FAIL_STACK_ERROR;
+ if (H5Rget_attr_name(obj, ref_ptr[i], attr_name, MAX_NAME) < 0) FAIL_STACK_ERROR;
printf("Attribute name: %s\n", attr_name);
- if (H5Aclose(obj) < 0) FAIL_STACK_ERROR;
- } else {
- if (H5Dclose(obj) < 0) FAIL_STACK_ERROR;
}
- if ((file == H5I_BADID) && (H5Fclose(loc) < 0)) FAIL_STACK_ERROR;
- ref_ptr = (href_t *)((uint8_t *) ref_ptr + ref_size);
+ if (H5Dclose(obj) < 0) FAIL_STACK_ERROR;
}
- if ((rtype == H5R_ATTR) || (rtype == H5R_REGION))
- if ((H5Dref_reclaim(ref_type, ref_space, H5P_DEFAULT, ref_buf)) < 0) FAIL_STACK_ERROR;
+ if ((H5Dref_reclaim(ref_type, ref_space, H5P_DEFAULT, ref_buf)) < 0) FAIL_STACK_ERROR;
if (H5Sclose(ref_space) < 0) FAIL_STACK_ERROR;
if (H5Tclose(ref_type) < 0) FAIL_STACK_ERROR;
if (H5Dclose(refs) < 0) FAIL_STACK_ERROR;
@@ -590,7 +593,7 @@ test_query_apply_view(const char *filename, hid_t fapl, unsigned idx_plugin)
+ ((float) (t2.tv_usec - t1.tv_usec)) / 1000.0f);
if (!(result & H5Q_REF_REG)) FAIL_STACK_ERROR;
- if (test_query_read_selection(file, fapl, view, H5R_REGION) < 0) FAIL_STACK_ERROR;
+ if (test_query_read_selection(1, &filename, &file, view, H5R_REGION) < 0) FAIL_STACK_ERROR;
if (H5Gclose(view) < 0) FAIL_STACK_ERROR;
if (test_query_close(query)) FAIL_STACK_ERROR;
@@ -603,7 +606,7 @@ test_query_apply_view(const char *filename, hid_t fapl, unsigned idx_plugin)
if ((view = H5Qapply(file, query, &result, H5P_DEFAULT)) < 0) FAIL_STACK_ERROR;
if (!(result & H5Q_REF_OBJ)) FAIL_STACK_ERROR;
- if (test_query_read_selection(file, fapl, view, H5R_OBJECT) < 0) FAIL_STACK_ERROR;
+ if (test_query_read_selection(1, &filename, &file, view, H5R_OBJECT) < 0) FAIL_STACK_ERROR;
if (H5Gclose(view) < 0) FAIL_STACK_ERROR;
if (test_query_close(query)) FAIL_STACK_ERROR;
@@ -616,7 +619,7 @@ test_query_apply_view(const char *filename, hid_t fapl, unsigned idx_plugin)
if ((view = H5Qapply(file, query, &result, H5P_DEFAULT)) < 0) FAIL_STACK_ERROR;
if (!(result & H5Q_REF_ATTR)) FAIL_STACK_ERROR;
- if (test_query_read_selection(file, fapl, view, H5R_ATTR) < 0) FAIL_STACK_ERROR;
+ if (test_query_read_selection(1, &filename, &file, view, H5R_ATTR) < 0) FAIL_STACK_ERROR;
if (H5Gclose(view) < 0) FAIL_STACK_ERROR;
if (test_query_close(query)) FAIL_STACK_ERROR;
@@ -638,23 +641,23 @@ error:
}
static herr_t
-test_query_apply_view_multi(const char filename[MULTI_NFILES][MAX_NAME], hid_t fapl, unsigned idx_plugin)
+test_query_apply_view_multi(const char *filenames[], hid_t fapl, unsigned idx_plugin)
{
hid_t files[MULTI_NFILES] = {H5I_BADID, H5I_BADID, H5I_BADID};
hid_t view = H5I_BADID;
hid_t query = H5I_BADID;
struct timeval t1, t2;
unsigned result = 0;
- int i;
+ unsigned int i;
printf(" ...\n---\n");
/* Create simple files for testing queries */
for (i = 0; i < MULTI_NFILES; i++) {
- printf("Creating test file \"%s\"\n", filename[i]);
- if ((test_query_create_simple_file(filename[i], fapl, idx_plugin)) < 0) FAIL_STACK_ERROR;
+ printf("Creating test file \"%s\"\n", filenames[i]);
+ if ((test_query_create_simple_file(filenames[i], fapl, idx_plugin)) < 0) FAIL_STACK_ERROR;
/* Open the file in read-only */
- if ((files[i] = H5Fopen(filename[i], H5F_ACC_RDONLY, fapl)) < 0) FAIL_STACK_ERROR;
+ if ((files[i] = H5Fopen(filenames[i], H5F_ACC_RDONLY, fapl)) < 0) FAIL_STACK_ERROR;
}
printf("\nRegion query\n");
@@ -674,7 +677,7 @@ test_query_apply_view_multi(const char filename[MULTI_NFILES][MAX_NAME], hid_t f
+ ((float) (t2.tv_usec - t1.tv_usec)) / 1000.0f);
if (!(result & H5Q_REF_REG)) FAIL_STACK_ERROR;
- if (test_query_read_selection(H5I_BADID, fapl, view, H5R_REGION) < 0) FAIL_STACK_ERROR;
+ if (test_query_read_selection(MULTI_NFILES, filenames, files, view, H5R_REGION) < 0) FAIL_STACK_ERROR;
if (H5Gclose(view) < 0) FAIL_STACK_ERROR;
if (test_query_close(query)) FAIL_STACK_ERROR;
@@ -687,7 +690,7 @@ test_query_apply_view_multi(const char filename[MULTI_NFILES][MAX_NAME], hid_t f
if ((view = H5Qapply_multi(MULTI_NFILES, files, query, &result, H5P_DEFAULT)) < 0) FAIL_STACK_ERROR;
if (!(result & H5Q_REF_OBJ)) FAIL_STACK_ERROR;
- if (test_query_read_selection(H5I_BADID, fapl, view, H5R_OBJECT) < 0) FAIL_STACK_ERROR;
+ if (test_query_read_selection(MULTI_NFILES, filenames, files, view, H5R_OBJECT) < 0) FAIL_STACK_ERROR;
if (H5Gclose(view) < 0) FAIL_STACK_ERROR;
if (test_query_close(query)) FAIL_STACK_ERROR;
@@ -700,7 +703,7 @@ test_query_apply_view_multi(const char filename[MULTI_NFILES][MAX_NAME], hid_t f
if ((view = H5Qapply_multi(MULTI_NFILES, files, query, &result, H5P_DEFAULT)) < 0) FAIL_STACK_ERROR;
if (!(result & H5Q_REF_ATTR)) FAIL_STACK_ERROR;
- if (test_query_read_selection(H5I_BADID, fapl, view, H5R_ATTR) < 0) FAIL_STACK_ERROR;
+ if (test_query_read_selection(MULTI_NFILES, filenames, files, view, H5R_ATTR) < 0) FAIL_STACK_ERROR;
if (H5Gclose(view) < 0) FAIL_STACK_ERROR;
if (test_query_close(query)) FAIL_STACK_ERROR;
@@ -730,8 +733,9 @@ main(void)
#ifdef H5_HAVE_FASTBIT
char filename_fastbit[MAX_NAME];
#endif
- char filename_multi[MULTI_NFILES][MAX_NAME];
+ char **filename_multi = NULL;
hid_t query = H5I_BADID, fapl = H5I_BADID;
+ int i;
/* Reset library */
h5_reset();
@@ -741,9 +745,14 @@ main(void)
#ifdef H5_HAVE_FASTBIT
h5_fixname(FILENAME[1], fapl, filename_fastbit, sizeof(filename_fastbit));
#endif
- h5_fixname(FILENAME[2], fapl, filename_multi[0], sizeof(filename_multi[0]));
- h5_fixname(FILENAME[3], fapl, filename_multi[1], sizeof(filename_multi[1]));
- h5_fixname(FILENAME[4], fapl, filename_multi[2], sizeof(filename_multi[2]));
+ filename_multi = (char **) HDmalloc(MULTI_NFILES * sizeof(char *));
+ for (i = 0; i < MULTI_NFILES; i++) {
+ filename_multi[i] = (char *) HDmalloc(MAX_NAME);
+ HDmemset(filename_multi[i], '\0', MAX_NAME);
+ }
+ h5_fixname(FILENAME[2], fapl, filename_multi[0], MAX_NAME);
+ h5_fixname(FILENAME[3], fapl, filename_multi[1], MAX_NAME);
+ h5_fixname(FILENAME[4], fapl, filename_multi[2], MAX_NAME);
/* Check that no object is left open */
H5Pset_fclose_degree(fapl, H5F_CLOSE_SEMI);
@@ -797,6 +806,9 @@ main(void)
puts("All query tests passed.");
h5_cleanup(FILENAME, fapl);
+ for (i = 0; i < MULTI_NFILES; i++)
+ HDfree(filename_multi[i]);
+ HDfree(filename_multi);
return 0;