• Home
  • Features
  • Pricing
  • Docs
  • Announcements
  • Sign In

openmc-dev / openmc / 17852533892

19 Sep 2025 08:10AM UTC coverage: 85.197% (-0.004%) from 85.201%
17852533892

push

github

web-flow
Combing for fission site sampling, and delayed neutron emission time (#2992)

Co-authored-by: Gavin Ridley <gavin.keith.ridley@gmail.com>
Co-authored-by: Paul Romano <paul.k.romano@gmail.com>

54 of 58 new or added lines in 7 files covered. (93.1%)

4 existing lines in 3 files now uncovered.

53173 of 62412 relevant lines covered (85.2%)

38753013.29 hits per line

Source File
Press 'n' to go to next uncovered line, 'b' for previous

88.91
/src/physics.cpp
1
#include "openmc/physics.h"
2

3
#include "openmc/bank.h"
4
#include "openmc/bremsstrahlung.h"
5
#include "openmc/chain.h"
6
#include "openmc/constants.h"
7
#include "openmc/distribution_multi.h"
8
#include "openmc/eigenvalue.h"
9
#include "openmc/endf.h"
10
#include "openmc/error.h"
11
#include "openmc/ifp.h"
12
#include "openmc/material.h"
13
#include "openmc/math_functions.h"
14
#include "openmc/message_passing.h"
15
#include "openmc/ncrystal_interface.h"
16
#include "openmc/nuclide.h"
17
#include "openmc/photon.h"
18
#include "openmc/physics_common.h"
19
#include "openmc/random_dist.h"
20
#include "openmc/random_lcg.h"
21
#include "openmc/reaction.h"
22
#include "openmc/search.h"
23
#include "openmc/secondary_uncorrelated.h"
24
#include "openmc/settings.h"
25
#include "openmc/simulation.h"
26
#include "openmc/string_utils.h"
27
#include "openmc/tallies/tally.h"
28
#include "openmc/thermal.h"
29
#include "openmc/weight_windows.h"
30

31
#include <fmt/core.h>
32

33
#include <algorithm> // for max, min, max_element
34
#include <cmath>     // for sqrt, exp, log, abs, copysign
35
#include <xtensor/xview.hpp>
36

37
namespace openmc {
38

39
//==============================================================================
40
// Non-member functions
41
//==============================================================================
42

43
void collision(Particle& p)
759,186,109✔
44
{
45
  // Add to collision counter for particle
46
  ++(p.n_collision());
759,186,109✔
47

48
  // Sample reaction for the material the particle is in
49
  switch (p.type()) {
759,186,109✔
50
  case ParticleType::neutron:
691,039,604✔
51
    sample_neutron_reaction(p);
691,039,604✔
52
    break;
691,039,604✔
53
  case ParticleType::photon:
16,213,293✔
54
    sample_photon_reaction(p);
16,213,293✔
55
    break;
16,213,293✔
56
  case ParticleType::electron:
51,847,521✔
57
    sample_electron_reaction(p);
51,847,521✔
58
    break;
51,847,521✔
59
  case ParticleType::positron:
85,691✔
60
    sample_positron_reaction(p);
85,691✔
61
    break;
85,691✔
62
  }
63

64
  if (settings::weight_window_checkpoint_collision)
759,186,109✔
65
    apply_weight_windows(p);
759,186,109✔
66

67
  // Kill particle if energy falls below cutoff
68
  int type = static_cast<int>(p.type());
759,186,109✔
69
  if (p.E() < settings::energy_cutoff[type]) {
759,186,109✔
70
    p.wgt() = 0.0;
5,007,669✔
71
  }
72

73
  // Display information about collision
74
  if (settings::verbosity >= 10 || p.trace()) {
759,186,109✔
75
    std::string msg;
66✔
76
    if (p.event() == TallyEvent::KILL) {
66✔
77
      msg = fmt::format("    Killed. Energy = {} eV.", p.E());
×
78
    } else if (p.type() == ParticleType::neutron) {
66✔
79
      msg = fmt::format("    {} with {}. Energy = {} eV.",
186✔
80
        reaction_name(p.event_mt()), data::nuclides[p.event_nuclide()]->name_,
132✔
81
        p.E());
66✔
82
    } else if (p.type() == ParticleType::photon) {
×
83
      msg = fmt::format("    {} with {}. Energy = {} eV.",
×
84
        reaction_name(p.event_mt()),
×
85
        to_element(data::nuclides[p.event_nuclide()]->name_), p.E());
×
86
    } else {
87
      msg = fmt::format("    Disappeared. Energy = {} eV.", p.E());
×
88
    }
89
    write_message(msg, 1);
66✔
90
  }
66✔
91
}
759,186,109✔
92

93
void sample_neutron_reaction(Particle& p)
691,039,604✔
94
{
95
  // Sample a nuclide within the material
96
  int i_nuclide = sample_nuclide(p);
691,039,604✔
97

98
  // Save which nuclide particle had collision with
99
  p.event_nuclide() = i_nuclide;
691,039,604✔
100

101
  // Create fission bank sites. Note that while a fission reaction is sampled,
102
  // it never actually "happens", i.e. the weight of the particle does not
103
  // change when sampling fission sites. The following block handles all
104
  // absorption (including fission)
105

106
  const auto& nuc {data::nuclides[i_nuclide]};
691,039,604✔
107

108
  if (nuc->fissionable_ && p.neutron_xs(i_nuclide).fission > 0.0) {
691,039,604✔
109
    auto& rx = sample_fission(i_nuclide, p);
88,035,708✔
110
    if (settings::run_mode == RunMode::EIGENVALUE) {
88,035,708✔
111
      create_fission_sites(p, i_nuclide, rx);
82,948,769✔
112
    } else if (settings::run_mode == RunMode::FIXED_SOURCE &&
5,086,939✔
113
               settings::create_fission_neutrons) {
114
      create_fission_sites(p, i_nuclide, rx);
5,073,035✔
115

116
      // Make sure particle population doesn't grow out of control for
117
      // subcritical multiplication problems.
118
      if (p.secondary_bank().size() >= settings::max_secondaries) {
5,073,035✔
119
        fatal_error(
×
120
          "The secondary particle bank appears to be growing without "
121
          "bound. You are likely running a subcritical multiplication problem "
122
          "with k-effective close to or greater than one.");
123
      }
124
    }
125
  }
126

127
  // Create secondary photons
128
  if (settings::photon_transport) {
691,039,604✔
129
    sample_secondary_photons(p, i_nuclide);
8,809,548✔
130
  }
131

132
  // If survival biasing is being used, the following subroutine adjusts the
133
  // weight of the particle. Otherwise, it checks to see if absorption occurs
134

135
  if (p.neutron_xs(i_nuclide).absorption > 0.0) {
691,039,604✔
136
    absorption(p, i_nuclide);
691,036,079✔
137
  }
138
  if (!p.alive())
691,039,604✔
139
    return;
17,151,039✔
140

141
  // Sample a scattering reaction and determine the secondary energy of the
142
  // exiting neutron
143
  const auto& ncrystal_mat = model::materials[p.material()]->ncrystal_mat();
673,888,565✔
144
  if (ncrystal_mat && p.E() < NCRYSTAL_MAX_ENERGY) {
673,888,565✔
145
    ncrystal_mat.scatter(p);
158,829✔
146
  } else {
147
    scatter(p, i_nuclide);
673,729,736✔
148
  }
149

150
  // Advance URR seed stream 'N' times after energy changes
151
  if (p.E() != p.E_last()) {
673,888,565✔
152
    advance_prn_seed(data::nuclides.size(), &p.seeds(STREAM_URR_PTABLE));
673,569,334✔
153
  }
154

155
  // Play russian roulette if survival biasing is turned on
156
  if (settings::survival_biasing) {
673,888,565✔
157
    // if survival normalization is on, use normalized weight cutoff and
158
    // normalized weight survive
159
    if (settings::survival_normalization) {
499,950✔
160
      if (p.wgt() < settings::weight_cutoff * p.wgt_born()) {
×
161
        russian_roulette(p, settings::weight_survive * p.wgt_born());
×
162
      }
163
    } else if (p.wgt() < settings::weight_cutoff) {
499,950✔
164
      russian_roulette(p, settings::weight_survive);
56,991✔
165
    }
166
  }
167
}
168

169
void create_fission_sites(Particle& p, int i_nuclide, const Reaction& rx)
88,021,804✔
170
{
171
  // If uniform fission source weighting is turned on, we increase or decrease
172
  // the expected number of fission sites produced
173
  double weight = settings::ufs_on ? ufs_get_weight(p) : 1.0;
88,021,804✔
174

175
  // Determine the expected number of neutrons produced
176
  double nu_t = p.wgt() / simulation::keff * weight *
88,021,804✔
177
                p.neutron_xs(i_nuclide).nu_fission /
88,021,804✔
178
                p.neutron_xs(i_nuclide).total;
88,021,804✔
179

180
  // Sample the number of neutrons produced
181
  int nu = static_cast<int>(nu_t);
88,021,804✔
182
  if (prn(p.current_seed()) <= (nu_t - nu))
88,021,804✔
183
    ++nu;
14,453,573✔
184

185
  // If no neutrons were produced then don't continue
186
  if (nu == 0)
88,021,804✔
187
    return;
69,681,269✔
188

189
  // Initialize the counter of delayed neutrons encountered for each delayed
190
  // group.
191
  double nu_d[MAX_DELAYED_GROUPS] = {0.};
18,340,535✔
192

193
  // Clear out particle's nu fission bank
194
  p.nu_bank().clear();
18,340,535✔
195

196
  p.fission() = true;
18,340,535✔
197

198
  // Determine whether to place fission sites into the shared fission bank
199
  // or the secondary particle bank.
200
  bool use_fission_bank = (settings::run_mode == RunMode::EIGENVALUE);
18,340,535✔
201

202
  // Counter for the number of fission sites successfully stored to the shared
203
  // fission bank or the secondary particle bank
204
  int n_sites_stored;
205

206
  for (n_sites_stored = 0; n_sites_stored < nu; n_sites_stored++) {
42,012,808✔
207
    // Initialize fission site object with particle data
208
    SourceSite site;
23,672,273✔
209
    site.r = p.r();
23,672,273✔
210
    site.particle = ParticleType::neutron;
23,672,273✔
211
    site.time = p.time();
23,672,273✔
212
    site.wgt = 1. / weight;
23,672,273✔
213
    site.surf_id = 0;
23,672,273✔
214

215
    // Sample delayed group and angle/energy for fission reaction
216
    sample_fission_neutron(i_nuclide, rx, &site, p);
23,672,273✔
217

218
    // Reject site if it exceeds time cutoff
219
    if (site.delayed_group > 0) {
23,672,273✔
220
      double t_cutoff = settings::time_cutoff[static_cast<int>(site.particle)];
147,911✔
221
      if (site.time > t_cutoff) {
147,911✔
NEW
222
        continue;
×
223
      }
224
    }
225

226
    // Set parent and progeny IDs
227
    site.parent_id = p.id();
23,672,273✔
228
    site.progeny_id = p.n_progeny()++;
23,672,273✔
229

230
    // Store fission site in bank
231
    if (use_fission_bank) {
23,672,273✔
232
      int64_t idx = simulation::fission_bank.thread_safe_append(site);
23,333,542✔
233
      if (idx == -1) {
23,333,542✔
234
        warning(
×
235
          "The shared fission bank is full. Additional fission sites created "
236
          "in this generation will not be banked. Results may be "
237
          "non-deterministic.");
238

239
        // Decrement number of particle progeny as storage was unsuccessful.
240
        // This step is needed so that the sum of all progeny is equal to the
241
        // size of the shared fission bank.
242
        p.n_progeny()--;
×
243

244
        // Break out of loop as no more sites can be added to fission bank
245
        break;
×
246
      }
247
      // Iterated Fission Probability (IFP) method
248
      if (settings::ifp_on) {
23,333,542✔
249
        ifp(p, site, idx);
218,295✔
250
      }
251
    } else {
252
      p.secondary_bank().push_back(site);
338,731✔
253
    }
254

255
    // Set the delayed group on the particle as well
256
    p.delayed_group() = site.delayed_group;
23,672,273✔
257

258
    // Increment the number of neutrons born delayed
259
    if (p.delayed_group() > 0) {
23,672,273✔
260
      nu_d[p.delayed_group() - 1]++;
147,911✔
261
    }
262

263
    // Write fission particles to nuBank
264
    NuBank& nu_bank_entry = p.nu_bank().emplace_back();
23,672,273✔
265
    nu_bank_entry.wgt = site.wgt;
23,672,273✔
266
    nu_bank_entry.E = site.E;
23,672,273✔
267
    nu_bank_entry.delayed_group = site.delayed_group;
23,672,273✔
268
  }
269

270
  // If shared fission bank was full, and no fissions could be added,
271
  // set the particle fission flag to false.
272
  if (n_sites_stored == 0) {
18,340,535✔
273
    p.fission() = false;
×
274
    return;
×
275
  }
276

277
  // Set nu to the number of fission sites successfully stored. If the fission
278
  // bank was not found to be full then these values are already equivalent.
279
  nu = n_sites_stored;
18,340,535✔
280

281
  // Store the total weight banked for analog fission tallies
282
  p.n_bank() = nu;
18,340,535✔
283
  p.wgt_bank() = nu / weight;
18,340,535✔
284
  for (size_t d = 0; d < MAX_DELAYED_GROUPS; d++) {
165,064,815✔
285
    p.n_delayed_bank(d) = nu_d[d];
146,724,280✔
286
  }
287
}
288

289
void sample_photon_reaction(Particle& p)
16,213,293✔
290
{
291
  // Kill photon if below energy cutoff -- an extra check is made here because
292
  // photons with energy below the cutoff may have been produced by neutrons
293
  // reactions or atomic relaxation
294
  int photon = static_cast<int>(ParticleType::photon);
16,213,293✔
295
  if (p.E() < settings::energy_cutoff[photon]) {
16,213,293✔
296
    p.E() = 0.0;
×
297
    p.wgt() = 0.0;
×
298
    return;
×
299
  }
300

301
  // Sample element within material
302
  int i_element = sample_element(p);
16,213,293✔
303
  const auto& micro {p.photon_xs(i_element)};
16,213,293✔
304
  const auto& element {*data::elements[i_element]};
16,213,293✔
305

306
  // Calculate photon energy over electron rest mass equivalent
307
  double alpha = p.E() / MASS_ELECTRON_EV;
16,213,293✔
308

309
  // For tallying purposes, this routine might be called directly. In that
310
  // case, we need to sample a reaction via the cutoff variable
311
  double prob = 0.0;
16,213,293✔
312
  double cutoff = prn(p.current_seed()) * micro.total;
16,213,293✔
313

314
  // Coherent (Rayleigh) scattering
315
  prob += micro.coherent;
16,213,293✔
316
  if (prob > cutoff) {
16,213,293✔
317
    p.mu() = element.rayleigh_scatter(alpha, p.current_seed());
894,955✔
318
    p.u() = rotate_angle(p.u(), p.mu(), nullptr, p.current_seed());
894,955✔
319
    p.event() = TallyEvent::SCATTER;
894,955✔
320
    p.event_mt() = COHERENT;
894,955✔
321
    return;
894,955✔
322
  }
323

324
  // Incoherent (Compton) scattering
325
  prob += micro.incoherent;
15,318,338✔
326
  if (prob > cutoff) {
15,318,338✔
327
    double alpha_out;
328
    int i_shell;
329
    element.compton_scatter(
20,643,184✔
330
      alpha, true, &alpha_out, &p.mu(), &i_shell, p.current_seed());
10,321,592✔
331

332
    // Determine binding energy of shell. The binding energy is 0.0 if
333
    // doppler broadening is not used.
334
    double e_b;
335
    if (i_shell == -1) {
10,321,592✔
336
      e_b = 0.0;
×
337
    } else {
338
      e_b = element.binding_energy_[i_shell];
10,321,592✔
339
    }
340

341
    // Create Compton electron
342
    double phi = uniform_distribution(0., 2.0 * PI, p.current_seed());
10,321,592✔
343
    double E_electron = (alpha - alpha_out) * MASS_ELECTRON_EV - e_b;
10,321,592✔
344
    int electron = static_cast<int>(ParticleType::electron);
10,321,592✔
345
    if (E_electron >= settings::energy_cutoff[electron]) {
10,321,592✔
346
      double mu_electron = (alpha - alpha_out * p.mu()) /
10,224,322✔
347
                           std::sqrt(alpha * alpha + alpha_out * alpha_out -
20,448,644✔
348
                                     2.0 * alpha * alpha_out * p.mu());
10,224,322✔
349
      Direction u = rotate_angle(p.u(), mu_electron, &phi, p.current_seed());
10,224,322✔
350
      p.create_secondary(p.wgt(), u, E_electron, ParticleType::electron);
10,224,322✔
351
    }
352

353
    // Allow electrons to fill orbital and produce Auger electrons and
354
    // fluorescent photons. Since Compton subshell data does not match atomic
355
    // relaxation data, use the mapping between the data to find the subshell
356
    if (i_shell >= 0 && element.subshell_map_[i_shell] >= 0) {
10,321,592✔
357
      element.atomic_relaxation(element.subshell_map_[i_shell], p);
10,321,592✔
358
    }
359

360
    phi += PI;
10,321,592✔
361
    p.E() = alpha_out * MASS_ELECTRON_EV;
10,321,592✔
362
    p.u() = rotate_angle(p.u(), p.mu(), &phi, p.current_seed());
10,321,592✔
363
    p.event() = TallyEvent::SCATTER;
10,321,592✔
364
    p.event_mt() = INCOHERENT;
10,321,592✔
365
    return;
10,321,592✔
366
  }
367

368
  // Photoelectric effect
369
  double prob_after = prob + micro.photoelectric;
4,996,746✔
370

371
  if (prob_after > cutoff) {
4,996,746✔
372
    // Get grid index, interpolation factor, and bounding subshell
373
    // cross sections
374
    int i_grid = micro.index_grid;
4,911,055✔
375
    double f = micro.interp_factor;
4,911,055✔
376
    const auto& xs_lower = xt::row(element.cross_sections_, i_grid);
4,911,055✔
377
    const auto& xs_upper = xt::row(element.cross_sections_, i_grid + 1);
4,911,055✔
378

379
    for (int i_shell = 0; i_shell < element.shells_.size(); ++i_shell) {
24,084,183✔
380
      const auto& shell {element.shells_[i_shell]};
24,084,183✔
381

382
      // Check threshold of reaction
383
      if (xs_lower(i_shell) == 0)
24,084,183✔
384
        continue;
10,260,132✔
385

386
      //  Evaluation subshell photoionization cross section
387
      prob += std::exp(
13,824,051✔
388
        xs_lower(i_shell) + f * (xs_upper(i_shell) - xs_lower(i_shell)));
13,824,051✔
389

390
      if (prob > cutoff) {
13,824,051✔
391
        // Determine binding energy based on whether atomic relaxation data is
392
        // present (if not, use value from Compton profile data)
393
        double binding_energy = element.has_atomic_relaxation_
4,911,055✔
394
                                  ? shell.binding_energy
4,911,055✔
395
                                  : element.binding_energy_[i_shell];
×
396

397
        // Determine energy of secondary electron
398
        double E_electron = p.E() - binding_energy;
4,911,055✔
399

400
        // Sample mu using non-relativistic Sauter distribution.
401
        // See Eqns 3.19 and 3.20 in "Implementing a photon physics
402
        // model in Serpent 2" by Toni Kaltiaisenaho
403
        double mu;
404
        while (true) {
405
          double r = prn(p.current_seed());
7,366,095✔
406
          if (4.0 * (1.0 - r) * r >= prn(p.current_seed())) {
7,366,095✔
407
            double rel_vel =
408
              std::sqrt(E_electron * (E_electron + 2.0 * MASS_ELECTRON_EV)) /
4,911,055✔
409
              (E_electron + MASS_ELECTRON_EV);
4,911,055✔
410
            mu =
4,911,055✔
411
              (2.0 * r + rel_vel - 1.0) / (2.0 * rel_vel * r - rel_vel + 1.0);
4,911,055✔
412
            break;
4,911,055✔
413
          }
414
        }
2,455,040✔
415

416
        double phi = uniform_distribution(0., 2.0 * PI, p.current_seed());
4,911,055✔
417
        Direction u;
4,911,055✔
418
        u.x = mu;
4,911,055✔
419
        u.y = std::sqrt(1.0 - mu * mu) * std::cos(phi);
4,911,055✔
420
        u.z = std::sqrt(1.0 - mu * mu) * std::sin(phi);
4,911,055✔
421

422
        // Create secondary electron
423
        p.create_secondary(p.wgt(), u, E_electron, ParticleType::electron);
4,911,055✔
424

425
        // Allow electrons to fill orbital and produce auger electrons
426
        // and fluorescent photons
427
        element.atomic_relaxation(i_shell, p);
4,911,055✔
428
        p.event() = TallyEvent::ABSORB;
4,911,055✔
429
        p.event_mt() = 533 + shell.index_subshell;
4,911,055✔
430
        p.wgt() = 0.0;
4,911,055✔
431
        p.E() = 0.0;
4,911,055✔
432
        return;
4,911,055✔
433
      }
434
    }
435
  }
9,822,110✔
436
  prob = prob_after;
85,691✔
437

438
  // Pair production
439
  prob += micro.pair_production;
85,691✔
440
  if (prob > cutoff) {
85,691✔
441
    double E_electron, E_positron;
442
    double mu_electron, mu_positron;
443
    element.pair_production(alpha, &E_electron, &E_positron, &mu_electron,
85,691✔
444
      &mu_positron, p.current_seed());
445

446
    // Create secondary electron
447
    Direction u = rotate_angle(p.u(), mu_electron, nullptr, p.current_seed());
85,691✔
448
    p.create_secondary(p.wgt(), u, E_electron, ParticleType::electron);
85,691✔
449

450
    // Create secondary positron
451
    u = rotate_angle(p.u(), mu_positron, nullptr, p.current_seed());
85,691✔
452
    p.create_secondary(p.wgt(), u, E_positron, ParticleType::positron);
85,691✔
453

454
    p.event() = TallyEvent::ABSORB;
85,691✔
455
    p.event_mt() = PAIR_PROD;
85,691✔
456
    p.wgt() = 0.0;
85,691✔
457
    p.E() = 0.0;
85,691✔
458
  }
459
}
460

461
void sample_electron_reaction(Particle& p)
51,847,521✔
462
{
463
  // TODO: create reaction types
464

465
  if (settings::electron_treatment == ElectronTreatment::TTB) {
51,847,521✔
466
    double E_lost;
467
    thick_target_bremsstrahlung(p, &E_lost);
51,844,859✔
468
  }
469

470
  p.E() = 0.0;
51,847,521✔
471
  p.wgt() = 0.0;
51,847,521✔
472
  p.event() = TallyEvent::ABSORB;
51,847,521✔
473
}
51,847,521✔
474

475
void sample_positron_reaction(Particle& p)
85,691✔
476
{
477
  // TODO: create reaction types
478

479
  if (settings::electron_treatment == ElectronTreatment::TTB) {
85,691✔
480
    double E_lost;
481
    thick_target_bremsstrahlung(p, &E_lost);
85,680✔
482
  }
483

484
  // Sample angle isotropically
485
  Direction u = isotropic_direction(p.current_seed());
85,691✔
486

487
  // Create annihilation photon pair traveling in opposite directions
488
  p.create_secondary(p.wgt(), u, MASS_ELECTRON_EV, ParticleType::photon);
85,691✔
489
  p.create_secondary(p.wgt(), -u, MASS_ELECTRON_EV, ParticleType::photon);
85,691✔
490

491
  p.E() = 0.0;
85,691✔
492
  p.wgt() = 0.0;
85,691✔
493
  p.event() = TallyEvent::ABSORB;
85,691✔
494
}
85,691✔
495

496
int sample_nuclide(Particle& p)
691,039,604✔
497
{
498
  // Sample cumulative distribution function
499
  double cutoff = prn(p.current_seed()) * p.macro_xs().total;
691,039,604✔
500

501
  // Get pointers to nuclide/density arrays
502
  const auto& mat {model::materials[p.material()]};
691,039,604✔
503
  int n = mat->nuclide_.size();
691,039,604✔
504

505
  double prob = 0.0;
691,039,604✔
506
  for (int i = 0; i < n; ++i) {
1,534,722,894✔
507
    // Get atom density
508
    int i_nuclide = mat->nuclide_[i];
1,534,722,894✔
509
    double atom_density = mat->atom_density(i, p.density_mult());
1,534,722,894✔
510

511
    // Increment probability to compare to cutoff
512
    prob += atom_density * p.neutron_xs(i_nuclide).total;
1,534,722,894✔
513
    if (prob >= cutoff)
1,534,722,894✔
514
      return i_nuclide;
691,039,604✔
515
  }
516

517
  // If we reach here, no nuclide was sampled
518
  p.write_restart();
×
519
  throw std::runtime_error {"Did not sample any nuclide during collision."};
×
520
}
521

522
int sample_element(Particle& p)
16,213,293✔
523
{
524
  // Sample cumulative distribution function
525
  double cutoff = prn(p.current_seed()) * p.macro_xs().total;
16,213,293✔
526

527
  // Get pointers to elements, densities
528
  const auto& mat {model::materials[p.material()]};
16,213,293✔
529

530
  double prob = 0.0;
16,213,293✔
531
  for (int i = 0; i < mat->element_.size(); ++i) {
39,362,126✔
532
    // Find atom density
533
    int i_element = mat->element_[i];
39,362,126✔
534
    double atom_density = mat->atom_density(i, p.density_mult());
39,362,126✔
535

536
    // Determine microscopic cross section
537
    double sigma = atom_density * p.photon_xs(i_element).total;
39,362,126✔
538

539
    // Increment probability to compare to cutoff
540
    prob += sigma;
39,362,126✔
541
    if (prob > cutoff) {
39,362,126✔
542
      // Save which nuclide particle had collision with for tally purpose
543
      p.event_nuclide() = mat->nuclide_[i];
16,213,293✔
544

545
      return i_element;
16,213,293✔
546
    }
547
  }
548

549
  // If we made it here, no element was sampled
550
  p.write_restart();
×
551
  fatal_error("Did not sample any element during collision.");
×
552
}
553

554
Reaction& sample_fission(int i_nuclide, Particle& p)
88,035,708✔
555
{
556
  // Get pointer to nuclide
557
  const auto& nuc {data::nuclides[i_nuclide]};
88,035,708✔
558

559
  // If we're in the URR, by default use the first fission reaction. We also
560
  // default to the first reaction if we know that there are no partial fission
561
  // reactions
562
  if (p.neutron_xs(i_nuclide).use_ptable || !nuc->has_partial_fission_) {
88,035,708✔
563
    return *nuc->fission_rx_[0];
88,010,071✔
564
  }
565

566
  // Check to see if we are in a windowed multipole range.  WMP only supports
567
  // the first fission reaction.
568
  if (nuc->multipole_) {
25,637✔
569
    if (p.E() >= nuc->multipole_->E_min_ && p.E() <= nuc->multipole_->E_max_) {
×
570
      return *nuc->fission_rx_[0];
×
571
    }
572
  }
573

574
  // Get grid index and interpolation factor and sample fission cdf
575
  const auto& micro = p.neutron_xs(i_nuclide);
25,637✔
576
  double cutoff = prn(p.current_seed()) * p.neutron_xs(i_nuclide).fission;
25,637✔
577
  double prob = 0.0;
25,637✔
578

579
  // Loop through each partial fission reaction type
580
  for (auto& rx : nuc->fission_rx_) {
25,688✔
581
    // add to cumulative probability
582
    prob += rx->xs(micro);
25,688✔
583

584
    // Create fission bank sites if fission occurs
585
    if (prob > cutoff)
25,688✔
586
      return *rx;
25,637✔
587
  }
588

589
  // If we reached here, no reaction was sampled
590
  throw std::runtime_error {
×
591
    "No fission reaction was sampled for " + nuc->name_};
×
592
}
593

594
void sample_photon_product(
1,337,644✔
595
  int i_nuclide, Particle& p, int* i_rx, int* i_product)
596
{
597
  // Get grid index and interpolation factor and sample photon production cdf
598
  const auto& micro = p.neutron_xs(i_nuclide);
1,337,644✔
599
  double cutoff = prn(p.current_seed()) * micro.photon_prod;
1,337,644✔
600
  double prob = 0.0;
1,337,644✔
601

602
  // Loop through each reaction type
603
  const auto& nuc {data::nuclides[i_nuclide]};
1,337,644✔
604
  for (int i = 0; i < nuc->reactions_.size(); ++i) {
21,995,831✔
605
    // Evaluate neutron cross section
606
    const auto& rx = nuc->reactions_[i];
21,995,831✔
607
    double xs = rx->xs(micro);
21,995,831✔
608

609
    // if cross section is zero for this reaction, skip it
610
    if (xs == 0.0)
21,995,831✔
611
      continue;
7,191,415✔
612

613
    for (int j = 0; j < rx->products_.size(); ++j) {
69,577,497✔
614
      if (rx->products_[j].particle_ == ParticleType::photon) {
56,110,725✔
615
        // For fission, artificially increase the photon yield to account
616
        // for delayed photons
617
        double f = 1.0;
43,656,910✔
618
        if (settings::delayed_photon_scaling) {
43,656,910✔
619
          if (is_fission(rx->mt_)) {
43,656,910✔
620
            if (nuc->prompt_photons_ && nuc->delayed_photons_) {
540,221✔
621
              double energy_prompt = (*nuc->prompt_photons_)(p.E());
540,221✔
622
              double energy_delayed = (*nuc->delayed_photons_)(p.E());
540,221✔
623
              f = (energy_prompt + energy_delayed) / (energy_prompt);
540,221✔
624
            }
625
          }
626
        }
627

628
        // add to cumulative probability
629
        prob += f * (*rx->products_[j].yield_)(p.E()) * xs;
43,656,910✔
630

631
        *i_rx = i;
43,656,910✔
632
        *i_product = j;
43,656,910✔
633
        if (prob > cutoff)
43,656,910✔
634
          return;
1,337,644✔
635
      }
636
    }
637
  }
638
}
639

640
void absorption(Particle& p, int i_nuclide)
691,036,079✔
641
{
642
  if (settings::survival_biasing) {
691,036,079✔
643
    // Determine weight absorbed in survival biasing
644
    const double wgt_absorb = p.wgt() * p.neutron_xs(i_nuclide).absorption /
499,950✔
645
                              p.neutron_xs(i_nuclide).total;
499,950✔
646

647
    // Adjust weight of particle by probability of absorption
648
    p.wgt() -= wgt_absorb;
499,950✔
649

650
    // Score implicit absorption estimate of keff
651
    if (settings::run_mode == RunMode::EIGENVALUE) {
499,950✔
652
      p.keff_tally_absorption() += wgt_absorb *
499,950✔
653
                                   p.neutron_xs(i_nuclide).nu_fission /
499,950✔
654
                                   p.neutron_xs(i_nuclide).absorption;
499,950✔
655
    }
656
  } else {
657
    // See if disappearance reaction happens
658
    if (p.neutron_xs(i_nuclide).absorption >
690,536,129✔
659
        prn(p.current_seed()) * p.neutron_xs(i_nuclide).total) {
690,536,129✔
660
      // Score absorption estimate of keff
661
      if (settings::run_mode == RunMode::EIGENVALUE) {
17,151,039✔
662
        p.keff_tally_absorption() += p.wgt() *
28,966,914✔
663
                                     p.neutron_xs(i_nuclide).nu_fission /
14,483,457✔
664
                                     p.neutron_xs(i_nuclide).absorption;
14,483,457✔
665
      }
666

667
      p.wgt() = 0.0;
17,151,039✔
668
      p.event() = TallyEvent::ABSORB;
17,151,039✔
669
      p.event_mt() = N_DISAPPEAR;
17,151,039✔
670
    }
671
  }
672
}
691,036,079✔
673

674
void scatter(Particle& p, int i_nuclide)
673,729,736✔
675
{
676
  // copy incoming direction
677
  Direction u_old {p.u()};
673,729,736✔
678

679
  // Get pointer to nuclide and grid index/interpolation factor
680
  const auto& nuc {data::nuclides[i_nuclide]};
673,729,736✔
681
  const auto& micro {p.neutron_xs(i_nuclide)};
673,729,736✔
682
  int i_temp = micro.index_temp;
673,729,736✔
683

684
  // For tallying purposes, this routine might be called directly. In that
685
  // case, we need to sample a reaction via the cutoff variable
686
  double cutoff = prn(p.current_seed()) * (micro.total - micro.absorption);
673,729,736✔
687
  bool sampled = false;
673,729,736✔
688

689
  // Calculate elastic cross section if it wasn't precalculated
690
  if (micro.elastic == CACHE_INVALID) {
673,729,736✔
691
    nuc->calculate_elastic_xs(p);
564,308,449✔
692
  }
693

694
  double prob = micro.elastic - micro.thermal;
673,729,736✔
695
  if (prob > cutoff) {
673,729,736✔
696
    // =======================================================================
697
    // NON-S(A,B) ELASTIC SCATTERING
698

699
    // Determine temperature
700
    double kT = nuc->multipole_ ? p.sqrtkT() * p.sqrtkT() : nuc->kTs_[i_temp];
584,147,375✔
701

702
    // Perform collision physics for elastic scattering
703
    elastic_scatter(i_nuclide, *nuc->reactions_[0], kT, p);
584,147,375✔
704

705
    p.event_mt() = ELASTIC;
584,147,375✔
706
    sampled = true;
584,147,375✔
707
  }
708

709
  prob = micro.elastic;
673,729,736✔
710
  if (prob > cutoff && !sampled) {
673,729,736✔
711
    // =======================================================================
712
    // S(A,B) SCATTERING
713

714
    sab_scatter(i_nuclide, micro.index_sab, p);
75,155,723✔
715

716
    p.event_mt() = ELASTIC;
75,155,723✔
717
    sampled = true;
75,155,723✔
718
  }
719

720
  if (!sampled) {
673,729,736✔
721
    // =======================================================================
722
    // INELASTIC SCATTERING
723

724
    int n = nuc->index_inelastic_scatter_.size();
14,426,638✔
725
    int i = 0;
14,426,638✔
726
    for (int j = 0; j < n && prob < cutoff; ++j) {
238,591,737✔
727
      i = nuc->index_inelastic_scatter_[j];
224,165,099✔
728

729
      // add to cumulative probability
730
      prob += nuc->reactions_[i]->xs(micro);
224,165,099✔
731
    }
732

733
    // Perform collision physics for inelastic scattering
734
    const auto& rx {nuc->reactions_[i]};
14,426,638✔
735
    inelastic_scatter(*nuc, *rx, p);
14,426,638✔
736
    p.event_mt() = rx->mt_;
14,426,638✔
737
  }
738

739
  // Set event component
740
  p.event() = TallyEvent::SCATTER;
673,729,736✔
741

742
  // Sample new outgoing angle for isotropic-in-lab scattering
743
  const auto& mat {model::materials[p.material()]};
673,729,736✔
744
  if (!mat->p0_.empty()) {
673,729,736✔
745
    int i_nuc_mat = mat->mat_nuclide_index_[i_nuclide];
326,370✔
746
    if (mat->p0_[i_nuc_mat]) {
326,370✔
747
      // Sample isotropic-in-lab outgoing direction
748
      p.u() = isotropic_direction(p.current_seed());
326,370✔
749
      p.mu() = u_old.dot(p.u());
326,370✔
750
    }
751
  }
752
}
673,729,736✔
753

754
void elastic_scatter(int i_nuclide, const Reaction& rx, double kT, Particle& p)
584,147,375✔
755
{
756
  // get pointer to nuclide
757
  const auto& nuc {data::nuclides[i_nuclide]};
584,147,375✔
758

759
  double vel = std::sqrt(p.E());
584,147,375✔
760
  double awr = nuc->awr_;
584,147,375✔
761

762
  // Neutron velocity in LAB
763
  Direction v_n = vel * p.u();
584,147,375✔
764

765
  // Sample velocity of target nucleus
766
  Direction v_t {};
584,147,375✔
767
  if (!p.neutron_xs(i_nuclide).use_ptable) {
584,147,375✔
768
    v_t = sample_target_velocity(*nuc, p.E(), p.u(), v_n,
1,131,493,304✔
769
      p.neutron_xs(i_nuclide).elastic, kT, p.current_seed());
565,746,652✔
770
  }
771

772
  // Velocity of center-of-mass
773
  Direction v_cm = (v_n + awr * v_t) / (awr + 1.0);
584,147,375✔
774

775
  // Transform to CM frame
776
  v_n -= v_cm;
584,147,375✔
777

778
  // Find speed of neutron in CM
779
  vel = v_n.norm();
584,147,375✔
780

781
  // Sample scattering angle, checking if angle distribution is present (assume
782
  // isotropic otherwise)
783
  double mu_cm;
784
  auto& d = rx.products_[0].distribution_[0];
584,147,375✔
785
  auto d_ = dynamic_cast<UncorrelatedAngleEnergy*>(d.get());
584,147,375✔
786
  if (!d_->angle().empty()) {
584,147,375✔
787
    mu_cm = d_->angle().sample(p.E(), p.current_seed());
584,147,375✔
788
  } else {
789
    mu_cm = uniform_distribution(-1., 1., p.current_seed());
×
790
  }
791

792
  // Determine direction cosines in CM
793
  Direction u_cm = v_n / vel;
584,147,375✔
794

795
  // Rotate neutron velocity vector to new angle -- note that the speed of the
796
  // neutron in CM does not change in elastic scattering. However, the speed
797
  // will change when we convert back to LAB
798
  v_n = vel * rotate_angle(u_cm, mu_cm, nullptr, p.current_seed());
584,147,375✔
799

800
  // Transform back to LAB frame
801
  v_n += v_cm;
584,147,375✔
802

803
  p.E() = v_n.dot(v_n);
584,147,375✔
804
  vel = std::sqrt(p.E());
584,147,375✔
805

806
  // compute cosine of scattering angle in LAB frame by taking dot product of
807
  // neutron's pre- and post-collision angle
808
  p.mu() = p.u().dot(v_n) / vel;
584,147,375✔
809

810
  // Set energy and direction of particle in LAB frame
811
  p.u() = v_n / vel;
584,147,375✔
812

813
  // Because of floating-point roundoff, it may be possible for mu_lab to be
814
  // outside of the range [-1,1). In these cases, we just set mu_lab to exactly
815
  // -1 or 1
816
  if (std::abs(p.mu()) > 1.0)
584,147,375✔
817
    p.mu() = std::copysign(1.0, p.mu());
×
818
}
584,147,375✔
819

820
void sab_scatter(int i_nuclide, int i_sab, Particle& p)
75,155,723✔
821
{
822
  // Determine temperature index
823
  const auto& micro {p.neutron_xs(i_nuclide)};
75,155,723✔
824
  int i_temp = micro.index_temp_sab;
75,155,723✔
825

826
  // Sample energy and angle
827
  double E_out;
828
  data::thermal_scatt[i_sab]->data_[i_temp].sample(
150,311,446✔
829
    micro, p.E(), &E_out, &p.mu(), p.current_seed());
75,155,723✔
830

831
  // Set energy to outgoing, change direction of particle
832
  p.E() = E_out;
75,155,723✔
833
  p.u() = rotate_angle(p.u(), p.mu(), nullptr, p.current_seed());
75,155,723✔
834
}
75,155,723✔
835

836
Direction sample_target_velocity(const Nuclide& nuc, double E, Direction u,
565,746,652✔
837
  Direction v_neut, double xs_eff, double kT, uint64_t* seed)
838
{
839
  // check if nuclide is a resonant scatterer
840
  ResScatMethod sampling_method;
841
  if (nuc.resonant_) {
565,746,652✔
842

843
    // sampling method to use
844
    sampling_method = settings::res_scat_method;
84,557✔
845

846
    // upper resonance scattering energy bound (target is at rest above this E)
847
    if (E > settings::res_scat_energy_max) {
84,557✔
848
      return {};
40,755✔
849

850
      // lower resonance scattering energy bound (should be no resonances below)
851
    } else if (E < settings::res_scat_energy_min) {
43,802✔
852
      sampling_method = ResScatMethod::cxs;
24,992✔
853
    }
854

855
    // otherwise, use free gas model
856
  } else {
857
    if (E >= FREE_GAS_THRESHOLD * kT && nuc.awr_ > 1.0) {
565,662,095✔
858
      return {};
227,656,947✔
859
    } else {
860
      sampling_method = ResScatMethod::cxs;
338,005,148✔
861
    }
862
  }
863

864
  // use appropriate target velocity sampling method
865
  switch (sampling_method) {
338,048,950✔
866
  case ResScatMethod::cxs:
338,030,140✔
867

868
    // sample target velocity with the constant cross section (cxs) approx.
869
    return sample_cxs_target_velocity(nuc.awr_, E, u, kT, seed);
338,030,140✔
870

871
  case ResScatMethod::dbrc:
18,810✔
872
  case ResScatMethod::rvs: {
873
    double E_red = std::sqrt(nuc.awr_ * E / kT);
18,810✔
874
    double E_low = std::pow(std::max(0.0, E_red - 4.0), 2) * kT / nuc.awr_;
18,810✔
875
    double E_up = (E_red + 4.0) * (E_red + 4.0) * kT / nuc.awr_;
18,810✔
876

877
    // find lower and upper energy bound indices
878
    // lower index
879
    int i_E_low;
880
    if (E_low < nuc.energy_0K_.front()) {
18,810✔
881
      i_E_low = 0;
×
882
    } else if (E_low > nuc.energy_0K_.back()) {
18,810✔
883
      i_E_low = nuc.energy_0K_.size() - 2;
×
884
    } else {
885
      i_E_low =
18,810✔
886
        lower_bound_index(nuc.energy_0K_.begin(), nuc.energy_0K_.end(), E_low);
18,810✔
887
    }
888

889
    // upper index
890
    int i_E_up;
891
    if (E_up < nuc.energy_0K_.front()) {
18,810✔
892
      i_E_up = 0;
×
893
    } else if (E_up > nuc.energy_0K_.back()) {
18,810✔
894
      i_E_up = nuc.energy_0K_.size() - 2;
×
895
    } else {
896
      i_E_up =
18,810✔
897
        lower_bound_index(nuc.energy_0K_.begin(), nuc.energy_0K_.end(), E_up);
18,810✔
898
    }
899

900
    if (i_E_up == i_E_low) {
18,810✔
901
      // Handle degenerate case -- if the upper/lower bounds occur for the same
902
      // index, then using cxs is probably a good approximation
903
      return sample_cxs_target_velocity(nuc.awr_, E, u, kT, seed);
18,810✔
904
    }
905

906
    if (sampling_method == ResScatMethod::dbrc) {
15,532✔
907
      // interpolate xs since we're not exactly at the energy indices
908
      double xs_low = nuc.elastic_0K_[i_E_low];
×
909
      double m = (nuc.elastic_0K_[i_E_low + 1] - xs_low) /
×
910
                 (nuc.energy_0K_[i_E_low + 1] - nuc.energy_0K_[i_E_low]);
×
911
      xs_low += m * (E_low - nuc.energy_0K_[i_E_low]);
×
912
      double xs_up = nuc.elastic_0K_[i_E_up];
×
913
      m = (nuc.elastic_0K_[i_E_up + 1] - xs_up) /
×
914
          (nuc.energy_0K_[i_E_up + 1] - nuc.energy_0K_[i_E_up]);
×
915
      xs_up += m * (E_up - nuc.energy_0K_[i_E_up]);
×
916

917
      // get max 0K xs value over range of practical relative energies
918
      double xs_max = *std::max_element(
×
919
        &nuc.elastic_0K_[i_E_low + 1], &nuc.elastic_0K_[i_E_up + 1]);
×
920
      xs_max = std::max({xs_low, xs_max, xs_up});
×
921

922
      while (true) {
923
        double E_rel;
924
        Direction v_target;
×
925
        while (true) {
926
          // sample target velocity with the constant cross section (cxs)
927
          // approx.
928
          v_target = sample_cxs_target_velocity(nuc.awr_, E, u, kT, seed);
×
929
          Direction v_rel = v_neut - v_target;
×
930
          E_rel = v_rel.dot(v_rel);
×
931
          if (E_rel < E_up)
×
932
            break;
×
933
        }
934

935
        // perform Doppler broadening rejection correction (dbrc)
936
        double xs_0K = nuc.elastic_xs_0K(E_rel);
×
937
        double R = xs_0K / xs_max;
×
938
        if (prn(seed) < R)
×
939
          return v_target;
×
940
      }
941

942
    } else if (sampling_method == ResScatMethod::rvs) {
15,532✔
943
      // interpolate xs CDF since we're not exactly at the energy indices
944
      // cdf value at lower bound attainable energy
945
      double cdf_low = 0.0;
15,532✔
946
      if (E_low > nuc.energy_0K_.front()) {
15,532✔
947
        double m = (nuc.xs_cdf_[i_E_low + 1] - nuc.xs_cdf_[i_E_low]) /
15,532✔
948
                   (nuc.energy_0K_[i_E_low + 1] - nuc.energy_0K_[i_E_low]);
15,532✔
949
        cdf_low = nuc.xs_cdf_[i_E_low] + m * (E_low - nuc.energy_0K_[i_E_low]);
15,532✔
950
      }
951

952
      // cdf value at upper bound attainable energy
953
      double m = (nuc.xs_cdf_[i_E_up + 1] - nuc.xs_cdf_[i_E_up]) /
15,532✔
954
                 (nuc.energy_0K_[i_E_up + 1] - nuc.energy_0K_[i_E_up]);
15,532✔
955
      double cdf_up = nuc.xs_cdf_[i_E_up] + m * (E_up - nuc.energy_0K_[i_E_up]);
15,532✔
956

957
      while (true) {
958
        // directly sample Maxwellian
959
        double E_t = -kT * std::log(prn(seed));
170,720✔
960

961
        // sample a relative energy using the xs cdf
962
        double cdf_rel = cdf_low + prn(seed) * (cdf_up - cdf_low);
170,720✔
963
        int i_E_rel = lower_bound_index(nuc.xs_cdf_.begin() + i_E_low,
170,720✔
964
          nuc.xs_cdf_.begin() + i_E_up + 2, cdf_rel);
170,720✔
965
        double E_rel = nuc.energy_0K_[i_E_low + i_E_rel];
170,720✔
966
        double m = (nuc.xs_cdf_[i_E_low + i_E_rel + 1] -
170,720✔
967
                     nuc.xs_cdf_[i_E_low + i_E_rel]) /
170,720✔
968
                   (nuc.energy_0K_[i_E_low + i_E_rel + 1] -
170,720✔
969
                     nuc.energy_0K_[i_E_low + i_E_rel]);
170,720✔
970
        E_rel += (cdf_rel - nuc.xs_cdf_[i_E_low + i_E_rel]) / m;
170,720✔
971

972
        // perform rejection sampling on cosine between
973
        // neutron and target velocities
974
        double mu = (E_t + nuc.awr_ * (E - E_rel)) /
170,720✔
975
                    (2.0 * std::sqrt(nuc.awr_ * E * E_t));
170,720✔
976

977
        if (std::abs(mu) < 1.0) {
170,720✔
978
          // set and accept target velocity
979
          E_t /= nuc.awr_;
15,532✔
980
          return std::sqrt(E_t) * rotate_angle(u, mu, nullptr, seed);
15,532✔
981
        }
982
      }
155,188✔
983
    }
984
  } // case RVS, DBRC
985
  } // switch (sampling_method)
986

987
  UNREACHABLE();
×
988
}
989

990
Direction sample_cxs_target_velocity(
338,033,418✔
991
  double awr, double E, Direction u, double kT, uint64_t* seed)
992
{
993
  double beta_vn = std::sqrt(awr * E / kT);
338,033,418✔
994
  double alpha = 1.0 / (1.0 + std::sqrt(PI) * beta_vn / 2.0);
338,033,418✔
995

996
  double beta_vt_sq;
997
  double mu;
998
  while (true) {
999
    // Sample two random numbers
1000
    double r1 = prn(seed);
400,271,398✔
1001
    double r2 = prn(seed);
400,271,398✔
1002

1003
    if (prn(seed) < alpha) {
400,271,398✔
1004
      // With probability alpha, we sample the distribution p(y) =
1005
      // y*e^(-y). This can be done with sampling scheme C45 from the Monte
1006
      // Carlo sampler
1007

1008
      beta_vt_sq = -std::log(r1 * r2);
101,144,688✔
1009

1010
    } else {
1011
      // With probability 1-alpha, we sample the distribution p(y) = y^2 *
1012
      // e^(-y^2). This can be done with sampling scheme C61 from the Monte
1013
      // Carlo sampler
1014

1015
      double c = std::cos(PI / 2.0 * prn(seed));
299,126,710✔
1016
      beta_vt_sq = -std::log(r1) - std::log(r2) * c * c;
299,126,710✔
1017
    }
1018

1019
    // Determine beta * vt
1020
    double beta_vt = std::sqrt(beta_vt_sq);
400,271,398✔
1021

1022
    // Sample cosine of angle between neutron and target velocity
1023
    mu = uniform_distribution(-1., 1., seed);
400,271,398✔
1024

1025
    // Determine rejection probability
1026
    double accept_prob =
1027
      std::sqrt(beta_vn * beta_vn + beta_vt_sq - 2 * beta_vn * beta_vt * mu) /
400,271,398✔
1028
      (beta_vn + beta_vt);
400,271,398✔
1029

1030
    // Perform rejection sampling on vt and mu
1031
    if (prn(seed) < accept_prob)
400,271,398✔
1032
      break;
338,033,418✔
1033
  }
62,237,980✔
1034

1035
  // Determine speed of target nucleus
1036
  double vt = std::sqrt(beta_vt_sq * kT / awr);
338,033,418✔
1037

1038
  // Determine velocity vector of target nucleus based on neutron's velocity
1039
  // and the sampled angle between them
1040
  return vt * rotate_angle(u, mu, nullptr, seed);
338,033,418✔
1041
}
1042

1043
void sample_fission_neutron(
23,672,273✔
1044
  int i_nuclide, const Reaction& rx, SourceSite* site, Particle& p)
1045
{
1046
  // Get attributes of particle
1047
  double E_in = p.E();
23,672,273✔
1048
  uint64_t* seed = p.current_seed();
23,672,273✔
1049

1050
  // Determine total nu, delayed nu, and delayed neutron fraction
1051
  const auto& nuc {data::nuclides[i_nuclide]};
23,672,273✔
1052
  double nu_t = nuc->nu(E_in, Nuclide::EmissionMode::total);
23,672,273✔
1053
  double nu_d = nuc->nu(E_in, Nuclide::EmissionMode::delayed);
23,672,273✔
1054
  double beta = nu_d / nu_t;
23,672,273✔
1055

1056
  if (prn(seed) < beta) {
23,672,273✔
1057
    // ====================================================================
1058
    // DELAYED NEUTRON SAMPLED
1059

1060
    // sampled delayed precursor group
1061
    double xi = prn(seed) * nu_d;
147,911✔
1062
    double prob = 0.0;
147,911✔
1063
    int group;
1064
    for (group = 1; group < nuc->n_precursor_; ++group) {
552,516✔
1065
      // determine delayed neutron precursor yield for group j
1066
      double yield = (*rx.products_[group].yield_)(E_in);
542,247✔
1067

1068
      // Check if this group is sampled
1069
      prob += yield;
542,247✔
1070
      if (xi < prob)
542,247✔
1071
        break;
137,642✔
1072
    }
1073

1074
    // if the sum of the probabilities is slightly less than one and the
1075
    // random number is greater, j will be greater than nuc %
1076
    // n_precursor -- check for this condition
1077
    group = std::min(group, nuc->n_precursor_);
147,911✔
1078

1079
    // set the delayed group for the particle born from fission
1080
    site->delayed_group = group;
147,911✔
1081

1082
    // Sample time of emission based on decay constant of precursor
1083
    double decay_rate = rx.products_[site->delayed_group].decay_rate_;
147,911✔
1084
    site->time -= std::log(prn(p.current_seed())) / decay_rate;
147,911✔
1085

1086
  } else {
1087
    // ====================================================================
1088
    // PROMPT NEUTRON SAMPLED
1089

1090
    // set the delayed group for the particle born from fission to 0
1091
    site->delayed_group = 0;
23,524,362✔
1092
  }
1093

1094
  // sample from prompt neutron energy distribution
1095
  int n_sample = 0;
23,672,273✔
1096
  double mu;
1097
  while (true) {
1098
    rx.products_[site->delayed_group].sample(E_in, site->E, mu, seed);
23,672,273✔
1099

1100
    // resample if energy is greater than maximum neutron energy
1101
    constexpr int neutron = static_cast<int>(ParticleType::neutron);
23,672,273✔
1102
    if (site->E < data::energy_max[neutron])
23,672,273✔
1103
      break;
23,672,273✔
1104

1105
    // check for large number of resamples
1106
    ++n_sample;
×
1107
    if (n_sample == MAX_SAMPLE) {
×
1108
      // particle_write_restart(p)
1109
      fatal_error("Resampled energy distribution maximum number of times "
×
1110
                  "for nuclide " +
×
1111
                  nuc->name_);
×
1112
    }
1113
  }
1114

1115
  // Sample azimuthal angle uniformly in [0, 2*pi) and assign angle
1116
  site->u = rotate_angle(p.u(), mu, nullptr, seed);
23,672,273✔
1117
}
23,672,273✔
1118

1119
void inelastic_scatter(const Nuclide& nuc, const Reaction& rx, Particle& p)
14,426,638✔
1120
{
1121
  // copy energy of neutron
1122
  double E_in = p.E();
14,426,638✔
1123

1124
  // sample outgoing energy and scattering cosine
1125
  double E;
1126
  double mu;
1127
  rx.products_[0].sample(E_in, E, mu, p.current_seed());
14,426,638✔
1128

1129
  // if scattering system is in center-of-mass, transfer cosine of scattering
1130
  // angle and outgoing energy from CM to LAB
1131
  if (rx.scatter_in_cm_) {
14,426,638✔
1132
    double E_cm = E;
14,376,726✔
1133

1134
    // determine outgoing energy in lab
1135
    double A = nuc.awr_;
14,376,726✔
1136
    E = E_cm + (E_in + 2.0 * mu * (A + 1.0) * std::sqrt(E_in * E_cm)) /
14,376,726✔
1137
                 ((A + 1.0) * (A + 1.0));
14,376,726✔
1138

1139
    // determine outgoing angle in lab
1140
    mu = mu * std::sqrt(E_cm / E) + 1.0 / (A + 1.0) * std::sqrt(E_in / E);
14,376,726✔
1141
  }
1142

1143
  // Because of floating-point roundoff, it may be possible for mu to be
1144
  // outside of the range [-1,1). In these cases, we just set mu to exactly -1
1145
  // or 1
1146
  if (std::abs(mu) > 1.0)
14,426,638✔
1147
    mu = std::copysign(1.0, mu);
×
1148

1149
  // Set outgoing energy and scattering angle
1150
  p.E() = E;
14,426,638✔
1151
  p.mu() = mu;
14,426,638✔
1152

1153
  // change direction of particle
1154
  p.u() = rotate_angle(p.u(), mu, nullptr, p.current_seed());
14,426,638✔
1155

1156
  // evaluate yield
1157
  double yield = (*rx.products_[0].yield_)(E_in);
14,426,638✔
1158
  if (std::floor(yield) == yield && yield > 0) {
14,426,638✔
1159
    // If yield is integral, create exactly that many secondary particles
1160
    for (int i = 0; i < static_cast<int>(std::round(yield)) - 1; ++i) {
14,493,876✔
1161
      p.create_secondary(p.wgt(), p.u(), p.E(), ParticleType::neutron);
67,314✔
1162
    }
1163
  } else {
14,426,562✔
1164
    // Otherwise, change weight of particle based on yield
1165
    p.wgt() *= yield;
76✔
1166
  }
1167
}
14,426,638✔
1168

1169
void sample_secondary_photons(Particle& p, int i_nuclide)
8,809,548✔
1170
{
1171
  // Sample the number of photons produced
1172
  double y_t =
1173
    p.neutron_xs(i_nuclide).photon_prod / p.neutron_xs(i_nuclide).total;
8,809,548✔
1174
  double photon_wgt = p.wgt();
8,809,548✔
1175
  int y = 1;
8,809,548✔
1176

1177
  if (settings::use_decay_photons) {
8,809,548✔
1178
    // For decay photons, sample a single photon and modify the weight
1179
    if (y_t <= 0.0)
72,006✔
1180
      return;
17,281✔
1181
    photon_wgt *= y_t;
54,725✔
1182
  } else {
1183
    // For prompt photons, sample an integral number of photons with weight
1184
    // equal to the neutron's weight
1185
    y = static_cast<int>(y_t);
8,737,542✔
1186
    if (prn(p.current_seed()) <= y_t - y)
8,737,542✔
1187
      ++y;
574,860✔
1188
  }
1189

1190
  // Sample each secondary photon
1191
  for (int i = 0; i < y; ++i) {
10,129,911✔
1192
    // Sample the reaction and product
1193
    int i_rx;
1194
    int i_product;
1195
    sample_photon_product(i_nuclide, p, &i_rx, &i_product);
1,337,644✔
1196

1197
    // Sample the outgoing energy and angle
1198
    auto& rx = data::nuclides[i_nuclide]->reactions_[i_rx];
1,337,644✔
1199
    double E;
1200
    double mu;
1201
    rx->products_[i_product].sample(p.E(), E, mu, p.current_seed());
1,337,644✔
1202

1203
    // Sample the new direction
1204
    Direction u = rotate_angle(p.u(), mu, nullptr, p.current_seed());
1,337,644✔
1205

1206
    // In a k-eigenvalue simulation, it's necessary to provide higher weight to
1207
    // secondary photons from non-fission reactions to properly balance energy
1208
    // release and deposition. See D. P. Griesheimer, S. J. Douglass, and M. H.
1209
    // Stedry, "Self-consistent energy normalization for quasistatic reactor
1210
    // calculations", Proc. PHYSOR, Cambridge, UK, Mar 29-Apr 2, 2020.
1211
    double wgt = photon_wgt;
1,337,644✔
1212
    if (settings::run_mode == RunMode::EIGENVALUE && !is_fission(rx->mt_)) {
1,337,644✔
1213
      wgt *= simulation::keff;
348,128✔
1214
    }
1215

1216
    // Create the secondary photon
1217
    bool created_photon = p.create_secondary(wgt, u, E, ParticleType::photon);
1,337,644✔
1218

1219
    // Tag secondary particle with parent nuclide
1220
    if (created_photon && settings::use_decay_photons) {
1,337,644✔
1221
      p.secondary_bank().back().parent_nuclide =
52,844✔
1222
        rx->products_[i_product].parent_nuclide_;
52,844✔
1223
    }
1224
  }
1225
}
1226

1227
} // namespace openmc
STATUS · Troubleshooting · Open an Issue · Sales · Support · CAREERS · ENTERPRISE · START FREE · SCHEDULE DEMO
ANNOUNCEMENTS · TWITTER · TOS & SLA · Supported CI Services · What's a CI service? · Automated Testing

© 2025 Coveralls, Inc