• Home
  • Features
  • Pricing
  • Docs
  • Announcements
  • Sign In

openmc-dev / openmc / 12776996362

14 Jan 2025 09:49PM UTC coverage: 84.938% (+0.2%) from 84.729%
12776996362

Pull #3133

github

web-flow
Merge 0495246d9 into 549cc0973
Pull Request #3133: Kinetics parameters using Iterated Fission Probability

318 of 330 new or added lines in 10 files covered. (96.36%)

1658 existing lines in 66 files now uncovered.

50402 of 59340 relevant lines covered (84.94%)

33987813.96 hits per line

Source File
Press 'n' to go to next uncovered line, 'b' for previous

56.1
/openmc/plotter.py
1
from __future__ import annotations
12✔
2
from itertools import chain
12✔
3
from numbers import Integral, Real
12✔
4
from typing import Dict, Iterable, List
12✔
5

6
import numpy as np
12✔
7

8
import openmc.checkvalue as cv
12✔
9
import openmc.data
12✔
10

11
# Supported keywords for continuous-energy cross section plotting
12
PLOT_TYPES = {'total', 'scatter', 'elastic', 'inelastic', 'fission',
12✔
13
              'absorption', 'capture', 'nu-fission', 'nu-scatter', 'unity',
14
              'slowing-down power', 'damage'}
15

16
# Supported keywords for multi-group cross section plotting
17
PLOT_TYPES_MGXS = {'total', 'absorption', 'scatter', 'fission',
12✔
18
                   'kappa-fission', 'nu-fission', 'prompt-nu-fission',
19
                   'deleyed-nu-fission', 'chi', 'chi-prompt', 'chi-delayed',
20
                   'inverse-velocity', 'beta', 'decay-rate', 'unity'}
21
# Create a dictionary which can be used to convert PLOT_TYPES_MGXS to the
22
# openmc.XSdata attribute name needed to access the data
23
_PLOT_MGXS_ATTR = {line: line.replace(' ', '_').replace('-', '_')
12✔
24
                   for line in PLOT_TYPES_MGXS}
25
_PLOT_MGXS_ATTR['scatter'] = 'scatter_matrix'
12✔
26

27
# Special MT values
28
UNITY_MT = -1
12✔
29
XI_MT = -2
12✔
30

31
# MTs to combine to generate associated plot_types
32
_INELASTIC = [mt for mt in openmc.data.SUM_RULES[3] if mt != 27]
12✔
33
PLOT_TYPES_MT = {
12✔
34
    'total': openmc.data.SUM_RULES[1],
35
    'scatter': [2] + _INELASTIC,
36
    'elastic': [2],
37
    'inelastic': _INELASTIC,
38
    'fission': [18],
39
    'absorption': [27],
40
    'capture': [101],
41
    'nu-fission': [18],
42
    'nu-scatter': [2] + _INELASTIC,
43
    'unity': [UNITY_MT],
44
    'slowing-down power': [2] + [XI_MT],
45
    'damage': [444]
46
}
47

48
# Types of plots to plot linearly in y
49
PLOT_TYPES_LINEAR = {'nu-fission / fission', 'nu-scatter / scatter',
12✔
50
                     'nu-fission / absorption', 'fission / absorption'}
51

52
# Minimum and maximum energies for plotting (units of eV)
53
_MIN_E = 1.e-5
12✔
54
_MAX_E = 20.e6
12✔
55

56

57
ELEMENT_NAMES = list(openmc.data.ELEMENT_SYMBOL.values())[1:]
12✔
58

59

60
def _get_legend_label(this, type):
12✔
61
    """Gets a label for the element or nuclide or material and reaction plotted"""
62
    if isinstance(this, str):
12✔
63
        if type in openmc.data.DADZ:
12✔
64
            if this in ELEMENT_NAMES:
12✔
65
                return f'{this} {type}'
12✔
66
            else:  # this is a nuclide so the legend can contain more information
67
                z, a, m = openmc.data.zam(this)
12✔
68
                da, dz = openmc.data.DADZ[type]
12✔
69
                gnds_name = openmc.data.gnds_name(z + dz, a + da, m)
12✔
70
                # makes a string with nuclide reaction and new nuclide
71
                # For example "Be9 (n,2n) Be8"
72
                return f'{this} {type} {gnds_name}'
12✔
73
        return f'{this} {type}'
12✔
74
    elif this.name == '':
12✔
75
        return f'Material {this.id} {type}'
12✔
76
    else:
UNCOV
77
        return f'{this.name} {type}'
×
78

79

80
def _get_yaxis_label(reactions, divisor_types):
12✔
81
    """Gets a y axis label for the type of data plotted"""
82

83
    heat_values = {"heating", "heating-local", "damage-energy"}
12✔
84

85
    # if all the types are heating a different stem and unit is needed
86
    if all(set(value).issubset(heat_values) for value in reactions.values()):
12✔
87
        stem = "Heating"
12✔
88
    elif all(isinstance(item, str) for item in reactions.keys()):
12✔
89
        for nuc_reactions in reactions.values():
12✔
90
            for reaction in nuc_reactions:
12✔
91
                if reaction in heat_values:
12✔
92
                    raise TypeError(
12✔
93
                        "Mixture of heating and Microscopic reactions. "
94
                        "Invalid type for plotting"
95
                    )
96
        stem = "Microscopic"
12✔
97
    elif all(isinstance(item, openmc.Material) for item in reactions.keys()):
12✔
98
        stem = 'Macroscopic'
12✔
99
    else:
100
        msg = "Mixture of openmc.Material and elements/nuclides. Invalid type for plotting"
12✔
101
        raise TypeError(msg)
12✔
102

103
    if divisor_types:
12✔
UNCOV
104
        mid, units = "Data", ""
×
105
    else:
106
        mid = "Cross Section"
