wire_mesh.cpp 7.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253
  1. #include "wire_mesh.h"
  2. #include "../../list_to_matrix.h"
  3. #include "../../PI.h"
  4. #include "convex_hull.h"
  5. #include "coplanar.h"
  6. #include "mesh_boolean.h"
  7. #include <Eigen/Geometry>
  8. #include <vector>
  9. template <
  10. typename DerivedWV,
  11. typename DerivedWE,
  12. typename Derivedth,
  13. typename DerivedV,
  14. typename DerivedF,
  15. typename DerivedJ>
  16. IGL_INLINE void igl::copyleft::cgal::wire_mesh(
  17. const Eigen::MatrixBase<DerivedWV> & WV,
  18. const Eigen::MatrixBase<DerivedWE> & WE,
  19. const Eigen::MatrixBase<Derivedth> & th,
  20. const int poly_size,
  21. const bool solid,
  22. Eigen::PlainObjectBase<DerivedV> & V,
  23. Eigen::PlainObjectBase<DerivedF> & F,
  24. Eigen::PlainObjectBase<DerivedJ> & J)
  25. {
  26. assert((th.size()==1 || th.size()==WE.rows()) &&
  27. "th should be scalar or size of WE");
  28. typedef typename DerivedWV::Scalar Scalar;
  29. // Canonical polygon to place at each endpoint
  30. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,3> MatrixX3S;
  31. MatrixX3S PV(poly_size,3);
  32. for(int p =0;p<PV.rows();p++)
  33. {
  34. const Scalar phi = (Scalar(p)/Scalar(PV.rows()))*2.*igl::PI;
  35. PV(p,0) = 0.5*cos(phi);
  36. PV(p,1) = 0.5*sin(phi);
  37. PV(p,2) = 0;
  38. }
  39. V.resize(WV.rows() + PV.rows() * 2 * WE.rows(),3);
  40. V.topLeftCorner(WV.rows(),3) = WV;
  41. // Signed adjacency list
  42. std::vector<std::vector<std::pair<int,int> > > A(WV.rows());
  43. // Inputs:
  44. // e index of edge
  45. // c index of endpoint [0,1]
  46. // p index of polygon vertex
  47. // Returns index of corresponding vertex in V
  48. const auto index =
  49. [&PV,&WV](const int e, const int c, const int p)->int
  50. {
  51. return WV.rows() + e*2*PV.rows() + PV.rows()*c + p;
  52. };
  53. //const auto unindex =
  54. // [&PV,&WV](int v, int & e, int & c, int & p)
  55. //{
  56. // assert(v>=WV.rows());
  57. // v = v-WV.rows();
  58. // e = v/(2*PV.rows());
  59. // v = v-e*(2*PV.rows());
  60. // c = v/(PV.rows());
  61. // v = v-c*(PV.rows());
  62. // p = v;
  63. //};
  64. // Count each vertex's indicident edges.
  65. std::vector<int> nedges(WV.rows(), 0);
  66. for(int e = 0;e<WE.rows();e++)
  67. {
  68. ++nedges[WE(e, 0)];
  69. ++nedges[WE(e, 1)];
  70. }
  71. // loop over all edges
  72. for(int e = 0;e<WE.rows();e++)
  73. {
  74. // Fill in adjacency list as we go
  75. A[WE(e,0)].emplace_back(e,0);
  76. A[WE(e,1)].emplace_back(e,1);
  77. typedef Eigen::Matrix<Scalar,1,3> RowVector3S;
  78. const RowVector3S ev = WV.row(WE(e,1))-WV.row(WE(e,0));
  79. const Scalar len = ev.norm();
  80. // Unit edge vector
  81. const RowVector3S uv = ev.normalized();
  82. Eigen::Quaternion<Scalar> q;
  83. q = q.FromTwoVectors(RowVector3S(0,0,1),uv);
  84. // loop over polygon vertices
  85. for(int p = 0;p<PV.rows();p++)
  86. {
  87. RowVector3S qp = q*(PV.row(p)*th(e%th.size()));
  88. // loop over endpoints
  89. for(int c = 0;c<2;c++)
  90. {
  91. // Direction moving along edge vector
  92. const Scalar dir = c==0?1:-1;
  93. // Amount (distance) to move along edge vector
  94. // Start with factor of thickness;
  95. // Max out amount at 1/3 of edge length so that there's always some
  96. // amount of edge
  97. // Zero out if vertex is incident on only one edge
  98. Scalar dist =
  99. std::min(1.*th(e%th.size()),len/3.0)*(nedges[WE(e,c)] > 1);
  100. // Move to endpoint, offset by amount
  101. V.row(index(e,c,p)) =
  102. qp+WV.row(WE(e,c)) + dist*dir*uv;
  103. }
  104. }
  105. }
  106. std::vector<std::vector<typename DerivedF::Index> > vF;
  107. std::vector<int> vJ;
  108. const auto append_hull =
  109. [&V,&vF,&vJ](const Eigen::VectorXi & I, const int j)
  110. {
  111. MatrixX3S Vv = V(I,Eigen::placeholders::all);
  112. if(coplanar(Vv))
  113. {
  114. return;
  115. }
  116. Eigen::MatrixXi Fv;
  117. convex_hull(Vv,Fv);
  118. for(int f = 0;f<Fv.rows();f++)
  119. {
  120. const Eigen::Array<int,1,3> face(I(Fv(f,0)), I(Fv(f,1)), I(Fv(f,2)));
  121. //const bool on_vertex = (face<WV.rows()).any();
  122. //if(!on_vertex)
  123. //{
  124. // // This correctly prunes fcaes on the "caps" of convex hulls around
  125. // // edges, but for convex hulls around vertices this will only work if
  126. // // the incoming edges are not overlapping.
  127. // //
  128. // // Q: For convex hulls around vertices, is the correct thing to do:
  129. // // check if all corners of face lie *on or _outside_* of plane of "cap"?
  130. // //
  131. // // H: Maybe, but if there's an intersection then the boundary of the
  132. // // incoming convex hulls around edges is still not going to match up
  133. // // with the boundary on the convex hull around the vertices.
  134. // //
  135. // // Might have to bite the bullet and always call self-union.
  136. // bool all_same = true;
  137. // int e0,c0,p0;
  138. // unindex(face(0),e0,c0,p0);
  139. // for(int i = 1;i<3;i++)
  140. // {
  141. // int ei,ci,pi;
  142. // unindex(face(i),ei,ci,pi);
  143. // all_same = all_same && (e0==ei && c0==ci);
  144. // }
  145. // if(all_same)
  146. // {
  147. // // don't add this face
  148. // continue;
  149. // }
  150. //}
  151. vF.push_back( { face(0),face(1),face(2)});
  152. vJ.push_back(j);
  153. }
  154. };
  155. // loop over each vertex
  156. for(int v = 0;v<WV.rows();v++)
  157. {
  158. // Gather together this vertex and the polygon vertices of all incident
  159. // edges
  160. Eigen::VectorXi I(1+A[v].size()*PV.rows());
  161. // This vertex
  162. I(0) = v;
  163. for(int n = 0;n<A[v].size();n++)
  164. {
  165. for(int p = 0;p<PV.rows();p++)
  166. {
  167. const int e = A[v][n].first;
  168. const int c = A[v][n].second;
  169. I(1+n*PV.rows()+p) = index(e,c,p);
  170. }
  171. }
  172. append_hull(I,v);
  173. }
  174. // loop over each edge
  175. for(int e = 0;e<WE.rows();e++)
  176. {
  177. // Gether together polygon vertices of both endpoints
  178. Eigen::VectorXi I(PV.rows()*2);
  179. for(int c = 0;c<2;c++)
  180. {
  181. for(int p = 0;p<PV.rows();p++)
  182. {
  183. I(c*PV.rows()+p) = index(e,c,p);
  184. }
  185. }
  186. append_hull(I,WV.rows()+e);
  187. }
  188. list_to_matrix(vF,F);
  189. if(solid)
  190. {
  191. // Self-union to clean up
  192. igl::copyleft::cgal::mesh_boolean(
  193. Eigen::MatrixXd(V),Eigen::MatrixXi(F),Eigen::MatrixXd(),Eigen::MatrixXi(),
  194. "union",
  195. V,F,J);
  196. for(int j=0;j<J.size();j++) J(j) = vJ[J(j)];
  197. }else
  198. {
  199. list_to_matrix(vJ,J);
  200. }
  201. }
  202. template <
  203. typename DerivedWV,
  204. typename DerivedWE,
  205. typename DerivedV,
  206. typename DerivedF,
  207. typename DerivedJ>
  208. IGL_INLINE void igl::copyleft::cgal::wire_mesh(
  209. const Eigen::MatrixBase<DerivedWV> & WV,
  210. const Eigen::MatrixBase<DerivedWE> & WE,
  211. const double th,
  212. const int poly_size,
  213. const bool solid,
  214. Eigen::PlainObjectBase<DerivedV> & V,
  215. Eigen::PlainObjectBase<DerivedF> & F,
  216. Eigen::PlainObjectBase<DerivedJ> & J)
  217. {
  218. return wire_mesh(
  219. WV,WE,(Eigen::VectorXd(1,1)<<th).finished(),poly_size,solid,V,F,J);
  220. }
  221. template <
  222. typename DerivedWV,
  223. typename DerivedWE,
  224. typename DerivedV,
  225. typename DerivedF,
  226. typename DerivedJ>
  227. IGL_INLINE void igl::copyleft::cgal::wire_mesh(
  228. const Eigen::MatrixBase<DerivedWV> & WV,
  229. const Eigen::MatrixBase<DerivedWE> & WE,
  230. const double th,
  231. const int poly_size,
  232. Eigen::PlainObjectBase<DerivedV> & V,
  233. Eigen::PlainObjectBase<DerivedF> & F,
  234. Eigen::PlainObjectBase<DerivedJ> & J)
  235. {
  236. return wire_mesh(WV,WE,th,poly_size,true,V,F,J);
  237. }
  238. #ifdef IGL_STATIC_LIBRARY
  239. // Explicit template instantiation
  240. template void igl::copyleft::cgal::wire_mesh<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, double, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  241. #endif