Skip to content

Commit dea887e

Browse files
committed
Bug fix in plotting utils. Now works for multifabs with multiple boxes.
1 parent b7c660a commit dea887e

File tree

2 files changed

+23
-6
lines changed

2 files changed

+23
-6
lines changed

swm_AMReX/plotting_utils/main.cpp

+23-6
Original file line numberDiff line numberDiff line change
@@ -42,6 +42,16 @@ void WriteMultiFabToHDF5(const std::string& plotfile, const std::string& hdf5fil
4242
H5Sclose(attr_space_id);
4343
}
4444

45+
// Add an attribute to specify the data layout of the following datasets (row-major or column-major)
46+
{
47+
hid_t attr_space_id = H5Screate(H5S_SCALAR);
48+
hid_t attr_id = H5Acreate(file_id, "data_layout", H5T_C_S1, attr_space_id, H5P_DEFAULT, H5P_DEFAULT);
49+
const char* layout = "row-major";
50+
H5Awrite(attr_id, H5T_C_S1, layout);
51+
H5Aclose(attr_id);
52+
H5Sclose(attr_space_id);
53+
}
54+
4555
// Get the dimensions of the MultiFab
4656
hsize_t dims[2] = {static_cast<hsize_t>(minimal_box.length(0)), static_cast<hsize_t>(minimal_box.length(1))};
4757

@@ -51,20 +61,27 @@ void WriteMultiFabToHDF5(const std::string& plotfile, const std::string& hdf5fil
5161
hid_t dataspace_id = H5Screate_simple(2, dims, NULL);
5262
hid_t dataset_id = H5Dcreate(file_id, varnames[component_idx].c_str(), H5T_NATIVE_DOUBLE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
5363

54-
// Write the data to the dataset
64+
std::vector<double> component_data(dims[0]*dims[1]);
65+
66+
// Loop over all sub-boxes in the MultiFab and fill component_data
5567
for (amrex::MFIter mfi(mf); mfi.isValid(); ++mfi) {
5668
const amrex::Box& bx = mfi.validbox();
69+
const amrex::Dim3 lo = amrex::lbound(bx);
70+
const amrex::Dim3 hi = amrex::ubound(bx);
5771
const amrex::Array4<const amrex::Real> &fab = mf.array(mfi);
5872
std::vector<double> data(bx.numPts());
59-
int idx = 0;
60-
for (int j = bx.smallEnd(1); j <= bx.bigEnd(1); ++j) {
61-
for (int i = bx.smallEnd(0); i <= bx.bigEnd(0); ++i) {
62-
data[idx++] = fab(i, j, 0, component_idx);
73+
74+
// Looping over the box in column-major order since that is how the multi-fab is stored
75+
for (int j = lo.y; j <= hi.y; ++j) {
76+
for (int i = lo.x; i <= hi.x; ++i) {
77+
const int idx = i*dims[1] + j; // Writing data to hdf5 file in row-major order
78+
component_data[idx] = fab(i, j, 0, component_idx);
6379
}
6480
}
65-
H5Dwrite(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, data.data());
6681
}
6782

83+
H5Dwrite(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, component_data.data());
84+
6885
H5Dclose(dataset_id);
6986
H5Sclose(dataspace_id);
7087
}
File renamed without changes.

0 commit comments

Comments
 (0)