-
Notifications
You must be signed in to change notification settings - Fork 647
Expand file tree
/
Copy pathtest_ww_list.py
More file actions
71 lines (61 loc) · 3.08 KB
/
test_ww_list.py
File metadata and controls
71 lines (61 loc) · 3.08 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
import h5py
import numpy as np
import pytest
import openmc
def test_ww_roundtrip(request, run_in_tmpdir):
# Load weight windows from a wwinp file
wwinp_file = request.path.with_name('wwinp_n')
wws = openmc.WeightWindowsList.from_wwinp(wwinp_file)
# Roundtrip them, writing to HDF5 and reading back in
wws.export_to_hdf5('ww.h5')
wws_new = openmc.WeightWindowsList.from_hdf5('ww.h5')
# Check that the new weight windows are the same as the original
assert len(wws) == len(wws_new)
for ww, ww_new in zip(wws, wws_new):
assert ww.particle_type == ww_new.particle_type
assert (ww.lower_ww_bounds == ww_new.lower_ww_bounds).all()
assert (ww.upper_ww_bounds == ww_new.upper_ww_bounds).all()
assert ww.survival_ratio == ww_new.survival_ratio
assert ww.num_energy_bins == ww_new.num_energy_bins
assert ww.max_split == ww_new.max_split
assert ww.weight_cutoff == ww_new.weight_cutoff
assert ww.mesh.id == ww_new.mesh.id
def test_export_hdf5_format(request, run_in_tmpdir):
# C++ openmc_weight_windows_import expects this on-disk layout.
wws = openmc.WeightWindowsList.from_wwinp(request.path.with_name('wwinp_n'))
wws.export_to_hdf5('ww.h5')
with h5py.File('ww.h5') as f:
assert f.attrs['filetype'] == b'weight_windows'
assert list(f.attrs['version']) == [1, 0]
wws_grp = f['weight_windows']
assert int(wws_grp.attrs['n_weight_windows']) == len(wws)
for ww in wws:
g = wws_grp[f'weight_windows_{ww.id}']
# 2D shape is required by the C++ tensor::Tensor<double> reader.
assert g['lower_ww_bounds'].ndim == 2
assert g['lower_ww_bounds'].shape[0] == ww.num_energy_bins
# max_lower_bound_ratio is read unconditionally by C++.
assert 'max_lower_bound_ratio' in g
m_grp = f['meshes']
for name in m_grp:
assert 'id' in m_grp[name].attrs
assert 'type' in m_grp[name]
@pytest.mark.parametrize('library', ('libmesh', 'moab'))
def test_export_hdf5_unstructured_mesh(request, run_in_tmpdir, library):
# UnstructuredMesh can't be serialized from pure Python; export_to_hdf5
# routes it through openmc.lib.export_weight_windows (a live session).
if library == 'libmesh' and not openmc.lib._libmesh_enabled():
pytest.skip('LibMesh not enabled in this build.')
if library == 'moab' and not openmc.lib._dagmc_enabled():
pytest.skip('DAGMC (and MOAB) not enabled in this build.')
mesh = openmc.UnstructuredMesh(
str(request.path.with_name('test_mesh_tets.exo')), library)
ww = openmc.WeightWindows(mesh, np.ones((12_000,)), upper_bound_ratio=5.0)
openmc.WeightWindowsList([ww]).export_to_hdf5('ww.h5')
with h5py.File('ww.h5') as f:
assert f.attrs['filetype'] == b'weight_windows'
assert list(f.attrs['version']) == [1, 0]
assert int(f['weight_windows'].attrs['n_weight_windows']) == 1
m_grp = f['meshes'][f'mesh {mesh.id}']
assert m_grp['type'][()] == b'unstructured'
assert 'filename' in m_grp