PolyFEM
Loading...
Searching...
No Matches
UpsampleMesh.hpp
Go to the documentation of this file.
1#pragma once
2
7
8#include <Eigen/Core>
9
10namespace polyfem::mesh
11{
18 void stitch_mesh(
19 const Eigen::MatrixXd &V,
20 const Eigen::MatrixXi &F,
21 Eigen::MatrixXd &V_out,
22 Eigen::MatrixXi &F_out,
23 const double epsilon = 1e-5);
24
36 void stitch_mesh(
37 const Eigen::MatrixXd &V,
38 const Eigen::MatrixXi &F,
39 const std::vector<Eigen::Triplet<double>> &W,
40 Eigen::MatrixXd &V_out,
41 Eigen::MatrixXi &F_out,
42 std::vector<Eigen::Triplet<double>> &W_out,
43 const double epsilon = 1e-5);
44
48 double max_edge_length(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F);
49
50 // Regular tessellation
51
57 const int n, Eigen::MatrixXd &V, Eigen::MatrixXi &F);
58
66 const Eigen::MatrixXd &V,
67 const Eigen::MatrixXi &F,
68 const double max_edge_length,
69 Eigen::MatrixXd &V_out,
70 Eigen::MatrixXi &F_out);
71
72 // Irregular tessellation
73
79 Eigen::MatrixXd
80 refine_edge(const VectorNd &a, const VectorNd &b, const double max_edge_length);
81
90 const VectorNd &a,
91 const VectorNd &b,
92 const VectorNd &c,
93 const double max_edge_len,
94 Eigen::MatrixXd &V,
95 Eigen::MatrixXi &E);
96
106 const Eigen::Vector3d &a,
107 const Eigen::Vector3d &b,
108 const Eigen::Vector3d &c,
109 const double max_edge_length,
110 Eigen::MatrixXd &UV,
111 Eigen::MatrixXi &F);
112
122 const Eigen::Vector3d &a,
123 const Eigen::Vector3d &b,
124 const Eigen::Vector3d &c,
125 const double max_edge_length,
126 Eigen::MatrixXd &V,
127 Eigen::MatrixXi &F);
128
137 const Eigen::MatrixXd &V,
138 const Eigen::MatrixXi &F,
139 const double max_edge_length,
140 Eigen::MatrixXd &V_out,
141 Eigen::MatrixXi &F_out);
142} // namespace polyfem::mesh
int V
Eigen::MatrixXd refine_edge(const VectorNd &a, const VectorNd &b, const double max_edge_length)
Refine an edge (a, b) so each refined edge has length at most max_edge_length.
void irregular_triangle(const Eigen::Vector3d &a, const Eigen::Vector3d &b, const Eigen::Vector3d &c, const double max_edge_length, Eigen::MatrixXd &V, Eigen::MatrixXi &F)
Refine a triangle (a, b, c) into a well shaped triangle mesh.
void irregular_tessellation(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const double max_edge_length, Eigen::MatrixXd &V_out, Eigen::MatrixXi &F_out)
Tessilate a triangle mesh (V, F) with well shaped triangles.
double max_edge_length(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
Compute the maximum edge length of a triangle mesh (V, F)
void regular_grid_triangle_barycentric_coordinates(const int n, Eigen::MatrixXd &V, Eigen::MatrixXi &F)
Compute the barycentric coordinates of a regular grid of triangles.
void irregular_triangle_barycentric_coordinates(const Eigen::Vector3d &a, const Eigen::Vector3d &b, const Eigen::Vector3d &c, const double max_edge_length, Eigen::MatrixXd &UV, Eigen::MatrixXi &F)
Refine a triangle (a, b, c) into a well shaped triangle mesh.
void regular_grid_tessellation(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const double out_max_edge_length, Eigen::MatrixXd &V_out, Eigen::MatrixXi &F_out)
Tessilate a triangle mesh (V, F) with regular grids of triangles of maximum edge length.
void refine_triangle_edges(const VectorNd &a, const VectorNd &b, const VectorNd &c, const double max_edge_len, Eigen::MatrixXd &V, Eigen::MatrixXi &E)
Refine the edges of a triangle (a, b, c) so each refined edge has length at most max_edge_length.
void stitch_mesh(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, Eigen::MatrixXd &V_out, Eigen::MatrixXi &F_out, const double epsilon)
Stitch a triangle mesh (V, F) together by removing duplicate vertices.
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:11