Skip to content

Commit 4eebdae

Browse files
committed
Bypass XML round-trip in WeightWindowsList.export_to_hdf5
The XML serialization raised MemoryError on bound arrays >~200M elements -- lxml's intermediate ASCII allocation fails before the text node can be built. Write HDF5 directly via h5py, mirroring the C++ WeightWindows::to_hdf5 writer. Critical details for C++ compatibility: - Bounds are 2D (ne, n_voxels) on disk (4D would segfault the C++ tensor::Tensor<double> reader). - max_lower_bound_ratio is written unconditionally (default 1.0). - Root attrs filetype and version are required by openmc_weight_windows_import. Mesh writing is handled by a private _write_mesh_group helper in weight_windows.py that dispatches by mesh type, matching the reference implementation. UnstructuredMesh raises NotImplementedError (wwinp cannot produce one). Also replaces the same ASCII serialization in WeightWindows.to_xml_element with np.char.mod('%.17g', ...), giving ~5x peak memory reduction for the XML path.
1 parent aa90915 commit 4eebdae

2 files changed

Lines changed: 133 additions & 20 deletions

File tree

openmc/weight_windows.py

Lines changed: 110 additions & 20 deletions
Original file line numberDiff line numberDiff line change
@@ -10,7 +10,10 @@
1010
import h5py
1111

1212
import openmc
13-
from openmc.mesh import MeshBase, RectilinearMesh, CylindricalMesh, SphericalMesh, UnstructuredMesh
13+
from openmc.mesh import (
14+
MeshBase, RegularMesh, RectilinearMesh, CylindricalMesh, SphericalMesh,
15+
UnstructuredMesh,
16+
)
1417
from openmc.tallies import Tallies
1518
import openmc.checkvalue as cv
1619
from openmc.checkvalue import PathLike
@@ -335,11 +338,16 @@ def to_xml_element(self) -> ET.Element:
335338
subelement = ET.SubElement(element, 'energy_bounds')
336339
subelement.text = ' '.join(str(e) for e in self.energy_bounds)
337340

341+
# '%.17g' round-trips IEEE 754 doubles through float() losslessly.
342+
def _arr_to_text(arr):
343+
return ' '.join(np.char.mod('%.17g',
344+
np.ascontiguousarray(arr.ravel('F'))))
345+
338346
subelement = ET.SubElement(element, 'lower_ww_bounds')
339-
subelement.text = ' '.join(str(b) for b in self.lower_ww_bounds.ravel('F'))
347+
subelement.text = _arr_to_text(self.lower_ww_bounds)
340348

341349
subelement = ET.SubElement(element, 'upper_ww_bounds')
342-
subelement.text = ' '.join(str(b) for b in self.upper_ww_bounds.ravel('F'))
350+
subelement.text = _arr_to_text(self.upper_ww_bounds)
343351

344352
subelement = ET.SubElement(element, 'survival_ratio')
345353
subelement.text = str(self.survival_ratio)
@@ -819,6 +827,51 @@ def hdf5_to_wws(path='weight_windows.h5') -> WeightWindowsList:
819827
return WeightWindowsList.from_hdf5(path)
820828

821829

