Skip to content
Snippets Groups Projects
Commit 496c4612 authored by Neil Vaytet's avatar Neil Vaytet
Browse files

Refs #24078: Re-introduced closing/opening of nexus file

Needed for instrument loading.
Also re-instated the findGroup() function in NexusGeometryParser as the refactored
version was disfunctional (segfaults)
parent 6f7a7f07
No related merge requests found
......@@ -646,12 +646,19 @@ void LoadEventNexus::loadEvents(API::Progress *const prog,
if (!m_ws->getInstrument() || !m_instrument_loaded_correctly) {
// Load the instrument (if not loaded before)
prog->report("Loading instrument");
// Note that closing an re-opening the file is needed here for loading
// instruments directly from the nexus file containing the event data.
// This may not be needed in the future if both LoadEventNexus and
// LoadInstrument are made to use the same Nexus/HDF5 library
m_file->close();
m_instrument_loaded_correctly =
loadInstrument(m_filename, m_ws, m_top_entry_name, this);
if (!m_instrument_loaded_correctly)
throw std::runtime_error(
"Instrument was not initialized correctly! Loading cannot continue.");
// reopen file
safeOpenFile(m_filename);
}
// top level file information
......
......@@ -154,8 +154,7 @@ std::string get1DStringDataset(const std::string &dataset, const Group &group) {
* size 1.
*/
std::vector<Group> openSubGroups(const Group &parentGroup,
const H5std_string &CLASS_TYPE,
const bool firstEntryOnly = false) {
const H5std_string &CLASS_TYPE) {
std::vector<Group> subGroups;
// Iterate over children, and determine if a group
for (hsize_t i = 0; i < parentGroup.getNumObjs(); ++i) {
......@@ -178,8 +177,6 @@ std::vector<Group> openSubGroups(const Group &parentGroup,
// If group of correct type, append to subGroup vector
if (classType == CLASS_TYPE) {
subGroups.push_back(childGroup);
if (firstEntryOnly)
return subGroups;
}
}
}
......@@ -188,11 +185,45 @@ std::vector<Group> openSubGroups(const Group &parentGroup,
return subGroups;
}
/// Find a single group inside parent (returns first match)
// TODO: refactor most of this and openSubGroups() to avoid duplication
Group findGroup(const Group &parentGroup, const H5std_string &CLASS_TYPE) {
Group childGroup;
// Iterate over children, and determine if a group
for (hsize_t i = 0; i < parentGroup.getNumObjs(); ++i) {
if (parentGroup.getObjTypeByIdx(i) == GROUP_TYPE) {
H5std_string childPath = parentGroup.getObjnameByIdx(i);
// Open the sub group
childGroup = parentGroup.openGroup(childPath);
// Iterate through attributes to find NX_class
for (uint32_t attribute_index = 0;
attribute_index < static_cast<uint32_t>(childGroup.getNumAttrs());
++attribute_index) {
// Test attribute at current index for NX_class
Attribute attribute = childGroup.openAttribute(attribute_index);
if (attribute.getName() == NX_CLASS) {
// Get attribute data type
DataType dataType = attribute.getDataType();
// Get the NX_class type
H5std_string classType;
attribute.read(dataType, classType);
// If group of correct type, return the childGroup
if (classType == CLASS_TYPE) {
return childGroup;
}
}
}
}
}
// TODO: would maybe be better to return boost::optional?
return childGroup;
}
// Get the instrument name
std::string instrumentName(const Group &root) {
std::vector<Group> instrumentGroup = openSubGroups(
openSubGroups(root, NX_ENTRY, true)[0], NX_INSTRUMENT, true);
return get1DStringDataset("name", instrumentGroup[0]);
Group entryGroup = findGroup(root, NX_ENTRY);
Group instrumentGroup = findGroup(entryGroup, NX_INSTRUMENT);
return get1DStringDataset("name", instrumentGroup);
}
// Open all detector groups into a vector
......@@ -572,12 +603,11 @@ parseNexusShape(const Group &detectorGroup, bool &searchTubes) {
void parseAndAddSource(const H5File &file, const Group &root,
InstrumentBuilder &builder) {
H5std_string sourcePath = "source";
std::vector<Group> sourceGroup =
openSubGroups(openSubGroups(openSubGroups(root, NX_ENTRY, true)[0],
NX_INSTRUMENT, true)[0],
sourcePath, true);
auto sourceName = get1DStringDataset("name", sourceGroup[0]);
auto sourceTransformations = getTransformations(file, sourceGroup[0]);
Group entryGroup = findGroup(root, NX_ENTRY);
Group instrumentGroup = findGroup(entryGroup, NX_INSTRUMENT);
Group sourceGroup = findGroup(instrumentGroup, sourcePath);
auto sourceName = get1DStringDataset("name", sourceGroup);
auto sourceTransformations = getTransformations(file, sourceGroup);
auto defaultPos = Eigen::Vector3d(0.0, 0.0, 0.0);
builder.addSource(sourceName, sourceTransformations * defaultPos);
}
......@@ -587,9 +617,9 @@ void parseAndAddSample(const H5File &file, const Group &root,
InstrumentBuilder &builder) {
std::string sampleName = "sample";
H5std_string samplePath = sampleName;
std::vector<Group> sampleGroup =
openSubGroups(openSubGroups(root, NX_ENTRY, true)[0], samplePath, true);
auto sampleTransforms = getTransformations(file, sampleGroup[0]);
Group entryGroup = findGroup(root, NX_ENTRY);
Group sampleGroup = findGroup(entryGroup, samplePath);
auto sampleTransforms = getTransformations(file, sampleGroup);
Eigen::Vector3d samplePos = sampleTransforms * Eigen::Vector3d(0.0, 0.0, 0.0);
builder.addSample(sampleName, samplePos);
}
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment