This notebook shows the how tallies can be combined (added, subtracted, multiplied, etc.) using the Python API in order to create derived tallies. Since no covariance information is obtained, it is assumed that tallies are completely independent of one another when propagating uncertainties. The target problem is a simple pin cell.


In [1]:
import glob

from IPython.display import Image
import numpy as np
import openmc

Generate Input Files

First we need to define materials that will be used in the problem. We'll create three materials for the fuel, water, and cladding of the fuel pin.


In [2]:
# 1.6 enriched fuel
fuel = openmc.Material(name='1.6% Fuel')
fuel.set_density('g/cm3', 10.31341)
fuel.add_nuclide('U235', 3.7503e-4)
fuel.add_nuclide('U238', 2.2625e-2)
fuel.add_nuclide('O16', 4.6007e-2)

# borated water
water = openmc.Material(name='Borated Water')
water.set_density('g/cm3', 0.740582)
water.add_nuclide('H1', 4.9457e-2)
water.add_nuclide('O16', 2.4732e-2)
water.add_nuclide('B10', 8.0042e-6)

# zircaloy
zircaloy = openmc.Material(name='Zircaloy')
zircaloy.set_density('g/cm3', 6.55)
zircaloy.add_nuclide('Zr90', 7.2758e-3)

With our three materials, we can now create a materials file object that can be exported to an actual XML file.


In [3]:
# Instantiate a Materials collection
materials_file = openmc.Materials([fuel, water, zircaloy])

# Export to "materials.xml"
materials_file.export_to_xml()

Now let's move on to the geometry. Our problem will have three regions for the fuel, the clad, and the surrounding coolant. The first step is to create the bounding surfaces -- in this case two cylinders and six planes.


In [4]:
# Create cylinders for the fuel and clad
fuel_outer_radius = openmc.ZCylinder(x0=0.0, y0=0.0, r=0.39218)
clad_outer_radius = openmc.ZCylinder(x0=0.0, y0=0.0, r=0.45720)

# Create boundary planes to surround the geometry
# Use both reflective and vacuum boundaries to make life interesting
min_x = openmc.XPlane(x0=-0.63, boundary_type='reflective')
max_x = openmc.XPlane(x0=+0.63, boundary_type='reflective')
min_y = openmc.YPlane(y0=-0.63, boundary_type='reflective')
max_y = openmc.YPlane(y0=+0.63, boundary_type='reflective')
min_z = openmc.ZPlane(z0=-100., boundary_type='vacuum')
max_z = openmc.ZPlane(z0=+100., boundary_type='vacuum')

With the surfaces defined, we can now create cells that are defined by intersections of half-spaces created by the surfaces.


In [5]:
# Create a Universe to encapsulate a fuel pin
pin_cell_universe = openmc.Universe(name='1.6% Fuel Pin')

# Create fuel Cell
fuel_cell = openmc.Cell(name='1.6% Fuel')
fuel_cell.fill = fuel
fuel_cell.region = -fuel_outer_radius
pin_cell_universe.add_cell(fuel_cell)

# Create a clad Cell
clad_cell = openmc.Cell(name='1.6% Clad')
clad_cell.fill = zircaloy
clad_cell.region = +fuel_outer_radius & -clad_outer_radius
pin_cell_universe.add_cell(clad_cell)

# Create a moderator Cell
moderator_cell = openmc.Cell(name='1.6% Moderator')
moderator_cell.fill = water
moderator_cell.region = +clad_outer_radius
pin_cell_universe.add_cell(moderator_cell)

OpenMC requires that there is a "root" universe. Let us create a root cell that is filled by the pin cell universe and then assign it to the root universe.


In [6]:
# Create root Cell
root_cell = openmc.Cell(name='root cell')
root_cell.fill = pin_cell_universe

# Add boundary planes
root_cell.region = +min_x & -max_x & +min_y & -max_y & +min_z & -max_z

# Create root Universe
root_universe = openmc.Universe(universe_id=0, name='root universe')
root_universe.add_cell(root_cell)

We now must create a geometry that is assigned a root universe, put the geometry into a geometry file, and export it to XML.


In [7]:
# Create Geometry and set root Universe
geometry = openmc.Geometry(root_universe)

In [8]:
# Export to "geometry.xml"
geometry.export_to_xml()

With the geometry and materials finished, we now just need to define simulation parameters. In this case, we will use 5 inactive batches and 15 active batches each with 2500 particles.


In [9]:
# OpenMC simulation parameters
batches = 20
inactive = 5
particles = 2500

