This IPython Notebook illustrates the use of the openmc.mgxs.Library class. The Library class is designed to automate the calculation of multi-group cross sections for use cases with one or more domains, cross section types, and/or nuclides. In particular, this Notebook illustrates the following features:

  • Calculation of multi-group cross sections for a fuel assembly
  • Automated creation, manipulation and storage of MGXS with openmc.mgxs.Library
  • Validation of multi-group cross sections with OpenMOC
  • Steady-state pin-by-pin fission rates comparison between OpenMC and OpenMOC

Note: This Notebook was created using OpenMOC to verify the multi-group cross-sections generated by OpenMC. In order to run this Notebook in its entirety, you must have OpenMOC installed on your system, along with OpenCG to convert the OpenMC geometries into OpenMOC geometries. In addition, this Notebook illustrates the use of Pandas DataFrames to containerize multi-group cross section data. We recommend using Pandas >v0.15.0 or later since OpenMC's Python API leverages the multi-indexing feature included in the most recent releases of Pandas.

Generate Input Files


In [1]:
import math
import pickle
from IPython.display import Image
import matplotlib.pylab as pylab
import numpy as np

import openmc
import openmc.mgxs
from openmc.statepoint import StatePoint
from openmc.summary import Summary
from openmc.source import Source
from openmc.stats import Box

import openmoc
import openmoc.process
from openmoc.compatible import get_openmoc_geometry
from openmoc.materialize import load_openmc_mgxs_lib

%matplotlib inline


/usr/local/lib/python2.7/dist-packages/matplotlib/__init__.py:1318: UserWarning:  This call to matplotlib.use() has no effect
because the backend has already been chosen;
matplotlib.use() must be called *before* pylab, matplotlib.pyplot,
or matplotlib.backends is imported for the first time.

  warnings.warn(_use_error_msg)

First we need to define materials that will be used in the problem. Before defining a material, we must create nuclides that are used in the material.


In [2]:
# Instantiate some Nuclides
h1 = openmc.Nuclide('H-1')
b10 = openmc.Nuclide('B-10')
o16 = openmc.Nuclide('O-16')
u235 = openmc.Nuclide('U-235')
u238 = openmc.Nuclide('U-238')
zr90 = openmc.Nuclide('Zr-90')

With the nuclides we defined, we will now create three materials for the fuel, water, and cladding of the fuel pins.


In [3]:
# 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 MaterialsFile object that can be exported to an actual XML file.


In [4]:
# Instantiate a MaterialsFile, add Materials
materials_file = openmc.MaterialsFile()
materials_file.add_material(fuel)
materials_file.add_material(water)
materials_file.add_material(zircaloy)
materials_file.default_xs = '71c'

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

Now let's move on to the geometry. This problem will be a square array of fuel pins and control rod guide tubes for which we can use OpenMC's lattice/universe feature. The basic universe will have three regions for the fuel, the clad, and the surrounding coolant. The first step is to create the bounding surfaces for fuel and clad, as well as the outer bounding surfaces of the problem.


In [5]:
# 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
min_x = openmc.XPlane(x0=-10.71, boundary_type='reflective')
max_x = openmc.XPlane(x0=+10.71, boundary_type='reflective')
min_y = openmc.YPlane(y0=-10.71, boundary_type='reflective')
max_y = openmc.YPlane(y0=+10.71, boundary_type='reflective')
min_z = openmc.ZPlane(z0=-10., boundary_type='reflective')
max_z = openmc.ZPlane(z0=+10., boundary_type='reflective')

With the surfaces defined, we can now construct a fuel pin cell from cells that are defined by intersections of half-spaces created by the surfaces.


In [6]:
# Create a Universe to encapsulate a fuel pin
fuel_pin_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
fuel_pin_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
fuel_pin_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
fuel_pin_universe.add_cell(moderator_cell)

Likewise, we can construct a control rod guide tube with the same surfaces.


In [7]:
# Create a Universe to encapsulate a control rod guide tube
guide_tube_universe = openmc.Universe(name='Guide Tube')

# Create guide tube Cell
guide_tube_cell = openmc.Cell(name='Guide Tube Water')
guide_tube_cell.fill = water
guide_tube_cell.region = -fuel_outer_radius
guide_tube_universe.add_cell(guide_tube_cell)

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

# Create a moderator Cell
moderator_cell = openmc.Cell(name='Guide Tube Moderator')
moderator_cell.fill = water
moderator_cell.region = +clad_outer_radius
guide_tube_universe.add_cell(moderator_cell)

Using the pin cell universe, we can construct a 17x17 rectangular lattice with a 1.26 cm pitch.


In [8]:
# Create fuel assembly Lattice
assembly = openmc.RectLattice(name='1.6% Fuel Assembly')
assembly.dimension = (17, 17)
assembly.pitch = (1.26, 1.26)
assembly.lower_left = [-1.26 * 17. / 2.0] * 2

Next, we create a NumPy array of fuel pin and guide tube universes for the lattice.