12✔
107
        units = {
12✔
108
            "Macroscopic": "[1/cm]",
109
            "Microscopic": "[b]",
110
            "Heating": "[eV-barn]",
111
        }[stem]
112

113
    return f'{stem} {mid} {units}'
12✔
114

115
def _get_title(reactions):
12✔
116
    """Gets a title for the type of data plotted"""
117
    if len(reactions) == 1:
12✔
118
        this, = reactions
12✔
119
        name = this.name if isinstance(this, openmc.Material) else this
12✔
120
        return f'Cross Section Plot For {name}'
12✔
121
    else:
122
        return 'Cross Section Plot'
12✔
123

124

125
def plot_xs(
12✔
126
    reactions: Dict[str | openmc.Material, List[str]],
127
    divisor_types: Iterable[str] | None = None,
128
    temperature: float = 294.0,
129
    axis: "plt.Axes" | None = None,
130
    sab_name: str | None = None,
131
    ce_cross_sections: str | None = None,
132
    mg_cross_sections: str | None = None,
133
    enrichment: float | None = None,
134
    plot_CE: bool = True,
135
    orders: Iterable[int] | None = None,
136
    divisor_orders: Iterable[int] | None = None,
137
    energy_axis_units: str = "eV",
138
    **kwargs,
139
) -> "plt.Figure" | None:
140
    """Creates a figure of continuous-energy cross sections for this item.
141

142
    Parameters
143
    ----------
144
    reactions : dict
145
        keys can be either a nuclide or element in string form or an
146
        openmc.Material object. Values are a list of the types of
147
        cross sections to include in the plot.
148
    divisor_types : Iterable of values of PLOT_TYPES, optional
149
        Cross section types which will divide those produced by types
150
        before plotting. A type of 'unity' can be used to effectively not
151
        divide some types.
152
    temperature : float, optional
153
        Temperature in Kelvin to plot. If not specified, a default
154
        temperature of 294K will be plotted. Note that the nearest
155
        temperature in the library for each nuclide will be used as opposed
156
        to using any interpolation.
157
    axis : matplotlib.axes, optional
158
        A previously generated axis to use for plotting. If not specified,
159
        a new axis and figure will be generated.
160
    sab_name : str, optional
161
        Name of S(a,b) library to apply to MT=2 data when applicable.
162
    ce_cross_sections : str, optional
163
        Location of cross_sections.xml file. Default is None.
164
    mg_cross_sections : str, optional
165
        Location of MGXS HDF5 Library file. Default is None.
166
    enrichment : float, optional
167
        Enrichment for U235 in weight percent. For example, input 4.95 for
168
        4.95 weight percent enriched U. Default is None.
169
    plot_CE : bool, optional
170
        Denotes whether or not continuous-energy will be plotted. Defaults to
171
        plotting the continuous-energy data.
172
    orders : Iterable of Integral, optional
173
        The scattering order or delayed group index to use for the
174
        corresponding entry in types. Defaults to the 0th order for scattering
175
        and the total delayed neutron data. This only applies to plots of
176
        multi-group data.
177
    divisor_orders : Iterable of Integral, optional
178
        Same as orders, but for divisor_types
179
    **kwargs :
180
        All keyword arguments are passed to
181
        :func:`matplotlib.pyplot.figure`.
182
    energy_axis_units : {'eV', 'keV', 'MeV'}
183
        Units used on the plot energy axis
184

185
        .. versionadded:: 0.15.0
186

187
    Returns
188
    -------
189
    fig : matplotlib.figure.Figure
190
        If axis is None, then a Matplotlib Figure of the generated
191
        cross section will be returned. Otherwise, a value of
192
        None will be returned as the figure and axes have already been
193
        generated.
194

195
    """
196
    import matplotlib.pyplot as plt
12✔
197

198
    cv.check_type("plot_CE", plot_CE, bool)
12✔
199
    cv.check_value("energy_axis_units", energy_axis_units, {"eV", "keV", "MeV"})
12✔
200

201
    axis_scaling_factor = {"eV": 1.0, "keV": 1e-3, "MeV": 1e-6}
12✔
202

203
    # Generate the plot
204
    if axis is None:
12✔
205
        fig, ax = plt.subplots(**kwargs)
12✔
206
    else:
UNCOV
207
        fig = None
×
208
        ax = axis
×
209

210
    all_types = []
12✔
211

212
    for this, types in reactions.items():
12✔
213
        all_types = all_types + types
12✔
214

215
        if plot_CE:
12✔
216
            cv.check_type("this", this, (str, openmc.Material))
12✔
217
            # Calculate for the CE cross sections
218
            E, data = calculate_cexs(this, types, temperature, sab_name,
12✔
219
                                    ce_cross_sections, enrichment)
220
            if divisor_types:
12✔
UNCOV
221
                cv.check_length('divisor types', divisor_types, len(types))
×
UNCOV
222
                Ediv, data_div = calculate_cexs(this, divisor_types, temperature,
×
223
                                                sab_name, ce_cross_sections,
224
                                                enrichment)
225

226
                # Create a new union grid, interpolate data and data_div on to that
227
                # grid, and then do the actual division
228
                Enum = E[:]
×
UNCOV
229
                E = np.union1d(Enum, Ediv)
×
UNCOV
230
                data_new = np.zeros((len(types), len(E)))
×
231

232
                for line in range(len(types)):
×
233
                    data_new[line, :] = \
×
234
                        np.divide(np.interp(E, Enum, data[line, :]),
235
                                np.interp(E, Ediv, data_div[line, :]))
UNCOV
236
                    if divisor_types[line] != 'unity':
×
UNCOV
237
                        types[line] = types[line] + ' / ' + divisor_types[line]
×
UNCOV
238
                data = data_new
×
239
        else:
240
            # Calculate for MG cross sections
241
            E, data = calculate_mgxs(this, types, orders, temperature,
×
242
                                    mg_cross_sections, ce_cross_sections,
243
                                    enrichment)
