InElementAABB.h 11 KB

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