/*@@ @file DumpUtils.c @date Fri Oct 6 2000 @author Thomas Radke @desc Utility routines for dumping data into HDF5 files. @enddesc @version $Id$ @@*/ #include #include "cctk.h" #include "cctk_Parameters.h" #include "CactusBase/IOUtil/src/ioGH.h" #include "CactusBase/IOUtil/src/ioutil_CheckpointRecovery.h" #include "ioHDF5UtilGH.h" /* the rcs ID and its dummy function to use it */ static char *rcsid = "$Id$"; CCTK_FILEVERSION(BetaThorns_IOHDF5Util_DumpUtils_c) /*@@ @routine IOHDF5Util_DumpGH @date Tue Oct 10 2000 @author Thomas Radke @desc Dump all variables of the given GH along with the global attributes and parameters to the (already opened) HDF5 checkpoint file. @enddesc @var GH @vdesc Pointer to CCTK grid hierarchy @vtype cGH * @vio in @endvar @var timers @vdesc array of timers for timing the checkpointing @vtype int [] @vio in @endvar @var file @vdesc HDF5 file handle @vtype hid_t @vio in @endvar @returntype int @returndesc always 0 (for success) @endreturndesc @@*/ int IOHDF5Util_DumpGH (cGH *GH, int timers[], hid_t file) { DECLARE_CCTK_PARAMETERS int index; int maxdim; int timelevel, current_timelevel; int *old_downsample; int old_out_single; ioGH *ioUtilGH; const char *thorn; const char *impl; ioHDF5Geo_t slab; /* Get the GH extension for IOUtil */ ioUtilGH = (ioGH *) CCTK_GHExtension (GH, "IO"); /* disable downsampling after saving original downsampling params */ maxdim = CCTK_MaxDim (); old_downsample = (int *) malloc (maxdim * sizeof (int)); for (index = 0; index < maxdim; index++) { old_downsample[index] = ioUtilGH->downsample[index]; ioUtilGH->downsample[index] = 1; slab.direction[index] = index; slab.downsample[index] = 1; } /* disable output in single precision */ old_out_single = ioUtilGH->out_single; ioUtilGH->out_single = 0; /* sync all active groups */ for (index = CCTK_NumGroups () - 1; index >= 0; index--) { if (CCTK_IsImplementationActive ( CCTK_ImpFromVarI (CCTK_FirstVarIndexI (index)))) { CCTK_SyncGroupI (GH, index); } } /* start CP_PARAMETERS_TIMER timer */ if (timers) { CCTK_TimerResetI (timers[CP_PARAMETERS_TIMER]); CCTK_TimerStartI (timers[CP_PARAMETERS_TIMER]); } /* Great; Now start dumping away! */ if (file >= 0) { /* all GH extension variables and parameter stuff which is not specific to any dataset goes into dedicated groups */ if (out3D_parameters) { IOHDF5Util_DumpParameters (GH, file); } IOHDF5Util_DumpGHExtensions (GH, file); } if (verbose) { CCTK_INFO ("Dumping variables ..."); } /* stop CP_PARAMETERS_TIMER timer and start CP_VARIABLES_TIMER */ if (timers) { CCTK_TimerStopI (timers[CP_PARAMETERS_TIMER]); CCTK_TimerResetI (timers[CP_VARIABLES_TIMER]); CCTK_TimerStartI (timers[CP_VARIABLES_TIMER]); } /* ... now the variables */ for (index = CCTK_NumVars () - 1; index >= 0; index--) { /* find out the thorn implementing variable with index */ impl = CCTK_ImpFromVarI (index); thorn = CCTK_ImplementationThorn (impl); if (! thorn) { thorn = impl; } /* let only variables pass which belong to an active thorn and have storage assigned */ if (! CCTK_IsThornActive (thorn) || ! CCTK_QueryGroupStorageI (GH, CCTK_GroupIndexFromVarI (index))) { continue; } if (verbose && file >= 0) { CCTK_VInfo (CCTK_THORNSTRING, " %s", CCTK_VarName (index)); } /* get the current timelevel */ current_timelevel = CCTK_NumTimeLevelsFromVarI (index) - 1; if (current_timelevel > 0) { current_timelevel--; } slab.sdim = slab.vdim = CCTK_GroupDimFromVarI (index); /* now dump all timelevels up to the current */ for (timelevel = 0; timelevel <= current_timelevel; timelevel++) { IOHDF5Util_DumpVar (GH, index, timelevel, &slab, file, 0); } } /* end of loop over all variables */ /* stop CP_VARIABLES_TIMER timer */ if (timers) { CCTK_TimerStopI (timers[CP_VARIABLES_TIMER]); } /* restore output precision flag */ ioUtilGH->out_single = old_out_single; /* restore original downsampling params */ memcpy (ioUtilGH->downsample, old_downsample, maxdim * sizeof (int)); free (old_downsample); return (0); } /*@@ @routine IOHDF5Util_DumpCommonAttributes @date May 21 1999 @author Thomas Radke @desc Add "Common" attributes to the given dataset, these are: @enddesc @var GH @vdesc Pointer to CCTK grid hierarchy @vtype cGH * @vio in @endvar @var index @vdesc CCTK index of the variable to dump @vtype int @vio in @endvar @var timelevel @vdesc timelevel of the variable to dump @vtype int @vio in @endvar @var global_shape @vdesc global shape of the grid @vtype CCTK_INT [] @vio in @endvar @var slab @vdesc pointer to hyperslab structure describing the dataset @vtype ioHDF5Geo_t * @vio in @endvar @var dataset @vdesc the dataset handle where the attributes should be attached to @vtype hid_t @vio in @endvar @@*/ void IOHDF5Util_DumpCommonAttributes (cGH *GH, int index, int timelevel, CCTK_INT global_shape[], ioHDF5Geo_t *slab, hid_t dataset) { DECLARE_CCTK_PARAMETERS int dim; char *groupname; CCTK_INT attr_int; CCTK_REAL *attr_real; ioHDF5UtilGH *myGH; /* Get the handle for IOHDF5Util extensions */ myGH = (ioHDF5UtilGH *) CCTK_GHExtension (GH, "IOHDF5Util"); /* attributes describing the variable */ dim = CCTK_GroupDimI (CCTK_GroupIndexFromVarI (index)); groupname = CCTK_GroupNameFromVarI (index); WRITE_ATTRIBUTE ("groupname", groupname, dataset, myGH->scalar_dataspace, 0, myGH->IOHDF5_STRING); free (groupname); attr_int = CCTK_GroupTypeFromVarI (index); WRITE_ATTRIBUTE ("grouptype", &attr_int, dataset, myGH->scalar_dataspace, 0, IOHDF5_INT); attr_int = CCTK_NumTimeLevelsFromVarI (index); WRITE_ATTRIBUTE ("ntimelevels", &attr_int, dataset, myGH->scalar_dataspace, 0, IOHDF5_INT); WRITE_ATTRIBUTE ("global_size", global_shape, dataset, myGH->array_dataspace, dim, IOHDF5_INT); WRITE_ATTRIBUTE ("time", &GH->cctk_time, dataset, myGH->scalar_dataspace, 0, IOHDF5_REAL); /* attributes describing the underlying grid */ dim = CCTK_MaxDim (); attr_real = (CCTK_REAL *) malloc (2 * dim * sizeof (CCTK_REAL)); CCTK_CoordRange (GH, &attr_real[0], &attr_real[dim + 0], -1, "x", "cart3d"); CCTK_CoordRange (GH, &attr_real[1], &attr_real[dim + 1], -1, "y", "cart3d"); CCTK_CoordRange (GH, &attr_real[2], &attr_real[dim + 2], -1, "z", "cart3d"); WRITE_ATTRIBUTE ("origin", attr_real, dataset, myGH->array_dataspace, dim, IOHDF5_REAL); WRITE_ATTRIBUTE ("min_ext", attr_real, dataset, myGH->array_dataspace, dim, IOHDF5_REAL); WRITE_ATTRIBUTE ("max_ext", attr_real + dim, dataset, myGH->array_dataspace, dim, IOHDF5_REAL); WRITE_ATTRIBUTE ("delta", GH->cctk_delta_space, dataset, myGH->array_dataspace, dim, IOHDF5_REAL); free (attr_real); /* attributes describing the hyperslab */ if (slab) { attr_real = (CCTK_REAL *) malloc (4 * slab->sdim * sizeof (CCTK_REAL)); for (dim = 0; dim < slab->sdim; dim++) { attr_real[dim + 0*slab->sdim] = slab->origin[slab->direction[dim]] * GH->cctk_delta_space[slab->direction[dim]]; attr_real[dim + 1*slab->sdim] = slab->origin[dim] * GH->cctk_delta_space[dim]; attr_real[dim + 2*slab->sdim] = (slab->origin[dim] + slab->actlen[dim]-1) * GH->cctk_delta_space[dim] * slab->downsample[dim]; attr_real[dim + 3*slab->sdim] = GH->cctk_delta_space[slab->direction[dim]] * slab->downsample[slab->direction[dim]]; } WRITE_ATTRIBUTE ("origin_slab", attr_real + 0*dim, dataset, myGH->array_dataspace, dim, IOHDF5_REAL); WRITE_ATTRIBUTE ("min_ext_slab", attr_real + 1*dim, dataset, myGH->array_dataspace, dim, IOHDF5_REAL); WRITE_ATTRIBUTE ("max_ext_slab", attr_real + 2*dim, dataset, myGH->array_dataspace, dim, IOHDF5_REAL); WRITE_ATTRIBUTE ("delta_slab", attr_real + 3*dim, dataset, myGH->array_dataspace, dim, IOHDF5_REAL); free (attr_real); } } /*@@ @routine IOHDF5Util_DumpParameters @date Thu Jan 20 2000 @author Thomas Radke @desc Collects the parameters of all active implementations into a single string and writes it as an attribute attached to the CACTUS_PARAMETERS_GROUP group in the HDF5 file. @enddesc @var GH @vdesc Pointer to CCTK grid hierarchy @vtype cGH * @vio in @endvar @var file @vdesc the HDF5 file handle @vtype hid_t @vio in @endvar @@*/ void IOHDF5Util_DumpParameters (cGH *GH, hid_t file) { DECLARE_CCTK_PARAMETERS ioHDF5UtilGH *myGH; char *parameters; hid_t group; if (verbose) { CCTK_INFO ("Dumping GH extensions ..."); } /* Get the parameter string buffer */ parameters = IOUtil_GetAllParameters (GH); if (parameters) { /* Get the handle for IOHDF5Util extensions */ myGH = (ioHDF5UtilGH *) CCTK_GHExtension (GH, "IOHDF5Util"); IOHDF5_ERROR (group = H5Gcreate (file, CACTUS_PARAMETERS_GROUP, 0)); WRITE_ATTRIBUTE (ALL_PARAMETERS, parameters, group, myGH->scalar_dataspace, 0, myGH->IOHDF5_STRING); IOHDF5_ERROR (H5Gclose (group)); free (parameters); } } /*@@ @routine IOHDF5Util_DumpGHExtensions @date Thu Jan 20 2000 @author Thomas Radke @desc Attaches the current values of important flesh and IO variables as attributes to the GLOBAL_ATTRIBUTES_GROUP group in the HDF5 file. @enddesc @var GH @vdesc Pointer to CCTK grid hierarchy @vtype cGH * @vio in @endvar @var file @vdesc the HDF5 file handle @vtype hid_t @vio in @endvar @@*/ void IOHDF5Util_DumpGHExtensions (cGH *GH, hid_t file) { DECLARE_CCTK_PARAMETERS int value; hid_t group; ioGH *ioUtilGH; ioHDF5UtilGH *myGH; if (verbose) { CCTK_INFO ("Dumping GH extensions ..."); } /* Get the handles for IOUtil and IOHDF5 extensions */ ioUtilGH = (ioGH *) CCTK_GHExtension (GH, "IO"); myGH = (ioHDF5UtilGH *) CCTK_GHExtension (GH, "IOHDF5Util"); IOHDF5_ERROR (group = H5Gcreate (file, GLOBAL_ATTRIBUTES_GROUP, 0)); value = CCTK_MainLoopIndex (); WRITE_ATTRIBUTE ("main_loop_index", &value, group, myGH->scalar_dataspace, 0, H5T_NATIVE_INT); value = CCTK_nProcs (GH); WRITE_ATTRIBUTE ("nprocs", &value, group, myGH->scalar_dataspace, 0, H5T_NATIVE_INT); WRITE_ATTRIBUTE ("ioproc_every", &ioUtilGH->ioproc_every, group, myGH->scalar_dataspace, 0, H5T_NATIVE_INT); WRITE_ATTRIBUTE ("unchunked", &ioUtilGH->unchunked, group, myGH->scalar_dataspace, 0, H5T_NATIVE_INT); WRITE_ATTRIBUTE ("cctk_time", &GH->cctk_time, group, myGH->scalar_dataspace, 0, IOHDF5_REAL); WRITE_ATTRIBUTE ("cctk_iteration", &GH->cctk_iteration, group, myGH->scalar_dataspace, 0, H5T_NATIVE_INT); IOHDF5_ERROR (H5Gclose (group)); }