Isis 3 Programmer Reference
EmbreeTargetShape.h
Go to the documentation of this file.
1 #ifndef EmbreeTargetShape_h
2 #define EmbreeTargetShape_h
3 
26 #include <QString>
27 
28 // Embree includes
29 #include <embree2/rtcore.h>
30 #include <embree2/rtcore_ray.h>
31 
32 // PointCloudLibrary includes
33 #include <pcl/io/auto_io.h>
34 #include <pcl/point_types.h>
35 #include <pcl/PolygonMesh.h>
36 
37 // WARNING
38 // The following undefes macros set in pcl. They conflict with ISIS's macros.
39 // If pcl tries to use these macros you will get compilation errors.
40 // - JAM
41 #undef DEG2RAD
42 #undef RAD2DEG
43 
44 #include "FileName.h"
45 #include "LinearAlgebra.h"
46 
47 namespace Isis {
48 
49 
50  class Pvl;
51 
62  RTCMultiHitRay(const std::vector<double> &origin,
63  const std::vector<double> &direction);
65  LinearAlgebra::Vector direction);
66 
67 // These are the inheritted members from RTCRay
68 // float org; //!< Ray origin
69 // float dir; //!< Ray direction
70 // float tnear; //!< Start of ray segment
71 // float tfar; //!< End of ray segment
72 // float time; //!< Time of this ray for motion blur.
73 // unsigned mask; //!< used to mask out objects during traversal
74 // float Ng; //!< Geometric normal.
75 // float u; //!< Barycentric u coordinate of hit
76 // float v; //!< Barycentric v coordinate of hit
77 // unsigned geomID; //!< geometry ID
78 // unsigned primID; //!< primitive ID
79 // unsigned instID; //!< instance ID
80 
81  // ray extensions
82  // we remember up to 16 hits to ignore duplicate hits
83  unsigned hitGeomIDs[16];
84  unsigned hitPrimIDs[16];
85  float hitUs[16];
86  float hitVs[16];
87  int lastHit;
88  };
89 
90 
100  RTCOcclusionRay();
101  RTCOcclusionRay(const std::vector<double> &origin,
102  const std::vector<double> &direction);
104  LinearAlgebra::Vector direction);
105 
106 // These are the inheritted members from RTCRay
107 // float org; //!< Ray origin
108 // float dir; //!< Ray direction
109 // float tnear; //!< Start of ray segment
110 // float tfar; //!< End of ray segment
111 // float time; //!< Time of this ray for motion blur.
112 // unsigned mask; //!< used to mask out objects during traversal
113 // float Ng; //!< Geometric normal.
114 // float u; //!< Barycentric u coordinate of hit
115 // float v; //!< Barycentric v coordinate of hit
116 // unsigned geomID; //!< geometry ID
117 // unsigned primID; //!< primitive ID
118 // unsigned instID; //!< instance ID
119 
120  // ray extensions
121  int lastHit;
122  unsigned ignorePrimID;
123  };
124 
125 
137 
140  int primID;
141  };
142 
143 
156  public:
157 
159  EmbreeTargetShape(pcl::PolygonMesh::Ptr mesh, const QString &name = "");
160  EmbreeTargetShape(const QString &dem, const Pvl *conf = 0);
161  virtual ~EmbreeTargetShape();
162 
163  QString name() const;
164  bool isValid() const;
165 
166  int numberOfPolygons() const;
167  int numberOfVertices() const;
168  RTCBounds sceneBounds() const;
169  double maximumSceneDistance() const;
170 
171  void intersectRay(RTCMultiHitRay &ray);
172  bool isOccluded(RTCOcclusionRay &ray);
173 
175 
176  static void multiHitFilter(void* userDataPtr, RTCMultiHitRay& ray);
177  static void occlusionFilter(void* userDataPtr, RTCOcclusionRay& ray);
178 
179  protected:
180  pcl::PolygonMesh::Ptr readDSK(FileName file);
181  pcl::PolygonMesh::Ptr readPC(FileName file);
182  void initMesh(pcl::PolygonMesh::Ptr mesh);
183  void addVertices(int geomID);
184  void addIndices(int geomID);
185 
186  private:
194  struct Vertex {
195  float x;
196  float y;
197  float z;
198  float a;
199  };
200 
209  struct Triangle {
210  int v0;
211  int v1;
212  int v2;
213  };
214 
215  QString m_name;
216  pcl::PolygonMesh::Ptr m_mesh;
219  pcl::PointCloud<pcl::PointXYZ> m_cloud;
224  RTCDevice m_device;
226  RTCScene m_scene;
232  };
233 
234 } // namespace Isis
235 
236 #endif
237 
Container that holds the body fixed intersection point and unit surface normal for a hit...
unsigned hitGeomIDs[16]
IDs of the geometries (bodies) hit.
int numberOfVertices() const
Return the number of vertices in the target shape.
File name manipulation and expansion.
Definition: FileName.h:116
RTCMultiHitRay()
Default constructor for RTCMultiHitRay.
RTCOcclusionRay()
Default constructor for RTCOcclussionRay.
int numberOfPolygons() const
Return the number of polygons in the target shape.
int lastHit
Index of the last hit in the hit containers.
double maximumSceneDistance() const
Return the maximum distance within the scene.
Struct for capturing multiple intersections when using embree::rtcintersectscene. ...
pcl::PolygonMesh::Ptr readPC(FileName file)
Read a PointCloudLibrary file into a PointCloudLibrary polygon mesh.
unsigned hitPrimIDs[16]
IDs of the primitives (trinagles) hit.
virtual ~EmbreeTargetShape()
Desctructor.
float hitUs[16]
Barycentric u coordinate of the hits.
void addIndices(int geomID)
Adds the polygon vertex indices from the internalized polygon mesh to the Embree scene.
Container for a tin, or triangular polygon.
LinearAlgebra::Vector surfaceNormal
The unit surface normal vector at the intersection.
int v1
The index of the second vertex in the tin.
boost::numeric::ublas::vector< double > Vector
Definition for an Isis::LinearAlgebra::Vector of doubles.
void initMesh(pcl::PolygonMesh::Ptr mesh)
Internalize a PointCloudLibrary polygon mesh in the target shape.
float y
Vertex y position.
RTCScene m_scene
!< The Embree device for rendering the scene.
RayHitInformation()
Default constructor for RayHitInformation.
unsigned ignorePrimID
IDs of the primitives (trinagles) which should be ignored.
int v0
The index of the first vertex in the tin.
Embree Target Shape for planetary bodies.
float x
Vertex x position.
void addVertices(int geomID)
Adds the vertices from the internalized vertex point cloud to the Embree scene.
void intersectRay(RTCMultiHitRay &ray)
Intersect a ray with the target shape.
pcl::PolygonMesh::Ptr m_mesh
!< The name of the target.
RTCBounds sceneBounds() const
Returns the bounds of the Embree scene.
Container for cube-like labels.
Definition: Pvl.h:135
float hitVs[16]
Barycentric v coordinate of the hits.
LinearAlgebra::Vector intersection
The (x, y, z) intersection location.
RayHitInformation getHitInformation(RTCMultiHitRay &ray, int hitIndex)
Extract the intersection point and unit surface normal from an RTCMultiHitRay that has been intersect...
QString name() const
Return the name of the target shape.
Struct for capturing occluded plates when using embree::rtcintersectscene.
bool isValid() const
Return if a valid mesh is internalized and ready for use.
Namespace for ISIS/Bullet specific routines.
Definition: Apollo.h:31
int primID
The primitive ID of the hit.
EmbreeTargetShape()
Default empty constructor.
float a
Extra float for aligning memory.
int lastHit
Index of the last hit in the hit containers.
bool isOccluded(RTCOcclusionRay &ray)
Check if a ray intersects the target body.
Container for a vertex.
RTCDevice m_device
!< The point cloud representation of the target.
pcl::PointCloud< pcl::PointXYZ > m_cloud
!< A boost shared pointer to the polygon mesh representation of the target.
static void occlusionFilter(void *userDataPtr, RTCOcclusionRay &ray)
Filter function for collecting multiple primitiveIDs This function is called by the Embree library du...
int v2
The index of the third vertex in the tin.
float z
Vertex z position.
static void multiHitFilter(void *userDataPtr, RTCMultiHitRay &ray)
Filter function for collecting multiple hits during ray intersection.
pcl::PolygonMesh::Ptr readDSK(FileName file)
Read a NAIF type 2 DSK file into a PointCloudLibrary polygon mesh.