In [9]:
# Create array indices for guide tube locations in lattice
template_x = np.array([5, 8, 11, 3, 13, 2, 5, 8, 11, 14, 2, 5, 8,
                       11, 14, 2, 5, 8, 11, 14, 3, 13, 5, 8, 11])
template_y = np.array([2, 2, 2, 3, 3, 5, 5, 5, 5, 5, 8, 8, 8, 8,
                       8, 11, 11, 11, 11, 11, 13, 13, 14, 14, 14])

# Initialize an empty 17x17 array of the lattice universes
universes = np.empty((17, 17), dtype=openmc.Universe)

# Fill the array with the fuel pin and guide tube universes
universes[:,:] = fuel_pin_universe
universes[template_x, template_y] = guide_tube_universe

# Store the array of universes in the lattice
assembly.universes = universes

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 [10]:
# Create root Cell
root_cell = openmc.Cell(name='root cell')
root_cell.fill = assembly

# 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 GeometryFile object, and export it to XML.


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

In [12]:
# Instantiate a GeometryFile
geometry_file = openmc.GeometryFile()
geometry_file.geometry = geometry

# Export to "geometry.xml"
geometry_file.export_to_xml()

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


In [13]:
# OpenMC simulation parameters
batches = 50
inactive = 10
particles = 2500

# Instantiate a SettingsFile
settings_file = openmc.SettingsFile()
settings_file.batches = batches
settings_file.inactive = inactive
settings_file.particles = particles
settings_file.output = {'tallies': False, 'summary': True}
source_bounds = [-10.71, -10.71, -10, 10.71, 10.71, 10.]
settings_file.source = Source(Box(
        source_bounds[:3], source_bounds[3:], only_fissionable=True))

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

Let us also create a PlotsFile that we can use to verify that our fuel assembly geometry was created successfully.


In [14]:
# Instantiate a Plot
plot = openmc.Plot(plot_id=1)
plot.filename = 'materials-xy'
plot.origin = [0, 0, 0]
plot.width = [21.5, 21.5]
plot.pixels = [250, 250]
plot.color = 'mat'

# Instantiate a PlotsFile, add Plot, and export to "plots.xml"
plot_file = openmc.PlotsFile()
plot_file.add_plot(plot)
plot_file.export_to_xml()

With the plots.xml file, we can now generate and view the plot. OpenMC outputs plots in .ppm format, which can be converted into a compressed format like .png with the convert utility.


In [15]:
# Run openmc in plotting mode
executor = openmc.Executor()
executor.plot_geometry(output=False)


Out[15]:
0

In [16]:
# Convert OpenMC's funky ppm to png
!convert materials-xy.ppm materials-xy.png

# Display the materials plot inline
Image(filename='materials-xy.png')


Out[16]:

As we can see from the plot, we have a nice array of fuel and guide tube pin cells with fuel, cladding, and water!

Create an MGXS Library

Now we are ready to generate multi-group cross sections! First, let's define a 2-group structure using the built-in EnergyGroups class.


In [17]:
# Instantiate a 2-group EnergyGroups object
groups = openmc.mgxs.EnergyGroups()
groups.group_edges = np.array([0., 0.625e-6, 20.])

Next, we will instantiate an openmc.mgxs.Library for the energy groups with our the fuel assembly geometry.


In [18]:
# Initialize an 2-group MGXS Library for OpenMOC
mgxs_lib = openmc.mgxs.Library(geometry)
mgxs_lib.energy_groups = groups

Now, we must specify to the Library which types of cross sections to compute. In particular, the following are the multi-group cross section MGXS subclasses that are mapped to string codes accepted by the Library class:

  • TotalXS ("total")
  • TransportXS ("transport")
  • AbsorptionXS ("absorption")
  • CaptureXS ("capture")
  • FissionXS ("fission")
  • NuFissionXS ("nu-fission")
  • ScatterXS ("scatter")
  • NuScatterXS ("nu-scatter")
  • ScatterMatrixXS ("scatter matrix")
  • NuScatterMatrixXS ("nu-scatter matrix")
  • Chi ("chi")

In this case, let's create the multi-group cross sections needed to run an OpenMOC simulation to verify the accuracy of our cross sections. In particular, we will define "transport", "nu-fission", "nu-scatter matrix" and "chi" cross sections for our Library.

Note: A variety of different approximate transport-corrected total multi-group cross sections (and corresponding scattering matrices) can be found in the literature. At the present time, the openmc.mgxs module only supports the "P0" transport correction. This correction can be turned on and off through the boolean Library.correction property which may take values of "P0" (default) or None.


In [19]:
# Specify multi-group cross section types to compute
mgxs_lib.mgxs_types = ['transport', 'nu-fission', 'nu-scatter matrix', 'chi']

Now we must specify the type of domain over which we would like the Library to compute multi-group cross sections. The domain type corresponds to the type of tally filter to be used in the tallies created to compute multi-group cross sections. At the present time, the Library supports "material," "cell," and "universe" domain types. We will use a "cell" domain type here to compute cross sections in each of the cells in the fuel assembly geometry.

