PolyFEM
Loading...
Searching...
No Matches
Refinement.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <geogram/mesh/mesh.h>
4#include <Eigen/Dense>
5#include <vector>
6
7namespace polyfem
8{
9 namespace mesh
10 {
11
12 // Compute a graph (V,E) where V are the edges of the input quad mesh, and E
13 // connects edges from opposite sides of the input quads.
14 //
15 // @param[in] Q #Q x 4 input quads
16 // @param[out] edge_index Map (f, lv) -> edge index for edge (lv, lv+1)
17 // @param[out] adj Adjacency graph
18 // @param[out] pairs_of_edges { List of mesh edges, corresponding to the
19 // vertices of the output graph }
20 // @param[out] pairs_of_quads List of adjacent quads
21 // @param[out] quad_index { Map (f, lv) -> index of the quad across edge
22 // (lv, lv+1) }
23 //
24 void edge_adjacency_graph(const Eigen::MatrixXi &Q, Eigen::MatrixXi &edge_index,
25 std::vector<std::vector<int>> &adj,
26 std::vector<std::pair<int, int>> *pairs_of_edges = nullptr,
27 std::vector<std::pair<int, int>> *pairs_of_quads = nullptr,
28 Eigen::MatrixXi *quad_index = nullptr);
29
30 typedef std::function<void(const Eigen::MatrixXd &, Eigen::MatrixXd &, int)> EvalParametersFunc;
31 typedef std::function<std::tuple<int, int, bool>(int, int)> GetAdjacentLocalEdge;
32
33 // Instantiate a periodic 2D pattern (triangle-mesh) on a given quad mesh
34 //
35 // @param[in] IV #IV x 3 input quad mesh vertices
36 // @param[in] IF #IF x 4 input quad mesh facets
37 // @param[in] PV #PV x (2|3) input pattern vertices in [0,1]^2
38 // @param[in] PF #PF x (3|4) input pattern facets
39 // @param[out] OV #OV x 3 output mesh vertices
40 // @param[out] OF #OF x 3 output mesh facets
41 // @param[in] SF #OV x 1 matrix of input source quad index, filled if pointer is non-zero
42 // @param[in] evalFunc { Evaluate the uv param of the pattern into a 2d or 3d
43 // position }
44 //
45 // @return { Return true in case of success. }
46 //
48 const Eigen::MatrixXd &IV, const Eigen::MatrixXi &IF,
49 const Eigen::MatrixXd &PV, const Eigen::MatrixXi &PF,
50 Eigen::MatrixXd &OV, Eigen::MatrixXi &OF,
51 Eigen::VectorXi *SF = nullptr,
52 EvalParametersFunc evalFunc = nullptr,
53 GetAdjacentLocalEdge getAdjLocalEdge = nullptr);
54
55 //
56 // Refine a quad-mesh by splitting each quad into 4 quads.
57 //
58 // @param[in] IV #IV x 3 input quad mesh vertices
59 // @param[in] IF #IF x 4 input quad mesh facets
60 // @param[out] OV #OV x 3 output mesh vertices
61 // @param[out] OF #OF x 4 output mesh facets
62 //
63 void refine_quad_mesh(const Eigen::MatrixXd &IV, const Eigen::MatrixXi &IF,
64 Eigen::MatrixXd &OV, Eigen::MatrixXi &OF);
65
66 namespace Polygons
67 {
68
69 typedef std::function<void(Eigen::MatrixXd &IV, Eigen::MatrixXd &OV, std::vector<std::vector<int>> &OF)> SplitFunction;
70
84 void polar_split(const Eigen::MatrixXd &IV, Eigen::MatrixXd &OV, std::vector<std::vector<int>> &OF, double t = 0.5);
85
88 {
89 return [t](const Eigen::MatrixXd &IV, Eigen::MatrixXd &OV, std::vector<std::vector<int>> &OF) { polar_split(IV, OV, OF, t); };
90 }
91
104 void catmul_clark_split(const Eigen::MatrixXd &IV, Eigen::MatrixXd &OV, std::vector<std::vector<int>> &OF);
105
107
111 void no_split(const Eigen::MatrixXd &IV, Eigen::MatrixXd &OV, std::vector<std::vector<int>> &OF);
112
114 } // namespace Polygons
115
126 void refine_polygonal_mesh(const GEO::Mesh &M_in, GEO::Mesh &M_out, Polygons::SplitFunction split_func);
127
134 void refine_triangle_mesh(const GEO::Mesh &M_in, GEO::Mesh &M_out);
135
136 // Kept for compatibility
137 [[deprecated]] inline void refine_polygonal_mesh(const GEO::Mesh &M_in, GEO::Mesh &M_out, bool refine_polygons = false, double t = 0.5)
138 {
139 if (refine_polygons == false)
140 {
142 }
143 else
144 {
146 }
147 }
148 } // namespace mesh
149} // namespace polyfem
vector< list< int > > adj
void polar_split(const Eigen::MatrixXd &IV, Eigen::MatrixXd &OV, std::vector< std::vector< int > > &OF, double t=0.5)
Split a polygon using polar refinement.
void catmul_clark_split(const Eigen::MatrixXd &IV, Eigen::MatrixXd &OV, std::vector< std::vector< int > > &OF)
Split a polygon using polar refinement.
std::function< void(Eigen::MatrixXd &IV, Eigen::MatrixXd &OV, std::vector< std::vector< int > > &OF)> SplitFunction
SplitFunction catmul_clark_split_func()
SplitFunction no_split_func()
SplitFunction polar_split_func(double t)
Helper function.
void no_split(const Eigen::MatrixXd &IV, Eigen::MatrixXd &OV, std::vector< std::vector< int > > &OF)
Don't split polygons.
void refine_triangle_mesh(const GEO::Mesh &M_in, GEO::Mesh &M_out)
Refine a triangle mesh.
void refine_quad_mesh(const Eigen::MatrixXd &IV, const Eigen::MatrixXi &IF, Eigen::MatrixXd &OV, Eigen::MatrixXi &OF)
void edge_adjacency_graph(const Eigen::MatrixXi &Q, Eigen::MatrixXi &edge_index, std::vector< std::vector< int > > &adj, std::vector< std::pair< int, int > > *pairs_of_edges=nullptr, std::vector< std::pair< int, int > > *pairs_of_quads=nullptr, Eigen::MatrixXi *quad_index=nullptr)
bool instantiate_pattern(const Eigen::MatrixXd &IV, const Eigen::MatrixXi &IF, const Eigen::MatrixXd &PV, const Eigen::MatrixXi &PF, Eigen::MatrixXd &OV, Eigen::MatrixXi &OF, Eigen::VectorXi *SF=nullptr, EvalParametersFunc evalFunc=nullptr, GetAdjacentLocalEdge getAdjLocalEdge=nullptr)
std::function< void(const Eigen::MatrixXd &, Eigen::MatrixXd &, int)> EvalParametersFunc
void refine_polygonal_mesh(const GEO::Mesh &M_in, GEO::Mesh &M_out, Polygons::SplitFunction split_func)
Refine a polygonal mesh.
std::function< std::tuple< int, int, bool >(int, int)> GetAdjacentLocalEdge