openmc.XSdata

class openmc.XSdata(name, energy_groups, temperatures=[294.0], representation='isotropic', num_delayed_groups=0)[source]

A multi-group cross section data set providing all the multi-group data necessary for a multi-group OpenMC calculation.

Parameters
  • name (str) – Name of the mgxs data set.

  • energy_groups (openmc.mgxs.EnergyGroups) – Energy group structure

  • representation ({'isotropic', 'angle'}, optional) – Method used in generating the MGXS (isotropic or angle-dependent flux weighting). Defaults to ‘isotropic’

  • temperatures (Iterable of float) – Temperatures (in units of Kelvin) of the provided datasets. Defaults to a single temperature at 294K.

  • num_delayed_groups (int) – Number of delayed groups

Variables
  • name (str) – Unique identifier for the xsdata object

  • atomic_weight_ratio (float) – Atomic weight ratio of an isotope. That is, the ratio of the mass of the isotope to the mass of a single neutron.

  • temperatures (numpy.ndarray) – Temperatures (in units of Kelvin) of the provided datasets. Defaults to a single temperature at 294K.

  • energy_groups (openmc.mgxs.EnergyGroups) – Energy group structure

  • num_delayed_groups (int) – Num delayed groups

  • fissionable (bool) – Whether or not this is a fissionable data set.

  • scatter_format ({'legendre', 'histogram', or 'tabular'}) – Angular distribution representation (legendre, histogram, or tabular)

  • order (int) – Either the Legendre order, number of bins, or number of points used to describe the angular distribution associated with each group-to-group transfer probability.

  • representation ({'isotropic', 'angle'}) – Method used in generating the MGXS (isotropic or angle-dependent flux weighting).

  • num_azimuthal (int) – Number of equal width angular bins that the azimuthal angular domain is subdivided into. This only applies when XSdata.representation is “angle”.

  • num_polar (int) – Number of equal width angular bins that the polar angular domain is subdivided into. This only applies when XSdata.representation is “angle”.

  • total (list of numpy.ndarray) – Group-wise total cross section.

  • absorption (list of numpy.ndarray) – Group-wise absorption cross section.

  • scatter_matrix (list of numpy.ndarray) – Scattering moment matrices presented with the columns representing incoming group and rows representing the outgoing group. That is, down-scatter will be above the diagonal of the resultant matrix.

  • multiplicity_matrix (list of numpy.ndarray) – Ratio of neutrons produced in scattering collisions to the neutrons which undergo scattering collisions; that is, the multiplicity provides the code with a scaling factor to account for neutrons produced in (n,xn) reactions.

  • fission (list of numpy.ndarray) – Group-wise fission cross section.

  • kappa_fission (list of numpy.ndarray) – Group-wise kappa_fission cross section.

  • chi (list of numpy.ndarray) – Group-wise fission spectra ordered by increasing group index (i.e., fast to thermal). This attribute should be used if making the common approximation that the fission spectra does not depend on incoming energy. If the user does not wish to make this approximation, then this should not be provided and this information included in the XSdata.nu_fission attribute instead.

  • chi_prompt (list of numpy.ndarray) – Group-wise prompt fission spectra ordered by increasing group index (i.e., fast to thermal). This attribute should be used if chi from prompt and delayed neutrons is being set separately.

  • chi_delayed (list of numpy.ndarray) – Group-wise delayed fission spectra ordered by increasing group index (i.e., fast to thermal). This attribute should be used if chi from prompt and delayed neutrons is being set separately.

  • nu_fission (list of numpy.ndarray) – Group-wise fission production cross section vector (i.e., if chi is provided), or is the group-wise fission production matrix.

  • prompt_nu_fission (list of numpy.ndarray) – Group-wise prompt fission production cross section vector.

  • delayed_nu_fission (list of numpy.ndarray) – Group-wise delayed fission production cross section vector.

  • beta (list of numpy.ndarray) – Delayed-group-wise delayed neutron fraction cross section vector.

  • decay_rate (list of numpy.ndarray) – Delayed-group-wise decay rate vector.

  • inverse_velocity (list of numpy.ndarray) – Inverse of velocity, in units of sec/cm.

  • xs_shapes (dict of iterable of int) – Dictionary with keys of _XS_SHAPES and iterable of int values with the corresponding shapes where “Order” corresponds to the pn scattering order, “G” corresponds to incoming energy group, “G’” corresponds to outgoing energy group, and “DG” corresponds to delayed group.

