openmc.data.IncidentPhoton

class openmc.data.IncidentPhoton(atomic_number)[source]

Photon interaction data.

This class stores photo-atomic, photo-nuclear, atomic relaxation, Compton profile, stopping power, and bremsstrahlung data assembled from different sources. To create an instance, the factory method IncidentPhoton.from_endf() can be used. To add atomic relaxation or Compton profile data, set the IncidentPhoton.atomic_relaxation and IncidentPhoton.compton_profiles attributes directly.

Parameters

atomic_number (int) – Number of protons in the target nucleus

Variables
  • atomic_number (int) – Number of protons in the target nucleus

  • atomic_relaxation (openmc.data.AtomicRelaxation or None) – Atomic relaxation data

  • bremsstrahlung (dict) – Dictionary of bremsstrahlung data with keys ‘I’ (mean excitation energy in [eV]), ‘num_electrons’ (number of electrons in each subshell), ‘ionization_energy’ (ionization potential of each subshell), ‘electron_energy’ (incident electron kinetic energy values in [eV]), ‘photon_energy’ (ratio of the energy of the emitted photon to the incident electron kinetic energy), and ‘dcs’ (cross section values in [b]). The cross sections are in scaled form: \((\beta^2/Z^2) E_k (d\sigma/dE_k)\), where \(E_k\) is the energy of the emitted photon. A negative number of electrons in a subshell indicates conduction electrons.

  • compton_profiles (dict) – Dictionary of Compton profile data with keys ‘num_electrons’ (number of electrons in each subshell), ‘binding_energy’ (ionization potential of each subshell), and ‘J’ (Hartree-Fock Compton profile as a function of the projection of the electron momentum on the scattering vector, \(p_z\) for each subshell). Note that subshell occupancies may not match the atomic relaxation data.

  • reactions (dict) – Contains the cross sections for each photon reaction. The keys are MT values and the values are instances of PhotonReaction.

export_to_hdf5(path, mode='a', libver='earliest')[source]

Export incident photon 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_ace(ace_or_filename)[source]

Generate incident photon data from an ACE table

Parameters

ace_or_filename (str or openmc.data.ace.Table) – ACE table to read from. If given as a string, it is assumed to be the filename for the ACE file.

Returns

Photon interaction data

Return type

openmc.data.IncidentPhoton

classmethod from_endf(photoatomic, relaxation=None)[source]

Generate incident photon data from an ENDF evaluation

Parameters
  • photoatomic (str or openmc.data.endf.Evaluation) – ENDF photoatomic data evaluation to read from. If given as a string, it is assumed to be the filename for the ENDF file.

  • relaxation (str or openmc.data.endf.Evaluation, optional) – ENDF atomic relaxation data evaluation to read from. If given as a string, it is assumed to be the filename for the ENDF file.

Returns

Photon interaction data

Return type

openmc.data.IncidentPhoton

classmethod from_hdf5(group_or_filename)[source]

Generate photon reaction from an 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

Photon interaction data

Return type

openmc.data.IncidentPhoton