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

openmc-dev / openmc / 25930573650

15 May 2026 05:01PM UTC coverage: 81.375% (+0.05%) from 81.326%
25930573650

Pull #3863

github

web-flow
Merge 95bd57fc1 into d56cda254
Pull Request #3863: Shared Secondary Particle Bank

17950 of 25871 branches covered (69.38%)

Branch coverage included in aggregate %.

407 of 417 new or added lines in 17 files covered. (97.6%)

1464 existing lines in 34 files now uncovered.

59095 of 68808 relevant lines covered (85.88%)

48517262.56 hits per line

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

72.33
/src/random_ray/flat_source_domain.cpp
1
#include "openmc/random_ray/flat_source_domain.h"
2

3
#include "openmc/cell.h"
4
#include "openmc/constants.h"
5
#include "openmc/eigenvalue.h"
6
#include "openmc/geometry.h"
7
#include "openmc/material.h"
8
#include "openmc/message_passing.h"
9
#include "openmc/mgxs_interface.h"
10
#include "openmc/output.h"
11
#include "openmc/plot.h"
12
#include "openmc/random_ray/random_ray.h"
13
#include "openmc/simulation.h"
14
#include "openmc/tallies/filter.h"
15
#include "openmc/tallies/tally.h"
16
#include "openmc/tallies/tally_scoring.h"
17
#include "openmc/timer.h"
18
#include "openmc/weight_windows.h"
19

20
#include <cstdio>
21
#include <numeric>
22

23
namespace openmc {
24

25
//==============================================================================
26
// FlatSourceDomain implementation
27
//==============================================================================
28

29
// Static Variable Declarations
30
RandomRayVolumeEstimator FlatSourceDomain::volume_estimator_ {
31
  RandomRayVolumeEstimator::HYBRID};
32
bool FlatSourceDomain::volume_normalized_flux_tallies_ {false};
33
bool FlatSourceDomain::adjoint_ {false};
34
bool FlatSourceDomain::fw_cadis_local_ {false};
35
double FlatSourceDomain::diagonal_stabilization_rho_ {1.0};
36
std::unordered_map<int, vector<std::pair<Source::DomainType, int>>>
37
  FlatSourceDomain::mesh_domain_map_;
38
std::vector<size_t> FlatSourceDomain::fw_cadis_local_targets_;
39

40
FlatSourceDomain::FlatSourceDomain() : negroups_(data::mg.num_energy_groups_)
792✔
41
{
42
  // Count the number of source regions, compute the cell offset
43
  // indices, and store the material type The reason for the offsets is that
44
  // some cell types may not have material fills, and therefore do not
45
  // produce FSRs. Thus, we cannot index into the global arrays directly
46
  int base_source_regions = 0;
792✔
47
  for (const auto& c : model::cells) {
6,285✔
48
    if (c->type_ != Fill::MATERIAL) {
5,493✔
49
      source_region_offsets_.push_back(-1);
2,425✔
50
    } else {
51
      source_region_offsets_.push_back(base_source_regions);
3,068✔
52
      base_source_regions += c->n_instances();
3,068✔
53
    }
54
  }
55

56
  // Initialize source regions.
57
  bool is_linear = RandomRay::source_shape_ != RandomRaySourceShape::FLAT;
792✔
58
  source_regions_ = SourceRegionContainer(negroups_, is_linear);
792✔
59

60
  // Initialize tally volumes
61
  if (volume_normalized_flux_tallies_) {
792✔
62
    tally_volumes_.resize(model::tallies.size());
506✔
63
    for (int i = 0; i < model::tallies.size(); i++) {
1,837✔
64
      //  Get the shape of the 3D result tensor
65
      auto shape = model::tallies[i]->results().shape();
1,331✔
66

67
      // Create a new 2D tensor with the same size as the first
68
      // two dimensions of the 3D tensor
69
      tally_volumes_[i] = tensor::Tensor<double>({shape[0], shape[1]});
2,662✔
70
    }
1,331✔
71
  }
72

73
  // Compute simulation domain volume based on ray source
74
  auto* is = dynamic_cast<IndependentSource*>(RandomRay::ray_source_.get());
792!
75
  SpatialDistribution* space_dist = is->space();
792!
76
  SpatialBox* sb = dynamic_cast<SpatialBox*>(space_dist);
792!
77
  Position dims = sb->upper_right() - sb->lower_left();
792✔
78
  simulation_volume_ = dims.x * dims.y * dims.z;
792✔
79
}
792✔
80

81
void FlatSourceDomain::batch_reset()
15,962✔
82
{
83
// Reset scalar fluxes and iteration volume tallies to zero
84
#pragma omp parallel for
8,712✔
85
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
32,862,035✔
86
    source_regions_.volume(sr) = 0.0;
32,854,785✔
87
    source_regions_.volume_sq(sr) = 0.0;
32,854,785✔
88
  }
89

90
#pragma omp parallel for
8,712✔
91
  for (int64_t se = 0; se < n_source_elements(); se++) {
37,526,975✔
92
    source_regions_.scalar_flux_new(se) = 0.0;
37,519,725✔
93
  }
94
}
15,962✔
95

96
void FlatSourceDomain::accumulate_iteration_flux()
7,486✔
97
{
98
#pragma omp parallel for
4,086✔
99
  for (int64_t se = 0; se < n_source_elements(); se++) {
17,236,600✔
100
    source_regions_.scalar_flux_final(se) +=
17,233,200✔
101
      source_regions_.scalar_flux_new(se);
17,233,200✔
102
  }
103
}
7,486✔
104

105
void FlatSourceDomain::update_single_neutron_source(SourceRegionHandle& srh)
45,065,426✔
106
{
107
  // Reset all source regions to zero (important for void regions)
108
  for (int g = 0; g < negroups_; g++) {
96,355,630✔
109
    srh.source(g) = 0.0;
51,290,204✔
110
  }
111

112
  // Add scattering + fission source
113
  int material = srh.material();
45,065,426✔
114
  int temp = srh.temperature_idx();
45,065,426✔
115
  double density_mult = srh.density_mult();
45,065,426✔
116
  if (material != MATERIAL_VOID) {
45,065,426✔
117
    double inverse_k_eff = 1.0 / k_eff_;
44,623,778✔
118
    const int material_offset = (material * ntemperature_ + temp) * negroups_;
44,623,778✔
119
    const int scatter_offset =
44,623,778✔
120
      (material * ntemperature_ + temp) * negroups_ * negroups_;
121
    for (int g_out = 0; g_out < negroups_; g_out++) {
95,471,566✔
122
      double sigma_t = sigma_t_[material_offset + g_out] * density_mult;
50,847,788✔
123
      double scatter_source = 0.0;
50,847,788✔
124
      double fission_source = 0.0;
50,847,788✔
125

126
      for (int g_in = 0; g_in < negroups_; g_in++) {
145,262,266✔
127
        double scalar_flux = srh.scalar_flux_old(g_in);
94,414,478!
128
        double sigma_s =
94,414,478✔
129
          sigma_s_[scatter_offset + g_out * negroups_ + g_in] * density_mult;
94,414,478!
130
        double nu_sigma_f = nu_sigma_f_[material_offset + g_in] * density_mult;
94,414,478✔
131
        double chi = chi_[material_offset + g_out];
94,414,478✔
132

133
        scatter_source += sigma_s * scalar_flux;
94,414,478✔
134
        if (settings::create_fission_neutrons) {
94,414,478!
135
          fission_source += nu_sigma_f * scalar_flux * chi;
94,414,478✔
136
        }
137
      }
138
      srh.source(g_out) =
50,847,788✔
139
        (scatter_source + fission_source * inverse_k_eff) / sigma_t;
50,847,788✔
140
    }
141
  }
142

143
  // Add external source if in fixed source mode
144
  if (settings::run_mode == RunMode::FIXED_SOURCE) {
45,065,426✔
145
    for (int g = 0; g < negroups_; g++) {
91,772,678✔
146
      srh.source(g) += srh.external_source(g);
47,322,361✔
147
    }
148
  }
149
}
45,065,426✔
150

151
// Compute new estimate of scattering + fission sources in each source region
152
// based on the flux estimate from the previous iteration.
153
void FlatSourceDomain::update_all_neutron_sources()
15,962✔
154
{
155
  simulation::time_update_src.start();
15,962✔
156

157
#pragma omp parallel for
8,712✔
158
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
32,862,035✔
159
    SourceRegionHandle srh = source_regions_.get_source_region_handle(sr);
32,854,785✔
160
    update_single_neutron_source(srh);
32,854,785✔
161
  }
162

163
  simulation::time_update_src.stop();
15,962✔
164
}
15,962✔
165

166
// Normalizes flux and updates simulation-averaged volume estimate
167
void FlatSourceDomain::normalize_scalar_flux_and_volumes(
8,427✔
168
  double total_active_distance_per_iteration)
169
{
170
  double normalization_factor = 1.0 / total_active_distance_per_iteration;
8,427✔
171
  double volume_normalization_factor =
8,427✔
172
    1.0 / (total_active_distance_per_iteration * simulation::current_batch);
8,427✔
173

174
// Normalize scalar flux to total distance travelled by all rays this
175
// iteration
176
#pragma omp parallel for
4,602✔
177
  for (int64_t se = 0; se < n_source_elements(); se++) {
23,285,535✔
178
    source_regions_.scalar_flux_new(se) *= normalization_factor;
23,281,710✔
179
  }
180

181
// Accumulate cell-wise ray length tallies collected this iteration, then
182
// update the simulation-averaged cell-wise volume estimates
183
#pragma omp parallel for
4,602✔
184
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
20,456,595✔
185
    source_regions_.volume_t(sr) += source_regions_.volume(sr);
20,452,770✔
186
    source_regions_.volume_sq_t(sr) += source_regions_.volume_sq(sr);
20,452,770✔
187
    source_regions_.volume_naive(sr) =
20,452,770✔
188
      source_regions_.volume(sr) * normalization_factor;
20,452,770✔
189
    source_regions_.volume_sq(sr) =
20,452,770✔
190
      source_regions_.volume_sq_t(sr) / source_regions_.volume_t(sr);
20,452,770✔
191
    source_regions_.volume(sr) =
20,452,770✔
192
      source_regions_.volume_t(sr) * volume_normalization_factor;
20,452,770✔
193
  }
194
}
8,427✔
195

196
void FlatSourceDomain::set_flux_to_flux_plus_source(
45,626,069✔
197
  int64_t sr, double volume, int g)
198
{
199
  int material = source_regions_.material(sr);
45,626,069✔
200
  int temp = source_regions_.temperature_idx(sr);
45,626,069✔
201
  if (material == MATERIAL_VOID) {
45,626,069✔
202
    source_regions_.scalar_flux_new(sr, g) /= volume;
882,416!
203
    if (settings::run_mode == RunMode::FIXED_SOURCE) {
882,416!
204
      source_regions_.scalar_flux_new(sr, g) +=
882,416✔
205
        0.5f * source_regions_.external_source(sr, g) *
882,416✔
206
        source_regions_.volume_sq(sr);
882,416✔
207
    }
208
  } else {
209
    double sigma_t =
44,743,653✔
210
      sigma_t_[(material * ntemperature_ + temp) * negroups_ + g] *
44,743,653✔
211
      source_regions_.density_mult(sr);
44,743,653✔
212
    source_regions_.scalar_flux_new(sr, g) /= (sigma_t * volume);
44,743,653✔
213
    source_regions_.scalar_flux_new(sr, g) += source_regions_.source(sr, g);
44,743,653✔
214
  }
215
}
45,626,069✔
216