Notes

The parameters containing cross section data have dimensionalities which depend upon the value of XSdata.representation as well as the number of Legendre or other angular dimensions as described by XSdata.order. The XSdata.xs_shapes are provided to obtain the dimensionality of the data for each temperature.

The following are cross sections which should use each of the properties. Note that some cross sections can be input in more than one shape so they are listed multiple times:

[G][G’][Order]: scatter_matrix

[G]: total, absorption, fission, kappa_fission, nu_fission,

prompt_nu_fission, delayed_nu_fission, inverse_velocity

[G’]: chi, chi_prompt, chi_delayed

[G][G’]: multiplicity_matrix, nu_fission, prompt_nu_fission

[DG]: beta, decay_rate

[DG][G]: delayed_nu_fission, beta, decay_rate

[DG][G’]: chi_delayed

[DG][G][G’]: delayed_nu_fission

add_temperature(temperature)[source]

This method re-sizes the attributes of this XSdata object so that it can accommodate an additional temperature. Note that the set_* methods will still need to be executed.

Parameters

temperature (float) – Temperature (in units of Kelvin) of the provided dataset.

convert_representation(target_representation, num_polar=None, num_azimuthal=None)[source]

Produce a new XSdata object with the same data, but converted to the new representation (isotropic or angle-dependent).

This method cannot be used to change the number of polar or azimuthal bins of an XSdata object that already uses an angular representation. Finally, this method simply uses an arithmetic mean to convert from an angular to isotropic representation; no flux-weighting is applied and therefore reaction rates will not be preserved.

Parameters
  • target_representation ({'isotropic', 'angle'}) – Representation of the MGXS (isotropic or angle-dependent flux weighting).

  • num_polar (int, optional) – Number of equal width angular bins that the polar angular domain is subdivided into. This is required when target_representation is “angle”.

  • num_azimuthal (int, optional) – Number of equal width angular bins that the azimuthal angular domain is subdivided into. This is required when target_representation is “angle”.

Returns

Multi-group cross section data with the same data as self, but represented as specified in target_representation.

Return type

openmc.XSdata

convert_scatter_format(target_format, target_order=None)[source]

Produce a new MGXSLibrary object with the same data, but converted to the new scatter format and order

Parameters
  • target_format ({'tabular', 'legendre', 'histogram'}) – Representation of the scattering angle distribution

  • target_order (int) – Either the Legendre target_order, number of bins, or number of points used to describe the angular distribution associated with each group-to-group transfer probability

Returns

Multi-group cross section data with the same data as in self, but represented as specified in target_format.

Return type

openmc.XSdata

classmethod from_hdf5(group, name, energy_groups, num_delayed_groups)[source]

Generate XSdata object from an HDF5 group

Parameters
  • group (h5py.Group) – HDF5 group to read from

  • name (str) – Name of the mgxs data set.

  • energy_groups (openmc.mgxs.EnergyGroups) – Energy group structure

  • num_delayed_groups (int) – Number of delayed groups

Returns

Multi-group cross section data

Return type

openmc.XSdata

