openmc.data.ThermalScattering¶
-
class
openmc.data.
ThermalScattering
(name, atomic_weight_ratio, energy_max, kTs)[source]¶ A ThermalScattering object contains thermal scattering data as represented by an S(alpha, beta) table.
Parameters: Variables: - atomic_weight_ratio (float) – Atomic mass ratio of the target nuclide.
- energy_max (float) – Maximum energy for thermal scattering data in [eV]
- elastic (openmc.data.ThermalScatteringReaction or None) – Elastic scattering derived in the coherent or incoherent approximation
- inelastic (openmc.data.ThermalScatteringReaction) – Inelastic scattering cross section derived in the incoherent approximation
- name (str) – Name of the material using GND convention, e.g. c_H_in_H2O
- temperatures (Iterable of str) – List of string representations the temperatures of the target nuclide in the data set. The temperatures are strings of the temperature, rounded to the nearest integer; e.g., ‘294K’
- kTs (Iterable of float) – List of temperatures of the target nuclide in the data set. The temperatures have units of eV.
- nuclides (Iterable of str) – Nuclide names that the thermal scattering data applies to
-
add_temperature_from_ace
(ace_or_filename, name=None)[source]¶ Add data to the ThermalScattering object from an ACE file at a different temperature.
Parameters: - ace_or_filename (openmc.data.ace.Table or str) – ACE table to read from. If given as a string, it is assumed to be the filename for the ACE file.
- name (str) – GND-conforming name of the material, e.g. c_H_in_H2O. If none is passed, the appropriate name is guessed based on the name of the ACE table.
Returns: Thermal scattering data
Return type:
-
export_to_hdf5
(path, mode='a', libver='earliest')[source]¶ Export table to an HDF5 file.
Parameters: - path (str) – Path to write HDF5 file to
- mode ({'r', r+', 'w', 'x', 'a'}) – Mode that is used to open the HDF5 file. This is the second argument
to the
h5py.File
constructor. - libver ({'earliest', 'latest'}) – Compatibility mode for the HDF5 file. ‘latest’ will produce files that are less backwards compatible but have performance benefits.
-
classmethod
from_ace
(ace_or_filename, name=None)[source]¶ Generate thermal scattering data from an ACE table
Parameters: - ace_or_filename (openmc.data.ace.Table or str) – ACE table to read from. If given as a string, it is assumed to be the filename for the ACE file.
- name (str) – GND-conforming name of the material, e.g. c_H_in_H2O. If none is passed, the appropriate name is guessed based on the name of the ACE table.
Returns: Thermal scattering data
Return type:
-
classmethod
from_endf
(ev_or_filename)[source]¶ 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. Returns: Thermal scattering data Return type: openmc.data.ThermalScattering
-
classmethod
from_hdf5
(group_or_filename)[source]¶ Generate thermal scattering data from HDF5 group
Parameters: group_or_filename (h5py.Group or str) – HDF5 group containing interaction data. If given as a string, it is assumed to be the filename for the HDF5 file, and the first group is used to read from. Returns: Neutron thermal scattering data Return type: openmc.data.ThermalScattering
-
classmethod
from_njoy
(filename, filename_thermal, temperatures=None, evaluation=None, evaluation_thermal=None, use_endf_data=True, **kwargs)[source]¶ Generate thermal scattering data by running NJOY.
Parameters: - filename (str) – Path to ENDF neutron sublibrary file
- filename_thermal (str) – Path to ENDF thermal scattering sublibrary file
- temperatures (iterable of float) – Temperatures in Kelvin to produce data at. If omitted, data is produced at all temperatures in the ENDF thermal scattering sublibrary.
- evaluation (openmc.data.endf.Evaluation, optional) – If the ENDF neutron sublibrary file contains multiple material evaluations, this argument indicates which evaluation to use.
- evaluation_thermal (openmc.data.endf.Evaluation, optional) – If the ENDF thermal scattering sublibrary file contains multiple material evaluations, this argument indicates which evaluation to use.
- use_endf_data (bool) – If the material has incoherent elastic scattering, the ENDF data will be used rather than the ACE data.
- **kwargs – Keyword arguments passed to
openmc.data.njoy.make_ace_thermal()
Returns: data – Thermal scattering data
Return type: