Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
State.cpp
Go to the documentation of this file.
1#include <polyfem/State.hpp>
2#include <polyfem/Common.hpp>
3
7
10
19
22
24
27
30
33
36
41
44
45#include <polysolve/linear/FEMSolver.hpp>
46
48
49#include <igl/edges.h>
50#include <igl/Timer.h>
51
52#include <iostream>
53#include <algorithm>
54#include <memory>
55#include <filesystem>
56
58
60
61using namespace Eigen;
62
63namespace polyfem
64{
65 using namespace assembler;
66 using namespace mesh;
67 using namespace io;
68 using namespace utils;
69
70 namespace
71 {
73 void build_in_node_to_in_primitive(const Mesh &mesh, const MeshNodes &mesh_nodes,
74 Eigen::VectorXi &in_node_to_in_primitive,
75 Eigen::VectorXi &in_node_offset)
76 {
77 const int num_vertex_nodes = mesh_nodes.num_vertex_nodes();
78 const int num_edge_nodes = mesh_nodes.num_edge_nodes();
79 const int num_face_nodes = mesh_nodes.num_face_nodes();
80 const int num_cell_nodes = mesh_nodes.num_cell_nodes();
81
82 const int num_nodes = num_vertex_nodes + num_edge_nodes + num_face_nodes + num_cell_nodes;
83
84 const long n_vertices = num_vertex_nodes;
85 const int num_in_primitives = n_vertices + mesh.n_edges() + mesh.n_faces() + mesh.n_cells();
86 const int num_primitives = mesh.n_vertices() + mesh.n_edges() + mesh.n_faces() + mesh.n_cells();
87
88 in_node_to_in_primitive.resize(num_nodes);
89 in_node_offset.resize(num_nodes);
90
91 // Only one node per vertex, so this is an identity map.
92 in_node_to_in_primitive.head(num_vertex_nodes).setLinSpaced(num_vertex_nodes, 0, num_vertex_nodes - 1); // vertex nodes
93 in_node_offset.head(num_vertex_nodes).setZero();
94
95 int prim_offset = n_vertices;
96 int node_offset = num_vertex_nodes;
97 auto foo = [&](const int num_prims, const int num_prim_nodes) {
98 if (num_prims <= 0 || num_prim_nodes <= 0)
99 return;
100 const Eigen::VectorXi range = Eigen::VectorXi::LinSpaced(num_prim_nodes, 0, num_prim_nodes - 1);
101 // TODO: This assumes isotropic degree of element.
102 const int node_per_prim = num_prim_nodes / num_prims;
103
104 in_node_to_in_primitive.segment(node_offset, num_prim_nodes) =
105 range.array() / node_per_prim + prim_offset;
106
107 in_node_offset.segment(node_offset, num_prim_nodes) =
108 range.unaryExpr([&](const int x) { return x % node_per_prim; });
109
110 prim_offset += num_prims;
111 node_offset += num_prim_nodes;
112 };
113
114 foo(mesh.n_edges(), num_edge_nodes);
115 foo(mesh.n_faces(), num_face_nodes);
116 foo(mesh.n_cells(), num_cell_nodes);
117 }
118
119 bool build_in_primitive_to_primitive(
120 const Mesh &mesh, const MeshNodes &mesh_nodes,
121 const Eigen::VectorXi &in_ordered_vertices,
122 const Eigen::MatrixXi &in_ordered_edges,
123 const Eigen::MatrixXi &in_ordered_faces,
124 Eigen::VectorXi &in_primitive_to_primitive)
125 {
126 // NOTE: Assume in_cells_to_cells is identity
127 const int num_vertex_nodes = mesh_nodes.num_vertex_nodes();
128 const int num_edge_nodes = mesh_nodes.num_edge_nodes();
129 const int num_face_nodes = mesh_nodes.num_face_nodes();
130 const int num_cell_nodes = mesh_nodes.num_cell_nodes();
131 const int num_nodes = num_vertex_nodes + num_edge_nodes + num_face_nodes + num_cell_nodes;
132
133 const long n_vertices = num_vertex_nodes;
134 const int num_in_primitives = n_vertices + mesh.n_edges() + mesh.n_faces() + mesh.n_cells();
135 const int num_primitives = mesh.n_vertices() + mesh.n_edges() + mesh.n_faces() + mesh.n_cells();
136
137 in_primitive_to_primitive.setLinSpaced(num_in_primitives, 0, num_in_primitives - 1);
138
139 igl::Timer timer;
140
141 // ------------
142 // Map vertices
143 // ------------
144
145 if (in_ordered_vertices.rows() != n_vertices)
146 {
147 logger().warn("Node ordering disabled, in_ordered_vertices != n_vertices, {} != {}", in_ordered_vertices.rows(), n_vertices);
148 return false;
149 }
150
151 in_primitive_to_primitive.head(n_vertices) = in_ordered_vertices;
152
153 int in_offset = n_vertices;
154 int offset = mesh.n_vertices();
155
156 // ---------
157 // Map edges
158 // ---------
159
160 logger().trace("Building Mesh edges to IDs...");
161 timer.start();
162 const auto edges_to_ids = mesh.edges_to_ids();
163 if (in_ordered_edges.rows() != edges_to_ids.size())
164 {
165 logger().warn("Node ordering disabled, in_ordered_edges != edges_to_ids, {} != {}", in_ordered_edges.rows(), edges_to_ids.size());
166 return false;
167 }
168 timer.stop();
169 logger().trace("Done (took {}s)", timer.getElapsedTime());
170
171 logger().trace("Building in-edge to edge mapping...");
172 timer.start();
173 for (int in_ei = 0; in_ei < in_ordered_edges.rows(); in_ei++)
174 {
175 const std::pair<int, int> in_edge(
176 in_ordered_edges.row(in_ei).minCoeff(),
177 in_ordered_edges.row(in_ei).maxCoeff());
178 in_primitive_to_primitive[in_offset + in_ei] =
179 offset + edges_to_ids.at(in_edge); // offset edge ids
180 }
181 timer.stop();
182 logger().trace("Done (took {}s)", timer.getElapsedTime());
183
184 in_offset += mesh.n_edges();
185 offset += mesh.n_edges();
186
187 // ---------
188 // Map faces
189 // ---------
190
191 if (mesh.is_volume())
192 {
193 logger().trace("Building Mesh faces to IDs...");
194 timer.start();
195 const auto faces_to_ids = mesh.faces_to_ids();
196 if (in_ordered_faces.rows() != faces_to_ids.size())
197 {
198 logger().warn("Node ordering disabled, in_ordered_faces != faces_to_ids, {} != {}", in_ordered_faces.rows(), faces_to_ids.size());
199 return false;
200 }
201 timer.stop();
202 logger().trace("Done (took {}s)", timer.getElapsedTime());
203
204 logger().trace("Building in-face to face mapping...");
205 timer.start();
206 for (int in_fi = 0; in_fi < in_ordered_faces.rows(); in_fi++)
207 {
208 std::vector<int> in_face(in_ordered_faces.cols());
209 for (int i = 0; i < in_face.size(); i++)
210 in_face[i] = in_ordered_faces(in_fi, i);
211 std::sort(in_face.begin(), in_face.end());
212
213 in_primitive_to_primitive[in_offset + in_fi] =
214 offset + faces_to_ids.at(in_face); // offset face ids
215 }
216 timer.stop();
217 logger().trace("Done (took {}s)", timer.getElapsedTime());
218
219 in_offset += mesh.n_faces();
220 offset += mesh.n_faces();
221 }
222
223 return true;
224 }
225 } // namespace
226
227 std::vector<int> State::primitive_to_node() const
228 {
229 auto indices = iso_parametric() ? mesh_nodes->primitive_to_node() : geom_mesh_nodes->primitive_to_node();
230 indices.resize(mesh->n_vertices());
231 return indices;
232 }
233
234 std::vector<int> State::node_to_primitive() const
235 {
236 auto p2n = primitive_to_node();
237 std::vector<int> indices;
238 indices.resize(n_geom_bases);
239 for (int i = 0; i < p2n.size(); i++)
240 indices[p2n[i]] = i;
241 return indices;
242 }
243
245 {
246 if (args["space"]["basis_type"] == "Spline")
247 {
248 logger().warn("Node ordering disabled, it dosent work for splines!");
249 return;
250 }
251
252 if (disc_orders.maxCoeff() >= 4 || disc_orders.maxCoeff() != disc_orders.minCoeff())
253 {
254 logger().warn("Node ordering disabled, it works only for p < 4 and uniform order!");
255 return;
256 }
257
258 if (!mesh->is_conforming())
259 {
260 logger().warn("Node ordering disabled, not supported for non-conforming meshes!");
261 return;
262 }
263
264 if (mesh->has_poly())
265 {
266 logger().warn("Node ordering disabled, not supported for polygonal meshes!");
267 return;
268 }
269
270 if (mesh->in_ordered_vertices().size() <= 0 || mesh->in_ordered_edges().size() <= 0 || (mesh->is_volume() && mesh->in_ordered_faces().size() <= 0))
271 {
272 logger().warn("Node ordering disabled, input vertices/edges/faces not computed!");
273 return;
274 }
275
276 const int num_vertex_nodes = mesh_nodes->num_vertex_nodes();
277 const int num_edge_nodes = mesh_nodes->num_edge_nodes();
278 const int num_face_nodes = mesh_nodes->num_face_nodes();
279 const int num_cell_nodes = mesh_nodes->num_cell_nodes();
280
281 const int num_nodes = num_vertex_nodes + num_edge_nodes + num_face_nodes + num_cell_nodes;
282
283 const long n_vertices = num_vertex_nodes;
284 const int num_in_primitives = n_vertices + mesh->n_edges() + mesh->n_faces() + mesh->n_cells();
285 const int num_primitives = mesh->n_vertices() + mesh->n_edges() + mesh->n_faces() + mesh->n_cells();
286
287 igl::Timer timer;
288
289 logger().trace("Building in-node to in-primitive mapping...");
290 timer.start();
291 Eigen::VectorXi in_node_to_in_primitive;
292 Eigen::VectorXi in_node_offset;
293 build_in_node_to_in_primitive(*mesh, *mesh_nodes, in_node_to_in_primitive, in_node_offset);
294 timer.stop();
295 logger().trace("Done (took {}s)", timer.getElapsedTime());
296
297 logger().trace("Building in-primitive to primitive mapping...");
298 timer.start();
299 bool ok = build_in_primitive_to_primitive(
300 *mesh, *mesh_nodes,
301 mesh->in_ordered_vertices(),
302 mesh->in_ordered_edges(),
303 mesh->in_ordered_faces(),
305 timer.stop();
306 logger().trace("Done (took {}s)", timer.getElapsedTime());
307
308 if (!ok)
309 {
310 in_node_to_node.resize(0);
312 return;
313 }
314 const auto &tmp = mesh_nodes->in_ordered_vertices();
315 int max_tmp = -1;
316 for (auto v : tmp)
317 max_tmp = std::max(max_tmp, v);
318
319 in_node_to_node.resize(max_tmp + 1);
320 for (int i = 0; i < tmp.size(); ++i)
321 {
322 if (tmp[i] >= 0)
323 in_node_to_node[tmp[i]] = i;
324 }
325 }
326
327 std::string State::formulation() const
328 {
329 if (args["materials"].is_null())
330 {
331 logger().error("specify some 'materials'");
332 assert(!args["materials"].is_null());
333 throw std::runtime_error("invalid input");
334 }
335
336 if (args["materials"].is_array())
337 {
338 std::string current = "";
339 for (const auto &m : args["materials"])
340 {
341 const std::string tmp = m["type"];
342 if (current.empty())
343 current = tmp;
344 else if (current != tmp)
345 {
347 {
349 current = "MultiModels";
350 else
351 {
352 logger().error("Current material is {}, new material is {}, multimaterial supported only for LinearElasticity and NeoHookean", current, tmp);
353 throw std::runtime_error("invalid input");
354 }
355 }
356 else
357 {
358 logger().error("Current material is {}, new material is {}, multimaterial supported only for LinearElasticity and NeoHookean", current, tmp);
359 throw std::runtime_error("invalid input");
360 }
361 }
362 }
363
364 return current;
365 }
366 else
367 return args["materials"]["type"];
368 }
369
370 void State::sol_to_pressure(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
371 {
372 if (n_pressure_bases <= 0)
373 {
374 logger().error("No pressure bases defined!");
375 return;
376 }
377
378 assert(mixed_assembler != nullptr);
379 Eigen::MatrixXd tmp = sol;
380
381 int fluid_offset = use_avg_pressure ? (assembler->is_fluid() ? 1 : 0) : 0;
382 sol = tmp.topRows(tmp.rows() - n_pressure_bases - fluid_offset);
383 assert(sol.size() == n_bases * (problem->is_scalar() ? 1 : mesh->dimension()));
384 pressure = tmp.middleRows(tmp.rows() - n_pressure_bases - fluid_offset, n_pressure_bases);
385 assert(pressure.size() == n_pressure_bases);
386 }
387
389 const Mesh3D &mesh,
390 const int n_bases,
391 const std::vector<basis::ElementBases> &bases,
392 const std::vector<basis::ElementBases> &gbases,
393 Eigen::MatrixXd &basis_integrals)
394 {
395 if (!mesh.is_volume())
396 {
397 logger().error("Works only on volumetric meshes!");
398 return;
399 }
400 assert(mesh.is_volume());
401
402 basis_integrals.resize(n_bases, 9);
403 basis_integrals.setZero();
404 Eigen::MatrixXd rhs(n_bases, 9);
405 rhs.setZero();
406
407 const int n_elements = mesh.n_elements();
408 for (int e = 0; e < n_elements; ++e)
409 {
410 // if (mesh.is_polytope(e)) {
411 // continue;
412 // }
413 // ElementAssemblyValues vals = values[e];
414 // const ElementAssemblyValues &gvals = gvalues[e];
416 vals.compute(e, mesh.is_volume(), bases[e], gbases[e]);
417
418 // Computes the discretized integral of the PDE over the element
419 const int n_local_bases = int(vals.basis_values.size());
420 for (int j = 0; j < n_local_bases; ++j)
421 {
422 const AssemblyValues &v = vals.basis_values[j];
423 const double integral_100 = (v.grad_t_m.col(0).array() * vals.det.array() * vals.quadrature.weights.array()).sum();
424 const double integral_010 = (v.grad_t_m.col(1).array() * vals.det.array() * vals.quadrature.weights.array()).sum();
425 const double integral_001 = (v.grad_t_m.col(2).array() * vals.det.array() * vals.quadrature.weights.array()).sum();
426
427 const double integral_110 = ((vals.val.col(1).array() * v.grad_t_m.col(0).array() + vals.val.col(0).array() * v.grad_t_m.col(1).array()) * vals.det.array() * vals.quadrature.weights.array()).sum();
428 const double integral_011 = ((vals.val.col(2).array() * v.grad_t_m.col(1).array() + vals.val.col(1).array() * v.grad_t_m.col(2).array()) * vals.det.array() * vals.quadrature.weights.array()).sum();
429 const double integral_101 = ((vals.val.col(0).array() * v.grad_t_m.col(2).array() + vals.val.col(2).array() * v.grad_t_m.col(0).array()) * vals.det.array() * vals.quadrature.weights.array()).sum();
430
431 const double integral_200 = 2 * (vals.val.col(0).array() * v.grad_t_m.col(0).array() * vals.det.array() * vals.quadrature.weights.array()).sum();
432 const double integral_020 = 2 * (vals.val.col(1).array() * v.grad_t_m.col(1).array() * vals.det.array() * vals.quadrature.weights.array()).sum();
433 const double integral_002 = 2 * (vals.val.col(2).array() * v.grad_t_m.col(2).array() * vals.det.array() * vals.quadrature.weights.array()).sum();
434
435 const double area = (v.val.array() * vals.det.array() * vals.quadrature.weights.array()).sum();
436
437 for (size_t ii = 0; ii < v.global.size(); ++ii)
438 {
439 basis_integrals(v.global[ii].index, 0) += integral_100 * v.global[ii].val;
440 basis_integrals(v.global[ii].index, 1) += integral_010 * v.global[ii].val;
441 basis_integrals(v.global[ii].index, 2) += integral_001 * v.global[ii].val;
442
443 basis_integrals(v.global[ii].index, 3) += integral_110 * v.global[ii].val;
444 basis_integrals(v.global[ii].index, 4) += integral_011 * v.global[ii].val;
445 basis_integrals(v.global[ii].index, 5) += integral_101 * v.global[ii].val;
446
447 basis_integrals(v.global[ii].index, 6) += integral_200 * v.global[ii].val;
448 basis_integrals(v.global[ii].index, 7) += integral_020 * v.global[ii].val;
449 basis_integrals(v.global[ii].index, 8) += integral_002 * v.global[ii].val;
450
451 rhs(v.global[ii].index, 6) += -2.0 * area * v.global[ii].val;
452 rhs(v.global[ii].index, 7) += -2.0 * area * v.global[ii].val;
453 rhs(v.global[ii].index, 8) += -2.0 * area * v.global[ii].val;
454 }
455 }
456 }
457
458 basis_integrals -= rhs;
459 }
460
462 {
463 if (mesh->has_poly())
464 return true;
465
466 if (args["space"]["basis_type"] == "Bernstein")
467 return false;
468
469 if (args["space"]["basis_type"] == "Spline")
470 return true;
471
472 if (mesh->is_rational())
473 return false;
474
475 if (args["space"]["use_p_ref"])
476 return false;
477
478 if (has_periodic_bc())
479 return false;
480
482 return false;
483
484 if (mesh->orders().size() <= 0)
485 {
486 if (args["space"]["discr_order"] == 1)
487 return true;
488 else
489 return args["space"]["advanced"]["isoparametric"];
490 }
491
492 if (mesh->orders().minCoeff() != mesh->orders().maxCoeff())
493 return false;
494
495 if (args["space"]["discr_order"] == mesh->orders().minCoeff())
496 return true;
497
498 // TODO:
499 // if (args["space"]["discr_order"] == 1 && args["force_linear_geometry"])
500 // return true;
501
502 return args["space"]["advanced"]["isoparametric"];
503 }
504
506 {
507 if (!mesh)
508 {
509 logger().error("Load the mesh first!");
510 return;
511 }
512
513 mesh->prepare_mesh();
514
515 bases.clear();
516 pressure_bases.clear();
517 geom_bases_.clear();
518 boundary_nodes.clear();
519 dirichlet_nodes.clear();
520 neumann_nodes.clear();
521 local_boundary.clear();
522 total_local_boundary.clear();
525 local_pressure_cavity.clear();
526 polys.clear();
527 poly_edge_to_data.clear();
528 rhs.resize(0, 0);
529 basis_nodes_to_gbasis_nodes.resize(0, 0);
530
531 if (assembler::MultiModel *mm = dynamic_cast<assembler::MultiModel *>(assembler.get()))
532 {
533 assert(args["materials"].is_array());
534
535 std::vector<std::string> materials(mesh->n_elements());
536
537 std::map<int, std::string> mats;
538
539 for (const auto &m : args["materials"])
540 mats[m["id"].get<int>()] = m["type"];
541
542 for (int i = 0; i < materials.size(); ++i)
543 materials[i] = mats.at(mesh->get_body_id(i));
544
545 mm->init_multimodels(materials);
546 }
547
548 n_bases = 0;
549 n_geom_bases = 0;
551
552 stats.reset();
553
554 disc_orders.resize(mesh->n_elements());
555
556 problem->init(*mesh);
557 logger().info("Building {} basis...", (iso_parametric() ? "isoparametric" : "not isoparametric"));
558 const bool has_polys = mesh->has_poly();
559
560 local_boundary.clear();
563 local_pressure_cavity.clear();
564 std::map<int, basis::InterfaceData> poly_edge_to_data_geom; // temp dummy variable
565
566 const auto &tmp_json = args["space"]["discr_order"];
567 if (tmp_json.is_number_integer())
568 {
569 disc_orders.setConstant(tmp_json);
570 }
571 else if (tmp_json.is_string())
572 {
573 const std::string discr_orders_path = utils::resolve_path(tmp_json, root_path());
574 Eigen::MatrixXi tmp;
575 read_matrix(discr_orders_path, tmp);
576 assert(tmp.size() == disc_orders.size());
577 assert(tmp.cols() == 1);
578 disc_orders = tmp;
579 }
580 else if (tmp_json.is_array())
581 {
582 const auto b_discr_orders = tmp_json;
583
584 std::map<int, int> b_orders;
585 for (size_t i = 0; i < b_discr_orders.size(); ++i)
586 {
587 assert(b_discr_orders[i]["id"].is_array() || b_discr_orders[i]["id"].is_number_integer());
588
589 const int order = b_discr_orders[i]["order"];
590 for (const int id : json_as_array<int>(b_discr_orders[i]["id"]))
591 {
592 b_orders[id] = order;
593 logger().trace("bid {}, discr {}", id, order);
594 }
595 }
596
597 for (int e = 0; e < mesh->n_elements(); ++e)
598 {
599 const int bid = mesh->get_body_id(e);
600 const auto order = b_orders.find(bid);
601 if (order == b_orders.end())
602 {
603 logger().debug("Missing discretization order for body {}; using 1", bid);
604 b_orders[bid] = 1;
605 disc_orders[e] = 1;
606 }
607 else
608 {
609 disc_orders[e] = order->second;
610 }
611 }
612 }
613 else
614 {
615 logger().error("space/discr_order must be either a number a path or an array");
616 throw std::runtime_error("invalid json");
617 }
618 // TODO: same for pressure!
619
620#ifdef POLYFEM_WITH_BEZIER
621 if (!mesh->is_simplicial())
622#else
623 if constexpr (true)
624#endif
625 {
626 args["space"]["advanced"]["count_flipped_els_continuous"] = false;
627 args["output"]["paraview"]["options"]["jacobian_validity"] = false;
628 args["solver"]["advanced"]["check_inversion"] = "Discrete";
629 }
630 else if (args["solver"]["advanced"]["check_inversion"] != "Discrete")
631 {
632 args["space"]["advanced"]["use_corner_quadrature"] = true;
633 }
634
635 Eigen::MatrixXi geom_disc_orders;
636 if (!iso_parametric())
637 {
638 if (mesh->orders().size() <= 0)
639 {
640 geom_disc_orders.resizeLike(disc_orders);
641 geom_disc_orders.setConstant(1);
642 }
643 else
644 geom_disc_orders = mesh->orders();
645 }
646
647 igl::Timer timer;
648 timer.start();
649 if (args["space"]["use_p_ref"])
650 {
652 *mesh,
653 args["space"]["advanced"]["B"],
654 args["space"]["advanced"]["h1_formula"],
655 args["space"]["discr_order"],
656 args["space"]["advanced"]["discr_order_max"],
657 stats,
659
660 logger().info("min p: {} max p: {}", disc_orders.minCoeff(), disc_orders.maxCoeff());
661 }
662
663 int quadrature_order = args["space"]["advanced"]["quadrature_order"].get<int>();
664 const int mass_quadrature_order = args["space"]["advanced"]["mass_quadrature_order"].get<int>();
665 if (mixed_assembler != nullptr)
666 {
667 const int disc_order = disc_orders.maxCoeff();
668 if (disc_order - disc_orders.minCoeff() != 0)
669 {
670 logger().error("p refinement not supported in mixed formulation!");
671 return;
672 }
673 }
674
675 // shape optimization needs continuous geometric basis
676 // const bool use_continuous_gbasis = optimization_enabled == solver::CacheLevel::Derivatives;
677 const bool use_continuous_gbasis = true;
678 const bool use_corner_quadrature = args["space"]["advanced"]["use_corner_quadrature"];
679
680 if (mesh->is_volume())
681 {
682 const Mesh3D &tmp_mesh = *dynamic_cast<Mesh3D *>(mesh.get());
683 if (args["space"]["basis_type"] == "Spline")
684 {
685 // if (!iso_parametric())
686 // {
687 // logger().error("Splines must be isoparametric, ignoring...");
688 // // LagrangeBasis3d::build_bases(tmp_mesh, quadrature_order, geom_disc_orders, has_polys, geom_bases_, local_boundary, poly_edge_to_data_geom, mesh_nodes);
689 // SplineBasis3d::build_bases(tmp_mesh, quadrature_order, geom_bases_, local_boundary, poly_edge_to_data);
690 // }
691
692 n_bases = basis::SplineBasis3d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, bases, local_boundary, poly_edge_to_data);
693
694 // if (iso_parametric() && args["fit_nodes"])
695 // SplineBasis3d::fit_nodes(tmp_mesh, n_bases, bases);
696 }
697 else
698 {
699 if (!iso_parametric())
700 n_geom_bases = basis::LagrangeBasis3d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, geom_disc_orders, false, false, has_polys, !use_continuous_gbasis, use_corner_quadrature, geom_bases_, local_boundary, poly_edge_to_data_geom, geom_mesh_nodes);
701
702 n_bases = basis::LagrangeBasis3d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, disc_orders, args["space"]["basis_type"] == "Bernstein", args["space"]["basis_type"] == "Serendipity", has_polys, false, use_corner_quadrature, bases, local_boundary, poly_edge_to_data, mesh_nodes);
703 }
704
705 // if(problem->is_mixed())
706 if (mixed_assembler != nullptr)
707 {
708 n_pressure_bases = basis::LagrangeBasis3d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, int(args["space"]["pressure_discr_order"]), args["space"]["basis_type"] == "Bernstein", false, has_polys, false, use_corner_quadrature, pressure_bases, local_boundary, poly_edge_to_data_geom, pressure_mesh_nodes);
709 }
710 }
711 else
712 {
713 const Mesh2D &tmp_mesh = *dynamic_cast<Mesh2D *>(mesh.get());
714 if (args["space"]["basis_type"] == "Spline")
715 {
716 // TODO:
717 // if (!iso_parametric())
718 // {
719 // logger().error("Splines must be isoparametric, ignoring...");
720 // // LagrangeBasis2d::build_bases(tmp_mesh, quadrature_order, disc_orders, has_polys, geom_bases_, local_boundary, poly_edge_to_data_geom, mesh_nodes);
721 // n_bases = SplineBasis2d::build_bases(tmp_mesh, quadrature_order, geom_bases_, local_boundary, poly_edge_to_data);
722 // }
723
724 n_bases = basis::SplineBasis2d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, bases, local_boundary, poly_edge_to_data);
725
726 // if (iso_parametric() && args["fit_nodes"])
727 // SplineBasis2d::fit_nodes(tmp_mesh, n_bases, bases);
728 }
729 else
730 {
731 if (!iso_parametric())
732 n_geom_bases = basis::LagrangeBasis2d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, geom_disc_orders, false, false, has_polys, !use_continuous_gbasis, use_corner_quadrature, geom_bases_, local_boundary, poly_edge_to_data_geom, geom_mesh_nodes);
733
734 n_bases = basis::LagrangeBasis2d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, disc_orders, args["space"]["basis_type"] == "Bernstein", args["space"]["basis_type"] == "Serendipity", has_polys, false, use_corner_quadrature, bases, local_boundary, poly_edge_to_data, mesh_nodes);
735 }
736
737 // if(problem->is_mixed())
738 if (mixed_assembler != nullptr)
739 {
740 n_pressure_bases = basis::LagrangeBasis2d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, int(args["space"]["pressure_discr_order"]), args["space"]["basis_type"] == "Bernstein", false, has_polys, false, use_corner_quadrature, pressure_bases, local_boundary, poly_edge_to_data_geom, pressure_mesh_nodes);
741 }
742 }
743
744 if (mixed_assembler != nullptr)
745 {
746 assert(bases.size() == pressure_bases.size());
747 for (int i = 0; i < pressure_bases.size(); ++i)
748 {
750 bases[i].compute_quadrature(b_quad);
751 pressure_bases[i].set_quadrature([b_quad](quadrature::Quadrature &quad) { quad = b_quad; });
752 }
753 }
754
755 timer.stop();
756
758
759 if (n_geom_bases == 0)
761
762 auto &gbases = geom_bases();
763
765 {
766 std::map<std::array<int, 2>, double> pairs;
767 for (int e = 0; e < gbases.size(); e++)
768 {
769 const auto &gbs = gbases[e].bases;
770 const auto &bs = bases[e].bases;
771
772 Eigen::MatrixXd local_pts;
773 const int order = bs.front().order();
774 if (mesh->is_volume())
775 {
776 if (mesh->is_simplex(e))
777 autogen::p_nodes_3d(order, local_pts);
778 else
779 autogen::q_nodes_3d(order, local_pts);
780 }
781 else
782 {
783 if (mesh->is_simplex(e))
784 autogen::p_nodes_2d(order, local_pts);
785 else
786 autogen::q_nodes_2d(order, local_pts);
787 }
788
790 vals.compute(e, mesh->is_volume(), local_pts, gbases[e], gbases[e]);
791
792 for (int i = 0; i < bs.size(); i++)
793 {
794 for (int j = 0; j < gbs.size(); j++)
795 {
796 if (std::abs(vals.basis_values[j].val(i)) > 1e-7)
797 {
798 std::array<int, 2> index = {{gbs[j].global()[0].index, bs[i].global()[0].index}};
799 pairs.insert({index, vals.basis_values[j].val(i)});
800 }
801 }
802 }
803 }
804
805 const int dim = mesh->dimension();
806 std::vector<Eigen::Triplet<double>> coeffs;
807 coeffs.clear();
808 for (const auto &iter : pairs)
809 for (int d = 0; d < dim; d++)
810 coeffs.emplace_back(iter.first[0] * dim + d, iter.first[1] * dim + d, iter.second);
811
812 basis_nodes_to_gbasis_nodes.resize(n_geom_bases * mesh->dimension(), n_bases * mesh->dimension());
813 basis_nodes_to_gbasis_nodes.setFromTriplets(coeffs.begin(), coeffs.end());
814 }
815
816 for (const auto &lb : local_boundary)
817 total_local_boundary.emplace_back(lb);
818
819 const int dim = mesh->dimension();
820 const int problem_dim = problem->is_scalar() ? 1 : dim;
821
822 // handle periodic bc
823 if (has_periodic_bc())
824 {
825 // collect periodic directions
826 json directions = args["boundary_conditions"]["periodic_boundary"]["correspondence"];
827 Eigen::MatrixXd tile_offset = Eigen::MatrixXd::Identity(dim, dim);
828
829 if (directions.size() > 0)
830 {
831 Eigen::VectorXd tmp;
832 for (int d = 0; d < dim; d++)
833 {
834 tmp = directions[d];
835 if (tmp.size() != dim)
836 log_and_throw_error("Invalid size of periodic directions!");
837 tile_offset.col(d) = tmp;
838 }
839 }
840
841 periodic_bc = std::make_shared<PeriodicBoundary>(problem->is_scalar(), n_bases, bases, mesh_nodes, tile_offset, args["boundary_conditions"]["periodic_boundary"]["tolerance"].get<double>());
842
843 macro_strain_constraint.init(dim, args["boundary_conditions"]["periodic_boundary"]);
844 }
845
846 if (args["space"]["advanced"]["count_flipped_els"])
848
849 const int prev_bases = n_bases;
851
852 {
853 igl::Timer timer2;
854 logger().debug("Building node mapping...");
855 timer2.start();
857 problem->update_nodes(in_node_to_node);
858 mesh->update_nodes(in_node_to_node);
859 timer2.stop();
860 logger().debug("Done (took {}s)", timer2.getElapsedTime());
861 }
862
863 logger().info("Building collision mesh...");
865 if (periodic_bc && args["contact"]["periodic"])
867 logger().info("Done!");
868
869 const int prev_b_size = local_boundary.size();
870 problem->setup_bc(*mesh, n_bases - obstacle.n_vertices(),
879
880 // setp nodal values
881 {
883 for (int n = 0; n < dirichlet_nodes.size(); ++n)
884 {
885 const int n_id = dirichlet_nodes[n];
886 bool found = false;
887 for (const auto &bs : bases)
888 {
889 for (const auto &b : bs.bases)
890 {
891 for (const auto &lg : b.global())
892 {
893 if (lg.index == n_id)
894 {
895 dirichlet_nodes_position[n] = lg.node;
896 found = true;
897 break;
898 }
899 }
900
901 if (found)
902 break;
903 }
904
905 if (found)
906 break;
907 }
908
909 assert(found);
910 }
911
913 for (int n = 0; n < neumann_nodes.size(); ++n)
914 {
915 const int n_id = neumann_nodes[n];
916 bool found = false;
917 for (const auto &bs : bases)
918 {
919 for (const auto &b : bs.bases)
920 {
921 for (const auto &lg : b.global())
922 {
923 if (lg.index == n_id)
924 {
925 neumann_nodes_position[n] = lg.node;
926 found = true;
927 break;
928 }
929 }
930
931 if (found)
932 break;
933 }
934
935 if (found)
936 break;
937 }
938
939 assert(found);
940 }
941 }
942
943 const bool has_neumann = local_neumann_boundary.size() > 0 || local_boundary.size() < prev_b_size;
944 use_avg_pressure = !has_neumann;
945
946 for (int i = prev_bases; i < n_bases; ++i)
947 {
948 for (int d = 0; d < problem_dim; ++d)
949 boundary_nodes.push_back(i * problem_dim + d);
950 }
951
952 std::sort(boundary_nodes.begin(), boundary_nodes.end());
953 auto it = std::unique(boundary_nodes.begin(), boundary_nodes.end());
954 boundary_nodes.resize(std::distance(boundary_nodes.begin(), it));
955
956 // for elastic pure periodic problem, find an internal node and force zero dirichlet
957 if ((!problem->is_time_dependent() || args["time"]["quasistatic"]) && boundary_nodes.size() == 0 && has_periodic_bc())
958 {
959 // find an internal node to force zero dirichlet
960 std::vector<bool> isboundary(n_bases, false);
961 for (const auto &lb : total_local_boundary)
962 {
963 const int e = lb.element_id();
964 for (int i = 0; i < lb.size(); ++i)
965 {
966 const auto nodes = bases[e].local_nodes_for_primitive(lb.global_primitive_id(i), *mesh);
967
968 for (int n : nodes)
969 isboundary[bases[e].bases[n].global()[0].index] = true;
970 }
971 }
972 int i = 0;
973 for (; i < n_bases; i++)
974 if (!isboundary[i]) // (!periodic_bc->is_periodic_dof(i))
975 break;
976 if (i >= n_bases)
977 log_and_throw_error("Failed to find a non-periodic node!");
978 const int actual_dim = problem->is_scalar() ? 1 : mesh->dimension();
979 for (int d = 0; d < actual_dim; d++)
980 {
981 boundary_nodes.push_back(i * actual_dim + d);
982 }
983 logger().info("Fix solution at node {} to remove singularity due to periodic BC", i);
984 }
985
986 const auto &curret_bases = geom_bases();
987 const int n_samples = 10;
988 stats.compute_mesh_size(*mesh, curret_bases, n_samples, args["output"]["advanced"]["curved_mesh_size"]);
990 {
992 }
994 {
996 }
997
998 if (is_contact_enabled())
999 {
1000 min_boundary_edge_length = std::numeric_limits<double>::max();
1001 for (const auto &edge : collision_mesh.edges().rowwise())
1002 {
1003 const VectorNd v0 = collision_mesh.rest_positions().row(edge(0));
1004 const VectorNd v1 = collision_mesh.rest_positions().row(edge(1));
1005 min_boundary_edge_length = std::min(min_boundary_edge_length, (v1 - v0).norm());
1006 }
1007
1008 double dhat = Units::convert(args["contact"]["dhat"], units.length());
1009 args["contact"]["epsv"] = Units::convert(args["contact"]["epsv"], units.velocity());
1010 args["contact"]["dhat"] = dhat;
1011
1012 if (!has_dhat && dhat > min_boundary_edge_length)
1013 {
1014 args["contact"]["dhat"] = double(args["contact"]["dhat_percentage"]) * min_boundary_edge_length;
1015 logger().info("dhat set to {}", double(args["contact"]["dhat"]));
1016 }
1017 else
1018 {
1019 if (dhat > min_boundary_edge_length)
1020 logger().warn("dhat larger than min boundary edge, {} > {}", dhat, min_boundary_edge_length);
1021 }
1022 }
1023
1024 logger().info("n_bases {}", n_bases);
1025
1026 timings.building_basis_time = timer.getElapsedTime();
1027 logger().info(" took {}s", timings.building_basis_time);
1028
1029 logger().info("flipped elements {}", stats.n_flipped);
1030 logger().info("h: {}", stats.mesh_size);
1031 logger().info("n bases: {}", n_bases);
1032 logger().info("n pressure bases: {}", n_pressure_bases);
1033
1036 if (n_bases <= args["solver"]["advanced"]["cache_size"])
1037 {
1038 timer.start();
1039 logger().info("Building cache...");
1040 ass_vals_cache.init(mesh->is_volume(), bases, curret_bases);
1041 mass_ass_vals_cache.init(mesh->is_volume(), bases, curret_bases, true);
1042 if (mixed_assembler != nullptr)
1043 pressure_ass_vals_cache.init(mesh->is_volume(), pressure_bases, curret_bases);
1044
1045 logger().info(" took {}s", timer.getElapsedTime());
1046 }
1047
1048 out_geom.build_grid(*mesh, args["output"]["advanced"]["sol_on_grid"]);
1049
1050 if ((!problem->is_time_dependent() || args["time"]["quasistatic"]) && boundary_nodes.empty())
1051 {
1052 if (has_periodic_bc())
1053 logger().warn("(Quasi-)Static problem without Dirichlet nodes, will fix solution at one node to find a unique solution!");
1054 else
1055 {
1056 if (args["constraints"]["hard"].empty())
1057 log_and_throw_error("Static problem need to have some Dirichlet nodes!");
1058 else
1059 logger().warn("Relying on hard constraints to avoid infinite solutions");
1060 }
1061 }
1062 }
1063
1065 {
1066 if (!mesh)
1067 {
1068 logger().error("Load the mesh first!");
1069 return;
1070 }
1071
1072 rhs.resize(0, 0);
1073
1074 if (poly_edge_to_data.empty() && polys.empty())
1075 {
1077 return;
1078 }
1079
1080 igl::Timer timer;
1081 timer.start();
1082 logger().info("Computing polygonal basis...");
1083
1084 // std::sort(boundary_nodes.begin(), boundary_nodes.end());
1085
1086 // mixed not supports polygonal bases
1087 assert(n_pressure_bases == 0 || poly_edge_to_data.size() == 0);
1088
1089 int new_bases = 0;
1090
1091 if (iso_parametric())
1092 {
1093 if (mesh->is_volume())
1094 {
1095 if (args["space"]["poly_basis_type"] == "MeanValue" || args["space"]["poly_basis_type"] == "Wachspress")
1096 logger().error("Barycentric bases not supported in 3D");
1097 assert(assembler->is_linear());
1099 *dynamic_cast<LinearAssembler *>(assembler.get()),
1100 args["space"]["advanced"]["n_harmonic_samples"],
1101 *dynamic_cast<Mesh3D *>(mesh.get()),
1102 n_bases,
1103 args["space"]["advanced"]["quadrature_order"],
1104 args["space"]["advanced"]["mass_quadrature_order"],
1105 args["space"]["advanced"]["integral_constraints"],
1106 bases,
1107 bases,
1109 polys_3d);
1110 }
1111 else
1112 {
1113 if (args["space"]["poly_basis_type"] == "MeanValue")
1114 {
1116 assembler->name(),
1117 assembler->is_tensor() ? 2 : 1,
1118 *dynamic_cast<Mesh2D *>(mesh.get()),
1119 n_bases,
1120 args["space"]["advanced"]["quadrature_order"],
1121 args["space"]["advanced"]["mass_quadrature_order"],
1123 }
1124 else if (args["space"]["poly_basis_type"] == "Wachspress")
1125 {
1127 assembler->name(),
1128 assembler->is_tensor() ? 2 : 1,
1129 *dynamic_cast<Mesh2D *>(mesh.get()),
1130 n_bases,
1131 args["space"]["advanced"]["quadrature_order"],
1132 args["space"]["advanced"]["mass_quadrature_order"],
1134 }
1135 else
1136 {
1137 assert(assembler->is_linear());
1139 *dynamic_cast<LinearAssembler *>(assembler.get()),
1140 args["space"]["advanced"]["n_harmonic_samples"],
1141 *dynamic_cast<Mesh2D *>(mesh.get()),
1142 n_bases,
1143 args["space"]["advanced"]["quadrature_order"],
1144 args["space"]["advanced"]["mass_quadrature_order"],
1145 args["space"]["advanced"]["integral_constraints"],
1146 bases,
1147 bases,
1149 polys);
1150 }
1151 }
1152 }
1153 else
1154 {
1155 if (mesh->is_volume())
1156 {
1157 if (args["space"]["poly_basis_type"] == "MeanValue" || args["space"]["poly_basis_type"] == "Wachspress")
1158 {
1159 logger().error("Barycentric bases not supported in 3D");
1160 throw std::runtime_error("not implemented");
1161 }
1162 else
1163 {
1164 assert(assembler->is_linear());
1166 *dynamic_cast<LinearAssembler *>(assembler.get()),
1167 args["space"]["advanced"]["n_harmonic_samples"],
1168 *dynamic_cast<Mesh3D *>(mesh.get()),
1169 n_bases,
1170 args["space"]["advanced"]["quadrature_order"],
1171 args["space"]["advanced"]["mass_quadrature_order"],
1172 args["space"]["advanced"]["integral_constraints"],
1173 bases,
1176 polys_3d);
1177 }
1178 }
1179 else
1180 {
1181 if (args["space"]["poly_basis_type"] == "MeanValue")
1182 {
1184 assembler->name(),
1185 assembler->is_tensor() ? 2 : 1,
1186 *dynamic_cast<Mesh2D *>(mesh.get()),
1187 n_bases, args["space"]["advanced"]["quadrature_order"],
1188 args["space"]["advanced"]["mass_quadrature_order"],
1190 }
1191 else if (args["space"]["poly_basis_type"] == "Wachspress")
1192 {
1194 assembler->name(),
1195 assembler->is_tensor() ? 2 : 1,
1196 *dynamic_cast<Mesh2D *>(mesh.get()),
1197 n_bases, args["space"]["advanced"]["quadrature_order"],
1198 args["space"]["advanced"]["mass_quadrature_order"],
1200 }
1201 else
1202 {
1203 assert(assembler->is_linear());
1205 *dynamic_cast<LinearAssembler *>(assembler.get()),
1206 args["space"]["advanced"]["n_harmonic_samples"],
1207 *dynamic_cast<Mesh2D *>(mesh.get()),
1208 n_bases,
1209 args["space"]["advanced"]["quadrature_order"],
1210 args["space"]["advanced"]["mass_quadrature_order"],
1211 args["space"]["advanced"]["integral_constraints"],
1212 bases,
1215 polys);
1216 }
1217 }
1218 }
1219
1220 timer.stop();
1221 timings.computing_poly_basis_time = timer.getElapsedTime();
1222 logger().info(" took {}s", timings.computing_poly_basis_time);
1223
1224 n_bases += new_bases;
1225 }
1226
1228 {
1229 assert(!mesh->is_volume());
1230 const int dim = mesh->dimension();
1231 const int n_tiles = 2;
1232
1233 if (mesh->dimension() != 2)
1234 log_and_throw_error("Periodic collision mesh is only implemented in 2D!");
1235
1236 Eigen::MatrixXd V(n_bases, dim);
1237 for (const auto &bs : bases)
1238 for (const auto &b : bs.bases)
1239 for (const auto &g : b.global())
1240 V.row(g.index) = g.node;
1241
1242 Eigen::MatrixXi E = collision_mesh.edges();
1243 for (int i = 0; i < E.size(); i++)
1244 E(i) = collision_mesh.to_full_vertex_id(E(i));
1245
1246 Eigen::MatrixXd bbox(V.cols(), 2);
1247 bbox.col(0) = V.colwise().minCoeff();
1248 bbox.col(1) = V.colwise().maxCoeff();
1249
1250 // remove boundary edges on periodic BC, buggy
1251 {
1252 std::vector<int> ind;
1253 for (int i = 0; i < E.rows(); i++)
1254 {
1255 if (!periodic_bc->is_periodic_dof(E(i, 0)) || !periodic_bc->is_periodic_dof(E(i, 1)))
1256 ind.push_back(i);
1257 }
1258
1259 E = E(ind, Eigen::all).eval();
1260 }
1261
1262 Eigen::MatrixXd Vtmp, Vnew;
1263 Eigen::MatrixXi Etmp, Enew;
1264 Vtmp.setZero(V.rows() * n_tiles * n_tiles, V.cols());
1265 Etmp.setZero(E.rows() * n_tiles * n_tiles, E.cols());
1266
1267 Eigen::MatrixXd tile_offset = periodic_bc->get_affine_matrix();
1268
1269 for (int i = 0, idx = 0; i < n_tiles; i++)
1270 {
1271 for (int j = 0; j < n_tiles; j++)
1272 {
1273 Eigen::Vector2d block_id;
1274 block_id << i, j;
1275
1276 Vtmp.middleRows(idx * V.rows(), V.rows()) = V;
1277 // Vtmp.block(idx * V.rows(), 0, V.rows(), 1).array() += tile_offset(0) * i;
1278 // Vtmp.block(idx * V.rows(), 1, V.rows(), 1).array() += tile_offset(1) * j;
1279 for (int vid = 0; vid < V.rows(); vid++)
1280 Vtmp.block(idx * V.rows() + vid, 0, 1, 2) += (tile_offset * block_id).transpose();
1281
1282 Etmp.middleRows(idx * E.rows(), E.rows()) = E.array() + idx * V.rows();
1283 idx += 1;
1284 }
1285 }
1286
1287 // clean duplicated vertices
1288 Eigen::VectorXi indices;
1289 {
1290 std::vector<int> tmp;
1291 for (int i = 0; i < V.rows(); i++)
1292 {
1293 if (periodic_bc->is_periodic_dof(i))
1294 tmp.push_back(i);
1295 }
1296
1297 indices.resize(tmp.size() * n_tiles * n_tiles);
1298 for (int i = 0; i < n_tiles * n_tiles; i++)
1299 {
1300 indices.segment(i * tmp.size(), tmp.size()) = Eigen::Map<Eigen::VectorXi, Eigen::Unaligned>(tmp.data(), tmp.size());
1301 indices.segment(i * tmp.size(), tmp.size()).array() += i * V.rows();
1302 }
1303 }
1304
1305 Eigen::VectorXi potentially_duplicate_mask(Vtmp.rows());
1306 potentially_duplicate_mask.setZero();
1307 potentially_duplicate_mask(indices).array() = 1;
1308
1309 Eigen::MatrixXd candidates = Vtmp(indices, Eigen::all);
1310
1311 Eigen::VectorXi SVI;
1312 std::vector<int> SVJ;
1313 SVI.setConstant(Vtmp.rows(), -1);
1314 int id = 0;
1315 const double eps = (bbox.col(1) - bbox.col(0)).maxCoeff() * args["boundary_conditions"]["periodic_boundary"]["tolerance"].get<double>();
1316 for (int i = 0; i < Vtmp.rows(); i++)
1317 {
1318 if (SVI[i] < 0)
1319 {
1320 SVI[i] = id;
1321 SVJ.push_back(i);
1322 if (potentially_duplicate_mask(i))
1323 {
1324 Eigen::VectorXd diffs = (candidates.rowwise() - Vtmp.row(i)).rowwise().norm();
1325 for (int j = 0; j < diffs.size(); j++)
1326 if (diffs(j) < eps)
1327 SVI[indices[j]] = id;
1328 }
1329 id++;
1330 }
1331 }
1332
1333 Vnew = Vtmp(SVJ, Eigen::all);
1334
1335 Enew.resizeLike(Etmp);
1336 for (int d = 0; d < Etmp.cols(); d++)
1337 Enew.col(d) = SVI(Etmp.col(d));
1338
1339 std::vector<bool> is_on_surface = ipc::CollisionMesh::construct_is_on_surface(Vnew.rows(), Enew);
1340
1341 Eigen::MatrixXi boundary_triangles;
1342 Eigen::SparseMatrix<double> displacement_map;
1343 periodic_collision_mesh = ipc::CollisionMesh(is_on_surface,
1344 std::vector<bool>(Vnew.rows(), false),
1345 Vnew,
1346 Enew,
1347 boundary_triangles,
1348 displacement_map);
1349
1350 periodic_collision_mesh.init_area_jacobians();
1351
1352 periodic_collision_mesh_to_basis.setConstant(Vnew.rows(), -1);
1353 for (int i = 0; i < V.rows(); i++)
1354 for (int j = 0; j < n_tiles * n_tiles; j++)
1355 periodic_collision_mesh_to_basis(SVI[j * V.rows() + i]) = i;
1356
1357 if (periodic_collision_mesh_to_basis.maxCoeff() + 1 != V.rows())
1358 log_and_throw_error("Failed to tile mesh!");
1359 }
1360
1362 {
1365 args, [this](const std::string &p) { return resolve_input_path(p); },
1367 }
1368
1370 const mesh::Mesh &mesh,
1371 const int n_bases,
1372 const std::vector<basis::ElementBases> &bases,
1373 const std::vector<basis::ElementBases> &geom_bases,
1374 const std::vector<mesh::LocalBoundary> &total_local_boundary,
1375 const mesh::Obstacle &obstacle,
1376 const json &args,
1377 const std::function<std::string(const std::string &)> &resolve_input_path,
1378 const Eigen::VectorXi &in_node_to_node,
1379 ipc::CollisionMesh &collision_mesh)
1380 {
1381 Eigen::MatrixXd collision_vertices;
1382 Eigen::VectorXi collision_codim_vids;
1383 Eigen::MatrixXi collision_edges, collision_triangles;
1384 std::vector<Eigen::Triplet<double>> displacement_map_entries;
1385
1386 if (args.contains("/contact/collision_mesh"_json_pointer)
1387 && args.at("/contact/collision_mesh/enabled"_json_pointer).get<bool>())
1388 {
1389 const json collision_mesh_args = args.at("/contact/collision_mesh"_json_pointer);
1390 if (collision_mesh_args.contains("linear_map"))
1391 {
1392 assert(displacement_map_entries.empty());
1393 assert(collision_mesh_args.contains("mesh"));
1394 const std::string root_path = utils::json_value<std::string>(args, "root_path", "");
1395 // TODO: handle transformation per geometry
1396 const json transformation = json_as_array(args["geometry"])[0]["transformation"];
1398 utils::resolve_path(collision_mesh_args["mesh"], root_path),
1399 utils::resolve_path(collision_mesh_args["linear_map"], root_path),
1400 in_node_to_node, transformation, collision_vertices, collision_codim_vids,
1401 collision_edges, collision_triangles, displacement_map_entries);
1402 }
1403 else if (collision_mesh_args.contains("max_edge_length"))
1404 {
1405 logger().debug(
1406 "Building collision proxy with max edge length={} ...",
1407 collision_mesh_args["max_edge_length"].get<double>());
1408 igl::Timer timer;
1409 timer.start();
1412 collision_mesh_args["max_edge_length"], collision_vertices,
1413 collision_triangles, displacement_map_entries,
1414 collision_mesh_args["tessellation_type"]);
1415 if (collision_triangles.size())
1416 igl::edges(collision_triangles, collision_edges);
1417 timer.stop();
1418 logger().debug(fmt::format(
1419 std::locale("en_US.UTF-8"),
1420 "Done (took {:g}s, {:L} vertices, {:L} triangles)",
1421 timer.getElapsedTime(),
1422 collision_vertices.rows(), collision_triangles.rows()));
1423 }
1424 else
1425 {
1428 collision_vertices, collision_edges, collision_triangles, displacement_map_entries);
1429 }
1430 }
1431 else
1432 {
1435 collision_vertices, collision_edges, collision_triangles, displacement_map_entries);
1436 }
1437
1438 std::vector<bool> is_orientable_vertex(collision_vertices.rows(), true);
1439
1440 // n_bases already contains the obstacle vertices
1441 const int num_fe_nodes = n_bases - obstacle.n_vertices();
1442 const int num_fe_collision_vertices = collision_vertices.rows();
1443 assert(collision_edges.size() == 0 || collision_edges.maxCoeff() < num_fe_collision_vertices);
1444 assert(collision_triangles.size() == 0 || collision_triangles.maxCoeff() < num_fe_collision_vertices);
1445
1446 // Append the obstacles to the collision mesh
1447 if (obstacle.n_vertices() > 0)
1448 {
1449 append_rows(collision_vertices, obstacle.v());
1450 append_rows(collision_codim_vids, obstacle.codim_v().array() + num_fe_collision_vertices);
1451 append_rows(collision_edges, obstacle.e().array() + num_fe_collision_vertices);
1452 append_rows(collision_triangles, obstacle.f().array() + num_fe_collision_vertices);
1453
1454 for (int i = 0; i < obstacle.n_vertices(); i++)
1455 {
1456 is_orientable_vertex.push_back(false);
1457 }
1458
1459 if (!displacement_map_entries.empty())
1460 {
1461 displacement_map_entries.reserve(displacement_map_entries.size() + obstacle.n_vertices());
1462 for (int i = 0; i < obstacle.n_vertices(); i++)
1463 {
1464 displacement_map_entries.emplace_back(num_fe_collision_vertices + i, num_fe_nodes + i, 1.0);
1465 }
1466 }
1467 }
1468
1469 std::vector<bool> is_on_surface = ipc::CollisionMesh::construct_is_on_surface(
1470 collision_vertices.rows(), collision_edges);
1471 for (const int vid : collision_codim_vids)
1472 {
1473 is_on_surface[vid] = true;
1474 }
1475
1476 Eigen::SparseMatrix<double> displacement_map;
1477 if (!displacement_map_entries.empty())
1478 {
1479 displacement_map.resize(collision_vertices.rows(), n_bases);
1480 displacement_map.setFromTriplets(displacement_map_entries.begin(), displacement_map_entries.end());
1481 }
1482
1483 collision_mesh = ipc::CollisionMesh(
1484 is_on_surface, is_orientable_vertex, collision_vertices, collision_edges, collision_triangles,
1485 displacement_map);
1486
1487 collision_mesh.can_collide = [&collision_mesh, num_fe_collision_vertices](size_t vi, size_t vj) {
1488 // obstacles do not collide with other obstacles
1489 return collision_mesh.to_full_vertex_id(vi) < num_fe_collision_vertices
1490 || collision_mesh.to_full_vertex_id(vj) < num_fe_collision_vertices;
1491 };
1492
1493 collision_mesh.init_area_jacobians();
1494 }
1495
1497 {
1498 if (!mesh)
1499 {
1500 logger().error("Load the mesh first!");
1501 return;
1502 }
1503 if (n_bases <= 0)
1504 {
1505 logger().error("Build the bases first!");
1506 return;
1507 }
1508 if (assembler->name() == "OperatorSplitting")
1509 {
1511 avg_mass = 1;
1512 return;
1513 }
1514
1515 if (!problem->is_time_dependent())
1516 {
1517 avg_mass = 1;
1519 return;
1520 }
1521
1522 mass.resize(0, 0);
1523
1524 igl::Timer timer;
1525 timer.start();
1526 logger().info("Assembling mass mat...");
1527
1528 if (mixed_assembler != nullptr)
1529 {
1530 StiffnessMatrix velocity_mass;
1531 mass_matrix_assembler->assemble(mesh->is_volume(), n_bases, bases, geom_bases(), mass_ass_vals_cache, 0, velocity_mass, true);
1532
1533 std::vector<Eigen::Triplet<double>> mass_blocks;
1534 mass_blocks.reserve(velocity_mass.nonZeros());
1535
1536 for (int k = 0; k < velocity_mass.outerSize(); ++k)
1537 {
1538 for (StiffnessMatrix::InnerIterator it(velocity_mass, k); it; ++it)
1539 {
1540 mass_blocks.emplace_back(it.row(), it.col(), it.value());
1541 }
1542 }
1543
1544 mass.resize(n_bases * assembler->size(), n_bases * assembler->size());
1545 mass.setFromTriplets(mass_blocks.begin(), mass_blocks.end());
1546 mass.makeCompressed();
1547 }
1548 else
1549 {
1550 mass_matrix_assembler->assemble(mesh->is_volume(), n_bases, bases, geom_bases(), mass_ass_vals_cache, 0, mass, true);
1551 }
1552
1553 assert(mass.size() > 0);
1554
1555 avg_mass = 0;
1556 for (int k = 0; k < mass.outerSize(); ++k)
1557 {
1558
1559 for (StiffnessMatrix::InnerIterator it(mass, k); it; ++it)
1560 {
1561 assert(it.col() == k);
1562 avg_mass += it.value();
1563 }
1564 }
1565
1566 avg_mass /= mass.rows();
1567 logger().info("average mass {}", avg_mass);
1568
1569 if (args["solver"]["advanced"]["lump_mass_matrix"])
1570 {
1572 }
1573
1574 timer.stop();
1575 timings.assembling_mass_mat_time = timer.getElapsedTime();
1576 logger().info(" took {}s", timings.assembling_mass_mat_time);
1577
1578 stats.nn_zero = mass.nonZeros();
1579 stats.num_dofs = mass.rows();
1580 stats.mat_size = (long long)mass.rows() * (long long)mass.cols();
1581 logger().info("sparsity: {}/{}", stats.nn_zero, stats.mat_size);
1582 }
1583
1584 std::shared_ptr<RhsAssembler> State::build_rhs_assembler(
1585 const int n_bases_,
1586 const std::vector<basis::ElementBases> &bases_,
1587 const assembler::AssemblyValsCache &ass_vals_cache_) const
1588 {
1589 json rhs_solver_params = args["solver"]["linear"];
1590 if (!rhs_solver_params.contains("Pardiso"))
1591 rhs_solver_params["Pardiso"] = {};
1592 rhs_solver_params["Pardiso"]["mtype"] = -2; // matrix type for Pardiso (2 = SPD)
1593
1594 const int size = problem->is_scalar() ? 1 : mesh->dimension();
1595
1596 return std::make_shared<RhsAssembler>(
1600 n_bases_, size, bases_, geom_bases(), ass_vals_cache_, *problem,
1601 args["space"]["advanced"]["bc_method"],
1602 rhs_solver_params);
1603 }
1604
1605 std::shared_ptr<PressureAssembler> State::build_pressure_assembler(
1606 const int n_bases_,
1607 const std::vector<basis::ElementBases> &bases_) const
1608 {
1609 const int size = problem->is_scalar() ? 1 : mesh->dimension();
1610
1611 return std::make_shared<PressureAssembler>(
1617 n_bases_, size, bases_, geom_bases(), *problem);
1618 }
1619
1621 {
1622 if (!mesh)
1623 {
1624 logger().error("Load the mesh first!");
1625 return;
1626 }
1627 if (n_bases <= 0)
1628 {
1629 logger().error("Build the bases first!");
1630 return;
1631 }
1632
1633 igl::Timer timer;
1634
1635 json p_params = {};
1636 p_params["formulation"] = assembler->name();
1637 p_params["root_path"] = root_path();
1638 {
1639 RowVectorNd min, max, delta;
1640 mesh->bounding_box(min, max);
1641 delta = (max - min) / 2. + min;
1642 if (mesh->is_volume())
1643 p_params["bbox_center"] = {delta(0), delta(1), delta(2)};
1644 else
1645 p_params["bbox_center"] = {delta(0), delta(1)};
1646 }
1647 problem->set_parameters(p_params);
1648
1649 rhs.resize(0, 0);
1650
1651 timer.start();
1652 logger().info("Assigning rhs...");
1653
1655 solve_data.rhs_assembler->assemble(mass_matrix_assembler->density(), rhs);
1656 rhs *= -1;
1657
1658 // if(problem->is_mixed())
1659 if (mixed_assembler != nullptr)
1660 {
1661 const int prev_size = rhs.size();
1662 const int n_larger = n_pressure_bases + (use_avg_pressure ? (assembler->is_fluid() ? 1 : 0) : 0);
1663 rhs.conservativeResize(prev_size + n_larger, rhs.cols());
1664 if (assembler->name() == "OperatorSplitting")
1665 {
1667 return;
1668 }
1669 // Divergence free rhs
1670 if (assembler->name() != "Bilaplacian" || local_neumann_boundary.empty())
1671 {
1672 rhs.block(prev_size, 0, n_larger, rhs.cols()).setZero();
1673 }
1674 else
1675 {
1676 Eigen::MatrixXd tmp(n_pressure_bases, 1);
1677 tmp.setZero();
1678
1679 std::shared_ptr<RhsAssembler> tmp_rhs_assembler = build_rhs_assembler(
1681
1682 tmp_rhs_assembler->set_bc(std::vector<LocalBoundary>(), std::vector<int>(), n_boundary_samples(), local_neumann_boundary, tmp);
1683 rhs.block(prev_size, 0, n_larger, rhs.cols()) = tmp;
1684 }
1685 }
1686
1687 timer.stop();
1688 timings.assigning_rhs_time = timer.getElapsedTime();
1689 logger().info(" took {}s", timings.assigning_rhs_time);
1690 }
1691
1692 void State::solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
1693 {
1694 if (!mesh)
1695 {
1696 logger().error("Load the mesh first!");
1697 return;
1698 }
1699 if (n_bases <= 0)
1700 {
1701 logger().error("Build the bases first!");
1702 return;
1703 }
1704
1705 // if (rhs.size() <= 0)
1706 // {
1707 // logger().error("Assemble the rhs first!");
1708 // return;
1709 // }
1710
1711 // sol.resize(0, 0);
1712 // pressure.resize(0, 0);
1713 stats.spectrum.setZero();
1714
1715 igl::Timer timer;
1716 timer.start();
1717 logger().info("Solving {}", assembler->name());
1718
1719 init_solve(sol, pressure);
1720
1721 if (problem->is_time_dependent())
1722 {
1723 const double t0 = args["time"]["t0"];
1724 const int time_steps = args["time"]["time_steps"];
1725 const double dt = args["time"]["dt"];
1726
1727 // Pre log the output path for easier watching
1728 if (args["output"]["advanced"]["save_time_sequence"])
1729 {
1730 logger().info("Time sequence of simulation will be written to: \"{}\"",
1731 resolve_output_path(args["output"]["paraview"]["file_name"]));
1732 }
1733
1734 if (assembler->name() == "NavierStokes")
1735 solve_transient_navier_stokes(time_steps, t0, dt, sol, pressure);
1736 else if (assembler->name() == "OperatorSplitting")
1737 solve_transient_navier_stokes_split(time_steps, dt, sol, pressure);
1738 else if (is_homogenization())
1739 solve_homogenization(time_steps, t0, dt, sol);
1740 else if (is_problem_linear())
1741 solve_transient_linear(time_steps, t0, dt, sol, pressure);
1742 else if (!assembler->is_linear() && problem->is_scalar())
1743 throw std::runtime_error("Nonlinear scalar problems are not supported yet!");
1744 else
1745 solve_transient_tensor_nonlinear(time_steps, t0, dt, sol);
1746 }
1747 else
1748 {
1749 if (assembler->name() == "NavierStokes")
1750 solve_navier_stokes(sol, pressure);
1751 else if (is_homogenization())
1752 solve_homogenization(/* time steps */ 0, /* t0 */ 0, /* dt */ 0, sol);
1753 else if (is_problem_linear())
1754 {
1755 init_linear_solve(sol);
1756 solve_linear(sol, pressure);
1758 cache_transient_adjoint_quantities(0, sol, Eigen::MatrixXd::Zero(mesh->dimension(), mesh->dimension()));
1759 }
1760 else if (!assembler->is_linear() && problem->is_scalar())
1761 throw std::runtime_error("Nonlinear scalar problems are not supported yet!");
1762 else
1763 {
1767 cache_transient_adjoint_quantities(0, sol, Eigen::MatrixXd::Zero(mesh->dimension(), mesh->dimension()));
1768
1769 const std::string state_path = resolve_output_path(args["output"]["data"]["state"]);
1770 if (!state_path.empty())
1771 write_matrix(state_path, "u", sol);
1772 }
1773 }
1774
1775 timer.stop();
1776 timings.solving_time = timer.getElapsedTime();
1777 logger().info(" took {}s", timings.solving_time);
1778 }
1779
1780} // namespace polyfem
int V
ElementAssemblyValues vals
Definition Assembler.cpp:22
int x
io::OutRuntimeData timings
runtime statistics
Definition State.hpp:590
assembler::MacroStrainValue macro_strain_constraint
Definition State.hpp:712
std::shared_ptr< utils::PeriodicBoundary > periodic_bc
periodic BC and periodic mesh utils
Definition State.hpp:386
int n_bases
number of bases
Definition State.hpp:178
void cache_transient_adjoint_quantities(const int current_step, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &disp_grad)
int n_pressure_bases
number of pressure bases
Definition State.hpp:180
std::string formulation() const
return the formulation (checks if the problem is scalar or not and deals with multiphysics)
Definition State.cpp:327
assembler::AssemblyValsCache ass_vals_cache
used to store assembly values for small problems
Definition State.hpp:196
int n_boundary_samples() const
quadrature used for projecting boundary conditions
Definition State.hpp:263
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:223
double starting_max_edge_length
Definition State.hpp:594
std::vector< mesh::LocalBoundary > local_pressure_boundary
mapping from elements to nodes for pressure boundary conditions
Definition State.hpp:442
void sol_to_pressure(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
splits the solution in solution and pressure for mixed problems
Definition State.cpp:370
std::string root_path() const
Get the root path for the state (e.g., args["root_path"] or ".")
void assemble_mass_mat()
assemble mass, step 4 of solve build mass matrix based on defined basis modifies mass (and maybe more...
Definition State.cpp:1496
Eigen::VectorXi in_node_to_node
Inpute nodes (including high-order) to polyfem nodes, only for isoparametric.
Definition State.hpp:455
void build_polygonal_basis()
builds bases for polygons, called inside build_basis
Definition State.cpp:1064
mesh::Obstacle obstacle
Obstacles used in collisions.
Definition State.hpp:473
std::shared_ptr< assembler::Assembler > assembler
assemblers
Definition State.hpp:155
ipc::CollisionMesh periodic_collision_mesh
IPC collision mesh under periodic BC.
Definition State.hpp:523
void build_periodic_collision_mesh()
Definition State.cpp:1227
io::OutGeometryData out_geom
visualization stuff
Definition State.hpp:588
bool use_avg_pressure
use average pressure for stokes problem to fix the additional dofs, true by default if false,...
Definition State.hpp:211
std::string resolve_output_path(const std::string &path) const
Resolve output path relative to output_dir if the path is not absolute.
std::vector< int > neumann_nodes
per node neumann
Definition State.hpp:451
ipc::CollisionMesh collision_mesh
IPC collision mesh.
Definition State.hpp:520
std::vector< basis::ElementBases > geom_bases_
Geometric mapping bases, if the elements are isoparametric, this list is empty.
Definition State.hpp:175
solver::CacheLevel optimization_enabled
Definition State.hpp:654
Eigen::VectorXi in_primitive_to_primitive
maps in vertices/edges/faces/cells to polyfem vertices/edges/faces/cells
Definition State.hpp:457
std::map< int, basis::InterfaceData > poly_edge_to_data
nodes on the boundary of polygonal elements, used for harmonic bases
Definition State.hpp:446
std::vector< int > pressure_boundary_nodes
list of neumann boundary nodes
Definition State.hpp:434
std::vector< basis::ElementBases > pressure_bases
FE pressure bases for mixed elements, the size is #elements.
Definition State.hpp:173
bool is_homogenization() const
Definition State.hpp:718
void solve_transient_tensor_nonlinear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol)
solves transient tensor nonlinear problem
double starting_min_edge_length
Definition State.hpp:593
std::shared_ptr< assembler::Mass > mass_matrix_assembler
Definition State.hpp:157
void build_node_mapping()
build the mapping from input nodes to polyfem nodes
Definition State.cpp:244
std::vector< int > primitive_to_node() const
Definition State.cpp:227
bool has_dhat
stores if input json contains dhat
Definition State.hpp:579
std::shared_ptr< polyfem::mesh::MeshNodes > pressure_mesh_nodes
Definition State.hpp:193
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:471
std::shared_ptr< polyfem::mesh::MeshNodes > mesh_nodes
Mapping from input nodes to FE nodes.
Definition State.hpp:193
StiffnessMatrix mass
Mass matrix, it is computed only for time dependent problems.
Definition State.hpp:202
void solve_homogenization(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol)
std::vector< int > dirichlet_nodes
per node dirichlet
Definition State.hpp:448
bool has_periodic_bc() const
Definition State.hpp:387
void solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves the problems
Definition State.cpp:1692
void build_basis()
builds the bases step 2 of solve modifies bases, pressure_bases, geom_bases_, boundary_nodes,...
Definition State.cpp:505
void solve_navier_stokes(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves a navier stokes
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:168
std::vector< RowVectorNd > dirichlet_nodes_position
Definition State.hpp:449
std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > polys_3d
polyhedra, used since poly have no geom mapping
Definition State.hpp:187
std::vector< int > node_to_primitive() const
Definition State.cpp:234
json args
main input arguments containing all defaults
Definition State.hpp:101
io::OutStatsData stats
Other statistics.
Definition State.hpp:592
std::vector< RowVectorNd > neumann_nodes_position
Definition State.hpp:452
void assemble_rhs()
compute rhs, step 3 of solve build rhs vector based on defined basis and given rhs of the problem mod...
Definition State.cpp:1620
double min_boundary_edge_length
Definition State.hpp:595
Eigen::VectorXi periodic_collision_mesh_to_basis
index mapping from periodic 2x2 collision mesh to FE periodic mesh
Definition State.hpp:525
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:171
void solve_transient_navier_stokes_split(const int time_steps, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves transient navier stokes with operator splitting
int n_geom_bases
number of geometric bases
Definition State.hpp:182
assembler::AssemblyValsCache pressure_ass_vals_cache
used to store assembly values for pressure for small problems
Definition State.hpp:199
std::vector< mesh::LocalBoundary > local_boundary
mapping from elements to nodes for dirichlet boundary conditions
Definition State.hpp:438
double avg_mass
average system mass, used for contact with IPC
Definition State.hpp:204
bool iso_parametric() const
check if using iso parametric bases
Definition State.cpp:461
std::shared_ptr< assembler::RhsAssembler > build_rhs_assembler() const
build a RhsAssembler for the problem
Definition State.hpp:248
std::vector< mesh::LocalBoundary > total_local_boundary
mapping from elements to nodes for all mesh
Definition State.hpp:436
void solve_transient_navier_stokes(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves transient navier stokes with FEM
void init_solve(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
initialize solver
void init_linear_solve(Eigen::MatrixXd &sol, const double t=1.0)
initialize the linear solve
assembler::AssemblyValsCache mass_ass_vals_cache
Definition State.hpp:197
void solve_tensor_nonlinear(Eigen::MatrixXd &sol, const int t=0, const bool init_lagging=true)
solves nonlinear problems
std::string resolve_input_path(const std::string &path, const bool only_if_exists=false) const
Resolve input path relative to root_path() if the path is not absolute.
std::vector< int > boundary_nodes
list of boundary nodes
Definition State.hpp:432
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:321
std::shared_ptr< polyfem::mesh::MeshNodes > geom_mesh_nodes
Definition State.hpp:193
void solve_linear(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves a linear problem
void init_nonlinear_tensor_solve(Eigen::MatrixXd &sol, const double t=1.0, const bool init_time_integrator=true)
initialize the nonlinear solver
bool is_contact_enabled() const
does the simulation have contact
Definition State.hpp:556
void build_collision_mesh()
extracts the boundary mesh for collision, called in build_basis
Definition State.cpp:1361
Eigen::VectorXi disc_orders
vector of discretization orders, used when not all elements have the same degree, one per element
Definition State.hpp:190
std::vector< mesh::LocalBoundary > local_neumann_boundary
mapping from elements to nodes for neumann boundary conditions
Definition State.hpp:440
std::shared_ptr< assembler::PressureAssembler > build_pressure_assembler() const
Definition State.hpp:256
void solve_transient_linear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves transient linear problem
bool is_problem_linear() const
Returns whether the system is linear. Collisions and pressure add nonlinearity to the problem.
Definition State.hpp:407
std::shared_ptr< assembler::MixedAssembler > mixed_assembler
Definition State.hpp:159
std::map< int, Eigen::MatrixXd > polys
polygons, used since poly have no geom mapping
Definition State.hpp:185
Eigen::MatrixXd rhs
System right-hand side.
Definition State.hpp:207
std::unordered_map< int, std::vector< mesh::LocalBoundary > > local_pressure_cavity
mapping from elements to nodes for pressure boundary conditions
Definition State.hpp:444
StiffnessMatrix basis_nodes_to_gbasis_nodes
Definition State.hpp:706
std::string velocity() const
Definition Units.hpp:29
static double convert(const json &val, const std::string &unit_type)
Definition Units.cpp:31
const std::string & length() const
Definition Units.hpp:19
static bool is_elastic_material(const std::string &material)
utility to check if material is one of the elastic materials
Caches basis evaluation and geometric mapping at every element.
void init(const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const bool is_mass=false)
computes the basis evaluation and geometric mapping for each of the given ElementBases in bases initi...
stores per local bases evaluations
std::vector< basis::Local2Global > global
stores per element basis values at given quadrature points and geometric mapping
void compute(const int el_index, const bool is_volume, const Eigen::MatrixXd &pts, const basis::ElementBases &basis, const basis::ElementBases &gbasis)
computes the per element values at the local (ref el) points (pts) sets basis_values,...
assemble matrix based on the local assembler local assembler is eg Laplace, LinearElasticity etc
void init(const int dim, const json &param)
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 bernstein, const bool serendipity, const bool has_polys, const bool is_geom_bases, const bool use_corner_quadrature, 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 bernstein, const bool serendipity, const bool has_polys, const bool is_geom_bases, const bool use_corner_quadrature, 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).
static int build_bases(const std::string &assembler_name, const int dim, const mesh::Mesh2D &mesh, const int n_bases, const int quadrature_order, const int mass_quadrature_order, std::vector< ElementBases > &bases, std::vector< mesh::LocalBoundary > &local_boundary, std::map< int, Eigen::MatrixXd > &mapped_boundary)
static int build_bases(const assembler::LinearAssembler &assembler, const int n_samples_per_edge, const mesh::Mesh2D &mesh, const int n_bases, const int quadrature_order, const int mass_quadrature_order, const int integral_constraints, std::vector< ElementBases > &bases, const std::vector< ElementBases > &gbases, const std::map< int, InterfaceData > &poly_edge_to_data, std::map< int, Eigen::MatrixXd > &mapped_boundary)
Build bases over the remaining polygons of a mesh.
static int build_bases(const assembler::LinearAssembler &assembler, const int n_samples_per_edge, const mesh::Mesh3D &mesh, const int n_bases, const int quadrature_order, const int mass_quadrature_order, const int integral_constraints, std::vector< ElementBases > &bases, const std::vector< ElementBases > &gbases, const std::map< int, InterfaceData > &poly_face_to_data, std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &mapped_boundary)
Build bases over the remaining polygons of a mesh.
static int build_bases(const mesh::Mesh2D &mesh, const std::string &assembler, const int quadrature_order, const int mass_quadrature_order, std::vector< ElementBases > &bases, std::vector< mesh::LocalBoundary > &local_boundary, std::map< int, InterfaceData > &poly_edge_to_data)
static int build_bases(const mesh::Mesh3D &mesh, const std::string &assembler, const int quadrature_order, const int mass_quadrature_order, std::vector< ElementBases > &bases, std::vector< mesh::LocalBoundary > &local_boundary, std::map< int, InterfaceData > &poly_face_to_data)
static int build_bases(const std::string &assembler_name, const int dim, const mesh::Mesh2D &mesh, const int n_bases, const int quadrature_order, const int mass_quadrature_order, std::vector< ElementBases > &bases, std::vector< mesh::LocalBoundary > &local_boundary, std::map< int, Eigen::MatrixXd > &mapped_boundary)
void build_grid(const polyfem::mesh::Mesh &mesh, const double spacing)
builds the grid to export the solution
Definition OutData.cpp:2419
static void extract_boundary_mesh(const mesh::Mesh &mesh, const int n_bases, const std::vector< basis::ElementBases > &bases, const std::vector< mesh::LocalBoundary > &total_local_boundary, Eigen::MatrixXd &node_positions, Eigen::MatrixXi &boundary_edges, Eigen::MatrixXi &boundary_triangles, std::vector< Eigen::Triplet< double > > &displacement_map_entries)
extracts the boundary mesh
Definition OutData.cpp:145
double assembling_stiffness_mat_time
time to assembly
Definition OutData.hpp:350
double assigning_rhs_time
time to computing the rhs
Definition OutData.hpp:354
double assembling_mass_mat_time
time to assembly mass
Definition OutData.hpp:352
double building_basis_time
time to construct the basis
Definition OutData.hpp:344
double solving_time
time to solve
Definition OutData.hpp:356
double computing_poly_basis_time
time to build the polygonal/polyhedral bases
Definition OutData.hpp:348
int n_flipped
number of flipped elements, compute only when using count_flipped_els (false by default)
Definition OutData.hpp:397
Eigen::Vector4d spectrum
spectrum of the stiffness matrix, enable only if POLYSOLVE_WITH_SPECTRA is ON (off by default)
Definition OutData.hpp:371
void count_flipped_elements(const polyfem::mesh::Mesh &mesh, const std::vector< polyfem::basis::ElementBases > &gbases)
counts the number of flipped elements
Definition OutData.cpp:2604
void compute_mesh_size(const polyfem::mesh::Mesh &mesh_in, const std::vector< polyfem::basis::ElementBases > &bases_in, const int n_samples, const bool use_curved_mesh_size)
computes the mesh size, it samples every edges n_samples times uses curved_mesh_size (false by defaul...
Definition OutData.cpp:2517
long long nn_zero
non zeros and sytem matrix size num dof is the total dof in the system
Definition OutData.hpp:389
double mesh_size
max edge lenght
Definition OutData.hpp:378
void reset()
clears all stats
Definition OutData.cpp:2595
double min_edge_length
min edge lenght
Definition OutData.hpp:380
bool is_volume() const override
checks if mesh is volume
Definition Mesh3D.hpp:28
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
int n_elements() const
utitlity to return the number of elements, cells or faces in 3d and 2d
Definition Mesh.hpp:161
const Eigen::MatrixXi & e() const
Definition Obstacle.hpp:44
const Eigen::MatrixXi & f() const
Definition Obstacle.hpp:43
const Eigen::MatrixXd & v() const
Definition Obstacle.hpp:41
const Eigen::VectorXi & codim_v() const
Definition Obstacle.hpp:42
static void p_refine(const mesh::Mesh &mesh, const double B, const bool h1_formula, const int base_p, const int discr_order_max, io::OutStatsData &stats, Eigen::VectorXi &disc_orders)
compute a priori prefinement
Definition APriori.cpp:242
std::shared_ptr< assembler::RhsAssembler > rhs_assembler
void q_nodes_2d(const int q, Eigen::MatrixXd &val)
void p_nodes_2d(const int p, Eigen::MatrixXd &val)
void p_nodes_3d(const int p, Eigen::MatrixXd &val)
void q_nodes_3d(const int q, Eigen::MatrixXd &val)
bool read_matrix(const std::string &path, Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &mat)
Reads a matrix from a file. Determines the file format based on the path's extension.
Definition MatrixIO.cpp:18
bool write_matrix(const std::string &path, const Mat &mat)
Writes a matrix to a file. Determines the file format based on the path's extension.
Definition MatrixIO.cpp:42
void load_collision_proxy(const std::string &mesh_filename, const std::string &weights_filename, const Eigen::VectorXi &in_node_to_node, const json &transformation, Eigen::MatrixXd &vertices, Eigen::VectorXi &codim_vertices, Eigen::MatrixXi &edges, Eigen::MatrixXi &faces, std::vector< Eigen::Triplet< double > > &displacement_map_entries)
Load a collision proxy mesh and displacement map from files.
void build_collision_proxy(const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &geom_bases, const std::vector< LocalBoundary > &total_local_boundary, const int n_bases, const int dim, const double max_edge_length, Eigen::MatrixXd &proxy_vertices, Eigen::MatrixXi &proxy_faces, std::vector< Eigen::Triplet< double > > &displacement_map_entries, const CollisionProxyTessellation tessellation)
Eigen::SparseMatrix< double > lump_matrix(const Eigen::SparseMatrix< double > &M)
Lump each row of a matrix into the diagonal.
std::string resolve_path(const std::string &path, const std::string &input_file_path, const bool only_if_exists=false)
std::vector< T > json_as_array(const json &j)
Return the value of a json object as an array.
Definition JSONUtils.hpp:38
void append_rows(DstMat &dst, const SrcMat &src)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:44
void compute_integral_constraints(const Mesh3D &mesh, const int n_bases, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, Eigen::MatrixXd &basis_integrals)
Definition State.cpp:388
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:11
nlohmann::json json
Definition Common.hpp:9
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:73
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22