PolyFEM
Loading...
Searching...
No Matches
Remesher.cpp
Go to the documentation of this file.
1#include "Remesher.hpp"
2
4
14#include <paraviewo/VTUWriter.hpp>
15
16#include <igl/boundary_facets.h>
17#include <igl/edges.h>
18
19namespace polyfem::mesh
20{
22 const Eigen::MatrixXd &obstacle_displacements,
23 const Eigen::MatrixXd &obstacle_quantities,
24 const double current_time,
25 const double starting_energy)
26 : state(state),
27 args(state.args["space"]["remesh"]),
28 m_obstacle_displacements(obstacle_displacements),
29 m_obstacle_quantities(obstacle_quantities),
30 current_time(current_time),
31 starting_energy(starting_energy)
32 {
33 }
34
36 const Eigen::MatrixXd &rest_positions,
37 const Eigen::MatrixXd &positions,
38 const Eigen::MatrixXi &elements,
39 const Eigen::MatrixXd &projection_quantities,
40 const BoundaryMap<int> &boundary_to_id,
41 const std::vector<int> &body_ids,
42 const EdgeMap<double> &elastic_energy,
43 const EdgeMap<double> &contact_energy)
44 {
45 assert(elements.size() > 0);
46
47 // --------------------------------------------------------------------
48 // Determine which vertices are on the boundary
49
50 // Partition mesh by body_ids
51 assert(body_ids.size() == elements.rows());
52 std::unordered_map<int, std::vector<int>> body_elements;
53 for (int i = 0; i < elements.rows(); ++i)
54 {
55 const int body_id = body_ids[i];
56 if (body_elements.find(body_id) == body_elements.end())
57 body_elements[body_id] = std::vector<int>();
58 body_elements[body_id].push_back(i);
59 }
60
61 // Determine boundary vertices
62 std::vector<bool> is_boundary_vertex(positions.rows(), false);
63 for (const auto &[body_id, rows] : body_elements)
64 {
65 Eigen::MatrixXi boundary_facets;
66 igl::boundary_facets(elements(rows, Eigen::all), boundary_facets);
67
68 for (int i = 0; i < boundary_facets.rows(); ++i)
69 for (int j = 0; j < boundary_facets.cols(); ++j)
70 is_boundary_vertex[boundary_facets(i, j)] = true;
71 }
72
73 // --------------------------------------------------------------------
74
75 // Initialize the mesh atributes and connectivity
77
78 // Save the vertex position in the vertex attributes
82 set_fixed(is_boundary_vertex);
83 set_boundary_ids(boundary_to_id);
85 }
86
93
95 {
96 POLYFEM_REMESHER_SCOPED_TIMER("Project quantities");
97
98 using namespace polyfem::assembler;
99 using namespace polyfem::basis;
100 using namespace polyfem::utils;
101
102 const std::unique_ptr<Mesh> from_mesh = Mesh::create(
104 std::vector<ElementBases> from_bases;
105 std::vector<LocalBoundary> _;
106 Eigen::VectorXi from_vertex_to_basis;
107 int n_from_basis = build_bases(*from_mesh, state.formulation(), from_bases, _, from_vertex_to_basis);
108 assert(from_vertex_to_basis.size() == global_projection_cache.rest_positions.rows());
109
110 // Old values of independent variables
111 Eigen::MatrixXd from_projection_quantities = reorder_matrix(
113 from_vertex_to_basis, n_from_basis, dim());
114 append_rows(from_projection_quantities, obstacle_quantities());
115 n_from_basis += obstacle().n_vertices();
116 assert(dim() * n_from_basis == from_projection_quantities.rows());
117
118 // --------------------------------------------------------------------
119
120 Eigen::MatrixXd rest_positions = this->rest_positions();
121 Eigen::MatrixXi elements = this->elements();
122
123 const std::unique_ptr<Mesh> to_mesh = Mesh::create(rest_positions, elements);
124 std::vector<ElementBases> to_bases;
125 Eigen::VectorXi to_vertex_to_basis;
126 int n_to_basis = build_bases(*to_mesh, state.formulation(), to_bases, _, to_vertex_to_basis);
127 assert(to_vertex_to_basis.size() == rest_positions.rows());
128
129 rest_positions = reorder_matrix(rest_positions, to_vertex_to_basis, n_to_basis);
130 elements = map_index_matrix(elements, to_vertex_to_basis);
131
132 // Interpolated values of independent variables
133 Eigen::MatrixXd to_projection_quantities = reorder_matrix(
134 projection_quantities(), to_vertex_to_basis, n_to_basis, dim());
135 append_rows(to_projection_quantities, obstacle_quantities());
136 n_to_basis += obstacle().n_vertices();
137 assert(dim() * n_to_basis == to_projection_quantities.rows());
138
139 // --------------------------------------------------------------------
140
141 // solve M x = A y for x where M is the mass matrix and A is the cross mass matrix.
142 Eigen::SparseMatrix<double> M, A;
143 {
144 MassMatrixAssembler assembler;
145 assembler::Mass mass_matrix_assembler;
147
148 mass_matrix_assembler.assemble(
149 is_volume(),
150 n_to_basis, to_bases, to_bases,
151 cache, 0, M, true);
152 assert(M.rows() == to_projection_quantities.rows());
153
154 // if (lump_mass_matrix)
155 // M = lump_matrix(M);
156
157 assembler.assemble_cross(
158 is_volume(), dim(),
159 n_from_basis, from_bases, from_bases,
160 n_to_basis, to_bases, to_bases,
161 cache, A);
162 assert(A.rows() == to_projection_quantities.rows());
163 assert(A.cols() == from_projection_quantities.rows());
164 }
165
166 // --------------------------------------------------------------------
167
168 ipc::CollisionMesh collision_mesh;
169 {
170 Eigen::MatrixXi boundary_edges, boundary_faces;
171 if (dim() == 2)
172 {
173 igl::boundary_facets(elements, boundary_edges);
174 }
175 else
176 {
177 igl::boundary_facets(elements, boundary_faces);
178 igl::edges(boundary_faces, boundary_edges);
179 }
180
181 if (obstacle().n_edges() > 0)
183 if (obstacle().n_faces() > 0)
185
187 assert(rest_positions.size() == to_projection_quantities.rows());
188
189 collision_mesh = ipc::CollisionMesh::build_from_full_mesh(
191 }
192
193 // --------------------------------------------------------------------
194
195 Eigen::MatrixXd projected_quantities(to_projection_quantities.rows(), n_quantities());
196 const int n_constrained_quantaties = n_quantities() / 3;
197 const int n_unconstrained_quantaties = n_quantities() - n_constrained_quantaties;
198
199 const std::vector<int> boundary_nodes = this->boundary_nodes(to_vertex_to_basis);
200 for (int i = 0; i < n_constrained_quantaties; ++i)
201 {
202 projected_quantities.col(i) = constrained_L2_projection(
203 state.make_nl_solver(/*for_al=*/false),
204 // L2 projection form
205 M, A, /*y=*/from_projection_quantities.col(i),
206 // Inversion-free form
208 // Contact form
209 collision_mesh, state.args["contact"]["dhat"],
211 ? state.solve_data.contact_form->barrier_stiffness()
212 : 1.0,
213 state.args["contact"]["use_convergent_formulation"],
214 state.args["solver"]["contact"]["CCD"]["broad_phase"],
215 state.args["solver"]["contact"]["CCD"]["tolerance"],
216 state.args["solver"]["contact"]["CCD"]["max_iterations"],
217 // Augmented lagrangian form
218 boundary_nodes, obstacle().ndof(), to_projection_quantities.col(i),
219 // Initial guess
220 to_projection_quantities.col(i));
221 }
222
223 // Set entry for obstacle to identity
224 // ┌ ┐ ┌ ┐ ┌ ┐
225 // │M 0│ │x_fem│ |y_fem|
226 // │ │ │ │ = | | ⟹ x_obs = y_obs
227 // │0 I│ │x_obs│ |y_obs|
228 // └ ┘ └ ┘ └ ┘
229
230 for (int i = 0; i < dim() * obstacle().n_vertices(); ++i)
231 {
232 const int row = M.rows() - i - 1; // fill from bottom
233 M.coeffRef(row, row) = 1.0;
234 }
235
236 for (int i = 0; i < dim() * obstacle().n_vertices(); ++i)
237 {
238 const int row = A.rows() - i - 1; // fill from bottom
239 const int col = A.cols() - i - 1; // fill from bottom
240 A.coeffRef(row, col) = 1.0;
241 }
242
243 // NOTE: no need for to_projection_quantities.rightCols(n_unconstrained_quantaties)
244 projected_quantities.rightCols(n_unconstrained_quantaties) = unconstrained_L2_projection(
245 M, A, from_projection_quantities.rightCols(n_unconstrained_quantaties));
246
247 // --------------------------------------------------------------------
248
249 assert(projected_quantities.rows() == dim() * n_to_basis);
251 projected_quantities, to_vertex_to_basis, to_vertex_to_basis.size(), dim()));
252 }
253
255 const Mesh &mesh,
256 const std::string &assembler_formulation,
257 std::vector<polyfem::basis::ElementBases> &bases,
258 std::vector<LocalBoundary> &local_boundary,
259 Eigen::VectorXi &vertex_to_basis)
260 {
261 using namespace polyfem::basis;
262
263 int n_bases;
264 std::map<int, basis::InterfaceData> poly_edge_to_data;
265 std::shared_ptr<mesh::MeshNodes> mesh_nodes;
266 if (mesh.dimension() == 2)
267 {
269 dynamic_cast<const Mesh2D &>(mesh),
270 assembler_formulation, /*quadrature_order=*/1,
271 /*mass_quadrature_order=*/2, /*discr_order=*/1,
272 /*serendipity=*/false, /*has_polys=*/false,
273 /*is_geom_bases=*/false, bases, local_boundary,
274 poly_edge_to_data, mesh_nodes);
275 }
276 else
277 {
278 assert(mesh.dimension() == 3);
280 dynamic_cast<const Mesh3D &>(mesh),
281 assembler_formulation, /*quadrature_order=*/1,
282 /*mass_quadrature_order=*/2, /*discr_order=*/1,
283 /*serendipity=*/false, /*has_polys=*/false,
284 /*is_geom_bases=*/false, bases, local_boundary,
285 poly_edge_to_data, mesh_nodes);
286 }
287
288 // TODO: use mesh_nodes to build vertex_to_basis
289 vertex_to_basis.setConstant(mesh.n_vertices(), -1);
290 for (const ElementBases &elm : bases)
291 {
292 for (const Basis &basis : elm.bases)
293 {
294 assert(basis.global().size() == 1);
295 const int basis_id = basis.global()[0].index;
296 const RowVectorNd v = basis.global()[0].node;
297
298 for (int i = 0; i < mesh.n_vertices(); i++)
299 {
300 if ((mesh.point(i).array() == v.array()).all())
301 {
302 if (vertex_to_basis[i] == -1)
303 vertex_to_basis[i] = basis_id;
304 assert(vertex_to_basis[i] == basis_id);
305 break;
306 }
307 }
308 }
309 }
310
311 return n_bases;
312 }
313
314 void Remesher::write_mesh(const std::string &path) const
315 {
316 assert(utils::StringUtils::endswith(path, ".vtu"));
317
318 paraviewo::VTUWriter writer;
319
320 Eigen::MatrixXd rest_positions = this->rest_positions();
321 Eigen::MatrixXd displacements = this->displacements();
322
323 const Eigen::MatrixXd projection_quantities = this->projection_quantities();
324 std::vector<Eigen::MatrixXd> unflattened_projection_quantities;
325 for (const auto &q : projection_quantities.colwise())
326 unflattened_projection_quantities.push_back(utils::unflatten(q, dim()));
327
328 Eigen::MatrixXi elements = this->elements();
329 std::vector<std::vector<int>> all_elements(elements.rows(), std::vector<int>(elements.cols()));
330 for (int i = 0; i < elements.rows(); i++)
331 for (int j = 0; j < elements.cols(); j++)
332 all_elements[i][j] = elements(i, j);
333
334 const int offset = rest_positions.rows();
335 if (obstacle().n_vertices() > 0)
336 {
339 for (int i = 0; i < unflattened_projection_quantities.size(); ++i)
341 unflattened_projection_quantities[i],
343
344 if (obstacle().n_edges() > 0)
345 {
346 const Eigen::MatrixXi obstacle_edges = obstacle().e().array() + offset;
347 all_elements.resize(all_elements.size() + obstacle_edges.rows());
348 for (int i = 0; i < obstacle().n_edges(); i++)
349 {
350 all_elements[i + elements.rows()] = std::vector<int>(obstacle_edges.cols());
351 for (int j = 0; j < obstacle_edges.cols(); j++)
352 all_elements[i + elements.rows()][j] = obstacle_edges(i, j);
353 }
354 }
355
356 if (obstacle().n_faces() > 0)
357 {
358 const Eigen::MatrixXi obstacle_faces = obstacle().f().array() + offset;
359 all_elements.resize(all_elements.size() + obstacle_faces.rows());
360 for (int i = 0; i < obstacle().n_edges(); i++)
361 {
362 all_elements[i + elements.rows()] = std::vector<int>(obstacle_faces.cols());
363 for (int j = 0; j < obstacle_faces.cols(); j++)
364 all_elements[i + elements.rows()][j] = obstacle_faces(i, j);
365 }
366 }
367 }
368
369 for (int i = 0; i < unflattened_projection_quantities.size(); ++i)
370 writer.add_field(fmt::format("projection_quantity_{:d}", i), unflattened_projection_quantities[i]);
371 writer.add_field("displacement", displacements);
372 writer.write_mesh(path, rest_positions, all_elements, /*is_simplicial=*/true, /*has_poly=*/false);
373 }
374
376 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> &time_integrator)
377 {
378 if (time_integrator == nullptr)
379 return Eigen::MatrixXd();
380
381 // not including current displacement as this will be handled as positions
382 Eigen::MatrixXd projection_quantities(
383 time_integrator->x_prev().size(), 3 * time_integrator->steps());
384 int i = 0;
385 for (const Eigen::VectorXd &x : time_integrator->x_prevs())
386 projection_quantities.col(i++) = x;
387 for (const Eigen::VectorXd &v : time_integrator->v_prevs())
388 projection_quantities.col(i++) = v;
389 for (const Eigen::VectorXd &a : time_integrator->a_prevs())
390 projection_quantities.col(i++) = a;
391 assert(i == projection_quantities.cols());
392
394 }
395
397 const Eigen::MatrixXd &quantities,
398 const int dim,
399 Eigen::MatrixXd &x_prevs,
400 Eigen::MatrixXd &v_prevs,
401 Eigen::MatrixXd &a_prevs)
402 {
403 if (quantities.size() == 0)
404 return;
405
406 assert(quantities.cols() % 3 == 0);
407 const int n_steps = quantities.cols() / 3;
408
409 x_prevs = quantities.leftCols(n_steps);
410 v_prevs = quantities.middleCols(n_steps, n_steps);
411 a_prevs = quantities.rightCols(n_steps);
412 }
413
415 {
416 if (!logger().should_log(spdlog::level::debug) || timings.empty())
417 return;
418
419 std::cout << "--------------------------------------------------------------------------------" << std::endl;
420
421 logger().debug("Total time: {:.3g}s", total_time);
422 double sum = 0;
423
424 // Copy key-value pair from Map to vector of pairs
425 std::vector<std::pair<std::string, utils::Timing>> sorted_timings(timings.begin(), timings.end());
426 // Sort timings by decreasing time
427 std::sort(sorted_timings.begin(), sorted_timings.end(), [](const auto &a, const auto &b) { return a.second > b.second; });
428 for (const auto &[name, time] : sorted_timings)
429 {
430 logger().debug("{:62s}: {:10.3g}s {:5.1f}% ({:6d} calls)", name, time, time / total_time * 100, time.count);
431 sum += time;
432 }
433
434 // logger().debug("Miscellaneous: {:.3g}s {:.1f}%", total_time - sum, (total_time - sum) / total_time * 100);
435 if (num_solves > 0)
436 logger().debug("Avg. # DOF per solve: {}", total_ndofs / double(num_solves));
437
438 std::cout << "--------------------------------------------------------------------------------" << std::endl;
439 }
440
441 // Static members must be initialized in the source file:
443 double Remesher::total_time = 0;
444 size_t Remesher::num_solves = 0;
445 size_t Remesher::total_ndofs = 0;
446
447} // namespace polyfem::mesh
std::unique_ptr< MatrixCache > cache
Definition Assembler.cpp:21
#define POLYFEM_REMESHER_SCOPED_TIMER(name)
Definition Remesher.hpp:15
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
std::string formulation() const
return the formulation (checks if the problem is scalar or not and deals with multiphysics)
Definition State.cpp:328
std::shared_ptr< polysolve::nonlinear::Solver > make_nl_solver(bool for_al) const
factory to create the nl solver depending on input
json args
main input arguments containing all defaults
Definition State.hpp:101
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:324
Caches basis evaluation and geometric mapping at every element.
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const override
computes and returns local stiffness matrix (1x1) for bases i,j (where i,j is passed in through data)...
Definition Mass.cpp:5
void assemble_cross(const bool is_volume, const int size, const int n_from_basis, const std::vector< basis::ElementBases > &from_bases, const std::vector< basis::ElementBases > &from_gbases, const int n_to_basis, const std::vector< basis::ElementBases > &to_bases, const std::vector< basis::ElementBases > &to_gbases, const AssemblyValsCache &cache, StiffnessMatrix &mass) const
Assembles the cross mass matrix between to function spaces.
Represents one basis function and its gradient.
Definition Basis.hpp:44
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d).
static int build_bases(const mesh::Mesh2D &mesh, const std::string &assembler, const int quadrature_order, const int mass_quadrature_order, const int discr_order, const bool serendipity, const bool has_polys, const bool is_geom_bases, std::vector< ElementBases > &bases, std::vector< mesh::LocalBoundary > &local_boundary, std::map< int, InterfaceData > &poly_edge_to_data, std::shared_ptr< mesh::MeshNodes > &mesh_nodes)
Builds FE basis functions over the entire mesh (P1, P2 over triangles, Q1, Q2 over quads).
static int build_bases(const mesh::Mesh3D &mesh, const std::string &assembler, const int quadrature_order, const int mass_quadrature_order, const int discr_order, const bool serendipity, const bool has_polys, const bool is_geom_bases, std::vector< ElementBases > &bases, std::vector< mesh::LocalBoundary > &local_boundary, std::map< int, InterfaceData > &poly_face_to_data, std::shared_ptr< mesh::MeshNodes > &mesh_nodes)
Builds FE basis functions over the entire mesh (P1, P2 over tets, Q1, Q2 over hes).
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
virtual int n_vertices() const =0
number of vertices
virtual RowVectorNd point(const int global_index) const =0
point coordinates
static std::unique_ptr< Mesh > create(const std::string &path, const bool non_conforming=false)
factory to build the proper mesh
Definition Mesh.cpp:173
int dimension() const
utily for dimension
Definition Mesh.hpp:151
const Eigen::MatrixXi & e() const
Definition Obstacle.hpp:44
const Eigen::MatrixXi & f() const
Definition Obstacle.hpp:43
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.
virtual void set_projection_quantities(const Eigen::MatrixXd &projection_quantities)=0
Set projected quantities of the stored mesh.
virtual std::vector< int > boundary_nodes(const Eigen::VectorXi &vertex_to_basis) const =0
Get the boundary nodes of the stored mesh.
Remesher(const State &state, const Eigen::MatrixXd &obstacle_displacements, const Eigen::MatrixXd &obstacle_quantities, const double current_time, const double starting_energy)
Construct a new Remesher object.
Definition Remesher.cpp:21
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 elements() const =0
Exports elements of the stored mesh.
const State & state
Reference to the simulation state.
Definition Remesher.hpp:191
static void log_timings()
Definition Remesher.cpp:414
virtual void set_positions(const Eigen::MatrixXd &positions)=0
Set deformed positions of the stored mesh.
virtual void set_boundary_ids(const BoundaryMap< int > &boundary_to_id)=0
Set the boundary IDs of all edges.
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
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 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)
std::shared_ptr< solver::ContactForm > contact_form
Used for test only.
Eigen::MatrixXd unconstrained_L2_projection(const Eigen::SparseMatrix< double > &M, const Eigen::SparseMatrix< double > &A, const Eigen::Ref< const Eigen::MatrixXd > &y)
Eigen::VectorXd constrained_L2_projection(std::shared_ptr< polysolve::nonlinear::Solver > nl_solver, const Eigen::SparseMatrix< double > &M, const Eigen::SparseMatrix< double > &A, const Eigen::VectorXd &y, const Eigen::MatrixXd &rest_positions, const Eigen::MatrixXi &elements, const int dim, const ipc::CollisionMesh &collision_mesh, const double dhat, const double barrier_stiffness, const bool use_convergent_formulation, const ipc::BroadPhaseMethod broad_phase_method, const double ccd_tolerance, const int ccd_max_iterations, const std::vector< int > &boundary_nodes, const size_t obstacle_ndof, const Eigen::VectorXd &target_x, const Eigen::VectorXd &x0)
bool endswith(const std::string &str, const std::string &suffix)
Eigen::MatrixXd reorder_matrix(const Eigen::MatrixXd &in, const Eigen::VectorXi &in_to_out, int out_blocks=-1, const int block_size=1)
Reorder row blocks in a matrix.
Eigen::MatrixXd unreorder_matrix(const Eigen::MatrixXd &out, const Eigen::VectorXi &in_to_out, int in_blocks=-1, const int block_size=1)
Undo the reordering of row blocks in a matrix.
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
void append_rows(DstMat &dst, const SrcMat &src)
Eigen::MatrixXi map_index_matrix(const Eigen::MatrixXi &in, const Eigen::VectorXi &index_mapping)
Map the entrys of an index matrix to new indices.
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
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