Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Correct normalization of thermal elastic in non standard ENDF-6 files #3234

Open
wants to merge 3 commits into
base: develop
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
115 changes: 62 additions & 53 deletions openmc/data/thermal.py
Original file line number Diff line number Diff line change
Expand Up @@ -804,7 +804,7 @@ def from_ace(cls, ace_or_filename, name=None):
@classmethod
def from_njoy(cls, filename, filename_thermal, temperatures=None,
evaluation=None, evaluation_thermal=None,
use_endf_data=True, **kwargs):
use_endf_data=True, nonstandard_endf=True, **kwargs):
"""Generate thermal scattering data by running NJOY.

Parameters
Expand All @@ -827,6 +827,10 @@ def from_njoy(cls, filename, filename_thermal, temperatures=None,
use_endf_data : bool
If the material has incoherent elastic scattering, the ENDF data
will be used rather than the ACE data.
nonstandard_endf: bool
Divide incoherent elastic cross section by
number of principal atoms. This is not part of the ENDF-6
standard but it is how it is processed by NJOY.
**kwargs
Keyword arguments passed to :func:`openmc.data.njoy.make_ace_thermal`

Expand All @@ -853,7 +857,7 @@ def from_njoy(cls, filename, filename_thermal, temperatures=None,

# Load ENDF data to replace incoherent elastic
if use_endf_data:
data_endf = cls.from_endf(filename_thermal)
data_endf = cls.from_endf(filename_thermal, nonstandard_endf)
if data_endf.elastic is not None:
# Get appropriate temperatures
if temperatures is None:
Expand All @@ -871,14 +875,18 @@ def from_njoy(cls, filename, filename_thermal, temperatures=None,
return data

@classmethod
def from_endf(cls, ev_or_filename):
def from_endf(cls, ev_or_filename, nonstandard_endf):
"""Generate thermal scattering data from an ENDF file

Parameters
----------
ev_or_filename : openmc.data.endf.Evaluation or str
ENDF evaluation to read from. If given as a string, it is assumed to
be the filename for the ENDF file.
nonstandard_endf: bool
Divide incoherent elastic cross section by
number of principal atoms. This is not part of the ENDF-6
standard but it is how it is processed by NJOY.

Returns
-------
Expand All @@ -891,56 +899,6 @@ def from_endf(cls, ev_or_filename):
else:
ev = endf.Evaluation(ev_or_filename)

# Read coherent/incoherent elastic data
elastic = None
if (7, 2) in ev.section:
# Define helper functions to avoid duplication
def get_coherent_elastic(file_obj):
# Get structure factor at first temperature
params, S = endf.get_tab1_record(file_obj)
strT = _temperature_str(params[0])
n_temps = params[2]
bragg_edges = S.x
xs = {strT: CoherentElastic(bragg_edges, S.y)}
distribution = {strT: CoherentElasticAE(xs[strT])}

# Get structure factor for subsequent temperatures
for _ in range(n_temps):
params, S = endf.get_list_record(file_obj)
strT = _temperature_str(params[0])
xs[strT] = CoherentElastic(bragg_edges, S)
distribution[strT] = CoherentElasticAE(xs[strT])
return xs, distribution

def get_incoherent_elastic(file_obj):
params, W = endf.get_tab1_record(file_obj)
bound_xs = params[0]
xs = {}
distribution = {}
for T, debye_waller in zip(W.x, W.y):
strT = _temperature_str(T)
xs[strT] = IncoherentElastic(bound_xs, debye_waller)
distribution[strT] = IncoherentElasticAE(debye_waller)
return xs, distribution

file_obj = StringIO(ev.section[7, 2])
lhtr = endf.get_head_record(file_obj)[2]
if lhtr == 1:
# coherent elastic
xs, distribution = get_coherent_elastic(file_obj)
elif lhtr == 2:
# incoherent elastic
xs, distribution = get_incoherent_elastic(file_obj)
elif lhtr == 3:
# mixed coherent / incoherent elastic
xs_c, dist_c = get_coherent_elastic(file_obj)
xs_i, dist_i = get_incoherent_elastic(file_obj)
assert sorted(xs_c) == sorted(xs_i)
xs = {T: Sum([xs_c[T], xs_i[T]]) for T in xs_c}
distribution = {T: MixedElasticAE(dist_c[T], dist_i[T]) for T in dist_c}