244

UNCOV
245
            if divisor_types:
×
UNCOV
246
                cv.check_length('divisor types', divisor_types, len(types))
×
UNCOV
247
                Ediv, data_div = calculate_mgxs(this, divisor_types,
×
248
                                                divisor_orders, temperature,
249
                                                mg_cross_sections,
250
                                                ce_cross_sections, enrichment)
251

252
                # Perform the division
UNCOV
253
                for line in range(len(types)):
×
UNCOV
254
                    data[line, :] /= data_div[line, :]
×
UNCOV
255
                    if divisor_types[line] != 'unity':
×
256
                        types[line] += ' / ' + divisor_types[line]
×
257

258
        E *= axis_scaling_factor[energy_axis_units]
12✔
259

260
        # Plot the data
261
        for i in range(len(data)):
12✔
262
            data[i, :] = np.nan_to_num(data[i, :])
12✔
263
            if np.sum(data[i, :]) > 0.:
12✔
264
                ax.plot(E, data[i, :], label=_get_legend_label(this, types[i]))
12✔
265

266
    # Set to loglog or semilogx depending on if we are plotting a data
267
    # type which we expect to vary linearly
268
    if set(all_types).issubset(PLOT_TYPES_LINEAR):
12✔
269
        ax.set_xscale('log')
×
UNCOV
270
        ax.set_yscale('linear')
×
271
    else:
272
        ax.set_xscale('log')
12✔
273
        ax.set_yscale('log')
12✔
274

275
    ax.set_xlabel(f"Energy [{energy_axis_units}]")
12✔
276
    if plot_CE:
12✔
277
        ax.set_xlim(
12✔
278
            _MIN_E * axis_scaling_factor[energy_axis_units],
279
            _MAX_E * axis_scaling_factor[energy_axis_units],
280
        )
281
    else:
UNCOV
282
        ax.set_xlim(E[-1], E[0])
×
283

284
    ax.set_ylabel(_get_yaxis_label(reactions, divisor_types))
12✔
285
    ax.legend(loc='best')
12✔
286
    ax.set_title(_get_title(reactions))
12✔
287

288
    return fig
12✔
289

290

291
def calculate_cexs(this, types, temperature=294., sab_name=None,
12✔
292
                   cross_sections=None, enrichment=None, ncrystal_cfg=None):
293
    """Calculates continuous-energy cross sections of a requested type.
294

295
    Parameters
296
    ----------
297
    this : str or openmc.Material
298
        Object to source data from. Nuclides and elements should be input as a
299
        str
300
    types : Iterable of values of PLOT_TYPES
301
        The type of cross sections to calculate
302
    temperature : float, optional
303
        Temperature in Kelvin to plot. If not specified, a default
304
        temperature of 294K will be plotted. Note that the nearest
305
        temperature in the library for each nuclide will be used as opposed
306
        to using any interpolation.
307
    sab_name : str, optional
308
        Name of S(a,b) library to apply to MT=2 data when applicable.
309
    cross_sections : str, optional
310
        Location of cross_sections.xml file. Default is None.
311
    enrichment : float, optional
312
        Enrichment for U235 in weight percent. For example, input 4.95 for
313
        4.95 weight percent enriched U. Default is None
314
        (natural composition).
315
    ncrystal_cfg : str, optional
316
        Configuration string for NCrystal material.
317

318
    Returns
319
    -------
320
    energy_grid : numpy.ndarray
321
        Energies at which cross sections are calculated, in units of eV
322
    data : numpy.ndarray
323
        Cross sections calculated at the energy grid described by energy_grid
324

325
    """
326

327
    # Check types
328
    cv.check_type('this', this, (str, openmc.Material))
12✔
329
    cv.check_type('temperature', temperature, Real)
12✔
330
    if sab_name:
12✔
331
        cv.check_type('sab_name', sab_name, str)
12✔
332
    if enrichment:
12✔
UNCOV
333
        cv.check_type('enrichment', enrichment, Real)
×
334

335
    if isinstance(this, str):
12✔
336
        if this in ELEMENT_NAMES:
12✔
337
            energy_grid, data = _calculate_cexs_elem_mat(
12✔
338
                this, types, temperature, cross_sections, sab_name, enrichment
339
            )
340
        else:
341
            energy_grid, xs = _calculate_cexs_nuclide(
12✔
342
                this, types, temperature, sab_name, cross_sections,
343
                ncrystal_cfg
344
            )
345

346
            # Convert xs (Iterable of Callable) to a grid of cross section values
347
            # calculated on the points in energy_grid for consistency with the
348
            # element and material functions.
349
            data = np.zeros((len(types), len(energy_grid)))
12✔
350
            for line in range(len(types)):
12✔
351
                data[line, :] = xs[line](energy_grid)
12✔
352
    else:
353
        energy_grid, data = _calculate_cexs_elem_mat(this, types, temperature,
12✔
354
                                                     cross_sections)
355

356
    return energy_grid, data
12✔
357

358

359
def _calculate_cexs_nuclide(this, types, temperature=294., sab_name=None,
12✔
360
                            cross_sections=None, ncrystal_cfg=None):
361
    """Calculates continuous-energy cross sections of a requested type.
362

363
    Parameters
364
    ----------
365
    this : str
366
        Nuclide object to source data from
367
    types : Iterable of str or Integral
368
        The type of cross sections to calculate; values can either be those
369
        in openmc.PLOT_TYPES or keys from openmc.data.REACTION_MT which
370
        correspond to a reaction description e.g '(n,2n)' or integers which
371
        correspond to reaction channel (MT) numbers.
372
    temperature : float, optional
373
        Temperature in Kelvin to plot. If not specified, a default
374
        temperature of 294K will be plotted. Note that the nearest
375
        temperature in the library for each nuclide will be used as opposed
376
        to using any interpolation.
377
    sab_name : str, optional
378
        Name of S(a,b) library to apply to MT=2 data when applicable.
379
    cross_sections : str, optional
380
        Location of cross_sections.xml file. Default is None.
381
    ncrystal_cfg : str, optional
382
        Configuration string for NCrystal material.
383

384
    Returns
385
    -------
386
    energy_grid : numpy.ndarray
387
        Energies at which cross sections are calculated, in units of eV
388
    data : Iterable of Callable
389
        Requested cross section functions
390

391
    """
