libigl v2.5.0
Loading...
Searching...
No Matches
igl::embree Namespace Reference

Classes

struct  EmbreeDevice
 keep track of embree device instance More...
 
class  EmbreeIntersector
 Simple class to wrap Embree's ray tracing functionality. More...
 

Functions

template<typename DerivedV , typename DerivedF , typename DerivedP , typename DerivedN , typename DerivedS >
void ambient_occlusion (const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedF > &F, const Eigen::MatrixBase< DerivedP > &P, const Eigen::MatrixBase< DerivedN > &N, const int num_samples, Eigen::PlainObjectBase< DerivedS > &S)
 Compute ambient occlusion per given point.
 
template<typename DerivedP , typename DerivedN , typename DerivedS >
void ambient_occlusion (const EmbreeIntersector &ei, const Eigen::MatrixBase< DerivedP > &P, const Eigen::MatrixBase< DerivedN > &N, const int num_samples, Eigen::PlainObjectBase< DerivedS > &S)
 
bool bone_heat (const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const Eigen::MatrixXd &C, const Eigen::VectorXi &P, const Eigen::MatrixXi &BE, const Eigen::MatrixXi &CE, Eigen::MatrixXd &W)
 Compute skinning weights W given a surface mesh (V,F) and an internal skeleton (C,BE) according to "Automatic Rigging" [Baran and Popovic 2007].
 
template<typename DerivedV , typename DerivedF , typename DerivedSD , typename Derivedflag >
void bone_visible (const Eigen::PlainObjectBase< DerivedV > &V, const Eigen::PlainObjectBase< DerivedF > &F, const Eigen::PlainObjectBase< DerivedSD > &s, const Eigen::PlainObjectBase< DerivedSD > &d, Eigen::PlainObjectBase< Derivedflag > &flag)
 Test whether vertices of mesh are "visible" to a given bone, where "visible" is defined as in [Baran & Popovic 07].
 
template<typename DerivedV , typename DerivedF , typename DerivedSD , typename Derivedflag >
void bone_visible (const Eigen::PlainObjectBase< DerivedV > &V, const Eigen::PlainObjectBase< DerivedF > &F, const EmbreeIntersector &ei, const Eigen::PlainObjectBase< DerivedSD > &s, const Eigen::PlainObjectBase< DerivedSD > &d, Eigen::PlainObjectBase< Derivedflag > &flag)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename ScalarMatrix , typename IndexMatrix >
ScalarMatrix line_mesh_intersection (const ScalarMatrix &V_source, const ScalarMatrix &N_source, const ScalarMatrix &V_target, const IndexMatrix &F_target)
 Project the point cloud V_source onto the triangle mesh V_target,F_target.
 
