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

openmc-dev / openmc / 13591584831

28 Feb 2025 03:46PM UTC coverage: 85.051% (+0.3%) from 84.722%
13591584831

Pull #3067

github

web-flow
Merge 08055e996 into c26fde666
Pull Request #3067: Implement user-configurable random number stride

36 of 44 new or added lines in 8 files covered. (81.82%)

3588 existing lines in 111 files now uncovered.

51062 of 60037 relevant lines covered (85.05%)

32650986.73 hits per line

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

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

3
#include "openmc/constants.h"
4
#include "openmc/geometry.h"
5
#include "openmc/message_passing.h"
6
#include "openmc/mgxs_interface.h"
7
#include "openmc/random_ray/flat_source_domain.h"
8
#include "openmc/random_ray/linear_source_domain.h"
9
#include "openmc/search.h"
10
#include "openmc/settings.h"
11
#include "openmc/simulation.h"
12

13
#include "openmc/distribution_spatial.h"
14
#include "openmc/random_dist.h"
15
#include "openmc/source.h"
16

17
namespace openmc {
18

19
//==============================================================================
20
// Non-method functions
21
//==============================================================================
22

23
// returns 1 - exp(-tau)
24
// Equivalent to -(_expm1f(-tau)), but faster
25
// Written by Colin Josey.
26
float cjosey_exponential(float tau)
580,664,664✔
27
{
28
  constexpr float c1n = -1.0000013559236386308f;
580,664,664✔
29
  constexpr float c2n = 0.23151368626911062025f;
580,664,664✔
30
  constexpr float c3n = -0.061481916409314966140f;
580,664,664✔
31
  constexpr float c4n = 0.0098619906458127653020f;
580,664,664✔
32
  constexpr float c5n = -0.0012629460503540849940f;
580,664,664✔
33
  constexpr float c6n = 0.00010360973791574984608f;
580,664,664✔
34
  constexpr float c7n = -0.000013276571933735820960f;
580,664,664✔
35

36
  constexpr float c0d = 1.0f;
580,664,664✔
37
  constexpr float c1d = -0.73151337729389001396f;
580,664,664✔
38
  constexpr float c2d = 0.26058381273536471371f;
580,664,664✔
39
  constexpr float c3d = -0.059892419041316836940f;
580,664,664✔
40
  constexpr float c4d = 0.0099070188241094279067f;
580,664,664✔
41
  constexpr float c5d = -0.0012623388962473160860f;
580,664,664✔
42
  constexpr float c6d = 0.00010361277635498731388f;
580,664,664✔
43
  constexpr float c7d = -0.000013276569500666698498f;
580,664,664✔
44

45
  float x = -tau;
580,664,664✔
46

47
  float den = c7d;
580,664,664✔
48
  den = den * x + c6d;
580,664,664✔
49
  den = den * x + c5d;
580,664,664✔
50
  den = den * x + c4d;
580,664,664✔
51
  den = den * x + c3d;
580,664,664✔
52
  den = den * x + c2d;
580,664,664✔
53
  den = den * x + c1d;
580,664,664✔
54
  den = den * x + c0d;
580,664,664✔
55

56
  float num = c7n;
580,664,664✔
57
  num = num * x + c6n;
580,664,664✔
58
  num = num * x + c5n;
580,664,664✔
59
  num = num * x + c4n;
580,664,664✔
60
  num = num * x + c3n;
580,664,664✔
61
  num = num * x + c2n;
580,664,664✔
62
  num = num * x + c1n;
580,664,664✔
63
  num = num * x;
580,664,664✔
64

65
  return num / den;
580,664,664✔
66
}
67

68
// The below two functions (exponentialG and exponentialG2) were developed
69
// by Colin Josey. The implementation of these functions is closely based
70
// on the OpenMOC versions of these functions. The OpenMOC license is given
71
// below:
72

73
// Copyright (C) 2012-2023 Massachusetts Institute of Technology and OpenMOC
74
// contributors
75
//
76
// Permission is hereby granted, free of charge, to any person obtaining a copy
77
// of this software and associated documentation files (the "Software"), to deal
78
// in the Software without restriction, including without limitation the rights
79
// to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
80
// copies of the Software, and to permit persons to whom the Software is
81
// furnished to do so, subject to the following conditions:
82
//
83
// The above copyright notice and this permission notice shall be included in
84
// all copies or substantial portions of the Software.
85
//
86
// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
87
// IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
88
// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
89
// AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
90
// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
91
// OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
92
// SOFTWARE.
93

94
// Computes y = 1/x-(1-exp(-x))/x**2 using a 5/6th order rational
95
// approximation. It is accurate to 2e-7 over [0, 1e5]. Developed by Colin
96
// Josey using Remez's algorithm, with original implementation in OpenMOC at:
97
// https://github.com/mit-crpg/OpenMOC/blob/develop/src/exponentials.h
98
float exponentialG(float tau)
867,898,140✔
99
{
100
  // Numerator coefficients in rational approximation for 1/x - (1 - exp(-x)) /
101
  // x^2
102
  constexpr float d0n = 0.5f;
867,898,140✔
103
  constexpr float d1n = 0.176558112351595f;
867,898,140✔
104
  constexpr float d2n = 0.04041584305811143f;
867,898,140✔
105
  constexpr float d3n = 0.006178333902037397f;
867,898,140✔
106
  constexpr float d4n = 0.0006429894635552992f;
867,898,140✔
107
  constexpr float d5n = 0.00006064409107557148f;
867,898,140✔
108

109
  // Denominator coefficients in rational approximation for 1/x - (1 - exp(-x))
110
  // / x^2
111
  constexpr float d0d = 1.0f;
867,898,140✔
112
  constexpr float d1d = 0.6864462055546078f;
867,898,140✔
113
  constexpr float d2d = 0.2263358514260129f;
867,898,140✔
114
  constexpr float d3d = 0.04721469893686252f;
867,898,140✔
115
  constexpr float d4d = 0.006883236664917246f;
867,898,140✔
116
  constexpr float d5d = 0.0007036272419147752f;
867,898,140✔
117
  constexpr float d6d = 0.00006064409107557148f;
867,898,140✔
118

119
  float x = tau;
867,898,140✔
120

121
  float num = d5n;
867,898,140✔
122
  num = num * x + d4n;
867,898,140✔
123
  num = num * x + d3n;
867,898,140✔
124
  num = num * x + d2n;
867,898,140✔
125
  num = num * x + d1n;
867,898,140✔
126
  num = num * x + d0n;
867,898,140✔
127

128
  float den = d6d;
867,898,140✔
129
  den = den * x + d5d;
867,898,140✔
130
  den = den * x + d4d;
867,898,140✔
131
  den = den * x + d3d;
867,898,140✔
132
  den = den * x + d2d;
867,898,140✔
133
  den = den * x + d1d;
867,898,140✔
134
  den = den * x + d0d;
867,898,140✔
135

136
  return num / den;
867,898,140✔
137
}
138

139
// Computes G2 : y = 2/3 - (1 + 2/x) * (1/x + 0.5 - (1 + 1/x) * (1-exp(-x)) /
140
// x) using a 5/5th order rational approximation. It is accurate to 1e-6 over
141
// [0, 1e6]. Developed by Colin Josey using Remez's algorithm, with original
142
// implementation in OpenMOC at:
143
// https://github.com/mit-crpg/OpenMOC/blob/develop/src/exponentials.h
144
float exponentialG2(float tau)
867,898,140✔
145
{
146

147
  // Coefficients for numerator in rational approximation
148
  constexpr float g1n = -0.08335775885589858f;
867,898,140✔
149
  constexpr float g2n = -0.003603942303847604f;
867,898,140✔
150
  constexpr float g3n = 0.0037673183263550827f;
867,898,140✔
151
  constexpr float g4n = 0.00001124183494990467f;
867,898,140✔
152
  constexpr float g5n = 0.00016837426505799449f;
867,898,140✔
153

154
  // Coefficients for denominator in rational approximation
155
  constexpr float g1d = 0.7454048371823628f;
867,898,140✔
156
  constexpr float g2d = 0.23794300531408347f;
867,898,140✔
157
  constexpr float g3d = 0.05367250964303789f;
867,898,140✔
158
  constexpr float g4d = 0.006125197988351906f;
867,898,140✔
159
  constexpr float g5d = 0.0010102514456857377f;
867,898,140✔
160

161
  float x = tau;
867,898,140✔
162

163
  float num = g5n;
867,898,140✔
164
  num = num * x + g4n;
867,898,140✔
165
  num = num * x + g3n;
867,898,140✔
166
  num = num * x + g2n;
867,898,140✔
167
  num = num * x + g1n;
867,898,140✔
168
  num = num * x;
867,898,140✔
169

170
  float den = g5d;
867,898,140✔
171
  den = den * x + g4d;
867,898,140✔
172
  den = den * x + g3d;
867,898,140✔
173
  den = den * x + g2d;
867,898,140✔
174
  den = den * x + g1d;
867,898,140✔
175
  den = den * x + 1.0f;
867,898,140✔
176

177
  return num / den;
867,898,140✔
178
}
179

180
// Implementation of the Fisher-Yates shuffle algorithm.
181
// Algorithm adapted from:
182
//    https://en.cppreference.com/w/cpp/algorithm/random_shuffle#Version_3
183
void fisher_yates_shuffle(vector<int64_t>& arr, uint64_t* seed)
1,728,000✔
184
{
185
  // Loop over the array from the last element down to the second
186
  for (int i = arr.size() - 1; i > 0; --i) {
7,452,000✔
187
    // Generate a random index in the range [0, i]
188
    int j = uniform_int_distribution(0, i, seed);
5,724,000✔
189
    std::swap(arr[i], arr[j]);
5,724,000✔
190
  }
191
}
1,728,000✔
192

193
// Function to generate randomized Halton sequence samples
194
//
195
// Algorithm adapted from:
196
//      A. B. Owen. A randomized halton algorithm in r. Arxiv, 6 2017.
197
//      URL https://arxiv.org/abs/1706.02808
198
vector<double> rhalton(int dim, uint64_t* seed, int64_t skip = 0)
12,000✔
199
{
200
  if (dim > 10) {
12,000✔
UNCOV
201
    fatal_error("Halton sampling dimension too large");
×
202
  }
203
  int64_t b, res, dig;
204
  double b2r, ans;
205
  const std::array<int64_t, 10> primes = {2, 3, 5, 7, 11, 13, 17, 19, 23, 29};
12,000✔
206
  vector<double> halton(dim, 0.0);
12,000✔
207

208
  vector<int64_t> perm;
12,000✔
209
  for (int D = 0; D < dim; ++D) {
72,000✔
210
    b = primes[D];
60,000✔
211
    perm.resize(b);
60,000✔
212
    b2r = 1.0 / b;
60,000✔
213
    res = skip;
60,000✔
214
    ans = 0.0;
60,000✔
215

216
    while ((1.0 - b2r) < 1.0) {
1,788,000✔
217
      std::iota(perm.begin(), perm.end(), 0);
1,728,000✔
218
      fisher_yates_shuffle(perm, seed);
1,728,000✔
219
      dig = res % b;
1,728,000✔
220
      ans += perm[dig] * b2r;
1,728,000✔
221
      res = (res - dig) / b;
1,728,000✔
222
      b2r /= b;
1,728,000✔
223
    }
224

225
    halton[D] = ans;
60,000✔
226
  }
227

228
  return halton;
24,000✔
229
}
12,000✔
230

231
//==============================================================================
232
// RandomRay implementation
233
//==============================================================================
234

235
// Static Variable Declarations
236
double RandomRay::distance_inactive_;
237
double RandomRay::distance_active_;
238
unique_ptr<Source> RandomRay::ray_source_;
239
RandomRaySourceShape RandomRay::source_shape_ {RandomRaySourceShape::FLAT};
240
RandomRaySampleMethod RandomRay::sample_method_ {RandomRaySampleMethod::PRNG};
241

242
RandomRay::RandomRay()
776,400✔
243
  : angular_flux_(data::mg.num_energy_groups_),
776,400✔
244
    delta_psi_(data::mg.num_energy_groups_),
776,400✔
245
    negroups_(data::mg.num_energy_groups_)
2,329,200✔
246
{
247
  if (source_shape_ == RandomRaySourceShape::LINEAR ||
776,400✔
248
      source_shape_ == RandomRaySourceShape::LINEAR_XY) {
512,400✔
249
    delta_moments_.resize(negroups_);
400,200✔
250
  }
251
}
776,400✔
252

253
RandomRay::RandomRay(uint64_t ray_id, FlatSourceDomain* domain) : RandomRay()
776,400✔
254
{
255
  initialize_ray(ray_id, domain);
776,400✔
256
}
776,400✔
257

258
// Transports ray until termination criteria are met
259
uint64_t RandomRay::transport_history_based_single_ray()
776,400✔
260
{
261
  using namespace openmc;
262
  while (alive()) {
373,969,572✔
263
    event_advance_ray();
373,969,572✔
264
    if (!alive())
373,969,572✔
265
      break;
776,400✔
266
    event_cross_surface();
373,193,172✔
267
  }
268

269
  return n_event();
776,400✔
270
}
271

272
// Transports ray across a single source region
273
void RandomRay::event_advance_ray()
373,969,572✔
274
{
275
  // Find the distance to the nearest boundary
276
  boundary() = distance_to_boundary(*this);
373,969,572✔
277
  double distance = boundary().distance;
373,969,572✔
278

279
  if (distance <= 0.0) {
373,969,572✔
UNCOV
280
    mark_as_lost("Negative transport distance detected for particle " +
×
UNCOV
281
                 std::to_string(id()));
×
UNCOV
282
    return;
×
283
  }
284

285
  if (is_active_) {
373,969,572✔
286
    // If the ray is in the active length, need to check if it has
287
    // reached its maximum termination distance. If so, reduce
288
    // the ray traced length so that the ray does not overrun the
289
    // maximum numerical length (so as to avoid numerical bias).
290
    if (distance_travelled_ + distance >= distance_active_) {
292,628,784✔
291
      distance = distance_active_ - distance_travelled_;
776,400✔
292
      wgt() = 0.0;
776,400✔
293
    }
294

295
    distance_travelled_ += distance;
292,628,784✔
296
    attenuate_flux(distance, true);
292,628,784✔
297
  } else {
298
    // If the ray is still in the dead zone, need to check if it
299
    // has entered the active phase. If so, split into two segments (one
300
    // representing the final part of the dead zone, the other representing the
301
    // first part of the active length) and attenuate each. Otherwise, if the
302
    // full length of the segment is within the dead zone, attenuate as normal.
303
    if (distance_travelled_ + distance >= distance_inactive_) {
81,340,788✔
304
      is_active_ = true;
776,400✔
305
      double distance_dead = distance_inactive_ - distance_travelled_;
776,400✔
306
      attenuate_flux(distance_dead, false);
776,400✔
307

308
      double distance_alive = distance - distance_dead;
776,400✔
309

310
      // Ensure we haven't travelled past the active phase as well
311
      if (distance_alive > distance_active_) {
776,400✔
UNCOV
312
        distance_alive = distance_active_;
×
UNCOV
313
        wgt() = 0.0;
×
314
      }
315

316
      attenuate_flux(distance_alive, true);
776,400✔
317
      distance_travelled_ = distance_alive;
776,400✔
318
    } else {
319
      distance_travelled_ += distance;
80,564,388✔
320
      attenuate_flux(distance, false);
80,564,388✔
321
    }
322
  }
323

324
  // Advance particle
325
  for (int j = 0; j < n_coord(); ++j) {
1,495,235,172✔
326
    coord(j).r += distance * coord(j).u;
1,121,265,600✔
327
  }
328
}
329

330
void RandomRay::attenuate_flux(double distance, bool is_active)
374,745,972✔
331
{
332

333
  switch (source_shape_) {
374,745,972✔
334
  case RandomRaySourceShape::FLAT:
168,937,356✔
335
    if (this->material() == MATERIAL_VOID) {
168,937,356✔
336
      attenuate_flux_flat_source_void(distance, is_active);
9,040,620✔
337
    } else {
338
      attenuate_flux_flat_source(distance, is_active);
159,896,736✔
339
    }
340
    break;
168,937,356✔
341
  case RandomRaySourceShape::LINEAR:
205,808,616✔
342
  case RandomRaySourceShape::LINEAR_XY:
343
    if (this->material() == MATERIAL_VOID) {
205,808,616✔
344
      attenuate_flux_linear_source_void(distance, is_active);
9,040,620✔
345
    } else {
346
      attenuate_flux_linear_source(distance, is_active);
196,767,996✔
347
    }
348
    break;
205,808,616✔
UNCOV
349
  default:
×
UNCOV
350
    fatal_error("Unknown source shape for random ray transport.");
×
351
  }
352
}
374,745,972✔
353

354
// This function forms the inner loop of the random ray transport process.
355
// It is responsible for several tasks. Based on the incoming angular flux
356
// of the ray and the source term in the region, the outgoing angular flux
357
// is computed. The delta psi between the incoming and outgoing fluxes is
358
// contributed to the estimate of the total scalar flux in the source region.
359
// Additionally, the contribution of the ray path to the stochastically
360
// estimated volume is also kept track of. All tasks involving writing
361
// to the data for the source region are done with a lock over the entire
362
// source region.  Locks are used instead of atomics as all energy groups
363
// must be written, such that locking once is typically much more efficient
364
// than use of many atomic operations corresponding to each energy group
365
// individually (at least on CPU). Several other bookkeeping tasks are also
366
// performed when inside the lock.
367
void RandomRay::attenuate_flux_flat_source(double distance, bool is_active)
159,896,736✔
368
{
369
  // The number of geometric intersections is counted for reporting purposes
370
  n_event()++;
159,896,736✔
371

372
  // Determine source region index etc.
373
  int i_cell = lowest_coord().cell;
159,896,736✔
374

375
  // The source region is the spatial region index
376
  int64_t sr = domain_->source_region_offsets_[i_cell] + cell_instance();
159,896,736✔
377

378
  // The source element is the energy-specific region index
379
  int material = this->material();
159,896,736✔
380

381
  // MOC incoming flux attenuation + source contribution/attenuation equation
382
  for (int g = 0; g < negroups_; g++) {
740,561,400✔
383
    float sigma_t = domain_->sigma_t_[material * negroups_ + g];
580,664,664✔
384
    float tau = sigma_t * distance;
580,664,664✔
385
    float exponential = cjosey_exponential(tau); // exponential = 1 - exp(-tau)
580,664,664✔
386
    float new_delta_psi =
387
      (angular_flux_[g] - domain_->source_regions_.source(sr, g)) * exponential;
580,664,664✔
388
    delta_psi_[g] = new_delta_psi;
580,664,664✔
389
    angular_flux_[g] -= new_delta_psi;
580,664,664✔
390
  }
391

392
  // If ray is in the active phase (not in dead zone), make contributions to
393
  // source region bookkeeping
394
  if (is_active) {
159,896,736✔
395

396
    // Aquire lock for source region
397
    domain_->source_regions_.lock(sr).lock();
123,858,444✔
398

399
    // Accumulate delta psi into new estimate of source region flux for
400
    // this iteration
401
    for (int g = 0; g < negroups_; g++) {
542,763,168✔
402
      domain_->source_regions_.scalar_flux_new(sr, g) += delta_psi_[g];
418,904,724✔
403
    }
404

405
    // Accomulate volume (ray distance) into this iteration's estimate
406
    // of the source region's volume
407
    domain_->source_regions_.volume(sr) += distance;
123,858,444✔
408

409
    // Tally valid position inside the source region (e.g., midpoint of
410
    // the ray) if not done already
411
    if (!domain_->source_regions_.position_recorded(sr)) {
123,858,444✔
412
      Position midpoint = r() + u() * (distance / 2.0);
396,100✔
413
      domain_->source_regions_.position(sr) = midpoint;
396,100✔
414
      domain_->source_regions_.position_recorded(sr) = 1;
396,100✔
415
    }
416

417
    // Release lock
418
    domain_->source_regions_.lock(sr).unlock();
123,858,444✔
419
  }
420
}
159,896,736✔
421

422
// Alternative flux attenuation function for true void regions.
423
void RandomRay::attenuate_flux_flat_source_void(double distance, bool is_active)
9,040,620✔
424
{
425
  // The number of geometric intersections is counted for reporting purposes
426
  n_event()++;
9,040,620✔
427

428
  // Determine source region index etc.
429
  int i_cell = lowest_coord().cell;
9,040,620✔
430

431
  // The source region is the spatial region index
432
  int64_t sr = domain_->source_region_offsets_[i_cell] + cell_instance();
9,040,620✔
433

434
  // If ray is in the active phase (not in dead zone), make contributions to
435
  // source region bookkeeping
436
  if (is_active) {
9,040,620✔
437

438
    // Aquire lock for source region
439
    domain_->source_regions_.lock(sr).lock();
7,518,984✔
440

441
    // Accumulate delta psi into new estimate of source region flux for
442
    // this iteration
443
    for (int g = 0; g < negroups_; g++) {
15,037,968✔
444
      domain_->source_regions_.scalar_flux_new(sr, g) +=
7,518,984✔
445
        angular_flux_[g] * distance;
7,518,984✔
446
    }
447

448
    // Accomulate volume (ray distance) into this iteration's estimate
449
    // of the source region's volume
450
    domain_->source_regions_.volume(sr) += distance;
7,518,984✔
451
    domain_->source_regions_.volume_sq(sr) += distance * distance;
7,518,984✔
452

453
    // Tally valid position inside the source region (e.g., midpoint of
454
    // the ray) if not done already
455
    if (!domain_->source_regions_.position_recorded(sr)) {
7,518,984✔
456
      Position midpoint = r() + u() * (distance / 2.0);
17,000✔
457
      domain_->source_regions_.position(sr) = midpoint;
17,000✔
458
      domain_->source_regions_.position_recorded(sr) = 1;
17,000✔
459
    }
460

461
    // Release lock
462
    domain_->source_regions_.lock(sr).unlock();
7,518,984✔
463
  }
464

465
  // Add source to incoming angular flux, assuming void region
466
  for (int g = 0; g < negroups_; g++) {
18,081,240✔
467
    angular_flux_[g] +=
9,040,620✔
468
      domain_->source_regions_.external_source(sr, g) * distance;
9,040,620✔
469
  }
470
}
9,040,620✔
471

472
void RandomRay::attenuate_flux_linear_source(double distance, bool is_active)
196,767,996✔
473
{
474
  // Cast domain to LinearSourceDomain
475
  LinearSourceDomain* domain = dynamic_cast<LinearSourceDomain*>(domain_);
196,767,996✔
476
  if (!domain) {
196,767,996✔
UNCOV
477
    fatal_error("RandomRay::attenuate_flux_linear_source() called with "
×
478
                "non-LinearSourceDomain domain.");
479
  }
480

481
  // The number of geometric intersections is counted for reporting purposes
482
  n_event()++;
196,767,996✔
483

484
  // Determine source region index etc.
485
  int i_cell = lowest_coord().cell;
196,767,996✔
486

487
  // The source region is the spatial region index
488
  int64_t sr = domain_->source_region_offsets_[i_cell] + cell_instance();
196,767,996✔
489

490
  // The source element is the energy-specific region index
491
  int material = this->material();
196,767,996✔
492

493
  Position& centroid = domain_->source_regions_.centroid(sr);
196,767,996✔
494
  Position midpoint = r() + u() * (distance / 2.0);
196,767,996✔
495

496
  // Determine the local position of the midpoint and the ray origin
497
  // relative to the source region's centroid
498
  Position rm_local;
196,767,996✔
499
  Position r0_local;
196,767,996✔
500

501
  // In the first few iterations of the simulation, the source region
502
  // may not yet have had any ray crossings, in which case there will
503
  // be no estimate of its centroid. We detect this by checking if it has
504
  // any accumulated volume. If its volume is zero, just use the midpoint
505
  // of the ray as the region's centroid.
506
  if (domain_->source_regions_.volume_t(sr)) {
196,767,996✔
507
    rm_local = midpoint - centroid;
192,311,652✔
508
    r0_local = r() - centroid;
192,311,652✔
509
  } else {
510
    rm_local = {0.0, 0.0, 0.0};
4,456,344✔
511
    r0_local = -u() * 0.5 * distance;
4,456,344✔
512
  }
513
  double distance_2 = distance * distance;
196,767,996✔
514

515
  // Linear Source MOC incoming flux attenuation + source
516
  // contribution/attenuation equation
517
  for (int g = 0; g < negroups_; g++) {
1,064,666,136✔
518

519
    // Compute tau, the optical thickness of the ray segment
520
    float sigma_t = domain_->sigma_t_[material * negroups_ + g];
867,898,140✔
521
    float tau = sigma_t * distance;
867,898,140✔
522

523
    // If tau is very small, set it to zero to avoid numerical issues.
524
    // The following computations will still work with tau = 0.
525
    if (tau < 1.0e-8f) {
867,898,140✔
526
      tau = 0.0f;
300✔
527
    }
528

529
    // Compute linear source terms, spatial and directional (dir),
530
    // calculated from the source gradients dot product with local centroid
531
    // and direction, respectively.
532
    float spatial_source =
533
      domain_->source_regions_.source(sr, g) +
867,898,140✔
534
      rm_local.dot(domain_->source_regions_.source_gradients(sr, g));
867,898,140✔
535
    float dir_source =
536
      u().dot(domain_->source_regions_.source_gradients(sr, g));
867,898,140✔
537

538
    float gn = exponentialG(tau);
867,898,140✔
539
    float f1 = 1.0f - tau * gn;
867,898,140✔
540
    float f2 = (2.0f * gn - f1) * distance_2;
867,898,140✔
541
    float new_delta_psi = (angular_flux_[g] - spatial_source) * f1 * distance -
867,898,140✔
542
                          0.5 * dir_source * f2;
867,898,140✔
543

544
    float h1 = f1 - gn;
867,898,140✔
545
    float g1 = 0.5f - h1;
867,898,140✔
546
    float g2 = exponentialG2(tau);
867,898,140✔
547
    g1 = g1 * spatial_source;
867,898,140✔
548
    g2 = g2 * dir_source * distance * 0.5f;
867,898,140✔
549
    h1 = h1 * angular_flux_[g];
867,898,140✔
550
    h1 = (g1 + g2 + h1) * distance_2;
867,898,140✔
551
    spatial_source = spatial_source * distance + new_delta_psi;
867,898,140✔
552

553
    // Store contributions for this group into arrays, so that they can
554
    // be accumulated into the source region's estimates inside of the locked
555
    // region.
556
    delta_psi_[g] = new_delta_psi;
867,898,140✔
557
    delta_moments_[g] = r0_local * spatial_source + u() * h1;
867,898,140✔
558

559
    // Update the angular flux for this group
560
    angular_flux_[g] -= new_delta_psi * sigma_t;
867,898,140✔
561

562
    // If 2D mode is enabled, the z-component of the flux moments is forced
563
    // to zero
564
    if (source_shape_ == RandomRaySourceShape::LINEAR_XY) {
867,898,140✔
565
      delta_moments_[g].z = 0.0;
504,226,956✔
566
    }
567
  }
568

569
  // If ray is in the active phase (not in dead zone), make contributions to
570
  // source region bookkeeping
571
  if (is_active) {
196,767,996✔
572
    // Compute an estimate of the spatial moments matrix for the source
573
    // region based on parameters from this ray's crossing
574
    MomentMatrix moment_matrix_estimate;
575
    moment_matrix_estimate.compute_spatial_moments_matrix(
154,508,772✔
576
      rm_local, u(), distance);
154,508,772✔
577

578
    // Aquire lock for source region
579
    domain_->source_regions_.lock(sr).lock();
154,508,772✔
580

581
    // Accumulate deltas into the new estimate of source region flux for this
582
    // iteration
583
    for (int g = 0; g < negroups_; g++) {
812,451,984✔
584
      domain_->source_regions_.scalar_flux_new(sr, g) += delta_psi_[g];
657,943,212✔
585
      domain_->source_regions_.flux_moments_new(sr, g) += delta_moments_[g];
657,943,212✔
586
    }
587

588
    // Accumulate the volume (ray segment distance), centroid, and spatial
589
    // momement estimates into the running totals for the iteration for this
590
    // source region. The centroid and spatial momements estimates are scaled by
591
    // the ray segment length as part of length averaging of the estimates.
592
    domain_->source_regions_.volume(sr) += distance;
154,508,772✔
593
    domain_->source_regions_.centroid_iteration(sr) += midpoint * distance;
154,508,772✔
594
    moment_matrix_estimate *= distance;
154,508,772✔
595
    domain_->source_regions_.mom_matrix(sr) += moment_matrix_estimate;
154,508,772✔
596

597
    // Tally valid position inside the source region (e.g., midpoint of
598
    // the ray) if not done already
599
    if (!domain_->source_regions_.position_recorded(sr)) {
154,508,772✔
600
      domain_->source_regions_.position(sr) = midpoint;
173,468✔
601
      domain_->source_regions_.position_recorded(sr) = 1;
173,468✔
602
    }
603

604
    // Release lock
605
    domain_->source_regions_.lock(sr).unlock();
154,508,772✔
606
  }
607
}
196,767,996✔
608

609
// If traveling through a void region, the source term is either zero
610
// or an external source. As all external sources are currently assumed
611
// to be flat, we don't really need this function and could instead just call
612
// the "attenuate_flux_flat_source_void" function and get the same numerical and
613
// tally results. However, computation of the flux moments in void regions is
614
// nonetheless useful as this information is still used by the plotter when
615
// estimating the flux at specific pixel coordinates. Thus, plots will look
616
// nicer/more accurate if we record flux moments, so this function is useful.
617
void RandomRay::attenuate_flux_linear_source_void(
9,040,620✔
618
  double distance, bool is_active)
619
{
620
  // Cast domain to LinearSourceDomain
621
  LinearSourceDomain* domain = dynamic_cast<LinearSourceDomain*>(domain_);
9,040,620✔
622
  if (!domain) {
9,040,620✔
UNCOV
623
    fatal_error("RandomRay::attenuate_flux_linear_source() called with "
×
624
                "non-LinearSourceDomain domain.");
625
  }
626

627
  // The number of geometric intersections is counted for reporting purposes
628
  n_event()++;
9,040,620✔
629

630
  // Determine source region index etc.
631
  int i_cell = lowest_coord().cell;
9,040,620✔
632

633
  // The source region is the spatial region index
634
  int64_t sr = domain_->source_region_offsets_[i_cell] + cell_instance();
9,040,620✔
635

636
  Position& centroid = domain_->source_regions_.centroid(sr);
9,040,620✔
637
  Position midpoint = r() + u() * (distance / 2.0);
9,040,620✔
638

639
  // Determine the local position of the midpoint and the ray origin
640
  // relative to the source region's centroid
641
  Position rm_local;
9,040,620✔
642
  Position r0_local;
9,040,620✔
643

644
  // In the first few iterations of the simulation, the source region
645
  // may not yet have had any ray crossings, in which case there will
646
  // be no estimate of its centroid. We detect this by checking if it has
647
  // any accumulated volume. If its volume is zero, just use the midpoint
648
  // of the ray as the region's centroid.
649
  if (domain_->source_regions_.volume_t(sr)) {
9,040,620✔
650
    rm_local = midpoint - centroid;
8,819,112✔
651
    r0_local = r() - centroid;
8,819,112✔
652
  } else {
653
    rm_local = {0.0, 0.0, 0.0};
221,508✔
654
    r0_local = -u() * 0.5 * distance;
221,508✔
655
  }
656
  double distance_2 = distance * distance;
9,040,620✔
657

658
  // Compared to linear flux attenuation through solid regions,
659
  // transport through a void region is greatly simplified. Here we
660
  // compute the updated flux moments.
661
  for (int g = 0; g < negroups_; g++) {
18,081,240✔
662
    float spatial_source = domain_->source_regions_.source(sr, g);
9,040,620✔
663
    float new_delta_psi = (angular_flux_[g] - spatial_source) * distance;
9,040,620✔
664
    float h1 = 0.5f;
9,040,620✔
665
    h1 = h1 * angular_flux_[g];
9,040,620✔
666
    h1 = h1 * distance_2;
9,040,620✔
667
    spatial_source = spatial_source * distance + new_delta_psi;
9,040,620✔
668

669
    // Store contributions for this group into arrays, so that they can
670
    // be accumulated into the source region's estimates inside of the locked
671
    // region.
672
    delta_moments_[g] = r0_local * spatial_source + u() * h1;
9,040,620✔
673

674
    // If 2D mode is enabled, the z-component of the flux moments is forced
675
    // to zero
676
    if (source_shape_ == RandomRaySourceShape::LINEAR_XY) {
9,040,620✔
UNCOV
677
      delta_moments_[g].z = 0.0;
×
678
    }
679
  }
680

681
  // If ray is in the active phase (not in dead zone), make contributions to
682
  // source region bookkeeping
683
  if (is_active) {
9,040,620✔
684
    // Compute an estimate of the spatial moments matrix for the source
685
    // region based on parameters from this ray's crossing
686
    MomentMatrix moment_matrix_estimate;
687
    moment_matrix_estimate.compute_spatial_moments_matrix(
7,518,984✔
688
      rm_local, u(), distance);
7,518,984✔
689

690
    // Aquire lock for source region
691
    domain_->source_regions_.lock(sr).lock();
7,518,984✔
692

693
    // Accumulate delta psi into new estimate of source region flux for
694
    // this iteration, and update flux momements
695
    for (int g = 0; g < negroups_; g++) {
15,037,968✔
696
      domain_->source_regions_.scalar_flux_new(sr, g) +=
15,037,968✔
697
        angular_flux_[g] * distance;
7,518,984✔
698
      domain_->source_regions_.flux_moments_new(sr, g) += delta_moments_[g];
7,518,984✔
699
    }
700

701
    // Accumulate the volume (ray segment distance), centroid, and spatial
702
    // momement estimates into the running totals for the iteration for this
703
    // source region. The centroid and spatial momements estimates are scaled by
704
    // the ray segment length as part of length averaging of the estimates.
705
    domain_->source_regions_.volume(sr) += distance;
7,518,984✔
706
    domain_->source_regions_.volume_sq(sr) += distance_2;
7,518,984✔
707
    domain_->source_regions_.centroid_iteration(sr) += midpoint * distance;
7,518,984✔
708
    moment_matrix_estimate *= distance;
7,518,984✔
709
    domain_->source_regions_.mom_matrix(sr) += moment_matrix_estimate;
7,518,984✔
710

711
    // Tally valid position inside the source region (e.g., midpoint of
712
    // the ray) if not done already
713
    if (!domain_->source_regions_.position_recorded(sr)) {
7,518,984✔
714
      domain_->source_regions_.position(sr) = midpoint;
17,000✔
715
      domain_->source_regions_.position_recorded(sr) = 1;
17,000✔
716
    }
717

718
    // Release lock
719
    domain_->source_regions_.lock(sr).unlock();
7,518,984✔
720
  }
721

722
  // Add source to incoming angular flux, assuming void region
723
  for (int g = 0; g < negroups_; g++) {
18,081,240✔
724
    angular_flux_[g] +=
9,040,620✔
725
      domain_->source_regions_.external_source(sr, g) * distance;
9,040,620✔
726
  }
727
}
9,040,620✔
728

729
void RandomRay::initialize_ray(uint64_t ray_id, FlatSourceDomain* domain)
776,400✔
730
{
731
  domain_ = domain;
776,400✔
732

733
  // Reset particle event counter
734
  n_event() = 0;
776,400✔
735

736
  is_active_ = (distance_inactive_ <= 0.0);
776,400✔
737

738
  wgt() = 1.0;
776,400✔
739

740
  // set identifier for particle
741
  id() = simulation::work_index[mpi::rank] + ray_id;
776,400✔
742

743
  // generate source site using sample method
744
  SourceSite site;
776,400✔
745
  switch (sample_method_) {
776,400✔
746
  case RandomRaySampleMethod::PRNG:
764,400✔
747
    site = sample_prng();
764,400✔
748
    break;
764,400✔
749
  case RandomRaySampleMethod::HALTON:
12,000✔
750
    site = sample_halton();
12,000✔
751
    break;
12,000✔
UNCOV
752
  default:
×
UNCOV
753
    fatal_error("Unknown sample method for random ray transport.");
×
754
  }
755

756
  site.E = lower_bound_index(
776,400✔
757
    data::mg.rev_energy_bins_.begin(), data::mg.rev_energy_bins_.end(), site.E);
758
  site.E = negroups_ - site.E - 1.;
776,400✔
759
  this->from_source(&site);
776,400✔
760

761
  // Locate ray
762
  if (lowest_coord().cell == C_NONE) {
776,400✔
763
    if (!exhaustive_find_cell(*this)) {
776,400✔
UNCOV
764
      this->mark_as_lost(
×
UNCOV
765
        "Could not find the cell containing particle " + std::to_string(id()));
×
766
    }
767

768
    // Set birth cell attribute
769
    if (cell_born() == C_NONE)
776,400✔
770
      cell_born() = lowest_coord().cell;
776,400✔
771
  }
772

773
  // Initialize ray's starting angular flux to starting location's isotropic
774
  // source
775
  int i_cell = lowest_coord().cell;
776,400✔
776
  int64_t sr = domain_->source_region_offsets_[i_cell] + cell_instance();
776,400✔
777

778
  for (int g = 0; g < negroups_; g++) {
2,956,800✔
779
    angular_flux_[g] = domain_->source_regions_.source(sr, g);
2,180,400✔
780
  }
781
}
776,400✔
782

783
SourceSite RandomRay::sample_prng()
764,400✔
784
{
785
  // set random number seed
786
  int64_t particle_seed =
787
    (simulation::current_batch - 1) * settings::n_particles + id();
764,400✔
788
  init_particle_seeds(particle_seed, seeds());
764,400✔
789
  stream() = STREAM_TRACKING;
764,400✔
790

791
  // Sample from ray source distribution
792
  SourceSite site {ray_source_->sample(current_seed())};
764,400✔
793

794
  return site;
764,400✔
795
}
796

797
SourceSite RandomRay::sample_halton()
12,000✔
798
{
799
  SourceSite site;
12,000✔
800

801
  // Set random number seed
802
  int64_t batch_seed = (simulation::current_batch - 1) * settings::n_particles;
12,000✔
803
  int64_t skip = id();
12,000✔
804
  init_particle_seeds(batch_seed, seeds());
12,000✔
805
  stream() = STREAM_TRACKING;
12,000✔
806

807
  // Calculate next samples in LDS across 5 dimensions
808
  vector<double> samples = rhalton(5, current_seed(), skip = skip);
12,000✔
809

810
  // Get spatial box of ray_source_
811
  SpatialBox* sb = dynamic_cast<SpatialBox*>(
12,000✔
812
    dynamic_cast<IndependentSource*>(RandomRay::ray_source_.get())->space());
12,000✔
813

814
  // Sample spatial distribution
815
  Position xi {samples[0], samples[1], samples[2]};
12,000✔
816
  // make a small shift in position to avoid geometry floating point issues
817
  Position shift {FP_COINCIDENT, FP_COINCIDENT, FP_COINCIDENT};
12,000✔
818
  site.r = (sb->lower_left() + shift) +
819
           xi * ((sb->upper_right() - shift) - (sb->lower_left() + shift));
12,000✔
820

821
  // Sample Polar cosine and azimuthal angles
822
  double mu = 2.0 * samples[3] - 1.0;
12,000✔
823
  double azi = 2.0 * PI * samples[4];
12,000✔
824
  // Convert to Cartesian coordinates
825
  double c = std::sqrt(1.0 - mu * mu);
12,000✔
826
  site.u.x = mu;
12,000✔
827
  site.u.y = std::cos(azi) * c;
12,000✔
828
  site.u.z = std::sin(azi) * c;
12,000✔
829

830
  return site;
24,000✔
831
}
12,000✔
832

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