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