template<typename DerivedV , typename DerivedF , typename DerivedI , typename DerivedC >
void reorient_facets_raycast (const Eigen::PlainObjectBase< DerivedV > &V, const Eigen::PlainObjectBase< DerivedF > &F, int rays_total, int rays_minimum, bool facet_wise, bool use_parity, bool is_verbose, Eigen::PlainObjectBase< DerivedI > &I, Eigen::PlainObjectBase< DerivedC > &C)
 Orient each component (identified by C) of a mesh (V,F) using ambient occlusion such that the front side is less occluded than back side, as described in "A Simple Method for Correcting Facet Orientations in Polygon Meshes Based on Ray Casting" [Takayama et al.
 
template<typename DerivedV , typename DerivedF , typename DerivedFF , typename DerivedI >
void reorient_facets_raycast (const Eigen::PlainObjectBase< DerivedV > &V, const Eigen::PlainObjectBase< DerivedF > &F, Eigen::PlainObjectBase< DerivedFF > &FF, Eigen::PlainObjectBase< DerivedI > &I)
 
template<typename DerivedP , typename DerivedN , typename DerivedS >
void shape_diameter_function (const EmbreeIntersector &ei, const Eigen::PlainObjectBase< DerivedP > &P, const Eigen::PlainObjectBase< DerivedN > &N, const int num_samples, Eigen::PlainObjectBase< DerivedS > &S)
 Compute shape diamter function per given point.
 
template<typename DerivedV , typename DerivedF , typename DerivedP , typename DerivedN , typename DerivedS >
void shape_diameter_function (const Eigen::PlainObjectBase< DerivedV > &V, const Eigen::PlainObjectBase< DerivedF > &F, const Eigen::PlainObjectBase< DerivedP > &P, const Eigen::PlainObjectBase< DerivedN > &N, const int num_samples, Eigen::PlainObjectBase< DerivedS > &S)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename Derivedobj >
int unproject_in_mesh (const Eigen::Vector2f &pos, const Eigen::Matrix4f &model, const Eigen::Matrix4f &proj, const Eigen::Vector4f &viewport, const EmbreeIntersector &ei, Eigen::PlainObjectBase< Derivedobj > &obj, std::vector< igl::Hit > &hits)
 Unproject a screen location (using current opengl viewport, projection, and model view) to a 3D position inside a given mesh.
 
template<typename Derivedobj >
int unproject_in_mesh (const Eigen::Vector2f &pos, const Eigen::Matrix4f &model, const Eigen::Matrix4f &proj, const Eigen::Vector4f &viewport, const EmbreeIntersector &ei, Eigen::PlainObjectBase< Derivedobj > &obj)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
bool unproject_onto_mesh (const Eigen::Vector2f &pos, const Eigen::Matrix4f &model, const Eigen::Matrix4f &proj, const Eigen::Vector4f &viewport, const EmbreeIntersector &ei, int &fid, Eigen::Vector3f &bc)
 Unproject a screen location (using the given model, proj and viewport) to find the first hit on a mesh.
 
bool unproject_onto_mesh (const Eigen::Vector2f &pos, const Eigen::MatrixXi &F, const Eigen::Matrix4f &model, const Eigen::Matrix4f &proj, const Eigen::Vector4f &viewport, const EmbreeIntersector &ei, int &fid, int &vid)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 

Function Documentation

◆ ambient_occlusion() [1/2]

template<typename DerivedV , typename DerivedF , typename DerivedP , typename DerivedN , typename DerivedS >
void igl::embree::ambient_occlusion ( const Eigen::MatrixBase< DerivedV > &  V,
const Eigen::MatrixBase< DerivedF > &  F,
const Eigen::MatrixBase< DerivedP > &  P,
const Eigen::MatrixBase< DerivedN > &  N,
const int  num_samples,
Eigen::PlainObjectBase< DerivedS > &  S 
)

Compute ambient occlusion per given point.

Parameters
[in]V#V by 3 list of mesh vertex positiosn
[in]F#F by 3 list of mesh triangle indices into rows of V
[in]P#P by 3 list of origin points
[in]N#P by 3 list of origin normals
[out]S#P list of ambient occlusion values between 1 (fully occluded) and 0 (not occluded)

◆ ambient_occlusion() [2/2]

template<typename DerivedP , typename DerivedN , typename DerivedS >
void igl::embree::ambient_occlusion ( const EmbreeIntersector ei,
const Eigen::MatrixBase< DerivedP > &  P,
const Eigen::MatrixBase< DerivedN > &  N,
const int  num_samples,
Eigen::PlainObjectBase< DerivedS > &  S 
)

◆ bone_heat()

bool igl::embree::bone_heat ( const Eigen::MatrixXd &  V,
const Eigen::MatrixXi &  F,
const Eigen::MatrixXd &  C,
const Eigen::VectorXi &  P,
const Eigen::MatrixXi &  BE,
const Eigen::MatrixXi &  CE,
Eigen::MatrixXd &  W 
)

Compute skinning weights W given a surface mesh (V,F) and an internal skeleton (C,BE) according to "Automatic Rigging" [Baran and Popovic 2007].

Parameters
[in]V#V by 3 list of mesh vertex positions
[in]F#F by 3 list of mesh corner indices into V
[in]C#C by 3 list of joint locations
[in]P#P list of point handle indices into C
[in]BE#BE by 2 list of bone edge indices into C
[in]CE#CE by 2 list of cage edge indices into P
[out]W#V by #P+#BE matrix of weights.
Returns
true only on success.

◆ bone_visible() [1/2]

template<typename DerivedV , typename DerivedF , typename DerivedSD , typename Derivedflag >
void igl::embree::bone_visible ( const Eigen::PlainObjectBase< DerivedV > &  V,
const Eigen::PlainObjectBase< DerivedF > &  F,
const Eigen::PlainObjectBase< DerivedSD > &  s,
const Eigen::PlainObjectBase< DerivedSD > &  d,
Eigen::PlainObjectBase< Derivedflag > &  flag 
)

Test whether vertices of mesh are "visible" to a given bone, where "visible" is defined as in [Baran & Popovic 07].

Instead of checking whether each point can see any of the bone, we just check if each point can see its own projection onto the bone segment. In other words, we project each vertex v onto the bone, projv. Then we check if there are any intersections between the line segment (projv-->v) and the mesh.

Parameters
[in]V#V by 3 list of vertex positions
[in]F#F by 3 list of triangle indices
[in]srow vector of position of start end point of bone
[in]drow vector of position of dest end point of bone
[out]flag#V by 1 list of bools (true) visible, (false) obstructed
Note
This checks for hits along the segment which are facing in any direction from the ray.

◆ bone_visible() [2/2]

template<typename DerivedV , typename DerivedF , typename DerivedSD , typename Derivedflag >
void igl::embree::bone_visible ( const Eigen::PlainObjectBase< DerivedV > &  V,
const Eigen::PlainObjectBase< DerivedF > &  F,
const EmbreeIntersector ei,
const Eigen::PlainObjectBase< DerivedSD > &  s,
const Eigen::PlainObjectBase< DerivedSD > &  d,
Eigen::PlainObjectBase< Derivedflag > &  flag 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Parameters
[in]eiEmbreeIntersector for mesh (V,F) should be double sided

◆ line_mesh_intersection()

template<typename ScalarMatrix , typename IndexMatrix >
ScalarMatrix igl::embree::line_mesh_intersection ( const ScalarMatrix &  V_source,
const ScalarMatrix &  N_source,
const ScalarMatrix &  V_target,
const IndexMatrix &  F_target 
)

Project the point cloud V_source onto the triangle mesh V_target,F_target.

A ray is casted for every vertex in the direction specified by N_source and its opposite.

Parameters
[in]V_source#Vx3 Vertices of the source mesh
[in]N_source#Vx3 Normals of the point cloud
[in]V_target#V2x3 Vertices of the target mesh
[in]F_target#F2x3 Faces of the target mesh
Returns
#Vx3 matrix of baricentric coordinate. Each row corresponds to a vertex of the projected mesh and it has the following format: id b1 b2. id is the id of a face of the source mesh. b1 and b2 are the barycentric coordinates wrt the first two edges of the triangle To convert to standard global coordinates, see barycentric_to_global.h

◆ reorient_facets_raycast() [1/2]

template<typename DerivedV , typename DerivedF , typename DerivedI , typename DerivedC >
void igl::embree::reorient_facets_raycast ( const Eigen::PlainObjectBase< DerivedV > &  V,
const Eigen::PlainObjectBase< DerivedF > &  F,
int  rays_total,
int  rays_minimum,
bool  facet_wise,
bool  use_parity,
bool  is_verbose,
Eigen::PlainObjectBase< DerivedI > &  I,
Eigen::PlainObjectBase< DerivedC > &  C 
)

Orient each component (identified by C) of a mesh (V,F) using ambient occlusion such that the front side is less occluded than back side, as described in "A Simple Method for Correcting Facet Orientations in Polygon Meshes Based on Ray Casting" [Takayama et al.

2014].

Parameters
[in]V#V by 3 list of vertex positions
[in]F#F by 3 list of triangle indices
[in]rays_totalTotal number of rays that will be shot
[in]rays_minimumMinimum number of rays that each patch should receive
[in]facet_wiseDecision made for each face independently, no use of patches (i.e., each face is treated as a patch)
[in]use_parityUse parity mode
[in]is_verboseVerbose output to cout
[out]I#F list of whether face has been flipped
[out]C#F list of patch ID (output of bfs_orient > manifold patches)

◆ reorient_facets_raycast() [2/2]

template<typename DerivedV , typename DerivedF , typename DerivedFF , typename DerivedI >
void igl::embree::reorient_facets_raycast ( const Eigen::PlainObjectBase< DerivedV > &  V,
const Eigen::PlainObjectBase< DerivedF > &  F,
Eigen::PlainObjectBase< DerivedFF > &  FF,
Eigen::PlainObjectBase< DerivedI > &  I 
)
Parameters
[out]FF#F by 3 list of reoriented faces

Defaults:

rays_total = F.rows()*100; rays_minimum = 10; facet_wise = false; use_parity = false; is_verbose = false;

◆ shape_diameter_function() [1/2]

template<typename DerivedP , typename DerivedN , typename DerivedS >
void igl::embree::shape_diameter_function ( const EmbreeIntersector ei,
const Eigen::PlainObjectBase< DerivedP > &  P,
const Eigen::PlainObjectBase< DerivedN > &  N,
const int  num_samples,
Eigen::PlainObjectBase< DerivedS > &  S 
)

Compute shape diamter function per given point.

Parameters
[in]eiEmbreeIntersector containing (V,F)
[in]P#P by 3 list of origin points
[in]N#P by 3 list of origin normals
[out]S#P list of shape diamater function values between bounding box diagonal (perfect sphere) and 0 (perfect needle hook)

◆ shape_diameter_function() [2/2]

template<typename DerivedV , typename DerivedF , typename DerivedP , typename DerivedN , typename DerivedS >
void igl::embree::shape_diameter_function ( const Eigen::PlainObjectBase< DerivedV > &  V,
const Eigen::PlainObjectBase< DerivedF > &  F,
const Eigen::PlainObjectBase< DerivedP > &  P,
const Eigen::PlainObjectBase< DerivedN > &  N,
const int  num_samples,
Eigen::PlainObjectBase< DerivedS > &  S 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Wrapper which builds new EmbreeIntersector for (V,F). That's expensive so avoid this if repeatedly calling.

Parameters
[in]V#V by 3 list of mesh vertex positiosn
[in]F#F by 3 list of mesh triangle indices into rows of V

◆ unproject_in_mesh() [1/2]

template<typename Derivedobj >
int igl::embree::unproject_in_mesh ( const Eigen::Vector2f &  pos,
const Eigen::Matrix4f &  model,
const Eigen::Matrix4f &  proj,
const Eigen::Vector4f &  viewport,
const EmbreeIntersector ei,
Eigen::PlainObjectBase< Derivedobj > &  obj,
std::vector< igl::Hit > &  hits 
)

Unproject a screen location (using current opengl viewport, projection, and model view) to a 3D position inside a given mesh.

If the ray through the given screen location (x,y) hits the mesh more than twice then the 3D midpoint between the first two hits is return. If it hits once, then that point is return. If it does not hit the mesh then obj is not set.

Parameters
[in]posscreen space coordinates
[in]modelmodel matrix
[in]projprojection matrix
[in]viewportvieweport vector
[in]eiEmbreeIntersector containing (V,F)
[out]obj3d unprojected mouse point in mesh
[out]hitsvector of embree hits
Returns
number of hits
Note
Previous prototype did not require model, proj, and viewport. This has been removed. Instead replace with:
 Eigen::Matrix4f model,proj;
 Eigen::Vector4f viewport;
 igl::opengl2::model_proj_viewport(model,proj,viewport);
 igl::embree::unproject_in_mesh(Vector2f(x,y),model,proj,viewport,ei,obj,hits);

◆ unproject_in_mesh() [2/2]

template<typename Derivedobj >
int igl::embree::unproject_in_mesh ( const Eigen::Vector2f &  pos,
const Eigen::Matrix4f &  model,
const Eigen::Matrix4f &  proj,
const Eigen::Vector4f &  viewport,
const EmbreeIntersector ei,
Eigen::PlainObjectBase< Derivedobj > &  obj 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ unproject_onto_mesh() [1/2]

bool igl::embree::unproject_onto_mesh ( const Eigen::Vector2f &  pos,
const Eigen::Matrix4f &  model,
const Eigen::Matrix4f &  proj,
const Eigen::Vector4f &  viewport,
const EmbreeIntersector ei,
int &  fid,
Eigen::Vector3f &  bc 
)

Unproject a screen location (using the given model, proj and viewport) to find the first hit on a mesh.

Parameters
[in]posscreen space coordinates
[in]modelmodel matrix
[in]projprojection matrix
[in]viewportvieweport vector
[in]eiEmbreeIntersector containing (V,F)
[out]fidid of the first face hit
[out]bcbarycentric coordinates of hit
Returns
true if there is a hit

◆ unproject_onto_mesh() [2/2]

bool igl::embree::unproject_onto_mesh ( const Eigen::Vector2f &  pos,
const Eigen::MatrixXi &  F,
const Eigen::Matrix4f &  model,
const Eigen::Matrix4f &  proj,
const Eigen::Vector4f &  viewport,
const EmbreeIntersector ei,
int &  fid,
int &  vid 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Parameters
[in]vidvertex id of the closest vertex hit