PolyFEM
Loading...
Searching...
No Matches
Remesher.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/State.hpp>
6
7#include <unordered_map>
8#include <variant>
9
11{
13} // namespace polyfem::time_integrator
14
15#define POLYFEM_REMESHER_SCOPED_TIMER(name) polyfem::utils::Timer __polyfem_timer(Remesher::timings[name])
16
17namespace polyfem::mesh
18{
20 {
21 // --------------------------------------------------------------------
22 // typedefs
23 public:
25 template <typename T>
26 using EdgeMap = std::unordered_map<
27 std::array<size_t, 2>, T,
31 template <typename T>
32 using FaceMap = std::unordered_map<
33 std::array<size_t, 3>, T,
36 template <typename T>
37 using TetMap = std::unordered_map<
38 std::array<size_t, 4>, T,
41 template <typename T>
42 using BoundaryMap = std::variant<EdgeMap<T>, FaceMap<T>>;
43
44 // --------------------------------------------------------------------
45 // constructors
46 public:
49 Remesher(const State &state,
50 const Eigen::MatrixXd &obstacle_displacements,
51 const Eigen::MatrixXd &obstacle_quantities,
52 const double current_time,
53 const double starting_energy);
54
55 virtual ~Remesher() = default;
56
64 virtual void init(
65 const Eigen::MatrixXd &rest_positions,
66 const Eigen::MatrixXd &positions,
67 const Eigen::MatrixXi &elements,
68 const Eigen::MatrixXd &projection_quantities,
69 const BoundaryMap<int> &boundary_to_id,
70 const std::vector<int> &body_ids,
71 const EdgeMap<double> &elastic_energy,
72 const EdgeMap<double> &contact_energy);
73
74 protected:
77 const size_t num_vertices,
78 const Eigen::MatrixXi &elements) = 0;
79
80 // --------------------------------------------------------------------
81 // main functions
82 public:
85 virtual bool execute() = 0;
86
87 protected:
89 void project_quantities();
90
92 void cache_before();
93
94 // --------------------------------------------------------------------
95 // getters
96 public:
98 virtual int dim() const = 0;
99
102 virtual bool is_volume() const { return dim() == 3; }
103
105 virtual Eigen::MatrixXd rest_positions() const = 0;
107 virtual Eigen::MatrixXd displacements() const = 0;
109 virtual Eigen::MatrixXd positions() const = 0;
111 virtual Eigen::MatrixXi edges() const = 0;
113 virtual Eigen::MatrixXi elements() const = 0;
115 virtual Eigen::MatrixXi boundary_edges() const = 0;
117 virtual Eigen::MatrixXi boundary_faces() const = 0;
119 virtual Eigen::MatrixXd projection_quantities() const = 0;
121 virtual BoundaryMap<int> boundary_ids() const = 0;
123 virtual std::vector<int> body_ids() const = 0;
125 virtual std::vector<int> boundary_nodes(const Eigen::VectorXi &vertex_to_basis) const = 0;
126
128 virtual int n_quantities() const = 0;
129
131 const Obstacle &obstacle() const { return state.obstacle; }
133 const Eigen::MatrixXd &obstacle_displacements() const { return m_obstacle_displacements; }
135 const Eigen::MatrixXd &obstacle_quantities() const { return m_obstacle_quantities; }
136
137 // --------------------------------------------------------------------
138 // setters
139 public:
141 virtual void set_rest_positions(const Eigen::MatrixXd &positions) = 0;
143 virtual void set_positions(const Eigen::MatrixXd &positions) = 0;
145 virtual void set_projection_quantities(const Eigen::MatrixXd &projection_quantities) = 0;
147 virtual void set_fixed(const std::vector<bool> &fixed) = 0;
149 virtual void set_boundary_ids(const BoundaryMap<int> &boundary_to_id) = 0;
151 virtual void set_body_ids(const std::vector<int> &body_ids) = 0;
152
153 // --------------------------------------------------------------------
154 // utilities
155 public:
158 void write_mesh(const std::string &path) const;
159
161 static Eigen::MatrixXd combine_time_integrator_quantities(
162 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> &time_integrator);
163
166 const Eigen::MatrixXd &quantities,
167 const int dim,
168 Eigen::MatrixXd &x_prevs,
169 Eigen::MatrixXd &v_prevs,
170 Eigen::MatrixXd &a_prevs);
171
175 void init_assembler(const std::vector<int> &body_ids) const;
176
183 static int build_bases(
184 const Mesh &mesh,
185 const std::string &assembler_formulation,
186 std::vector<polyfem::basis::ElementBases> &bases,
187 std::vector<LocalBoundary> &local_boundary,
188 Eigen::VectorXi &vertex_to_basis);
189
191 const State &state;
192
193 // --------------------------------------------------------------------
194 // members
195 public:
197
198 protected:
199 // TODO: Drop this and only use a local EdgeOperationCache
201 {
203 Eigen::MatrixXd rest_positions;
205 Eigen::MatrixXi elements;
207 Eigen::MatrixXd projection_quantities;
208 };
209
211
213 const json args;
215 const Eigen::MatrixXd m_obstacle_displacements;
217 Eigen::MatrixXd m_obstacle_quantities;
219 const double current_time;
221 const double starting_energy;
222
223 // --------------------------------------------------------------------
224 // statistics
225 public:
226 static void log_timings();
227
229 static std::unordered_map<std::string, utils::Timing> timings;
230 static double total_time; // = 0;
231 static size_t num_solves; // = 0;
232 static size_t total_ndofs; // = 0;
233 };
234
235} // namespace polyfem::mesh
main class that contains the polyfem solver and all its state
Definition State.hpp:79
mesh::Obstacle obstacle
Obstacles used in collisions.
Definition State.hpp:468
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
const Eigen::MatrixXd m_obstacle_displacements
Collision obstacles' displacements.
Definition Remesher.hpp:215
const Eigen::MatrixXd & obstacle_displacements() const
Get a reference to the collision obstacles' displacements.
Definition Remesher.hpp:133
static int build_bases(const Mesh &mesh, const std::string &assembler_formulation, std::vector< polyfem::basis::ElementBases > &bases, std::vector< LocalBoundary > &local_boundary, Eigen::VectorXi &vertex_to_basis)
Build bases for a given mesh (V, F)
Definition Remesher.cpp:254
std::unordered_map< std::array< size_t, 2 >, T, polyfem::utils::HashUnorderedArray< size_t, 2 >, polyfem::utils::EqualUnorderedArray< size_t, 2 > > EdgeMap
Map from a (sorted) edge to an integer (ID)
Definition Remesher.hpp:29
virtual void set_body_ids(const std::vector< int > &body_ids)=0
Set the body IDs of all elements.
virtual void init_attributes_and_connectivity(const size_t num_vertices, const Eigen::MatrixXi &elements)=0
Create an internal mesh representation and associate attributes.
const double current_time
Current time.
Definition Remesher.hpp:219
virtual void set_projection_quantities(const Eigen::MatrixXd &projection_quantities)=0
Set projected quantities of the stored mesh.
virtual bool execute()=0
Execute the remeshing.
Eigen::MatrixXd m_obstacle_quantities
Collision obstacles' extra quantities.
Definition Remesher.hpp:217
const double starting_energy
Starting energy.
Definition Remesher.hpp:221
virtual std::vector< int > boundary_nodes(const Eigen::VectorXi &vertex_to_basis) const =0
Get the boundary nodes of the stored mesh.
const Obstacle & obstacle() const
Get a reference to the collision obstacles.
Definition Remesher.hpp:131
virtual Eigen::MatrixXd positions() const =0
Exports displacements of the stored mesh.
virtual Eigen::MatrixXd rest_positions() const =0
Exports rest positions of the stored mesh.
std::variant< EdgeMap< T >, FaceMap< T > > BoundaryMap
Definition Remesher.hpp:42
virtual Eigen::MatrixXi boundary_edges() const =0
Exports boundary edges of the stored mesh.
GlobalProjectionCache global_projection_cache
Definition Remesher.hpp:210
static double total_time
Definition Remesher.hpp:230
virtual void set_fixed(const std::vector< bool > &fixed)=0
Set if a vertex is fixed.
void project_quantities()
Update the mesh positions and other projection quantities.
Definition Remesher.cpp:94
const Eigen::MatrixXd & obstacle_quantities() const
Get a reference to the collision obstacles' extra quantities.
Definition Remesher.hpp:135
static Eigen::MatrixXd combine_time_integrator_quantities(const std::shared_ptr< time_integrator::ImplicitTimeIntegrator > &time_integrator)
Combine the quantities of a time integrator into a single matrix (one column per quantity)
Definition Remesher.cpp:375
void write_mesh(const std::string &path) const
Writes a VTU mesh file.
Definition Remesher.cpp:314
virtual bool is_volume() const
Is the mesh a volumetric mesh.
Definition Remesher.hpp:102
virtual Eigen::MatrixXi boundary_faces() const =0
Exports boundary faces of the stored mesh.
void cache_before()
Cache quantities before applying an operation.
Definition Remesher.cpp:87
virtual std::vector< int > body_ids() const =0
Exports body ids of the stored mesh.
virtual Eigen::MatrixXd displacements() const =0
Exports positions of the stored mesh.
virtual Eigen::MatrixXi edges() const =0
Exports edges of the stored mesh.
virtual Eigen::MatrixXi elements() const =0
Exports elements of the stored mesh.
const State & state
Reference to the simulation state.
Definition Remesher.hpp:191
std::unordered_map< std::array< size_t, 3 >, T, polyfem::utils::HashUnorderedArray< size_t, 3 >, polyfem::utils::EqualUnorderedArray< size_t, 3 > > FaceMap
Map from a (sorted) edge to an integer (ID)
Definition Remesher.hpp:35
static void log_timings()
Definition Remesher.cpp:414
virtual BoundaryMap< int > boundary_ids() const =0
Exports boundary ids of the stored mesh.
virtual void set_positions(const Eigen::MatrixXd &positions)=0
Set deformed positions of the stored mesh.
const json args
Copy of remesh args.
Definition Remesher.hpp:213
virtual void set_boundary_ids(const BoundaryMap< int > &boundary_to_id)=0
Set the boundary IDs of all edges.
void init_assembler(const std::vector< int > &body_ids) const
Create an assembler object.
static size_t total_ndofs
Definition Remesher.hpp:232
virtual void set_rest_positions(const Eigen::MatrixXd &positions)=0
Set rest positions of the stored mesh.
static void split_time_integrator_quantities(const Eigen::MatrixXd &quantities, const int dim, Eigen::MatrixXd &x_prevs, Eigen::MatrixXd &v_prevs, Eigen::MatrixXd &a_prevs)
Split the quantities of a time integrator into separate vectors.
Definition Remesher.cpp:396
static std::unordered_map< std::string, utils::Timing > timings
Timings for the remeshing operations.
Definition Remesher.hpp:229
std::unordered_map< std::array< size_t, 4 >, T, polyfem::utils::HashUnorderedArray< size_t, 4 >, polyfem::utils::EqualUnorderedArray< size_t, 4 > > TetMap
Definition Remesher.hpp:40
virtual int dim() const =0
Dimension of the mesh.
static size_t num_solves
Definition Remesher.hpp:231
virtual Eigen::MatrixXd projection_quantities() const =0
Exports projected quantities of the stored mesh.
virtual ~Remesher()=default
virtual void init(const Eigen::MatrixXd &rest_positions, const Eigen::MatrixXd &positions, const Eigen::MatrixXi &elements, const Eigen::MatrixXd &projection_quantities, const BoundaryMap< int > &boundary_to_id, const std::vector< int > &body_ids, const EdgeMap< double > &elastic_energy, const EdgeMap< double > &contact_energy)
Initialize the mesh.
Definition Remesher.cpp:35
virtual int n_quantities() const =0
Number of projection quantities (not including the position)
Implicit time integrator of a second order ODE (equivently a system of coupled first order ODEs).
nlohmann::json json
Definition Common.hpp:9
Eigen::MatrixXi elements
Elements before an operation.
Definition Remesher.hpp:205
Eigen::MatrixXd rest_positions
Rest positions of the mesh before an operation.
Definition Remesher.hpp:203
Eigen::MatrixXd projection_quantities
dim rows per vertex and 1 column per quantity
Definition Remesher.hpp:207
Hash function for an array where the order does not matter.
Definition HashUtils.hpp:35