principal_curvature.cpp 28 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154
  1. #include "principal_curvature.h"
  2. #include <iostream>
  3. #include <fstream>
  4. #include <iomanip>
  5. #include <queue>
  6. #include <list>
  7. #include <cmath>
  8. #include <limits>
  9. #include <Eigen/SparseCholesky>
  10. // Lib IGL includes
  11. #include <igl/adjacency_list.h>
  12. #include <igl/per_face_normals.h>
  13. #include <igl/per_vertex_normals.h>
  14. #include <igl/vf.h>
  15. using std::setfill;
  16. using std::setw;
  17. using std::cout;
  18. using std::endl;
  19. typedef enum
  20. {
  21. SPHERE_SEARCH,
  22. K_RING_SEARCH
  23. } searchType;
  24. typedef enum
  25. {
  26. AVERAGE,
  27. PROJ_PLANE
  28. } normalType;
  29. class CurvatureCalculator
  30. {
  31. public:
  32. /* Row number i represents the i-th vertex, whose columns are:
  33. curv[i][0] : K2
  34. curv[i][1] : K1
  35. curvDir[i][0] : PD1
  36. curvDir[i][1] : PD2
  37. */
  38. std::vector< std::vector<double> > curv;
  39. std::vector< std::vector<Eigen::Vector3d> > curvDir;
  40. bool curvatureComputed;
  41. class Quadric
  42. {
  43. public:
  44. IGL_INLINE Quadric ()
  45. {
  46. a() = b() = c() = d() = e() = 1.0;
  47. }
  48. IGL_INLINE Quadric(double av, double bv, double cv, double dv, double ev)
  49. {
  50. a() = av;
  51. b() = bv;
  52. c() = cv;
  53. d() = dv;
  54. e() = ev;
  55. }
  56. IGL_INLINE double& a() { return data[0];}
  57. IGL_INLINE double& b() { return data[1];}
  58. IGL_INLINE double& c() { return data[2];}
  59. IGL_INLINE double& d() { return data[3];}
  60. IGL_INLINE double& e() { return data[4];}
  61. double data[5];
  62. IGL_INLINE double evaluate(double u, double v)
  63. {
  64. return a()*u*u + b()*u*v + c()*v*v + d()*u + e()*v;
  65. }
  66. IGL_INLINE double du(double u, double v)
  67. {
  68. return 2.0*a()*u + b()*v + d();
  69. }
  70. IGL_INLINE double dv(double u, double v)
  71. {
  72. return 2.0*c()*v + b()*u + e();
  73. }
  74. IGL_INLINE double duv(double u, double v)
  75. {
  76. return b();
  77. }
  78. IGL_INLINE double duu(double u, double v)
  79. {
  80. return 2.0*a();
  81. }
  82. IGL_INLINE double dvv(double u, double v)
  83. {
  84. return 2.0*c();
  85. }
  86. IGL_INLINE static Quadric fit(vector<Eigen::Vector3d> &VV, bool zeroDetCheck, bool svd)
  87. {
  88. assert(VV.size() >= 5);
  89. Eigen::MatrixXd A(VV.size(),5);
  90. Eigen::MatrixXd b(VV.size(),1);
  91. Eigen::MatrixXd sol(5,1);
  92. for(unsigned int c=0; c < VV.size(); ++c)
  93. {
  94. double u = VV[c][0];
  95. double v = VV[c][1];
  96. double n = VV[c][2];
  97. A(c,0) = u*u;
  98. A(c,1) = u*v;
  99. A(c,2) = v*v;
  100. A(c,3) = u;
  101. A(c,4) = v;
  102. b(c) = n;
  103. }
  104. static int count = 0, index = 0;
  105. double min = 0.000000000001; //1.0e-12
  106. //Devo fare solving svd del sistema A*sol=b
  107. if (zeroDetCheck && ((A.transpose()*A).determinant() < min && (A.transpose()*A).determinant() > -min))
  108. {
  109. printf("Quadric: unsolvable vertex %d %d\n", count, ++index);
  110. sol=A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
  111. return Quadric(sol(0),sol(1),sol(2),sol(3),sol(4));
  112. }
  113. count++;
  114. //for (int i = 0; i < 100; i++)
  115. {
  116. if (svd)
  117. sol=A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
  118. else
  119. sol = ((A.transpose()*A).inverse()*A.transpose())*b;
  120. }
  121. return Quadric(sol(0),sol(1),sol(2),sol(3),sol(4));
  122. }
  123. };
  124. class Quadric6
  125. {
  126. public:
  127. IGL_INLINE Quadric6 ()
  128. {
  129. a() = b() = c() = d() = e() = f() = 1.0;
  130. }
  131. IGL_INLINE Quadric6 (double av, double bv, double cv, double dv, double ev, double fv)
  132. {
  133. a() = av;
  134. b() = bv;
  135. c() = cv;
  136. d() = dv;
  137. e() = ev;
  138. f() = fv;
  139. }
  140. IGL_INLINE double& a() { return data[0];}
  141. IGL_INLINE double& b() { return data[1];}
  142. IGL_INLINE double& c() { return data[2];}
  143. IGL_INLINE double& d() { return data[3];}
  144. IGL_INLINE double& e() { return data[4];}
  145. IGL_INLINE double& f() { return data[5];}
  146. double data[6];
  147. IGL_INLINE double evaluate(double u, double v)
  148. {
  149. return a()*u*u + b()*u*v + c()*v*v + d()*u + e()*v + f();
  150. }
  151. IGL_INLINE double du(double u, double v)
  152. {
  153. return 2.0*a()*u + b()*v + d();
  154. }
  155. IGL_INLINE double dv(double u, double v)
  156. {
  157. return 2.0*c()*v + b()*u + e();
  158. }
  159. IGL_INLINE double duv(double u, double v)
  160. {
  161. return b();
  162. }
  163. IGL_INLINE double duu(double u, double v)
  164. {
  165. return 2.0*a();
  166. }
  167. IGL_INLINE double dvv(double u, double v)
  168. {
  169. return 2.0*c();
  170. }
  171. IGL_INLINE static Quadric6 fit(vector<Eigen::Vector3d> &VV, bool zeroDetCheck, bool svd)
  172. {
  173. assert(VV.size() >= 6);
  174. Eigen::MatrixXd A(VV.size(),6);
  175. Eigen::MatrixXd b(VV.size(),1);
  176. Eigen::MatrixXd sol(VV.size(),1);
  177. for(unsigned int c=0; c < VV.size(); ++c)
  178. {
  179. double u = VV[c][0];
  180. double v = VV[c][1];
  181. double n = VV[c][2];
  182. A(c,0) = u*u;
  183. A(c,1) = u*v;
  184. A(c,2) = v*v;
  185. A(c,3) = u;
  186. A(c,4) = v;
  187. A(c,5) = 1;
  188. b(c) = n;
  189. }
  190. static int count = 0, index = 0;
  191. double min = 0.000000000001; //1.0e-12
  192. /*
  193. if (!count)
  194. printf("GNE %e\n", min);
  195. */
  196. if (zeroDetCheck && ((A.transpose()*A).determinant() < min && (A.transpose()*A).determinant() > -min))
  197. {
  198. //A.svd().solve(b, &sol); A.svd().solve(b, &sol);
  199. //cout << sol << endl;
  200. printf("Quadric: unsolvable vertex %d %d\n", count, ++index);
  201. //return Quadric6 (1, 1, 1, 1, 1, 1);
  202. sol=A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
  203. return Quadric6(sol(0),sol(1),sol(2),sol(3),sol(4),sol(5));
  204. }
  205. count++;
  206. if (svd)
  207. sol=A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
  208. else
  209. sol = ((A.transpose()*A).inverse()*A.transpose())*b;
  210. return Quadric6(sol(0),sol(1),sol(2),sol(3),sol(4),sol(5));
  211. }
  212. };
  213. class Linear
  214. {
  215. public:
  216. IGL_INLINE Linear ()
  217. {
  218. a() = b() = c() = 1.0;
  219. }
  220. IGL_INLINE Linear (double av, double bv, double cv)
  221. {
  222. a() = av;
  223. b() = bv;
  224. c() = cv;
  225. }
  226. IGL_INLINE double& a() { return data[0];}
  227. IGL_INLINE double& b() { return data[1];}
  228. IGL_INLINE double& c() { return data[2];}
  229. double data[3];
  230. IGL_INLINE double evaluate(double u, double v)
  231. {
  232. return a()*u + b()*v + c();
  233. }
  234. IGL_INLINE double du()
  235. {
  236. return a();
  237. }
  238. IGL_INLINE double dv()
  239. {
  240. return b();
  241. }
  242. IGL_INLINE static Linear fit(vector<Eigen::Vector3d> &VV, bool svdRes, bool zeroDetCheck)
  243. {
  244. assert(VV.size() >= 3);
  245. Eigen::MatrixXd A(VV.size(),3);
  246. Eigen::MatrixXd b(VV.size(),1);
  247. Eigen::MatrixXd sol(VV.size(),1);
  248. for(unsigned int c=0; c < VV.size(); ++c)
  249. {
  250. double u = VV[c][0];
  251. double v = VV[c][1];
  252. double n = VV[c][2];
  253. A(c,0) = u;
  254. A(c,1) = v;
  255. A(c,2) = 1;
  256. b(c) = n;
  257. }
  258. static int count = 0, index = 0;
  259. double min = 0.000000000001; //1.0e-12
  260. /*
  261. if (!count)
  262. printf("GNE %e\n", min);
  263. */
  264. if (zeroDetCheck && ((A.transpose()*A).determinant() < min && (A.transpose()*A).determinant() > -min))
  265. {
  266. //A.svd().solve(b, &sol); A.svd().solve(b, &sol);
  267. //cout << sol << endl;
  268. printf("Quadric: unsolvable vertex %d %d\n", count, ++index);
  269. //return Linear (1, 1, 1);
  270. sol=A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
  271. return Linear(sol(0),sol(1),sol(2));
  272. }
  273. count++;
  274. if (svdRes)
  275. sol=A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
  276. else
  277. sol = ((A.transpose()*A).inverse()*A.transpose())*b;
  278. return Linear(sol(0),sol(1),sol(2));
  279. }
  280. };
  281. class Discrete
  282. {
  283. public:
  284. IGL_INLINE Discrete ()
  285. {
  286. a() = b() = 1.0;
  287. }
  288. IGL_INLINE Discrete (double av, double bv)
  289. {
  290. a() = av;
  291. b() = bv;
  292. }
  293. IGL_INLINE double& a() { return data[0];}
  294. IGL_INLINE double& b() { return data[1];}
  295. double data[2];
  296. IGL_INLINE static Discrete fit(vector<Eigen::Vector3d> &VV, bool svdRes, bool detCheck)
  297. {
  298. assert(VV.size() >= 2);
  299. Eigen::MatrixXd A(VV.size(),2);
  300. Eigen::MatrixXd b(VV.size(),1);
  301. Eigen::MatrixXd sol(2/*VV.size()*/,1);
  302. for(unsigned int c=0; c < VV.size(); ++c)
  303. {
  304. double u = VV[c][0];
  305. double v = VV[c][1];
  306. double n = VV[c][2];
  307. A(c,0) = u;
  308. A(c,1) = v;
  309. b(c) = n;
  310. }
  311. static int count = 0, index = 0;
  312. double min = 0.000000000001; //1.0e-12
  313. /*
  314. if (!count)
  315. printf("GNE %e\n", min);
  316. */
  317. if (detCheck && ((A.transpose()*A).determinant() < min && (A.transpose()*A).determinant() > -min))
  318. {
  319. //A.svd().solve(b, &sol); A.svd().solve(b, &sol);
  320. //cout << sol << endl;
  321. printf("Discrete: unsolvable vertex %d %d\n", count, ++index);
  322. //return Discrete (1, 1);
  323. sol=A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
  324. return Discrete (sol(0),sol(1));
  325. }
  326. count++;
  327. if (svdRes)
  328. {
  329. //printf ("SSSVVVDDDDDDIIIIIIII!!!!!1!!!!!!1111!1!!111\n");
  330. sol=A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
  331. }
  332. else
  333. sol = ((A.transpose()*A).inverse()*A.transpose())*b;
  334. return Discrete (sol(0),sol(1));
  335. }
  336. };
  337. class Bicubic
  338. {
  339. public:
  340. IGL_INLINE Bicubic ()
  341. {
  342. a() = b() = c() = d() = 1.0;
  343. }
  344. IGL_INLINE Bicubic (double av, double bv, double cv, double dv)
  345. {
  346. a() = av;
  347. b() = bv;
  348. c() = cv;
  349. d() = dv;
  350. }
  351. IGL_INLINE double& a() { return data[0];}
  352. IGL_INLINE double& b() { return data[1];}
  353. IGL_INLINE double& c() { return data[2];}
  354. IGL_INLINE double& d() { return data[3];}
  355. double data[4];
  356. IGL_INLINE double du()
  357. {
  358. return 6.0 * a();
  359. }
  360. IGL_INLINE double dv()
  361. {
  362. return 6.0 * b();
  363. }
  364. IGL_INLINE static Bicubic fit(vector<Eigen::Vector3d> &VV, bool svdRes, bool detCheck)
  365. {
  366. assert(VV.size() >= 4);
  367. Eigen::MatrixXd A(VV.size(),4);
  368. Eigen::MatrixXd b(VV.size(),1);
  369. Eigen::MatrixXd sol(4,1);
  370. for(unsigned int c=0; c < VV.size(); ++c)
  371. {
  372. double u = VV[c][0];
  373. double v = VV[c][1];
  374. double n = VV[c][2];
  375. A(c,0) = u*u*u;
  376. A(c,1) = v*v*v;
  377. A(c,2) = u*u*v;
  378. A(c,3) = u*v*v;
  379. b(c) = n;
  380. }
  381. static int count = 0, index = 0;
  382. double min = 0.000000000001; //1.0e-12
  383. /*
  384. if (!count)
  385. printf("GNE %e\n", min);
  386. */
  387. if (detCheck && ((A.transpose()*A).determinant() < min && (A.transpose()*A).determinant() > -min))
  388. {
  389. //A.svd().solve(b, &sol); A.svd().solve(b, &sol);
  390. //cout << sol << endl;
  391. printf("Bicubic: unsolvable vertex %d %d\n", count, ++index);
  392. //return Bicubic (1, 1, 1, 1, 1, 1);
  393. sol=A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
  394. return Bicubic(sol(0),sol(1),sol(2),sol(3));
  395. }
  396. count++;
  397. if (svdRes)
  398. sol=A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
  399. else
  400. sol = ((A.transpose()*A).inverse()*A.transpose())*b;
  401. return Bicubic(sol(0),sol(1),sol(2),sol(3));
  402. }
  403. };
  404. public:
  405. Eigen::MatrixXd vertices;
  406. // Face list of current mesh (#F x 3) or (#F x 4)
  407. // The i-th row contains the indices of the vertices that forms the i-th face in ccw order
  408. Eigen::MatrixXi faces;
  409. std::vector<std::vector<int> > vertex_to_vertices;
  410. std::vector<std::vector<int> > vertex_to_faces;
  411. std::vector<std::vector<int> > vertex_to_faces_index;
  412. Eigen::MatrixXd face_normals;
  413. Eigen::MatrixXd vertex_normals;
  414. /* Size of the neighborhood */
  415. double sphereRadius;
  416. int kRing;
  417. bool localMode; /* Use local mode */
  418. bool projectionPlaneCheck; /* Check collected vertices on tangent plane */
  419. bool montecarlo;
  420. bool svd; /* Use svd calculation instead of pseudoinverse */
  421. bool zeroDetCheck; /* Check if the determinant is close to zero */
  422. unsigned int montecarloN;
  423. searchType st; /* Use either a sphere search or a k-ring search */
  424. normalType nt;
  425. double lastRadius;
  426. double scaledRadius;
  427. std::string lastMeshName;
  428. /* Benchmark related variables */
  429. bool expStep; /* True if we want the radius to increase exponentially */
  430. int step; /* If expStep==false, by how much rhe radius increases on every step */
  431. int maxSize; /* The maximum limit of the radius in the benchmark */
  432. IGL_INLINE CurvatureCalculator();
  433. IGL_INLINE void init(const Eigen::MatrixXd& V, const Eigen::MatrixXi& F);
  434. IGL_INLINE void finalEigenStuff (int, vector<Eigen::Vector3d>, Quadric );
  435. IGL_INLINE void fitQuadric (Eigen::Vector3d, std::vector<Eigen::Vector3d> ref, const std::vector<int>& , Quadric *);
  436. IGL_INLINE void applyProjOnPlane(Eigen::Vector3d, std::vector<int>, std::vector<int>&);
  437. IGL_INLINE void getSphere(const int, const double, std::vector<int>&, int min);
  438. IGL_INLINE void getKRing(const int, const double,std::vector<int>&);
  439. IGL_INLINE Eigen::Vector3d project(Eigen::Vector3d, Eigen::Vector3d, Eigen::Vector3d);
  440. IGL_INLINE void computeReferenceFrame(int, Eigen::Vector3d, std::vector<Eigen::Vector3d>&);
  441. IGL_INLINE void getAverageNormal(int, std::vector<int>, Eigen::Vector3d&);
  442. IGL_INLINE void getProjPlane(int, std::vector<int>, Eigen::Vector3d&);
  443. IGL_INLINE void applyMontecarlo(std::vector<int>&,std::vector<int>*);
  444. IGL_INLINE void computeCurvature();
  445. IGL_INLINE void printCurvature(std::string outpath);
  446. IGL_INLINE double getAverageEdge();
  447. IGL_INLINE static int rotateForward (float *v0, float *v1, float *v2)
  448. {
  449. float t;
  450. if (abs(*v2) >= abs(*v1) && abs(*v2) >= abs(*v0))
  451. return 0;
  452. t = *v0;
  453. *v0 = *v2;
  454. *v2 = *v1;
  455. *v1 = t;
  456. return 1 + rotateForward (v0, v1, v2);
  457. }
  458. IGL_INLINE static void rotateBackward (int nr, float *v0, float *v1, float *v2)
  459. {
  460. float t;
  461. if (nr == 0)
  462. return;
  463. t = *v2;
  464. *v2 = *v0;
  465. *v0 = *v1;
  466. *v1 = t;
  467. rotateBackward (nr - 1, v0, v1, v2);
  468. }
  469. IGL_INLINE static Eigen::Vector3d chooseMax (Eigen::Vector3d n, Eigen::Vector3d abc, float ab)
  470. {
  471. int i, max_i;
  472. float max_sp;
  473. Eigen::Vector3d nt[8];
  474. n.normalize ();
  475. abc.normalize ();
  476. max_sp = - std::numeric_limits<float>::max();
  477. for (i = 0; i < 4; i++)
  478. {
  479. nt[i] = n;
  480. if (ab > 0)
  481. {
  482. switch (i)
  483. {
  484. case 0:
  485. break;
  486. case 1:
  487. nt[i][2] = -n[2];
  488. break;
  489. case 2:
  490. nt[i][0] = -n[0];
  491. nt[i][1] = -n[1];
  492. break;
  493. case 3:
  494. nt[i][0] = -n[0];
  495. nt[i][1] = -n[1];
  496. nt[i][2] = -n[2];
  497. break;
  498. }
  499. }
  500. else
  501. {
  502. switch (i)
  503. {
  504. case 0:
  505. nt[i][0] = -n[0];
  506. break;
  507. case 1:
  508. nt[i][1] = -n[1];
  509. break;
  510. case 2:
  511. nt[i][0] = -n[0];
  512. nt[i][2] = -n[2];
  513. break;
  514. case 3:
  515. nt[i][1] = -n[1];
  516. nt[i][2] = -n[2];
  517. break;
  518. }
  519. }
  520. if (nt[i].dot(abc) > max_sp)
  521. {
  522. max_sp = nt[i].dot(abc);
  523. max_i = i;
  524. }
  525. }
  526. return nt[max_i];
  527. }
  528. };
  529. class comparer
  530. {
  531. public:
  532. IGL_INLINE bool operator() (const std::pair<int, double>& lhs, const std::pair<int, double>&rhs) const
  533. {
  534. return lhs.second>rhs.second;
  535. }
  536. };
  537. IGL_INLINE CurvatureCalculator::CurvatureCalculator()
  538. {
  539. this->localMode=false;
  540. this->projectionPlaneCheck=false;
  541. this->sphereRadius=5;
  542. this->st=SPHERE_SEARCH;
  543. this->nt=PROJ_PLANE;
  544. this->montecarlo=false;
  545. this->montecarloN=0;
  546. this->kRing=3;
  547. this->svd=true;
  548. this->zeroDetCheck=true;
  549. this->curvatureComputed=false;
  550. this->expStep=true;
  551. }
  552. IGL_INLINE void CurvatureCalculator::init(const Eigen::MatrixXd& V, const Eigen::MatrixXi& F)
  553. {
  554. vertices = V;
  555. faces = F;
  556. igl::adjacency_list(F, vertex_to_vertices);
  557. igl::vf(V, F, vertex_to_faces, vertex_to_faces_index);
  558. igl::per_face_normals(V, F, face_normals);
  559. igl::per_vertex_normals(V, F, face_normals, vertex_normals);
  560. }
  561. IGL_INLINE void CurvatureCalculator::fitQuadric (Eigen::Vector3d v, std::vector<Eigen::Vector3d> ref, const std::vector<int>& vv, Quadric *q)
  562. {
  563. vector<Eigen::Vector3d> points;
  564. points.reserve (vv.size());
  565. for (unsigned int i = 0; i < vv.size(); ++i) {
  566. Eigen::Vector3d cp = vertices.row(vv[i]);
  567. // vtang non e` il v tangente!!!
  568. Eigen::Vector3d vTang = cp - v;
  569. double x = vTang.dot(ref[0]);
  570. double y = vTang.dot(ref[1]);
  571. double z = vTang.dot(ref[2]);
  572. points.push_back(Eigen::Vector3d (x,y,z));
  573. }
  574. *q = Quadric::fit (points, zeroDetCheck, svd);
  575. }
  576. IGL_INLINE void CurvatureCalculator::finalEigenStuff (int i, vector<Eigen::Vector3d> ref, Quadric q)
  577. {
  578. double a = q.a();
  579. double b = q.b();
  580. double c = q.c();
  581. double d = q.d();
  582. double e = q.e();
  583. double E = 1.0 + d*d;
  584. double F = d*e;
  585. double G = 1.0 + e*e;
  586. Eigen::Vector3d n = Eigen::Vector3d(-d,-e,1.0).normalized();
  587. double L = 2.0 * a * n[2];
  588. double M = b * n[2];
  589. double N = 2 * c * n[2];
  590. // ----------------- Eigen stuff
  591. Eigen::Matrix2d m;
  592. m << L*G - M*F, M*E-L*F, M*E-L*F, N*E-M*F;
  593. m = m / (E*G-F*F);
  594. Eigen::SelfAdjointEigenSolver<Eigen::Matrix2d> eig(m);
  595. Eigen::Vector2d c_val = eig.eigenvalues();
  596. Eigen::Matrix2d c_vec = eig.eigenvectors();
  597. c_val = -c_val;
  598. Eigen::Vector3d v1, v2;
  599. v1[0] = c_vec(0);
  600. v1[1] = c_vec(1);
  601. v1[2] = d * v1[0] + e * v1[1];
  602. v2[0] = c_vec(2);
  603. v2[1] = c_vec(3);
  604. v2[2] = d * v2[0] + e * v2[1];
  605. v1 = v1.normalized();
  606. v2 = v2.normalized();
  607. Eigen::Vector3d v1global = ref[0] * v1[0] + ref[1] * v1[1] + ref[2] * v1[2];
  608. Eigen::Vector3d v2global = ref[0] * v2[0] + ref[1] * v2[1] + ref[2] * v2[2];
  609. v1global.normalize();
  610. v2global.normalize();
  611. v1global *= c_val(0);
  612. v2global *= c_val(1);
  613. if (c_val[0] > c_val[1])
  614. {
  615. curv[i]=std::vector<double>(2);
  616. curv[i][0]=c_val(1);
  617. curv[i][1]=c_val(0);
  618. curvDir[i]=std::vector<Eigen::Vector3d>(2);
  619. curvDir[i][0]=v2global;
  620. curvDir[i][1]=v1global;
  621. }
  622. else
  623. {
  624. curv[i]=std::vector<double>(2);
  625. curv[i][0]=c_val(0);
  626. curv[i][1]=c_val(1);
  627. curvDir[i]=std::vector<Eigen::Vector3d>(2);
  628. curvDir[i][0]=v1global;
  629. curvDir[i][1]=v2global;
  630. }
  631. // ---- end Eigen stuff
  632. }
  633. IGL_INLINE void CurvatureCalculator::getKRing(const int start, const double r, std::vector<int>&vv)
  634. {
  635. int bufsize=vertices.rows();
  636. vv.reserve(bufsize);
  637. std::list<std::pair<int,int> > queue;
  638. bool* visited = (bool*)calloc(bufsize,sizeof(bool));
  639. queue.push_back(std::pair<int,int>(start,0));
  640. visited[start]=true;
  641. while (!queue.empty())
  642. {
  643. int toVisit=queue.front().first;
  644. int distance=queue.front().second;
  645. queue.pop_front();
  646. vv.push_back(toVisit);
  647. if (distance<(int)r)
  648. {
  649. for (unsigned int i=0; i<vertex_to_vertices[toVisit].size(); i++)
  650. {
  651. int neighbor=vertex_to_vertices[toVisit][i];
  652. if (!visited[neighbor])
  653. {
  654. queue.push_back(std::pair<int,int> (neighbor,distance+1));
  655. visited[neighbor]=true;
  656. }
  657. }
  658. }
  659. }
  660. free(visited);
  661. return;
  662. }
  663. IGL_INLINE void CurvatureCalculator::getSphere(const int start, const double r, std::vector<int> &vv, int min)
  664. {
  665. int bufsize=vertices.rows();
  666. vv.reserve(bufsize);
  667. std::list<int>* queue= new std::list<int>();
  668. bool* visited = (bool*)calloc(bufsize,sizeof(bool));
  669. queue->push_back(start);
  670. visited[start]=true;
  671. Eigen::Vector3d me=vertices.row(start);
  672. std::priority_queue<std::pair<int, double>, vector<std::pair<int, double> >, comparer >* extra_candidates= new std::priority_queue<std::pair<int, double>, vector<std::pair<int, double> >, comparer >();
  673. while (!queue->empty())
  674. {
  675. int toVisit=queue->front();
  676. queue->pop_front();
  677. vv.push_back(toVisit);
  678. for (unsigned int i=0; i<vertex_to_vertices[toVisit].size(); i++)
  679. {
  680. int neighbor=vertex_to_vertices[toVisit][i];
  681. if (!visited[neighbor])
  682. {
  683. Eigen::Vector3d neigh=vertices.row(neighbor);
  684. float distance=(me-neigh).norm();
  685. if (distance<r)
  686. queue->push_back(neighbor);
  687. else if (vv.size()<min)
  688. extra_candidates->push(std::pair<int,double>(neighbor,distance));
  689. visited[neighbor]=true;
  690. }
  691. }
  692. }
  693. while (!extra_candidates->empty() && vv.size()<min)
  694. {
  695. std::pair<int, double> cand=extra_candidates->top();
  696. extra_candidates->pop();
  697. vv.push_back(cand.first);
  698. for (unsigned int i=0; i<vertex_to_vertices[cand.first].size(); i++)
  699. {
  700. int neighbor=vertex_to_vertices[cand.first][i];
  701. if (!visited[neighbor])
  702. {
  703. Eigen::Vector3d neigh=vertices.row(neighbor);
  704. float distance=(me-neigh).norm();
  705. extra_candidates->push(std::pair<int,double>(neighbor,distance));
  706. visited[neighbor]=true;
  707. }
  708. }
  709. }
  710. free(extra_candidates);
  711. free(queue);
  712. free(visited);
  713. }
  714. IGL_INLINE inline Eigen::Vector3d CurvatureCalculator::project(Eigen::Vector3d v, Eigen::Vector3d vp, Eigen::Vector3d ppn)
  715. {
  716. return (vp - (ppn * ((vp - v).dot(ppn))));
  717. }
  718. IGL_INLINE void CurvatureCalculator::computeReferenceFrame(int i, Eigen::Vector3d normal, std::vector<Eigen::Vector3d>& ref )
  719. {
  720. Eigen::Vector3d longest_v=Eigen::Vector3d::Zero();
  721. longest_v=Eigen::Vector3d(vertices.row(vertex_to_vertices[i][0]));
  722. longest_v=(project(vertices.row(i),longest_v,normal)-Eigen::Vector3d(vertices.row(i))).normalized();
  723. /* L'ultimo asse si ottiene come prodotto vettoriale tra i due
  724. * calcolati */
  725. Eigen::Vector3d y_axis=(normal.cross(longest_v)).normalized();
  726. ref[0]=longest_v;
  727. ref[1]=y_axis;
  728. ref[2]=normal;
  729. }
  730. IGL_INLINE void CurvatureCalculator::getAverageNormal(int j, std::vector<int> vv, Eigen::Vector3d& normal)
  731. {
  732. normal=(vertex_normals.row(j)).normalized();
  733. if (localMode)
  734. return;
  735. for (unsigned int i=0; i<vv.size(); i++)
  736. {
  737. normal+=vertex_normals.row(vv[i]).normalized();
  738. }
  739. normal.normalize();
  740. }
  741. IGL_INLINE void CurvatureCalculator::getProjPlane(int j, std::vector<int> vv, Eigen::Vector3d& ppn)
  742. {
  743. int nr;
  744. float a, b, c;
  745. float nx, ny, nz;
  746. float abcq;
  747. a = b = c = 0;
  748. if (localMode)
  749. {
  750. for (unsigned int i=0; i<vertex_to_faces.at(j).size(); ++i)
  751. {
  752. Eigen::Vector3d faceNormal=face_normals.row(vertex_to_faces.at(j).at(i));
  753. a += faceNormal[0];
  754. b += faceNormal[1];
  755. c += faceNormal[2];
  756. }
  757. }
  758. else
  759. {
  760. for (unsigned int i=0; i<vv.size(); ++i)
  761. {
  762. a+= vertex_normals.row(vv[i])[0];
  763. b+= vertex_normals.row(vv[i])[1];
  764. c+= vertex_normals.row(vv[i])[2];
  765. }
  766. }
  767. nr = rotateForward (&a, &b, &c);
  768. abcq = a*a + b*b + c*c;
  769. nx = sqrt (a*a / abcq);
  770. ny = sqrt (b*b / abcq);
  771. nz = sqrt (1 - nx*nx - ny*ny);
  772. rotateBackward (nr, &a, &b, &c);
  773. rotateBackward (nr, &nx, &ny, &nz);
  774. ppn = chooseMax (Eigen::Vector3d(nx, ny, nz), Eigen::Vector3d (a, b, c), a * b);
  775. ppn.normalize();
  776. }
  777. IGL_INLINE double CurvatureCalculator::getAverageEdge()
  778. {
  779. double sum = 0;
  780. int count = 0;
  781. for (int i = 0; i<faces.rows(); i++)
  782. {
  783. for (short unsigned j=0; j<3; j++)
  784. {
  785. Eigen::Vector3d p1=vertices.row(faces.row(i)[j]);
  786. Eigen::Vector3d p2=vertices.row(faces.row(i)[(j+1)%3]);
  787. double l = (p1-p2).norm();
  788. sum+=l;
  789. ++count;
  790. }
  791. }
  792. return (sum/(double)count);
  793. }
  794. IGL_INLINE void CurvatureCalculator::applyProjOnPlane(Eigen::Vector3d ppn, std::vector<int> vin, std::vector<int> &vout)
  795. {
  796. for (vector<int>::iterator vpi = vin.begin(); vpi != vin.end(); ++vpi)
  797. if (vertex_normals.row(*vpi) * ppn > 0.0f)
  798. vout.push_back (*vpi);
  799. }
  800. IGL_INLINE void CurvatureCalculator::applyMontecarlo(std::vector<int>& vin, std::vector<int> *vout)
  801. {
  802. if (montecarloN >= vin.size ())
  803. {
  804. *vout = vin;
  805. return;
  806. }
  807. float p = ((float) montecarloN) / (float) vin.size();
  808. for (std::vector<int>::iterator vpi = vin.begin(); vpi != vin.end(); ++vpi)
  809. {
  810. float r;
  811. if ((r = ((float)rand () / RAND_MAX)) < p)
  812. {
  813. vout->push_back (*vpi);
  814. }
  815. }
  816. }
  817. IGL_INLINE void CurvatureCalculator::computeCurvature()
  818. {
  819. //CHECK che esista la mesh
  820. size_t vertices_count=vertices.rows() ;
  821. if (vertices_count <=0)
  822. return;
  823. curvDir=std::vector< std::vector<Eigen::Vector3d> >(vertices_count);
  824. curv=std::vector<std::vector<double> >(vertices_count);
  825. scaledRadius=getAverageEdge()*sphereRadius;
  826. std::vector<int> vv;
  827. std::vector<int> vvtmp;
  828. Eigen::Vector3d normal;
  829. double time_spent;
  830. double searchtime=0, ref_time=0, fit_time=0, final_time=0;
  831. for (size_t i=0; i<vertices_count; ++i)
  832. {
  833. vv.clear();
  834. vvtmp.clear();
  835. Eigen::Vector3d me=vertices.row(i);
  836. switch (st)
  837. {
  838. case SPHERE_SEARCH:
  839. getSphere(i,scaledRadius,vv,6);
  840. break;
  841. case K_RING_SEARCH:
  842. getKRing(i,kRing,vv);
  843. break;
  844. default:
  845. fprintf(stderr,"Error: search type not recognized");
  846. return;
  847. }
  848. std::vector<Eigen::Vector3d> ref(3);
  849. if (vv.size()<6)
  850. {
  851. std::cerr << "Could not compute curvature of radius " << scaledRadius << endl;
  852. return;
  853. }
  854. switch (nt)
  855. {
  856. case AVERAGE:
  857. getAverageNormal(i,vv,normal);
  858. break;
  859. case PROJ_PLANE:
  860. getProjPlane(i,vv,normal);
  861. break;
  862. default:
  863. fprintf(stderr,"Error: normal type not recognized");
  864. return;
  865. }
  866. if (projectionPlaneCheck)
  867. {
  868. vvtmp.reserve (vv.size ());
  869. applyProjOnPlane (normal, vv, vvtmp);
  870. if (vvtmp.size() >= 6)
  871. vv = vvtmp;
  872. }
  873. if (vv.size()<6)
  874. {
  875. std::cerr << "Could not compute curvature of radius " << scaledRadius << endl;
  876. return;
  877. }
  878. if (montecarlo)
  879. {
  880. if(montecarloN<6)
  881. break;
  882. vvtmp.reserve(vv.size());
  883. applyMontecarlo(vv,&vvtmp);
  884. vv=vvtmp;
  885. }
  886. if (vv.size()<6)
  887. return;
  888. computeReferenceFrame(i,normal,ref);
  889. Quadric q;
  890. fitQuadric (me, ref, vv, &q);
  891. finalEigenStuff(i,ref,q);
  892. }
  893. lastRadius=sphereRadius;
  894. curvatureComputed=true;
  895. }
  896. IGL_INLINE void CurvatureCalculator::printCurvature(std::string outpath)
  897. {
  898. if (!curvatureComputed)
  899. return;
  900. std::ofstream of;
  901. of.open(outpath.c_str());
  902. if (!of)
  903. {
  904. fprintf(stderr, "Error: could not open output file %s\n", outpath.c_str());
  905. return;
  906. }
  907. int vertices_count=vertices.rows();
  908. of << vertices_count << endl;
  909. for (int i=0; i<vertices_count; i++)
  910. {
  911. of << curv[i][0] << " " << curv[i][1] << " " << curvDir[i][0][0] << " " << curvDir[i][0][1] << " " << curvDir[i][0][2] << " " <<
  912. curvDir[i][1][0] << " " << curvDir[i][1][1] << " " << curvDir[i][1][2] << endl;
  913. }
  914. of.close();
  915. }
  916. template <typename DerivedV, typename DerivedF>
  917. IGL_INLINE void igl::principal_curvature(
  918. const Eigen::PlainObjectBase<DerivedV>& V,
  919. const Eigen::PlainObjectBase<DerivedF>& F,
  920. Eigen::PlainObjectBase<DerivedV>& PD1,
  921. Eigen::PlainObjectBase<DerivedV>& PD2,
  922. unsigned radius
  923. )
  924. {
  925. // Preallocate memory
  926. PD1.resize(V.rows(),3);
  927. PD2.resize(V.rows(),3);
  928. // Precomputation
  929. CurvatureCalculator cc;
  930. cc.init(V.template cast<double>(),F.template cast<int>());
  931. cc.sphereRadius = radius;
  932. // Compute
  933. cc.computeCurvature();
  934. // Copy it back
  935. for (unsigned i=0; i<V.rows(); i++)
  936. {
  937. Eigen::Vector3d d1;
  938. Eigen::Vector3d d2;
  939. d1 << cc.curvDir[i][0][0], cc.curvDir[i][0][1], cc.curvDir[i][0][2];
  940. d2 << cc.curvDir[i][1][0], cc.curvDir[i][1][1], cc.curvDir[i][1][2];
  941. d1.normalize();
  942. d2.normalize();
  943. d1 *= cc.curv[i][0];
  944. d2 *= cc.curv[i][1];
  945. PD1.row(i) = d1;
  946. PD2.row(i) = d2;
  947. }
  948. }