392

393
    # Load the library
394
    library = openmc.data.DataLibrary.from_xml(cross_sections)
12✔
395

396
    # Convert temperature to format needed for access in the library
397
    strT = f"{int(round(temperature))}K"
12✔
398
    T = temperature
12✔
399

400
    # Now we can create the data sets to be plotted
401
    energy_grid = []
12✔
402
    xs = []
12✔
403
    lib = library.get_by_material(this)
12✔
404
    if lib is not None:
12✔
405
        nuc = openmc.data.IncidentNeutron.from_hdf5(lib['path'])
12✔
406
        # Obtain the nearest temperature
407
        if strT in nuc.temperatures:
12✔
408
            nucT = strT
12✔
409
        else:
UNCOV
410
            delta_T = np.array(nuc.kTs) - T * openmc.data.K_BOLTZMANN
×
UNCOV
411
            closest_index = np.argmin(np.abs(delta_T))
×
UNCOV
412
            nucT = nuc.temperatures[closest_index]
×
413

414
        # Prep S(a,b) data if needed
415
        if sab_name:
12✔
416
            sab = openmc.data.ThermalScattering.from_hdf5(sab_name)
12✔
417
            # Obtain the nearest temperature
418
            if strT in sab.temperatures:
12✔
419
                sabT = strT
×
420
            else:
421
                delta_T = np.array(sab.kTs) - T * openmc.data.K_BOLTZMANN
12✔
422
                closest_index = np.argmin(np.abs(delta_T))
12✔
423
                sabT = sab.temperatures[closest_index]
12✔
424

425
            # Create an energy grid composed the S(a,b) and the nuclide's grid
426
            grid = nuc.energy[nucT]
12✔
427
            sab_Emax = 0.
12✔
428
            sab_funcs = []
12✔
429
            if sab.elastic is not None:
12✔
UNCOV
430
                elastic = sab.elastic.xs[sabT]
×
UNCOV
431
                if isinstance(elastic, openmc.data.CoherentElastic):
×
UNCOV
432
                    grid = np.union1d(grid, elastic.bragg_edges)
×
UNCOV
433
                    if elastic.bragg_edges[-1] > sab_Emax:
×
UNCOV
434
                        sab_Emax = elastic.bragg_edges[-1]
×
UNCOV
435
                elif isinstance(elastic, openmc.data.Tabulated1D):
×
UNCOV
436
                    grid = np.union1d(grid, elastic.x)
×
UNCOV
437
                    if elastic.x[-1] > sab_Emax:
×
UNCOV
438
                        sab_Emax = elastic.x[-1]
×
UNCOV
439
                sab_funcs.append(elastic)
×
440
            if sab.inelastic is not None:
12✔
441
                inelastic = sab.inelastic.xs[sabT]
12✔
442
                grid = np.union1d(grid, inelastic.x)
12✔
443
                if inelastic.x[-1] > sab_Emax:
12✔
444
                        sab_Emax = inelastic.x[-1]
12✔
445
                sab_funcs.append(inelastic)
12✔
446
            energy_grid = grid
12✔
447
        else:
448
            energy_grid = nuc.energy[nucT]
12✔
449

450
        # Parse the types
451
        mts = []
12✔
452
        ops = []
12✔
453
        yields = []
12✔
454
        for line in types:
12✔
455
            if line in PLOT_TYPES:
12✔
456
                tmp_mts = [mtj for mti in PLOT_TYPES_MT[line] for mtj in
12✔
457
                           nuc.get_reaction_components(mti)]
458
                mts.append(tmp_mts)
12✔
459
                if line.startswith('nu'):
12✔
UNCOV
460
                    yields.append(True)
×
461
                else:
462
                    yields.append(False)
12✔
463
                if XI_MT in tmp_mts:
12✔
UNCOV
464
                    ops.append((np.add,) * (len(tmp_mts) - 2) + (np.multiply,))
×
465
                else:
466
                    ops.append((np.add,) * (len(tmp_mts) - 1))
12✔
467
            elif line in openmc.data.REACTION_MT:
12✔
468
                mt_number = openmc.data.REACTION_MT[line]
12✔
469
                cv.check_type('MT in types', mt_number, Integral)
12✔
470
                cv.check_greater_than('MT in types', mt_number, 0)
12✔
471
                tmp_mts = nuc.get_reaction_components(mt_number)
12✔
472
                mts.append(tmp_mts)
12✔
473
                ops.append((np.add,) * (len(tmp_mts) - 1))
12✔
474
                yields.append(False)
12✔
475
            elif isinstance(line, int):
12✔
476
                # Not a built-in type, we have to parse it ourselves
477
                cv.check_type('MT in types', line, Integral)
12✔
478
                cv.check_greater_than('MT in types', line, 0)
12✔
479
                tmp_mts = nuc.get_reaction_components(line)
12✔
480
                mts.append(tmp_mts)
12✔
481
                ops.append((np.add,) * (len(tmp_mts) - 1))
12✔
482
                yields.append(False)
12✔
483
            else:
UNCOV
484
                raise TypeError("Invalid type", line)
×
485

486
        for i, mt_set in enumerate(mts):
12✔
487
            # Get the reaction xs data from the nuclide
488
            funcs = []
12✔
489
            op = ops[i]
12✔
490
            for mt in mt_set:
12✔
491
                if mt == 2:
12✔
492
                    if sab_name:
12✔
493
                        # Then we need to do a piece-wise function of
