aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorrhaas <rhaas@7842ec3a-9562-4be5-9c5b-06ba18f2b668>2014-02-13 22:09:01 +0000
committerrhaas <rhaas@7842ec3a-9562-4be5-9c5b-06ba18f2b668>2014-02-13 22:09:01 +0000
commit2e45c05ece97c145ae10b6fa9d94da3bfa7333d8 (patch)
tree04d31ed2790ab9670e62c31a734ab15115911e92
parent0bc3455b35874c7348a0a2302d1e572b52ac6aa7 (diff)
eckpoint cctk_delta_timesvn
the PUGH driver currently did not checkpoint the stepsize which makes using it with adaptive timestepping hard. Carpet on the other hand checkpoints the stepsize. This patch writes cckt_delta_time into checkpoints. When reading a checkpoint file that does not set the timestep, the user will see level 1 warning that the attribute cctk_delta_time was not found but the code will not abort and will behave as it does right now otherwise (which is to use cctk_delta_time based on the values from parameters/set in basegrid). git-svn-id: http://svn.cactuscode.org/arrangements/CactusPUGHIO/IOHDF5Util/trunk@174 7842ec3a-9562-4be5-9c5b-06ba18f2b668
-rw-r--r--src/DumpUtils.c2
-rw-r--r--src/RecoverVar.c15
2 files changed, 13 insertions, 4 deletions
diff --git a/src/DumpUtils.c b/src/DumpUtils.c
index ad99e00..d559db7 100644
--- a/src/DumpUtils.c
+++ b/src/DumpUtils.c
@@ -429,6 +429,8 @@ void IOHDF5Util_DumpGHExtensions (const cGH *GH, hid_t file)
WRITE_ATTRIBUTE ("unchunked", &ioUtilGH->unchunked, group, myGH, 0,
H5T_NATIVE_INT);
WRITE_ATTRIBUTE ("cctk_time", &GH->cctk_time, group, myGH, 0, HDF5_REAL);
+ WRITE_ATTRIBUTE ("cctk_delta_time", &GH->cctk_delta_time, group, myGH, 0,
+ HDF5_REAL);
WRITE_ATTRIBUTE ("cctk_iteration", &GH->cctk_iteration, group, myGH, 0,
H5T_NATIVE_INT);
version = CCTK_FullVersion ();
diff --git a/src/RecoverVar.c b/src/RecoverVar.c
index 19a6cb3..ace099c 100644
--- a/src/RecoverVar.c
+++ b/src/RecoverVar.c
@@ -193,7 +193,7 @@ int IOHDF5Util_RecoverVariables (cGH *GH, const fileinfo_t *fileinfo)
int IOHDF5Util_RecoverGHextensions (cGH *GH, const fileinfo_t *fileinfo)
{
hid_t group;
- CCTK_REAL realBuffer;
+ CCTK_REAL realBuffer[2];
CCTK_INT4 int4Buffer[3];
@@ -208,7 +208,13 @@ int IOHDF5Util_RecoverGHextensions (cGH *GH, const fileinfo_t *fileinfo)
{
READ_ATTRIBUTE (group, "cctk_iteration", HDF5_INT4, &int4Buffer[1]);
READ_ATTRIBUTE (group, "main_loop_index", HDF5_INT4, &int4Buffer[2]);
- READ_ATTRIBUTE (group, "cctk_time", HDF5_REAL, &realBuffer);
+ READ_ATTRIBUTE (group, "cctk_time", HDF5_REAL, &realBuffer[0]);
+ // pre rev 173 checkpoints do not contain cctk_delta_time. When reading
+ // such a checkpoint, we output a warning and continue to use the old
+ // cctk_delta_time. For this we need to initialize realBuffer[1] since we
+ // cannot tell if READ_ATTRIBUTE succeeded or not.
+ realBuffer[1] = GH->cctk_delta_time;
+ READ_ATTRIBUTE (group, "cctk_delta_time", HDF5_REAL, &realBuffer[1]);
HDF5_ERROR (H5Gclose (group));
}
@@ -227,14 +233,15 @@ int IOHDF5Util_RecoverGHextensions (cGH *GH, const fileinfo_t *fileinfo)
CACTUS_MPI_ERROR (MPI_Bcast (int4Buffer, 3, PUGH_MPI_INT4, 0,MPI_COMM_WORLD));
if (int4Buffer[0])
{
- CACTUS_MPI_ERROR (MPI_Bcast (&realBuffer, 1, PUGH_MPI_REAL, 0,
+ CACTUS_MPI_ERROR (MPI_Bcast (realBuffer, 2, PUGH_MPI_REAL, 0,
MPI_COMM_WORLD));
}
#endif
if (int4Buffer[0])
{
- GH->cctk_time = realBuffer;
+ GH->cctk_time = realBuffer[0];
+ GH->cctk_delta_time = realBuffer[1];
GH->cctk_iteration = int4Buffer[1];
CCTK_SetMainLoopIndex ((int) int4Buffer[2]);
}