py_doc.cpp 65 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396
  1. const char *__doc_igl_active_set = R"igl_Qu8mg5v7(// Known Bugs: rows of [Aeq;Aieq] **must** be linearly independent. Should be
  2. // using QR decomposition otherwise:
  3. // http://www.okstate.edu/sas/v8/sashtml/ormp/chap5/sect32.htm
  4. //
  5. // ACTIVE_SET Minimize quadratic energy
  6. //
  7. // 0.5*Z'*A*Z + Z'*B + C with constraints
  8. //
  9. // that Z(known) = Y, optionally also subject to the constraints Aeq*Z = Beq,
  10. // and further optionally subject to the linear inequality constraints that
  11. // Aieq*Z <= Bieq and constant inequality constraints lx <= x <= ux
  12. //
  13. // Inputs:
  14. // A n by n matrix of quadratic coefficients
  15. // B n by 1 column of linear coefficients
  16. // known list of indices to known rows in Z
  17. // Y list of fixed values corresponding to known rows in Z
  18. // Aeq meq by n list of linear equality constraint coefficients
  19. // Beq meq by 1 list of linear equality constraint constant values
  20. // Aieq mieq by n list of linear inequality constraint coefficients
  21. // Bieq mieq by 1 list of linear inequality constraint constant values
  22. // lx n by 1 list of lower bounds [] implies -Inf
  23. // ux n by 1 list of upper bounds [] implies Inf
  24. // params struct of additional parameters (see below)
  25. // Z if not empty, is taken to be an n by 1 list of initial guess values
  26. // (see output)
  27. // Outputs:
  28. // Z n by 1 list of solution values
  29. // Returns true on success, false on error
  30. //
  31. // Benchmark: For a harmonic solve on a mesh with 325K facets, matlab 2.2
  32. // secs, igl/min_quad_with_fixed.h 7.1 secs
  33. //)igl_Qu8mg5v7";
  34. const char *__doc_igl_arap_precomputation = R"igl_Qu8mg5v7(// Compute necessary information to start using an ARAP deformation
  35. //
  36. // Inputs:
  37. // V #V by dim list of mesh positions
  38. // F #F by simplex-size list of triangle|tet indices into V
  39. // dim dimension being used at solve time. For deformation usually dim =
  40. // V.cols(), for surface parameterization V.cols() = 3 and dim = 2
  41. // b #b list of "boundary" fixed vertex indices into V
  42. // Outputs:
  43. // data struct containing necessary precomputation)igl_Qu8mg5v7";
  44. const char *__doc_igl_arap_solve = R"igl_Qu8mg5v7(// Inputs:
  45. // bc #b by dim list of boundary conditions
  46. // data struct containing necessary precomputation and parameters
  47. // U #V by dim initial guess)igl_Qu8mg5v7";
  48. const char *__doc_igl_avg_edge_length = R"igl_Qu8mg5v7(// Compute the average edge length for the given triangle mesh
  49. // Templates:
  50. // DerivedV derived from vertex positions matrix type: i.e. MatrixXd
  51. // DerivedF derived from face indices matrix type: i.e. MatrixXi
  52. // DerivedL derived from edge lengths matrix type: i.e. MatrixXd
  53. // Inputs:
  54. // V eigen matrix #V by 3
  55. // F #F by simplex-size list of mesh faces (must be simplex)
  56. // Outputs:
  57. // l average edge length
  58. //
  59. // See also: adjacency_matrix)igl_Qu8mg5v7";
  60. const char *__doc_igl_barycenter = R"igl_Qu8mg5v7(// Computes the barycenter of every simplex
  61. //
  62. // Inputs:
  63. // V #V x dim matrix of vertex coordinates
  64. // F #F x simplex_size matrix of indices of simplex corners into V
  65. // Output:
  66. // BC #F x dim matrix of 3d vertices
  67. //)igl_Qu8mg5v7";
  68. const char *__doc_igl_barycentric_coordinates = R"igl_Qu8mg5v7(// Compute barycentric coordinates in a tet
  69. //
  70. // Inputs:
  71. // P #P by 3 Query points in 3d
  72. // A #P by 3 Tet corners in 3d
  73. // B #P by 3 Tet corners in 3d
  74. // C #P by 3 Tet corners in 3d
  75. // D #P by 3 Tet corners in 3d
  76. // Outputs:
  77. // L #P by 4 list of barycentric coordinates
  78. // )igl_Qu8mg5v7";
  79. const char *__doc_igl_barycentric_to_global = R"igl_Qu8mg5v7(// Converts barycentric coordinates in the embree form to 3D coordinates
  80. // Embree stores barycentric coordinates as triples: fid, bc1, bc2
  81. // fid is the id of a face, bc1 is the displacement of the point wrt the
  82. // first vertex v0 and the edge v1-v0. Similarly, bc2 is the displacement
  83. // wrt v2-v0.
  84. //
  85. // Input:
  86. // V: #Vx3 Vertices of the mesh
  87. // F: #Fxe Faces of the mesh
  88. // bc: #Xx3 Barycentric coordinates, one row per point
  89. //
  90. // Output:
  91. // #X: #Xx3 3D coordinates of all points in bc
  92. // )igl_Qu8mg5v7";
  93. const char *__doc_igl_bbw_bbw = R"igl_Qu8mg5v7(// Compute Bounded Biharmonic Weights on a given domain (V,Ele) with a given
  94. // set of boundary conditions
  95. //
  96. // Templates
  97. // DerivedV derived type of eigen matrix for V (e.g. MatrixXd)
  98. // DerivedF derived type of eigen matrix for F (e.g. MatrixXi)
  99. // Derivedb derived type of eigen matrix for b (e.g. VectorXi)
  100. // Derivedbc derived type of eigen matrix for bc (e.g. MatrixXd)
  101. // DerivedW derived type of eigen matrix for W (e.g. MatrixXd)
  102. // Inputs:
  103. // V #V by dim vertex positions
  104. // Ele #Elements by simplex-size list of element indices
  105. // b #b boundary indices into V
  106. // bc #b by #W list of boundary values
  107. // data object containing options, intial guess --> solution and results
  108. // Outputs:
  109. // W #V by #W list of *unnormalized* weights to normalize use
  110. // igl::normalize_row_sums(W,W);
  111. // Returns true on success, false on failure)igl_Qu8mg5v7";
  112. const char *__doc_igl_boundary_conditions = R"igl_Qu8mg5v7(// Compute boundary conditions for automatic weights computation. This
  113. // function expects that the given mesh (V,Ele) has sufficient samples
  114. // (vertices) exactly at point handle locations and exactly along bone and
  115. // cage edges.
  116. //
  117. // Inputs:
  118. // V #V by dim list of domain vertices
  119. // Ele #Ele by simplex-size list of simplex indices
  120. // C #C by dim list of handle positions
  121. // P #P by 1 list of point handle indices into C
  122. // BE #BE by 2 list of bone edge indices into C
  123. // CE #CE by 2 list of cage edge indices into *P*
  124. // Outputs:
  125. // b #b list of boundary indices (indices into V of vertices which have
  126. // known, fixed values)
  127. // bc #b by #weights list of known/fixed values for boundary vertices
  128. // (notice the #b != #weights in general because #b will include all the
  129. // intermediary samples along each bone, etc.. The ordering of the
  130. // weights corresponds to [P;BE]
  131. // Returns true if boundary conditions make sense)igl_Qu8mg5v7";
  132. const char *__doc_igl_boundary_facets = R"igl_Qu8mg5v7(// BOUNDARY_FACETS Determine boundary faces (edges) of tetrahedra (triangles)
  133. // stored in T (analogous to qptoolbox's `outline` and `boundary_faces`).
  134. //
  135. // Templates:
  136. // IntegerT integer-value: e.g. int
  137. // IntegerF integer-value: e.g. int
  138. // Input:
  139. // T tetrahedron (triangle) index list, m by 4 (3), where m is the number of tetrahedra
  140. // Output:
  141. // F list of boundary faces, n by 3 (2), where n is the number of boundary faces
  142. //
  143. //)igl_Qu8mg5v7";
  144. const char *__doc_igl_boundary_loop = R"igl_Qu8mg5v7(// Compute list of ordered boundary loops for a manifold mesh.
  145. //
  146. // Templates:
  147. // Index index type
  148. // Inputs:
  149. // F #V by dim list of mesh faces
  150. // Outputs:
  151. // L list of loops where L[i] = ordered list of boundary vertices in loop i
  152. //)igl_Qu8mg5v7";
  153. const char *__doc_igl_cat = R"igl_Qu8mg5v7(// Perform concatenation of a two matrices along a single dimension
  154. // If dim == 1, then C = [A;B]. If dim == 2 then C = [A B]
  155. //
  156. // Template:
  157. // Scalar scalar data type for sparse matrices like double or int
  158. // Mat matrix type for all matrices (e.g. MatrixXd, SparseMatrix)
  159. // MatC matrix type for ouput matrix (e.g. MatrixXd) needs to support
  160. // resize
  161. // Inputs:
  162. // A first input matrix
  163. // B second input matrix
  164. // dim dimension along which to concatenate, 0 or 1
  165. // Outputs:
  166. // C output matrix
  167. // )igl_Qu8mg5v7";
  168. const char *__doc_igl_collapse_edge = R"igl_Qu8mg5v7(See collapse_edge for the documentation.)igl_Qu8mg5v7";
  169. const char *__doc_igl_colon = R"igl_Qu8mg5v7(// Colon operator like matlab's colon operator. Enumerats values between low
  170. // and hi with step step.
  171. // Templates:
  172. // L should be a eigen matrix primitive type like int or double
  173. // S should be a eigen matrix primitive type like int or double
  174. // H should be a eigen matrix primitive type like int or double
  175. // T should be a eigen matrix primitive type like int or double
  176. // Inputs:
  177. // low starting value if step is valid then this is *always* the first
  178. // element of I
  179. // step step difference between sequential elements returned in I,
  180. // remember this will be cast to template T at compile time. If low<hi
  181. // then step must be positive. If low>hi then step must be negative.
  182. // Otherwise I will be set to empty.
  183. // hi ending value, if (hi-low)%step is zero then this will be the last
  184. // element in I. If step is positive there will be no elements greater
  185. // than hi, vice versa if hi<low
  186. // Output:
  187. // I list of values from low to hi with step size step)igl_Qu8mg5v7";
  188. const char *__doc_igl_column_to_quats = R"igl_Qu8mg5v7(// "Columnize" a list of quaternions (q1x,q1y,q1z,q1w,q2x,q2y,q2z,q2w,...)
  189. //
  190. // Inputs:
  191. // Q n*4-long list of coefficients
  192. // Outputs:
  193. // vQ n-long list of quaternions
  194. // Returns false if n%4!=0)igl_Qu8mg5v7";
  195. const char *__doc_igl_comb_cross_field = R"igl_Qu8mg5v7(// Inputs:
  196. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  197. // F #F by 4 eigen Matrix of face (quad) indices
  198. // PD1in #F by 3 eigen Matrix of the first per face cross field vector
  199. // PD2in #F by 3 eigen Matrix of the second per face cross field vector
  200. // Output:
  201. // PD1out #F by 3 eigen Matrix of the first combed cross field vector
  202. // PD2out #F by 3 eigen Matrix of the second combed cross field vector
  203. //)igl_Qu8mg5v7";
  204. const char *__doc_igl_comb_frame_field = R"igl_Qu8mg5v7(// Inputs:
  205. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  206. // F #F by 4 eigen Matrix of face (quad) indices
  207. // PD1 #F by 3 eigen Matrix of the first per face cross field vector
  208. // PD2 #F by 3 eigen Matrix of the second per face cross field vector
  209. // BIS1_combed #F by 3 eigen Matrix of the first combed bisector field vector
  210. // BIS2_combed #F by 3 eigen Matrix of the second combed bisector field vector
  211. // Output:
  212. // PD1_combed #F by 3 eigen Matrix of the first combed cross field vector
  213. // PD2_combed #F by 3 eigen Matrix of the second combed cross field vector
  214. //)igl_Qu8mg5v7";
  215. const char *__doc_igl_compute_frame_field_bisectors = R"igl_Qu8mg5v7(// Compute bisectors of a frame field defined on mesh faces
  216. // Inputs:
  217. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  218. // F #F by 3 eigen Matrix of face (triangle) indices
  219. // B1 #F by 3 eigen Matrix of face (triangle) base vector 1
  220. // B2 #F by 3 eigen Matrix of face (triangle) base vector 2
  221. // PD1 #F by 3 eigen Matrix of the first per face frame field vector
  222. // PD2 #F by 3 eigen Matrix of the second per face frame field vector
  223. // Output:
  224. // BIS1 #F by 3 eigen Matrix of the first per face frame field bisector
  225. // BIS2 #F by 3 eigen Matrix of the second per face frame field bisector
  226. //)igl_Qu8mg5v7";
  227. const char *__doc_igl_copyleft_cgal_mesh_boolean = R"igl_Qu8mg5v7(// MESH_BOOLEAN Compute boolean csg operations on "solid", consistently
  228. // oriented meshes.
  229. //
  230. // Inputs:
  231. // VA #VA by 3 list of vertex positions of first mesh
  232. // FA #FA by 3 list of triangle indices into VA
  233. // VB #VB by 3 list of vertex positions of second mesh
  234. // FB #FB by 3 list of triangle indices into VB
  235. // type type of boolean operation
  236. // Outputs:
  237. // VC #VC by 3 list of vertex positions of boolean result mesh
  238. // FC #FC by 3 list of triangle indices into VC
  239. // J #FC list of indices into [FA;FA.rows()+FB] revealing "birth" facet
  240. // Returns true if inputs induce a piecewise constant winding number
  241. // field and type is valid
  242. //
  243. // See also: mesh_boolean_cork, intersect_other,
  244. // remesh_self_intersections)igl_Qu8mg5v7";
  245. const char *__doc_igl_copyleft_cgal_remesh_self_intersections = R"igl_Qu8mg5v7(// Given a triangle mesh (V,F) compute a new mesh (VV,FF) which is the same
  246. // as (V,F) except that any self-intersecting triangles in (V,F) have been
  247. // subdivided (new vertices and face created) so that the self-intersection
  248. // contour lies exactly on edges in (VV,FF). New vertices will appear in
  249. // original faces or on original edges. New vertices on edges are "merged"
  250. // only across original faces sharing that edge. This means that if the input
  251. // triangle mesh is a closed manifold the output will be too.
  252. //
  253. // Inputs:
  254. // V #V by 3 list of vertex positions
  255. // F #F by 3 list of triangle indices into V
  256. // params struct of optional parameters
  257. // Outputs:
  258. // VV #VV by 3 list of vertex positions
  259. // FF #FF by 3 list of triangle indices into VV
  260. // IF #intersecting face pairs by 2 list of intersecting face pairs,
  261. // indexing F
  262. // J #FF list of indices into F denoting birth triangle
  263. // IM #VV list of indices into VV of unique vertices.
  264. //
  265. // Known bugs: If an existing edge in (V,F) lies exactly on another face then
  266. // any resulting additional vertices along that edge may not get properly
  267. // connected so that the output mesh has the same global topology. This is
  268. // because
  269. //
  270. // Example:
  271. // // resolve intersections
  272. // igl::copyleft::cgal::remesh_self_intersections(V,F,params,VV,FF,IF,J,IM);
  273. // // _apply_ duplicate vertex mapping IM to FF
  274. // for_each(FF.data(),FF.data()+FF.size(),[&IM](int & a){a=IM(a);});
  275. // // remove any vertices now unreferenced after duplicate mapping.
  276. // igl::remove_unreferenced(VV,FF,SV,SF,UIM);
  277. // // Now (SV,SF) is ready to extract outer hull
  278. // igl::copyleft::cgal::outer_hull(SV,SF,G,J,flip);
  279. //)igl_Qu8mg5v7";
  280. const char *__doc_igl_copyleft_comiso_miq = R"igl_Qu8mg5v7(// Inputs:
  281. // V #V by 3 list of mesh vertex 3D positions
  282. // F #F by 3 list of faces indices in V
  283. // PD1 #V by 3 first line of the Jacobian per triangle
  284. // PD2 #V by 3 second line of the Jacobian per triangle
  285. // (optional, if empty it will be a vector in the tangent plane orthogonal to PD1)
  286. // scale global scaling for the gradient (controls the quads resolution)
  287. // stiffness weight for the stiffness iterations
  288. // direct_round greedily round all integer variables at once (greatly improves optimization speed but lowers quality)
  289. // iter stiffness iterations (0 = no stiffness)
  290. // local_iter number of local iterations for the integer rounding
  291. // do_round enables the integer rounding (disabling it could be useful for debugging)
  292. // round_vertices id of additional vertices that should be snapped to integer coordinates
  293. // hard_features #H by 2 list of pairs of vertices that belongs to edges that should be snapped to integer coordinates
  294. //
  295. // Output:
  296. // UV #UV by 2 list of vertices in 2D
  297. // FUV #FUV by 3 list of face indices in UV
  298. //
  299. // TODO: rename the parameters name in the cpp consistenly
  300. // improve the handling of hard_features, right now it might fail in difficult cases)igl_Qu8mg5v7";
  301. const char *__doc_igl_copyleft_comiso_nrosy = R"igl_Qu8mg5v7(// Generate a N-RoSy field from a sparse set of constraints
  302. //
  303. // Inputs:
  304. // V #V by 3 list of mesh vertex coordinates
  305. // F #F by 3 list of mesh faces (must be triangles)
  306. // b #B by 1 list of constrained face indices
  307. // bc #B by 3 list of representative vectors for the constrained
  308. // faces
  309. // b_soft #S by 1 b for soft constraints
  310. // w_soft #S by 1 weight for the soft constraints (0-1)
  311. // bc_soft #S by 3 bc for soft constraints
  312. // N the degree of the N-RoSy vector field
  313. // soft the strenght of the soft contraints w.r.t. smoothness
  314. // (0 -> smoothness only, 1->constraints only)
  315. // Outputs:
  316. // R #F by 3 the representative vectors of the interpolated field
  317. // S #V by 1 the singularity index for each vertex (0 = regular))igl_Qu8mg5v7";
  318. const char *__doc_igl_copyleft_marching_cubes = R"igl_Qu8mg5v7(// marching_cubes( values, points, x_res, y_res, z_res, vertices, faces )
  319. //
  320. // performs marching cubes reconstruction on the grid defined by values, and
  321. // points, and generates vertices and faces
  322. //
  323. // Input:
  324. // values #number_of_grid_points x 1 array -- the scalar values of an
  325. // implicit function defined on the grid points (<0 in the inside of the
  326. // surface, 0 on the border, >0 outside)
  327. // points #number_of_grid_points x 3 array -- 3-D positions of the grid
  328. // points, ordered in x,y,z order:
  329. // points[index] = the point at (x,y,z) where :
  330. // x = (index % (xres -1),
  331. // y = (index / (xres-1)) %(yres-1),
  332. // z = index / (xres -1) / (yres -1) ).
  333. // where x,y,z index x, y, z dimensions
  334. // i.e. index = x + y*xres + z*xres*yres
  335. // xres resolutions of the grid in x dimension
  336. // yres resolutions of the grid in y dimension
  337. // zres resolutions of the grid in z dimension
  338. // Output:
  339. // vertices #V by 3 list of mesh vertex positions
  340. // faces #F by 3 list of mesh triangle indices
  341. //)igl_Qu8mg5v7";
  342. const char *__doc_igl_copyleft_swept_volume = R"igl_Qu8mg5v7(// Compute the surface of the swept volume of a solid object with surface
  343. // (V,F) mesh under going rigid motion.
  344. //
  345. // Inputs:
  346. // V #V by 3 list of mesh positions in reference pose
  347. // F #F by 3 list of mesh indices into V
  348. // transform function handle so that transform(t) returns the rigid
  349. // transformation at time t∈[0,1]
  350. // steps number of time steps: steps=3 --> t∈{0,0.5,1}
  351. // grid_res number of grid cells on the longest side containing the
  352. // motion (isolevel+1 cells will also be added on each side as padding)
  353. // isolevel distance level to be contoured as swept volume
  354. // Outputs:
  355. // SV #SV by 3 list of mesh positions of the swept surface
  356. // SF #SF by 3 list of mesh faces into SV)igl_Qu8mg5v7";
  357. const char *__doc_igl_copyleft_tetgen_tetrahedralize = R"igl_Qu8mg5v7(// Mesh the interior of a surface mesh (V,F) using tetgen
  358. //
  359. // Inputs:
  360. // V #V by 3 vertex position list
  361. // F #F list of polygon face indices into V (0-indexed)
  362. // switches string of tetgen options (See tetgen documentation) e.g.
  363. // "pq1.414a0.01" tries to mesh the interior of a given surface with
  364. // quality and area constraints
  365. // "" will mesh the convex hull constrained to pass through V (ignores F)
  366. // Outputs:
  367. // TV #V by 3 vertex position list
  368. // TT #T by 4 list of tet face indices
  369. // TF #F by 3 list of triangle face indices
  370. // Returns status:
  371. // 0 success
  372. // 1 tetgen threw exception
  373. // 2 tetgen did not crash but could not create any tets (probably there are
  374. // holes, duplicate faces etc.)
  375. // -1 other error)igl_Qu8mg5v7";
  376. const char *__doc_igl_cotmatrix = R"igl_Qu8mg5v7(// Constructs the cotangent stiffness matrix (discrete laplacian) for a given
  377. // mesh (V,F).
  378. //
  379. // Templates:
  380. // DerivedV derived type of eigen matrix for V (e.g. derived from
  381. // MatrixXd)
  382. // DerivedF derived type of eigen matrix for F (e.g. derived from
  383. // MatrixXi)
  384. // Scalar scalar type for eigen sparse matrix (e.g. double)
  385. // Inputs:
  386. // V #V by dim list of mesh vertex positions
  387. // F #F by simplex_size list of mesh faces (must be triangles)
  388. // Outputs:
  389. // L #V by #V cotangent matrix, each row i corresponding to V(i,:)
  390. //
  391. // See also: adjacency_matrix
  392. //
  393. // Note: This Laplacian uses the convention that diagonal entries are
  394. // **minus** the sum of off-diagonal entries. The diagonal entries are
  395. // therefore in general negative and the matrix is **negative** semi-definite
  396. // (immediately, -L is **positive** semi-definite)
  397. //
  398. // Known bugs: off by 1e-16 on regular grid. I think its a problem of
  399. // arithmetic order in cotmatrix_entries.h: C(i,e) = (arithmetic)/dblA/4)igl_Qu8mg5v7";
  400. const char *__doc_igl_covariance_scatter_matrix = R"igl_Qu8mg5v7(// Construct the covariance scatter matrix for a given arap energy
  401. // Inputs:
  402. // V #V by Vdim list of initial domain positions
  403. // F #F by 3 list of triangle indices into V
  404. // energy ARAPEnergyType enum value defining which energy is being used.
  405. // See ARAPEnergyType.h for valid options and explanations.
  406. // Outputs:
  407. // CSM dim*#V/#F by dim*#V sparse matrix containing special laplacians along
  408. // the diagonal so that when multiplied by V gives covariance matrix
  409. // elements, can be used to speed up covariance matrix computation)igl_Qu8mg5v7";
  410. const char *__doc_igl_cross_field_missmatch = R"igl_Qu8mg5v7(// Inputs:
  411. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  412. // F #F by 3 eigen Matrix of face (quad) indices
  413. // PD1 #F by 3 eigen Matrix of the first per face cross field vector
  414. // PD2 #F by 3 eigen Matrix of the second per face cross field vector
  415. // isCombed boolean, specifying whether the field is combed (i.e. matching has been precomputed.
  416. // If not, the field is combed first.
  417. // Output:
  418. // Handle_MMatch #F by 3 eigen Matrix containing the integer missmatch of the cross field
  419. // across all face edges
  420. //)igl_Qu8mg5v7";
  421. const char *__doc_igl_cut_mesh_from_singularities = R"igl_Qu8mg5v7(// Given a mesh (V,F) and the integer mismatch of a cross field per edge
  422. // (MMatch), finds the cut_graph connecting the singularities (seams) and the
  423. // degree of the singularities singularity_index
  424. //
  425. // Input:
  426. // V #V by 3 list of mesh vertex positions
  427. // F #F by 3 list of faces
  428. // MMatch #F by 3 list of per corner integer mismatch
  429. // Outputs:
  430. // seams #F by 3 list of per corner booleans that denotes if an edge is a
  431. // seam or not
  432. //)igl_Qu8mg5v7";
  433. const char *__doc_igl_deform_skeleton = R"igl_Qu8mg5v7(// Deform a skeleton.
  434. //
  435. // Inputs:
  436. // C #C by 3 list of joint positions
  437. // BE #BE by 2 list of bone edge indices
  438. // vA #BE list of bone transformations
  439. // Outputs
  440. // CT #BE*2 by 3 list of deformed joint positions
  441. // BET #BE by 2 list of bone edge indices (maintains order)
  442. //)igl_Qu8mg5v7";
  443. const char *__doc_igl_directed_edge_orientations = R"igl_Qu8mg5v7(// Determine rotations that take each edge from the x-axis to its given rest
  444. // orientation.
  445. //
  446. // Inputs:
  447. // C #C by 3 list of edge vertex positions
  448. // E #E by 2 list of directed edges
  449. // Outputs:
  450. // Q #E list of quaternions
  451. //)igl_Qu8mg5v7";
  452. const char *__doc_igl_directed_edge_parents = R"igl_Qu8mg5v7(// Recover "parents" (preceeding edges) in a tree given just directed edges.
  453. //
  454. // Inputs:
  455. // E #E by 2 list of directed edges
  456. // Outputs:
  457. // P #E list of parent indices into E (-1) means root
  458. //)igl_Qu8mg5v7";
  459. const char *__doc_igl_doublearea = R"igl_Qu8mg5v7(// DOUBLEAREA computes twice the area for each input triangle[quad]
  460. //
  461. // Templates:
  462. // DerivedV derived type of eigen matrix for V (e.g. derived from
  463. // MatrixXd)
  464. // DerivedF derived type of eigen matrix for F (e.g. derived from
  465. // MatrixXi)
  466. // DeriveddblA derived type of eigen matrix for dblA (e.g. derived from
  467. // MatrixXd)
  468. // Inputs:
  469. // V #V by dim list of mesh vertex positions
  470. // F #F by simplex_size list of mesh faces (must be triangles or quads)
  471. // Outputs:
  472. // dblA #F list of triangle[quad] double areas (SIGNED only for 2D input)
  473. //
  474. // Known bug: For dim==3 complexity is O(#V + #F)!! Not just O(#F). This is a big deal
  475. // if you have 1million unreferenced vertices and 1 face)igl_Qu8mg5v7";
  476. const char *__doc_igl_doublearea_single = R"igl_Qu8mg5v7(// Single triangle in 2D!
  477. //
  478. // This should handle streams of corners not just single corners)igl_Qu8mg5v7";
  479. const char *__doc_igl_doublearea_quad = R"igl_Qu8mg5v7(// DOUBLEAREA_QUAD computes twice the area for each input quadrilateral
  480. //
  481. // Inputs:
  482. // V #V by dim list of mesh vertex positions
  483. // F #F by simplex_size list of mesh faces (must be quadrilaterals)
  484. // Outputs:
  485. // dblA #F list of quadrilateral double areas
  486. //)igl_Qu8mg5v7";
  487. const char *__doc_igl_dqs = R"igl_Qu8mg5v7(// Dual quaternion skinning
  488. //
  489. // Inputs:
  490. // V #V by 3 list of rest positions
  491. // W #W by #C list of weights
  492. // vQ #C list of rotation quaternions
  493. // vT #C list of translation vectors
  494. // Outputs:
  495. // U #V by 3 list of new positions)igl_Qu8mg5v7";
  496. const char *__doc_igl_edge_lengths = R"igl_Qu8mg5v7(// Constructs a list of lengths of edges opposite each index in a face
  497. // (triangle/tet) list
  498. //
  499. // Templates:
  500. // DerivedV derived from vertex positions matrix type: i.e. MatrixXd
  501. // DerivedF derived from face indices matrix type: i.e. MatrixXi
  502. // DerivedL derived from edge lengths matrix type: i.e. MatrixXd
  503. // Inputs:
  504. // V eigen matrix #V by 3
  505. // F #F by 2 list of mesh edges
  506. // or
  507. // F #F by 3 list of mesh faces (must be triangles)
  508. // or
  509. // T #T by 4 list of mesh elements (must be tets)
  510. // Outputs:
  511. // L #F by {1|3|6} list of edge lengths
  512. // for edges, column of lengths
  513. // for triangles, columns correspond to edges [1,2],[2,0],[0,1]
  514. // for tets, columns correspond to edges
  515. // [3 0],[3 1],[3 2],[1 2],[2 0],[0 1]
  516. //)igl_Qu8mg5v7";
  517. const char *__doc_igl_edge_topology = R"igl_Qu8mg5v7(// Initialize Edges and their topological relations
  518. //
  519. // Output:
  520. // EV : #Ex2, Stores the edge description as pair of indices to vertices
  521. // FE : #Fx3, Stores the Triangle-Edge relation
  522. // EF : #Ex2: Stores the Edge-Triangle relation
  523. //
  524. // TODO: This seems to be a duplicate of edge_flaps.h)igl_Qu8mg5v7";
  525. const char *__doc_igl_eigs = R"igl_Qu8mg5v7(See eigs for the documentation.)igl_Qu8mg5v7";
  526. const char *__doc_igl_embree_ambient_occlusion = R"igl_Qu8mg5v7(// Compute ambient occlusion per given point
  527. //
  528. // Inputs:
  529. // ei EmbreeIntersector containing (V,F)
  530. // P #P by 3 list of origin points
  531. // N #P by 3 list of origin normals
  532. // Outputs:
  533. // S #P list of ambient occlusion values between 1 (fully occluded) and
  534. // 0 (not occluded)
  535. //)igl_Qu8mg5v7";
  536. const char *__doc_igl_embree_reorient_facets_raycast = R"igl_Qu8mg5v7(// Orient each component (identified by C) of a mesh (V,F) using ambient
  537. // occlusion such that the front side is less occluded than back side, as
  538. // described in "A Simple Method for Correcting Facet Orientations in
  539. // Polygon Meshes Based on Ray Casting" [Takayama et al. 2014].
  540. //
  541. // Inputs:
  542. // V #V by 3 list of vertex positions
  543. // F #F by 3 list of triangle indices
  544. // rays_total Total number of rays that will be shot
  545. // rays_minimum Minimum number of rays that each patch should receive
  546. // facet_wise Decision made for each face independently, no use of patches
  547. // (i.e., each face is treated as a patch)
  548. // use_parity Use parity mode
  549. // is_verbose Verbose output to cout
  550. // Outputs:
  551. // I #F list of whether face has been flipped
  552. // C #F list of patch ID (output of bfs_orient > manifold patches))igl_Qu8mg5v7";
  553. const char *__doc_igl_embree_line_mesh_intersection = R"igl_Qu8mg5v7(// Project the point cloud V_source onto the triangle mesh
  554. // V_target,F_target.
  555. // A ray is casted for every vertex in the direction specified by
  556. // N_source and its opposite.
  557. //
  558. // Input:
  559. // V_source: #Vx3 Vertices of the source mesh
  560. // N_source: #Vx3 Normals of the point cloud
  561. // V_target: #V2x3 Vertices of the target mesh
  562. // F_target: #F2x3 Faces of the target mesh
  563. //
  564. // Output:
  565. // #Vx3 matrix of baricentric coordinate. Each row corresponds to
  566. // a vertex of the projected mesh and it has the following format:
  567. // id b1 b2. id is the id of a face of the source mesh. b1 and b2 are
  568. // the barycentric coordinates wrt the first two edges of the triangle
  569. // To convert to standard global coordinates, see barycentric_to_global.h)igl_Qu8mg5v7";
  570. const char *__doc_igl_find_cross_field_singularities = R"igl_Qu8mg5v7(// Inputs:
  571. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  572. // F #F by 3 eigen Matrix of face (quad) indices
  573. // Handle_MMatch #F by 3 eigen Matrix containing the integer missmatch of the cross field
  574. // across all face edges
  575. // Output:
  576. // isSingularity #V by 1 boolean eigen Vector indicating the presence of a singularity on a vertex
  577. // singularityIndex #V by 1 integer eigen Vector containing the singularity indices
  578. //)igl_Qu8mg5v7";
  579. const char *__doc_igl_fit_rotations = R"igl_Qu8mg5v7(// Known issues: This seems to be implemented in Eigen/Geometry:
  580. // Eigen::umeyama
  581. //
  582. // FIT_ROTATIONS Given an input mesh and new positions find rotations for
  583. // every covariance matrix in a stack of covariance matrices
  584. //
  585. // Inputs:
  586. // S nr*dim by dim stack of covariance matrices
  587. // single_precision whether to use single precision (faster)
  588. // Outputs:
  589. // R dim by dim * nr list of rotations
  590. //)igl_Qu8mg5v7";
  591. const char *__doc_igl_fit_rotations_planar = R"igl_Qu8mg5v7(// FIT_ROTATIONS Given an input mesh and new positions find 2D rotations for
  592. // every vertex that best maps its one ring to the new one ring
  593. //
  594. // Inputs:
  595. // S nr*dim by dim stack of covariance matrices, third column and every
  596. // third row will be ignored
  597. // Outputs:
  598. // R dim by dim * nr list of rotations, third row and third column of each
  599. // rotation will just be identity
  600. //)igl_Qu8mg5v7";
  601. const char *__doc_igl_fit_rotations_SSE = R"igl_Qu8mg5v7(See fit_rotations_SSE for the documentation.)igl_Qu8mg5v7";
  602. const char *__doc_igl_floor = R"igl_Qu8mg5v7(// Floor a given matrix to nearest integers
  603. //
  604. // Inputs:
  605. // X m by n matrix of scalars
  606. // Outputs:
  607. // Y m by n matrix of floored integers)igl_Qu8mg5v7";
  608. const char *__doc_igl_forward_kinematics = R"igl_Qu8mg5v7(// Given a skeleton and a set of relative bone rotations compute absolute
  609. // rigid transformations for each bone.
  610. //
  611. // Inputs:
  612. // C #C by dim list of joint positions
  613. // BE #BE by 2 list of bone edge indices
  614. // P #BE list of parent indices into BE
  615. // dQ #BE list of relative rotations
  616. // dT #BE list of relative translations
  617. // Outputs:
  618. // vQ #BE list of absolute rotations
  619. // vT #BE list of absolute translations)igl_Qu8mg5v7";
  620. const char *__doc_igl_gaussian_curvature = R"igl_Qu8mg5v7(// Compute discrete local integral gaussian curvature (angle deficit, without
  621. // averaging by local area).
  622. //
  623. // Inputs:
  624. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  625. // F #F by 3 eigen Matrix of face (triangle) indices
  626. // Output:
  627. // K #V by 1 eigen Matrix of discrete gaussian curvature values
  628. //)igl_Qu8mg5v7";
  629. const char *__doc_igl_get_seconds = R"igl_Qu8mg5v7(// Return the current time in seconds since program start
  630. //
  631. // Example:
  632. // const auto & tictoc = []()
  633. // {
  634. // static double t_start = igl::get_seconds();
  635. // double diff = igl::get_seconds()-t_start;
  636. // t_start += diff;
  637. // return diff;
  638. // };
  639. // tictoc();
  640. // ... // part 1
  641. // cout<<"part 1: "<<tictoc()<<endl;
  642. // ... // part 2
  643. // cout<<"part 2: "<<tictoc()<<endl;
  644. // ... // etc)igl_Qu8mg5v7";
  645. const char *__doc_igl_grad = R"igl_Qu8mg5v7(// Gradient of a scalar function defined on piecewise linear elements (mesh)
  646. // is constant on each triangle i,j,k:
  647. // grad(Xijk) = (Xj-Xi) * (Vi - Vk)^R90 / 2A + (Xk-Xi) * (Vj - Vi)^R90 / 2A
  648. // where Xi is the scalar value at vertex i, Vi is the 3D position of vertex
  649. // i, and A is the area of triangle (i,j,k). ^R90 represent a rotation of
  650. // 90 degrees
  651. //)igl_Qu8mg5v7";
  652. const char *__doc_igl_harmonic = R"igl_Qu8mg5v7(// Compute k-harmonic weight functions "coordinates".
  653. //
  654. //
  655. // Inputs:
  656. // V #V by dim vertex positions
  657. // F #F by simplex-size list of element indices
  658. // b #b boundary indices into V
  659. // bc #b by #W list of boundary values
  660. // k power of harmonic operation (1: harmonic, 2: biharmonic, etc)
  661. // Outputs:
  662. // W #V by #W list of weights
  663. //)igl_Qu8mg5v7";
  664. const char *__doc_igl_hsv_to_rgb = R"igl_Qu8mg5v7(// Convert RGB to HSV
  665. //
  666. // Inputs:
  667. // h hue value (degrees: [0,360])
  668. // s saturation value ([0,1])
  669. // v value value ([0,1])
  670. // Outputs:
  671. // r red value ([0,1])
  672. // g green value ([0,1])
  673. // b blue value ([0,1]))igl_Qu8mg5v7";
  674. const char *__doc_igl_internal_angles = R"igl_Qu8mg5v7(// Compute internal angles for a triangle mesh
  675. //
  676. // Inputs:
  677. // V #V by dim eigen Matrix of mesh vertex nD positions
  678. // F #F by poly-size eigen Matrix of face (triangle) indices
  679. // Output:
  680. // K #F by poly-size eigen Matrix of internal angles
  681. // for triangles, columns correspond to edges [1,2],[2,0],[0,1]
  682. //
  683. // Known Issues:
  684. // if poly-size ≠ 3 then dim must equal 3.)igl_Qu8mg5v7";
  685. const char *__doc_igl_invert_diag = R"igl_Qu8mg5v7(// Templates:
  686. // T should be a eigen sparse matrix primitive type like int or double
  687. // Inputs:
  688. // X an m by n sparse matrix
  689. // Outputs:
  690. // Y an m by n sparse matrix)igl_Qu8mg5v7";
  691. const char *__doc_igl_is_irregular_vertex = R"igl_Qu8mg5v7(// Determine if a vertex is irregular, i.e. it has more than 6 (triangles)
  692. // or 4 (quads) incident edges. Vertices on the boundary are ignored.
  693. //
  694. // Inputs:
  695. // V #V by dim list of vertex positions
  696. // F #F by 3[4] list of triangle[quads] indices
  697. // Returns #V vector of bools revealing whether vertices are singular
  698. //)igl_Qu8mg5v7";
  699. const char *__doc_igl_jet = R"igl_Qu8mg5v7(// JET like MATLAB's jet
  700. //
  701. // Inputs:
  702. // m number of colors
  703. // Outputs:
  704. // J m by list of RGB colors between 0 and 1
  705. //
  706. //#ifndef IGL_NO_EIGEN
  707. // void jet(const int m, Eigen::MatrixXd & J);
  708. //#endif
  709. // Wrapper for directly computing [r,g,b] values for a given factor f between
  710. // 0 and 1
  711. //
  712. // Inputs:
  713. // f factor determining color value as if 0 was min and 1 was max
  714. // Outputs:
  715. // r red value
  716. // g green value
  717. // b blue value)igl_Qu8mg5v7";
  718. const char *__doc_igl_lbs_matrix = R"igl_Qu8mg5v7(// LBS_MATRIX Linear blend skinning can be expressed by V' = M * T where V' is
  719. // a #V by dim matrix of deformed vertex positions (one vertex per row), M is a
  720. // #V by (dim+1)*#T (composed of weights and rest positions) and T is a
  721. // #T*(dim+1) by dim matrix of #T stacked transposed transformation matrices.
  722. // See equations (1) and (2) in "Fast Automatic Skinning Transformations"
  723. // [Jacobson et al 2012]
  724. //
  725. // Inputs:
  726. // V #V by dim list of rest positions
  727. // W #V+ by #T list of weights
  728. // Outputs:
  729. // M #V by #T*(dim+1)
  730. //
  731. // In MATLAB:
  732. // kron(ones(1,size(W,2)),[V ones(size(V,1),1)]).*kron(W,ones(1,size(V,2)+1)))igl_Qu8mg5v7";
  733. const char *__doc_igl_lbs_matrix_column = R"igl_Qu8mg5v7(// LBS_MATRIX construct a matrix that when multiplied against a column of
  734. // affine transformation entries computes new coordinates of the vertices
  735. //
  736. // I'm not sure it makes since that the result is stored as a sparse matrix.
  737. // The number of non-zeros per row *is* dependent on the number of mesh
  738. // vertices and handles.
  739. //
  740. // Inputs:
  741. // V #V by dim list of vertex rest positions
  742. // W #V by #handles list of correspondence weights
  743. // Output:
  744. // M #V * dim by #handles * dim * (dim+1) matrix such that
  745. // new_V(:) = LBS(V,W,A) = reshape(M * A,size(V)), where A is a column
  746. // vectors formed by the entries in each handle's dim by dim+1
  747. // transformation matrix. Specifcally, A =
  748. // reshape(permute(Astack,[3 1 2]),n*dim*(dim+1),1)
  749. // or A = [Lxx;Lyx;Lxy;Lyy;tx;ty], and likewise for other dim
  750. // if Astack(:,:,i) is the dim by (dim+1) transformation at handle i)igl_Qu8mg5v7";
  751. const char *__doc_igl_local_basis = R"igl_Qu8mg5v7(// Compute a local orthogonal reference system for each triangle in the given mesh
  752. // Templates:
  753. // DerivedV derived from vertex positions matrix type: i.e. MatrixXd
  754. // DerivedF derived from face indices matrix type: i.e. MatrixXi
  755. // Inputs:
  756. // V eigen matrix #V by 3
  757. // F #F by 3 list of mesh faces (must be triangles)
  758. // Outputs:
  759. // B1 eigen matrix #F by 3, each vector is tangent to the triangle
  760. // B2 eigen matrix #F by 3, each vector is tangent to the triangle and perpendicular to B1
  761. // B3 eigen matrix #F by 3, normal of the triangle
  762. //
  763. // See also: adjacency_matrix)igl_Qu8mg5v7";
  764. const char *__doc_igl_lscm = R"igl_Qu8mg5v7(// Compute a Least-squares conformal map parametrization (equivalently
  765. // derived in "Intrinsic Parameterizations of Surface Meshes" [Desbrun et al.
  766. // 2002] and "Least Squares Conformal Maps for Automatic Texture Atlas
  767. // Generation" [Lévy et al. 2002]), though this implementation follows the
  768. // derivation in: "Spectral Conformal Parameterization" [Mullen et al. 2008]
  769. // (note, this does **not** implement the Eigen-decomposition based method in
  770. // [Mullen et al. 2008], which is not equivalent). Input should be a manifold
  771. // mesh (also no unreferenced vertices) and "boundary" (fixed vertices) `b`
  772. // should contain at least two vertices per connected component.
  773. //
  774. // Inputs:
  775. // V #V by 3 list of mesh vertex positions
  776. // F #F by 3 list of mesh faces (must be triangles)
  777. // b #b boundary indices into V
  778. // bc #b by 3 list of boundary values
  779. // Outputs:
  780. // UV #V by 2 list of 2D mesh vertex positions in UV space
  781. // Returns true only on solver success.
  782. //)igl_Qu8mg5v7";
  783. const char *__doc_igl_map_vertices_to_circle = R"igl_Qu8mg5v7(// Map the vertices whose indices are in a given boundary loop (bnd) on the
  784. // unit circle with spacing proportional to the original boundary edge
  785. // lengths.
  786. //
  787. // Inputs:
  788. // V #V by dim list of mesh vertex positions
  789. // b #W list of vertex ids
  790. // Outputs:
  791. // UV #W by 2 list of 2D position on the unit circle for the vertices in b)igl_Qu8mg5v7";
  792. const char *__doc_igl_massmatrix = R"igl_Qu8mg5v7(// Constructs the mass (area) matrix for a given mesh (V,F).
  793. //
  794. // Templates:
  795. // DerivedV derived type of eigen matrix for V (e.g. derived from
  796. // MatrixXd)
  797. // DerivedF derived type of eigen matrix for F (e.g. derived from
  798. // MatrixXi)
  799. // Scalar scalar type for eigen sparse matrix (e.g. double)
  800. // Inputs:
  801. // V #V by dim list of mesh vertex positions
  802. // F #F by simplex_size list of mesh faces (must be triangles)
  803. // type one of the following ints:
  804. // MASSMATRIX_TYPE_BARYCENTRIC barycentric
  805. // MASSMATRIX_TYPE_VORONOI voronoi-hybrid {default}
  806. // MASSMATRIX_TYPE_FULL full {not implemented}
  807. // Outputs:
  808. // M #V by #V mass matrix
  809. //
  810. // See also: adjacency_matrix
  811. //)igl_Qu8mg5v7";
  812. const char *__doc_igl_min_quad_with_fixed_precompute = R"igl_Qu8mg5v7(// Known Bugs: rows of Aeq **should probably** be linearly independent.
  813. // During precomputation, the rows of a Aeq are checked via QR. But in case
  814. // they're not then resulting probably will no longer be sparse: it will be
  815. // slow.
  816. //
  817. // MIN_QUAD_WITH_FIXED Minimize quadratic energy
  818. //
  819. // 0.5*Z'*A*Z + Z'*B + C with
  820. //
  821. // constraints that Z(known) = Y, optionally also subject to the constraints
  822. // Aeq*Z = Beq
  823. //
  824. // Templates:
  825. // T should be a eigen matrix primitive type like int or double
  826. // Inputs:
  827. // A n by n matrix of quadratic coefficients
  828. // known list of indices to known rows in Z
  829. // Y list of fixed values corresponding to known rows in Z
  830. // Aeq m by n list of linear equality constraint coefficients
  831. // pd flag specifying whether A(unknown,unknown) is positive definite
  832. // Outputs:
  833. // data factorization struct with all necessary information to solve
  834. // using min_quad_with_fixed_solve
  835. // Returns true on success, false on error
  836. //
  837. // Benchmark: For a harmonic solve on a mesh with 325K facets, matlab 2.2
  838. // secs, igl/min_quad_with_fixed.h 7.1 secs
  839. //)igl_Qu8mg5v7";
  840. const char *__doc_igl_min_quad_with_fixed_solve = R"igl_Qu8mg5v7(// Solves a system previously factored using min_quad_with_fixed_precompute
  841. //
  842. // Template:
  843. // T type of sparse matrix (e.g. double)
  844. // DerivedY type of Y (e.g. derived from VectorXd or MatrixXd)
  845. // DerivedZ type of Z (e.g. derived from VectorXd or MatrixXd)
  846. // Inputs:
  847. // data factorization struct with all necessary precomputation to solve
  848. // B n by 1 column of linear coefficients
  849. // Y b by 1 list of constant fixed values
  850. // Beq m by 1 list of linear equality constraint constant values
  851. // Outputs:
  852. // Z n by cols solution
  853. // sol #unknowns+#lagrange by cols solution to linear system
  854. // Returns true on success, false on error)igl_Qu8mg5v7";
  855. const char *__doc_igl_min_quad_with_fixed = R"igl_Qu8mg5v7(See min_quad_with_fixed for the documentation.)igl_Qu8mg5v7";
  856. const char *__doc_igl_n_polyvector = R"igl_Qu8mg5v7(// Inputs:
  857. // v0, v1 the two #3 by 1 vectors
  858. // normalized boolean, if false, then the vectors are normalized prior to the calculation
  859. // Output:
  860. // 3 by 3 rotation matrix that takes v0 to v1
  861. //)igl_Qu8mg5v7";
  862. const char *__doc_igl_normalize_row_lengths = R"igl_Qu8mg5v7(// Obsolete: just use A.rowwise().normalize() or B=A.rowwise().normalized();
  863. //
  864. // Normalize the rows in A so that their lengths are each 1 and place the new
  865. // entries in B
  866. // Inputs:
  867. // A #rows by k input matrix
  868. // Outputs:
  869. // B #rows by k input matrix, can be the same as A)igl_Qu8mg5v7";
  870. const char *__doc_igl_normalize_row_sums = R"igl_Qu8mg5v7(// Normalize the rows in A so that their sums are each 1 and place the new
  871. // entries in B
  872. // Inputs:
  873. // A #rows by k input matrix
  874. // Outputs:
  875. // B #rows by k input matrix, can be the same as A
  876. //
  877. // Note: This is just calling an Eigen one-liner.)igl_Qu8mg5v7";
  878. const char *__doc_igl_parula = R"igl_Qu8mg5v7(// PARULA like MATLAB's parula
  879. //
  880. // Inputs:
  881. // m number of colors
  882. // Outputs:
  883. // J m by list of RGB colors between 0 and 1
  884. //
  885. // Wrapper for directly computing [r,g,b] values for a given factor f between
  886. // 0 and 1
  887. //
  888. // Inputs:
  889. // f factor determining color value as if 0 was min and 1 was max
  890. // Outputs:
  891. // r red value
  892. // g green value
  893. // b blue value)igl_Qu8mg5v7";
  894. const char *__doc_igl_per_corner_normals = R"igl_Qu8mg5v7(// Compute vertex normals via vertex position list, face list
  895. // Inputs:
  896. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  897. // F #F by 3 eigne Matrix of face (triangle) indices
  898. // corner_threshold threshold in degrees on sharp angles
  899. // Output:
  900. // CN #F*3 by 3 eigen Matrix of mesh vertex 3D normals, where the normal
  901. // for corner F(i,j) is at CN(i*3+j,:) )igl_Qu8mg5v7";
  902. const char *__doc_igl_per_edge_normals = R"igl_Qu8mg5v7(// Compute face normals via vertex position list, face list
  903. // Inputs:
  904. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  905. // F #F by 3 eigen Matrix of face (triangle) indices
  906. // weight weighting type
  907. // FN #F by 3 matrix of 3D face normals per face
  908. // Output:
  909. // N #2 by 3 matrix of mesh edge 3D normals per row
  910. // E #E by 2 matrix of edge indices per row
  911. // EMAP #E by 1 matrix of indices from all edges to E
  912. //)igl_Qu8mg5v7";
  913. const char *__doc_igl_per_face_normals = R"igl_Qu8mg5v7(// Compute face normals via vertex position list, face list
  914. // Inputs:
  915. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  916. // F #F by 3 eigen Matrix of face (triangle) indices
  917. // Z 3 vector normal given to faces with degenerate normal.
  918. // Output:
  919. // N #F by 3 eigen Matrix of mesh face (triangle) 3D normals
  920. //
  921. // Example:
  922. // // Give degenerate faces (1/3,1/3,1/3)^0.5
  923. // per_face_normals(V,F,Vector3d(1,1,1).normalized(),N);)igl_Qu8mg5v7";
  924. const char *__doc_igl_per_face_normals_stable = R"igl_Qu8mg5v7(// Special version where order of face indices is guaranteed not to effect
  925. // output.)igl_Qu8mg5v7";
  926. const char *__doc_igl_per_vertex_normals = R"igl_Qu8mg5v7(// Compute vertex normals via vertex position list, face list
  927. // Inputs:
  928. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  929. // F #F by 3 eigne Matrix of face (triangle) indices
  930. // weighting Weighting type
  931. // Output:
  932. // N #V by 3 eigen Matrix of mesh vertex 3D normals)igl_Qu8mg5v7";
  933. const char *__doc_igl_planarize_quad_mesh = R"igl_Qu8mg5v7(// Inputs:
  934. // Vin #V by 3 eigen Matrix of mesh vertex 3D positions
  935. // F #F by 4 eigen Matrix of face (quad) indices
  936. // maxIter maximum numbers of iterations
  937. // threshold minimum allowed threshold for non-planarity
  938. // Output:
  939. // Vout #V by 3 eigen Matrix of planar mesh vertex 3D positions
  940. //)igl_Qu8mg5v7";
  941. const char *__doc_igl_png_readPNG = R"igl_Qu8mg5v7(// Read an image from a .png file into 4 memory buffers
  942. //
  943. // Input:
  944. // png_file path to .png file
  945. // Output:
  946. // R,G,B,A texture channels
  947. // Returns true on success, false on failure
  948. //)igl_Qu8mg5v7";
  949. const char *__doc_igl_png_writePNG = R"igl_Qu8mg5v7(// Writes an image to a png file
  950. //
  951. // Input:
  952. // R,G,B,A texture channels
  953. // Output:
  954. // png_file path to .png file
  955. // Returns true on success, false on failure
  956. //)igl_Qu8mg5v7";
  957. const char *__doc_igl_point_mesh_squared_distance = R"igl_Qu8mg5v7(// Compute distances from a set of points P to a triangle mesh (V,F)
  958. //
  959. // Inputs:
  960. // P #P by 3 list of query point positions
  961. // V #V by 3 list of vertex positions
  962. // Ele #Ele by (3|2|1) list of (triangle|edge|point) indices
  963. // Outputs:
  964. // sqrD #P list of smallest squared distances
  965. // I #P list of primitive indices corresponding to smallest distances
  966. // C #P by 3 list of closest points
  967. //
  968. // Known bugs: This only computes distances to given primitivess. So
  969. // unreferenced vertices are ignored. However, degenerate primitives are
  970. // handled correctly: triangle [1 2 2] is treated as a segment [1 2], and
  971. // triangle [1 1 1] is treated as a point. So one _could_ add extra
  972. // combinatorially degenerate rows to Ele for all unreferenced vertices to
  973. // also get distances to points.)igl_Qu8mg5v7";
  974. const char *__doc_igl_polar_svd = R"igl_Qu8mg5v7(// Computes the polar decomposition (R,T) of a matrix A using SVD singular
  975. // value decomposition
  976. //
  977. // Inputs:
  978. // A 3 by 3 matrix to be decomposed
  979. // Outputs:
  980. // R 3 by 3 rotation matrix part of decomposition (**always rotataion**)
  981. // T 3 by 3 stretch matrix part of decomposition
  982. // U 3 by 3 left-singular vectors
  983. // S 3 by 1 singular values
  984. // V 3 by 3 right-singular vectors
  985. //
  986. //)igl_Qu8mg5v7";
  987. const char *__doc_igl_principal_curvature = R"igl_Qu8mg5v7(// Compute the principal curvature directions and magnitude of the given triangle mesh
  988. // DerivedV derived from vertex positions matrix type: i.e. MatrixXd
  989. // DerivedF derived from face indices matrix type: i.e. MatrixXi
  990. // Inputs:
  991. // V eigen matrix #V by 3
  992. // F #F by 3 list of mesh faces (must be triangles)
  993. // radius controls the size of the neighbourhood used, 1 = average edge lenght
  994. //
  995. // Outputs:
  996. // PD1 #V by 3 maximal curvature direction for each vertex.
  997. // PD2 #V by 3 minimal curvature direction for each vertex.
  998. // PV1 #V by 1 maximal curvature value for each vertex.
  999. // PV2 #V by 1 minimal curvature value for each vertex.
  1000. //
  1001. // See also: average_onto_faces, average_onto_vertices
  1002. //
  1003. // This function has been developed by: Nikolas De Giorgis, Luigi Rocca and Enrico Puppo.
  1004. // The algorithm is based on:
  1005. // Efficient Multi-scale Curvature and Crease Estimation
  1006. // Daniele Panozzo, Enrico Puppo, Luigi Rocca
  1007. // GraVisMa, 2010)igl_Qu8mg5v7";
  1008. const char *__doc_igl_quad_planarity = R"igl_Qu8mg5v7(// Compute planarity of the faces of a quad mesh
  1009. // Inputs:
  1010. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  1011. // F #F by 4 eigen Matrix of face (quad) indices
  1012. // Output:
  1013. // P #F by 1 eigen Matrix of mesh face (quad) planarities
  1014. //)igl_Qu8mg5v7";
  1015. const char *__doc_igl_randperm = R"igl_Qu8mg5v7(// Like matlab's randperm(n) but minus 1
  1016. //
  1017. // Inputs:
  1018. // n number of elements
  1019. // Outputs:
  1020. // I n list of rand permutation of 0:n-1)igl_Qu8mg5v7";
  1021. const char *__doc_igl_readDMAT = R"igl_Qu8mg5v7(See readDMAT for the documentation.)igl_Qu8mg5v7";
  1022. const char *__doc_igl_readMESH = R"igl_Qu8mg5v7(// load a tetrahedral volume mesh from a .mesh file
  1023. //
  1024. // Templates:
  1025. // Scalar type for positions and vectors (will be read as double and cast
  1026. // to Scalar)
  1027. // Index type for indices (will be read as int and cast to Index)
  1028. // Input:
  1029. // mesh_file_name path of .mesh file
  1030. // Outputs:
  1031. // V double matrix of vertex positions #V by 3
  1032. // T #T list of tet indices into vertex positions
  1033. // F #F list of face indices into vertex positions
  1034. //
  1035. // Known bugs: Holes and regions are not supported)igl_Qu8mg5v7";
  1036. const char *__doc_igl_readOBJ = R"igl_Qu8mg5v7(// Read a mesh from an ascii obj file, filling in vertex positions, normals
  1037. // and texture coordinates. Mesh may have faces of any number of degree
  1038. //
  1039. // Templates:
  1040. // Scalar type for positions and vectors (will be read as double and cast
  1041. // to Scalar)
  1042. // Index type for indices (will be read as int and cast to Index)
  1043. // Inputs:
  1044. // str path to .obj file
  1045. // Outputs:
  1046. // V double matrix of vertex positions #V by 3
  1047. // TC double matrix of texture coordinats #TC by 2
  1048. // N double matrix of corner normals #N by 3
  1049. // F #F list of face indices into vertex positions
  1050. // FTC #F list of face indices into vertex texture coordinates
  1051. // FN #F list of face indices into vertex normals
  1052. // Returns true on success, false on errors)igl_Qu8mg5v7";
  1053. const char *__doc_igl_readOFF = R"igl_Qu8mg5v7(// Read a mesh from an ascii obj file, filling in vertex positions, normals
  1054. // and texture coordinates. Mesh may have faces of any number of degree
  1055. //
  1056. // Templates:
  1057. // Scalar type for positions and vectors (will be read as double and cast
  1058. // to Scalar)
  1059. // Index type for indices (will be read as int and cast to Index)
  1060. // Inputs:
  1061. // str path to .obj file
  1062. // Outputs:
  1063. // V double matrix of vertex positions #V by 3
  1064. // F #F list of face indices into vertex positions
  1065. // TC double matrix of texture coordinats #TC by 2
  1066. // FTC #F list of face indices into vertex texture coordinates
  1067. // N double matrix of corner normals #N by 3
  1068. // FN #F list of face indices into vertex normals
  1069. // Returns true on success, false on errors)igl_Qu8mg5v7";
  1070. const char *__doc_igl_readTGF = R"igl_Qu8mg5v7(// READTGF
  1071. //
  1072. // [V,E,P,BE,CE,PE] = readTGF(filename)
  1073. //
  1074. // Read a graph from a .tgf file
  1075. //
  1076. // Input:
  1077. // filename .tgf file name
  1078. // Ouput:
  1079. // V # vertices by 3 list of vertex positions
  1080. // E # edges by 2 list of edge indices
  1081. // P # point-handles list of point handle indices
  1082. // BE # bone-edges by 2 list of bone-edge indices
  1083. // CE # cage-edges by 2 list of cage-edge indices
  1084. // PE # pseudo-edges by 2 list of pseudo-edge indices
  1085. //
  1086. // Assumes that graph vertices are 3 dimensional)igl_Qu8mg5v7";
  1087. const char *__doc_igl_read_triangle_mesh = R"igl_Qu8mg5v7(// read mesh from an ascii file with automatic detection of file format.
  1088. // supported: obj, off, stl, wrl, ply, mesh)
  1089. //
  1090. // Templates:
  1091. // Scalar type for positions and vectors (will be read as double and cast
  1092. // to Scalar)
  1093. // Index type for indices (will be read as int and cast to Index)
  1094. // Inputs:
  1095. // str path to file
  1096. // Outputs:
  1097. // V eigen double matrix #V by 3
  1098. // F eigen int matrix #F by 3
  1099. // Returns true iff success)igl_Qu8mg5v7";
  1100. const char *__doc_igl_remove_duplicate_vertices = R"igl_Qu8mg5v7(// REMOVE_DUPLICATE_VERTICES Remove duplicate vertices upto a uniqueness
  1101. // tolerance (epsilon)
  1102. //
  1103. // Inputs:
  1104. // V #V by dim list of vertex positions
  1105. // epsilon uniqueness tolerance (significant digit), can probably think of
  1106. // this as a tolerance on L1 distance
  1107. // Outputs:
  1108. // SV #SV by dim new list of vertex positions
  1109. // SVI #V by 1 list of indices so SV = V(SVI,:)
  1110. // SVJ #SV by 1 list of indices so V = SV(SVJ,:)
  1111. //
  1112. // Example:
  1113. // % Mesh in (V,F)
  1114. // [SV,SVI,SVJ] = remove_duplicate_vertices(V,1e-7);
  1115. // % remap faces
  1116. // SF = SVJ(F);
  1117. //)igl_Qu8mg5v7";
  1118. const char *__doc_igl_rotate_vectors = R"igl_Qu8mg5v7(// Rotate the vectors V by A radiants on the tangent plane spanned by B1 and
  1119. // B2
  1120. //
  1121. // Inputs:
  1122. // V #V by 3 eigen Matrix of vectors
  1123. // A #V eigen vector of rotation angles or a single angle to be applied
  1124. // to all vectors
  1125. // B1 #V by 3 eigen Matrix of base vector 1
  1126. // B2 #V by 3 eigen Matrix of base vector 2
  1127. //
  1128. // Output:
  1129. // Returns the rotated vectors
  1130. //)igl_Qu8mg5v7";
  1131. const char *__doc_igl_setdiff = R"igl_Qu8mg5v7(// Set difference of elements of matrices
  1132. //
  1133. // Inputs:
  1134. // A m-long vector of indices
  1135. // B n-long vector of indices
  1136. // Outputs:
  1137. // C (k<=m)-long vector of unique elements appearing in A but not in B
  1138. // IA (k<=m)-long list of indices into A so that C = A(IA)
  1139. //)igl_Qu8mg5v7";
  1140. const char *__doc_igl_signed_distance = R"igl_Qu8mg5v7(// Computes signed distance to a mesh
  1141. //
  1142. // Inputs:
  1143. // P #P by 3 list of query point positions
  1144. // V #V by 3 list of vertex positions
  1145. // F #F by ss list of triangle indices, ss should be 3 unless sign_type ==
  1146. // SIGNED_DISTANCE_TYPE_UNSIGNED
  1147. // sign_type method for computing distance _sign_ S
  1148. // Outputs:
  1149. // S #P list of smallest signed distances
  1150. // I #P list of facet indices corresponding to smallest distances
  1151. // C #P by 3 list of closest points
  1152. // N #P by 3 list of closest normals (only set if
  1153. // sign_type=SIGNED_DISTANCE_TYPE_PSEUDONORMAL)
  1154. //
  1155. // Known bugs: This only computes distances to triangles. So unreferenced
  1156. // vertices and degenerate triangles are ignored.)igl_Qu8mg5v7";
  1157. const char *__doc_igl_signed_distance_pseudonormal = R"igl_Qu8mg5v7(// Computes signed distance to mesh
  1158. //
  1159. // Inputs:
  1160. // tree AABB acceleration tree (see AABB.h)
  1161. // F #F by 3 list of triangle indices
  1162. // FN #F by 3 list of triangle normals
  1163. // VN #V by 3 list of vertex normals (ANGLE WEIGHTING)
  1164. // EN #E by 3 list of edge normals (UNIFORM WEIGHTING)
  1165. // EMAP #F*3 mapping edges in F to E
  1166. // q Query point
  1167. // Returns signed distance to mesh
  1168. //)igl_Qu8mg5v7";
  1169. const char *__doc_igl_signed_distance_winding_number = R"igl_Qu8mg5v7(// Inputs:
  1170. // tree AABB acceleration tree (see cgal/point_mesh_squared_distance.h)
  1171. // hier Winding number evaluation hierarchy
  1172. // q Query point
  1173. // Returns signed distance to mesh)igl_Qu8mg5v7";
  1174. const char *__doc_igl_slice = R"igl_Qu8mg5v7(// Act like the matlab X(row_indices,col_indices) operator, where
  1175. // row_indices, col_indices are non-negative integer indices.
  1176. //
  1177. // Inputs:
  1178. // X m by n matrix
  1179. // R list of row indices
  1180. // C list of column indices
  1181. // Output:
  1182. // Y #R by #C matrix
  1183. //
  1184. // See also: slice_mask)igl_Qu8mg5v7";
  1185. const char *__doc_igl_slice_into = R"igl_Qu8mg5v7(// Act like the matlab Y(row_indices,col_indices) = X
  1186. //
  1187. // Inputs:
  1188. // X xm by xn rhs matrix
  1189. // R list of row indices
  1190. // C list of column indices
  1191. // Y ym by yn lhs matrix
  1192. // Output:
  1193. // Y ym by yn lhs matrix, same as input but Y(R,C) = X)igl_Qu8mg5v7";
  1194. const char *__doc_igl_slice_mask = R"igl_Qu8mg5v7(// Act like the matlab X(row_mask,col_mask) operator, where
  1195. // row_mask, col_mask are non-negative integer indices.
  1196. //
  1197. // Inputs:
  1198. // X m by n matrix
  1199. // R m list of row bools
  1200. // C n list of column bools
  1201. // Output:
  1202. // Y #trues-in-R by #trues-in-C matrix
  1203. //
  1204. // See also: slice_mask)igl_Qu8mg5v7";
  1205. const char *__doc_igl_slice_tets = R"igl_Qu8mg5v7(// SLICE_TETS Slice through a tet mesh (V,T) along a given plane (via its
  1206. // implicit equation).
  1207. //
  1208. // Inputs:
  1209. // V #V by 3 list of tet mesh vertices
  1210. // T #T by 4 list of tet indices into V
  1211. // plane list of 4 coefficients in the plane equation: [x y z 1]'*plane = 0
  1212. // Optional:
  1213. // 'Manifold' followed by whether to stitch together triangles into a
  1214. // manifold mesh {true}: results in more compact U but slightly slower.
  1215. // Outputs:
  1216. // U #U by 3 list of triangle mesh vertices along slice
  1217. // G #G by 3 list of triangles indices into U
  1218. // J #G list of indices into T revealing from which tet each faces comes
  1219. // BC #U by #V list of barycentric coordinates (or more generally: linear
  1220. // interpolation coordinates) so that U = BC*V
  1221. // )igl_Qu8mg5v7";
  1222. const char *__doc_igl_sortrows = R"igl_Qu8mg5v7(// Act like matlab's [Y,I] = sortrows(X)
  1223. //
  1224. // Templates:
  1225. // DerivedX derived scalar type, e.g. MatrixXi or MatrixXd
  1226. // DerivedI derived integer type, e.g. MatrixXi
  1227. // Inputs:
  1228. // X m by n matrix whose entries are to be sorted
  1229. // ascending sort ascending (true, matlab default) or descending (false)
  1230. // Outputs:
  1231. // Y m by n matrix whose entries are sorted (**should not** be same
  1232. // reference as X)
  1233. // I m list of indices so that
  1234. // Y = X(I,:);)igl_Qu8mg5v7";
  1235. const char *__doc_igl_streamlines_init = R"igl_Qu8mg5v7(// Given a mesh and a field the function computes the /data/ necessary for tracing the field'
  1236. // streamlines, and creates the initial /state/ for the tracing.
  1237. // Inputs:
  1238. // V #V by 3 list of mesh vertex coordinates
  1239. // F #F by 3 list of mesh faces
  1240. // temp_field #F by 3n list of the 3D coordinates of the per-face vectors
  1241. // (n-degrees stacked horizontally for each triangle)
  1242. // treat_as_symmetric
  1243. // if true, adds n symmetry directions to the field (N = 2n). Else N = n
  1244. // percentage [0-1] percentage of faces sampled
  1245. // Outputs:
  1246. // data struct containing topology information of the mesh and field
  1247. // state struct containing the state of the tracing)igl_Qu8mg5v7";
  1248. const char *__doc_igl_streamlines_next = R"igl_Qu8mg5v7(// The function computes the next state for each point in the sample
  1249. // V #V by 3 list of mesh vertex coordinates
  1250. // F #F by 3 list of mesh faces
  1251. // data struct containing topology information
  1252. // state struct containing the state of the tracing)igl_Qu8mg5v7";
  1253. const char *__doc_igl_triangle_triangle_adjacency = R"igl_Qu8mg5v7(// Constructs the triangle-triangle adjacency matrix for a given
  1254. // mesh (V,F).
  1255. //
  1256. // Templates:
  1257. // Scalar derived type of eigen matrix for V (e.g. derived from
  1258. // MatrixXd)
  1259. // Index derived type of eigen matrix for F (e.g. derived from
  1260. // MatrixXi)
  1261. // Inputs:
  1262. // F #F by simplex_size list of mesh faces (must be triangles)
  1263. // Outputs:
  1264. // TT #F by #3 adjacent matrix, the element i,j is the id of the triangle adjacent to the j edge of triangle i
  1265. // TTi #F by #3 adjacent matrix, the element i,j is the id of edge of the triangle TT(i,j) that is adjacent with triangle i
  1266. // NOTE: the first edge of a triangle is [0,1] the second [1,2] and the third [2,3].
  1267. // this convention is DIFFERENT from cotmatrix_entries.h
  1268. // Known bug: this should not need to take V as input.)igl_Qu8mg5v7";
  1269. const char *__doc_igl_triangle_triangle_adjacency_preprocess = R"igl_Qu8mg5v7(// Preprocessing)igl_Qu8mg5v7";
  1270. const char *__doc_igl_triangle_triangle_adjacency_extractTT = R"igl_Qu8mg5v7(// Extract the face adjacencies)igl_Qu8mg5v7";
  1271. const char *__doc_igl_triangle_triangle_adjacency_extractTTi = R"igl_Qu8mg5v7(// Extract the face adjacencies indices (needed for fast traversal))igl_Qu8mg5v7";
  1272. const char *__doc_igl_triangle_triangulate = R"igl_Qu8mg5v7(// Triangulate the interior of a polygon using the triangle library.
  1273. //
  1274. // Inputs:
  1275. // V #V by 2 list of 2D vertex positions
  1276. // E #E by 2 list of vertex ids forming unoriented edges of the boundary of the polygon
  1277. // H #H by 2 coordinates of points contained inside holes of the polygon
  1278. // flags string of options pass to triangle (see triangle documentation)
  1279. // Outputs:
  1280. // V2 #V2 by 2 coordinates of the vertives of the generated triangulation
  1281. // F2 #F2 by 3 list of indices forming the faces of the generated triangulation
  1282. //
  1283. // TODO: expose the option to prevent Steiner points on the boundary
  1284. //)igl_Qu8mg5v7";
  1285. const char *__doc_igl_unique = R"igl_Qu8mg5v7(// Act like matlab's [C,IA,IC] = unique(X)
  1286. //
  1287. // Templates:
  1288. // T comparable type T
  1289. // Inputs:
  1290. // A #A vector of type T
  1291. // Outputs:
  1292. // C #C vector of unique entries in A
  1293. // IA #C index vector so that C = A(IA);
  1294. // IC #A index vector so that A = C(IC);)igl_Qu8mg5v7";
  1295. const char *__doc_igl_unique_rows = R"igl_Qu8mg5v7(// Act like matlab's [C,IA,IC] = unique(X,'rows')
  1296. //
  1297. // Templates:
  1298. // DerivedA derived scalar type, e.g. MatrixXi or MatrixXd
  1299. // DerivedIA derived integer type, e.g. MatrixXi
  1300. // DerivedIC derived integer type, e.g. MatrixXi
  1301. // Inputs:
  1302. // A m by n matrix whose entries are to unique'd according to rows
  1303. // Outputs:
  1304. // C #C vector of unique rows in A
  1305. // IA #C index vector so that C = A(IA,:);
  1306. // IC #A index vector so that A = C(IC,:);)igl_Qu8mg5v7";
  1307. const char *__doc_igl_unproject_onto_mesh = R"igl_Qu8mg5v7(// Unproject a screen location (using current opengl viewport, projection, and
  1308. // model view) to a 3D position _onto_ a given mesh, if the ray through the
  1309. // given screen location (x,y) _hits_ the mesh.
  1310. //
  1311. // Inputs:
  1312. // pos screen space coordinates
  1313. // model model matrix
  1314. // proj projection matrix
  1315. // viewport vieweport vector
  1316. // V #V by 3 list of mesh vertex positions
  1317. // F #F by 3 list of mesh triangle indices into V
  1318. // Outputs:
  1319. // fid id of the first face hit
  1320. // bc barycentric coordinates of hit
  1321. // Returns true if there's a hit)igl_Qu8mg5v7";
  1322. const char *__doc_igl_upsample = R"igl_Qu8mg5v7(// Subdivide a mesh without moving vertices: loop subdivision but odd
  1323. // vertices stay put and even vertices are just edge midpoints
  1324. //
  1325. // Templates:
  1326. // MatV matrix for vertex positions, e.g. MatrixXd
  1327. // MatF matrix for vertex positions, e.g. MatrixXi
  1328. // Inputs:
  1329. // V #V by dim mesh vertices
  1330. // F #F by 3 mesh triangles
  1331. // Outputs:
  1332. // NV new vertex positions, V is guaranteed to be at top
  1333. // NF new list of face indices
  1334. //
  1335. // NOTE: V should not be the same as NV,
  1336. // NOTE: F should not be the same as NF, use other proto
  1337. //
  1338. // Known issues:
  1339. // - assumes (V,F) is edge-manifold.)igl_Qu8mg5v7";
  1340. const char *__doc_igl_winding_number = R"igl_Qu8mg5v7(// WINDING_NUMBER Compute the sum of solid angles of a triangle/tetrahedron
  1341. // described by points (vectors) V
  1342. //
  1343. // Templates:
  1344. // dim dimension of input
  1345. // Inputs:
  1346. // V n by 3 list of vertex positions
  1347. // F #F by 3 list of triangle indices, minimum index is 0
  1348. // O no by 3 list of origin positions
  1349. // Outputs:
  1350. // S no by 1 list of winding numbers
  1351. //
  1352. // 3d)igl_Qu8mg5v7";
  1353. const char *__doc_igl_winding_number_3 = R"igl_Qu8mg5v7(// Inputs:
  1354. // V pointer to array containing #V by 3 vertex positions along rows,
  1355. // given in column major order
  1356. // n number of mesh vertices
  1357. // F pointer to array containing #F by 3 face indices along rows,
  1358. // given in column major order
  1359. // m number of faces
  1360. // O pointer to array containing #O by 3 query positions along rows,
  1361. // given in column major order
  1362. // no number of origins
  1363. // Outputs:
  1364. // S no by 1 list of winding numbers)igl_Qu8mg5v7";
  1365. const char *__doc_igl_winding_number_2 = R"igl_Qu8mg5v7(//// Only one evaluation origin
  1366. //template <typename DerivedF>
  1367. //IGL_INLINE void winding_number_3(
  1368. // const double * V,
  1369. // const int n,
  1370. // const DerivedF * F,
  1371. // const int m,
  1372. // const double * O,
  1373. // double * S);
  1374. // 2d)igl_Qu8mg5v7";
  1375. const char *__doc_igl_writeMESH = R"igl_Qu8mg5v7(// save a tetrahedral volume mesh to a .mesh file
  1376. //
  1377. // Templates:
  1378. // Scalar type for positions and vectors (will be cast as double)
  1379. // Index type for indices (will be cast to int)
  1380. // Input:
  1381. // mesh_file_name path of .mesh file
  1382. // V double matrix of vertex positions #V by 3
  1383. // T #T list of tet indices into vertex positions
  1384. // F #F list of face indices into vertex positions
  1385. //
  1386. // Known bugs: Holes and regions are not supported)igl_Qu8mg5v7";
  1387. const char *__doc_igl_writeOBJ = R"igl_Qu8mg5v7(// Write a mesh in an ascii obj file
  1388. // Inputs:
  1389. // str path to outputfile
  1390. // V #V by 3 mesh vertex positions
  1391. // F #F by 3|4 mesh indices into V
  1392. // CN #CN by 3 normal vectors
  1393. // FN #F by 3|4 corner normal indices into CN
  1394. // TC #TC by 2|3 texture coordinates
  1395. // FTC #F by 3|4 corner texture coord indices into TC
  1396. // Returns true on success, false on error)igl_Qu8mg5v7";