494
                        # The S(a,b) and non-thermal data
495
                        sab_sum = openmc.data.Sum(sab_funcs)
12✔
496
                        pw_funcs = openmc.data.Regions1D(
12✔
497
                            [sab_sum, nuc[mt].xs[nucT]],
498
                            [sab_Emax])
499
                        funcs.append(pw_funcs)
12✔
500
                    elif ncrystal_cfg:
12✔
UNCOV
501
                        import NCrystal
×
UNCOV
502
                        nc_scatter = NCrystal.createScatter(ncrystal_cfg)
×
UNCOV
503
                        nc_func = nc_scatter.crossSectionNonOriented
×
UNCOV
504
                        nc_emax = 5 # eV # this should be obtained from NCRYSTAL_MAX_ENERGY
×
UNCOV
505
                        energy_grid = np.union1d(np.geomspace(min(energy_grid),
×
506
                                                              1.1*nc_emax,
507
                                                              1000),energy_grid) # NCrystal does not have
508
                                                                                 # an intrinsic energy grid
509
                        pw_funcs = openmc.data.Regions1D(
×
510
                            [nc_func, nuc[mt].xs[nucT]],
511
                            [nc_emax])
UNCOV
512
                        funcs.append(pw_funcs)
×
513
                    else:
514
                        funcs.append(nuc[mt].xs[nucT])
12✔
515
                elif mt in nuc:
12✔
516
                    if yields[i]:
12✔
517
                        # Get the total yield first if available. This will be
518
                        # used primarily for fission.
UNCOV
519
                        for prod in chain(nuc[mt].products,
×
520
                                          nuc[mt].derived_products):
521
                            if prod.particle == 'neutron' and \
×
522
                                prod.emission_mode == 'total':
523
                                func = openmc.data.Combination(
×
524
                                    [nuc[mt].xs[nucT], prod.yield_],
525
                                    [np.multiply])
526
                                funcs.append(func)
×
UNCOV
527
                                break
×
528
                        else:
529
                            # Total doesn't exist so we have to create from
530
                            # prompt and delayed. This is used for scatter
531
                            # multiplication.
UNCOV
532
                            func = None
×
UNCOV
533
                            for prod in chain(nuc[mt].products,
×
534
                                              nuc[mt].derived_products):
UNCOV
535
                                if prod.particle == 'neutron' and \
×
536
                                    prod.emission_mode != 'total':
UNCOV
537
                                    if func:
×
538
                                        func = openmc.data.Combination(
×
539
                                            [prod.yield_, func], [np.add])
540
                                    else:
541
                                        func = prod.yield_
×
542
                            if func:
×
543
                                funcs.append(openmc.data.Combination(
×
544
                                    [func, nuc[mt].xs[nucT]], [np.multiply]))
545
                            else:
546
                                # If func is still None, then there were no
547
                                # products. In that case, assume the yield is
548
                                # one as its not provided for some summed
549
                                # reactions like MT=4
UNCOV
550
                                funcs.append(nuc[mt].xs[nucT])
×
551
                    else:
552
                        funcs.append(nuc[mt].xs[nucT])
12✔
553
                elif mt == UNITY_MT:
×
UNCOV
554
                    funcs.append(lambda x: 1.)
×
UNCOV
555
                elif mt == XI_MT:
×
UNCOV
556
                    awr = nuc.atomic_weight_ratio
×
UNCOV
557
                    alpha = ((awr - 1.) / (awr + 1.))**2
×
UNCOV
558
                    xi = 1. + alpha * np.log(alpha) / (1. - alpha)
×
UNCOV
559
                    funcs.append(lambda x: xi)
×
560
                else:
UNCOV
561
                    funcs.append(lambda x: 0.)
×
562
            funcs = funcs if funcs else [lambda x: 0.]
12✔
563
            xs.append(openmc.data.Combination(funcs, op))
12✔
564
    else:
UNCOV
565
        raise ValueError(this + " not in library")
×
566

567
    return energy_grid, xs
12✔
568

569

570
def _calculate_cexs_elem_mat(this, types, temperature=294.,
12✔
571
                             cross_sections=None, sab_name=None,
572
                             enrichment=None):
573
    """Calculates continuous-energy cross sections of a requested type.
574

575
    Parameters
576
    ----------
577
    this : openmc.Material or str
578
        Object to source data from. Element can be input as str
579
    types : Iterable of values of PLOT_TYPES
580
        The type of cross sections to calculate
581
    temperature : float, optional
582
        Temperature in Kelvin to plot. If not specified, a default
583
        temperature of 294K will be plotted. Note that the nearest
584
        temperature in the library for each nuclide will be used as opposed
585
        to using any interpolation.
586
    cross_sections : str, optional
587
        Location of cross_sections.xml file. Default is None.
588
    sab_name : str, optional
589
        Name of S(a,b) library to apply to MT=2 data when applicable.
590
    enrichment : float, optional
591
        Enrichment for U235 in weight percent. For example, input 4.95 for
592
        4.95 weight percent enriched U. Default is None
593
        (natural composition).
594

595
    Returns
596
    -------
597
    energy_grid : numpy.ndarray
598
        Energies at which cross sections are calculated, in units of eV
599
    data : numpy.ndarray
600
        Cross sections calculated at the energy grid described by energy_grid
601

602
    """
603

604
    if isinstance(this, openmc.Material):
12✔
605
        if this.temperature is not None:
12✔
UNCOV
606
            T = this.temperature
×
607
        else:
608
            T = temperature
12✔
609
    else:
610
        T = temperature
12✔
611

612
    # Load the library
613
    library = openmc.data.DataLibrary.from_xml(cross_sections)
12✔
614

615
    ncrystal_cfg = None
12✔
616
    if isinstance(this, openmc.Material):
12✔
617
        # Expand elements in to nuclides with atomic densities
618
        nuc_fractions = this.get_nuclide_atom_densities()