Note: By default, the Library class will instantiate MGXS objects for each and every domain (material, cell or universe) in the geometry of interest. However, one may specify a subset of these domains to the Library.domains property. In our case, we wish to compute multi-group cross sections in each and every cell since they will be needed in our downstream OpenMOC calculation on the identical combinatorial geometry mesh.


In [20]:
# Specify a "cell" domain type for the cross section tally filters
mgxs_lib.domain_type = "cell"

# Specify the cell domains over which to compute multi-group cross sections
mgxs_lib.domains = geometry.get_all_material_cells()

We can easily instruct the Library to compute multi-group cross sections on a nuclide-by-nuclide basis with the boolean Library.by_nuclide property. By default, by_nuclide is set to False, but we will set it to True here.


In [21]:
# Compute cross sections on a nuclide-by-nuclide basis
mgxs_lib.by_nuclide = True

Lastly, we use the Library to construct the tallies needed to compute all of the requested multi-group cross sections in each domain and nuclide.


In [22]:
# Construct all tallies needed for the multi-group cross section library
mgxs_lib.build_library()

The tallies can now be export to a "tallies.xml" input file for OpenMC.

NOTE: At this point the Library has constructed nearly 100 distinct Tally objects. The overhead to tally in OpenMC scales as $O(N)$ for $N$ tallies, which can become a bottleneck for large tally datasets. To compensate for this, the Python API's Tally, Filter and TalliesFile classes allow for the smart merging of tallies when possible. The Library class supports this runtime optimization with the use of the optional merge paramter (False by default) for the Library.add_to_tallies_file(...) method, as shown below.


In [23]:
# Create a "tallies.xml" file for the MGXS Library
tallies_file = openmc.TalliesFile()
mgxs_lib.add_to_tallies_file(tallies_file, merge=True)

In addition, we instantiate a fission rate mesh tally to compare with OpenMOC.


In [24]:
# Instantiate a tally Mesh
mesh = openmc.Mesh(mesh_id=1)
mesh.type = 'regular'
mesh.dimension = [17, 17]
mesh.lower_left = [-10.71, -10.71]
mesh.width = [1.26, 1.26]

# Instantiate tally Filter
mesh_filter = openmc.Filter()
mesh_filter.mesh = mesh

# Instantiate the Tally
tally = openmc.Tally(name='mesh tally')
tally.add_filter(mesh_filter)
tally.add_score('fission')
tally.add_score('nu-fission')

# Add mesh and Tally to TalliesFile
tallies_file.add_mesh(mesh)
tallies_file.add_tally(tally)

In [25]:
# Export all tallies to a "tallies.xml" file
tallies_file.export_to_xml()

In [26]:
# Run OpenMC
executor.run_simulation()


       .d88888b.                             888b     d888  .d8888b.
      d88P" "Y88b                            8888b   d8888 d88P  Y88b
      888     888                            88888b.d88888 888    888
      888     888 88888b.   .d88b.  88888b.  888Y88888P888 888       
      888     888 888 "88b d8P  Y8b 888 "88b 888 Y888P 888 888       
      888     888 888  888 88888888 888  888 888  Y8P  888 888    888
      Y88b. .d88P 888 d88P Y8b.     888  888 888   "   888 Y88b  d88P
       "Y88888P"  88888P"   "Y8888  888  888 888       888  "Y8888P"
