-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy path_ray_tracing_8cc_source.html
493 lines (491 loc) · 79.6 KB
/
_ray_tracing_8cc_source.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "https://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=9"/>
<meta name="generator" content="Doxygen 1.9.1"/>
<meta name="viewport" content="width=device-width, initial-scale=1"/>
<title>adamantine: /home/bt2/Documents/adamantine_doc/doxygen/adamantine-master/source/RayTracing.cc Source File</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="search/search.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="search/searchdata.js"></script>
<script type="text/javascript" src="search/search.js"></script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
<tbody>
<tr style="height: 56px;">
<td id="projectalign" style="padding-left: 0.5em;">
<div id="projectname">adamantine
</div>
</td>
</tr>
</tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.9.1 -->
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&dn=gpl-2.0.txt GPL-v2 */
var searchBox = new SearchBox("searchBox", "search",false,'Search','.html');
/* @license-end */
</script>
<script type="text/javascript" src="menudata.js"></script>
<script type="text/javascript" src="menu.js"></script>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&dn=gpl-2.0.txt GPL-v2 */
$(function() {
initMenu('',true,false,'search.php','Search');
$(document).ready(function() { init_search(); });
});
/* @license-end */</script>
<div id="main-nav"></div>
<!-- window showing the filter options -->
<div id="MSearchSelectWindow"
onmouseover="return searchBox.OnSearchSelectShow()"
onmouseout="return searchBox.OnSearchSelectHide()"
onkeydown="return searchBox.OnSearchSelectKey(event)">
</div>
<!-- iframe showing the search results (closed by default) -->
<div id="MSearchResultsWindow">
<iframe src="javascript:void(0)" frameborder="0"
name="MSearchResults" id="MSearchResults">
</iframe>
</div>
<div id="nav-path" class="navpath">
<ul>
<li class="navelem"><a class="el" href="dir_b2f33c71d4aa5e7af42a1ca61ff5af1b.html">source</a></li> </ul>
</div>
</div><!-- top -->
<div class="header">
<div class="headertitle">
<div class="title">RayTracing.cc</div> </div>
</div><!--header-->
<div class="contents">
<a href="_ray_tracing_8cc.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <span class="comment">/* Copyright (c) 2023 - 2024, the adamantine authors.</span></div>
<div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="comment"> * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception</span></div>
<div class="line"><a name="l00003"></a><span class="lineno"> 3</span> <span class="comment"> */</span></div>
<div class="line"><a name="l00004"></a><span class="lineno"> 4</span>  </div>
<div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="preprocessor">#include <<a class="code" href="_ray_tracing_8hh.html">RayTracing.hh</a>></span></div>
<div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="preprocessor">#include <<a class="code" href="utils_8hh.html">utils.hh</a>></span></div>
<div class="line"><a name="l00007"></a><span class="lineno"> 7</span>  </div>
<div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="preprocessor">#include <deal.II/arborx/distributed_tree.h></span></div>
<div class="line"><a name="l00009"></a><span class="lineno"> 9</span> <span class="preprocessor">#include <deal.II/grid/filtered_iterator.h></span></div>
<div class="line"><a name="l00010"></a><span class="lineno"> 10</span>  </div>
<div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="preprocessor">#include <Kokkos_Core.hpp></span></div>
<div class="line"><a name="l00012"></a><span class="lineno"> 12</span>  </div>
<div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="preprocessor">#include <fstream></span></div>
<div class="line"><a name="l00014"></a><span class="lineno"> 14</span> <span class="preprocessor">#include <regex></span></div>
<div class="line"><a name="l00015"></a><span class="lineno"> 15</span>  </div>
<div class="line"><a name="l00016"></a><span class="lineno"> 16</span> <span class="preprocessor">#include <ArborX_Ray.hpp></span></div>
<div class="line"><a name="l00017"></a><span class="lineno"> 17</span>  </div>
<div class="line"><a name="l00018"></a><span class="lineno"> 18</span> <span class="keyword">namespace </span><a class="code" href="namespaceadamantine.html">adamantine</a></div>
<div class="line"><a name="l00019"></a><span class="lineno"> 19</span> {</div>
<div class="line"><a name="l00025"></a><span class="lineno"><a class="line" href="classadamantine_1_1_ray_nearest_predicate.html"> 25</a></span> <span class="keyword">class </span><a class="code" href="classadamantine_1_1_ray_nearest_predicate.html">RayNearestPredicate</a></div>
<div class="line"><a name="l00026"></a><span class="lineno"> 26</span> {</div>
<div class="line"><a name="l00027"></a><span class="lineno"> 27</span> <span class="keyword">public</span>:</div>
<div class="line"><a name="l00032"></a><span class="lineno"><a class="line" href="classadamantine_1_1_ray_nearest_predicate.html#a69a8378532d6cc9bf4bdeaebbc764ff5"> 32</a></span>  <a class="code" href="classadamantine_1_1_ray_nearest_predicate.html#a69a8378532d6cc9bf4bdeaebbc764ff5">RayNearestPredicate</a>(std::vector<<a class="code" href="structadamantine_1_1_ray.html">Ray<3></a>> <span class="keyword">const</span> &rays) : <a class="code" href="classadamantine_1_1_ray_nearest_predicate.html#a1206cdcc666feb98cacb292ebacb41d5">_rays</a>(rays) {}</div>
<div class="line"><a name="l00033"></a><span class="lineno"> 33</span>  </div>
<div class="line"><a name="l00037"></a><span class="lineno"><a class="line" href="classadamantine_1_1_ray_nearest_predicate.html#a2021c9b114a0dd6c22874e305622e062"> 37</a></span>  std::size_t <a class="code" href="classadamantine_1_1_ray_nearest_predicate.html#a2021c9b114a0dd6c22874e305622e062">size</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> <a class="code" href="classadamantine_1_1_ray_nearest_predicate.html#a1206cdcc666feb98cacb292ebacb41d5">_rays</a>.size(); }</div>
<div class="line"><a name="l00038"></a><span class="lineno"> 38</span>  </div>
<div class="line"><a name="l00042"></a><span class="lineno"><a class="line" href="classadamantine_1_1_ray_nearest_predicate.html#a2d52437d01690e56d8eaabf8354d37a4"> 42</a></span>  <a class="code" href="structadamantine_1_1_ray.html">Ray<3></a> <span class="keyword">const</span> &<a class="code" href="classadamantine_1_1_ray_nearest_predicate.html#a2d52437d01690e56d8eaabf8354d37a4">get</a>(<span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> i)<span class="keyword"> const </span>{ <span class="keywordflow">return</span> <a class="code" href="classadamantine_1_1_ray_nearest_predicate.html#a1206cdcc666feb98cacb292ebacb41d5">_rays</a>[i]; }</div>
<div class="line"><a name="l00043"></a><span class="lineno"> 43</span>  </div>
<div class="line"><a name="l00044"></a><span class="lineno"> 44</span> <span class="keyword">private</span>:</div>
<div class="line"><a name="l00045"></a><span class="lineno"><a class="line" href="classadamantine_1_1_ray_nearest_predicate.html#a1206cdcc666feb98cacb292ebacb41d5"> 45</a></span>  std::vector<Ray<3>> <a class="code" href="classadamantine_1_1_ray_nearest_predicate.html#a1206cdcc666feb98cacb292ebacb41d5">_rays</a>;</div>
<div class="line"><a name="l00046"></a><span class="lineno"> 46</span> };</div>
<div class="line"><a name="l00047"></a><span class="lineno"> 47</span> } <span class="comment">// namespace adamantine</span></div>
<div class="line"><a name="l00048"></a><span class="lineno"> 48</span>  </div>
<div class="line"><a name="l00049"></a><span class="lineno"><a class="line" href="namespace_arbor_x.html"> 49</a></span> <span class="keyword">namespace </span><a class="code" href="namespace_arbor_x.html">ArborX</a></div>
<div class="line"><a name="l00050"></a><span class="lineno"> 50</span> {</div>
<div class="line"><a name="l00051"></a><span class="lineno"> 51</span> <span class="keyword">template</span> <></div>
<div class="line"><a name="l00052"></a><span class="lineno"><a class="line" href="struct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4.html"> 52</a></span> <span class="keyword">struct </span>AccessTraits<<a class="code" href="namespaceadamantine.html">adamantine</a>::RayNearestPredicate, PredicatesTag></div>
<div class="line"><a name="l00053"></a><span class="lineno"> 53</span> {</div>
<div class="line"><a name="l00054"></a><span class="lineno"><a class="line" href="struct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4.html#a4e1b6669289c70e41ddf3a71092362da"> 54</a></span>  <span class="keyword">using</span> <a class="code" href="struct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4.html#a4e1b6669289c70e41ddf3a71092362da">memory_space</a> = Kokkos::HostSpace;</div>
<div class="line"><a name="l00055"></a><span class="lineno"> 55</span>  </div>
<div class="line"><a name="l00056"></a><span class="lineno"><a class="line" href="struct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4.html#a66932a60bd94e9171a4057e8dcbddd06"> 56</a></span>  <span class="keyword">static</span> std::size_t <a class="code" href="struct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4.html#a66932a60bd94e9171a4057e8dcbddd06">size</a>(<a class="code" href="classadamantine_1_1_ray_nearest_predicate.html">adamantine::RayNearestPredicate</a> <span class="keyword">const</span> &ray_nearest)</div>
<div class="line"><a name="l00057"></a><span class="lineno"> 57</span>  {</div>
<div class="line"><a name="l00058"></a><span class="lineno"> 58</span>  <span class="keywordflow">return</span> ray_nearest.<a class="code" href="classadamantine_1_1_ray_nearest_predicate.html#a2021c9b114a0dd6c22874e305622e062">size</a>();</div>
<div class="line"><a name="l00059"></a><span class="lineno"> 59</span>  }</div>
<div class="line"><a name="l00060"></a><span class="lineno"> 60</span>  </div>
<div class="line"><a name="l00061"></a><span class="lineno"><a class="line" href="struct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4.html#ac0ee537fc30b13dfde371f038ff44f63"> 61</a></span>  <span class="keyword">static</span> <span class="keyword">auto</span> <a class="code" href="struct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4.html#ac0ee537fc30b13dfde371f038ff44f63">get</a>(<a class="code" href="classadamantine_1_1_ray_nearest_predicate.html">adamantine::RayNearestPredicate</a> <span class="keyword">const</span> &ray_nearest,</div>
<div class="line"><a name="l00062"></a><span class="lineno"> 62</span>  std::size_t i)</div>
<div class="line"><a name="l00063"></a><span class="lineno"> 63</span>  {</div>
<div class="line"><a name="l00064"></a><span class="lineno"> 64</span>  <span class="keyword">auto</span> <span class="keyword">const</span> &ray = ray_nearest.<a class="code" href="classadamantine_1_1_ray_nearest_predicate.html#a2d52437d01690e56d8eaabf8354d37a4">get</a>(i);</div>
<div class="line"><a name="l00065"></a><span class="lineno"> 65</span>  <span class="keyword">auto</span> <span class="keyword">const</span> &origin = ray.origin;</div>
<div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  <span class="keyword">auto</span> <span class="keyword">const</span> &direction = ray.direction;</div>
<div class="line"><a name="l00067"></a><span class="lineno"> 67</span>  ArborX::Experimental::Ray arborx_ray = {</div>
<div class="line"><a name="l00068"></a><span class="lineno"> 68</span>  ArborX::Point{(float)origin[0], (<span class="keywordtype">float</span>)origin[1], (float)origin[2]},</div>
<div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  ArborX::Experimental::Vector{(float)direction[0], (<span class="keywordtype">float</span>)direction[1],</div>
<div class="line"><a name="l00070"></a><span class="lineno"> 70</span>  (float)direction[2]}};</div>
<div class="line"><a name="l00071"></a><span class="lineno"> 71</span>  <span class="comment">// When the mesh is unstructured, bounding boxes do not tightly bound the</span></div>
<div class="line"><a name="l00072"></a><span class="lineno"> 72</span>  <span class="comment">// cells and so many of them overlap. A ray may hit a bounding box but may</span></div>
<div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  <span class="comment">// missed the cell. We need to ask for more bounding boxes to increase the</span></div>
<div class="line"><a name="l00074"></a><span class="lineno"> 74</span>  <span class="comment">// chance that we don't miss a ray-cell interaction. Asking for eight</span></div>
<div class="line"><a name="l00075"></a><span class="lineno"> 75</span>  <span class="comment">// bounding boxes seems to work pretty well.</span></div>
<div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  <span class="keywordflow">return</span> nearest(arborx_ray, 8);</div>
<div class="line"><a name="l00077"></a><span class="lineno"> 77</span>  }</div>
<div class="line"><a name="l00078"></a><span class="lineno"> 78</span> };</div>
<div class="line"><a name="l00079"></a><span class="lineno"> 79</span> } <span class="comment">// namespace ArborX</span></div>
<div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  </div>
<div class="line"><a name="l00081"></a><span class="lineno"> 81</span> <span class="keyword">namespace</span></div>
<div class="line"><a name="l00082"></a><span class="lineno"> 82</span> {</div>
<div class="line"><a name="l00083"></a><span class="lineno"> 83</span> <span class="keyword">template</span> <<span class="keywordtype">int</span> dim></div>
<div class="line"><a name="l00084"></a><span class="lineno"> 84</span> <span class="keywordtype">bool</span> point_in_triangle(dealii::Point<dim> <span class="keyword">const</span> &a, dealii::Point<dim> <span class="keyword">const</span> &b,</div>
<div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  dealii::Point<dim> <span class="keyword">const</span> &c, dealii::Point<dim> <span class="keyword">const</span> &p)</div>
<div class="line"><a name="l00086"></a><span class="lineno"> 86</span> {</div>
<div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  <span class="comment">// To decide if a point is inside the triangle, we compute the barycentric</span></div>
<div class="line"><a name="l00088"></a><span class="lineno"> 88</span>  <span class="comment">// coordinate of the points. If they are positive, the point is in the</span></div>
<div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  <span class="comment">// triangle. See:</span></div>
<div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  <span class="comment">// https://math.stackexchange.com/questions/544946/determine-if-projection-of-3d-point-onto-plane-is-within-a-triangle/544947</span></div>
<div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  <span class="comment">// and</span></div>
<div class="line"><a name="l00092"></a><span class="lineno"> 92</span>  <span class="comment">// https://gamedev.stackexchange.com/questions/23743/whats-the-most-efficient-way-to-find-barycentric-coordinates</span></div>
<div class="line"><a name="l00093"></a><span class="lineno"> 93</span>  dealii::Tensor<1, dim> ab({b[0] - a[0], b[1] - a[1], b[2] - a[2]});</div>
<div class="line"><a name="l00094"></a><span class="lineno"> 94</span>  dealii::Tensor<1, dim> ac({c[0] - a[0], c[1] - a[1], c[2] - a[2]});</div>
<div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  dealii::Tensor<1, dim> ap({p[0] - a[0], p[1] - a[1], p[2] - a[2]});</div>
<div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> ab_ab = ab * ab;</div>
<div class="line"><a name="l00097"></a><span class="lineno"> 97</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> ab_ac = ab * ac;</div>
<div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> ac_ac = ac * ac;</div>
<div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> ap_ab = ap * ab;</div>
<div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> ap_ac = ap * ac;</div>
<div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> denom = ab_ab * ac_ac - ab_ac * ab_ac;</div>
<div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> v = (ac_ac * ap_ab - ab_ac * ap_ac) / denom;</div>
<div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> w = (ab_ab * ap_ac - ab_ac * ap_ab) / denom;</div>
<div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> u = 1. - v - w;</div>
<div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  <span class="keywordflow">if</span> ((u >= 0.) && (v >= 0.) && (w >= 0.))</div>
<div class="line"><a name="l00106"></a><span class="lineno"> 106</span>  {</div>
<div class="line"><a name="l00107"></a><span class="lineno"> 107</span>  <span class="keywordflow">return</span> <span class="keyword">true</span>;</div>
<div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  }</div>
<div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  </div>
<div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  <span class="keywordflow">return</span> <span class="keyword">false</span>;</div>
<div class="line"><a name="l00111"></a><span class="lineno"> 111</span> }</div>
<div class="line"><a name="l00112"></a><span class="lineno"> 112</span> } <span class="comment">// namespace</span></div>
<div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  </div>
<div class="line"><a name="l00114"></a><span class="lineno"> 114</span> <span class="keyword">namespace </span><a class="code" href="namespaceadamantine.html">adamantine</a></div>
<div class="line"><a name="l00115"></a><span class="lineno"> 115</span> {</div>
<div class="line"><a name="l00116"></a><span class="lineno"><a class="line" href="classadamantine_1_1_ray_tracing.html#a311145cc3e9dccf9674939baa8205f0e"> 116</a></span> <a class="code" href="classadamantine_1_1_ray_tracing.html#a311145cc3e9dccf9674939baa8205f0e">RayTracing::RayTracing</a>(boost::property_tree::ptree <span class="keyword">const</span> &experiment_database,</div>
<div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  dealii::DoFHandler<3> <span class="keyword">const</span> &dof_handler)</div>
<div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  : _dof_handler(dof_handler)</div>
<div class="line"><a name="l00119"></a><span class="lineno"> 119</span> {</div>
<div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  </div>
<div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  <span class="comment">// Format of the file names: the format is pretty arbitrary, #frame and</span></div>
<div class="line"><a name="l00122"></a><span class="lineno"> 122</span>  <span class="comment">// #camera are replaced by the frame and the camera number.</span></div>
<div class="line"><a name="l00123"></a><span class="lineno"> 123</span>  <span class="comment">// PropertyTreeInput experiment.file</span></div>
<div class="line"><a name="l00124"></a><span class="lineno"> 124</span>  <a class="code" href="classadamantine_1_1_ray_tracing.html#a24a8f437b8e3e1e73f880e52a66fa456">_data_filename</a> = experiment_database.get<std::string>(<span class="stringliteral">"file"</span>);</div>
<div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  <span class="comment">// PropertyTreeInput experiment.first_frame</span></div>
<div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  <a class="code" href="classadamantine_1_1_ray_tracing.html#aa58ad004836a0c29a24bde18d3ba4947">_next_frame</a> = experiment_database.get(<span class="stringliteral">"first_frame"</span>, 0);</div>
<div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  <span class="comment">// PropertyTreeInput experiment.first_camera_id</span></div>
<div class="line"><a name="l00128"></a><span class="lineno"> 128</span>  <a class="code" href="classadamantine_1_1_ray_tracing.html#a5e1ffea3480cec7e73cbe0a5c354d72d">_first_camera_id</a> = experiment_database.get<<span class="keywordtype">unsigned</span> <span class="keywordtype">int</span>>(<span class="stringliteral">"first_camera_id"</span>);</div>
<div class="line"><a name="l00129"></a><span class="lineno"> 129</span>  <span class="comment">// PropertyTreeInput experiment.last_camera_id</span></div>
<div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  <a class="code" href="classadamantine_1_1_ray_tracing.html#ac6abef2dbdc651660422877350428b20">_last_camera_id</a> = experiment_database.get<<span class="keywordtype">int</span>>(<span class="stringliteral">"last_camera_id"</span>);</div>
<div class="line"><a name="l00131"></a><span class="lineno"> 131</span> }</div>
<div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  </div>
<div class="line"><a name="l00133"></a><span class="lineno"><a class="line" href="classadamantine_1_1_ray_tracing.html#a3cb1e406035092d2f0825455202c4521"> 133</a></span> <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> <a class="code" href="classadamantine_1_1_ray_tracing.html#a3cb1e406035092d2f0825455202c4521">RayTracing::read_next_frame</a>()</div>
<div class="line"><a name="l00134"></a><span class="lineno"> 134</span> {</div>
<div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  <a class="code" href="classadamantine_1_1_ray_tracing.html#af24254cda95c0dc2084f74411a7661fe">_rays_current_frame</a>.clear();</div>
<div class="line"><a name="l00136"></a><span class="lineno"> 136</span>  <a class="code" href="classadamantine_1_1_ray_tracing.html#a14cdd076d9145977877d8d595a4b70d4">_values_current_frame</a>.clear();</div>
<div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  <span class="keywordflow">for</span> (<span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> camera_id = <a class="code" href="classadamantine_1_1_ray_tracing.html#a5e1ffea3480cec7e73cbe0a5c354d72d">_first_camera_id</a>;</div>
<div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  camera_id < <a class="code" href="classadamantine_1_1_ray_tracing.html#ac6abef2dbdc651660422877350428b20">_last_camera_id</a> + 1; ++camera_id)</div>
<div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  {</div>
<div class="line"><a name="l00140"></a><span class="lineno"> 140</span>  <span class="comment">// Use regex to get the next file to read</span></div>
<div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  std::regex camera_regex(<span class="stringliteral">"#camera"</span>);</div>
<div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  std::regex frame_regex(<span class="stringliteral">"#frame"</span>);</div>
<div class="line"><a name="l00143"></a><span class="lineno"> 143</span>  <span class="keyword">auto</span> filename =</div>
<div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  std::regex_replace((std::regex_replace(<a class="code" href="classadamantine_1_1_ray_tracing.html#a24a8f437b8e3e1e73f880e52a66fa456">_data_filename</a>, camera_regex,</div>
<div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  std::to_string(camera_id))),</div>
<div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  frame_regex, std::to_string(<a class="code" href="classadamantine_1_1_ray_tracing.html#aa58ad004836a0c29a24bde18d3ba4947">_next_frame</a>));</div>
<div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  <a class="code" href="namespaceadamantine.html#a737967bbc854d35e1281ad917ba24814">wait_for_file</a>(filename, <span class="stringliteral">"Waiting for the next frame: "</span> + filename);</div>
<div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  </div>
<div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  <span class="comment">// Read and parse the file</span></div>
<div class="line"><a name="l00150"></a><span class="lineno"> 150</span>  std::ifstream file;</div>
<div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  file.open(filename);</div>
<div class="line"><a name="l00152"></a><span class="lineno"> 152</span>  std::string <a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474">line</a>;</div>
<div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  std::getline(file, <a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474">line</a>); <span class="comment">// skip the header</span></div>
<div class="line"><a name="l00154"></a><span class="lineno"> 154</span>  <span class="keywordflow">while</span> (std::getline(file, <a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474">line</a>))</div>
<div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  {</div>
<div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  std::size_t pos = 0;</div>
<div class="line"><a name="l00157"></a><span class="lineno"> 157</span>  std::size_t last_pos = 0;</div>
<div class="line"><a name="l00158"></a><span class="lineno"> 158</span>  std::size_t line_length = <a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474">line</a>.length();</div>
<div class="line"><a name="l00159"></a><span class="lineno"> 159</span>  <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> i = 0;</div>
<div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  dealii::Point<dim> <a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a78ee54aa8f813885fe2fe20d232518b9">point</a>;</div>
<div class="line"><a name="l00161"></a><span class="lineno"> 161</span>  dealii::Tensor<1, dim> direction;</div>
<div class="line"><a name="l00162"></a><span class="lineno"> 162</span>  <span class="keywordtype">double</span> value = 0.;</div>
<div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  <span class="keywordflow">while</span> (last_pos < line_length + 1)</div>
<div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  {</div>
<div class="line"><a name="l00165"></a><span class="lineno"> 165</span>  pos = <a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474">line</a>.find_first_of(<span class="charliteral">','</span>, last_pos);</div>
<div class="line"><a name="l00166"></a><span class="lineno"> 166</span>  <span class="comment">// If no comma was found that we read until the end of the file</span></div>
<div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  <span class="keywordflow">if</span> (pos == std::string::npos)</div>
<div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  {</div>
<div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  pos = line_length;</div>
<div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  }</div>
<div class="line"><a name="l00171"></a><span class="lineno"> 171</span>  </div>
<div class="line"><a name="l00172"></a><span class="lineno"> 172</span>  <span class="keywordflow">if</span> (pos != last_pos)</div>
<div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  {</div>
<div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  <span class="keywordtype">char</span> *end = <a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474">line</a>.data() + pos;</div>
<div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  <span class="keywordflow">if</span> (i < <a class="code" href="classadamantine_1_1_ray_tracing.html#a11ad68edb125a2f2aca6538a1f4c33f4">dim</a>)</div>
<div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  {</div>
<div class="line"><a name="l00177"></a><span class="lineno"> 177</span>  <a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a78ee54aa8f813885fe2fe20d232518b9">point</a>[i] = std::strtod(<a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474">line</a>.data() + last_pos, &end);</div>
<div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  }</div>
<div class="line"><a name="l00179"></a><span class="lineno"> 179</span>  <span class="keywordflow">else</span> <span class="keywordflow">if</span> (i < 2 * <a class="code" href="classadamantine_1_1_ray_tracing.html#a11ad68edb125a2f2aca6538a1f4c33f4">dim</a>)</div>
<div class="line"><a name="l00180"></a><span class="lineno"> 180</span>  {</div>
<div class="line"><a name="l00181"></a><span class="lineno"> 181</span>  <span class="comment">// Calculate the direction from the first and second points in the</span></div>
<div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  <span class="comment">// file</span></div>
<div class="line"><a name="l00183"></a><span class="lineno"> 183</span>  direction[i - <a class="code" href="classadamantine_1_1_ray_tracing.html#a11ad68edb125a2f2aca6538a1f4c33f4">dim</a>] =</div>
<div class="line"><a name="l00184"></a><span class="lineno"> 184</span>  std::strtod(<a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474">line</a>.data() + last_pos, &end) - <a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a78ee54aa8f813885fe2fe20d232518b9">point</a>[i - <a class="code" href="classadamantine_1_1_ray_tracing.html#a11ad68edb125a2f2aca6538a1f4c33f4">dim</a>];</div>
<div class="line"><a name="l00185"></a><span class="lineno"> 185</span>  }</div>
<div class="line"><a name="l00186"></a><span class="lineno"> 186</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00187"></a><span class="lineno"> 187</span>  {</div>
<div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  value = std::strtod(<a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474">line</a>.data() + last_pos, &end);</div>
<div class="line"><a name="l00189"></a><span class="lineno"> 189</span>  }</div>
<div class="line"><a name="l00190"></a><span class="lineno"> 190</span>  </div>
<div class="line"><a name="l00191"></a><span class="lineno"> 191</span>  ++i;</div>
<div class="line"><a name="l00192"></a><span class="lineno"> 192</span>  }</div>
<div class="line"><a name="l00193"></a><span class="lineno"> 193</span>  </div>
<div class="line"><a name="l00194"></a><span class="lineno"> 194</span>  last_pos = pos + 1;</div>
<div class="line"><a name="l00195"></a><span class="lineno"> 195</span>  }</div>
<div class="line"><a name="l00196"></a><span class="lineno"> 196</span>  </div>
<div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  <a class="code" href="structadamantine_1_1_ray.html">Ray<dim></a> ray{<a class="code" href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a78ee54aa8f813885fe2fe20d232518b9">point</a>, direction};</div>
<div class="line"><a name="l00198"></a><span class="lineno"> 198</span>  <a class="code" href="classadamantine_1_1_ray_tracing.html#af24254cda95c0dc2084f74411a7661fe">_rays_current_frame</a>.push_back(ray);</div>
<div class="line"><a name="l00199"></a><span class="lineno"> 199</span>  <a class="code" href="classadamantine_1_1_ray_tracing.html#a14cdd076d9145977877d8d595a4b70d4">_values_current_frame</a>.push_back(value);</div>
<div class="line"><a name="l00200"></a><span class="lineno"> 200</span>  }</div>
<div class="line"><a name="l00201"></a><span class="lineno"> 201</span>  }</div>
<div class="line"><a name="l00202"></a><span class="lineno"> 202</span>  </div>
<div class="line"><a name="l00203"></a><span class="lineno"> 203</span>  <span class="keywordflow">return</span> <a class="code" href="classadamantine_1_1_ray_tracing.html#aa58ad004836a0c29a24bde18d3ba4947">_next_frame</a>++;</div>
<div class="line"><a name="l00204"></a><span class="lineno"> 204</span> }</div>
<div class="line"><a name="l00205"></a><span class="lineno"> 205</span>  </div>
<div class="line"><a name="l00206"></a><span class="lineno"><a class="line" href="classadamantine_1_1_ray_tracing.html#a7e9becfe37db0db6a01ca884f52b3335"> 206</a></span> <a class="code" href="structadamantine_1_1_points_values.html">PointsValues<3></a> <a class="code" href="classadamantine_1_1_ray_tracing.html#a7e9becfe37db0db6a01ca884f52b3335">RayTracing::get_points_values</a>()</div>
<div class="line"><a name="l00207"></a><span class="lineno"> 207</span> {</div>
<div class="line"><a name="l00208"></a><span class="lineno"> 208</span>  <span class="comment">// Perform the ray tracing to get the cells that are intersected by rays</span></div>
<div class="line"><a name="l00209"></a><span class="lineno"> 209</span>  </div>
<div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  <span class="comment">// Create the bounding boxes associated to the locally owned cells with FE</span></div>
<div class="line"><a name="l00211"></a><span class="lineno"> 211</span>  <span class="comment">// index = 0</span></div>
<div class="line"><a name="l00212"></a><span class="lineno"> 212</span>  std::vector<dealii::BoundingBox<dim>> bounding_boxes;</div>
<div class="line"><a name="l00213"></a><span class="lineno"> 213</span>  std::vector<typename dealii::DoFHandler<dim>::active_cell_iterator></div>
<div class="line"><a name="l00214"></a><span class="lineno"> 214</span>  cell_iterators;</div>
<div class="line"><a name="l00215"></a><span class="lineno"> 215</span>  <span class="keywordflow">for</span> (<span class="keyword">auto</span> <span class="keyword">const</span> &cell : dealii::filter_iterators(</div>
<div class="line"><a name="l00216"></a><span class="lineno"> 216</span>  <a class="code" href="classadamantine_1_1_ray_tracing.html#a74976dc4d69bf1b3f23c0c25f1824384">_dof_handler</a>.active_cell_iterators(),</div>
<div class="line"><a name="l00217"></a><span class="lineno"> 217</span>  dealii::IteratorFilters::LocallyOwnedCell(),</div>
<div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  dealii::IteratorFilters::ActiveFEIndexEqualTo(0)))</div>
<div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  {</div>
<div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  bounding_boxes.push_back(cell->bounding_box());</div>
<div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  cell_iterators.push_back(cell);</div>
<div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  }</div>
<div class="line"><a name="l00223"></a><span class="lineno"> 223</span>  </div>
<div class="line"><a name="l00224"></a><span class="lineno"> 224</span>  <span class="comment">// Use ArborX to find where the rays intersect the activated cells. All the</span></div>
<div class="line"><a name="l00225"></a><span class="lineno"> 225</span>  <span class="comment">// processors have access to all the rays but we still need to use</span></div>
<div class="line"><a name="l00226"></a><span class="lineno"> 226</span>  <span class="comment">// DistributedTree because some rays can be stopped by activated cells on a</span></div>
<div class="line"><a name="l00227"></a><span class="lineno"> 227</span>  <span class="comment">// different processors. Since the rays are on all the processors, we don't</span></div>
<div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  <span class="comment">// need to communicate the results to other processors. Note that we may lose</span></div>
<div class="line"><a name="l00229"></a><span class="lineno"> 229</span>  <span class="comment">// some rays because the bounding boxes are larger than the cells and so a ray</span></div>
<div class="line"><a name="l00230"></a><span class="lineno"> 230</span>  <span class="comment">// can be stopped by a bounding box and then miss the associated cell. We</span></div>
<div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  <span class="comment">// could get around this by asking ArborX for all the bounding boxes that the</span></div>
<div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  <span class="comment">// rays intersect but then we would need to keep track of the order of the</span></div>
<div class="line"><a name="l00233"></a><span class="lineno"> 233</span>  <span class="comment">// intersections ourselves.</span></div>
<div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  <span class="comment">// TODO The current code can be simplified a lot once we can use ArborX 2.0</span></div>
<div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  <span class="comment">// This version of ArborX supports distributed ray tracing on triangle.</span></div>
<div class="line"><a name="l00236"></a><span class="lineno"> 236</span>  <span class="comment">// Currently we create a bounding box for each cell, use ArborX to perform a</span></div>
<div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  <span class="comment">// coarse search, and finally perform a fine search. With ArborX 2.0, we will</span></div>
<div class="line"><a name="l00238"></a><span class="lineno"> 238</span>  <span class="comment">// need to create triangles for each face of each cell and then, ArborX will</span></div>
<div class="line"><a name="l00239"></a><span class="lineno"> 239</span>  <span class="comment">// perform the ray tracing on these triangles. This will simplify the current</span></div>
<div class="line"><a name="l00240"></a><span class="lineno"> 240</span>  <span class="comment">// code and will solve the problem of missing ray-cell intersections.</span></div>
<div class="line"><a name="l00241"></a><span class="lineno"> 241</span>  <span class="keyword">auto</span> communicator = <a class="code" href="classadamantine_1_1_ray_tracing.html#a74976dc4d69bf1b3f23c0c25f1824384">_dof_handler</a>.get_communicator();</div>
<div class="line"><a name="l00242"></a><span class="lineno"> 242</span>  dealii::ArborXWrappers::DistributedTree distributed_tree(communicator,</div>
<div class="line"><a name="l00243"></a><span class="lineno"> 243</span>  bounding_boxes);</div>
<div class="line"><a name="l00244"></a><span class="lineno"> 244</span>  <a class="code" href="classadamantine_1_1_ray_nearest_predicate.html">RayNearestPredicate</a> ray_nearest(<a class="code" href="classadamantine_1_1_ray_tracing.html#af24254cda95c0dc2084f74411a7661fe">_rays_current_frame</a>);</div>
<div class="line"><a name="l00245"></a><span class="lineno"> 245</span>  <span class="keyword">auto</span> [indices_ranks, offset] = distributed_tree.query(ray_nearest);</div>
<div class="line"><a name="l00246"></a><span class="lineno"> 246</span>  </div>
<div class="line"><a name="l00247"></a><span class="lineno"> 247</span>  <span class="comment">// Find the exact intersections points</span></div>
<div class="line"><a name="l00248"></a><span class="lineno"> 248</span>  <span class="comment">// See https://en.wikipedia.org/wiki/Line%E2%80%93plane_intersection</span></div>
<div class="line"><a name="l00249"></a><span class="lineno"> 249</span>  <span class="comment">// NOTE that we assume that the faces are flat. If the faces are curved,</span></div>
<div class="line"><a name="l00250"></a><span class="lineno"> 250</span>  <span class="comment">// this is wrong.</span></div>
<div class="line"><a name="l00251"></a><span class="lineno"> 251</span>  <span class="keywordtype">int</span> <span class="keyword">const</span> my_rank = dealii::Utilities::MPI::this_mpi_process(communicator);</div>
<div class="line"><a name="l00252"></a><span class="lineno"> 252</span>  <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> <span class="keyword">const</span> n_rays = <a class="code" href="classadamantine_1_1_ray_tracing.html#af24254cda95c0dc2084f74411a7661fe">_rays_current_frame</a>.size();</div>
<div class="line"><a name="l00253"></a><span class="lineno"> 253</span>  <span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> n_intersections = 0;</div>
<div class="line"><a name="l00254"></a><span class="lineno"> 254</span>  <span class="keywordflow">for</span> (<span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> i = 0; i < n_rays; ++i)</div>
<div class="line"><a name="l00255"></a><span class="lineno"> 255</span>  {</div>
<div class="line"><a name="l00256"></a><span class="lineno"> 256</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> j = offset[i]; j < offset[i + 1]; ++j)</div>
<div class="line"><a name="l00257"></a><span class="lineno"> 257</span>  {</div>
<div class="line"><a name="l00258"></a><span class="lineno"> 258</span>  <span class="keywordflow">if</span> (indices_ranks[j].second == my_rank)</div>
<div class="line"><a name="l00259"></a><span class="lineno"> 259</span>  {</div>
<div class="line"><a name="l00260"></a><span class="lineno"> 260</span>  ++n_intersections;</div>
<div class="line"><a name="l00261"></a><span class="lineno"> 261</span>  }</div>
<div class="line"><a name="l00262"></a><span class="lineno"> 262</span>  }</div>
<div class="line"><a name="l00263"></a><span class="lineno"> 263</span>  }</div>
<div class="line"><a name="l00264"></a><span class="lineno"> 264</span>  </div>
<div class="line"><a name="l00265"></a><span class="lineno"> 265</span>  std::vector<dealii::Point<dim>> points;</div>
<div class="line"><a name="l00266"></a><span class="lineno"> 266</span>  std::vector<double> values;</div>
<div class="line"><a name="l00267"></a><span class="lineno"> 267</span>  <span class="keywordflow">if</span> (n_intersections != 0)</div>
<div class="line"><a name="l00268"></a><span class="lineno"> 268</span>  {</div>
<div class="line"><a name="l00269"></a><span class="lineno"> 269</span>  points.reserve(n_intersections);</div>
<div class="line"><a name="l00270"></a><span class="lineno"> 270</span>  values.reserve(n_intersections);</div>
<div class="line"><a name="l00271"></a><span class="lineno"> 271</span>  <span class="keyword">auto</span> constexpr reference_cell =</div>
<div class="line"><a name="l00272"></a><span class="lineno"> 272</span>  dealii::ReferenceCells::get_hypercube<dim>();</div>
<div class="line"><a name="l00273"></a><span class="lineno"> 273</span>  <span class="keywordtype">double</span> constexpr tol = 1e-10;</div>
<div class="line"><a name="l00274"></a><span class="lineno"> 274</span>  <span class="keywordflow">for</span> (<span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> i = 0; i < n_rays; ++i)</div>
<div class="line"><a name="l00275"></a><span class="lineno"> 275</span>  {</div>
<div class="line"><a name="l00276"></a><span class="lineno"> 276</span>  dealii::Point<dim> intersection;</div>
<div class="line"><a name="l00277"></a><span class="lineno"> 277</span>  <span class="keywordtype">double</span> distance = std::numeric_limits<double>::max();</div>
<div class="line"><a name="l00278"></a><span class="lineno"> 278</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> j = offset[i]; j < offset[i + 1]; ++j)</div>
<div class="line"><a name="l00279"></a><span class="lineno"> 279</span>  {</div>
<div class="line"><a name="l00280"></a><span class="lineno"> 280</span>  <span class="comment">// The bounding boxes found can be on different processors. Instead of</span></div>
<div class="line"><a name="l00281"></a><span class="lineno"> 281</span>  <span class="comment">// performing the fine search and then comparing results between</span></div>
<div class="line"><a name="l00282"></a><span class="lineno"> 282</span>  <span class="comment">// processors, we just filter out the bounding boxes found on a</span></div>
<div class="line"><a name="l00283"></a><span class="lineno"> 283</span>  <span class="comment">// processor different than indices_ranks[offet[i]].second. Some rays</span></div>
<div class="line"><a name="l00284"></a><span class="lineno"> 284</span>  <span class="comment">// will be lost but the loss should be minimal.</span></div>
<div class="line"><a name="l00285"></a><span class="lineno"> 285</span>  <span class="keywordflow">if</span> ((indices_ranks[j].second == indices_ranks[offset[i]].second) &&</div>
<div class="line"><a name="l00286"></a><span class="lineno"> 286</span>  (indices_ranks[j].second == my_rank))</div>
<div class="line"><a name="l00287"></a><span class="lineno"> 287</span>  {</div>
<div class="line"><a name="l00288"></a><span class="lineno"> 288</span>  <span class="keyword">auto</span> <span class="keyword">const</span> &cell = cell_iterators[indices_ranks[j].first];</div>
<div class="line"><a name="l00289"></a><span class="lineno"> 289</span>  </div>
<div class="line"><a name="l00290"></a><span class="lineno"> 290</span>  <span class="comment">// We know that the ray intersects the bounding box but we don't know</span></div>
<div class="line"><a name="l00291"></a><span class="lineno"> 291</span>  <span class="comment">// where it intersects the cells. We need to check the intersection of</span></div>
<div class="line"><a name="l00292"></a><span class="lineno"> 292</span>  <span class="comment">// the ray with each face of the cell. To do that, we split each face</span></div>
<div class="line"><a name="l00293"></a><span class="lineno"> 293</span>  <span class="comment">// in two triangles and check if the point belongs either triangles.</span></div>
<div class="line"><a name="l00294"></a><span class="lineno"> 294</span>  <span class="comment">// We then use the barycentric coordinate method which involves</span></div>
<div class="line"><a name="l00295"></a><span class="lineno"> 295</span>  <span class="comment">// calculating the "weights" of the point relative to each vertex of</span></div>
<div class="line"><a name="l00296"></a><span class="lineno"> 296</span>  <span class="comment">// the triangle; if all weights are non-negative, then the point lies</span></div>
<div class="line"><a name="l00297"></a><span class="lineno"> 297</span>  <span class="comment">// inside the triangle.</span></div>
<div class="line"><a name="l00298"></a><span class="lineno"> 298</span>  <span class="keywordflow">for</span> (<span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> f = 0; f < reference_cell.n_faces(); ++f)</div>
<div class="line"><a name="l00299"></a><span class="lineno"> 299</span>  {</div>
<div class="line"><a name="l00300"></a><span class="lineno"> 300</span>  <span class="keyword">auto</span> <span class="keyword">const</span> point_0 = cell->face(f)->vertex(0);</div>
<div class="line"><a name="l00301"></a><span class="lineno"> 301</span>  <span class="keyword">auto</span> <span class="keyword">const</span> point_1 = cell->face(f)->vertex(1);</div>
<div class="line"><a name="l00302"></a><span class="lineno"> 302</span>  <span class="keyword">auto</span> <span class="keyword">const</span> point_2 = cell->face(f)->vertex(2);</div>
<div class="line"><a name="l00303"></a><span class="lineno"> 303</span>  <span class="keyword">auto</span> <span class="keyword">const</span> point_3 = cell->face(f)->vertex(3);</div>
<div class="line"><a name="l00304"></a><span class="lineno"> 304</span>  <span class="comment">// First we check if the ray is parallel to the face. If this is the</span></div>
<div class="line"><a name="l00305"></a><span class="lineno"> 305</span>  <span class="comment">// case, either the ray misses the face or the ray hits the edge of</span></div>
<div class="line"><a name="l00306"></a><span class="lineno"> 306</span>  <span class="comment">// the face. In that last case, the ray is also orthogonal to</span></div>
<div class="line"><a name="l00307"></a><span class="lineno"> 307</span>  <span class="comment">// another face and it is safe to discard all rays parallel to a</span></div>
<div class="line"><a name="l00308"></a><span class="lineno"> 308</span>  <span class="comment">// face.</span></div>
<div class="line"><a name="l00309"></a><span class="lineno"> 309</span>  dealii::Tensor<1, dim> <span class="keyword">const</span> edge_01({point_1[0] - point_0[0],</div>
<div class="line"><a name="l00310"></a><span class="lineno"> 310</span>  point_1[1] - point_0[1],</div>
<div class="line"><a name="l00311"></a><span class="lineno"> 311</span>  point_1[2] - point_0[2]});</div>
<div class="line"><a name="l00312"></a><span class="lineno"> 312</span>  dealii::Tensor<1, dim> <span class="keyword">const</span> edge_02({point_2[0] - point_0[0],</div>
<div class="line"><a name="l00313"></a><span class="lineno"> 313</span>  point_2[1] - point_0[1],</div>
<div class="line"><a name="l00314"></a><span class="lineno"> 314</span>  point_2[2] - point_0[2]});</div>
<div class="line"><a name="l00315"></a><span class="lineno"> 315</span>  dealii::Tensor<1, dim> <span class="keyword">const</span> ray_direction(</div>
<div class="line"><a name="l00316"></a><span class="lineno"> 316</span>  {<span class="keyword">static_cast<</span><span class="keywordtype">double</span><span class="keyword">></span>(<a class="code" href="classadamantine_1_1_ray_tracing.html#af24254cda95c0dc2084f74411a7661fe">_rays_current_frame</a>[i].direction[0]),</div>
<div class="line"><a name="l00317"></a><span class="lineno"> 317</span>  <span class="keyword">static_cast<</span><span class="keywordtype">double</span><span class="keyword">></span>(<a class="code" href="classadamantine_1_1_ray_tracing.html#af24254cda95c0dc2084f74411a7661fe">_rays_current_frame</a>[i].direction[1]),</div>
<div class="line"><a name="l00318"></a><span class="lineno"> 318</span>  <span class="keyword">static_cast<</span><span class="keywordtype">double</span><span class="keyword">></span>(<a class="code" href="classadamantine_1_1_ray_tracing.html#af24254cda95c0dc2084f74411a7661fe">_rays_current_frame</a>[i].direction[2])});</div>
<div class="line"><a name="l00319"></a><span class="lineno"> 319</span>  dealii::Tensor<2, dim> matrix(</div>
<div class="line"><a name="l00320"></a><span class="lineno"> 320</span>  {{-ray_direction[0], -ray_direction[1], -ray_direction[2]},</div>
<div class="line"><a name="l00321"></a><span class="lineno"> 321</span>  {edge_01[0], edge_01[1], edge_01[2]},</div>
<div class="line"><a name="l00322"></a><span class="lineno"> 322</span>  {edge_02[0], edge_02[1], edge_02[2]}});</div>
<div class="line"><a name="l00323"></a><span class="lineno"> 323</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> det = dealii::determinant(matrix);</div>
<div class="line"><a name="l00324"></a><span class="lineno"> 324</span>  <span class="keywordtype">double</span> face_area = 0.;</div>
<div class="line"><a name="l00325"></a><span class="lineno"> 325</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> e = 0; e < <a class="code" href="classadamantine_1_1_ray_tracing.html#a11ad68edb125a2f2aca6538a1f4c33f4">dim</a>; ++e)</div>
<div class="line"><a name="l00326"></a><span class="lineno"> 326</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> k = 0; k < <a class="code" href="classadamantine_1_1_ray_tracing.html#a11ad68edb125a2f2aca6538a1f4c33f4">dim</a>; ++k)</div>
<div class="line"><a name="l00327"></a><span class="lineno"> 327</span>  face_area += std::abs(edge_01[e] * edge_02[k]);</div>
<div class="line"><a name="l00328"></a><span class="lineno"> 328</span>  <span class="comment">// If determinant is close to zero, the ray is parallel to the face</span></div>
<div class="line"><a name="l00329"></a><span class="lineno"> 329</span>  <span class="comment">// and we go to the next face.</span></div>
<div class="line"><a name="l00330"></a><span class="lineno"> 330</span>  <span class="keywordflow">if</span> (std::abs(det) < tol * face_area)</div>
<div class="line"><a name="l00331"></a><span class="lineno"> 331</span>  <span class="keywordflow">continue</span>;</div>
<div class="line"><a name="l00332"></a><span class="lineno"> 332</span>  <span class="comment">// Compute the distance along the ray direction between the origin</span></div>
<div class="line"><a name="l00333"></a><span class="lineno"> 333</span>  <span class="comment">// of the ray and the intersection point.</span></div>
<div class="line"><a name="l00334"></a><span class="lineno"> 334</span>  <span class="keyword">auto</span> <span class="keyword">const</span> cross_product =</div>
<div class="line"><a name="l00335"></a><span class="lineno"> 335</span>  dealii::cross_product_3d(edge_01, edge_02);</div>
<div class="line"><a name="l00336"></a><span class="lineno"> 336</span>  <span class="keyword">auto</span> <span class="keyword">const</span> &ray_origin = <a class="code" href="classadamantine_1_1_ray_tracing.html#af24254cda95c0dc2084f74411a7661fe">_rays_current_frame</a>[i].origin;</div>
<div class="line"><a name="l00337"></a><span class="lineno"> 337</span>  dealii::Tensor<1, dim> p0_ray({ray_origin[0] - point_0[0],</div>
<div class="line"><a name="l00338"></a><span class="lineno"> 338</span>  ray_origin[1] - point_0[1],</div>
<div class="line"><a name="l00339"></a><span class="lineno"> 339</span>  ray_origin[2] - point_0[2]});</div>
<div class="line"><a name="l00340"></a><span class="lineno"> 340</span>  <span class="keywordtype">double</span> <span class="keyword">const</span> d = cross_product * p0_ray / det;</div>
<div class="line"><a name="l00341"></a><span class="lineno"> 341</span>  <span class="comment">// If d is less than 0, this means that the ray intersects the plane</span></div>
<div class="line"><a name="l00342"></a><span class="lineno"> 342</span>  <span class="comment">// of the face but not the face itself. We can exit early.</span></div>
<div class="line"><a name="l00343"></a><span class="lineno"> 343</span>  <span class="keywordflow">if</span> (d < 0)</div>
<div class="line"><a name="l00344"></a><span class="lineno"> 344</span>  <span class="keywordflow">continue</span>;</div>
<div class="line"><a name="l00345"></a><span class="lineno"> 345</span>  <span class="comment">// We can finally compute the intersection point. It is possible</span></div>
<div class="line"><a name="l00346"></a><span class="lineno"> 346</span>  <span class="comment">// that a ray intersects multiple faces. For instance if the mesh is</span></div>
<div class="line"><a name="l00347"></a><span class="lineno"> 347</span>  <span class="comment">// a cube the ray will get into the cube from one face and it will</span></div>
<div class="line"><a name="l00348"></a><span class="lineno"> 348</span>  <span class="comment">// get out of the cube by the opposite face. The correct</span></div>
<div class="line"><a name="l00349"></a><span class="lineno"> 349</span>  <span class="comment">// intersection point is the one with the smallest distance.</span></div>
<div class="line"><a name="l00350"></a><span class="lineno"> 350</span>  <span class="keywordflow">if</span> (d < distance)</div>
<div class="line"><a name="l00351"></a><span class="lineno"> 351</span>  {</div>
<div class="line"><a name="l00352"></a><span class="lineno"> 352</span>  dealii::Point<dim> face_intersection =</div>
<div class="line"><a name="l00353"></a><span class="lineno"> 353</span>  ray_origin + d * ray_direction;</div>
<div class="line"><a name="l00354"></a><span class="lineno"> 354</span>  <span class="comment">// Check if the point intersect the first triangle</span></div>
<div class="line"><a name="l00355"></a><span class="lineno"> 355</span>  <span class="keywordflow">if</span> (point_in_triangle(point_0, point_1, point_2,</div>
<div class="line"><a name="l00356"></a><span class="lineno"> 356</span>  face_intersection))</div>
<div class="line"><a name="l00357"></a><span class="lineno"> 357</span>  {</div>
<div class="line"><a name="l00358"></a><span class="lineno"> 358</span>  distance = d;</div>
<div class="line"><a name="l00359"></a><span class="lineno"> 359</span>  intersection = face_intersection;</div>
<div class="line"><a name="l00360"></a><span class="lineno"> 360</span>  }</div>
<div class="line"><a name="l00361"></a><span class="lineno"> 361</span>  <span class="keywordflow">else</span> <span class="keywordflow">if</span> (point_in_triangle(point_1, point_2, point_3,</div>
<div class="line"><a name="l00362"></a><span class="lineno"> 362</span>  face_intersection))</div>
<div class="line"><a name="l00363"></a><span class="lineno"> 363</span>  {</div>
<div class="line"><a name="l00364"></a><span class="lineno"> 364</span>  distance = d;</div>
<div class="line"><a name="l00365"></a><span class="lineno"> 365</span>  intersection = face_intersection;</div>
<div class="line"><a name="l00366"></a><span class="lineno"> 366</span>  }</div>
<div class="line"><a name="l00367"></a><span class="lineno"> 367</span>  }</div>
<div class="line"><a name="l00368"></a><span class="lineno"> 368</span>  }</div>
<div class="line"><a name="l00369"></a><span class="lineno"> 369</span>  }</div>
<div class="line"><a name="l00370"></a><span class="lineno"> 370</span>  }</div>
<div class="line"><a name="l00371"></a><span class="lineno"> 371</span>  <span class="keywordflow">if</span> (distance < std::numeric_limits<double>::max())</div>
<div class="line"><a name="l00372"></a><span class="lineno"> 372</span>  {</div>
<div class="line"><a name="l00373"></a><span class="lineno"> 373</span>  points.push_back(intersection);</div>
<div class="line"><a name="l00374"></a><span class="lineno"> 374</span>  values.push_back(<a class="code" href="classadamantine_1_1_ray_tracing.html#a14cdd076d9145977877d8d595a4b70d4">_values_current_frame</a>[i]);</div>
<div class="line"><a name="l00375"></a><span class="lineno"> 375</span>  }</div>
<div class="line"><a name="l00376"></a><span class="lineno"> 376</span>  }</div>
<div class="line"><a name="l00377"></a><span class="lineno"> 377</span>  <span class="comment">// We know that there are at most n_intersections between the rays and the</span></div>
<div class="line"><a name="l00378"></a><span class="lineno"> 378</span>  <span class="comment">// mesh</span></div>
<div class="line"><a name="l00379"></a><span class="lineno"> 379</span>  <a class="code" href="utils_8hh.html#aa06eedd6f738a415870e97a375337d51">ASSERT</a>(points.size() <= n_intersections,</div>
<div class="line"><a name="l00380"></a><span class="lineno"> 380</span>  <span class="stringliteral">"Error when computing the intersection of rays with the mesh."</span>);</div>
<div class="line"><a name="l00381"></a><span class="lineno"> 381</span>  }</div>
<div class="line"><a name="l00382"></a><span class="lineno"> 382</span>  </div>
<div class="line"><a name="l00383"></a><span class="lineno"> 383</span>  <span class="comment">// We need all the processors to know the intersection points of the rays with</span></div>
<div class="line"><a name="l00384"></a><span class="lineno"> 384</span>  <span class="comment">// the mesh and the values at these points.</span></div>
<div class="line"><a name="l00385"></a><span class="lineno"> 385</span>  <a class="code" href="structadamantine_1_1_points_values.html">PointsValues<dim></a> points_values;</div>
<div class="line"><a name="l00386"></a><span class="lineno"> 386</span>  <span class="keyword">auto</span> all_points = dealii::Utilities::MPI::all_gather(communicator, points);</div>
<div class="line"><a name="l00387"></a><span class="lineno"> 387</span>  <span class="keyword">auto</span> all_values = dealii::Utilities::MPI::all_gather(communicator, values);</div>
<div class="line"><a name="l00388"></a><span class="lineno"> 388</span>  <span class="keywordflow">for</span> (<span class="keywordtype">unsigned</span> <span class="keywordtype">int</span> i = 0; i < all_points.size(); ++i)</div>
<div class="line"><a name="l00389"></a><span class="lineno"> 389</span>  {</div>
<div class="line"><a name="l00390"></a><span class="lineno"> 390</span>  points_values.<a class="code" href="structadamantine_1_1_points_values.html#a9938ef730772cc6779992aa1e2caf9e2">points</a>.insert(points_values.<a class="code" href="structadamantine_1_1_points_values.html#a9938ef730772cc6779992aa1e2caf9e2">points</a>.end(),</div>
<div class="line"><a name="l00391"></a><span class="lineno"> 391</span>  all_points[i].begin(), all_points[i].end());</div>
<div class="line"><a name="l00392"></a><span class="lineno"> 392</span>  points_values.<a class="code" href="structadamantine_1_1_points_values.html#a9d7c8d9774a8a4e04d285ef624e42c9f">values</a>.insert(points_values.<a class="code" href="structadamantine_1_1_points_values.html#a9d7c8d9774a8a4e04d285ef624e42c9f">values</a>.end(),</div>
<div class="line"><a name="l00393"></a><span class="lineno"> 393</span>  all_values[i].begin(), all_values[i].end());</div>
<div class="line"><a name="l00394"></a><span class="lineno"> 394</span>  }</div>
<div class="line"><a name="l00395"></a><span class="lineno"> 395</span>  </div>
<div class="line"><a name="l00396"></a><span class="lineno"> 396</span>  <span class="keywordflow">return</span> points_values;</div>
<div class="line"><a name="l00397"></a><span class="lineno"> 397</span> }</div>
<div class="line"><a name="l00398"></a><span class="lineno"> 398</span>  </div>
<div class="line"><a name="l00399"></a><span class="lineno"> 399</span> } <span class="comment">// namespace adamantine</span></div>
<div class="ttc" id="a_ray_tracing_8hh_html"><div class="ttname"><a href="_ray_tracing_8hh.html">RayTracing.hh</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_nearest_predicate_html"><div class="ttname"><a href="classadamantine_1_1_ray_nearest_predicate.html">adamantine::RayNearestPredicate</a></div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00025">RayTracing.cc:26</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_nearest_predicate_html_a1206cdcc666feb98cacb292ebacb41d5"><div class="ttname"><a href="classadamantine_1_1_ray_nearest_predicate.html#a1206cdcc666feb98cacb292ebacb41d5">adamantine::RayNearestPredicate::_rays</a></div><div class="ttdeci">std::vector< Ray< 3 > > _rays</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00045">RayTracing.cc:45</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_nearest_predicate_html_a2021c9b114a0dd6c22874e305622e062"><div class="ttname"><a href="classadamantine_1_1_ray_nearest_predicate.html#a2021c9b114a0dd6c22874e305622e062">adamantine::RayNearestPredicate::size</a></div><div class="ttdeci">std::size_t size() const</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00037">RayTracing.cc:37</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_nearest_predicate_html_a2d52437d01690e56d8eaabf8354d37a4"><div class="ttname"><a href="classadamantine_1_1_ray_nearest_predicate.html#a2d52437d01690e56d8eaabf8354d37a4">adamantine::RayNearestPredicate::get</a></div><div class="ttdeci">Ray< 3 > const & get(unsigned int i) const</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00042">RayTracing.cc:42</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_nearest_predicate_html_a69a8378532d6cc9bf4bdeaebbc764ff5"><div class="ttname"><a href="classadamantine_1_1_ray_nearest_predicate.html#a69a8378532d6cc9bf4bdeaebbc764ff5">adamantine::RayNearestPredicate::RayNearestPredicate</a></div><div class="ttdeci">RayNearestPredicate(std::vector< Ray< 3 >> const &rays)</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00032">RayTracing.cc:32</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_a11ad68edb125a2f2aca6538a1f4c33f4"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#a11ad68edb125a2f2aca6538a1f4c33f4">adamantine::RayTracing::dim</a></div><div class="ttdeci">static constexpr int dim</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8hh_source.html#l00036">RayTracing.hh:36</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_a14cdd076d9145977877d8d595a4b70d4"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#a14cdd076d9145977877d8d595a4b70d4">adamantine::RayTracing::_values_current_frame</a></div><div class="ttdeci">std::vector< double > _values_current_frame</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8hh_source.html#l00076">RayTracing.hh:76</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_a24a8f437b8e3e1e73f880e52a66fa456"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#a24a8f437b8e3e1e73f880e52a66fa456">adamantine::RayTracing::_data_filename</a></div><div class="ttdeci">std::string _data_filename</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8hh_source.html#l00064">RayTracing.hh:64</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_a311145cc3e9dccf9674939baa8205f0e"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#a311145cc3e9dccf9674939baa8205f0e">adamantine::RayTracing::RayTracing</a></div><div class="ttdeci">RayTracing(boost::property_tree::ptree const &experiment_database, dealii::DoFHandler< dim > const &dof_handler)</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00116">RayTracing.cc:116</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_a3cb1e406035092d2f0825455202c4521"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#a3cb1e406035092d2f0825455202c4521">adamantine::RayTracing::read_next_frame</a></div><div class="ttdeci">unsigned int read_next_frame() override</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00133">RayTracing.cc:133</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_a5e1ffea3480cec7e73cbe0a5c354d72d"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#a5e1ffea3480cec7e73cbe0a5c354d72d">adamantine::RayTracing::_first_camera_id</a></div><div class="ttdeci">unsigned int _first_camera_id</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8hh_source.html#l00056">RayTracing.hh:56</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_a74976dc4d69bf1b3f23c0c25f1824384"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#a74976dc4d69bf1b3f23c0c25f1824384">adamantine::RayTracing::_dof_handler</a></div><div class="ttdeci">dealii::DoFHandler< dim > const & _dof_handler</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8hh_source.html#l00068">RayTracing.hh:68</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_a7e9becfe37db0db6a01ca884f52b3335"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#a7e9becfe37db0db6a01ca884f52b3335">adamantine::RayTracing::get_points_values</a></div><div class="ttdeci">PointsValues< dim > get_points_values() override</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00206">RayTracing.cc:206</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_aa58ad004836a0c29a24bde18d3ba4947"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#aa58ad004836a0c29a24bde18d3ba4947">adamantine::RayTracing::_next_frame</a></div><div class="ttdeci">unsigned int _next_frame</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8hh_source.html#l00052">RayTracing.hh:52</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_ac6abef2dbdc651660422877350428b20"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#ac6abef2dbdc651660422877350428b20">adamantine::RayTracing::_last_camera_id</a></div><div class="ttdeci">unsigned int _last_camera_id</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8hh_source.html#l00060">RayTracing.hh:60</a></div></div>
<div class="ttc" id="aclassadamantine_1_1_ray_tracing_html_af24254cda95c0dc2084f74411a7661fe"><div class="ttname"><a href="classadamantine_1_1_ray_tracing.html#af24254cda95c0dc2084f74411a7661fe">adamantine::RayTracing::_rays_current_frame</a></div><div class="ttdeci">std::vector< Ray< dim > > _rays_current_frame</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8hh_source.html#l00072">RayTracing.hh:72</a></div></div>
<div class="ttc" id="anamespace_arbor_x_html"><div class="ttname"><a href="namespace_arbor_x.html">ArborX</a></div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00049">RayTracing.cc:50</a></div></div>
<div class="ttc" id="anamespaceadamantine_html"><div class="ttname"><a href="namespaceadamantine.html">adamantine</a></div><div class="ttdef"><b>Definition:</b> <a href="_beam_heat_source_properties_8hh_source.html#l00014">BeamHeatSourceProperties.hh:15</a></div></div>
<div class="ttc" id="anamespaceadamantine_html_a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474"><div class="ttname"><a href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a6438c669e0d0de98e6929c2cc0fac474">adamantine::ScanPathSegmentType::line</a></div><div class="ttdeci">@ line</div></div>
<div class="ttc" id="anamespaceadamantine_html_a15248f4195a1e0e7e349886de9657265a78ee54aa8f813885fe2fe20d232518b9"><div class="ttname"><a href="namespaceadamantine.html#a15248f4195a1e0e7e349886de9657265a78ee54aa8f813885fe2fe20d232518b9">adamantine::ScanPathSegmentType::point</a></div><div class="ttdeci">@ point</div></div>
<div class="ttc" id="anamespaceadamantine_html_a737967bbc854d35e1281ad917ba24814"><div class="ttname"><a href="namespaceadamantine.html#a737967bbc854d35e1281ad917ba24814">adamantine::wait_for_file</a></div><div class="ttdeci">void wait_for_file(std::string const &filename, std::string const &message)</div><div class="ttdef"><b>Definition:</b> <a href="utils_8hh_source.html#l00024">utils.hh:24</a></div></div>
<div class="ttc" id="astruct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4_html_a4e1b6669289c70e41ddf3a71092362da"><div class="ttname"><a href="struct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4.html#a4e1b6669289c70e41ddf3a71092362da">ArborX::AccessTraits< adamantine::RayNearestPredicate, PredicatesTag >::memory_space</a></div><div class="ttdeci">Kokkos::HostSpace memory_space</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00054">RayTracing.cc:54</a></div></div>
<div class="ttc" id="astruct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4_html_a66932a60bd94e9171a4057e8dcbddd06"><div class="ttname"><a href="struct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4.html#a66932a60bd94e9171a4057e8dcbddd06">ArborX::AccessTraits< adamantine::RayNearestPredicate, PredicatesTag >::size</a></div><div class="ttdeci">static std::size_t size(adamantine::RayNearestPredicate const &ray_nearest)</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00056">RayTracing.cc:56</a></div></div>
<div class="ttc" id="astruct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4_html_ac0ee537fc30b13dfde371f038ff44f63"><div class="ttname"><a href="struct_arbor_x_1_1_access_traits_3_01adamantine_1_1_ray_nearest_predicate_00_01_predicates_tag_01_4.html#ac0ee537fc30b13dfde371f038ff44f63">ArborX::AccessTraits< adamantine::RayNearestPredicate, PredicatesTag >::get</a></div><div class="ttdeci">static auto get(adamantine::RayNearestPredicate const &ray_nearest, std::size_t i)</div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8cc_source.html#l00061">RayTracing.cc:61</a></div></div>
<div class="ttc" id="astructadamantine_1_1_points_values_html"><div class="ttname"><a href="structadamantine_1_1_points_values.html">adamantine::PointsValues</a></div><div class="ttdef"><b>Definition:</b> <a href="experimental__data__utils_8hh_source.html#l00021">experimental_data_utils.hh:22</a></div></div>
<div class="ttc" id="astructadamantine_1_1_points_values_html_a9938ef730772cc6779992aa1e2caf9e2"><div class="ttname"><a href="structadamantine_1_1_points_values.html#a9938ef730772cc6779992aa1e2caf9e2">adamantine::PointsValues::points</a></div><div class="ttdeci">std::vector< dealii::Point< dim > > points</div><div class="ttdef"><b>Definition:</b> <a href="experimental__data__utils_8hh_source.html#l00023">experimental_data_utils.hh:23</a></div></div>
<div class="ttc" id="astructadamantine_1_1_points_values_html_a9d7c8d9774a8a4e04d285ef624e42c9f"><div class="ttname"><a href="structadamantine_1_1_points_values.html#a9d7c8d9774a8a4e04d285ef624e42c9f">adamantine::PointsValues::values</a></div><div class="ttdeci">std::vector< double > values</div><div class="ttdef"><b>Definition:</b> <a href="experimental__data__utils_8hh_source.html#l00024">experimental_data_utils.hh:24</a></div></div>
<div class="ttc" id="astructadamantine_1_1_ray_html"><div class="ttname"><a href="structadamantine_1_1_ray.html">adamantine::Ray</a></div><div class="ttdef"><b>Definition:</b> <a href="_ray_tracing_8hh_source.html#l00018">RayTracing.hh:19</a></div></div>
<div class="ttc" id="autils_8hh_html"><div class="ttname"><a href="utils_8hh.html">utils.hh</a></div></div>
<div class="ttc" id="autils_8hh_html_aa06eedd6f738a415870e97a375337d51"><div class="ttname"><a href="utils_8hh.html#aa06eedd6f738a415870e97a375337d51">ASSERT</a></div><div class="ttdeci">#define ASSERT(condition, message)</div><div class="ttdef"><b>Definition:</b> <a href="utils_8hh_source.html#l00068">utils.hh:68</a></div></div>
</div><!-- fragment --></div><!-- contents -->
<!-- start footer part -->
<hr class="footer"/><address class="footer"><small>
Generated by <a href="https://www.doxygen.org/index.html"><img class="footer" src="doxygen.svg" width="104" height="31" alt="doxygen"/></a> 1.9.1
</small></address>
</body>
</html>