12✔
619
        # Create a dict of [nuclide name] = nuclide object to carry forward
620
        # with a common nuclides format between openmc.Material and Elements
621
        nuclides = {nuclide: nuclide for nuclide in nuc_fractions}
12✔
622
        # Add NCrystal cfg string if it exists
623
        ncrystal_cfg = this.ncrystal_cfg
12✔
624
    else:
625
        # Expand elements in to nuclides with atomic densities
626
        nuclides = openmc.Element(this).expand(1., 'ao', enrichment=enrichment,
12✔
627
                               cross_sections=cross_sections)
628
        # For ease of processing split out the nuclide and its fraction
629
        nuc_fractions = {nuclide[0]: nuclide[1] for nuclide in nuclides}
12✔
630
        # Create a dict of [nuclide name] = nuclide object to carry forward
631
        # with a common nuclides format between openmc.Material and Elements
632
        nuclides = {nuclide[0]: nuclide[0] for nuclide in nuclides}
12✔
633

634
    # Identify the nuclides which have S(a,b) data
635
    sabs = {}
12✔
636
    for nuclide in nuclides.items():
12✔
637
        sabs[nuclide[0]] = None
12✔
638
    if isinstance(this, openmc.Material):
12✔
639
        for sab_name, _ in this._sab:
12✔
640
            sab = openmc.data.ThermalScattering.from_hdf5(
12✔
641
                library.get_by_material(sab_name, data_type='thermal')['path'])
642
            for nuc in sab.nuclides:
12✔
643
                sabs[nuc] = library.get_by_material(sab_name,
12✔
644
                        data_type='thermal')['path']
645
    else:
646
        if sab_name:
12✔
UNCOV
647
            sab = openmc.data.ThermalScattering.from_hdf5(sab_name)
×
UNCOV
648
            for nuc in sab.nuclides:
×
UNCOV
649
                sabs[nuc] = library.get_by_material(sab_name,
×
650
                        data_type='thermal')['path']
651

652
    # Now we can create the data sets to be plotted
653
    xs = {}
12✔
654
    E = []
12✔
655
    for nuclide in nuclides.items():
12✔
656
        name = nuclide[0]
12✔
657
        nuc = nuclide[1]
12✔
658
        sab_tab = sabs[name]
12✔
659
        temp_E, temp_xs = calculate_cexs(nuc, types, T, sab_tab, cross_sections,
12✔
660
                                         ncrystal_cfg=ncrystal_cfg
661
                                         )
662
        E.append(temp_E)
12✔
663
        # Since the energy grids are different, store the cross sections as
664
        # a tabulated function so they can be calculated on any grid needed.
665
        xs[name] = [openmc.data.Tabulated1D(temp_E, temp_xs[line])
12✔
666
                    for line in range(len(types))]
667

668
    # Condense the data for every nuclide
669
    # First create a union energy grid
670
    energy_grid = E[0]
12✔
671
    for grid in E[1:]:
12✔
672
        energy_grid = np.union1d(energy_grid, grid)
12✔
673

674
    # Now we can combine all the nuclidic data
675
    data = np.zeros((len(types), len(energy_grid)))
12✔
676
    for line in range(len(types)):
12✔
677
        if types[line] == 'unity':
12✔
UNCOV
678
            data[line, :] = 1.
×
679
        else:
680
            for nuclide in nuclides.items():
12✔
681
                name = nuclide[0]
12✔
682
                data[line, :] += (nuc_fractions[name] *
12✔
683
                                  xs[name][line](energy_grid))
684

685
    return energy_grid, data
12✔
686

687

688
def calculate_mgxs(this, types, orders=None, temperature=294.,
12✔
689
                   cross_sections=None, ce_cross_sections=None,
690
                   enrichment=None):
691
    """Calculates multi-group cross sections of a requested type.
692

693
    If the data for the nuclide or macroscopic object in the library is
694
    represented as angle-dependent data then this method will return the
695
    geometric average cross section over all angles.
696

697
    Parameters
698
    ----------
699
    this : str or openmc.Material
700
        Object to source data from. Nuclides and elements can be input as a str
701
    types : Iterable of values of PLOT_TYPES_MGXS
702
        The type of cross sections to calculate
703
    orders : Iterable of Integral, optional
704
        The scattering order or delayed group index to use for the
705
        corresponding entry in types. Defaults to the 0th order for scattering
706
        and the total delayed neutron data.
707
    temperature : float, optional
708
        Temperature in Kelvin to plot. If not specified, a default
709
        temperature of 294K will be plotted. Note that the nearest
710
        temperature in the library for each nuclide will be used as opposed
711
        to using any interpolation.
712
    cross_sections : str, optional
713
        Location of MGXS HDF5 Library file. Default is None.
714
    ce_cross_sections : str, optional
715
        Location of continuous-energy cross_sections.xml file. Default is None.
716
    enrichment : float, optional
717
        Enrichment for U235 in weight percent. For example, input 4.95 for
718
        4.95 weight percent enriched U. Default is None
719
        (natural composition).
720

721
    Returns
722
    -------
723
    energy_grid : numpy.ndarray
724
        Energies at which cross sections are calculated, in units of eV
725
    data : numpy.ndarray
726
        Cross sections calculated at the energy grid described by energy_grid
727

728
    """
729

730
    # Check types
731
    cv.check_type('temperature', temperature, Real)
×
732
    if enrichment:
×
UNCOV
733
        cv.check_type('enrichment', enrichment, Real)
×
UNCOV
734
    cv.check_iterable_type('types', types, str)
×
735

UNCOV
736
    cv.check_type("cross_sections", cross_sections, str)
×
UNCOV
737
    library = openmc.MGXSLibrary.from_hdf5(cross_sections)
×
738

739
    if this in ELEMENT_NAMES or isinstance(this, openmc.Material):
