Skip to content

Commit ace73ab

Browse files
GuyStenshimwellpaulromano
authored
Fixed a bug in charged particle energy deposition. (#3416)
Co-authored-by: Jonathan Shimwell <drshimwell@gmail.com> Co-authored-by: Paul Romano <paul.k.romano@gmail.com>
1 parent cb95c78 commit ace73ab

7 files changed

Lines changed: 166 additions & 47 deletions

File tree

docs/source/methods/energy_deposition.rst

Lines changed: 70 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -25,19 +25,36 @@ KERMA (Kinetic Energy Release in Materials) [Mack97]_ coefficients for reaction
2525
:math:`\times` cross-section (e.g., eV-barn) and can be used much like a reaction
2626
cross section for the purpose of tallying energy deposition.
2727

28-
KERMA coefficients can be computed using the energy-balance method with
29-
a nuclear data processing code like NJOY, which performs the following
30-
iteration over all reactions :math:`r` for all isotopes :math:`i`
31-
requested
28+
KERMA coefficients can be computed using the energy-balance method with a
29+
nuclear data processing code like NJOY, which estimates the KERMA coefficients
30+
using the following equation:
3231

3332
.. math::
3433
35-
k_{i, r}(E) = \left(E + Q_{i, r} - \bar{E}_{i, r, n}
34+
k_{i, r}(E) = \left(E + Q_{i, r} - \sum\limits_x \bar{E}_{i, r, x}
35+
\right)\sigma_{i, r}(E),
36+
37+
where the summation is over each secondary particle type :math:`x`. This
38+
equation states that the energy deposited is equal to the energy of the incident
39+
particle plus the reaction :math:`Q` value less the energy of secondary
40+
particles that are transported away from the reaction site. For neutron
41+
interactions, the energy-balance KERMA coefficient is
42+
43+
.. math::
44+
45+
k_{i, r}(E) = \left(E + Q_{i, r} - \sum\limits_x \bar{E}_{i, r, n}
3646
- \bar{E}_{i, r, \gamma}\right)\sigma_{i, r}(E),
3747
38-
removing the energy of neutral particles (neutrons and photons) that are
39-
transported away from the reaction site :math:`\bar{E}`, and the reaction
40-
:math:`Q` value.
48+
where :math:`\bar{E}_{i, r, n}` is the average energy of secondary neutrons and
49+
:math:`\bar{E}_{i, r, \gamma}` is the average energy of secondary photons. For
50+
photon and charged particle interactions, the :math:`Q` value is zero and thus
51+
the KERMA coefficient is
52+
53+
.. math::
54+
:label: energy-balance-photon
55+
56+
k_{i, r}(E) = \left(E - \sum\limits_x \bar{E}_{i, r, x}
57+
\right)\sigma_{i, r}(E).
4158
4259
-------
4360
Fission
@@ -120,7 +137,7 @@ run with :math:`N918` reflecting fission heating computed from NJOY.
120137
This modified heating data is stored as the MT=901 reaction and will be scored
121138
if ``heating-local`` is included in :attr:`openmc.Tally.scores`.
122139

123-
Coupled neutron-photon transport
140+
Coupled Neutron-Photon Transport
124141
--------------------------------
125142

126143
Here, OpenMC instructs ``heatr`` to assume that energy from photons is not
@@ -138,6 +155,50 @@ Let :math:`N301` represent the total heating number returned from this
138155
This modified heating data is stored as the MT=301 reaction and will be scored
139156
if ``heating`` is included in :attr:`openmc.Tally.scores`.
140157

158+
Photons and Charged Particles
159+
-----------------------------
160+
161+
In OpenMC, energy deposition from photons or charged particles is scored using
162+
the energy balance method based on Equation :eq:`energy-balance-photon`. Special
163+
consideration is given to electrons and positrons as described below.
164+
165+
+++++++++++++++++
166+
Charged Particles
167+
+++++++++++++++++
168+
169+
OpenMC tracks photons interaction by interaction so the energy deposited in each
170+
collision is easily attributed back to the nuclide and reaction for which the
171+
photon interacted with. Charged particles (electrons and photons) aren't tracked
172+
in the same way. For charged particles, OpenMC assumes that all their energy
173+
(less the energy of bremsstrahlung radiation) is deposited in the material in
174+
which they were born. In this way it is harder to trace how much energy should
175+
be attributed in each nuclide.
176+
177+
According to the CSDA approximation (see :ref:`ttb`) the energy deposited by a
178+
charged particle with kinetic energy :math:`T` in the :math:`i`-th element can
179+
be calculated as:
180+
181+
.. math::
182+
183+
E_{i} = \int_{0}^{R(T)} w_{i}S_{\text{col,i}} dx
184+
185+
where :math:`R(T)` is the CSDA range of the charged particle,
186+
:math:`S_{\text{col},i}` is the collision stopping power of the charged particle
187+
in the :math:`i`-th element and :math:`w_i` is the mass fraction of the
188+
:math:`i`-th element. According to the Bethe formula the collision stopping
189+
power of the :math:`i`-th element is proportional to :math:`Z_i/A_i`, so the
190+
fractional collision stopping power from the :math:`i`-th element is:
191+
192+
.. math::
193+
194+
\frac{w_{i}S_{\text{col},i}(T)}{S_{\text{col}}(T)} =
195+
\frac{\frac{w_{i}Z_{i}}{A_{i}}}{\sum_{i}\frac{w_{i}Z_{i}}{A_{i}}} =
196+
\frac{\gamma_i Z_{i}}{\sum_{i}\gamma_i Z_{i}}.
197+
198+
where :math:`\gamma_i` is the atomic fraction of the :math:`i`-th element.
199+
Therefore, the energy deposited by charged particles should be attributed to
200+
a given element according to its fractional charge density.
201+
141202
----------
142203
References
143204
----------

include/openmc/material.h

Lines changed: 5 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -107,6 +107,10 @@ class Material {
107107
//! \return Density in [g/cm^3]
108108
double density_gpcc() const { return density_gpcc_; }
109109

110+
//! Get charge density in [e/b-cm]
111+
//! \return Charge density in [e/b-cm]
112+
double charge_density() const { return charge_density_; };
113+
110114
//! Get name
111115
//! \return Material name
112116
const std::string& name() const { return name_; }
@@ -177,6 +181,7 @@ class Material {
177181
xt::xtensor<double, 1> atom_density_; //!< Nuclide atom density in [atom/b-cm]
178182
double density_; //!< Total atom density in [atom/b-cm]
179183
double density_gpcc_; //!< Total atom density in [g/cm^3]
184+
double charge_density_; //!< Total charge density in [e/b-cm]
180185
double volume_ {-1.0}; //!< Volume in [cm^3]
181186
vector<bool> p0_; //!< Indicate which nuclides are to be treated with
182187
//!< iso-in-lab scattering

src/material.cpp

Lines changed: 9 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -454,12 +454,15 @@ void Material::normalize_density()
454454
// Calculate nuclide atom densities
455455
atom_density_ *= density_;
456456

457-
// Calculate density in g/cm^3.
457+
// Calculate density in [g/cm^3] and charge density in [e/b-cm]
458458
density_gpcc_ = 0.0;
459+
charge_density_ = 0.0;
459460
for (int i = 0; i < nuclide_.size(); ++i) {
460461
int i_nuc = nuclide_[i];
461462
double awr = settings::run_CE ? data::nuclides[i_nuc]->awr_ : 1.0;
463+
int z = settings::run_CE ? data::nuclides[i_nuc]->Z_ : 0.0;
462464
density_gpcc_ += atom_density_(i) * awr * MASS_NEUTRON / N_AVOGADRO;
465+
charge_density_ += atom_density_(i) * z;
463466
}
464467
}
465468

@@ -982,12 +985,15 @@ void Material::set_density(double density, const std::string& units)
982985
// Recalculate nuclide atom densities based on given density
983986
atom_density_ *= density;
984987

985-
// Calculate density in g/cm^3.
988+
// Calculate density in g/cm^3 and charge density in [e/b-cm]
986989
density_gpcc_ = 0.0;
990+
charge_density_ = 0.0;
987991
for (int i = 0; i < nuclide_.size(); ++i) {
988992
int i_nuc = nuclide_[i];
989993
double awr = data::nuclides[i_nuc]->awr_;
994+
int z = settings::run_CE ? data::nuclides[i_nuc]->Z_ : 0.0;
990995
density_gpcc_ += atom_density_(i) * awr * MASS_NEUTRON / N_AVOGADRO;
996+
charge_density_ += atom_density_(i) * z;
991997
}
992998
} else if (units == "g/cm3" || units == "g/cc") {
993999
// Determine factor by which to change densities
@@ -998,6 +1004,7 @@ void Material::set_density(double density, const std::string& units)
9981004
density_gpcc_ = density;
9991005
density_ *= f;
10001006
atom_density_ *= f;
1007+
charge_density_ *= f;
10011008
} else {
10021009
throw std::invalid_argument {
10031010
"Invalid units '" + std::string(units.data()) + "' specified."};

src/tallies/tally_scoring.cpp

Lines changed: 37 additions & 30 deletions
Original file line numberDiff line numberDiff line change
@@ -325,6 +325,39 @@ double score_neutron_heating(const Particle& p, const Tally& tally, double flux,
325325
return score;
326326
}
327327

328+
//! Helper function to obtain particle heating [eV]
329+
330+
double score_particle_heating(const Particle& p, const Tally& tally,
331+
double flux, int rxn_bin, int i_nuclide, double atom_density)
332+
{
333+
if (p.type() == ParticleType::neutron)
334+
return score_neutron_heating(
335+
p, tally, flux, rxn_bin, i_nuclide, atom_density);
336+
if (i_nuclide == -1 || i_nuclide == p.event_nuclide() ||
337+
p.event_nuclide() == -1) {
338+
// Get the pre-collision energy of the particle.
339+
auto E = p.E_last();
340+
// The energy deposited is the difference between the pre-collision
341+
// and post-collision energy...
342+
double score = E - p.E();
343+
// ...less the energy of any secondary particles since they will be
344+
// transported individually later
345+
score -= p.bank_second_E();
346+
score *= p.wgt_last();
347+
348+
// if no event_nuclide (charged particle) scale energy deposition by
349+
// fractional charge density
350+
if (i_nuclide != -1 && p.event_nuclide() == -1) {
351+
const auto& mat {model::materials[p.material()]};
352+
int z = data::nuclides[i_nuclide]->Z_;
353+
auto i = mat->mat_nuclide_index_[i_nuclide];
354+
score *= (z * mat->atom_density_[i] / mat->charge_density());
355+
}
356+
return score;
357+
}
358+
return 0.0;
359+
}
360+
328361
//! Helper function for nu-fission tallies with energyout filters.
329362
//
330363
//! In this case, we may need to score to multiple bins if there were multiple
@@ -1007,23 +1040,8 @@ void score_general_ce_nonanalog(Particle& p, int i_tally, int start_index,
10071040
break;
10081041

10091042
case HEATING:
1010-
if (p.type() == Type::neutron) {
1011-
score = score_neutron_heating(
1012-
p, tally, flux, HEATING, i_nuclide, atom_density);
1013-
} else {
1014-
if (i_nuclide == -1 || i_nuclide == p.event_nuclide()) {
1015-
// The energy deposited is the difference between the pre-collision
1016-
// and post-collision energy...
1017-
score = E - p.E();
1018-
// ...less the energy of any secondary particles since they will be
1019-
// transported individually later
1020-
score -= p.bank_second_E();
1021-
1022-
score *= p.wgt_last();
1023-
} else {
1024-
score = 0.0;
1025-
}
1026-
}
1043+
score = score_particle_heating(
1044+
p, tally, flux, HEATING, i_nuclide, atom_density);
10271045
break;
10281046

10291047
default:
@@ -1539,19 +1557,8 @@ void score_general_ce_analog(Particle& p, int i_tally, int start_index,
15391557
break;
15401558

15411559
case HEATING:
1542-
if (p.type() == Type::neutron) {
1543-
score = score_neutron_heating(
1544-
p, tally, flux, HEATING, i_nuclide, atom_density);
1545-
} else {
1546-
// The energy deposited is the difference between the pre-collision and
1547-
// post-collision energy...
1548-
score = E - p.E();
1549-
// ...less the energy of any secondary particles since they will be
1550-
// transported individually later
1551-
score -= p.bank_second_E();
1552-
1553-
score *= p.wgt_last();
1554-
}
1560+
score = score_particle_heating(
1561+
p, tally, flux, HEATING, i_nuclide, atom_density);
15551562
break;
15561563

15571564
default:

tests/regression_tests/photon_production/results_true.dat

Lines changed: 4 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -67,8 +67,8 @@ tally 3:
6767
0.000000E+00
6868
0.000000E+00
6969
0.000000E+00
70-
0.000000E+00
71-
0.000000E+00
70+
1.774484E+05
71+
3.148794E+10
7272
0.000000E+00
7373
0.000000E+00
7474
0.000000E+00
@@ -79,8 +79,8 @@ tally 3:
7979
0.000000E+00
8080
0.000000E+00
8181
0.000000E+00
82-
0.000000E+00
83-
0.000000E+00
82+
7.692488E+03
83+
5.917437E+07
8484
0.000000E+00
8585
0.000000E+00
8686
0.000000E+00

tests/unit_tests/test_data_multipole.py

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -52,13 +52,13 @@ def test_from_endf():
5252
pytest.importorskip('vectfit')
5353
endf_data = os.environ['OPENMC_ENDF_DATA']
5454
endf_file = os.path.join(endf_data, 'neutrons', 'n-001_H_001.endf')
55-
return openmc.data.WindowedMultipole.from_endf(
55+
assert openmc.data.WindowedMultipole.from_endf(
5656
endf_file, log=True, wmp_options={"n_win": 400, "n_cf": 3})
5757

5858

5959
def test_from_endf_search():
6060
pytest.importorskip('vectfit')
6161
endf_data = os.environ['OPENMC_ENDF_DATA']
6262
endf_file = os.path.join(endf_data, 'neutrons', 'n-095_Am_244.endf')
63-
return openmc.data.WindowedMultipole.from_endf(
63+
assert openmc.data.WindowedMultipole.from_endf(
6464
endf_file, log=True, wmp_options={"search": True, 'rtol':1e-2})
Lines changed: 39 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,39 @@
1+
import openmc
2+
from pytest import approx
3+
4+
5+
def test_nuclide_heating(run_in_tmpdir):
6+
mat = openmc.Material()
7+
mat.add_nuclide("Li6", 0.5)
8+
mat.add_nuclide("Li7", 0.5)
9+
mat.set_density("g/cm3", 1.0)
10+
11+
sphere = openmc.Sphere(r=20, boundary_type="reflective")
12+
inside_sphere = openmc.Cell(fill=mat, region=-sphere)
13+
model = openmc.Model()
14+
model.geometry = openmc.Geometry([inside_sphere])
15+
16+
model.settings.particles = 1000
17+
model.settings.batches = 1
18+
model.settings.photon_transport = True
19+
model.settings.electron_treatment = "ttb"
20+
model.settings.cutoff = {"energy_photon": 1000}
21+
model.settings.run_mode = "fixed source"
22+
model.settings.source = openmc.IndependentSource(
23+
energy=openmc.stats.delta_function(10.0e6),
24+
particle="photon"
25+
)
26+
27+
# Create two tallies, one with heating by nuclide and one with total heating
28+
tally1 = openmc.Tally()
29+
tally1.scores = ["heating"]
30+
tally1.nuclides = mat.get_nuclides()
31+
tally2 = openmc.Tally()
32+
tally2.scores = ["heating"]
33+
model.tallies = [tally1, tally2]
34+
35+
# Run the model
36+
model.run(apply_tally_results=True)
37+
38+
# Make sure the heating results are consistent
39+
assert tally1.mean.sum() == approx(tally2.mean.sum())

0 commit comments

Comments
 (0)