217
void FlatSourceDomain::set_flux_to_old_flux(int64_t sr, int g)
1,542✔
218
{
219
  source_regions_.scalar_flux_new(sr, g) =
1,542✔
220
    source_regions_.scalar_flux_old(sr, g);
1,542✔
221
}
1,542✔
222

223
void FlatSourceDomain::set_flux_to_source(int64_t sr, int g)
8,927,547✔
224
{
225
  source_regions_.scalar_flux_new(sr, g) = source_regions_.source(sr, g);
8,927,547✔
226
}
8,927,547✔
227

228
// Combine transport flux contributions and flat source contributions from the
229
// previous iteration to generate this iteration's estimate of scalar flux.
230
int64_t FlatSourceDomain::add_source_to_scalar_flux()
15,962✔
231
{
232
  int64_t n_hits = 0;
15,962✔
233
  double inverse_batch = 1.0 / simulation::current_batch;
15,962✔
234

235
#pragma omp parallel for reduction(+ : n_hits)
8,712✔
236
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
33,919,185✔
237

238
    double volume_simulation_avg = source_regions_.volume(sr);
33,911,935✔
239
    double volume_iteration = source_regions_.volume_naive(sr);
33,911,935✔
240

241
    // Increment the number of hits if cell was hit this iteration
242
    if (volume_iteration) {
33,911,935✔
243
      n_hits++;
29,852,780✔
244
    }
245

246
    // Set the SR to small status if its expected number of hits
247
    // per iteration is less than 1.5
248
    if (source_regions_.n_hits(sr) * inverse_batch < MIN_HITS_PER_BATCH) {
33,911,935✔
249
      source_regions_.is_small(sr) = 1;
6,812,705✔
250
    } else {
251
      source_regions_.is_small(sr) = 0;
27,099,230✔
252
    }
253

254
    // The volume treatment depends on the volume estimator type
255
    // and whether or not an external source is present in the cell.
256
    double volume;
33,911,935✔
257
    switch (volume_estimator_) {
33,911,935!
258
    case RandomRayVolumeEstimator::NAIVE:
259
      volume = volume_iteration;
260
      break;
261
    case RandomRayVolumeEstimator::SIMULATION_AVERAGED:
432,000✔
262
      volume = volume_simulation_avg;
432,000✔
263
      break;
432,000✔
264
    case RandomRayVolumeEstimator::HYBRID:
32,290,735✔
265
      if (source_regions_.external_source_present(sr) ||
32,290,735✔
266
          source_regions_.is_small(sr)) {
28,100,835✔
267
        volume = volume_iteration;
268
      } else {
269
        volume = volume_simulation_avg;
270
      }
271
      break;
272
    default:
273
      fatal_error("Invalid volume estimator type");
274
    }
275

276
    for (int g = 0; g < negroups_; g++) {
72,648,410✔
277
      // There are three scenarios we need to consider:
278
      if (volume_iteration > 0.0) {
38,736,475✔
279
        // 1. If the FSR was hit this iteration, then the new flux is equal to
280
        // the flat source from the previous iteration plus the contributions
281
        // from rays passing through the source region (computed during the
282
        // transport sweep)
283
        set_flux_to_flux_plus_source(sr, volume, g);
34,677,170✔
284
      } else if (volume_simulation_avg > 0.0) {
4,059,305✔
285
        // 2. If the FSR was not hit this iteration, but has been hit some
286
        // previous iteration, then we need to make a choice about what
287
        // to do. Naively we will usually want to set the flux to be equal
288
        // to the reduced source. However, in fixed source problems where
289
        // there is a strong external source present in the cell, and where
290
        // the cell has a very low cross section, this approximation will
291
        // cause a huge upward bias in the flux estimate of the cell (in these
292
        // conditions, the flux estimate can be orders of magnitude too large).
293
        // Thus, to avoid this bias, if any external source is present
294
        // in the cell we will use the previous iteration's flux estimate. This
295
        // injects a small degree of correlation into the simulation, but this
296
        // is going to be trivial when the miss rate is a few percent or less.
297
        if (source_regions_.external_source_present(sr)) {
4,059,295✔
298
          set_flux_to_old_flux(sr, g);
1,320✔
299
        } else {
300
          set_flux_to_source(sr, g);
4,057,975✔
301
        }
302
      }
303
      // Halt if NaN implosion is detected
304
      if (!std::isfinite(source_regions_.scalar_flux_new(sr, g))) {
38,736,475!
305
        fatal_error("A source region scalar flux is not finite. "
306
                    "This indicates a numerical instability in the "
307
                    "simulation. Consider increasing ray density or adjusting "
308
                    "the source region mesh.");
309
      }
310
    }
311
  }
312

313
  // Return the number of source regions that were hit this iteration
314
  return n_hits;
15,962✔
315
}
316

317
// Generates new estimate of k_eff based on the differences between this
318
// iteration's estimate of the scalar flux and the last iteration's estimate.
319
void FlatSourceDomain::compute_k_eff()
5,610✔
320
{
321
  double fission_rate_old = 0;
5,610✔
322
  double fission_rate_new = 0;
5,610✔
323

324
  // Vector for gathering fission source terms for Shannon entropy calculation
325
  vector<float> p(n_source_regions(), 0.0f);
5,610✔
326

327
#pragma omp parallel for reduction(+ : fission_rate_old, fission_rate_new)
3,060✔
328
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
388,740✔
329

330
    // If simulation averaged volume is zero, don't include this cell
331
    double volume = source_regions_.volume(sr);
386,190!
332
    if (volume == 0.0) {
386,190!
333
      continue;
334
    }
335

336
    int material = source_regions_.material(sr);
386,190!
337
    int temp = source_regions_.temperature_idx(sr);
386,190!
338
    if (material == MATERIAL_VOID) {
386,190!
339
      continue;
340
    }
341

342
    double sr_fission_source_old = 0;
343
    double sr_fission_source_new = 0;
344

345
    for (int g = 0; g < negroups_; g++) {
2,986,920✔
346
      double nu_sigma_f =
2,600,730✔
347
        nu_sigma_f_[(material * ntemperature_ + temp) * negroups_ + g] *
2,600,730✔
348
        source_regions_.density_mult(sr);
2,600,730✔
349
      sr_fission_source_old +=
2,600,730✔
350
        nu_sigma_f * source_regions_.scalar_flux_old(sr, g);
2,600,730✔
351
      sr_fission_source_new +=
2,600,730✔
352
        nu_sigma_f * source_regions_.scalar_flux_new(sr, g);
2,600,730✔
353
    }
354

355
    // Compute total fission rates in FSR
356
    sr_fission_source_old *= volume;
386,190✔
357
    sr_fission_source_new *= volume;
386,190✔
358

359
    // Accumulate totals
360
    fission_rate_old += sr_fission_source_old;
386,190✔
361
    fission_rate_new += sr_fission_source_new;
386,190✔
362

363
    // Store total fission rate in the FSR for Shannon calculation
364
    p[sr] = sr_fission_source_new;
386,190✔
365
  }
366

367
  double k_eff_new = k_eff_ * (fission_rate_new / fission_rate_old);
5,610✔
368

369
  double H = 0.0;
5,610✔
370
  // defining an inverse sum for better performance
371
  double inverse_sum = 1 / fission_rate_new;
5,610✔
372

373
#pragma omp parallel for reduction(+ : H)
3,060✔
374
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
388,740✔
375
    // Only if FSR has non-negative and non-zero fission source
376
    if (p[sr] > 0.0f) {
386,190✔
377
      // Normalize to total weight of bank sites. p_i for better performance
378
      float p_i = p[sr] * inverse_sum;
164,195✔
379
      // Sum values to obtain Shannon entropy.
380
      H -= p_i * std::log2(p_i);
164,195✔
381
    }
382
  }
383

384
  // Adds entropy value to shared entropy vector in openmc namespace.
385
  simulation::entropy.push_back(H);
5,610✔
386

387
  fission_rate_ = fission_rate_new;
5,610✔
388
  k_eff_ = k_eff_new;
5,610✔
389
}
5,610✔
390

391
// This function is responsible for generating a mapping between random
392
// ray flat source regions (cell instances) and tally bins. The mapping
393
// takes the form of a "TallyTask" object, which accounts for one single
394
// score being applied to a single tally. Thus, a single source region
395
// may have anywhere from zero to many tally tasks associated with it ---
396
// meaning that the global "tally_task" data structure is in 2D. The outer
397
// dimension corresponds to the source element (i.e., each entry corresponds
398
// to a specific energy group within a specific source region), and the
399
// inner dimension corresponds to the tallying task itself. Mechanically,
400
// the mapping between FSRs and spatial filters is done by considering
401
// the location of a single known ray midpoint that passed through the
402
// FSR. I.e., during transport, the first ray to pass through a given FSR
403
// will write down its midpoint for use with this function. This is a cheap
404
// and easy way of mapping FSRs to spatial tally filters, but comes with
405
// the downside of adding the restriction that spatial tally filters must
406
// share boundaries with the physical geometry of the simulation (so as
407
// not to subdivide any FSR). It is acceptable for a spatial tally region
408
// to contain multiple FSRs, but not the other way around.
409

410
// TODO: In future work, it would be preferable to offer a more general
411
// (but perhaps slightly more expensive) option for handling arbitrary
412
// spatial tallies that would be allowed to subdivide FSRs.
413

414
// Besides generating the mapping structure, this function also keeps track
415
// of whether or not all flat source regions have been hit yet. This is
416
// required, as there is no guarantee that all flat source regions will
417
// be hit every iteration, such that in the first few iterations some FSRs
418
// may not have a known position within them yet to facilitate mapping to
419
// spatial tally filters. However, after several iterations, if all FSRs
420
// have been hit and have had a tally map generated, then this status will
421
// be passed back to the caller to alert them that this function doesn't
422
// need to be called for the remainder of the simulation.
423

424
// It takes as an argument the starting index in the source region array,
425
// and it will operate from that index until the end of the array. This
426
// is useful as it can be called for both explicit user source regions or
427
// when a source region mesh is overlaid.
428