×
740
        mgxs = _calculate_mgxs_elem_mat(this, types, library, orders,
×
741
                                        temperature, ce_cross_sections,
742
                                        enrichment)
UNCOV
743
    elif isinstance(this, str):
×
UNCOV
744
        mgxs = _calculate_mgxs_nuc_macro(this, types, library, orders,
×
745
                                         temperature)
746
    else:
747
        raise TypeError("Invalid type")
×
748

749
    # Convert the data to the format needed
UNCOV
750
    data = np.zeros((len(types), 2 * library.energy_groups.num_groups))
×
751
    energy_grid = np.zeros(2 * library.energy_groups.num_groups)
×
UNCOV
752
    for g in range(library.energy_groups.num_groups):
×
UNCOV
753
        energy_grid[g * 2: g * 2 + 2] = \
×
754
            library.energy_groups.group_edges[g: g + 2]
755
    # Ensure the energy will show on a log-axis by replacing 0s with a
756
    # sufficiently small number
UNCOV
757
    energy_grid[0] = max(energy_grid[0], _MIN_E)
×
758

UNCOV
759
    for line in range(len(types)):
×
UNCOV
760
        for g in range(library.energy_groups.num_groups):
×
UNCOV
761
            data[line, g * 2: g * 2 + 2] = mgxs[line, g]
×
762

UNCOV
763
    return energy_grid[::-1], data
×
764

765

766
def _calculate_mgxs_nuc_macro(this, types, library, orders=None,
12✔
767
                              temperature=294.):
768
    """Determines the multi-group cross sections of a nuclide or macroscopic
769
    object.
770

771
    If the data for the nuclide or macroscopic object in the library is
772
    represented as angle-dependent data then this method will return the
773
    geometric average cross section over all angles.
774

775
    Parameters
776
    ----------
777
    this : str
778
        Object to source data from
779
    types : Iterable of str
780
        The type of cross sections to calculate; values can either be those
781
        in openmc.PLOT_TYPES_MGXS
782
    library : openmc.MGXSLibrary
783
        MGXS Library containing the data of interest
784
    orders : Iterable of Integral, optional
785
        The scattering order or delayed group index to use for the
786
        corresponding entry in types. Defaults to the 0th order for scattering
787
        and the total delayed neutron data.
788
    temperature : float, optional
789
        Temperature in Kelvin to plot. If not specified, a default
790
        temperature of 294K will be plotted. Note that the nearest
791
        temperature in the library for each nuclide will be used as opposed
792
        to using any interpolation.
793

794
    Returns
795
    -------
796
    data : numpy.ndarray
797
        Cross sections calculated at the energy grid described by energy_grid
798

799
    """
800

801
    # Check the parameters and grab order/delayed groups
UNCOV
802
    if orders:
×
803
        cv.check_iterable_type('orders', orders, Integral,
×
804
                               min_depth=len(types), max_depth=len(types))
805
    else:
UNCOV
806
        orders = [None] * len(types)
×
UNCOV
807
    for i, line in enumerate(types):
×
808
        cv.check_type("line", line, str)
×
809
        cv.check_value("line", line, PLOT_TYPES_MGXS)
×
810
        if orders[i]:
×
811
            cv.check_greater_than("order value", orders[i], 0, equality=True)
×
812

813
    xsdata = library.get_by_name(this)
×
814

UNCOV
815
    if xsdata is not None:
×
816
        # Obtain the nearest temperature
UNCOV
817
        t = np.abs(xsdata.temperatures - temperature).argmin()
×
818

819
        # Get the data
UNCOV
820
        data = np.zeros((len(types), library.energy_groups.num_groups))
×
UNCOV
821
        for i, line in enumerate(types):
×
UNCOV
822
            if 'fission' in line and not xsdata.fissionable:
×
823
                continue
×
824
            elif line == 'unity':
×
UNCOV
825
                data[i, :] = 1.
×
826
            else:
827
                # Now we have to get the cross section data and properly
828
                # treat it depending on the requested type.
829
                # First get the data in a generic fashion
UNCOV
830
                temp_data = getattr(xsdata, _PLOT_MGXS_ATTR[line])[t]
×
UNCOV
831
                shape = temp_data.shape[:]
×
832
                # If we have angular data, then want the geometric
833
                # average over all provided angles.  Since the angles are
834
                # equi-distant, un-weighted averaging will suffice
UNCOV
835
                if xsdata.representation == 'angle':
×
UNCOV
836
                    temp_data = np.mean(temp_data, axis=(0, 1))
×
837

838
                # Now we can look at the shape of the data to identify how
839
                # it should be modified to produce an array of values
840
                # with groups.
UNCOV
841
                if shape in (xsdata.xs_shapes["[G']"],
×
842
                             xsdata.xs_shapes["[G]"]):
843
                    # Then the data is already an array vs groups so copy
844
                    # and move along
845
                    data[i, :] = temp_data
×
UNCOV
846
                elif shape == xsdata.xs_shapes["[G][G']"]:
×
847
                    # Sum the data over outgoing groups to create our array vs
848
                    # groups
UNCOV
849
                    data[i, :] = np.sum(temp_data, axis=1)
×
UNCOV
850
                elif shape == xsdata.xs_shapes["[DG]"]:
×
851
                    # Then we have a constant vs groups with a value for each
852
                    # delayed group. The user-provided value of orders tells us
853
                    # which delayed group we want. If none are provided, then
854
                    # we sum all the delayed groups together.
855
                    if orders[i]:
×
856
                        if orders[i] < len(shape[0]):
×
UNCOV
857
                            data[i, :] = temp_data[orders[i]]
×
858
                    else:
859
                        data[i, :] = np.sum(temp_data[:])
×
UNCOV
860
                elif shape in (xsdata.xs_shapes["[DG][G']"],
×
861
                               xsdata.xs_shapes["[DG][G]"]):
862
                    # Then we have an array vs groups with values for each
863
                    # delayed group. The user-provided value of orders tells us
