aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authortradke <tradke@7842ec3a-9562-4be5-9c5b-06ba18f2b668>2002-07-05 11:33:14 +0000
committertradke <tradke@7842ec3a-9562-4be5-9c5b-06ba18f2b668>2002-07-05 11:33:14 +0000
commit48bf661bfedaf8b7dbf6f46378874495c860e80f (patch)
tree7bb1b1429fb3b1b04c46af46f516d3b2033da3bc
parenta4587ec295e073276c1f4ef6b27b6664b712adf5 (diff)
Fixed a resource leak in the recovery routines.
Also print the iteration number of a recovered variable. git-svn-id: http://svn.cactuscode.org/arrangements/CactusPUGHIO/IOHDF5Util/trunk@74 7842ec3a-9562-4be5-9c5b-06ba18f2b668
-rw-r--r--src/RecoverVar.c118
1 files changed, 62 insertions, 56 deletions
diff --git a/src/RecoverVar.c b/src/RecoverVar.c
index b0d80bc..29a3390 100644
--- a/src/RecoverVar.c
+++ b/src/RecoverVar.c
@@ -721,7 +721,7 @@ static herr_t processDataset (hid_t group, const char *datasetname, void *arg)
{
ioGH *ioUtilGH;
ioHDF5UtilGH *myGH;
- int vindex, vtype, gtype, timelevel, iteration, is_group;
+ int vindex, vtype, gtype, timelevel, iteration, is_group, retval;
iterate_info_t *it_info = (iterate_info_t *) arg;
recover_info_t rec_info;
hid_t dataset;
@@ -737,6 +737,8 @@ static herr_t processDataset (hid_t group, const char *datasetname, void *arg)
return (0);
}
+ retval = 0;
+
HDF5_ERROR (H5Gget_objinfo (group, datasetname, 0, &object_info));
is_group = object_info.type == H5G_GROUP;
if (is_group)
@@ -755,66 +757,73 @@ static herr_t processDataset (hid_t group, const char *datasetname, void *arg)
{
CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
"Ignoring dataset '%s'", datasetname);
- return (0);
}
-
- ioUtilGH = (ioGH *) CCTK_GHExtension (it_info->GH, "IO");
- myGH = (ioHDF5UtilGH *) CCTK_GHExtension (it_info->GH, "IOHDF5Util");
-
- /* if we read in initial data via the file reader interface
- check whether the user wants to have this variable with a specific
- iteration number restored */
- if (ioUtilGH->do_inVars && ioUtilGH->do_inVars[vindex] >= 0 &&
- ioUtilGH->do_inVars[vindex] != iteration + 1)
+ else
{
- if (CCTK_Equals (verbose, "full"))
+ ioUtilGH = (ioGH *) CCTK_GHExtension (it_info->GH, "IO");
+ myGH = (ioHDF5UtilGH *) CCTK_GHExtension (it_info->GH, "IOHDF5Util");
+
+ /* if we read in initial data via the file reader interface
+ check whether the user wants to have this variable with a specific
+ iteration number restored */
+ if (ioUtilGH->do_inVars && ioUtilGH->do_inVars[vindex] >= 0 &&
+ ioUtilGH->do_inVars[vindex] != iteration + 1)
{
- CCTK_VInfo (CCTK_THORNSTRING, "Ignoring dataset '%s' for file reader "
- "recovery", datasetname);
+ if (CCTK_Equals (verbose, "full"))
+ {
+ CCTK_VInfo (CCTK_THORNSTRING, "Ignoring dataset '%s' for file reader "
+ "recovery", datasetname);
+ }
}
- return (0);
- }
-
- if (CCTK_Equals (verbose, "full"))
- {
- fullname = CCTK_FullName (vindex);
- CCTK_VInfo (CCTK_THORNSTRING, "Restoring variable '%s' (timelevel %d)",
- fullname, timelevel);
- free (fullname);
- }
+ else
+ {
+ if (CCTK_Equals (verbose, "full"))
+ {
+ fullname = CCTK_FullName (vindex);
+ CCTK_VInfo (CCTK_THORNSTRING, "Restoring variable '%s' (timelevel %d, "
+ "cctk_iteration %d)", fullname, timelevel, iteration);
+ free (fullname);
+ }
- vtype = CCTK_VarTypeI (vindex);
- rec_info.it_info = it_info;
- rec_info.element_size = CCTK_VarTypeSize (vtype);
+ vtype = CCTK_VarTypeI (vindex);
+ rec_info.it_info = it_info;
+ rec_info.element_size = CCTK_VarTypeSize (vtype);
#ifdef CCTK_MPI
- rec_info.mpi_type = PUGH_MPIDataType (PUGH_pGH (it_info->GH), vtype);
+ rec_info.mpi_type = PUGH_MPIDataType (PUGH_pGH (it_info->GH), vtype);
#endif
- rec_info.hdf5type = IOHDF5Util_DataType (myGH, vtype);
- if (rec_info.element_size <= 0 ||
+ rec_info.hdf5type = IOHDF5Util_DataType (myGH, vtype);
+ if (rec_info.element_size <= 0 ||
#ifdef CCTK_MPI
- rec_info.mpi_type == MPI_DATATYPE_NULL ||
+ rec_info.mpi_type == MPI_DATATYPE_NULL ||
#endif
- rec_info.hdf5type < 0)
- {
- CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
- "Unsupported variable datatype %d", vtype);
- return (1);
- }
-
- /* Read in the data */
- switch (gtype)
- {
- case CCTK_SCALAR:
- IOHDF5Util_RestoreGS (dataset, vindex, timelevel, &rec_info);
- break;
- case CCTK_GF:
- case CCTK_ARRAY:
- IOHDF5Util_RestoreGA (dataset, vindex, timelevel, &rec_info);
- break;
- default:
- CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
- "Unsupported group type %d", gtype);
- return (1);
+ rec_info.hdf5type < 0)
+ {
+ CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
+ "Unsupported variable datatype %d", vtype);
+ retval = 1;
+ }
+ else if (gtype != CCTK_SCALAR && gtype != CCTK_GF && gtype != CCTK_ARRAY)
+ {
+ CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
+ "Unsupported group type %d", gtype);
+ retval = 1;
+ }
+ else
+ {
+ /* Read in the data */
+ if (gtype == CCTK_SCALAR)
+ {
+ IOHDF5Util_RestoreGS (dataset, vindex, timelevel, &rec_info);
+ }
+ else
+ {
+ IOHDF5Util_RestoreGA (dataset, vindex, timelevel, &rec_info);
+ }
+
+ /* increment counter for total number of recovered variables */
+ it_info->num_recovered++;
+ }
+ }
}
if (is_group)
@@ -826,8 +835,5 @@ static herr_t processDataset (hid_t group, const char *datasetname, void *arg)
HDF5_ERROR (H5Dclose (dataset));
}
- /* increment counter for total number of recovered variables */
- it_info->num_recovered++;
-
- return (0);
+ return (retval);
}