Browse Source

* Updated embree interface to version 2.0

Former-commit-id: 167f46d80443f3ea1990f849ba26f2eb393e8154
schuellc 11 years ago
parent
commit
a077cad5f3

+ 12 - 11
examples/embree/example.cpp

@@ -1,3 +1,4 @@
+#define IGL_HEADER_ONLY
 #include <igl/OpenGL_convenience.h>
 #include <igl/per_face_normals.h>
 #include <igl/read.h>
@@ -37,9 +38,9 @@ double bbd;
 // Faces
 Eigen::MatrixXi F;
 // Embree intersection structure
-igl::EmbreeIntersector<Eigen::MatrixXd,Eigen::MatrixXi,Eigen::Vector3d> ei;
+igl::EmbreeIntersector<double,int> ei;
 // Hits collected
-std::vector<embree::Hit > hits;
+std::vector<igl::Hit > hits;
 // Ray information, "projection screen" corners
 Eigen::Vector3d win_s,s,d,dir,NW,NE,SE,SW;
 // Textures and framebuffers for "projection screen"
@@ -213,7 +214,7 @@ void display()
   // Draw all hits
   glBegin(GL_POINTS);
   glColor3f(0,0.2,0.2);
-  for(vector<embree::Hit>::iterator hit = hits.begin();
+  for(vector<igl::Hit>::iterator hit = hits.begin();
       hit != hits.end();
       hit++)
   {
@@ -221,9 +222,9 @@ void display()
     const double w1 = hit->u;
     const double w2 = hit->v;
     VectorXd hitP = 
-      w0 * V.row(F(hit->id0,0)) + 
-      w1 * V.row(F(hit->id0,1)) + 
-      w2 * V.row(F(hit->id0,2));
+      w0 * V.row(F(hit->id,0)) + 
+      w1 * V.row(F(hit->id,1)) + 
+      w2 * V.row(F(hit->id,2));
     glVertex3dv(hitP.data());
   }
   glEnd();
@@ -327,14 +328,14 @@ void mouse_move(int mouse_x, int mouse_y)
   dir = d-s;
   int num_rays_shot;
   ei.intersectRay(s,dir,hits,num_rays_shot);
-  for(vector<embree::Hit>::iterator hit = hits.begin();
+  for(vector<igl::Hit>::iterator hit = hits.begin();
       hit != hits.end();
       hit++)
   {
     // Change color of hit faces
-    C(hit->id0,0) = 1;
-    C(hit->id0,1) = 0.4;
-    C(hit->id0,2) = 0.4;
+    C(hit->id,0) = 1;
+    C(hit->id,1) = 0.4;
+    C(hit->id,2) = 0.4;
   }
 }
 
@@ -463,7 +464,7 @@ int main(int argc, char * argv[])
     V.colwise().minCoeff()).maxCoeff();
 
   // Init embree
-  ei = EmbreeIntersector<MatrixXd,MatrixXi,Vector3d>(V,F);
+  ei = EmbreeIntersector<double,int>(V,F);
 
   // Init glut
   glutInit(&argc,argv);

+ 0 - 314
include/igl/embree/Embree2Intersector.h