# Instantiate a Settings object
settings_file = openmc.Settings()
settings_file.batches = batches
settings_file.inactive = inactive
settings_file.particles = particles
settings_file.output = {'tallies': True}

# Create an initial uniform spatial source distribution over fissionable zones
bounds = [-0.63, -0.63, -100., 0.63, 0.63, 100.]
uniform_dist = openmc.stats.Box(bounds[:3], bounds[3:], only_fissionable=True)
settings_file.source = openmc.Source(space=uniform_dist)

# Export to "settings.xml"
settings_file.export_to_xml()

Let us also create a plot file that we can use to verify that our pin cell geometry was created successfully.


In [10]:
# Instantiate a Plot
plot = openmc.Plot(plot_id=1)
plot.filename = 'materials-xy'
plot.origin = [0, 0, 0]
plot.width = [1.26, 1.26]
plot.pixels = [250, 250]
plot.color_by = 'material'

# Show plot
openmc.plot_inline(plot)


As we can see from the plot, we have a nice pin cell with fuel, cladding, and water! Before we run our simulation, we need to tell the code what we want to tally. The following code shows how to create a variety of tallies.


In [11]:
# Instantiate an empty Tallies object
tallies_file = openmc.Tallies()

In [12]:
# Create Tallies to compute microscopic multi-group cross-sections

# Instantiate energy filter for multi-group cross-section Tallies
energy_filter = openmc.EnergyFilter([0., 0.625, 20.0e6])

# Instantiate flux Tally in moderator and fuel
tally = openmc.Tally(name='flux')
tally.filters = [openmc.CellFilter([fuel_cell, moderator_cell])]
tally.filters.append(energy_filter)
tally.scores = ['flux']
tallies_file.append(tally)

# Instantiate reaction rate Tally in fuel
tally = openmc.Tally(name='fuel rxn rates')
tally.filters = [openmc.CellFilter(fuel_cell)]
tally.filters.append(energy_filter)
tally.scores = ['nu-fission', 'scatter']
tally.nuclides = ['U238', 'U235']
tallies_file.append(tally)

# Instantiate reaction rate Tally in moderator
tally = openmc.Tally(name='moderator rxn rates')
tally.filters = [openmc.CellFilter(moderator_cell)]
tally.filters.append(energy_filter)
tally.scores = ['absorption', 'total']
tally.nuclides = ['O16', 'H1']
tallies_file.append(tally)

# Instantiate a tally mesh
mesh = openmc.RegularMesh(mesh_id=1)
mesh.dimension = [1, 1, 1]
mesh.lower_left = [-0.63, -0.63, -100.]
mesh.width = [1.26, 1.26, 200.]
meshsurface_filter = openmc.MeshSurfaceFilter(mesh)

# Instantiate thermal, fast, and total leakage tallies
leak = openmc.Tally(name='leakage')
leak.filters = [meshsurface_filter]
leak.scores = ['current']
tallies_file.append(leak)

thermal_leak = openmc.Tally(name='thermal leakage')
thermal_leak.filters = [meshsurface_filter, openmc.EnergyFilter([0., 0.625])]
thermal_leak.scores = ['current']
tallies_file.append(thermal_leak)

fast_leak = openmc.Tally(name='fast leakage')
fast_leak.filters = [meshsurface_filter, openmc.EnergyFilter([0.625, 20.0e6])]
fast_leak.scores = ['current']
tallies_file.append(fast_leak)

In [13]:
# K-Eigenvalue (infinity) tallies
fiss_rate = openmc.Tally(name='fiss. rate')
abs_rate = openmc.Tally(name='abs. rate')
fiss_rate.scores = ['nu-fission']
abs_rate.scores = ['absorption']
tallies_file += (fiss_rate, abs_rate)

In [14]:
# Resonance Escape Probability tallies
therm_abs_rate = openmc.Tally(name='therm. abs. rate')
therm_abs_rate.scores = ['absorption']
therm_abs_rate.filters = [openmc.EnergyFilter([0., 0.625])]
tallies_file.append(therm_abs_rate)

In [15]:
# Thermal Flux Utilization tallies
fuel_therm_abs_rate = openmc.Tally(name='fuel therm. abs. rate')
fuel_therm_abs_rate.scores = ['absorption']
fuel_therm_abs_rate.filters = [openmc.EnergyFilter([0., 0.625]),
                               openmc.CellFilter([fuel_cell])]
tallies_file.append(fuel_therm_abs_rate)

