10. Depletion¶
OpenMC supports coupled depletion, or burnup, calculations through the
openmc.deplete
Python module. OpenMC solves the transport equation to
obtain transmutation reaction rates, and then the reaction rates are used to
solve a set of transmutation equations that determine the evolution of nuclide
densities within a material. The nuclide densities predicted as some future time
are then used to determine updated reaction rates, and the process is repeated
for as many timesteps as are requested.
The depletion module is designed such that the flux/reaction rate solution (the
transport “operator”) is completely isolated from the solution of the
transmutation equations and the method used for advancing time. At present, the
openmc.deplete
module offers a single transport operator,
openmc.deplete.Operator
(which uses the OpenMC transport solver), but
in principle additional operator classes based on other transport codes could be
implemented and no changes to the depletion solver itself would be needed. The
operator class requires a openmc.Geometry
instance and a
openmc.Settings
instance:
geom = openmc.Geometry()
settings = openmc.Settings()
...
op = openmc.deplete.Operator(geom, settings)
openmc.deplete
supports multiple time-integration methods for determining
material compositions over time. Each method appears as a different class.
For example, openmc.deplete.CECMIntegrator
runs a depletion calculation
using the CE/CM algorithm (deplete over a timestep using the middle-of-step
reaction rates). An instance of openmc.deplete.Operator
is passed to
one of these functions along with the power level and timesteps:
power = 1200.0e6
days = 24*60*60
timesteps = [10.0*days, 10.0*days, 10.0*days]
openmc.deplete.CECMIntegrator(op, power, timesteps).integrate()
The coupled transport-depletion problem is executed, and once it is done a
depletion_results.h5
file is written. The results can be analyzed using the
openmc.deplete.ResultsList
class. This class has methods that allow for
easy retrieval of k-effective, nuclide concentrations, and reaction rates over
time:
results = openmc.deplete.ResultsList.from_hdf5("depletion_results.h5")
time, keff = results.get_eigenvalue()
Note that the coupling between the transport solver and the transmutation solver happens in-memory rather than by reading/writing files on disk.
10.1. Caveats¶
10.1.1. Energy Deposition¶
The default energy deposition mode, "fission-q"
, instructs the
openmc.deplete.Operator
to normalize reaction rates using the product
of fission reaction rates and fission Q values taken from the depletion chain.
This approach does not consider indirect contributions to energy deposition,
such as neutron heating and energy from secondary photons. In doing this,
the energy deposited during a transport calculation will be lower than expected.
This causes the reaction rates to be over-adjusted to hit the user-specific power,
or power density, leading to an over-depletion of burnable materials.
There are some remedies. First, the fission Q values can be directly set in a variety of ways. This requires knowing what the total fission energy release should be, including indirect components. Some examples are provided below:
# use a dictionary of fission_q values
fission_q = {"U235": 202} # energy in MeV
# create a modified chain and write it to a new file
chain = openmc.deplete.Chain.from_xml("chain.xml", fission_q)
chain.export_to_xml("chain_mod_q.xml")
op = openmc.deplete.Operator(geometry, setting, "chain_mod_q.xml")
# alternatively, pass the modified fission Q directly to the operator
op = openmc.deplete.Operator(geometry, setting, "chain.xml",
fission_q=fission_q)
A more complete way to model the energy deposition is to use the modified heating reactions described in Heating and Energy Deposition. These values can be used to normalize reaction rates instead of using the fission reaction rates with:
op = openmc.deplete.Operator(geometry, settings, "chain.xml",
energy_mode="energy-deposition")
These modified heating libraries can be generated by running the latest version
of openmc.data.IncidentNeutron.from_njoy()
, and will eventually be bundled into
the distributed libraries.
10.1.2. Local Spectra and Repeated Materials¶
It is not uncommon to explicitly create a single burnable material across many locations.
From a pure transport perspective, there is nothing wrong with creating a single
3.5 wt.% enriched fuel fuel_3
, and placing that fuel in every fuel pin in an assembly
or even full core problem. This certainly expedites the model making process, but can pose
issues with depletion.
Under this setup, openmc.deplete
will deplete a single fuel_3
material using
a single set of reaction rates, and produce a single new composition for the next time
step. This can be problematic if the same fuel_3
is used in very different regions
of the problem.
As an example, consider a full-scale power reactor core with vacuum boundary
conditions, and with fuel pins solely composed of the same fuel_3
material.
The fuel pins towards the center of the problem will surely experience a more intense
neutron flux and greater reaction rates than those towards the edge of the domain.
This indicates that the fuel in the center should be at a more depleted state than
periphery pins, at least for the fist depletion step.
However, without any other instructions, OpenMC will deplete fuel_3
as a single
material, and all of the fuel pins will have an identical composition at the next
transport step.
This can be countered by instructing the operator to treat repeated instances of the same material as a unique material definition with:
op = openmc.deplete.Operator(geometry, settings, chain_file,
diff_burnable_mats=True)
For our example problem, this would deplete fuel on the outer region of the problem with different reaction rates than those in the center. Materials will be depleted corresponding to their local neutron spectra, and have unique compositions at each transport step.
Note
This will increase the total memory usage and run time due to an increased number of tallies and material definitions.