830+
def _write_mesh_group(meshes_grp: h5py.Group, mesh: MeshBase):
831+
"""Write *mesh* into *meshes_grp* as a ``mesh <id>`` subgroup.
832+
833+
Matches the layout produced by ``Mesh::to_hdf5`` in src/mesh.cpp so the
834+
resulting file is readable by ``openmc_weight_windows_import``.
835+
"""
836+
grp = meshes_grp.create_group(f'mesh {mesh.id}')
837+
grp.attrs['id'] = np.int32(mesh.id)
838+
grp.create_dataset('name', data=np.bytes_(mesh.name or ''))
839+
840+
if isinstance(mesh, RegularMesh):
841+
grp.create_dataset('type', data=np.bytes_('regular'))
842+
grp.create_dataset('dimension',
843+
data=np.asarray(mesh.dimension, dtype=np.int32))
844+
grp.create_dataset('lower_left',
845+
data=np.asarray(mesh.lower_left, dtype=float))
846+
grp.create_dataset('upper_right',
847+
data=np.asarray(mesh.upper_right, dtype=float))
848+
grp.create_dataset('width', data=np.asarray(mesh.width, dtype=float))
849+
elif isinstance(mesh, RectilinearMesh):
850+
grp.create_dataset('type', data=np.bytes_('rectilinear'))
851+
grp.create_dataset('x_grid', data=np.asarray(mesh.x_grid, dtype=float))
852+
grp.create_dataset('y_grid', data=np.asarray(mesh.y_grid, dtype=float))
853+
grp.create_dataset('z_grid', data=np.asarray(mesh.z_grid, dtype=float))
854+
elif isinstance(mesh, CylindricalMesh):
855+
grp.create_dataset('type', data=np.bytes_('cylindrical'))
856+
grp.create_dataset('r_grid', data=np.asarray(mesh.r_grid, dtype=float))
857+
grp.create_dataset('phi_grid', data=np.asarray(mesh.phi_grid, dtype=float))
858+
grp.create_dataset('z_grid', data=np.asarray(mesh.z_grid, dtype=float))
859+
grp.create_dataset('origin', data=np.asarray(mesh.origin, dtype=float))
860+
elif isinstance(mesh, SphericalMesh):
861+
grp.create_dataset('type', data=np.bytes_('spherical'))
862+
grp.create_dataset('r_grid', data=np.asarray(mesh.r_grid, dtype=float))
863+
grp.create_dataset('theta_grid', data=np.asarray(mesh.theta_grid, dtype=float))
864+
grp.create_dataset('phi_grid', data=np.asarray(mesh.phi_grid, dtype=float))
865+
grp.create_dataset('origin', data=np.asarray(mesh.origin, dtype=float))
866+
elif isinstance(mesh, UnstructuredMesh):
867+
raise NotImplementedError(
868+
"Exporting weight windows on an UnstructuredMesh is not yet "
869+
"implemented in Python; use openmc.lib.export_weight_windows()."
870+
)
871+
else:
872+
raise TypeError(f"Unknown mesh type: {type(mesh).__name__}")
873+
874+
822875
class WeightWindowsList(list):
823876
"""A list of WeightWindows objects.
824877
@@ -1061,29 +1114,66 @@ def from_wwinp(cls, path: PathLike) -> Self:
10611114
def export_to_hdf5(self, path: PathLike = 'weight_windows.h5', **init_kwargs):
10621115
"""Write weight windows to an HDF5 file.
10631116
1117+
Writes the file directly via :mod:`h5py`, mirroring the layout
1118+
produced by :func:`openmc.lib.export_weight_windows`. The previous
1119+
XML round-trip raised :class:`MemoryError` on multi-GB bound arrays
1120+
because of the intermediate ASCII allocation inside lxml.
1121+
10641122
Parameters
10651123
----------
10661124
path : PathLike
10671125
Path to the file to write weight windows to
10681126
**init_kwargs
1069-
Keyword arguments passed to :func:`openmc.lib.init`
1070-
1127+
Unused. Retained for backward compatibility (previously forwarded
1128+
to :func:`openmc.lib.init`).
10711129
"""
1072-
import openmc.lib
10731130
cv.check_type('path', path, PathLike)
1074-
1075-
# Create a temporary model with the weight windows
1076-
model = openmc.Model()
1077-
sph = openmc.Sphere(boundary_type='vacuum')
1078-
cell = openmc.Cell(region=-sph)
1079-
model.geometry = openmc.Geometry([cell])
1080-
model.settings.weight_windows = self
1081-
model.settings.particles = 100
1082-
model.settings.batches = 1
1083-
1084-
# Get absolute path before moving to temporary directory
10851131
path = Path(path).resolve()
10861132

1087-
# Load the model with openmc.lib and then export it to an HDF5 file
1088-
with openmc.lib.TemporarySession(model, **init_kwargs):
1089-
openmc.lib.export_weight_windows(path)
1133+
with h5py.File(path, 'w') as f:
1134+
f.attrs['filetype'] = np.bytes_('weight_windows')
1135+
f.attrs['version'] = np.asarray([1, 0], dtype=np.int32)
1136+
1137+
meshes_grp = f.create_group('meshes')
1138+
wws_grp = f.create_group('weight_windows')
1139+
1140+
seen_mesh_ids = set()
1141+
ww_ids = []
1142+
for ww in self:
1143+
if ww.mesh.id not in seen_mesh_ids:
1144+
_write_mesh_group(meshes_grp, ww.mesh)
1145+
seen_mesh_ids.add(ww.mesh.id)
1146+
1147+
g = wws_grp.create_group(f'weight_windows_{ww.id}')
1148+
ww_ids.append(int(ww.id))
1149+
1150+
g.create_dataset('mesh', data=np.int32(ww.mesh.id))
1151+
g.create_dataset('particle_type',
1152+
data=np.bytes_(str(ww.particle_type)))
1153+
g.create_dataset('energy_bounds',
1154+
data=np.asarray(ww.energy_bounds, dtype=float))
1155+
1156+
# 2D (ne, n_voxels) C-contiguous: the C++ reader expects this layout.
1157+
ne = ww.lower_ww_bounds.shape[-1]
1158+
lo = np.ascontiguousarray(ww.lower_ww_bounds.T).reshape(ne, -1)
1159+
up = np.ascontiguousarray(ww.upper_ww_bounds.T).reshape(ne, -1)
1160+
g.create_dataset('lower_ww_bounds', data=lo)
1161+
g.create_dataset('upper_ww_bounds', data=up)
1162+
1163+
g.create_dataset('survival_ratio',
1164+
data=float(ww.survival_ratio))
1165+
# max_lower_bound_ratio is read unconditionally by C++; default 1.0 when unset.
1166+
mlbr = ww.max_lower_bound_ratio
1167+
g.create_dataset(
1168+
'max_lower_bound_ratio',
1169+
data=float(mlbr if mlbr is not None else 1.0),
1170+
)
1171+
g.create_dataset('max_split', data=np.int32(ww.max_split))
1172+
g.create_dataset('weight_cutoff',
1173+
data=float(ww.weight_cutoff))
1174+
1175+
wws_grp.attrs['ids'] = np.asarray(ww_ids, dtype=np.int32)
1176+
wws_grp.attrs['n_weight_windows'] = np.int32(len(ww_ids))
1177+
meshes_grp.attrs['ids'] = np.asarray(sorted(seen_mesh_ids),
1178+
dtype=np.int32)
1179+
meshes_grp.attrs['n_meshes'] = np.int32(len(seen_mesh_ids))

tests/unit_tests/weightwindows/test_ww_list.py

Lines changed: 23 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,3 +1,4 @@
1+
import h5py
12
import openmc
23

34

@@ -21,3 +22,25 @@ def test_ww_roundtrip(request, run_in_tmpdir):
2122
assert ww.max_split == ww_new.max_split
2223
assert ww.weight_cutoff == ww_new.weight_cutoff
2324
assert ww.mesh.id == ww_new.mesh.id
25+
26+
27+
def test_export_hdf5_format(request, run_in_tmpdir):
28+
# C++ openmc_weight_windows_import expects this on-disk layout.
29+
wws = openmc.WeightWindowsList.from_wwinp(request.path.with_name('wwinp_n'))
30+
wws.export_to_hdf5('ww.h5')
31+
with h5py.File('ww.h5') as f:
32+
assert f.attrs['filetype'] == b'weight_windows'
33+
assert list(f.attrs['version']) == [1, 0]
34+
wws_grp = f['weight_windows']
35+
assert int(wws_grp.attrs['n_weight_windows']) == len(wws)
36+
for ww in wws:
37+
g = wws_grp[f'weight_windows_{ww.id}']
38+
# 2D shape is required by the C++ tensor::Tensor<double> reader.
39+
assert g['lower_ww_bounds'].ndim == 2
40+
assert g['lower_ww_bounds'].shape[0] == ww.num_energy_bins
41+
# max_lower_bound_ratio is read unconditionally by C++.
42+
assert 'max_lower_bound_ratio' in g
43+
m_grp = f['meshes']
44+
for name in m_grp:
45+
assert 'id' in m_grp[name].attrs
46+
assert 'type' in m_grp[name]

0 commit comments

Comments
 (0)