In [16]:
# Fast Fission Factor tallies
therm_fiss_rate = openmc.Tally(name='therm. fiss. rate')
therm_fiss_rate.scores = ['nu-fission']
therm_fiss_rate.filters = [openmc.EnergyFilter([0., 0.625])]
tallies_file.append(therm_fiss_rate)

In [17]:
# Instantiate energy filter to illustrate Tally slicing
fine_energy_filter = openmc.EnergyFilter(np.logspace(np.log10(1e-2), np.log10(20.0e6), 10))

# Instantiate flux Tally in moderator and fuel
tally = openmc.Tally(name='need-to-slice')
tally.filters = [openmc.CellFilter([fuel_cell, moderator_cell])]
tally.filters.append(fine_energy_filter)
tally.scores = ['nu-fission', 'scatter']
tally.nuclides = ['H1', 'U238']
tallies_file.append(tally)

In [18]:
# Export to "tallies.xml"
tallies_file.export_to_xml()


/home/romano/openmc/openmc/mixin.py:71: IDWarning: Another Filter instance already exists with id=6.
  warn(msg, IDWarning)
/home/romano/openmc/openmc/mixin.py:71: IDWarning: Another Filter instance already exists with id=3.
  warn(msg, IDWarning)
/home/romano/openmc/openmc/mixin.py:71: IDWarning: Another Filter instance already exists with id=2.
  warn(msg, IDWarning)

Now we a have a complete set of inputs, so we can go ahead and run our simulation.


In [19]:
# Run OpenMC!
openmc.run()


                                %%%%%%%%%%%%%%%
                           %%%%%%%%%%%%%%%%%%%%%%%%
                        %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
                      %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
                    %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
                   %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
                                    %%%%%%%%%%%%%%%%%%%%%%%%
                                     %%%%%%%%%%%%%%%%%%%%%%%%
                 ###############      %%%%%%%%%%%%%%%%%%%%%%%%
                ##################     %%%%%%%%%%%%%%%%%%%%%%%
                ###################     %%%%%%%%%%%%%%%%%%%%%%%
                ####################     %%%%%%%%%%%%%%%%%%%%%%
                #####################     %%%%%%%%%%%%%%%%%%%%%
                ######################     %%%%%%%%%%%%%%%%%%%%
                #######################     %%%%%%%%%%%%%%%%%%
                 #######################     %%%%%%%%%%%%%%%%%
                 ######################     %%%%%%%%%%%%%%%%%
                  ####################     %%%%%%%%%%%%%%%%%
                    #################     %%%%%%%%%%%%%%%%%
                     ###############     %%%%%%%%%%%%%%%%
                       ############     %%%%%%%%%%%%%%%
                          ########     %%%%%%%%%%%%%%
                                      %%%%%%%%%%%

                   | The OpenMC Monte Carlo Code
         Copyright | 2011-2019 MIT and OpenMC contributors
           License | http://openmc.readthedocs.io/en/latest/license.html
           Version | 0.11.0-dev
          Git SHA1 | 61c911cffdae2406f9f4bc667a9a6954748bb70c
         Date/Time | 2019-07-18 22:51:02
    OpenMP Threads | 4

 Reading settings XML file...
 Reading cross sections XML file...
 Reading materials XML file...
 Reading geometry XML file...
 Reading U235 from /opt/data/hdf5/nndc_hdf5_v15/U235.h5
 Reading U238 from /opt/data/hdf5/nndc_hdf5_v15/U238.h5
 Reading O16 from /opt/data/hdf5/nndc_hdf5_v15/O16.h5
 Reading H1 from /opt/data/hdf5/nndc_hdf5_v15/H1.h5
 Reading B10 from /opt/data/hdf5/nndc_hdf5_v15/B10.h5
 Reading Zr90 from /opt/data/hdf5/nndc_hdf5_v15/Zr90.h5
 Maximum neutron transport energy: 20000000.000000 eV for U235
 Reading tallies XML file...
 Writing summary.h5 file...
 Initializing source particles...

 ====================>     K EIGENVALUE SIMULATION     <====================

  Bat./Gen.      k            Average k
  =========   ========   ====================
        1/1    0.96168
        2/1    0.96651
        3/1    1.00678
        4/1    0.98773
        5/1    1.01883
        6/1    1.02959
        7/1    0.99859    1.01409 +/- 0.01550
        8/1    1.03441    1.02086 +/- 0.01123
        9/1    1.06097    1.03089 +/- 0.01279
       10/1    1.06132    1.03698 +/- 0.01163
       11/1    1.04687    1.03863 +/- 0.00964
       12/1    1.02982    1.03737 +/- 0.00824
       13/1    1.03520    1.03710 +/- 0.00714
       14/1    0.99508    1.03243 +/- 0.00784
       15/1    1.03987    1.03317 +/- 0.00705
       16/1    1.02743    1.03265 +/- 0.00640
       17/1    1.02975    1.03241 +/- 0.00585
       18/1    0.99671    1.02966 +/- 0.00604
       19/1    1.02040    1.02900 +/- 0.00563
       20/1    1.02024    1.02842 +/- 0.00527
 Creating state point statepoint.20.h5...

 =======================>     TIMING STATISTICS     <=======================

 Total time for initialization     = 3.4427e-01 seconds
   Reading cross sections          = 3.1628e-01 seconds
 Total time in simulation          = 3.7319e+00 seconds
   Time in transport only          = 3.6302e+00 seconds
   Time in inactive batches        = 4.9601e-01 seconds
   Time in active batches          = 3.2359e+00 seconds
   Time synchronizing fission bank = 2.8100e-03 seconds
     Sampling source sites         = 2.4682e-03 seconds
     SEND/RECV source sites        = 3.2484e-04 seconds
   Time accumulating tallies       = 4.4538e-05 seconds
 Total time for finalization       = 9.3656e-04 seconds
 Total time elapsed                = 4.0859e+00 seconds
 Calculation Rate (inactive)       = 25201.2 particles/second
 Calculation Rate (active)         = 11588.7 particles/second

 ============================>     RESULTS     <============================

 k-effective (Collision)     = 1.02889 +/- 0.00492
 k-effective (Track-length)  = 1.02842 +/- 0.00527
 k-effective (Absorption)    = 1.02637 +/- 0.00349
 Combined k-effective        = 1.02700 +/- 0.00291
 Leakage Fraction            = 0.01717 +/- 0.00107

