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