py_doc.cpp 47 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019
  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_boundary_facets = R"igl_Qu8mg5v7(// BOUNDARY_FACETS Determine boundary faces (edges) of tetrahedra (triangles)
  80. // stored in T (analogous to qptoolbox's `outline` and `boundary_faces`).
  81. //
  82. // Templates:
  83. // IntegerT integer-value: e.g. int
  84. // IntegerF integer-value: e.g. int
  85. // Input:
  86. // T tetrahedron (triangle) index list, m by 4 (3), where m is the number of tetrahedra
  87. // Output:
  88. // F list of boundary faces, n by 3 (2), where n is the number of boundary faces
  89. //
  90. //)igl_Qu8mg5v7";
  91. const char *__doc_igl_boundary_loop = R"igl_Qu8mg5v7(// Compute list of ordered boundary loops for a manifold mesh.
  92. //
  93. // Templates:
  94. // Index index type
  95. // Inputs:
  96. // F #V by dim list of mesh faces
  97. // Outputs:
  98. // L list of loops where L[i] = ordered list of boundary vertices in loop i
  99. //)igl_Qu8mg5v7";
  100. const char *__doc_igl_cat = R"igl_Qu8mg5v7(// Perform concatenation of a two matrices along a single dimension
  101. // If dim == 1, then C = [A;B]. If dim == 2 then C = [A B]
  102. //
  103. // Template:
  104. // Scalar scalar data type for sparse matrices like double or int
  105. // Mat matrix type for all matrices (e.g. MatrixXd, SparseMatrix)
  106. // MatC matrix type for ouput matrix (e.g. MatrixXd) needs to support
  107. // resize
  108. // Inputs:
  109. // A first input matrix
  110. // B second input matrix
  111. // dim dimension along which to concatenate, 0 or 1
  112. // Outputs:
  113. // C output matrix
  114. // )igl_Qu8mg5v7";
  115. const char *__doc_igl_colon = R"igl_Qu8mg5v7(// Colon operator like matlab's colon operator. Enumerats values between low
  116. // and hi with step step.
  117. // Templates:
  118. // L should be a eigen matrix primitive type like int or double
  119. // S should be a eigen matrix primitive type like int or double
  120. // H should be a eigen matrix primitive type like int or double
  121. // T should be a eigen matrix primitive type like int or double
  122. // Inputs:
  123. // low starting value if step is valid then this is *always* the first
  124. // element of I
  125. // step step difference between sequential elements returned in I,
  126. // remember this will be cast to template T at compile time. If low<hi
  127. // then step must be positive. If low>hi then step must be negative.
  128. // Otherwise I will be set to empty.
  129. // hi ending value, if (hi-low)%step is zero then this will be the last
  130. // element in I. If step is positive there will be no elements greater
  131. // than hi, vice versa if hi<low
  132. // Output:
  133. // I list of values from low to hi with step size step)igl_Qu8mg5v7";
  134. const char *__doc_igl_comb_cross_field = R"igl_Qu8mg5v7(// Inputs:
  135. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  136. // F #F by 4 eigen Matrix of face (quad) indices
  137. // PD1in #F by 3 eigen Matrix of the first per face cross field vector
  138. // PD2in #F by 3 eigen Matrix of the second per face cross field vector
  139. // Output:
  140. // PD1out #F by 3 eigen Matrix of the first combed cross field vector
  141. // PD2out #F by 3 eigen Matrix of the second combed cross field vector
  142. //)igl_Qu8mg5v7";
  143. const char *__doc_igl_comb_frame_field = R"igl_Qu8mg5v7(// Inputs:
  144. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  145. // F #F by 4 eigen Matrix of face (quad) indices
  146. // PD1 #F by 3 eigen Matrix of the first per face cross field vector
  147. // PD2 #F by 3 eigen Matrix of the second per face cross field vector
  148. // BIS1_combed #F by 3 eigen Matrix of the first combed bisector field vector
  149. // BIS2_combed #F by 3 eigen Matrix of the second combed bisector field vector
  150. // Output:
  151. // PD1_combed #F by 3 eigen Matrix of the first combed cross field vector
  152. // PD2_combed #F by 3 eigen Matrix of the second combed cross field vector
  153. //)igl_Qu8mg5v7";
  154. const char *__doc_igl_compute_frame_field_bisectors = R"igl_Qu8mg5v7(// Compute bisectors of a frame field defined on mesh faces
  155. // Inputs:
  156. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  157. // F #F by 3 eigen Matrix of face (triangle) indices
  158. // B1 #F by 3 eigen Matrix of face (triangle) base vector 1
  159. // B2 #F by 3 eigen Matrix of face (triangle) base vector 2
  160. // PD1 #F by 3 eigen Matrix of the first per face frame field vector
  161. // PD2 #F by 3 eigen Matrix of the second per face frame field vector
  162. // Output:
  163. // BIS1 #F by 3 eigen Matrix of the first per face frame field bisector
  164. // BIS2 #F by 3 eigen Matrix of the second per face frame field bisector
  165. //)igl_Qu8mg5v7";
  166. const char *__doc_igl_copyleft_cgal_mesh_boolean = R"igl_Qu8mg5v7(// MESH_BOOLEAN Compute boolean csg operations on "solid", consistently
  167. // oriented meshes.
  168. //
  169. // Inputs:
  170. // VA #VA by 3 list of vertex positions of first mesh
  171. // FA #FA by 3 list of triangle indices into VA
  172. // VB #VB by 3 list of vertex positions of second mesh
  173. // FB #FB by 3 list of triangle indices into VB
  174. // type type of boolean operation
  175. // Outputs:
  176. // VC #VC by 3 list of vertex positions of boolean result mesh
  177. // FC #FC by 3 list of triangle indices into VC
  178. // J #FC list of indices into [FA;FA.rows()+FB] revealing "birth" facet
  179. // Returns true if inputs induce a piecewise constant winding number
  180. // field and type is valid
  181. //
  182. // See also: mesh_boolean_cork, intersect_other,
  183. // remesh_self_intersections)igl_Qu8mg5v7";
  184. const char *__doc_igl_copyleft_comiso_miq = R"igl_Qu8mg5v7(// Inputs:
  185. // V #V by 3 list of mesh vertex 3D positions
  186. // F #F by 3 list of faces indices in V
  187. // PD1 #V by 3 first line of the Jacobian per triangle
  188. // PD2 #V by 3 second line of the Jacobian per triangle
  189. // (optional, if empty it will be a vector in the tangent plane orthogonal to PD1)
  190. // scale global scaling for the gradient (controls the quads resolution)
  191. // stiffness weight for the stiffness iterations
  192. // direct_round greedily round all integer variables at once (greatly improves optimization speed but lowers quality)
  193. // iter stiffness iterations (0 = no stiffness)
  194. // local_iter number of local iterations for the integer rounding
  195. // do_round enables the integer rounding (disabling it could be useful for debugging)
  196. // round_vertices id of additional vertices that should be snapped to integer coordinates
  197. // hard_features #H by 2 list of pairs of vertices that belongs to edges that should be snapped to integer coordinates
  198. //
  199. // Output:
  200. // UV #UV by 2 list of vertices in 2D
  201. // FUV #FUV by 3 list of face indices in UV
  202. //
  203. // TODO: rename the parameters name in the cpp consistenly
  204. // improve the handling of hard_features, right now it might fail in difficult cases)igl_Qu8mg5v7";
  205. const char *__doc_igl_copyleft_comiso_nrosy = R"igl_Qu8mg5v7(// Generate a N-RoSy field from a sparse set of constraints
  206. //
  207. // Inputs:
  208. // V #V by 3 list of mesh vertex coordinates
  209. // F #F by 3 list of mesh faces (must be triangles)
  210. // b #B by 1 list of constrained face indices
  211. // bc #B by 3 list of representative vectors for the constrained
  212. // faces
  213. // b_soft #S by 1 b for soft constraints
  214. // w_soft #S by 1 weight for the soft constraints (0-1)
  215. // bc_soft #S by 3 bc for soft constraints
  216. // N the degree of the N-RoSy vector field
  217. // soft the strenght of the soft contraints w.r.t. smoothness
  218. // (0 -> smoothness only, 1->constraints only)
  219. // Outputs:
  220. // R #F by 3 the representative vectors of the interpolated field
  221. // S #V by 1 the singularity index for each vertex (0 = regular))igl_Qu8mg5v7";
  222. const char *__doc_igl_copyleft_tetgen_tetrahedralize = R"igl_Qu8mg5v7(// Mesh the interior of a surface mesh (V,F) using tetgen
  223. //
  224. // Inputs:
  225. // V #V by 3 vertex position list
  226. // F #F list of polygon face indices into V (0-indexed)
  227. // switches string of tetgen options (See tetgen documentation) e.g.
  228. // "pq1.414a0.01" tries to mesh the interior of a given surface with
  229. // quality and area constraints
  230. // "" will mesh the convex hull constrained to pass through V (ignores F)
  231. // Outputs:
  232. // TV #V by 3 vertex position list
  233. // TT #T by 4 list of tet face indices
  234. // TF #F by 3 list of triangle face indices
  235. // Returns status:
  236. // 0 success
  237. // 1 tetgen threw exception
  238. // 2 tetgen did not crash but could not create any tets (probably there are
  239. // holes, duplicate faces etc.)
  240. // -1 other error)igl_Qu8mg5v7";
  241. const char *__doc_igl_cotmatrix = R"igl_Qu8mg5v7(// Constructs the cotangent stiffness matrix (discrete laplacian) for a given
  242. // mesh (V,F).
  243. //
  244. // Templates:
  245. // DerivedV derived type of eigen matrix for V (e.g. derived from
  246. // MatrixXd)
  247. // DerivedF derived type of eigen matrix for F (e.g. derived from
  248. // MatrixXi)
  249. // Scalar scalar type for eigen sparse matrix (e.g. double)
  250. // Inputs:
  251. // V #V by dim list of mesh vertex positions
  252. // F #F by simplex_size list of mesh faces (must be triangles)
  253. // Outputs:
  254. // L #V by #V cotangent matrix, each row i corresponding to V(i,:)
  255. //
  256. // See also: adjacency_matrix
  257. //
  258. // Note: This Laplacian uses the convention that diagonal entries are
  259. // **minus** the sum of off-diagonal entries. The diagonal entries are
  260. // therefore in general negative and the matrix is **negative** semi-definite
  261. // (immediately, -L is **positive** semi-definite)
  262. //
  263. // Known bugs: off by 1e-16 on regular grid. I think its a problem of
  264. // arithmetic order in cotmatrix_entries.h: C(i,e) = (arithmetic)/dblA/4)igl_Qu8mg5v7";
  265. const char *__doc_igl_covariance_scatter_matrix = R"igl_Qu8mg5v7(// Construct the covariance scatter matrix for a given arap energy
  266. // Inputs:
  267. // V #V by Vdim list of initial domain positions
  268. // F #F by 3 list of triangle indices into V
  269. // energy ARAPEnergyType enum value defining which energy is being used.
  270. // See ARAPEnergyType.h for valid options and explanations.
  271. // Outputs:
  272. // CSM dim*#V/#F by dim*#V sparse matrix containing special laplacians along
  273. // the diagonal so that when multiplied by V gives covariance matrix
  274. // elements, can be used to speed up covariance matrix computation)igl_Qu8mg5v7";
  275. const char *__doc_igl_cross_field_missmatch = R"igl_Qu8mg5v7(// Inputs:
  276. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  277. // F #F by 3 eigen Matrix of face (quad) indices
  278. // PD1 #F by 3 eigen Matrix of the first per face cross field vector
  279. // PD2 #F by 3 eigen Matrix of the second per face cross field vector
  280. // isCombed boolean, specifying whether the field is combed (i.e. matching has been precomputed.
  281. // If not, the field is combed first.
  282. // Output:
  283. // Handle_MMatch #F by 3 eigen Matrix containing the integer missmatch of the cross field
  284. // across all face edges
  285. //)igl_Qu8mg5v7";
  286. 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
  287. // (MMatch), finds the cut_graph connecting the singularities (seams) and the
  288. // degree of the singularities singularity_index
  289. //
  290. // Input:
  291. // V #V by 3 list of mesh vertex positions
  292. // F #F by 3 list of faces
  293. // MMatch #F by 3 list of per corner integer mismatch
  294. // Outputs:
  295. // seams #F by 3 list of per corner booleans that denotes if an edge is a
  296. // seam or not
  297. //)igl_Qu8mg5v7";
  298. const char *__doc_igl_doublearea = R"igl_Qu8mg5v7(// DOUBLEAREA computes twice the area for each input triangle[quad]
  299. //
  300. // Templates:
  301. // DerivedV derived type of eigen matrix for V (e.g. derived from
  302. // MatrixXd)
  303. // DerivedF derived type of eigen matrix for F (e.g. derived from
  304. // MatrixXi)
  305. // DeriveddblA derived type of eigen matrix for dblA (e.g. derived from
  306. // MatrixXd)
  307. // Inputs:
  308. // V #V by dim list of mesh vertex positions
  309. // F #F by simplex_size list of mesh faces (must be triangles or quads)
  310. // Outputs:
  311. // dblA #F list of triangle[quad] double areas (SIGNED only for 2D input)
  312. //
  313. // Known bug: For dim==3 complexity is O(#V + #F)!! Not just O(#F). This is a big deal
  314. // if you have 1million unreferenced vertices and 1 face)igl_Qu8mg5v7";
  315. const char *__doc_igl_doublearea_single = R"igl_Qu8mg5v7(// Single triangle in 2D!
  316. //
  317. // This should handle streams of corners not just single corners)igl_Qu8mg5v7";
  318. const char *__doc_igl_doublearea_quad = R"igl_Qu8mg5v7(// DOUBLEAREA_QUAD computes twice the area for each input quadrilateral
  319. //
  320. // Inputs:
  321. // V #V by dim list of mesh vertex positions
  322. // F #F by simplex_size list of mesh faces (must be quadrilaterals)
  323. // Outputs:
  324. // dblA #F list of quadrilateral double areas
  325. //)igl_Qu8mg5v7";
  326. const char *__doc_igl_edge_lengths = R"igl_Qu8mg5v7(// Constructs a list of lengths of edges opposite each index in a face
  327. // (triangle/tet) list
  328. //
  329. // Templates:
  330. // DerivedV derived from vertex positions matrix type: i.e. MatrixXd
  331. // DerivedF derived from face indices matrix type: i.e. MatrixXi
  332. // DerivedL derived from edge lengths matrix type: i.e. MatrixXd
  333. // Inputs:
  334. // V eigen matrix #V by 3
  335. // F #F by 2 list of mesh edges
  336. // or
  337. // F #F by 3 list of mesh faces (must be triangles)
  338. // or
  339. // T #T by 4 list of mesh elements (must be tets)
  340. // Outputs:
  341. // L #F by {1|3|6} list of edge lengths
  342. // for edges, column of lengths
  343. // for triangles, columns correspond to edges [1,2],[2,0],[0,1]
  344. // for tets, columns correspond to edges
  345. // [3 0],[3 1],[3 2],[1 2],[2 0],[0 1]
  346. //)igl_Qu8mg5v7";
  347. const char *__doc_igl_eigs = R"igl_Qu8mg5v7(See eigs for the documentation.)igl_Qu8mg5v7";
  348. const char *__doc_igl_embree_ambient_occlusion = R"igl_Qu8mg5v7(// Compute ambient occlusion per given point
  349. //
  350. // Inputs:
  351. // ei EmbreeIntersector containing (V,F)
  352. // P #P by 3 list of origin points
  353. // N #P by 3 list of origin normals
  354. // Outputs:
  355. // S #P list of ambient occlusion values between 1 (fully occluded) and
  356. // 0 (not occluded)
  357. //)igl_Qu8mg5v7";
  358. const char *__doc_igl_find_cross_field_singularities = R"igl_Qu8mg5v7(// Inputs:
  359. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  360. // F #F by 3 eigen Matrix of face (quad) indices
  361. // Handle_MMatch #F by 3 eigen Matrix containing the integer missmatch of the cross field
  362. // across all face edges
  363. // Output:
  364. // isSingularity #V by 1 boolean eigen Vector indicating the presence of a singularity on a vertex
  365. // singularityIndex #V by 1 integer eigen Vector containing the singularity indices
  366. //)igl_Qu8mg5v7";
  367. const char *__doc_igl_fit_rotations = R"igl_Qu8mg5v7(// Known issues: This seems to be implemented in Eigen/Geometry:
  368. // Eigen::umeyama
  369. //
  370. // FIT_ROTATIONS Given an input mesh and new positions find rotations for
  371. // every covariance matrix in a stack of covariance matrices
  372. //
  373. // Inputs:
  374. // S nr*dim by dim stack of covariance matrices
  375. // single_precision whether to use single precision (faster)
  376. // Outputs:
  377. // R dim by dim * nr list of rotations
  378. //)igl_Qu8mg5v7";
  379. const char *__doc_igl_fit_rotations_planar = R"igl_Qu8mg5v7(// FIT_ROTATIONS Given an input mesh and new positions find 2D rotations for
  380. // every vertex that best maps its one ring to the new one ring
  381. //
  382. // Inputs:
  383. // S nr*dim by dim stack of covariance matrices, third column and every
  384. // third row will be ignored
  385. // Outputs:
  386. // R dim by dim * nr list of rotations, third row and third column of each
  387. // rotation will just be identity
  388. //)igl_Qu8mg5v7";
  389. const char *__doc_igl_fit_rotations_SSE = R"igl_Qu8mg5v7(See fit_rotations_SSE for the documentation.)igl_Qu8mg5v7";
  390. const char *__doc_igl_floor = R"igl_Qu8mg5v7(// Floor a given matrix to nearest integers
  391. //
  392. // Inputs:
  393. // X m by n matrix of scalars
  394. // Outputs:
  395. // Y m by n matrix of floored integers)igl_Qu8mg5v7";
  396. const char *__doc_igl_gaussian_curvature = R"igl_Qu8mg5v7(// Compute discrete local integral gaussian curvature (angle deficit, without
  397. // averaging by local area).
  398. //
  399. // Inputs:
  400. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  401. // F #F by 3 eigen Matrix of face (triangle) indices
  402. // Output:
  403. // K #V by 1 eigen Matrix of discrete gaussian curvature values
  404. //)igl_Qu8mg5v7";
  405. const char *__doc_igl_grad = R"igl_Qu8mg5v7(// Gradient of a scalar function defined on piecewise linear elements (mesh)
  406. // is constant on each triangle i,j,k:
  407. // grad(Xijk) = (Xj-Xi) * (Vi - Vk)^R90 / 2A + (Xk-Xi) * (Vj - Vi)^R90 / 2A
  408. // where Xi is the scalar value at vertex i, Vi is the 3D position of vertex
  409. // i, and A is the area of triangle (i,j,k). ^R90 represent a rotation of
  410. // 90 degrees
  411. //)igl_Qu8mg5v7";
  412. const char *__doc_igl_harmonic = R"igl_Qu8mg5v7(// Compute k-harmonic weight functions "coordinates".
  413. //
  414. //
  415. // Inputs:
  416. // V #V by dim vertex positions
  417. // F #F by simplex-size list of element indices
  418. // b #b boundary indices into V
  419. // bc #b by #W list of boundary values
  420. // k power of harmonic operation (1: harmonic, 2: biharmonic, etc)
  421. // Outputs:
  422. // W #V by #W list of weights
  423. //)igl_Qu8mg5v7";
  424. const char *__doc_igl_internal_angles = R"igl_Qu8mg5v7(// Compute internal angles for a triangle mesh
  425. //
  426. // Inputs:
  427. // V #V by dim eigen Matrix of mesh vertex nD positions
  428. // F #F by poly-size eigen Matrix of face (triangle) indices
  429. // Output:
  430. // K #F by poly-size eigen Matrix of internal angles
  431. // for triangles, columns correspond to edges [1,2],[2,0],[0,1]
  432. //
  433. // Known Issues:
  434. // if poly-size ≠ 3 then dim must equal 3.)igl_Qu8mg5v7";
  435. const char *__doc_igl_invert_diag = R"igl_Qu8mg5v7(// Templates:
  436. // T should be a eigen sparse matrix primitive type like int or double
  437. // Inputs:
  438. // X an m by n sparse matrix
  439. // Outputs:
  440. // Y an m by n sparse matrix)igl_Qu8mg5v7";
  441. const char *__doc_igl_is_irregular_vertex = R"igl_Qu8mg5v7(// Determine if a vertex is irregular, i.e. it has more than 6 (triangles)
  442. // or 4 (quads) incident edges. Vertices on the boundary are ignored.
  443. //
  444. // Inputs:
  445. // V #V by dim list of vertex positions
  446. // F #F by 3[4] list of triangle[quads] indices
  447. // Returns #V vector of bools revealing whether vertices are singular
  448. //)igl_Qu8mg5v7";
  449. const char *__doc_igl_jet = R"igl_Qu8mg5v7(// JET like MATLAB's jet
  450. //
  451. // Inputs:
  452. // m number of colors
  453. // Outputs:
  454. // J m by list of RGB colors between 0 and 1
  455. //
  456. //#ifndef IGL_NO_EIGEN
  457. // void jet(const int m, Eigen::MatrixXd & J);
  458. //#endif
  459. // Wrapper for directly computing [r,g,b] values for a given factor f between
  460. // 0 and 1
  461. //
  462. // Inputs:
  463. // f factor determining color value as if 0 was min and 1 was max
  464. // Outputs:
  465. // r red value
  466. // g green value
  467. // b blue value)igl_Qu8mg5v7";
  468. const char *__doc_igl_local_basis = R"igl_Qu8mg5v7(// Compute a local orthogonal reference system for each triangle in the given mesh
  469. // Templates:
  470. // DerivedV derived from vertex positions matrix type: i.e. MatrixXd
  471. // DerivedF derived from face indices matrix type: i.e. MatrixXi
  472. // Inputs:
  473. // V eigen matrix #V by 3
  474. // F #F by 3 list of mesh faces (must be triangles)
  475. // Outputs:
  476. // B1 eigen matrix #F by 3, each vector is tangent to the triangle
  477. // B2 eigen matrix #F by 3, each vector is tangent to the triangle and perpendicular to B1
  478. // B3 eigen matrix #F by 3, normal of the triangle
  479. //
  480. // See also: adjacency_matrix)igl_Qu8mg5v7";
  481. const char *__doc_igl_lscm = R"igl_Qu8mg5v7(// Compute a Least-squares conformal map parametrization (equivalently
  482. // derived in "Intrinsic Parameterizations of Surface Meshes" [Desbrun et al.
  483. // 2002] and "Least Squares Conformal Maps for Automatic Texture Atlas
  484. // Generation" [Lévy et al. 2002]), though this implementation follows the
  485. // derivation in: "Spectral Conformal Parameterization" [Mullen et al. 2008]
  486. // (note, this does **not** implement the Eigen-decomposition based method in
  487. // [Mullen et al. 2008], which is not equivalent). Input should be a manifold
  488. // mesh (also no unreferenced vertices) and "boundary" (fixed vertices) `b`
  489. // should contain at least two vertices per connected component.
  490. //
  491. // Inputs:
  492. // V #V by 3 list of mesh vertex positions
  493. // F #F by 3 list of mesh faces (must be triangles)
  494. // b #b boundary indices into V
  495. // bc #b by 3 list of boundary values
  496. // Outputs:
  497. // UV #V by 2 list of 2D mesh vertex positions in UV space
  498. // Returns true only on solver success.
  499. //)igl_Qu8mg5v7";
  500. 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
  501. // unit circle with spacing proportional to the original boundary edge
  502. // lengths.
  503. //
  504. // Inputs:
  505. // V #V by dim list of mesh vertex positions
  506. // b #W list of vertex ids
  507. // Outputs:
  508. // UV #W by 2 list of 2D position on the unit circle for the vertices in b)igl_Qu8mg5v7";
  509. const char *__doc_igl_massmatrix = R"igl_Qu8mg5v7(// Constructs the mass (area) matrix for a given mesh (V,F).
  510. //
  511. // Templates:
  512. // DerivedV derived type of eigen matrix for V (e.g. derived from
  513. // MatrixXd)
  514. // DerivedF derived type of eigen matrix for F (e.g. derived from
  515. // MatrixXi)
  516. // Scalar scalar type for eigen sparse matrix (e.g. double)
  517. // Inputs:
  518. // V #V by dim list of mesh vertex positions
  519. // F #F by simplex_size list of mesh faces (must be triangles)
  520. // type one of the following ints:
  521. // MASSMATRIX_TYPE_BARYCENTRIC barycentric
  522. // MASSMATRIX_TYPE_VORONOI voronoi-hybrid {default}
  523. // MASSMATRIX_TYPE_FULL full {not implemented}
  524. // Outputs:
  525. // M #V by #V mass matrix
  526. //
  527. // See also: adjacency_matrix
  528. //)igl_Qu8mg5v7";
  529. const char *__doc_igl_min_quad_with_fixed_precompute = R"igl_Qu8mg5v7(// Known Bugs: rows of Aeq **should probably** be linearly independent.
  530. // During precomputation, the rows of a Aeq are checked via QR. But in case
  531. // they're not then resulting probably will no longer be sparse: it will be
  532. // slow.
  533. //
  534. // MIN_QUAD_WITH_FIXED Minimize quadratic energy
  535. //
  536. // 0.5*Z'*A*Z + Z'*B + C with
  537. //
  538. // constraints that Z(known) = Y, optionally also subject to the constraints
  539. // Aeq*Z = Beq
  540. //
  541. // Templates:
  542. // T should be a eigen matrix primitive type like int or double
  543. // Inputs:
  544. // A n by n matrix of quadratic coefficients
  545. // known list of indices to known rows in Z
  546. // Y list of fixed values corresponding to known rows in Z
  547. // Aeq m by n list of linear equality constraint coefficients
  548. // pd flag specifying whether A(unknown,unknown) is positive definite
  549. // Outputs:
  550. // data factorization struct with all necessary information to solve
  551. // using min_quad_with_fixed_solve
  552. // Returns true on success, false on error
  553. //
  554. // Benchmark: For a harmonic solve on a mesh with 325K facets, matlab 2.2
  555. // secs, igl/min_quad_with_fixed.h 7.1 secs
  556. //)igl_Qu8mg5v7";
  557. const char *__doc_igl_min_quad_with_fixed_solve = R"igl_Qu8mg5v7(// Solves a system previously factored using min_quad_with_fixed_precompute
  558. //
  559. // Template:
  560. // T type of sparse matrix (e.g. double)
  561. // DerivedY type of Y (e.g. derived from VectorXd or MatrixXd)
  562. // DerivedZ type of Z (e.g. derived from VectorXd or MatrixXd)
  563. // Inputs:
  564. // data factorization struct with all necessary precomputation to solve
  565. // B n by 1 column of linear coefficients
  566. // Y b by 1 list of constant fixed values
  567. // Beq m by 1 list of linear equality constraint constant values
  568. // Outputs:
  569. // Z n by cols solution
  570. // sol #unknowns+#lagrange by cols solution to linear system
  571. // Returns true on success, false on error)igl_Qu8mg5v7";
  572. const char *__doc_igl_min_quad_with_fixed = R"igl_Qu8mg5v7(See min_quad_with_fixed for the documentation.)igl_Qu8mg5v7";
  573. const char *__doc_igl_n_polyvector = R"igl_Qu8mg5v7(// Inputs:
  574. // v0, v1 the two #3 by 1 vectors
  575. // normalized boolean, if false, then the vectors are normalized prior to the calculation
  576. // Output:
  577. // 3 by 3 rotation matrix that takes v0 to v1
  578. //)igl_Qu8mg5v7";
  579. const char *__doc_igl_parula = R"igl_Qu8mg5v7(// PARULA like MATLAB's parula
  580. //
  581. // Inputs:
  582. // m number of colors
  583. // Outputs:
  584. // J m by list of RGB colors between 0 and 1
  585. //
  586. // Wrapper for directly computing [r,g,b] values for a given factor f between
  587. // 0 and 1
  588. //
  589. // Inputs:
  590. // f factor determining color value as if 0 was min and 1 was max
  591. // Outputs:
  592. // r red value
  593. // g green value
  594. // b blue value)igl_Qu8mg5v7";
  595. const char *__doc_igl_per_corner_normals = R"igl_Qu8mg5v7(// Compute vertex normals via vertex position list, face list
  596. // Inputs:
  597. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  598. // F #F by 3 eigne Matrix of face (triangle) indices
  599. // corner_threshold threshold in degrees on sharp angles
  600. // Output:
  601. // CN #F*3 by 3 eigen Matrix of mesh vertex 3D normals, where the normal
  602. // for corner F(i,j) is at CN(i*3+j,:) )igl_Qu8mg5v7";
  603. const char *__doc_igl_per_edge_normals = R"igl_Qu8mg5v7(// Compute face normals via vertex position list, face list
  604. // Inputs:
  605. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  606. // F #F by 3 eigen Matrix of face (triangle) indices
  607. // weight weighting type
  608. // FN #F by 3 matrix of 3D face normals per face
  609. // Output:
  610. // N #2 by 3 matrix of mesh edge 3D normals per row
  611. // E #E by 2 matrix of edge indices per row
  612. // EMAP #E by 1 matrix of indices from all edges to E
  613. //)igl_Qu8mg5v7";
  614. const char *__doc_igl_per_face_normals = R"igl_Qu8mg5v7(// Compute face normals via vertex position list, face list
  615. // Inputs:
  616. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  617. // F #F by 3 eigen Matrix of face (triangle) indices
  618. // Z 3 vector normal given to faces with degenerate normal.
  619. // Output:
  620. // N #F by 3 eigen Matrix of mesh face (triangle) 3D normals
  621. //
  622. // Example:
  623. // // Give degenerate faces (1/3,1/3,1/3)^0.5
  624. // per_face_normals(V,F,Vector3d(1,1,1).normalized(),N);)igl_Qu8mg5v7";
  625. const char *__doc_igl_per_face_normals_stable = R"igl_Qu8mg5v7(// Special version where order of face indices is guaranteed not to effect
  626. // output.)igl_Qu8mg5v7";
  627. const char *__doc_igl_per_vertex_normals = R"igl_Qu8mg5v7(// Compute vertex normals via vertex position list, face list
  628. // Inputs:
  629. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  630. // F #F by 3 eigne Matrix of face (triangle) indices
  631. // weighting Weighting type
  632. // Output:
  633. // N #V by 3 eigen Matrix of mesh vertex 3D normals)igl_Qu8mg5v7";
  634. const char *__doc_igl_planarize_quad_mesh = R"igl_Qu8mg5v7(// Inputs:
  635. // Vin #V by 3 eigen Matrix of mesh vertex 3D positions
  636. // F #F by 4 eigen Matrix of face (quad) indices
  637. // maxIter maximum numbers of iterations
  638. // threshold minimum allowed threshold for non-planarity
  639. // Output:
  640. // Vout #V by 3 eigen Matrix of planar mesh vertex 3D positions
  641. //)igl_Qu8mg5v7";
  642. const char *__doc_igl_png_readPNG = R"igl_Qu8mg5v7(// Read an image from a .png file into 4 memory buffers
  643. //
  644. // Input:
  645. // png_file path to .png file
  646. // Output:
  647. // R,G,B,A texture channels
  648. // Returns true on success, false on failure
  649. //)igl_Qu8mg5v7";
  650. const char *__doc_igl_png_writePNG = R"igl_Qu8mg5v7(// Writes an image to a png file
  651. //
  652. // Input:
  653. // R,G,B,A texture channels
  654. // Output:
  655. // png_file path to .png file
  656. // Returns true on success, false on failure
  657. //)igl_Qu8mg5v7";
  658. 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)
  659. //
  660. // Inputs:
  661. // P #P by 3 list of query point positions
  662. // V #V by 3 list of vertex positions
  663. // Ele #Ele by (3|2|1) list of (triangle|edge|point) indices
  664. // Outputs:
  665. // sqrD #P list of smallest squared distances
  666. // I #P list of primitive indices corresponding to smallest distances
  667. // C #P by 3 list of closest points
  668. //
  669. // Known bugs: This only computes distances to given primitivess. So
  670. // unreferenced vertices are ignored. However, degenerate primitives are
  671. // handled correctly: triangle [1 2 2] is treated as a segment [1 2], and
  672. // triangle [1 1 1] is treated as a point. So one _could_ add extra
  673. // combinatorially degenerate rows to Ele for all unreferenced vertices to
  674. // also get distances to points.)igl_Qu8mg5v7";
  675. const char *__doc_igl_polar_svd = R"igl_Qu8mg5v7(// Computes the polar decomposition (R,T) of a matrix A using SVD singular
  676. // value decomposition
  677. //
  678. // Inputs:
  679. // A 3 by 3 matrix to be decomposed
  680. // Outputs:
  681. // R 3 by 3 rotation matrix part of decomposition (**always rotataion**)
  682. // T 3 by 3 stretch matrix part of decomposition
  683. // U 3 by 3 left-singular vectors
  684. // S 3 by 1 singular values
  685. // V 3 by 3 right-singular vectors
  686. //
  687. //)igl_Qu8mg5v7";
  688. const char *__doc_igl_principal_curvature = R"igl_Qu8mg5v7(// Compute the principal curvature directions and magnitude of the given triangle mesh
  689. // DerivedV derived from vertex positions matrix type: i.e. MatrixXd
  690. // DerivedF derived from face indices matrix type: i.e. MatrixXi
  691. // Inputs:
  692. // V eigen matrix #V by 3
  693. // F #F by 3 list of mesh faces (must be triangles)
  694. // radius controls the size of the neighbourhood used, 1 = average edge lenght
  695. //
  696. // Outputs:
  697. // PD1 #V by 3 maximal curvature direction for each vertex.
  698. // PD2 #V by 3 minimal curvature direction for each vertex.
  699. // PV1 #V by 1 maximal curvature value for each vertex.
  700. // PV2 #V by 1 minimal curvature value for each vertex.
  701. //
  702. // See also: average_onto_faces, average_onto_vertices
  703. //
  704. // This function has been developed by: Nikolas De Giorgis, Luigi Rocca and Enrico Puppo.
  705. // The algorithm is based on:
  706. // Efficient Multi-scale Curvature and Crease Estimation
  707. // Daniele Panozzo, Enrico Puppo, Luigi Rocca
  708. // GraVisMa, 2010)igl_Qu8mg5v7";
  709. const char *__doc_igl_quad_planarity = R"igl_Qu8mg5v7(// Compute planarity of the faces of a quad mesh
  710. // Inputs:
  711. // V #V by 3 eigen Matrix of mesh vertex 3D positions
  712. // F #F by 4 eigen Matrix of face (quad) indices
  713. // Output:
  714. // P #F by 1 eigen Matrix of mesh face (quad) planarities
  715. //)igl_Qu8mg5v7";
  716. const char *__doc_igl_readDMAT = R"igl_Qu8mg5v7(See readDMAT for the documentation.)igl_Qu8mg5v7";
  717. const char *__doc_igl_readMESH = R"igl_Qu8mg5v7(// load a tetrahedral volume mesh from a .mesh file
  718. //
  719. // Templates:
  720. // Scalar type for positions and vectors (will be read as double and cast
  721. // to Scalar)
  722. // Index type for indices (will be read as int and cast to Index)
  723. // Input:
  724. // mesh_file_name path of .mesh file
  725. // Outputs:
  726. // V double matrix of vertex positions #V by 3
  727. // T #T list of tet indices into vertex positions
  728. // F #F list of face indices into vertex positions
  729. //
  730. // Known bugs: Holes and regions are not supported)igl_Qu8mg5v7";
  731. const char *__doc_igl_readOBJ = R"igl_Qu8mg5v7(// Read a mesh from an ascii obj file, filling in vertex positions, normals
  732. // and texture coordinates. Mesh may have faces of any number of degree
  733. //
  734. // Templates:
  735. // Scalar type for positions and vectors (will be read as double and cast
  736. // to Scalar)
  737. // Index type for indices (will be read as int and cast to Index)
  738. // Inputs:
  739. // str path to .obj file
  740. // Outputs:
  741. // V double matrix of vertex positions #V by 3
  742. // TC double matrix of texture coordinats #TC by 2
  743. // N double matrix of corner normals #N by 3
  744. // F #F list of face indices into vertex positions
  745. // FTC #F list of face indices into vertex texture coordinates
  746. // FN #F list of face indices into vertex normals
  747. // Returns true on success, false on errors)igl_Qu8mg5v7";
  748. const char *__doc_igl_readOFF = R"igl_Qu8mg5v7(// Read a mesh from an ascii obj file, filling in vertex positions, normals
  749. // and texture coordinates. Mesh may have faces of any number of degree
  750. //
  751. // Templates:
  752. // Scalar type for positions and vectors (will be read as double and cast
  753. // to Scalar)
  754. // Index type for indices (will be read as int and cast to Index)
  755. // Inputs:
  756. // str path to .obj file
  757. // Outputs:
  758. // V double matrix of vertex positions #V by 3
  759. // F #F list of face indices into vertex positions
  760. // TC double matrix of texture coordinats #TC by 2
  761. // FTC #F list of face indices into vertex texture coordinates
  762. // N double matrix of corner normals #N by 3
  763. // FN #F list of face indices into vertex normals
  764. // Returns true on success, false on errors)igl_Qu8mg5v7";
  765. const char *__doc_igl_read_triangle_mesh = R"igl_Qu8mg5v7(// read mesh from an ascii file with automatic detection of file format.
  766. // supported: obj, off, stl, wrl, ply, mesh)
  767. //
  768. // Templates:
  769. // Scalar type for positions and vectors (will be read as double and cast
  770. // to Scalar)
  771. // Index type for indices (will be read as int and cast to Index)
  772. // Inputs:
  773. // str path to file
  774. // Outputs:
  775. // V eigen double matrix #V by 3
  776. // F eigen int matrix #F by 3
  777. // Returns true iff success)igl_Qu8mg5v7";
  778. const char *__doc_igl_rotate_vectors = R"igl_Qu8mg5v7(// Rotate the vectors V by A radiants on the tangent plane spanned by B1 and
  779. // B2
  780. //
  781. // Inputs:
  782. // V #V by 3 eigen Matrix of vectors
  783. // A #V eigen vector of rotation angles or a single angle to be applied
  784. // to all vectors
  785. // B1 #V by 3 eigen Matrix of base vector 1
  786. // B2 #V by 3 eigen Matrix of base vector 2
  787. //
  788. // Output:
  789. // Returns the rotated vectors
  790. //)igl_Qu8mg5v7";
  791. const char *__doc_igl_setdiff = R"igl_Qu8mg5v7(// Set difference of elements of matrices
  792. //
  793. // Inputs:
  794. // A m-long vector of indices
  795. // B n-long vector of indices
  796. // Outputs:
  797. // C (k<=m)-long vector of unique elements appearing in A but not in B
  798. // IA (k<=m)-long list of indices into A so that C = A(IA)
  799. //)igl_Qu8mg5v7";
  800. const char *__doc_igl_signed_distance = R"igl_Qu8mg5v7(// Computes signed distance to a mesh
  801. //
  802. // Inputs:
  803. // P #P by 3 list of query point positions
  804. // V #V by 3 list of vertex positions
  805. // F #F by ss list of triangle indices, ss should be 3 unless sign_type ==
  806. // SIGNED_DISTANCE_TYPE_UNSIGNED
  807. // sign_type method for computing distance _sign_ S
  808. // Outputs:
  809. // S #P list of smallest signed distances
  810. // I #P list of facet indices corresponding to smallest distances
  811. // C #P by 3 list of closest points
  812. // N #P by 3 list of closest normals (only set if
  813. // sign_type=SIGNED_DISTANCE_TYPE_PSEUDONORMAL)
  814. //
  815. // Known bugs: This only computes distances to triangles. So unreferenced
  816. // vertices and degenerate triangles are ignored.)igl_Qu8mg5v7";
  817. const char *__doc_igl_signed_distance_pseudonormal = R"igl_Qu8mg5v7(// Computes signed distance to mesh
  818. //
  819. // Inputs:
  820. // tree AABB acceleration tree (see AABB.h)
  821. // F #F by 3 list of triangle indices
  822. // FN #F by 3 list of triangle normals
  823. // VN #V by 3 list of vertex normals (ANGLE WEIGHTING)
  824. // EN #E by 3 list of edge normals (UNIFORM WEIGHTING)
  825. // EMAP #F*3 mapping edges in F to E
  826. // q Query point
  827. // Returns signed distance to mesh
  828. //)igl_Qu8mg5v7";
  829. const char *__doc_igl_signed_distance_winding_number = R"igl_Qu8mg5v7(// Inputs:
  830. // tree AABB acceleration tree (see cgal/point_mesh_squared_distance.h)
  831. // hier Winding number evaluation hierarchy
  832. // q Query point
  833. // Returns signed distance to mesh)igl_Qu8mg5v7";
  834. const char *__doc_igl_slice = R"igl_Qu8mg5v7(// Act like the matlab X(row_indices,col_indices) operator, where
  835. // row_indices, col_indices are non-negative integer indices.
  836. //
  837. // Inputs:
  838. // X m by n matrix
  839. // R list of row indices
  840. // C list of column indices
  841. // Output:
  842. // Y #R by #C matrix
  843. //
  844. // See also: slice_mask)igl_Qu8mg5v7";
  845. const char *__doc_igl_slice_into = R"igl_Qu8mg5v7(// Act like the matlab Y(row_indices,col_indices) = X
  846. //
  847. // Inputs:
  848. // X xm by xn rhs matrix
  849. // R list of row indices
  850. // C list of column indices
  851. // Y ym by yn lhs matrix
  852. // Output:
  853. // Y ym by yn lhs matrix, same as input but Y(R,C) = X)igl_Qu8mg5v7";
  854. const char *__doc_igl_slice_mask = R"igl_Qu8mg5v7(// Act like the matlab X(row_mask,col_mask) operator, where
  855. // row_mask, col_mask are non-negative integer indices.
  856. //
  857. // Inputs:
  858. // X m by n matrix
  859. // R m list of row bools
  860. // C n list of column bools
  861. // Output:
  862. // Y #trues-in-R by #trues-in-C matrix
  863. //
  864. // See also: slice_mask)igl_Qu8mg5v7";
  865. const char *__doc_igl_slice_tets = R"igl_Qu8mg5v7(// SLICE_TETS Slice through a tet mesh (V,T) along a given plane (via its
  866. // implicit equation).
  867. //
  868. // Inputs:
  869. // V #V by 3 list of tet mesh vertices
  870. // T #T by 4 list of tet indices into V
  871. // plane list of 4 coefficients in the plane equation: [x y z 1]'*plane = 0
  872. // Optional:
  873. // 'Manifold' followed by whether to stitch together triangles into a
  874. // manifold mesh {true}: results in more compact U but slightly slower.
  875. // Outputs:
  876. // U #U by 3 list of triangle mesh vertices along slice
  877. // G #G by 3 list of triangles indices into U
  878. // J #G list of indices into T revealing from which tet each faces comes
  879. // BC #U by #V list of barycentric coordinates (or more generally: linear
  880. // interpolation coordinates) so that U = BC*V
  881. // )igl_Qu8mg5v7";
  882. const char *__doc_igl_sortrows = R"igl_Qu8mg5v7(// Act like matlab's [Y,I] = sortrows(X)
  883. //
  884. // Templates:
  885. // DerivedX derived scalar type, e.g. MatrixXi or MatrixXd
  886. // DerivedI derived integer type, e.g. MatrixXi
  887. // Inputs:
  888. // X m by n matrix whose entries are to be sorted
  889. // ascending sort ascending (true, matlab default) or descending (false)
  890. // Outputs:
  891. // Y m by n matrix whose entries are sorted (**should not** be same
  892. // reference as X)
  893. // I m list of indices so that
  894. // Y = X(I,:);)igl_Qu8mg5v7";
  895. const char *__doc_igl_triangle_triangulate = R"igl_Qu8mg5v7(// Triangulate the interior of a polygon using the triangle library.
  896. //
  897. // Inputs:
  898. // V #V by 2 list of 2D vertex positions
  899. // E #E by 2 list of vertex ids forming unoriented edges of the boundary of the polygon
  900. // H #H by 2 coordinates of points contained inside holes of the polygon
  901. // flags string of options pass to triangle (see triangle documentation)
  902. // Outputs:
  903. // V2 #V2 by 2 coordinates of the vertives of the generated triangulation
  904. // F2 #F2 by 3 list of indices forming the faces of the generated triangulation
  905. //
  906. // TODO: expose the option to prevent Steiner points on the boundary
  907. //)igl_Qu8mg5v7";
  908. const char *__doc_igl_unique = R"igl_Qu8mg5v7(// Act like matlab's [C,IA,IC] = unique(X)
  909. //
  910. // Templates:
  911. // T comparable type T
  912. // Inputs:
  913. // A #A vector of type T
  914. // Outputs:
  915. // C #C vector of unique entries in A
  916. // IA #C index vector so that C = A(IA);
  917. // IC #A index vector so that A = C(IC);)igl_Qu8mg5v7";
  918. const char *__doc_igl_unique_rows = R"igl_Qu8mg5v7(// Act like matlab's [C,IA,IC] = unique(X,'rows')
  919. //
  920. // Templates:
  921. // DerivedA derived scalar type, e.g. MatrixXi or MatrixXd
  922. // DerivedIA derived integer type, e.g. MatrixXi
  923. // DerivedIC derived integer type, e.g. MatrixXi
  924. // Inputs:
  925. // A m by n matrix whose entries are to unique'd according to rows
  926. // Outputs:
  927. // C #C vector of unique rows in A
  928. // IA #C index vector so that C = A(IA,:);
  929. // IC #A index vector so that A = C(IC,:);)igl_Qu8mg5v7";
  930. const char *__doc_igl_unproject_onto_mesh = R"igl_Qu8mg5v7(// Unproject a screen location (using current opengl viewport, projection, and
  931. // model view) to a 3D position _onto_ a given mesh, if the ray through the
  932. // given screen location (x,y) _hits_ the mesh.
  933. //
  934. // Inputs:
  935. // pos screen space coordinates
  936. // model model matrix
  937. // proj projection matrix
  938. // viewport vieweport vector
  939. // V #V by 3 list of mesh vertex positions
  940. // F #F by 3 list of mesh triangle indices into V
  941. // Outputs:
  942. // fid id of the first face hit
  943. // bc barycentric coordinates of hit
  944. // Returns true if there's a hit)igl_Qu8mg5v7";
  945. const char *__doc_igl_upsample = R"igl_Qu8mg5v7(// Subdivide a mesh without moving vertices: loop subdivision but odd
  946. // vertices stay put and even vertices are just edge midpoints
  947. //
  948. // Templates:
  949. // MatV matrix for vertex positions, e.g. MatrixXd
  950. // MatF matrix for vertex positions, e.g. MatrixXi
  951. // Inputs:
  952. // V #V by dim mesh vertices
  953. // F #F by 3 mesh triangles
  954. // Outputs:
  955. // NV new vertex positions, V is guaranteed to be at top
  956. // NF new list of face indices
  957. //
  958. // NOTE: V should not be the same as NV,
  959. // NOTE: F should not be the same as NF, use other proto
  960. //
  961. // Known issues:
  962. // - assumes (V,F) is edge-manifold.)igl_Qu8mg5v7";
  963. const char *__doc_igl_winding_number = R"igl_Qu8mg5v7(// WINDING_NUMBER Compute the sum of solid angles of a triangle/tetrahedron
  964. // described by points (vectors) V
  965. //
  966. // Templates:
  967. // dim dimension of input
  968. // Inputs:
  969. // V n by 3 list of vertex positions
  970. // F #F by 3 list of triangle indices, minimum index is 0
  971. // O no by 3 list of origin positions
  972. // Outputs:
  973. // S no by 1 list of winding numbers
  974. //
  975. // 3d)igl_Qu8mg5v7";
  976. const char *__doc_igl_winding_number_3 = R"igl_Qu8mg5v7(// Inputs:
  977. // V pointer to array containing #V by 3 vertex positions along rows,
  978. // given in column major order
  979. // n number of mesh vertices
  980. // F pointer to array containing #F by 3 face indices along rows,
  981. // given in column major order
  982. // m number of faces
  983. // O pointer to array containing #O by 3 query positions along rows,
  984. // given in column major order
  985. // no number of origins
  986. // Outputs:
  987. // S no by 1 list of winding numbers)igl_Qu8mg5v7";
  988. const char *__doc_igl_winding_number_2 = R"igl_Qu8mg5v7(//// Only one evaluation origin
  989. //template <typename DerivedF>
  990. //IGL_INLINE void winding_number_3(
  991. // const double * V,
  992. // const int n,
  993. // const DerivedF * F,
  994. // const int m,
  995. // const double * O,
  996. // double * S);
  997. // 2d)igl_Qu8mg5v7";
  998. const char *__doc_igl_writeMESH = R"igl_Qu8mg5v7(// save a tetrahedral volume mesh to a .mesh file
  999. //
  1000. // Templates:
  1001. // Scalar type for positions and vectors (will be cast as double)
  1002. // Index type for indices (will be cast to int)
  1003. // Input:
  1004. // mesh_file_name path of .mesh file
  1005. // V double matrix of vertex positions #V by 3
  1006. // T #T list of tet indices into vertex positions
  1007. // F #F list of face indices into vertex positions
  1008. //
  1009. // Known bugs: Holes and regions are not supported)igl_Qu8mg5v7";
  1010. const char *__doc_igl_writeOBJ = R"igl_Qu8mg5v7(// Write a mesh in an ascii obj file
  1011. // Inputs:
  1012. // str path to outputfile
  1013. // V #V by 3 mesh vertex positions
  1014. // F #F by 3|4 mesh indices into V
  1015. // CN #CN by 3 normal vectors
  1016. // FN #F by 3|4 corner normal indices into CN
  1017. // TC #TC by 2|3 texture coordinates
  1018. // FTC #F by 3|4 corner texture coord indices into TC
  1019. // Returns true on success, false on error)igl_Qu8mg5v7";