Tally Data Processing

Our simulation ran successfully and created a statepoint file with all the tally data in it. We begin our analysis here loading the statepoint file and 'reading' the results. By default, the tally results are not read into memory because they might be large, even large enough to exceed the available memory on a computer.


In [20]:
# Load the statepoint file
sp = openmc.StatePoint('statepoint.20.h5')

We have a tally of the total fission rate and the total absorption rate, so we can calculate k-eff as: $$k_{eff} = \frac{\langle \nu \Sigma_f \phi \rangle}{\langle \Sigma_a \phi \rangle + \langle L \rangle}$$ In this notation, $\langle \cdot \rangle^a_b$ represents an OpenMC that is integrated over region $a$ and energy range $b$. If $a$ or $b$ is not reported, it means the value represents an integral over all space or all energy, respectively.


In [21]:
# Get the fission and absorption rate tallies
fiss_rate = sp.get_tally(name='fiss. rate')
abs_rate = sp.get_tally(name='abs. rate')

# Get the leakage tally
leak = sp.get_tally(name='leakage')
leak = leak.summation(filter_type=openmc.MeshSurfaceFilter, remove_filter=True)

# Compute k-infinity using tally arithmetic
keff = fiss_rate / (abs_rate + leak)
keff.get_pandas_dataframe()


Out[21]:
nuclide score mean std. dev.
0 total (nu-fission / (absorption + current)) 1.023002 0.006647

Notice that even though the neutron production rate, absorption rate, and current are separate tallies, we still get a first-order estimate of the uncertainty on the quotient of them automatically!

Often in textbooks you'll see k-eff represented using the six-factor formula $$k_{eff} = p \epsilon f \eta P_{FNL} P_{TNL}.$$ Let's analyze each of these factors, starting with the resonance escape probability which is defined as $$p=\frac{\langle\Sigma_a\phi\rangle_T + \langle L \rangle_T}{\langle\Sigma_a\phi\rangle + \langle L \rangle_T}$$ where the subscript $T$ means thermal energies.


In [22]:
# Compute resonance escape probability using tally arithmetic
therm_abs_rate = sp.get_tally(name='therm. abs. rate')
thermal_leak = sp.get_tally(name='thermal leakage')
thermal_leak = thermal_leak.summation(filter_type=openmc.MeshSurfaceFilter, remove_filter=True)
res_esc = (therm_abs_rate + thermal_leak) / (abs_rate + thermal_leak)
res_esc.get_pandas_dataframe()


Out[22]:
energy low [eV] energy high [eV] nuclide score mean std. dev.
0 0.0 0.625 total ((absorption + current) / (absorption + current)) 0.694368 0.004606