429
void FlatSourceDomain::convert_source_regions_to_tallies(int64_t start_sr_id)
914✔
430
{
431
  openmc::simulation::time_tallies.start();
914✔
432

433
  // Tracks if we've generated a mapping yet for all source regions.
434
  bool all_source_regions_mapped = true;
914✔
435

436
// Attempt to generate mapping for all source regions
437
#pragma omp parallel for
499✔
438
  for (int64_t sr = start_sr_id; sr < n_source_regions(); sr++) {
1,057,565✔
439

440
    // If this source region has not been hit by a ray yet, then
441
    // we aren't going to be able to map it, so skip it.
442
    if (!source_regions_.position_recorded(sr)) {
1,057,150!
443
      all_source_regions_mapped = false;
444
      continue;
445
    }
446

447
    // A particle located at the recorded midpoint of a ray
448
    // crossing through this source region is used to estabilish
449
    // the spatial location of the source region
450
    Particle p;
1,057,150✔
451
    p.r() = source_regions_.position(sr);
1,057,150✔
452
    p.r_last() = source_regions_.position(sr);
1,057,150✔
453
    p.u() = {1.0, 0.0, 0.0};
1,057,150✔
454
    bool found = exhaustive_find_cell(p);
1,057,150✔
455

456
    // Loop over energy groups (so as to support energy filters)
457
    for (int g = 0; g < negroups_; g++) {
2,273,900✔
458

459
      // Set particle to the current energy
460
      p.g() = g;
1,216,750✔
461
      p.g_last() = g;
1,216,750✔
462
      p.E() = data::mg.energy_bin_avg_[p.g()];
1,216,750!
463
      p.E_last() = p.E();
1,216,750✔
464

465
      int64_t source_element = sr * negroups_ + g;
1,216,750✔
466

467
      // If this task has already been populated, we don't need to do
468
      // it again.
469
      if (source_regions_.tally_task(sr, g).size() > 0) {
1,216,750!
470
        continue;
471
      }
472

473
      // Loop over all active tallies. This logic is essentially identical
474
      // to what happens when scanning for applicable tallies during
475
      // MC transport.
476
      for (int i_tally = 0; i_tally < model::tallies.size(); i_tally++) {
4,654,740✔
477
        Tally& tally {*model::tallies[i_tally]};
3,437,990✔
478

479
        // Initialize an iterator over valid filter bin combinations.
480
        // If there are no valid combinations, use a continue statement
481
        // to ensure we skip the assume_separate break below.
482
        auto filter_iter = FilterBinIter(tally, p);
3,437,990✔
483
        auto end = FilterBinIter(tally, true, &p.filter_matches());
3,437,990✔
484
        if (filter_iter == end)
3,437,990✔
485
          continue;
2,083,840✔
486

487
        // Loop over filter bins.
488
        for (; filter_iter != end; ++filter_iter) {
4,062,450✔
489
          auto filter_index = filter_iter.index_;
490
          auto filter_weight = filter_iter.weight_;
3,065,860✔
491

492
          // Loop over scores
493
          for (int score = 0; score < tally.scores_.size(); score++) {
3,065,860✔
494
            auto score_bin = tally.scores_[score];
1,711,710✔
495
            // If a valid tally, filter, and score combination has been found,
496
            // then add it to the list of tally tasks for this source element.
497
            TallyTask task(i_tally, filter_index, score, score_bin);
1,711,710✔
498
            source_regions_.tally_task(sr, g).push_back(task);
1,711,710✔
499

500
            // Also add this task to the list of volume tasks for this source
501
            // region.
502
            source_regions_.volume_task(sr).insert(task);
1,711,710✔
503
          }
504
        }
505
      }
506
      // Reset all the filter matches for the next tally event.
507
      for (auto& match : p.filter_matches())
5,151,360✔
508
        match.bins_present_ = false;
3,934,610✔
509
    }
510
  }
1,057,150✔
511
  openmc::simulation::time_tallies.stop();
914✔
512

513
  mapped_all_tallies_ = all_source_regions_mapped;
914✔
514
}
914✔
515

516
// Set the volume accumulators to zero for all tallies
517
void FlatSourceDomain::reset_tally_volumes()
7,486✔
518
{
519
  if (volume_normalized_flux_tallies_) {
7,486✔
520
#pragma omp parallel for
3,180✔
521
    for (int i = 0; i < tally_volumes_.size(); i++) {
8,650✔
522
      auto& tensor = tally_volumes_[i];
6,000✔
523
      tensor.fill(0.0); // Set all elements of the tensor to 0.0
12,000✔
524
    }
525
  }
526
}
7,486✔
527

528
// In fixed source mode, due to the way that volumetric fixed sources are
529
// converted and applied as volumetric sources in one or more source regions,
530
// we need to perform an additional normalization step to ensure that the
531
// reported scalar fluxes are in units per source neutron. This allows for
532
// direct comparison of reported tallies to Monte Carlo flux results.
533
// This factor needs to be computed at each iteration, as it is based on the
534
// volume estimate of each FSR, which improves over the course of the
535
// simulation
536
double FlatSourceDomain::compute_fixed_source_normalization_factor() const
8,369✔
537
{
538
  // Eigenvalue mode normalization
539
  if (settings::run_mode == RunMode::EIGENVALUE) {
8,369✔
540
    // Normalize fluxes by total number of fission neutrons produced. This
541
    // ensures consistent scaling of the eigenvector such that its magnitude is
542
    // comparable to the eigenvector produced by the Monte Carlo solver.
543
    // Multiplying by the eigenvalue is unintuitive, but it is necessary.
544
    // If the eigenvalue is 1.2, per starting source neutron, you will
545
    // generate 1.2 neutrons. Thus if we normalize to generating only ONE
546
    // neutron in total for the whole domain, then we don't actually have enough
547
    // flux to generate the required 1.2 neutrons. We only know the flux
548
    // required to generate 1 neutron (which would have required less than one
549
    // starting neutron). Thus, you have to scale the flux up by the eigenvalue
550
    // such that 1.2 neutrons are generated, so as to be consistent with the
551
    // bookkeeping in MC which is all done per starting source neutron (not per
552
    // neutron produced).
553
    return k_eff_ / (fission_rate_ * simulation_volume_);
3,671✔
554
  }
555

556
  // If we are in adjoint mode of a fixed source problem, the external
557
  // source is already normalized, such that all resulting fluxes are
558
  // also normalized.
559
  if (adjoint_) {
4,698✔
560
    return 1.0;
561
  }
562

563
  // Fixed source mode normalization
564

565
  // Step 1 is to sum over all source regions and energy groups to get the
566
  // total external source strength in the simulation.
567
  double simulation_external_source_strength = 0.0;
2,296✔
568
#pragma omp parallel for reduction(+ : simulation_external_source_strength)
2,296✔
569
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
14,761,638✔
570
    int material = source_regions_.material(sr);
14,759,770✔
571
    int temp = source_regions_.temperature_idx(sr);
14,759,770✔
572
    double volume = source_regions_.volume(sr) * simulation_volume_;
14,759,770✔
573
    for (int g = 0; g < negroups_; g++) {
30,062,420✔
574
      // For non-void regions, we store the external source pre-divided by
575
      // sigma_t. We need to multiply non-void regions back up by sigma_t
576
      // to get the total source strength in the expected units.
577
      double sigma_t = 1.0;
15,302,650✔
578
      if (material != MATERIAL_VOID) {
15,302,650✔
579
        sigma_t = sigma_t_[(material * ntemperature_ + temp) * negroups_ + g] *
15,092,410✔
580
                  source_regions_.density_mult(sr);
15,092,410✔
581
      }
582
      simulation_external_source_strength +=
15,302,650✔
583
        source_regions_.external_source(sr, g) * sigma_t * volume;
15,302,650✔
584
    }
585
  }
586

587
  // Step 2 is to determine the total user-specified external source strength
588
  double user_external_source_strength = 0.0;
4,164✔
589
  for (auto& ext_source : model::external_sources) {
8,328✔
590
    user_external_source_strength += ext_source->strength();
4,164✔
591
  }
592

593
  // The correction factor is the ratio of the user-specified external source
594
  // strength to the simulation external source strength.
595
  double source_normalization_factor =
4,164✔
596
    user_external_source_strength / simulation_external_source_strength;
597

598
  return source_normalization_factor;
4,164✔
599
}
600

601
// Tallying in random ray is not done directly during transport, rather,
602
// it is done only once after each power iteration. This is made possible
603
// by way of a mapping data structure that relates spatial source regions
604
// (FSRs) to tally/filter/score combinations. The mechanism by which the
605
// mapping is done (and the limitations incurred) is documented in the
606
// "convert_source_regions_to_tallies()" function comments above. The present
607
// tally function simply traverses the mapping data structure and executes
608
// the scoring operations to OpenMC's native tally result arrays.
609

610
void FlatSourceDomain::random_ray_tally()
7,486✔
611
{
612
  openmc::simulation::time_tallies.start();
7,486✔
613

614
  // Reset our tally volumes to zero
615
  reset_tally_volumes();
7,486✔
616

617
  double source_normalization_factor =
7,486✔
618
    compute_fixed_source_normalization_factor();
7,486✔
619

620
// We loop over all source regions and energy groups. For each
621
// element, we check if there are any scores needed and apply
622
// them.
623
#pragma omp parallel for
4,086✔
624
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
15,556,000✔
625
    // The fsr.volume_ is the unitless fractional simulation averaged volume
626
    // (i.e., it is the FSR's fraction of the overall simulation volume). The
627
    // simulation_volume_ is the total 3D physical volume in cm^3 of the
628
    // entire global simulation domain (as defined by the ray source box).
629
    // Thus, the FSR's true 3D spatial volume in cm^3 is found by multiplying
630
    // its fraction of the total volume by the total volume. Not important in
631
    // eigenvalue solves, but useful in fixed source solves for returning the
632
    // flux shape with a magnitude that makes sense relative to the fixed
633
    // source strength.
634
    double volume = source_regions_.volume(sr) * simulation_volume_;
15,552,600✔
635

636
    double material = source_regions_.material(sr);
15,552,600✔
637
    int temp = source_regions_.temperature_idx(sr);
15,552,600✔
638
    double density_mult = source_regions_.density_mult(sr);
15,552,600✔
639
    for (int g = 0; g < negroups_; g++) {
32,785,800✔
640
      double flux =
17,233,200✔
641
        source_regions_.scalar_flux_new(sr, g) * source_normalization_factor;
17,233,200✔
642

643
      // Determine numerical score value
644
      for (auto& task : source_regions_.tally_task(sr, g)) {
41,117,800✔
645
        double score = 0.0;
23,884,600✔
646
        switch (task.score_type) {
23,884,600!
647

648
        case SCORE_FLUX:
20,036,800✔
649
          score = flux * volume;
20,036,800✔
650
          break;
20,036,800✔
651

652
        case SCORE_TOTAL:
653
          if (material != MATERIAL_VOID) {
×
654
            score =
655
              flux * volume *
656
              sigma_t_[(material * ntemperature_ + temp) * negroups_ + g] *
657
              density_mult;
658
          }
659
          break;
660

661
        case SCORE_FISSION:
1,923,600✔
662
          if (material != MATERIAL_VOID) {
1,923,600!
663
            score =
1,923,600✔
664
              flux * volume *
1,923,600✔
665
              sigma_f_[(material * ntemperature_ + temp) * negroups_ + g] *
1,923,600✔
666
              density_mult;
667
          }
668
          break;
669

670
        case SCORE_NU_FISSION:
1,923,600✔
671
          if (material != MATERIAL_VOID) {
1,923,600!
672
            score =
1,923,600✔
673
              flux * volume *
1,923,600✔
674
              nu_sigma_f_[(material * ntemperature_ + temp) * negroups_ + g] *
1,923,600✔
675
              density_mult;
676
          }
677
          break;
678

679
        case SCORE_EVENTS:
680
          score = 1.0;
681
          break;
682

683
        case SCORE_KAPPA_FISSION:
600✔
684
          score =
600✔
685
            flux * volume *
600✔
686
            kappa_fission_[(material * ntemperature_ + temp) * negroups_ + g] *
600✔
687
            density_mult;
688
          break;
600✔
689

690
        default:
691
          fatal_error("Invalid score specified in tallies.xml. Only flux, "
692
                      "total, fission, nu-fission, kappa-fission, and events "
693
                      "are supported in random ray mode.");
694
          break;
23,884,600✔
695
        }
696
        // Apply score to the appropriate tally bin
697
        Tally& tally {*model::tallies[task.tally_idx]};
23,884,600✔
698
#pragma omp atomic
699
        tally.results_(task.filter_idx, task.score_idx, TallyResult::VALUE) +=
23,884,600✔
700
          score;
701
      }
702
    }
703

704
    // For flux tallies, the total volume of the spatial region is needed
705
    // for normalizing the flux. We store this volume in a separate tensor.
706
    // We only contribute to each volume tally bin once per FSR.
707
    if (volume_normalized_flux_tallies_) {
15,552,600✔
708
      for (const auto& task : source_regions_.volume_task(sr)) {
36,769,200✔
709
        if (task.score_type == SCORE_FLUX) {
21,364,800✔
710
#pragma omp atomic
711
          tally_volumes_[task.tally_idx](task.filter_idx, task.score_idx) +=
19,167,600✔
712
            volume;
713
        }
714
      }
715
    }
716
  } // end FSR loop