elastic = ThermalScatteringReaction(xs, distribution)

# Read incoherent inelastic data
assert (7, 4) in ev.section, 'No MF=7, MT=4 found in thermal scattering'
file_obj = StringIO(ev.section[7, 4])
Expand Down Expand Up @@ -999,6 +957,57 @@ def get_incoherent_elastic(file_obj):
_, Teff = endf.get_tab1_record(file_obj)
data['effective_temperature'].append(Teff)

# Read coherent/incoherent elastic data
elastic = None
if (7, 2) in ev.section:
# Define helper functions to avoid duplication
def get_coherent_elastic(file_obj):
# Get structure factor at first temperature
params, S = endf.get_tab1_record(file_obj)
strT = _temperature_str(params[0])
n_temps = params[2]
bragg_edges = S.x
xs = {strT: CoherentElastic(bragg_edges, S.y)}
distribution = {strT: CoherentElasticAE(xs[strT])}

# Get structure factor for subsequent temperatures
for _ in range(n_temps):
params, S = endf.get_list_record(file_obj)
strT = _temperature_str(params[0])
xs[strT] = CoherentElastic(bragg_edges, S)
distribution[strT] = CoherentElasticAE(xs[strT])
return xs, distribution

def get_incoherent_elastic(file_obj, natom):
params, W = endf.get_tab1_record(file_obj)
bound_xs = params[0]/natom
xs = {}
distribution = {}
for T, debye_waller in zip(W.x, W.y):
strT = _temperature_str(T)
xs[strT] = IncoherentElastic(bound_xs, debye_waller)
distribution[strT] = IncoherentElasticAE(debye_waller)
return xs, distribution

file_obj = StringIO(ev.section[7, 2])
lhtr = endf.get_head_record(file_obj)[2]
natom = data['M0'] if nonstandard_endf else 1
if lhtr == 1:
# coherent elastic
xs, distribution = get_coherent_elastic(file_obj)
elif lhtr == 2:
# incoherent elastic
xs, distribution = get_incoherent_elastic(file_obj, natom)
elif lhtr == 3:
# mixed coherent / incoherent elastic
xs_c, dist_c = get_coherent_elastic(file_obj)
xs_i, dist_i = get_incoherent_elastic(file_obj, natom)
assert sorted(xs_c) == sorted(xs_i)
xs = {T: Sum([xs_c[T], xs_i[T]]) for T in xs_c}
distribution = {T: MixedElasticAE(dist_c[T], dist_i[T]) for T in dist_c}

elastic = ThermalScatteringReaction(xs, distribution)

name = ev.target['zsymam'].strip()
instance = cls(name, awr, energy_max, kTs)
if elastic is not None:
Expand Down
6 changes: 3 additions & 3 deletions tests/unit_tests/test_data_thermal.py
Original file line number Diff line number Diff line change
Expand Up @@ -41,7 +41,7 @@ def hzrh():
"""H in ZrH thermal scattering data."""
endf_data = os.environ['OPENMC_ENDF_DATA']
filename = os.path.join(endf_data, 'thermal_scatt', 'tsl-HinZrH.endf')
return openmc.data.ThermalScattering.from_endf(filename)
return openmc.data.ThermalScattering.from_endf(filename, nonstandard_endf=True)


@pytest.fixture(scope='module')
Expand All @@ -64,7 +64,7 @@ def sio2():
"""SiO2 thermal scattering data."""
endf_data = os.environ['OPENMC_ENDF_DATA']
filename = os.path.join(endf_data, 'thermal_scatt', 'tsl-SiO2.endf')
return openmc.data.ThermalScattering.from_endf(filename)
return openmc.data.ThermalScattering.from_endf(filename, nonstandard_endf=True)


def test_h2o_attributes(h2o):
Expand Down Expand Up @@ -144,7 +144,7 @@ def test_continuous_dist(h2o_njoy):
def test_h2o_endf():
endf_data = os.environ['OPENMC_ENDF_DATA']
filename = os.path.join(endf_data, 'thermal_scatt', 'tsl-HinH2O.endf')
h2o = openmc.data.ThermalScattering.from_endf(filename)
h2o = openmc.data.ThermalScattering.from_endf(filename, nonstandard_endf=True)
assert not h2o.elastic
assert h2o.atomic_weight_ratio == pytest.approx(0.99917)
assert h2o.energy_max == pytest.approx(3.99993)
Expand Down
Loading