The fast fission factor can be calculated as $$\epsilon=\frac{\langle\nu\Sigma_f\phi\rangle}{\langle\nu\Sigma_f\phi\rangle_T}$$


In [23]:
# Compute fast fission factor factor using tally arithmetic
therm_fiss_rate = sp.get_tally(name='therm. fiss. rate')
fast_fiss = fiss_rate / therm_fiss_rate
fast_fiss.get_pandas_dataframe()


Out[23]:
energy low [eV] energy high [eV] nuclide score mean std. dev.
0 0.0 0.625 total (nu-fission / nu-fission) 1.203099 0.009615

The thermal flux utilization is calculated as $$f=\frac{\langle\Sigma_a\phi\rangle^F_T}{\langle\Sigma_a\phi\rangle_T}$$ where the superscript $F$ denotes fuel.


In [24]:
# Compute thermal flux utilization factor using tally arithmetic
fuel_therm_abs_rate = sp.get_tally(name='fuel therm. abs. rate')
therm_util = fuel_therm_abs_rate / therm_abs_rate
therm_util.get_pandas_dataframe()


Out[24]:
energy low [eV] energy high [eV] cell nuclide score mean std. dev.
0 0.0 0.625 1 total (absorption / absorption) 0.749423 0.006089

The next factor is the number of fission neutrons produced per absorption in fuel, calculated as $$\eta = \frac{\langle \nu\Sigma_f\phi \rangle_T}{\langle \Sigma_a \phi \rangle^F_T}$$


In [25]:
# Compute neutrons produced per absorption (eta) using tally arithmetic
eta = therm_fiss_rate / fuel_therm_abs_rate
eta.get_pandas_dataframe()


Out[25]:
energy low [eV] energy high [eV] cell nuclide score mean std. dev.
0 0.0 0.625 1 total (nu-fission / absorption) 1.663727 0.014403

There are two leakage factors to account for fast and thermal leakage. The fast non-leakage probability is computed as $$P_{FNL} = \frac{\langle \Sigma_a\phi \rangle + \langle L \rangle_T}{\langle \Sigma_a \phi \rangle + \langle L \rangle}$$


In [26]:
p_fnl = (abs_rate + thermal_leak) / (abs_rate + leak)
p_fnl.get_pandas_dataframe()


Out[26]:
energy low [eV] energy high [eV] nuclide score mean std. dev.
0 0.0 0.625 total ((absorption + current) / (absorption + current)) 0.984668 0.005509

The final factor is the thermal non-leakage probability and is computed as $$P_{TNL} = \frac{\langle \Sigma_a\phi \rangle_T}{\langle \Sigma_a \phi \rangle_T + \langle L \rangle_T}$$


In [27]:
p_tnl = therm_abs_rate / (therm_abs_rate + thermal_leak)
p_tnl.get_pandas_dataframe()


Out[27]:
energy low [eV] energy high [eV] nuclide score mean std. dev.
0 0.0 0.625 total (absorption / (absorption + current)) 0.997439 0.007548

Now we can calculate $k_{eff}$ using the product of the factors form the four-factor formula.


In [28]:
keff = res_esc * fast_fiss * therm_util * eta * p_fnl * p_tnl
keff.get_pandas_dataframe()


