https://mooseframework.inl.gov
LotsOfRaysRayStudy.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 "LotsOfRaysRayStudy.h"
11 
12 // Local includes
14 
15 // libMesh includes
16 #include "libmesh/parallel_algebra.h"
17 
18 registerMooseObject("RayTracingTestApp", LotsOfRaysRayStudy);
19 
22 {
24 
25  params.addParam<bool>(
26  "vertex_to_vertex",
27  true,
28  "Enable generation of rays from vertices on boundary sides in the direction of "
29  "the vertices on the other side of the same elem");
30  params.addParam<bool>("centroid_to_vertex",
31  true,
32  "Enable generation of rays from centroids on boundary sides in the "
33  "direction of all other nodes in the same elem");
34  params.addParam<bool>("centroid_to_centroid",
35  true,
36  "Enable generation of rays from centroids on boundary sides to centroids "
37  "of all boundary elements");
38  params.addParam<bool>("edge_to_edge",
39  false,
40  "Enable generation of rays from centroids of boundary edges in the "
41  "direction of all other edge centroids in the same elem");
42  params.addParam<bool>("side_aq",
43  false,
44  "Enable generation of rays from boundary side centroids in the direction "
45  "of angular quadrature");
46  params.addParam<bool>("centroid_aq",
47  false,
48  "Enable generation of rays from boundary element centroids in the "
49  "direction of angular quadrature.");
50 
51  params.addRangeCheckedParam<unsigned int>(
52  "polar_quad_order",
53  5,
54  "polar_quad_order % 2",
55  "Order of the polar quadrature angular quadrature generation. Polar angle is between ray and "
56  "the normal direction. Must be odd.");
57  params.addRangeCheckedParam<unsigned int>(
58  "azimuthal_quad_order",
59  2,
60  "azimuthal_quad_order > 0",
61  "Order of the azimuthal quadrature per quadrant for angular quadrature generation. The "
62  "azimuthal angle is measured in a plane perpendicular to the normal. Not needed in 2D.");
63 
64  params.addParam<bool>(
65  "compute_expected_distance",
66  false,
67  "Whether or not to compute the expected distance for all of the Rays generated");
68 
69  params.addParam<bool>(
70  "use_unsized_rays", false, "Whether or not to size Ray data on defineRays()");
71 
72  params.addParam<bool>("set_incoming_side", true, "Whether or not to set the incoming side");
73 
74  // When we define rays, they are already localized on the processors that start them
75  // and the starting elem/incoming side is set
76  params.set<bool>("_claim_after_define_rays") = false;
77  // Can't have replicated Rays if the Rays are on their starting processors only
78  params.set<bool>("_define_rays_replicated") = false;
79  // Don't need to use registration here. For the purposes of testing and this study,
80  // we'll really only ever access the Rays by ID (if at all)
81  params.set<bool>("_use_ray_registration") = false;
82 
83  return params;
84 }
85 
87  : RepeatableRayStudyBase(parameters),
88  _vertex_to_vertex(getParam<bool>("vertex_to_vertex")),
89  _centroid_to_vertex(getParam<bool>("centroid_to_vertex")),
90  _centroid_to_centroid(getParam<bool>("centroid_to_centroid")),
91  _edge_to_edge(getParam<bool>("edge_to_edge")),
92  _side_aq(getParam<bool>("side_aq")),
93  _centroid_aq(getParam<bool>("centroid_aq")),
94  _compute_expected_distance(getParam<bool>("compute_expected_distance")),
95  _polar_quad_order(getParam<unsigned int>("polar_quad_order")),
96  _azimuthal_quad_order(getParam<unsigned int>("azimuthal_quad_order")),
97  _use_unsized_rays(getParam<bool>("use_unsized_rays")),
98  _set_incoming_side(getParam<bool>("set_incoming_side")),
99  _expected_distance(declareRestartableData<Real>("expected_distance"))
100 {
101  // Computing the expected distance only works with rectanglar domains because we have to compute
102  // the intersection with the Ray and the domain bounding box (which only represents the domain if
103  // it is rectangular)
105  paramError("compute_expected_distance", "Not supported when the domain is non-rectangular.");
106 
107  const bool polar_set = parameters.isParamSetByUser("polar_quad_order");
108  const bool azimuthal_set = parameters.isParamSetByUser("azimuthal_quad_order");
109  if (_mesh.dimension() == 1)
110  {
111  if (polar_set)
112  paramError("polar_quad_order", "Not needed for 1D");
113  if (azimuthal_set)
114  paramError("azimuthal_quad_order", "Not needed for 1D");
115  }
116  if (!_side_aq && !_centroid_aq)
117  {
118  if (polar_set)
119  paramError("polar_quad_order", "Not needed without side_aq or centroid_aq enabled");
120  if (azimuthal_set)
121  paramError("azimuthal_quad_order", "Not needed without side_aq or centroid_aq enabled");
122  }
123 }
124 
125 void
127 {
129  {
131  std::make_unique<BoundingBoxIntersectionHelper>(boundingBox(), _mesh.dimension());
132  _expected_distance = 0;
133  }
134 
135  // Gather the boundary elem centroids
136  std::vector<Point> boundary_centroids;
138  {
139  for (const auto & elem : *_mesh.getActiveLocalElementRange())
140  if (elem->on_boundary())
141  boundary_centroids.push_back(elem->vertex_average());
142  _comm.allgather(boundary_centroids);
143  }
144 
145  const std::unique_ptr<RayTracingAngularQuadrature> half_aq =
146  _side_aq ? std::make_unique<RayTracingAngularQuadrature>(
148  : nullptr;
149 
150  for (const auto & bnd_elem : *_mesh.getBoundaryElementRange())
151  {
152  const auto side = bnd_elem->_side;
153  const Elem * elem = bnd_elem->_elem;
154 
155  // Starting from local elems only
156  if (elem->processor_id() != _pid)
157  continue;
158 
159  // Centroid on this boundary face
160  const auto & side_centroid = elemSide(*elem, side).vertex_average();
161 
162  // Vertices on this boundary side -> in direction of all other vertices not on said boundary
163  // side on elem
164  if (_vertex_to_vertex)
165  for (const auto n : elem->nodes_on_side(side))
166  if (elem->is_vertex(n))
167  for (unsigned int v_to = 0; v_to < elem->n_vertices(); ++v_to)
168  if (!elem->is_node_on_side(v_to, side))
169  defineRay(elem, side, elem->point(n), elem->point(v_to), false);
170 
171  // Centroid on this boundary side -> in direction of all other vertices not on said boundary
172  // side on elem
174  for (unsigned int v_to = 0; v_to < elem->n_vertices(); ++v_to)
175  if (!elem->is_node_on_side(v_to, side))
176  defineRay(elem, side, side_centroid, elem->point(v_to), false);
177 
178  // Centroid on this boundary side -> centroids of all boundary elements
180  for (const auto & other_centroid : boundary_centroids)
181  defineRay(elem, RayTracingCommon::invalid_side, side_centroid, other_centroid, true);
182 
183  // Centroids of edges on this boundary side -> in direction of all other edge centroids
184  if (_edge_to_edge)
185  for (const auto edge : elem->edge_index_range())
186  if (elem->is_edge_on_side(edge, side))
187  {
188  const Point edge_centroid = elem->build_edge_ptr(edge)->vertex_average();
189  for (const auto edge_to : elem->edge_index_range())
190  if (edge != edge_to)
191  defineRay(
192  elem, side, edge_centroid, elem->build_edge_ptr(edge_to)->vertex_average(), false);
193  }
194 
195  if (_side_aq)
196  {
197  const auto inward_normal = -1.0 * getSideNormal(elem, side, /* tid = */ 0);
198  half_aq->rotate(inward_normal);
199 
200  for (std::size_t l = 0; l < half_aq->numDirections(); ++l)
201  defineRay(elem, side, side_centroid, side_centroid + half_aq->getDirection(l), false);
202  }
203  }
204 
205  if (_centroid_aq)
206  {
209 
210  for (const Elem * elem : *_mesh.getActiveLocalElementRange())
211  if (elem->processor_id() == _pid && elem->on_boundary())
212  {
213  const auto centroid = elem->vertex_average();
214  for (std::size_t l = 0; l < full_aq.numDirections(); ++l)
215  defineRay(elem,
217  centroid,
218  centroid + full_aq.getDirection(l),
219  false);
220  }
221  }
222 
225 
226  // Insertion point for other derived test objects to modify the Rays
227  modifyRays();
228 }
229 
230 void
231 LotsOfRaysRayStudy::defineRay(const Elem * starting_elem,
232  const unsigned short incoming_side,
233  const Point & p1,
234  const Point & p2,
235  const bool ends_within_mesh)
236 {
237  std::shared_ptr<Ray> ray = _use_unsized_rays ? acquireUnsizedRay() : acquireRay();
238 
239  ray->setStart(
240  p1, starting_elem, _set_incoming_side ? incoming_side : RayTracingCommon::invalid_side);
241 
242  if (ends_within_mesh)
243  {
244  ray->setStartingEndPoint(p2);
245 
247  _expected_distance += (p2 - p1).norm();
248  }
249  else
250  {
251  ray->setStartingDirection(p2 - p1);
252 
254  {
255  const Point end = _bbox_intersection_helper->intersection(p1, ray->direction());
257  mooseError("Expected distance end intersection not found\n\n", ray->getInfo());
258 
259  _expected_distance += (end - p1).norm();
260  }
261  }
262 
263  _rays.emplace_back(std::move(ray));
264 }
265 
266 void
268 {
269 }
static const unsigned short invalid_side
Identifier for an invalid side index.
const bool _centroid_to_centroid
libMesh::ConstElemRange * getActiveLocalElementRange()
const bool _vertex_to_vertex
const unsigned int _polar_quad_order
Polar angular quadrature order for aq tests.
MooseMesh & _mesh
The Mesh.
const bool _centroid_to_vertex
std::shared_ptr< Ray > acquireRay()
User APIs for constructing Rays within the RayTracingStudy.
const unsigned int _azimuthal_quad_order
Azimuthal angular quadrature order for aq tests.
virtual void defineRays() override
Entry point for the user to create Rays.
const Parallel::Communicator & _communicator
static InputParameters validParams()
const BoundingBox & boundingBox() const
Get the nodal bounding box for the domain.
Real & _expected_distance
The expected total distance Rays should travel.
bool isRectangularDomain() const
Whether or not the domain is rectangular (if it is prefectly encompassed by its bounding box) ...
std::unique_ptr< BoundingBoxIntersectionHelper > _bbox_intersection_helper
Helper for computing the end point for Rays that don&#39;t end within mesh.
void defineRay(const Elem *starting_elem, const unsigned short incoming_side, const Point &p1, const Point &p2, const bool ends_within_mesh)
virtual unsigned int dimension() const
const bool _compute_expected_distance
Whether or not to compute the expected distance for generated rays.
void paramError(const std::string &param, Args... args) const
auto norm(const T &a) -> decltype(std::abs(a))
const bool _use_unsized_rays
virtual void modifyRays()
Insertion point for after _rays is defined for other derived test studies to modify the Rays...
bool isParamSetByUser(const std::string &name) 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.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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 mooseError(Args &&... args) const
std::vector< std::shared_ptr< Ray > > & _rays
Vector of Rays that the user will fill into in defineRays() (restartable)
const processor_id_type _pid
The rank of this processor (this actually takes time to lookup - so just do it once) ...
const InputParameters & parameters() const
const Parallel::Communicator & _comm
The Communicator.
libMesh::StoredRange< MooseMesh::const_bnd_elem_iterator, const BndElement *> * getBoundaryElementRange()
static const libMesh::Point invalid_point(invalid_distance, invalid_distance, invalid_distance)
Identifier for an invalid point.
const bool _set_incoming_side
LotsOfRaysRayStudy(const InputParameters &parameters)
A RayTracingStudy used for generating a lot of rays for testing purposes.
static InputParameters validParams()
A RayTracingStudy that generates and traces Rays repeatedly that a user defines only once...
std::shared_ptr< Ray > acquireUnsizedRay()
Acquire a Ray from the pool of Rays within generateRays(), without resizing the data (sizes the data ...
void ErrorVector unsigned int
registerMooseObject("RayTracingTestApp", LotsOfRaysRayStudy)
Point vertex_average() const