// This file is part of libigl, a simple c++ geometry processing library. // // Copyright (C) 2013 Alec Jacobson // // This Source Code Form is subject to the terms of the Mozilla Public License // v. 2.0. If a copy of the MPL was not distributed with this file, You can // obtain one at http://mozilla.org/MPL/2.0/. #ifndef IGL_COPYLEFT_TETGEN_TETRAHEDRALIZE_H #define IGL_COPYLEFT_TETGEN_TETRAHEDRALIZE_H #include "../../igl_inline.h" #include #include #include #ifndef TETLIBRARY #define TETLIBRARY #endif #include // Defined REAL namespace igl { namespace copyleft { namespace tetgen { /// Mesh the interior of a surface mesh (V,F) using tetgen /// /// @param[in] V #V by 3 vertex position list /// @param[in] F #F list of polygon face indices into V (0-indexed) /// @param[in] H #H by 3 list of seed points inside holes /// @param[in] VM #VM list of vertex markers /// @param[in] FM #FM list of face markers /// @param[in] R #R by 5 list of region attributes /// @param[in] switches string of tetgen options (See tetgen documentation) e.g. /// "pq1.414a0.01" tries to mesh the interior of a given surface with /// quality and area constraints /// "" will mesh the convex hull constrained to pass through V (ignores F) /// @param[out] TV #TV by 3 vertex position list /// @param[out] TT #TT by 4 list of tet face indices /// @param[out] TF #TF by 3 list of triangle face indices ('f', else /// `boundary_facets` is called on TT) /// @param[out] TR #TT list of region ID for each tetrahedron /// @param[out] TN #TT by 4 list of indices neighbors for each tetrahedron ('n') /// @param[out] PT #TV list of incident tetrahedron for a vertex ('m') /// @param[out] FT #TF by 2 list of tetrahedrons sharing a triface ('nn') /// @param[out] num_regions Number of regions in output mesh /// @return status: /// 0 success /// 1 tetgen threw exception /// 2 tetgen did not crash but could not create any tets (probably there are /// holes, duplicate faces etc.) /// -1 other error /// /// \note Tetgen mixes integer region ids in with other region data `attr /// = (int) in->regionlist[i + 3];`. So it's declared safe to use integer /// types for `TR` since this also assumes that there's a single tet /// attribute and that it's the region id. /// /// #### Example /// /// ```cpp /// Eigen::MatrixXd V; /// Eigen::MatrixXi F; /// … /// Eigen::VectorXi VM,FM; /// Eigen::MatrixXd H,R; /// Eigen::VectorXi TM,TR,PT; /// Eigen::MatrixXi FT,TN; /// int numRegions; /// tetrahedralize(V,F,H,VM,FM,R,switches,TV,TT,TF,TM,TR,TN,PT,FT,numRegions); /// ``` template < typename DerivedV, typename DerivedF, typename DerivedH, typename DerivedVM, typename DerivedFM, typename DerivedR, typename DerivedTV, typename DerivedTT, typename DerivedTF, typename DerivedTM, typename DerivedTR, typename DerivedTN, typename DerivedPT, typename DerivedFT> IGL_INLINE int tetrahedralize( const Eigen::MatrixBase& V, const Eigen::MatrixBase& F, const Eigen::MatrixBase& H, const Eigen::MatrixBase& VM, const Eigen::MatrixBase& FM, const Eigen::MatrixBase& R, const std::string switches, Eigen::PlainObjectBase& TV, Eigen::PlainObjectBase& TT, Eigen::PlainObjectBase& TF, Eigen::PlainObjectBase& TM, Eigen::PlainObjectBase& TR, Eigen::PlainObjectBase& TN, Eigen::PlainObjectBase& PT, Eigen::PlainObjectBase& FT, int & num_regions); /// \overload template < typename DerivedV, typename DerivedF, typename DerivedTV, typename DerivedTT, typename DerivedTF> IGL_INLINE int tetrahedralize( const Eigen::MatrixBase& V, const Eigen::MatrixBase& F, const std::string switches, Eigen::PlainObjectBase& TV, Eigen::PlainObjectBase& TT, Eigen::PlainObjectBase& TF); } } } #ifndef IGL_STATIC_LIBRARY # include "tetrahedralize.cpp" #endif #endif