InElementAABB.h 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407
  1. #ifndef IGL_IN_ELEMENT_AABB_H
  2. #define IGL_IN_ELEMENT_AABB_H
  3. #include <Eigen/Core>
  4. #include <memory>
  5. #include <vector>
  6. namespace igl
  7. {
  8. class InElementAABB
  9. {
  10. public:
  11. std::shared_ptr<InElementAABB> m_left, m_right;
  12. Eigen::RowVectorXd m_bb_min,m_bb_max;
  13. // -1 non-leaf
  14. int m_element;
  15. InElementAABB():
  16. m_left(NULL), m_right(NULL),
  17. m_bb_min(), m_bb_max(), m_element(-1)
  18. {}
  19. // Build an Axis-Aligned Bounding Box tree for a given mesh and given
  20. // serialization of a previous AABB tree.
  21. //
  22. // Inputs:
  23. // V #V by dim list of mesh vertex positions.
  24. // Ele #Ele by dim+1 list of mesh indices into #V.
  25. // bb_mins max_tree by dim list of bounding box min corner positions
  26. // bb_maxs max_tree by dim list of bounding box max corner positions
  27. // elements max_tree list of element or (not leaf id) indices into Ele
  28. // i recursive call index {0}
  29. inline void init(
  30. const Eigen::MatrixXd & V,
  31. const Eigen::MatrixXi & Ele,
  32. const Eigen::MatrixXd & bb_mins,
  33. const Eigen::MatrixXd & bb_maxs,
  34. const Eigen::VectorXi & elements,
  35. const int i = 0);
  36. // Wrapper for root with empty serialization
  37. inline void init(
  38. const Eigen::MatrixXd & V,
  39. const Eigen::MatrixXi & Ele);
  40. // Build an Axis-Aligned Bounding Box tree for a given mesh.
  41. //
  42. // Inputs:
  43. // V #V by dim list of mesh vertex positions.
  44. // Ele #Ele by dim+1 list of mesh indices into #V.
  45. // SI #Ele by dim list revealing for each coordinate where Ele's
  46. // barycenters would be sorted: SI(e,d) = i --> the dth coordinate of
  47. // the barycenter of the eth element would be placed at position i in a
  48. // sorted list.
  49. // I #I list of indices into Ele of elements to include (for recursive
  50. // calls)
  51. //
  52. inline void init(
  53. const Eigen::MatrixXd & V,
  54. const Eigen::MatrixXi & Ele,
  55. const Eigen::MatrixXi & SI,
  56. const Eigen::VectorXi & I);
  57. // Find the indices of elements containing given point.
  58. //
  59. // Inputs:
  60. // V #V by dim list of mesh vertex positions. **Should be same as used to
  61. // construct mesh.**
  62. // Ele #Ele by dim+1 list of mesh indices into #V. **Should be same as used to
  63. // construct mesh.**
  64. // q dim row-vector query position
  65. // first whether to only return first element containing q
  66. // Returns:
  67. // list of indices of elements containing q
  68. inline std::vector<int> find(
  69. const Eigen::MatrixXd & V,
  70. const Eigen::MatrixXi & Ele,
  71. const Eigen::RowVectorXd & q,
  72. const bool first=false) const;
  73. // If number of elements m then total tree size should be 2*h where h is
  74. // the deepest depth 2^ceil(log(#Ele*2-1))
  75. inline int subtree_size();
  76. // Serialize this class into 3 arrays (so we can pass it pack to matlab)
  77. //
  78. // Outputs:
  79. // bb_mins max_tree by dim list of bounding box min corner positions
  80. // bb_maxs max_tree by dim list of bounding box max corner positions
  81. // elements max_tree list of element or (not leaf id) indices into Ele
  82. // i recursive call index into these arrays {0}
  83. inline void serialize(
  84. Eigen::MatrixXd & bb_mins,
  85. Eigen::MatrixXd & bb_maxs,
  86. Eigen::VectorXi & elements,
  87. const int i = 0);
  88. };
  89. }
  90. // Implementation
  91. #include <igl/volume.h>
  92. #include <igl/colon.h>
  93. #include <igl/doublearea.h>
  94. #include <igl/matlab_format.h>
  95. #include <igl/colon.h>
  96. #include <igl/sort.h>
  97. #include <igl/barycenter.h>
  98. #include <iostream>
  99. inline void igl::InElementAABB::init(
  100. const Eigen::MatrixXd & V,
  101. const Eigen::MatrixXi & Ele,
  102. const Eigen::MatrixXd & bb_mins,
  103. const Eigen::MatrixXd & bb_maxs,
  104. const Eigen::VectorXi & elements,
  105. const int i)
  106. {
  107. using namespace std;
  108. using namespace Eigen;
  109. using namespace igl;
  110. if(bb_mins.size() > 0)
  111. {
  112. assert(bb_mins.rows() == bb_maxs.rows() && "Serial tree arrays must match");
  113. assert(bb_mins.cols() == V.cols() && "Serial tree array dim must match V");
  114. assert(bb_mins.cols() == bb_maxs.cols() && "Serial tree arrays must match");
  115. assert(bb_mins.rows() == elements.rows() &&
  116. "Serial tree arrays must match");
  117. // construct from serialization
  118. m_bb_min = bb_mins.row(i);
  119. m_bb_max = bb_maxs.row(i);
  120. m_element = elements(i);
  121. // Not leaf then recurse
  122. if(m_element == -1)
  123. {
  124. m_left = make_shared<InElementAABB>();
  125. m_left->init( V,Ele,bb_mins,bb_maxs,elements,2*i+1);
  126. m_right = make_shared<InElementAABB>();
  127. m_right->init( V,Ele,bb_mins,bb_maxs,elements,2*i+2);
  128. }
  129. }else
  130. {
  131. VectorXi allI = colon<int>(0,Ele.rows()-1);
  132. MatrixXd BC;
  133. barycenter(V,Ele,BC);
  134. MatrixXi SI(BC.rows(),BC.cols());
  135. {
  136. MatrixXd _;
  137. MatrixXi IS;
  138. igl::sort(BC,1,true,_,IS);
  139. // Need SI(i) to tell which place i would be sorted into
  140. const int dim = IS.cols();
  141. for(int i = 0;i<IS.rows();i++)
  142. {
  143. for(int d = 0;d<dim;d++)
  144. {
  145. SI(IS(i,d),d) = i;
  146. }
  147. }
  148. }
  149. init(V,Ele,SI,allI);
  150. }
  151. }
  152. inline void igl::InElementAABB::init(
  153. const Eigen::MatrixXd & V,
  154. const Eigen::MatrixXi & Ele)
  155. {
  156. using namespace Eigen;
  157. return init(V,Ele,MatrixXd(),MatrixXd(),VectorXi(),0);
  158. }
  159. inline void igl::InElementAABB::init(
  160. const Eigen::MatrixXd & V,
  161. const Eigen::MatrixXi & Ele,
  162. const Eigen::MatrixXi & SI,
  163. const Eigen::VectorXi & I)
  164. {
  165. using namespace Eigen;
  166. using namespace std;
  167. using namespace igl;
  168. const int dim = V.cols();
  169. const double inf = numeric_limits<double>::infinity();
  170. m_bb_min.setConstant(1,dim, inf);
  171. m_bb_max.setConstant(1,dim,-inf);
  172. // Compute bounding box
  173. for(int i = 0;i<I.rows();i++)
  174. {
  175. for(int c = 0;c<Ele.cols();c++)
  176. {
  177. for(int d = 0;d<dim;d++)
  178. {
  179. m_bb_min(d) = min(m_bb_min(d),V(Ele(I(i),c),d));
  180. m_bb_max(d) = max(m_bb_max(d),V(Ele(I(i),c),d));
  181. }
  182. }
  183. }
  184. switch(I.size())
  185. {
  186. case 0:
  187. {
  188. assert(false);
  189. }
  190. case 1:
  191. {
  192. m_element = I(0);
  193. break;
  194. }
  195. default:
  196. {
  197. // Compute longest direction
  198. int max_d = -1;
  199. double max_len = -inf;
  200. for(int d = 0;d<dim;d++)
  201. {
  202. const auto diff = (m_bb_max[d] - m_bb_min[d]);
  203. if( diff > max_len )
  204. {
  205. max_len = diff;
  206. max_d = d;
  207. }
  208. }
  209. // Can't use median on BC directly because many may have same value,
  210. // but can use median on sorted BC indices
  211. VectorXi SIdI(I.rows());
  212. for(int i = 0;i<I.rows();i++)
  213. {
  214. SIdI(i) = SI(I(i),max_d);
  215. }
  216. // Since later I use <= I think I don't need to worry about odd/even
  217. // Pass by copy to avoid changing input
  218. const auto median = [](VectorXi A)->double
  219. {
  220. size_t n = A.size()/2;
  221. nth_element(A.data(),A.data()+n,A.data()+A.size());
  222. if(A.rows() % 2 == 1)
  223. {
  224. return A(n);
  225. }else
  226. {
  227. nth_element(A.data(),A.data()+n-1,A.data()+A.size());
  228. return 0.5*(A(n)+A(n-1));
  229. }
  230. };
  231. const double med = median(SIdI);
  232. VectorXi LI((I.rows()+1)/2),RI(I.rows()/2);
  233. assert(LI.rows()+RI.rows() == I.rows());
  234. // Distribute left and right
  235. {
  236. int li = 0;
  237. int ri = 0;
  238. for(int i = 0;i<I.rows();i++)
  239. {
  240. if(SIdI(i)<=med)
  241. {
  242. LI(li++) = I(i);
  243. }else
  244. {
  245. RI(ri++) = I(i);
  246. }
  247. }
  248. }
  249. if(LI.rows()>0)
  250. {
  251. m_left = make_shared<InElementAABB>();
  252. m_left->init(V,Ele,SI,LI);
  253. }
  254. if(RI.rows()>0)
  255. {
  256. m_right = make_shared<InElementAABB>();
  257. m_right->init(V,Ele,SI,RI);
  258. }
  259. }
  260. }
  261. }
  262. inline std::vector<int> igl::InElementAABB::find(
  263. const Eigen::MatrixXd & V,
  264. const Eigen::MatrixXi & Ele,
  265. const Eigen::RowVectorXd & q,
  266. const bool first) const
  267. {
  268. using namespace std;
  269. using namespace igl;
  270. using namespace Eigen;
  271. bool inside = true;
  272. const int dim = m_bb_max.size();
  273. assert(q.size() == m_bb_max.size());
  274. const double epsilon = 1e-14;
  275. // Check if outside bounding box
  276. for(int d = 0;d<q.size()&&inside;d++)
  277. {
  278. inside &= (q(d)-m_bb_min(d))>=epsilon;
  279. inside &= (m_bb_max(d)-q(d))>=epsilon;
  280. }
  281. cout<<"searching..."<<endl;
  282. if(!inside)
  283. {
  284. cout<<"not in bb"<<endl;
  285. return std::vector<int>();
  286. }
  287. if(m_element != -1)
  288. {
  289. // Initialize to some value > -epsilon
  290. double a1=1,a2=1,a3=1,a4=1;
  291. switch(dim)
  292. {
  293. case 3:
  294. {
  295. // Barycentric coordinates
  296. const RowVector3d V1 = V.row(Ele(m_element,0));
  297. const RowVector3d V2 = V.row(Ele(m_element,1));
  298. const RowVector3d V3 = V.row(Ele(m_element,2));
  299. const RowVector3d V4 = V.row(Ele(m_element,3));
  300. a1 = volume_single(V2,V4,V3,(RowVector3d)q);
  301. a2 = volume_single(V1,V3,V4,(RowVector3d)q);
  302. a3 = volume_single(V1,V4,V2,(RowVector3d)q);
  303. a4 = volume_single(V1,V2,V3,(RowVector3d)q);
  304. break;
  305. }
  306. case 2:
  307. {
  308. // Barycentric coordinates
  309. const Vector2d V1 = V.row(Ele(m_element,0));
  310. const Vector2d V2 = V.row(Ele(m_element,1));
  311. const Vector2d V3 = V.row(Ele(m_element,2));
  312. double a0 = doublearea_single(V1,V2,V3);
  313. a1 = doublearea_single(V1,V2,(Vector2d)q);
  314. a2 = doublearea_single(V2,V3,(Vector2d)q);
  315. a3 = doublearea_single(V3,V1,(Vector2d)q);
  316. cout<<
  317. a0<<" "<<
  318. a1<<" "<<
  319. a2<<" "<<
  320. a3<<" "<<
  321. endl;
  322. break;
  323. }
  324. default:assert(false);
  325. }
  326. if(
  327. a1>=-epsilon &&
  328. a2>=-epsilon &&
  329. a3>=-epsilon &&
  330. a4>=-epsilon)
  331. {
  332. return std::vector<int>(1,m_element);
  333. }else
  334. {
  335. return std::vector<int>();
  336. }
  337. }
  338. std::vector<int> left = m_left->find(V,Ele,q,first);
  339. if(first && !left.empty())
  340. {
  341. return left;
  342. }
  343. std::vector<int> right = m_right->find(V,Ele,q,first);
  344. if(first)
  345. {
  346. return right;
  347. }
  348. left.insert(left.end(),right.begin(),right.end());
  349. return left;
  350. }
  351. inline int igl::InElementAABB::subtree_size()
  352. {
  353. // 1 for self
  354. int n = 1;
  355. int n_left = 0,n_right = 0;
  356. if(m_left != NULL)
  357. {
  358. n_left = m_left->subtree_size();
  359. }
  360. if(m_right != NULL)
  361. {
  362. n_right = m_right->subtree_size();
  363. }
  364. n += 2*std::max(n_left,n_right);
  365. return n;
  366. }
  367. inline void igl::InElementAABB::serialize(
  368. Eigen::MatrixXd & bb_mins,
  369. Eigen::MatrixXd & bb_maxs,
  370. Eigen::VectorXi & elements,
  371. const int i)
  372. {
  373. using namespace std;
  374. // Calling for root then resize output
  375. if(i==0)
  376. {
  377. const int m = subtree_size();
  378. //cout<<"m: "<<m<<endl;
  379. bb_mins.resize(m,m_bb_min.size());
  380. bb_maxs.resize(m,m_bb_max.size());
  381. elements.resize(m,1);
  382. }
  383. //cout<<i<<" ";
  384. bb_mins.row(i) = m_bb_min;
  385. bb_maxs.row(i) = m_bb_max;
  386. elements(i) = m_element;
  387. if(m_left != NULL)
  388. {
  389. m_left->serialize(bb_mins,bb_maxs,elements,2*i+1);
  390. }
  391. if(m_right != NULL)
  392. {
  393. m_right->serialize(bb_mins,bb_maxs,elements,2*i+2);
  394. }
  395. }
  396. #endif