openmc.mgxs.ArbitraryMatrixXS¶
- class openmc.mgxs.ArbitraryMatrixXS(rxn_type, domain=None, domain_type=None, energy_groups=None, by_nuclide=False, name='', num_polar=1, num_azimuthal=1)[source]¶
A multi-group matrix cross section for an arbitrary reaction type.
This class can be used for both OpenMC input generation and tally data post-processing to compute spatially-homogenized and energy-integrated multi-group cross sections for multi-group neutronics calculations. At a minimum, one needs to set the
ArbitraryMatrixXS.energy_groups
andArbitraryMatrixXS.domain
properties. Tallies for the flux and appropriate reaction rates over the specified domain are generated automatically via theArbitraryMatrixXS.tallies
property, which can then be appended to aopenmc.Tallies
instance.For post-processing, the
MGXS.load_from_statepoint()
will pull in the necessary data to compute multi-group cross sections from aopenmc.StatePoint
instance. The derived multi-group cross section can then be obtained from theArbitraryMatrixXS.xs_tally
property.For a spatial domain \(V\), incoming energy group \([E_{g'},E_{g'-1}]\), and outgoing energy group \([E_g,E_{g-1}]\), the fission production is calculated as:
\[\begin{aligned} \langle \sigma_{X,g'\rightarrow g} \phi \rangle &= \int_{r \in V} dr \int_{4\pi} d\Omega' \int_{E_{g'}}^{E_{g'-1}} dE' \int_{E_g}^{E_{g-1}} dE \; \chi(E) \sigma_X (r, E') \psi(r, E', \Omega')\\ \langle \phi \rangle &= \int_{r \in V} dr \int_{4\pi} d\Omega \int_{E_g}^{E_{g-1}} dE \; \psi (r, E, \Omega) \\ \sigma_{X,g'\rightarrow g} &= \frac{\langle \sigma_{X,g'\rightarrow g} \phi \rangle}{\langle \phi \rangle} \end{aligned}\]where \(\sigma_X\) is the requested reaction type of interest.
- Parameters
domain (openmc.Material or openmc.Cell or openmc.Universe or openmc.RegularMesh) – The domain for spatial homogenization
domain_type ({'material', 'cell', 'distribcell', 'universe', 'mesh'}) – The domain type for spatial homogenization
energy_groups (openmc.mgxs.EnergyGroups) – The energy group structure for energy condensation
by_nuclide (bool) – If true, computes cross sections for each nuclide in domain
name (str, optional) – Name of the multi-group cross section. Used as a label to identify tallies in OpenMC ‘tallies.xml’ file.
num_polar (Integral, optional) – Number of equi-width polar angle bins for angle discretization; defaults to one bin
num_azimuthal (Integral, optional) – Number of equi-width azimuthal angle bins for angle discretization; defaults to one bin
- Variables
name (str, optional) – Name of the multi-group cross section
rxn_type (str) – Reaction type (e.g., ‘total’, ‘nu-fission’, etc.)
by_nuclide (bool) – If true, computes cross sections for each nuclide in domain
domain (openmc.Material or openmc.Cell or openmc.Universe or openmc.RegularMesh) – Domain for spatial homogenization
domain_type ({'material', 'cell', 'distribcell', 'universe', 'mesh'}) – Domain type for spatial homogenization
energy_groups (openmc.mgxs.EnergyGroups) – Energy group structure for energy condensation
num_polar (Integral) – Number of equi-width polar angle bins for angle discretization
num_azimuthal (Integral) – Number of equi-width azimuthal angle bins for angle discretization
tally_trigger (openmc.Trigger) – An (optional) tally precision trigger given to each tally used to compute the cross section
scores (list of str) – The scores in each tally used to compute the multi-group cross section
filters (list of openmc.Filter) – The filters in each tally used to compute the multi-group cross section
tally_keys (list of str) – The keys into the tallies dictionary for each tally used to compute the multi-group cross section
estimator ({'tracklength', 'collision', 'analog'}) – The tally estimator used to compute the multi-group cross section
tallies (dict) – OpenMC tallies needed to compute the multi-group cross section
rxn_rate_tally (openmc.Tally) – Derived tally for the reaction rate tally used in the numerator to compute the multi-group cross section. This attribute is None unless the multi-group cross section has been computed.
xs_tally (openmc.Tally) – Derived tally for the multi-group cross section. This attribute is None unless the multi-group cross section has been computed.
num_subdomains (int) – The number of subdomains is unity for ‘material’, ‘cell’ and ‘universe’ domain types. This is equal to the number of cell instances for ‘distribcell’ domain types (it is equal to unity prior to loading tally data from a statepoint file) and the number of mesh cells for ‘mesh’ domain types.
num_nuclides (int) – The number of nuclides for which the multi-group cross section is being tracked. This is unity if the by_nuclide attribute is False.
nuclides (Iterable of str or 'sum') – The optional user-specified nuclides for which to compute cross sections (e.g., ‘U238’, ‘O16’). If by_nuclide is True but nuclides are not specified by the user, all nuclides in the spatial domain are included. This attribute is ‘sum’ if by_nuclide is false.
sparse (bool) – Whether or not the MGXS’ tallies use SciPy’s LIL sparse matrix format for compressed data storage
loaded_sp (bool) – Whether or not a statepoint file has been loaded with tally data
derived (bool) – Whether or not the MGXS is merged from one or more other MGXS
mgxs_type (str) –
The name of this MGXS type, to be used when printing and indexing in an HDF5 data store
New in version 0.13.1.
- build_hdf5_store(filename='mgxs.h5', directory='mgxs', subdomains='all', nuclides='all', xs_type='macro', row_column='inout', append=True, libver='earliest')¶
Export the multi-group cross section data to an HDF5 binary file.
This method constructs an HDF5 file which stores the multi-group cross section data. The data is stored in a hierarchy of HDF5 groups from the domain type, domain id, subdomain id (for distribcell domains), nuclides and cross section type. Two datasets for the mean and standard deviation are stored for each subdomain entry in the HDF5 file.
Note
This requires the h5py Python package.
- Parameters
filename (str) – Filename for the HDF5 file. Defaults to ‘mgxs.h5’.
directory (str) – Directory for the HDF5 file. Defaults to ‘mgxs’.
subdomains (Iterable of Integral or 'all') – The subdomain IDs of the cross sections to include in the report. Defaults to ‘all’.
nuclides (Iterable of str or 'all' or 'sum') – The nuclides of the cross-sections to include in the report. This may be a list of nuclide name strings (e.g., [‘U235’, ‘U238’]). The special string ‘all’ will report the cross sections for all nuclides in the spatial domain. The special string ‘sum’ will report the cross sections summed over all nuclides. Defaults to ‘all’.
xs_type ({'macro', 'micro'}) – Store the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.
row_column ({'inout', 'outin'}) – Store scattering matrices indexed first by incoming group and second by outgoing group (‘inout’), or vice versa (‘outin’). Defaults to ‘inout’.
append (bool) – If true, appends to an existing HDF5 file with the same filename directory (if one exists). Defaults to True.
libver ({'earliest', 'latest'}) – Compatibility mode for the HDF5 file. ‘latest’ will produce files that are less backwards compatible but have performance benefits.
- Raises
ValueError – When this method is called before the multi-group cross section is computed from tally data.
- can_merge(other)¶
Determine if another MGXS can be merged with this one
If results have been loaded from a statepoint, then MGXS are only mergeable along one and only one of enegy groups or nuclides.
- Parameters
other (openmc.mgxs.MGXS) – MGXS to check for merging
- export_xs_data(filename='mgxs', directory='mgxs', format='csv', groups='all', xs_type='macro')¶
Export the multi-group cross section data to a file.
This method leverages the functionality in the Pandas library to export the multi-group cross section data in a variety of output file formats for storage and/or post-processing.
- Parameters
filename (str) – Filename for the exported file. Defaults to ‘mgxs’.
directory (str) – Directory for the exported file. Defaults to ‘mgxs’.
format ({'csv', 'excel', 'pickle', 'latex'}) – The format for the exported data file. Defaults to ‘csv’.
groups (Iterable of Integral or 'all') – Energy groups of interest. Defaults to ‘all’.
xs_type ({'macro', 'micro'}) – Store the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.
- get_condensed_xs(coarse_groups)¶
Construct an energy-condensed version of this cross section.
- Parameters
coarse_groups (openmc.mgxs.EnergyGroups) – The coarse energy group structure of interest
- Returns
A new MGXS condensed to the group structure of interest
- Return type
- get_flux(groups='all', subdomains='all', order_groups='increasing', value='mean', squeeze=True, **kwargs)¶
Returns an array of the fluxes used to weight the MGXS.
This method constructs a 2D NumPy array for the requested weighting flux for one or more subdomains (1st dimension), and energy groups (2nd dimension).
- Parameters
groups (Iterable of Integral or 'all') – Energy groups of interest. Defaults to ‘all’.
subdomains (Iterable of Integral or 'all') – Subdomain IDs of interest. Defaults to ‘all’.
order_groups ({'increasing', 'decreasing'}) – Return the cross section indexed according to increasing or decreasing energy groups (decreasing or increasing energies). Defaults to ‘increasing’.
value ({'mean', 'std_dev', 'rel_err'}) – A string for the type of value to return. Defaults to ‘mean’.
squeeze (bool) – A boolean representing whether to eliminate the extra dimensions of the multi-dimensional array to be returned. Defaults to True.
- Returns
A NumPy array of the flux indexed in the order each group and subdomain is listed in the parameters.
- Return type
- Raises
ValueError – When this method is called before the data is available from tally data, or, when this is used on an MGXS type without a flux score.
- get_homogenized_mgxs(other_mgxs)¶
Construct a homogenized mgxs with other MGXS objects.
- Parameters
other_mgxs (openmc.mgxs.MGXS or Iterable of openmc.mgxs.MGXS) – The MGXS to homogenize with this one.
- Returns
A new homogenized MGXS
- Return type
- Raises
ValueError – If the other_mgxs is of a different type.
- static get_mgxs(mgxs_type, domain=None, domain_type=None, energy_groups=None, by_nuclide=False, name='', num_polar=1, num_azimuthal=1)¶
Return a MGXS subclass object for some energy group structure within some spatial domain for some reaction type.
This is a factory method which can be used to quickly create MGXS subclass objects for various reaction types.
- Parameters
mgxs_type (str or Integral) – The type of multi-group cross section object to return; valid values are members of MGXS_TYPES, or the reaction types that are the keys of REACTION_MT. Note that if a reaction type from REACTION_MT is used, it can be appended with ‘ matrix’ to obtain a multigroup matrix (from incoming to outgoing energy groups) for reactions with a neutron in an outgoing channel.
domain (openmc.Material or openmc.Cell or openmc.Universe or openmc.RegularMesh) – The domain for spatial homogenization
domain_type ({'material', 'cell', 'distribcell', 'universe', 'mesh'}) – The domain type for spatial homogenization
energy_groups (openmc.mgxs.EnergyGroups) – The energy group structure for energy condensation
by_nuclide (bool) – If true, computes cross sections for each nuclide in domain. Defaults to False
name (str, optional) – Name of the multi-group cross section. Used as a label to identify tallies in OpenMC ‘tallies.xml’ file. Defaults to the empty string.
num_polar (Integral, optional) – Number of equi-width polar angles for angle discretization; defaults to no discretization
num_azimuthal (Integral, optional) – Number of equi-width azimuthal angles for angle discretization; defaults to no discretization
- Returns
A subclass of the abstract MGXS class for the multi-group cross section type requested by the user
- Return type
- get_nuclide_densities(nuclides='all')¶
Get an array of atomic number densities in units of atom/b-cm for all nuclides in the cross section’s spatial domain.
- Parameters
nuclides (Iterable of str or 'all' or 'sum') – A list of nuclide name strings (e.g., [‘U235’, ‘U238’]). The special string ‘all’ will return the atom densities for all nuclides in the spatial domain. The special string ‘sum’ will return the atom density summed across all nuclides in the spatial domain. Defaults to ‘all’.
- Returns
An array of the atomic number densities (atom/b-cm) for each of the nuclides in the spatial domain
- Return type
numpy.ndarray of float
- Raises
ValueError – When this method is called before the spatial domain has been set.
- get_nuclide_density(nuclide)¶
Get the atomic number density in units of atoms/b-cm for a nuclide in the cross section’s spatial domain.
- get_nuclides()¶
Get all nuclides in the cross section’s spatial domain.
- Returns
A list of the string names for each nuclide in the spatial domain (e.g., [‘U235’, ‘U238’, ‘O16’])
- Return type
list of str
- Raises
ValueError – When this method is called before the spatial domain has been set.
- get_pandas_dataframe(groups='all', nuclides='all', xs_type='macro', paths=True)¶
Build a Pandas DataFrame for the MGXS data.
This method leverages
openmc.Tally.get_pandas_dataframe()
, but renames the columns with terminology appropriate for cross section data.- Parameters
groups (Iterable of Integral or 'all') – Energy groups of interest. Defaults to ‘all’.
nuclides (Iterable of str or 'all' or 'sum') – The nuclides of the cross-sections to include in the dataframe. This may be a list of nuclide name strings (e.g., [‘U235’, ‘U238’]). The special string ‘all’ will include the cross sections for all nuclides in the spatial domain. The special string ‘sum’ will include the cross sections summed over all nuclides. Defaults to ‘all’.
xs_type ({'macro', 'micro'}) – Return macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.
paths (bool, optional) – Construct columns for distribcell tally filters (default is True). The geometric information in the Summary object is embedded into a Multi-index column with a geometric “path” to each distribcell instance.
- Returns
A Pandas DataFrame for the cross section data.
- Return type
- Raises
ValueError – When this method is called before the multi-group cross section is computed from tally data.
- get_slice(nuclides=[], in_groups=[], out_groups=[])¶
Build a sliced MatrixMGXS object for the specified nuclides and energy groups.
This method constructs a new MGXS to encapsulate a subset of the data represented by this MGXS. The subset of data to include in the tally slice is determined by the nuclides and energy groups specified in the input parameters.
- Parameters
nuclides (list of str) – A list of nuclide name strings (e.g., [‘U235’, ‘U238’]; default is [])
in_groups (list of int) – A list of incoming energy group indices starting at 1 for the high energies (e.g., [1, 2, 3]; default is [])
out_groups (list of int) – A list of outgoing energy group indices starting at 1 for the high energies (e.g., [1, 2, 3]; default is [])
- Returns
A new MatrixMGXS object which encapsulates the subset of data requested for the nuclide(s) and/or energy group(s) requested in the parameters.
- Return type
- get_subdomain_avg_xs(subdomains='all')¶
Construct a subdomain-averaged version of this cross section.
This method is useful for averaging cross sections across distribcell instances. The method performs spatial homogenization to compute the scalar flux-weighted average cross section across the subdomains.
- Parameters
subdomains (Iterable of Integral or 'all') – The subdomain IDs to average across. Defaults to ‘all’.
- Returns
A new MGXS averaged across the subdomains of interest
- Return type
- Raises
ValueError – When this method is called before the multi-group cross section is computed from tally data.
- get_units(xs_type='macro')¶
This method returns the units of a MGXS based on a desired xs_type.
- Parameters
xs_type ({'macro', 'micro'}) – Return the macro or micro cross section units. Defaults to ‘macro’.
- Returns
A string representing the units of the MGXS.
- Return type
- get_xs(in_groups='all', out_groups='all', subdomains='all', nuclides='all', xs_type='macro', order_groups='increasing', row_column='inout', value='mean', squeeze=True, **kwargs)¶
Returns an array of multi-group cross sections.
This method constructs a 4D NumPy array for the requested multi-group cross section data for one or more subdomains (1st dimension), energy groups in (2nd dimension), energy groups out (3rd dimension), and nuclides (4th dimension).
- Parameters
in_groups (Iterable of Integral or 'all') – Incoming energy groups of interest. Defaults to ‘all’.
out_groups (Iterable of Integral or 'all') – Outgoing energy groups of interest. Defaults to ‘all’.
subdomains (Iterable of Integral or 'all') – Subdomain IDs of interest. Defaults to ‘all’.
nuclides (Iterable of str or 'all' or 'sum') – A list of nuclide name strings (e.g., [‘U235’, ‘U238’]). The special string ‘all’ will return the cross sections for all nuclides in the spatial domain. The special string ‘sum’ will return the cross section summed over all nuclides. Defaults to ‘all’.
xs_type ({'macro', 'micro'}) – Return the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.
order_groups ({'increasing', 'decreasing'}) – Return the cross section indexed according to increasing or decreasing energy groups (decreasing or increasing energies). Defaults to ‘increasing’.
row_column ({'inout', 'outin'}) – Return the cross section indexed first by incoming group and second by outgoing group (‘inout’), or vice versa (‘outin’). Defaults to ‘inout’.
value ({'mean', 'std_dev', 'rel_err'}) – A string for the type of value to return. Defaults to ‘mean’.
squeeze (bool) – A boolean representing whether to eliminate the extra dimensions of the multi-dimensional array to be returned. Defaults to True.
- Returns
A NumPy array of the multi-group cross section indexed in the order each group and subdomain is listed in the parameters.
- Return type
- Raises
ValueError – When this method is called before the multi-group cross section is computed from tally data.
- load_from_statepoint(statepoint)¶
Extracts tallies in an OpenMC StatePoint with the data needed to compute multi-group cross sections.
This method is needed to compute cross section data from tallies in an OpenMC StatePoint object.
Note
The statepoint must be linked with an OpenMC Summary object.
- Parameters
statepoint (openmc.StatePoint) – An OpenMC StatePoint object with tally data
- Raises
ValueError – When this method is called with a statepoint that has not been linked with a summary object.
- merge(other)¶
Merge another MGXS with this one
MGXS are only mergeable if their energy groups and nuclides are either identical or mutually exclusive. If results have been loaded from a statepoint, then MGXS are only mergeable along one and only one of energy groups or nuclides.
- Parameters
other (openmc.mgxs.MGXS) – MGXS to merge with this one
- Returns
merged_mgxs – Merged MGXS
- Return type
- print_xs(subdomains='all', nuclides='all', xs_type='macro')¶
Prints a string representation for the multi-group cross section.
- Parameters
subdomains (Iterable of Integral or 'all') – The subdomain IDs of the cross sections to include in the report. Defaults to ‘all’.
nuclides (Iterable of str or 'all' or 'sum') – The nuclides of the cross-sections to include in the report. This may be a list of nuclide name strings (e.g., [‘U235’, ‘U238’]). The special string ‘all’ will report the cross sections for all nuclides in the spatial domain. The special string ‘sum’ will report the cross sections summed over all nuclides. Defaults to ‘all’.
xs_type ({'macro', 'micro'}) – Return the macro or micro cross section in units of cm^-1 or barns. Defaults to ‘macro’.