@@ -1,314 +0,0 @@
-#ifndef IGL_EMBREE_INTERSECTOR_H
-#define IGL_EMBREE_INTERSECTOR_H
-
-#include "include/embree.h"
-#include "include/intersector1.h"
-#include "common/ray.h"
-#include <vector>
-
-namespace igl
-{
-  struct Hit
-  {
-    int id; // primitive id
-    float u,v; // barycentric coordinates
-    float t; // distance = direction*t to intersection
-  };
-
-  template <
-  typename Scalar,
-  typename Index>
-  class EmbreeIntersector
-  {
-    typedef Eigen::Matrix<Scalar,Eigen::Dynamic,Eigen::Dynamic> PointMatrixType;
-    typedef Eigen::Matrix<Index,Eigen::Dynamic,Eigen::Dynamic>  FaceMatrixType;
-    typedef Eigen::Matrix<Scalar,1,3> RowVector3;
-  public:
-    // V  #V by 3 list of vertex positions
-    // F  #F by 3 list of Oriented triangles
-    EmbreeIntersector();
-    EmbreeIntersector(
-      const PointMatrixType & V,
-      const FaceMatrixType & F,
-      const char* structure = "default",
-      const char* builder = "default",
-      const char* traverser = "default");
-    virtual ~EmbreeIntersector();
-  
-    // Given a ray find the first hit
-    // 
-    // Inputs:
-    //   origin     3d origin point of ray
-    //   direction  3d (not necessarily normalized) direction vector of ray
-    // Output:
-    //   hit        information about hit
-    // Returns true if and only if there was a hit
-    bool intersectRay(
-      const RowVector3& origin, 
-      const RowVector3& direction,
-      Hit& hit,
-      float tnear = 0,
-      float tfar = embree::inf) const;
-
-    // Given a ray find the all hits in order
-    // 
-    // Inputs:
-    //   origin     3d origin point of ray
-    //   direction  3d (not necessarily normalized) direction vector of ray
-    // Output:
-    //   hit        information about hit
-    //   num_rays   number of rays shot (at least one)
-    // Returns true if and only if there was a hit
-    bool intersectRay(
-      const RowVector3& origin,
-      const RowVector3& direction,
-      std::vector<Hit > &hits,
-      int& num_rays,
-      float tnear = 0,
-      float tfar = embree::inf) const;
-
-    // Given a ray find the first hit
-    // 
-    // Inputs:
-    //   a    3d first end point of segment
-    //   ab   3d vector from a to other endpoint b
-    // Output:
-    //   hit  information about hit
-    // Returns true if and only if there was a hit
-    bool intersectSegment(const RowVector3& a, const RowVector3& ab, Hit &hit) const;
-    
-  private:
-    embree::RTCGeometry* mesh;
-    embree::RTCTriangle* triangles;
-    embree::RTCVertex *vertices;
-    embree::Intersector1 *intersector;
-  };
-}
-
-// Implementation
-#include <igl/EPS.h>
-
-template <typename RowVector3>
-inline embree::Vector3f toVector3f(const RowVector3 &p) { return embree::Vector3f((float)p[0], (float)p[1], (float)p[2]); }
-
-template <
-typename Scalar,
-typename Index>
-igl::EmbreeIntersector < Scalar, Index>
-::EmbreeIntersector()
-{
-  static bool inited = false;
-  if(!inited)
-  {
-    embree::rtcInit();
-    embree::rtcStartThreads();
-    inited = true;
-  }
-}
-
-template <
-typename Scalar,
-typename Index>
-igl::EmbreeIntersector < Scalar, Index>
-::EmbreeIntersector(const PointMatrixType & V,
-                    const FaceMatrixType & F,
-                    const char* structure,
-                    const char* builder,
-                    const char* traverser)
-{
-  static bool inited = false;
-  if(!inited)
-  {
-    embree::rtcInit();
-#ifdef VERBOSE
-    embree::rtcSetVerbose(3);
-#endif
-    embree::rtcStartThreads();
-    inited = true;
-  }
-  
-  mesh = embree::rtcNewTriangleMesh(F.rows(),V.rows(),structure);
-
-  // fill vertex buffer
-  vertices = embree::rtcMapPositionBuffer(mesh);
-  for(int i=0;i<(int)V.rows();i++)
-  {
-    vertices[i] = embree::RTCVertex((float)V(i,0),(float)V(i,1),(float)V(i,2));
-  }
-  embree::rtcUnmapPositionBuffer(mesh);
-
-  // fill triangle buffer
-  triangles = embree::rtcMapTriangleBuffer(mesh);
-  for(int i=0;i<(int)F.rows();i++)
-  {
-    triangles[i] = embree::RTCTriangle((int)F(i,0),(int)F(i,1),(int)F(i,2),i);
-  }
-  embree::rtcUnmapTriangleBuffer(mesh);
-
-  embree::rtcBuildAccel(mesh,builder);
-  embree::rtcCleanupGeometry(mesh);
-  
-  intersector = embree::rtcQueryIntersector1(mesh,traverser);
-}
-
-template <
-typename Scalar,
-typename Index>
-igl::EmbreeIntersector < Scalar, Index>
-::~EmbreeIntersector()
-{
-  embree::rtcDeleteIntersector1(intersector);
-  embree::rtcDeleteGeometry(mesh);
-//  embree::rtcStopThreads();
-//  embree::rtcExit();
-//  embree::rtcFreeMemory();
-}
-
-template <
-typename Scalar,
-typename Index>
-bool 
-igl::EmbreeIntersector< Scalar, Index>
-::intersectRay(
-  const RowVector3& origin,
-  const RowVector3& direction,
-  Hit& hit,
-  float tnear,
-  float tfar) const
-{
-  embree::Ray ray(toVector3f(origin), toVector3f(direction), tnear, tfar);
-  intersector->intersect(ray);
-  
-  if(ray)
-  {
-    hit.id = ray.id0;
-    hit.u = ray.u;
-    hit.v = ray.v;
-    hit.t = ray.tfar;
-    return true;
-  }
-
-  return false;
-}
-
-template <
-typename Scalar,
-typename Index>
-bool 
-igl::EmbreeIntersector < Scalar, Index>
-::intersectRay(
-  const RowVector3& origin, 
-  const RowVector3& direction,
-  std::vector<Hit > &hits,
-  int& num_rays,
-  float tnear,
-  float tfar) const
-{
-  using namespace std;
-  num_rays = 0;
-  hits.clear();
-  int last_id0 = -1;
-  double self_hits = 0;
-  // This epsilon is directly correleated to the number of missed hits, smaller
-  // means more accurate and slower
-  //const double eps = DOUBLE_EPS;
-  const double eps = FLOAT_EPS;
-  double min_t = tnear;
-  bool large_hits_warned = false;
-  embree::Ray ray(toVector3f(origin),toVector3f(direction));
-
-  while(true)
-  {
-    ray.tnear = min_t;
-    ray.tfar = tfar;
-    ray.id0 = -1;
-    num_rays++;
-    intersector->intersect(ray);
-    if(ray)
-    {
-      // Hit self again, progressively advance
-      if(ray.id0 == last_id0 || ray.tfar <= min_t)
-      {
-        // push min_t a bit more
-        //double t_push = pow(2.0,self_hits-4)*(hit.t<eps?eps:hit.t);
-        double t_push = pow(2.0,self_hits)*eps;
-#ifdef VERBOSE
-        cout<<"  t_push: "<<t_push<<endl;
-#endif
-        //o = o+t_push*d;
-        min_t += t_push;
-        self_hits++;
-      }
-      else
-      {
-        Hit hit;
-        hit.id = ray.id0;
-        hit.u = ray.u;
-        hit.v = ray.v;
-        hit.t = ray.tfar;
-        hits.push_back(hit);
-#ifdef VERBOSE
-        cout<<"  t: "<<hit.t<<endl;
-#endif
-        // Instead of moving origin, just change min_t. That way calculations
-        // all use exactly same origin values
-        min_t = ray.tfar;
-
-        // reset t_scale
-        self_hits = 0;
-      }
-      last_id0 = ray.id0;
-    }
-    else
-      break; // no more hits
-    
-    if(hits.size()>1000 && !large_hits_warned)
-    {
-      cerr<<"Warning: Large number of hits..."<<endl;
-      cerr<<"[ ";
-      for(vector<Hit>::iterator hit = hits.begin(); hit != hits.end();hit++)
-      {
-        cerr<<(hit->id+1)<<" ";
-      }
-      
-      cerr.precision(std::numeric_limits< double >::digits10);
-      cerr<<"[ ";
-      
-      for(vector<Hit>::iterator hit = hits.begin(); hit != hits.end(); hit++)
-      {
-        cerr<<(hit->t)<<endl;;
-      }
-
-      cerr<<"]"<<endl;
-      large_hits_warned = true;
-
-      return hits.empty();
-    }
-  }
-
-  return hits.empty();
-}
-
-template <
-typename Scalar,
-typename Index>
-bool 
-igl::EmbreeIntersector < Scalar, Index>
-::intersectSegment(const RowVector3& a, const RowVector3& ab, Hit &hit) const
-{
-  embree::Ray ray(toVector3f(a), toVector3f(ab), embree::zero, embree::one);
-  intersector->intersect(ray);
-
-  if(ray)
-  {
-    hit.id = ray.id0;
-    hit.u = ray.u;
-    hit.v = ray.v;
-    hit.t = ray.tfar;
-    return true;
-  }
-
-  return false;
-}
-
-#endif //EMBREE_INTERSECTOR_H