__________________888______________________________________________________
                  888
                  888

      Copyright:      2011-2015 Massachusetts Institute of Technology
      License:        http://mit-crpg.github.io/openmc/license.html
      Version:        0.7.1
      Git SHA1:       ea9fb637f63f9374c7436456141afa850b84acf9
      Date/Time:      2016-01-14 08:12:09

 ===========================================================================
 ========================>     INITIALIZATION     <=========================
 ===========================================================================

 Reading settings XML file...
 Reading cross sections XML file...
 Reading geometry XML file...
 Reading materials XML file...
 Reading tallies XML file...
 Building neighboring cells lists for each surface...
 Loading ACE cross section table: 92235.71c
 Loading ACE cross section table: 92238.71c
 Loading ACE cross section table: 8016.71c
 Loading ACE cross section table: 1001.71c
 Loading ACE cross section table: 5010.71c
 Loading ACE cross section table: 40090.71c
 Maximum neutron transport energy: 20.0000 MeV for 92235.71c
 Initializing source particles...

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

  Bat./Gen.      k            Average k         
  =========   ========   ====================   
        1/1    1.02650                       
        2/1    1.01386                       
        3/1    1.01045                       
        4/1    1.05511                       
        5/1    1.04873                       
        6/1    1.04558                       
        7/1    1.03840                       
        8/1    1.02086                       
        9/1    1.08845                       
       10/1    1.03932                       
       11/1    1.01271                       
       12/1    1.03448    1.02360 +/- 0.01088
       13/1    1.04395    1.03038 +/- 0.00925
       14/1    1.05477    1.03648 +/- 0.00894
       15/1    1.00485    1.03015 +/- 0.00938
       16/1    1.04523    1.03267 +/- 0.00806
       17/1    1.01328    1.02990 +/- 0.00735
       18/1    1.01476    1.02800 +/- 0.00664
       19/1    1.01490    1.02655 +/- 0.00604
       20/1    1.00926    1.02482 +/- 0.00567
       21/1    0.98504    1.02120 +/- 0.00627
       22/1    1.00397    1.01977 +/- 0.00591
       23/1    1.02556    1.02021 +/- 0.00545
       24/1    0.99808    1.01863 +/- 0.00529
       25/1    0.99638    1.01715 +/- 0.00514
       26/1    0.99615    1.01584 +/- 0.00499
       27/1    1.01843    1.01599 +/- 0.00469
       28/1    1.00315    1.01528 +/- 0.00447
       29/1    1.00633    1.01480 +/- 0.00426
       30/1    1.02159    1.01514 +/- 0.00405
       31/1    1.03395    1.01604 +/- 0.00396
       32/1    1.02672    1.01652 +/- 0.00381
       33/1    1.03778    1.01745 +/- 0.00375
       34/1    1.03807    1.01831 +/- 0.00369
       35/1    1.07854    1.02072 +/- 0.00428
       36/1    1.03524    1.02128 +/- 0.00415
       37/1    1.03100    1.02164 +/- 0.00401
       38/1    1.03853    1.02224 +/- 0.00391
       39/1    1.04089    1.02288 +/- 0.00383
       40/1    1.02150    1.02284 +/- 0.00370
       41/1    0.98470    1.02161 +/- 0.00379
       42/1    1.00658    1.02114 +/- 0.00370
       43/1    0.98652    1.02009 +/- 0.00373
       44/1    1.02787    1.02032 +/- 0.00363
       45/1    0.98800    1.01939 +/- 0.00364
       46/1    1.00286    1.01893 +/- 0.00357
       47/1    1.02559    1.01911 +/- 0.00348
       48/1    1.03729    1.01959 +/- 0.00342
       49/1    1.02538    1.01974 +/- 0.00333
       50/1    1.01478    1.01962 +/- 0.00325
 Creating state point statepoint.50.h5...

 ===========================================================================
 ======================>     SIMULATION FINISHED     <======================
 ===========================================================================


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

 Total time for initialization     =  4.1800E-01 seconds
   Reading cross sections          =  1.4300E-01 seconds
 Total time in simulation          =  4.1206E+01 seconds
   Time in transport only          =  4.1193E+01 seconds
   Time in inactive batches        =  4.1760E+00 seconds
   Time in active batches          =  3.7030E+01 seconds
   Time synchronizing fission bank =  3.0000E-03 seconds
     Sampling source sites         =  2.0000E-03 seconds
     SEND/RECV source sites        =  1.0000E-03 seconds
   Time accumulating tallies       =  0.0000E+00 seconds
 Total time for finalization       =  0.0000E+00 seconds
 Total time elapsed                =  4.1648E+01 seconds
 Calculation Rate (inactive)       =  5986.59 neutrons/second
 Calculation Rate (active)         =  2700.51 neutrons/second

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

 k-effective (Collision)     =  1.01805 +/-  0.00261
 k-effective (Track-length)  =  1.01962 +/-  0.00325
 k-effective (Absorption)    =  1.01554 +/-  0.00339
 Combined k-effective        =  1.01711 +/-  0.00235
 Leakage Fraction            =  0.00000 +/-  0.00000

Out[26]:
0

Tally Data Processing

Our simulation ran successfully and created statepoint and summary output files. We begin our analysis by instantiating a StatePoint object.


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

In addition to the statepoint file, our simulation also created a summary file which encapsulates information about the materials and geometry. This is necessary for the openmc.mgxs module to properly process the tally data. We first create a Summary object and link it with the statepoint.


In [28]:
su = openmc.Summary('summary.h5')
sp.link_with_summary(su)

The statepoint is now ready to be analyzed by the Library. We simply have to load the tallies from the statepoint into the Library and our MGXS objects will compute the cross sections for us under-the-hood.


In [29]:
# Initialize MGXS Library with OpenMC statepoint data
mgxs_lib.load_from_statepoint(sp)

Voila! Our multi-group cross sections are now ready to rock 'n roll!

Extracting and Storing MGXS Data

The Library supports a rich API to automate a variety of tasks, including multi-group cross section data retrieval and storage. We will highlight a few of these features here. First, the Library.get_mgxs(...) method allows one to extract an MGXS object from the Library for a particular domain and cross section type. The following cell illustrates how one may extract the NuFissionXS object for the fuel cell.

Note: The MGXS.get_mgxs(...) method will accept either the domain or the integer domain ID of interest.


In [30]:
# Retrieve the NuFissionXS object for the fuel cell from the library
fuel_mgxs = mgxs_lib.get_mgxs(fuel_cell, 'nu-fission')

