boundary_conditions.cpp 4.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172
  1. #include "boundary_conditions.h"
  2. #include "verbose.h"
  3. #include "EPS.h"
  4. #include "project_to_line.h"
  5. #include <vector>
  6. #include <map>
  7. #include <iostream>
  8. IGL_INLINE bool igl::boundary_conditions(
  9. const Eigen::MatrixXd & V ,
  10. const Eigen::MatrixXi & /*Ele*/,
  11. const Eigen::MatrixXd & C ,
  12. const Eigen::VectorXi & P ,
  13. const Eigen::MatrixXi & BE ,
  14. const Eigen::MatrixXi & CE ,
  15. Eigen::VectorXi & b ,
  16. Eigen::MatrixXd & bc )
  17. {
  18. using namespace Eigen;
  19. using namespace igl;
  20. using namespace std;
  21. if(P.size()+BE.rows() == 0)
  22. {
  23. verbose("^%s: Error: no handles found\n",__FUNCTION__);
  24. return false;
  25. }
  26. vector<int> bci;
  27. vector<int> bcj;
  28. vector<int> bcv;
  29. // loop over points
  30. for(int p = 0;p<P.size();p++)
  31. {
  32. VectorXd pos = C.row(P(p));
  33. // loop over domain vertices
  34. for(int i = 0;i<V.rows();i++)
  35. {
  36. // Find samples just on pos
  37. //Vec3 vi(V(i,0),V(i,1),V(i,2));
  38. // EIGEN GOTCHA:
  39. // double sqrd = (V.row(i)-pos).array().pow(2).sum();
  40. // Must first store in temporary
  41. VectorXd vi = V.row(i);
  42. double sqrd = (vi-pos).array().pow(2).sum();
  43. if(sqrd <= FLOAT_EPS)
  44. {
  45. cout<<"sum((["<<
  46. V(i,0)<<" "<<
  47. V(i,1)<<" "<<
  48. V(i,2)<<"] - ["<<
  49. pos(0)<<" "<<
  50. pos(1)<<" "<<
  51. pos(2)<<"]).^2) = "<<sqrd<<endl;
  52. bci.push_back(i);
  53. bcj.push_back(p);
  54. bcv.push_back(1.0);
  55. }
  56. }
  57. }
  58. // loop over bone edges
  59. for(int e = 0;e<BE.rows();e++)
  60. {
  61. // loop over domain vertices
  62. for(int i = 0;i<V.rows();i++)
  63. {
  64. // Find samples from tip up to tail
  65. VectorXd tip = C.row(BE(e,0));
  66. VectorXd tail = C.row(BE(e,1));
  67. // Compute parameter along bone and squared distance
  68. double t,sqrd;
  69. project_to_line(
  70. V(i,0),V(i,1),V(i,2),
  71. tip(0),tip(1),tip(2),
  72. tail(0),tail(1),tail(2),
  73. t,sqrd);
  74. if(t>=-FLOAT_EPS && t<=(1.0f+FLOAT_EPS) && sqrd<=FLOAT_EPS)
  75. {
  76. bci.push_back(i);
  77. bcj.push_back(P.size()+e);
  78. bcv.push_back(1.0);
  79. }
  80. }
  81. }
  82. // Cage edges are not considered yet
  83. // loop over cage edges
  84. for(int e = 0;e<CE.rows();e++)
  85. {
  86. // loop over domain vertices
  87. for(int i = 0;i<V.rows();i++)
  88. {
  89. // Find samples from tip up to tail
  90. VectorXd tip = C.row(P(CE(e,0)));
  91. VectorXd tail = C.row(P(CE(e,1)));
  92. // Compute parameter along bone and squared distance
  93. double t,sqrd;
  94. project_to_line(
  95. V(i,0),V(i,1),V(i,2),
  96. tip(0),tip(1),tip(2),
  97. tail(0),tail(1),tail(2),
  98. t,sqrd);
  99. if(t>=-FLOAT_EPS && t<=(1.0f+FLOAT_EPS) && sqrd<=FLOAT_EPS)
  100. {
  101. bci.push_back(i);
  102. bcj.push_back(CE(e,0));
  103. bcv.push_back(1.0-t);
  104. bci.push_back(i);
  105. bcj.push_back(CE(e,1));
  106. bcv.push_back(t);
  107. }
  108. }
  109. }
  110. // find unique boundary indices
  111. vector<int> vb = bci;
  112. sort(vb.begin(),vb.end());
  113. vb.erase(unique(vb.begin(), vb.end()), vb.end());
  114. b.resize(vb.size());
  115. bc = MatrixXd::Zero(vb.size(),P.size()+BE.rows());
  116. // Map from boundary index to index in boundary
  117. map<int,int> bim;
  118. int i = 0;
  119. // Also fill in b
  120. for(vector<int>::iterator bit = vb.begin();bit != vb.end();bit++)
  121. {
  122. b(i) = *bit;
  123. bim[*bit] = i;
  124. i++;
  125. }
  126. // Build BC
  127. for(i = 0;i < (int)bci.size();i++)
  128. {
  129. assert(bim.find(bci[i]) != bim.end());
  130. bc(bim[bci[i]],bcj[i]) = bcv[i];
  131. }
  132. // Normalize accross rows so that conditions sum to one
  133. for(i = 0;i<bc.rows();i++)
  134. {
  135. double sum = bc.row(i).sum();
  136. assert(sum != 0);
  137. bc.row(i).array() /= sum;
  138. }
  139. // Check that every Weight function has at least one boundary value of 1 and
  140. // one value of 0
  141. for(i = 0;i<bc.cols();i++)
  142. {
  143. double min_abs_c = bc.col(i).array().abs().minCoeff();
  144. double max_c = bc.col(i).maxCoeff();
  145. if(min_abs_c > FLOAT_EPS)
  146. {
  147. verbose("^%s: Error: handle %d does not receive 0 weight\n",__FUNCTION__,i);
  148. return false;
  149. }
  150. if(max_c< (1-FLOAT_EPS))
  151. {
  152. verbose("^%s: Error: handle %d does not receive 1 weight\n",__FUNCTION__,i);
  153. return false;
  154. }
  155. }
  156. return true;
  157. }