Skip to content

Commit

Permalink
For readability, move bulk of write code into IOFamily (#244)
Browse files Browse the repository at this point in the history
Now the methods IOFamily::writePartitioned and IOFamily::writeSerial
handle the data write, making the logic in IODataOrganizer::write
easier to follow.
  • Loading branch information
trevilo committed Feb 2, 2024
1 parent 38afd0b commit 6a19a93
Show file tree
Hide file tree
Showing 2 changed files with 60 additions and 39 deletions.
94 changes: 56 additions & 38 deletions src/io.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -585,6 +585,58 @@ void IOFamily::serializeForWrite() {
}
}

void IOFamily::writePartitioned(hid_t file) {
hsize_t dims[1];
hid_t group = -1;
hid_t dataspace = -1;

// use all ranks to write data from this->.pfunc_
dims[0] = pfunc_->ParFESpace()->GetNDofs();
dataspace = H5Screate_simple(1, dims, NULL);
assert(dataspace >= 0);
group = H5Gcreate(file, group_.c_str(), H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
assert(group >= 0);

// get pointer to raw data
const double *data = pfunc_->HostRead();

// save raw data
for (auto var : vars_) write_soln_data(group, var.varName_, dataspace, data + var.index_ * dims[0], rank0_);

if (group >= 0) H5Gclose(group);
if (dataspace >= 0) H5Sclose(dataspace);
}

void IOFamily::writeSerial(hid_t file) {
// Only get here if need to serialize
serializeForWrite();

// Only rank 0 writes data (which has just been collected into this->serial_sol)
if (rank0_) {
hsize_t dims[1];
hid_t group = -1;
hid_t dataspace = -1;

assert(serial_fes != NULL);
dims[0] = serial_fes->GetNDofs();

dataspace = H5Screate_simple(1, dims, NULL);
assert(dataspace >= 0);
group = H5Gcreate(file, group_.c_str(), H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
assert(group >= 0);

// get pointer to raw data
assert(serial_sol != NULL);
const double *data = serial_sol->HostRead();

// save raw data
for (auto var : vars_) write_soln_data(group, var.varName_, dataspace, data + var.index_ * dims[0], rank0_);

if (group >= 0) H5Gclose(group);
if (dataspace >= 0) H5Sclose(dataspace);
}
}

void IOFamily::readPartitioned(hid_t file) {
hid_t dataspace;
hsize_t numInSoln = 0;
Expand Down Expand Up @@ -801,55 +853,20 @@ void IODataOrganizer::write(hid_t file, bool serial) {

// Loop over defined IO families to save desired output
for (auto fam : families_) {
if (require_serialization) fam.serializeForWrite();

// determine if rank 0
const int rank = fam.pfunc_->ParFESpace()->GetMyRank();
const bool rank0 = (rank == 0);

hsize_t dims[1];
hid_t group = -1;
hid_t dataspace = -1;

if (rank0) {
grvy_printf(ginfo, "\nCreating HDF5 group for defined IO families\n");
grvy_printf(ginfo, "--> %s : %s\n", fam.group_.c_str(), fam.description_.c_str());
}

// Write handled by appropriate method from IOFamily
if (require_serialization) {
// if require_serialization, then we just populated fam.serial_sol above, and now we use only rank0 to write
if (rank0) {
assert(fam.serial_fes != NULL);
dims[0] = fam.serial_fes->GetNDofs();

dataspace = H5Screate_simple(1, dims, NULL);
assert(dataspace >= 0);
group = H5Gcreate(file, fam.group_.c_str(), H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
assert(group >= 0);

// get pointer to raw data
assert(fam.serial_sol != NULL);
const double *data = fam.serial_sol->HostRead();

// save raw data
for (auto var : fam.vars_) write_soln_data(group, var.varName_, dataspace, data + var.index_ * dims[0], rank0);
}
fam.writeSerial(file);
} else {
// otherwise, use all ranks to write data from fam.pfunc_
dims[0] = fam.pfunc_->ParFESpace()->GetNDofs();
dataspace = H5Screate_simple(1, dims, NULL);
assert(dataspace >= 0);
group = H5Gcreate(file, fam.group_.c_str(), H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
assert(group >= 0);

// get pointer to raw data
const double *data = fam.pfunc_->HostRead();

// save raw data
for (auto var : fam.vars_) write_soln_data(group, var.varName_, dataspace, data + var.index_ * dims[0], rank0);
fam.writePartitioned(file);
}
if (group >= 0) H5Gclose(group);
if (dataspace >= 0) H5Sclose(dataspace);
}
}

Expand Down Expand Up @@ -885,6 +902,7 @@ void IODataOrganizer::read(hid_t file, bool serial, int read_order) {
change_order = (fam_order != read_order);
}

// Read handled by appropriate method from IOFamily
if (change_order) {
if (serial) {
if (rank0) grvy_printf(gerror, "[ERROR]: Serial read does not support order change.\n");
Expand Down
5 changes: 4 additions & 1 deletion src/io.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -52,6 +52,8 @@ class IOFamily {
protected:
bool rank0_;

void serializeForWrite();

public:
std::string description_; // family description
std::string group_; // HDF5 group name
Expand Down Expand Up @@ -80,7 +82,8 @@ class IOFamily {

IOFamily(std::string desc, std::string grp, mfem::ParGridFunction *pf);

void serializeForWrite();
void writePartitioned(hid_t file);
void writeSerial(hid_t file);
void readPartitioned(hid_t file);
void readSerial(hid_t file);
void readChangeOrder(hid_t file, int read_order);
Expand Down

0 comments on commit 6a19a93

Please sign in to comment.