Skip to content

Commit 219d827

Browse files
shimwelljon-proximafusionpaulromano
authored
Allow tracklength estimator for neutron heating (#3915)
Co-authored-by: Jon Shimwell <jon@proximafusion.com> Co-authored-by: Paul Romano <paul.k.romano@gmail.com>
1 parent 2dd5322 commit 219d827

2 files changed

Lines changed: 62 additions & 2 deletions

File tree

src/tallies/tally.cpp

Lines changed: 18 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -644,8 +644,24 @@ void Tally::set_scores(const vector<std::string>& scores)
644644
break;
645645

646646
case HEATING:
647-
if (settings::photon_transport)
648-
estimator_ = TallyEstimator::COLLISION;
647+
if (settings::photon_transport) {
648+
// Photon heating requires a collision estimator (analog energy
649+
// balance). However, if the tally only scores neutrons, we can keep the
650+
// tracklength estimator since neutron heating uses kerma coefficients
651+
// that support tracklength scoring.
652+
bool neutron_only = false;
653+
for (auto i_filt : filters_) {
654+
auto pf =
655+
dynamic_cast<ParticleFilter*>(model::tally_filters[i_filt].get());
656+
if (pf && pf->particles().size() == 1 &&
657+
pf->particles()[0].is_neutron()) {
658+
neutron_only = true;
659+
break;
660+
}
661+
}
662+
if (!neutron_only)
663+
estimator_ = TallyEstimator::COLLISION;
664+
}
649665
break;
650666

651667
case SCORE_PULSE_HEIGHT: {

tests/unit_tests/test_tally.py

Lines changed: 44 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -17,3 +17,47 @@ def test_tally_init_args():
1717
assert tally.filters == [filter]
1818
assert tally.nuclides == ['U235']
1919
assert tally.estimator == 'tracklength'
20+
21+
22+
def test_heating_estimator_with_photon_transport(run_in_tmpdir):
23+
"""Test that a neutron-only heating tally uses the tracklength estimator
24+
even when photon transport is enabled.
25+
26+
Without a neutron particle filter, photon heating requires the collision
27+
estimator (analog energy balance). But neutron heating has kerma
28+
coefficients that support tracklength scoring, so a neutron-only tally
29+
should keep the tracklength estimator.
30+
"""
31+
mat = openmc.Material()
32+
mat.add_nuclide('Si28', 1.0)
33+
mat.set_density('g/cm3', 2.3)
34+
sph = openmc.Sphere(r=10.0, boundary_type='vacuum')
35+
cell = openmc.Cell(fill=mat, region=-sph)
36+
37+
model = openmc.Model()
38+
model.geometry = openmc.Geometry([cell])
39+
model.settings.run_mode = 'fixed source'
40+
model.settings.particles = 100
41+
model.settings.batches = 1
42+
model.settings.photon_transport = True
43+
model.settings.source = openmc.IndependentSource(
44+
space=openmc.stats.Point()
45+
)
46+
47+
neutron_filter = openmc.ParticleFilter(['neutron'])
48+
49+
# Neutron-only heating — should get tracklength estimator
50+
t_neutron = openmc.Tally(name='heating_neutron')
51+
t_neutron.filters = [neutron_filter]
52+
t_neutron.scores = ['heating']
53+
54+
# No particle filter — should get collision estimator
55+
t_all = openmc.Tally(name='heating_all')
56+
t_all.scores = ['heating']
57+
58+
model.tallies = [t_neutron, t_all]
59+
60+
sp_filename = model.run()
61+
with openmc.StatePoint(sp_filename) as sp:
62+
assert sp.tallies[t_neutron.id].estimator == 'tracklength'
63+
assert sp.tallies[t_all.id].estimator == 'collision'

0 commit comments

Comments
 (0)