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