openmc.data.WindowedMultipole¶
- class openmc.data.WindowedMultipole(name)[source]¶
Resonant cross sections represented in the windowed multipole format.
- Parameters
name (str) – Name of the nuclide using the GNDS naming convention
- Variables
name (str) – Name of the nuclide using the GNDS naming convention
spacing (float) – The width of each window in sqrt(E)-space. For example, the frst window will end at (sqrt(E_min) + spacing)**2 and the second window at (sqrt(E_min) + 2*spacing)**2.
sqrtAWR (float) – Square root of the atomic weight ratio of the target nuclide.
E_min (float) – Lowest energy in eV the library is valid for.
E_max (float) – Highest energy in eV the library is valid for.
data (np.ndarray) – A 2D array of complex poles and residues. data[i, 0] gives the energy at which pole i is located. data[i, 1:] gives the residues associated with the i-th pole. There are 3 residues, one each for the scattering, absorption, and fission channels.
windows (np.ndarray) – A 2D array of Integral values. windows[i, 0] - 1 is the index of the first pole in window i. windows[i, 1] - 1 is the index of the last pole in window i.
broaden_poly (np.ndarray) – A 1D array of boolean values indicating whether or not the polynomial curvefit in that window should be Doppler broadened.
curvefit (np.ndarray) – A 3D array of Real curvefit polynomial coefficients. curvefit[i, 0, :] gives coefficients for the scattering cross section in window i. curvefit[i, 1, :] gives absorption coefficients and curvefit[i, 2, :] gives fission coefficients. The polynomial terms are increasing powers of sqrt(E) starting with 1/E e.g: a/E + b/sqrt(E) + c + d sqrt(E) + …
- export_to_hdf5(path, mode='a', libver='earliest')[source]¶
Export windowed multipole data to an HDF5 file.
- Parameters
path (str) – Path to write HDF5 file to
mode ({'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_endf(endf_file, log=False, vf_options=None, wmp_options=None)[source]¶
Generate windowed multipole neutron data from an ENDF evaluation.
New in version 0.12.1.
- Parameters
endf_file (str) – Path to ENDF evaluation
log (bool or int, optional) – Whether to print running logs (use int for verbosity control)
vf_options (dict, optional) – Dictionary of keyword arguments, e.g. {‘njoy_error’: 0.001}, passed to
openmc.data.multipole.vectfit_nuclide()
wmp_options (dict, optional) – Dictionary of keyword arguments, e.g. {‘search’: True, ‘rtol’: 0.01}, passed to
openmc.data.WindowedMultipole.from_multipole()
- Returns
Resonant cross sections represented in the windowed multipole format.
- Return type
- classmethod from_hdf5(group_or_filename)[source]¶
Construct a WindowedMultipole object from an HDF5 group or file.
- Parameters
group_or_filename (h5py.Group or str) – HDF5 group containing multipole 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
Resonant cross sections represented in the windowed multipole format.
- Return type
- classmethod from_multipole(mp_data, search=None, log=False, **kwargs)[source]¶
Generate windowed multipole neutron data from multipole data.
- Parameters
mp_data (dictionary or str) – Dictionary or Path to the multipole data stored in a pickle file
search (bool, optional) – Whether to search for optimal window size and curvefit order. Defaults to True if no windowing parameters are specified.
log (bool or int, optional) – Whether to print running logs (use int for verbosity control)
**kwargs – Keyword arguments passed to
openmc.data.multipole._windowing()
- Returns
Resonant cross sections represented in the windowed multipole format.
- Return type