wire_mesh.cpp 7.6 KB

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