boundary_conditions.cpp 7.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <[email protected]>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include "boundary_conditions.h"
  9. #include "verbose.h"
  10. #include "EPS.h"
  11. #include "project_to_line.h"
  12. #include <vector>
  13. #include <map>
  14. #include <iostream>
  15. template <
  16. typename DerivedV,
  17. typename DerivedEle,
  18. typename DerivedC,
  19. typename DerivedP,
  20. typename DerivedBE,
  21. typename DerivedCE,
  22. typename DerivedCF,
  23. typename Derivedb,
  24. typename Derivedbc>
  25. IGL_INLINE bool igl::boundary_conditions(
  26. const Eigen::MatrixBase<DerivedV> & V,
  27. const Eigen::MatrixBase<DerivedEle> & Ele,
  28. const Eigen::MatrixBase<DerivedC> & C,
  29. const Eigen::MatrixBase<DerivedP> & P,
  30. const Eigen::MatrixBase<DerivedBE> & BE,
  31. const Eigen::MatrixBase<DerivedCE> & CE,
  32. const Eigen::MatrixBase<DerivedCF> & CF,
  33. Eigen::PlainObjectBase<Derivedb> & b,
  34. Eigen::PlainObjectBase<Derivedbc> & bc)
  35. {
  36. if(P.size()+BE.rows() == 0)
  37. {
  38. verbose("^%s: Error: no handles found\n",__FUNCTION__);
  39. return false;
  40. }
  41. std::vector<int> bci;
  42. std::vector<int> bcj;
  43. std::vector<double> bcv;
  44. // loop over points
  45. for(int p = 0;p<P.size();p++)
  46. {
  47. Eigen::VectorXd pos = C.row(P(p));
  48. // loop over domain vertices
  49. for(int i = 0;i<V.rows();i++)
  50. {
  51. // Find samples just on pos
  52. //Vec3 vi(V(i,0),V(i,1),V(i,2));
  53. // EIGEN GOTCHA:
  54. // double sqrd = (V.row(i)-pos).array().pow(2).sum();
  55. // Must first store in temporary
  56. Eigen::VectorXd vi = V.row(i);
  57. double sqrd = (vi-pos).squaredNorm();
  58. if(sqrd <= FLOAT_EPS)
  59. {
  60. //cout<<"sum((["<<
  61. // V(i,0)<<" "<<
  62. // V(i,1)<<" "<<
  63. // V(i,2)<<"] - ["<<
  64. // pos(0)<<" "<<
  65. // pos(1)<<" "<<
  66. // pos(2)<<"]).^2) = "<<sqrd<<endl;
  67. bci.push_back(i);
  68. bcj.push_back(p);
  69. bcv.push_back(1.0);
  70. }
  71. }
  72. }
  73. // loop over bone edges
  74. for(int e = 0;e<BE.rows();e++)
  75. {
  76. // loop over domain vertices
  77. for(int i = 0;i<V.rows();i++)
  78. {
  79. // Find samples from tip up to tail
  80. Eigen::VectorXd tip = C.row(BE(e,0));
  81. Eigen::VectorXd tail = C.row(BE(e,1));
  82. // Compute parameter along bone and squared distance
  83. double t,sqrd;
  84. project_to_line(
  85. V(i,0),V(i,1),V(i,2),
  86. tip(0),tip(1),tip(2),
  87. tail(0),tail(1),tail(2),
  88. t,sqrd);
  89. if(t>=-FLOAT_EPS && t<=(1.0f+FLOAT_EPS) && sqrd<=FLOAT_EPS)
  90. {
  91. bci.push_back(i);
  92. bcj.push_back(P.size()+e);
  93. bcv.push_back(1.0);
  94. }
  95. }
  96. }
  97. // loop over cage edges
  98. for(int e = 0;e<CE.rows();e++)
  99. {
  100. // loop over domain vertices
  101. for(int i = 0;i<V.rows();i++)
  102. {
  103. // Find samples from tip up to tail
  104. Eigen::VectorXd tip = C.row(P(CE(e,0)));
  105. Eigen::VectorXd tail = C.row(P(CE(e,1)));
  106. // Compute parameter along bone and squared distance
  107. double t,sqrd;
  108. project_to_line(
  109. V(i,0),V(i,1),V(i,2),
  110. tip(0),tip(1),tip(2),
  111. tail(0),tail(1),tail(2),
  112. t,sqrd);
  113. if(t>=-FLOAT_EPS && t<=(1.0f+FLOAT_EPS) && sqrd<=FLOAT_EPS)
  114. {
  115. bci.push_back(i);
  116. bcj.push_back(CE(e,0));
  117. bcv.push_back(1.0-t);
  118. bci.push_back(i);
  119. bcj.push_back(CE(e,1));
  120. bcv.push_back(t);
  121. }
  122. }
  123. }
  124. std::vector<bool> vertices_marked(V.rows(), false);
  125. // loop over cage faces
  126. for(int f = 0;f<CF.rows();f++)
  127. {
  128. Eigen::Vector3d v_0 = C.row(P(CF(f, 0)));
  129. Eigen::Vector3d v_1 = C.row(P(CF(f, 1)));
  130. Eigen::Vector3d v_2 = C.row(P(CF(f, 2)));
  131. Eigen::Vector3d n = (v_1 - v_0).cross(v_2 - v_1);
  132. n.normalize();
  133. // loop over domain vertices
  134. for (int i = 0;i<V.rows();i++)
  135. {
  136. // ensure each vertex is associated with only one face
  137. if (vertices_marked[i])
  138. {
  139. continue;
  140. }
  141. Eigen::Vector3d point = V.row(i);
  142. Eigen::Vector3d v = point - v_0;
  143. double dist = abs(v.dot(n));
  144. if (dist <= 1.e-1f)
  145. {
  146. //barycentric coordinates
  147. Eigen::Vector3d vec_0 = v_1 - v_0, vec_1 = v_2 - v_0, vec_2 = point - v_0;
  148. double d00 = vec_0.dot(vec_0);
  149. double d01 = vec_0.dot(vec_1);
  150. double d11 = vec_1.dot(vec_1);
  151. double d20 = vec_2.dot(vec_0);
  152. double d21 = vec_2.dot(vec_1);
  153. double denom = d00 * d11 - d01 * d01;
  154. double v = (d11 * d20 - d01 * d21) / denom;
  155. double w = (d00 * d21 - d01 * d20) / denom;
  156. double u = 1.0 - v - w;
  157. if (u>=0. && u<=1.0 && v>=0. && v<=1.0 && w >=0. && w<=1.0)
  158. {
  159. vertices_marked[i] = true;
  160. bci.push_back(i);
  161. bcj.push_back(CF(f, 0));
  162. bcv.push_back(u);
  163. bci.push_back(i);
  164. bcj.push_back(CF(f, 1));
  165. bcv.push_back(v);
  166. bci.push_back(i);
  167. bcj.push_back(CF(f, 2));
  168. bcv.push_back(w);
  169. }
  170. }
  171. }
  172. }
  173. // find unique boundary indices
  174. std::vector<int> vb = bci;
  175. std::sort(vb.begin(),vb.end());
  176. vb.erase(unique(vb.begin(), vb.end()), vb.end());
  177. b.resize(vb.size());
  178. bc = Eigen::MatrixXd::Zero(vb.size(),P.size()+BE.rows());
  179. // Map from boundary index to index in boundary
  180. std::map<int,int> bim;
  181. int i = 0;
  182. // Also fill in b
  183. for(std::vector<int>::iterator bit = vb.begin();bit != vb.end();bit++)
  184. {
  185. b(i) = *bit;
  186. bim[*bit] = i;
  187. i++;
  188. }
  189. // Build BC
  190. for(i = 0;i < (int)bci.size();i++)
  191. {
  192. assert(bim.find(bci[i]) != bim.end());
  193. bc(bim[bci[i]],bcj[i]) = bcv[i];
  194. }
  195. // Normalize across rows so that conditions sum to one
  196. for(i = 0;i<bc.rows();i++)
  197. {
  198. double sum = bc.row(i).sum();
  199. assert(sum != 0 && "Some boundary vertex getting all zero BCs");
  200. bc.row(i).array() /= sum;
  201. }
  202. if(bc.size() == 0)
  203. {
  204. verbose("^%s: Error: boundary conditions are empty.\n",__FUNCTION__);
  205. return false;
  206. }
  207. // If there's only a single boundary condition, the following tests
  208. // are overzealous.
  209. if(bc.cols() == 1)
  210. {
  211. // If there is only one weight function,
  212. // then we expect that there is only one handle.
  213. assert(P.rows() + BE.rows() == 1);
  214. return true;
  215. }
  216. // Check that every Weight function has at least one boundary value of 1 and
  217. // one value of 0
  218. for(i = 0;i<bc.cols();i++)
  219. {
  220. double min_abs_c = bc.col(i).array().abs().minCoeff();
  221. double max_c = bc.col(i).maxCoeff();
  222. if(min_abs_c > FLOAT_EPS)
  223. {
  224. verbose("^%s: Error: handle %d does not receive 0 weight\n",__FUNCTION__,i);
  225. return false;
  226. }
  227. if(max_c< (1-FLOAT_EPS))
  228. {
  229. verbose("^%s: Error: handle %d does not receive 1 weight\n",__FUNCTION__,i);
  230. return false;
  231. }
  232. }
  233. return true;
  234. }
  235. #ifdef IGL_STATIC_LIBRARY
  236. // Explicit template instantiation
  237. template bool igl::boundary_conditions<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::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -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&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1>> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1>> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1>> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1>> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1>> const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1>>&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1>>&);
  238. #endif