The NuFissionXS object supports all of the methods described previously the openmc.mgxs tutorials, such as Pandas DataFrames:


In [31]:
df = fuel_mgxs.get_pandas_dataframe()
df


/home/romano/openmc/openmc/tallies.py:1642: RuntimeWarning: invalid value encountered in true_divide
  self_rel_err = data['self']['std. dev.'] / data['self']['mean']
Out[31]:
cell group in nuclide mean std. dev.
3 10000 1 U-235 0.008064 4.062984e-05
4 10000 1 U-238 0.007336 4.459335e-05
5 10000 1 O-16 0.000000 0.000000e+00
0 10000 2 U-235 0.361327 1.902492e-03
1 10000 2 U-238 0.000001 3.536787e-09
2 10000 2 O-16 0.000000 0.000000e+00

Similarly, we can use the MGXS.print_xs(...) method to view a string representation of the multi-group cross section data.


In [32]:
fuel_mgxs.print_xs()


Multi-Group XS
	Reaction Type  =	nu-fission
	Domain Type    =	cell
	Domain ID      =	10000
	Nuclide        =	U-235
	Cross Sections [cm^-1]:
            Group 1 [6.25e-07   - 20.0      MeV]:	8.06e-03 +/- 5.04e-01%
            Group 2 [0.0        - 6.25e-07  MeV]:	3.61e-01 +/- 5.27e-01%

	Nuclide        =	U-238
	Cross Sections [cm^-1]:
            Group 1 [6.25e-07   - 20.0      MeV]:	7.34e-03 +/- 6.08e-01%
            Group 2 [0.0        - 6.25e-07  MeV]:	6.74e-07 +/- 5.25e-01%

	Nuclide        =	O-16
	Cross Sections [cm^-1]:
            Group 1 [6.25e-07   - 20.0      MeV]:	0.00e+00 +/- nan%
            Group 2 [0.0        - 6.25e-07  MeV]:	0.00e+00 +/- nan%



One can export the entire Library to HDF5 with the Library.build_hdf5_store(...) method as follows:


In [33]:
# Store the cross section data in an "mgxs/mgxs.h5" HDF5 binary file
mgxs_lib.build_hdf5_store(filename='mgxs.h5', directory='mgxs')

The HDF5 store will contain the numerical multi-group cross section data indexed by domain, nuclide and cross section type. Some data workflows may be optimized by storing and retrieving binary representations of the MGXS objects in the Library. This feature is supported through the Library.dump_to_file(...) and Library.load_from_file(...) routines which use Python's pickle module. This is illustrated as follows.


In [34]:
# Store a Library and its MGXS objects in a pickled binary file "mgxs/mgxs.pkl"
mgxs_lib.dump_to_file(filename='mgxs', directory='mgxs')

In [35]:
# Instantiate a new MGXS Library from the pickled binary file "mgxs/mgxs.pkl"
mgxs_lib = openmc.mgxs.Library.load_from_file(filename='mgxs', directory='mgxs')

The Library class may be used to leverage the energy condensation features supported by the MGXS class. In particular, one can use the Library.get_condensed_library(...) with a coarse group structure which is a subset of the original "fine" group structure as shown below.


In [36]:
# Create a 1-group structure
coarse_groups = openmc.mgxs.EnergyGroups(group_edges=[0., 20.])

# Create a new MGXS Library on the coarse 1-group structure
coarse_mgxs_lib = mgxs_lib.get_condensed_library(coarse_groups)

In [37]:
# Retrieve the NuFissionXS object for the fuel cell from the 1-group library
coarse_fuel_mgxs = coarse_mgxs_lib.get_mgxs(fuel_cell, 'nu-fission')

# Show the Pandas DataFrame for the 1-group MGXS
coarse_fuel_mgxs.get_pandas_dataframe()


Out[37]:
cell group in nuclide mean std. dev.
0 10000 1 U-235 0.074383 0.000280
1 10000 1 U-238 0.005959 0.000036
2 10000 1 O-16 0.000000 0.000000

Verification with OpenMOC

Of course it is always a good idea to verify that one's cross sections are accurate. We can easily do so here with the deterministic transport code OpenMOC. We will extract an OpenCG geometry from the summary file and convert it into an equivalent OpenMOC geometry.


In [38]:
# Create an OpenMOC Geometry from the OpenCG Geometry
openmoc_geometry = get_openmoc_geometry(mgxs_lib.opencg_geometry)

Now, we can inject the multi-group cross sections into the equivalent fuel assembly OpenMOC geometry. The openmoc.materialize module supports the loading of Library objects from OpenMC as illustrated below.


In [39]:
# Load the library into the OpenMOC geometry
materials = load_openmc_mgxs_lib(mgxs_lib, openmoc_geometry)

We are now ready to run OpenMOC to verify our cross-sections from OpenMC.


In [40]:
# Generate tracks for OpenMOC
openmoc_geometry.initializeFlatSourceRegions()
track_generator = openmoc.TrackGenerator(openmoc_geometry, num_azim=32, spacing=0.1)
track_generator.generateTracks()