717

718
  // Normalize any flux scores by the total volume of the FSRs scoring to that
719
  // bin. To do this, we loop over all tallies, and then all filter bins,
720
  // and then scores. For each score, we check the tally data structure to
721
  // see what index that score corresponds to. If that score is a flux score,
722
  // then we divide it by volume.
723
  if (volume_normalized_flux_tallies_) {
7,486✔
724
    for (int i = 0; i < model::tallies.size(); i++) {
19,030✔
725
      Tally& tally {*model::tallies[i]};
13,200✔
726
#pragma omp parallel for
7,200✔
727
      for (int bin = 0; bin < tally.n_filter_bins(); bin++) {
724,775✔
728
        for (int score_idx = 0; score_idx < tally.n_scores(); score_idx++) {
1,460,050✔
729
          auto score_type = tally.scores_[score_idx];
741,275✔
730
          if (score_type == SCORE_FLUX) {
741,275✔
731
            double vol = tally_volumes_[i](bin, score_idx);
718,775!
732
            if (vol > 0.0) {
718,775!
733
              tally.results_(bin, score_idx, TallyResult::VALUE) /= vol;
718,775✔
734
            }
735
          }
736
        }
737
      }
738
    }
739
  }
740

741
  openmc::simulation::time_tallies.stop();
7,486✔
742
}
7,486✔
743

UNCOV
744
double FlatSourceDomain::evaluate_flux_at_point(
×
745
  Position r, int64_t sr, int g) const
746
{
UNCOV
747
  return source_regions_.scalar_flux_final(sr, g) /
×
UNCOV
748
         (settings::n_batches - settings::n_inactive);
×
749
}
750

751
// Outputs all basic material, FSR ID, multigroup flux, and
752
// fission source data to .vtk file that can be directly
753
// loaded and displayed by Paraview. Note that .vtk binary
754
// files require big endian byte ordering, so endianness
755
// is checked and flipped if necessary.
UNCOV
756
void FlatSourceDomain::output_to_vtk() const
×
757
{
758
  // Rename .h5 plot filename(s) to .vtk filenames
759
  for (int p = 0; p < model::plots.size(); p++) {
×
760
    PlottableInterface* plot = model::plots[p].get();
×
UNCOV
761
    plot->path_plot() =
×
UNCOV
762
      plot->path_plot().substr(0, plot->path_plot().find_last_of('.')) + ".vtk";
×
763
  }
764

765
  // Print header information
UNCOV
766
  print_plot();
×
767

768
  // Outer loop over plots
UNCOV
769
  for (int plt = 0; plt < model::plots.size(); plt++) {
×
770

771
    // Get handle to OpenMC plot object and extract params
UNCOV
772
    Plot* openmc_plot = dynamic_cast<Plot*>(model::plots[plt].get());
×
773

774
    // Random ray plots only support voxel plots
UNCOV
775
    if (!openmc_plot) {
×
UNCOV
776
      warning(fmt::format("Plot {} is invalid plot type -- only voxel plotting "
×
777
                          "is allowed in random ray mode.",
778
        plt));
779
      continue;
×
UNCOV
780
    } else if (openmc_plot->type_ != Plot::PlotType::voxel) {
×
UNCOV
781
      warning(fmt::format("Plot {} is invalid plot type -- only voxel plotting "
×
782
                          "is allowed in random ray mode.",
783
        plt));
UNCOV
784
      continue;
×
785
    }
786

787
    int Nx = openmc_plot->pixels_[0];
×
788
    int Ny = openmc_plot->pixels_[1];
×
789
    int Nz = openmc_plot->pixels_[2];
×
790
    Position origin = openmc_plot->origin_;
×
791
    Position width = openmc_plot->width_;
×
792
    Position ll = origin - width / 2.0;
×
793
    double x_delta = width.x / Nx;
×
794
    double y_delta = width.y / Ny;
×
UNCOV
795
    double z_delta = width.z / Nz;
×
UNCOV
796
    std::string filename = openmc_plot->path_plot();
×
797

798
    // Perform sanity checks on file size
799
    uint64_t bytes = Nx * Ny * Nz * (negroups_ + 1 + 1 + 1) * sizeof(float);
×
800
    write_message(5, "Processing plot {}: {}... (Estimated size is {} MB)",
×
801
      openmc_plot->id(), filename, bytes / 1.0e6);
×
UNCOV
802
    if (bytes / 1.0e9 > 1.0) {
×
803
      warning("Voxel plot specification is very large (>1 GB). Plotting may be "
×
804
              "slow.");
UNCOV
805
    } else if (bytes / 1.0e9 > 100.0) {
×
UNCOV
806
      fatal_error("Voxel plot specification is too large (>100 GB). Exiting.");
×
807
    }
808

809
    // Relate voxel spatial locations to random ray source regions
810
    vector<int> voxel_indices(Nx * Ny * Nz);
×
811
    vector<Position> voxel_positions(Nx * Ny * Nz);
×
UNCOV
812
    vector<double> weight_windows(Nx * Ny * Nz);
×
UNCOV
813
    float min_weight = 1e20;
×
814
#pragma omp parallel for collapse(3) reduction(min : min_weight)
815
    for (int z = 0; z < Nz; z++) {
×
816
      for (int y = 0; y < Ny; y++) {
×
817
        for (int x = 0; x < Nx; x++) {
×
818
          Position sample;
819
          sample.z = ll.z + z_delta / 2.0 + z * z_delta;
820
          sample.y = ll.y + y_delta / 2.0 + y * y_delta;
821
          sample.x = ll.x + x_delta / 2.0 + x * x_delta;
822
          Particle p;
×
823
          p.r() = sample;
×
824
          p.r_last() = sample;
825
          p.E() = 1.0;
826
          p.E_last() = 1.0;
827
          p.u() = {1.0, 0.0, 0.0};
×
828

829
          bool found = exhaustive_find_cell(p);
×
830
          if (!found) {
×
831
            voxel_indices[z * Ny * Nx + y * Nx + x] = -1;
832
            voxel_positions[z * Ny * Nx + y * Nx + x] = sample;
833
            weight_windows[z * Ny * Nx + y * Nx + x] = 0.0;
834
            continue;
835
          }
836

837
          SourceRegionKey sr_key = lookup_source_region_key(p);
×
838
          int64_t sr = -1;
839
          auto it = source_region_map_.find(sr_key);
×
840
          if (it != source_region_map_.end()) {
×
841
            sr = it->second;
842
          }
843

844
          voxel_indices[z * Ny * Nx + y * Nx + x] = sr;
×
845
          voxel_positions[z * Ny * Nx + y * Nx + x] = sample;
846

847
          if (variance_reduction::weight_windows.size() == 1) {
×
848
            auto [ww_found, ww] =
849
              variance_reduction::weight_windows[0]->get_weight_window(p);
×
850
            float weight = ww.lower_weight;
851
            weight_windows[z * Ny * Nx + y * Nx + x] = weight;
×
852
            if (weight < min_weight)
×
853
              min_weight = weight;
854
          }
855
        }
856
      }
857
    }
858

UNCOV
859
    double source_normalization_factor =
×
UNCOV
860
      compute_fixed_source_normalization_factor();
×
861

862
    // Open file for writing
UNCOV
863
    std::FILE* plot = std::fopen(filename.c_str(), "wb");
×
864

865
    // Write vtk metadata
866
    std::fprintf(plot, "# vtk DataFile Version 2.0\n");
×
867
    std::fprintf(plot, "Dataset File\n");
×
868
    std::fprintf(plot, "BINARY\n");
×
869
    std::fprintf(plot, "DATASET STRUCTURED_POINTS\n");
×
870
    std::fprintf(plot, "DIMENSIONS %d %d %d\n", Nx, Ny, Nz);
×
871
    std::fprintf(plot, "ORIGIN %lf %lf %lf\n", ll.x, ll.y, ll.z);
×
UNCOV
872
    std::fprintf(plot, "SPACING %lf %lf %lf\n", x_delta, y_delta, z_delta);
×
UNCOV
873
    std::fprintf(plot, "POINT_DATA %d\n", Nx * Ny * Nz);
×
874

875
    int64_t num_neg = 0;
876
    int64_t num_samples = 0;
877
    float min_flux = 0.0;
878
    float max_flux = -1.0e20;
879
    // Plot multigroup flux data
880
    for (int g = 0; g < negroups_; g++) {
×
881
      std::fprintf(plot, "SCALARS flux_group_%d float\n", g);
×
882
      std::fprintf(plot, "LOOKUP_TABLE default\n");
×
883
      for (int i = 0; i < Nx * Ny * Nz; i++) {
×
884
        int64_t fsr = voxel_indices[i];
×
885
        int64_t source_element = fsr * negroups_ + g;
×
886
        float flux = 0;
×
887
        if (fsr >= 0) {
×
888
          flux = evaluate_flux_at_point(voxel_positions[i], fsr, g);
×
889
          if (flux < 0.0)
×
UNCOV
890
            flux = FlatSourceDomain::evaluate_flux_at_point(
×
891
              voxel_positions[i], fsr, g);
×
892
        }
893
        if (flux < 0.0) {
×
894
          num_neg++;
×
UNCOV
895
          if (flux < min_flux) {
×
UNCOV
896
            min_flux = flux;
×
897
          }
898
        }
899
        if (flux > max_flux)
×
900
          max_flux = flux;
×
901
        num_samples++;
×
UNCOV
902
        flux = convert_to_big_endian<float>(flux);
×
UNCOV
903
        std::fwrite(&flux, sizeof(float), 1, plot);
×
904
      }
905
    }
906

907
    // Slightly negative fluxes can be normal when sampling corners of linear
908
    // source regions. However, very common and high magnitude negative fluxes
909
    // may indicate numerical instability.
UNCOV
910
    if (num_neg > 0) {
×
911
      warning(fmt::format("{} plot samples ({:.4f}%) contained negative fluxes "
×
912
                          "(minumum found = {:.2e} maximum_found = {:.2e})",
UNCOV
913
        num_neg, (100.0 * num_neg) / num_samples, min_flux, max_flux));
×
914
    }
915

916
    // Plot FSRs
917
    std::fprintf(plot, "SCALARS FSRs float\n");
×
918
    std::fprintf(plot, "LOOKUP_TABLE default\n");
×
919
    for (int fsr : voxel_indices) {
×
920
      float value = future_prn(10, fsr);
×
UNCOV
921
      value = convert_to_big_endian<float>(value);
×
UNCOV
922
      std::fwrite(&value, sizeof(float), 1, plot);
×
923
    }
924

925
    // Plot Materials
926
    std::fprintf(plot, "SCALARS Materials int\n");
×
927
    std::fprintf(plot, "LOOKUP_TABLE default\n");
×
928
    for (int fsr : voxel_indices) {
×
929
      int mat = -1;
×
930
      if (fsr >= 0)
×
931
        mat = source_regions_.material(fsr);
×
UNCOV
932
      mat = convert_to_big_endian<int>(mat);
×
UNCOV
933
      std::fwrite(&mat, sizeof(int), 1, plot);
×
934
    }
935

936
    // Plot fission source
937
    if (settings::run_mode == RunMode::EIGENVALUE) {
×
938
      std::fprintf(plot, "SCALARS total_fission_source float\n");
×
939
      std::fprintf(plot, "LOOKUP_TABLE default\n");
×
940
      for (int i = 0; i < Nx * Ny * Nz; i++) {
×
941
        int64_t fsr = voxel_indices[i];
×
942
        float total_fission = 0.0;
×
943
        if (fsr >= 0) {
×
944
          int mat = source_regions_.material(fsr);
×
945
          int temp = source_regions_.temperature_idx(fsr);
×
946
          if (mat != MATERIAL_VOID) {
×
947
            for (int g = 0; g < negroups_; g++) {
×
948
              int64_t source_element = fsr * negroups_ + g;
×
949
              float flux = evaluate_flux_at_point(voxel_positions[i], fsr, g);
×
950
              double sigma_f =
×
951
                sigma_f_[(mat * ntemperature_ + temp) * negroups_ + g] *
×
UNCOV
952
                source_regions_.density_mult(fsr);
×
UNCOV
953
              total_fission += sigma_f * flux;
×
954
            }
955
          }
956
        }
UNCOV
957
        total_fission = convert_to_big_endian<float>(total_fission);
×
UNCOV
958
        std::fwrite(&total_fission, sizeof(float), 1, plot);
×
959
      }
960
    } else {
961
      std::fprintf(plot, "SCALARS external_source float\n");
×
962
      std::fprintf(plot, "LOOKUP_TABLE default\n");
×
963
      for (int i = 0; i < Nx * Ny * Nz; i++) {
×
964
        int64_t fsr = voxel_indices[i];
×
965
        int mat = source_regions_.material(fsr);
×
966
        int temp = source_regions_.temperature_idx(fsr);
×
967
        float total_external = 0.0f;
×
UNCOV
968
        if (fsr >= 0) {
×
UNCOV
969
          for (int g = 0; g < negroups_; g++) {
×
970
            // External sources are already divided by sigma_t, so we need to
971
            // multiply it back to get the true external source.
972
            double sigma_t = 1.0;
×
973
            if (mat != MATERIAL_VOID) {
×
UNCOV
974
              sigma_t = sigma_t_[(mat * ntemperature_ + temp) * negroups_ + g] *
×
975
                        source_regions_.density_mult(fsr);
×
976
            }
UNCOV
977
            total_external += source_regions_.external_source(fsr, g) * sigma_t;
×
978
          }
979
        }
UNCOV
980
        total_external = convert_to_big_endian<float>(total_external);
×
UNCOV
981
        std::fwrite(&total_external, sizeof(float), 1, plot);
×
982
      }
983
    }
984

985
    // Plot weight window data
986
    if (variance_reduction::weight_windows.size() == 1) {
×
987
      std::fprintf(plot, "SCALARS weight_window_lower float\n");
×
988
      std::fprintf(plot, "LOOKUP_TABLE default\n");
×
989
      for (int i = 0; i < Nx * Ny * Nz; i++) {
×
990
        float weight = weight_windows[i];
×
991
        if (weight == 0.0)
×
992
          weight = min_weight;
×
UNCOV
993
        weight = convert_to_big_endian<float>(weight);
×
UNCOV
994
        std::fwrite(&weight, sizeof(float), 1, plot);
×
995
      }
996
    }
997

998
    std::fclose(plot);
×
UNCOV
999
  }
×
UNCOV
1000
}
×
1001