+ 184 - 147
include/igl/embree/EmbreeIntersector.h

@@ -1,196 +1,243 @@
+// igl function interface for Embree2.0
+//
+// Necessary changes to switch from previous Embree versions:
+// * Use igl:Hit instead of embree:Hit (where id0 -> id)
+// * Embree2.0 finds now also back face intersections
+
 #ifndef IGL_EMBREE_INTERSECTOR_H
 #define IGL_EMBREE_INTERSECTOR_H
-#include "Embree_convenience.h"
 
+#include "include/embree.h"
+#include "include/intersector1.h"
+#include "common/ray.h"
 #include <vector>
-//#include "types.h"
 
 namespace igl
 {
+  struct Hit
+  {
+    int id; // primitive id
+    float u,v; // barycentric coordinates
+    float t; // distance = direction*t to intersection
+  };
+
   template <
-  typename PointMatrixType,
-  typename FaceMatrixType,
-  typename RowVector3>
+  typename Scalar,
+  typename Index>
   class EmbreeIntersector
   {
+    typedef Eigen::Matrix<Scalar,Eigen::Dynamic,Eigen::Dynamic> PointMatrixType;
+    typedef Eigen::Matrix<Index,Eigen::Dynamic,Eigen::Dynamic>  FaceMatrixType;
+    typedef Eigen::Matrix<Scalar,1,3> RowVector3;
   public:
     // V  #V by 3 list of vertex positions
     // F  #F by 3 list of Oriented triangles
-    //
-    // Note: this will only find front-facing hits. To consider all hits then
-    // pass [F;fliplr(F)]
-    EmbreeIntersector(const PointMatrixType & V = PointMatrixType(), const FaceMatrixType &  F = FaceMatrixType());
+    EmbreeIntersector();
+    EmbreeIntersector(
+      const PointMatrixType & V,
+      const FaceMatrixType & F,
+      const char* structure = "default",
+      const char* builder = "default",
+      const char* traverser = "default");
     virtual ~EmbreeIntersector();
   
-    // Given a ray find the first *front-facing* hit
+    // Given a ray find the first hit
     // 
     // Inputs:
-    //   origin  3d origin point of ray
+    //   origin     3d origin point of ray
     //   direction  3d (not necessarily normalized) direction vector of ray
     // Output:
-    //   hit  embree information about hit
+    //   hit        information about hit
     // Returns true if and only if there was a hit
     bool intersectRay(
       const RowVector3& origin, 
       const RowVector3& direction,
-      embree::Hit &hit) const;
-    // Given a ray find the all *front-facing* hits in order
+      Hit& hit,
+      float tnear = 0,
+      float tfar = embree::inf) const;
+
+    // Given a ray find the all hits in order
     // 
     // Inputs:
-    //   origin  3d origin point of ray
+    //   origin     3d origin point of ray
     //   direction  3d (not necessarily normalized) direction vector of ray
     // Output:
-    //   hit  embree information about hit
-    //   num_rays  number of rays shot (at least one)
+    //   hit        information about hit
+    //   num_rays   number of rays shot (at least one)
     // Returns true if and only if there was a hit
     bool intersectRay(
-      const RowVector3& origin, 
-      const RowVector3& direction, 
-      std::vector<embree::Hit > &hits,
-      int & num_rays) const;
-  
-    // Given a ray find the first *front-facing* hit
+      const RowVector3& origin,
+      const RowVector3& direction,
+      std::vector<Hit > &hits,
+      int& num_rays,
+      float tnear = 0,
+      float tfar = embree::inf) const;
+
+    // Given a ray find the first hit
     // 
     // Inputs:
