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