1002
void FlatSourceDomain::apply_external_source_to_source_region(
4,962✔
1003
  int src_idx, SourceRegionHandle& srh)
1004
{
1005
  auto s = (adjoint_ && !model::adjoint_sources.empty())
176!
1006
             ? model::adjoint_sources[src_idx].get()
5,138✔
1007
             : model::external_sources[src_idx].get();
4,786✔
1008
  auto is = dynamic_cast<IndependentSource*>(s);
4,962!
1009
  auto discrete = dynamic_cast<Discrete*>(is->energy());
4,962!
1010
  double strength_factor = is->strength();
4,962✔
1011
  const auto& discrete_energies = discrete->x();
4,962✔
1012
  const auto& discrete_probs = discrete->prob();
4,962✔
1013

1014
  srh.external_source_present() = 1;
4,962✔
1015

1016
  for (int i = 0; i < discrete_energies.size(); i++) {
10,056✔
1017
    int g = data::mg.get_group_index(discrete_energies[i]);
5,094✔
1018
    srh.external_source(g) += discrete_probs[i] * strength_factor;
5,094✔
1019
  }
1020
}
4,962✔
1021

1022
void FlatSourceDomain::apply_external_source_to_cell_instances(int32_t i_cell,
465✔
1023
  int src_idx, int target_material_id, const vector<int32_t>& instances)
1024
{
1025
  Cell& cell = *model::cells[i_cell];
465!
1026

1027
  if (cell.type_ != Fill::MATERIAL)
465!
1028
    return;
1029

1030
  for (int j : instances) {
29,220✔
1031
    int cell_material_idx = cell.material(j);
28,755!
1032
    int cell_material_id;
28,755✔
1033
    if (cell_material_idx == MATERIAL_VOID) {
28,755✔
1034
      cell_material_id = MATERIAL_VOID;
1035
    } else {
1036
      cell_material_id = model::materials[cell_material_idx]->id();
28,515✔
1037
    }
1038
    if (target_material_id == C_NONE ||
28,755✔
1039
        cell_material_id == target_material_id) {
28,755✔
1040
      int64_t source_region = source_region_offsets_[i_cell] + j;
2,955✔
1041
      external_volumetric_source_map_[source_region].push_back(src_idx);
2,955✔
1042
    }
1043
  }
1044
}
1045

1046
void FlatSourceDomain::apply_external_source_to_cell_and_children(
495✔
1047
  int32_t i_cell, int src_idx, int32_t target_material_id)
1048
{
1049
  Cell& cell = *model::cells[i_cell];
495✔
1050

1051
  if (cell.type_ == Fill::MATERIAL) {
495✔
1052
    vector<int> instances(cell.n_instances());
465✔
1053
    std::iota(instances.begin(), instances.end(), 0);
465✔
1054
    apply_external_source_to_cell_instances(
465✔
1055
      i_cell, src_idx, target_material_id, instances);
1056
  } else if (target_material_id == C_NONE) {
495!
1057
    std::unordered_map<int32_t, vector<int32_t>> cell_instance_list =
×
1058
      cell.get_contained_cells(0, nullptr);
×
UNCOV
1059
    for (const auto& pair : cell_instance_list) {
×
1060
      int32_t i_child_cell = pair.first;
×
UNCOV
1061
      apply_external_source_to_cell_instances(
×
UNCOV
1062
        i_child_cell, src_idx, target_material_id, pair.second);
×
1063
    }
UNCOV
1064
  }
×
1065
}
495✔
1066

1067
void FlatSourceDomain::count_external_source_regions()
1,319✔
1068
{
1069
  n_external_source_regions_ = 0;
1,319✔
1070
#pragma omp parallel for reduction(+ : n_external_source_regions_)
792✔
1071
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
1,228,137✔
1072
    if (source_regions_.external_source_present(sr)) {
1,227,610✔
1073
      n_external_source_regions_++;
157,500✔
1074
    }
1075
  }
1076
}
1,319✔
1077

1078
void FlatSourceDomain::convert_external_sources(bool use_adjoint_sources)
436✔
1079
{
1080
  // Determine whether forward or (local) adjoint sources are desired
1081
  const auto& sources =
436✔
1082
    use_adjoint_sources ? model::adjoint_sources : model::external_sources;
1083

1084
  // Loop over external sources
1085
  for (int es = 0; es < sources.size(); es++) {
872✔
1086

1087
    // Extract source information
1088
    Source* s = sources[es].get();
436!
1089
    IndependentSource* is = dynamic_cast<IndependentSource*>(s);
436!
1090
    Discrete* energy = dynamic_cast<Discrete*>(is->energy());
436✔
1091
    const std::unordered_set<int32_t>& domain_ids = is->domain_ids();
436✔
1092
    double strength_factor = is->strength();
436✔
1093

1094
    // If there is no domain constraint specified, then this must be a point
1095
    // source. In this case, we need to find the source region that contains the
1096
    // point source and apply or relate it to the external source.
1097
    if (is->domain_ids().size() == 0) {
436✔
1098

1099
      // Extract the point source coordinate and find the base source region at
1100
      // that point
1101
      auto sp = dynamic_cast<SpatialPoint*>(is->space());
16!
1102
      GeometryState gs;
16✔
1103
      gs.r() = sp->r();
16✔
1104
      gs.r_last() = sp->r();
16✔
1105
      gs.u() = {1.0, 0.0, 0.0};
16✔
1106
      bool found = exhaustive_find_cell(gs);
16✔
1107
      if (!found) {
16!
UNCOV
1108
        fatal_error(fmt::format("Could not find cell containing external "
×
1109
                                "point source at {}",
UNCOV
1110
          sp->r()));
×
1111
      }
1112
      SourceRegionKey key = lookup_source_region_key(gs);
16✔
1113

1114
      // With the source region and mesh bin known, we can use the
1115
      // accompanying SourceRegionKey as a key into a map that stores the
1116
      // corresponding external source index for the point source. Notably, we
1117
      // do not actually apply the external source to any source regions here,
1118
      // as if mesh subdivision is enabled, they haven't actually been
1119
      // discovered & initilized yet. When discovered, they will read from the
1120
      // external_source_map to determine if there are any external source
1121
      // terms that should be applied.
1122
      external_point_source_map_[key].push_back(es);
16✔
1123

1124
    } else {
16✔
1125
      // If not a point source, then use the volumetric domain constraints to
1126
      // determine which source regions to apply the external source to.
1127
      if (is->domain_type() == Source::DomainType::MATERIAL) {
420✔
1128
        for (int32_t material_id : domain_ids) {
60✔
1129
          for (int i_cell = 0; i_cell < model::cells.size(); i_cell++) {
120✔
1130
            apply_external_source_to_cell_and_children(i_cell, es, material_id);
90✔
1131
          }
1132
        }
1133
      } else if (is->domain_type() == Source::DomainType::CELL) {
390✔
1134
        for (int32_t cell_id : domain_ids) {
75✔
1135
          int32_t i_cell = model::cell_map[cell_id];
45✔
1136
          apply_external_source_to_cell_and_children(i_cell, es, C_NONE);
45✔
1137
        }
1138
      } else if (is->domain_type() == Source::DomainType::UNIVERSE) {
360!
1139
        for (int32_t universe_id : domain_ids) {
720✔
1140
          int32_t i_universe = model::universe_map[universe_id];
360✔
1141
          Universe& universe = *model::universes[i_universe];
360✔
1142
          for (int32_t i_cell : universe.cells_) {
720✔
1143
            apply_external_source_to_cell_and_children(i_cell, es, C_NONE);
360✔
1144
          }
1145
        }
1146
      }
1147
    }
1148
  } // End loop over external sources
1149
}
436✔
1150

1151
void FlatSourceDomain::flux_swap()
15,962✔
1152
{
1153
  source_regions_.flux_swap();
15,962✔
1154
}
15,962✔
1155

