1import copy
2from numbers import Real, Integral
3import os
4
5import h5py
6import numpy as np
7from scipy.integrate import simps
8from scipy.interpolate import interp1d
9from scipy.special import eval_legendre
10
11import openmc
12import openmc.mgxs
13from openmc.mgxs import SCATTER_TABULAR, SCATTER_LEGENDRE, SCATTER_HISTOGRAM
14from .checkvalue import check_type, check_value, check_greater_than, \
15    check_iterable_type, check_less_than, check_filetype_version
16
17ROOM_TEMPERATURE_KELVIN = 294.0
18
19# Supported incoming particle MGXS angular treatment representations
20REPRESENTATION_ISOTROPIC = 'isotropic'
21REPRESENTATION_ANGLE = 'angle'
22_REPRESENTATIONS = [
23    REPRESENTATION_ISOTROPIC,
24    REPRESENTATION_ANGLE
25]
26
27# Supported scattering angular distribution representations
28_SCATTER_TYPES = [
29    SCATTER_TABULAR,
30    SCATTER_LEGENDRE,
31    SCATTER_HISTOGRAM
32]
33
34# List of MGXS indexing schemes
35_XS_SHAPES = ["[G][G'][Order]", "[G]", "[G']", "[G][G']", "[DG]", "[DG][G]",
36              "[DG][G']", "[DG][G][G']"]
37
38# Number of mu points for conversion between scattering formats
39_NMU = 257
40
41# Filetype name of the MGXS Library
42_FILETYPE_MGXS_LIBRARY = 'mgxs'
43
44# Current version of the MGXS Library Format
45_VERSION_MGXS_LIBRARY = 1
46
47
48class XSdata:
49    """A multi-group cross section data set providing all the
50    multi-group data necessary for a multi-group OpenMC calculation.
51
52    Parameters
53    ----------
54    name : str
55        Name of the mgxs data set.
56    energy_groups : openmc.mgxs.EnergyGroups
57        Energy group structure
58    representation : {'isotropic', 'angle'}, optional
59        Method used in generating the MGXS (isotropic or angle-dependent flux
60        weighting). Defaults to 'isotropic'
61    temperatures : Iterable of float
62        Temperatures (in units of Kelvin) of the provided datasets.  Defaults
63        to a single temperature at 294K.
64    num_delayed_groups : int
65        Number of delayed groups
66
67    Attributes
68    ----------
69    name : str
70        Unique identifier for the xsdata object
71    atomic_weight_ratio : float
72        Atomic weight ratio of an isotope.  That is, the ratio of the mass
73        of the isotope to the mass of a single neutron.
74    temperatures : numpy.ndarray
75        Temperatures (in units of Kelvin) of the provided datasets.  Defaults
76        to a single temperature at 294K.
77    energy_groups : openmc.mgxs.EnergyGroups
78        Energy group structure
79    num_delayed_groups : int
80        Num delayed groups
81    fissionable : bool
82        Whether or not this is a fissionable data set.
83    scatter_format : {'legendre', 'histogram', or 'tabular'}
84        Angular distribution representation (legendre, histogram, or tabular)
85    order : int
86        Either the Legendre order, number of bins, or number of points used to
87        describe the angular distribution associated with each group-to-group
88        transfer probability.
89    representation : {'isotropic', 'angle'}
90        Method used in generating the MGXS (isotropic or angle-dependent flux
91        weighting).
92    num_azimuthal : int
93        Number of equal width angular bins that the azimuthal angular domain is
94        subdivided into. This only applies when :attr:`XSdata.representation`
95        is "angle".
96    num_polar : int
97        Number of equal width angular bins that the polar angular domain is
98        subdivided into. This only applies when :attr:`XSdata.representation`
99        is "angle".
100    total : list of numpy.ndarray
101        Group-wise total cross section.
102    absorption : list of numpy.ndarray
103        Group-wise absorption cross section.
104    scatter_matrix : list of numpy.ndarray
105        Scattering moment matrices presented with the columns representing
106        incoming group and rows representing the outgoing group.  That is,
107        down-scatter will be above the diagonal of the resultant matrix.
108    multiplicity_matrix : list of numpy.ndarray
109        Ratio of neutrons produced in scattering collisions to the neutrons
110        which undergo scattering collisions; that is, the multiplicity provides
111        the code with a scaling factor to account for neutrons produced in
112        (n,xn) reactions.
113    fission : list of numpy.ndarray
114        Group-wise fission cross section.
115    kappa_fission : list of numpy.ndarray
116        Group-wise kappa_fission cross section.
117    chi : list of numpy.ndarray
118        Group-wise fission spectra ordered by increasing group index (i.e.,
119        fast to thermal). This attribute should be used if making the common
120        approximation that the fission spectra does not depend on incoming
121        energy. If the user does not wish to make this approximation, then
122        this should not be provided and this information included in the
123        :attr:`XSdata.nu_fission` attribute instead.
124    chi_prompt : list of numpy.ndarray
125        Group-wise prompt fission spectra ordered by increasing group index
126        (i.e., fast to thermal). This attribute should be used if chi from
127        prompt and delayed neutrons is being set separately.
128    chi_delayed : list of numpy.ndarray
129        Group-wise delayed fission spectra ordered by increasing group index
130        (i.e., fast to thermal). This attribute should be used if chi from
131        prompt and delayed neutrons is being set separately.
132    nu_fission : list of numpy.ndarray
133        Group-wise fission production cross section vector (i.e., if ``chi`` is
134        provided), or is the group-wise fission production matrix.
135    prompt_nu_fission : list of numpy.ndarray
136        Group-wise prompt fission production cross section vector.
137    delayed_nu_fission : list of numpy.ndarray
138        Group-wise delayed fission production cross section vector.
139    beta : list of numpy.ndarray
140        Delayed-group-wise delayed neutron fraction cross section vector.
141    decay_rate : list of numpy.ndarray
142        Delayed-group-wise decay rate vector.
143    inverse_velocity : list of numpy.ndarray
144        Inverse of velocity, in units of sec/cm.
145    xs_shapes : dict of iterable of int
146        Dictionary with keys of _XS_SHAPES and iterable of int values with the
147        corresponding shapes where "Order" corresponds to the pn scattering
148        order, "G" corresponds to incoming energy group, "G'" corresponds to
149        outgoing energy group, and "DG" corresponds to delayed group.
150
151    Notes
152    -----
153    The parameters containing cross section data have dimensionalities which
154    depend upon the value of :attr:`XSdata.representation` as well as the
155    number of Legendre or other angular dimensions as described by
156    :attr:`XSdata.order`. The :attr:`XSdata.xs_shapes` are provided to obtain
157    the dimensionality of the data for each temperature.
158
159    The following are cross sections which should use each of the properties.
160    Note that some cross sections can be input in more than one shape so they
161    are listed multiple times:
162
163    [G][G'][Order]: scatter_matrix
164
165    [G]: total, absorption, fission, kappa_fission, nu_fission,
166         prompt_nu_fission, delayed_nu_fission, inverse_velocity
167
168    [G']: chi, chi_prompt, chi_delayed
169
170    [G][G']: multiplicity_matrix, nu_fission, prompt_nu_fission
171
172    [DG]: beta, decay_rate
173
174    [DG][G]: delayed_nu_fission, beta, decay_rate
175
176    [DG][G']: chi_delayed
177
178    [DG][G][G']: delayed_nu_fission
179
180    """
181
182    def __init__(self, name, energy_groups, temperatures=[ROOM_TEMPERATURE_KELVIN],
183                 representation=REPRESENTATION_ISOTROPIC, num_delayed_groups=0):
184
185        # Initialize class attributes
186        self.name = name
187        self.energy_groups = energy_groups
188        self.num_delayed_groups = num_delayed_groups
189        self.temperatures = temperatures
190        self.representation = representation
191        self._atomic_weight_ratio = None
192        self._fissionable = False
193        self._scatter_format = SCATTER_LEGENDRE
194        self._order = None
195        self._num_polar = None
196        self._num_azimuthal = None
197        self._total = len(temperatures) * [None]
198        self._absorption = len(temperatures) * [None]
199        self._scatter_matrix = len(temperatures) * [None]
200        self._multiplicity_matrix = len(temperatures) * [None]
201        self._fission = len(temperatures) * [None]
202        self._nu_fission = len(temperatures) * [None]
203        self._prompt_nu_fission = len(temperatures) * [None]
204        self._delayed_nu_fission = len(temperatures) * [None]
205        self._kappa_fission = len(temperatures) * [None]
206        self._chi = len(temperatures) * [None]
207        self._chi_prompt = len(temperatures) * [None]
208        self._chi_delayed = len(temperatures) * [None]
209        self._beta = len(temperatures) * [None]
210        self._decay_rate = len(temperatures) * [None]
211        self._inverse_velocity = len(temperatures) * [None]
212        self._xs_shapes = None
213
214    def __deepcopy__(self, memo):
215        existing = memo.get(id(self))
216
217        # If this is the first time we have tried to copy this object, copy it
218        if existing is None:
219            clone = type(self).__new__(type(self))
220            clone._name = self.name
221            clone._energy_groups = copy.deepcopy(self.energy_groups, memo)
222            clone._num_delayed_groups = self.num_delayed_groups
223            clone._temperatures = copy.deepcopy(self.temperatures, memo)
224            clone._representation = self.representation
225            clone._atomic_weight_ratio = self._atomic_weight_ratio
226            clone._fissionable = self._fissionable
227            clone._scatter_format = self._scatter_format
228            clone._order = self._order
229            clone._num_polar = self._num_polar
230            clone._num_azimuthal = self._num_azimuthal
231            clone._total = copy.deepcopy(self._total, memo)
232            clone._absorption = copy.deepcopy(self._absorption, memo)
233            clone._scatter_matrix = copy.deepcopy(self._scatter_matrix, memo)
234            clone._multiplicity_matrix = \
235                copy.deepcopy(self._multiplicity_matrix, memo)
236            clone._fission = copy.deepcopy(self._fission, memo)
237            clone._nu_fission = copy.deepcopy(self._nu_fission, memo)
238            clone._prompt_nu_fission = \
239                copy.deepcopy(self._prompt_nu_fission, memo)
240            clone._delayed_nu_fission = \
241                copy.deepcopy(self._delayed_nu_fission, memo)
242            clone._kappa_fission = copy.deepcopy(self._kappa_fission, memo)
243            clone._chi = copy.deepcopy(self._chi, memo)
244            clone._chi_prompt = copy.deepcopy(self._chi_prompt, memo)
245            clone._chi_delayed = copy.deepcopy(self._chi_delayed, memo)
246            clone._beta = copy.deepcopy(self._beta, memo)
247            clone._decay_rate = copy.deepcopy(self._decay_rate, memo)
248            clone._inverse_velocity = \
249                copy.deepcopy(self._inverse_velocity, memo)
250            clone._xs_shapes = copy.deepcopy(self._xs_shapes, memo)
251
252            memo[id(self)] = clone
253
254            return clone
255
256        # If this object has been copied before, return the first copy made
257        else:
258            return existing
259
260    @property
261    def name(self):
262        return self._name
263
264    @property
265    def energy_groups(self):
266        return self._energy_groups
267
268    @property
269    def num_delayed_groups(self):
270        return self._num_delayed_groups
271
272    @property
273    def representation(self):
274        return self._representation
275
276    @property
277    def atomic_weight_ratio(self):
278        return self._atomic_weight_ratio
279
280    @property
281    def fissionable(self):
282        return self._fissionable
283
284    @property
285    def temperatures(self):
286        return self._temperatures
287
288    @property
289    def scatter_format(self):
290        return self._scatter_format
291
292    @property
293    def order(self):
294        return self._order
295
296    @property
297    def num_polar(self):
298        return self._num_polar
299
300    @property
301    def num_azimuthal(self):
302        return self._num_azimuthal
303
304    @property
305    def total(self):
306        return self._total
307
308    @property
309    def absorption(self):
310        return self._absorption
311
312    @property
313    def scatter_matrix(self):
314        return self._scatter_matrix
315
316    @property
317    def multiplicity_matrix(self):
318        return self._multiplicity_matrix
319
320    @property
321    def fission(self):
322        return self._fission
323
324    @property
325    def nu_fission(self):
326        return self._nu_fission
327
328    @property
329    def prompt_nu_fission(self):
330        return self._prompt_nu_fission
331
332    @property
333    def delayed_nu_fission(self):
334        return self._delayed_nu_fission
335
336    @property
337    def kappa_fission(self):
338        return self._kappa_fission
339
340    @property
341    def chi(self):
342        return self._chi
343
344    @property
345    def chi_prompt(self):
346        return self._chi_prompt
347
348    @property
349    def chi_delayed(self):
350        return self._chi_delayed
351
352    @property
353    def beta(self):
354        return self._beta
355
356    @property
357    def decay_rate(self):
358        return self._decay_rate
359
360    @property
361    def inverse_velocity(self):
362        return self._inverse_velocity
363
364    @property
365    def num_orders(self):
366        if self._order is None:
367            raise ValueError('Order has not been set.')
368
369        if self._scatter_format in (None, SCATTER_LEGENDRE):
370            return self._order + 1
371        else:
372            return self._order
373
374    @property
375    def xs_shapes(self):
376
377        if self._xs_shapes is None:
378
379            self._xs_shapes = {}
380            self._xs_shapes["[G]"] = (self.energy_groups.num_groups,)
381            self._xs_shapes["[G']"] = (self.energy_groups.num_groups,)
382            self._xs_shapes["[G][G']"] = (self.energy_groups.num_groups,
383                                          self.energy_groups.num_groups)
384            self._xs_shapes["[DG]"] = (self.num_delayed_groups,)
385            self._xs_shapes["[DG][G]"] = (self.num_delayed_groups,
386                                          self.energy_groups.num_groups)
387            self._xs_shapes["[DG][G']"] = (self.num_delayed_groups,
388                                           self.energy_groups.num_groups)
389            self._xs_shapes["[DG][G][G']"] = (self.num_delayed_groups,
390                                              self.energy_groups.num_groups,
391                                              self.energy_groups.num_groups)
392
393            self._xs_shapes["[G][G'][Order]"] \
394                = (self.energy_groups.num_groups,
395                   self.energy_groups.num_groups, self.num_orders)
396
397            # If representation is by angle prepend num polar and num azim
398            if self.representation == REPRESENTATION_ANGLE:
399                for key, shapes in self._xs_shapes.items():
400                    self._xs_shapes[key] \
401                        = (self.num_polar, self.num_azimuthal) + shapes
402
403        return self._xs_shapes
404
405    @name.setter
406    def name(self, name):
407
408        check_type('name for XSdata', name, str)
409        self._name = name
410
411    @energy_groups.setter
412    def energy_groups(self, energy_groups):
413
414        check_type('energy_groups', energy_groups, openmc.mgxs.EnergyGroups)
415        if energy_groups.group_edges is None:
416            msg = 'Unable to assign an EnergyGroups object ' \
417                  'with uninitialized group edges'
418            raise ValueError(msg)
419
420        self._energy_groups = energy_groups
421
422    @num_delayed_groups.setter
423    def num_delayed_groups(self, num_delayed_groups):
424
425        check_type('num_delayed_groups', num_delayed_groups, Integral)
426        check_less_than('num_delayed_groups', num_delayed_groups,
427                        openmc.mgxs.MAX_DELAYED_GROUPS, equality=True)
428        check_greater_than('num_delayed_groups', num_delayed_groups, 0,
429                           equality=True)
430        self._num_delayed_groups = num_delayed_groups
431
432    @representation.setter
433    def representation(self, representation):
434
435        check_value('representation', representation, _REPRESENTATIONS)
436        self._representation = representation
437
438    @atomic_weight_ratio.setter
439    def atomic_weight_ratio(self, atomic_weight_ratio):
440
441        check_type('atomic_weight_ratio', atomic_weight_ratio, Real)
442        check_greater_than('atomic_weight_ratio', atomic_weight_ratio, 0.0)
443        self._atomic_weight_ratio = atomic_weight_ratio
444
445    @temperatures.setter
446    def temperatures(self, temperatures):
447
448        check_iterable_type('temperatures', temperatures, Real)
449        self._temperatures = np.array(temperatures)
450
451    @scatter_format.setter
452    def scatter_format(self, scatter_format):
453
454        check_value('scatter_format', scatter_format, _SCATTER_TYPES)
455        self._scatter_format = scatter_format
456
457    @order.setter
458    def order(self, order):
459
460        check_type('order', order, Integral)
461        check_greater_than('order', order, 0, equality=True)
462        self._order = order
463
464    @num_polar.setter
465    def num_polar(self, num_polar):
466
467        check_type('num_polar', num_polar, Integral)
468        check_greater_than('num_polar', num_polar, 0)
469        self._num_polar = num_polar
470
471    @num_azimuthal.setter
472    def num_azimuthal(self, num_azimuthal):
473
474        check_type('num_azimuthal', num_azimuthal, Integral)
475        check_greater_than('num_azimuthal', num_azimuthal, 0)
476        self._num_azimuthal = num_azimuthal
477
478    def add_temperature(self, temperature):
479        """This method re-sizes the attributes of this XSdata object so that it
480        can accomodate an additional temperature.  Note that the set_* methods
481        will still need to be executed.
482
483        Parameters
484        ----------
485        temperature : float
486            Temperature (in units of Kelvin) of the provided dataset.
487
488        """
489
490        check_type('temperature', temperature, Real)
491
492        temp_store = self.temperatures.tolist().append(temperature)
493        self.temperatures = temp_store
494
495        self._total.append(None)
496        self._absorption.append(None)
497        self._scatter_matrix.append(None)
498        self._multiplicity_matrix.append(None)
499        self._fission.append(None)
500        self._nu_fission.append(None)
501        self._prompt_nu_fission.append(None)
502        self._delayed_nu_fission.append(None)
503        self._kappa_fission.append(None)
504        self._chi.append(None)
505        self._chi_prompt.append(None)
506        self._chi_delayed.append(None)
507        self._beta.append(None)
508        self._decay_rate.append(None)
509        self._inverse_velocity.append(None)
510
511    def _check_temperature(self, temperature):
512        check_type('temperature', temperature, Real)
513        check_value('temperature', temperature, self.temperatures)
514
515    def _temperature_index(self, temperature):
516        return np.where(self.temperatures == temperature)[0][0]
517
518    def _set_fissionable(self, array):
519        if np.sum(array) > 0:
520            self._fissionable = True
521
522    def set_total(self, total, temperature=ROOM_TEMPERATURE_KELVIN):
523        """This method sets the cross section for this XSdata object at the
524        provided temperature.
525
526        Parameters
527        ----------
528        total: np.ndarray
529            Total Cross Section
530        temperature : float
531            Temperature (in Kelvin) of the data. Defaults to room temperature
532            (294K).
533
534        See also
535        --------
536        openmc.mgxs_library.set_total_mgxs()
537
538        """
539
540        # Get the accepted shapes for this xs
541        shapes = [self.xs_shapes["[G]"]]
542
543        # Convert to a numpy array so we can easily get the shape for checking
544        total = np.asarray(total)
545        check_value('total shape', total.shape, shapes)
546        self._check_temperature(temperature)
547
548        i = self._temperature_index(temperature)
549        self._total[i] = total
550
551    def set_absorption(self, absorption, temperature=ROOM_TEMPERATURE_KELVIN):
552        """This method sets the cross section for this XSdata object at the
553        provided temperature.
554
555        Parameters
556        ----------
557        absorption: np.ndarray
558            Absorption Cross Section
559        temperature : float
560            Temperature (in Kelvin) of the data. Defaults to room temperature
561            (294K).
562
563        See also
564        --------
565        openmc.mgxs_library.set_absorption_mgxs()
566
567        """
568
569        # Get the accepted shapes for this xs
570        shapes = [self.xs_shapes["[G]"]]
571
572        # Convert to a numpy array so we can easily get the shape for checking
573        absorption = np.asarray(absorption)
574        check_value('absorption shape', absorption.shape, shapes)
575        self._check_temperature(temperature)
576
577        i = self._temperature_index(temperature)
578        self._absorption[i] = absorption
579
580    def set_fission(self, fission, temperature=ROOM_TEMPERATURE_KELVIN):
581        """This method sets the cross section for this XSdata object at the
582        provided temperature.
583
584        Parameters
585        ----------
586        fission: np.ndarray
587            Fission Cross Section
588        temperature : float
589            Temperature (in Kelvin) of the data. Defaults to room temperature
590            (294K).
591
592        See also
593        --------
594        openmc.mgxs_library.set_fission_mgxs()
595
596        """
597
598        # Get the accepted shapes for this xs
599        shapes = [self.xs_shapes["[G]"]]
600
601        # Convert to a numpy array so we can easily get the shape for checking
602        fission = np.asarray(fission)
603        check_value('fission shape', fission.shape, shapes)
604        self._check_temperature(temperature)
605
606        i = self._temperature_index(temperature)
607        self._fission[i] = fission
608
609        self._set_fissionable(fission)
610
611    def set_kappa_fission(self, kappa_fission, temperature=ROOM_TEMPERATURE_KELVIN):
612        """This method sets the cross section for this XSdata object at the
613        provided temperature.
614
615        Parameters
616        ----------
617        kappa_fission: np.ndarray
618            Kappa-Fission Cross Section
619        temperature : float
620            Temperature (in Kelvin) of the data. Defaults to room temperature
621            (294K).
622
623        See also
624        --------
625        openmc.mgxs_library.set_kappa_fission_mgxs()
626
627        """
628
629        # Get the accepted shapes for this xs
630        shapes = [self.xs_shapes["[G]"]]
631
632        # Convert to a numpy array so we can easily get the shape for checking
633        kappa_fission = np.asarray(kappa_fission)
634        check_value('kappa fission shape', kappa_fission.shape, shapes)
635        self._check_temperature(temperature)
636
637        i = self._temperature_index(temperature)
638        self._kappa_fission[i] = kappa_fission
639
640        self._set_fissionable(kappa_fission)
641
642    def set_chi(self, chi, temperature=ROOM_TEMPERATURE_KELVIN):
643        """This method sets the cross section for this XSdata object at the
644        provided temperature.
645
646        Parameters
647        ----------
648        chi: np.ndarray
649            Fission Spectrum
650        temperature : float
651            Temperature (in Kelvin) of the data. Defaults to room temperature
652            (294K).
653
654        See also
655        --------
656        openmc.mgxs_library.set_chi_mgxs()
657
658        """
659
660        # Get the accepted shapes for this xs
661        shapes = [self.xs_shapes["[G']"]]
662
663        # Convert to a numpy array so we can easily get the shape for checking
664        chi = np.asarray(chi)
665        check_value('chi shape', chi.shape, shapes)
666        self._check_temperature(temperature)
667
668        i = self._temperature_index(temperature)
669        self._chi[i] = chi
670
671    def set_chi_prompt(self, chi_prompt, temperature=ROOM_TEMPERATURE_KELVIN):
672        """This method sets the cross section for this XSdata object at the
673        provided temperature.
674
675        Parameters
676        ----------
677        chi_prompt : np.ndarray
678            Prompt fission Spectrum
679        temperature : float
680            Temperature (in units of Kelvin) of the provided dataset. Defaults
681            to 294K
682
683        See also
684        --------
685        openmc.mgxs_library.set_chi_prompt_mgxs()
686
687        """
688
689        # Get the accepted shapes for this xs
690        shapes = [self.xs_shapes["[G']"]]
691
692        # Convert to a numpy array so we can easily get the shape for checking
693        chi_prompt = np.asarray(chi_prompt)
694        check_value('chi prompt shape', chi_prompt.shape, shapes)
695        self._check_temperature(temperature)
696
697        i = self._temperature_index(temperature)
698        self._chi_prompt[i] = chi_prompt
699
700    def set_chi_delayed(self, chi_delayed, temperature=ROOM_TEMPERATURE_KELVIN):
701        """This method sets the cross section for this XSdata object at the
702        provided temperature.
703
704        Parameters
705        ----------
706        chi_delayed : np.ndarray
707            Delayed fission Spectrum
708        temperature : float
709            Temperature (in units of Kelvin) of the provided dataset. Defaults
710            to 294K
711
712        See also
713        --------
714        openmc.mgxs_library.set_chi_delayed_mgxs()
715
716        """
717
718        # Get the accepted shapes for this xs
719        shapes = [self.xs_shapes["[G']"], self.xs_shapes["[DG][G']"]]
720
721        # Convert to a numpy array so we can easily get the shape for checking
722        chi_delayed = np.asarray(chi_delayed)
723        check_value('chi delayed shape', chi_delayed.shape, shapes)
724        self._check_temperature(temperature)
725        check_type('temperature', temperature, Real)
726        check_value('temperature', temperature, self.temperatures)
727
728        i = self._temperature_index(temperature)
729        self._chi_delayed[i] = chi_delayed
730
731    def set_beta(self, beta, temperature=ROOM_TEMPERATURE_KELVIN):
732        """This method sets the cross section for this XSdata object at the
733        provided temperature.
734
735        Parameters
736        ----------
737        beta : np.ndarray
738            Delayed fission spectrum
739        temperature : float
740            Temperature (in units of Kelvin) of the provided dataset. Defaults
741            to 294K
742
743        See also
744        --------
745        openmc.mgxs_library.set_beta_mgxs()
746
747        """
748
749        # Get the accepted shapes for this xs
750        shapes = [self.xs_shapes["[DG]"], self.xs_shapes["[DG][G]"]]
751
752        # Convert to a numpy array so we can easily get the shape for checking
753        beta = np.asarray(beta)
754        check_value('beta shape', beta.shape, shapes)
755        self._check_temperature(temperature)
756
757        i = self._temperature_index(temperature)
758        self._beta[i] = beta
759
760    def set_decay_rate(self, decay_rate, temperature=ROOM_TEMPERATURE_KELVIN):
761        """This method sets the cross section for this XSdata object at the
762        provided temperature.
763
764        Parameters
765        ----------
766        decay_rate : np.ndarray
767            Delayed neutron precursor decay rate
768        temperature : float
769            Temperature (in units of Kelvin) of the provided dataset. Defaults
770            to 294K
771
772        See also
773        --------
774        openmc.mgxs_library.set_decay_rate_mgxs()
775
776        """
777
778        # Get the accepted shapes for this xs
779        shapes = [self.xs_shapes["[DG]"]]
780
781        # Convert to a numpy array so we can easily get the shape for checking
782        decay_rate = np.asarray(decay_rate)
783        check_value('decay rate shape', decay_rate.shape, shapes)
784        self._check_temperature(temperature)
785
786        i = self._temperature_index(temperature)
787        self._decay_rate[i] = decay_rate
788
789    def set_scatter_matrix(self, scatter, temperature=ROOM_TEMPERATURE_KELVIN):
790        """This method sets the cross section for this XSdata object at the
791        provided temperature.
792
793        Parameters
794        ----------
795        scatter: np.ndarray
796            Scattering Matrix Cross Section
797        temperature : float
798            Temperature (in Kelvin) of the data. Defaults to room temperature
799            (294K).
800
801        See also
802        --------
803        openmc.mgxs_library.set_scatter_matrix_mgxs()
804
805        """
806
807        # Get the accepted shapes for this xs
808        shapes = [self.xs_shapes["[G][G'][Order]"]]
809
810        # Convert to a numpy array so we can easily get the shape for checking
811        scatter = np.asarray(scatter)
812        check_iterable_type('scatter', scatter, Real,
813                            max_depth=len(scatter.shape))
814        check_value('scatter shape', scatter.shape, shapes)
815        self._check_temperature(temperature)
816
817        i = self._temperature_index(temperature)
818        self._scatter_matrix[i] = scatter
819
820    def set_multiplicity_matrix(self, multiplicity, temperature=ROOM_TEMPERATURE_KELVIN):
821        """This method sets the cross section for this XSdata object at the
822        provided temperature.
823
824        Parameters
825        ----------
826        multiplicity: np.ndarray
827            Multiplicity Matrix Cross Section
828        temperature : float
829            Temperature (in Kelvin) of the data. Defaults to room temperature
830            (294K).
831
832        See also
833        --------
834        openmc.mgxs_library.set_multiplicity_matrix_mgxs()
835
836        """
837
838        # Get the accepted shapes for this xs
839        shapes = [self.xs_shapes["[G][G']"]]
840
841        # Convert to a numpy array so we can easily get the shape for checking
842        multiplicity = np.asarray(multiplicity)
843        check_iterable_type('multiplicity', multiplicity, Real,
844                            max_depth=len(multiplicity.shape))
845        check_value('multiplicity shape', multiplicity.shape, shapes)
846        self._check_temperature(temperature)
847
848        i = self._temperature_index(temperature)
849        self._multiplicity_matrix[i] = multiplicity
850
851    def set_nu_fission(self, nu_fission, temperature=ROOM_TEMPERATURE_KELVIN):
852        """This method sets the cross section for this XSdata object at the
853        provided temperature.
854
855        Parameters
856        ----------
857        nu_fission: np.ndarray
858            Nu-fission Cross Section
859        temperature : float
860            Temperature (in Kelvin) of the data. Defaults to room temperature
861            (294K).
862
863        See also
864        --------
865        openmc.mgxs_library.set_nu_fission_mgxs()
866
867        """
868
869        # Get the accepted shapes for this xs
870        shapes = [self.xs_shapes["[G]"], self.xs_shapes["[G][G']"]]
871
872        # Convert to a numpy array so we can easily get the shape for checking
873        nu_fission = np.asarray(nu_fission)
874        check_value('nu_fission shape', nu_fission.shape, shapes)
875        check_iterable_type('nu_fission', nu_fission, Real,
876                            max_depth=len(nu_fission.shape))
877        self._check_temperature(temperature)
878
879        i = self._temperature_index(temperature)
880        self._nu_fission[i] = nu_fission
881        self._set_fissionable(nu_fission)
882
883    def set_prompt_nu_fission(self, prompt_nu_fission, temperature=ROOM_TEMPERATURE_KELVIN):
884        """This method sets the cross section for this XSdata object at the
885        provided temperature.
886
887        Parameters
888        ----------
889        prompt_nu_fission: np.ndarray
890            Prompt-nu-fission Cross Section
891        temperature : float
892            Temperature (in units of Kelvin) of the provided dataset. Defaults
893            to 294K
894
895        See also
896        --------
897        openmc.mgxs_library.set_prompt_nu_fission_mgxs()
898
899        """
900
901        # Get the accepted shapes for this xs
902        shapes = [self.xs_shapes["[G]"], self.xs_shapes["[G][G']"]]
903
904        # Convert to a numpy array so we can easily get the shape for checking
905        prompt_nu_fission = np.asarray(prompt_nu_fission)
906        check_value('prompt_nu_fission shape', prompt_nu_fission.shape, shapes)
907        check_iterable_type('prompt_nu_fission', prompt_nu_fission, Real,
908                            max_depth=len(prompt_nu_fission.shape))
909        self._check_temperature(temperature)
910
911        i = self._temperature_index(temperature)
912        self._prompt_nu_fission[i] = prompt_nu_fission
913        self._set_fissionable(prompt_nu_fission)
914
915    def set_delayed_nu_fission(self, delayed_nu_fission, temperature=ROOM_TEMPERATURE_KELVIN):
916        """This method sets the cross section for this XSdata object at the
917        provided temperature.
918
919        Parameters
920        ----------
921        delayed_nu_fission: np.ndarray
922            Delayed-nu-fission Cross Section
923        temperature : float
924            Temperature (in units of Kelvin) of the provided dataset. Defaults
925            to 294K
926
927        See also
928        --------
929        openmc.mgxs_library.set_delayed_nu_fission_mgxs()
930
931        """
932
933        # Get the accepted shapes for this xs
934        shapes = [self.xs_shapes["[DG][G]"], self.xs_shapes["[DG][G][G']"]]
935
936        # Convert to a numpy array so we can easily get the shape for checking
937        delayed_nu_fission = np.asarray(delayed_nu_fission)
938        check_value('delayed_nu_fission shape', delayed_nu_fission.shape,
939                    shapes)
940        check_iterable_type('delayed_nu_fission', delayed_nu_fission, Real,
941                            max_depth=len(delayed_nu_fission.shape))
942        self._check_temperature(temperature)
943
944        i = self._temperature_index(temperature)
945        self._delayed_nu_fission[i] = delayed_nu_fission
946        self._set_fissionable(delayed_nu_fission)
947
948    def set_inverse_velocity(self, inv_vel, temperature=ROOM_TEMPERATURE_KELVIN):
949        """This method sets the inverse velocity for this XSdata object at the
950        provided temperature.
951
952        Parameters
953        ----------
954        inv_vel: np.ndarray
955            Inverse velocity in units of sec/cm.
956        temperature : float
957            Temperature (in Kelvin) of the data. Defaults to room temperature
958            (294K).
959
960        """
961
962        # Get the accepted shapes for this xs
963        shapes = [self.xs_shapes["[G]"]]
964
965        # Convert to a numpy array so we can easily get the shape for checking
966        inv_vel = np.asarray(inv_vel)
967        check_value('inverse_velocity shape', inv_vel.shape, shapes)
968        self._check_temperature(temperature)
969
970        i = self._temperature_index(temperature)
971        self._inverse_velocity[i] = inv_vel
972
973    def set_total_mgxs(self, total, temperature=ROOM_TEMPERATURE_KELVIN, nuclide='total',
974                       xs_type='macro', subdomain=None):
975        """This method allows for an openmc.mgxs.TotalXS or
976        openmc.mgxs.TransportXS to be used to set the total cross section for
977        this XSdata object.
978
979        Parameters
980        ----------
981        total: openmc.mgxs.TotalXS or openmc.mgxs.TransportXS
982            MGXS Object containing the total, transport or nu-transport cross
983            section for the domain of interest.
984        temperature : float
985            Temperature (in Kelvin) of the data. Defaults to room temperature
986            (294K).
987        nuclide : str
988            Individual nuclide (or 'total' if obtaining material-wise data)
989            to gather data for.  Defaults to 'total'.
990        xs_type: {'macro', 'micro'}
991            Provide the macro or micro cross section in units of cm^-1 or
992            barns. Defaults to 'macro'.
993        subdomain : iterable of int
994            If the MGXS contains a mesh domain type, the subdomain parameter
995            specifies which mesh cell (i.e., [i, j, k] index) to use.
996
997        See also
998        --------
999        openmc.mgxs.Library.create_mg_library()
1000        openmc.mgxs.Library.get_xsdata()
1001
1002        """
1003
1004        check_type('total', total, (openmc.mgxs.TotalXS,
1005                                    openmc.mgxs.TransportXS))
1006        check_value('energy_groups', total.energy_groups, [self.energy_groups])
1007        check_value('domain_type', total.domain_type, openmc.mgxs.DOMAIN_TYPES)
1008        self._check_temperature(temperature)
1009
1010        i = self._temperature_index(temperature)
1011        self._total[i] = total.get_xs(nuclides=nuclide, xs_type=xs_type,
1012                                      subdomains=subdomain)
1013
1014    def set_absorption_mgxs(self, absorption, temperature=ROOM_TEMPERATURE_KELVIN,
1015                            nuclide='total', xs_type='macro', subdomain=None):
1016        """This method allows for an openmc.mgxs.AbsorptionXS
1017        to be used to set the absorption cross section for this XSdata object.
1018
1019        Parameters
1020        ----------
1021        absorption: openmc.mgxs.AbsorptionXS
1022            MGXS Object containing the absorption cross section
1023            for the domain of interest.
1024        temperature : float
1025            Temperature (in Kelvin) of the data. Defaults to room temperature
1026            (294K).
1027        nuclide : str
1028            Individual nuclide (or 'total' if obtaining material-wise data)
1029            to gather data for.  Defaults to 'total'.
1030        xs_type: {'macro', 'micro'}
1031            Provide the macro or micro cross section in units of cm^-1 or
1032            barns. Defaults to 'macro'.
1033        subdomain : iterable of int
1034            If the MGXS contains a mesh domain type, the subdomain parameter
1035            specifies which mesh cell (i.e., [i, j, k] index) to use.
1036
1037        See also
1038        --------
1039        openmc.mgxs.Library.create_mg_library()
1040        openmc.mgxs.Library.get_xsdata()
1041
1042        """
1043
1044        check_type('absorption', absorption, openmc.mgxs.AbsorptionXS)
1045        check_value('energy_groups', absorption.energy_groups,
1046                    [self.energy_groups])
1047        check_value('domain_type', absorption.domain_type,
1048                    openmc.mgxs.DOMAIN_TYPES)
1049        self._check_temperature(temperature)
1050
1051        i = self._temperature_index(temperature)
1052        self._absorption[i] = absorption.get_xs(nuclides=nuclide,
1053                                                xs_type=xs_type,
1054                                                subdomains=subdomain)
1055
1056    def set_fission_mgxs(self, fission, temperature=ROOM_TEMPERATURE_KELVIN, nuclide='total',
1057                         xs_type='macro', subdomain=None):
1058        """This method allows for an openmc.mgxs.FissionXS
1059        to be used to set the fission cross section for this XSdata object.
1060
1061        Parameters
1062        ----------
1063        fission: openmc.mgxs.FissionXS
1064            MGXS Object containing the fission cross section
1065            for the domain of interest.
1066        temperature : float
1067            Temperature (in Kelvin) of the data. Defaults to room temperature
1068            (294K).
1069        nuclide : str
1070            Individual nuclide (or 'total' if obtaining material-wise data)
1071            to gather data for.  Defaults to 'total'.
1072        xs_type: {'macro', 'micro'}
1073            Provide the macro or micro cross section in units of cm^-1 or
1074            barns. Defaults to 'macro'.
1075        subdomain : iterable of int
1076            If the MGXS contains a mesh domain type, the subdomain parameter
1077            specifies which mesh cell (i.e., [i, j, k] index) to use.
1078
1079        See also
1080        --------
1081        openmc.mgxs.Library.create_mg_library()
1082        openmc.mgxs.Library.get_xsdata()
1083
1084        """
1085
1086        check_type('fission', fission, openmc.mgxs.FissionXS)
1087        check_value('energy_groups', fission.energy_groups,
1088                    [self.energy_groups])
1089        check_value('domain_type', fission.domain_type,
1090                    openmc.mgxs.DOMAIN_TYPES)
1091        self._check_temperature(temperature)
1092
1093        i = self._temperature_index(temperature)
1094        self._fission[i] = fission.get_xs(nuclides=nuclide,
1095                                          xs_type=xs_type,
1096                                          subdomains=subdomain)
1097
1098    def set_nu_fission_mgxs(self, nu_fission, temperature=ROOM_TEMPERATURE_KELVIN,
1099                            nuclide='total', xs_type='macro', subdomain=None):
1100        """This method allows for an openmc.mgxs.FissionXS
1101        to be used to set the nu-fission cross section for this XSdata object.
1102
1103        Parameters
1104        ----------
1105        nu_fission: openmc.mgxs.FissionXS
1106            MGXS Object containing the nu-fission cross section
1107            for the domain of interest.
1108        temperature : float
1109            Temperature (in Kelvin) of the data. Defaults to room temperature
1110            (294K).
1111        nuclide : str
1112            Individual nuclide (or 'total' if obtaining material-wise data)
1113            to gather data for.  Defaults to 'total'.
1114        xs_type: {'macro', 'micro'}
1115            Provide the macro or micro cross section in units of cm^-1 or
1116            barns. Defaults to 'macro'.
1117        subdomain : iterable of int
1118            If the MGXS contains a mesh domain type, the subdomain parameter
1119            specifies which mesh cell (i.e., [i, j, k] index) to use.
1120
1121        See also
1122        --------
1123        openmc.mgxs.Library.create_mg_library()
1124        openmc.mgxs.Library.get_xsdata()
1125
1126        """
1127
1128        check_type('nu_fission', nu_fission, (openmc.mgxs.FissionXS,
1129                                              openmc.mgxs.NuFissionMatrixXS))
1130        if isinstance(nu_fission, openmc.mgxs.FissionXS):
1131            check_value('nu', nu_fission.nu, [True])
1132        check_value('energy_groups', nu_fission.energy_groups,
1133                    [self.energy_groups])
1134        check_value('domain_type', nu_fission.domain_type,
1135                    openmc.mgxs.DOMAIN_TYPES)
1136        self._check_temperature(temperature)
1137
1138        i = self._temperature_index(temperature)
1139        self._nu_fission[i] = nu_fission.get_xs(nuclides=nuclide,
1140                                                xs_type=xs_type,
1141                                                subdomains=subdomain)
1142
1143        self._set_fissionable(self._nu_fission)
1144
1145    def set_prompt_nu_fission_mgxs(self, prompt_nu_fission, temperature=ROOM_TEMPERATURE_KELVIN,
1146                                   nuclide='total', xs_type='macro',
1147                                   subdomain=None):
1148        """Sets the prompt-nu-fission cross section.
1149
1150        This method allows for an openmc.mgxs.FissionXS or
1151        openmc.mgxs.NuFissionMatrixXS to be used to set the prompt-nu-fission
1152        cross section for this XSdata object.
1153
1154        Parameters
1155        ----------
1156        prompt_nu_fission: openmc.mgxs.FissionXS or openmc.mgxs.NuFissionMatrixXS
1157            MGXS Object containing the prompt-nu-fission cross section
1158            for the domain of interest.
1159        temperature : float
1160            Temperature (in units of Kelvin) of the provided dataset. Defaults
1161            to 294K
1162        nuclide : str
1163            Individual nuclide (or 'total' if obtaining material-wise data)
1164            to gather data for.  Defaults to 'total'.
1165        xs_type: {'macro', 'micro'}
1166            Provide the macro or micro cross section in units of cm^-1 or
1167            barns. Defaults to 'macro'.
1168        subdomain : iterable of int
1169            If the MGXS contains a mesh domain type, the subdomain parameter
1170            specifies which mesh cell (i.e., [i, j, k] index) to use.
1171
1172        See also
1173        --------
1174        openmc.mgxs.Library.create_mg_library()
1175        openmc.mgxs.Library.get_xsdata()
1176
1177        """
1178
1179        check_type('prompt_nu_fission', prompt_nu_fission,
1180                   (openmc.mgxs.FissionXS, openmc.mgxs.NuFissionMatrixXS))
1181        check_value('prompt', prompt_nu_fission.prompt, [True])
1182        check_value('energy_groups', prompt_nu_fission.energy_groups,
1183                    [self.energy_groups])
1184        check_value('domain_type', prompt_nu_fission.domain_type,
1185                    openmc.mgxs.DOMAIN_TYPES)
1186        self._check_temperature(temperature)
1187
1188        i = self._temperature_index(temperature)
1189        self._prompt_nu_fission[i] = prompt_nu_fission.get_xs(
1190            nuclides=nuclide, xs_type=xs_type, subdomains=subdomain)
1191
1192        self._set_fissionable(self._prompt_nu_fission)
1193
1194    def set_delayed_nu_fission_mgxs(self, delayed_nu_fission, temperature=ROOM_TEMPERATURE_KELVIN,
1195                                    nuclide='total', xs_type='macro',
1196                                    subdomain=None):
1197        """This method allows for an openmc.mgxs.DelayedNuFissionXS or
1198        openmc.mgxs.DelayedNuFissionMatrixXS to be used to set the
1199        delayed-nu-fission cross section for this XSdata object.
1200
1201        Parameters
1202        ----------
1203        delayed_nu_fission: openmc.mgxs.DelayedNuFissionXS or openmc.mgxs.DelayedNuFissionMatrixXS
1204            MGXS Object containing the delayed-nu-fission cross section
1205            for the domain of interest.
1206        temperature : float
1207            Temperature (in units of Kelvin) of the provided dataset. Defaults
1208            to 294K
1209        nuclide : str
1210            Individual nuclide (or 'total' if obtaining material-wise data)
1211            to gather data for.  Defaults to 'total'.
1212        xs_type: {'macro', 'micro'}
1213            Provide the macro or micro cross section in units of cm^-1 or
1214            barns. Defaults to 'macro'.
1215        subdomain : iterable of int
1216            If the MGXS contains a mesh domain type, the subdomain parameter
1217            specifies which mesh cell (i.e., [i, j, k] index) to use.
1218
1219        See also
1220        --------
1221        openmc.mgxs.Library.create_mg_library()
1222        openmc.mgxs.Library.get_xsdata()
1223
1224        """
1225
1226        check_type('delayed_nu_fission', delayed_nu_fission,
1227                   (openmc.mgxs.DelayedNuFissionXS,
1228                    openmc.mgxs.DelayedNuFissionMatrixXS))
1229        check_value('energy_groups', delayed_nu_fission.energy_groups,
1230                    [self.energy_groups])
1231        check_value('num_delayed_groups', delayed_nu_fission.num_delayed_groups,
1232                    [self.num_delayed_groups])
1233        check_value('domain_type', delayed_nu_fission.domain_type,
1234                    openmc.mgxs.DOMAIN_TYPES)
1235        self._check_temperature(temperature)
1236
1237        i = self._temperature_index(temperature)
1238        self._delayed_nu_fission[i] = delayed_nu_fission.get_xs(
1239            nuclides=nuclide, xs_type=xs_type, subdomains=subdomain)
1240
1241        self._set_fissionable(self._delayed_nu_fission)
1242
1243    def set_kappa_fission_mgxs(self, k_fission, temperature=ROOM_TEMPERATURE_KELVIN,
1244                               nuclide='total', xs_type='macro',
1245                               subdomain=None):
1246        """This method allows for an openmc.mgxs.KappaFissionXS
1247        to be used to set the kappa-fission cross section for this XSdata
1248        object.
1249
1250        Parameters
1251        ----------
1252        kappa_fission: openmc.mgxs.KappaFissionXS
1253            MGXS Object containing the kappa-fission cross section
1254            for the domain of interest.
1255        temperature : float
1256            Temperature (in Kelvin) of the data. Defaults to room temperature
1257            (294K).
1258        nuclide : str
1259            Individual nuclide (or 'total' if obtaining material-wise data)
1260            to gather data for.  Defaults to 'total'.
1261        xs_type: {'macro', 'micro'}
1262            Provide the macro or micro cross section in units of cm^-1 or
1263            barns. Defaults to 'macro'.
1264        subdomain : iterable of int
1265            If the MGXS contains a mesh domain type, the subdomain parameter
1266            specifies which mesh cell (i.e., [i, j, k] index) to use.
1267
1268        See also
1269        --------
1270        openmc.mgxs.Library.create_mg_library()
1271        openmc.mgxs.Library.get_xsdata()
1272
1273        """
1274
1275        check_type('kappa_fission', k_fission, openmc.mgxs.KappaFissionXS)
1276        check_value('energy_groups', k_fission.energy_groups,
1277                    [self.energy_groups])
1278        check_value('domain_type', k_fission.domain_type,
1279                    openmc.mgxs.DOMAIN_TYPES)
1280        self._check_temperature(temperature)
1281
1282        i = self._temperature_index(temperature)
1283        self._kappa_fission[i] = k_fission.get_xs(nuclides=nuclide,
1284                                                  xs_type=xs_type,
1285                                                  subdomains=subdomain)
1286
1287    def set_chi_mgxs(self, chi, temperature=ROOM_TEMPERATURE_KELVIN, nuclide='total',
1288                     xs_type='macro', subdomain=None):
1289        """This method allows for an openmc.mgxs.Chi
1290        to be used to set chi for this XSdata object.
1291
1292        Parameters
1293        ----------
1294        chi: openmc.mgxs.Chi
1295            MGXS Object containing chi for the domain of interest.
1296        temperature : float
1297            Temperature (in Kelvin) of the data. Defaults to room temperature
1298            (294K).
1299        nuclide : str
1300            Individual nuclide (or 'total' if obtaining material-wise data)
1301            to gather data for.  Defaults to 'total'.
1302        xs_type: {'macro', 'micro'}
1303            Provide the macro or micro cross section in units of cm^-1 or
1304            barns. Defaults to 'macro'.
1305        subdomain : iterable of int
1306            If the MGXS contains a mesh domain type, the subdomain parameter
1307            specifies which mesh cell (i.e., [i, j, k] index) to use.
1308
1309        See also
1310        --------
1311        openmc.mgxs.Library.create_mg_library()
1312        openmc.mgxs.Library.get_xsdata()
1313
1314        """
1315
1316        check_type('chi', chi, openmc.mgxs.Chi)
1317        check_value('energy_groups', chi.energy_groups, [self.energy_groups])
1318        check_value('domain_type', chi.domain_type, openmc.mgxs.DOMAIN_TYPES)
1319        self._check_temperature(temperature)
1320
1321        i = self._temperature_index(temperature)
1322        self._chi[i] = chi.get_xs(nuclides=nuclide, xs_type=xs_type,
1323                                  subdomains=subdomain)
1324
1325    def set_chi_prompt_mgxs(self, chi_prompt, temperature=ROOM_TEMPERATURE_KELVIN,
1326                            nuclide='total', xs_type='macro', subdomain=None):
1327        """This method allows for an openmc.mgxs.Chi to be used to set
1328        chi-prompt for this XSdata object.
1329
1330        Parameters
1331        ----------
1332        chi_prompt: openmc.mgxs.Chi
1333            MGXS Object containing chi-prompt for the domain of interest.
1334        temperature : float
1335            Temperature (in units of Kelvin) of the provided dataset. Defaults
1336            to 294K
1337        nuclide : str
1338            Individual nuclide (or 'total' if obtaining material-wise data)
1339            to gather data for.  Defaults to 'total'.
1340        xs_type: {'macro', 'micro'}
1341            Provide the macro or micro cross section in units of cm^-1 or
1342            barns. Defaults to 'macro'.
1343        subdomain : iterable of int
1344            If the MGXS contains a mesh domain type, the subdomain parameter
1345            specifies which mesh cell (i.e., [i, j, k] index) to use.
1346
1347        See also
1348        --------
1349        openmc.mgxs.Library.create_mg_library()
1350        openmc.mgxs.Library.get_xsdata()
1351
1352        """
1353
1354        check_type('chi_prompt', chi_prompt, openmc.mgxs.Chi)
1355        check_value('prompt', chi_prompt.prompt, [True])
1356        check_value('energy_groups', chi_prompt.energy_groups,
1357                    [self.energy_groups])
1358        check_value('domain_type', chi_prompt.domain_type,
1359                    openmc.mgxs.DOMAIN_TYPES)
1360        self._check_temperature(temperature)
1361
1362        i = self._temperature_index(temperature)
1363        self._chi_prompt[i] = chi_prompt.get_xs(nuclides=nuclide,
1364                                                xs_type=xs_type,
1365                                                subdomains=subdomain)
1366
1367    def set_chi_delayed_mgxs(self, chi_delayed, temperature=ROOM_TEMPERATURE_KELVIN,
1368                             nuclide='total', xs_type='macro', subdomain=None):
1369        """This method allows for an openmc.mgxs.ChiDelayed
1370        to be used to set chi-delayed for this XSdata object.
1371
1372        Parameters
1373        ----------
1374        chi_delayed: openmc.mgxs.ChiDelayed
1375            MGXS Object containing chi-delayed for the domain of interest.
1376        temperature : float
1377            Temperature (in units of Kelvin) of the provided dataset. Defaults
1378            to 294K
1379        nuclide : str
1380            Individual nuclide (or 'total' if obtaining material-wise data)
1381            to gather data for.  Defaults to 'total'.
1382        xs_type: {'macro', 'micro'}
1383            Provide the macro or micro cross section in units of cm^-1 or
1384            barns. Defaults to 'macro'.
1385        subdomain : iterable of int
1386            If the MGXS contains a mesh domain type, the subdomain parameter
1387            specifies which mesh cell (i.e., [i, j, k] index) to use.
1388
1389        See also
1390        --------
1391        openmc.mgxs.Library.create_mg_library()
1392        openmc.mgxs.Library.get_xsdata()
1393
1394        """
1395
1396        check_type('chi_delayed', chi_delayed, openmc.mgxs.ChiDelayed)
1397        check_value('energy_groups', chi_delayed.energy_groups,
1398                    [self.energy_groups])
1399        check_value('num_delayed_groups', chi_delayed.num_delayed_groups,
1400                    [self.num_delayed_groups])
1401        check_value('domain_type', chi_delayed.domain_type,
1402                    openmc.mgxs.DOMAIN_TYPES)
1403        self._check_temperature(temperature)
1404
1405        i = self._temperature_index(temperature)
1406        self._chi_delayed[i] = chi_delayed.get_xs(nuclides=nuclide,
1407                                                  xs_type=xs_type,
1408                                                  subdomains=subdomain)
1409
1410    def set_beta_mgxs(self, beta, temperature=ROOM_TEMPERATURE_KELVIN,
1411                      nuclide='total', xs_type='macro', subdomain=None):
1412        """This method allows for an openmc.mgxs.Beta
1413        to be used to set beta for this XSdata object.
1414
1415        Parameters
1416        ----------
1417        beta : openmc.mgxs.Beta
1418            MGXS Object containing beta for the domain of interest.
1419        temperature : float
1420            Temperature (in units of Kelvin) of the provided dataset. Defaults
1421            to 294K
1422        nuclide : str
1423            Individual nuclide (or 'total' if obtaining material-wise data)
1424            to gather data for.  Defaults to 'total'.
1425        xs_type: {'macro', 'micro'}
1426            Provide the macro or micro cross section in units of cm^-1 or
1427            barns. Defaults to 'macro'.
1428        subdomain : iterable of int
1429            If the MGXS contains a mesh domain type, the subdomain parameter
1430            specifies which mesh cell (i.e., [i, j, k] index) to use.
1431
1432        See also
1433        --------
1434        openmc.mgxs.Library.create_mg_library()
1435        openmc.mgxs.Library.get_xsdata()
1436
1437        """
1438
1439        check_type('beta', beta, openmc.mgxs.Beta)
1440        check_value('num_delayed_groups', beta.num_delayed_groups,
1441                    [self.num_delayed_groups])
1442        check_value('domain_type', beta.domain_type, openmc.mgxs.DOMAIN_TYPES)
1443        self._check_temperature(temperature)
1444
1445        i = self._temperature_index(temperature)
1446        self._beta[i] = beta.get_xs(nuclides=nuclide,
1447                                    xs_type=xs_type,
1448                                    subdomains=subdomain)
1449
1450    def set_decay_rate_mgxs(self, decay_rate, temperature=ROOM_TEMPERATURE_KELVIN,
1451                            nuclide='total', xs_type='macro', subdomain=None):
1452        """This method allows for an openmc.mgxs.DecayRate
1453        to be used to set decay rate for this XSdata object.
1454
1455        Parameters
1456        ----------
1457        decay_rate : openmc.mgxs.DecayRate
1458            MGXS Object containing decay rate for the domain of interest.
1459        temperature : float
1460            Temperature (in units of Kelvin) of the provided dataset. Defaults
1461            to 294K
1462        nuclide : str
1463            Individual nuclide (or 'total' if obtaining material-wise data)
1464            to gather data for.  Defaults to 'total'.
1465        xs_type: {'macro', 'micro'}
1466            Provide the macro or micro cross section in units of cm^-1 or
1467            barns. Defaults to 'macro'.
1468        subdomain : iterable of int
1469            If the MGXS contains a mesh domain type, the subdomain parameter
1470            specifies which mesh cell (i.e., [i, j, k] index) to use.
1471
1472        See also
1473        --------
1474        openmc.mgxs.Library.create_mg_library()
1475        openmc.mgxs.Library.get_xsdata()
1476
1477        """
1478
1479        check_type('decay_rate', decay_rate, openmc.mgxs.DecayRate)
1480        check_value('num_delayed_groups', decay_rate.num_delayed_groups,
1481                    [self.num_delayed_groups])
1482        check_value('domain_type', decay_rate.domain_type,
1483                    openmc.mgxs.DOMAIN_TYPES)
1484        self._check_temperature(temperature)
1485
1486        i = self._temperature_index(temperature)
1487        self._decay_rate[i] = decay_rate.get_xs(nuclides=nuclide,
1488                                                xs_type=xs_type,
1489                                                subdomains=subdomain)
1490
1491    def set_scatter_matrix_mgxs(self, scatter, temperature=ROOM_TEMPERATURE_KELVIN,
1492                                nuclide='total', xs_type='macro',
1493                                subdomain=None):
1494        """This method allows for an openmc.mgxs.ScatterMatrixXS
1495        to be used to set the scatter matrix cross section for this XSdata
1496        object.  If the XSdata.order attribute has not yet been set, then
1497        it will be set based on the properties of scatter.
1498
1499        Parameters
1500        ----------
1501        scatter: openmc.mgxs.ScatterMatrixXS
1502            MGXS Object containing the scatter matrix cross section
1503            for the domain of interest.
1504        temperature : float
1505            Temperature (in Kelvin) of the data. Defaults to room temperature
1506            (294K).
1507        nuclide : str
1508            Individual nuclide (or 'total' if obtaining material-wise data)
1509            to gather data for.  Defaults to 'total'.
1510        xs_type: {'macro', 'micro'}
1511            Provide the macro or micro cross section in units of cm^-1 or
1512            barns. Defaults to 'macro'.
1513        subdomain : iterable of int
1514            If the MGXS contains a mesh domain type, the subdomain parameter
1515            specifies which mesh cell (i.e., [i, j, k] index) to use.
1516
1517        See also
1518        --------
1519        openmc.mgxs.Library.create_mg_library()
1520        openmc.mgxs.Library.get_xsdata()
1521
1522        """
1523
1524        check_type('scatter', scatter, openmc.mgxs.ScatterMatrixXS)
1525        check_value('energy_groups', scatter.energy_groups,
1526                    [self.energy_groups])
1527        check_value('domain_type', scatter.domain_type,
1528                    openmc.mgxs.DOMAIN_TYPES)
1529        self._check_temperature(temperature)
1530
1531        # Set the value of scatter_format based on the same value within
1532        # scatter
1533        self.scatter_format = scatter.scatter_format
1534
1535        # If the user has not defined XSdata.order, then we will set
1536        # the order based on the data within scatter.
1537        # Otherwise, we will check to see that XSdata.order matches
1538        # the order of scatter
1539        if self.scatter_format == SCATTER_LEGENDRE:
1540            if self.order is None:
1541                self.order = scatter.legendre_order
1542            else:
1543                check_value('legendre_order', scatter.legendre_order,
1544                            [self.order])
1545        elif self.scatter_format == SCATTER_HISTOGRAM:
1546            if self.order is None:
1547                self.order = scatter.histogram_bins
1548            else:
1549                check_value('histogram_bins', scatter.histogram_bins,
1550                            [self.order])
1551
1552        i = self._temperature_index(temperature)
1553        if self.scatter_format == SCATTER_LEGENDRE:
1554            self._scatter_matrix[i] = \
1555                np.zeros(self.xs_shapes["[G][G'][Order]"])
1556            # Get the scattering orders in the outermost dimension
1557            if self.representation == REPRESENTATION_ISOTROPIC:
1558                for moment in range(self.num_orders):
1559                    self._scatter_matrix[i][:, :, moment] = \
1560                        scatter.get_xs(nuclides=nuclide, xs_type=xs_type,
1561                                       moment=moment, subdomains=subdomain)
1562            elif self.representation == REPRESENTATION_ANGLE:
1563                for moment in range(self.num_orders):
1564                    self._scatter_matrix[i][:, :, :, :, moment] = \
1565                        scatter.get_xs(nuclides=nuclide, xs_type=xs_type,
1566                                       moment=moment, subdomains=subdomain)
1567        else:
1568            self._scatter_matrix[i] = \
1569                scatter.get_xs(nuclides=nuclide, xs_type=xs_type,
1570                               subdomains=subdomain)
1571
1572    def set_multiplicity_matrix_mgxs(self, nuscatter, scatter=None,
1573                                     temperature=ROOM_TEMPERATURE_KELVIN, nuclide='total',
1574                                     xs_type='macro', subdomain=None):
1575        """This method allows for either the direct use of only an
1576        openmc.mgxs.MultiplicityMatrixXS or an openmc.mgxs.ScatterMatrixXS and
1577        openmc.mgxs.ScatterMatrixXS to be used to set the scattering
1578        multiplicity for this XSdata object. Multiplicity, in OpenMC parlance,
1579        is a factor used to account for the production of neutrons introduced by
1580        scattering multiplication reactions, i.e., (n,xn) events. In this sense,
1581        the multiplication matrix is simply defined as the ratio of the
1582        nu-scatter and scatter matrices.
1583
1584        Parameters
1585        ----------
1586        nuscatter: openmc.mgxs.ScatterMatrixXS or openmc.mgxs.MultiplicityMatrixXS
1587            MGXS Object containing the matrix cross section for the domain
1588            of interest.
1589        scatter: openmc.mgxs.ScatterMatrixXS
1590            MGXS Object containing the scattering matrix cross section
1591            for the domain of interest.
1592        temperature : float
1593            Temperature (in Kelvin) of the data. Defaults to room temperature
1594            (294K).
1595        nuclide : str
1596            Individual nuclide (or 'total' if obtaining material-wise data)
1597            to gather data for.  Defaults to 'total'.
1598        xs_type: {'macro', 'micro'}
1599            Provide the macro or micro cross section in units of cm^-1 or
1600            barns. Defaults to 'macro'.
1601        subdomain : iterable of int
1602            If the MGXS contains a mesh domain type, the subdomain parameter
1603            specifies which mesh cell (i.e., [i, j, k] index) to use.
1604
1605        See also
1606        --------
1607        openmc.mgxs.Library.create_mg_library()
1608        openmc.mgxs.Library.get_xsdata()
1609
1610        """
1611
1612        check_type('nuscatter', nuscatter, (openmc.mgxs.ScatterMatrixXS,
1613                                            openmc.mgxs.MultiplicityMatrixXS))
1614        check_value('energy_groups', nuscatter.energy_groups,
1615                    [self.energy_groups])
1616        check_value('domain_type', nuscatter.domain_type,
1617                    openmc.mgxs.DOMAIN_TYPES)
1618        self._check_temperature(temperature)
1619
1620        if scatter is not None:
1621            check_type('scatter', scatter, openmc.mgxs.ScatterMatrixXS)
1622            if isinstance(nuscatter, openmc.mgxs.MultiplicityMatrixXS):
1623                msg = 'Either an MultiplicityMatrixXS object must be passed ' \
1624                      'for "nuscatter" or the "scatter" argument must be ' \
1625                      'provided.'
1626                raise ValueError(msg)
1627            check_value('energy_groups', scatter.energy_groups,
1628                        [self.energy_groups])
1629            check_value('domain_type', scatter.domain_type,
1630                        openmc.mgxs.DOMAIN_TYPES)
1631        i = self._temperature_index(temperature)
1632        nuscatt = nuscatter.get_xs(nuclides=nuclide,
1633                                   xs_type=xs_type, moment=0,
1634                                   subdomains=subdomain)
1635        if isinstance(nuscatter, openmc.mgxs.MultiplicityMatrixXS):
1636            self._multiplicity_matrix[i] = nuscatt
1637        else:
1638            scatt = scatter.get_xs(nuclides=nuclide,
1639                                   xs_type=xs_type, moment=0,
1640                                   subdomains=subdomain)
1641            if scatter.scatter_format == SCATTER_HISTOGRAM:
1642                scatt = np.sum(scatt, axis=2)
1643            if nuscatter.scatter_format == SCATTER_HISTOGRAM:
1644                nuscatt = np.sum(nuscatt, axis=2)
1645            self._multiplicity_matrix[i] = np.divide(nuscatt, scatt)
1646
1647        self._multiplicity_matrix[i] = \
1648            np.nan_to_num(self._multiplicity_matrix[i])
1649
1650    def set_inverse_velocity_mgxs(self, inverse_velocity, temperature=ROOM_TEMPERATURE_KELVIN,
1651                                  nuclide='total', xs_type='macro',
1652                                  subdomain=None):
1653        """This method allows for an openmc.mgxs.InverseVelocity
1654        to be used to set the inverse velocity for this XSdata object.
1655
1656        Parameters
1657        ----------
1658        inverse_velocity : openmc.mgxs.InverseVelocity
1659            MGXS object containing the inverse velocity for the domain of
1660            interest.
1661        temperature : float
1662            Temperature (in Kelvin) of the data. Defaults to room temperature
1663            (294K).
1664        nuclide : str
1665            Individual nuclide (or 'total' if obtaining material-wise data)
1666            to gather data for.  Defaults to 'total'.
1667        xs_type: {'macro', 'micro'}
1668            Provide the macro or micro cross section in units of cm^-1 or
1669            barns. Defaults to 'macro'.
1670        subdomain : iterable of int
1671            If the MGXS contains a mesh domain type, the subdomain parameter
1672            specifies which mesh cell (i.e., [i, j, k] index) to use.
1673
1674        See also
1675        --------
1676        openmc.mgxs.Library.create_mg_library()
1677        openmc.mgxs.Library.get_xsdata()
1678
1679        """
1680
1681        check_type('inverse_velocity', inverse_velocity,
1682                   openmc.mgxs.InverseVelocity)
1683        check_value('energy_groups', inverse_velocity.energy_groups,
1684                    [self.energy_groups])
1685        check_value('domain_type', inverse_velocity.domain_type,
1686                    openmc.mgxs.DOMAIN_TYPES)
1687        self._check_temperature(temperature)
1688
1689        i = self._temperature_index(temperature)
1690        self._inverse_velocity[i] = inverse_velocity.get_xs(
1691            nuclides=nuclide, xs_type=xs_type, subdomains=subdomain)
1692
1693    def convert_representation(self, target_representation, num_polar=None,
1694                               num_azimuthal=None):
1695        """Produce a new XSdata object with the same data, but converted to the
1696        new representation (isotropic or angle-dependent).
1697
1698        This method cannot be used to change the number of polar or
1699        azimuthal bins of an XSdata object that already uses an angular
1700        representation. Finally, this method simply uses an arithmetic mean to
1701        convert from an angular to isotropic representation; no flux-weighting
1702        is applied and therefore reaction rates will not be preserved.
1703
1704        Parameters
1705        ----------
1706        target_representation : {'isotropic', 'angle'}
1707            Representation of the MGXS (isotropic or angle-dependent flux
1708            weighting).
1709        num_polar : int, optional
1710            Number of equal width angular bins that the polar angular domain is
1711            subdivided into. This is required when `target_representation` is
1712            "angle".
1713
1714        num_azimuthal : int, optional
1715            Number of equal width angular bins that the azimuthal angular domain
1716            is subdivided into. This is required when `target_representation` is
1717            "angle".
1718
1719        Returns
1720        -------
1721        openmc.XSdata
1722
1723            Multi-group cross section data with the same data as self, but
1724            represented as specified in `target_representation`.
1725
1726        """
1727
1728        check_value('target_representation', target_representation,
1729                    _REPRESENTATIONS)
1730        if target_representation == REPRESENTATION_ANGLE:
1731            check_type('num_polar', num_polar, Integral)
1732            check_type('num_azimuthal', num_azimuthal, Integral)
1733            check_greater_than('num_polar', num_polar, 0)
1734            check_greater_than('num_azimuthal', num_azimuthal, 0)
1735
1736        xsdata = copy.deepcopy(self)
1737
1738        # First handle the case where the current and requested
1739        # representations are the same
1740        if target_representation == self.representation:
1741            # Check to make sure the num_polar and num_azimuthal values match
1742            if target_representation == REPRESENTATION_ANGLE:
1743                if num_polar != self.num_polar or num_azimuthal != self.num_azimuthal:
1744                    raise ValueError("Cannot translate between `angle`"
1745                                     " representations with different angle"
1746                                     " bin structures")
1747            # Nothing to do as the same structure was requested
1748            return xsdata
1749
1750        xsdata.representation = target_representation
1751        # We have different actions depending on the representation conversion
1752        if target_representation == REPRESENTATION_ISOTROPIC:
1753            # This is not needed for the correct functionality, but these
1754            # values are changed back to None for clarity
1755            xsdata._num_polar = None
1756            xsdata._num_azimuthal = None
1757
1758        elif target_representation == REPRESENTATION_ANGLE:
1759            xsdata.num_polar = num_polar
1760            xsdata.num_azimuthal = num_azimuthal
1761
1762        # Reset xs_shapes so it is recalculated the next time it is needed
1763        xsdata._xs_shapes = None
1764
1765        for i, temp in enumerate(xsdata.temperatures):
1766            for xs in ['total', 'absorption', 'fission', 'nu_fission',
1767                       'scatter_matrix', 'multiplicity_matrix',
1768                       'prompt_nu_fission', 'delayed_nu_fission',
1769                       'kappa_fission', 'chi', 'chi_prompt', 'chi_delayed',
1770                       'beta', 'decay_rate', 'inverse_velocity']:
1771                # Get the original data
1772                orig_data = getattr(self, '_' + xs)[i]
1773                if orig_data is not None:
1774
1775                    if target_representation == 'isotropic':
1776                        # Since we are going from angle to isotropic, the
1777                        # current data is just the average over the angle bins
1778                        new_data = orig_data.mean(axis=(0, 1))
1779
1780                    elif target_representation == REPRESENTATION_ANGLE:
1781                        # Since we are going from isotropic to angle, the
1782                        # current data is just copied for every angle bin
1783                        new_shape = (num_polar, num_azimuthal) + \
1784                            orig_data.shape
1785                        new_data = np.resize(orig_data, new_shape)
1786
1787                    setter = getattr(xsdata, 'set_' + xs)
1788                    setter(new_data, temp)
1789
1790        return xsdata
1791
1792    def convert_scatter_format(self, target_format, target_order=None):
1793        """Produce a new MGXSLibrary object with the same data, but converted
1794        to the new scatter format and order
1795
1796        Parameters
1797        ----------
1798        target_format : {'tabular', 'legendre', 'histogram'}
1799            Representation of the scattering angle distribution
1800        target_order : int
1801            Either the Legendre target_order, number of bins, or number of
1802            points used to describe the angular distribution associated with
1803            each group-to-group transfer probability
1804
1805        Returns
1806        -------
1807        openmc.XSdata
1808            Multi-group cross section data with the same data as in self, but
1809            represented as specified in `target_format`.
1810
1811        """
1812
1813        check_value('target_format', target_format, _SCATTER_TYPES)
1814        check_type('target_order', target_order, Integral)
1815        if target_format == SCATTER_LEGENDRE:
1816            check_greater_than('target_order', target_order, 0, equality=True)
1817        else:
1818            check_greater_than('target_order', target_order, 0)
1819
1820        xsdata = copy.deepcopy(self)
1821        xsdata.scatter_format = target_format
1822        xsdata.order = target_order
1823
1824        # Reset and re-generate XSdata.xs_shapes with the new scattering format
1825        xsdata._xs_shapes = None
1826
1827        for i, temp in enumerate(xsdata.temperatures):
1828            orig_data = self._scatter_matrix[i]
1829            new_shape = orig_data.shape[:-1] + (xsdata.num_orders,)
1830            new_data = np.zeros(new_shape)
1831
1832            if self.scatter_format == SCATTER_LEGENDRE:
1833                if target_format == SCATTER_LEGENDRE:
1834                    # Then we are changing orders and only need to change
1835                    # dimensionality of the mu data and pad/truncate as needed
1836                    order = min(xsdata.num_orders, self.num_orders)
1837                    new_data[..., :order] = orig_data[..., :order]
1838
1839                elif target_format == SCATTER_TABULAR:
1840                    mu = np.linspace(-1, 1, xsdata.num_orders)
1841                    # Evaluate the legendre on the mu grid
1842                    for imu in range(len(mu)):
1843                        for l in range(self.num_orders):
1844                            new_data[..., imu] += (
1845                                 (l + 0.5) * eval_legendre(l, mu[imu]) *
1846                                 orig_data[..., l])
1847
1848                elif target_format == SCATTER_HISTOGRAM:
1849                    # This code uses the vectorized integration capabilities
1850                    # instead of having an isotropic and angle representation
1851                    # path.
1852                    # Set the histogram mu grid
1853                    mu = np.linspace(-1, 1, xsdata.num_orders + 1)
1854                    # For every bin perform simpson integration of a finely
1855                    # sampled orig_data
1856                    for h_bin in range(xsdata.num_orders):
1857                        mu_fine = np.linspace(mu[h_bin], mu[h_bin + 1], _NMU)
1858                        table_fine = np.zeros(new_data.shape[:-1] + (_NMU,))
1859                        for imu in range(len(mu_fine)):
1860                            for l in range(self.num_orders):
1861                                table_fine[..., imu] += ((l + 0.5)
1862                                     * eval_legendre(l, mu_fine[imu]) *
1863                                     orig_data[..., l])
1864                        new_data[..., h_bin] = simps(table_fine, mu_fine)
1865
1866            elif self.scatter_format == SCATTER_TABULAR:
1867                # Calculate the mu points of the current data
1868                mu_self = np.linspace(-1, 1, self.num_orders)
1869
1870                if target_format == SCATTER_LEGENDRE:
1871                    # Find the Legendre coefficients via integration. To best
1872                    # use the vectorized integration capabilities of scipy,
1873                    # this is done with fixed sample integration routines.
1874                    mu_fine = np.linspace(-1, 1, _NMU)
1875                    for l in range(xsdata.num_orders):
1876                        y = (interp1d(mu_self, orig_data)(mu_fine) *
1877                             eval_legendre(l, mu_fine))
1878                        new_data[..., l] = simps(y, mu_fine)
1879
1880                elif target_format == SCATTER_TABULAR:
1881                    # Simply use an interpolating function to get the new data
1882                    mu = np.linspace(-1, 1, xsdata.num_orders)
1883                    new_data[..., :] = interp1d(mu_self, orig_data)(mu)
1884
1885                elif target_format == SCATTER_HISTOGRAM:
1886                    # Use an interpolating function to do the bin-wise
1887                    # integrals
1888                    mu = np.linspace(-1, 1, xsdata.num_orders + 1)
1889
1890                    # Like the tabular -> legendre path above, this code will
1891                    # be written to utilize the vectorized integration
1892                    # capabilities instead of having an isotropic and
1893                    # angle representation path.
1894                    interp = interp1d(mu_self, orig_data)
1895                    for h_bin in range(xsdata.num_orders):
1896                        mu_fine = np.linspace(mu[h_bin], mu[h_bin + 1], _NMU)
1897                        new_data[..., h_bin] = simps(interp(mu_fine), mu_fine)
1898
1899            elif self.scatter_format == SCATTER_HISTOGRAM:
1900                # The histogram format does not have enough information to
1901                # convert to the other forms without inducing some amount of
1902                # error. We will make the assumption that the center of the bin
1903                # has the value of the bin. The mu=-1 and 1 points will be
1904                # extrapolated from the shape.
1905                mu_midpoint = np.linspace(-1, 1, self.num_orders,
1906                                          endpoint=False)
1907                mu_midpoint += (mu_midpoint[1] - mu_midpoint[0]) * 0.5
1908                interp = interp1d(mu_midpoint, orig_data,
1909                                  fill_value='extrapolate')
1910                # Now get the distribution normalization factor to take from
1911                # an integral quantity to a point-wise quantity
1912                norm = float(self.num_orders) / 2.0
1913
1914                # We now have a tabular distribution in tab_data on mu_self.
1915                # We now proceed just like the tabular branch above.
1916                if target_format == SCATTER_LEGENDRE:
1917                    # find the legendre coefficients via integration. To best
1918                    # use the vectorized integration capabilities of scipy,
1919                    # this will be done with fixed sample integration routines.
1920                    mu_fine = np.linspace(-1, 1, _NMU)
1921                    for l in range(xsdata.num_orders):
1922                        y = interp(mu_fine) * norm * eval_legendre(l, mu_fine)
1923                        new_data[..., l] = simps(y, mu_fine)
1924
1925                elif target_format == SCATTER_TABULAR:
1926                    # Simply use an interpolating function to get the new data
1927                    mu = np.linspace(-1, 1, xsdata.num_orders)
1928                    new_data[..., :] = interp(mu) * norm
1929
1930                elif target_format == SCATTER_HISTOGRAM:
1931                    # Use an interpolating function to do the bin-wise
1932                    # integrals
1933                    mu = np.linspace(-1, 1, xsdata.num_orders + 1)
1934
1935                    # Like the tabular -> legendre path above, this code will
1936                    # be written to utilize the vectorized integration
1937                    # capabilities instead of having an isotropic and
1938                    # angle representation path.
1939                    for h_bin in range(xsdata.num_orders):
1940                        mu_fine = np.linspace(mu[h_bin], mu[h_bin + 1], _NMU)
1941                        new_data[..., h_bin] = \
1942                            norm * simps(interp(mu_fine), mu_fine)
1943
1944            # Remove small values resulting from numerical precision issues
1945            new_data[..., np.abs(new_data) < 1.E-10] = 0.
1946
1947            xsdata.set_scatter_matrix(new_data, temp)
1948
1949        return xsdata
1950
1951    def to_hdf5(self, file):
1952        """Write XSdata to an HDF5 file
1953
1954        Parameters
1955        ----------
1956        file : h5py.File
1957            HDF5 File (a root Group) to write to
1958
1959        """
1960
1961        grp = file.create_group(self.name)
1962        if self.atomic_weight_ratio is not None:
1963            grp.attrs['atomic_weight_ratio'] = self.atomic_weight_ratio
1964        if self.fissionable is not None:
1965            grp.attrs['fissionable'] = self.fissionable
1966
1967        if self.representation is not None:
1968            grp.attrs['representation'] = np.string_(self.representation)
1969            if self.representation == REPRESENTATION_ANGLE:
1970                if self.num_azimuthal is not None:
1971                    grp.attrs['num_azimuthal'] = self.num_azimuthal
1972
1973                if self.num_polar is not None:
1974                    grp.attrs['num_polar'] = self.num_polar
1975
1976        grp.attrs['scatter_shape'] = np.string_("[G][G'][Order]")
1977        if self.scatter_format is not None:
1978            grp.attrs['scatter_format'] = np.string_(self.scatter_format)
1979        if self.order is not None:
1980            grp.attrs['order'] = self.order
1981
1982        ktg = grp.create_group('kTs')
1983        for temperature in self.temperatures:
1984            temp_label = str(int(np.round(temperature))) + "K"
1985            kT = temperature * openmc.data.K_BOLTZMANN
1986            ktg.create_dataset(temp_label, data=kT)
1987
1988        # Create the temperature datasets
1989        for i, temperature in enumerate(self.temperatures):
1990
1991            xs_grp = grp.create_group(str(int(np.round(temperature))) + "K")
1992
1993            if self._total[i] is None:
1994                raise ValueError('total data must be provided when writing '
1995                                 'the HDF5 library')
1996
1997            xs_grp.create_dataset("total", data=self._total[i])
1998
1999            if self._absorption[i] is None:
2000                raise ValueError('absorption data must be provided when '
2001                                 'writing the HDF5 library')
2002
2003            xs_grp.create_dataset("absorption", data=self._absorption[i])
2004
2005            if self.fissionable:
2006                if self._fission[i] is not None:
2007                    xs_grp.create_dataset("fission", data=self._fission[i])
2008
2009                if self._kappa_fission[i] is not None:
2010                    xs_grp.create_dataset("kappa-fission",
2011                                          data=self._kappa_fission[i])
2012
2013                if self._chi[i] is not None:
2014                    xs_grp.create_dataset("chi", data=self._chi[i])
2015
2016                if self._chi_prompt[i] is not None:
2017                    xs_grp.create_dataset("chi-prompt",
2018                                          data=self._chi_prompt[i])
2019
2020                if self._chi_delayed[i] is not None:
2021                    xs_grp.create_dataset("chi-delayed",
2022                                          data=self._chi_delayed[i])
2023
2024                if self._nu_fission[i] is None and \
2025                   (self._delayed_nu_fission[i] is None or \
2026                    self._prompt_nu_fission[i] is None):
2027                    raise ValueError('nu-fission or prompt-nu-fission and '
2028                                     'delayed-nu-fission data must be '
2029                                     'provided when writing the HDF5 library')
2030
2031                if self._nu_fission[i] is not None:
2032                    xs_grp.create_dataset("nu-fission",
2033                                          data=self._nu_fission[i])
2034
2035                if self._prompt_nu_fission[i] is not None:
2036                    xs_grp.create_dataset("prompt-nu-fission",
2037                                          data=self._prompt_nu_fission[i])
2038
2039                if self._delayed_nu_fission[i] is not None:
2040                    xs_grp.create_dataset("delayed-nu-fission",
2041                                          data=self._delayed_nu_fission[i])
2042
2043                if self._beta[i] is not None:
2044                    xs_grp.create_dataset("beta", data=self._beta[i])
2045
2046                if self._decay_rate[i] is not None:
2047                    xs_grp.create_dataset("decay-rate",
2048                                          data=self._decay_rate[i])
2049
2050            if self._scatter_matrix[i] is None:
2051                raise ValueError('Scatter matrix must be provided when '
2052                                 'writing the HDF5 library')
2053
2054            # Get the sparse scattering data to print to the library
2055            G = self.energy_groups.num_groups
2056            if self.representation == REPRESENTATION_ISOTROPIC:
2057                Np = 1
2058                Na = 1
2059            elif self.representation == REPRESENTATION_ANGLE:
2060                Np = self.num_polar
2061                Na = self.num_azimuthal
2062
2063            g_out_bounds = np.zeros((Np, Na, G, 2), dtype=int)
2064            for p in range(Np):
2065                for a in range(Na):
2066                    for g_in in range(G):
2067                        if self.scatter_format == SCATTER_LEGENDRE:
2068                            if self.representation == REPRESENTATION_ISOTROPIC:
2069                                matrix = \
2070                                    self._scatter_matrix[i][g_in, :, 0]
2071                            elif self.representation == REPRESENTATION_ANGLE:
2072                                matrix = \
2073                                    self._scatter_matrix[i][p, a, g_in, :, 0]
2074                        else:
2075                            if self.representation == REPRESENTATION_ISOTROPIC:
2076                                matrix = \
2077                                    np.sum(self._scatter_matrix[i][g_in, :, :],
2078                                           axis=1)
2079                            elif self.representation == REPRESENTATION_ANGLE:
2080                                matrix = \
2081                                    np.sum(self._scatter_matrix[i][p, a, g_in, :, :],
2082                                           axis=1)
2083                        nz = np.nonzero(matrix)
2084                        # It is possible that there only zeros in matrix
2085                        # and therefore nz will be empty, in that case set
2086                        # g_out_bounds to 0s
2087                        if len(nz[0]) == 0:
2088                            g_out_bounds[p, a, g_in, :] = 0
2089                        else:
2090                            g_out_bounds[p, a, g_in, 0] = nz[0][0]
2091                            g_out_bounds[p, a, g_in, 1] = nz[0][-1]
2092
2093            # Now create the flattened scatter matrix array
2094            flat_scatt = []
2095            for p in range(Np):
2096                for a in range(Na):
2097                    if self.representation == REPRESENTATION_ISOTROPIC:
2098                        matrix = self._scatter_matrix[i][:, :, :]
2099                    elif self.representation == REPRESENTATION_ANGLE:
2100                        matrix = self._scatter_matrix[i][p, a, :, :, :]
2101                    for g_in in range(G):
2102                        for g_out in range(g_out_bounds[p, a, g_in, 0],
2103                                           g_out_bounds[p, a, g_in, 1] + 1):
2104                            for l in range(len(matrix[g_in, g_out, :])):
2105                                flat_scatt.append(matrix[g_in, g_out, l])
2106
2107            # And write it.
2108            scatt_grp = xs_grp.create_group('scatter_data')
2109            scatt_grp.create_dataset("scatter_matrix",
2110                                     data=np.array(flat_scatt))
2111
2112            # Repeat for multiplicity
2113            if self._multiplicity_matrix[i] is not None:
2114
2115                # Now create the flattened scatter matrix array
2116                flat_mult = []
2117                for p in range(Np):
2118                    for a in range(Na):
2119                        if self.representation == REPRESENTATION_ISOTROPIC:
2120                            matrix = self._multiplicity_matrix[i][:, :]
2121                        elif self.representation == REPRESENTATION_ANGLE:
2122                            matrix = self._multiplicity_matrix[i][p, a, :, :]
2123                        for g_in in range(G):
2124                            for g_out in range(g_out_bounds[p, a, g_in, 0],
2125                                               g_out_bounds[p, a, g_in, 1] + 1):
2126                                flat_mult.append(matrix[g_in, g_out])
2127
2128                # And write it.
2129                scatt_grp.create_dataset("multiplicity_matrix",
2130                                         data=np.array(flat_mult))
2131
2132            # And finally, adjust g_out_bounds for 1-based group counting
2133            # and write it.
2134            g_out_bounds[:, :, :, :] += 1
2135            if self.representation == REPRESENTATION_ISOTROPIC:
2136                scatt_grp.create_dataset("g_min", data=g_out_bounds[0, 0, :, 0])
2137                scatt_grp.create_dataset("g_max", data=g_out_bounds[0, 0, :, 1])
2138            elif self.representation == REPRESENTATION_ANGLE:
2139                scatt_grp.create_dataset("g_min", data=g_out_bounds[:, :, :, 0])
2140                scatt_grp.create_dataset("g_max", data=g_out_bounds[:, :, :, 1])
2141
2142            # Add the kinetics data
2143            if self._inverse_velocity[i] is not None:
2144                xs_grp.create_dataset("inverse-velocity",
2145                                      data=self._inverse_velocity[i])
2146
2147    @classmethod
2148    def from_hdf5(cls, group, name, energy_groups, num_delayed_groups):
2149        """Generate XSdata object from an HDF5 group
2150
2151        Parameters
2152        ----------
2153        group : h5py.Group
2154            HDF5 group to read from
2155        name : str
2156            Name of the mgxs data set.
2157        energy_groups : openmc.mgxs.EnergyGroups
2158            Energy group structure
2159        num_delayed_groups : int
2160            Number of delayed groups
2161
2162        Returns
2163        -------
2164        openmc.XSdata
2165            Multi-group cross section data
2166
2167        """
2168
2169        # Get a list of all the subgroups which will contain our temperature
2170        # strings
2171        subgroups = group.keys()
2172        temperatures = []
2173        for subgroup in subgroups:
2174            if subgroup != 'kTs':
2175                temperatures.append(subgroup)
2176
2177        # To ensure the actual floating point temperature used when creating
2178        # the new library is consistent with that used when originally creating
2179        # the file, get the floating point temperatures straight from the kTs
2180        # group.
2181        kTs_group = group['kTs']
2182        float_temperatures = []
2183        for temperature in temperatures:
2184            kT = kTs_group[temperature][()]
2185            float_temperatures.append(kT / openmc.data.K_BOLTZMANN)
2186
2187        attrs = group.attrs.keys()
2188        if 'representation' in attrs:
2189            representation = group.attrs['representation'].decode()
2190        else:
2191            representation = REPRESENTATION_ISOTROPIC
2192
2193        data = cls(name, energy_groups, float_temperatures, representation,
2194                   num_delayed_groups)
2195
2196        if 'scatter_format' in attrs:
2197            data.scatter_format = group.attrs['scatter_format'].decode()
2198
2199        # Get the remaining optional attributes
2200        if 'atomic_weight_ratio' in attrs:
2201            data.atomic_weight_ratio = group.attrs['atomic_weight_ratio']
2202        if 'order' in attrs:
2203            data.order = group.attrs['order']
2204        if data.representation == REPRESENTATION_ANGLE:
2205            data.num_azimuthal = group.attrs['num_azimuthal']
2206            data.num_polar = group.attrs['num_polar']
2207
2208        # Read the temperature-dependent datasets
2209        for temp, float_temp in zip(temperatures, float_temperatures):
2210            xs_types = ['total', 'absorption', 'fission', 'kappa-fission',
2211                        'chi', 'chi-prompt', 'chi-delayed', 'nu-fission',
2212                        'prompt-nu-fission', 'delayed-nu-fission', 'beta',
2213                        'decay-rate', 'inverse-velocity']
2214
2215            temperature_group = group[temp]
2216
2217            for xs_type in xs_types:
2218                set_func = 'set_' + xs_type.replace(' ', '_').replace('-', '_')
2219                if xs_type in temperature_group:
2220                    getattr(data, set_func)(temperature_group[xs_type][()],
2221                                            float_temp)
2222
2223            scatt_group = temperature_group['scatter_data']
2224
2225            # Get scatter matrix and 'un-flatten' it
2226            g_max = scatt_group['g_max']
2227            g_min = scatt_group['g_min']
2228            flat_scatter = scatt_group['scatter_matrix'][()]
2229            scatter_matrix = np.zeros(data.xs_shapes["[G][G'][Order]"])
2230            G = data.energy_groups.num_groups
2231            if data.representation == REPRESENTATION_ISOTROPIC:
2232                Np = 1
2233                Na = 1
2234            elif data.representation == REPRESENTATION_ANGLE:
2235                Np = data.num_polar
2236                Na = data.num_azimuthal
2237            flat_index = 0
2238            for p in range(Np):
2239                for a in range(Na):
2240                    for g_in in range(G):
2241                        if data.representation == REPRESENTATION_ISOTROPIC:
2242                            g_mins = g_min[g_in]
2243                            g_maxs = g_max[g_in]
2244                        elif data.representation == REPRESENTATION_ANGLE:
2245                            g_mins = g_min[p, a, g_in]
2246                            g_maxs = g_max[p, a, g_in]
2247                        for g_out in range(g_mins - 1, g_maxs):
2248                            for ang in range(data.num_orders):
2249                                if data.representation == REPRESENTATION_ISOTROPIC:
2250                                    scatter_matrix[g_in, g_out, ang] = \
2251                                        flat_scatter[flat_index]
2252                                elif data.representation == REPRESENTATION_ANGLE:
2253                                    scatter_matrix[p, a, g_in, g_out, ang] = \
2254                                        flat_scatter[flat_index]
2255                                flat_index += 1
2256            data.set_scatter_matrix(scatter_matrix, float_temp)
2257
2258            # Repeat for multiplicity
2259            if 'multiplicity_matrix' in scatt_group:
2260                flat_mult = scatt_group['multiplicity_matrix'][()]
2261                mult_matrix = np.zeros(data.xs_shapes["[G][G']"])
2262                flat_index = 0
2263                for p in range(Np):
2264                    for a in range(Na):
2265                        for g_in in range(G):
2266                            if data.representation == REPRESENTATION_ISOTROPIC:
2267                                g_mins = g_min[g_in]
2268                                g_maxs = g_max[g_in]
2269                            elif data.representation == REPRESENTATION_ANGLE:
2270                                g_mins = g_min[p, a, g_in]
2271                                g_maxs = g_max[p, a, g_in]
2272                            for g_out in range(g_mins - 1, g_maxs):
2273                                if data.representation == REPRESENTATION_ISOTROPIC:
2274                                    mult_matrix[g_in, g_out] = \
2275                                        flat_mult[flat_index]
2276                                elif data.representation == REPRESENTATION_ANGLE:
2277                                    mult_matrix[p, a, g_in, g_out] = \
2278                                        flat_mult[flat_index]
2279                                flat_index += 1
2280                data.set_multiplicity_matrix(mult_matrix, float_temp)
2281
2282        return data
2283
2284
2285class MGXSLibrary:
2286    """Multi-Group Cross Sections file used for an OpenMC simulation.
2287    Corresponds directly to the MG version of the cross_sections.xml input
2288    file.
2289
2290    Parameters
2291    ----------
2292    energy_groups : openmc.mgxs.EnergyGroups
2293        Energy group structure
2294    num_delayed_groups : int
2295        Num delayed groups
2296
2297    Attributes
2298    ----------
2299    energy_groups : openmc.mgxs.EnergyGroups
2300        Energy group structure.
2301    num_delayed_groups : int
2302        Num delayed groups
2303    xsdatas : Iterable of openmc.XSdata
2304        Iterable of multi-Group cross section data objects
2305    """
2306
2307    def __init__(self, energy_groups, num_delayed_groups=0):
2308        self.energy_groups = energy_groups
2309        self.num_delayed_groups = num_delayed_groups
2310        self._xsdatas = []
2311
2312    def __deepcopy__(self, memo):
2313        existing = memo.get(id(self))
2314
2315        # If this is the first time we have tried to copy this object, copy it
2316        if existing is None:
2317            clone = type(self).__new__(type(self))
2318            clone._energy_groups = copy.deepcopy(self.energy_groups, memo)
2319            clone._num_delayed_groups = self.num_delayed_groups
2320            clone._xsdatas = copy.deepcopy(self.xsdatas, memo)
2321
2322            memo[id(self)] = clone
2323
2324            return clone
2325
2326        # If this object has been copied before, return the first copy made
2327        else:
2328            return existing
2329
2330    @property
2331    def energy_groups(self):
2332        return self._energy_groups
2333
2334    @property
2335    def num_delayed_groups(self):
2336        return self._num_delayed_groups
2337
2338    @property
2339    def xsdatas(self):
2340        return self._xsdatas
2341
2342    @property
2343    def names(self):
2344        return [xsdata.name for xsdata in self.xsdatas]
2345
2346    @energy_groups.setter
2347    def energy_groups(self, energy_groups):
2348        check_type('energy groups', energy_groups, openmc.mgxs.EnergyGroups)
2349        self._energy_groups = energy_groups
2350
2351    @num_delayed_groups.setter
2352    def num_delayed_groups(self, num_delayed_groups):
2353        check_type('num_delayed_groups', num_delayed_groups, Integral)
2354        check_greater_than('num_delayed_groups', num_delayed_groups, 0,
2355                           equality=True)
2356        check_less_than('num_delayed_groups', num_delayed_groups,
2357                        openmc.mgxs.MAX_DELAYED_GROUPS, equality=True)
2358        self._num_delayed_groups = num_delayed_groups
2359
2360    def add_xsdata(self, xsdata):
2361        """Add an XSdata entry to the file.
2362
2363        Parameters
2364        ----------
2365        xsdata : openmc.XSdata
2366            MGXS information to add
2367
2368        """
2369
2370        if not isinstance(xsdata, XSdata):
2371            msg = 'Unable to add a non-XSdata "{0}" to the ' \
2372                  'MGXSLibrary instance'.format(xsdata)
2373            raise ValueError(msg)
2374
2375        if xsdata.energy_groups != self._energy_groups:
2376            msg = 'Energy groups of XSdata do not match that of MGXSLibrary.'
2377            raise ValueError(msg)
2378
2379        self._xsdatas.append(xsdata)
2380
2381    def add_xsdatas(self, xsdatas):
2382        """Add multiple XSdatas to the file.
2383
2384        Parameters
2385        ----------
2386        xsdatas : tuple or list of openmc.XSdata
2387            XSdatas to add
2388
2389        """
2390
2391        check_iterable_type('xsdatas', xsdatas, XSdata)
2392
2393        for xsdata in xsdatas:
2394            self.add_xsdata(xsdata)
2395
2396    def remove_xsdata(self, xsdata):
2397        """Remove a xsdata from the file
2398
2399        Parameters
2400        ----------
2401        xsdata : openmc.XSdata
2402            XSdata to remove
2403
2404        """
2405
2406        if not isinstance(xsdata, XSdata):
2407            msg = 'Unable to remove a non-XSdata "{0}" from the ' \
2408                  'MGXSLibrary instance'.format(xsdata)
2409            raise ValueError(msg)
2410
2411        self._xsdatas.remove(xsdata)
2412
2413    def get_by_name(self, name):
2414        """Access the XSdata objects by name
2415
2416        Parameters
2417        ----------
2418        name : str
2419            Name of openmc.XSdata object to obtain
2420
2421        Returns
2422        -------
2423        result : openmc.XSdata or None
2424            Provides the matching XSdata object or None, if not found
2425
2426        """
2427        check_type("name", name, str)
2428        result = None
2429        for xsdata in self.xsdatas:
2430            if name == xsdata.name:
2431                result = xsdata
2432        return result
2433
2434    def convert_representation(self, target_representation, num_polar=None,
2435                               num_azimuthal=None):
2436        """Produce a new XSdata object with the same data, but converted to the
2437        new representation (isotropic or angle-dependent).
2438
2439        This method cannot be used to change the number of polar or
2440        azimuthal bins of an XSdata object that already uses an angular
2441        representation. Finally, this method simply uses an arithmetic mean to
2442        convert from an angular to isotropic representation; no flux-weighting
2443        is applied and therefore the reaction rates will not be preserved.
2444
2445        Parameters
2446        ----------
2447        target_representation : {'isotropic', 'angle'}
2448            Representation of the MGXS (isotropic or angle-dependent flux
2449            weighting).
2450        num_polar : int, optional
2451            Number of equal width angular bins that the polar angular domain is
2452            subdivided into. This is required when `target_representation` is
2453            "angle".
2454        num_azimuthal : int, optional
2455            Number of equal width angular bins that the azimuthal angular domain
2456            is subdivided into. This is required when `target_representation` is
2457            "angle".
2458
2459        Returns
2460        -------
2461        openmc.MGXSLibrary
2462            Multi-group Library with the same data as self, but represented as
2463            specified in `target_representation`.
2464
2465        """
2466
2467        library = copy.deepcopy(self)
2468        for i, xsdata in enumerate(self.xsdatas):
2469            library.xsdatas[i] = \
2470                xsdata.convert_representation(target_representation,
2471                                              num_polar, num_azimuthal)
2472        return library
2473
2474    def convert_scatter_format(self, target_format, target_order):
2475        """Produce a new MGXSLibrary object with the same data, but converted
2476        to the new scatter format and order
2477
2478        Parameters
2479        ----------
2480        target_format : {'tabular', 'legendre', 'histogram'}
2481            Representation of the scattering angle distribution
2482        target_order : int
2483            Either the Legendre target_order, number of bins, or number of
2484            points used to describe the angular distribution associated with
2485            each group-to-group transfer probability
2486
2487        Returns
2488        -------
2489        openmc.MGXSLibrary
2490            Multi-group Library with the same data as self, but with the scatter
2491            format represented as specified in `target_format` and
2492            `target_order`.
2493
2494        """
2495
2496        library = copy.deepcopy(self)
2497        for i, xsdata in enumerate(self.xsdatas):
2498            library.xsdatas[i] = \
2499                xsdata.convert_scatter_format(target_format, target_order)
2500
2501        return library
2502
2503    def export_to_hdf5(self, filename='mgxs.h5', libver='earliest'):
2504        """Create an hdf5 file that can be used for a simulation.
2505
2506        Parameters
2507        ----------
2508        filename : str
2509            Filename of file, default is mgxs.h5.
2510        libver : {'earliest', 'latest'}
2511            Compatibility mode for the HDF5 file. 'latest' will produce files
2512            that are less backwards compatible but have performance benefits.
2513
2514        """
2515
2516        check_type('filename', filename, str)
2517
2518        # Create and write to the HDF5 file
2519        file = h5py.File(filename, "w", libver=libver)
2520        file.attrs['filetype'] = np.string_(_FILETYPE_MGXS_LIBRARY)
2521        file.attrs['version'] = [_VERSION_MGXS_LIBRARY, 0]
2522        file.attrs['energy_groups'] = self.energy_groups.num_groups
2523        file.attrs['delayed_groups'] = self.num_delayed_groups
2524        file.attrs['group structure'] = self.energy_groups.group_edges
2525
2526        for xsdata in self._xsdatas:
2527            xsdata.to_hdf5(file)
2528
2529        file.close()
2530
2531    @classmethod
2532    def from_hdf5(cls, filename=None):
2533        """Generate an MGXS Library from an HDF5 group or file
2534
2535        Parameters
2536        ----------
2537        filename : str, optional
2538            Name of HDF5 file containing MGXS data. Default is None.
2539            If not provided, the value of the OPENMC_MG_CROSS_SECTIONS
2540            environmental variable will be used
2541
2542        Returns
2543        -------
2544        openmc.MGXSLibrary
2545            Multi-group cross section data object.
2546
2547        """
2548        # If filename is None, get the cross sections from the
2549        # OPENMC_CROSS_SECTIONS environment variable
2550        if filename is None:
2551            filename = os.environ.get('OPENMC_MG_CROSS_SECTIONS')
2552
2553        # Check to make sure there was an environmental variable.
2554        if filename is None:
2555            raise ValueError("Either path or OPENMC_MG_CROSS_SECTIONS "
2556                             "environmental variable must be set")
2557
2558        check_type('filename', filename, str)
2559        file = h5py.File(filename, 'r')
2560
2561        # Check filetype and version
2562        check_filetype_version(file, _FILETYPE_MGXS_LIBRARY,
2563                               _VERSION_MGXS_LIBRARY)
2564
2565        group_structure = file.attrs['group structure']
2566        num_delayed_groups = file.attrs['delayed_groups']
2567        energy_groups = openmc.mgxs.EnergyGroups(group_structure)
2568        data = cls(energy_groups, num_delayed_groups)
2569
2570        for group_name, group in file.items():
2571            data.add_xsdata(openmc.XSdata.from_hdf5(group, group_name,
2572                                                    energy_groups,
2573                                                    num_delayed_groups))
2574
2575        return data
2576