aboutsummaryrefslogtreecommitdiff
path: root/Carpet
diff options
context:
space:
mode:
authorRoland Haas <roland.haas@physics.gatech.edu>2010-10-02 18:58:55 -0400
committerBarry Wardell <barry.wardell@gmail.com>2011-12-14 18:25:33 +0000
commit116af5710a4ff435bb7365d21778c82d2b756047 (patch)
tree4e6889d7f87ea190cbad3322b75c381eced37e9c /Carpet
parentd03220ad193e13be7e3fea74be259114a1ef067f (diff)
CarpetIOHDF5: make sure metadata matches old code
* add some attributes to metadata to match what the old code wrote * only tag object names with components of there are more than one
Diffstat (limited to 'Carpet')
-rw-r--r--Carpet/CarpetIOHDF5/src/CarpetIOHDF5.cc21
-rw-r--r--Carpet/CarpetIOHDF5/src/CarpetIOHDF5.hh5
-rw-r--r--Carpet/CarpetIOHDF5/src/OutputSlice.cc35
3 files changed, 49 insertions, 12 deletions
diff --git a/Carpet/CarpetIOHDF5/src/CarpetIOHDF5.cc b/Carpet/CarpetIOHDF5/src/CarpetIOHDF5.cc
index 669d57e2d..13410d47c 100644
--- a/Carpet/CarpetIOHDF5/src/CarpetIOHDF5.cc
+++ b/Carpet/CarpetIOHDF5/src/CarpetIOHDF5.cc
@@ -293,9 +293,14 @@ hid_t CCTKtoHDF5_Datatype (const cGH* const cctkGH,
int AddSliceAttributes(const cGH* const cctkGH,
const char* const fullname,
const int refinementlevel,
+ const int multigridlevel,
+ const int map,
+ const int timelevel,
const vector<double>& origin,
const vector<double>& delta,
const vector<int>& iorigin,
+ const vector<int>& bbox,
+ const vector<int>& nghostzones,
hid_t& dataset)
{
int error_count = 0;
@@ -304,9 +309,23 @@ int AddSliceAttributes(const cGH* const cctkGH,
error_count += WriteAttribute(dataset, "timestep", cctkGH->cctk_iteration);
error_count += WriteAttribute(dataset, "name", fullname);
error_count += WriteAttribute(dataset, "level", refinementlevel);
+ error_count += WriteAttribute(dataset, "carpet_mglevel", multigridlevel);
+ error_count += WriteAttribute(dataset, "group_timelevel", timelevel);
error_count += WriteAttribute(dataset, "origin", &origin[0], origin.size());
error_count += WriteAttribute(dataset, "delta", &delta[0], delta.size());
error_count += WriteAttribute(dataset, "iorigin", &iorigin[0], iorigin.size());
+ // bbox and nghostzones are only used for grid functions and grid arrays
+ if (bbox.size() > 0) {
+ error_count += WriteAttribute(dataset, "cctk_bbox", &bbox[0], bbox.size());
+ }
+ if (nghostzones.size() > 0) {
+ error_count += WriteAttribute(dataset, "cctk_nghostzones", &nghostzones[0], nghostzones.size());
+ }
+ // Specify whether the coordinate system is Cartesian or not
+ if (CCTK_IsFunctionAliased ("MultiPatch_MapIsCartesian")) {
+ int const map_is_cartesian = MultiPatch_MapIsCartesian (map);
+ error_count += WriteAttribute(dataset, "MapIsCartesian", map_is_cartesian);
+ }
return error_count;
}
@@ -1321,7 +1340,7 @@ static int WriteAttribute (hid_t const group,
int error_count = 0;
HDF5_ERROR (datatype = H5Tcopy (H5T_C_S1));
- HDF5_ERROR (H5Tset_size (datatype, strlen (svalue) + 1));
+ HDF5_ERROR (H5Tset_size (datatype, strlen (svalue)));
HDF5_ERROR (dataspace = H5Screate (H5S_SCALAR));
HDF5_ERROR (attribute = H5Acreate (group, name, datatype,
dataspace, H5P_DEFAULT));
diff --git a/Carpet/CarpetIOHDF5/src/CarpetIOHDF5.hh b/Carpet/CarpetIOHDF5/src/CarpetIOHDF5.hh
index f42c8f23b..f22605430 100644
--- a/Carpet/CarpetIOHDF5/src/CarpetIOHDF5.hh
+++ b/Carpet/CarpetIOHDF5/src/CarpetIOHDF5.hh
@@ -134,9 +134,14 @@ namespace CarpetIOHDF5
int AddSliceAttributes(const cGH* const cctkGH,
const char* const fullname,
const int refinementlevel,
+ const int multigridlevel,
+ const int map,
+ const int timelevel,
const vector<double>& origin,
const vector<double>& delta,
const vector<int>& iorigin,
+ const vector<int>& bbox,
+ const vector<int>& nghostzones,
hid_t& dataset);
// returns an HDF5 datatype corresponding to the given CCTK datatype
diff --git a/Carpet/CarpetIOHDF5/src/OutputSlice.cc b/Carpet/CarpetIOHDF5/src/OutputSlice.cc
index 5225e8aee..c503b95f0 100644
--- a/Carpet/CarpetIOHDF5/src/OutputSlice.cc
+++ b/Carpet/CarpetIOHDF5/src/OutputSlice.cc
@@ -539,7 +539,8 @@ namespace CarpetIOHDF5 {
const int num_tl = CCTK_NumTimeLevelsFromVarI (vindex);
assert (num_tl >= 1);
- const int numvars = CCTK_NumVarsInGroupI(group);
+ const int numvars = one_file_per_group ?
+ CCTK_NumVarsInGroupI(group) : 1;
@@ -1237,12 +1238,10 @@ namespace CarpetIOHDF5 {
assert (outdim<=dim);
cGroup groupdata;
- {
- int const gi = CCTK_GroupIndexFromVarI (vi);
- assert (gi >= 0);
- int const ierr = CCTK_GroupData (gi, & groupdata);
- assert (not ierr);
- }
+ int const gi = CCTK_GroupIndexFromVarI (vi);
+ assert (gi >= 0);
+ int const ierr = CCTK_GroupData (gi, & groupdata);
+ assert (not ierr);
// boolean that says if we are doing 1D-diagonal output
// This is not beautiful, but works for the moment
@@ -1296,8 +1295,9 @@ namespace CarpetIOHDF5 {
if (maps > 1) datasetname_suffix << " m=" << m;
datasetname_suffix << " rl=" << rl;
}
- if (groupdata.grouptype == CCTK_GF or
- groupdata.disttype != CCTK_DISTRIB_CONSTANT) {
+ if (arrdata.at(gi).at(m).dd->
+ light_boxes.at(ml).at(rl).size () > 1 and
+ groupdata.disttype != CCTK_DISTRIB_CONSTANT) {
datasetname_suffix << " c=" << output_component;
}
@@ -1371,6 +1371,7 @@ namespace CarpetIOHDF5 {
vector<int> iorigin(rank, 0);
vector<double> delta(rank, 0), origin(rank, 0);
+ vector<int> bbox(2*rank, 0), nghostzones(rank, 0);
for (int d = 0; d < outdim; d++) {
assert(gfext.upper()[dirs[d]] - gfext.lower()[dirs[d]] >= 0);
iorigin[d] = ext.lower()[d];
@@ -1383,6 +1384,18 @@ namespace CarpetIOHDF5 {
iorigin[d] /= gfext.stride()[dirs[d]];
}
}
+ // store cctk_bbox and cctk_nghostzones (for grid arrays only)
+ if (groupdata.grouptype != CCTK_SCALAR) {
+ const b2vect obnds = vhh.at(m)->outer_boundaries(rl,c);
+ const i2vect ghost_width = arrdata.at(gi).at(m).dd->ghost_widths.AT(rl);
+ for (int d = 0; d < outdim; d++) {
+ nghostzones[d] = ghost_width[0][dirs[d]];
+ assert (all (ghost_width[0] == ghost_width[1]));
+
+ bbox[2*d] = obnds[0][dirs[d]];
+ bbox[2*d+1] = obnds[1][dirs[d]];
+ }
+ }
// now loop over all variables
for (size_t n = 0; n < gfdatas.size(); n++) {
@@ -1407,8 +1420,8 @@ namespace CarpetIOHDF5 {
HDF5_ERROR(H5Dwrite (dataset, mem_type, mem_space, H5S_ALL,
H5P_DEFAULT, gfdatas[n]->storage()));
error_count +=
- AddSliceAttributes (cctkGH, fullname, rl, origin, delta, iorigin,
- dataset);
+ AddSliceAttributes (cctkGH, fullname, rl, ml, m, tl, origin, delta,
+ iorigin, bbox, nghostzones, dataset);
HDF5_ERROR(H5Dclose (dataset));
free (fullname);