1156
void FlatSourceDomain::flatten_xs()
792✔
1157
{
1158
  // Temperature and angle indices, if using multiple temperature
1159
  // data sets and/or anisotropic data sets.
1160
  // TODO: Currently assumes we are only using single angle data.
1161
  const int a = 0;
792✔
1162

1163
  n_materials_ = data::mg.macro_xs_.size();
792✔
1164
  ntemperature_ = 1;
792✔
1165
  for (int i = 0; i < n_materials_; i++) {
2,945✔
1166
    ntemperature_ =
4,306✔
1167
      std::max(ntemperature_, data::mg.macro_xs_[i].n_temperature_points());
2,168✔
1168
  }
1169

1170
  for (int i = 0; i < n_materials_; i++) {
2,945✔
1171
    auto& m = data::mg.macro_xs_[i];
2,153✔
1172
    for (int t = 0; t < ntemperature_; t++) {
4,336✔
1173
      for (int g_out = 0; g_out < negroups_; g_out++) {
11,249✔
1174
        if (m.exists_in_model && t < m.n_temperature_points()) {
9,066✔
1175
          double sigma_t =
8,901✔
1176
            m.get_xs(MgxsType::TOTAL, g_out, NULL, NULL, NULL, t, a);
8,901✔
1177
          sigma_t_.push_back(sigma_t);
8,901✔
1178

1179
          if (sigma_t < MINIMUM_MACRO_XS) {
8,901✔
1180
            Material* mat = model::materials[i].get();
15✔
1181
            warning(fmt::format(
21✔
1182
              "Material \"{}\" (id: {}) has a group {} total cross section "
1183
              "({:.3e}) below the minimum threshold "
1184
              "({:.3e}). Material will be treated as pure void.",
1185
              mat->name(), mat->id(), g_out, sigma_t, MINIMUM_MACRO_XS));
27✔
1186
          }
1187

1188
          double nu_sigma_f =
8,901✔
1189
            m.get_xs(MgxsType::NU_FISSION, g_out, NULL, NULL, NULL, t, a);
8,901✔
1190
          nu_sigma_f_.push_back(nu_sigma_f);
8,901✔
1191

1192
          double sigma_f =
8,901✔
1193
            m.get_xs(MgxsType::FISSION, g_out, NULL, NULL, NULL, t, a);
8,901✔
1194
          sigma_f_.push_back(sigma_f);
8,901✔
1195

1196
          double chi =
8,901✔
1197
            m.get_xs(MgxsType::CHI_PROMPT, g_out, &g_out, NULL, NULL, t, a);
8,901✔
1198
          if (!std::isfinite(chi)) {
8,901✔
1199
            // MGXS interface may return NaN in some cases, such as when
1200
            // material is fissionable but has very small sigma_f.
1201
            chi = 0.0;
660✔
1202
          }
1203
          chi_.push_back(chi);
8,901✔
1204

1205
          double kappa_fission =
8,901✔
1206
            m.get_xs(MgxsType::KAPPA_FISSION, g_out, NULL, NULL, NULL, t, a);
8,901✔
1207
          kappa_fission_.push_back(kappa_fission);
8,901✔
1208

1209
          for (int g_in = 0; g_in < negroups_; g_in++) {
260,198✔
1210
            double sigma_s =
251,297✔
1211
              m.get_xs(MgxsType::NU_SCATTER, g_in, &g_out, NULL, NULL, t, a);
251,297✔
1212
            sigma_s_.push_back(sigma_s);
251,297✔
1213
            // For transport corrected XS data, diagonal elements may be
1214
            // negative. In this case, set a flag to enable transport
1215
            // stabilization for the simulation.
1216
            if (g_out == g_in && sigma_s < 0.0)
251,297✔
1217
              is_transport_stabilization_needed_ = true;
825✔
1218
          }
1219
        } else {
1220
          sigma_t_.push_back(0);
165✔
1221
          nu_sigma_f_.push_back(0);
165✔
1222
          sigma_f_.push_back(0);
165✔
1223
          chi_.push_back(0);
165✔
1224
          kappa_fission_.push_back(0);
165✔
1225
          for (int g_in = 0; g_in < negroups_; g_in++) {
960✔
1226
            sigma_s_.push_back(0);
795✔
1227
          }
1228
        }
1229
      }
1230
    }
1231
  }
1232
}
792✔
1233

1234
void FlatSourceDomain::set_fw_adjoint_sources()
76✔
1235
{
1236
  // Set the adjoint external source to 1/forward_flux. If the forward flux is
1237
  // negative, zero, or extremely close to zero, set the adjoint source to zero,
1238
  // as this is likely a very small source region that we don't need to bother
1239
  // trying to vector particles towards. In the case of flux "being extremely
1240
  // close to zero", we define this as being a fixed fraction of the maximum
1241
  // forward flux, below which we assume the flux would be physically
1242
  // undetectable.
1243

1244
  // First, find the maximum forward flux value
1245
  double max_flux = 0.0;
76✔
1246
#pragma omp parallel for reduction(max : max_flux)
46✔
1247
  for (int64_t se = 0; se < n_source_elements(); se++) {
169,270✔
1248
    double flux = source_regions_.scalar_flux_final(se);
169,240✔
1249
    if (flux > max_flux) {
169,240✔
1250
      max_flux = flux;
30✔
1251
    }
1252
  }
1253

1254
  // Then, compute the adjoint source for each source region
1255
#pragma omp parallel for
46✔
1256
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
169,270✔
1257
    for (int g = 0; g < negroups_; g++) {
338,480✔
1258
      double flux = source_regions_.scalar_flux_final(sr, g);
169,240✔
1259
      if (flux <= ZERO_FLUX_CUTOFF * max_flux) {
169,240✔
1260
        source_regions_.external_source(sr, g) = 0.0;
325✔
1261
      } else {
1262
        source_regions_.external_source(sr, g) = 1.0 / flux;
168,915!
1263
        if (!std::isfinite(source_regions_.external_source(sr, g))) {
168,915!
1264
          // If the flux is NaN or Inf, set the adjoint source to zero
1265
          source_regions_.external_source(sr, g) = 0.0;
1266
        }
1267
      }
1268
      if (flux > 0.0) {
169,240✔
1269
        source_regions_.external_source_present(sr) = 1;
168,915✔
1270
      }
1271
      source_regions_.scalar_flux_final(sr, g) = 0.0;
169,240✔
1272
    }
1273
  }
1274

1275
  // "Small" source regions in OpenMC are defined as those that are hit by
1276
  // MIN_HITS_PER_BATCH rays or fewer each batch. These regions typically have
1277
  // very small volumes combined with a low aspect ratio, and are often
1278
  // generated when applying a source region mesh that clips the edge of a
1279
  // curved surface. As perhaps only a few rays will visit these regions over
1280
  // the entire forward simulation, the forward flux estimates are extremely
1281
  // noisy and unreliable. In some cases, the noise may make the forward fluxes
1282
  // extremely low, leading to unphysically large adjoint source terms,
1283
  // resulting in weight windows that aggressively try to drive particles
1284
  // towards these regions. To fix this, we simply filter out any "small" source
1285
  // regions from consideration. If a source region is "small", we
1286
  // set its adjoint source to zero. This adds negligible bias to the adjoint
1287
  // flux solution, as the true total adjoint source contribution from small
1288
  // regions is likely to be negligible.
1289
#pragma omp parallel for
46✔
1290
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
169,270✔
1291
    if (source_regions_.is_small(sr)) {
169,240!
1292
      for (int g = 0; g < negroups_; g++) {
×
1293
        source_regions_.external_source(sr, g) = 0.0;
1294
      }
1295
      source_regions_.external_source_present(sr) = 0;
1296
    }
1297
  }
1298

1299
  // Divide the fixed source term by sigma t (to save time when applying each
1300
  // iteration)
1301
#pragma omp parallel for
46✔
1302
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
169,270✔
1303
    int material = source_regions_.material(sr);
169,240!
1304
    int temp = source_regions_.temperature_idx(sr);
169,240!
1305
    if (material == MATERIAL_VOID) {
169,240!
1306
      continue;
1307
    }
1308
    for (int g = 0; g < negroups_; g++) {
338,480✔
1309
      double sigma_t =
169,240✔
1310
        sigma_t_[(material * ntemperature_ + temp) * negroups_ + g] *
169,240!
1311
        source_regions_.density_mult(sr);
169,240!
1312
      source_regions_.external_source(sr, g) /= sigma_t;
169,240!
1313
      if (!std::isfinite(source_regions_.external_source(sr, g))) {
169,240!
1314
        // If the flux is NaN or Inf, set the adjoint source to zero
1315
        source_regions_.external_source(sr, g) = 0.0;
1316
      }
1317
    }
1318
  }
1319

1320
  if (fw_cadis_local_) {
76✔
1321
// Only external sources that have a non-mesh type tally task should remain
1322
// non-zero. Everything else gets zero'd out.
1323
#pragma omp parallel for
9✔
1324
    for (int64_t sr = 0; sr < n_source_regions(); sr++) {
13,726✔
1325

1326
      // If there is already no external source, don't need to do anything
1327
      if (source_regions_.external_source_present(sr) == 0) {
13,720!
1328
        continue;
1329
      }
1330

1331
      // If there is an adjoint source term here, then we need to check it.
1332

1333
      // We will track if ANY group has a valid local FW-CADIS source term
1334
      bool has_any_sources = false;
1335

1336
      // Now, loop over groups
1337
      for (int g = 0; g < negroups_; g++) {
27,440✔
1338

1339
        // If there are no tally tasks associated with this source element
1340
        // then it is not a local FW-CADIS source, so we continue to the next
1341
        // group
1342
        if (source_regions_.tally_task(sr, g).empty()) {
13,720!
1343
          source_regions_.external_source(sr, g) = 0.0;
1344
          continue;
1345
        }
1346

1347
        // If there are tally tasks, we can through them and check if
1348
        // any of them are local FW-CADIS targets.
1349

1350
        // We track if ANY of the tasks are local FW-CADIS target tallies
1351
        bool local_fw_cadis_target_region = false;
27,360✔
1352

1353
        // Now we loop through
1354
        for (const auto& task : source_regions_.tally_task(sr, g)) {
27,360✔
1355
          Tally& tally {*model::tallies[task.tally_idx]};
13,720✔
1356
          const auto t_id = tally.id();
13,720✔
1357

1358
          // Search for target tallies
1359
          if (std::find(fw_cadis_local_targets_.begin(),
13,720✔
1360
                fw_cadis_local_targets_.end(),
1361
                t_id) != fw_cadis_local_targets_.end()) {
13,720✔
1362
            local_fw_cadis_target_region = true;
80✔
1363
            break;
80✔
1364
          }
1365
        }
1366

1367
        // If ANY of the tasks is a local FW-CADIS target,
1368
        // Then we keep the source term and set that this
1369
        // source region has a valid FW-CADIS source term.
1370
        // Otherwise, we zero out the source term.
1371
        if (local_fw_cadis_target_region) {
80✔
1372
          has_any_sources = true;
1373
        } else {
1374
          source_regions_.external_source(sr, g) = 0.0;
13,640✔
1375
        }
1376
      } // End loop over groups
1377

1378
      // If there were any valid FW-CADIS source terms for any
1379
      // of the groups, then the SR as a whole counts as a source
1380
      if (has_any_sources) {
13,720✔
1381
        source_regions_.external_source_present(sr) = 1;
80✔
1382
      } else {
1383
        source_regions_.external_source_present(sr) = 0;
13,640✔
1384
      }
1385
    } // End loop over source regions
1386
  } // End local FW-CADIS logic
1387
}
76✔
1388

1389
void FlatSourceDomain::set_local_adjoint_sources()
15✔
1390
{
1391
  // Set the external source to user-specified adjoint sources.
1392
  convert_external_sources(true);
15✔
1393
}
15✔
1394

1395
void FlatSourceDomain::transpose_scattering_matrix()
106✔
1396
{
1397
  // Transpose the inner two dimensions for each material
1398
#pragma omp parallel for
64✔
1399
  for (int m = 0; m < n_materials_; ++m) {
162✔
1400
    for (int t = 0; t < ntemperature_; t++) {
240✔
1401
      int material_offset = (m * ntemperature_ + t) * negroups_ * negroups_;
120✔
1402
      for (int i = 0; i < negroups_; ++i) {
312✔
1403
        for (int j = i + 1; j < negroups_; ++j) {
444✔
1404
          // Calculate indices of the elements to swap
1405
          int idx1 = material_offset + i * negroups_ + j;
252✔
1406
          int idx2 = material_offset + j * negroups_ + i;
252✔
1407

1408
          // Swap the elements to transpose the matrix
1409
          std::swap(sigma_s_[idx1], sigma_s_[idx2]);
252✔
1410
        }
1411
      }
1412
    }
1413
  }
1414
}
106✔
1415