Out[28]:
energy low [eV] energy high [eV] cell nuclide score mean std. dev.
0 0.0 0.625 1 total (((((((absorption + current) / (absorption + c... 1.023002 0.018791

We see that the value we've obtained here has exactly the same mean as before. However, because of the way it was calculated, the standard deviation appears to be larger.

Let's move on to a more complicated example now. Before we set up tallies to get reaction rates in the fuel and moderator in two energy groups for two different nuclides. We can use tally arithmetic to divide each of these reaction rates by the flux to get microscopic multi-group cross sections.


In [29]:
# Compute microscopic multi-group cross-sections
flux = sp.get_tally(name='flux')
flux = flux.get_slice(filters=[openmc.CellFilter], filter_bins=[(fuel_cell.id,)])
fuel_rxn_rates = sp.get_tally(name='fuel rxn rates')
mod_rxn_rates = sp.get_tally(name='moderator rxn rates')

In [30]:
fuel_xs = fuel_rxn_rates / flux
fuel_xs.get_pandas_dataframe()


Out[30]:
cell energy low [eV] energy high [eV] nuclide score mean std. dev.
0 1 0.000 6.250000e-01 (U238 / total) (nu-fission / flux) 6.659486e-07 5.627975e-09
1 1 0.000 6.250000e-01 (U238 / total) (scatter / flux) 2.099901e-01 1.748379e-03
2 1 0.000 6.250000e-01 (U235 / total) (nu-fission / flux) 3.566329e-01 3.030782e-03
3 1 0.000 6.250000e-01 (U235 / total) (scatter / flux) 5.555466e-03 4.635318e-05
4 1 0.625 2.000000e+07 (U238 / total) (nu-fission / flux) 7.251304e-03 5.161998e-05
5 1 0.625 2.000000e+07 (U238 / total) (scatter / flux) 2.272661e-01 9.576939e-04
6 1 0.625 2.000000e+07 (U235 / total) (nu-fission / flux) 7.920169e-03 5.751231e-05
7 1 0.625 2.000000e+07 (U235 / total) (scatter / flux) 3.358280e-03 1.341281e-05

We see that when the two tallies with multiple bins were divided, the derived tally contains the outer product of the combinations. If the filters/scores are the same, no outer product is needed. The get_values(...) method allows us to obtain a subset of tally scores. In the following example, we obtain just the neutron production microscopic cross sections.


In [31]:
# Show how to use Tally.get_values(...) with a CrossScore
nu_fiss_xs = fuel_xs.get_values(scores=['(nu-fission / flux)'])
print(nu_fiss_xs)


[[[6.65948580e-07]
  [3.56632881e-01]]

 [[7.25130446e-03]
  [7.92016892e-03]]]

The same idea can be used not only for scores but also for filters and nuclides.


In [32]:
# Show how to use Tally.get_values(...) with a CrossScore and CrossNuclide
u235_scatter_xs = fuel_xs.get_values(nuclides=['(U235 / total)'], 
                                scores=['(scatter / flux)'])
print(u235_scatter_xs)


[[[0.00555547]]

 [[0.00335828]]]

In [33]:
# Show how to use Tally.get_values(...) with a CrossFilter and CrossScore
fast_scatter_xs = fuel_xs.get_values(filters=[openmc.EnergyFilter], 
                                     filter_bins=[((0.625, 20.0e6),)], 
                                     scores=['(scatter / flux)'])
print(fast_scatter_xs)


[[[0.22726611]
  [0.00335828]]]

A more advanced method is to use get_slice(...) to create a new derived tally that is a subset of an existing tally. This has the benefit that we can use get_pandas_dataframe() to see the tallies in a more human-readable format.


In [34]:
# "Slice" the nu-fission data into a new derived Tally
nu_fission_rates = fuel_rxn_rates.get_slice(scores=['nu-fission'])
nu_fission_rates.get_pandas_dataframe()


Out[34]:
cell energy low [eV] energy high [eV] nuclide score mean std. dev.
0 1 0.000 6.250000e-01 U238 nu-fission 0.000002 9.679304e-09
1 1 0.000 6.250000e-01 U235 nu-fission 0.854805 5.239673e-03
2 1 0.625 2.000000e+07 U238 nu-fission 0.082978 5.346135e-04
3 1 0.625 2.000000e+07 U235 nu-fission 0.090632 5.981942e-04

In [35]:
# "Slice" the H-1 scatter data in the moderator Cell into a new derived Tally
need_to_slice = sp.get_tally(name='need-to-slice')
slice_test = need_to_slice.get_slice(scores=['scatter'], nuclides=['H1'],
                                     filters=[openmc.CellFilter], filter_bins=[(moderator_cell.id,)])
slice_test.get_pandas_dataframe()


Out[35]:
cell energy low [eV] energy high [eV] nuclide score mean std. dev.
0 3 1.000000e-02 1.080060e-01 H1 scatter 4.541188 0.025230
1 3 1.080060e-01 1.166529e+00 H1 scatter 2.001332 0.006754
2 3 1.166529e+00 1.259921e+01 H1 scatter 1.639292 0.011374
3 3 1.259921e+01 1.360790e+02 H1 scatter 1.821633 0.009590
4 3 1.360790e+02 1.469734e+03 H1 scatter 2.032395 0.009953
5 3 1.469734e+03 1.587401e+04 H1 scatter 2.120745 0.011090
6 3 1.587401e+04 1.714488e+05 H1 scatter 2.181709 0.013602
7 3 1.714488e+05 1.851749e+06 H1 scatter 2.013644 0.009219
8 3 1.851749e+06 2.000000e+07 H1 scatter 0.372640 0.002903