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
6
9
18
21
23
26
29
32
35
40
43
44#include <polysolve/linear/FEMSolver.hpp>
45
47
48#include <igl/edges.h>
49#include <igl/Timer.h>
50
51#include <iostream>
52#include <algorithm>
53#include <memory>
54#include <filesystem>
55
57
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
315 const auto primitive_offset = [&](int node) {
316 if (mesh_nodes->is_vertex_node(node))
317 return 0;
318 else if (mesh_nodes->is_edge_node(node))
319 return mesh->n_vertices();
320 else if (mesh_nodes->is_face_node(node))
321 return mesh->n_vertices() + mesh->n_edges();
322 else if (mesh_nodes->is_cell_node(node))
323 return mesh->n_vertices() + mesh->n_edges() + mesh->n_faces();
324 throw std::runtime_error("Invalid node ID!");
325 };
326
327 logger().trace("Building primitive to node mapping...");
328 timer.start();
329 std::vector<std::vector<int>> primitive_to_nodes(num_primitives);
330 const std::vector<int> &grouped_nodes = mesh_nodes->primitive_to_node();
331 int node_count = 0;
332 for (int i = 0; i < grouped_nodes.size(); i++)
333 {
334 int node = grouped_nodes[i];
335 assert(node < num_nodes);
336 if (node >= 0)
337 {
338 int primitive = mesh_nodes->node_to_primitive_gid().at(node) + primitive_offset(i);
339 assert(primitive < num_primitives);
340 primitive_to_nodes[primitive].push_back(node);
341 node_count++;
342 }
343 }
344 assert(node_count == num_nodes);
345 timer.stop();
346 logger().trace("Done (took {}s)", timer.getElapsedTime());
347
348 logger().trace("Combining mappings...");
349 timer.start();
350 in_node_to_node.setConstant(num_nodes, -1);
351 for (int i = 0; i < num_nodes; i++)
352 {
353 // input node id -> input primitive -> primitive -> node(s)
354 const std::vector<int> &possible_nodes =
355 primitive_to_nodes[in_primitive_to_primitive[in_node_to_in_primitive[i]]];
356
357 if (possible_nodes.size() == 1)
358 in_node_to_node[i] = possible_nodes[0];
359 else
360 {
361 assert(possible_nodes.size() > 1);
362
363 // TODO: The following code assumes multiple nodes must come from an edge.
364 // This only true for P3. P4+ has multiple face nodes and P5+ have multiple cell nodes.
365
366 int e_id = in_primitive_to_primitive[in_node_to_in_primitive[i]] - mesh->n_vertices();
367 assert(e_id < mesh->n_edges());
368 assert(in_node_to_node[mesh->edge_vertex(e_id, 0)] >= 0); // Vertex nodes should be computed first
369 RowVectorNd v0 = mesh_nodes->node_position(in_node_to_node[mesh->edge_vertex(e_id, 0)]);
370 RowVectorNd a = mesh_nodes->node_position(possible_nodes[0]);
371 RowVectorNd b = mesh_nodes->node_position(possible_nodes[1]);
372 // Assume possible nodes are ordered, so only need to check order of two nodes
373
374 // Input edges are sorted, so if a is closer to v0 then the order is correct
375 // otherwise the nodes are flipped.
376 assert(mesh->edge_vertex(e_id, 0) < mesh->edge_vertex(e_id, 1));
377 int offset = (a - v0).squaredNorm() < (b - v0).squaredNorm()
378 ? in_node_offset[i]
379 : (possible_nodes.size() - in_node_offset[i] - 1);
380 in_node_to_node[i] = possible_nodes[offset];
381 }
382 }
383 timer.stop();
384 logger().trace("Done (took {}s)", timer.getElapsedTime());
385 }
386
387 std::string State::formulation() const
388 {
389 if (args["materials"].is_null())
390 {
391 logger().error("specify some 'materials'");
392 assert(!args["materials"].is_null());
393 throw "invalid input";
394 }
395
396 if (args["materials"].is_array())
397 {
398 std::string current = "";
399 for (const auto &m : args["materials"])
400 {
401 const std::string tmp = m["type"];
402 if (current.empty())
403 current = tmp;
404 else if (current != tmp)
405 {
406 if (current == "LinearElasticity"
407 || current == "NeoHookean"
408 || current == "SaintVenant"
409 || current == "HookeLinearElasticity"
410 || current == "MooneyRivlin"
411 || current == "MooneyRivlin3Param"
412 || current == "UnconstrainedOgden"
413 || current == "IncompressibleOgden"
414 || current == "MultiModels")
415 {
416 if (tmp == "LinearElasticity"
417 || tmp == "NeoHookean"
418 || tmp == "SaintVenant"
419 || tmp == "HookeLinearElasticity"
420 || tmp == "MooneyRivlin"
421 || current == "MooneyRivlin3Param"
422 || tmp == "UnconstrainedOgden"
423 || tmp == "IncompressibleOgden")
424 current = "MultiModels";
425 else
426 {
427 logger().error("Current material is {}, new material is {}, multimaterial supported only for LinearElasticity and NeoHookean", current, tmp);
428 throw "invalid input";
429 }
430 }
431 else
432 {
433 logger().error("Current material is {}, new material is {}, multimaterial supported only for LinearElasticity and NeoHookean", current, tmp);
434 throw "invalid input";
435 }
436 }
437 }
438
439 return current;
440 }
441 else
442 return args["materials"]["type"];
443 }
444
445 void State::sol_to_pressure(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
446 {
447 if (n_pressure_bases <= 0)
448 {
449 logger().error("No pressure bases defined!");
450 return;
451 }
452
453 assert(mixed_assembler != nullptr);
454 Eigen::MatrixXd tmp = sol;
455
456 int fluid_offset = use_avg_pressure ? (assembler->is_fluid() ? 1 : 0) : 0;
457 sol = tmp.topRows(tmp.rows() - n_pressure_bases - fluid_offset);
458 assert(sol.size() == n_bases * (problem->is_scalar() ? 1 : mesh->dimension()));
459 pressure = tmp.middleRows(tmp.rows() - n_pressure_bases - fluid_offset, n_pressure_bases);
460 assert(pressure.size() == n_pressure_bases);
461 }
462
464 const Mesh3D &mesh,
465 const int n_bases,
466 const std::vector<basis::ElementBases> &bases,
467 const std::vector<basis::ElementBases> &gbases,
468 Eigen::MatrixXd &basis_integrals)
469 {
470 if (!mesh.is_volume())
471 {
472 logger().error("Works only on volumetric meshes!");
473 return;
474 }
475 assert(mesh.is_volume());
476
477 basis_integrals.resize(n_bases, 9);
478 basis_integrals.setZero();
479 Eigen::MatrixXd rhs(n_bases, 9);
480 rhs.setZero();
481
482 const int n_elements = mesh.n_elements();
483 for (int e = 0; e < n_elements; ++e)
484 {
485 // if (mesh.is_polytope(e)) {
486 // continue;
487 // }
488 // ElementAssemblyValues vals = values[e];
489 // const ElementAssemblyValues &gvals = gvalues[e];
491 vals.compute(e, mesh.is_volume(), bases[e], gbases[e]);
492
493 // Computes the discretized integral of the PDE over the element
494 const int n_local_bases = int(vals.basis_values.size());
495 for (int j = 0; j < n_local_bases; ++j)
496 {
497 const AssemblyValues &v = vals.basis_values[j];
498 const double integral_100 = (v.grad_t_m.col(0).array() * vals.det.array() * vals.quadrature.weights.array()).sum();
499 const double integral_010 = (v.grad_t_m.col(1).array() * vals.det.array() * vals.quadrature.weights.array()).sum();
500 const double integral_001 = (v.grad_t_m.col(2).array() * vals.det.array() * vals.quadrature.weights.array()).sum();
501
502 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();
503 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();
504 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();
505
506 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();
507 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();
508 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();
509
510 const double area = (v.val.array() * vals.det.array() * vals.quadrature.weights.array()).sum();
511
512 for (size_t ii = 0; ii < v.global.size(); ++ii)
513 {
514 basis_integrals(v.global[ii].index, 0) += integral_100 * v.global[ii].val;
515 basis_integrals(v.global[ii].index, 1) += integral_010 * v.global[ii].val;
516 basis_integrals(v.global[ii].index, 2) += integral_001 * v.global[ii].val;
517
518 basis_integrals(v.global[ii].index, 3) += integral_110 * v.global[ii].val;
519 basis_integrals(v.global[ii].index, 4) += integral_011 * v.global[ii].val;
520 basis_integrals(v.global[ii].index, 5) += integral_101 * v.global[ii].val;
521
522 basis_integrals(v.global[ii].index, 6) += integral_200 * v.global[ii].val;
523 basis_integrals(v.global[ii].index, 7) += integral_020 * v.global[ii].val;
524 basis_integrals(v.global[ii].index, 8) += integral_002 * v.global[ii].val;
525
526 rhs(v.global[ii].index, 6) += -2.0 * area * v.global[ii].val;
527 rhs(v.global[ii].index, 7) += -2.0 * area * v.global[ii].val;
528 rhs(v.global[ii].index, 8) += -2.0 * area * v.global[ii].val;
529 }
530 }
531 }
532
533 basis_integrals -= rhs;
534 }
535
537 {
538 if (mesh->has_poly())
539 return true;
540
541 if (args["space"]["basis_type"] == "Spline")
542 return true;
543
544 if (mesh->is_rational())
545 return false;
546
547 if (args["space"]["use_p_ref"])
548 return false;
549
550 if (has_periodic_bc())
551 return false;
552
554 return false;
555
556 if (mesh->orders().size() <= 0)
557 {
558 if (args["space"]["discr_order"] == 1)
559 return true;
560 else
561 return args["space"]["advanced"]["isoparametric"];
562 }
563
564 if (mesh->orders().minCoeff() != mesh->orders().maxCoeff())
565 return false;
566
567 if (args["space"]["discr_order"] == mesh->orders().minCoeff())
568 return true;
569
570 // TODO:
571 // if (args["space"]["discr_order"] == 1 && args["force_linear_geometry"])
572 // return true;
573
574 return args["space"]["advanced"]["isoparametric"];
575 }
576
578 {
579 if (!mesh)
580 {
581 logger().error("Load the mesh first!");
582 return;
583 }
584
585 mesh->prepare_mesh();
586
587 bases.clear();
588 pressure_bases.clear();
589 geom_bases_.clear();
590 boundary_nodes.clear();
591 dirichlet_nodes.clear();
592 neumann_nodes.clear();
593 local_boundary.clear();
594 total_local_boundary.clear();
597 local_pressure_cavity.clear();
598 polys.clear();
599 poly_edge_to_data.clear();
600 rhs.resize(0, 0);
601
602 if (assembler::MultiModel *mm = dynamic_cast<assembler::MultiModel *>(assembler.get()))
603 {
604 assert(args["materials"].is_array());
605
606 std::vector<std::string> materials(mesh->n_elements());
607
608 std::map<int, std::string> mats;
609
610 for (const auto &m : args["materials"])
611 mats[m["id"].get<int>()] = m["type"];
612
613 for (int i = 0; i < materials.size(); ++i)
614 materials[i] = mats.at(mesh->get_body_id(i));
615
616 mm->init_multimodels(materials);
617 }
618
619 n_bases = 0;
620 n_geom_bases = 0;
622
623 stats.reset();
624
625 disc_orders.resize(mesh->n_elements());
626
627 problem->init(*mesh);
628 logger().info("Building {} basis...", (iso_parametric() ? "isoparametric" : "not isoparametric"));
629 const bool has_polys = mesh->has_poly();
630
631 local_boundary.clear();
634 local_pressure_cavity.clear();
635 std::map<int, basis::InterfaceData> poly_edge_to_data_geom; // temp dummy variable
636
637 const auto &tmp_json = args["space"]["discr_order"];
638 if (tmp_json.is_number_integer())
639 {
640 disc_orders.setConstant(tmp_json);
641 }
642 else if (tmp_json.is_string())
643 {
644 const std::string discr_orders_path = tmp_json;
645 Eigen::MatrixXi tmp;
646 read_matrix(discr_orders_path, tmp);
647 assert(tmp.size() == disc_orders.size());
648 assert(tmp.cols() == 1);
649 disc_orders = tmp;
650 }
651 else if (tmp_json.is_array())
652 {
653 const auto b_discr_orders = tmp_json;
654
655 std::map<int, int> b_orders;
656 for (size_t i = 0; i < b_discr_orders.size(); ++i)
657 {
658 assert(b_discr_orders[i]["id"].is_array() || b_discr_orders[i]["id"].is_number_integer());
659
660 const int order = b_discr_orders[i]["order"];
661 for (const int id : json_as_array<int>(b_discr_orders[i]["id"]))
662 {
663 b_orders[id] = order;
664 logger().trace("bid {}, discr {}", id, order);
665 }
666 }
667
668 for (int e = 0; e < mesh->n_elements(); ++e)
669 {
670 const int bid = mesh->get_body_id(e);
671 const auto order = b_orders.find(bid);
672 if (order == b_orders.end())
673 {
674 logger().debug("Missing discretization order for body {}; using 1", bid);
675 b_orders[bid] = 1;
676 disc_orders[e] = 1;
677 }
678 else
679 {
680 disc_orders[e] = order->second;
681 }
682 }
683 }
684 else
685 {
686 logger().error("space/discr_order must be either a number a path or an array");
687 throw "invalid json";
688 }
689 // TODO: same for pressure!
690
691 Eigen::MatrixXi geom_disc_orders;
692 if (!iso_parametric())
693 {
694 if (mesh->orders().size() <= 0)
695 {
696 geom_disc_orders.resizeLike(disc_orders);
697 geom_disc_orders.setConstant(1);
698 }
699 else
700 geom_disc_orders = mesh->orders();
701 }
702
703 igl::Timer timer;
704 timer.start();
705 if (args["space"]["use_p_ref"])
706 {
708 *mesh,
709 args["space"]["advanced"]["B"],
710 args["space"]["advanced"]["h1_formula"],
711 args["space"]["discr_order"],
712 args["space"]["advanced"]["discr_order_max"],
713 stats,
715
716 logger().info("min p: {} max p: {}", disc_orders.minCoeff(), disc_orders.maxCoeff());
717 }
718
719 int quadrature_order = args["space"]["advanced"]["quadrature_order"].get<int>();
720 const int mass_quadrature_order = args["space"]["advanced"]["mass_quadrature_order"].get<int>();
721 if (mixed_assembler != nullptr)
722 {
723 const int disc_order = disc_orders.maxCoeff();
724 if (disc_order - disc_orders.minCoeff() != 0)
725 {
726 logger().error("p refinement not supported in mixed formulation!");
727 return;
728 }
729 }
730
731 // shape optimization needs continuous geometric basis
732 // const bool use_continuous_gbasis = optimization_enabled == solver::CacheLevel::Derivatives;
733 const bool use_continuous_gbasis = true;
734
735 if (mesh->is_volume())
736 {
737 const Mesh3D &tmp_mesh = *dynamic_cast<Mesh3D *>(mesh.get());
738 if (args["space"]["basis_type"] == "Spline")
739 {
740 // if (!iso_parametric())
741 // {
742 // logger().error("Splines must be isoparametric, ignoring...");
743 // // LagrangeBasis3d::build_bases(tmp_mesh, quadrature_order, geom_disc_orders, has_polys, geom_bases_, local_boundary, poly_edge_to_data_geom, mesh_nodes);
744 // SplineBasis3d::build_bases(tmp_mesh, quadrature_order, geom_bases_, local_boundary, poly_edge_to_data);
745 // }
746
747 n_bases = basis::SplineBasis3d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, bases, local_boundary, poly_edge_to_data);
748
749 // if (iso_parametric() && args["fit_nodes"])
750 // SplineBasis3d::fit_nodes(tmp_mesh, n_bases, bases);
751 }
752 else
753 {
754 if (!iso_parametric())
755 n_geom_bases = basis::LagrangeBasis3d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, geom_disc_orders, false, has_polys, !use_continuous_gbasis, geom_bases_, local_boundary, poly_edge_to_data_geom, geom_mesh_nodes);
756
757 n_bases = basis::LagrangeBasis3d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, disc_orders, args["space"]["basis_type"] == "Serendipity", has_polys, false, bases, local_boundary, poly_edge_to_data, mesh_nodes);
758 }
759
760 // if(problem->is_mixed())
761 if (mixed_assembler != nullptr)
762 {
763 n_pressure_bases = basis::LagrangeBasis3d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, int(args["space"]["pressure_discr_order"]), false, has_polys, false, pressure_bases, local_boundary, poly_edge_to_data_geom, pressure_mesh_nodes);
764 }
765 }
766 else
767 {
768 const Mesh2D &tmp_mesh = *dynamic_cast<Mesh2D *>(mesh.get());
769 if (args["space"]["basis_type"] == "Spline")
770 {
771 // TODO:
772 // if (!iso_parametric())
773 // {
774 // logger().error("Splines must be isoparametric, ignoring...");
775 // // LagrangeBasis2d::build_bases(tmp_mesh, quadrature_order, disc_orders, has_polys, geom_bases_, local_boundary, poly_edge_to_data_geom, mesh_nodes);
776 // n_bases = SplineBasis2d::build_bases(tmp_mesh, quadrature_order, geom_bases_, local_boundary, poly_edge_to_data);
777 // }
778
779 n_bases = basis::SplineBasis2d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, bases, local_boundary, poly_edge_to_data);
780
781 // if (iso_parametric() && args["fit_nodes"])
782 // SplineBasis2d::fit_nodes(tmp_mesh, n_bases, bases);
783 }
784 else
785 {
786 if (!iso_parametric())
787 n_geom_bases = basis::LagrangeBasis2d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, geom_disc_orders, false, has_polys, !use_continuous_gbasis, geom_bases_, local_boundary, poly_edge_to_data_geom, geom_mesh_nodes);
788
789 n_bases = basis::LagrangeBasis2d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, disc_orders, args["space"]["basis_type"] == "Serendipity", has_polys, false, bases, local_boundary, poly_edge_to_data, mesh_nodes);
790 }
791
792 // if(problem->is_mixed())
793 if (mixed_assembler != nullptr)
794 {
795 n_pressure_bases = basis::LagrangeBasis2d::build_bases(tmp_mesh, assembler->name(), quadrature_order, mass_quadrature_order, int(args["space"]["pressure_discr_order"]), false, has_polys, false, pressure_bases, local_boundary, poly_edge_to_data_geom, pressure_mesh_nodes);
796 }
797 }
798
799 if (mixed_assembler != nullptr)
800 {
801 assert(bases.size() == pressure_bases.size());
802 for (int i = 0; i < pressure_bases.size(); ++i)
803 {
805 bases[i].compute_quadrature(b_quad);
806 pressure_bases[i].set_quadrature([b_quad](quadrature::Quadrature &quad) { quad = b_quad; });
807 }
808 }
809
810 timer.stop();
811
813
814 if (n_geom_bases == 0)
816
817 auto &gbases = geom_bases();
818
820 {
821 std::map<std::array<int, 2>, double> pairs;
822 for (int e = 0; e < gbases.size(); e++)
823 {
824 const auto &gbs = gbases[e].bases;
825 const auto &bs = bases[e].bases;
826
827 Eigen::MatrixXd local_pts;
828 const int order = bs.front().order();
829 if (mesh->is_volume())
830 {
831 if (mesh->is_simplex(e))
832 autogen::p_nodes_3d(order, local_pts);
833 else
834 autogen::q_nodes_3d(order, local_pts);
835 }
836 else
837 {
838 if (mesh->is_simplex(e))
839 autogen::p_nodes_2d(order, local_pts);
840 else
841 autogen::q_nodes_2d(order, local_pts);
842 }
843
845 vals.compute(e, mesh->is_volume(), local_pts, gbases[e], gbases[e]);
846
847 for (int i = 0; i < bs.size(); i++)
848 {
849 for (int j = 0; j < gbs.size(); j++)
850 {
851 if (std::abs(vals.basis_values[j].val(i)) > 1e-7)
852 {
853 std::array<int, 2> index = {{gbs[j].global()[0].index, bs[i].global()[0].index}};
854 pairs.insert({index, vals.basis_values[j].val(i)});
855 }
856 }
857 }
858 }
859
860 const int dim = mesh->dimension();
861 std::vector<Eigen::Triplet<double>> coeffs;
862 coeffs.clear();
863 for (const auto &iter : pairs)
864 for (int d = 0; d < dim; d++)
865 coeffs.emplace_back(iter.first[0] * dim + d, iter.first[1] * dim + d, iter.second);
866
867 gbasis_nodes_to_basis_nodes.resize(n_geom_bases * mesh->dimension(), n_bases * mesh->dimension());
868 gbasis_nodes_to_basis_nodes.setFromTriplets(coeffs.begin(), coeffs.end());
869 }
870
871 for (const auto &lb : local_boundary)
872 total_local_boundary.emplace_back(lb);
873
874 const int dim = mesh->dimension();
875 const int problem_dim = problem->is_scalar() ? 1 : dim;
876
877 // handle periodic bc
878 if (has_periodic_bc())
879 {
880 // collect periodic directions
881 json directions = args["boundary_conditions"]["periodic_boundary"]["correspondence"];
882 Eigen::MatrixXd tile_offset = Eigen::MatrixXd::Identity(dim, dim);
883
884 if (directions.size() > 0)
885 {
886 Eigen::VectorXd tmp;
887 for (int d = 0; d < dim; d++)
888 {
889 tmp = directions[d];
890 if (tmp.size() != dim)
891 log_and_throw_error("Invalid size of periodic directions!");
892 tile_offset.col(d) = tmp;
893 }
894 }
895
896 periodic_bc = std::make_shared<PeriodicBoundary>(problem->is_scalar(), n_bases, bases, mesh_nodes, tile_offset, args["boundary_conditions"]["periodic_boundary"]["tolerance"].get<double>());
897
898 macro_strain_constraint.init(dim, args["boundary_conditions"]["periodic_boundary"]);
899 }
900
901 if (args["space"]["advanced"]["count_flipped_els"])
903
904 const int prev_bases = n_bases;
906
907 {
908 igl::Timer timer2;
909 logger().debug("Building node mapping...");
910 timer2.start();
912 problem->update_nodes(in_node_to_node);
913 mesh->update_nodes(in_node_to_node);
914 timer2.stop();
915 logger().debug("Done (took {}s)", timer2.getElapsedTime());
916 }
917
918 logger().info("Building collision mesh...");
920 if (periodic_bc && args["contact"]["periodic"])
922 logger().info("Done!");
923
924 const int prev_b_size = local_boundary.size();
925 problem->setup_bc(*mesh, n_bases - obstacle.n_vertices(),
934
935 // setp nodal values
936 {
938 for (int n = 0; n < dirichlet_nodes.size(); ++n)
939 {
940 const int n_id = dirichlet_nodes[n];
941 bool found = false;
942 for (const auto &bs : bases)
943 {
944 for (const auto &b : bs.bases)
945 {
946 for (const auto &lg : b.global())
947 {
948 if (lg.index == n_id)
949 {
950 dirichlet_nodes_position[n] = lg.node;
951 found = true;
952 break;
953 }
954 }
955
956 if (found)
957 break;
958 }
959
960 if (found)
961 break;
962 }
963
964 assert(found);
965 }
966
968 for (int n = 0; n < neumann_nodes.size(); ++n)
969 {
970 const int n_id = neumann_nodes[n];
971 bool found = false;
972 for (const auto &bs : bases)
973 {
974 for (const auto &b : bs.bases)
975 {
976 for (const auto &lg : b.global())
977 {
978 if (lg.index == n_id)
979 {
980 neumann_nodes_position[n] = lg.node;
981 found = true;
982 break;
983 }
984 }
985
986 if (found)
987 break;
988 }
989
990 if (found)
991 break;
992 }
993
994 assert(found);
995 }
996 }
997
998 const bool has_neumann = local_neumann_boundary.size() > 0 || local_boundary.size() < prev_b_size;
999 use_avg_pressure = !has_neumann;
1000
1001 for (int i = prev_bases; i < n_bases; ++i)
1002 {
1003 for (int d = 0; d < problem_dim; ++d)
1004 boundary_nodes.push_back(i * problem_dim + d);
1005 }
1006
1007 std::sort(boundary_nodes.begin(), boundary_nodes.end());
1008 auto it = std::unique(boundary_nodes.begin(), boundary_nodes.end());
1009 boundary_nodes.resize(std::distance(boundary_nodes.begin(), it));
1010
1011 // for elastic pure periodic problem, find an internal node and force zero dirichlet
1012 if ((!problem->is_time_dependent() || args["time"]["quasistatic"]) && boundary_nodes.size() == 0 && has_periodic_bc())
1013 {
1014 // find an internal node to force zero dirichlet
1015 std::vector<bool> isboundary(n_bases, false);
1016 for (const auto &lb : total_local_boundary)
1017 {
1018 const int e = lb.element_id();
1019 for (int i = 0; i < lb.size(); ++i)
1020 {
1021 const auto nodes = bases[e].local_nodes_for_primitive(lb.global_primitive_id(i), *mesh);
1022
1023 for (int n : nodes)
1024 isboundary[bases[e].bases[n].global()[0].index] = true;
1025 }
1026 }
1027 int i = 0;
1028 for (; i < n_bases; i++)
1029 if (!isboundary[i]) // (!periodic_bc->is_periodic_dof(i))
1030 break;
1031 if (i >= n_bases)
1032 log_and_throw_error("Failed to find a non-periodic node!");
1033 const int actual_dim = problem->is_scalar() ? 1 : mesh->dimension();
1034 for (int d = 0; d < actual_dim; d++)
1035 {
1036 boundary_nodes.push_back(i * actual_dim + d);
1037 }
1038 logger().info("Fix solution at node {} to remove singularity due to periodic BC", i);
1039 }
1040
1041 const auto &curret_bases = geom_bases();
1042 const int n_samples = 10;
1043 stats.compute_mesh_size(*mesh, curret_bases, n_samples, args["output"]["advanced"]["curved_mesh_size"]);
1045 {
1047 }
1049 {
1051 }
1052
1053 if (is_contact_enabled())
1054 {
1055 min_boundary_edge_length = std::numeric_limits<double>::max();
1056 for (const auto &edge : collision_mesh.edges().rowwise())
1057 {
1058 const VectorNd v0 = collision_mesh.rest_positions().row(edge(0));
1059 const VectorNd v1 = collision_mesh.rest_positions().row(edge(1));
1060 min_boundary_edge_length = std::min(min_boundary_edge_length, (v1 - v0).norm());
1061 }
1062
1063 double dhat = Units::convert(args["contact"]["dhat"], units.length());
1064 args["contact"]["epsv"] = Units::convert(args["contact"]["epsv"], units.velocity());
1065 args["contact"]["dhat"] = dhat;
1066
1067 if (!has_dhat && dhat > min_boundary_edge_length)
1068 {
1069 args["contact"]["dhat"] = double(args["contact"]["dhat_percentage"]) * min_boundary_edge_length;
1070 logger().info("dhat set to {}", double(args["contact"]["dhat"]));
1071 }
1072 else
1073 {
1074 if (dhat > min_boundary_edge_length)
1075 logger().warn("dhat larger than min boundary edge, {} > {}", dhat, min_boundary_edge_length);
1076 }
1077 }
1078
1079 logger().info("n_bases {}", n_bases);
1080
1081 timings.building_basis_time = timer.getElapsedTime();
1082 logger().info(" took {}s", timings.building_basis_time);
1083
1084 logger().info("flipped elements {}", stats.n_flipped);
1085 logger().info("h: {}", stats.mesh_size);
1086 logger().info("n bases: {}", n_bases);
1087 logger().info("n pressure bases: {}", n_pressure_bases);
1088
1091 if (n_bases <= args["solver"]["advanced"]["cache_size"])
1092 {
1093 timer.start();
1094 logger().info("Building cache...");
1095 ass_vals_cache.init(mesh->is_volume(), bases, curret_bases);
1096 mass_ass_vals_cache.init(mesh->is_volume(), bases, curret_bases, true);
1097 if (mixed_assembler != nullptr)
1098 pressure_ass_vals_cache.init(mesh->is_volume(), pressure_bases, curret_bases);
1099
1100 logger().info(" took {}s", timer.getElapsedTime());
1101 }
1102
1103 out_geom.build_grid(*mesh, args["output"]["advanced"]["sol_on_grid"]);
1104
1105 if ((!problem->is_time_dependent() || args["time"]["quasistatic"]) && boundary_nodes.empty())
1106 {
1107 if (has_periodic_bc())
1108 logger().warn("(Quasi-)Static problem without Dirichlet nodes, will fix solution at one node to find a unique solution!");
1109 else
1110 log_and_throw_error("Static problem need to have some Dirichlet nodes!");
1111 }
1112 }
1113
1115 {
1116 if (!mesh)
1117 {
1118 logger().error("Load the mesh first!");
1119 return;
1120 }
1121
1122 rhs.resize(0, 0);
1123
1124 if (poly_edge_to_data.empty() && polys.empty())
1125 {
1127 return;
1128 }
1129
1130 igl::Timer timer;
1131 timer.start();
1132 logger().info("Computing polygonal basis...");
1133
1134 // std::sort(boundary_nodes.begin(), boundary_nodes.end());
1135
1136 // mixed not supports polygonal bases
1137 assert(n_pressure_bases == 0 || poly_edge_to_data.size() == 0);
1138
1139 int new_bases = 0;
1140
1141 if (iso_parametric())
1142 {
1143 if (mesh->is_volume())
1144 {
1145 if (args["space"]["poly_basis_type"] == "MeanValue" || args["space"]["poly_basis_type"] == "Wachspress")
1146 logger().error("Barycentric bases not supported in 3D");
1147 assert(assembler->is_linear());
1149 *dynamic_cast<LinearAssembler *>(assembler.get()),
1150 args["space"]["advanced"]["n_harmonic_samples"],
1151 *dynamic_cast<Mesh3D *>(mesh.get()),
1152 n_bases,
1153 args["space"]["advanced"]["quadrature_order"],
1154 args["space"]["advanced"]["mass_quadrature_order"],
1155 args["space"]["advanced"]["integral_constraints"],
1156 bases,
1157 bases,
1159 polys_3d);
1160 }
1161 else
1162 {
1163 if (args["space"]["poly_basis_type"] == "MeanValue")
1164 {
1166 assembler->name(),
1167 assembler->is_tensor() ? 2 : 1,
1168 *dynamic_cast<Mesh2D *>(mesh.get()),
1169 n_bases,
1170 args["space"]["advanced"]["quadrature_order"],
1171 args["space"]["advanced"]["mass_quadrature_order"],
1173 }
1174 else if (args["space"]["poly_basis_type"] == "Wachspress")
1175 {
1177 assembler->name(),
1178 assembler->is_tensor() ? 2 : 1,
1179 *dynamic_cast<Mesh2D *>(mesh.get()),
1180 n_bases,
1181 args["space"]["advanced"]["quadrature_order"],
1182 args["space"]["advanced"]["mass_quadrature_order"],
1184 }
1185 else
1186 {
1187 assert(assembler->is_linear());
1189 *dynamic_cast<LinearAssembler *>(assembler.get()),
1190 args["space"]["advanced"]["n_harmonic_samples"],
1191 *dynamic_cast<Mesh2D *>(mesh.get()),
1192 n_bases,
1193 args["space"]["advanced"]["quadrature_order"],
1194 args["space"]["advanced"]["mass_quadrature_order"],
1195 args["space"]["advanced"]["integral_constraints"],
1196 bases,
1197 bases,
1199 polys);
1200 }
1201 }
1202 }
1203 else
1204 {
1205 if (mesh->is_volume())
1206 {
1207 if (args["space"]["poly_basis_type"] == "MeanValue" || args["space"]["poly_basis_type"] == "Wachspress")
1208 {
1209 logger().error("Barycentric bases not supported in 3D");
1210 throw "not implemented";
1211 }
1212 else
1213 {
1214 assert(assembler->is_linear());
1216 *dynamic_cast<LinearAssembler *>(assembler.get()),
1217 args["space"]["advanced"]["n_harmonic_samples"],
1218 *dynamic_cast<Mesh3D *>(mesh.get()),
1219 n_bases,
1220 args["space"]["advanced"]["quadrature_order"],
1221 args["space"]["advanced"]["mass_quadrature_order"],
1222 args["space"]["advanced"]["integral_constraints"],
1223 bases,
1226 polys_3d);
1227 }
1228 }
1229 else
1230 {
1231 if (args["space"]["poly_basis_type"] == "MeanValue")
1232 {
1234 assembler->name(),
1235 assembler->is_tensor() ? 2 : 1,
1236 *dynamic_cast<Mesh2D *>(mesh.get()),
1237 n_bases, args["space"]["advanced"]["quadrature_order"],
1238 args["space"]["advanced"]["mass_quadrature_order"],
1240 }
1241 else if (args["space"]["poly_basis_type"] == "Wachspress")
1242 {
1244 assembler->name(),
1245 assembler->is_tensor() ? 2 : 1,
1246 *dynamic_cast<Mesh2D *>(mesh.get()),
1247 n_bases, args["space"]["advanced"]["quadrature_order"],
1248 args["space"]["advanced"]["mass_quadrature_order"],
1250 }
1251 else
1252 {
1253 assert(assembler->is_linear());
1255 *dynamic_cast<LinearAssembler *>(assembler.get()),
1256 args["space"]["advanced"]["n_harmonic_samples"],
1257 *dynamic_cast<Mesh2D *>(mesh.get()),
1258 n_bases,
1259 args["space"]["advanced"]["quadrature_order"],
1260 args["space"]["advanced"]["mass_quadrature_order"],
1261 args["space"]["advanced"]["integral_constraints"],
1262 bases,
1265 polys);
1266 }
1267 }
1268 }
1269
1270 timer.stop();
1271 timings.computing_poly_basis_time = timer.getElapsedTime();
1272 logger().info(" took {}s", timings.computing_poly_basis_time);
1273
1274 n_bases += new_bases;
1275 }
1276
1278 {
1279 assert(!mesh->is_volume());
1280 const int dim = mesh->dimension();
1281 const int n_tiles = 2;
1282
1283 if (mesh->dimension() != 2)
1284 log_and_throw_error("Periodic collision mesh is only implemented in 2D!");
1285
1286 Eigen::MatrixXd V(n_bases, dim);
1287 for (const auto &bs : bases)
1288 for (const auto &b : bs.bases)
1289 for (const auto &g : b.global())
1290 V.row(g.index) = g.node;
1291
1292 Eigen::MatrixXi E = collision_mesh.edges();
1293 for (int i = 0; i < E.size(); i++)
1294 E(i) = collision_mesh.to_full_vertex_id(E(i));
1295
1296 Eigen::MatrixXd bbox(V.cols(), 2);
1297 bbox.col(0) = V.colwise().minCoeff();
1298 bbox.col(1) = V.colwise().maxCoeff();
1299
1300 // remove boundary edges on periodic BC, buggy
1301 {
1302 std::vector<int> ind;
1303 for (int i = 0; i < E.rows(); i++)
1304 {
1305 if (!periodic_bc->is_periodic_dof(E(i, 0)) || !periodic_bc->is_periodic_dof(E(i, 1)))
1306 ind.push_back(i);
1307 }
1308 E = E(ind, Eigen::all).eval();
1309 }
1310
1311 Eigen::MatrixXd Vtmp, Vnew;
1312 Eigen::MatrixXi Etmp, Enew;
1313 Vtmp.setZero(V.rows() * n_tiles * n_tiles, V.cols());
1314 Etmp.setZero(E.rows() * n_tiles * n_tiles, E.cols());
1315
1316 Eigen::MatrixXd tile_offset = periodic_bc->get_affine_matrix();
1317
1318 for (int i = 0, idx = 0; i < n_tiles; i++)
1319 {
1320 for (int j = 0; j < n_tiles; j++)
1321 {
1322 Eigen::Vector2d block_id;
1323 block_id << i, j;
1324
1325 Vtmp.middleRows(idx * V.rows(), V.rows()) = V;
1326 // Vtmp.block(idx * V.rows(), 0, V.rows(), 1).array() += tile_offset(0) * i;
1327 // Vtmp.block(idx * V.rows(), 1, V.rows(), 1).array() += tile_offset(1) * j;
1328 for (int vid = 0; vid < V.rows(); vid++)
1329 Vtmp.block(idx * V.rows() + vid, 0, 1, 2) += (tile_offset * block_id).transpose();
1330
1331 Etmp.middleRows(idx * E.rows(), E.rows()) = E.array() + idx * V.rows();
1332 idx += 1;
1333 }
1334 }
1335
1336 // clean duplicated vertices
1337 Eigen::VectorXi indices;
1338 {
1339 std::vector<int> tmp;
1340 for (int i = 0; i < V.rows(); i++)
1341 {
1342 if (periodic_bc->is_periodic_dof(i))
1343 tmp.push_back(i);
1344 }
1345
1346 indices.resize(tmp.size() * n_tiles * n_tiles);
1347 for (int i = 0; i < n_tiles * n_tiles; i++)
1348 {
1349 indices.segment(i * tmp.size(), tmp.size()) = Eigen::Map<Eigen::VectorXi, Eigen::Unaligned>(tmp.data(), tmp.size());
1350 indices.segment(i * tmp.size(), tmp.size()).array() += i * V.rows();
1351 }
1352 }
1353
1354 Eigen::VectorXi potentially_duplicate_mask(Vtmp.rows());
1355 potentially_duplicate_mask.setZero();
1356 potentially_duplicate_mask(indices).array() = 1;
1357 Eigen::MatrixXd candidates = Vtmp(indices, Eigen::all);
1358
1359 Eigen::VectorXi SVI;
1360 std::vector<int> SVJ;
1361 SVI.setConstant(Vtmp.rows(), -1);
1362 int id = 0;
1363 const double eps = (bbox.col(1) - bbox.col(0)).maxCoeff() * args["boundary_conditions"]["periodic_boundary"]["tolerance"].get<double>();
1364 for (int i = 0; i < Vtmp.rows(); i++)
1365 {
1366 if (SVI[i] < 0)
1367 {
1368 SVI[i] = id;
1369 SVJ.push_back(i);
1370 if (potentially_duplicate_mask(i))
1371 {
1372 Eigen::VectorXd diffs = (candidates.rowwise() - Vtmp.row(i)).rowwise().norm();
1373 for (int j = 0; j < diffs.size(); j++)
1374 if (diffs(j) < eps)
1375 SVI[indices[j]] = id;
1376 }
1377 id++;
1378 }
1379 }
1380 Vnew = Vtmp(SVJ, Eigen::all);
1381
1382 Enew.resizeLike(Etmp);
1383 for (int d = 0; d < Etmp.cols(); d++)
1384 Enew.col(d) = SVI(Etmp.col(d));
1385
1386 std::vector<bool> is_on_surface = ipc::CollisionMesh::construct_is_on_surface(Vnew.rows(), Enew);
1387
1388 Eigen::MatrixXi boundary_triangles;
1389 Eigen::SparseMatrix<double> displacement_map;
1390 periodic_collision_mesh = ipc::CollisionMesh(is_on_surface,
1391 Vnew,
1392 Enew,
1393 boundary_triangles,
1394 displacement_map);
1395
1396 periodic_collision_mesh.init_area_jacobians();
1397
1398 periodic_collision_mesh_to_basis.setConstant(Vnew.rows(), -1);
1399 for (int i = 0; i < V.rows(); i++)
1400 for (int j = 0; j < n_tiles * n_tiles; j++)
1401 periodic_collision_mesh_to_basis(SVI[j * V.rows() + i]) = i;
1402
1403 if (periodic_collision_mesh_to_basis.maxCoeff() + 1 != V.rows())
1404 log_and_throw_error("Failed to tile mesh!");
1405 }
1406
1408 {
1411 args, [this](const std::string &p) { return resolve_input_path(p); },
1413 }
1414
1416 const mesh::Mesh &mesh,
1417 const int n_bases,
1418 const std::vector<basis::ElementBases> &bases,
1419 const std::vector<basis::ElementBases> &geom_bases,
1420 const std::vector<mesh::LocalBoundary> &total_local_boundary,
1421 const mesh::Obstacle &obstacle,
1422 const json &args,
1423 const std::function<std::string(const std::string &)> &resolve_input_path,
1424 const Eigen::VectorXi &in_node_to_node,
1425 ipc::CollisionMesh &collision_mesh)
1426 {
1427 Eigen::MatrixXd collision_vertices;
1428 Eigen::VectorXi collision_codim_vids;
1429 Eigen::MatrixXi collision_edges, collision_triangles;
1430 std::vector<Eigen::Triplet<double>> displacement_map_entries;
1431
1432 if (args.contains("/contact/collision_mesh"_json_pointer)
1433 && args.at("/contact/collision_mesh/enabled"_json_pointer).get<bool>())
1434 {
1435 const json collision_mesh_args = args.at("/contact/collision_mesh"_json_pointer);
1436 if (collision_mesh_args.contains("linear_map"))
1437 {
1438 assert(displacement_map_entries.empty());
1439 assert(collision_mesh_args.contains("mesh"));
1440 const std::string root_path = utils::json_value<std::string>(args, "root_path", "");
1441 // TODO: handle transformation per geometry
1442 const json transformation = json_as_array(args["geometry"])[0]["transformation"];
1444 utils::resolve_path(collision_mesh_args["mesh"], root_path),
1445 utils::resolve_path(collision_mesh_args["linear_map"], root_path),
1446 in_node_to_node, transformation, collision_vertices, collision_codim_vids,
1447 collision_edges, collision_triangles, displacement_map_entries);
1448 }
1449 else
1450 {
1451 assert(collision_mesh_args.contains("max_edge_length"));
1452 logger().debug(
1453 "Building collision proxy with max edge length={} ...",
1454 collision_mesh_args["max_edge_length"].get<double>());
1455 igl::Timer timer;
1456 timer.start();
1459 collision_mesh_args["max_edge_length"], collision_vertices,
1460 collision_triangles, displacement_map_entries,
1461 collision_mesh_args["tessellation_type"]);
1462 if (collision_triangles.size())
1463 igl::edges(collision_triangles, collision_edges);
1464 timer.stop();
1465 logger().debug(fmt::format(
1466 std::locale("en_US.UTF-8"),
1467 "Done (took {:g}s, {:L} vertices, {:L} triangles)",
1468 timer.getElapsedTime(),
1469 collision_vertices.rows(), collision_triangles.rows()));
1470 }
1471 }
1472 else
1473 {
1476 collision_vertices, collision_edges, collision_triangles, displacement_map_entries);
1477 }
1478
1479 // n_bases already contains the obstacle vertices
1480 const int num_fe_nodes = n_bases - obstacle.n_vertices();
1481 const int num_fe_collision_vertices = collision_vertices.rows();
1482 assert(collision_edges.size() == 0 || collision_edges.maxCoeff() < num_fe_collision_vertices);
1483 assert(collision_triangles.size() == 0 || collision_triangles.maxCoeff() < num_fe_collision_vertices);
1484
1485 // Append the obstacles to the collision mesh
1486 if (obstacle.n_vertices() > 0)
1487 {
1488 append_rows(collision_vertices, obstacle.v());
1489 append_rows(collision_codim_vids, obstacle.codim_v().array() + num_fe_collision_vertices);
1490 append_rows(collision_edges, obstacle.e().array() + num_fe_collision_vertices);
1491 append_rows(collision_triangles, obstacle.f().array() + num_fe_collision_vertices);
1492
1493 if (!displacement_map_entries.empty())
1494 {
1495 displacement_map_entries.reserve(displacement_map_entries.size() + obstacle.n_vertices());
1496 for (int i = 0; i < obstacle.n_vertices(); i++)
1497 {
1498 displacement_map_entries.emplace_back(num_fe_collision_vertices + i, num_fe_nodes + i, 1.0);
1499 }
1500 }
1501 }
1502
1503 std::vector<bool> is_on_surface = ipc::CollisionMesh::construct_is_on_surface(
1504 collision_vertices.rows(), collision_edges);
1505 for (const int vid : collision_codim_vids)
1506 {
1507 is_on_surface[vid] = true;
1508 }
1509
1510 Eigen::SparseMatrix<double> displacement_map;
1511 if (!displacement_map_entries.empty())
1512 {
1513 displacement_map.resize(collision_vertices.rows(), n_bases);
1514 displacement_map.setFromTriplets(displacement_map_entries.begin(), displacement_map_entries.end());
1515 }
1516
1517 collision_mesh = ipc::CollisionMesh(
1518 is_on_surface, collision_vertices, collision_edges, collision_triangles,
1519 displacement_map);
1520
1521 collision_mesh.can_collide = [&collision_mesh, num_fe_collision_vertices](size_t vi, size_t vj) {
1522 // obstacles do not collide with other obstacles
1523 return collision_mesh.to_full_vertex_id(vi) < num_fe_collision_vertices
1524 || collision_mesh.to_full_vertex_id(vj) < num_fe_collision_vertices;
1525 };
1526
1527 collision_mesh.init_area_jacobians();
1528 }
1529
1531 {
1532 if (!mesh)
1533 {
1534 logger().error("Load the mesh first!");
1535 return;
1536 }
1537 if (n_bases <= 0)
1538 {
1539 logger().error("Build the bases first!");
1540 return;
1541 }
1542 if (assembler->name() == "OperatorSplitting")
1543 {
1545 avg_mass = 1;
1546 return;
1547 }
1548
1549 if (!problem->is_time_dependent())
1550 {
1551 avg_mass = 1;
1553 return;
1554 }
1555
1556 mass.resize(0, 0);
1557
1558 igl::Timer timer;
1559 timer.start();
1560 logger().info("Assembling mass mat...");
1561
1562 if (mixed_assembler != nullptr)
1563 {
1564 StiffnessMatrix velocity_mass;
1565 mass_matrix_assembler->assemble(mesh->is_volume(), n_bases, bases, geom_bases(), mass_ass_vals_cache, 0, velocity_mass, true);
1566
1567 std::vector<Eigen::Triplet<double>> mass_blocks;
1568 mass_blocks.reserve(velocity_mass.nonZeros());
1569
1570 for (int k = 0; k < velocity_mass.outerSize(); ++k)
1571 {
1572 for (StiffnessMatrix::InnerIterator it(velocity_mass, k); it; ++it)
1573 {
1574 mass_blocks.emplace_back(it.row(), it.col(), it.value());
1575 }
1576 }
1577
1578 mass.resize(n_bases * assembler->size(), n_bases * assembler->size());
1579 mass.setFromTriplets(mass_blocks.begin(), mass_blocks.end());
1580 mass.makeCompressed();
1581 }
1582 else
1583 {
1584 mass_matrix_assembler->assemble(mesh->is_volume(), n_bases, bases, geom_bases(), mass_ass_vals_cache, 0, mass, true);
1585 }
1586
1587 assert(mass.size() > 0);
1588
1589 avg_mass = 0;
1590 for (int k = 0; k < mass.outerSize(); ++k)
1591 {
1592
1593 for (StiffnessMatrix::InnerIterator it(mass, k); it; ++it)
1594 {
1595 assert(it.col() == k);
1596 avg_mass += it.value();
1597 }
1598 }
1599
1600 avg_mass /= mass.rows();
1601 logger().info("average mass {}", avg_mass);
1602
1603 if (args["solver"]["advanced"]["lump_mass_matrix"])
1604 {
1606 }
1607
1608 timer.stop();
1609 timings.assembling_mass_mat_time = timer.getElapsedTime();
1610 logger().info(" took {}s", timings.assembling_mass_mat_time);
1611
1612 stats.nn_zero = mass.nonZeros();
1613 stats.num_dofs = mass.rows();
1614 stats.mat_size = (long long)mass.rows() * (long long)mass.cols();
1615 logger().info("sparsity: {}/{}", stats.nn_zero, stats.mat_size);
1616 }
1617
1618 std::shared_ptr<RhsAssembler> State::build_rhs_assembler(
1619 const int n_bases_,
1620 const std::vector<basis::ElementBases> &bases_,
1621 const assembler::AssemblyValsCache &ass_vals_cache_) const
1622 {
1623 json rhs_solver_params = args["solver"]["linear"];
1624 if (!rhs_solver_params.contains("Pardiso"))
1625 rhs_solver_params["Pardiso"] = {};
1626 rhs_solver_params["Pardiso"]["mtype"] = -2; // matrix type for Pardiso (2 = SPD)
1627
1628 const int size = problem->is_scalar() ? 1 : mesh->dimension();
1629
1630 return std::make_shared<RhsAssembler>(
1634 n_bases_, size, bases_, geom_bases(), ass_vals_cache_, *problem,
1635 args["space"]["advanced"]["bc_method"],
1636 rhs_solver_params);
1637 }
1638
1639 std::shared_ptr<PressureAssembler> State::build_pressure_assembler(
1640 const int n_bases_,
1641 const std::vector<basis::ElementBases> &bases_) const
1642 {
1643 const int size = problem->is_scalar() ? 1 : mesh->dimension();
1644
1645 return std::make_shared<PressureAssembler>(
1651 n_bases_, size, bases_, geom_bases(), *problem);
1652 }
1653
1655 {
1656 if (!mesh)
1657 {
1658 logger().error("Load the mesh first!");
1659 return;
1660 }
1661 if (n_bases <= 0)
1662 {
1663 logger().error("Build the bases first!");
1664 return;
1665 }
1666
1667 igl::Timer timer;
1668
1669 json p_params = {};
1670 p_params["formulation"] = assembler->name();
1671 p_params["root_path"] = root_path();
1672 {
1673 RowVectorNd min, max, delta;
1674 mesh->bounding_box(min, max);
1675 delta = (max - min) / 2. + min;
1676 if (mesh->is_volume())
1677 p_params["bbox_center"] = {delta(0), delta(1), delta(2)};
1678 else
1679 p_params["bbox_center"] = {delta(0), delta(1)};
1680 }
1681 problem->set_parameters(p_params);
1682
1683 rhs.resize(0, 0);
1684
1685 timer.start();
1686 logger().info("Assigning rhs...");
1687
1689 solve_data.rhs_assembler->assemble(mass_matrix_assembler->density(), rhs);
1690 rhs *= -1;
1691
1692 // if(problem->is_mixed())
1693 if (mixed_assembler != nullptr)
1694 {
1695 const int prev_size = rhs.size();
1696 const int n_larger = n_pressure_bases + (use_avg_pressure ? (assembler->is_fluid() ? 1 : 0) : 0);
1697 rhs.conservativeResize(prev_size + n_larger, rhs.cols());
1698 if (assembler->name() == "OperatorSplitting")
1699 {
1701 return;
1702 }
1703 // Divergence free rhs
1704 if (assembler->name() != "Bilaplacian" || local_neumann_boundary.empty())
1705 {
1706 rhs.block(prev_size, 0, n_larger, rhs.cols()).setZero();
1707 }
1708 else
1709 {
1710 Eigen::MatrixXd tmp(n_pressure_bases, 1);
1711 tmp.setZero();
1712
1713 std::shared_ptr<RhsAssembler> tmp_rhs_assembler = build_rhs_assembler(
1715
1716 tmp_rhs_assembler->set_bc(std::vector<LocalBoundary>(), std::vector<int>(), n_boundary_samples(), local_neumann_boundary, tmp);
1717 rhs.block(prev_size, 0, n_larger, rhs.cols()) = tmp;
1718 }
1719 }
1720
1721 timer.stop();
1722 timings.assigning_rhs_time = timer.getElapsedTime();
1723 logger().info(" took {}s", timings.assigning_rhs_time);
1724 }
1725
1726 void State::solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
1727 {
1728 if (!mesh)
1729 {
1730 logger().error("Load the mesh first!");
1731 return;
1732 }
1733 if (n_bases <= 0)
1734 {
1735 logger().error("Build the bases first!");
1736 return;
1737 }
1738
1739 // if (rhs.size() <= 0)
1740 // {
1741 // logger().error("Assemble the rhs first!");
1742 // return;
1743 // }
1744
1745 // sol.resize(0, 0);
1746 // pressure.resize(0, 0);
1747 stats.spectrum.setZero();
1748
1749 igl::Timer timer;
1750 timer.start();
1751 logger().info("Solving {}", assembler->name());
1752
1753 init_solve(sol, pressure);
1754
1755 if (problem->is_time_dependent())
1756 {
1757 const double t0 = args["time"]["t0"];
1758 const int time_steps = args["time"]["time_steps"];
1759 const double dt = args["time"]["dt"];
1760
1761 // Pre log the output path for easier watching
1762 if (args["output"]["advanced"]["save_time_sequence"])
1763 {
1764 logger().info("Time sequence of simulation will be written to: \"{}\"",
1765 resolve_output_path(args["output"]["paraview"]["file_name"]));
1766 }
1767
1768 if (assembler->name() == "NavierStokes")
1769 solve_transient_navier_stokes(time_steps, t0, dt, sol, pressure);
1770 else if (assembler->name() == "OperatorSplitting")
1771 solve_transient_navier_stokes_split(time_steps, dt, sol, pressure);
1772 else if (is_homogenization())
1773 solve_homogenization(time_steps, t0, dt, sol);
1774 else if (is_problem_linear())
1775 solve_transient_linear(time_steps, t0, dt, sol, pressure);
1776 else if (!assembler->is_linear() && problem->is_scalar())
1777 throw std::runtime_error("Nonlinear scalar problems are not supported yet!");
1778 else
1779 solve_transient_tensor_nonlinear(time_steps, t0, dt, sol);
1780 }
1781 else
1782 {
1783 if (assembler->name() == "NavierStokes")
1784 solve_navier_stokes(sol, pressure);
1785 else if (is_homogenization())
1786 solve_homogenization(/* time steps */ 0, /* t0 */ 0, /* dt */ 0, sol);
1787 else if (is_problem_linear())
1788 {
1789 init_linear_solve(sol);
1790 solve_linear(sol, pressure);
1792 cache_transient_adjoint_quantities(0, sol, Eigen::MatrixXd::Zero(mesh->dimension(), mesh->dimension()));
1793 }
1794 else if (!assembler->is_linear() && problem->is_scalar())
1795 throw std::runtime_error("Nonlinear scalar problems are not supported yet!");
1796 else
1797 {
1801 cache_transient_adjoint_quantities(0, sol, Eigen::MatrixXd::Zero(mesh->dimension(), mesh->dimension()));
1802
1803 const std::string state_path = resolve_output_path(args["output"]["data"]["state"]);
1804 if (!state_path.empty())
1805 write_matrix(state_path, "u", sol);
1806 }
1807 }
1808
1809 timer.stop();
1810 timings.solving_time = timer.getElapsedTime();
1811 logger().info(" took {}s", timings.solving_time);
1812 }
1813
1814} // namespace polyfem
int V
ElementAssemblyValues vals
Definition Assembler.cpp:22
int x
DECLARE_DIFFSCALAR_BASE()
io::OutRuntimeData timings
runtime statistics
Definition State.hpp:603
assembler::MacroStrainValue macro_strain_constraint
Definition State.hpp:725
std::shared_ptr< utils::PeriodicBoundary > periodic_bc
periodic BC and periodic mesh utils
Definition State.hpp:389
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:387
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:607
std::vector< mesh::LocalBoundary > local_pressure_boundary
mapping from elements to nodes for pressure boundary conditions
Definition State.hpp:437
void sol_to_pressure(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
splits the solution in solution and pressure for mixed problems
Definition State.cpp:445
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:1530
Eigen::VectorXi in_node_to_node
Inpute nodes (including high-order) to polyfem nodes, only for isoparametric.
Definition State.hpp:450
void build_polygonal_basis()
builds bases for polygons, called inside build_basis
Definition State.cpp:1114
mesh::Obstacle obstacle
Obstacles used in collisions.
Definition State.hpp:468
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:541
void build_periodic_collision_mesh()
Definition State.cpp:1277
io::OutGeometryData out_geom
visualization stuff
Definition State.hpp:601
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:446
ipc::CollisionMesh collision_mesh
IPC collision mesh.
Definition State.hpp:538
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:667
Eigen::VectorXi in_primitive_to_primitive
maps in vertices/edges/faces/cells to polyfem vertices/edges/faces/cells
Definition State.hpp:452
std::map< int, basis::InterfaceData > poly_edge_to_data
nodes on the boundary of polygonal elements, used for harmonic bases
Definition State.hpp:441
std::vector< int > pressure_boundary_nodes
list of neumann boundary nodes
Definition State.hpp:429
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:731
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:606
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:586
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:466
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:443
StiffnessMatrix gbasis_nodes_to_basis_nodes
Definition State.hpp:719
bool has_periodic_bc() const
Definition State.hpp:390
void solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves the problems
Definition State.cpp:1726
void build_basis()
builds the bases step 2 of solve modifies bases, pressure_bases, geom_bases_, boundary_nodes,...
Definition State.cpp:577
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:444
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:605
std::vector< RowVectorNd > neumann_nodes_position
Definition State.hpp:447
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:1654
double min_boundary_edge_length
Definition State.hpp:608
Eigen::VectorXi periodic_collision_mesh_to_basis
index mapping from periodic 2x2 collision mesh to FE periodic mesh
Definition State.hpp:543
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:433
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:536
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:431
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:427
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:324
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 has contact
Definition State.hpp:574
void build_collision_mesh()
extracts the boundary mesh for collision, called in build_basis
Definition State.cpp:1407
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:435
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:410
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:439
std::string velocity() const
Definition Units.hpp:27
static double convert(const json &val, const std::string &unit_type)
Definition Units.cpp:31
const std::string & length() const
Definition Units.hpp:19
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 serendipity, const bool has_polys, const bool is_geom_bases, std::vector< ElementBases > &bases, std::vector< mesh::LocalBoundary > &local_boundary, std::map< int, InterfaceData > &poly_edge_to_data, std::shared_ptr< mesh::MeshNodes > &mesh_nodes)
Builds FE basis functions over the entire mesh (P1, P2 over triangles, Q1, Q2 over quads).
static int build_bases(const mesh::Mesh3D &mesh, const std::string &assembler, const int quadrature_order, const int mass_quadrature_order, const int discr_order, const bool serendipity, const bool has_polys, const bool is_geom_bases, std::vector< ElementBases > &bases, std::vector< mesh::LocalBoundary > &local_boundary, std::map< int, InterfaceData > &poly_face_to_data, std::shared_ptr< mesh::MeshNodes > &mesh_nodes)
Builds FE basis functions over the entire mesh (P1, P2 over tets, Q1, Q2 over hes).
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:2314
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:146
double assembling_stiffness_mat_time
time to assembly
Definition OutData.hpp:373
double assigning_rhs_time
time to computing the rhs
Definition OutData.hpp:377
double assembling_mass_mat_time
time to assembly mass
Definition OutData.hpp:375
double building_basis_time
time to construct the basis
Definition OutData.hpp:367
double solving_time
time to solve
Definition OutData.hpp:379
double computing_poly_basis_time
time to build the polygonal/polyhedral bases
Definition OutData.hpp:371
int n_flipped
number of flipped elements, compute only when using count_flipped_els (false by default)
Definition OutData.hpp:420
Eigen::Vector4d spectrum
spectrum of the stiffness matrix, enable only if POLYSOLVE_WITH_SPECTRA is ON (off by default)
Definition OutData.hpp:394
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:2499
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:2412
long long nn_zero
non zeros and sytem matrix size num dof is the total dof in the system
Definition OutData.hpp:412
double mesh_size
max edge lenght
Definition OutData.hpp:401
void reset()
clears all stats
Definition OutData.cpp:2490
double min_edge_length
min edge lenght
Definition OutData.hpp:403
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:158
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:42
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:463
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:9
nlohmann::json json
Definition Common.hpp:9
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:11
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:20