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