OPALX (Object Oriented Parallel Accelerator Library for Exascal) MINIorX
OPALX
BoundaryGeometry.cpp
Go to the documentation of this file.
1//
2// Declaration of the BoundaryGeometry class
3//
4// Copyright (c) 200x - 2020, Achim Gsell,
5// Paul Scherrer Institut, Villigen PSI, Switzerland
6// All rights reserved.
7//
8// This file is part of OPAL.
9//
10// OPAL is free software: you can redistribute it and/or modify
11// it under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// You should have received a copy of the GNU General Public License
16// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
17//
18
19// #define ENABLE_DEBUG
20
22
23#include <algorithm>
24#include <cmath>
25#include <ctime>
26#include <fstream>
27#include <string>
28
29#include <cfloat>
30#include "H5hut.h"
31
33#include "PartBunch/PartBunch.h"
36#include "Physics/Physics.h"
38#include "Utilities/Options.h"
39
40#include <boost/filesystem.hpp>
41
42#include <gsl/gsl_sys.h>
43
44extern Inform* gmsg;
45
46#define SQR(x) ((x) * (x))
47#define PointID(triangle_id, vertex_id) Triangles_m[triangle_id][vertex_id]
48#define Point(triangle_id, vertex_id) Points_m[Triangles_m[triangle_id][vertex_id]]
49
50/*
51 In the following namespaces various approximately floating point
52 comparisons are implemented. The used implementation is selected
53 via
54
55 namespaces cmp = IMPLEMENTATION;
56
57*/
58
59/*
60 First we define some macros for function common in all namespaces.
61*/
62#define FUNC_EQ(x, y) \
63 inline bool eq(double x, double y) { return almost_eq(x, y); }
64
65#define FUNC_EQ_ZERO(x) \
66 inline bool eq_zero(double x) { return almost_eq_zero(x); }
67
68#define FUNC_LE(x, y) \
69 inline bool le(double x, double y) { \
70 if (almost_eq(x, y)) { \
71 return true; \
72 } \
73 return x < y; \
74 }
75
76#define FUNC_LE_ZERO(x) \
77 inline bool le_zero(double x) { \
78 if (almost_eq_zero(x)) { \
79 return true; \
80 } \
81 return x < 0.0; \
82 }
83
84#define FUNC_LT(x, y) \
85 inline bool lt(double x, double y) { \
86 if (almost_eq(x, y)) { \
87 return false; \
88 } \
89 return x < y; \
90 }
91
92#define FUNC_LT_ZERO(x) \
93 inline bool lt_zero(double x) { \
94 if (almost_eq_zero(x)) { \
95 return false; \
96 } \
97 return x < 0.0; \
98 }
99
100#define FUNC_GE(x, y) \
101 inline bool ge(double x, double y) { \
102 if (almost_eq(x, y)) { \
103 return true; \
104 } \
105 return x > y; \
106 }
107
108#define FUNC_GE_ZERO(x) \
109 inline bool ge_zero(double x) { \
110 if (almost_eq_zero(x)) { \
111 return true; \
112 } \
113 return x > 0.0; \
114 }
115
116#define FUNC_GT(x, y) \
117 inline bool gt(double x, double y) { \
118 if (almost_eq(x, y)) { \
119 return false; \
120 } \
121 return x > y; \
122 }
123
124#define FUNC_GT_ZERO(x) \
125 inline bool gt_zero(double x) { \
126 if (almost_eq_zero(x)) { \
127 return false; \
128 } \
129 return x > 0.0; \
130 }
131
132namespace cmp_diff {
133
134 /*
135 Link:
136 https://randomascii.wordpress.com/2012/02/25/comparing-floating-point-numbers-2012-edition/
137 */
138 inline bool almost_eq(
139 double A, double B, double maxDiff = 1e-15, double maxRelDiff = DBL_EPSILON) {
140 // Check if the numbers are really close -- needed
141 // when comparing numbers near zero.
142 const double diff = std::abs(A - B);
143 if (diff <= maxDiff)
144 return true;
145
146 A = std::abs(A);
147 B = std::abs(B);
148 const double largest = (B > A) ? B : A;
149
150 if (diff <= largest * maxRelDiff)
151 return true;
152 return false;
153 }
154
155 inline bool almost_eq_zero(double A, double maxDiff = 1e-15) {
156 const double diff = std::abs(A);
157 return (diff <= maxDiff);
158 }
159
160 FUNC_EQ(x, y);
162 FUNC_LE(x, y);
164 FUNC_LT(x, y);
166 FUNC_GE(x, y);
168 FUNC_GT(x, y);
170} // namespace cmp_diff
171
173 /*
174 See:
175 https://www.cygnus-software.com/papers/comparingfloats/comparing_floating_point_numbers_obsolete.htm
176 */
177 inline bool almost_eq(double A, double B, double maxDiff = 1e-20, int maxUlps = 1000) {
178 // Make sure maxUlps is non-negative and small enough that the
179 // default NAN won't compare as equal to anything.
180 // assert(maxUlps > 0 && maxUlps < 4 * 1024 * 1024);
181
182 // handle NaN's
183 // Note: comparing something with a NaN is always false!
184 if (std::isnan(A) || std::isnan(B)) {
185 return false;
186 }
187
188 if (std::abs(A - B) <= maxDiff) {
189 return true;
190 }
191
192#pragma GCC diagnostic push
193#pragma GCC diagnostic ignored "-Wstrict-aliasing"
194 auto aInt = *(int64_t*)&A;
195#pragma GCC diagnostic pop
196 // Make aInt lexicographically ordered as a twos-complement int
197 if (aInt < 0) {
198 aInt = 0x8000000000000000 - aInt;
199 }
200
201#pragma GCC diagnostic push
202#pragma GCC diagnostic ignored "-Wstrict-aliasing"
203 auto bInt = *(int64_t*)&B;
204#pragma GCC diagnostic pop
205 // Make bInt lexicographically ordered as a twos-complement int
206 if (bInt < 0) {
207 bInt = 0x8000000000000000 - bInt;
208 }
209
210 if (std::abs(aInt - bInt) <= maxUlps) {
211 return true;
212 }
213 return false;
214 }
215
216 inline bool almost_eq_zero(double A, double maxDiff = 1e-15) {
217 // no need to handle NaN's!
218 return (std::abs(A) <= maxDiff);
219 }
220 FUNC_EQ(x, y);
222 FUNC_LE(x, y);
224 FUNC_LT(x, y);
226 FUNC_GE(x, y);
228 FUNC_GT(x, y);
230} // namespace cmp_ulp_obsolete
231
232namespace cmp_ulp {
233 /*
234 See:
235 https://randomascii.wordpress.com/2012/02/25/comparing-floating-point-numbers-2012-edition/
236 */
237
238 inline bool almost_eq(double A, double B, double maxDiff = 1e-20, int maxUlps = 1000) {
239 // handle NaN's
240 if (std::isnan(A) || std::isnan(B)) {
241 return false;
242 }
243
244 // Check if the numbers are really close -- needed
245 // when comparing numbers near zero.
246 if (std::abs(A - B) <= maxDiff)
247 return true;
248
249#pragma GCC diagnostic push
250#pragma GCC diagnostic ignored "-Wstrict-aliasing"
251 auto aInt = *(int64_t*)&A;
252 auto bInt = *(int64_t*)&B;
253#pragma GCC diagnostic pop
254
255 // Different signs means they do not match.
256 // Note: a negative floating point number is also negative as integer.
257 if (std::signbit(aInt) != std::signbit(bInt))
258 return false;
259
260 // Find the difference in ULPs.
261 return (std::abs(aInt - bInt) <= maxUlps);
262 }
263
264 inline bool almost_eq_zero(double A, double maxDiff = 1e-15) {
265 return (std::abs(A) <= maxDiff);
266 }
267 FUNC_EQ(x, y);
269 FUNC_LE(x, y);
271 FUNC_LT(x, y);
273 FUNC_GE(x, y);
275 FUNC_GT(x, y);
277} // namespace cmp_ulp
278
279namespace cmp = cmp_ulp;
280/*
281
282 Some
283 _ _ _
284 | | | | ___| |_ __ ___ _ __
285 | |_| |/ _ \ | '_ \ / _ \ '__|
286 | _ | __/ | |_) | __/ |
287 |_| |_|\___|_| .__/ \___|_|
288 |_|
289
290 functions
291 */
292namespace {
293 struct VectorLessX {
294 bool operator()(Vector_t<double, 3> x1, Vector_t<double, 3> x2) {
295 return cmp::lt(x1(0), x2(0));
296 }
297 };
298
299 struct VectorLessY {
300 bool operator()(Vector_t<double, 3> x1, Vector_t<double, 3> x2) {
301 return cmp::lt(x1(1), x2(1));
302 }
303 };
304
305 struct VectorLessZ {
306 bool operator()(Vector_t<double, 3> x1, Vector_t<double, 3> x2) {
307 return cmp::lt(x1(2), x2(2));
308 }
309 };
310
314 Vector_t<double, 3> get_max_extent(std::vector<Vector_t<double, 3>>& coords) {
315 const Vector_t<double, 3> x = *max_element(coords.begin(), coords.end(), VectorLessX());
316 const Vector_t<double, 3> y = *max_element(coords.begin(), coords.end(), VectorLessY());
317 const Vector_t<double, 3> z = *max_element(coords.begin(), coords.end(), VectorLessZ());
318 return Vector_t<double, 3>(x(0), y(1), z(2));
319 }
320
321 /*
322 Compute the minimum of coordinates of geometry, i.e the minimum of X,Y,Z
323 */
324 Vector_t<double, 3> get_min_extent(std::vector<Vector_t<double, 3>>& coords) {
325 const Vector_t<double, 3> x = *min_element(coords.begin(), coords.end(), VectorLessX());
326 const Vector_t<double, 3> y = *min_element(coords.begin(), coords.end(), VectorLessY());
327 const Vector_t<double, 3> z = *min_element(coords.begin(), coords.end(), VectorLessZ());
328 return Vector_t<double, 3>(x(0), y(1), z(2));
329 }
330
331 /*
332 write legacy VTK file of voxel mesh
333 */
334 static void write_voxel_mesh(
335 std::string fname, const std::unordered_map<int, std::unordered_set<int>>& ids,
336 const Vector_t<double, 3>& hr_m, const Vector_t<int, 3>& nr,
337 const Vector_t<double, 3>& origin) {
338 /*----------------------------------------------------------------------*/
339 const size_t numpoints = 8 * ids.size();
340 std::ofstream of;
341
342 *gmsg << level2 << "* Writing VTK file of voxel mesh '" << fname << "'" << endl;
343 of.open(fname);
344 PAssert(of.is_open());
345 of.precision(6);
346
347 of << "# vtk DataFile Version 2.0" << std::endl;
348 of << "generated using BoundaryGeometry::computeMeshVoxelization" << std::endl;
349 of << "ASCII" << std::endl << std::endl;
350 of << "DATASET UNSTRUCTURED_GRID" << std::endl;
351 of << "POINTS " << numpoints << " float" << std::endl;
352
353 const auto nr0_times_nr1 = nr[0] * nr[1];
354 for (auto& elem : ids) {
355 auto id = elem.first;
356 int k = (id - 1) / nr0_times_nr1;
357 int rest = (id - 1) % nr0_times_nr1;
358 int j = rest / nr[0];
359 int i = rest % nr[0];
360
362 P[0] = i * hr_m[0] + origin[0];
363 P[1] = j * hr_m[1] + origin[1];
364 P[2] = k * hr_m[2] + origin[2];
365
366 of << P[0] << " " << P[1] << " " << P[2] << std::endl;
367 of << P[0] + hr_m[0] << " " << P[1] << " " << P[2] << std::endl;
368 of << P[0] << " " << P[1] + hr_m[1] << " " << P[2] << std::endl;
369 of << P[0] + hr_m[0] << " " << P[1] + hr_m[1] << " " << P[2] << std::endl;
370 of << P[0] << " " << P[1] << " " << P[2] + hr_m[2] << std::endl;
371 of << P[0] + hr_m[0] << " " << P[1] << " " << P[2] + hr_m[2] << std::endl;
372 of << P[0] << " " << P[1] + hr_m[1] << " " << P[2] + hr_m[2] << std::endl;
373 of << P[0] + hr_m[0] << " " << P[1] + hr_m[1] << " " << P[2] + hr_m[2] << std::endl;
374 }
375 of << std::endl;
376 const auto num_cells = ids.size();
377 of << "CELLS " << num_cells << " " << 9 * num_cells << std::endl;
378 for (size_t i = 0; i < num_cells; i++)
379 of << "8 " << 8 * i << " " << 8 * i + 1 << " " << 8 * i + 2 << " " << 8 * i + 3 << " "
380 << 8 * i + 4 << " " << 8 * i + 5 << " " << 8 * i + 6 << " " << 8 * i + 7
381 << std::endl;
382 of << "CELL_TYPES " << num_cells << std::endl;
383 for (size_t i = 0; i < num_cells; i++)
384 of << "11" << std::endl;
385 of << "CELL_DATA " << num_cells << std::endl;
386 of << "SCALARS "
387 << "cell_attribute_data"
388 << " float "
389 << "1" << std::endl;
390 of << "LOOKUP_TABLE "
391 << "default" << std::endl;
392 for (size_t i = 0; i < num_cells; i++)
393 of << (float)(i) << std::endl;
394 of << std::endl;
395 of << "COLOR_SCALARS "
396 << "BBoxColor " << 4 << std::endl;
397 for (size_t i = 0; i < num_cells; i++) {
398 of << "1.0"
399 << " 1.0 "
400 << "0.0 "
401 << "1.0" << std::endl;
402 }
403 of << std::endl;
404 }
405} // namespace
406
407/*___________________________________________________________________________
408
409 Triangle-cube intersection test.
410
411 See:
412 http://tog.acm.org/resources/GraphicsGems/gemsiii/triangleCube.c
413
414 */
415
416#define INSIDE 0
417#define OUTSIDE 1
418
419class Triangle {
420public:
422 }
425 const Vector_t<double, 3>& v3) {
426 pts[0] = v1;
427 pts[1] = v2;
428 pts[2] = v3;
429 }
430
431 inline const Vector_t<double, 3>& v1() const {
432 return pts[0];
433 }
434 inline double v1(int i) const {
435 return pts[0][i];
436 }
437 inline const Vector_t<double, 3>& v2() const {
438 return pts[1];
439 }
440 inline double v2(int i) const {
441 return pts[1][i];
442 }
443 inline const Vector_t<double, 3>& v3() const {
444 return pts[2];
445 }
446 inline double v3(int i) const {
447 return pts[2][i];
448 }
449
450 inline void scale(const Vector_t<double, 3>& scaleby, const Vector_t<double, 3>& shiftby) {
451 pts[0][0] *= scaleby[0];
452 pts[0][1] *= scaleby[1];
453 pts[0][2] *= scaleby[2];
454 pts[1][0] *= scaleby[0];
455 pts[1][1] *= scaleby[1];
456 pts[1][2] *= scaleby[2];
457 pts[2][0] *= scaleby[0];
458 pts[2][1] *= scaleby[1];
459 pts[2][2] *= scaleby[2];
460 pts[0] -= shiftby;
461 pts[1] -= shiftby;
462 pts[2] -= shiftby;
463 }
464
466};
467
468/*___________________________________________________________________________*/
469
470/* Which of the six face-plane(s) is point P outside of? */
471
472static inline int face_plane(const Vector_t<double, 3>& p) {
473 int outcode_fcmp = 0;
474
475 if (cmp::gt(p[0], 0.5))
476 outcode_fcmp |= 0x01;
477 if (cmp::lt(p[0], -0.5))
478 outcode_fcmp |= 0x02;
479 if (cmp::gt(p[1], 0.5))
480 outcode_fcmp |= 0x04;
481 if (cmp::lt(p[1], -0.5))
482 outcode_fcmp |= 0x08;
483 if (cmp::gt(p[2], 0.5))
484 outcode_fcmp |= 0x10;
485 if (cmp::lt(p[2], -0.5))
486 outcode_fcmp |= 0x20;
487
488 return (outcode_fcmp);
489}
490
491/*. . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . */
492
493/* Which of the twelve edge plane(s) is point P outside of? */
494
495static inline int bevel_2d(const Vector_t<double, 3>& p) {
496 int outcode_fcmp = 0;
497
498 if (cmp::gt(p[0] + p[1], 1.0))
499 outcode_fcmp |= 0x001;
500 if (cmp::gt(p[0] - p[1], 1.0))
501 outcode_fcmp |= 0x002;
502 if (cmp::gt(-p[0] + p[1], 1.0))
503 outcode_fcmp |= 0x004;
504 if (cmp::gt(-p[0] - p[1], 1.0))
505 outcode_fcmp |= 0x008;
506 if (cmp::gt(p[0] + p[2], 1.0))
507 outcode_fcmp |= 0x010;
508 if (cmp::gt(p[0] - p[2], 1.0))
509 outcode_fcmp |= 0x020;
510 if (cmp::gt(-p[0] + p[2], 1.0))
511 outcode_fcmp |= 0x040;
512 if (cmp::gt(-p[0] - p[2], 1.0))
513 outcode_fcmp |= 0x080;
514 if (cmp::gt(p[1] + p[2], 1.0))
515 outcode_fcmp |= 0x100;
516 if (cmp::gt(p[1] - p[2], 1.0))
517 outcode_fcmp |= 0x200;
518 if (cmp::gt(-p[1] + p[2], 1.0))
519 outcode_fcmp |= 0x400;
520 if (cmp::gt(-p[1] - p[2], 1.0))
521 outcode_fcmp |= 0x800;
522
523 return (outcode_fcmp);
524}
525
526/*. . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
527
528 Which of the eight corner plane(s) is point P outside of?
529*/
530static inline int bevel_3d(const Vector_t<double, 3>& p) {
531 int outcode_fcmp = 0;
532
533 if (cmp::gt(p[0] + p[1] + p[2], 1.5))
534 outcode_fcmp |= 0x01;
535 if (cmp::gt(p[0] + p[1] - p[2], 1.5))
536 outcode_fcmp |= 0x02;
537 if (cmp::gt(p[0] - p[1] + p[2], 1.5))
538 outcode_fcmp |= 0x04;
539 if (cmp::gt(p[0] - p[1] - p[2], 1.5))
540 outcode_fcmp |= 0x08;
541 if (cmp::gt(-p[0] + p[1] + p[2], 1.5))
542 outcode_fcmp |= 0x10;
543 if (cmp::gt(-p[0] + p[1] - p[2], 1.5))
544 outcode_fcmp |= 0x20;
545 if (cmp::gt(-p[0] - p[1] + p[2], 1.5))
546 outcode_fcmp |= 0x40;
547 if (cmp::gt(-p[0] - p[1] - p[2], 1.5))
548 outcode_fcmp |= 0x80;
549
550 return (outcode_fcmp);
551}
552
553/*. . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
554
555 Test the point "alpha" of the way from P1 to P2
556 See if it is on a face of the cube
557 Consider only faces in "mask"
558*/
559
560static inline int check_point(
561 const Vector_t<double, 3>& p1, const Vector_t<double, 3>& p2, const double alpha,
562 const int mask) {
563 Vector_t<double, 3> plane_point;
564
565#define LERP(a, b, t) (a + t * (b - a))
566 // with C++20 we can use: std::lerp(a, b, t)
567 plane_point[0] = LERP(p1[0], p2[0], alpha);
568 plane_point[1] = LERP(p1[1], p2[1], alpha);
569 plane_point[2] = LERP(p1[2], p2[2], alpha);
570#undef LERP
571 return (face_plane(plane_point) & mask);
572}
573
574/*. . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
575
576 Compute intersection of P1 --> P2 line segment with face planes
577 Then test intersection point to see if it is on cube face
578 Consider only face planes in "outcode_diff"
579 Note: Zero bits in "outcode_diff" means face line is outside of
580*/
581static inline int check_line(
582 const Vector_t<double, 3>& p1, const Vector_t<double, 3>& p2, const int outcode_diff) {
583 if ((0x01 & outcode_diff) != 0)
584 if (check_point(p1, p2, (.5 - p1[0]) / (p2[0] - p1[0]), 0x3e) == INSIDE)
585 return (INSIDE);
586 if ((0x02 & outcode_diff) != 0)
587 if (check_point(p1, p2, (-.5 - p1[0]) / (p2[0] - p1[0]), 0x3d) == INSIDE)
588 return (INSIDE);
589 if ((0x04 & outcode_diff) != 0)
590 if (check_point(p1, p2, (.5 - p1[1]) / (p2[1] - p1[1]), 0x3b) == INSIDE)
591 return (INSIDE);
592 if ((0x08 & outcode_diff) != 0)
593 if (check_point(p1, p2, (-.5 - p1[1]) / (p2[1] - p1[1]), 0x37) == INSIDE)
594 return (INSIDE);
595 if ((0x10 & outcode_diff) != 0)
596 if (check_point(p1, p2, (.5 - p1[2]) / (p2[2] - p1[2]), 0x2f) == INSIDE)
597 return (INSIDE);
598 if ((0x20 & outcode_diff) != 0)
599 if (check_point(p1, p2, (-.5 - p1[2]) / (p2[2] - p1[2]), 0x1f) == INSIDE)
600 return (INSIDE);
601 return (OUTSIDE);
602}
603
604/*. . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
605
606 Test if 3D point is inside 3D triangle
607*/
608constexpr double EPS = 10e-15;
609static inline int SIGN3(Vector_t<double, 3> A) {
610 return (
611 ((A[0] < EPS) ? 4 : 0) | ((A[0] > -EPS) ? 32 : 0) | ((A[1] < EPS) ? 2 : 0)
612 | ((A[1] > -EPS) ? 16 : 0) | ((A[2] < EPS) ? 1 : 0) | ((A[2] > -EPS) ? 8 : 0));
613}
614
615static int point_triangle_intersection(const Vector_t<double, 3>& p, const Triangle& t) {
616 /*
617 First, a quick bounding-box test:
618 If P is outside triangle bbox, there cannot be an intersection.
619 */
620 if (cmp::gt(p[0], std::max({t.v1(0), t.v2(0), t.v3(0)})))
621 return (OUTSIDE);
622 if (cmp::gt(p[1], std::max({t.v1(1), t.v2(1), t.v3(1)})))
623 return (OUTSIDE);
624 if (cmp::gt(p[2], std::max({t.v1(2), t.v2(2), t.v3(2)})))
625 return (OUTSIDE);
626 if (cmp::lt(p[0], std::min({t.v1(0), t.v2(0), t.v3(0)})))
627 return (OUTSIDE);
628 if (cmp::lt(p[1], std::min({t.v1(1), t.v2(1), t.v3(1)})))
629 return (OUTSIDE);
630 if (cmp::lt(p[2], std::min({t.v1(2), t.v2(2), t.v3(2)})))
631 return (OUTSIDE);
632
633 /*
634 For each triangle side, make a vector out of it by subtracting vertexes;
635 make another vector from one vertex to point P.
636 The crossproduct of these two vectors is orthogonal to both and the
637 signs of its X,Y,Z components indicate whether P was to the inside or
638 to the outside of this triangle side.
639 */
640 const Vector_t<double, 3> vect12 = t.v1() - t.v2();
641 const Vector_t<double, 3> vect1h = t.v1() - p;
642 const Vector_t<double, 3> cross12_1p = cross(vect12, vect1h);
643 const int sign12 = SIGN3(cross12_1p); /* Extract X,Y,Z signs as 0..7 or 0...63 integer */
644
645 const Vector_t<double, 3> vect23 = t.v2() - t.v3();
646 const Vector_t<double, 3> vect2h = t.v2() - p;
647 const Vector_t<double, 3> cross23_2p = cross(vect23, vect2h);
648 const int sign23 = SIGN3(cross23_2p);
649
650 const Vector_t<double, 3> vect31 = t.v3() - t.v1();
651 const Vector_t<double, 3> vect3h = t.v3() - p;
652 const Vector_t<double, 3> cross31_3p = cross(vect31, vect3h);
653 const int sign31 = SIGN3(cross31_3p);
654
655 /*
656 If all three crossproduct vectors agree in their component signs,
657 then the point must be inside all three.
658 P cannot be OUTSIDE all three sides simultaneously.
659 */
660 return ((sign12 & sign23 & sign31) == 0) ? OUTSIDE : INSIDE;
661}
662
663/*. . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
664
665 This is the main algorithm procedure.
666 Triangle t is compared with a unit cube,
667 centered on the origin.
668 It returns INSIDE (0) or OUTSIDE(1) if t
669 intersects or does not intersect the cube.
670*/
671static int triangle_intersects_cube(const Triangle& t) {
672 int v1_test;
673 int v2_test;
674 int v3_test;
675
676 /*
677 First compare all three vertexes with all six face-planes
678 If any vertex is inside the cube, return immediately!
679 */
680 if ((v1_test = face_plane(t.v1())) == INSIDE)
681 return (INSIDE);
682 if ((v2_test = face_plane(t.v2())) == INSIDE)
683 return (INSIDE);
684 if ((v3_test = face_plane(t.v3())) == INSIDE)
685 return (INSIDE);
686
687 /*
688 If all three vertexes were outside of one or more face-planes,
689 return immediately with a trivial rejection!
690 */
691 if ((v1_test & v2_test & v3_test) != 0)
692 return (OUTSIDE);
693
694 /*
695 Now do the same trivial rejection test for the 12 edge planes
696 */
697 v1_test |= bevel_2d(t.v1()) << 8;
698 v2_test |= bevel_2d(t.v2()) << 8;
699 v3_test |= bevel_2d(t.v3()) << 8;
700 if ((v1_test & v2_test & v3_test) != 0)
701 return (OUTSIDE);
702
703 /*
704 Now do the same trivial rejection test for the 8 corner planes
705 */
706 v1_test |= bevel_3d(t.v1()) << 24;
707 v2_test |= bevel_3d(t.v2()) << 24;
708 v3_test |= bevel_3d(t.v3()) << 24;
709 if ((v1_test & v2_test & v3_test) != 0)
710 return (OUTSIDE);
711
712 /*
713 If vertex 1 and 2, as a pair, cannot be trivially rejected
714 by the above tests, then see if the v1-->v2 triangle edge
715 intersects the cube. Do the same for v1-->v3 and v2-->v3./
716 Pass to the intersection algorithm the "OR" of the outcode
717 bits, so that only those cube faces which are spanned by
718 each triangle edge need be tested.
719 */
720 if ((v1_test & v2_test) == 0)
721 if (check_line(t.v1(), t.v2(), v1_test | v2_test) == INSIDE)
722 return (INSIDE);
723 if ((v1_test & v3_test) == 0)
724 if (check_line(t.v1(), t.v3(), v1_test | v3_test) == INSIDE)
725 return (INSIDE);
726 if ((v2_test & v3_test) == 0)
727 if (check_line(t.v2(), t.v3(), v2_test | v3_test) == INSIDE)
728 return (INSIDE);
729
730 /*
731 By now, we know that the triangle is not off to any side,
732 and that its sides do not penetrate the cube. We must now
733 test for the cube intersecting the interior of the triangle.
734 We do this by looking for intersections between the cube
735 diagonals and the triangle...first finding the intersection
736 of the four diagonals with the plane of the triangle, and
737 then if that intersection is inside the cube, pursuing
738 whether the intersection point is inside the triangle itself.
739
740 To find plane of the triangle, first perform crossproduct on
741 two triangle side vectors to compute the normal vector.
742 */
743 Vector_t<double, 3> vect12 = t.v1() - t.v2();
744 Vector_t<double, 3> vect13 = t.v1() - t.v3();
745 Vector_t<double, 3> norm = cross(vect12, vect13);
746
747 /*
748 The normal vector "norm" X,Y,Z components are the coefficients
749 of the triangles AX + BY + CZ + D = 0 plane equation. If we
750 solve the plane equation for X=Y=Z (a diagonal), we get
751 -D/(A+B+C) as a metric of the distance from cube center to the
752 diagonal/plane intersection. If this is between -0.5 and 0.5,
753 the intersection is inside the cube. If so, we continue by
754 doing a point/triangle intersection.
755 Do this for all four diagonals.
756 */
757 double d = norm[0] * t.v1(0) + norm[1] * t.v1(1) + norm[2] * t.v1(2);
758
759 /*
760 if one of the diagonals is parallel to the plane, the other will
761 intersect the plane
762 */
763 double denom = norm[0] + norm[1] + norm[2];
764 if (cmp::eq_zero(std::abs(denom)) == false) {
765 /* skip parallel diagonals to the plane; division by 0 can occure */
766 Vector_t<double, 3> hitpp = d / denom;
767 if (cmp::le(std::abs(hitpp[0]), 0.5))
768 if (point_triangle_intersection(hitpp, t) == INSIDE)
769 return (INSIDE);
770 }
771 denom = norm[0] + norm[1] - norm[2];
772 if (cmp::eq_zero(std::abs(denom)) == false) {
774 hitpn[2] = -(hitpn[0] = hitpn[1] = d / denom);
775 if (cmp::le(std::abs(hitpn[0]), 0.5))
776 if (point_triangle_intersection(hitpn, t) == INSIDE)
777 return (INSIDE);
778 }
779 denom = norm[0] - norm[1] + norm[2];
780 if (cmp::eq_zero(std::abs(denom)) == false) {
782 hitnp[1] = -(hitnp[0] = hitnp[2] = d / denom);
783 if (cmp::le(std::abs(hitnp[0]), 0.5))
784 if (point_triangle_intersection(hitnp, t) == INSIDE)
785 return (INSIDE);
786 }
787 denom = norm[0] - norm[1] - norm[2];
788 if (cmp::eq_zero(std::abs(denom)) == false) {
790 hitnn[1] = hitnn[2] = -(hitnn[0] = d / denom);
791 if (cmp::le(std::abs(hitnn[0]), 0.5))
792 if (point_triangle_intersection(hitnn, t) == INSIDE)
793 return (INSIDE);
794 }
795
796 /*
797 No edge touched the cube; no cube diagonal touched the triangle.
798 We're done...there was no intersection.
799 */
800 return (OUTSIDE);
801}
802
803/*
804 * Ray class, for use with the optimized ray-box intersection test
805 * described in:
806 *
807 * Amy Williams, Steve Barrus, R. Keith Morley, and Peter Shirley
808 * "An Efficient and Robust Ray-Box Intersection Algorithm"
809 * Journal of graphics tools, 10(1):49-54, 2005
810 *
811 */
812
813class Ray {
814public:
815 Ray() {
816 }
818 origin = o;
819 direction = d;
820 inv_direction = Vector_t<double, 3>(1 / d[0], 1 / d[1], 1 / d[2]);
821 sign[0] = (inv_direction[0] < 0);
822 sign[1] = (inv_direction[1] < 0);
823 sign[2] = (inv_direction[2] < 0);
824 }
825 Ray(const Ray& r) {
826 origin = r.origin;
829 sign[0] = r.sign[0];
830 sign[1] = r.sign[1];
831 sign[2] = r.sign[2];
832 }
833 const Ray& operator=(const Ray& a) = delete;
834
838 int sign[3];
839};
840
841/*
842 * Axis-aligned bounding box class, for use with the optimized ray-box
843 * intersection test described in:
844 *
845 * Amy Williams, Steve Barrus, R. Keith Morley, and Peter Shirley
846 * "An Efficient and Robust Ray-Box Intersection Algorithm"
847 * Journal of graphics tools, 10(1):49-54, 2005
848 *
849 */
850
851class Voxel {
852public:
854 }
856 pts[0] = min;
857 pts[1] = max;
858 }
859 inline void scale(const Vector_t<double, 3>& scale) {
860 pts[0][0] *= scale[0];
861 pts[0][1] *= scale[1];
862 pts[0][2] *= scale[2];
863 pts[1][0] *= scale[0];
864 pts[1][1] *= scale[1];
865 pts[1][2] *= scale[2];
866 }
867
868 // (t0, t1) is the interval for valid hits
870 const Ray& r,
871 double& tmin, // tmin and tmax are unchanged, if there is
872 double& tmax // no intersection
873 ) const {
874 double tmin_ = (pts[r.sign[0]][0] - r.origin[0]) * r.inv_direction[0];
875 double tmax_ = (pts[1 - r.sign[0]][0] - r.origin[0]) * r.inv_direction[0];
876 const double tymin = (pts[r.sign[1]][1] - r.origin[1]) * r.inv_direction[1];
877 const double tymax = (pts[1 - r.sign[1]][1] - r.origin[1]) * r.inv_direction[1];
878 if (cmp::gt(tmin_, tymax) || cmp::gt(tymin, tmax_))
879 return 0; // no intersection
880 if (cmp::gt(tymin, tmin_))
881 tmin_ = tymin;
882 if (cmp::lt(tymax, tmax_))
883 tmax_ = tymax;
884 const double tzmin = (pts[r.sign[2]][2] - r.origin[2]) * r.inv_direction[2];
885 const double tzmax = (pts[1 - r.sign[2]][2] - r.origin[2]) * r.inv_direction[2];
886 if (cmp::gt(tmin_, tzmax) || cmp::gt(tzmin, tmax_))
887 return 0; // no intersection
888 if (cmp::gt(tzmin, tmin_))
889 tmin_ = tzmin;
890 tmin = tmin_;
891 if (cmp::lt(tzmax, tmax_))
892 tmax_ = tzmax;
893 tmax = tmax_;
894 return cmp::ge_zero(tmax);
895 }
896
897 inline bool intersect(const Ray& r) const {
898 double tmin = 0.0;
899 double tmax = 0.0;
900 return intersect(r, tmin, tmax);
901 }
902
903 inline int intersect(const Triangle& t) const {
904 Voxel v_ = *this;
905 Triangle t_ = t;
906 const Vector_t<double, 3> scaleby = 1.0 / v_.extent();
907 v_.scale(scaleby);
908 t_.scale(scaleby, v_.pts[0] + 0.5);
909 return triangle_intersects_cube(t_);
910 }
911
913 return (pts[1] - pts[0]);
914 }
915
916 inline bool isInside(const Vector_t<double, 3>& P) const {
917 return (
918 cmp::ge(P[0], pts[0][0]) && cmp::ge(P[1], pts[0][1]) && cmp::ge(P[2], pts[0][2])
919 && cmp::le(P[0], pts[1][0]) && cmp::le(P[1], pts[1][1]) && cmp::le(P[2], pts[1][2]));
920 }
921
923};
924
925static inline Vector_t<double, 3> normalVector(
926 const Vector_t<double, 3>& A, const Vector_t<double, 3>& B, const Vector_t<double, 3>& C) {
927 const Vector_t<double, 3> N = cross(B - A, C - A);
928 const double magnitude = std::sqrt(SQR(N(0)) + SQR(N(1)) + SQR(N(2)));
929 PAssert(cmp::gt_zero(magnitude)); // in case we have degenerated triangles
930 return N / magnitude;
931}
932
933// Calculate the area of triangle given by id.
934static inline double computeArea(
935 const Vector_t<double, 3>& A, const Vector_t<double, 3>& B, const Vector_t<double, 3>& C) {
936 const Vector_t<double, 3> AB = A - B;
937 const Vector_t<double, 3> AC = C - A;
938 return (0.5 * std::sqrt(dot(AB, AB) * dot(AC, AC) - dot(AB, AC) * dot(AB, AC)));
939}
940
941/*
942 ____ _
943 / ___| ___ ___ _ __ ___ ___| |_ _ __ _ _
944| | _ / _ \/ _ \| '_ ` _ \ / _ \ __| '__| | | |
945| |_| | __/ (_) | | | | | | __/ |_| | | |_| |
946 \____|\___|\___/|_| |_| |_|\___|\__|_| \__, |
947 |___/
948*/
949
951 : Definition(SIZE, "GEOMETRY", "The \"GEOMETRY\" statement defines the beam pipe geometry.") {
952 itsAttr[FGEOM] = Attributes::makeString("FGEOM", "Specifies the geometry file [H5hut]", "");
953
955 "TOPO", "If FGEOM is selected topo is over-written. ",
956 {"RECTANGULAR", "BOXCORNER", "ELLIPTIC"}, "ELLIPTIC");
957
959 "LENGTH", "Specifies the length of a tube shaped elliptic beam pipe [m]", 1.0);
960
962 "S", "Specifies the start of a tube shaped elliptic beam pipe [m]", 0.0);
963
965 "A", "Specifies the major semi-axis of a tube shaped elliptic beam pipe [m]", 0.025);
966
968 "B", "Specifies the major semi-axis of a tube shaped elliptic beam pipe [m]", 0.025);
969
971 "L1", "In case of BOXCORNER Specifies first part with height == B [m]", 0.5);
972
974 "L2", "In case of BOXCORNER Specifies first second with height == B-C [m]", 0.2);
975
976 itsAttr[C] =
977 Attributes::makeReal("C", "In case of BOXCORNER Specifies height of corner C [m]", 0.01);
978
980 Attributes::makeReal("XYZSCALE", "Multiplicative scaling factor for coordinates ", 1.0);
981
982 itsAttr[XSCALE] =
983 Attributes::makeReal("XSCALE", "Multiplicative scaling factor for X coordinates ", 1.0);
984
985 itsAttr[YSCALE] =
986 Attributes::makeReal("YSCALE", "Multiplicative scaling factor for Y coordinates ", 1.0);
987
988 itsAttr[ZSCALE] =
989 Attributes::makeReal("ZSCALE", "Multiplicative scaling factor for Z coordinates ", 1.0);
990
991 itsAttr[ZSHIFT] = Attributes::makeReal("ZSHIFT", "Shift in z direction", 0.0);
992
993 itsAttr[INSIDEPOINT] = Attributes::makeRealArray("INSIDEPOINT", "A point inside the geometry");
994
996
997 BoundaryGeometry* defGeometry = clone("UNNAMED_GEOMETRY");
998 defGeometry->builtin = true;
999
1000 Tinitialize_m = IpplTimings::getTimer("Initialize geometry");
1001 TisInside_m = IpplTimings::getTimer("Inside test");
1002 TfastIsInside_m = IpplTimings::getTimer("Fast inside test");
1003 TRayTrace_m = IpplTimings::getTimer("Ray tracing");
1004 TPartInside_m = IpplTimings::getTimer("Particle Inside");
1005
1007
1008 try {
1009 defGeometry->update();
1010 OpalData::getInstance()->define(defGeometry);
1011 } catch (...) {
1012 delete defGeometry;
1013 }
1014 gsl_rng_env_setup();
1015 randGen_m = gsl_rng_alloc(gsl_rng_default);
1016
1017 if (!h5FileName_m.empty())
1018 initialize();
1019}
1020
1022 : Definition(name, parent) {
1023 gsl_rng_env_setup();
1024 randGen_m = gsl_rng_alloc(gsl_rng_default);
1025
1026 Tinitialize_m = IpplTimings::getTimer("Initialize geometry");
1027 TisInside_m = IpplTimings::getTimer("Inside test");
1028 TfastIsInside_m = IpplTimings::getTimer("Fast inside test");
1029 TRayTrace_m = IpplTimings::getTimer("Ray tracing");
1030 TPartInside_m = IpplTimings::getTimer("Particle Inside");
1031
1033 if (!h5FileName_m.empty())
1034 initialize();
1035}
1036
1040
1042 // Can replace only by another GEOMETRY.
1043 return dynamic_cast<BGeometryBase*>(object) != 0;
1044}
1045
1047 return new BoundaryGeometry(name, this);
1048}
1049
1051 if (getOpalName().empty())
1052 setOpalName("UNNAMED_GEOMETRY");
1053}
1054
1056 update();
1057 Tinitialize_m = IpplTimings::getTimer("Initialize geometry");
1058 TisInside_m = IpplTimings::getTimer("Inside test");
1059 TfastIsInside_m = IpplTimings::getTimer("Fast inside test");
1060 TRayTrace_m = IpplTimings::getTimer("Ray tracing");
1061 TPartInside_m = IpplTimings::getTimer("Particle Inside");
1062}
1063
1064BoundaryGeometry* BoundaryGeometry::find(const std::string& name) {
1065 BoundaryGeometry* geom = dynamic_cast<BoundaryGeometry*>(OpalData::getInstance()->find(name));
1066
1067 if (geom == 0)
1068 throw OpalException("BoundaryGeometry::find()", "Geometry \"" + name + "\" not found.");
1069 return geom;
1070}
1071
1074
1076 const int triangle_id, const int i, const int j, const int k) {
1077 const Triangle t(getPoint(triangle_id, 1), getPoint(triangle_id, 2), getPoint(triangle_id, 3));
1078
1079 const Vector_t<double, 3> P(
1080 i * voxelMesh_m.sizeOfVoxel[0] + voxelMesh_m.minExtent[0],
1081 j * voxelMesh_m.sizeOfVoxel[1] + voxelMesh_m.minExtent[1],
1082 k * voxelMesh_m.sizeOfVoxel[2] + voxelMesh_m.minExtent[2]);
1083
1084 Voxel v(P, P + voxelMesh_m.sizeOfVoxel);
1085
1086 return v.intersect(t);
1087}
1088
1089/*
1090 Find the 3D intersection of a line segment, ray or line with a triangle.
1091
1092 Input:
1093 kind: type of test: SEGMENT, RAY or LINE
1094 P0, P0: defining
1095 a line segment from P0 to P1 or
1096 a ray starting at P0 with directional vector P1-P0 or
1097 a line through P0 and P1
1098 V0, V1, V2: the triangle vertices
1099
1100 Output:
1101 I: intersection point (when it exists)
1102
1103 Return values for line segment and ray test :
1104 -1 = triangle is degenerated (a segment or point)
1105 0 = disjoint (no intersect)
1106 1 = are in the same plane
1107 2 = intersect in unique point I1
1108
1109 Return values for line intersection test :
1110 -1: triangle is degenerated (a segment or point)
1111 0: disjoint (no intersect)
1112 1: are in the same plane
1113 2: intersect in unique point I1, with r < 0.0
1114 3: intersect in unique point I1, with 0.0 <= r <= 1.0
1115 4: intersect in unique point I1, with 1.0 < r
1116
1117 For algorithm and implementation see:
1118 http://geomalgorithms.com/a06-_intersect-2.html
1119
1120 Copyright 2001 softSurfer, 2012 Dan Sunday
1121 This code may be freely used and modified for any purpose
1122 providing that this copyright notice is included with it.
1123 SoftSurfer makes no warranty for this code, and cannot be held
1124 liable for any real or imagined damage resulting from its use.
1125 Users of this code must verify correctness for their application.
1126 */
1127
1129 const enum INTERSECTION_TESTS kind, const Vector_t<double, 3>& P0,
1130 const Vector_t<double, 3>& P1, const int triangle_id, Vector_t<double, 3>& I) {
1131 const Vector_t<double, 3> V0 = getPoint(triangle_id, 1);
1132 const Vector_t<double, 3> V1 = getPoint(triangle_id, 2);
1133 const Vector_t<double, 3> V2 = getPoint(triangle_id, 3);
1134
1135 // get triangle edge vectors and plane normal
1136 const Vector_t<double, 3> u = V1 - V0; // triangle vectors
1137 const Vector_t<double, 3> v = V2 - V0;
1138 const Vector_t<double, 3> n = cross(u, v);
1139 /* ADA
1140 if (n == (Vector_t<double, 3>(0))) // triangle is degenerate
1141 return -1; // do not deal with this case
1142 */
1143 const Vector_t<double, 3> dir = P1 - P0; // ray direction vector
1144 const Vector_t<double, 3> w0 = P0 - V0;
1145 const double a = -dot(n, w0);
1146 const double b = dot(n, dir);
1147 if (cmp::eq_zero(b)) { // ray is parallel to triangle plane
1148 if (cmp::eq_zero(a)) { // ray lies in triangle plane
1149 return 1;
1150 } else { // ray disjoint from plane
1151 return 0;
1152 }
1153 }
1154
1155 // get intersect point of ray with triangle plane
1156 const double r = a / b;
1157 switch (kind) {
1158 case RAY:
1159 if (cmp::lt_zero(r)) { // ray goes away from triangle
1160 return 0; // => no intersect
1161 }
1162 break;
1163 case SEGMENT:
1164 if (cmp::lt_zero(r) || cmp::lt(1.0, r)) { // intersection on extended
1165 return 0; // segment
1166 }
1167 break;
1168 case LINE:
1169 break;
1170 };
1171 I = P0 + r * dir; // intersect point of ray and plane
1172
1173 // is I inside T?
1174 const double uu = dot(u, u);
1175 const double uv = dot(u, v);
1176 const double vv = dot(v, v);
1177 const Vector_t<double, 3> w = I - V0;
1178 const double wu = dot(w, u);
1179 const double wv = dot(w, v);
1180 const double D = uv * uv - uu * vv;
1181
1182 // get and test parametric coords
1183 const double s = (uv * wv - vv * wu) / D;
1184 if (cmp::lt_zero(s) || cmp::gt(s, 1.0)) { // I is outside T
1185 return 0;
1186 }
1187 const double t = (uv * wu - uu * wv) / D;
1188 if (cmp::lt_zero(t) || cmp::gt((s + t), 1.0)) { // I is outside T
1189 return 0;
1190 }
1191 // intersection point is in triangle
1192 if (cmp::lt_zero(r)) { // in extended segment in opposite
1193 return 2; // direction of ray
1194 } else if (cmp::ge_zero(r) && cmp::le(r, 1.0)) { // in segment
1195 return 3;
1196 } else { // in extended segment in
1197 return 4; // direction of ray
1198 }
1199}
1200
1201static inline double magnitude(const Vector_t<double, 3>& v) {
1202 return std::sqrt(dot(v, v));
1203}
1204
1205bool BoundaryGeometry::isInside(const Vector_t<double, 3>& P // [in] pt to test
1206) {
1207 /*
1208 select a "close" reference pt outside the bounding box
1209 */
1210 // right boundary of bounding box (x direction)
1211 double x = minExtent_m[0] - 0.01;
1212 double distance = P[0] - x;
1213 Vector_t<double, 3> ref_pt{x, P[1], P[2]};
1214
1215 // left boundary of bounding box (x direction)
1216 x = maxExtent_m[0] + 0.01;
1217 if (cmp::lt(x - P[0], distance)) {
1218 distance = x - P[0];
1219 ref_pt = {x, P[1], P[2]};
1220 }
1221
1222 // lower boundary of bounding box (y direction)
1223 double y = minExtent_m[1] - 0.01;
1224 if (cmp::lt(P[1] - y, distance)) {
1225 distance = P[1] - y;
1226 ref_pt = {P[0], y, P[1]};
1227 }
1228
1229 // upper boundary of bounding box (y direction)
1230 y = maxExtent_m[1] + 0.01;
1231 if (cmp::lt(y - P[1], distance)) {
1232 distance = y - P[1];
1233 ref_pt = {P[0], y, P[2]};
1234 }
1235 // front boundary of bounding box (z direction)
1236 double z = minExtent_m[2] - 0.01;
1237 if (cmp::lt(P[2] - z, distance)) {
1238 distance = P[2] - z;
1239 ref_pt = {P[0], P[1], z};
1240 }
1241 // back boundary of bounding box (z direction)
1242 z = maxExtent_m[2] + 0.01;
1243 if (cmp::lt(z - P[2], distance)) {
1244 ref_pt = {P[0], P[1], z};
1245 }
1246
1247 /*
1248 the test returns the number of intersections =>
1249 since the reference point is outside, P is inside
1250 if the result is odd.
1251 */
1252 int k = fastIsInside(ref_pt, P);
1253 return (k % 2) == 1;
1254}
1255
1256/*
1257 searching a point inside the geometry.
1258
1259 sketch of the algorithm:
1260 In a first step, we try to find a line segment defined by one
1261 point outside the bounding box and a point somewhere inside the
1262 bounding box which has intersects with the geometry.
1263
1264 If the number of intersections is odd, the center point is inside
1265 the geometry and we are already done.
1266
1267 If the number of intersections is even, there must be points on
1268 this line segment which are inside the geometry. In the next step
1269 we have to find one if these points.
1270
1271
1272 A bit more in detail:
1273
1274 1. Finding a line segment intersecting the geometry
1275 For the fast isInside test it is of advantage to choose line segments
1276 parallel to the X, Y or Z axis. In this implementation we choose as
1277 point outside the bounding box a point on an axis but close to the
1278 bounding box and the center of the bounding box. This gives us six
1279 line segments to test. This covers not all possible geometries but
1280 most likely almost all. If not, it's easy to extend.
1281
1282 2. Searching for a point inside the geometry
1283 In the first step we get a line segment from which we know, that one
1284 point is ouside the geometry (P_out) and the other inside the bounding
1285 box (Q). We also know the number of intersections n_i of this line
1286 segment with the geometry.
1287
1288 If n_i is odd, Q is inside the boundary!
1289
1290 while (true); do
1291 bisect the line segment [P_out, Q], let B the bisecting point.
1292
1293 compute number of intersections of the line segment [P_out, B]
1294 and the geometry.
1295
1296 If the number of intersections is odd, then B is inside the geometry
1297 and we are done. Set P_in = B and exit loop.
1298
1299 Otherwise we have either no or an even number of intersections.
1300 In both cases this implies that B is a point outside the geometry.
1301
1302 If the number of intersection of [P_out, B] is even but not equal zero,
1303 it might be that *all* intersections are in this line segment and none in
1304 [B, Q].
1305 In this case we continue with the line segment [P_out, Q] = [P_out, B],
1306 otherwise with the line segment [P_out, Q] = [B, Q].
1307*/
1309 *gmsg << level2 << "* Searching for a point inside the geometry..." << endl;
1310 /*
1311 find line segment
1312 */
1314 std::vector<Vector_t<double, 3>> P_outs{
1315 {minExtent_m[0] - 0.01, Q[1], Q[2]}, {maxExtent_m[0] + 0.01, Q[1], Q[2]},
1316 {Q[0], minExtent_m[1] - 0.01, Q[2]}, {Q[0], maxExtent_m[1] + 0.01, Q[2]},
1317 {Q[0], Q[1], minExtent_m[2] - 0.01}, {Q[0], Q[1], maxExtent_m[2] + 0.01}};
1318 int n_i = 0;
1319 Vector_t<double, 3> P_out;
1320 for (const auto& P : P_outs) {
1321 n_i = fastIsInside(P, Q);
1322 if (n_i != 0) {
1323 P_out = P;
1324 break;
1325 }
1326 }
1327 if (n_i == 0) {
1328 // this is possible with some obscure geometries.
1329 return false;
1330 }
1331
1332 /*
1333 if the number of intersections is odd, Q is inside the geometry
1334 */
1335 if (n_i % 2 == 1) {
1336 insidePoint_m = Q;
1337 return true;
1338 }
1339 while (true) {
1340 Vector_t<double, 3> B{(P_out + Q) / 2};
1341 int n = fastIsInside(P_out, B);
1342 if (n % 2 == 1) {
1343 insidePoint_m = B;
1344 return true;
1345 } else if (n == n_i) {
1346 Q = B;
1347 } else {
1348 P_out = B;
1349 }
1350 n_i = n;
1351 }
1352 // never reached
1353 return false;
1354}
1355
1356/*
1357 Game plan:
1358 Count number of intersection of the line segment defined by P and a reference
1359 pt with the boundary. If the reference pt is inside the boundary and the number
1360 of intersections is even, then P is inside the geometry. Otherwise P is outside.
1361 To count the number of intersection, we divide the line segment in N segments
1362 and run the line-segment boundary intersection test for all these segments.
1363 N must be choosen carefully. It shouldn't be to large to avoid needless test.
1364 */
1366 const Vector_t<double, 3>& reference_pt, // [in] reference pt inside the boundary
1367 const Vector_t<double, 3>& P // [in] pt to test
1368) {
1369 const Voxel c(minExtent_m, maxExtent_m);
1370 if (!c.isInside(P))
1371 return 1;
1372 IpplTimings::startTimer(TfastIsInside_m);
1373#ifdef ENABLE_DEBUG
1374 int saved_flags = debugFlags_m;
1376 *gmsg << "* " << __func__ << ": "
1377 << "reference_pt=" << reference_pt << ", P=" << P << endl;
1379 }
1380#endif
1381 const Vector_t<double, 3> v = reference_pt - P;
1382 const int N = std::ceil(
1383 magnitude(v)
1384 / std::min(
1385 {voxelMesh_m.sizeOfVoxel[0], voxelMesh_m.sizeOfVoxel[1], voxelMesh_m.sizeOfVoxel[2]}));
1386 const Vector_t<double, 3> v_ = v / N;
1387 Vector_t<double, 3> P0 = P;
1388 Vector_t<double, 3> P1 = P + v_;
1390 int triangle_id = -1;
1391 int result = 0;
1392 for (int i = 0; i < N; i++) {
1393 result += intersectTinyLineSegmentBoundary(P0, P1, I, triangle_id);
1394 P0 = P1;
1395 P1 += v_;
1396 }
1397#ifdef ENABLE_DEBUG
1399 *gmsg << "* " << __func__ << ": "
1400 << "result: " << result << endl;
1401 debugFlags_m = saved_flags;
1402 }
1403#endif
1404 IpplTimings::stopTimer(TfastIsInside_m);
1405 return result;
1406}
1407
1408/*
1409 P must be *inside* the boundary geometry!
1410
1411 return value:
1412 0 no intersection
1413 1 intersection found, I is set to the first intersection coordinates in
1414 ray direction
1415 */
1418 IpplTimings::startTimer(TRayTrace_m);
1419#ifdef ENABLE_DEBUG
1420 int saved_flags = debugFlags_m;
1422 *gmsg << "* " << __func__ << ": "
1423 << " ray: "
1424 << " origin=" << P << " dir=" << v << endl;
1426 }
1427#endif
1428 /*
1429 set P1 to intersection of ray with bbox of voxel mesh
1430 run line segment boundary intersection test with P and P1
1431 */
1432 Ray r = Ray(P, v);
1433 Voxel c = Voxel(
1434 voxelMesh_m.minExtent + 0.25 * voxelMesh_m.sizeOfVoxel,
1435 voxelMesh_m.maxExtent - 0.25 * voxelMesh_m.sizeOfVoxel);
1436 double tmin = 0.0;
1437 double tmax = 0.0;
1438 c.intersect(r, tmin, tmax);
1439 int triangle_id = -1;
1440 int result = (intersectLineSegmentBoundary(P, P + (tmax * v), I, triangle_id) > 0) ? 1 : 0;
1441#ifdef ENABLE_DEBUG
1443 *gmsg << "* " << __func__ << ": "
1444 << " result=" << result << " I=" << I << endl;
1445 debugFlags_m = saved_flags;
1446 }
1447#endif
1448 IpplTimings::stopTimer(TRayTrace_m);
1449 return result;
1450}
1451
1452/*
1453 Map point to unique voxel ID.
1454
1455 Remember:
1456 * hr_m: is the mesh size
1457 * nr_m: number of mesh points
1458 */
1459inline int BoundaryGeometry::mapVoxelIndices2ID(const int i, const int j, const int k) {
1460 if (i < 0 || i >= voxelMesh_m.nr_m[0] || j < 0 || j >= voxelMesh_m.nr_m[1] || k < 0
1461 || k >= voxelMesh_m.nr_m[2]) {
1462 return 0;
1463 }
1464 return 1 + k * voxelMesh_m.nr_m[0] * voxelMesh_m.nr_m[1] + j * voxelMesh_m.nr_m[0] + i;
1465}
1466
1467#define mapPoint2VoxelIndices(pt, i, j, k) \
1468 { \
1469 i = floor((pt[0] - voxelMesh_m.minExtent[0]) / voxelMesh_m.sizeOfVoxel[0]); \
1470 j = floor((pt[1] - voxelMesh_m.minExtent[1]) / voxelMesh_m.sizeOfVoxel[1]); \
1471 k = floor((pt[2] - voxelMesh_m.minExtent[2]) / voxelMesh_m.sizeOfVoxel[2]); \
1472 if (!(0 <= i && i < voxelMesh_m.nr_m[0] && 0 <= j && j < voxelMesh_m.nr_m[1] && 0 <= k \
1473 && k < voxelMesh_m.nr_m[2])) { \
1474 *gmsg << level2 << "* " << __func__ << ":" \
1475 << " WARNING: pt=" << pt << " is outside the bbox" \
1476 << " i=" << i << " j=" << j << " k=" << k << endl; \
1477 } \
1478 }
1479
1481 const int i, const int j, const int k) {
1482 return Vector_t<double, 3>(
1483 i * voxelMesh_m.sizeOfVoxel[0] + voxelMesh_m.minExtent[0],
1484 j * voxelMesh_m.sizeOfVoxel[1] + voxelMesh_m.minExtent[1],
1485 k * voxelMesh_m.sizeOfVoxel[2] + voxelMesh_m.minExtent[2]);
1486}
1487
1489 const int i = std::floor((pt[0] - voxelMesh_m.minExtent[0]) / voxelMesh_m.sizeOfVoxel[0]);
1490 const int j = std::floor((pt[1] - voxelMesh_m.minExtent[1]) / voxelMesh_m.sizeOfVoxel[1]);
1491 const int k = std::floor((pt[2] - voxelMesh_m.minExtent[2]) / voxelMesh_m.sizeOfVoxel[2]);
1492
1493 return mapIndices2Voxel(i, j, k);
1494}
1495
1497 for (unsigned int triangle_id = 0; triangle_id < Triangles_m.size(); triangle_id++) {
1498 Vector_t<double, 3> v1 = getPoint(triangle_id, 1);
1499 Vector_t<double, 3> v2 = getPoint(triangle_id, 2);
1500 Vector_t<double, 3> v3 = getPoint(triangle_id, 3);
1501 Vector_t<double, 3> bbox_min = {
1502 std::min({v1[0], v2[0], v3[0]}), std::min({v1[1], v2[1], v3[1]}),
1503 std::min({v1[2], v2[2], v3[2]})};
1504 Vector_t<double, 3> bbox_max = {
1505 std::max({v1[0], v2[0], v3[0]}), std::max({v1[1], v2[1], v3[1]}),
1506 std::max({v1[2], v2[2], v3[2]})};
1507 int i_min, j_min, k_min;
1508 int i_max, j_max, k_max;
1509 mapPoint2VoxelIndices(bbox_min, i_min, j_min, k_min);
1510 mapPoint2VoxelIndices(bbox_max, i_max, j_max, k_max);
1511
1512 for (int i = i_min; i <= i_max; i++) {
1513 for (int j = j_min; j <= j_max; j++) {
1514 for (int k = k_min; k <= k_max; k++) {
1515 // test if voxel (i,j,k) has an intersection with triangle
1516 if (intersectTriangleVoxel(triangle_id, i, j, k) == INSIDE) {
1517 auto id = mapVoxelIndices2ID(i, j, k);
1518 voxelMesh_m.ids[id].insert(triangle_id);
1519 }
1520 }
1521 }
1522 }
1523 } // for_each triangle
1524 *gmsg << level2 << "* Mesh voxelization done" << endl;
1525
1526 // write voxel mesh into VTK file
1527 if (ippl::Comm->rank() == 0 && Options::enableVTK) {
1528 std::string vtkFileName = Util::combineFilePath(
1529 {OpalData::getInstance()->getAuxiliaryOutputDirectory(), "testBBox.vtk"});
1530 bool writeVTK = false;
1531
1532 if (!boost::filesystem::exists(vtkFileName)) {
1533 writeVTK = true;
1534 } else {
1535 std::time_t t_geom = boost::filesystem::last_write_time(h5FileName_m);
1536 std::time_t t_vtk = boost::filesystem::last_write_time(vtkFileName);
1537 if (std::difftime(t_geom, t_vtk) > 0)
1538 writeVTK = true;
1539 }
1540
1541 if (writeVTK) {
1542 write_voxel_mesh(
1543 vtkFileName, voxelMesh_m.ids, voxelMesh_m.sizeOfVoxel, voxelMesh_m.nr_m,
1544 voxelMesh_m.minExtent);
1545 }
1546 }
1547}
1548
1550 class Local {
1551 public:
1552 static void computeGeometryInterval(BoundaryGeometry* bg) {
1553 bg->minExtent_m = get_min_extent(bg->Points_m);
1554 bg->maxExtent_m = get_max_extent(bg->Points_m);
1555
1556 /*
1557 Calculate the maximum size of triangles. This value will be used to
1558 define the voxel size
1559 */
1560 double longest_edge_max_m = 0.0;
1561 for (unsigned int i = 0; i < bg->Triangles_m.size(); i++) {
1562 // compute length of longest edge
1563 const Vector_t<double, 3> x1 = bg->getPoint(i, 1);
1564 const Vector_t<double, 3> x2 = bg->getPoint(i, 2);
1565 const Vector_t<double, 3> x3 = bg->getPoint(i, 3);
1566 const double length_edge1 =
1567 std::sqrt(SQR(x1[0] - x2[0]) + SQR(x1[1] - x2[1]) + SQR(x1[2] - x2[2]));
1568 const double length_edge2 =
1569 std::sqrt(SQR(x3[0] - x2[0]) + SQR(x3[1] - x2[1]) + SQR(x3[2] - x2[2]));
1570 const double length_edge3 =
1571 std::sqrt(SQR(x3[0] - x1[0]) + SQR(x3[1] - x1[1]) + SQR(x3[2] - x1[2]));
1572
1573 double max = length_edge1;
1574 if (length_edge2 > max)
1575 max = length_edge2;
1576 if (length_edge3 > max)
1577 max = length_edge3;
1578
1579 // save min and max of length of longest edge
1580 if (longest_edge_max_m < max)
1581 longest_edge_max_m = max;
1582 }
1583
1584 /*
1585 In principal the number of discretization nr_m is the extent of
1586 the geometry divided by the extent of the largest triangle. Whereby
1587 the extent of a triangle is defined as the lenght of its longest
1588 edge. Thus the largest triangle is the triangle with the longest edge.
1589
1590 But if the hot spot, i.e., the multipacting/field emission zone is
1591 too small that the normal bounding box covers the whole hot spot, the
1592 expensive triangle-line intersection tests will be frequently called.
1593 In these cases, we have to use smaller bounding box size to speed up
1594 simulation.
1595
1596 Todo:
1597 The relation between bounding box size and simulation time step &
1598 geometry shape maybe need to be summarized and modeled in a more
1599 flexible manner and could be adjusted in input file.
1600 */
1601 Vector_t<double, 3> extent = bg->maxExtent_m - bg->minExtent_m;
1602 bg->voxelMesh_m.nr_m(0) = 16 * (int)std::floor(extent[0] / longest_edge_max_m);
1603 bg->voxelMesh_m.nr_m(1) = 16 * (int)std::floor(extent[1] / longest_edge_max_m);
1604 bg->voxelMesh_m.nr_m(2) = 16 * (int)std::floor(extent[2] / longest_edge_max_m);
1605
1606 bg->voxelMesh_m.sizeOfVoxel = extent / bg->voxelMesh_m.nr_m;
1607 bg->voxelMesh_m.minExtent = bg->minExtent_m - 0.5 * bg->voxelMesh_m.sizeOfVoxel;
1608 bg->voxelMesh_m.maxExtent = bg->maxExtent_m + 0.5 * bg->voxelMesh_m.sizeOfVoxel;
1609 bg->voxelMesh_m.nr_m += 1;
1610 }
1611
1612 /*
1613 To speed up ray-triangle intersection tests, the normal vector of
1614 all triangles are pointing inward. Since this is clearly not
1615 guaranteed for the triangles in the H5hut file, this must be checked
1616 for each triangle and - if necessary changed - after reading the mesh.
1617
1618 To test whether the normal of a triangle is pointing inward or outward,
1619 we choose a random point P close to the center of the triangle and test
1620 whether this point is inside or outside the geometry. The way we choose
1621 P guarantees that the vector spanned by P and a vertex of the triangle
1622 points into the same direction as the normal vector. From this it
1623 follows that if P is inside the geometry the normal vector is pointing
1624 to the inside and vise versa.
1625
1626 Since the inside-test is computational expensive we perform this test
1627 for one reference triangle T (per sub-mesh) only. Knowing the adjacent
1628 triangles for all three edges of a triangle for all triangles of the
1629 mesh facilitates another approach using the orientation of the
1630 reference triangle T. Assuming that the normal vector of T points to
1631 the inside of the geometry an adjacent triangle of T has an inward
1632 pointing normal vector if and only if it has the same orientation as
1633 T.
1634
1635 Starting with the reference triangle T we can change the orientation
1636 of the adjancent triangle of T and so on.
1637
1638 NOTE: For the time being we do not make use of the inward pointing
1639 normals.
1640
1641 FIXME: Describe the basic ideas behind the following comment! Without
1642 it is completely unclear.
1643
1644 Following combinations are possible:
1645 1,1 && 2,2 1,2 && 2,1 1,3 && 2,1
1646 1,1 && 2,3 1,2 && 2,3 1,3 && 2,2
1647 1,1 && 3,2 1,2 && 3,1 1,3 && 3,1
1648 1,1 && 3,3 1,2 && 3,3 1,3 && 3,2
1649
1650 (2,1 && 1,2) (2,2 && 1,1) (2,3 && 1,1)
1651 (2,1 && 1,3) (2,2 && 1,3) (2,3 && 1,2)
1652 2,1 && 3,2 2,2 && 3,1 2,3 && 3,1
1653 2,1 && 3,3 2,2 && 3,3 2,3 && 3,2
1654
1655 (3,1 && 1,2) (3,2 && 1,1) (3,3 && 1,1)
1656 (3,1 && 1,3) (3,2 && 1,3) (3,3 && 1,2)
1657 (3,1 && 2,2) (3,2 && 2,1) (3,3 && 2,1)
1658 (3,1 && 2,3) (3,2 && 2,3) (3,3 && 2,2)
1659
1660 Note:
1661 Since we find vertices with lower enumeration first, we
1662 can ignore combinations in ()
1663
1664 2 2 2 3 3 2 3 3
1665 * * * *
1666 /|\ /|\ /|\ /|\
1667 / | \ / | \ / | \ / | \
1668 / | \ / | \ / | \ / | \
1669 / | \ / | \ / | \ / | \
1670 *----*----* *----*----* *----*----* *----*----*
1671 3 1 1 3 3 1 1 2 2 1 1 3 2 1 1 2
1672diff: (1,1) (1,2) (2,1) (2,2)
1673change orient.: yes no no yes
1674
1675
1676 2 1 2 3 3 1 3 3
1677 * * * *
1678 /|\ /|\ /|\ /|\
1679 / | \ / | \ / | \ / | \
1680 / | \ / | \ / | \ / | \
1681 / | \ / | \ / | \ / | \
1682 *----*----* *----*----* *----*----* *----*----*
1683 3 1 2 3 3 1 2 1 2 1 2 3 2 1 2 1
1684diff: (1,-1) (1,1) (2,-1) (2,1)
1685change orient.: no yes yes no
1686
1687
1688 2 1 2 2 3 1 3 2
1689 * * * *
1690 /|\ /|\ /|\ /|\
1691 / | \ / | \ / | \ / | \
1692 / | \ / | \ / | \ / | \
1693 / | \ / | \ / | \ / | \
1694 *----*----* *----*----* *----*----* *----*----*
1695 3 1 3 2 3 1 3 1 2 1 3 2 2 1 3 1
1696diff: (1,-2) (1,-1) (2,-2) (2,-1)
1697change orient.: yes no no yes
1698
1699 3 2 3 3
1700 * *
1701 /|\ /|\
1702 / | \ / | \
1703 / | \ / | \
1704 / | \ / | \
1705 *----*----* *----*----*
1706 1 2 1 3 1 2 1 2
1707diff: (1,1) (1,2)
1708change orient.: yes no
1709
1710 3 1 3 3
1711 * *
1712 /|\ /|\
1713 / | \ / | \
1714 / | \ / | \
1715 / | \ / | \
1716 *----*----* *----*----*
1717 1 2 2 3 1 2 2 1
1718diff: (1,-1) (1,1)
1719change orient.: no yes
1720
1721 3 1 3 2
1722 * *
1723 /|\ /|\
1724 / | \ / | \
1725 / | \ / | \
1726 / | \ / | \
1727 *----*----* *----*----*
1728 1 2 3 2 1 2 3 1
1729diff: (1,-2) (1,-1)
1730change orient.: yes no
1731
1732
1733Change orientation if diff is:
1734(1,1) || (1,-2) || (2,2) || (2,-1) || (2,-1)
1735
1736 */
1737
1738 static void computeTriangleNeighbors(
1739 BoundaryGeometry* bg, std::vector<std::set<unsigned int>>& neighbors) {
1740 std::vector<std::set<unsigned int>> adjacencies_to_pt(bg->Points_m.size());
1741
1742 // for each triangles find adjacent triangles for each vertex
1743 for (unsigned int triangle_id = 0; triangle_id < bg->Triangles_m.size();
1744 triangle_id++) {
1745 for (unsigned int j = 1; j <= 3; j++) {
1746 auto pt_id = bg->PointID(triangle_id, j);
1747 PAssert(pt_id < bg->Points_m.size());
1748 adjacencies_to_pt[pt_id].insert(triangle_id);
1749 }
1750 }
1751
1752 for (unsigned int triangle_id = 0; triangle_id < bg->Triangles_m.size();
1753 triangle_id++) {
1754 std::set<unsigned int> to_A = adjacencies_to_pt[bg->PointID(triangle_id, 1)];
1755 std::set<unsigned int> to_B = adjacencies_to_pt[bg->PointID(triangle_id, 2)];
1756 std::set<unsigned int> to_C = adjacencies_to_pt[bg->PointID(triangle_id, 3)];
1757
1758 std::set<unsigned int> intersect;
1759 std::set_intersection(
1760 to_A.begin(), to_A.end(), to_B.begin(), to_B.end(),
1761 std::inserter(intersect, intersect.begin()));
1762 std::set_intersection(
1763 to_B.begin(), to_B.end(), to_C.begin(), to_C.end(),
1764 std::inserter(intersect, intersect.begin()));
1765 std::set_intersection(
1766 to_C.begin(), to_C.end(), to_A.begin(), to_A.end(),
1767 std::inserter(intersect, intersect.begin()));
1768 intersect.erase(triangle_id);
1769
1770 neighbors[triangle_id] = intersect;
1771 }
1772 *gmsg << level2 << "* " << __func__ << ": Computing neighbors done" << endl;
1773 }
1774
1775 /*
1776 Helper function for hasInwardPointingNormal()
1777
1778 Determine if a point x is outside or inside the geometry or just on
1779 the boundary. Return true if point is inside geometry or on the
1780 boundary, false otherwise
1781
1782 The basic idea here is:
1783 If a line segment from the point to test to a random point outside
1784 the geometry has has an even number of intersections with the
1785 boundary, the point is outside the geometry.
1786
1787 Note:
1788 If the point is on the boundary, the number of intersections is 1.
1789 Points on the boundary are handled as inside.
1790
1791 A random selection of the reference point outside the boundary avoids
1792 some specific issues, like line parallel to boundary.
1793 */
1794 static inline bool isInside(BoundaryGeometry* bg, const Vector_t<double, 3> x) {
1795 IpplTimings::startTimer(bg->TisInside_m);
1796
1798 bg->maxExtent_m[0] * (1.1 + gsl_rng_uniform(bg->randGen_m)),
1799 bg->maxExtent_m[1] * (1.1 + gsl_rng_uniform(bg->randGen_m)),
1800 bg->maxExtent_m[2] * (1.1 + gsl_rng_uniform(bg->randGen_m)));
1801
1802 std::vector<Vector_t<double, 3>> intersection_points;
1803 // int num_intersections = 0;
1804
1805 for (unsigned int triangle_id = 0; triangle_id < bg->Triangles_m.size();
1806 triangle_id++) {
1807 Vector_t<double, 3> result;
1808 if (bg->intersectLineTriangle(SEGMENT, x, y, triangle_id, result)) {
1809 intersection_points.push_back(result);
1810 // num_intersections++;
1811 }
1812 }
1813 IpplTimings::stopTimer(bg->TisInside_m);
1814 return ((intersection_points.size() % 2) == 1);
1815 }
1816
1817 // helper for function makeTriangleNormalInwardPointing()
1818 static bool hasInwardPointingNormal(BoundaryGeometry* const bg, const int triangle_id) {
1819 const Vector_t<double, 3>& A = bg->getPoint(triangle_id, 1);
1820 const Vector_t<double, 3>& B = bg->getPoint(triangle_id, 2);
1821 const Vector_t<double, 3>& C = bg->getPoint(triangle_id, 3);
1822 const Vector_t<double, 3> triNormal = normalVector(A, B, C);
1823
1824 // choose a point P close to the center of the triangle
1825 // const Vector_t<double, 3> P = (A+B+C)/3 + triNormal * 0.1;
1826 double minvoxelmesh = bg->voxelMesh_m.sizeOfVoxel[0];
1827 if (minvoxelmesh > bg->voxelMesh_m.sizeOfVoxel[1])
1828 minvoxelmesh = bg->voxelMesh_m.sizeOfVoxel[1];
1829 if (minvoxelmesh > bg->voxelMesh_m.sizeOfVoxel[2])
1830 minvoxelmesh = bg->voxelMesh_m.sizeOfVoxel[2];
1831 const Vector_t<double, 3> P = (A + B + C) / 3 + triNormal * minvoxelmesh;
1832 /*
1833 The triangle normal points inward, if P is
1834 - outside the geometry and the dot product is negativ
1835 - or inside the geometry and the dot product is positiv
1836
1837 Remember:
1838 The dot product is positiv only if both vectors are
1839 pointing in the same direction.
1840 */
1841 const bool is_inside = isInside(bg, P);
1842 const Vector_t<double, 3> d = P - A;
1843 const double dotPA_N = dot(d, triNormal);
1844 return (is_inside && dotPA_N >= 0) || (!is_inside && dotPA_N < 0);
1845 }
1846
1847 // helper for function makeTriangleNormalInwardPointing()
1848 static void orientTriangle(BoundaryGeometry* bg, int ref_id, int triangle_id) {
1849 // find pts of common edge
1850 int ic[2];
1851 int id[2];
1852 int n = 0;
1853 for (int i = 1; i <= 3; i++) {
1854 for (int j = 1; j <= 3; j++) {
1855 if (bg->PointID(triangle_id, j) == bg->PointID(ref_id, i)) {
1856 id[n] = j;
1857 ic[n] = i;
1858 n++;
1859 if (n == 2)
1860 goto edge_found;
1861 }
1862 }
1863 }
1864 PAssert(n == 2);
1865edge_found:
1866 int diff = id[1] - id[0];
1867 if ((((ic[1] - ic[0]) == 1) && ((diff == 1) || (diff == -2)))
1868 || (((ic[1] - ic[0]) == 2) && ((diff == -1) || (diff == 2)))) {
1869 std::swap(bg->PointID(triangle_id, id[0]), bg->PointID(triangle_id, id[1]));
1870 }
1871 }
1872
1873 static void makeTriangleNormalInwardPointing(BoundaryGeometry* bg) {
1874 std::vector<std::set<unsigned int>> neighbors(bg->Triangles_m.size());
1875
1876 computeTriangleNeighbors(bg, neighbors);
1877
1878 // loop over all sub-meshes
1879 int triangle_id = 0;
1880 int parts = 0;
1881 std::vector<unsigned int> triangles(bg->Triangles_m.size());
1882 std::vector<unsigned int>::size_type queue_cursor = 0;
1883 std::vector<unsigned int>::size_type queue_end = 0;
1884 std::vector<bool> isOriented(bg->Triangles_m.size(), false);
1885 do {
1886 parts++;
1887 /*
1888 Find next untested triangle, trivial for the first sub-mesh.
1889 There is a least one not yet tested triangle!
1890 */
1891 while (isOriented[triangle_id])
1892 triangle_id++;
1893
1894 // ensure that normal of this triangle is inward pointing
1895 if (!hasInwardPointingNormal(bg, triangle_id)) {
1896 std::swap(bg->PointID(triangle_id, 2), bg->PointID(triangle_id, 3));
1897 }
1898 isOriented[triangle_id] = true;
1899
1900 // loop over all triangles in sub-mesh
1901 triangles[queue_end++] = triangle_id;
1902 do {
1903 for (auto neighbor_id : neighbors[triangle_id]) {
1904 if (isOriented[neighbor_id])
1905 continue;
1906 orientTriangle(bg, triangle_id, neighbor_id);
1907 isOriented[neighbor_id] = true;
1908 triangles[queue_end++] = neighbor_id;
1909 }
1910 queue_cursor++;
1911 } while (queue_cursor < queue_end && (triangle_id = triangles[queue_cursor], true));
1912 } while (queue_end < bg->Triangles_m.size());
1913
1914 if (parts == 1) {
1915 *gmsg << level2 << "* " << __func__ << ": mesh is contiguous" << endl;
1916 } else {
1917 *gmsg << level2 << "* " << __func__ << ": mesh is discontiguous (" << parts
1918 << ") parts" << endl;
1919 }
1920 *gmsg << level2 << "* Triangle Normal built done" << endl;
1921 }
1922 };
1923
1924 debugFlags_m = 0;
1925 *gmsg << level2 << "* Initializing Boundary Geometry..." << endl;
1926 IpplTimings::startTimer(Tinitialize_m);
1927
1928 if (!boost::filesystem::exists(h5FileName_m)) {
1929 throw OpalException(
1930 "BoundaryGeometry::initialize",
1931 "Failed to open file '" + h5FileName_m + "', please check if it exists");
1932 }
1933
1934 double xscale = Attributes::getReal(itsAttr[XSCALE]);
1935 double yscale = Attributes::getReal(itsAttr[YSCALE]);
1936 double zscale = Attributes::getReal(itsAttr[ZSCALE]);
1937 double xyzscale = Attributes::getReal(itsAttr[XYZSCALE]);
1938 double zshift = (double)(Attributes::getReal(itsAttr[ZSHIFT]));
1939
1940 h5_int64_t rc;
1941#if defined(NDEBUG)
1942 (void)rc;
1943#endif
1944 rc = H5SetErrorHandler(H5AbortErrorhandler);
1945 PAssert(rc != H5_ERR);
1946 H5SetVerbosityLevel(1);
1947
1948 h5_prop_t props = H5CreateFileProp();
1949 MPI_Comm comm = ippl::Comm->getCommunicator();
1950 H5SetPropFileMPIOCollective(props, &comm);
1951 h5_file_t f = H5OpenFile(h5FileName_m.c_str(), H5_O_RDONLY, props);
1952 H5CloseProp(props);
1953
1954 h5t_mesh_t* m = nullptr;
1955 H5FedOpenTriangleMesh(f, "0", &m);
1956 H5FedSetLevel(m, 0);
1957
1958 auto numTriangles = H5FedGetNumElementsTotal(m);
1959 Triangles_m.resize(numTriangles);
1960
1961 // iterate over all co-dim 0 entities, i.e. elements
1962 h5_loc_id_t local_id;
1963 int i = 0;
1964 h5t_iterator_t* iter = H5FedBeginTraverseEntities(m, 0);
1965 while ((local_id = H5FedTraverseEntities(iter)) >= 0) {
1966 h5_loc_id_t local_vids[4];
1967 H5FedGetVertexIndicesOfEntity(m, local_id, local_vids);
1968 PointID(i, 0) = 0;
1969 PointID(i, 1) = local_vids[0];
1970 PointID(i, 2) = local_vids[1];
1971 PointID(i, 3) = local_vids[2];
1972 i++;
1973 }
1974 H5FedEndTraverseEntities(iter);
1975
1976 // loop over all vertices
1977 int num_points = H5FedGetNumVerticesTotal(m);
1978 Points_m.reserve(num_points);
1979 for (i = 0; i < num_points; i++) {
1980 h5_float64_t P[3];
1981 H5FedGetVertexCoordsByIndex(m, i, P);
1982 Points_m.push_back(Vector_t<double, 3>(
1983 P[0] * xyzscale * xscale, P[1] * xyzscale * yscale, P[2] * xyzscale * zscale + zshift));
1984 }
1985 H5FedCloseMesh(m);
1986 H5CloseFile(f);
1987 *gmsg << level2 << "* Reading mesh done" << endl;
1988
1989 Local::computeGeometryInterval(this);
1991 haveInsidePoint_m = false;
1992 std::vector<double> pt = Attributes::getRealArray(itsAttr[INSIDEPOINT]);
1993 if (!pt.empty()) {
1994 if (pt.size() != 3) {
1995 throw OpalException(
1996 "BoundaryGeometry::initialize()", "Dimension of INSIDEPOINT must be 3");
1997 }
1998 /* test whether this point is inside */
1999 insidePoint_m = {pt[0], pt[1], pt[2]};
2000 bool is_inside = isInside(insidePoint_m);
2001 if (is_inside == false) {
2002 throw OpalException(
2003 "BoundaryGeometry::initialize()", "INSIDEPOINT is not inside the geometry");
2004 }
2005 haveInsidePoint_m = true;
2006 } else {
2008 }
2009 if (haveInsidePoint_m == true) {
2010 *gmsg << level2 << "* using as point inside the geometry: (" << insidePoint_m[0] << ", "
2011 << insidePoint_m[1] << ", " << insidePoint_m[2] << ")" << endl;
2012 } else {
2013 *gmsg << level2 << "* no point inside the geometry found!" << endl;
2014 }
2015
2016 Local::makeTriangleNormalInwardPointing(this);
2017
2018 TriNormals_m.resize(Triangles_m.size());
2019 TriAreas_m.resize(Triangles_m.size());
2020
2021 for (size_t i = 0; i < Triangles_m.size(); i++) {
2022 const Vector_t<double, 3>& A = getPoint(i, 1);
2023 const Vector_t<double, 3>& B = getPoint(i, 2);
2024 const Vector_t<double, 3>& C = getPoint(i, 3);
2025
2026 TriAreas_m[i] = computeArea(A, B, C);
2027 TriNormals_m[i] = normalVector(A, B, C);
2028 }
2029 *gmsg << level2 << "* Triangle barycent built done" << endl;
2030
2031 *gmsg << *this << endl;
2032 ippl::Comm->barrier();
2033 IpplTimings::stopTimer(Tinitialize_m);
2034}
2035
2036/*
2037 Line segment triangle intersection test. This method should be used only
2038 for "tiny" line segments or, to be more exact, if the number of
2039 voxels covering the bounding box of the line segment is small (<<100).
2040
2041 Actually the method can be used for any line segment, but may not perform
2042 well. Performace depends on the size of the bounding box of the line
2043 segment.
2044
2045 The method returns the number of intersections of the line segment defined
2046 by the points P and Q with the boundary. If there are multiple intersections,
2047 the nearest intersection point with respect to P wil be returned.
2048 */
2050 const Vector_t<double, 3>& P, // [i] starting point of ray
2051 const Vector_t<double, 3>& Q, // [i] end point of ray
2052 Vector_t<double, 3>& intersect_pt, // [o] intersection with boundary
2053 int& triangle_id // [o] intersected triangle
2054) {
2055#ifdef ENABLE_DEBUG
2057 *gmsg << "* " << __func__ << ": "
2058 << " P = " << P << ", Q = " << Q << endl;
2059 }
2060#endif
2061 const Vector_t<double, 3> v_ = Q - P;
2062 const Ray r = Ray(P, v_);
2063 const Vector_t<double, 3> bbox_min = {
2064 std::min(P[0], Q[0]), std::min(P[1], Q[1]), std::min(P[2], Q[2])};
2065 const Vector_t<double, 3> bbox_max = {
2066 std::max(P[0], Q[0]), std::max(P[1], Q[1]), std::max(P[2], Q[2])};
2067 int i_min, i_max;
2068 int j_min, j_max;
2069 int k_min, k_max;
2070 mapPoint2VoxelIndices(bbox_min, i_min, j_min, k_min);
2071 mapPoint2VoxelIndices(bbox_max, i_max, j_max, k_max);
2072
2073 Vector_t<double, 3> tmp_intersect_pt = Q;
2074 double tmin = 1.1;
2075
2076 /*
2077 Triangles can - and in many cases do - intersect with more than one
2078 voxel. If we loop over all voxels intersecting with the line segment
2079 spaned by the points P and Q, we might perform the same line-triangle
2080 intersection test more than once. We must this into account when
2081 counting the intersections with the boundary.
2082
2083 To avoid multiple counting we can either
2084 - build a set of all relevant triangles and loop over this set
2085 - or we loop over all voxels and remember the intersecting triangles.
2086
2087 The first solution is implemented here.
2088 */
2089 std::unordered_set<int> triangle_ids;
2090 for (int i = i_min; i <= i_max; i++) {
2091 for (int j = j_min; j <= j_max; j++) {
2092 for (int k = k_min; k <= k_max; k++) {
2093 const Vector_t<double, 3> bmin = mapIndices2Voxel(i, j, k);
2094 const Voxel v(bmin, bmin + voxelMesh_m.sizeOfVoxel);
2095#ifdef ENABLE_DEBUG
2097 *gmsg << "* " << __func__ << ": "
2098 << " Test voxel: (" << i << ", " << j << ", " << k << "), " << v.pts[0]
2099 << v.pts[1] << endl;
2100 }
2101#endif
2102 /*
2103 do line segment and voxel intersect? continue if not
2104 */
2105 if (!v.intersect(r)) {
2106 continue;
2107 }
2108
2109 /*
2110 get triangles intersecting with this voxel and add them to
2111 the to be tested triangles.
2112 */
2113 const int voxel_id = mapVoxelIndices2ID(i, j, k);
2114 const auto triangles_intersecting_voxel = voxelMesh_m.ids.find(voxel_id);
2115 if (triangles_intersecting_voxel != voxelMesh_m.ids.end()) {
2116 triangle_ids.insert(
2117 triangles_intersecting_voxel->second.begin(),
2118 triangles_intersecting_voxel->second.end());
2119 }
2120 }
2121 }
2122 }
2123 /*
2124 test all triangles intersecting with one of the above voxels
2125 if there is more than one intersection, return closest
2126 */
2127 int num_intersections = 0;
2128 int tmp_intersect_result = 0;
2129
2130 for (auto it = triangle_ids.begin(); it != triangle_ids.end(); it++) {
2131 tmp_intersect_result = intersectLineTriangle(LINE, P, Q, *it, tmp_intersect_pt);
2132#ifdef ENABLE_DEBUG
2134 *gmsg << "* " << __func__ << ": "
2135 << " Test triangle: " << *it << " intersect: " << tmp_intersect_result
2136 << getPoint(*it, 1) << getPoint(*it, 2) << getPoint(*it, 3) << endl;
2137 }
2138#endif
2139 switch (tmp_intersect_result) {
2140 case 0: // no intersection
2141 case 2: // both points are outside
2142 case 4: // both points are inside
2143 break;
2144 case 1: // line and triangle are in same plane
2145 case 3: // unique intersection in segment
2146 double t;
2147 if (cmp::eq_zero(Q[0] - P[0]) == false) {
2148 t = (tmp_intersect_pt[0] - P[0]) / (Q[0] - P[0]);
2149 } else if (cmp::eq_zero(Q[1] - P[1]) == false) {
2150 t = (tmp_intersect_pt[1] - P[1]) / (Q[1] - P[1]);
2151 } else {
2152 t = (tmp_intersect_pt[2] - P[2]) / (Q[2] - P[2]);
2153 }
2154 num_intersections++;
2155 if (t < tmin) {
2156#ifdef ENABLE_DEBUG
2158 *gmsg << "* " << __func__ << ": "
2159 << " set triangle" << endl;
2160 }
2161#endif
2162 tmin = t;
2163 intersect_pt = tmp_intersect_pt;
2164 triangle_id = (*it);
2165 }
2166 break;
2167 case -1: // triangle is degenerated
2168 PAssert(tmp_intersect_result != -1);
2169 exit(42); // terminate even if NDEBUG is set
2170 }
2171 } // end for all triangles
2172 return num_intersections;
2173}
2174
2175/*
2176 General purpose line segment boundary intersection test.
2177
2178 The method returns with a value > 0 if an intersection was found.
2179 */
2181 const Vector_t<double, 3>& P0, // [in] starting point of ray
2182 const Vector_t<double, 3>& P1, // [in] end point of ray
2183 Vector_t<double, 3>& intersect_pt, // [out] intersection with boundary
2184 int& triangle_id // [out] triangle the line segment intersects with
2185) {
2186#ifdef ENABLE_DEBUG
2187 int saved_flags = debugFlags_m;
2189 *gmsg << "* " << __func__ << ": "
2190 << " P0 = " << P0 << " P1 = " << P1 << endl;
2192 }
2193#endif
2194 triangle_id = -1;
2195
2196 const Vector_t<double, 3> v = P1 - P0;
2197 int intersect_result = 0;
2198 int n = 0;
2199 int i_min, j_min, k_min;
2200 int i_max, j_max, k_max;
2201 do {
2202 n++;
2203 Vector_t<double, 3> Q = P0 + v / n;
2204 Vector_t<double, 3> bbox_min = {
2205 std::min(P0[0], Q[0]), std::min(P0[1], Q[1]), std::min(P0[2], Q[2])};
2206 Vector_t<double, 3> bbox_max = {
2207 std::max(P0[0], Q[0]), std::max(P0[1], Q[1]), std::max(P0[2], Q[2])};
2208 mapPoint2VoxelIndices(bbox_min, i_min, j_min, k_min);
2209 mapPoint2VoxelIndices(bbox_max, i_max, j_max, k_max);
2210 } while (((i_max - i_min + 1) * (j_max - j_min + 1) * (k_max - k_min + 1)) > 27);
2211 Vector_t<double, 3> P = P0;
2213 const Vector_t<double, 3> v_ = v / n;
2214
2215 for (int l = 1; l <= n; l++, P = Q) {
2216 Q = P0 + l * v_;
2217 intersect_result = intersectTinyLineSegmentBoundary(P, Q, intersect_pt, triangle_id);
2218 if (triangle_id != -1) {
2219 break;
2220 }
2221 }
2222#ifdef ENABLE_DEBUG
2224 *gmsg << "* " << __func__ << ": "
2225 << " result=" << intersect_result << " intersection pt: " << intersect_pt << endl;
2226 debugFlags_m = saved_flags;
2227 }
2228#endif
2229 return intersect_result;
2230}
2231
2241 const Vector_t<double, 3>& r, // [in] particle position
2242 const Vector_t<double, 3>& v, // [in] momentum
2243 const double dt, // [in]
2244 Vector_t<double, 3>& intersect_pt, // [out] intersection with boundary
2245 int& triangle_id // [out] intersected triangle
2246) {
2247#ifdef ENABLE_DEBUG
2248 int saved_flags = debugFlags_m;
2250 *gmsg << "* " << __func__ << ": "
2251 << " r=" << r << " v=" << v << " dt=" << dt << endl;
2253 }
2254#endif
2255 int ret = -1; // result defaults to no collision
2256
2257 // nothing to do if momenta == 0
2258 /* ADA
2259 if (v == (Vector_t<double, 3>)0)
2260 return ret;
2261 */
2262
2263 IpplTimings::startTimer(TPartInside_m);
2264
2265 // P0, P1: particle position in time steps n and n+1
2266 const Vector_t<double, 3> P0 = r;
2267 const Vector_t<double, 3> P1 = r + (Physics::c * v * dt / std::sqrt(1.0 + dot(v, v)));
2268
2269 Vector_t<double, 3> tmp_intersect_pt = 0.0;
2270 int tmp_triangle_id = -1;
2271 intersectTinyLineSegmentBoundary(P0, P1, tmp_intersect_pt, tmp_triangle_id);
2272 if (tmp_triangle_id >= 0) {
2273 intersect_pt = tmp_intersect_pt;
2274 triangle_id = tmp_triangle_id;
2275 ret = 0;
2276 }
2277#ifdef ENABLE_DEBUG
2279 *gmsg << "* " << __func__ << ":"
2280 << " result=" << ret;
2281 if (ret == 0) {
2282 *gmsg << " intersetion=" << intersect_pt;
2283 }
2284 *gmsg << endl;
2285 debugFlags_m = saved_flags;
2286 }
2287#endif
2288 IpplTimings::stopTimer(TPartInside_m);
2289 return ret;
2290}
2291
2293 std::ofstream of;
2294 of.open(fn.c_str());
2295 PAssert(of.is_open());
2296 of.precision(6);
2297 of << "# vtk DataFile Version 2.0" << std::endl;
2298 of << "generated using DataSink::writeGeoToVtk" << std::endl;
2299 of << "ASCII" << std::endl << std::endl;
2300 of << "DATASET UNSTRUCTURED_GRID" << std::endl;
2301 of << "POINTS " << Points_m.size() << " float" << std::endl;
2302 for (unsigned int i = 0; i < Points_m.size(); i++)
2303 of << Points_m[i](0) << " " << Points_m[i](1) << " " << Points_m[i](2) << std::endl;
2304 of << std::endl;
2305
2306 of << "CELLS " << Triangles_m.size() << " " << 4 * Triangles_m.size() << std::endl;
2307 for (size_t i = 0; i < Triangles_m.size(); i++)
2308 of << "3 " << PointID(i, 1) << " " << PointID(i, 2) << " " << PointID(i, 3) << std::endl;
2309 of << "CELL_TYPES " << Triangles_m.size() << std::endl;
2310 for (size_t i = 0; i < Triangles_m.size(); i++)
2311 of << "5" << std::endl;
2312 of << "CELL_DATA " << Triangles_m.size() << std::endl;
2313 of << "SCALARS "
2314 << "cell_attribute_data"
2315 << " float "
2316 << "1" << std::endl;
2317 of << "LOOKUP_TABLE "
2318 << "default" << std::endl;
2319 for (size_t i = 0; i < Triangles_m.size(); i++)
2320 of << (float)(i) << std::endl;
2321 of << std::endl;
2322}
2323
2324Inform& BoundaryGeometry::printInfo(Inform& os) const {
2325 os << endl;
2326 os << "* ************* B O U N D A R Y G E O M E T R Y *********************************** "
2327 << endl;
2328 os << "* GEOMETRY " << getOpalName() << '\n'
2329 << "* FGEOM '" << Attributes::getString(itsAttr[FGEOM]) << "'\n"
2330 << "* TOPO " << Attributes::getString(itsAttr[TOPO]) << '\n'
2331 << "* XSCALE " << Attributes::getReal(itsAttr[XSCALE]) << '\n'
2332 << "* YSCALE " << Attributes::getReal(itsAttr[YSCALE]) << '\n'
2333 << "* ZSCALE " << Attributes::getReal(itsAttr[ZSCALE]) << '\n'
2334 << "* XYZSCALE " << Attributes::getReal(itsAttr[XYZSCALE]) << '\n'
2335 << "* LENGTH " << Attributes::getReal(itsAttr[LENGTH]) << '\n'
2336 << "* S " << Attributes::getReal(itsAttr[S]) << '\n'
2337 << "* A " << Attributes::getReal(itsAttr[A]) << '\n'
2338 << "* B " << Attributes::getReal(itsAttr[B]) << '\n';
2340 os << "* C " << Attributes::getReal(itsAttr[C]) << '\n'
2341 << "* L1 " << Attributes::getReal(itsAttr[L1]) << '\n'
2342 << "* L2 " << Attributes::getReal(itsAttr[L2]) << '\n';
2343 }
2344 /* ADA os << "* Total triangle num " << Triangles_m.size() << '\n'
2345 << "* Total points num " << Points_m.size() << '\n'
2346 << "* Geometry bounds(m) Max = " << maxExtent_m << '\n'
2347 << "* Min = " << minExtent_m << '\n'
2348 << "* Geometry length(m) " << maxExtent_m - minExtent_m << '\n'
2349 << "* Resolution of voxel mesh " << voxelMesh_m.nr_m << '\n'
2350 << "* Size of voxel " << voxelMesh_m.sizeOfVoxel << '\n'
2351 << "* Number of voxels in mesh " << voxelMesh_m.ids.size() << endl;
2352 */
2353 os << "* ********************************************************************************** "
2354 << endl;
2355 return os;
2356}
Vector3D cross(const Vector3D &lhs, const Vector3D &rhs)
Vector cross product.
Definition Vector3D.cpp:111
double dot(const Vector3D &lhs, const Vector3D &rhs)
Vector dot product.
Definition Vector3D.cpp:118
Inform * gmsg
Definition changes.cpp:7
ippl::Vector< T, Dim > Vector_t
const int nr
#define INSIDE
constexpr double EPS
#define PointID(triangle_id, vertex_id)
#define LERP(a, b, t)
#define FUNC_LE(x, y)
#define FUNC_GE(x, y)
#define FUNC_LT_ZERO(x)
#define FUNC_EQ(x, y)
#define OUTSIDE
#define FUNC_GT(x, y)
#define FUNC_EQ_ZERO(x)
#define SQR(x)
#define mapPoint2VoxelIndices(pt, i, j, k)
#define FUNC_LE_ZERO(x)
#define FUNC_GT_ZERO(x)
#define FUNC_GE_ZERO(x)
#define FUNC_LT(x, y)
double getReal(const Attribute &attr)
Return real value.
Attribute makePredefinedString(const std::string &name, const std::string &help, const std::initializer_list< std::string > &predefinedStrings)
Make predefined string attribute.
Attribute makeReal(const std::string &name, const std::string &help)
Make real attribute.
Attribute makeRealArray(const std::string &name, const std::string &help)
Create real array attribute.
std::vector< double > getRealArray(const Attribute &attr)
Get array value.
std::string getString(const Attribute &attr)
Get string value.
Attribute makeString(const std::string &name, const std::string &help)
Make string attribute.
constexpr double c
The velocity of light in m/s.
Definition Physics.h:45
bool almost_eq(double A, double B, double maxDiff=1e-15, double maxRelDiff=DBL_EPSILON)
bool almost_eq_zero(double A, double maxDiff=1e-15)
bool almost_eq(double A, double B, double maxDiff=1e-20, int maxUlps=1000)
bool almost_eq_zero(double A, double maxDiff=1e-15)
bool almost_eq(double A, double B, double maxDiff=1e-20, int maxUlps=1000)
bool ge_zero(double x)
bool ge(double x, double y)
bool gt(double x, double y)
bool gt_zero(double x)
bool lt_zero(double x)
bool eq_zero(double x)
bool almost_eq_zero(double A, double maxDiff=1e-15)
bool lt(double x, double y)
bool le(double x, double y)
bool enableVTK
If true VTK files are written.
Definition Options.cpp:83
std::string combineFilePath(std::initializer_list< std::string > ilist)
Definition Util.cpp:197
Definition(int size, const char *name, const char *help)
Constructor for exemplars.
void registerOwnership(const AttributeHandler::OwnerType &itsClass) const
Definition Object.cpp:189
const std::string & getOpalName() const
Return object name.
Definition Object.cpp:308
Object(int size, const char *name, const char *help)
Constructor for exemplars.
Definition Object.cpp:354
void setOpalName(const std::string &name)
Set object name.
Definition Object.cpp:329
std::vector< Attribute > itsAttr
The object attributes.
Definition Object.h:216
bool builtin
Built-in flag.
Definition Object.h:233
Object * find(const std::string &name)
Find entry.
Definition OpalData.cpp:563
static OpalData * getInstance()
Definition OpalData.cpp:195
void define(Object *newObject)
Define a new object.
Definition OpalData.cpp:486
std::string getAuxiliaryOutputDirectory() const
get the name of the the additional data directory
Definition OpalData.cpp:677
Abstract base class for accelerator geometry classes.
Definition Geometry.h:43
const Vector_t< double, 3 > & v1() const
double v1(int i) const
Vector_t< double, 3 > pts[3]
Triangle(const Vector_t< double, 3 > &v1, const Vector_t< double, 3 > &v2, const Vector_t< double, 3 > &v3)
double v2(int i) const
const Vector_t< double, 3 > & v2() const
const Vector_t< double, 3 > & v3() const
void scale(const Vector_t< double, 3 > &scaleby, const Vector_t< double, 3 > &shiftby)
double v3(int i) const
Ray(const Ray &r)
Vector_t< double, 3 > inv_direction
Vector_t< double, 3 > direction
Ray(Vector_t< double, 3 > o, Vector_t< double, 3 > d)
const Ray & operator=(const Ray &a)=delete
int sign[3]
Vector_t< double, 3 > origin
bool intersect(const Ray &r, double &tmin, double &tmax) const
bool isInside(const Vector_t< double, 3 > &P) const
Vector_t< double, 3 > extent() const
Voxel(const Vector_t< double, 3 > &min, const Vector_t< double, 3 > &max)
void scale(const Vector_t< double, 3 > &scale)
int intersect(const Triangle &t) const
Vector_t< double, 3 > pts[2]
bool intersect(const Ray &r) const
IpplTimings::TimerRef TisInside_m
int fastIsInside(const Vector_t< double, 3 > &reference_pt, const Vector_t< double, 3 > &P)
Vector_t< double, 3 > mapIndices2Voxel(const int, const int, const int)
std::vector< double > TriAreas_m
virtual void execute()
Execute the command.
virtual void update()
Update this object.
IpplTimings::TimerRef TRayTrace_m
struct BoundaryGeometry::@117056140100230112133004050263065356254101277345 voxelMesh_m
virtual bool canReplaceBy(Object *object)
Test if replacement is allowed.
Vector_t< double, 3 > maxExtent_m
const Vector_t< double, 3 > & getPoint(const int triangle_id, const int vertex_id)
std::vector< std::array< unsigned int, 4 > > Triangles_m
std::string h5FileName_m
int intersectLineTriangle(const enum INTERSECTION_TESTS kind, const Vector_t< double, 3 > &P0, const Vector_t< double, 3 > &P1, const int triangle_id, Vector_t< double, 3 > &I)
static BoundaryGeometry * find(const std::string &name)
IpplTimings::TimerRef Tinitialize_m
int mapVoxelIndices2ID(const int i, const int j, const int k)
int intersectRayBoundary(const Vector_t< double, 3 > &P, const Vector_t< double, 3 > &v, Vector_t< double, 3 > &I)
IpplTimings::TimerRef TPartInside_m
IpplTimings::TimerRef TfastIsInside_m
int intersectTriangleVoxel(const int triangle_id, const int i, const int j, const int k)
int intersectLineSegmentBoundary(const Vector_t< double, 3 > &P0, const Vector_t< double, 3 > &P1, Vector_t< double, 3 > &intersection_pt, int &triangle_id)
Vector_t< double, 3 > insidePoint_m
int partInside(const Vector_t< double, 3 > &r, const Vector_t< double, 3 > &v, const double dt, Vector_t< double, 3 > &intecoords, int &triId)
int intersectTinyLineSegmentBoundary(const Vector_t< double, 3 > &, const Vector_t< double, 3 > &, Vector_t< double, 3 > &, int &)
void writeGeomToVtk(std::string fn)
virtual BoundaryGeometry * clone(const std::string &name)
Return a clone.
bool isInside(const Vector_t< double, 3 > &P)
Topology getTopology() const
Inform & printInfo(Inform &os) const
Vector_t< double, 3 > mapPoint2Voxel(const Vector_t< double, 3 > &)
void computeMeshVoxelization(void)
std::vector< Vector_t< double, 3 > > TriNormals_m
std::vector< Vector_t< double, 3 > > Points_m
Vector_t< double, 3 > minExtent_m
void updateElement(ElementBase *element)
The base class for all OPAL exceptions.