-    //   a  3d first end point of segment
-    //   ab  3d vector from a to other endpoint b
+    //   a    3d first end point of segment
+    //   ab   3d vector from a to other endpoint b
     // Output:
-    //   hit  embree information about hit
+    //   hit  information about hit
     // Returns true if and only if there was a hit
-    bool intersectSegment(const RowVector3& a, const RowVector3& ab, embree::Hit &hit) const;
-    //RowVector3 hit_position(const embree::Hit & hit) const;
+    bool intersectSegment(const RowVector3& a, const RowVector3& ab, Hit &hit) const;
     
   private:
-    embree::Ref<embree::Accel> _accel;
-    embree::Ref<embree::Intersector> _intersector;
+    embree::RTCGeometry* mesh;
+    embree::RTCTriangle* triangles;
+    embree::RTCVertex *vertices;
+    embree::Intersector1 *intersector;
   };
 }
 
 // Implementation
 #include <igl/EPS.h>
-#include <iostream>
-#include <sstream>
 
 template <typename RowVector3>
-inline embree::Vec3f toVec3f(const RowVector3 &p) { return embree::Vec3f((float)p[0], (float)p[1], (float)p[2]); }
+inline embree::Vector3f toVector3f(const RowVector3 &p) { return embree::Vector3f((float)p[0], (float)p[1], (float)p[2]); }
 
 template <