set_absorption(absorption, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • absorption (np.ndarray) – Absorption Cross Section

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

See also

openmc.mgxs_library.set_absorption_mgxs

set_absorption_mgxs(absorption, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.AbsorptionXS to be used to set the absorption cross section for this XSdata object.

Parameters
  • absorption (openmc.mgxs.AbsorptionXS) – MGXS Object containing the absorption cross section for the domain of interest.

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_beta(beta, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • beta (np.ndarray) – Delayed fission spectrum

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

See also

openmc.mgxs_library.set_beta_mgxs

set_beta_mgxs(beta, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.Beta to be used to set beta for this XSdata object.

Parameters
  • beta (openmc.mgxs.Beta) – MGXS Object containing beta for the domain of interest.

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_chi(chi, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • chi (np.ndarray) – Fission Spectrum

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

See also

openmc.mgxs_library.set_chi_mgxs

set_chi_delayed(chi_delayed, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • chi_delayed (np.ndarray) – Delayed fission Spectrum

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

See also

openmc.mgxs_library.set_chi_delayed_mgxs

set_chi_delayed_mgxs(chi_delayed, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.ChiDelayed to be used to set chi-delayed for this XSdata object.

Parameters
  • chi_delayed (openmc.mgxs.ChiDelayed) – MGXS Object containing chi-delayed for the domain of interest.

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_chi_mgxs(chi, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.Chi to be used to set chi for this XSdata object.

Parameters
  • chi (openmc.mgxs.Chi) – MGXS Object containing chi for the domain of interest.

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_chi_prompt(chi_prompt, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • chi_prompt (np.ndarray) – Prompt fission Spectrum

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

See also

openmc.mgxs_library.set_chi_prompt_mgxs

set_chi_prompt_mgxs(chi_prompt, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.Chi to be used to set chi-prompt for this XSdata object.

Parameters
  • chi_prompt (openmc.mgxs.Chi) – MGXS Object containing chi-prompt for the domain of interest.

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_decay_rate(decay_rate, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • decay_rate (np.ndarray) – Delayed neutron precursor decay rate

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

See also

openmc.mgxs_library.set_decay_rate_mgxs

set_decay_rate_mgxs(decay_rate, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.DecayRate to be used to set decay rate for this XSdata object.

Parameters
  • decay_rate (openmc.mgxs.DecayRate) – MGXS Object containing decay rate for the domain of interest.

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_delayed_nu_fission(delayed_nu_fission, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • delayed_nu_fission (np.ndarray) – Delayed-nu-fission Cross Section

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

See also

openmc.mgxs_library.set_delayed_nu_fission_mgxs

set_delayed_nu_fission_mgxs(delayed_nu_fission, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.DelayedNuFissionXS or openmc.mgxs.DelayedNuFissionMatrixXS to be used to set the delayed-nu-fission cross section for this XSdata object.

Parameters
  • delayed_nu_fission (openmc.mgxs.DelayedNuFissionXS or openmc.mgxs.DelayedNuFissionMatrixXS) – MGXS Object containing the delayed-nu-fission cross section for the domain of interest.

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_fission(fission, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • fission (np.ndarray) – Fission Cross Section

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

See also

openmc.mgxs_library.set_fission_mgxs

set_fission_mgxs(fission, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.FissionXS to be used to set the fission cross section for this XSdata object.

Parameters
  • fission (openmc.mgxs.FissionXS) – MGXS Object containing the fission cross section for the domain of interest.

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_inverse_velocity(inv_vel, temperature=294.0)[source]

This method sets the inverse velocity for this XSdata object at the provided temperature.

Parameters
  • inv_vel (np.ndarray) – Inverse velocity in units of sec/cm.

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

set_inverse_velocity_mgxs(inverse_velocity, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.InverseVelocity to be used to set the inverse velocity for this XSdata object.

Parameters
  • inverse_velocity (openmc.mgxs.InverseVelocity) – MGXS object containing the inverse velocity for the domain of interest.

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_kappa_fission(kappa_fission, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • kappa_fission (np.ndarray) – Kappa-Fission Cross Section

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

See also

openmc.mgxs_library.set_kappa_fission_mgxs

set_kappa_fission_mgxs(k_fission, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.KappaFissionXS to be used to set the kappa-fission cross section for this XSdata object.

Parameters
  • kappa_fission (openmc.mgxs.KappaFissionXS) – MGXS Object containing the kappa-fission cross section for the domain of interest.

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_multiplicity_matrix(multiplicity, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • multiplicity (np.ndarray) – Multiplicity Matrix Cross Section

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

See also

openmc.mgxs_library.set_multiplicity_matrix_mgxs

set_multiplicity_matrix_mgxs(nuscatter, scatter=None, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for either the direct use of only an openmc.mgxs.MultiplicityMatrixXS or an openmc.mgxs.ScatterMatrixXS and openmc.mgxs.ScatterMatrixXS to be used to set the scattering multiplicity for this XSdata object. Multiplicity, in OpenMC parlance, is a factor used to account for the production of neutrons introduced by scattering multiplication reactions, i.e., (n,xn) events. In this sense, the multiplication matrix is simply defined as the ratio of the nu-scatter and scatter matrices.

Parameters
  • nuscatter (openmc.mgxs.ScatterMatrixXS or openmc.mgxs.MultiplicityMatrixXS) – MGXS Object containing the matrix cross section for the domain of interest.

  • scatter (openmc.mgxs.ScatterMatrixXS) – MGXS Object containing the scattering matrix cross section for the domain of interest.

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_nu_fission(nu_fission, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • nu_fission (np.ndarray) – Nu-fission Cross Section

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

See also

openmc.mgxs_library.set_nu_fission_mgxs

set_nu_fission_mgxs(nu_fission, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.FissionXS to be used to set the nu-fission cross section for this XSdata object.

Parameters
  • nu_fission (openmc.mgxs.FissionXS) – MGXS Object containing the nu-fission cross section for the domain of interest.

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_prompt_nu_fission(prompt_nu_fission, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • prompt_nu_fission (np.ndarray) – Prompt-nu-fission Cross Section

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

See also

openmc.mgxs_library.set_prompt_nu_fission_mgxs

set_prompt_nu_fission_mgxs(prompt_nu_fission, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

Sets the prompt-nu-fission cross section.

This method allows for an openmc.mgxs.FissionXS or openmc.mgxs.NuFissionMatrixXS to be used to set the prompt-nu-fission cross section for this XSdata object.

Parameters
  • prompt_nu_fission (openmc.mgxs.FissionXS or openmc.mgxs.NuFissionMatrixXS) – MGXS Object containing the prompt-nu-fission cross section for the domain of interest.

  • temperature (float) – Temperature (in units of Kelvin) of the provided dataset. Defaults to 294K

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_scatter_matrix(scatter, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • scatter (np.ndarray) – Scattering Matrix Cross Section

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

See also

openmc.mgxs_library.set_scatter_matrix_mgxs

set_scatter_matrix_mgxs(scatter, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.ScatterMatrixXS to be used to set the scatter matrix cross section for this XSdata object. If the XSdata.order attribute has not yet been set, then it will be set based on the properties of scatter.

Parameters
  • scatter (openmc.mgxs.ScatterMatrixXS) – MGXS Object containing the scatter matrix cross section for the domain of interest.

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

set_total(total, temperature=294.0)[source]

This method sets the cross section for this XSdata object at the provided temperature.

Parameters
  • total (np.ndarray) – Total Cross Section

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

See also

openmc.mgxs_library.set_total_mgxs

set_total_mgxs(total, temperature=294.0, nuclide='total', xs_type='macro', subdomain=None)[source]

This method allows for an openmc.mgxs.TotalXS or openmc.mgxs.TransportXS to be used to set the total cross section for this XSdata object.

Parameters
  • total (openmc.mgxs.TotalXS or openmc.mgxs.TransportXS) – MGXS Object containing the total, transport or nu-transport cross section for the domain of interest.

  • temperature (float) – Temperature (in Kelvin) of the data. Defaults to room temperature (294K).

  • nuclide (str) – Individual nuclide (or ‘total’ if obtaining material-wise data) to gather data for. Defaults to ‘total’.

  • xs_type ({'macro', 'micro'}) – Provide the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.

  • subdomain (iterable of int) – If the MGXS contains a mesh domain type, the subdomain parameter specifies which mesh cell (i.e., [i, j, k] index) to use.

to_hdf5(file)[source]

Write XSdata to an HDF5 file

Parameters

file (h5py.File) – HDF5 File (a root Group) to write to