# Run OpenMOC
solver = openmoc.CPUSolver(track_generator)
solver.computeEigenvalue()


[  NORMAL ]  Ray tracing for track segmentation...
[  NORMAL ]  Dumping tracks to file...
[  NORMAL ]  Computing the eigenvalue...
[  NORMAL ]  Iteration 0:	k_eff = 0.854316	res = 0.000E+00
[  NORMAL ]  Iteration 1:	k_eff = 0.801593	res = 1.522E-01
[  NORMAL ]  Iteration 2:	k_eff = 0.761131	res = 6.380E-02
[  NORMAL ]  Iteration 3:	k_eff = 0.731467	res = 5.066E-02
[  NORMAL ]  Iteration 4:	k_eff = 0.709897	res = 3.910E-02
[  NORMAL ]  Iteration 5:	k_eff = 0.695111	res = 2.954E-02
[  NORMAL ]  Iteration 6:	k_eff = 0.685967	res = 2.085E-02
[  NORMAL ]  Iteration 7:	k_eff = 0.681511	res = 1.317E-02
[  NORMAL ]  Iteration 8:	k_eff = 0.680926	res = 6.520E-03
[  NORMAL ]  Iteration 9:	k_eff = 0.683509	res = 1.046E-03
[  NORMAL ]  Iteration 10:	k_eff = 0.688659	res = 3.848E-03
[  NORMAL ]  Iteration 11:	k_eff = 0.695861	res = 7.565E-03
[  NORMAL ]  Iteration 12:	k_eff = 0.704674	res = 1.048E-02
[  NORMAL ]  Iteration 13:	k_eff = 0.714726	res = 1.269E-02
[  NORMAL ]  Iteration 14:	k_eff = 0.725701	res = 1.428E-02
[  NORMAL ]  Iteration 15:	k_eff = 0.737329	res = 1.537E-02
[  NORMAL ]  Iteration 16:	k_eff = 0.749388	res = 1.604E-02
[  NORMAL ]  Iteration 17:	k_eff = 0.761691	res = 1.637E-02
[  NORMAL ]  Iteration 18:	k_eff = 0.774081	res = 1.643E-02
[  NORMAL ]  Iteration 19:	k_eff = 0.786431	res = 1.628E-02
[  NORMAL ]  Iteration 20:	k_eff = 0.798638	res = 1.597E-02
[  NORMAL ]  Iteration 21:	k_eff = 0.810618	res = 1.553E-02
[  NORMAL ]  Iteration 22:	k_eff = 0.822303	res = 1.501E-02
[  NORMAL ]  Iteration 23:	k_eff = 0.833643	res = 1.443E-02
[  NORMAL ]  Iteration 24:	k_eff = 0.844598	res = 1.380E-02
[  NORMAL ]  Iteration 25:	k_eff = 0.855140	res = 1.315E-02
[  NORMAL ]  Iteration 26:	k_eff = 0.865249	res = 1.249E-02
[  NORMAL ]  Iteration 27:	k_eff = 0.874914	res = 1.183E-02
[  NORMAL ]  Iteration 28:	k_eff = 0.884128	res = 1.118E-02
[  NORMAL ]  Iteration 29:	k_eff = 0.892891	res = 1.054E-02
[  NORMAL ]  Iteration 30:	k_eff = 0.901206	res = 9.920E-03
[  NORMAL ]  Iteration 31:	k_eff = 0.909080	res = 9.320E-03
[  NORMAL ]  Iteration 32:	k_eff = 0.916522	res = 8.745E-03
[  NORMAL ]  Iteration 33:	k_eff = 0.923545	res = 8.194E-03
[  NORMAL ]  Iteration 34:	k_eff = 0.930162	res = 7.669E-03
[  NORMAL ]  Iteration 35:	k_eff = 0.936387	res = 7.171E-03
[  NORMAL ]  Iteration 36:	k_eff = 0.942236	res = 6.698E-03
[  NORMAL ]  Iteration 37:	k_eff = 0.947725	res = 6.252E-03
[  NORMAL ]  Iteration 38:	k_eff = 0.952869	res = 5.830E-03
[  NORMAL ]  Iteration 39:	k_eff = 0.957686	res = 5.433E-03
[  NORMAL ]  Iteration 40:	k_eff = 0.962192	res = 5.060E-03
[  NORMAL ]  Iteration 41:	k_eff = 0.966404	res = 4.710E-03
[  NORMAL ]  Iteration 42:	k_eff = 0.970336	res = 4.381E-03
[  NORMAL ]  Iteration 43:	k_eff = 0.974005	res = 4.073E-03
[  NORMAL ]  Iteration 44:	k_eff = 0.977426	res = 3.785E-03
[  NORMAL ]  Iteration 45:	k_eff = 0.980613	res = 3.515E-03
[  NORMAL ]  Iteration 46:	k_eff = 0.983580	res = 3.264E-03
[  NORMAL ]  Iteration 47:	k_eff = 0.986340	res = 3.029E-03
[  NORMAL ]  Iteration 48:	k_eff = 0.988907	res = 2.809E-03
[  NORMAL ]  Iteration 49:	k_eff = 0.991293	res = 2.605E-03
[  NORMAL ]  Iteration 50:	k_eff = 0.993509	res = 2.415E-03
[  NORMAL ]  Iteration 51:	k_eff = 0.995566	res = 2.238E-03
[  NORMAL ]  Iteration 52:	k_eff = 0.997475	res = 2.073E-03
[  NORMAL ]  Iteration 53:	k_eff = 0.999246	res = 1.920E-03
[  NORMAL ]  Iteration 54:	k_eff = 1.000887	res = 1.777E-03
[  NORMAL ]  Iteration 55:	k_eff = 1.002408	res = 1.645E-03
[  NORMAL ]  Iteration 56:	k_eff = 1.003818	res = 1.522E-03
[  NORMAL ]  Iteration 57:	k_eff = 1.005123	res = 1.408E-03
[  NORMAL ]  Iteration 58:	k_eff = 1.006331	res = 1.302E-03
[  NORMAL ]  Iteration 59:	k_eff = 1.007449	res = 1.203E-03
[  NORMAL ]  Iteration 60:	k_eff = 1.008483	res = 1.112E-03
[  NORMAL ]  Iteration 61:	k_eff = 1.009439	res = 1.028E-03
[  NORMAL ]  Iteration 62:	k_eff = 1.010324	res = 9.496E-04
[  NORMAL ]  Iteration 63:	k_eff = 1.011141	res = 8.771E-04
[  NORMAL ]  Iteration 64:	k_eff = 1.011896	res = 8.100E-04
[  NORMAL ]  Iteration 65:	k_eff = 1.012594	res = 7.479E-04
[  NORMAL ]  Iteration 66:	k_eff = 1.013238	res = 6.903E-04
[  NORMAL ]  Iteration 67:	k_eff = 1.013833	res = 6.372E-04
[  NORMAL ]  Iteration 68:	k_eff = 1.014382	res = 5.880E-04
[  NORMAL ]  Iteration 69:	k_eff = 1.014889	res = 5.425E-04
[  NORMAL ]  Iteration 70:	k_eff = 1.015357	res = 5.004E-04
[  NORMAL ]  Iteration 71:	k_eff = 1.015788	res = 4.615E-04
[  NORMAL ]  Iteration 72:	k_eff = 1.016186	res = 4.256E-04
[  NORMAL ]  Iteration 73:	k_eff = 1.016553	res = 3.924E-04
[  NORMAL ]  Iteration 74:	k_eff = 1.016892	res = 3.617E-04
[  NORMAL ]  Iteration 75:	k_eff = 1.017204	res = 3.334E-04
[  NORMAL ]  Iteration 76:	k_eff = 1.017491	res = 3.073E-04
[  NORMAL ]  Iteration 77:	k_eff = 1.017756	res = 2.831E-04
[  NORMAL ]  Iteration 78:	k_eff = 1.018000	res = 2.608E-04
[  NORMAL ]  Iteration 79:	k_eff = 1.018225	res = 2.403E-04
[  NORMAL ]  Iteration 80:	k_eff = 1.018432	res = 2.213E-04
[  NORMAL ]  Iteration 81:	k_eff = 1.018623	res = 2.038E-04
[  NORMAL ]  Iteration 82:	k_eff = 1.018799	res = 1.877E-04
[  NORMAL ]  Iteration 83:	k_eff = 1.018961	res = 1.729E-04
[  NORMAL ]  Iteration 84:	k_eff = 1.019110	res = 1.591E-04
[  NORMAL ]  Iteration 85:	k_eff = 1.019247	res = 1.465E-04
[  NORMAL ]  Iteration 86:	k_eff = 1.019373	res = 1.349E-04
[  NORMAL ]  Iteration 87:	k_eff = 1.019489	res = 1.241E-04
[  NORMAL ]  Iteration 88:	k_eff = 1.019596	res = 1.142E-04
[  NORMAL ]  Iteration 89:	k_eff = 1.019695	res = 1.051E-04
[  NORMAL ]  Iteration 90:	k_eff = 1.019785	res = 9.673E-05
[  NORMAL ]  Iteration 91:	k_eff = 1.019869	res = 8.899E-05
[  NORMAL ]  Iteration 92:	k_eff = 1.019945	res = 8.187E-05
[  NORMAL ]  Iteration 93:	k_eff = 1.020016	res = 7.532E-05
[  NORMAL ]  Iteration 94:	k_eff = 1.020081	res = 6.925E-05
[  NORMAL ]  Iteration 95:	k_eff = 1.020141	res = 6.372E-05
[  NORMAL ]  Iteration 96:	k_eff = 1.020195	res = 5.861E-05
[  NORMAL ]  Iteration 97:	k_eff = 1.020246	res = 5.388E-05
[  NORMAL ]  Iteration 98:	k_eff = 1.020292	res = 4.956E-05
[  NORMAL ]  Iteration 99:	k_eff = 1.020335	res = 4.555E-05
[  NORMAL ]  Iteration 100:	k_eff = 1.020374	res = 4.187E-05
[  NORMAL ]  Iteration 101:	k_eff = 1.020410	res = 3.850E-05
[  NORMAL ]  Iteration 102:	k_eff = 1.020443	res = 3.539E-05
[  NORMAL ]  Iteration 103:	k_eff = 1.020473	res = 3.254E-05
[  NORMAL ]  Iteration 104:	k_eff = 1.020501	res = 2.991E-05
[  NORMAL ]  Iteration 105:	k_eff = 1.020527	res = 2.747E-05
[  NORMAL ]  Iteration 106:	k_eff = 1.020551	res = 2.526E-05
[  NORMAL ]  Iteration 107:	k_eff = 1.020573	res = 2.321E-05
[  NORMAL ]  Iteration 108:	k_eff = 1.020592	res = 2.134E-05
[  NORMAL ]  Iteration 109:	k_eff = 1.020611	res = 1.961E-05
[  NORMAL ]  Iteration 110:	k_eff = 1.020628	res = 1.802E-05
[  NORMAL ]  Iteration 111:	k_eff = 1.020643	res = 1.653E-05
[  NORMAL ]  Iteration 112:	k_eff = 1.020657	res = 1.519E-05
[  NORMAL ]  Iteration 113:	k_eff = 1.020671	res = 1.398E-05
[  NORMAL ]  Iteration 114:	k_eff = 1.020683	res = 1.283E-05
[  NORMAL ]  Iteration 115:	k_eff = 1.020694	res = 1.178E-05
[  NORMAL ]  Iteration 116:	k_eff = 1.020704	res = 1.082E-05

