Modeling a Pin-Cell¶
This notebook is intended to demonstrate the basic features of the Python API for constructing input files and running OpenMC. In it, we will show how to create a basic reflective pin-cell model that is equivalent to modeling an infinite array of fuel pins. If you have never used OpenMC, this can serve as a good starting point to learn the Python API. We highly recommend having a copy of the Python API reference documentation open in another browser tab that you can refer to.
%matplotlib inline
import openmc
Defining Materials¶
Materials in OpenMC are defined as a set of nuclides with specified atom/weight fractions. To begin, we will create a material by making an instance of the Material class. In OpenMC, many objects, including materials, are identified by a "unique ID" that is simply just a positive integer. These IDs are used when exporting XML files that the solver reads in. They also appear in the output and can be used for identification. Since an integer ID is not very useful by itself, you can also give a material a name as well.
uo2 = openmc.Material(1, "uo2")
print(uo2)
Material ID = 1 Name = uo2 Temperature = None Density = None [sum] S(a,b) Tables Nuclides
On the XML side, you have no choice but to supply an ID. However, in the Python API, if you don't give an ID, one will be automatically generated for you:
mat = openmc.Material()
print(mat)
Material ID = 2 Name = Temperature = None Density = None [sum] S(a,b) Tables Nuclides
We see that an ID of 2 was automatically assigned. Let's now move on to adding nuclides to our uo2 material. The Material object has a method add_nuclide() whose first argument is the name of the nuclide and second argument is the atom or weight fraction.
help(uo2.add_nuclide)
Help on method add_nuclide in module openmc.material:
add_nuclide(nuclide: str, percent: float, percent_type: str = 'ao') method of openmc.material.Material instance
Add a nuclide to the material
Parameters
----------
nuclide : str
Nuclide to add, e.g., 'Mo95'
percent : float
Atom or weight percent
percent_type : {'ao', 'wo'}
'ao' for atom percent and 'wo' for weight percent
We see that by default it assumes we want an atom fraction.
# Add nuclides to uo2
uo2.add_nuclide('U235', 0.03)
uo2.add_nuclide('U238', 0.97)
uo2.add_nuclide('O16', 2.0)
Now we need to assign a total density to the material. We'll use the set_density for this.
uo2.set_density('g/cm3', 10.0)
You may sometimes be given a material specification where all the nuclide densities are in units of atom/b-cm. In this case, you just want the density to be the sum of the constituents. In that case, you can simply run mat.set_density('sum').
With UO2 finished, let's now create materials for the clad and coolant. Note the use of add_element() for zirconium.
zirconium = openmc.Material(name="zirconium")
zirconium.add_element('Zr', 1.0)
zirconium.set_density('g/cm3', 6.6)
water = openmc.Material(name="h2o")
water.add_nuclide('H1', 2.0)
water.add_nuclide('O16', 1.0)
water.set_density('g/cm3', 1.0)
An astute observer might now point out that this water material we just created will only use free-atom cross sections. We need to tell it to use an $S(\alpha,\beta)$ table so that the bound atom cross section is used at thermal energies. To do this, there's an add_s_alpha_beta() method. Note the use of the GND-style name "c_H_in_H2O".
water.add_s_alpha_beta('c_H_in_H2O')
When we go to run the transport solver in OpenMC, it is going to look for a materials.xml file. Thus far, we have only created objects in memory. To actually create a materials.xml file, we need to instantiate a Materials collection and export it to XML.
materials = openmc.Materials([uo2, zirconium, water])
Note that Materials is actually a subclass of Python's built-in list, so we can use methods like append(), insert(), pop(), etc.
materials = openmc.Materials()
materials.append(uo2)
materials += [zirconium, water]
isinstance(materials, list)
True
Finally, we can create the XML file with the export_to_xml() method. In a Jupyter notebook, we can run a shell command by putting ! before it, so in this case we are going to display the materials.xml file that we created.
materials.export_to_xml()
!cat materials.xml
<?xml version='1.0' encoding='utf-8'?>
<materials>
<material depletable="true" id="1" name="uo2">
<density units="g/cm3" value="10.0" />
<nuclide ao="0.03" name="U235" />
<nuclide ao="0.97" name="U238" />
<nuclide ao="2.0" name="O16" />
</material>
<material id="3" name="zirconium">
<density units="g/cm3" value="6.6" />
<nuclide ao="0.5145" name="Zr90" />
<nuclide ao="0.1122" name="Zr91" />
<nuclide ao="0.1715" name="Zr92" />
<nuclide ao="0.1738" name="Zr94" />
<nuclide ao="0.028" name="Zr96" />
</material>
<material id="4" name="h2o">
<density units="g/cm3" value="1.0" />
<nuclide ao="2.0" name="H1" />
<nuclide ao="1.0" name="O16" />
<sab name="c_H_in_H2O" />
</material>
</materials>
Element Expansion¶
Did you notice something really cool that happened to our Zr element? OpenMC automatically turned it into a list of nuclides when it exported it! The way this feature works is as follows:
- First, it checks whether
Materials.cross_sectionshas been set, indicating the path to across_sections.xmlfile. - If
Materials.cross_sectionsisn't set, it looks for theOPENMC_CROSS_SECTIONSenvironment variable. - If either of these are found, it scans the file to see what nuclides are actually available and will expand elements accordingly.
Let's see what happens if we change O16 in water to elemental O.
water.remove_nuclide('O16')
water.add_element('O', 1.0)
materials.export_to_xml()
!cat materials.xml
<?xml version='1.0' encoding='utf-8'?>
<materials>
<material depletable="true" id="1" name="uo2">
<density units="g/cm3" value="10.0" />
<nuclide ao="0.03" name="U235" />
<nuclide ao="0.97" name="U238" />
<nuclide ao="2.0" name="O16" />
</material>
<material id="3" name="zirconium">
<density units="g/cm3" value="6.6" />
<nuclide ao="0.5145" name="Zr90" />
<nuclide ao="0.1122" name="Zr91" />
<nuclide ao="0.1715" name="Zr92" />
<nuclide ao="0.1738" name="Zr94" />
<nuclide ao="0.028" name="Zr96" />
</material>
<material id="4" name="h2o">
<density units="g/cm3" value="1.0" />
<nuclide ao="2.0" name="H1" />
<nuclide ao="0.999621" name="O16" />
<nuclide ao="0.000379" name="O17" />
<sab name="c_H_in_H2O" />
</material>
</materials>
We see that now O16 and O17 were automatically added. O18 is missing because our cross sections file (which is based on ENDF/B-VII.1) doesn't have O18. If OpenMC didn't know about the cross sections file, it would have assumed that all isotopes exist.
The cross_sections.xml file¶
The cross_sections.xml tells OpenMC where it can find nuclide cross sections and $S(\alpha,\beta)$ tables. It serves the same purpose as MCNP's xsdir file and Serpent's xsdata file. As we mentioned, this can be set either by the OPENMC_CROSS_SECTIONS environment variable or the Materials.cross_sections attribute.
Let's have a look at what's inside this file:
!cat $OPENMC_CROSS_SECTIONS | head -n 10
print(' ...')
!cat $OPENMC_CROSS_SECTIONS | tail -n 10
<?xml version='1.0' encoding='utf-8'?>
<cross_sections>
<library materials="Ac225" path="Ac225.h5" type="neutron" />
<library materials="Ac226" path="Ac226.h5" type="neutron" />
<library materials="Ac227" path="Ac227.h5" type="neutron" />
<library materials="Ag107" path="Ag107.h5" type="neutron" />
<library materials="Ag109" path="Ag109.h5" type="neutron" />
<library materials="Ag110_m1" path="Ag110_m1.h5" type="neutron" />
<library materials="Ag111" path="Ag111.h5" type="neutron" />
<library materials="Al27" path="Al27.h5" type="neutron" />
...
<library materials="Cf253" path="wmp/098253.h5" type="wmp" />
<library materials="Cf254" path="wmp/098254.h5" type="wmp" />
<library materials="Es251" path="wmp/099251.h5" type="wmp" />
<library materials="Es252" path="wmp/099252.h5" type="wmp" />
<library materials="Es253" path="wmp/099253.h5" type="wmp" />
<library materials="Es254" path="wmp/099254.h5" type="wmp" />
<library materials="Es254_m1" path="wmp/099254m1.h5" type="wmp" />
<library materials="Es255" path="wmp/099255.h5" type="wmp" />
<library materials="Fm255" path="wmp/100255.h5" type="wmp" />
</cross_sections>
Enrichment¶
Note that the add_element() method has a special argument enrichment that can be used for Uranium. For example, if we know that we want to create 3% enriched UO2, the following would work:
uo2_three = openmc.Material()
uo2_three.add_element('U', 1.0, enrichment=3.0)
uo2_three.add_element('O', 2.0)
uo2_three.set_density('g/cc', 10.0)
Mixtures¶
In OpenMC it is also possible to define materials by mixing existing materials. For example, if we wanted to create MOX fuel out of a mixture of UO2 (97 wt%) and PuO2 (3 wt%) we could do the following:
# Create PuO2 material
puo2 = openmc.Material()
puo2.add_nuclide('Pu239', 0.94)
puo2.add_nuclide('Pu240', 0.06)
puo2.add_nuclide('O16', 2.0)
puo2.set_density('g/cm3', 11.5)
# Create the mixture
mox = openmc.Material.mix_materials([uo2, puo2], [0.97, 0.03], 'wo')
The 'wo' argument in the mix_materials() method specifies that the fractions are weight fractions. Materials can also be mixed by atomic and volume fractions with 'ao' and 'vo', respectively. For 'ao' and 'wo' the fractions must sum to one. For 'vo', if fractions do not sum to one, the remaining fraction is set as void.
Defining Geometry¶
At this point, we have three materials defined, exported to XML, and ready to be used in our model. To finish our model, we need to define the geometric arrangement of materials. OpenMC represents physical volumes using constructive solid geometry (CSG), also known as combinatorial geometry. The object that allows us to assign a material to a region of space is called a Cell (same concept in MCNP, for those familiar). In order to define a region that we can assign to a cell, we must first define surfaces which bound the region. A surface is a locus of zeros of a function of Cartesian coordinates $x$, $y$, and $z$, e.g.
- A plane perpendicular to the x axis: $x - x_0 = 0$
- A cylinder parallel to the z axis: $(x - x_0)^2 + (y - y_0)^2 - R^2 = 0$
- A sphere: $(x - x_0)^2 + (y - y_0)^2 + (z - z_0)^2 - R^2 = 0$
Between those three classes of surfaces (planes, cylinders, spheres), one can construct a wide variety of models. It is also possible to define cones and general second-order surfaces (tori are not currently supported).
Note that defining a surface is not sufficient to specify a volume -- in order to define an actual volume, one must reference the half-space of a surface. A surface half-space is the region whose points satisfy a positive or negative inequality of the surface equation. For example, for a sphere of radius one centered at the origin, the surface equation is $f(x,y,z) = x^2 + y^2 + z^2 - 1 = 0$. Thus, we say that the negative half-space of the sphere, is defined as the collection of points satisfying $f(x,y,z) < 0$, which one can reason is the inside of the sphere. Conversely, the positive half-space of the sphere would correspond to all points outside of the sphere.
Let's go ahead and create a sphere and confirm that what we've told you is true.
sphere = openmc.Sphere(r=1.0)
Note that by default the sphere is centered at the origin so we didn't have to supply x0, y0, or z0 arguments. Strictly speaking, we could have omitted r as well since it defaults to one. To get the negative or positive half-space, we simply need to apply the - or + unary operators, respectively.
(NOTE: Those unary operators are defined by special methods: __pos__ and __neg__ in this case).
inside_sphere = -sphere
outside_sphere = +sphere
Now let's see if inside_sphere actually contains points inside the sphere:
print((0,0,0) in inside_sphere, (0,0,2) in inside_sphere)
print((0,0,0) in outside_sphere, (0,0,2) in outside_sphere)
True False False True
Everything works as expected! Now that we understand how to create half-spaces, we can create more complex volumes by combining half-spaces using Boolean operators: & (intersection), | (union), and ~ (complement). For example, let's say we want to define a region that is the top part of the sphere (all points inside the sphere that have $z > 0$.
z_plane = openmc.ZPlane(0)
northern_hemisphere = -sphere & +z_plane
For many regions, OpenMC can automatically determine a bounding box. To get the bounding box, we use the bounding_box property of a region, which returns a tuple of the lower-left and upper-right Cartesian coordinates for the bounding box:
northern_hemisphere.bounding_box
(array([-1., -1., 0.]), array([1., 1., 1.]))
Now that we see how to create volumes, we can use them to create a cell.
cell = openmc.Cell()
cell.region = northern_hemisphere
# or...
cell = openmc.Cell(region=northern_hemisphere)
By default, the cell is not filled by any material (void). In order to assign a material, we set the fill property of a Cell.
cell.fill = water
Universes and in-line plotting¶
A collection of cells is known as a universe (again, this will be familiar to MCNP/Serpent users) and can be used as a repeatable unit when creating a model. Although we don't need it yet, the benefit of creating a universe is that we can visualize our geometry while we're creating it.
universe = openmc.Universe()
universe.add_cell(cell)
# this also works
universe = openmc.Universe(cells=[cell])
The Universe object has a plot method that will display our the universe as current constructed:
universe.plot(width=(2.0, 2.0), origin=(0.0, 0.0, 0.1))
<matplotlib.image.AxesImage at 0x7fe7dfdc5490>
By default, the plot will appear in the $x$-$y$ plane. We can change that with the basis argument.
universe.plot(width=(2.0, 2.0), basis='xz')
<matplotlib.image.AxesImage at 0x7fe7d7edbeb0>
If we have particular fondness for, say, fuchsia, we can tell the plot() method to make our cell that color.
universe.plot(width=(2.0, 2.0), basis='xz',
colors={cell: 'fuchsia'})
<matplotlib.image.AxesImage at 0x7fe7d550f130>
Pin cell geometry¶
We now have enough knowledge to create our pin-cell. We need three surfaces to define the fuel and clad:
- The outer surface of the fuel -- a cylinder parallel to the z axis
- The inner surface of the clad -- same as above
- The outer surface of the clad -- same as above
These three surfaces will all be instances of openmc.ZCylinder, each with a different radius according to the specification.
fuel_outer_radius = openmc.ZCylinder(r=0.39)
clad_inner_radius = openmc.ZCylinder(r=0.40)
clad_outer_radius = openmc.ZCylinder(r=0.46)
With the surfaces created, we can now take advantage of the built-in operators on surfaces to create regions for the fuel, the gap, and the clad:
fuel_region = -fuel_outer_radius
gap_region = +fuel_outer_radius & -clad_inner_radius
clad_region = +clad_inner_radius & -clad_outer_radius
Now we can create corresponding cells that assign materials to these regions. As with materials, cells have unique IDs that are assigned either manually or automatically. Note that the gap cell doesn't have any material assigned (it is void by default).
fuel = openmc.Cell(name='fuel')
fuel.fill = uo2
fuel.region = fuel_region
gap = openmc.Cell(name='air gap')
gap.region = gap_region
clad = openmc.Cell(name='clad')
clad.fill = zirconium
clad.region = clad_region
Finally, we need to handle the coolant outside of our fuel pin. To do this, we create x- and y-planes that bound the geometry.
pitch = 1.26
left = openmc.XPlane(-pitch/2, boundary_type='reflective')
right = openmc.XPlane(pitch/2, boundary_type='reflective')
bottom = openmc.YPlane(-pitch/2, boundary_type='reflective')
top = openmc.YPlane(pitch/2, boundary_type='reflective')
The water region is going to be everything outside of the clad outer radius and within the box formed as the intersection of four half-spaces.
water_region = +left & -right & +bottom & -top & +clad_outer_radius
moderator = openmc.Cell(name='moderator')
moderator.fill = water
moderator.region = water_region
OpenMC also includes a factory function that generates a rectangular prism that could have made our lives easier.
box = openmc.rectangular_prism(width=pitch, height=pitch,
boundary_type='reflective')
type(box)
openmc.region.Intersection
Pay attention here -- the object that was returned is NOT a surface. It is actually the intersection of four surface half-spaces, just like we created manually before. Thus, we don't need to apply the unary operator (-box). Instead, we can directly combine it with +clad_or.
water_region = box & +clad_outer_radius
The final step is to assign the cells we created to a universe and tell OpenMC that this universe is the "root" universe in our geometry. The Geometry is the final object that is actually exported to XML.
root_universe = openmc.Universe(cells=(fuel, gap, clad, moderator))
geometry = openmc.Geometry()
geometry.root_universe = root_universe
# or...
geometry = openmc.Geometry(root_universe)
geometry.export_to_xml()
!cat geometry.xml
<?xml version='1.0' encoding='utf-8'?> <geometry> <cell id="3" material="1" name="fuel" region="-3" universe="3" /> <cell id="4" material="void" name="air gap" region="3 -4" universe="3" /> <cell id="5" material="3" name="clad" region="4 -5" universe="3" /> <cell id="6" material="4" name="moderator" region="6 -7 8 -9 5" universe="3" /> <surface coeffs="0.0 0.0 0.39" id="3" type="z-cylinder" /> <surface coeffs="0.0 0.0 0.4" id="4" type="z-cylinder" /> <surface coeffs="0.0 0.0 0.46" id="5" type="z-cylinder" /> <surface boundary="reflective" coeffs="-0.63" id="6" type="x-plane" /> <surface boundary="reflective" coeffs="0.63" id="7" type="x-plane" /> <surface boundary="reflective" coeffs="-0.63" id="8" type="y-plane" /> <surface boundary="reflective" coeffs="0.63" id="9" type="y-plane" /> </geometry>
Starting source and settings¶
The Python API has a module openmc.stats with various univariate and multivariate probability distributions. We can use these distributions to create a starting source using the openmc.Source object.
# Create a point source
point = openmc.stats.Point((0, 0, 0))
source = openmc.Source(space=point)
Now let's create a Settings object and give it the source we created along with specifying how many batches and particles we want to run.
settings = openmc.Settings()
settings.source = source
settings.batches = 100
settings.inactive = 10
settings.particles = 1000
settings.export_to_xml()
!cat settings.xml
<?xml version='1.0' encoding='utf-8'?>
<settings>
<run_mode>eigenvalue</run_mode>
<particles>1000</particles>
<batches>100</batches>
<inactive>10</inactive>
<source strength="1.0">
<space type="point">
<parameters>0 0 0</parameters>
</space>
</source>
</settings>
User-defined tallies¶
We actually have all the required files needed to run a simulation. Before we do that though, let's give a quick example of how to create tallies. We will show how one would tally the total, fission, absorption, and (n,$\gamma$) reaction rates for $^{235}$U in the cell containing fuel. Recall that filters allow us to specify where in phase-space we want events to be tallied and scores tell us what we want to tally:
$$X = \underbrace{\int d\mathbf{r} \int d\mathbf{\Omega} \int dE}_{\text{filters}} \; \underbrace{f(\mathbf{r},\mathbf{\Omega},E)}_{\text{scores}} \psi (\mathbf{r},\mathbf{\Omega},E)$$
In this case, the where is "the fuel cell". So, we will create a cell filter specifying the fuel cell.
cell_filter = openmc.CellFilter(fuel)
tally = openmc.Tally(1)
tally.filters = [cell_filter]
The what is the total, fission, absorption, and (n,$\gamma$) reaction rates in $^{235}$U. By default, if we only specify what reactions, it will gives us tallies over all nuclides. We can use the nuclides attribute to name specific nuclides we're interested in.
tally.nuclides = ['U235']
tally.scores = ['total', 'fission', 'absorption', '(n,gamma)']
Similar to the other files, we need to create a Tallies collection and export it to XML.
tallies = openmc.Tallies([tally])
tallies.export_to_xml()
!cat tallies.xml
<?xml version='1.0' encoding='utf-8'?>
<tallies>
<filter id="1" type="cell">
<bins>3</bins>
</filter>
<tally id="1">
<filters>1</filters>
<nuclides>U235</nuclides>
<scores>total fission absorption (n,gamma)</scores>
</tally>
</tallies>
Running OpenMC¶
Running OpenMC from Python can be done using the openmc.run() function. This function allows you to set the number of MPI processes and OpenMP threads, if need be.
openmc.run()
%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%
############### %%%%%%%%%%%%%%%%%%%%%%%%
################## %%%%%%%%%%%%%%%%%%%%%%%
################### %%%%%%%%%%%%%%%%%%%%%%%
#################### %%%%%%%%%%%%%%%%%%%%%%
##################### %%%%%%%%%%%%%%%%%%%%%
###################### %%%%%%%%%%%%%%%%%%%%
####################### %%%%%%%%%%%%%%%%%%
####################### %%%%%%%%%%%%%%%%%
###################### %%%%%%%%%%%%%%%%%
#################### %%%%%%%%%%%%%%%%%
################# %%%%%%%%%%%%%%%%%
############### %%%%%%%%%%%%%%%%
############ %%%%%%%%%%%%%%%
######## %%%%%%%%%%%%%%
%%%%%%%%%%%
| The OpenMC Monte Carlo Code
Copyright | 2011-2022 MIT, UChicago Argonne LLC, and contributors
License | https://docs.openmc.org/en/latest/license.html
Version | 0.13.1
Git SHA1 | 33bc948f4b855c037975f16d16091fe4ecd12de3
Date/Time | 2022-10-04 13:41:58
OpenMP Threads | 2
Reading settings XML file...
Reading cross sections XML file...
Reading materials XML file...
Reading geometry XML file...
Reading U235 from /home/pshriwise/data/xs/openmc/nndc_hdf5/U235.h5
Reading U238 from /home/pshriwise/data/xs/openmc/nndc_hdf5/U238.h5
Reading O16 from /home/pshriwise/data/xs/openmc/nndc_hdf5/O16.h5
Reading Zr90 from /home/pshriwise/data/xs/openmc/nndc_hdf5/Zr90.h5
Reading Zr91 from /home/pshriwise/data/xs/openmc/nndc_hdf5/Zr91.h5
Reading Zr92 from /home/pshriwise/data/xs/openmc/nndc_hdf5/Zr92.h5
Reading Zr94 from /home/pshriwise/data/xs/openmc/nndc_hdf5/Zr94.h5
Reading Zr96 from /home/pshriwise/data/xs/openmc/nndc_hdf5/Zr96.h5
Reading H1 from /home/pshriwise/data/xs/openmc/nndc_hdf5/H1.h5
Reading O17 from /home/pshriwise/data/xs/openmc/nndc_hdf5/O17.h5
Reading c_H_in_H2O from /home/pshriwise/data/xs/openmc/nndc_hdf5/c_H_in_H2O.h5
Minimum neutron data temperature: 294 K
Maximum neutron data temperature: 294 K
Reading tallies XML file...
Preparing distributed cell instances...
Reading plot XML file...
Writing summary.h5 file...
Maximum neutron transport energy: 20000000 eV for U235
Initializing source particles...
====================> K EIGENVALUE SIMULATION <====================
Bat./Gen. k Average k
========= ======== ====================
1/1 1.40026
2/1 1.45321
3/1 1.34098
4/1 1.36425
5/1 1.42225
6/1 1.37529
7/1 1.38197
8/1 1.44688
9/1 1.40427
10/1 1.40456
11/1 1.36640
12/1 1.39125 1.37882 +/- 0.01243
13/1 1.45403 1.40389 +/- 0.02608
14/1 1.39302 1.40117 +/- 0.01864
15/1 1.34153 1.38924 +/- 0.01873
16/1 1.39034 1.38943 +/- 0.01529
17/1 1.36353 1.38573 +/- 0.01344
18/1 1.30980 1.37624 +/- 0.01502
19/1 1.34877 1.37318 +/- 0.01359
20/1 1.41439 1.37730 +/- 0.01284
21/1 1.41695 1.38091 +/- 0.01216
22/1 1.33334 1.37694 +/- 0.01179
23/1 1.38863 1.37784 +/- 0.01088
24/1 1.37182 1.37741 +/- 0.01008
25/1 1.44857 1.38216 +/- 0.01052
26/1 1.37250 1.38155 +/- 0.00986
27/1 1.38933 1.38201 +/- 0.00927
28/1 1.38796 1.38234 +/- 0.00874
29/1 1.41205 1.38390 +/- 0.00842
30/1 1.29854 1.37964 +/- 0.00906
31/1 1.33041 1.37729 +/- 0.00893
32/1 1.33484 1.37536 +/- 0.00873
33/1 1.41818 1.37722 +/- 0.00854
34/1 1.34015 1.37568 +/- 0.00833
35/1 1.37602 1.37569 +/- 0.00799
36/1 1.34516 1.37452 +/- 0.00776
37/1 1.44145 1.37700 +/- 0.00787
38/1 1.47618 1.38054 +/- 0.00837
39/1 1.40590 1.38141 +/- 0.00812
40/1 1.44112 1.38340 +/- 0.00810
41/1 1.47602 1.38639 +/- 0.00838
42/1 1.34668 1.38515 +/- 0.00821
43/1 1.40106 1.38563 +/- 0.00797
44/1 1.35355 1.38469 +/- 0.00779
45/1 1.35236 1.38377 +/- 0.00762
46/1 1.37841 1.38362 +/- 0.00741
47/1 1.48521 1.38636 +/- 0.00771
48/1 1.34912 1.38538 +/- 0.00757
49/1 1.38579 1.38539 +/- 0.00737
50/1 1.37628 1.38516 +/- 0.00719
51/1 1.45863 1.38696 +/- 0.00724
52/1 1.46987 1.38893 +/- 0.00733
53/1 1.39011 1.38896 +/- 0.00716
54/1 1.46214 1.39062 +/- 0.00719
55/1 1.39540 1.39073 +/- 0.00703
56/1 1.42968 1.39157 +/- 0.00693
57/1 1.40848 1.39193 +/- 0.00679
58/1 1.46630 1.39348 +/- 0.00682
59/1 1.42602 1.39415 +/- 0.00672
60/1 1.43031 1.39487 +/- 0.00662
61/1 1.37934 1.39457 +/- 0.00650
62/1 1.41366 1.39493 +/- 0.00638
63/1 1.33133 1.39373 +/- 0.00637
64/1 1.32888 1.39253 +/- 0.00637
65/1 1.41176 1.39288 +/- 0.00626
66/1 1.39418 1.39290 +/- 0.00615
67/1 1.44874 1.39388 +/- 0.00612
68/1 1.51719 1.39601 +/- 0.00638
69/1 1.30328 1.39444 +/- 0.00646
70/1 1.49383 1.39610 +/- 0.00657
71/1 1.48548 1.39756 +/- 0.00662
72/1 1.35073 1.39681 +/- 0.00656
73/1 1.26807 1.39476 +/- 0.00677
74/1 1.40339 1.39490 +/- 0.00666
75/1 1.42471 1.39536 +/- 0.00658
76/1 1.31608 1.39415 +/- 0.00659
77/1 1.31429 1.39296 +/- 0.00659
78/1 1.35579 1.39242 +/- 0.00652
79/1 1.40558 1.39261 +/- 0.00643
80/1 1.45121 1.39344 +/- 0.00639
81/1 1.48841 1.39478 +/- 0.00644
82/1 1.42474 1.39520 +/- 0.00636
83/1 1.42245 1.39557 +/- 0.00629
84/1 1.45461 1.39637 +/- 0.00625
85/1 1.40745 1.39652 +/- 0.00617
86/1 1.42791 1.39693 +/- 0.00610
87/1 1.47831 1.39799 +/- 0.00611
88/1 1.36950 1.39762 +/- 0.00605
89/1 1.40043 1.39766 +/- 0.00597
90/1 1.35277 1.39710 +/- 0.00592
91/1 1.46321 1.39791 +/- 0.00590
92/1 1.49950 1.39915 +/- 0.00596
93/1 1.46649 1.39996 +/- 0.00595
94/1 1.36392 1.39953 +/- 0.00589
95/1 1.35214 1.39898 +/- 0.00585
96/1 1.39164 1.39889 +/- 0.00578
97/1 1.40895 1.39901 +/- 0.00571
98/1 1.43429 1.39941 +/- 0.00566
99/1 1.47610 1.40027 +/- 0.00566
100/1 1.44518 1.40077 +/- 0.00562
Creating state point statepoint.100.h5...
=======================> TIMING STATISTICS <=======================
Total time for initialization = 1.0530e+00 seconds
Reading cross sections = 1.0453e+00 seconds
Total time in simulation = 2.0599e+01 seconds
Time in transport only = 2.0580e+01 seconds
Time in inactive batches = 1.6682e+00 seconds
Time in active batches = 1.8931e+01 seconds
Time synchronizing fission bank = 8.3413e-03 seconds
Sampling source sites = 7.6383e-03 seconds
SEND/RECV source sites = 6.4672e-04 seconds
Time accumulating tallies = 4.0019e-04 seconds
Time writing statepoints = 1.8284e-03 seconds
Total time for finalization = 1.0948e-04 seconds
Total time elapsed = 2.1662e+01 seconds
Calculation Rate (inactive) = 5994.49 particles/second
Calculation Rate (active) = 4754.23 particles/second
============================> RESULTS <============================
k-effective (Collision) = 1.40349 +/- 0.00437
k-effective (Track-length) = 1.40077 +/- 0.00562
k-effective (Absorption) = 1.40672 +/- 0.00326
Combined k-effective = 1.40544 +/- 0.00302
Leakage Fraction = 0.00000 +/- 0.00000
Great! OpenMC already told us our k-effective. It also spit out a file called tallies.out that shows our tallies. This is a very basic method to look at tally data; for more sophisticated methods, see other example notebooks.
!cat tallies.out
============================> TALLY 1 <============================
Cell 3
U235
Total Reaction Rate 0.730682 +/- 0.00279166
Fission Rate 0.5477 +/- 0.00230053
Absorption Rate 0.657455 +/- 0.00272209
(n,gamma) 0.109755 +/- 0.000435716
Geometry plotting¶
We saw before that we could call the Universe.plot() method to show a universe while we were creating our geometry. There is also a built-in plotter in the codebase that is much faster than the Python plotter and has more options. The interface looks somewhat similar to the Universe.plot() method. Instead though, we create Plot instances, assign them to a Plots collection, export it to XML, and then run OpenMC in geometry plotting mode. As an example, let's specify that we want the plot to be colored by material (rather than by cell) and we assign yellow to fuel and blue to water.
plot = openmc.Plot()
plot.filename = 'pinplot'
plot.width = (pitch, pitch)
plot.pixels = (200, 200)
plot.color_by = 'material'
plot.colors = {uo2: 'yellow', water: 'blue'}
With our plot created, we need to add it to a Plots collection which can be exported to XML.
plots = openmc.Plots([plot])
plots.export_to_xml()
!cat plots.xml
<?xml version='1.0' encoding='utf-8'?>
<plots>
<plot basis="xy" color_by="material" filename="pinplot" id="4" type="slice">
<origin>0.0 0.0 0.0</origin>
<width>1.26 1.26</width>
<pixels>200 200</pixels>
<color id="1" rgb="255 255 0" />
<color id="4" rgb="0 0 255" />
</plot>
</plots>
Now we can run OpenMC in plotting mode by calling the plot_geometry() function. Under the hood this is calling openmc --plot.
openmc.plot_geometry()
%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%
############### %%%%%%%%%%%%%%%%%%%%%%%%
################## %%%%%%%%%%%%%%%%%%%%%%%
################### %%%%%%%%%%%%%%%%%%%%%%%
#################### %%%%%%%%%%%%%%%%%%%%%%
##################### %%%%%%%%%%%%%%%%%%%%%
###################### %%%%%%%%%%%%%%%%%%%%
####################### %%%%%%%%%%%%%%%%%%
####################### %%%%%%%%%%%%%%%%%
###################### %%%%%%%%%%%%%%%%%
#################### %%%%%%%%%%%%%%%%%
################# %%%%%%%%%%%%%%%%%
############### %%%%%%%%%%%%%%%%
############ %%%%%%%%%%%%%%%
######## %%%%%%%%%%%%%%
%%%%%%%%%%%
| The OpenMC Monte Carlo Code
Copyright | 2011-2022 MIT, UChicago Argonne LLC, and contributors
License | https://docs.openmc.org/en/latest/license.html
Version | 0.13.1
Git SHA1 | 33bc948f4b855c037975f16d16091fe4ecd12de3
Date/Time | 2022-10-04 13:42:20
OpenMP Threads | 2
Reading settings XML file...
Reading cross sections XML file...
Reading materials XML file...
Reading geometry XML file...
Reading tallies XML file...
Preparing distributed cell instances...
Reading plot XML file...
=======================> PLOTTING SUMMARY <========================
Plot ID: 4
Plot file: pinplot.png
Universe depth: -1
Plot Type: Slice
Origin: 0 0 0
Width: 1.26 1.26
Coloring: Materials
Basis: XY
Pixels: 200 200
Processing plot 4: pinplot.png...
Now, we can use functionality from IPython to display the .png image inline in our notebook:
from IPython.display import Image
Image("pinplot.png")
OpenMC also provides us with a method on the Plot class that simplifies the workflow.
plot.to_ipython_image()