-typename PointMatrixType,
-typename FaceMatrixType,
-typename RowVector3>
-igl::EmbreeIntersector < PointMatrixType, FaceMatrixType, RowVector3>
-::EmbreeIntersector(const PointMatrixType & V, const FaceMatrixType & F)
-  :
-    _accel(),
-    _intersector()
+typename Scalar,
+typename Index>
+igl::EmbreeIntersector < Scalar, Index>
+::EmbreeIntersector()
 {
-  using namespace std;
   static bool inited = false;
   if(!inited)
   {
-	//embree::TaskScheduler::start();//init();
+    embree::rtcInit();
+#ifdef VERBOSE
+  embree::rtcSetVerbose(3);
+#endif
+    embree::rtcStartThreads();
     inited = true;
   }
-  
-  if(V.size() == 0 || F.size() == 0)
+}
+
+template <
+typename Scalar,
+typename Index>
+igl::EmbreeIntersector < Scalar, Index>
+::EmbreeIntersector(const PointMatrixType & V,
+                    const FaceMatrixType & F,
+                    const char* structure,
+                    const char* builder,
+                    const char* traverser)
+{
+  static bool inited = false;
+  if(!inited)
   {
-    return;
+    embree::rtcInit();
+#ifdef VERBOSE
+    embree::rtcSetVerbose(3);
+#endif
+    embree::rtcStartThreads();
+    inited = true;
   }
-  size_t numVertices = 0;
-  size_t numTriangles = 0;
   
-  embree::BuildTriangle *triangles;
-  embree::BuildVertex *vertices;
-  triangles = (embree::BuildTriangle*) embree::rtcMalloc(sizeof(embree::BuildTriangle) * F.rows());
-  vertices = (embree::BuildVertex*) embree::rtcMalloc(sizeof(embree::BuildVertex) * V.rows());
+  mesh = embree::rtcNewTriangleMesh(F.rows(),V.rows(),structure);
 
-  for(int i = 0; i < (int)V.rows(); ++i)
+  // fill vertex buffer
+  vertices = embree::rtcMapPositionBuffer(mesh);
+  for(int i=0;i<(int)V.rows();i++)
   {
-    vertices[numVertices++] = embree::BuildVertex((float)V(i,0),(float)V(i,1),(float)V(i,2));
+    vertices[i] = embree::RTCVertex((float)V(i,0),(float)V(i,1),(float)V(i,2));
   }
-  
-  for(int i = 0; i < (int)F.rows(); ++i)
+  embree::rtcUnmapPositionBuffer(mesh);
+
+  // fill triangle buffer
+  triangles = embree::rtcMapTriangleBuffer(mesh);
+  for(int i=0;i<(int)F.rows();i++)
   {
-    triangles[numTriangles++] = embree::BuildTriangle((int)F(i,0),(int)F(i,1),(int)F(i,2),i);
+    triangles[i] = embree::RTCTriangle((int)F(i,0),(int)F(i,1),(int)F(i,2),i);
   }
+  embree::rtcUnmapTriangleBuffer(mesh);
+
+  embree::rtcBuildAccel(mesh,builder);
+  embree::rtcCleanupGeometry(mesh);
   
-  // turn off verbose output by disabling cout
-  // http://stackoverflow.com/a/8246536/148668
-  streambuf *old = cout.rdbuf(); // <-- save        
-  stringstream ss;
-  cout.rdbuf (ss.rdbuf());       // <-- redirect
-  _accel = embree::rtcCreateAccel("default", "default", triangles, numTriangles, vertices, numVertices);
-  _intersector = _accel->queryInterface<embree::Intersector>();
-  cout.rdbuf (old);  
+  intersector = embree::rtcQueryIntersector1(mesh,traverser);
 }
 
 template <
