aboutsummaryrefslogtreecommitdiff
path: root/CarpetAttic
diff options
context:
space:
mode:
authorschnetter <>2003-12-10 17:08:00 +0000
committerschnetter <>2003-12-10 17:08:00 +0000
commit30d2aeeeca6bb0c92def1bd16f5e8b7be0f04e3a (patch)
tree28e561794d6931ff060279bb6194fd3be4a55154 /CarpetAttic
parenta665cf5aa6f067c2f1b475d9f9aa537bc367ac6d (diff)
Add the content of the cGH structure to each dataset. This allows to
Add the content of the cGH structure to each dataset. This allows to extract all the information independent of FlexIO. There is also a version number "carpet_flexio_version" that will increase for every incompatible change. darcs-hash:20031210170853-07bb3-31f1ce9c51d6dfb67d265f1ee806522c7391aac5.gz
Diffstat (limited to 'CarpetAttic')
-rw-r--r--CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc66
1 files changed, 66 insertions, 0 deletions
diff --git a/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc b/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc
index 2a00ebe5c..ccd0af66e 100644
--- a/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc
+++ b/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc
@@ -72,6 +72,15 @@ namespace CarpetIOFlexIO {
const char* const varlist, const int vindex);
static void SetFlag (int index, const char* optstring, void* arg);
+ static void WriteAttribute (IObase* writer, const char* name,
+ int value);
+ static void WriteAttribute (IObase* writer, const char* name,
+ const int* values, int nvalues);
+ static void WriteAttribute (IObase* writer, const char* name,
+ CCTK_REAL value);
+ static void WriteAttribute (IObase* writer, const char* name,
+ const CCTK_REAL* values, int nvalues);
+
int CarpetIOFlexIOStartup ()
@@ -317,6 +326,25 @@ namespace CarpetIOFlexIO {
dims[d] = (ext.shape() / ext.stride())[d];
}
amrwriter->write (origin, dims, (void*)tmp->storage());
+
+ // Write some additional attributes
+ WriteAttribute (writer, "carpet_flexio_version", 1);
+ WriteAttribute (writer, "cctk_dim", cgh->cctk_dim);
+ WriteAttribute (writer, "cctk_iteration", cgh->cctk_iteration);
+ WriteAttribute (writer, "cctk_gsh", cgh->cctk_gsh, dim);
+ WriteAttribute (writer, "cctk_lsh", cgh->cctk_lsh, dim);
+ WriteAttribute (writer, "cctk_lbnd", cgh->cctk_lbnd, dim);
+ WriteAttribute (writer, "cctk_delta_time", cgh->cctk_delta_time);
+ WriteAttribute (writer, "cctk_delta_space", cgh->cctk_delta_space, dim);
+ WriteAttribute (writer, "cctk_origin_space", cgh->cctk_origin_space, dim);
+ WriteAttribute (writer, "cctk_bbox", cgh->cctk_bbox, 2*dim);
+ WriteAttribute (writer, "cctk_levfac", cgh->cctk_levfac, dim);
+ WriteAttribute (writer, "cctk_levoff", cgh->cctk_levoff, dim);
+ WriteAttribute (writer, "cctk_levoffdenom", cgh->cctk_levoffdenom, dim);
+ WriteAttribute (writer, "cctk_timefac", cgh->cctk_timefac);
+ WriteAttribute (writer, "cctk_convlevel", cgh->cctk_convlevel);
+ WriteAttribute (writer, "cctk_nghostzones", cgh->cctk_nghostzones, dim);
+ WriteAttribute (writer, "cctk_time", cgh->cctk_time);
}
// Delete temporary copy
@@ -698,4 +726,42 @@ namespace CarpetIOFlexIO {
+ void WriteAttribute (IObase* writer, const char* name, int value)
+ {
+ WriteAttribute (writer, name, &value, 1);
+ }
+
+ void WriteAttribute (IObase* writer, const char* name,
+ const int* values, int nvalues)
+ {
+ assert (writer);
+ assert (name);
+ assert (values);
+ vector<CCTK_INT4> values1(nvalues);
+ for (int i=0; i<nvalues; ++i) {
+ values1[i] = values[i];
+ }
+ writer->writeAttribute (name, IObase::Int32, nvalues, &values1[0]);
+ }
+
+ void WriteAttribute (IObase* writer, const char* name, CCTK_REAL value)
+ {
+ WriteAttribute (writer, name, &value, 1);
+ }
+
+ void WriteAttribute (IObase* writer, const char* name,
+ const CCTK_REAL* values, int nvalues)
+ {
+ assert (writer);
+ assert (name);
+ assert (values);
+ vector<CCTK_REAL8> values1(nvalues);
+ for (int i=0; i<nvalues; ++i) {
+ values1[i] = values[i];
+ }
+ writer->writeAttribute (name, IObase::Float64, nvalues, &values1[0]);
+ }
+
+
+
} // namespace CarpetIOFlexIO