UNCOV
1416
void FlatSourceDomain::serialize_final_fluxes(vector<double>& flux)
×
1417
{
1418
  // Ensure array is correct size
UNCOV
1419
  flux.resize(n_source_regions() * negroups_);
×
1420
// Serialize the final fluxes for output
1421
#pragma omp parallel for
1422
  for (int64_t se = 0; se < n_source_elements(); se++) {
×
1423
    flux[se] = source_regions_.scalar_flux_final(se);
1424
  }
UNCOV
1425
}
×
1426

1427
void FlatSourceDomain::apply_mesh_to_cell_instances(int32_t i_cell,
1,217✔
1428
  int32_t mesh_idx, int target_material_id, const vector<int32_t>& instances,
1429
  bool is_target_void)
1430
{
1431
  Cell& cell = *model::cells[i_cell];
1,217!
1432
  if (cell.type_ != Fill::MATERIAL)
1,217!
1433
    return;
1434
  for (int32_t j : instances) {
217,054✔
1435
    int cell_material_idx = cell.material(j);
215,837!
1436
    int cell_material_id = (cell_material_idx == C_NONE)
215,837✔
1437
                             ? C_NONE
215,837✔
1438
                             : model::materials[cell_material_idx]->id();
215,836✔
1439

1440
    if ((target_material_id == C_NONE && !is_target_void) ||
215,837✔
1441
        cell_material_id == target_material_id) {
1442
      int64_t sr = source_region_offsets_[i_cell] + j;
185,837!
1443
      // Check if the key is already present in the mesh_map_
1444
      if (mesh_map_.find(sr) != mesh_map_.end()) {
185,837!
UNCOV
1445
        fatal_error(fmt::format("Source region {} already has mesh idx {} "
×
1446
                                "applied, but trying to apply mesh idx {}",
UNCOV
1447
          sr, mesh_map_[sr], mesh_idx));
×
1448
      }
1449
      // If the SR has not already been assigned, then we can write to it
1450
      mesh_map_[sr] = mesh_idx;
185,837✔
1451
    }
1452
  }
1453
}
1454

1455
void FlatSourceDomain::apply_mesh_to_cell_and_children(int32_t i_cell,
1,036✔
1456
  int32_t mesh_idx, int32_t target_material_id, bool is_target_void)
1457
{
1458
  Cell& cell = *model::cells[i_cell];
1,036✔
1459

1460
  if (cell.type_ == Fill::MATERIAL) {
1,036✔
1461
    vector<int> instances(cell.n_instances());
885✔
1462
    std::iota(instances.begin(), instances.end(), 0);
885✔
1463
    apply_mesh_to_cell_instances(
885✔
1464
      i_cell, mesh_idx, target_material_id, instances, is_target_void);
1465
  } else if (target_material_id == C_NONE && !is_target_void) {
1,036✔
1466
    for (int j = 0; j < cell.n_instances(); j++) {
182✔
1467
      std::unordered_map<int32_t, vector<int32_t>> cell_instance_list =
91✔
1468
        cell.get_contained_cells(j, nullptr);
91✔
1469
      for (const auto& pair : cell_instance_list) {
423✔
1470
        int32_t i_child_cell = pair.first;
332✔
1471
        apply_mesh_to_cell_instances(i_child_cell, mesh_idx, target_material_id,
332✔
1472
          pair.second, is_target_void);
332✔
1473
      }
1474
    }
91✔
1475
  }
1476
}
1,036✔
1477

1478
void FlatSourceDomain::apply_meshes()
792✔
1479
{
1480
  // Skip if there are no mappings between mesh IDs and domains
1481
  if (mesh_domain_map_.empty())
792✔
1482
    return;
1483

1484
  // Loop over meshes
1485
  for (int mesh_idx = 0; mesh_idx < model::meshes.size(); mesh_idx++) {
782✔
1486
    Mesh* mesh = model::meshes[mesh_idx].get();
436✔
1487
    int mesh_id = mesh->id();
436✔
1488

1489
    // Skip if mesh id is not present in the map
1490
    if (mesh_domain_map_.find(mesh_id) == mesh_domain_map_.end())
436✔
1491
      continue;
30✔
1492

1493
    // Loop over domains associated with the mesh
1494
    for (auto& domain : mesh_domain_map_[mesh_id]) {
812✔
1495
      Source::DomainType domain_type = domain.first;
406✔
1496
      int domain_id = domain.second;
406✔
1497

1498
      if (domain_type == Source::DomainType::MATERIAL) {
406✔
1499
        for (int i_cell = 0; i_cell < model::cells.size(); i_cell++) {
180✔
1500
          if (domain_id == C_NONE) {
150!
UNCOV
1501
            apply_mesh_to_cell_and_children(i_cell, mesh_idx, domain_id, true);
×
1502
          } else {
1503
            apply_mesh_to_cell_and_children(i_cell, mesh_idx, domain_id, false);
150✔
1504
          }
1505
        }
1506
      } else if (domain_type == Source::DomainType::CELL) {
376✔
1507
        int32_t i_cell = model::cell_map[domain_id];
30✔
1508
        apply_mesh_to_cell_and_children(i_cell, mesh_idx, C_NONE, false);
30✔
1509
      } else if (domain_type == Source::DomainType::UNIVERSE) {
346!
1510
        int32_t i_universe = model::universe_map[domain_id];
346✔
1511
        Universe& universe = *model::universes[i_universe];
346✔
1512
        for (int32_t i_cell : universe.cells_) {
1,202✔
1513
          apply_mesh_to_cell_and_children(i_cell, mesh_idx, C_NONE, false);
856✔
1514
        }
1515
      }
1516
    }
1517
  }
1518
}
1519

1520
SourceRegionHandle FlatSourceDomain::get_subdivided_source_region_handle(
1,495,487,076✔
1521
  SourceRegionKey sr_key, Position r, Direction u)
1522
{
1523
  // Case 1: Check if the source region key is already present in the permanent
1524
  // map. This is the most common condition, as any source region visited in a
1525
  // previous power iteration will already be present in the permanent map. If
1526
  // the source region key is found, we translate the key into a specific 1D
1527
  // source region index and return a handle its position in the
1528
  // source_regions_ vector.
1529
  auto it = source_region_map_.find(sr_key);
1,495,487,076✔
1530
  if (it != source_region_map_.end()) {
1,495,487,076✔
1531
    int64_t sr = it->second;
1,445,661,016✔
1532
    return source_regions_.get_source_region_handle(sr);
1,445,661,016✔
1533
  }
1534

1535
  // Case 2: Check if the source region key is present in the temporary (thread
1536
  // safe) map. This is a common occurrence in the first power iteration when
1537
  // the source region has already been visited already by some other ray. We
1538
  // begin by locking the temporary map before any operations are performed. The
1539
  // lock is not global over the full data structure -- it will be dependent on
1540
  // which key is used.
1541
  discovered_source_regions_.lock(sr_key);
49,826,060✔
1542

1543
  // If the key is found in the temporary map, then we return a handle to the
1544
  // source region that is stored in the temporary map.
1545
  if (discovered_source_regions_.contains(sr_key)) {
49,826,060✔
1546
    SourceRegionHandle handle {discovered_source_regions_[sr_key]};
47,398,680✔
1547
    discovered_source_regions_.unlock(sr_key);
47,398,680✔
1548
    return handle;
47,398,680✔
1549
  }
1550

1551
  // Case 3: The source region key is not present anywhere, but it is only due
1552
  // to floating point artifacts. These artifacts occur when the overlaid mesh
1553
  // overlaps with actual geometry surfaces. In these cases, roundoff error may
1554
  // result in the ray tracer detecting an additional (very short) segment
1555
  // though a mesh bin that is actually past the physical source region
1556
  // boundary. This is a result of the the multi-level ray tracing treatment in
1557
  // OpenMC, which depending on the number of universes in the hierarchy etc can
1558
  // result in the wrong surface being selected as the nearest. This can happen
1559
  // in a lattice when there are two directions that both are very close in
1560
  // distance, within the tolerance of FP_REL_PRECISION, and the are thus
1561
  // treated as being equivalent so alternative logic is used. However, when we
1562
  // go and ray trace on this with the mesh tracer we may go past the surface
1563
  // bounding the current source region.
1564
  //
1565
  // To filter out this case, before we create the new source region, we double
1566
  // check that the actual starting point of this segment (r) is still in the
1567
  // same geometry source region that we started in. If an artifact is detected,
1568
  // we discard the segment (and attenuation through it) as it is not really a
1569
  // valid source region and will have only an infinitessimally small cell
1570
  // combined with the mesh bin. Thankfully, this is a fairly rare condition,
1571
  // and only triggers for very short ray lengths. It can be fixed by decreasing
1572
  // the value of FP_REL_PRECISION in constants.h, but this may have unknown
1573
  // consequences for the general ray tracer, so for now we do the below sanity
1574
  // checks before generating phantom source regions. A significant extra cost
1575
  // is incurred in instantiating the GeometryState object and doing a cell
1576
  // lookup, but again, this is going to be an extremely rare thing to check
1577
  // after the first power iteration has completed.
1578

1579
  // Sanity check on source region id
1580
  GeometryState gs;
2,427,380✔
1581
  gs.r() = r + TINY_BIT * u;
2,427,380✔
1582
  gs.u() = {1.0, 0.0, 0.0};
2,427,380✔
1583
  exhaustive_find_cell(gs);
2,427,380✔
1584
  int64_t sr_found = lookup_base_source_region_idx(gs);
2,427,380✔
1585
  if (sr_found != sr_key.base_source_region_id) {
2,427,380✔
1586
    discovered_source_regions_.unlock(sr_key);
121✔
1587
    SourceRegionHandle handle;
121✔
1588
    handle.is_numerical_fp_artifact_ = true;
121✔
1589
    return handle;
121✔
1590
  }
1591

1592
  // Sanity check on mesh bin
1593
  int mesh_idx = lookup_mesh_idx(sr_key.base_source_region_id);
2,427,259✔
1594
  if (mesh_idx == C_NONE) {
2,427,259✔
1595
    if (sr_key.mesh_bin != 0) {
393,756!
UNCOV
1596
      discovered_source_regions_.unlock(sr_key);
×
UNCOV
1597
      SourceRegionHandle handle;
×
UNCOV
1598
      handle.is_numerical_fp_artifact_ = true;
×
UNCOV
1599
      return handle;
×
1600
    }
1601
  } else {
1602
    Mesh* mesh = model::meshes[mesh_idx].get();
2,033,503✔
1603
    int bin_found = mesh->get_bin(r + TINY_BIT * u);
2,033,503✔
1604
    if (bin_found != sr_key.mesh_bin) {
2,033,503!
UNCOV
1605
      discovered_source_regions_.unlock(sr_key);
×
UNCOV
1606
      SourceRegionHandle handle;
×
UNCOV
1607
      handle.is_numerical_fp_artifact_ = true;
×
UNCOV
1608
      return handle;
×
1609
    }
1610
  }
1611

1612
  // Case 4: The source region key is valid, but is not present anywhere. This
1613
  // condition only occurs the first time the source region is discovered
1614
  // (typically in the first power iteration). In this case, we need to handle
1615
  // creation of the new source region and its storage into the parallel map.
1616
  // Additionally, we need to determine the source region's material, initialize
1617
  // the starting scalar flux guess, and apply any known external sources.
1618

1619
  // Call the basic constructor for the source region and store in the parallel
1620
  // map.
1621
  bool is_linear = RandomRay::source_shape_ != RandomRaySourceShape::FLAT;
2,427,259✔
1622
  SourceRegion* sr_ptr =
2,427,259✔
1623
    discovered_source_regions_.emplace(sr_key, {negroups_, is_linear});
2,427,259✔
1624
  SourceRegionHandle handle {*sr_ptr};
2,427,259✔
1625

1626
  // Determine the material
1627
  int gs_i_cell = gs.lowest_coord().cell();
2,427,259!
1628
  Cell& cell = *model::cells[gs_i_cell];
2,427,259!
1629
  int material = cell.material(gs.cell_instance());
2,427,259!
1630
  int temp = 0;
2,427,259✔
1631

1632
  // If material total XS is extremely low, just set it to void to avoid
1633
  // problems with 1/Sigma_t
1634
  if (material != MATERIAL_VOID) {
2,427,259✔
1635
    temp = data::mg.macro_xs_[material].get_temperature_index(
4,810,390✔
1636
      cell.sqrtkT(gs.cell_instance()));
2,405,195✔
1637
    for (int g = 0; g < negroups_; g++) {
5,161,511✔
1638
      double sigma_t =
2,756,404✔
1639
        sigma_t_[(material * ntemperature_ + temp) * negroups_ + g];
2,756,404✔
1640
      if (sigma_t < MINIMUM_MACRO_XS) {
2,756,404✔
1641
        material = MATERIAL_VOID;
1642
        temp = 0;
1643
        break;
1644
      }
1645
    }
1646
  }
1647

1648
  handle.material() = material;
2,427,259✔
1649
  handle.temperature_idx() = temp;
2,427,259✔
1650

1651
  handle.density_mult() = cell.density_mult(gs.cell_instance());
2,427,259✔
1652

1653
  // Store the mesh index (if any) assigned to this source region
1654
  handle.mesh() = mesh_idx;
2,427,259✔
1655

1656
  if (settings::run_mode == RunMode::FIXED_SOURCE) {
2,427,259✔
1657
    // Determine if there are any volumetric sources, and apply them.
1658
    // Volumetric sources are specifc only to the base SR idx.
1659
    auto it_vol =
2,371,236✔
1660
      external_volumetric_source_map_.find(sr_key.base_source_region_id);
2,371,236✔
1661
    if (it_vol != external_volumetric_source_map_.end()) {
2,371,236✔
1662
      const vector<int>& vol_sources = it_vol->second;
4,950✔
1663
      for (int src_idx : vol_sources) {
9,900✔
1664
        apply_external_source_to_source_region(src_idx, handle);
4,950✔
1665
      }
1666
    }
1667

1668
    // Determine if there are any point sources, and apply them.
1669
    // Point sources are specific to the source region key.
1670
    auto it_point = external_point_source_map_.find(sr_key);
2,371,236✔
1671
    if (it_point != external_point_source_map_.end()) {
2,371,236✔
1672
      const vector<int>& point_sources = it_point->second;
12✔
1673
      for (int src_idx : point_sources) {
24✔
1674
        apply_external_source_to_source_region(src_idx, handle);
12✔
1675
      }
1676
    }
1677

1678
    // Divide external source term by sigma_t
1679
    if (material != C_NONE) {
2,371,236✔
1680
      for (int g = 0; g < negroups_; g++) {
4,744,127✔
1681
        double sigma_t =
2,395,043✔
1682
          sigma_t_[(material * ntemperature_ + temp) * negroups_ + g] *
2,395,043✔
1683
          handle.density_mult();
2,395,043✔
1684
        handle.external_source(g) /= sigma_t;
2,395,043✔
1685
      }
1686
    }
1687
  }
1688

1689
  // Compute the combined source term
1690
  update_single_neutron_source(handle);
2,427,259✔
1691

1692
  // Unlock the parallel map. Note: we may be tempted to release
1693
  // this lock earlier, and then just use the source region's lock to protect
1694
  // the flux/source initialization stages above. However, the rest of the code
1695
  // only protects updates to the new flux and volume fields, and assumes that
1696
  // the source is constant for the duration of transport. Thus, using just the
1697
  // source region's lock by itself would result in other threads potentially
1698
  // reading from the source before it is computed, as they won't use the lock
1699
  // when only reading from the SR's source. It would be expensive to protect
1700
  // those operations, whereas generating the SR is only done once, so we just
1701
  // hold the map's bucket lock until the source region is fully initialized.
1702
  discovered_source_regions_.unlock(sr_key);
2,427,259✔
1703

1704
  return handle;
2,427,259✔
1705
}
2,427,380✔
1706

