aboutsummaryrefslogtreecommitdiff
path: root/CarpetAttic
diff options
context:
space:
mode:
authorcott <>2004-01-09 14:43:00 +0000
committercott <>2004-01-09 14:43:00 +0000
commit73de628ef31f00b9a42476e9fdb6168e09c60c45 (patch)
tree31768081c5ca11d8ef64d01efb7d2c9d1ed91d2a /CarpetAttic
parentfdabf8e0b555306b25527e4d685e315c2ef80784 (diff)
a few changes that help with recovery
darcs-hash:20040109144346-19929-67096355895264e0861c3792c7bd3893ee932801.gz
Diffstat (limited to 'CarpetAttic')
-rw-r--r--CarpetAttic/CarpetIOFlexIOCheckpoint/src/checkpointrestart.cc113
-rw-r--r--CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexio.cc18
-rw-r--r--CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexioUtil.cc8
3 files changed, 121 insertions, 18 deletions
diff --git a/CarpetAttic/CarpetIOFlexIOCheckpoint/src/checkpointrestart.cc b/CarpetAttic/CarpetIOFlexIOCheckpoint/src/checkpointrestart.cc
index d58289594..f92cd526e 100644
--- a/CarpetAttic/CarpetIOFlexIOCheckpoint/src/checkpointrestart.cc
+++ b/CarpetAttic/CarpetIOFlexIOCheckpoint/src/checkpointrestart.cc
@@ -42,13 +42,14 @@
#include "ggf.hh"
#include "vect.hh"
+
//#include "StoreNamedData.h"
#include "carpet.hh"
#include "ioflexio.hh"
extern "C" {
- static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/CarpetAttic/CarpetIOFlexIOCheckpoint/src/checkpointrestart.cc,v 1.20 2004/01/08 19:43:33 cott Exp $";
+ static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/CarpetAttic/CarpetIOFlexIOCheckpoint/src/checkpointrestart.cc,v 1.21 2004/01/09 15:43:46 cott Exp $";
CCTK_FILEVERSION(Carpet_CarpetIOFlexIO_checkpointrestart_cc);
}
@@ -234,10 +235,46 @@ int CarpetIOFlexIO_Recover (cGH* cgh, const char *basefilename, int called_from)
if(myproc==0)
amrreader = new AmrGridReader(*reader);
+
+
+ CCTK_INT4 numberoftimes;
+ IObase::DataType datatype;
+ int i,dim;
+ if (myproc == 0) {
+
+ /* we need all the times on the individual levels */
+ i = reader->readAttributeInfo ("numberoftimes", datatype, dim);
+ if(i >=0 && datatype == FLEXIO_INT && dim == 1) {
+ char buffer[100];
+ reader->readAttribute (i, &numberoftimes);
+ assert(numberoftimes==refleveltimes.size());
+ for(int lcv=0;lcv<numberoftimes;lcv++) {
+ sprintf(buffer,"refleveltime%d",lcv);
+ i = reader->readAttributeInfo (buffer, datatype, dim);
+ if(i >=0 && datatype == FLEXIO_REAL && dim == 1) {
+ reader->readAttribute (i, &refleveltimes[lcv]);
+ }
+ else {
+ CCTK_WARN(0,"BAD BAD BAD! Can't read refleveltime!!");
+ }
+ }
+ }
+ else
+ {
+ CCTK_WARN (0, "Unable to restore reflevel times!");
+ }
+ }
+
+ CACTUS_MPI_ERROR (MPI_Bcast (&numberoftimes, 1, CARPET_MPI_INT4, 0,MPI_COMM_WORLD));
+ CACTUS_MPI_ERROR (MPI_Bcast (&refleveltimes[0], numberoftimes, CARPET_MPI_REAL, 0, MPI_COMM_WORLD));
+
+
+
BEGIN_REFLEVEL_LOOP(cgh) {
BEGIN_MGLEVEL_LOOP(cgh) {
+
/* make sure we are looking at the first dataset where
all the good stuff ist stored! */
if(myproc==0)
@@ -246,6 +283,11 @@ int CarpetIOFlexIO_Recover (cGH* cgh, const char *basefilename, int called_from)
/* Recover GH extentions */
CCTK_INFO ("Recovering GH extensions");
result += RecoverGHextensions (cgh, reader);
+
+ cout << refleveltimes[reflevel]<<endl;
+ tt->set_time(reflevel,mglevel,(CCTK_REAL) cgh->cctk_iteration/maxreflevelfact);
+ cout << "tt " << tt->time(0,reflevel,mglevel) << endl;
+
} END_MGLEVEL_LOOP;
} END_REFLEVEL_LOOP;
@@ -259,6 +301,8 @@ int CarpetIOFlexIO_Recover (cGH* cgh, const char *basefilename, int called_from)
delete reader;
delete amrreader;
}
+
+
}
@@ -371,7 +415,8 @@ int CarpetIOFlexIO_Recover (cGH* cgh, const char *basefilename, int called_from)
int i, type,dim;
CCTK_REAL realBuffer;
CCTK_INT4 int4Buffer[2];
-
+ CCTK_INT4 intbuffer,numberoftimes;
+
IObase::DataType datatype;
if (CCTK_MyProc (GH) == 0)
@@ -414,6 +459,29 @@ int CarpetIOFlexIO_Recover (cGH* cgh, const char *basefilename, int called_from)
CCTK_WARN (1, "Unable to restore GH->cctk_time, defaulting to 0.0");
realBuffer = 0.0;
}
+
+ /* finally, we need all the times on the individual levels */
+ i = reader->readAttributeInfo ("numberoftimes", datatype, dim);
+ if(i >=0 && datatype == FLEXIO_INT && dim == 1) {
+ char buffer[100];
+ reader->readAttribute (i, &numberoftimes);
+ assert(numberoftimes==refleveltimes.size());
+ for(int lcv=0;lcv<numberoftimes;lcv++) {
+ sprintf(buffer,"refleveltime%d",lcv);
+ i = reader->readAttributeInfo (buffer, datatype, dim);
+ if(i >=0 && datatype == FLEXIO_REAL && dim == 1) {
+ reader->readAttribute (i, &refleveltimes[lcv]);
+ }
+ else {
+ CCTK_WARN(0,"BAD BAD BAD! Can't read refleveltime!!");
+ }
+ }
+ }
+ else
+ {
+ CCTK_WARN (0, "Unable to restore reflevel times!");
+ }
+
}
#ifdef CCTK_MPI
@@ -421,14 +489,22 @@ int CarpetIOFlexIO_Recover (cGH* cgh, const char *basefilename, int called_from)
/* NOTE: We have to use MPI_COMM_WORLD here
because PUGH_COMM_WORLD is not yet set up at parameter recovery time.
We also assume that PUGH_MPI_INT4 is a compile-time defined datatype. */
+ CACTUS_MPI_ERROR (MPI_Bcast (&numberoftimes, 1, CARPET_MPI_INT4, 0,MPI_COMM_WORLD));
+ CACTUS_MPI_ERROR (MPI_Bcast (&refleveltimes[0], numberoftimes, CARPET_MPI_REAL, 0, MPI_COMM_WORLD));
+ CACTUS_MPI_ERROR (MPI_Bcast (int4Buffer, 2, CARPET_MPI_INT4, 0,MPI_COMM_WORLD));
CACTUS_MPI_ERROR (MPI_Bcast (int4Buffer, 2, CARPET_MPI_INT4, 0,MPI_COMM_WORLD));
CACTUS_MPI_ERROR (MPI_Bcast (&realBuffer, 1, CARPET_MPI_REAL,0,MPI_COMM_WORLD));
#endif
- GH->cctk_time = realBuffer;
+ GH->cctk_time = refleveltimes[reflevel];
+
+ cout << "cctk_time " << realBuffer << endl;
+
GH->cctk_iteration = (int) int4Buffer[0];
CCTK_SetMainLoopIndex ((int) int4Buffer[1]);
+ cout << "refleveltimes" << refleveltimes[0] << " " << refleveltimes[1] << " " << refleveltimes[2] << endl;
+
return (0);
}
@@ -607,9 +683,17 @@ int CarpetIOFlexIO_Recover (cGH* cgh, const char *basefilename, int called_from)
DumpGHExtensions(cgh,writer);
+ /* finally, we need all the times on the individual levels */
+ const int numberoftimes=refleveltimes.size();
+ writer->writeAttribute("numberoftimes",FlexIODataType(CCTK_VARIABLE_INT),1,&numberoftimes);
+ for(int i=0;i < numberoftimes;i++) {
+ char buffer[100];
+ sprintf(buffer,"refleveltime%d",i);
+ writer->writeAttribute(buffer,FlexIODataType(CCTK_VARIABLE_REAL),1,&refleveltimes[i]);
+ }
}
-
+ // cout << "refleveltimes: " << refleveltimes[0] << "," << refleveltimes[1] << endl;
// now dump the grid varibles for all reflevels and components, sorted by groups
BEGIN_REFLEVEL_LOOP(cgh) {
@@ -661,13 +745,22 @@ int CarpetIOFlexIO_Recover (cGH* cgh, const char *basefilename, int called_from)
gdata.numtimelevels = CCTK_GroupStorageIncrease (cgh, 1, &group,
&gdata.numtimelevels,NULL);
- /* dump all timelevels except the oldest (for multi-level groups) */
+
+
+
+
+
CCTK_GroupData (group, &gdata);
+
+ /* dump all timelevels except the oldest (for multi-level groups) */
+ /* COMMENTED OUT!! We _need_ all timelevels!!!
if (gdata.numtimelevels > 1)
{
gdata.numtimelevels--;
}
+ */
+
int first_vindex = CCTK_FirstVarIndexI (group);
@@ -761,10 +854,10 @@ int CarpetIOFlexIO_Recover (cGH* cgh, const char *basefilename, int called_from)
{
if (myGH->cp_filename_list[myGH->cp_filename_index])
{
- // if (checkpoint_keep > 0)
- //{
- // remove (myGH->cp_filename_list[myGH->cp_filename_index]);
- //}
+ if (checkpoint_keep > 0)
+ {
+ remove (myGH->cp_filename_list[myGH->cp_filename_index]);
+ }
free (myGH->cp_filename_list[myGH->cp_filename_index]);
}
myGH->cp_filename_list[myGH->cp_filename_index] = strdup (cp_filename);
@@ -781,8 +874,10 @@ int CarpetIOFlexIO_Recover (cGH* cgh, const char *basefilename, int called_from)
int retval = 0;
int myproc = CCTK_MyProc (cgh);
int currdataset,ndatasets;
+
CCTK_VInfo(CCTK_THORNSTRING,"Starting to recover data!!!");
+
if(myproc==0) {
ndatasets = reader->nDatasets();
//CCTK_VInfo (CCTK_THORNSTRING, "ndatasets=%d", ndatasets);
diff --git a/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexio.cc b/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexio.cc
index a4eb98298..336b0b53a 100644
--- a/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexio.cc
+++ b/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexio.cc
@@ -45,7 +45,7 @@
extern "C" {
- static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexio.cc,v 1.17 2004/01/08 19:43:33 cott Exp $";
+ static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexio.cc,v 1.18 2004/01/09 15:43:46 cott Exp $";
CCTK_FILEVERSION(Carpet_CarpetIOFlexIO_ioflexio_cc);
}
@@ -286,8 +286,8 @@ namespace CarpetIOFlexIO {
// Ignore ghost zones if desired
#warning "need to check size of gdyndata.bbox"
for (int d=0; d<dim; ++d) {
- const int max_lower_ghosts = (gdyndata.bbox[2*d ] && !out3D_output_outer_boundary) ? -1 : out3D_max_num_lower_ghosts;
- const int max_upper_ghosts = (gdyndata.bbox[2*d+1] && !out3D_output_outer_boundary) ? -1 : out3D_max_num_upper_ghosts;
+ const int max_lower_ghosts = (gdyndata.bbox[2*d ] && out3D_output_outer_boundary) ? -1 : out3D_max_num_lower_ghosts;
+ const int max_upper_ghosts = (gdyndata.bbox[2*d+1] && out3D_output_outer_boundary) ? -1 : out3D_max_num_upper_ghosts;
const int num_lower_ghosts = max_lower_ghosts == -1 ? gdyndata.nghostzones[d] : min(out3D_max_num_lower_ghosts, gdyndata.nghostzones[d]);
const int num_upper_ghosts = max_upper_ghosts == -1 ? gdyndata.nghostzones[d] : min(out3D_max_num_upper_ghosts, gdyndata.nghostzones[d]);
@@ -539,7 +539,7 @@ namespace CarpetIOFlexIO {
int asize,i;
IObase::DataType datatype;
int group,varindex;
-
+ CCTK_REAL cctk_time;
if(myproc==0) {
// read the name of the variable
@@ -600,6 +600,16 @@ namespace CarpetIOFlexIO {
{
CCTK_WARN (0, "Something is wrong! Can't read multi group level!!!");
}
+
+ i = reader->readAttributeInfo ("cctk_time", datatype, asize);
+ if (i >= 0 && datatype == FLEXIO_REAL && asize > 0)
+ {
+ reader->readAttribute (i, &cctk_time);
+ }
+ else
+ {
+ CCTK_WARN (0, "Something is wrong! Can't read coordinate time!!!");
+ }
diff --git a/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexioUtil.cc b/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexioUtil.cc
index 97ee800cf..17d2c699a 100644
--- a/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexioUtil.cc
+++ b/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexioUtil.cc
@@ -145,12 +145,10 @@ void DumpCommonAttributes (const cGH *cgh, IObase* writer, ioRequest* request)
-#if 0
- // already dumped by amrwriter if we are dealing with a grid function or a grid array:
- if (CCTK_GroupTypeFromVarI (request->vindex) == CCTK_SCALAR)
- writer->writeAttribute("time",FlexIODataType(CCTK_VARIABLE_REAL),1,&cgh->cctk_time);
-#endif
+ // already dumped by amrwriter (obsolete!)
+ writer->writeAttribute("cctk_time",FlexIODataType(CCTK_VARIABLE_REAL),1,&cgh->cctk_time);
+
/* attributes describing the underlying grid
These are only stored for CCTK_GF variables if they are associated