We report the eigenvalues computed by OpenMC and OpenMOC here together to summarize our results.


In [41]:
# Print report of keff and bias with OpenMC
openmoc_keff = solver.getKeff()
openmc_keff = sp.k_combined[0]
bias = (openmoc_keff - openmc_keff) * 1e5

print('openmc keff = {0:1.6f}'.format(openmc_keff))
print('openmoc keff = {0:1.6f}'.format(openmoc_keff))
print('bias [pcm]: {0:1.1f}'.format(bias))


openmc keff = 1.017105
openmoc keff = 1.020704
bias [pcm]: 359.8

There is a non-trivial bias between the eigenvalues computed by OpenMC and OpenMOC. One can show that these biases do not converge to <100 pcm with more particle histories. For heterogeneous geometries, additional measures must be taken to address the following three sources of bias:

  • Appropriate transport-corrected cross sections
  • Spatial discretization of OpenMOC's mesh
  • Constant-in-angle multi-group cross sections

Flux and Pin Power Visualizations

We will conclude this tutorial by illustrating how to visualize the fission rates computed by OpenMOC and OpenMC. First, we extract volume-integrated fission rates from OpenMC's mesh fission rate tally for each pin cell in the fuel assembly.


In [42]:
# Get the OpenMC fission rate mesh tally data
mesh_tally = sp.get_tally(name='mesh tally')
openmc_fission_rates = mesh_tally.get_values(scores=['nu-fission'])

