Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Rocon/src/ray-triangle.C
Go to the documentation of this file.
1 /* Ray-Triangle Intersection Test Routines */
2 /* Different optimizations of my and Ben Trumbore's */
3 /* code from journals of graphics tools (JGT) */
4 /* http://www.acm.org/jgt/ */
5 /* by Tomas Moller, May 2000 */
6 
7 using namespace std;
8 
9 #include <cmath>
10 
11 #include "NVec.h"
12 using namespace nvc;
13 
14 #define EPSILON 1e-24
15 
16 #define CROSS(dest,v1,v2) \
17  dest[0]=v1[1]*v2[2]-v1[2]*v2[1]; \
18  dest[1]=v1[2]*v2[0]-v1[0]*v2[2]; \
19  dest[2]=v1[0]*v2[1]-v1[1]*v2[0];
20 #define DOT(v1,v2) (v1[0]*v2[0]+v1[1]*v2[1]+v1[2]*v2[2])
21 #define SUB(dest,v1,v2) \
22  dest[0]=v1[0]-v2[0]; \
23  dest[1]=v1[1]-v2[1]; \
24  dest[2]=v1[2]-v2[2];
25 #define PLUCKER(pline, pt1, pt2) \
26  pline[0] = pt1[0]*pt2[1] - pt2[0]*pt1[1]; \
27  pline[1] = pt1[0]*pt2[2] - pt2[0]*pt1[2]; \
28  pline[2] = pt1[0] - pt2[0]; \
29  pline[3] = pt1[1]*pt2[2] - pt2[1]*pt1[2]; \
30  pline[4] = pt1[2] - pt2[2]; \
31  pline[5] = pt2[1] - pt1[1];
32 #define SIDE(p,q) (p[0]*q[4]+p[1]*q[5]+p[2]*q[3]+p[3]*q[2]+p[4]*q[0]+p[5]*q[1])
33 #define RAYPOINT(result,start,end,dist) { result[0]=start[0]+d*(end[0]-start[0]); result[1]=start[1]+d*(end[1]-start[1]); result[2]=start[2]+d*(end[2]-start[2]); }
34 
35 /* Sample implementation of the line segment—triangle intersection test
36  * presented in the paper:
37  *
38  * Fast 3D Line Segment—Triangle Intersection Test
39  * Nick Chirkov
40  * journal of graphics tools 10(3):13-18, 2005
41  *
42  */
43 
44 struct RAYTRI
45 {
46  double org[3];
47  double end[3];
48  double dir[3];
49  double v0[3],v1[3],v2[3];
50 
51  struct PLANE
52  {
53  double x, y, z, d;
54  enum MAIN_AXIS { X, Y, Z };
56  };
58 };
59 
60 //Plucker Co-ordinates
61 int P_triRay(double *p, double *d, double *v0,
62  double *v1, double *v2, Vec3D & I)
63 {
64  int flag = 0;
65  double *q,*L,*e1,*e2,*e3, s1, s2,s3;
66  q = new double[3];
67  for(int i=0;i<3;++i)q[i] = p[i]+d[i];
68 
69  L = new double[6];
70  e1 = new double[6];
71  e2 = new double[6];
72  e3 = new double[6];
73 
74  PLUCKER(L,p,q);
75  PLUCKER(e1,v0,v1);
76  PLUCKER(e2,v1,v2);
77  PLUCKER(e3,v2,v0);
78 
79  s1 = SIDE(L,e1);
80  s2 = SIDE(L,e2);
81  s3 = SIDE(L,e3);
82 
83  if(s1==0.0 && s2==0.0 && s3==0.0)
84  { // co-planar
85  int pts[4][3] = {{0,0,0},{0,1,0},{0,0,1},{1,0,0}};
86  int edge[3] = {0,0,0};
87  int vertex[3] = {0,0,0};
88  int *x,i,intersect=0;
89  double cs1,cs2,cs3;
90  for(i=0;i<4;++i)
91  {
92  PLUCKER(e1,v0,v1);
93  PLUCKER(e2,v1,pts[i]);
94  PLUCKER(e3,pts[i],v0);
95  cs1 = SIDE(L,e1);
96  cs2 = SIDE(L,e2);
97  cs3 = SIDE(L,e3);
98  if(!(cs1==0 && cs2==0 && cs3==0)) break;
99  }
100 
101  x = pts[i];
102 
103  for(i=0;i<3;++i)
104  {
105  if(i==0){
106  PLUCKER(e2,v1,x);
107  PLUCKER(e3,x,v0);
108  cs2 = SIDE(L,e2);
109  cs3 = SIDE(L,e3);
110  } else if(i==1){
111  PLUCKER(e2,v2,x);
112  PLUCKER(e3,x,v1);
113  cs2 = SIDE(L,e2);
114  cs3 = SIDE(L,e3);
115  } else {
116  PLUCKER(e2,v0,x);
117  PLUCKER(e3,x,v2);
118  cs2 = SIDE(L,e2);
119  cs3 = SIDE(L,e3);
120  }
121 
122  if(cs2*cs3<0);
123  else if(cs2*cs3>0)
124  { intersect++;
125  edge[i]++;
126  }
127  else if(cs2==0 && cs3==0)
128  { intersect+=2;
129  edge[i]+=2;
130  break;}
131  else if (cs2==0)
132  {intersect++;
133  vertex[(i+1)-3*((i+1)/3)]++;
134  }
135  else if (cs3==0)
136  {intersect++;
137  vertex[i]++;
138  }
139  }
140  // implement line-line intersect
141  // implement vertex-line intersect
142  if(intersect!=0)
143  {
144  if(intersect<=3){
145  if(vertex[0]!=0 && !flag)
146  {
147  I[0] = v0[0]; I[1]=v0[1]; I[2]=v0[2];
148  flag = 1;
149  }
150  if(vertex[1]!=0 && !flag)
151  {
152  I[0] = v1[0]; I[1]=v1[1]; I[2]=v1[2];
153  flag = 1;
154  }
155  if(vertex[2]!=0 && !flag)
156  {
157  I[0] = v2[0]; I[1]=v2[1]; I[2]=v2[2];
158  flag = 1;
159  }
160  }
161  if(!flag && (edge[0]==1||edge[0]==2) )
162  {
163  I[0]=(v0[0]+v1[0])/2;I[1]=(v0[1]+v1[1])/2;I[2]=(v0[2]+v1[2])/2;
164  flag = 1;
165  }
166  if(!flag && (edge[1]==1||edge[1]==2) )
167  {
168  I[0]=(v1[0]+v2[0])/2;I[1]=(v1[1]+v2[1])/2;I[2]=(v1[2]+v2[2])/2;
169  flag = 1;
170  }
171  if(!flag && (edge[2]==1||edge[2]==2) )
172  {
173  I[0]=(v2[0]+v0[0])/2;I[1]=(v0[1]+v2[1])/2;I[2]=(v0[2]+v2[2])/2;
174  flag = 1;
175  }
176  }
177  } else if ((s1>0 && s2>0 && s3>0) || (s1<0 && s2<0 && s3<0)){
178  double al,be,ga;
179  al = s1/(s1+s2+s3);
180  be = s2/(s1+s2+s3);
181  ga = s3/(s1+s2+s3);
182  for(int i=0;i<3;++i)
183  I[i] = al*v2[i] + be*v0[i] + ga*v1[i];
184  flag = 1;
185  } else if ( s1==0 && s2*s3>0){
186  I[0]=(v0[0]+v1[0])/2;I[1]=(v0[1]+v1[1])/2;I[2]=(v0[2]+v1[2])/2;
187  flag = 1;
188  } else if ( s2==0 && s1*s3>0){
189  I[0]=(v1[0]+v2[0])/2;I[1]=(v1[1]+v2[1])/2;I[2]=(v1[2]+v2[2])/2;
190  flag = 1;
191  } else if ( s3==0 && s1*s2>0){
192  I[0]=(v2[0]+v0[0])/2;I[1]=(v0[1]+v2[1])/2;I[2]=(v0[2]+v2[2])/2;
193  flag = 1;
194  } else if ( s1==0 && s2==0){
195  I[0]=v1[0];I[1]=v1[1];I[2]=v1[2];
196  flag = 1;
197  } else if ( s1==0 && s3==0){
198  I[0]=v0[0];I[1]=v0[1];I[2]=v0[2];
199  flag = 1;
200  } else if ( s2==0 && s3==0){
201  I[0]=v2[0];I[1]=v2[1];I[2]=v2[2];
202  flag = 1;
203  } else
204  flag = 0;
205 
206  delete [] q;
207  delete [] L;
208  delete [] e1;
209  delete [] e2;
210  delete [] e3;
211 
212  return flag;
213 }
214 
215 int c2005(const RAYTRI* rt)
216 {
217  double e0[3],e1[3],e2[3],norm[3],point[3],v[3],av[3],vb[3],vc[3];
218  SUB(e0, rt->v1, rt->v0);
219  SUB(e1, rt->v2, rt->v0);
220  CROSS(norm,e0,e1);
221 
222  double pd = DOT(norm, rt->v0);
223 
224  double signSrc = DOT(norm, rt->org) -pd;
225  double signDst = DOT(norm, rt->end) -pd;
226 
227  if(signSrc*signDst > 0.0) return 0;
228 
229  double d = signSrc/(signSrc - signDst);
230 
231  RAYPOINT(point, rt->org, rt->end,d);
232  SUB(v, point, rt->v0);
233  CROSS(av,e0,v);
234  CROSS(vb,v,e1);
235 
236  if(DOT(av,vb) >= 0.0)
237  {
238  SUB(e2, rt->v1, rt->v2);
239  SUB(v, point, rt->v1);
240  CROSS(vc,v,e2);
241  if(DOT(av,vc) >= 0.0) return 1;
242  }
243  return 0;
244 }
245 
246 // Copyright 2001, softSurfer (www.softsurfer.com)
247 // This code may be freely used and modified for any purpose
248 // providing that this copyright notice is included with it.
249 // SoftSurfer makes no warranty for this code, and cannot be held
250 // liable for any real or imagined damage resulting from its use.
251 // Users of this code must verify correctness for their application.
252 
253 // Assume that classes are already given for the objects:
254 // Point and Vector with
255 // coordinates {double x, y, z;}
256 // operators for:
257 // == to test equality
258 // != to test inequality
259 // (Vector)0 = (0,0,0) (null vector)
260 // Point = Point ± Vector
261 // Vector = Point - Point
262 // Vector = Scalar * Vector (scalar product)
263 // Vector = Vector * Vector (cross product)
264 // Line and Ray and Segment with defining points {Point P0, P1;}
265 // (a Line is infinite, Rays and Segments start at P0)
266 // (a Ray extends beyond P1, but a Segment ends at P1)
267 // Plane with a point and a normal {Point V0; Vector n;}
268 // Triangle with defining vertices {Point V0, V1, V2;}
269 // Polyline and Polygon with n vertices {int n; Point *V;}
270 // (a Polygon has V[n]=V[0])
271 //===================================================================
272 
273 #define SMALL_NUM 1e-24 // anything that avoids division overflow
274 // dot product (3D) which allows vector operations in arguments
275 
276 // intersect_RayTriangle(): intersect a ray with a 3D triangle
277 // Input: a ray R, and a triangle T
278 // Output: *I = intersection point (when it exists)
279 // Return: -1 = triangle is degenerate (a segment or point)
280 // 0 = disjoint (no intersect)
281 // 1 = intersect in unique point I1
282 // 2 = are in the same plane
283 int intersect_ray_triangle( const Vec3D &RP0, const Vec3D & dir,Vec3D & TV0,
284  Vec3D & TV1,Vec3D & TV2, Vec3D & I )
285 {
286  Vec3D u, v, n; // triangle vectors
287  Vec3D w0, w; // ray vectors
288  double r, a, b; // params to calc ray-plane intersect
289 
290  // get triangle edge vectors and plane normal
291  u = TV1 - TV0;
292  v = TV2 - TV0;
293  n = cross(u,v); // cross product
294  if (n == Vec3D(0,0,0)) // triangle is degenerate
295  return 0; // do not deal with this case
296 
297  //dir = RP1 - RP0; // ray direction vector
298  w0 = RP0 - TV0;
299  a = -1.0*(n*w0);
300  b = n*dir;
301  if (fabs(b) < SMALL_NUM) { // ray is parallel to triangle plane
302  if (a != 0) // ray lies in triangle plane
303  //return 1;
304  //else
305  return 0; // ray disjoint from plane
306  }
307 
308  // get intersect point of ray with triangle plane
309  r = a / b;
310  if (r < 0.0) // ray goes away from triangle
311  return 0; // => no intersect
312  if (r > 1.0) // ray goes away from triangle
313  return 0; //change by prateek, 0 instead of 3
314  // for a segment, also test if (r > 1.0) => no intersect
315 
316  I = RP0 + (r * dir); // intersect point of ray and plane
317 
318  // is I inside T?
319  double uu, uv, vv, wu, wv, D;
320  uu = (u*u);
321  uv = (u*v);
322  vv = (v*v);
323  w = I - TV0;
324  wu = w*u;
325  wv = w*v;
326  D = (uv * uv) -( uu * vv);
327 
328  // get and test parametric coords
329  double s, t;
330  s = (uv * wv - vv * wu) / D;
331  if (s < 0.0 || s > 1.0) // I is outside T
332  return 0;
333  t = (uv * wu - uu * wv) / D;
334  if (t < 0.0 || (s + t) > 1.0) // I is outside T
335  return 0;
336 
337  return 1; // I is in T
338 }
339 
340 /*
341  *the original jgt code
342  *modified by prateek -- checking if 0<=t<=1 for line segment
343 */
344 int intersect_triangle(double orig[3], double dir[3],
345  double vert0[3], double vert1[3], double vert2[3],
346  double *t, double *u, double *v)
347 {
348  double edge1[3], edge2[3], tvec[3], pvec[3], qvec[3];
349  double det,inv_det;
350 
351  /* find vectors for two edges sharing vert0 */
352  SUB(edge1, vert1, vert0);
353  SUB(edge2, vert2, vert0);
354 
355  /* begin calculating determinant - also used to calculate U parameter */
356  CROSS(pvec, dir, edge2);
357 
358  /* if determinant is near zero, ray lies in plane of triangle */
359  det = DOT(edge1, pvec);
360 
361  if (det > -EPSILON && det < EPSILON)
362  return 0;
363  inv_det = 1.0 / det;
364 
365  /* calculate distance from vert0 to ray origin */
366  SUB(tvec, orig, vert0);
367 
368  /* calculate U parameter and test bounds */
369  *u = DOT(tvec, pvec) * inv_det;
370  if (*u < 0.0 || *u > 1.0)
371  return 0;
372 
373  /* prepare to test V parameter */
374  CROSS(qvec, tvec, edge1);
375 
376  /* calculate V parameter and test bounds */
377  *v = DOT(dir, qvec) * inv_det;
378  if (*v < 0.0 || *u + *v > 1.0)
379  return 0;
380 
381  /* calculate t, ray intersects triangle */
382  *t = DOT(edge2, qvec) * inv_det;
383 
384  if(*t<0 || *t>1)
385  return 0; //enforce that ray is line segment
386 
387  return 1;
388 }
389 
390 
391 /* code rewritten to do tests on the sign of the determinant */
392 /* the division is at the end in the code */
393 int intersect_triangle1(double orig[3], double dir[3],
394  double vert0[3], double vert1[3], double vert2[3],
395  double *t, double *u, double *v)
396 {
397  double edge1[3], edge2[3], tvec[3], pvec[3], qvec[3];
398  double det,inv_det;
399 
400  /* find vectors for two edges sharing vert0 */
401  SUB(edge1, vert1, vert0);
402  SUB(edge2, vert2, vert0);
403 
404  /* begin calculating determinant - also used to calculate U parameter */
405  CROSS(pvec, dir, edge2);
406 
407  /* if determinant is near zero, ray lies in plane of triangle */
408  det = DOT(edge1, pvec);
409 
410  if (det > EPSILON)
411  {
412  /* calculate distance from vert0 to ray origin */
413  SUB(tvec, orig, vert0);
414 
415  /* calculate U parameter and test bounds */
416  *u = DOT(tvec, pvec);
417  if (*u < 0.0 || *u > det)
418  return 0;
419 
420  /* prepare to test V parameter */
421  CROSS(qvec, tvec, edge1);
422 
423  /* calculate V parameter and test bounds */
424  *v = DOT(dir, qvec);
425  if (*v < 0.0 || *u + *v > det)
426  return 0;
427 
428  }
429  else if(det < -EPSILON)
430  {
431  /* calculate distance from vert0 to ray origin */
432  SUB(tvec, orig, vert0);
433 
434  /* calculate U parameter and test bounds */
435  *u = DOT(tvec, pvec);
436 /* printf("*u=%f\n",(double)*u); */
437 /* printf("det=%f\n",det); */
438  if (*u > 0.0 || *u < det)
439  return 0;
440 
441  /* prepare to test V parameter */
442  CROSS(qvec, tvec, edge1);
443 
444  /* calculate V parameter and test bounds */
445  *v = DOT(dir, qvec) ;
446  if (*v > 0.0 || *u + *v < det)
447  return 0;
448  }
449  else return 0; /* ray is parallell to the plane of the triangle */
450 
451 
452  inv_det = 1.0 / det;
453 
454  /* calculate t, ray intersects triangle */
455  *t = DOT(edge2, qvec) * inv_det;
456  (*u) *= inv_det;
457  (*v) *= inv_det;
458 
459  return 1;
460 }
461 
462 /* code rewritten to do tests on the sign of the determinant */
463 /* the division is before the test of the sign of the det */
464 int intersect_triangle2(double orig[3], double dir[3],
465  double vert0[3], double vert1[3], double vert2[3],
466  double *t, double *u, double *v)
467 {
468  double edge1[3], edge2[3], tvec[3], pvec[3], qvec[3];
469  double det,inv_det;
470 
471  /* find vectors for two edges sharing vert0 */
472  SUB(edge1, vert1, vert0);
473  SUB(edge2, vert2, vert0);
474 
475  /* begin calculating determinant - also used to calculate U parameter */
476  CROSS(pvec, dir, edge2);
477 
478  /* if determinant is near zero, ray lies in plane of triangle */
479  det = DOT(edge1, pvec);
480 
481  /* calculate distance from vert0 to ray origin */
482  SUB(tvec, orig, vert0);
483  inv_det = 1.0 / det;
484 
485  if (det > EPSILON)
486  {
487  /* calculate U parameter and test bounds */
488  *u = DOT(tvec, pvec);
489  if (*u < 0.0 || *u > det)
490  return 0;
491 
492  /* prepare to test V parameter */
493  CROSS(qvec, tvec, edge1);
494 
495  /* calculate V parameter and test bounds */
496  *v = DOT(dir, qvec);
497  if (*v < 0.0 || *u + *v > det)
498  return 0;
499 
500  }
501  else if(det < -EPSILON)
502  {
503  /* calculate U parameter and test bounds */
504  *u = DOT(tvec, pvec);
505  if (*u > 0.0 || *u < det)
506  return 0;
507 
508  /* prepare to test V parameter */
509  CROSS(qvec, tvec, edge1);
510 
511  /* calculate V parameter and test bounds */
512  *v = DOT(dir, qvec) ;
513  if (*v > 0.0 || *u + *v < det)
514  return 0;
515  }
516  else return 0; /* ray is parallell to the plane of the triangle */
517 
518  /* calculate t, ray intersects triangle */
519  *t = DOT(edge2, qvec) * inv_det;
520  (*u) *= inv_det;
521  (*v) *= inv_det;
522 
523  return 1;
524 }
525 
526 /* code rewritten to do tests on the sign of the determinant */
527 /* the division is before the test of the sign of the det */
528 /* and one CROSS has been moved out from the if-else if-else */
529 int intersect_triangle3(double orig[3], double dir[3],
530  double vert0[3], double vert1[3], double vert2[3],
531  double *t, double *u, double *v)
532 {
533  double edge1[3], edge2[3], tvec[3], pvec[3], qvec[3];
534  double det,inv_det;
535 
536  /* find vectors for two edges sharing vert0 */
537  SUB(edge1, vert1, vert0);
538  SUB(edge2, vert2, vert0);
539 
540  /* begin calculating determinant - also used to calculate U parameter */
541  CROSS(pvec, dir, edge2);
542 
543  /* if determinant is near zero, ray lies in plane of triangle */
544  det = DOT(edge1, pvec);
545 
546  /* calculate distance from vert0 to ray origin */
547  SUB(tvec, orig, vert0);
548  inv_det = 1.0 / det;
549 
550  CROSS(qvec, tvec, edge1);
551 
552  if (det > EPSILON)
553  {
554  *u = DOT(tvec, pvec);
555  if (*u < 0.0 || *u > det)
556  return 0;
557 
558  /* calculate V parameter and test bounds */
559  *v = DOT(dir, qvec);
560  if (*v < 0.0 || *u + *v > det)
561  return 0;
562 
563  }
564  else if(det < -EPSILON)
565  {
566  /* calculate U parameter and test bounds */
567  *u = DOT(tvec, pvec);
568  if (*u > 0.0 || *u < det)
569  return 0;
570 
571  /* calculate V parameter and test bounds */
572  *v = DOT(dir, qvec) ;
573  if (*v > 0.0 || *u + *v < det)
574  return 0;
575  }
576  else return 0; /* ray is parallell to the plane of the triangle */
577 
578  *t = DOT(edge2, qvec) * inv_det;
579  (*u) *= inv_det;
580  (*v) *= inv_det;
581 
582  return 1;
583 }
584 
585 /* Wrapper for all other methods ? */
586 bool rayIntersectsTriangle(double *p, double *d,
587  double *v0, double *v1, double *v2, Vec3D & I)
588 {
589  /*Vec3D org(p[0],p[1],p[2]);
590  Vec3D dir(d[0],d[1],d[2]);
591  Vec3D ver0(v0[0],v0[1],v0[2]);
592  Vec3D ver1(v1[0],v1[1],v1[2]);
593  Vec3D ver2(v2[0],v2[1],v2[2]);
594  int result = intersect_ray_triangle(org,dir,ver0,ver1,ver2,I);*/
595 
596  /*RAYTRI * input = new RAYTRI;
597  for(int i =0;i<3;++i)
598  {
599  input->org[i] = p[i];
600  input->dir[i] = d[i];
601  input->v0[i] = v0[i];
602  input->v1[i] = v1[i];
603  input->v2[i] = v2[i];
604  input->end[i] = p[i]+d[i];
605  }
606  int result = c2005(input);*/
607 
608  /*double *t,*u,*v;
609  t = new double;
610  u = new double;
611  v = new double;
612  *t = *u = *v = 0;
613  int result = intersect_triangle(p,d,v0,v1,v2,t,u,v);
614  for(int j=0;j<3;j++)
615  {
616  I[j]=p[j]+ (*t)*d[j];//I[j]=v0[j]+ (*u)*v1[j] + (*v)*v2[j];
617  }*/
618 
619  int result = P_triRay(p,d,v0,v1,v2,I);
620 
621  if(result)
622  { double *t,*u,*v;
623  t = new double;
624  u = new double;
625  v = new double;
626  *t = *u = *v = 0;
627  result = intersect_triangle(p,d,v0,v1,v2,t,u,v);
628  for(int j=0;j<3;j++)
629  {
630  I[j]=p[j]+ (*t)*d[j];
631  }
632  delete t;
633  delete u;
634  delete v;
635  }
636  return result;
637 }
#define SUB(dest, v1, v2)
NVec< 3, T > cross(const NVec< 3, T > &u, const NVec< 3, T > &v)
const NT & d
void int int REAL REAL * y
Definition: read.cpp:74
double s
Definition: blastest.C:80
void intersect(T_VertexSet s1, T_VertexSet s2, T_VertexSet *inter)
int intersect_triangle3(double orig[3], double dir[3], double vert0[3], double vert1[3], double vert2[3], double *t, double *u, double *v)
int intersect_ray_triangle(const Vec3D &RP0, const Vec3D &dir, Vec3D &TV0, Vec3D &TV1, Vec3D &TV2, Vec3D &I)
int intersect_triangle2(double orig[3], double dir[3], double vert0[3], double vert1[3], double vert2[3], double *t, double *u, double *v)
NVec< 3, double > Vec3D
T norm(const NVec< DIM, T > &v)
*********************************************************************Illinois Open Source License ****University of Illinois NCSA **Open Source License University of Illinois All rights reserved ****Developed free of to any person **obtaining a copy of this software and associated documentation to deal with the Software without including without limitation the rights to and or **sell copies of the and to permit persons to whom the **Software is furnished to do subject to the following this list of conditions and the following disclaimers ****Redistributions in binary form must reproduce the above **copyright this list of conditions and the following **disclaimers in the documentation and or other materials **provided with the distribution ****Neither the names of the Center for Simulation of Advanced the University of nor the names of its **contributors may be used to endorse or promote products derived **from this Software without specific prior written permission ****THE SOFTWARE IS PROVIDED AS WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **ARISING OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE v
Definition: roccomf90.h:20
#define DOT(v1, v2)
blockLoc i
Definition: read.cpp:79
void int int REAL * x
Definition: read.cpp:74
const NT & n
bool rayIntersectsTriangle(double *p, double *d, double *v0, double *v1, double *v2, Vec3D &I)
double det(const Matrix3D &A)
#define PLUCKER(pline, pt1, pt2)
j indices j
Definition: Indexing.h:6
NT q
#define RAYPOINT(result, start, end, dist)
#define SIDE(p, q)
#define SMALL_NUM
#define EPSILON
#define CROSS(dest, v1, v2)
int P_triRay(double *p, double *d, double *v0, double *v1, double *v2, Vec3D &I)
int intersect_triangle1(double orig[3], double dir[3], double vert0[3], double vert1[3], double vert2[3], double *t, double *u, double *v)
int intersect_triangle(double orig[3], double dir[3], double vert0[3], double vert1[3], double vert2[3], double *t, double *u, double *v)
int c2005(const RAYTRI *rt)