7#include <unordered_map>
15#define POLYFEM_REMESHER_SCOPED_TIMER(name) polyfem::utils::Timer __polyfem_timer(Remesher::timings[name])
27 std::array<size_t, 2>, T,
33 std::array<size_t, 3>, T,
38 std::array<size_t, 4>, T,
77 const size_t num_vertices,
78 const Eigen::MatrixXi &
elements) = 0;
98 virtual int dim()
const = 0;
111 virtual Eigen::MatrixXi
edges()
const = 0;
125 virtual std::vector<int>
boundary_nodes(
const Eigen::VectorXi &vertex_to_basis)
const = 0;
147 virtual void set_fixed(
const std::vector<bool> &fixed) = 0;
158 void write_mesh(
const std::string &path)
const;
162 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> &time_integrator);
166 const Eigen::MatrixXd &quantities,
168 Eigen::MatrixXd &x_prevs,
169 Eigen::MatrixXd &v_prevs,
170 Eigen::MatrixXd &a_prevs);
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);
229 static std::unordered_map<std::string, utils::Timing>
timings;
main class that contains the polyfem solver and all its state
mesh::Obstacle obstacle
Obstacles used in collisions.
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
const Eigen::MatrixXd m_obstacle_displacements
Collision obstacles' displacements.
const Eigen::MatrixXd & obstacle_displacements() const
Get a reference to the collision obstacles' displacements.
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)
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)
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.
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.
const double starting_energy
Starting energy.
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.
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
virtual Eigen::MatrixXi boundary_edges() const =0
Exports boundary edges of the stored mesh.
GlobalProjectionCache global_projection_cache
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.
const Eigen::MatrixXd & obstacle_quantities() const
Get a reference to the collision obstacles' extra quantities.
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)
void write_mesh(const std::string &path) const
Writes a VTU mesh file.
virtual bool is_volume() const
Is the mesh a volumetric mesh.
virtual Eigen::MatrixXi boundary_faces() const =0
Exports boundary faces of the stored mesh.
void cache_before()
Cache quantities before applying an operation.
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.
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)
static void log_timings()
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.
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
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.
static std::unordered_map< std::string, utils::Timing > timings
Timings for the remeshing operations.
std::unordered_map< std::array< size_t, 4 >, T, polyfem::utils::HashUnorderedArray< size_t, 4 >, polyfem::utils::EqualUnorderedArray< size_t, 4 > > TetMap
virtual int dim() const =0
Dimension of the mesh.
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.
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).
Eigen::MatrixXi elements
Elements before an operation.
Eigen::MatrixXd rest_positions
Rest positions of the mesh before an operation.
Eigen::MatrixXd projection_quantities
dim rows per vertex and 1 column per quantity
Hash function for an array where the order does not matter.