AABB.h 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 Alec Jacobson <alecjacobson@gmail.com>
  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. #ifndef IGL_AABB_H
  9. #define IGL_AABB_H
  10. // Implementation of semi-general purpose axis-aligned bounding box hierarchy.
  11. // The mesh (V,Ele) is stored and managed by the caller and each routine here
  12. // simply takes it as references (it better not change between calls).
  13. //
  14. #include <Eigen/Core>
  15. #include <Eigen/Geometry>
  16. #include <vector>
  17. namespace igl
  18. {
  19. template <typename DerivedV, int DIM>
  20. class AABB
  21. {
  22. public:
  23. typedef typename DerivedV::Scalar Scalar;
  24. typedef Eigen::Matrix<Scalar,1,DIM> RowVectorDIMS;
  25. typedef Eigen::Matrix<Scalar,DIM,1> VectorDIMS;
  26. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,DIM> MatrixXDIMS;
  27. // Shared pointers are slower...
  28. AABB * m_left;
  29. AABB * m_right;
  30. Eigen::AlignedBox<Scalar,DIM> m_box;
  31. // -1 non-leaf
  32. int m_primitive;
  33. AABB():
  34. m_left(NULL), m_right(NULL),
  35. m_box(), m_primitive(-1)
  36. {}
  37. // http://stackoverflow.com/a/3279550/148668
  38. AABB(const AABB& other):
  39. m_left(other.m_left ? new AABB(*other.m_left) : NULL),
  40. m_right(other.m_right ? new AABB(*other.m_right) : NULL),
  41. m_box(other.m_box),
  42. m_primitive(other.m_primitive)
  43. {
  44. }
  45. // copy-swap idiom
  46. friend void swap(AABB& first, AABB& second)
  47. {
  48. // Enable ADL
  49. using std::swap;
  50. swap(first.m_left,second.m_left);
  51. swap(first.m_right,second.m_right);
  52. swap(first.m_box,second.m_box);
  53. swap(first.m_primitive,second.m_primitive);
  54. }
  55. // Pass-by-value (aka copy)
  56. AABB& operator=(AABB other)
  57. {
  58. swap(*this,other);
  59. return *this;
  60. }
  61. AABB(AABB&& other):
  62. // initialize via default constructor
  63. AABB()
  64. {
  65. swap(*this,other);
  66. }
  67. // Seems like there should have been an elegant solution to this using
  68. // the copy-swap idiom above:
  69. inline void deinit()
  70. {
  71. m_primitive = -1;
  72. m_box = Eigen::AlignedBox<Scalar,DIM>();
  73. delete m_left;
  74. m_left = NULL;
  75. delete m_right;
  76. m_right = NULL;
  77. }
  78. ~AABB()
  79. {
  80. deinit();
  81. }
  82. // Build an Axis-Aligned Bounding Box tree for a given mesh and given
  83. // serialization of a previous AABB tree.
  84. //
  85. // Inputs:
  86. // V #V by dim list of mesh vertex positions.
  87. // Ele #Ele by dim+1 list of mesh indices into #V.
  88. // bb_mins max_tree by dim list of bounding box min corner positions
  89. // bb_maxs max_tree by dim list of bounding box max corner positions
  90. // elements max_tree list of element or (not leaf id) indices into Ele
  91. // i recursive call index {0}
  92. inline void init(
  93. const Eigen::PlainObjectBase<DerivedV> & V,
  94. const Eigen::MatrixXi & Ele,
  95. const MatrixXDIMS & bb_mins,
  96. const MatrixXDIMS & bb_maxs,
  97. const Eigen::VectorXi & elements,
  98. const int i = 0);
  99. // Wrapper for root with empty serialization
  100. inline void init(
  101. const Eigen::PlainObjectBase<DerivedV> & V,
  102. const Eigen::MatrixXi & Ele);
  103. // Build an Axis-Aligned Bounding Box tree for a given mesh.
  104. //
  105. // Inputs:
  106. // V #V by dim list of mesh vertex positions.
  107. // Ele #Ele by dim+1 list of mesh indices into #V.
  108. // SI #Ele by dim list revealing for each coordinate where Ele's
  109. // barycenters would be sorted: SI(e,d) = i --> the dth coordinate of
  110. // the barycenter of the eth element would be placed at position i in a
  111. // sorted list.
  112. // I #I list of indices into Ele of elements to include (for recursive
  113. // calls)
  114. //
  115. inline void init(
  116. const Eigen::PlainObjectBase<DerivedV> & V,
  117. const Eigen::MatrixXi & Ele,
  118. const Eigen::MatrixXi & SI,
  119. const Eigen::VectorXi & I);
  120. // Return whether at leaf node
  121. inline bool is_leaf() const;
  122. // Find the indices of elements containing given point.
  123. //
  124. // Inputs:
  125. // V #V by dim list of mesh vertex positions. **Should be same as used to
  126. // construct mesh.**
  127. // Ele #Ele by dim+1 list of mesh indices into #V. **Should be same as used to
  128. // construct mesh.**
  129. // q dim row-vector query position
  130. // first whether to only return first element containing q
  131. // Returns:
  132. // list of indices of elements containing q
  133. inline std::vector<int> find(
  134. const Eigen::PlainObjectBase<DerivedV> & V,
  135. const Eigen::MatrixXi & Ele,
  136. const RowVectorDIMS & q,
  137. const bool first=false) const;
  138. // If number of elements m then total tree size should be 2*h where h is
  139. // the deepest depth 2^ceil(log(#Ele*2-1))
  140. inline int subtree_size() const;
  141. // Serialize this class into 3 arrays (so we can pass it pack to matlab)
  142. //
  143. // Outputs:
  144. // bb_mins max_tree by dim list of bounding box min corner positions
  145. // bb_maxs max_tree by dim list of bounding box max corner positions
  146. // elements max_tree list of element or (not leaf id) indices into Ele
  147. // i recursive call index into these arrays {0}
  148. inline void serialize(
  149. MatrixXDIMS & bb_mins,
  150. MatrixXDIMS & bb_maxs,
  151. Eigen::VectorXi & elements,
  152. const int i = 0) const;
  153. // Compute squared distance to a query point
  154. //
  155. // Inputs:
  156. // V #V by dim list of vertex positions
  157. // Ele #Ele by dim list of simplex indices
  158. // P 3 list of query point coordinates
  159. // min_sqr_d current minimum squared distance (only find distances
  160. // less than this)
  161. // Outputs:
  162. // I #P list of facet indices corresponding to smallest distances
  163. // C #P by 3 list of closest points
  164. // Returns squared distance
  165. //
  166. // Known bugs: currently assumes Elements are triangles regardless of
  167. // dimension.
  168. inline Scalar squared_distance(
  169. const Eigen::PlainObjectBase<DerivedV> & V,
  170. const Eigen::MatrixXi & Ele,
  171. const RowVectorDIMS & p,
  172. int & i,
  173. RowVectorDIMS & c) const;
  174. private:
  175. inline Scalar squared_distance(
  176. const Eigen::PlainObjectBase<DerivedV> & V,
  177. const Eigen::MatrixXi & Ele,
  178. const RowVectorDIMS & p,
  179. const Scalar min_sqr_d,
  180. int & i,
  181. RowVectorDIMS & c) const;
  182. public:
  183. template <
  184. typename DerivedP,
  185. typename DerivedsqrD,
  186. typename DerivedI,
  187. typename DerivedC>
  188. inline void squared_distance(
  189. const Eigen::PlainObjectBase<DerivedV> & V,
  190. const Eigen::MatrixXi & Ele,
  191. const Eigen::PlainObjectBase<DerivedP> & P,
  192. Eigen::PlainObjectBase<DerivedsqrD> & sqrD,
  193. Eigen::PlainObjectBase<DerivedI> & I,
  194. Eigen::PlainObjectBase<DerivedC> & C) const;
  195. private:
  196. // Helper function for leaves: works in-place on sqr_d
  197. inline void leaf_squared_distance(
  198. const Eigen::PlainObjectBase<DerivedV> & V,
  199. const Eigen::MatrixXi & Ele,
  200. const RowVectorDIMS & p,
  201. Scalar & sqr_d,
  202. int & i,
  203. RowVectorDIMS & c) const;
  204. inline void set_min(
  205. const RowVectorDIMS & p,
  206. const Scalar sqr_d_candidate,
  207. const int i_candidate,
  208. const RowVectorDIMS & c_candidate,
  209. Scalar & sqr_d,
  210. int & i,
  211. RowVectorDIMS & c) const;
  212. public:
  213. template <int SS>
  214. static
  215. void barycentric_coordinates(
  216. const RowVectorDIMS & p,
  217. const RowVectorDIMS & a,
  218. const RowVectorDIMS & b,
  219. const RowVectorDIMS & c,
  220. Eigen::Matrix<Scalar,1,SS> & bary);
  221. };
  222. }
  223. // Implementation
  224. #include "EPS.h"
  225. #include "barycenter.h"
  226. #include "colon.h"
  227. #include "colon.h"
  228. #include "doublearea.h"
  229. #include "matlab_format.h"
  230. #include "project_to_line_segment.h"
  231. #include "sort.h"
  232. #include "volume.h"
  233. #include <iostream>
  234. #include <limits>
  235. template <typename DerivedV, int DIM>
  236. inline void igl::AABB<DerivedV,DIM>::init(
  237. const Eigen::PlainObjectBase<DerivedV> & V,
  238. const Eigen::MatrixXi & Ele,
  239. const MatrixXDIMS & bb_mins,
  240. const MatrixXDIMS & bb_maxs,
  241. const Eigen::VectorXi & elements,
  242. const int i)
  243. {
  244. using namespace std;
  245. using namespace Eigen;
  246. if(bb_mins.size() > 0)
  247. {
  248. assert(bb_mins.rows() == bb_maxs.rows() && "Serial tree arrays must match");
  249. assert(bb_mins.cols() == V.cols() && "Serial tree array dim must match V");
  250. assert(bb_mins.cols() == bb_maxs.cols() && "Serial tree arrays must match");
  251. assert(bb_mins.rows() == elements.rows() &&
  252. "Serial tree arrays must match");
  253. // construct from serialization
  254. m_box.extend(bb_mins.row(i).transpose());
  255. m_box.extend(bb_maxs.row(i).transpose());
  256. m_primitive = elements(i);
  257. // Not leaf then recurse
  258. if(m_primitive == -1)
  259. {
  260. m_left = new AABB();
  261. m_left->init( V,Ele,bb_mins,bb_maxs,elements,2*i+1);
  262. m_right = new AABB();
  263. m_right->init( V,Ele,bb_mins,bb_maxs,elements,2*i+2);
  264. }
  265. }else
  266. {
  267. VectorXi allI = colon<int>(0,Ele.rows()-1);
  268. MatrixXDIMS BC;
  269. barycenter(V,Ele,BC);
  270. MatrixXi SI(BC.rows(),BC.cols());
  271. {
  272. MatrixXDIMS _;
  273. MatrixXi IS;
  274. igl::sort(BC,1,true,_,IS);
  275. // Need SI(i) to tell which place i would be sorted into
  276. const int dim = IS.cols();
  277. for(int i = 0;i<IS.rows();i++)
  278. {
  279. for(int d = 0;d<dim;d++)
  280. {
  281. SI(IS(i,d),d) = i;
  282. }
  283. }
  284. }
  285. init(V,Ele,SI,allI);
  286. }
  287. }
  288. template <typename DerivedV, int DIM>
  289. inline void igl::AABB<DerivedV,DIM>::init(
  290. const Eigen::PlainObjectBase<DerivedV> & V,
  291. const Eigen::MatrixXi & Ele)
  292. {
  293. using namespace Eigen;
  294. return init(V,Ele,MatrixXDIMS(),MatrixXDIMS(),VectorXi(),0);
  295. }
  296. template <typename DerivedV, int DIM>
  297. inline void igl::AABB<DerivedV,DIM>::init(
  298. const Eigen::PlainObjectBase<DerivedV> & V,
  299. const Eigen::MatrixXi & Ele,
  300. const Eigen::MatrixXi & SI,
  301. const Eigen::VectorXi & I)
  302. {
  303. using namespace Eigen;
  304. using namespace std;
  305. assert(DIM == V.cols() && "V.cols() should matched declared dimension");
  306. const Scalar inf = numeric_limits<Scalar>::infinity();
  307. m_box = AlignedBox<Scalar,DIM>();
  308. // Compute bounding box
  309. for(int i = 0;i<I.rows();i++)
  310. {
  311. for(int c = 0;c<Ele.cols();c++)
  312. {
  313. m_box.extend(V.row(Ele(I(i),c)).transpose());
  314. m_box.extend(V.row(Ele(I(i),c)).transpose());
  315. }
  316. }
  317. switch(I.size())
  318. {
  319. case 0:
  320. {
  321. assert(false);
  322. }
  323. case 1:
  324. {
  325. m_primitive = I(0);
  326. break;
  327. }
  328. default:
  329. {
  330. // Compute longest direction
  331. int max_d = -1;
  332. m_box.diagonal().maxCoeff(&max_d);
  333. // Can't use median on BC directly because many may have same value,
  334. // but can use median on sorted BC indices
  335. VectorXi SIdI(I.rows());
  336. for(int i = 0;i<I.rows();i++)
  337. {
  338. SIdI(i) = SI(I(i),max_d);
  339. }
  340. // Since later I use <= I think I don't need to worry about odd/even
  341. // Pass by copy to avoid changing input
  342. const auto median = [](VectorXi A)->Scalar
  343. {
  344. size_t n = A.size()/2;
  345. nth_element(A.data(),A.data()+n,A.data()+A.size());
  346. if(A.rows() % 2 == 1)
  347. {
  348. return A(n);
  349. }else
  350. {
  351. nth_element(A.data(),A.data()+n-1,A.data()+A.size());
  352. return 0.5*(A(n)+A(n-1));
  353. }
  354. };
  355. const Scalar med = median(SIdI);
  356. VectorXi LI((I.rows()+1)/2),RI(I.rows()/2);
  357. assert(LI.rows()+RI.rows() == I.rows());
  358. // Distribute left and right
  359. {
  360. int li = 0;
  361. int ri = 0;
  362. for(int i = 0;i<I.rows();i++)
  363. {
  364. if(SIdI(i)<=med)
  365. {
  366. LI(li++) = I(i);
  367. }else
  368. {
  369. RI(ri++) = I(i);
  370. }
  371. }
  372. }
  373. if(LI.rows()>0)
  374. {
  375. m_left = new AABB();
  376. m_left->init(V,Ele,SI,LI);
  377. }
  378. if(RI.rows()>0)
  379. {
  380. m_right = new AABB();
  381. m_right->init(V,Ele,SI,RI);
  382. }
  383. }
  384. }
  385. }
  386. template <typename DerivedV, int DIM>
  387. inline bool igl::AABB<DerivedV,DIM>::is_leaf() const
  388. {
  389. return m_primitive != -1;
  390. }
  391. template <typename DerivedV, int DIM>
  392. inline std::vector<int> igl::AABB<DerivedV,DIM>::find(
  393. const Eigen::PlainObjectBase<DerivedV> & V,
  394. const Eigen::MatrixXi & Ele,
  395. const RowVectorDIMS & q,
  396. const bool first) const
  397. {
  398. using namespace std;
  399. using namespace Eigen;
  400. assert(q.size() == DIM &&
  401. "Query dimension should match aabb dimension");
  402. assert(Ele.cols() == V.cols()+1 &&
  403. "AABB::find only makes sense for (d+1)-simplices");
  404. const Scalar epsilon = igl::EPS<Scalar>();
  405. // Check if outside bounding box
  406. bool inside = m_box.contains(q.transpose());
  407. if(!inside)
  408. {
  409. return std::vector<int>();
  410. }
  411. assert(m_primitive==-1 || (m_left == NULL && m_right == NULL));
  412. if(is_leaf())
  413. {
  414. // Initialize to some value > -epsilon
  415. Scalar a1=1,a2=1,a3=1,a4=1;
  416. switch(DIM)
  417. {
  418. case 3:
  419. {
  420. // Barycentric coordinates
  421. typedef Eigen::Matrix<Scalar,1,3> RowVector3S;
  422. const RowVector3S V1 = V.row(Ele(m_primitive,0));
  423. const RowVector3S V2 = V.row(Ele(m_primitive,1));
  424. const RowVector3S V3 = V.row(Ele(m_primitive,2));
  425. const RowVector3S V4 = V.row(Ele(m_primitive,3));
  426. a1 = volume_single(V2,V4,V3,(RowVector3S)q);
  427. a2 = volume_single(V1,V3,V4,(RowVector3S)q);
  428. a3 = volume_single(V1,V4,V2,(RowVector3S)q);
  429. a4 = volume_single(V1,V2,V3,(RowVector3S)q);
  430. break;
  431. }
  432. case 2:
  433. {
  434. // Barycentric coordinates
  435. typedef Eigen::Matrix<Scalar,2,1> Vector2S;
  436. const Vector2S V1 = V.row(Ele(m_primitive,0));
  437. const Vector2S V2 = V.row(Ele(m_primitive,1));
  438. const Vector2S V3 = V.row(Ele(m_primitive,2));
  439. // Hack for now to keep templates simple. If becomes bottleneck
  440. // consider using std::enable_if_t
  441. const Vector2S q2 = q.head(2);
  442. Scalar a0 = doublearea_single(V1,V2,V3);
  443. a1 = doublearea_single(V1,V2,q2);
  444. a2 = doublearea_single(V2,V3,q2);
  445. a3 = doublearea_single(V3,V1,q2);
  446. break;
  447. }
  448. default:assert(false);
  449. }
  450. if(
  451. a1>=-epsilon &&
  452. a2>=-epsilon &&
  453. a3>=-epsilon &&
  454. a4>=-epsilon)
  455. {
  456. return std::vector<int>(1,m_primitive);
  457. }else
  458. {
  459. return std::vector<int>();
  460. }
  461. }
  462. std::vector<int> left = m_left->find(V,Ele,q,first);
  463. if(first && !left.empty())
  464. {
  465. return left;
  466. }
  467. std::vector<int> right = m_right->find(V,Ele,q,first);
  468. if(first)
  469. {
  470. return right;
  471. }
  472. left.insert(left.end(),right.begin(),right.end());
  473. return left;
  474. }
  475. template <typename DerivedV, int DIM>
  476. inline int igl::AABB<DerivedV,DIM>::subtree_size() const
  477. {
  478. // 1 for self
  479. int n = 1;
  480. int n_left = 0,n_right = 0;
  481. if(m_left != NULL)
  482. {
  483. n_left = m_left->subtree_size();
  484. }
  485. if(m_right != NULL)
  486. {
  487. n_right = m_right->subtree_size();
  488. }
  489. n += 2*std::max(n_left,n_right);
  490. return n;
  491. }
  492. template <typename DerivedV, int DIM>
  493. inline void igl::AABB<DerivedV,DIM>::serialize(
  494. MatrixXDIMS & bb_mins,
  495. MatrixXDIMS & bb_maxs,
  496. Eigen::VectorXi & elements,
  497. const int i) const
  498. {
  499. using namespace std;
  500. using namespace Eigen;
  501. // Calling for root then resize output
  502. if(i==0)
  503. {
  504. const int m = subtree_size();
  505. //cout<<"m: "<<m<<endl;
  506. bb_mins.resize(m,DIM);
  507. bb_maxs.resize(m,DIM);
  508. elements.resize(m,1);
  509. }
  510. //cout<<i<<" ";
  511. bb_mins.row(i) = m_box.min();
  512. bb_maxs.row(i) = m_box.max();
  513. elements(i) = m_primitive;
  514. if(m_left != NULL)
  515. {
  516. m_left->serialize(bb_mins,bb_maxs,elements,2*i+1);
  517. }
  518. if(m_right != NULL)
  519. {
  520. m_right->serialize(bb_mins,bb_maxs,elements,2*i+2);
  521. }
  522. }
  523. template <typename DerivedV, int DIM>
  524. inline typename igl::AABB<DerivedV,DIM>::Scalar
  525. igl::AABB<DerivedV,DIM>::squared_distance(
  526. const Eigen::PlainObjectBase<DerivedV> & V,
  527. const Eigen::MatrixXi & Ele,
  528. const RowVectorDIMS & p,
  529. int & i,
  530. RowVectorDIMS & c) const
  531. {
  532. return squared_distance(V,Ele,p,std::numeric_limits<Scalar>::infinity(),i,c);
  533. }
  534. template <typename DerivedV, int DIM>
  535. inline typename igl::AABB<DerivedV,DIM>::Scalar
  536. igl::AABB<DerivedV,DIM>::squared_distance(
  537. const Eigen::PlainObjectBase<DerivedV> & V,
  538. const Eigen::MatrixXi & Ele,
  539. const RowVectorDIMS & p,
  540. Scalar min_sqr_d,
  541. int & i,
  542. RowVectorDIMS & c) const
  543. {
  544. using namespace Eigen;
  545. using namespace std;
  546. using namespace igl;
  547. Scalar sqr_d = min_sqr_d;
  548. assert(DIM == 3 && "Code has only been tested for DIM == 3");
  549. assert((Ele.cols() == 3 || Ele.cols() == 2 || Ele.cols() == 1)
  550. && "Code has only been tested for simplex sizes 3,2,1");
  551. assert(m_primitive==-1 || (m_left == NULL && m_right == NULL));
  552. if(is_leaf())
  553. {
  554. leaf_squared_distance(V,Ele,p,sqr_d,i,c);
  555. }else
  556. {
  557. bool looked_left = false;
  558. bool looked_right = false;
  559. const auto & look_left = [&]()
  560. {
  561. int i_left;
  562. RowVectorDIMS c_left;
  563. Scalar sqr_d_left = m_left->squared_distance(V,Ele,p,sqr_d,i_left,c_left);
  564. set_min(p,sqr_d_left,i_left,c_left,sqr_d,i,c);
  565. looked_left = true;
  566. };
  567. const auto & look_right = [&]()
  568. {
  569. int i_right;
  570. RowVectorDIMS c_right;
  571. Scalar sqr_d_right = m_right->squared_distance(V,Ele,p,sqr_d,i_right,c_right);
  572. set_min(p,sqr_d_right,i_right,c_right,sqr_d,i,c);
  573. looked_right = true;
  574. };
  575. // must look left or right if in box
  576. if(m_left->m_box.contains(p.transpose()))
  577. {
  578. look_left();
  579. }
  580. if(m_right->m_box.contains(p.transpose()))
  581. {
  582. look_right();
  583. }
  584. // if haven't looked left and could be less than current min, then look
  585. Scalar left_min_sqr_d = m_left->m_box.squaredExteriorDistance(p.transpose());
  586. Scalar right_min_sqr_d = m_right->m_box.squaredExteriorDistance(p.transpose());
  587. if(left_min_sqr_d < right_min_sqr_d)
  588. {
  589. if(!looked_left && left_min_sqr_d<sqr_d)
  590. {
  591. look_left();
  592. }
  593. if( !looked_right && right_min_sqr_d<sqr_d)
  594. {
  595. look_right();
  596. }
  597. }else
  598. {
  599. if( !looked_right && right_min_sqr_d<sqr_d)
  600. {
  601. look_right();
  602. }
  603. if(!looked_left && left_min_sqr_d<sqr_d)
  604. {
  605. look_left();
  606. }
  607. }
  608. }
  609. return sqr_d;
  610. }
  611. template <typename DerivedV, int DIM>
  612. template <
  613. typename DerivedP,
  614. typename DerivedsqrD,
  615. typename DerivedI,
  616. typename DerivedC>
  617. inline void igl::AABB<DerivedV,DIM>::squared_distance(
  618. const Eigen::PlainObjectBase<DerivedV> & V,
  619. const Eigen::MatrixXi & Ele,
  620. const Eigen::PlainObjectBase<DerivedP> & P,
  621. Eigen::PlainObjectBase<DerivedsqrD> & sqrD,
  622. Eigen::PlainObjectBase<DerivedI> & I,
  623. Eigen::PlainObjectBase<DerivedC> & C) const
  624. {
  625. assert(P.cols() == V.cols() && "cols in P should match dim of cols in V");
  626. sqrD.resize(P.rows(),1);
  627. I.resize(P.rows(),1);
  628. C.resize(P.rows(),P.cols());
  629. for(int p = 0;p<P.rows();p++)
  630. {
  631. RowVectorDIMS Pp = P.row(p), c;
  632. int Ip;
  633. sqrD(p) = squared_distance(V,Ele,Pp,Ip,c);
  634. I(p) = Ip;
  635. C.row(p) = c;
  636. }
  637. }
  638. template <typename DerivedV, int DIM>
  639. inline void igl::AABB<DerivedV,DIM>::leaf_squared_distance(
  640. const Eigen::PlainObjectBase<DerivedV> & V,
  641. const Eigen::MatrixXi & Ele,
  642. const RowVectorDIMS & p,
  643. Scalar & sqr_d,
  644. int & i,
  645. RowVectorDIMS & c) const
  646. {
  647. using namespace Eigen;
  648. using namespace igl;
  649. using namespace std;
  650. // Simplex size
  651. const size_t ss = Ele.cols();
  652. // Only one element per node
  653. // plane unit normal
  654. bool inside_triangle = false;
  655. Scalar d_j = std::numeric_limits<Scalar>::infinity();
  656. RowVectorDIMS pp;
  657. // Only consider triangles, and non-degenerate triangles at that
  658. if(ss == 3 &&
  659. Ele(m_primitive,0) != Ele(m_primitive,1) &&
  660. Ele(m_primitive,1) != Ele(m_primitive,2) &&
  661. Ele(m_primitive,2) != Ele(m_primitive,0))
  662. {
  663. const RowVectorDIMS v10 = (V.row(Ele(m_primitive,1))- V.row(Ele(m_primitive,0)));
  664. const RowVectorDIMS v20 = (V.row(Ele(m_primitive,2))- V.row(Ele(m_primitive,0)));
  665. const RowVectorDIMS n = v10.cross(v20);
  666. Scalar n_norm = n.norm();
  667. if(n_norm > 0)
  668. {
  669. const RowVectorDIMS un = n/n.norm();
  670. // vector to plane
  671. const RowVectorDIMS bc =
  672. 1./3.*
  673. ( V.row(Ele(m_primitive,0))+
  674. V.row(Ele(m_primitive,1))+
  675. V.row(Ele(m_primitive,2)));
  676. const auto & v = p-bc;
  677. // projected point on plane
  678. d_j = v.dot(un);
  679. pp = p - d_j*un;
  680. // determine if pp is inside triangle
  681. Eigen::Matrix<Scalar,1,3> b;
  682. barycentric_coordinates(
  683. pp,
  684. V.row(Ele(m_primitive,0)),
  685. V.row(Ele(m_primitive,1)),
  686. V.row(Ele(m_primitive,2)),
  687. b);
  688. inside_triangle = fabs(fabs(b(0)) + fabs(b(1)) + fabs(b(2)) - 1.) <= 1e-10;
  689. }
  690. }
  691. const auto & point_point_squared_distance = [&](const RowVectorDIMS & s)
  692. {
  693. const Scalar sqr_d_s = (p-s).squaredNorm();
  694. set_min(p,sqr_d_s,m_primitive,s,sqr_d,i,c);
  695. };
  696. if(inside_triangle)
  697. {
  698. // point-triangle squared distance
  699. const Scalar sqr_d_j = d_j*d_j;
  700. //cout<<"point-triangle..."<<endl;
  701. set_min(p,sqr_d_j,m_primitive,pp,sqr_d,i,c);
  702. }else
  703. {
  704. if(ss >= 2)
  705. {
  706. // point-segment distance
  707. // number of edges
  708. size_t ne = ss==3?3:1;
  709. for(size_t x = 0;x<ne;x++)
  710. {
  711. const size_t e1 = Ele(m_primitive,(x+1)%ss);
  712. const size_t e2 = Ele(m_primitive,(x+2)%ss);
  713. const RowVectorDIMS & s = V.row(e1);
  714. const RowVectorDIMS & d = V.row(e2);
  715. // Degenerate edge
  716. if(e1 == e2 || (s-d).squaredNorm()==0)
  717. {
  718. // only consider once
  719. if(e1 < e2)
  720. {
  721. point_point_squared_distance(s);
  722. }
  723. continue;
  724. }
  725. Matrix<Scalar,1,1> sqr_d_j_x(1,1);
  726. Matrix<Scalar,1,1> t(1,1);
  727. project_to_line_segment(p,s,d,t,sqr_d_j_x);
  728. const RowVectorDIMS q = s+t(0)*(d-s);
  729. set_min(p,sqr_d_j_x(0),m_primitive,q,sqr_d,i,c);
  730. }
  731. }else
  732. {
  733. // So then Ele is just a list of points...
  734. assert(ss == 1);
  735. const RowVectorDIMS & s = V.row(Ele(m_primitive,0));
  736. point_point_squared_distance(s);
  737. }
  738. }
  739. }
  740. template <typename DerivedV, int DIM>
  741. inline void igl::AABB<DerivedV,DIM>::set_min(
  742. const RowVectorDIMS & p,
  743. const Scalar sqr_d_candidate,
  744. const int i_candidate,
  745. const RowVectorDIMS & c_candidate,
  746. Scalar & sqr_d,
  747. int & i,
  748. RowVectorDIMS & c) const
  749. {
  750. #ifndef NDEBUG
  751. const Scalar diff = fabs(sqr_d_candidate - (p-c_candidate).squaredNorm());
  752. assert(diff<=1e-10 && "distance should match norm of difference");
  753. #endif
  754. if(sqr_d_candidate < sqr_d)
  755. {
  756. i = i_candidate;
  757. c = c_candidate;
  758. sqr_d = sqr_d_candidate;
  759. }
  760. }
  761. template <typename DerivedV, int DIM>
  762. template <int SS>
  763. inline void
  764. igl::AABB<DerivedV,DIM>::barycentric_coordinates(
  765. const RowVectorDIMS & p,
  766. const RowVectorDIMS & a,
  767. const RowVectorDIMS & b,
  768. const RowVectorDIMS & c,
  769. Eigen::Matrix<Scalar,1,SS> & bary)
  770. {
  771. assert(SS==3);
  772. // http://gamedev.stackexchange.com/a/23745
  773. const RowVectorDIMS v0 = b - a;
  774. const RowVectorDIMS v1 = c - a;
  775. const RowVectorDIMS v2 = p - a;
  776. Scalar d00 = v0.dot(v0);
  777. Scalar d01 = v0.dot(v1);
  778. Scalar d11 = v1.dot(v1);
  779. Scalar d20 = v2.dot(v0);
  780. Scalar d21 = v2.dot(v1);
  781. Scalar denom = d00 * d11 - d01 * d01;
  782. bary(1) = (d11 * d20 - d01 * d21) / denom;
  783. bary(2) = (d00 * d21 - d01 * d20) / denom;
  784. bary(0) = 1.0f - bary(1) - bary(2);
  785. };
  786. #endif