# Reshape array to 2D for plotting
openmc_fission_rates.shape = (17,17)

# Normalize to the average pin power
openmc_fission_rates /= np.mean(openmc_fission_rates)

Next, we extract OpenMOC's volume-averaged fission rates into a 2D 17x17 NumPy array.


In [43]:
# Export OpenMOC's fission rates for each pin cell instance in the fuel assembly
openmoc.process.compute_fission_rates(solver)

# Open the pickle file with the fission rates
fission_rates = pickle.load(open('fission-rates/fission-rates.pkl', 'rb' ))

# Allocate array for fission rates in each fuel pin
openmoc_fission_rates = np.zeros((17, 17))

# Extract fission rates for each fuel pin
for key, value in fission_rates.items():
    lat_x = int(key.split(':')[1].split()[3][1:-1])
    lat_y = int(key.split(':')[1].split()[4][:-1])    
    openmoc_fission_rates[lat_x, lat_y] = value

# Normalize to the average pin fission rate
openmoc_fission_rates /= np.mean(openmoc_fission_rates)

Now we can easily use Matplotlib to visualize the fission rates from OpenMC and OpenMOC side-by-side.


In [44]:
# Plot OpenMC's fission rates in the left subplot
fig = pylab.subplot(121)
pylab.imshow(openmc_fission_rates, interpolation='none', cmap='jet')
pylab.title('OpenMC Fission Rates')

# Plot OpenMOC's fission rates in the right subplot
fig2 = pylab.subplot(122)
pylab.imshow(openmoc_fission_rates, interpolation='none', cmap='jet')
pylab.title('OpenMOC Fission Rates')


Out[44]:
<matplotlib.text.Text at 0x7fa72eda50d0>