1707
void FlatSourceDomain::finalize_discovered_source_regions()
15,962✔
1708
{
1709
  // Extract keys for entries with a valid volume.
1710
  vector<SourceRegionKey> keys;
15,962✔
1711
  for (const auto& pair : discovered_source_regions_) {
4,886,442✔
1712
    if (pair.second.volume_ > 0.0) {
2,427,259✔
1713
      keys.push_back(pair.first);
2,325,817✔
1714
    }
1715
  }
1716

1717
  if (!keys.empty()) {
15,962✔
1718
    // Sort the keys, so as to ensure reproducible ordering given that source
1719
    // regions may have been added to discovered_source_regions_ in an arbitrary
1720
    // order due to shared memory threading.
1721
    std::sort(keys.begin(), keys.end());
914✔
1722

1723
    // Remember the index of the first new source region
1724
    int64_t start_sr_id = source_regions_.n_source_regions();
914✔
1725

1726
    // Append the source regions in the sorted key order.
1727
    for (const auto& key : keys) {
2,326,731✔
1728
      const SourceRegion& sr = discovered_source_regions_[key];
2,325,817✔
1729
      source_region_map_[key] = source_regions_.n_source_regions();
2,325,817✔
1730
      source_regions_.push_back(sr);
2,325,817✔
1731
    }
1732

1733
    // Map all new source regions to tallies
1734
    convert_source_regions_to_tallies(start_sr_id);
914✔
1735
  }
1736

1737
  discovered_source_regions_.clear();
15,962✔
1738
}
15,962✔
1739

1740
// This is the "diagonal stabilization" technique developed by Gunow et al. in:
1741
//
1742
// Geoffrey Gunow, Benoit Forget, Kord Smith, Stabilization of multi-group
1743
// neutron transport with transport-corrected cross-sections, Annals of Nuclear
1744
// Energy, Volume 126, 2019, Pages 211-219, ISSN 0306-4549,
1745
// https://doi.org/10.1016/j.anucene.2018.10.036.
1746
void FlatSourceDomain::apply_transport_stabilization()
15,962✔
1747
{
1748
  // Don't do anything if all in-group scattering
1749
  // cross sections are positive
1750
  if (!is_transport_stabilization_needed_) {
15,962✔
1751
    return;
1752
  }
1753

1754
  // Apply the stabilization factor to all source elements
1755
#pragma omp parallel for
120✔
1756
  for (int64_t sr = 0; sr < n_source_regions(); sr++) {
1,300✔
1757
    int material = source_regions_.material(sr);
1,200!
1758
    int temp = source_regions_.temperature_idx(sr);
1,200!
1759
    double density_mult = source_regions_.density_mult(sr);
1,200!
1760
    if (material == MATERIAL_VOID) {
1,200!
1761
      continue;
1762
    }
1763
    for (int g = 0; g < negroups_; g++) {
85,200✔
1764
      // Only apply stabilization if the diagonal (in-group) scattering XS is
1765
      // negative
1766
      double sigma_s =
84,000✔
1767
        sigma_s_[((material * ntemperature_ + temp) * negroups_ + g) *
84,000✔
1768
                   negroups_ +
84,000✔
1769
                 g] *
84,000✔
1770
        density_mult;
84,000✔
1771
      if (sigma_s < 0.0) {
84,000✔
1772
        double sigma_t =
22,000✔
1773
          sigma_t_[(material * ntemperature_ + temp) * negroups_ + g] *
22,000✔
1774
          density_mult;
22,000✔
1775
        double phi_new = source_regions_.scalar_flux_new(sr, g);
22,000✔
1776
        double phi_old = source_regions_.scalar_flux_old(sr, g);
22,000✔
1777

1778
        // Equation 18 in the above Gunow et al. 2019 paper. For a default
1779
        // rho of 1.0, this ensures there are no negative diagonal elements
1780
        // in the iteration matrix. A lesser rho could be used (or exposed
1781
        // as a user input parameter) to reduce the negative impact on
1782
        // convergence rate though would need to be experimentally tested to see
1783
        // if it doesn't become unstable. rho = 1.0 is good as it gives the
1784
        // highest assurance of stability, and the impacts on convergence rate
1785
        // are pretty mild.
1786
        double D = diagonal_stabilization_rho_ * sigma_s / sigma_t;
22,000✔
1787

1788
        // Equation 16 in the above Gunow et al. 2019 paper
1789
        source_regions_.scalar_flux_new(sr, g) =
22,000✔
1790
          (phi_new - D * phi_old) / (1.0 - D);
22,000✔
1791
      }
1792
    }
1793
  }
1794
}
1795

1796
// Determines the base source region index (i.e., a material filled cell
1797
// instance) that corresponds to a particular location in the geometry. Requires
1798
// that the "gs" object passed in has already been initialized and has called
1799
// find_cell etc.
1800
int64_t FlatSourceDomain::lookup_base_source_region_idx(
1,009,291,276✔
1801
  const GeometryState& gs) const
1802
{
1803
  int i_cell = gs.lowest_coord().cell();
1,009,291,276✔
1804
  int64_t sr = source_region_offsets_[i_cell] + gs.cell_instance();
1,009,291,276✔
1805
  return sr;
1,009,291,276✔
1806
}
1807

1808
// Determines the index of the mesh (if any) that has been applied
1809
// to a particular base source region index.
1810
int FlatSourceDomain::lookup_mesh_idx(int64_t sr) const
1,009,291,155✔
1811
{
1812
  int mesh_idx = C_NONE;
1,009,291,155✔
1813
  auto mesh_it = mesh_map_.find(sr);
1,009,291,155✔
1814
  if (mesh_it != mesh_map_.end()) {
1,009,291,155✔
1815
    mesh_idx = mesh_it->second;
543,833,690✔
1816
  }
1817
  return mesh_idx;
1,009,291,155✔
1818
}
1819

1820
// Determines the source region key that corresponds to a particular location in
1821
// the geometry. This takes into account both the base source region index as
1822
// well as the mesh bin if a mesh is applied to this source region for
1823
// subdivision.
1824
SourceRegionKey FlatSourceDomain::lookup_source_region_key(
2,317,336✔
1825
  const GeometryState& gs) const
1826
{
1827
  int64_t sr = lookup_base_source_region_idx(gs);
2,317,336✔
1828
  int64_t mesh_bin = lookup_mesh_bin(sr, gs.r());
2,317,336✔
1829
  return SourceRegionKey {sr, mesh_bin};
2,317,336✔
1830
}
1831

1832
// Determines the mesh bin that corresponds to a particular base source region
1833
// index and position.
1834
int64_t FlatSourceDomain::lookup_mesh_bin(int64_t sr, Position r) const
2,317,336✔
1835
{
1836
  int mesh_idx = lookup_mesh_idx(sr);
2,317,336✔
1837
  int mesh_bin = 0;
2,317,336✔
1838
  if (mesh_idx != C_NONE) {
2,317,336✔
1839
    mesh_bin = model::meshes[mesh_idx]->get_bin(r);
1,453,836✔
1840
  }
1841
  return mesh_bin;
2,317,336✔
1842
}
1843

1844
} // 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

© 2026 Coveralls, Inc