-typename PointMatrixType,
-typename FaceMatrixType,
-typename RowVector3>
-igl::EmbreeIntersector < PointMatrixType, FaceMatrixType, RowVector3>
+typename Scalar,
+typename Index>
+igl::EmbreeIntersector < Scalar, Index>
 ::~EmbreeIntersector()
 {
-	embree::rtcFreeMemory();
+  embree::rtcDeleteIntersector1(intersector);
+  embree::rtcDeleteGeometry(mesh);
+//  embree::rtcStopThreads();
+//  embree::rtcExit();
+//  embree::rtcFreeMemory();
 }
 
 template <
-typename PointMatrixType,
-typename FaceMatrixType,
-typename RowVector3>
+typename Scalar,
+typename Index>
 bool 
-igl::EmbreeIntersector < PointMatrixType, FaceMatrixType, RowVector3>
-::intersectRay(const RowVector3& origin, const RowVector3& direction, embree::Hit &hit) const
+igl::EmbreeIntersector< Scalar, Index>
+::intersectRay(
+  const RowVector3& origin,
+  const RowVector3& direction,
+  Hit& hit,
+  float tnear,
+  float tfar) const
 {
-  embree::Ray ray(toVec3f(origin), toVec3f(direction), 1e-4f);
-  _intersector->intersect(ray, hit);
-  return hit ; 
+  embree::Ray ray(toVector3f(origin), toVector3f(direction), tnear, tfar);
+  intersector->intersect(ray);
+  
+  if(ray)
+  {
+    hit.id = ray.id0;
+    hit.u = ray.u;
+    hit.v = ray.v;
+    hit.t = ray.tfar;
+    return true;
+  }
+
+  return false;
 }
 
 template <
-typename PointMatrixType,
-typename FaceMatrixType,
-typename RowVector3>
+typename Scalar,
+typename Index>
 bool 
-igl::EmbreeIntersector < PointMatrixType, FaceMatrixType, RowVector3>
+igl::EmbreeIntersector < Scalar, Index>
 ::intersectRay(
   const RowVector3& origin, 
-  const RowVector3& direction, 
-  std::vector<embree::Hit > &hits,
-  int & num_rays) const
+  const RowVector3& direction,
+  std::vector<Hit > &hits,
+  int& num_rays,
+  float tnear,
+  float tfar) const
 {
   using namespace std;
   num_rays = 0;
   hits.clear();
-  embree::Vec3f o = toVec3f(origin);
-  embree::Vec3f d = toVec3f(direction);
   int last_id0 = -1;
   double self_hits = 0;
   // This epsilon is directly correleated to the number of missed hits, smaller
   // means more accurate and slower
   //const double eps = DOUBLE_EPS;
   const double eps = FLOAT_EPS;
-  double min_t = embree::zero;
+  double min_t = tnear;
   bool large_hits_warned = false;
+  embree::Ray ray(toVector3f(origin),toVector3f(direction));
+
   while(true)
   {
-#ifdef VERBOSE
-    cout<<
-      o[0]<<" "<<o[1]<<" "<<o[2]<<" + t*"<<
-      d[0]<<" "<<d[1]<<" "<<d[2]<<" ---> "<<
-      endl;
-#endif
-    embree::Hit hit;
-    embree::Ray ray(o,d,min_t);
+    ray.tnear = min_t;
+    ray.tfar = tfar;
+    ray.id0 = -1;
     num_rays++;
-    _intersector->intersect(ray, hit);
-    if(hit)
+    intersector->intersect(ray);
+    if(ray)
     {
       // Hit self again, progressively advance
-      if(hit.id0 == last_id0 || hit.t <= min_t)
+      if(ray.id0 == last_id0 || ray.tfar <= min_t)
       {
-        // sanity check
-        assert(hit.t<1);
         // push min_t a bit more
         //double t_push = pow(2.0,self_hits-4)*(hit.t<eps?eps:hit.t);
         double t_push = pow(2.0,self_hits)*eps;
@@ -200,87 +247,77 @@ igl::EmbreeIntersector < PointMatrixType, FaceMatrixType, RowVector3>
         //o = o+t_push*d;
         min_t += t_push;
         self_hits++;
-      }else
+      }
+      else
       {
+        Hit hit;
+        hit.id = ray.id0;
+        hit.u = ray.u;
+        hit.v = ray.v;
+        hit.t = ray.tfar;
         hits.push_back(hit);
 #ifdef VERBOSE
         cout<<"  t: "<<hit.t<<endl;
 #endif
         // Instead of moving origin, just change min_t. That way calculations
         // all use exactly same origin values
-        min_t = hit.t;
+        min_t = ray.tfar;
 
         // reset t_scale
         self_hits = 0;
       }
-      last_id0 = hit.id0;
-      //cout<<"  id0: "<<hit.id0<<endl;
-    }else
-    {
-      break;
+      last_id0 = ray.id0;
     }
+    else
+      break; // no more hits
+    
     if(hits.size()>1000 && !large_hits_warned)
     {
       cerr<<"Warning: Large number of hits..."<<endl;
       cerr<<"[ ";
-      for(vector<embree::Hit>::iterator hit = hits.begin();
-          hit != hits.end();
-          hit++)
+      for(vector<Hit>::iterator hit = hits.begin(); hit != hits.end();hit++)
       {
-        cerr<<(hit->id0+1)<<" ";
+        cerr<<(hit->id+1)<<" ";
       }
+      
       cerr.precision(std::numeric_limits< double >::digits10);
       cerr<<"[ ";
-      for(vector<embree::Hit>::iterator hit = hits.begin();
-          hit != hits.end();
-          hit++)
+      
+      for(vector<Hit>::iterator hit = hits.begin(); hit != hits.end(); hit++)
       {
         cerr<<(hit->t)<<endl;;
       }
 
       cerr<<"]"<<endl;
       large_hits_warned = true;
+
       return hits.empty();
     }
   }
