main.cpp 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347
  1. #include <igl/pathinfo.h>
  2. #include <igl/readOBJ.h>
  3. #include <igl/readOFF.h>
  4. #include <igl/readMESH.h>
  5. #include <igl/sample_edges.h>
  6. #include <igl/cat.h>
  7. #include <igl/faces_first.h>
  8. #include <igl/readTGF.h>
  9. #include <igl/tetgen/tetrahedralize.h>
  10. #include <igl/launch_medit.h>
  11. #include <igl/boundary_conditions.h>
  12. #include <igl/mosek/bbw.h>
  13. #include <igl/writeDMAT.h>
  14. #include <igl/writeMESH.h>
  15. #include <Eigen/Dense>
  16. #include <Eigen/Dense>
  17. #include <iostream>
  18. #include <string>
  19. // Whether medit program is install
  20. const bool WITH_MEDIT = true;
  21. const char * USAGE=
  22. "Usage:\n"
  23. " ./bbw_demo shape{.obj|.off|.mesh} skeleton{.tgf|.bf}\n"
  24. ;
  25. // Read a surface mesh from a {.obj|.off|.mesh} files
  26. // Inputs:
  27. // mesh_filename path to {.obj|.off|.mesh} file
  28. // Outputs:
  29. // V #V by 3 list of mesh vertex positions
  30. // F #F by 3 list of triangle indices
  31. // Returns true only if successfuly able to read file
  32. bool load_mesh_from_file(
  33. const std::string mesh_filename,
  34. Eigen::MatrixXd & V,
  35. Eigen::MatrixXi & F)
  36. {
  37. using namespace std;
  38. using namespace igl;
  39. using namespace Eigen;
  40. string dirname, basename, extension, filename;
  41. pathinfo(mesh_filename,dirname,basename,extension,filename);
  42. transform(extension.begin(), extension.end(), extension.begin(), ::tolower);
  43. bool success = false;
  44. if(extension == "obj")
  45. {
  46. success = readOBJ(mesh_filename,V,F);
  47. }else if(extension == "off")
  48. {
  49. success = readOFF(mesh_filename,V,F);
  50. }else if(extension == "mesh")
  51. {
  52. // Unused Tets read from .mesh file
  53. MatrixXi Tets;
  54. success = readMESH(mesh_filename,V,Tets,F);
  55. // We're not going to use any input tets. Only the surface
  56. if(Tets.size() > 0 && F.size() == 0)
  57. {
  58. // If Tets read, but no faces then use surface of tet volume
  59. }else
  60. {
  61. // Rearrange vertices so that faces come first
  62. VectorXi IM;
  63. faces_first(V,F,IM);
  64. // Dont' bother reordering Tets, but this is how one would:
  65. //Tets =
  66. // Tets.unaryExpr(bind1st(mem_fun( static_cast<VectorXi::Scalar&
  67. // (VectorXi::*)(VectorXi::Index)>(&VectorXi::operator())),
  68. // &IM)).eval();
  69. // Don't throw away any interior vertices, since user may want weights
  70. // there
  71. }
  72. }else
  73. {
  74. cerr<<"Error: Unknown shape file format extension: ."<<extension<<endl;
  75. return false;
  76. }
  77. return success;
  78. }
  79. // Load a skeleton (bones, points and cage edges) from a {.bf|.tgf} file
  80. //
  81. // Inputs:
  82. // skel_filename path to skeleton {.bf|.tgf} file
  83. // Outputs:
  84. // C # vertices by 3 list of vertex positions
  85. // P # point-handles list of point handle indices
  86. // BE # bone-edges by 2 list of bone-edge indices
  87. // CE # cage-edges by 2 list of cage-edge indices
  88. bool load_skeleton_from_file(
  89. const std::string skel_filename,
  90. Eigen::MatrixXd & C,
  91. Eigen::VectorXi & P,
  92. Eigen::MatrixXi & BE,
  93. Eigen::MatrixXi & CE)
  94. {
  95. using namespace std;
  96. using namespace igl;
  97. using namespace Eigen;
  98. string dirname, basename, extension, filename;
  99. pathinfo(skel_filename,dirname,basename,extension,filename);
  100. transform(extension.begin(), extension.end(), extension.begin(), ::tolower);
  101. bool success = false;
  102. if(extension == "tgf")
  103. {
  104. // Phony space for unused all edges and pseudo edges
  105. MatrixXi E;
  106. MatrixXi PE;
  107. success = readTGF(skel_filename,C,E,P,BE,CE,PE);
  108. }else
  109. {
  110. cerr<<"Error: Unknown skeleton file format extension: ."<<extension<<endl;
  111. return false;
  112. }
  113. return success;
  114. }
  115. // Mesh the interior of a given surface with tetrahedra which are graded (tend
  116. // to be small near the surface and large inside) and conform to the given
  117. // handles and samplings thereof.
  118. //
  119. // Inputs:
  120. // V #V by 3 list of mesh vertex positions
  121. // F #F by 3 list of triangle indices
  122. // C #C by 3 list of vertex positions
  123. // P #P list of point handle indices
  124. // BE #BE by 2 list of bone-edge indices
  125. // CE #CE by 2 list of cage-edge indices
  126. // Outputs:
  127. // VV #VV by 3 list of tet-mesh vertex positions
  128. // TT #TT by 4 list of tetrahedra indices
  129. // FF #FF by 3 list of surface triangle indices
  130. // Returns true only on success
  131. bool mesh_with_skeleton(
  132. const Eigen::MatrixXd & V,
  133. const Eigen::MatrixXi & F,
  134. const Eigen::MatrixXd & C,
  135. const Eigen::VectorXi & /*P*/,
  136. const Eigen::MatrixXi & BE,
  137. const Eigen::MatrixXi & CE,
  138. Eigen::MatrixXd & VV,
  139. Eigen::MatrixXi & TT,
  140. Eigen::MatrixXi & FF)
  141. {
  142. using namespace Eigen;
  143. using namespace igl;
  144. using namespace std;
  145. // Collect all edges that need samples:
  146. MatrixXi BECE = cat(1,BE,CE);
  147. MatrixXd S;
  148. // Sample each edge with 10 samples. (Choice of 10 doesn't seem to matter so
  149. // much, but could under some circumstances)
  150. sample_edges(C,BECE,10,S);
  151. // Vertices we'll constrain tet mesh to meet
  152. MatrixXd VS = cat(1,V,S);
  153. // Boundary faces
  154. MatrixXi BF;
  155. // Use tetgen to mesh the interior of surface, this assumes surface:
  156. // * has no holes
  157. // * has no non-manifold edges or vertices
  158. // * has consistent orientation
  159. // * has no self-intersections
  160. // * has no 0-volume pieces
  161. // Default settings pq100 tell tetgen to mesh interior of triangle mesh and
  162. // to produce a graded tet mesh
  163. cerr<<"tetgen begin()"<<endl;
  164. int status = tetrahedralize( VS,F,"pq100",VV,TT,FF);
  165. cerr<<"tetgen end()"<<endl;
  166. if(FF.rows() != F.rows())
  167. {
  168. // Issue a warning if the surface has changed
  169. cerr<<"mesh_with_skeleton: Warning: boundary faces != input faces"<<endl;
  170. }
  171. if(status != 0)
  172. {
  173. cerr<<
  174. "***************************************************************"<<endl<<
  175. "***************************************************************"<<endl<<
  176. "***************************************************************"<<endl<<
  177. "***************************************************************"<<endl<<
  178. "* mesh_with_skeleton: tetgen failed. Just meshing convex hull *"<<endl<<
  179. "***************************************************************"<<endl<<
  180. "***************************************************************"<<endl<<
  181. "***************************************************************"<<endl<<
  182. "***************************************************************"<<endl;
  183. // If meshing convex hull then use more regular mesh
  184. status = tetrahedralize(VS,F,"q1.414",VV,TT,FF);
  185. // I suppose this will fail if the skeleton is outside the mesh
  186. assert(FF.maxCoeff() < VV.rows());
  187. if(status != 0)
  188. {
  189. cerr<<"mesh_with_skeleton: tetgen failed again."<<endl;
  190. return false;
  191. }
  192. }
  193. // If you have medit installed then it's convenient to visualize the tet mesh
  194. // at this point
  195. if(WITH_MEDIT)
  196. {
  197. launch_medit(VV,TT,FF,false);
  198. }
  199. return true;
  200. }
  201. // Writes output files to /path/to/input/mesh-skeleton.dmat,
  202. // mesh-volume.dmat, mesh-volume.mesh if input mesh was
  203. // located at /path/to/input/mesh.obj and input skeleton was at
  204. // /other/path/to/input/skel.tgf
  205. //
  206. // Writes:
  207. //// mesh.dmat dense weights matrix corresponding to original input
  208. //// vertices V
  209. // mesh-volume.dmat dense weights matrix corresponding to all
  210. // vertices in tet mesh used for computation VV
  211. // mesh-volume.mesh Tet mesh used for computation
  212. //
  213. // Inputs:
  214. // mesh_filename path to {.obj|.off|.mesh} file
  215. // skel_filename path to skeleton {.bf|.tgf} file
  216. // V #V by 3 list of original mesh vertex positions
  217. // F #F by 3 list of original triangle indices
  218. // VV #VV by 3 list of tet-mesh vertex positions
  219. // TT #TT by 4 list of tetrahedra indices
  220. // FF #FF by 3 list of surface triangle indices
  221. // W #VV by #W weights matrix
  222. // Returns true on success
  223. bool save_output(
  224. const std::string mesh_filename,
  225. const std::string /*skel_filename*/,
  226. const Eigen::MatrixXd & V,
  227. const Eigen::MatrixXi & /*F*/,
  228. const Eigen::MatrixXd & VV,
  229. const Eigen::MatrixXi & TT,
  230. const Eigen::MatrixXi & FF,
  231. const Eigen::MatrixXd & W)
  232. {
  233. using namespace std;
  234. using namespace igl;
  235. using namespace Eigen;
  236. // build filename prefix out of input base names
  237. string prefix = "";
  238. {
  239. string dirname, basename, extension, filename;
  240. pathinfo(mesh_filename,dirname,basename,extension,filename);
  241. transform(extension.begin(), extension.end(), extension.begin(), ::tolower);
  242. prefix += dirname + "/" + filename;
  243. }
  244. //{
  245. // string dirname, basename, extension, filename;
  246. // pathinfo(skel_filename,dirname,basename,extension,filename);
  247. // transform(extension.begin(), extension.end(), extension.begin(), ::tolower);
  248. // prefix += "-" + filename;
  249. //}
  250. // Keep track if any fail
  251. bool success = true;
  252. //// Weights matrix for just V. Assumes V prefaces VV
  253. //MatrixXd WV = W.block(0,0,V.rows(),W.cols());
  254. //// write dmat
  255. //success &= writeDMAT(prefix + ".dmat",WV);
  256. // write volume weights dmat
  257. success &= writeDMAT(prefix + "-volume.dmat",W);
  258. // write volume mesh
  259. success &= writeMESH(prefix + "-volume.mesh",VV,TT,FF);
  260. //// write surface OBJ with pseudocolor
  261. return success;
  262. }
  263. int main(int argc, char * argv[])
  264. {
  265. using namespace std;
  266. using namespace Eigen;
  267. using namespace igl;
  268. if(argc<3)
  269. {
  270. cerr<<USAGE<<endl;
  271. return 1;
  272. }
  273. // #V by 3 list of mesh vertex positions
  274. MatrixXd V;
  275. // #F by 3 list of triangle indices
  276. MatrixXi F;
  277. // load mesh from .obj, .off or .mesh
  278. if(!load_mesh_from_file(argv[1],V,F))
  279. {
  280. return 1;
  281. }
  282. // "Skeleton" (handles) descriptors:
  283. // List of control and joint (bone endpoint) positions
  284. MatrixXd C;
  285. // List of point handles indexing C
  286. VectorXi P;
  287. // List of bone edges indexing C
  288. MatrixXi BE;
  289. // List of cage edges indexing *P*
  290. MatrixXi CE;
  291. // load skeleton (.tgf or .bf)
  292. if(!load_skeleton_from_file(argv[2],C,P,BE,CE))
  293. {
  294. return 1;
  295. }
  296. // Mesh with samples on skeleton
  297. // New vertices of tet mesh, V prefaces VV
  298. MatrixXd VV;
  299. // Tetrahedra
  300. MatrixXi TT;
  301. // New surface faces FF
  302. MatrixXi FF;
  303. if(!mesh_with_skeleton(V,F,C,P,BE,CE,VV,TT,FF))
  304. {
  305. return 1;
  306. }
  307. // Compute boundary conditions (aka fixed value constraints)
  308. // List of boundary indices (aka fixed value indices into VV)
  309. VectorXi b;
  310. // List of boundary conditions of each weight function
  311. MatrixXd bc;
  312. if(!boundary_conditions(VV,TT,C,P,BE,CE,b,bc))
  313. {
  314. return 1;
  315. }
  316. cout<<"b=["<<b<<"];"<<endl;
  317. cout<<"bc=["<<bc<<"];"<<endl;
  318. // compute BBW
  319. // Default bbw data and flags
  320. BBWData bbw_data;
  321. // Weights matrix
  322. MatrixXd W;
  323. if(!bbw(VV,TT,b,bc,bbw_data,W))
  324. {
  325. return 1;
  326. }
  327. // Save output
  328. save_output(argv[1],argv[2],V,F,VV,TT,FF,W);
  329. return 0;
  330. }