864
                    # which delayed group we want. If none are provided, then
865
                    # we sum all the delayed groups together.
866
                    if orders[i]:
×
867
                        if orders[i] < len(shape[0]):
×
UNCOV
868
                            data[i, :] = temp_data[orders[i], :]
×
869
                    else:
870
                        data[i, :] = np.sum(temp_data[:, :], axis=0)
×
UNCOV
871
                elif shape == xsdata.xs_shapes["[DG][G][G']"]:
×
872
                    # Then we have a delayed group matrix. We will first
873
                    # remove the outgoing group dependency
UNCOV
874
                    temp_data = np.sum(temp_data, axis=-1)
×
875
                    # And then proceed in exactly the same manner as the
876
                    # "[DG][G']" or "[DG][G]" shapes in the previous block.
877
                    if orders[i]:
×
UNCOV
878
                        if orders[i] < len(shape[0]):
×
879
                            data[i, :] = temp_data[orders[i], :]
×
880
                    else:
UNCOV
881
                        data[i, :] = np.sum(temp_data[:, :], axis=0)
×
UNCOV
882
                elif shape == xsdata.xs_shapes["[G][G'][Order]"]:
×
883
                    # This is a scattering matrix with angular data
884
                    # First remove the outgoing group dependence
UNCOV
885
                    temp_data = np.sum(temp_data, axis=1)
×
886
                    # The user either provided a specific order or we resort
887
                    # to the default 0th order
888
                    if orders[i]:
×
UNCOV
889
                        order = orders[i]
×
890
                    else:
UNCOV
891
                        order = 0
×
892
                    # If the order is available, store the data for that order
893
                    # if it is not available, then the expansion coefficient
894
                    # is zero and thus we already have the correct value.
UNCOV
895
                    if order < shape[1]:
×
UNCOV
896
                        data[i, :] = temp_data[:, order]
×
897
    else:
UNCOV
898
        raise ValueError(f"{this} not present in provided MGXS library")
×
899

UNCOV
900
    return data
×
901

902

903
def _calculate_mgxs_elem_mat(this, types, library, orders=None,
12✔
904
                             temperature=294., ce_cross_sections=None,
905
                             enrichment=None):
906
    """Determines the multi-group cross sections of an element or material
907
    object.
908

909
    If the data for the nuclide or macroscopic object in the library is
910
    represented as angle-dependent data then this method will return the
911
    geometric average cross section over all angles.
912

913
    Parameters
914
    ----------
915
    this : str or openmc.Material
916
        Object to source data from. Elements can be input as a str
917
    types : Iterable of str
918
        The type of cross sections to calculate; values can either be those
919
        in openmc.PLOT_TYPES_MGXS
920
    library : openmc.MGXSLibrary
921
        MGXS Library containing the data of interest
922
    orders : Iterable of Integral, optional
923
        The scattering order or delayed group index to use for the
924
        corresponding entry in types. Defaults to the 0th order for scattering
925
        and the total delayed neutron data.
926
    temperature : float, optional
927
        Temperature in Kelvin to plot. If not specified, a default
928
        temperature of 294K will be plotted. Note that the nearest
929
        temperature in the library for each nuclide will be used as opposed
930
        to using any interpolation.
931
    ce_cross_sections : str, optional
932
        Location of continuous-energy cross_sections.xml file. Default is None.
933
        This is used only for expanding the elements
934
    enrichment : float, optional
935
        Enrichment for U235 in weight percent. For example, input 4.95 for
936
        4.95 weight percent enriched U. Default is None
937
        (natural composition).
938

939
    Returns
940
    -------
941
    data : numpy.ndarray
942
        Cross sections calculated at the energy grid described by energy_grid
943

944
    """
945

946
    if isinstance(this, openmc.Material):
×
UNCOV
947
        if this.temperature is not None:
×
UNCOV
948
            T = this.temperature
×
949
        else:
UNCOV
950
            T = temperature
×
951

952
        # Check to see if we have nuclides/elements or a macroscopic object
953
        if this._macroscopic is not None:
×
954
            # We have macroscopics
UNCOV
955
            nuclides = {this._macroscopic: this.density}
×
956
        else:
957
            # Expand elements in to nuclides with atomic densities
UNCOV
958
            nuclides = this.get_nuclide_atom_densities()
×
959

960
        # For ease of processing split out nuc and nuc_density
961
        nuc_fraction = list(nuclides.values())
×
962
    else:
UNCOV
963
        T = temperature
×
964
        # Expand elements in to nuclides with atomic densities
965
        nuclides = openmc.Element(this).expand(100., 'ao', enrichment=enrichment,
×
966
                               cross_sections=ce_cross_sections)
967

968
        # For ease of processing split out nuc and nuc_fractions
UNCOV
969
        nuc_fraction = [nuclide[1] for nuclide in nuclides]
×
970

971
    nuc_data = []
×
UNCOV
972
    for nuclide in nuclides.items():
×
973
        nuc_data.append(_calculate_mgxs_nuc_macro(nuclide[0], types, library,
×
974
                                                  orders, T))
975

976
    # Combine across the nuclides
UNCOV
977
    data = np.zeros((len(types), library.energy_groups.num_groups))
×
UNCOV
978
    for line in range(len(types)):
×
UNCOV
979
        if types[line] == 'unity':
×
UNCOV
980
            data[line, :] = 1.
×
981
        else:
UNCOV
982
            for n in range(len(nuclides)):
×
UNCOV
983
                data[line, :] += nuc_fraction[n] * nuc_data[n][line, :]
×
984

UNCOV
985
    return data
×
STATUS · Troubleshooting · Open an Issue · Sales · Support · CAREERS · ENTERPRISE · START FREE · SCHEDULE DEMO
ANNOUNCEMENTS · TWITTER · TOS & SLA · Supported CI Services · What's a CI service? · Automated Testing

© 2025 Coveralls, Inc