summaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorAnthony LaTorre <telatorre@gmail.com>2011-05-09 13:58:12 -0400
committerAnthony LaTorre <telatorre@gmail.com>2011-05-09 13:58:12 -0400
commit7e61bfbe7df445ff43abfc802df9471cf66b55ca (patch)
treebe4dc6797a87b64db60c4d9d53c9e63e62e7145c /src
parentbcaef46bb56feb2f92c4feae1bff9e041a4f84cf (diff)
downloadchroma-7e61bfbe7df445ff43abfc802df9471cf66b55ca.tar.gz
chroma-7e61bfbe7df445ff43abfc802df9471cf66b55ca.tar.bz2
chroma-7e61bfbe7df445ff43abfc802df9471cf66b55ca.zip
improve triangle intersection algorithm by allowing it to terminate early
Diffstat (limited to 'src')
-rw-r--r--src/intersect.cu38
1 files changed, 19 insertions, 19 deletions
diff --git a/src/intersect.cu b/src/intersect.cu
index 1c157dc..f3c5657 100644
--- a/src/intersect.cu
+++ b/src/intersect.cu
@@ -1,36 +1,36 @@
//-*-c-*-
-__device__ Matrix inv(const Matrix&m, float& determinant)
-{
- determinant = det(m);
-
- return make_matrix(m.a11*m.a22 - m.a12*m.a21,
- m.a02*m.a21 - m.a01*m.a22,
- m.a01*m.a12 - m.a02*m.a11,
- m.a12*m.a20 - m.a10*m.a22,
- m.a00*m.a22 - m.a02*m.a20,
- m.a02*m.a10 - m.a00*m.a12,
- m.a10*m.a21 - m.a11*m.a20,
- m.a01*m.a20 - m.a00*m.a21,
- m.a00*m.a11 - m.a01*m.a10)/determinant;
-}
-
__device__ bool intersect_triangle(const float3 &x, const float3 &p, float3 *triangle, float3 &intersection)
{
float3 v0 = triangle[0];
float3 v1 = triangle[1];
float3 v2 = triangle[2];
- float determinant;
- float3 u = inv(make_matrix(v1-v0,v2-v0,-p), determinant)*(x-v0);
+ Matrix m = make_matrix(v1-v0, v2-v0, -p);
+
+ float determinant = det(m);
if (determinant == 0.0)
return false;
- if (u.x < 0.0 || u.y < 0.0 || u.z < 0.0 || (1-u.x-u.y) < 0.0)
+ float3 b = x-v0;
+
+ float u1 = ((m.a11*m.a22 - m.a12*m.a21)*b.x + (m.a02*m.a21 - m.a01*m.a22)*b.y + (m.a01*m.a12 - m.a02*m.a11)*b.z)/determinant;
+
+ if (u1 < 0.0)
+ return false;
+
+ float u2 = ((m.a12*m.a20 - m.a10*m.a22)*b.x + (m.a00*m.a22 - m.a02*m.a20)*b.y + (m.a02*m.a10 - m.a00*m.a12)*b.z)/determinant;
+
+ if (u2 < 0.0)
+ return false;
+
+ float u3 = ((m.a10*m.a21 - m.a11*m.a20)*b.x + (m.a01*m.a20 - m.a00*m.a21)*b.y + (m.a00*m.a11 - m.a01*m.a10)*b.z)/determinant;
+
+ if (u3 < 0.0 || (1-u1-u2) < 0.0)
return false;
- intersection = x + p*u.z;
+ intersection = x + p*u3;
return true;
}