https://mooseframework.inl.gov
RayTracingStudy.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "RayTracingStudy.h"
11 
12 // Local Includes
13 #include "AuxRayKernel.h"
15 #include "RayKernel.h"
16 #include "TraceRay.h"
17 #include "TraceRayTools.h"
18 
19 // MOOSE Includes
20 #include "AuxiliarySystem.h"
21 #include "Assembly.h"
22 #include "NonlinearSystemBase.h"
23 
24 // libMesh Includes
25 #include "libmesh/enum_to_string.h"
26 #include "libmesh/mesh_tools.h"
27 #include "libmesh/parallel_sync.h"
28 #include "libmesh/remote_elem.h"
29 
30 using namespace libMesh;
31 
34 {
35  auto params = GeneralUserObject::validParams();
36 
37  // Parameters for the execution of the rays
39 
40  params.addRangeCheckedParam<Real>("ray_distance",
41  std::numeric_limits<Real>::max(),
42  "ray_distance > 0",
43  "The maximum distance all Rays can travel");
44 
45  params.addParam<bool>(
46  "tolerate_failure", false, "Whether or not to tolerate a ray tracing failure");
47 
48  MooseEnum work_buffers("lifo circular", "circular");
49  params.addParam<MooseEnum>("work_buffer_type", work_buffers, "The work buffer type to use");
50 
51  params.addParam<bool>(
52  "ray_kernel_coverage_check", true, "Whether or not to perform coverage checks on RayKernels");
53  params.addParam<bool>("warn_non_planar",
54  true,
55  "Whether or not to produce a warning if any element faces are non-planar.");
56 
57  params.addParam<bool>(
58  "always_cache_traces",
59  false,
60  "Whether or not to cache the Ray traces on every execution, primarily for use in output. "
61  "Warning: this can get expensive very quick with a large number of rays!");
62  params.addParam<bool>("data_on_cache_traces",
63  false,
64  "Whether or not to also cache the Ray's data when caching its traces");
65  params.addParam<bool>("aux_data_on_cache_traces",
66  false,
67  "Whether or not to also cache the Ray's aux data when caching its traces");
68  params.addParam<bool>(
69  "segments_on_cache_traces",
70  true,
71  "Whether or not to cache individual segments when trace caching is enabled. If false, we "
72  "will instead cache a segment for each part of the trace where the direction is the same. "
73  "This minimizes the number of segments requied to represent the Ray's path, but removes the "
74  "ability to show Ray field data on each segment through an element.");
75 
76  params.addParam<bool>("use_internal_sidesets",
77  false,
78  "Whether or not to use internal sidesets for RayBCs in ray tracing");
79 
80  params.addParam<bool>("warn_subdomain_hmax",
81  true,
82  "Whether or not to warn if the approximated hmax (constant on subdomain) "
83  "varies significantly for an element");
84 
85  params.addParam<bool>(
86  "verify_rays",
87  true,
88  "Whether or not to verify the generated Rays. This includes checking their "
89  "starting information and the uniqueness of Rays before and after execution. This is also "
90  "used by derived studies for more specific verification.");
91  params.addParam<bool>("verify_trace_intersections",
92  true,
93  "Whether or not to verify the trace intersections in devel and dbg modes. "
94  "Trace intersections are not verified regardless of this parameter in "
95  "optimized modes (opt, oprof).");
96 
97  params.addParam<bool>("allow_other_flags_with_prekernels",
98  false,
99  "Whether or not to allow the list of execution flags to have PRE_KERNELS "
100  "mixed with other flags. If this parameter is not set then if PRE_KERNELS "
101  "is provided it must be the only execution flag.");
102 
103  ExecFlagEnum & exec_enum = params.set<ExecFlagEnum>("execute_on", true);
105 
106  params.addParamNamesToGroup(
107  "always_cache_traces data_on_cache_traces aux_data_on_cache_traces segments_on_cache_traces",
108  "Trace cache");
109  params.addParamNamesToGroup("warn_non_planar warn_subdomain_hmax", "Tracing Warnings");
110  params.addParamNamesToGroup("ray_kernel_coverage_check verify_rays verify_trace_intersections",
111  "Checks and verifications");
112 
113  // Whether or not each Ray must be registered using the registerRay() API
114  params.addPrivateParam<bool>("_use_ray_registration", true);
115  // Whether or not to bank Rays on completion
116  params.addPrivateParam<bool>("_bank_rays_on_completion", true);
118  params.addPrivateParam<bool>("_ray_dependent_subdomain_setup", true);
119 
120  // Add a point neighbor relationship manager
121  params.addRelationshipManager("ElementPointNeighborLayers",
124  [](const InputParameters &, InputParameters & rm_params)
125  { rm_params.set<unsigned short>("layers") = 1; });
126 
127  return params;
128 }
129 
131  : GeneralUserObject(parameters),
132  _mesh(_fe_problem.mesh()),
133  _comm(_mesh.comm()),
134  _pid(_comm.rank()),
135 
136  _ray_kernel_coverage_check(getParam<bool>("ray_kernel_coverage_check")),
137  _warn_non_planar(getParam<bool>("warn_non_planar")),
138  _use_ray_registration(getParam<bool>("_use_ray_registration")),
139  _use_internal_sidesets(getParam<bool>("use_internal_sidesets")),
140  _tolerate_failure(getParam<bool>("tolerate_failure")),
141  _bank_rays_on_completion(getParam<bool>("_bank_rays_on_completion")),
142  _ray_dependent_subdomain_setup(getParam<bool>("_ray_dependent_subdomain_setup")),
143 
144  _always_cache_traces(getParam<bool>("always_cache_traces")),
145  _data_on_cache_traces(getParam<bool>("data_on_cache_traces")),
146  _aux_data_on_cache_traces(getParam<bool>("aux_data_on_cache_traces")),
147  _segments_on_cache_traces(getParam<bool>("segments_on_cache_traces")),
148  _ray_max_distance(getParam<Real>("ray_distance")),
149  _verify_rays(getParam<bool>("verify_rays")),
150 #ifndef NDEBUG
151  _verify_trace_intersections(getParam<bool>("verify_trace_intersections")),
152 #endif
153 
154  _threaded_elem_side_builders(libMesh::n_threads()),
155 
156  _registered_ray_map(
157  declareRestartableData<std::unordered_map<std::string, RayID>>("registered_ray_map")),
158  _reverse_registered_ray_map(
159  declareRestartableData<std::vector<std::string>>("reverse_registered_ray_map")),
160 
161  _threaded_cached_traces(libMesh::n_threads()),
162 
163  _num_cached(libMesh::n_threads(), 0),
164 
165  _has_non_planar_sides(true),
166  _has_same_level_active_elems(sameLevelActiveElems()),
167 
168  _b_box(MeshTools::create_nodal_bounding_box(_mesh.getMesh())),
169  _domain_max_length(1.01 * (_b_box.max() - _b_box.min()).norm()),
170  _total_volume(computeTotalVolume()),
171 
172  _threaded_cache_ray_kernel(libMesh::n_threads()),
173  _threaded_cache_ray_bc(libMesh::n_threads()),
174  _threaded_ray_object_registration(libMesh::n_threads()),
175  _threaded_current_ray_kernels(libMesh::n_threads()),
176  _threaded_trace_ray(libMesh::n_threads()),
177  _threaded_fe_face(libMesh::n_threads()),
178  _threaded_q_face(libMesh::n_threads()),
179  _threaded_cached_normals(libMesh::n_threads()),
180  _threaded_next_ray_id(libMesh::n_threads()),
181 
182  _parallel_ray_study(std::make_unique<ParallelRayStudy>(*this, _threaded_trace_ray)),
183 
184  _local_trace_ray_results(TraceRay::FAILED_TRACES + 1, 0),
185 
186  _called_initial_setup(false),
187 
188  _elem_index_helper(_mesh.getMesh(), name() + "_elem_index")
189 {
190  // Initialize a tracing object for each thread
191  for (THREAD_ID tid = 0; tid < libMesh::n_threads(); ++tid)
192  {
193  // Initialize a tracing object for each thread
194  _threaded_trace_ray[tid] = std::make_shared<TraceRay>(*this, tid);
195 
196  // Setup the face FEs for normal computation on the fly
197  _threaded_fe_face[tid] = FEBase::build(_mesh.dimension(), FEType(CONSTANT, MONOMIAL));
198  _threaded_q_face[tid] = QBase::build(QGAUSS, _mesh.dimension() - 1, CONSTANT);
199  _threaded_fe_face[tid]->attach_quadrature_rule(_threaded_q_face[tid].get());
200  _threaded_fe_face[tid]->get_normals();
201  }
202 
203  // Evaluating on residual and Jacobian evaluation
205  {
206  if (!getParam<bool>("allow_other_flags_with_prekernels") && _execute_enum.size() > 1)
207  paramError("execute_on",
208  "PRE_KERNELS cannot be mixed with any other execution flag.\nThat is, you cannot "
209  "currently "
210  "mix RayKernels that contribute to the Jacobian/residual with those that do not.");
211 
212  if (_app.useEigenvalue())
213  mooseError("Execution on residual and Jacobian evaluation (execute_on = PRE_KERNELS)\n",
214  "is not supported for an eigenvalue solve.");
215  }
216 
219 
220  // Scale the bounding box for loose checking
222  _loose_b_box.scale(TOLERANCE * TOLERANCE);
223 }
224 
225 void
227 {
228  // Keep track of initialSetup call to avoid registration of various things
229  _called_initial_setup = true;
230 
231  // Sets up a local index for each elem this proc knows about
233 
234  // Check for RayKernel coverage
235  coverageChecks();
236 
237  // Make sure the dependencies exist, if any
239 
240  // Check for traceable element types
242 
243  // Setup for internal sidesets
245 
247 
248  // Setup approximate hmax for each subdomain
250 
251  // Call initial setup on all of the objects
252  for (auto & rto : getRayTracingObjects())
253  rto->initialSetup();
254 
255  // Check for proper exec flags with RayKernels
256  std::vector<RayKernelBase *> ray_kernels;
257  getRayKernels(ray_kernels, 0);
258  for (const auto & rkb : ray_kernels)
259  if (dynamic_cast<RayKernel *>(rkb) && !_execute_enum.isValueSet(EXEC_PRE_KERNELS))
260  mooseError("This study has RayKernel objects that contribute to residuals and Jacobians.",
261  "\nIn this case, the study must use the execute_on = PRE_KERNELS");
262 
263  // Build 1D quadrature rule for along a segment
265  QBase::build(QGAUSS, 1, _fe_problem.getSystemBase(_sys.number()).getMinQuadratureOrder());
266 }
267 
268 void
270 {
271  for (THREAD_ID tid = 0; tid < libMesh::n_threads(); ++tid)
272  mooseAssert(_num_cached[tid] == 0, "Cached residuals/Jacobians not empty");
273 
274  for (auto & rto : getRayTracingObjects())
275  rto->residualSetup();
276 }
277 
278 void
280 {
281  for (THREAD_ID tid = 0; tid < libMesh::n_threads(); ++tid)
282  mooseAssert(_num_cached[tid] == 0, "Cached residuals/Jacobians not empty");
283 
284  for (auto & rto : getRayTracingObjects())
285  rto->jacobianSetup();
286 }
287 
288 void
290 {
291  for (auto & rto : getRayTracingObjects())
292  rto->timestepSetup();
293 }
294 
295 void
297 {
302 
304 
305  for (const auto & trace_ray : _threaded_trace_ray)
306  trace_ray->meshChanged();
307 }
308 
309 void
311 {
312  executeStudy();
313 }
314 
315 void
317 {
318  // Check for coverage of RayKernels on domain
320  {
321  std::vector<RayKernelBase *> ray_kernels;
322  getRayKernels(ray_kernels, 0);
323 
324  std::set<SubdomainID> ray_kernel_blocks;
325  for (const auto & rk : ray_kernels)
326  ray_kernel_blocks.insert(rk->blockIDs().begin(), rk->blockIDs().end());
327 
328  std::set<SubdomainID> missing;
329  std::set_difference(_mesh.meshSubdomains().begin(),
330  _mesh.meshSubdomains().end(),
331  ray_kernel_blocks.begin(),
332  ray_kernel_blocks.end(),
333  std::inserter(missing, missing.begin()));
334 
335  if (!missing.empty() && !ray_kernel_blocks.count(Moose::ANY_BLOCK_ID))
336  {
337  std::ostringstream error;
338  error << "Subdomains { ";
339  std::copy(missing.begin(), missing.end(), std::ostream_iterator<SubdomainID>(error, " "));
340  error << "} do not have RayKernels defined!";
341 
342  mooseError(error.str());
343  }
344  }
345 }
346 
347 void
349 {
350  std::vector<RayTracingObject *> ray_tracing_objects;
351 
352  getRayKernels(ray_tracing_objects, 0);
353  verifyDependenciesExist(ray_tracing_objects);
354 
355  getRayBCs(ray_tracing_objects, 0);
356  verifyDependenciesExist(ray_tracing_objects);
357 }
358 
359 void
360 RayTracingStudy::verifyDependenciesExist(const std::vector<RayTracingObject *> & rtos)
361 {
362  for (const auto & rto : rtos)
363  for (const auto & dep_name : rto->getRequestedItems())
364  {
365  bool found = false;
366  for (const auto & rto_search : rtos)
367  if (rto_search->name() == dep_name)
368  {
369  found = true;
370  break;
371  }
372 
373  if (!found)
374  rto->paramError("depends_on",
375  "The ",
376  rto->parameters().get<std::string>("_moose_base"),
377  " '",
378  dep_name,
379  "' does not exist");
380  }
381 }
382 
383 void
385 {
386 
387  for (const auto & elem : *_mesh.getActiveLocalElementRange())
388  {
389  if (_fe_problem.adaptivity().isOn())
390  {
392  mooseError("Element type ",
393  Utility::enum_to_string(elem->type()),
394  " is not supported in ray tracing with adaptivity");
395  }
396  else if (!TraceRayTools::isTraceableElem(elem))
397  mooseError("Element type ",
398  Utility::enum_to_string(elem->type()),
399  " is not supported in ray tracing");
400  }
401 }
402 
403 void
405 {
406  // TODO: We could probably minimize this to local active elements followed by
407  // boundary point neighbors, but if using distibuted mesh it really shouldn't matter
408  _elem_index_helper.initialize(_mesh.getMesh().active_element_ptr_range());
409 }
410 
411 void
413 {
414  // Even if we have _use_internal_sidesets == false, we will make sure the user didn't add RayBCs
415  // on internal boundaries
416 
417  // Clear the data structures and size the map based on the elements that we know about
418  _internal_sidesets.clear();
419  _internal_sidesets_map.clear();
421 
422  // First, we are going to store all elements with internal sidesets (if any) that have active
423  // RayBCs on them as elem -> vector of (side, vector of boundary ids)
424  for (const auto & bnd_elem : *_mesh.getBoundaryElementRange())
425  {
426  Elem * elem = bnd_elem->_elem;
427  const unsigned int side = bnd_elem->_side;
428  const auto bnd_id = bnd_elem->_bnd_id;
429 
430  // Not internal
431  const Elem * const neighbor = elem->neighbor_ptr(side);
432  if (!neighbor || neighbor == remote_elem)
433  continue;
434 
435  // No RayBCs on this sideset
436  std::vector<RayBoundaryConditionBase *> result;
437  getRayBCs(result, bnd_id, 0);
438  if (result.empty())
439  continue;
440 
441  if (neighbor->subdomain_id() == elem->subdomain_id())
442  mooseError("RayBCs exist on internal sidesets that are not bounded by a different",
443  "\nsubdomain on each side.",
444  "\n\nIn order to use RayBCs on internal sidesets, said sidesets must have",
445  "\na different subdomain on each side.");
446 
447  // Mark that this boundary is an internal sideset with RayBC(s)
448  _internal_sidesets.insert(bnd_id);
449 
450  // Get elem's entry in the internal sidset data structure
451  const auto index = _elem_index_helper.getIndex(elem);
452  auto & entry = _internal_sidesets_map[index];
453 
454  // Initialize this elem's sides if they have not been already
455  if (entry.empty())
456  entry.resize(elem->n_sides(), std::vector<BoundaryID>());
457 
458  // Add the internal boundary to the side entry
459  entry[side].push_back(bnd_id);
460  }
461 
463  mooseError("RayBCs are defined on internal sidesets, but the study is not set to use ",
464  "internal sidesets during tracing.",
465  "\n\nSet the parameter use_internal_sidesets = true to enable this capability.");
466 }
467 
468 void
470 {
471  _has_non_planar_sides = false;
472  bool warned = !_warn_non_planar;
473 
474  // Nothing to do here for 2D or 1D
475  if (_mesh.dimension() != 3)
476  return;
477 
478  // Clear the data structure and size it based on the elements that we know about
479  _non_planar_sides.clear();
481 
482  for (const Elem * elem : _mesh.getMesh().active_element_ptr_range())
483  {
484  const auto index = _elem_index_helper.getIndex(elem);
485  auto & entry = _non_planar_sides[index];
486  entry.resize(elem->n_sides(), 0);
487 
488  for (const auto s : elem->side_index_range())
489  {
490  const auto & side = elemSide(*elem, s);
491  if (side.n_vertices() < 4)
492  continue;
493 
494  if (!side.has_affine_map())
495  {
496  entry[s] = 1;
497  _has_non_planar_sides = true;
498 
499  if (!warned)
500  {
501  mooseWarning("The mesh contains non-planar faces.\n\n",
502  "Ray tracing on non-planar faces is an approximation and may fail.\n\n",
503  "Use at your own risk! You can disable this warning by setting the\n",
504  "parameter 'warn_non_planar' to false.");
505  warned = true;
506  }
507  }
508  }
509  }
510 }
511 
512 void
514 {
515  // Setup map with subdomain keys
516  _subdomain_hmax.clear();
517  for (const auto subdomain_id : _mesh.meshSubdomains())
518  _subdomain_hmax[subdomain_id] = std::numeric_limits<Real>::min();
519 
520  // Set local max for each subdomain
521  for (const auto & elem : *_mesh.getActiveLocalElementRange())
522  {
523  auto & entry = _subdomain_hmax.at(elem->subdomain_id());
524  entry = std::max(entry, elem->hmax());
525  }
526 
527  // Accumulate global max for each subdomain
529 
530  if (getParam<bool>("warn_subdomain_hmax"))
531  {
532  const auto warn_prefix = type() + " '" + name() + "': ";
533  const auto warn_suffix =
534  "\n\nRay tracing uses an approximate element size for each subdomain to scale the\n"
535  "tolerances used in computing ray intersections. This warning suggests that the\n"
536  "approximate element size is not a good approximation. This is likely due to poor\n"
537  "element aspect ratios.\n\n"
538  "This warning is only output for the first element affected.\n"
539  "To disable this warning, set warn_subdomain_hmax = false.\n";
540 
541  for (const auto & elem : *_mesh.getActiveLocalElementRange())
542  {
543  const auto hmin = elem->hmin();
544  const auto hmax = elem->hmax();
545  const auto max_hmax = subdomainHmax(elem->subdomain_id());
546 
547  const auto hmax_rel = hmax / max_hmax;
548  if (hmax_rel < 1.e-2 || hmax_rel > 1.e2)
549  mooseDoOnce(mooseWarning(warn_prefix,
550  "Element hmax varies significantly from subdomain hmax.\n",
551  warn_suffix,
552  "First element affected:\n",
553  Moose::stringify(*elem)););
554 
555  const auto h_rel = max_hmax / hmin;
556  if (h_rel > 1.e2)
557  mooseDoOnce(mooseWarning(warn_prefix,
558  "Element hmin varies significantly from subdomain hmax.\n",
559  warn_suffix,
560  "First element affected:\n",
561  Moose::stringify(*elem)););
562  }
563  }
564 }
565 
566 void
568 {
569  // First, clear the objects associated with each Ray on each thread
570  const auto num_rays = _registered_ray_map.size();
571  for (auto & entry : _threaded_ray_object_registration)
572  {
573  entry.clear();
574  entry.resize(num_rays);
575  }
576 
577  const auto rtos = getRayTracingObjects();
578 
580  {
581  // All of the registered ray names - used when a RayTracingObject did not specify
582  // any Rays so it should be associated with all Rays.
583  std::vector<std::string> all_ray_names;
584  all_ray_names.reserve(_registered_ray_map.size());
585  for (const auto & pair : _registered_ray_map)
586  all_ray_names.push_back(pair.first);
587 
588  for (auto & rto : rtos)
589  {
590  // The Ray names associated with this RayTracingObject
591  const auto & ray_names = rto->parameters().get<std::vector<std::string>>("rays");
592  // The registration for RayTracingObjects for the thread rto is on
593  const auto tid = rto->parameters().get<THREAD_ID>("_tid");
594  auto & registration = _threaded_ray_object_registration[tid];
595 
596  // Register each Ray for this object in the registration
597  for (const auto & ray_name : (ray_names.empty() ? all_ray_names : ray_names))
598  {
599  const auto id = registeredRayID(ray_name, /* graceful = */ true);
600  if (ray_names.size() && id == Ray::INVALID_RAY_ID)
601  rto->paramError(
602  "rays", "Supplied ray '", ray_name, "' is not a registered Ray in ", typeAndName());
603  registration[id].insert(rto);
604  }
605  }
606  }
607  // Not using Ray registration
608  else
609  {
610  for (const auto & rto : rtos)
611  if (rto->parameters().get<std::vector<std::string>>("rays").size())
612  rto->paramError(
613  "rays",
614  "Rays cannot be supplied when the study does not require Ray registration.\n\n",
615  type(),
616  " does not require Ray registration.");
617  }
618 }
619 
620 void
622 {
623  std::set<std::string> vars_to_be_zeroed;
624  std::vector<RayKernelBase *> ray_kernels;
625  getRayKernels(ray_kernels, 0);
626  for (auto & rk : ray_kernels)
627  {
628  AuxRayKernel * aux_rk = dynamic_cast<AuxRayKernel *>(rk);
629  if (aux_rk)
630  vars_to_be_zeroed.insert(aux_rk->variable().name());
631  }
632 
633  std::vector<std::string> vars_to_be_zeroed_vec(vars_to_be_zeroed.begin(),
634  vars_to_be_zeroed.end());
635  _fe_problem.getAuxiliarySystem().zeroVariables(vars_to_be_zeroed_vec);
636 }
637 
638 void
640  const THREAD_ID tid,
641  const RayID ray_id)
642 {
643  mooseAssert(currentlyPropagating(), "Should not call while not propagating");
644 
645  // Call subdomain setup on FE
646  _fe_problem.subdomainSetup(subdomain, tid);
647 
648  std::set<MooseVariableFEBase *> needed_moose_vars;
649  std::unordered_set<unsigned int> needed_mat_props;
650 
651  // Get RayKernels and their dependencies and call subdomain setup
652  getRayKernels(_threaded_current_ray_kernels[tid], subdomain, tid, ray_id);
653  for (auto & rkb : _threaded_current_ray_kernels[tid])
654  {
655  rkb->subdomainSetup();
656 
657  const auto & mv_deps = rkb->getMooseVariableDependencies();
658  needed_moose_vars.insert(mv_deps.begin(), mv_deps.end());
659 
660  const auto & mp_deps = rkb->getMatPropDependencies();
661  needed_mat_props.insert(mp_deps.begin(), mp_deps.end());
662  }
663 
664  // Prepare aux vars
665  for (auto & var : needed_moose_vars)
666  if (var->kind() == Moose::VarKindType::VAR_AUXILIARY)
667  var->prepareAux();
668 
669  _fe_problem.setActiveElementalMooseVariables(needed_moose_vars, tid);
670  _fe_problem.prepareMaterials(needed_mat_props, subdomain, tid);
671 }
672 
673 void
675  const Elem * elem, const Point & start, const Point & end, const Real length, THREAD_ID tid)
676 {
677  mooseAssert(MooseUtils::absoluteFuzzyEqual((start - end).norm(), length), "Invalid length");
678  mooseAssert(currentlyPropagating(), "Should not call while not propagating");
679 
681 
682  // If we have any variables or material properties that are active, we definitely need to reinit
683  bool reinit = _fe_problem.hasActiveElementalMooseVariables(tid) ||
685  // If not, make sure that the RayKernels have not requested a reinit (this could happen when a
686  // RayKernel doesn't have variables or materials but still does an integration and needs qps)
687  if (!reinit)
688  for (const RayKernelBase * rk : currentRayKernels(tid))
689  if (rk->needSegmentReinit())
690  {
691  reinit = true;
692  break;
693  }
694 
695  if (reinit)
696  {
697  _fe_problem.prepare(elem, tid);
698 
699  std::vector<Point> points;
700  std::vector<Real> weights;
701  buildSegmentQuadrature(start, end, length, points, weights);
702  _fe_problem.reinitElemPhys(elem, points, tid);
703  _fe_problem.assembly(tid, _sys.number()).modifyArbitraryWeights(weights);
704 
706  }
707 }
708 
709 void
711  const Point & end,
712  const Real length,
713  std::vector<Point> & points,
714  std::vector<Real> & weights) const
715 {
716  points.resize(_segment_qrule->n_points());
717  weights.resize(_segment_qrule->n_points());
718 
719  const Point diff = end - start;
720  const Point sum = end + start;
721  mooseAssert(MooseUtils::absoluteFuzzyEqual(length, diff.norm()), "Invalid length");
722 
723  // The standard quadrature rule should be on x = [-1, 1]
724  // To scale the points, you...
725  // - Scale to size of the segment in 3D
726  // initial_scaled_qp = x_qp * 0.5 * (end - start) = 0.5 * x_qp * diff
727  // - Shift quadrature midpoint to segment midpoint
728  // final_qp = initial_scaled_qp + 0.5 * (end - start) = initial_scaled_qp + 0.5 * sum
729  // = 0.5 * (x_qp * diff + sum)
730  for (unsigned int qp = 0; qp < _segment_qrule->n_points(); ++qp)
731  {
732  points[qp] = 0.5 * (_segment_qrule->qp(qp)(0) * diff + sum);
733  weights[qp] = 0.5 * _segment_qrule->w(qp) * length;
734  }
735 }
736 
737 void
738 RayTracingStudy::postOnSegment(const THREAD_ID tid, const std::shared_ptr<Ray> & /* ray */)
739 {
740  mooseAssert(currentlyPropagating(), "Should not call while not propagating");
742  mooseAssert(_num_cached[tid] == 0,
743  "Values should only be cached when computing Jacobian/residual");
744 
745  // Fill into cached Jacobian/residuals if necessary
747  {
749 
750  if (++_num_cached[tid] == 20)
751  {
752  Threads::spin_mutex::scoped_lock lock(_spin_mutex);
754  _num_cached[tid] = 0;
755  }
756  }
758  {
760 
761  if (++_num_cached[tid] == 20)
762  {
763  Threads::spin_mutex::scoped_lock lock(_spin_mutex);
765  _num_cached[tid] = 0;
766  }
767  }
768 }
769 
770 void
772 {
773  TIME_SECTION("executeStudy", 2, "Executing Study");
774 
775  mooseAssert(_called_initial_setup, "Initial setup not called");
776 
777  // Reset ray start/complete timers
779  _max_intersections = 0;
781 
782  // Reset physical tracing stats
783  for (auto & val : _local_trace_ray_results)
784  val = 0;
785 
786  // Reset crossing and intersection
793  _ending_distance = 0;
794  _total_distance = 0;
795 
796  // Zero the AuxVariables that our AuxRayKernels contribute to before they accumulate
798 
799  preExecuteStudy();
800  for (auto & rto : getRayTracingObjects())
801  rto->preExecuteStudy();
802 
803  _ray_bank.clear();
804 
805  for (THREAD_ID tid = 0; tid < libMesh::n_threads(); ++tid)
806  {
807  _threaded_trace_ray[tid]->preExecute();
808  _threaded_cached_normals[tid].clear();
809  }
810 
812  _execution_start_time = std::chrono::steady_clock::now();
813 
814  _parallel_ray_study->preExecute();
815 
816  {
817  {
818  auto generation_start_time = std::chrono::steady_clock::now();
819 
820  TIME_SECTION("generateRays", 2, "Generating Rays");
821 
822  generateRays();
823 
824  _generation_time = std::chrono::steady_clock::now() - generation_start_time;
825  }
826 
827  // At this point, nobody is working so this is good time to make sure
828  // Rays are unique across all processors in the working buffer
829  if (verifyRays())
830  {
831  verifyUniqueRays(_parallel_ray_study->workBuffer().begin(),
832  _parallel_ray_study->workBuffer().end(),
833  /* error_suffix = */ "after generateRays()");
834 
835  verifyUniqueRayIDs(_parallel_ray_study->workBuffer().begin(),
836  _parallel_ray_study->workBuffer().end(),
837  /* global = */ true,
838  /* error_suffix = */ "after generateRays()");
839  }
840 
842 
843  {
844  TIME_SECTION("propagateRays", 2, "Propagating Rays");
845 
846  const auto propagation_start_time = std::chrono::steady_clock::now();
847 
848  _parallel_ray_study->execute();
849 
850  _propagation_time = std::chrono::steady_clock::now() - propagation_start_time;
851  }
852  }
853 
854  _execution_time = std::chrono::steady_clock::now() - _execution_start_time;
855 
856  if (verifyRays())
857  {
858  verifyUniqueRays(_parallel_ray_study->workBuffer().begin(),
859  _parallel_ray_study->workBuffer().end(),
860  /* error_suffix = */ "after tracing completed");
861 
862 #ifndef NDEBUG
863  // Outside of debug, _ray_bank always holds all of the Rays that have ended on this processor
864  // We can use this as a global point to check for unique IDs for every Ray that has traced
866  _ray_bank.end(),
867  /* global = */ true,
868  /* error_suffix = */ "after tracing completed");
869 #endif
870  }
871 
872  // Update counters from the threaded trace objects
873  for (const auto & tr : _threaded_trace_ray)
874  for (std::size_t i = 0; i < _local_trace_ray_results.size(); ++i)
875  _local_trace_ray_results[i] += tr->results()[i];
876 
877  // Update local ending counters
884  // ...and communicate the global values
891 
892  // Throw a warning with the number of failed (tolerated) traces
893  if (_tolerate_failure)
894  {
896  _communicator.sum(failures);
897  if (failures)
898  mooseWarning(
899  type(), " '", name(), "': ", failures, " ray tracing failures were tolerated.\n");
900  }
901 
902  // Clear the current RayKernels
903  for (THREAD_ID tid = 0; tid < libMesh::n_threads(); ++tid)
905 
906  // Move the threaded cache trace information into the full cached trace vector
907  // Here, we only clear the cached vectors so that we might not have to
908  // reallocate on future traces
909  std::size_t num_entries = 0;
910  for (THREAD_ID tid = 0; tid < libMesh::n_threads(); ++tid)
911  num_entries += _threaded_cached_traces[tid].size();
912  _cached_traces.clear();
913  _cached_traces.reserve(num_entries);
914  for (THREAD_ID tid = 0; tid < libMesh::n_threads(); ++tid)
915  {
916  for (const auto & entry : _threaded_cached_traces[tid])
917  _cached_traces.emplace_back(std::move(entry));
918  _threaded_cached_traces[tid].clear();
919  }
920 
921  // Add any stragglers that contribute to the Jacobian or residual
922  for (THREAD_ID tid = 0; tid < libMesh::n_threads(); ++tid)
923  if (_num_cached[tid] != 0)
924  {
925  mooseAssert(_fe_problem.currentlyComputingJacobian() ||
927  "Should not have cached values without Jacobian/residual computation");
928 
931  else
933 
934  _num_cached[tid] = 0;
935  }
936 
937  // AuxRayKernels may have modified AuxVariables
940 
941  // Clear FE
942  for (THREAD_ID tid = 0; tid < libMesh::n_threads(); ++tid)
943  {
946  }
947 
949  for (auto & rto : getRayTracingObjects())
950  rto->postExecuteStudy();
951 }
952 
953 void
954 RayTracingStudy::onCompleteRay(const std::shared_ptr<Ray> & ray)
955 {
956  mooseAssert(currentlyPropagating(), "Should only be called during Ray propagation");
957 
958  _ending_processor_crossings += ray->processorCrossings();
960  std::max(_ending_max_processor_crossings, ray->processorCrossings());
961  _ending_intersections += ray->intersections();
962  _ending_max_intersections = std::max(_ending_max_intersections, ray->intersections());
964  std::max(_ending_max_trajectory_changes, ray->trajectoryChanges());
965  _ending_distance += ray->distance();
966 
967 #ifdef NDEBUG
968  // In non-opt modes, we will always bank the Rays for debugging
970 #endif
971  _ray_bank.emplace_back(ray);
972 }
973 
975 RayTracingStudy::registerRayDataInternal(const std::string & name, const bool aux)
976 {
977  Threads::spin_mutex::scoped_lock lock(_spin_mutex);
978 
980  mooseError("Cannot register Ray ", (aux ? "aux " : ""), "data after initialSetup()");
981 
982  auto & map = aux ? _ray_aux_data_map : _ray_data_map;
983  const auto find = map.find(name);
984  if (find != map.end())
985  return find->second;
986 
987  auto & other_map = aux ? _ray_data_map : _ray_aux_data_map;
988  if (other_map.find(name) != other_map.end())
989  mooseError("Cannot register Ray aux data with name ",
990  name,
991  " because Ray ",
992  (aux ? "(non-aux)" : "aux"),
993  " data already exists with said name.");
994 
995  // Add into the name -> index map
996  map.emplace(name, map.size());
997 
998  // Add into the index -> names vector
999  auto & vector = aux ? _ray_aux_data_names : _ray_data_names;
1000  vector.push_back(name);
1001 
1002  return map.size() - 1;
1003 }
1004 
1005 std::vector<RayDataIndex>
1006 RayTracingStudy::registerRayDataInternal(const std::vector<std::string> & names, const bool aux)
1007 {
1008  std::vector<RayDataIndex> indices(names.size());
1009  for (std::size_t i = 0; i < names.size(); ++i)
1010  indices[i] = registerRayDataInternal(names[i], aux);
1011  return indices;
1012 }
1013 
1016  const bool aux,
1017  const bool graceful) const
1018 {
1019  Threads::spin_mutex::scoped_lock lock(_spin_mutex);
1020 
1021  const auto & map = aux ? _ray_aux_data_map : _ray_data_map;
1022  const auto find = map.find(name);
1023  if (find != map.end())
1024  return find->second;
1025 
1026  if (graceful)
1028 
1029  const auto & other_map = aux ? _ray_data_map : _ray_aux_data_map;
1030  if (other_map.find(name) != other_map.end())
1031  mooseError("Ray data with name '",
1032  name,
1033  "' was not found.\n\n",
1034  "However, Ray ",
1035  (aux ? "non-aux" : "aux"),
1036  " data with said name was found.\n",
1037  "Did you mean to use ",
1038  (aux ? "getRayDataIndex()/getRayDataIndices()?"
1039  : "getRayAuxDataIndex()/getRayAuxDataIndices()"),
1040  "?");
1041 
1042  mooseError("Unknown Ray ", (aux ? "aux " : ""), "data with name ", name);
1043 }
1044 
1045 std::vector<RayDataIndex>
1046 RayTracingStudy::getRayDataIndicesInternal(const std::vector<std::string> & names,
1047  const bool aux,
1048  const bool graceful) const
1049 {
1050  std::vector<RayDataIndex> indices(names.size());
1051  for (std::size_t i = 0; i < names.size(); ++i)
1052  indices[i] = getRayDataIndexInternal(names[i], aux, graceful);
1053  return indices;
1054 }
1055 
1056 const std::string &
1057 RayTracingStudy::getRayDataNameInternal(const RayDataIndex index, const bool aux) const
1058 {
1059  Threads::spin_mutex::scoped_lock lock(_spin_mutex);
1060 
1061  if ((aux ? rayAuxDataSize() : rayDataSize()) < index)
1062  mooseError("Unknown Ray ", aux ? "aux " : "", "data with index ", index);
1063  return aux ? _ray_aux_data_names[index] : _ray_data_names[index];
1064 }
1065 
1068 {
1069  return registerRayDataInternal(name, /* aux = */ false);
1070 }
1071 
1072 std::vector<RayDataIndex>
1073 RayTracingStudy::registerRayData(const std::vector<std::string> & names)
1074 {
1075  return registerRayDataInternal(names, /* aux = */ false);
1076 }
1077 
1079 RayTracingStudy::getRayDataIndex(const std::string & name, const bool graceful /* = false */) const
1080 {
1081  return getRayDataIndexInternal(name, /* aux = */ false, graceful);
1082 }
1083 
1084 std::vector<RayDataIndex>
1085 RayTracingStudy::getRayDataIndices(const std::vector<std::string> & names,
1086  const bool graceful /* = false */) const
1087 {
1088  return getRayDataIndicesInternal(names, /* aux = */ false, graceful);
1089 }
1090 
1091 const std::string &
1093 {
1094  return getRayDataNameInternal(index, /* aux = */ false);
1095 }
1096 
1099 {
1100  return registerRayDataInternal(name, /* aux = */ true);
1101 }
1102 
1103 std::vector<RayDataIndex>
1104 RayTracingStudy::registerRayAuxData(const std::vector<std::string> & names)
1105 {
1106  return registerRayDataInternal(names, /* aux = */ true);
1107 }
1108 
1111  const bool graceful /* = false */) const
1112 {
1113  return getRayDataIndexInternal(name, /* aux = */ true, graceful);
1114 }
1115 
1116 std::vector<RayDataIndex>
1117 RayTracingStudy::getRayAuxDataIndices(const std::vector<std::string> & names,
1118  const bool graceful /* = false */) const
1119 {
1120  return getRayDataIndicesInternal(names, /* aux = */ true, graceful);
1121 }
1122 
1123 const std::string &
1125 {
1126  return getRayDataNameInternal(index, /* aux = */ true);
1127 }
1128 
1129 bool
1131 {
1132  std::vector<RayKernelBase *> result;
1133  getRayKernels(result, tid);
1134  return result.size();
1135 }
1136 
1137 void
1138 RayTracingStudy::getRayKernels(std::vector<RayKernelBase *> & result, SubdomainID id, THREAD_ID tid)
1139 {
1140  // If the cache doesn't have any attributes yet, it means that we haven't set
1141  // the conditions yet. We do this so that it can be generated on the fly on first use.
1142  if (!_threaded_cache_ray_kernel[tid].numAttribs())
1143  {
1144  if (!_called_initial_setup)
1145  mooseError("Should not call getRayKernels() before initialSetup()");
1146 
1147  auto query = _fe_problem.theWarehouse()
1148  .query()
1149  .condition<AttribRayTracingStudy>(this)
1150  .condition<AttribSystem>("RayKernel")
1151  .condition<AttribThread>(tid);
1152  _threaded_cache_ray_kernel[tid] = query.clone();
1153  }
1154 
1155  _threaded_cache_ray_kernel[tid].queryInto(result, id);
1156 }
1157 
1158 void
1159 RayTracingStudy::getRayKernels(std::vector<RayKernelBase *> & result,
1160  SubdomainID id,
1161  THREAD_ID tid,
1162  RayID ray_id)
1163 {
1164  // No Ray registration: no need to sift through objects
1165  if (!_use_ray_registration)
1166  {
1167  getRayKernels(result, id, tid);
1168  }
1169  // Has Ray registration: only pick the objects associated with ray_id
1170  else
1171  {
1172  // Get all of the kernels on this block
1173  std::vector<RayKernelBase *> rkbs;
1174  getRayKernels(rkbs, id, tid);
1175 
1176  // The RayTracingObjects associated with this ray
1177  const auto & ray_id_rtos = _threaded_ray_object_registration[tid][ray_id];
1178 
1179  // The result is the union of all of the kernels and the objects associated with this Ray
1180  result.clear();
1181  for (auto rkb : rkbs)
1182  if (ray_id_rtos.count(rkb))
1183  result.push_back(rkb);
1184  }
1185 }
1186 
1187 void
1188 RayTracingStudy::getRayBCs(std::vector<RayBoundaryConditionBase *> & result,
1189  BoundaryID id,
1190  THREAD_ID tid)
1191 {
1192  // If the cache doesn't have any attributes yet, it means that we haven't set
1193  // the conditions yet. We do this so that it can be generated on the fly on first use.
1194  if (!_threaded_cache_ray_bc[tid].numAttribs())
1195  {
1196  if (!_called_initial_setup)
1197  mooseError("Should not call getRayBCs() before initialSetup()");
1198 
1199  auto query = _fe_problem.theWarehouse()
1200  .query()
1201  .condition<AttribRayTracingStudy>(this)
1202  .condition<AttribSystem>("RayBoundaryCondition")
1203  .condition<AttribThread>(tid);
1204  _threaded_cache_ray_bc[tid] = query.clone();
1205  }
1206 
1207  _threaded_cache_ray_bc[tid].queryInto(result, std::make_tuple(id, false));
1208 }
1209 
1210 void
1211 RayTracingStudy::getRayBCs(std::vector<RayBoundaryConditionBase *> & result,
1212  const std::vector<TraceRayBndElement> & bnd_elems,
1213  THREAD_ID tid,
1214  RayID ray_id)
1215 {
1216  // No Ray registration: no need to sift through objects
1217  if (!_use_ray_registration)
1218  {
1219  if (bnd_elems.size() == 1)
1220  getRayBCs(result, bnd_elems[0].bnd_id, tid);
1221  else
1222  {
1223  std::vector<BoundaryID> bnd_ids(bnd_elems.size());
1224  for (MooseIndex(bnd_elems.size()) i = 0; i < bnd_elems.size(); ++i)
1225  bnd_ids[i] = bnd_elems[i].bnd_id;
1226  getRayBCs(result, bnd_ids, tid);
1227  }
1228  }
1229  // Has Ray registration: only pick the objects associated with ray_id
1230  else
1231  {
1232  // Get all of the RayBCs on these boundaries
1233  std::vector<RayBoundaryConditionBase *> rbcs;
1234  if (bnd_elems.size() == 1)
1235  getRayBCs(rbcs, bnd_elems[0].bnd_id, tid);
1236  else
1237  {
1238  std::vector<BoundaryID> bnd_ids(bnd_elems.size());
1239  for (MooseIndex(bnd_elems.size()) i = 0; i < bnd_elems.size(); ++i)
1240  bnd_ids[i] = bnd_elems[i].bnd_id;
1241  getRayBCs(rbcs, bnd_ids, tid);
1242  }
1243 
1244  // The RayTracingObjects associated with this ray
1245  mooseAssert(ray_id < _threaded_ray_object_registration[tid].size(), "Not in registration");
1246  const auto & ray_id_rtos = _threaded_ray_object_registration[tid][ray_id];
1247 
1248  // The result is the union of all of the kernels and the objects associated with this Ray
1249  result.clear();
1250  for (auto rbc : rbcs)
1251  if (ray_id_rtos.count(rbc))
1252  result.push_back(rbc);
1253  }
1254 }
1255 
1256 std::vector<RayTracingObject *>
1258 {
1259  std::vector<RayTracingObject *> result;
1260  _fe_problem.theWarehouse().query().condition<AttribRayTracingStudy>(this).queryInto(result);
1261  return result;
1262 }
1263 
1264 const std::vector<std::shared_ptr<Ray>> &
1266 {
1268  mooseError("The Ray bank is not available because the private parameter "
1269  "'_bank_rays_on_completion' is set to false.");
1271  mooseError("Cannot get the Ray bank during generation or propagation.");
1272 
1273  return _ray_bank;
1274 }
1275 
1276 std::shared_ptr<Ray>
1278 {
1279  // This is only a linear search - can be improved on with a map in the future
1280  // if this is used on a larger scale
1281  std::shared_ptr<Ray> ray;
1282  for (const std::shared_ptr<Ray> & possible_ray : rayBank())
1283  if (possible_ray->id() == ray_id)
1284  {
1285  ray = possible_ray;
1286  break;
1287  }
1288 
1289  // Make sure one and only one processor has the Ray
1290  unsigned int have_ray = ray ? 1 : 0;
1291  _communicator.sum(have_ray);
1292  if (have_ray == 0)
1293  mooseError("Could not find a Ray with the ID ", ray_id, " in the Ray banks.");
1294 
1295  // This should never happen... but let's make sure
1296  mooseAssert(have_ray == 1, "Multiple rays with the same ID were found in the Ray banks");
1297 
1298  return ray;
1299 }
1300 
1301 RayData
1303  const RayDataIndex index,
1304  const bool aux) const
1305 {
1306  // Will be a nullptr shared_ptr if this processor doesn't own the Ray
1307  const std::shared_ptr<Ray> ray = getBankedRay(ray_id);
1308 
1309  Real value = ray ? (aux ? ray->auxData(index) : ray->data(index)) : 0;
1311  return value;
1312 }
1313 
1314 RayData
1315 RayTracingStudy::getBankedRayData(const RayID ray_id, const RayDataIndex index) const
1316 {
1317  return getBankedRayDataInternal(ray_id, index, /* aux = */ false);
1318 }
1319 
1320 RayData
1322 {
1323  return getBankedRayDataInternal(ray_id, index, /* aux = */ true);
1324 }
1325 
1326 RayID
1328 {
1329  libmesh_parallel_only(comm());
1330 
1331  Threads::spin_mutex::scoped_lock lock(_spin_mutex);
1332 
1333  if (!_use_ray_registration)
1334  mooseError("Cannot use registerRay() with Ray registration disabled");
1335 
1336  // This is parallel only for now. We could likely stagger the ID building like we do with
1337  // the unique IDs, but it would require a sync point which isn't there right now
1338  libmesh_parallel_only(comm());
1339 
1340  const auto & it = _registered_ray_map.find(name);
1341  if (it != _registered_ray_map.end())
1342  return it->second;
1343 
1344  const auto id = _reverse_registered_ray_map.size();
1345  _registered_ray_map.emplace(name, id);
1346  _reverse_registered_ray_map.push_back(name);
1347  return id;
1348 }
1349 
1350 RayID
1351 RayTracingStudy::registeredRayID(const std::string & name, const bool graceful /* = false */) const
1352 {
1353  Threads::spin_mutex::scoped_lock lock(_spin_mutex);
1354 
1355  if (!_use_ray_registration)
1356  mooseError("Should not use registeredRayID() with Ray registration disabled");
1357 
1358  const auto search = _registered_ray_map.find(name);
1359  if (search != _registered_ray_map.end())
1360  return search->second;
1361 
1362  if (graceful)
1363  return Ray::INVALID_RAY_ID;
1364 
1365  mooseError("Attempted to obtain ID of registered Ray ",
1366  name,
1367  ", but a Ray with said name is not registered.");
1368 }
1369 
1370 const std::string &
1372 {
1373  Threads::spin_mutex::scoped_lock lock(_spin_mutex);
1374 
1375  if (!_use_ray_registration)
1376  mooseError("Should not use registeredRayName() with Ray registration disabled");
1377 
1378  if (_reverse_registered_ray_map.size() > ray_id)
1379  return _reverse_registered_ray_map[ray_id];
1380 
1381  mooseError("Attempted to obtain name of registered Ray with ID ",
1382  ray_id,
1383  ", but a Ray with said ID is not registered.");
1384 }
1385 
1386 Real
1388 {
1389  Real volume = 0;
1390  for (const auto & elem : *_mesh.getActiveLocalElementRange())
1391  volume += elem->volume();
1393  return volume;
1394 }
1395 
1396 const std::vector<std::vector<BoundaryID>> &
1397 RayTracingStudy::getInternalSidesets(const Elem * elem) const
1398 {
1399  mooseAssert(_use_internal_sidesets, "Not using internal sidesets");
1400  mooseAssert(hasInternalSidesets(), "Processor does not have internal sidesets");
1401  mooseAssert(_internal_sidesets_map.size() > _elem_index_helper.getIndex(elem),
1402  "Internal sideset map not initialized");
1403 
1404  const auto index = _elem_index_helper.getIndex(elem);
1405  return _internal_sidesets_map[index];
1406 }
1407 
1408 TraceData &
1409 RayTracingStudy::initThreadedCachedTrace(const std::shared_ptr<Ray> & ray, THREAD_ID tid)
1410 {
1411  mooseAssert(shouldCacheTrace(ray), "Not caching trace");
1412  mooseAssert(currentlyPropagating(), "Should only use while tracing");
1413 
1414  _threaded_cached_traces[tid].emplace_back(ray);
1415  return _threaded_cached_traces[tid].back();
1416 }
1417 
1418 void
1419 RayTracingStudy::verifyUniqueRayIDs(const std::vector<std::shared_ptr<Ray>>::const_iterator begin,
1420  const std::vector<std::shared_ptr<Ray>>::const_iterator end,
1421  const bool global,
1422  const std::string & error_suffix) const
1423 {
1424  // Determine the unique set of Ray IDs on this processor,
1425  // and if not locally unique throw an error. Once we build this set,
1426  // we will send it to rank 0 to verify globally
1427  std::set<RayID> local_rays;
1428  for (const std::shared_ptr<Ray> & ray : as_range(begin, end))
1429  {
1430  mooseAssert(ray, "Null ray");
1431 
1432  // Try to insert into the set; the second entry in the pair
1433  // will be false if it was not inserted
1434  if (!local_rays.insert(ray->id()).second)
1435  {
1436  for (const std::shared_ptr<Ray> & other_ray : as_range(begin, end))
1437  if (ray.get() != other_ray.get() && ray->id() == other_ray->id())
1438  mooseError("Multiple Rays exist with ID ",
1439  ray->id(),
1440  " on processor ",
1441  _pid,
1442  " ",
1443  error_suffix,
1444  "\n\nOffending Ray information:\n\n",
1445  ray->getInfo(),
1446  "\n",
1447  other_ray->getInfo());
1448  }
1449  }
1450 
1451  // Send IDs from all procs to rank 0 and verify on rank 0
1452  if (global)
1453  {
1454  // Package our local IDs and send to rank 0
1455  std::map<processor_id_type, std::vector<RayID>> send_ids;
1456  if (local_rays.size())
1457  send_ids.emplace(std::piecewise_construct,
1458  std::forward_as_tuple(0),
1459  std::forward_as_tuple(local_rays.begin(), local_rays.end()));
1460  local_rays.clear();
1461 
1462  // Mapping on rank 0 from ID -> processor ID
1463  std::map<RayID, processor_id_type> global_map;
1464 
1465  // Verify another processor's IDs against the global map on rank 0
1466  const auto check_ids =
1467  [this, &global_map, &error_suffix](processor_id_type pid, const std::vector<RayID> & ids)
1468  {
1469  for (const RayID id : ids)
1470  {
1471  const auto emplace_pair = global_map.emplace(id, pid);
1472 
1473  // Means that this ID already exists in the map
1474  if (!emplace_pair.second)
1475  mooseError("Ray with ID ",
1476  id,
1477  " exists on ranks ",
1478  emplace_pair.first->second,
1479  " and ",
1480  pid,
1481  "\n",
1482  error_suffix);
1483  }
1484  };
1485 
1486  Parallel::push_parallel_vector_data(_communicator, send_ids, check_ids);
1487  }
1488 }
1489 
1490 void
1491 RayTracingStudy::verifyUniqueRays(const std::vector<std::shared_ptr<Ray>>::const_iterator begin,
1492  const std::vector<std::shared_ptr<Ray>>::const_iterator end,
1493  const std::string & error_suffix)
1494 {
1495  std::set<const Ray *> rays;
1496  for (const std::shared_ptr<Ray> & ray : as_range(begin, end))
1497  if (!rays.insert(ray.get()).second) // false if not inserted into rays
1498  mooseError("Multiple shared_ptrs were found that point to the same Ray ",
1499  error_suffix,
1500  "\n\nOffending Ray:\n",
1501  ray->getInfo());
1502 }
1503 
1504 void
1505 RayTracingStudy::moveRayToBuffer(std::shared_ptr<Ray> & ray)
1506 {
1507  mooseAssert(currentlyGenerating(), "Can only use while generating");
1508  mooseAssert(ray, "Null ray");
1509  mooseAssert(ray->shouldContinue(), "Ray is not continuing");
1510 
1511  _parallel_ray_study->moveWorkToBuffer(ray, /* tid = */ 0);
1512 }
1513 
1514 void
1515 RayTracingStudy::moveRaysToBuffer(std::vector<std::shared_ptr<Ray>> & rays)
1516 {
1517  mooseAssert(currentlyGenerating(), "Can only use while generating");
1518 #ifndef NDEBUG
1519  for (const std::shared_ptr<Ray> & ray : rays)
1520  {
1521  mooseAssert(ray, "Null ray");
1522  mooseAssert(ray->shouldContinue(), "Ray is not continuing");
1523  }
1524 #endif
1525 
1526  _parallel_ray_study->moveWorkToBuffer(rays, /* tid = */ 0);
1527 }
1528 
1529 void
1531  const THREAD_ID tid,
1532  const AcquireMoveDuringTraceKey &)
1533 {
1534  mooseAssert(ray, "Null ray");
1535  mooseAssert(currentlyPropagating(), "Can only use while tracing");
1536 
1537  _parallel_ray_study->moveWorkToBuffer(ray, tid);
1538 }
1539 
1540 void
1541 RayTracingStudy::reserveRayBuffer(const std::size_t size)
1542 {
1543  if (!currentlyGenerating())
1544  mooseError("Can only reserve in Ray buffer during generateRays()");
1545 
1546  _parallel_ray_study->reserveBuffer(size);
1547 }
1548 
1549 const Point &
1550 RayTracingStudy::getSideNormal(const Elem * elem, unsigned short side, const THREAD_ID tid)
1551 {
1552  std::unordered_map<std::pair<const Elem *, unsigned short>, Point> & cache =
1554 
1555  // See if we've already cached this side normal
1556  const auto elem_side_pair = std::make_pair(elem, side);
1557  const auto search = cache.find(elem_side_pair);
1558 
1559  // Haven't cached this side normal: compute it and then cache it
1560  if (search == cache.end())
1561  {
1562  _threaded_fe_face[tid]->reinit(elem, side);
1563  const auto & normal = _threaded_fe_face[tid]->get_normals()[0];
1564  cache.emplace(elem_side_pair, normal);
1565  return normal;
1566  }
1567 
1568  // Have cached this side normal: simply return it
1569  return search->second;
1570 }
1571 
1572 bool
1574 {
1575  unsigned int min_level = std::numeric_limits<unsigned int>::max();
1576  unsigned int max_level = std::numeric_limits<unsigned int>::min();
1577 
1578  for (const auto & elem : *_mesh.getActiveLocalElementRange())
1579  {
1580  const auto level = elem->level();
1581  min_level = std::min(level, min_level);
1582  max_level = std::max(level, max_level);
1583  }
1584 
1585  _communicator.min(min_level);
1586  _communicator.max(max_level);
1587 
1588  return min_level == max_level;
1589 }
1590 
1591 Real
1593 {
1594  const auto find = _subdomain_hmax.find(subdomain_id);
1595  if (find == _subdomain_hmax.end())
1596  mooseError("Subdomain ", subdomain_id, " not found in subdomain hmax map");
1597  return find->second;
1598 }
1599 
1600 bool
1602 {
1603  Real bbox_volume = 1;
1604  for (unsigned int d = 0; d < _mesh.dimension(); ++d)
1605  bbox_volume *= std::abs(_b_box.max()(d) - _b_box.min()(d));
1606 
1607  return MooseUtils::absoluteFuzzyEqual(bbox_volume, totalVolume(), TOLERANCE);
1608 }
1609 
1610 void
1612 {
1613  libmesh_parallel_only(comm());
1614 
1615  Threads::spin_mutex::scoped_lock lock(_spin_mutex);
1616 
1617  mooseAssert(!currentlyGenerating() && !currentlyPropagating(),
1618  "Cannot be reset during generation or propagation");
1619 
1620  for (THREAD_ID tid = 0; tid < libMesh::n_threads(); ++tid)
1622 }
1623 
1624 RayID
1626 {
1627  // Get the current ID to return
1628  const auto id = _threaded_next_ray_id[tid];
1629 
1630  // Advance so that the next call has the correct ID
1632 
1633  return id;
1634 }
1635 
1636 void
1638 {
1639  libmesh_parallel_only(comm());
1640 
1641  Threads::spin_mutex::scoped_lock lock(_spin_mutex);
1642 
1643  mooseAssert(!currentlyGenerating() && !currentlyPropagating(),
1644  "Cannot be reset during generation or propagation");
1645 
1647 }
1648 
1649 RayID
1651 {
1652  return _replicated_next_ray_id++;
1653 }
1654 
1655 bool
1657  const unsigned short side,
1658  const Point & direction,
1659  const THREAD_ID tid)
1660 {
1661  const auto & normal = getSideNormal(elem, side, tid);
1662  const auto dot = normal * direction;
1663  return dot < TraceRayTools::TRACE_TOLERANCE;
1664 }
1665 
1666 std::shared_ptr<Ray>
1668 {
1669  mooseAssert(currentlyGenerating(), "Can only use during generateRays()");
1670 
1671  return _parallel_ray_study->acquireParallelData(
1672  /* tid = */ 0,
1673  this,
1674  generateUniqueRayID(/* tid = */ 0),
1675  rayDataSize(),
1676  rayAuxDataSize(),
1677  /* reset = */ true,
1679 }
1680 
1681 std::shared_ptr<Ray>
1683 {
1684  mooseAssert(currentlyGenerating(), "Can only use during generateRays()");
1685 
1686  return _parallel_ray_study->acquireParallelData(/* tid = */ 0,
1687  this,
1688  generateUniqueRayID(/* tid = */ 0),
1689  /* data_size = */ 0,
1690  /* aux_data_size = */ 0,
1691  /* reset = */ true,
1693 }
1694 
1695 std::shared_ptr<Ray>
1697 {
1698  mooseAssert(currentlyGenerating(), "Can only use during generateRays()");
1699  libmesh_parallel_only(comm());
1700 
1701  return _parallel_ray_study->acquireParallelData(
1702  /* tid = */ 0,
1703  this,
1705  rayDataSize(),
1706  rayAuxDataSize(),
1707  /* reset = */ true,
1709 }
1710 
1711 std::shared_ptr<Ray>
1713 {
1714  mooseAssert(currentlyGenerating(), "Can only use during generateRays()");
1715 
1716  // Either register a Ray or get an already registered Ray id
1717  const RayID id = registerRay(name);
1718 
1719  // Acquire a Ray with the properly sized data initialized to zero
1720  return _parallel_ray_study->acquireParallelData(
1721  /* tid = */ 0,
1722  this,
1723  id,
1724  rayDataSize(),
1725  rayAuxDataSize(),
1726  /* reset = */ true,
1728 }
1729 
1730 std::shared_ptr<Ray>
1732 {
1733  mooseAssert(currentlyGenerating(), "Can only use during generateRays()");
1734  return _parallel_ray_study->acquireParallelData(
1735  /* tid = */ 0, &ray, Ray::ConstructRayKey());
1736 }
1737 
1738 std::shared_ptr<Ray>
1740 {
1741  mooseAssert(currentlyPropagating(), "Can only use during propagation");
1742  return _parallel_ray_study->acquireParallelData(tid,
1743  this,
1744  generateUniqueRayID(tid),
1745  rayDataSize(),
1746  rayAuxDataSize(),
1747  /* reset = */ true,
1749 }
bool hasRayKernels(const THREAD_ID tid)
Whether or not there are currently any active RayKernel objects.
bool currentlyGenerating() const
Whether or not the study is generating.
void getRayKernels(std::vector< RayKernelBase *> &result, SubdomainID id, THREAD_ID tid)
Fills the active RayKernels associated with this study and a block into result.
void internalSidesetSetup()
Does the setup for internal sidesets.
RayDataIndex registerRayData(const std::string &name)
Register a value to be filled in the data on a Ray with a given name.
RayData getBankedRayDataInternal(const RayID ray_id, const RayDataIndex index, const bool aux) const
Internal method for getting the value (replicated across all processors) in a Ray&#39;s data or aux data ...
virtual void postOnSegment(const THREAD_ID tid, const std::shared_ptr< Ray > &ray)
Called at the end of a Ray segment.
const bool _tolerate_failure
Whether or not to tolerate a Ray Tracing failure.
Real totalVolume() const
Get the current total volume of the domain.
void traceableMeshChecks()
Check for if all of the element types in the mesh are supported by ray tracing.
virtual void onCompleteRay(const std::shared_ptr< Ray > &ray)
Entry point for acting on a ray when it is completed (shouldContinue() == false)
libMesh::ConstElemRange * getActiveLocalElementRange()
const std::string & getRayDataNameInternal(const RayDataIndex index, const bool aux) const
Internal method for getting the name of Ray data or Ray aux data.
const std::vector< std::shared_ptr< Ray > > & rayBank() const
Get the Ray bank.
static const RayDataIndex INVALID_RAY_DATA_INDEX
Invalid index into a Ray&#39;s data.
Definition: Ray.h:200
std::vector< std::vector< std::set< const RayTracingObject * > > > _threaded_ray_object_registration
Threaded storage for all of the RayTracingObjects associated with a single Ray.
std::vector< std::string > & _reverse_registered_ray_map
Map from registered Ray ID to name.
unsigned long long int _total_intersections
Total number of Ray/element intersections.
RayID registerRay(const std::string &name)
Registers a Ray with a given name.
unsigned long long int _ending_intersections
Total number of Ray/element intersections for Rays that finished on this processor.
unsigned int _max_trajectory_changes
Max number of trajectory changes for a single Ray.
unsigned long int RayID
Type for a Ray&#39;s ID.
Definition: Ray.h:43
void verifyUniqueRayIDs(const std::vector< std::shared_ptr< Ray >>::const_iterator begin, const std::vector< std::shared_ptr< Ray >>::const_iterator end, const bool global, const std::string &error_suffix) const
Verifies that the Rays in the given range have unique Ray IDs.
unsigned int n_threads()
bool _has_same_level_active_elems
Whether or not the mesh has active elements of the same level.
auto norm() const -> decltype(std::norm(Real()))
bool sameLevelActiveElems() const
Determine whether or not the mesh currently has active elements that are all the same level...
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
virtual void initialSetup() override
void moveRayToBufferDuringTrace(std::shared_ptr< Ray > &ray, const THREAD_ID tid, const AcquireMoveDuringTraceKey &)
INTERNAL method for moving a Ray into the buffer during tracing.
RayDataIndex registerRayAuxData(const std::string &name)
Register a value to be filled in the aux data on a Ray with a given name.
Class that is used as a parameter to the public constructors/reset methods.
Definition: Ray.h:87
MooseMesh & _mesh
The Mesh.
T & getMesh(MooseMesh &mesh)
function to cast mesh
Definition: SCM.h:35
Attribute for the RayTracingStudy a RayTracingObject is associated with.
static InputParameters validParams()
virtual void prepare(const Elem *elem, const THREAD_ID tid) override
NumericVector< Number > & solution()
std::vector< std::vector< TraceData > > _threaded_cached_traces
The threaded storage for cached traces.
const bool _use_internal_sidesets
Whether or not to use the internal sidesets in ray tracing.
Data structure that stores information for output of a partial trace of a Ray on a processor...
Definition: TraceData.h:42
QGAUSS
bool isTraceableElem(const Elem *elem)
bool isAdaptivityTraceableElem(const Elem *elem)
std::shared_ptr< Ray > acquireRay()
User APIs for constructing Rays within the RayTracingStudy.
std::size_t rayAuxDataSize() const
The registered size of values in the Ray aux data.
std::unordered_map< SubdomainID, Real > _subdomain_hmax
The cached hmax for all elements in a subdomain.
libMesh::dof_id_type maxIndex() const
Gets the maximum index generated using this object.
libMesh::BoundingBox _loose_b_box
Loose nodal bounding box for the domain.
T & set(const std::string &name, bool quiet_mode=false)
const std::set< BoundaryID > & getInternalSidesets() const
Gets the internal sidesets (that have RayBCs) within the local domain.
std::vector< std::vector< unsigned short > > _non_planar_sides
Non planar side data, which is for quick checking if an elem side is non-planar We use unsigned short...
RayData getBankedRayData(const RayID ray_id, const RayDataIndex index) const
Gets the data value for a banked ray with a given ID.
Base object for the RayKernel syntax.
Definition: RayKernelBase.h:27
if(subdm)
MeshBase & mesh
void initialize(const libMesh::SimpleRange< libMesh::MeshBase::element_iterator > elems)
Initializes the indices in a contiguous manner for the given element range.
Real _ending_distance
Total distance traveled by Rays that end on this processor.
virtual void reinitSegment(const Elem *elem, const Point &start, const Point &end, const Real length, THREAD_ID tid)
Reinitialize objects for a Ray segment for ray tracing.
const std::string & name() const override
unsigned int size() const
RayDataIndex registerRayDataInternal(const std::string &name, const bool aux)
Internal method for registering Ray data or Ray aux data with a name.
const std::vector< RayKernelBase * > & currentRayKernels(THREAD_ID tid) const
Gets the current RayKernels for a thread, which are set in segmentSubdomainSetup() ...
const Parallel::Communicator & comm() const
Key that is used for restricting access to moveRayToBufferDuringTrace() and acquireRayDuringTrace().
void addAvailableFlags(const ExecFlagType &flag, Args... flags)
std::shared_ptr< Ray > acquireReplicatedRay()
Acquire a Ray from the pool of Rays within generateRays() in a replicated fashion.
std::shared_ptr< Ray > acquireRayDuringTrace(const THREAD_ID tid, const AcquireMoveDuringTraceKey &)
INTERNAL methods for acquiring a Ray during a trace in RayKernels and RayBCs.
const Parallel::Communicator & _communicator
std::shared_ptr< Ray > acquireRegisteredRay(const std::string &name)
Acquires a Ray with a given name within generateRays().
void clearActiveMaterialProperties(const THREAD_ID tid)
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
RayID _replicated_next_ray_id
Storage for the next available replicated RayID, obtained via generateReplicatedRayID() ...
std::set< BoundaryID > _internal_sidesets
The BoundaryIDs on the local mesh that have internal RayBCs.
virtual RayID generateUniqueRayID(const THREAD_ID tid)
Generates a unique RayID to be used for a Ray.
bool isOn()
std::vector< std::vector< RayKernelBase * > > _threaded_current_ray_kernels
The current RayKernel objects for each thread.
virtual const std::string & name() const
void mooseWarning(Args &&... args) const
auto max(const L &left, const R &right)
std::vector< TraceData > _cached_traces
Storage for the cached traces.
RayDataIndex getRayAuxDataIndex(const std::string &name, const bool graceful=false) const
Gets the index associated with a registered value in the Ray aux data.
std::vector< std::shared_ptr< TraceRay > > _threaded_trace_ray
The TraceRay objects for each thread (they do the physical tracing)
virtual Assembly & assembly(const THREAD_ID tid, const unsigned int sys_num) override
Real subdomainHmax(const SubdomainID subdomain_id) const
Get the cached hmax for all elements in a subdomain.
const Real TRACE_TOLERANCE
The standard tolerance to use in tracing.
Definition: TraceRayTools.h:48
std::vector< RayID > _threaded_next_ray_id
Storage for the next available unique RayID, obtained via generateUniqueRayID()
bool isRectangularDomain() const
Whether or not the domain is rectangular (if it is prefectly encompassed by its bounding box) ...
std::vector< std::string > _ray_aux_data_names
The names for each Ray aux data entry.
virtual void setActiveElementalMooseVariables(const std::set< MooseVariableFEBase * > &moose_vars, const THREAD_ID tid) override
Threads::spin_mutex _spin_mutex
Spin mutex object for locks.
uint8_t processor_id_type
std::shared_ptr< Ray > getBankedRay(const RayID ray_id) const
Gets the Ray with the ID ray_id from the Ray bank.
unsigned long long int _total_processor_crossings
Total number of processor crossings.
processor_id_type n_processors() const
CONSTANT
void nonPlanarSideSetup()
Sets up the caching of whether or not each element side is non-planar, which is stored in _non_planar...
bool verifyRays() const
Whether or not to verify if Rays have valid information before being traced.
unsigned int RayDataIndex
Type for the index into the data and aux data on a Ray.
Definition: Ray.h:51
void coverageChecks()
Perform coverage checks (coverage of RayMaterials and RayKernels, if enabled)
bool & useEigenvalue()
void resetReplicatedRayIDs()
Resets the generation of unique replicated RayIDs accessed via generateReplicatedRayID().
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
const bool _bank_rays_on_completion
Whether or not to bank rays on completion.
MeshBase & getMesh()
void min(const T &r, T &o, Request &req) const
unsigned int _ending_max_intersections
Max number of intersections for Rays that finished on this processor.
sideset clear()
TheWarehouse & theWarehouse() const
virtual void postExecuteStudy()
Entry point after study execution.
MooseVariableFE< Real > & variable()
Gets the variable this AuxRayKernel contributes to.
Definition: AuxRayKernel.h:38
unsigned int _max_processor_crossings
Max number of processor crossings for all Rays.
const std::unique_ptr< ParallelRayStudy > _parallel_ray_study
The study that used is to actually execute (trace) the Rays.
virtual unsigned int dimension() const
std::vector< RayTracingObject * > getRayTracingObjects()
Gets all of the currently active RayTracingObjects.
virtual void residualSetup() override
const Point & min() const
const bool & currentlyComputingResidual() const
boundary_id_type BoundaryID
const std::string name
Definition: Setup.h:20
std::vector< TheWarehouse::QueryCache< AttribSubdomains > > _threaded_cache_ray_kernel
Threaded cached subdomain query for RayKernelBase objects pertaining to this study.
virtual void execute() override
Executes the study (generates and propagates Rays)
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
unsigned int _max_intersections
Max number of intersections for a single Ray.
unsigned int _ending_max_trajectory_changes
Max number of trajectory changes for Rays that finished on this processor.
unsigned int _ending_max_processor_crossings
Max number of total processor crossings for Rays that finished on this processor. ...
const std::string & getRayDataName(const RayDataIndex index) const
Gets the name associated with a registered value in the Ray data.
void reserveRayBuffer(const std::size_t size)
Reserve size entires in the Ray buffer.
void zeroAuxVariables()
Zero the AuxVariables that the registered AuxRayKernels contribute to.
std::vector< std::unique_ptr< libMesh::FEBase > > _threaded_fe_face
Face FE used for computing face normals for each thread.
const std::string & type() const
void dependencyChecks()
Perform checks to see if the listed dependencies in the RayTracingObjects exist.
static InputParameters validParams()
const std::string & getRayAuxDataName(const RayDataIndex index) const
Gets the name associated with a registered value in the Ray aux data.
virtual void cacheResidual(const THREAD_ID tid) override
Basic datastructure for a ray that will traverse the mesh.
Definition: Ray.h:56
std::chrono::steady_clock::time_point _execution_start_time
Timing.
std::vector< TheWarehouse::QueryCache< AttribBoundaries > > _threaded_cache_ray_bc
Threaded cached boundary query for RayBC objects pertaining to this study.
SystemBase & _sys
virtual void jacobianSetup() override
std::string typeAndName() const
void paramError(const std::string &param, Args... args) const
bool isValueSet(const std::string &value) const
unsigned int number() const
auto norm(const T &a) -> decltype(std::abs(a))
std::string stringify(const T &t)
virtual void preExecuteStudy()
Entry point before study execution.
std::unordered_map< std::string, RayDataIndex > _ray_aux_data_map
The map from Ray aux data names to index.
AuxiliarySystem & getAuxiliarySystem()
Real volume(const MeshBase &mesh, unsigned int dim=libMesh::invalid_uint)
const bool _warn_non_planar
Whether not to warn if non-planar faces are found.
libMesh::BoundingBox _b_box
Nodal bounding box for the domain.
virtual void meshChanged() override
libMesh::dof_id_type getIndex(const libMesh::Elem *elem) const
Get the index associated with the element elem.
virtual void setCurrentSubdomainID(const Elem *elem, const THREAD_ID tid) override
virtual void generateRays()=0
Subclasses should override this to determine how to generate Rays.
RayID generateReplicatedRayID()
Generates a Ray ID that is replicated across all processors.
virtual void close()=0
MONOMIAL
const ExecFlagEnum & _execute_enum
RayDataIndex getRayDataIndexInternal(const std::string &name, const bool aux, const bool graceful) const
Internal method for getting the index of Ray data or Ray aux data.
void verifyDependenciesExist(const std::vector< RayTracingObject *> &rtos)
Verifies that the dependencies exist for a set of RayTracingObjects.
unsigned long long int _ending_processor_crossings
Total number of processor crossings for Rays that finished on this processor.
query_obj query
ElemIndexHelper _elem_index_helper
Helper for defining a local contiguous index for each element.
RayTracingStudy(const InputParameters &parameters)
virtual unsigned int n_sides() const=0
std::vector< std::unordered_map< std::pair< const Elem *, unsigned short >, Point > > _threaded_cached_normals
Threaded cache for side normals that have been computed already during tracing.
TraceData & initThreadedCachedTrace(const std::shared_ptr< Ray > &ray, THREAD_ID tid)
Initialize a Ray in the threaded cached trace map to be filled with segments.
virtual const SystemBase & getSystemBase(const unsigned int sys_num) const
const Elem * neighbor_ptr(unsigned int i) const
virtual const Point & getSideNormal(const Elem *elem, const unsigned short side, const THREAD_ID tid)
Get the outward normal for a given element side.
const SubdomainID ANY_BLOCK_ID
unsigned int level() const
std::vector< RayDataIndex > getRayAuxDataIndices(const std::vector< std::string > &names, const bool graceful=false) const
Gets the indices associated with registered values in the Ray aux data.
std::unordered_map< std::string, RayDataIndex > _ray_data_map
The map from Ray data names to index.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< RayDataIndex > getRayDataIndices(const std::vector< std::string > &names, const bool graceful=false) const
Gets the indices associated with registered values in the Ray data.
bool currentlyPropagating() const
Whether or not the study is propagating (tracing Rays)
subdomain_id_type subdomain_id() const
bool _has_non_planar_sides
Whether or not the local mesh has elements with non-planar sides.
virtual void subdomainSetup(SubdomainID subdomain, const THREAD_ID tid)
const libMesh::Elem & elemSide(const libMesh::Elem &elem, const unsigned int s, const THREAD_ID tid=0)
Get an element&#39;s side pointer without excessive memory allocation.
void max(const T &r, T &o, Request &req) const
void reinitMaterials(SubdomainID blk_id, const THREAD_ID tid, bool swap_stateful=true)
FEProblemBase & _fe_problem
void getRayBCs(std::vector< RayBoundaryConditionBase *> &result, BoundaryID id, THREAD_ID tid)
Fills the active RayBCs associated with this study and a boundary into result.
std::vector< std::size_t > _num_cached
Number of currently cached objects for Jacobian/residual for each thread.
std::chrono::steady_clock::duration _propagation_time
Query query()
virtual void segmentSubdomainSetup(const SubdomainID subdomain, const THREAD_ID tid, const RayID ray_id)
Setup for on subdomain change or subdomain AND ray change during ray tracing.
RayData getBankedRayAuxData(const RayID ray_id, const RayDataIndex index) const
Gets the data value for a banked ray with a given ID.
bool _called_initial_setup
Whether or not we&#39;ve called initial setup - used to stop from late registration.
void verifyUniqueRays(const std::vector< std::shared_ptr< Ray >>::const_iterator begin, const std::vector< std::shared_ptr< Ray >>::const_iterator end, const std::string &error_suffix)
Verifies that the Rays in the given range are unique.
std::vector< RayDataIndex > getRayDataIndicesInternal(const std::vector< std::string > &names, const bool aux, const bool graceful) const
Internal method for getting the indicies of Ray data or Ray aux data.
void mooseError(Args &&... args) const
Real _total_distance
Total distance traveled by all Rays.
const processor_id_type _pid
The rank of this processor (this actually takes time to lookup - so just do it once) ...
std::unordered_map< std::string, RayID > & _registered_ray_map
Map from registered Ray name to ID.
std::vector< std::string > _ray_data_names
The names for each Ray data entry.
float RayData
Type for a Ray&#39;s data.
Definition: Ray.h:46
const std::string & registeredRayName(const RayID ray_id) const
Gets the name of a registered ray.
void moveRayToBuffer(std::shared_ptr< Ray > &ray)
Moves a ray to the buffer to be traced during generateRays().
const ExecFlagType EXEC_PRE_KERNELS
RayID registeredRayID(const std::string &name, const bool graceful=false) const
Gets the ID of a registered ray.
virtual void buildSegmentQuadrature(const Point &start, const Point &end, const Real length, std::vector< Point > &points, std::vector< Real > &weights) const
Builds quadrature points for a given segment using the _segment_qrule.
libMesh::StoredRange< MooseMesh::const_bnd_elem_iterator, const BndElement *> * getBoundaryElementRange()
const bool _use_ray_registration
Whether or not to use Ray registration.
virtual bool hasActiveElementalMooseVariables(const THREAD_ID tid) const
const bool _ray_kernel_coverage_check
Whether or not to perform coverage checks on RayKernels.
void localElemIndexSetup()
Sets up the _elem_index_helper, which is used for obtaining a contiguous index for all elements that ...
const bool & currentlyComputingJacobian() const
std::vector< std::shared_ptr< Ray > > _ray_bank
Cumulative Ray bank - stored only when _bank_rays_on_completion.
libMesh::BoundingBox create_nodal_bounding_box(const MeshBase &mesh)
const Point & max() const
std::vector< unsigned long long int > _local_trace_ray_results
Cumulative results on this processor from the threaded TraceRay objects.
std::chrono::steady_clock::duration _generation_time
void registeredRaySetup()
Sets up the maps from Ray to associated RayTracingObjects if _use_ray_registration.
void scale(const Real factor)
void moveRaysToBuffer(std::vector< std::shared_ptr< Ray >> &rays)
Moves rays to the buffer to be traced during generateRays().
void executeStudy()
Method for executing the study so that it can be called out of the standard UO execute() ...
void prepareMaterials(const std::unordered_set< unsigned int > &consumer_needed_mat_props, const SubdomainID blk_id, const THREAD_ID tid)
Adaptivity & adaptivity()
bool hasActiveMaterialProperties(const THREAD_ID tid) const
virtual void cacheJacobian(const THREAD_ID tid) override
auto min(const L &left, const R &right)
virtual void reinitElemPhys(const Elem *elem, const std::vector< Point > &phys_points_in_elem, const THREAD_ID tid) override
std::shared_ptr< Ray > acquireUnsizedRay()
Acquire a Ray from the pool of Rays within generateRays(), without resizing the data (sizes the data ...
RayDataIndex getRayDataIndex(const std::string &name, const bool graceful=false) const
Gets the index associated with a registered value in the Ray data.
void resetUniqueRayIDs()
Resets the generation of unique RayIDs via generateUniqueRayID() to the beginning of the range...
virtual void addCachedResidual(const THREAD_ID tid) override
std::chrono::steady_clock::duration _execution_time
virtual void zeroVariables(std::vector< std::string > &vars_to_be_zeroed)
virtual void clearActiveElementalMooseVariables(const THREAD_ID tid) override
std::vector< std::unique_ptr< libMesh::QBase > > _threaded_q_face
Face quadrature used for computing face normals for each thread.
std::shared_ptr< Ray > acquireCopiedRay(const Ray &ray)
Acquires a Ray that that is copied from another Ray within generateRays().
Traces Rays through the mesh on a single processor.
Definition: TraceRay.h:46
std::vector< std::vector< std::vector< BoundaryID > > > _internal_sidesets_map
Internal sideset data, if internal sidesets exist (indexed with getLocalElemIndex()) ...
void subdomainHMaxSetup()
Caches the hmax for all elements in each subdomain.
bool sideIsIncoming(const Elem *const elem, const unsigned short side, const Point &direction, const THREAD_ID tid)
Whether or not side is incoming on element elem in direction direction.
static const RayID INVALID_RAY_ID
Invalid Ray ID.
Definition: Ray.h:202
std::unique_ptr< libMesh::QBase > _segment_qrule
Quadrature rule for laying points across a 1D ray segment.
bool hasInternalSidesets() const
Whether or not the local mesh has internal sidesets that have RayBCs on them.
std::size_t rayDataSize() const
The registered size of values in the Ray data.
unsigned int THREAD_ID
virtual void timestepSetup() override
const std::set< SubdomainID > & meshSubdomains() const
virtual void addCachedJacobian(const THREAD_ID tid) override
virtual bool shouldCacheTrace(const std::shared_ptr< Ray > &) const
Virtual that allows for selection in if a Ray should be cached or not (only used when _cache_traces)...
const RemoteElem * remote_elem
Real computeTotalVolume()
Helper function for computing the total domain volume.