+
   return hits.empty();
 }
 
 template <
-typename PointMatrixType,
-typename FaceMatrixType,
-typename RowVector3>
+typename Scalar,
+typename Index>
 bool 
-igl::EmbreeIntersector < PointMatrixType, FaceMatrixType, RowVector3>
-::intersectSegment(const RowVector3& a, const RowVector3& ab, embree::Hit &hit) const
+igl::EmbreeIntersector < Scalar, Index>
+::intersectSegment(const RowVector3& a, const RowVector3& ab, Hit &hit) const
 {
-  embree::Ray ray(toVec3f(a), toVec3f(ab), embree::zero, embree::one);
-  _intersector->intersect(ray, hit);
-  return hit ; 
-}
+  embree::Ray ray(toVector3f(a), toVector3f(ab), embree::zero, embree::one);
+  intersector->intersect(ray);
 
-//template <
-//typename PointMatrixType,
-//typename FaceMatrixType,
-//typename RowVector3>
-//RowVector3 
-//igl::EmbreeIntersector < PointMatrixType, FaceMatrixType, RowVector3>
-//::hit_position(const embree::Hit & hit) const
-//{
-//  // Barycenteric coordinates
-//  const double w0 = (1.0-hit.u-hit.v);
-//  const double w1 = hit.u;
-//  const double w2 = hit.v;
-//  RowVector3 hitP;
-//  hitP.resize(3);
-//  for(int d = 0;d<3;d++)
-//  {
-//    hitP(d) = 
-//    w0 * vertices[triangles[hit.id0](0)](d) + 
-//    w1 * vertices[triangles[hit.id0](1)](d) + 
-//    w2 * vertices[triangles[hit.id0](2)](d);
-//  }
-//  return hitP;
-//}
+  if(ray)
+  {
+    hit.id = ray.id0;
+    hit.u = ray.u;
+    hit.v = ray.v;
+    hit.t = ray.tfar;
+    return true;
+  }
+
+  return false;
+}
 
 #endif //EMBREE_INTERSECTOR_H

+ 1 - 2
include/igl/embree/Embree_convenience.h

@@ -16,8 +16,7 @@
 // This is a hack
 #  pragma GCC system_header
 #endif
-#include <common/intersector.h>
-#include <common/accel.h>
+#include <common/ray.h>
 #ifdef __GNUC__
 #  if __GNUC__ >= 4
 #    if __GNUC_MINOR__ >= 6