PolyFEM
Loading...
Searching...
No Matches
OutData.cpp
Go to the documentation of this file.
1#include "OutData.hpp"
2
3#include "Evaluator.hpp"
4
5#include <polyfem/State.hpp>
6
11
13
17
19
32
40
43
44#include <paraviewo/VTMWriter.hpp>
45#include <paraviewo/PVDWriter.hpp>
46
47#include <SimpleBVH/BVH.hpp>
48
49#include <igl/write_triangle_mesh.h>
50#include <igl/edges.h>
51#include <igl/facet_adjacency_matrix.h>
52#include <igl/connected_components.h>
53
54#include <ipc/ipc.hpp>
55
56#include <filesystem>
57
58namespace polyfem::io
59{
60 namespace
61 {
62 void compute_traction_forces(const State &state, const Eigen::MatrixXd &solution, const double t, Eigen::MatrixXd &traction_forces, bool skip_dirichlet = true)
63 {
64 int actual_dim = 1;
65 if (!state.problem->is_scalar())
66 actual_dim = state.mesh->dimension();
67 else
68 return;
69
70 const std::vector<basis::ElementBases> &bases = state.bases;
71 const std::vector<basis::ElementBases> &gbases = state.geom_bases();
72
73 Eigen::MatrixXd uv, samples, gtmp, rhs_fun, deform_mat, trafo;
74 Eigen::VectorXi global_primitive_ids;
75 Eigen::MatrixXd points, normals;
76 Eigen::VectorXd weights;
78
79 traction_forces.setZero(state.n_bases * actual_dim, 1);
80
81 for (const auto &lb : state.total_local_boundary)
82 {
83 const int e = lb.element_id();
84 bool has_samples = utils::BoundarySampler::boundary_quadrature(lb, state.n_boundary_samples(), *state.mesh, false, uv, points, normals, weights, global_primitive_ids);
85
86 if (!has_samples)
87 continue;
88
89 const basis::ElementBases &gbs = gbases[e];
90 const basis::ElementBases &bs = bases[e];
91
92 vals.compute(e, state.mesh->is_volume(), points, bs, gbs);
93
94 for (int n = 0; n < normals.rows(); ++n)
95 {
96 trafo = vals.jac_it[n].inverse();
97
98 if (solution.size() > 0)
99 {
100 assert(actual_dim == 2 || actual_dim == 3);
101 deform_mat.resize(actual_dim, actual_dim);
102 deform_mat.setZero();
103 for (const auto &b : vals.basis_values)
104 {
105 for (const auto &g : b.global)
106 {
107 for (int d = 0; d < actual_dim; ++d)
108 {
109 deform_mat.row(d) += solution(g.index * actual_dim + d) * b.grad.row(n);
110 }
111 }
112 }
113
114 trafo += deform_mat;
115 }
116
117 normals.row(n) = normals.row(n) * trafo.inverse();
118 normals.row(n).normalize();
119 }
120
121 std::vector<assembler::Assembler::NamedMatrix> tensor_flat;
122 state.assembler->compute_tensor_value(assembler::OutputData(t, e, bs, gbs, points, solution), tensor_flat);
123
124 for (long n = 0; n < vals.basis_values.size(); ++n)
125 {
127
128 const int g_index = v.global[0].index * actual_dim;
129
130 for (int q = 0; q < points.rows(); ++q)
131 {
132 // TF computed only from cauchy stress
133 assert(tensor_flat[0].first == "cauchy_stess");
134 assert(tensor_flat[0].second.row(q).size() == actual_dim * actual_dim);
135
136 Eigen::MatrixXd stress_tensor = utils::unflatten(tensor_flat[0].second.row(q), actual_dim);
137
138 traction_forces.block(g_index, 0, actual_dim, 1) += stress_tensor * normals.row(q).transpose() * v.val(q) * weights(q);
139 }
140 }
141 }
142 }
143 } // namespace
144
146 const mesh::Mesh &mesh,
147 const int n_bases,
148 const std::vector<basis::ElementBases> &bases,
149 const std::vector<mesh::LocalBoundary> &total_local_boundary,
150 Eigen::MatrixXd &node_positions,
151 Eigen::MatrixXi &boundary_edges,
152 Eigen::MatrixXi &boundary_triangles,
153 std::vector<Eigen::Triplet<double>> &displacement_map_entries)
154 {
155 using namespace polyfem::mesh;
156
157 displacement_map_entries.clear();
158
159 if (mesh.is_volume())
160 {
161 if (mesh.has_poly())
162 {
163 logger().warn("Skipping as the mesh has polygons");
164 return;
165 }
166
167 const bool is_simplicial = mesh.is_simplicial();
168
169 node_positions.resize(n_bases + (is_simplicial ? 0 : mesh.n_faces()), 3);
170 node_positions.setZero();
171 const Mesh3D &mesh3d = dynamic_cast<const Mesh3D &>(mesh);
172
173 std::vector<std::tuple<int, int, int>> tris;
174
175 std::vector<bool> visited_node(n_bases, false);
176
177 std::stringstream print_warning;
178
179 for (const LocalBoundary &lb : total_local_boundary)
180 {
181 const basis::ElementBases &b = bases[lb.element_id()];
182
183 for (int j = 0; j < lb.size(); ++j)
184 {
185 const int eid = lb.global_primitive_id(j);
186 const int lid = lb[j];
187 const Eigen::VectorXi nodes = b.local_nodes_for_primitive(eid, mesh3d);
188
189 if (mesh.is_cube(lb.element_id()))
190 {
191 assert(!is_simplicial);
192 assert(!mesh.has_poly());
193 std::vector<int> loc_nodes;
194 RowVectorNd bary = RowVectorNd::Zero(3);
195
196 for (long n = 0; n < nodes.size(); ++n)
197 {
198 auto &bs = b.bases[nodes(n)];
199 const auto &glob = bs.global();
200 if (glob.size() != 1)
201 continue;
202
203 int gindex = glob.front().index;
204 node_positions.row(gindex) = glob.front().node;
205 bary += glob.front().node;
206 loc_nodes.push_back(gindex);
207 }
208
209 if (loc_nodes.size() != 4)
210 {
211 logger().trace("skipping element {} since it is not Q1", eid);
212 continue;
213 }
214
215 bary /= 4;
216
217 const int new_node = n_bases + eid;
218 node_positions.row(new_node) = bary;
219 tris.emplace_back(loc_nodes[1], loc_nodes[0], new_node);
220 tris.emplace_back(loc_nodes[2], loc_nodes[1], new_node);
221 tris.emplace_back(loc_nodes[3], loc_nodes[2], new_node);
222 tris.emplace_back(loc_nodes[0], loc_nodes[3], new_node);
223
224 for (int q = 0; q < 4; ++q)
225 {
226 if (!visited_node[loc_nodes[q]])
227 displacement_map_entries.emplace_back(loc_nodes[q], loc_nodes[q], 1);
228
229 visited_node[loc_nodes[q]] = true;
230 displacement_map_entries.emplace_back(new_node, loc_nodes[q], 0.25);
231 }
232
233 continue;
234 }
235
236 if (!mesh.is_simplex(lb.element_id()))
237 {
238 logger().trace("skipping element {} since it is not a simplex or hex", eid);
239 continue;
240 }
241
242 assert(mesh.is_simplex(lb.element_id()));
243
244 std::vector<int> loc_nodes;
245
246 bool is_follower = false;
247 if (!mesh3d.is_conforming())
248 {
249 for (long n = 0; n < nodes.size(); ++n)
250 {
251 auto &bs = b.bases[nodes(n)];
252 const auto &glob = bs.global();
253 if (glob.size() != 1)
254 {
255 is_follower = true;
256 break;
257 }
258 }
259 }
260
261 if (is_follower)
262 continue;
263
264 for (long n = 0; n < nodes.size(); ++n)
265 {
266 const basis::Basis &bs = b.bases[nodes(n)];
267 const std::vector<basis::Local2Global> &glob = bs.global();
268 if (glob.size() != 1)
269 continue;
270
271 int gindex = glob.front().index;
272 node_positions.row(gindex) = glob.front().node;
273 loc_nodes.push_back(gindex);
274 }
275
276 if (loc_nodes.size() == 3)
277 {
278 tris.emplace_back(loc_nodes[0], loc_nodes[1], loc_nodes[2]);
279 }
280 else if (loc_nodes.size() == 6)
281 {
282 tris.emplace_back(loc_nodes[0], loc_nodes[3], loc_nodes[5]);
283 tris.emplace_back(loc_nodes[3], loc_nodes[1], loc_nodes[4]);
284 tris.emplace_back(loc_nodes[4], loc_nodes[2], loc_nodes[5]);
285 tris.emplace_back(loc_nodes[3], loc_nodes[4], loc_nodes[5]);
286 }
287 else if (loc_nodes.size() == 10)
288 {
289 tris.emplace_back(loc_nodes[0], loc_nodes[3], loc_nodes[8]);
290 tris.emplace_back(loc_nodes[3], loc_nodes[4], loc_nodes[9]);
291 tris.emplace_back(loc_nodes[4], loc_nodes[1], loc_nodes[5]);
292 tris.emplace_back(loc_nodes[5], loc_nodes[6], loc_nodes[9]);
293 tris.emplace_back(loc_nodes[6], loc_nodes[2], loc_nodes[7]);
294 tris.emplace_back(loc_nodes[7], loc_nodes[8], loc_nodes[9]);
295 tris.emplace_back(loc_nodes[8], loc_nodes[3], loc_nodes[9]);
296 tris.emplace_back(loc_nodes[9], loc_nodes[4], loc_nodes[5]);
297 tris.emplace_back(loc_nodes[6], loc_nodes[7], loc_nodes[9]);
298 }
299 else if (loc_nodes.size() == 15)
300 {
301 tris.emplace_back(loc_nodes[0], loc_nodes[3], loc_nodes[11]);
302 tris.emplace_back(loc_nodes[3], loc_nodes[4], loc_nodes[12]);
303 tris.emplace_back(loc_nodes[3], loc_nodes[12], loc_nodes[11]);
304 tris.emplace_back(loc_nodes[12], loc_nodes[10], loc_nodes[11]);
305 tris.emplace_back(loc_nodes[4], loc_nodes[5], loc_nodes[13]);
306 tris.emplace_back(loc_nodes[4], loc_nodes[13], loc_nodes[12]);
307 tris.emplace_back(loc_nodes[12], loc_nodes[13], loc_nodes[14]);
308 tris.emplace_back(loc_nodes[12], loc_nodes[14], loc_nodes[10]);
309 tris.emplace_back(loc_nodes[14], loc_nodes[9], loc_nodes[10]);
310 tris.emplace_back(loc_nodes[5], loc_nodes[1], loc_nodes[6]);
311 tris.emplace_back(loc_nodes[5], loc_nodes[6], loc_nodes[13]);
312 tris.emplace_back(loc_nodes[6], loc_nodes[7], loc_nodes[13]);
313 tris.emplace_back(loc_nodes[13], loc_nodes[7], loc_nodes[14]);
314 tris.emplace_back(loc_nodes[7], loc_nodes[8], loc_nodes[14]);
315 tris.emplace_back(loc_nodes[14], loc_nodes[8], loc_nodes[9]);
316 tris.emplace_back(loc_nodes[8], loc_nodes[2], loc_nodes[9]);
317 }
318 else
319 {
320 print_warning << loc_nodes.size() << " ";
321 // assert(false);
322 }
323
324 if (!is_simplicial)
325 {
326 for (int k = 0; k < loc_nodes.size(); ++k)
327 {
328 if (!visited_node[loc_nodes[k]])
329 displacement_map_entries.emplace_back(loc_nodes[k], loc_nodes[k], 1);
330
331 visited_node[loc_nodes[k]] = true;
332 }
333 }
334 }
335 }
336
337 if (print_warning.str().size() > 0)
338 logger().warn("Skipping faces as theys have {} nodes, boundary export supported up to p4", print_warning.str());
339
340 boundary_triangles.resize(tris.size(), 3);
341 for (int i = 0; i < tris.size(); ++i)
342 {
343 boundary_triangles.row(i) << std::get<0>(tris[i]), std::get<2>(tris[i]), std::get<1>(tris[i]);
344 }
345
346 if (boundary_triangles.rows() > 0)
347 {
348 igl::edges(boundary_triangles, boundary_edges);
349 }
350 }
351 else
352 {
353 node_positions.resize(n_bases, 2);
354 node_positions.setZero();
355 const Mesh2D &mesh2d = dynamic_cast<const Mesh2D &>(mesh);
356
357 std::vector<std::pair<int, int>> edges;
358
359 for (const LocalBoundary &lb : total_local_boundary)
360 {
361 const basis::ElementBases &b = bases[lb.element_id()];
362
363 for (int j = 0; j < lb.size(); ++j)
364 {
365 const int eid = lb.global_primitive_id(j);
366 const int lid = lb[j];
367 const Eigen::VectorXi nodes = b.local_nodes_for_primitive(eid, mesh2d);
368
369 int prev_node = -1;
370
371 for (long n = 0; n < nodes.size(); ++n)
372 {
373 const basis::Basis &bs = b.bases[nodes(n)];
374 const std::vector<basis::Local2Global> &glob = bs.global();
375 if (glob.size() != 1)
376 continue;
377
378 int gindex = glob.front().index;
379 node_positions.row(gindex) = glob.front().node.head<2>();
380
381 if (prev_node >= 0)
382 edges.emplace_back(prev_node, gindex);
383 prev_node = gindex;
384 }
385 }
386 }
387
388 boundary_triangles.resize(0, 0);
389 boundary_edges.resize(edges.size(), 2);
390 for (int i = 0; i < edges.size(); ++i)
391 {
392 boundary_edges.row(i) << edges[i].first, edges[i].second;
393 }
394 }
395 }
396
398 const mesh::Mesh &mesh,
399 const std::vector<basis::ElementBases> &bases,
400 const std::vector<basis::ElementBases> &gbases,
401 const std::vector<mesh::LocalBoundary> &total_local_boundary,
402 const Eigen::MatrixXd &solution,
403 const int problem_dim,
404 Eigen::MatrixXd &boundary_vis_vertices,
405 Eigen::MatrixXd &boundary_vis_local_vertices,
406 Eigen::MatrixXi &boundary_vis_elements,
407 Eigen::MatrixXi &boundary_vis_elements_ids,
408 Eigen::MatrixXi &boundary_vis_primitive_ids,
409 Eigen::MatrixXd &boundary_vis_normals,
410 Eigen::MatrixXd &displaced_boundary_vis_normals) const
411 {
412 using namespace polyfem::mesh;
413
414 std::vector<Eigen::MatrixXd> lv, vertices, allnormals, displaced_allnormals;
415 std::vector<int> el_ids, global_primitive_ids;
416 Eigen::MatrixXd uv, local_pts, tmp_n, normals, displaced_normals, trafo, deform_mat;
418 const auto &sampler = ref_element_sampler;
419 const int n_samples = sampler.num_samples();
420 int size = 0;
421
422 std::vector<std::pair<int, int>> edges;
423 std::vector<std::tuple<int, int, int>> tris;
424
425 for (auto it = total_local_boundary.begin(); it != total_local_boundary.end(); ++it)
426 {
427 const auto &lb = *it;
428 const auto &gbs = gbases[lb.element_id()];
429 const auto &bs = bases[lb.element_id()];
430
431 for (int k = 0; k < lb.size(); ++k)
432 {
433 switch (lb.type())
434 {
435 case BoundaryType::TRI_LINE:
437 utils::BoundarySampler::sample_parametric_tri_edge(lb[k], n_samples, uv, local_pts);
438 break;
439 case BoundaryType::QUAD_LINE:
441 utils::BoundarySampler::sample_parametric_quad_edge(lb[k], n_samples, uv, local_pts);
442 break;
443 case BoundaryType::QUAD:
445 utils::BoundarySampler::sample_parametric_quad_face(lb[k], n_samples, uv, local_pts);
446 break;
447 case BoundaryType::TRI:
449 utils::BoundarySampler::sample_parametric_tri_face(lb[k], n_samples, uv, local_pts);
450 break;
451 case BoundaryType::POLYGON:
452 utils::BoundarySampler::normal_for_polygon_edge(lb.element_id(), lb.global_primitive_id(k), mesh, tmp_n);
453 utils::BoundarySampler::sample_polygon_edge(lb.element_id(), lb.global_primitive_id(k), n_samples, mesh, uv, local_pts);
454 break;
455 case BoundaryType::POLYHEDRON:
456 assert(false);
457 break;
458 case BoundaryType::INVALID:
459 assert(false);
460 break;
461 default:
462 assert(false);
463 }
464
465 vertices.emplace_back();
466 lv.emplace_back(local_pts);
467 el_ids.push_back(lb.element_id());
468 global_primitive_ids.push_back(lb.global_primitive_id(k));
469 gbs.eval_geom_mapping(local_pts, vertices.back());
470 vals.compute(lb.element_id(), mesh.is_volume(), local_pts, bs, gbs);
471 const int tris_start = tris.size();
472
473 if (mesh.is_volume())
474 {
475 if (lb.type() == BoundaryType::QUAD)
476 {
477 const auto map = [n_samples, size](int i, int j) { return j * n_samples + i + size; };
478
479 for (int j = 0; j < n_samples - 1; ++j)
480 {
481 for (int i = 0; i < n_samples - 1; ++i)
482 {
483 tris.emplace_back(map(i, j), map(i + 1, j), map(i, j + 1));
484 tris.emplace_back(map(i + 1, j + 1), map(i, j + 1), map(i + 1, j));
485 }
486 }
487 }
488 else if (lb.type() == BoundaryType::TRI)
489 {
490 int index = 0;
491 std::vector<int> mapp(n_samples * n_samples, -1);
492 for (int j = 0; j < n_samples; ++j)
493 {
494 for (int i = 0; i < n_samples - j; ++i)
495 {
496 mapp[j * n_samples + i] = index;
497 ++index;
498 }
499 }
500 const auto map = [mapp, n_samples](int i, int j) {
501 if (j * n_samples + i >= mapp.size())
502 return -1;
503 return mapp[j * n_samples + i];
504 };
505
506 for (int j = 0; j < n_samples - 1; ++j)
507 {
508 for (int i = 0; i < n_samples - j; ++i)
509 {
510 if (map(i, j) >= 0 && map(i + 1, j) >= 0 && map(i, j + 1) >= 0)
511 tris.emplace_back(map(i, j) + size, map(i + 1, j) + size, map(i, j + 1) + size);
512
513 if (map(i + 1, j + 1) >= 0 && map(i, j + 1) >= 0 && map(i + 1, j) >= 0)
514 tris.emplace_back(map(i + 1, j + 1) + size, map(i, j + 1) + size, map(i + 1, j) + size);
515 }
516 }
517 }
518 else
519 {
520 assert(false);
521 }
522 }
523 else
524 {
525 for (int i = 0; i < vertices.back().rows() - 1; ++i)
526 edges.emplace_back(i + size, i + size + 1);
527 }
528
529 normals.resize(vals.jac_it.size(), tmp_n.cols());
530 displaced_normals.resize(vals.jac_it.size(), tmp_n.cols());
531
532 for (int n = 0; n < vals.jac_it.size(); ++n)
533 {
534 trafo = vals.jac_it[n].inverse();
535
536 if (problem_dim == 2 || problem_dim == 3)
537 {
538
539 if (solution.size() > 0)
540 {
541 deform_mat.resize(problem_dim, problem_dim);
542 deform_mat.setZero();
543 for (const auto &b : vals.basis_values)
544 for (const auto &g : b.global)
545 for (int d = 0; d < problem_dim; ++d)
546 deform_mat.row(d) += solution(g.index * problem_dim + d) * b.grad.row(n);
547
548 trafo += deform_mat;
549 }
550 }
551
552 normals.row(n) = tmp_n * vals.jac_it[n];
553 normals.row(n).normalize();
554
555 displaced_normals.row(n) = tmp_n * trafo.inverse();
556 displaced_normals.row(n).normalize();
557 }
558
559 allnormals.push_back(normals);
560 displaced_allnormals.push_back(displaced_normals);
561
562 tmp_n.setZero();
563 for (int n = 0; n < vals.jac_it.size(); ++n)
564 {
565 tmp_n += normals.row(n);
566 }
567
568 if (mesh.is_volume())
569 {
570 Eigen::Vector3d e1 = vertices.back().row(std::get<1>(tris.back()) - size) - vertices.back().row(std::get<0>(tris.back()) - size);
571 Eigen::Vector3d e2 = vertices.back().row(std::get<2>(tris.back()) - size) - vertices.back().row(std::get<0>(tris.back()) - size);
572
573 Eigen::Vector3d n = e1.cross(e2);
574 Eigen::Vector3d nn = tmp_n.transpose();
575
576 if (n.dot(nn) < 0)
577 {
578 for (int i = tris_start; i < tris.size(); ++i)
579 {
580 tris[i] = std::tuple<int, int, int>(std::get<0>(tris[i]), std::get<2>(tris[i]), std::get<1>(tris[i]));
581 }
582 }
583 }
584
585 size += vertices.back().rows();
586 }
587 }
588
589 boundary_vis_vertices.resize(size, vertices.front().cols());
590 boundary_vis_local_vertices.resize(size, vertices.front().cols());
591 boundary_vis_elements_ids.resize(size, 1);
592 boundary_vis_primitive_ids.resize(size, 1);
593 boundary_vis_normals.resize(size, vertices.front().cols());
594 displaced_boundary_vis_normals.resize(size, vertices.front().cols());
595
596 if (mesh.is_volume())
597 boundary_vis_elements.resize(tris.size(), 3);
598 else
599 boundary_vis_elements.resize(edges.size(), 2);
600
601 int index = 0;
602 int ii = 0;
603 for (const auto &v : vertices)
604 {
605 boundary_vis_vertices.block(index, 0, v.rows(), v.cols()) = v;
606 boundary_vis_local_vertices.block(index, 0, v.rows(), v.cols()) = lv[ii];
607 boundary_vis_elements_ids.block(index, 0, v.rows(), 1).setConstant(el_ids[ii]);
608 boundary_vis_primitive_ids.block(index, 0, v.rows(), 1).setConstant(global_primitive_ids[ii++]);
609 index += v.rows();
610 }
611
612 index = 0;
613 for (const auto &n : allnormals)
614 {
615 boundary_vis_normals.block(index, 0, n.rows(), n.cols()) = n;
616 index += n.rows();
617 }
618
619 index = 0;
620 for (const auto &n : displaced_allnormals)
621 {
622 displaced_boundary_vis_normals.block(index, 0, n.rows(), n.cols()) = n;
623 index += n.rows();
624 }
625
626 index = 0;
627 if (mesh.is_volume())
628 {
629 for (const auto &t : tris)
630 {
631 boundary_vis_elements.row(index) << std::get<0>(t), std::get<1>(t), std::get<2>(t);
632 ++index;
633 }
634 }
635 else
636 {
637 for (const auto &e : edges)
638 {
639 boundary_vis_elements.row(index) << e.first, e.second;
640 ++index;
641 }
642 }
643 }
644
646 const mesh::Mesh &mesh,
647 const Eigen::VectorXi &disc_orders,
648 const std::vector<basis::ElementBases> &gbases,
649 const std::map<int, Eigen::MatrixXd> &polys,
650 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
651 const bool boundary_only,
652 Eigen::MatrixXd &points,
653 Eigen::MatrixXi &tets,
654 Eigen::MatrixXi &el_id,
655 Eigen::MatrixXd &discr) const
656 {
657 // if (!mesh)
658 // {
659 // logger().error("Load the mesh first!");
660 // return;
661 // }
662 // if (n_bases <= 0)
663 // {
664 // logger().error("Build the bases first!");
665 // return;
666 // }
667
668 const auto &sampler = ref_element_sampler;
669
670 const auto &current_bases = gbases;
671 int tet_total_size = 0;
672 int pts_total_size = 0;
673
674 Eigen::MatrixXd vis_pts_poly;
675 Eigen::MatrixXi vis_faces_poly, vis_edges_poly;
676
677 for (size_t i = 0; i < current_bases.size(); ++i)
678 {
679 const auto &bs = current_bases[i];
680
681 if (boundary_only && mesh.is_volume() && !mesh.is_boundary_element(i))
682 continue;
683
684 if (mesh.is_simplex(i))
685 {
686 tet_total_size += sampler.simplex_volume().rows();
687 pts_total_size += sampler.simplex_points().rows();
688 }
689 else if (mesh.is_cube(i))
690 {
691 tet_total_size += sampler.cube_volume().rows();
692 pts_total_size += sampler.cube_points().rows();
693 }
694 else
695 {
696 if (mesh.is_volume())
697 {
698 sampler.sample_polyhedron(polys_3d.at(i).first, polys_3d.at(i).second, vis_pts_poly, vis_faces_poly, vis_edges_poly);
699
700 tet_total_size += vis_faces_poly.rows();
701 pts_total_size += vis_pts_poly.rows();
702 }
703 else
704 {
705 sampler.sample_polygon(polys.at(i), vis_pts_poly, vis_faces_poly, vis_edges_poly);
706
707 tet_total_size += vis_faces_poly.rows();
708 pts_total_size += vis_pts_poly.rows();
709 }
710 }
711 }
712
713 points.resize(pts_total_size, mesh.dimension());
714 tets.resize(tet_total_size, mesh.is_volume() ? 4 : 3);
715
716 el_id.resize(pts_total_size, 1);
717 discr.resize(pts_total_size, 1);
718
719 Eigen::MatrixXd mapped, tmp;
720 int tet_index = 0, pts_index = 0;
721
722 for (size_t i = 0; i < current_bases.size(); ++i)
723 {
724 const auto &bs = current_bases[i];
725
726 if (boundary_only && mesh.is_volume() && !mesh.is_boundary_element(i))
727 continue;
728
729 if (mesh.is_simplex(i))
730 {
731 bs.eval_geom_mapping(sampler.simplex_points(), mapped);
732
733 tets.block(tet_index, 0, sampler.simplex_volume().rows(), tets.cols()) = sampler.simplex_volume().array() + pts_index;
734 tet_index += sampler.simplex_volume().rows();
735
736 points.block(pts_index, 0, mapped.rows(), points.cols()) = mapped;
737 discr.block(pts_index, 0, mapped.rows(), 1).setConstant(disc_orders(i));
738 el_id.block(pts_index, 0, mapped.rows(), 1).setConstant(i);
739 pts_index += mapped.rows();
740 }
741 else if (mesh.is_cube(i))
742 {
743 bs.eval_geom_mapping(sampler.cube_points(), mapped);
744
745 tets.block(tet_index, 0, sampler.cube_volume().rows(), tets.cols()) = sampler.cube_volume().array() + pts_index;
746 tet_index += sampler.cube_volume().rows();
747
748 points.block(pts_index, 0, mapped.rows(), points.cols()) = mapped;
749 discr.block(pts_index, 0, mapped.rows(), 1).setConstant(disc_orders(i));
750 el_id.block(pts_index, 0, mapped.rows(), 1).setConstant(i);
751 pts_index += mapped.rows();
752 }
753 else
754 {
755 if (mesh.is_volume())
756 {
757 sampler.sample_polyhedron(polys_3d.at(i).first, polys_3d.at(i).second, vis_pts_poly, vis_faces_poly, vis_edges_poly);
758 bs.eval_geom_mapping(vis_pts_poly, mapped);
759
760 tets.block(tet_index, 0, vis_faces_poly.rows(), tets.cols()) = vis_faces_poly.array() + pts_index;
761 tet_index += vis_faces_poly.rows();
762
763 points.block(pts_index, 0, mapped.rows(), points.cols()) = mapped;
764 discr.block(pts_index, 0, mapped.rows(), 1).setConstant(-1);
765 el_id.block(pts_index, 0, mapped.rows(), 1).setConstant(i);
766 pts_index += mapped.rows();
767 }
768 else
769 {
770 sampler.sample_polygon(polys.at(i), vis_pts_poly, vis_faces_poly, vis_edges_poly);
771 bs.eval_geom_mapping(vis_pts_poly, mapped);
772
773 tets.block(tet_index, 0, vis_faces_poly.rows(), tets.cols()) = vis_faces_poly.array() + pts_index;
774 tet_index += vis_faces_poly.rows();
775
776 points.block(pts_index, 0, mapped.rows(), points.cols()) = mapped;
777 discr.block(pts_index, 0, mapped.rows(), 1).setConstant(-1);
778 el_id.block(pts_index, 0, mapped.rows(), 1).setConstant(i);
779 pts_index += mapped.rows();
780 }
781 }
782 }
783
784 assert(pts_index == points.rows());
785 assert(tet_index == tets.rows());
786 }
787
789 const mesh::Mesh &mesh,
790 const Eigen::VectorXi &disc_orders,
791 const std::vector<basis::ElementBases> &bases,
792 Eigen::MatrixXd &points,
793 std::vector<std::vector<int>> &elements,
794 Eigen::MatrixXi &el_id,
795 Eigen::MatrixXd &discr) const
796 {
797 // if (!mesh)
798 // {
799 // logger().error("Load the mesh first!");
800 // return;
801 // }
802 // if (n_bases <= 0)
803 // {
804 // logger().error("Build the bases first!");
805 // return;
806 // }
807 assert(mesh.is_linear());
808
809 std::vector<RowVectorNd> nodes;
810 int pts_total_size = 0;
811 elements.resize(bases.size());
812 Eigen::MatrixXd ref_pts;
813
814 for (size_t i = 0; i < bases.size(); ++i)
815 {
816 const auto &bs = bases[i];
817 if (mesh.is_volume())
818 {
819 if (mesh.is_simplex(i))
820 autogen::p_nodes_3d(disc_orders(i), ref_pts);
821 else if (mesh.is_cube(i))
822 autogen::q_nodes_3d(disc_orders(i), ref_pts);
823 else
824 continue;
825 }
826 else
827 {
828 if (mesh.is_simplex(i))
829 autogen::p_nodes_2d(disc_orders(i), ref_pts);
830 else if (mesh.is_cube(i))
831 autogen::q_nodes_2d(disc_orders(i), ref_pts);
832 else
833 {
834 const int n_v = static_cast<const mesh::Mesh2D &>(mesh).n_face_vertices(i);
835 ref_pts.resize(n_v, 2);
836 }
837 }
838
839 pts_total_size += ref_pts.rows();
840 }
841
842 points.resize(pts_total_size, mesh.dimension());
843
844 el_id.resize(pts_total_size, 1);
845 discr.resize(pts_total_size, 1);
846
847 Eigen::MatrixXd mapped;
848 int pts_index = 0;
849
850 std::string error_msg = "";
851
852 for (size_t i = 0; i < bases.size(); ++i)
853 {
854 const auto &bs = bases[i];
855 if (mesh.is_volume())
856 {
857 if (mesh.is_simplex(i))
858 autogen::p_nodes_3d(disc_orders(i), ref_pts);
859 else if (mesh.is_cube(i))
860 autogen::q_nodes_3d(disc_orders(i), ref_pts);
861 else
862 continue;
863 }
864 else
865 {
866 if (mesh.is_simplex(i))
867 autogen::p_nodes_2d(disc_orders(i), ref_pts);
868 else if (mesh.is_cube(i))
869 autogen::q_nodes_2d(disc_orders(i), ref_pts);
870 else
871 continue;
872 }
873
874 bs.eval_geom_mapping(ref_pts, mapped);
875
876 for (int j = 0; j < mapped.rows(); ++j)
877 {
878 points.row(pts_index) = mapped.row(j);
879 el_id(pts_index) = i;
880 discr(pts_index) = disc_orders(i);
881 elements[i].push_back(pts_index);
882
883 pts_index++;
884 }
885
886 if (mesh.is_simplex(i))
887 {
888 if (mesh.is_volume())
889 {
890 const int n_nodes = elements[i].size();
891 if (disc_orders(i) >= 3)
892 {
893 std::swap(elements[i][16], elements[i][17]);
894 std::swap(elements[i][17], elements[i][18]);
895 std::swap(elements[i][18], elements[i][19]);
896 }
897 if (disc_orders(i) > 4)
898 error_msg = "Saving high-order meshes not implemented for P5+ elements!";
899 }
900 else
901 {
902 if (disc_orders(i) == 4)
903 {
904 const int n_nodes = elements[i].size();
905 std::swap(elements[i][n_nodes - 1], elements[i][n_nodes - 2]);
906 }
907 if (disc_orders(i) > 4)
908 error_msg = "Saving high-order meshes not implemented for P5+ elements!";
909 }
910 }
911 else if (disc_orders(i) > 1)
912 error_msg = "Saving high-order meshes not implemented for Q2+ elements!";
913 }
914
915 if (!error_msg.empty())
916 logger().warn(error_msg);
917
918 for (size_t i = 0; i < bases.size(); ++i)
919 {
920 if (mesh.is_volume() || !mesh.is_polytope(i))
921 continue;
922
923 const auto &mesh2d = static_cast<const mesh::Mesh2D &>(mesh);
924 const int n_v = mesh2d.n_face_vertices(i);
925
926 for (int j = 0; j < n_v; ++j)
927 {
928 points.row(pts_index) = mesh2d.point(mesh2d.face_vertex(i, j));
929 el_id(pts_index) = i;
930 discr(pts_index) = disc_orders(i);
931 elements[i].push_back(pts_index);
932
933 pts_index++;
934 }
935 }
936
937 assert(pts_index == points.rows());
938 }
939
941 const State &state,
942 const Eigen::MatrixXd &sol,
943 const Eigen::MatrixXd &pressure,
944 const bool is_time_dependent,
945 const double tend_in,
946 const double dt,
947 const ExportOptions &opts,
948 const std::string &vis_mesh_path,
949 const std::string &nodes_path,
950 const std::string &solution_path,
951 const std::string &stress_path,
952 const std::string &mises_path,
953 const bool is_contact_enabled,
954 std::vector<SolutionFrame> &solution_frames) const
955 {
956 if (!state.mesh)
957 {
958 logger().error("Load the mesh first!");
959 return;
960 }
961 const int n_bases = state.n_bases;
962 const std::vector<basis::ElementBases> &bases = state.bases;
963 const std::vector<basis::ElementBases> &gbases = state.geom_bases();
964 const mesh::Mesh &mesh = *state.mesh;
965 const Eigen::VectorXi &in_node_to_node = state.in_node_to_node;
966 const Eigen::MatrixXd &rhs = state.rhs;
967 const assembler::Problem &problem = *state.problem;
968
969 if (n_bases <= 0)
970 {
971 logger().error("Build the bases first!");
972 return;
973 }
974 // if (rhs.size() <= 0)
975 // {
976 // logger().error("Assemble the rhs first!");
977 // return;
978 // }
979 if (sol.size() <= 0)
980 {
981 logger().error("Solve the problem first!");
982 return;
983 }
984
985 if (!solution_path.empty())
986 {
987 std::ofstream out(solution_path);
988 out.precision(100);
989 out << std::scientific;
990 if (opts.reorder_output)
991 {
992 int problem_dim = (problem.is_scalar() ? 1 : mesh.dimension());
993 Eigen::VectorXi reordering(n_bases);
994 reordering.setConstant(-1);
995
996 for (int i = 0; i < in_node_to_node.size(); ++i)
997 {
998 reordering[in_node_to_node[i]] = i;
999 }
1000 Eigen::MatrixXd tmp_sol = utils::unflatten(sol, problem_dim);
1001 Eigen::MatrixXd tmp(tmp_sol.rows(), tmp_sol.cols());
1002
1003 for (int i = 0; i < reordering.size(); ++i)
1004 {
1005 if (reordering[i] < 0)
1006 continue;
1007
1008 tmp.row(reordering[i]) = tmp_sol.row(i);
1009 }
1010
1011 for (int i = 0; i < tmp.rows(); ++i)
1012 {
1013 for (int j = 0; j < tmp.cols(); ++j)
1014 out << tmp(i, j) << " ";
1015
1016 out << std::endl;
1017 }
1018 }
1019 else
1020 out << sol << std::endl;
1021 out.close();
1022 }
1023
1024 double tend = tend_in;
1025 if (tend <= 0)
1026 tend = 1;
1027
1028 if (!vis_mesh_path.empty() && !is_time_dependent)
1029 {
1030 save_vtu(
1031 vis_mesh_path, state, sol, pressure,
1032 tend, dt, opts,
1033 is_contact_enabled, solution_frames);
1034 }
1035 if (!nodes_path.empty())
1036 {
1037 Eigen::MatrixXd nodes(n_bases, mesh.dimension());
1038 for (const basis::ElementBases &eb : bases)
1039 {
1040 for (const basis::Basis &b : eb.bases)
1041 {
1042 // for(const auto &lg : b.global())
1043 for (size_t ii = 0; ii < b.global().size(); ++ii)
1044 {
1045 const auto &lg = b.global()[ii];
1046 nodes.row(lg.index) = lg.node;
1047 }
1048 }
1049 }
1050 std::ofstream out(nodes_path);
1051 out.precision(100);
1052 out << nodes;
1053 out.close();
1054 }
1055 if (!stress_path.empty())
1056 {
1057 Eigen::MatrixXd result;
1058 Eigen::VectorXd mises;
1060 mesh, problem.is_scalar(),
1061 bases, gbases, state.disc_orders, *state.assembler,
1062 sol, tend, result, mises);
1063 std::ofstream out(stress_path);
1064 out.precision(20);
1065 out << result;
1066 }
1067 if (!mises_path.empty())
1068 {
1069 Eigen::MatrixXd result;
1070 Eigen::VectorXd mises;
1072 mesh, problem.is_scalar(),
1073 bases, gbases, state.disc_orders, *state.assembler,
1074 sol, tend, result, mises);
1075 std::ofstream out(mises_path);
1076 out.precision(20);
1077 out << mises;
1078 }
1079 }
1080
1081 OutGeometryData::ExportOptions::ExportOptions(const json &args, const bool is_mesh_linear, const bool is_problem_scalar, const bool solve_export_to_file)
1082 {
1083 volume = args["output"]["paraview"]["volume"];
1084 surface = args["output"]["paraview"]["surface"];
1085 wire = args["output"]["paraview"]["wireframe"];
1086 points = args["output"]["paraview"]["points"];
1087 contact_forces = args["output"]["paraview"]["options"]["contact_forces"] && !is_problem_scalar;
1088 friction_forces = args["output"]["paraview"]["options"]["friction_forces"] && !is_problem_scalar;
1089
1090 use_sampler = !(is_mesh_linear && solve_export_to_file && args["output"]["paraview"]["high_order_mesh"]);
1091 boundary_only = use_sampler && args["output"]["advanced"]["vis_boundary_only"];
1092 material_params = args["output"]["paraview"]["options"]["material"];
1093 body_ids = args["output"]["paraview"]["options"]["body_ids"];
1094 sol_on_grid = args["output"]["advanced"]["sol_on_grid"] > 0;
1095 velocity = args["output"]["paraview"]["options"]["velocity"];
1096 acceleration = args["output"]["paraview"]["options"]["acceleration"];
1097 forces = args["output"]["paraview"]["options"]["forces"] && !is_problem_scalar;
1098 jacobian_validity = args["output"]["paraview"]["options"]["jacobian_validity"] && !is_problem_scalar;
1099
1100 scalar_values = args["output"]["paraview"]["options"]["scalar_values"];
1101 tensor_values = args["output"]["paraview"]["options"]["tensor_values"] && !is_problem_scalar;
1102 discretization_order = args["output"]["paraview"]["options"]["discretization_order"];
1103 nodes = args["output"]["paraview"]["options"]["nodes"] && !is_problem_scalar;
1104
1105 use_spline = args["space"]["basis_type"] == "Spline";
1106
1107 reorder_output = args["output"]["data"]["advanced"]["reorder_nodes"];
1108
1109 use_hdf5 = args["output"]["paraview"]["options"]["use_hdf5"];
1110
1111 this->solve_export_to_file = solve_export_to_file;
1112 }
1113
1115 const std::string &path,
1116 const State &state,
1117 const Eigen::MatrixXd &sol,
1118 const Eigen::MatrixXd &pressure,
1119 const double t,
1120 const double dt,
1121 const ExportOptions &opts,
1122 const bool is_contact_enabled,
1123 std::vector<SolutionFrame> &solution_frames) const
1124 {
1125 if (!state.mesh)
1126 {
1127 logger().error("Load the mesh first!");
1128 return;
1129 }
1130 const mesh::Mesh &mesh = *state.mesh;
1131 const Eigen::MatrixXd &rhs = state.rhs;
1132
1133 if (state.n_bases <= 0)
1134 {
1135 logger().error("Build the bases first!");
1136 return;
1137 }
1138 // if (rhs.size() <= 0)
1139 // {
1140 // logger().error("Assemble the rhs first!");
1141 // return;
1142 // }
1143 if (sol.size() <= 0)
1144 {
1145 logger().error("Solve the problem first!");
1146 return;
1147 }
1148
1149 const std::filesystem::path fs_path(path);
1150 const std::string path_stem = fs_path.stem().string();
1151 const std::string base_path = (fs_path.parent_path() / path_stem).string();
1152
1153 if (opts.volume)
1154 {
1155 save_volume(base_path + opts.file_extension(), state, sol, pressure, t, dt, opts, solution_frames);
1156 }
1157
1158 if (opts.surface)
1159 {
1160 save_surface(base_path + "_surf" + opts.file_extension(), state, sol, pressure, t, dt, opts,
1161 is_contact_enabled, solution_frames);
1162 }
1163
1164 if (is_contact_enabled && (opts.contact_forces || opts.friction_forces))
1165 {
1166 save_contact_surface(base_path + "_surf" + opts.file_extension(), state, sol, pressure, t, dt, opts,
1167 is_contact_enabled, solution_frames);
1168 }
1169
1170 if (opts.wire)
1171 {
1172 save_wire(base_path + "_wire" + opts.file_extension(), state, sol, t, opts, solution_frames);
1173 }
1174
1175 if (opts.points)
1176 {
1177 save_points(base_path + "_points" + opts.file_extension(), state, sol, opts, solution_frames);
1178 }
1179
1180 if (!opts.solve_export_to_file)
1181 return;
1182
1183 paraviewo::VTMWriter vtm(t);
1184 if (opts.volume)
1185 vtm.add_dataset("Volume", "data", path_stem + opts.file_extension());
1186 if (opts.surface)
1187 vtm.add_dataset("Surface", "data", path_stem + "_surf" + opts.file_extension());
1188 if (is_contact_enabled && (opts.contact_forces || opts.friction_forces))
1189 vtm.add_dataset("Contact", "data", path_stem + "_surf_contact" + opts.file_extension());
1190 if (opts.wire)
1191 vtm.add_dataset("Wireframe", "data", path_stem + "_wire" + opts.file_extension());
1192 if (opts.points)
1193 vtm.add_dataset("Points", "data", path_stem + "_points" + opts.file_extension());
1194 vtm.save(base_path + ".vtm");
1195 }
1196
1198 const std::string &path,
1199 const State &state,
1200 const Eigen::MatrixXd &sol,
1201 const Eigen::MatrixXd &pressure,
1202 const double t,
1203 const double dt,
1204 const ExportOptions &opts,
1205 std::vector<SolutionFrame> &solution_frames) const
1206 {
1207 const Eigen::VectorXi &disc_orders = state.disc_orders;
1208 const auto &density = state.mass_matrix_assembler->density();
1209 const std::vector<basis::ElementBases> &bases = state.bases;
1210 const std::vector<basis::ElementBases> &pressure_bases = state.pressure_bases;
1211 const std::vector<basis::ElementBases> &gbases = state.geom_bases();
1212 const std::map<int, Eigen::MatrixXd> &polys = state.polys;
1213 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d = state.polys_3d;
1214 const assembler::Assembler &assembler = *state.assembler;
1215 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> time_integrator = state.solve_data.time_integrator;
1216 const mesh::Mesh &mesh = *state.mesh;
1217 const mesh::Obstacle &obstacle = state.obstacle;
1218 const assembler::Problem &problem = *state.problem;
1219
1220 Eigen::MatrixXd points;
1221 Eigen::MatrixXi tets;
1222 Eigen::MatrixXi el_id;
1223 Eigen::MatrixXd discr;
1224 std::vector<std::vector<int>> elements;
1225
1226 if (opts.use_sampler)
1227 build_vis_mesh(mesh, disc_orders, gbases,
1228 state.polys, state.polys_3d, opts.boundary_only,
1229 points, tets, el_id, discr);
1230 else
1231 build_high_order_vis_mesh(mesh, disc_orders, bases,
1232 points, elements, el_id, discr);
1233
1234 Eigen::MatrixXd fun, exact_fun, err, node_fun;
1235
1236 if (opts.sol_on_grid)
1237 {
1238 const int problem_dim = problem.is_scalar() ? 1 : mesh.dimension();
1239 Eigen::MatrixXd tmp, tmp_grad;
1240 Eigen::MatrixXd tmp_p, tmp_grad_p;
1241 Eigen::MatrixXd res(grid_points_to_elements.size(), problem_dim);
1242 res.setConstant(std::numeric_limits<double>::quiet_NaN());
1243 Eigen::MatrixXd res_grad(grid_points_to_elements.size(), problem_dim * problem_dim);
1244 res_grad.setConstant(std::numeric_limits<double>::quiet_NaN());
1245
1246 Eigen::MatrixXd res_p(grid_points_to_elements.size(), 1);
1247 res_p.setConstant(std::numeric_limits<double>::quiet_NaN());
1248 Eigen::MatrixXd res_grad_p(grid_points_to_elements.size(), problem_dim);
1249 res_grad_p.setConstant(std::numeric_limits<double>::quiet_NaN());
1250
1251 for (int i = 0; i < grid_points_to_elements.size(); ++i)
1252 {
1253 const int el_id = grid_points_to_elements(i);
1254 if (el_id < 0)
1255 continue;
1256 assert(mesh.is_simplex(el_id));
1257 const Eigen::MatrixXd bc = grid_points_bc.row(i);
1258 Eigen::MatrixXd pt(1, bc.cols() - 1);
1259 for (int d = 1; d < bc.cols(); ++d)
1260 pt(d - 1) = bc(d);
1262 mesh, problem.is_scalar(), bases, gbases,
1263 el_id, pt, sol, tmp, tmp_grad);
1264
1265 res.row(i) = tmp;
1266 res_grad.row(i) = tmp_grad;
1267
1268 if (state.mixed_assembler != nullptr)
1269 {
1271 mesh, 1, pressure_bases, gbases,
1272 el_id, pt, pressure, tmp_p, tmp_grad_p);
1273 res_p.row(i) = tmp_p;
1274 res_grad_p.row(i) = tmp_grad_p;
1275 }
1276 }
1277
1278 std::ofstream os(path + "_sol.txt");
1279 os << res;
1280
1281 std::ofstream osg(path + "_grad.txt");
1282 osg << res_grad;
1283
1284 std::ofstream osgg(path + "_grid.txt");
1285 osgg << grid_points;
1286
1287 if (state.mixed_assembler != nullptr)
1288 {
1289 std::ofstream osp(path + "_p_sol.txt");
1290 osp << res_p;
1291
1292 std::ofstream osgp(path + "_p_grad.txt");
1293 osgp << res_grad_p;
1294 }
1295 }
1296
1297 Eigen::Vector<bool, -1> validity;
1298 if (opts.jacobian_validity)
1300 mesh, gbases, bases, state.disc_orders,
1301 state.polys, state.polys_3d, ref_element_sampler,
1302 points.rows(), sol, validity, opts.use_sampler, opts.boundary_only);
1303
1305 mesh, problem.is_scalar(), bases, state.disc_orders,
1306 state.polys, state.polys_3d, ref_element_sampler,
1307 points.rows(), sol, fun, opts.use_sampler, opts.boundary_only);
1308
1309 {
1310 Eigen::MatrixXd tmp = Eigen::VectorXd::LinSpaced(sol.size(), 0, sol.size() - 1);
1311
1313 mesh, problem.is_scalar(), bases, state.disc_orders,
1314 state.polys, state.polys_3d, ref_element_sampler,
1315 points.rows(), tmp, node_fun, opts.use_sampler, opts.boundary_only);
1316 }
1317
1318 if (obstacle.n_vertices() > 0)
1319 {
1320 fun.conservativeResize(fun.rows() + obstacle.n_vertices(), fun.cols());
1321 node_fun.conservativeResize(node_fun.rows() + obstacle.n_vertices(), node_fun.cols());
1322 node_fun.bottomRows(obstacle.n_vertices()).setZero();
1323 // obstacle.update_displacement(t, fun);
1324 // NOTE: Assuming the obstacle displacement is the last part of the solution
1325 fun.bottomRows(obstacle.n_vertices()) = utils::unflatten(sol.bottomRows(obstacle.ndof()), fun.cols());
1326 }
1327
1328 if (problem.has_exact_sol())
1329 {
1330 problem.exact(points, t, exact_fun);
1331 err = (fun - exact_fun).eval().rowwise().norm();
1332
1333 if (obstacle.n_vertices() > 0)
1334 {
1335 exact_fun.conservativeResize(exact_fun.rows() + obstacle.n_vertices(), exact_fun.cols());
1336 // obstacle.update_displacement(t, exact_fun);
1337 exact_fun.bottomRows(obstacle.n_vertices()) = utils::unflatten(sol.bottomRows(obstacle.ndof()), fun.cols());
1338
1339 err.conservativeResize(err.rows() + obstacle.n_vertices(), 1);
1340 err.bottomRows(obstacle.n_vertices()).setZero();
1341 }
1342 }
1343
1344 std::shared_ptr<paraviewo::ParaviewWriter> tmpw;
1345 if (opts.use_hdf5)
1346 tmpw = std::make_shared<paraviewo::HDF5VTUWriter>();
1347 else
1348 tmpw = std::make_shared<paraviewo::VTUWriter>();
1349 paraviewo::ParaviewWriter &writer = *tmpw;
1350
1351 if (validity.size())
1352 writer.add_field("validity", validity.cast<double>());
1353
1354 if (opts.solve_export_to_file && opts.nodes)
1355 writer.add_field("nodes", node_fun);
1356
1357 if (problem.is_time_dependent())
1358 {
1359 bool is_time_integrator_valid = time_integrator != nullptr;
1360
1361 if (opts.velocity)
1362 {
1363 const Eigen::VectorXd velocity =
1364 is_time_integrator_valid ? (time_integrator->v_prev()) : Eigen::VectorXd::Zero(sol.size());
1365 save_volume_vector_field(state, points, opts, "velocity", velocity, writer);
1366 }
1367
1368 if (opts.acceleration)
1369 {
1370 const Eigen::VectorXd acceleration =
1371 is_time_integrator_valid ? (time_integrator->a_prev()) : Eigen::VectorXd::Zero(sol.size());
1372 save_volume_vector_field(state, points, opts, "acceleration", acceleration, writer);
1373 }
1374 }
1375
1376 if (opts.forces)
1377 {
1378 const double s = state.solve_data.time_integrator
1379 ? state.solve_data.time_integrator->acceleration_scaling()
1380 : 1;
1381
1382 for (const auto &[name, form] : state.solve_data.named_forms())
1383 {
1384 // NOTE: Assumes this form will be null for the entire sim
1385 if (form == nullptr)
1386 continue;
1387
1388 Eigen::VectorXd force;
1389 if (form->enabled())
1390 {
1391 form->first_derivative(sol, force);
1392 force *= -1.0 / s; // Divide by acceleration scaling to get units of force
1393 }
1394 else
1395 {
1396 force.setZero(sol.size());
1397 }
1398
1399 save_volume_vector_field(state, points, opts, name + "_forces", force, writer);
1400 }
1401 }
1402
1403 // if(problem->is_mixed())
1404 if (state.mixed_assembler != nullptr)
1405 {
1406 Eigen::MatrixXd interp_p;
1408 mesh, 1, // FIXME: state.disc_orders should use pressure discr orders, works only with sampler
1409 pressure_bases, state.disc_orders, state.polys, state.polys_3d, ref_element_sampler,
1410 points.rows(), pressure, interp_p, opts.use_sampler, opts.boundary_only);
1411
1412 if (obstacle.n_vertices() > 0)
1413 {
1414 interp_p.conservativeResize(interp_p.size() + obstacle.n_vertices(), 1);
1415 interp_p.bottomRows(obstacle.n_vertices()).setZero();
1416 }
1417
1418 if (opts.solve_export_to_file)
1419 writer.add_field("pressure", interp_p);
1420 else
1421 solution_frames.back().pressure = interp_p;
1422 }
1423
1424 if (obstacle.n_vertices() > 0)
1425 {
1426 discr.conservativeResize(discr.size() + obstacle.n_vertices(), 1);
1427 discr.bottomRows(obstacle.n_vertices()).setZero();
1428 }
1429
1431 writer.add_field("discr", discr);
1432
1433 if (problem.has_exact_sol())
1434 {
1435 if (opts.solve_export_to_file)
1436 {
1437 writer.add_field("exact", exact_fun);
1438 writer.add_field("error", err);
1439 }
1440 else
1441 {
1442 solution_frames.back().exact = exact_fun;
1443 solution_frames.back().error = err;
1444 }
1445 }
1446
1447 if (fun.cols() != 1)
1448 {
1449 std::vector<assembler::Assembler::NamedMatrix> vals, tvals;
1451 mesh, problem.is_scalar(), bases, gbases,
1452 state.disc_orders, state.polys, state.polys_3d,
1453 *state.assembler,
1454 ref_element_sampler, points.rows(), sol, t, vals, opts.use_sampler, opts.boundary_only);
1455
1456 for (auto &[_, v] : vals)
1458
1459 if (opts.scalar_values)
1460 {
1461 if (opts.solve_export_to_file)
1462 {
1463 for (const auto &[name, v] : vals)
1464 writer.add_field(name, v);
1465 }
1466 else if (vals.size() > 0)
1467 solution_frames.back().scalar_value = vals[0].second;
1468 }
1469
1470 if (opts.solve_export_to_file && opts.tensor_values)
1471 {
1473 mesh, problem.is_scalar(), bases, gbases, state.disc_orders,
1474 state.polys, state.polys_3d, *state.assembler, ref_element_sampler,
1475 points.rows(), sol, t, tvals, opts.use_sampler, opts.boundary_only);
1476
1477 for (auto &[_, v] : tvals)
1479
1480 for (const auto &[name, v] : tvals)
1481 {
1482 const int stride = mesh.dimension();
1483 assert(v.cols() % stride == 0);
1484
1485 for (int i = 0; i < v.cols(); i += stride)
1486 {
1487 const Eigen::MatrixXd tmp = v.middleCols(i, stride);
1488 assert(tmp.cols() == stride);
1489
1490 const int ii = (i / stride) + 1;
1491 writer.add_field(fmt::format("{:s}_{:d}", name, ii), tmp);
1492 }
1493 }
1494 }
1495
1496 if (!opts.use_spline)
1497 {
1499 mesh, problem.is_scalar(), state.n_bases, bases, gbases,
1500 state.disc_orders, state.polys, state.polys_3d, *state.assembler,
1501 ref_element_sampler, t, points.rows(), sol, vals, tvals,
1502 opts.use_sampler, opts.boundary_only);
1503
1504 if (obstacle.n_vertices() > 0)
1505 {
1506 for (auto &v : vals)
1507 {
1508 v.second.conservativeResize(v.second.size() + obstacle.n_vertices(), 1);
1509 v.second.bottomRows(obstacle.n_vertices()).setZero();
1510 }
1511 }
1512
1513 if (opts.scalar_values)
1514 {
1515 if (opts.solve_export_to_file)
1516 {
1517 for (const auto &v : vals)
1518 writer.add_field(fmt::format("{:s}_avg", v.first), v.second);
1519 }
1520 else if (vals.size() > 0)
1521 solution_frames.back().scalar_value_avg = vals[0].second;
1522 }
1523 // for(int i = 0; i < tvals.cols(); ++i){
1524 // const int ii = (i / mesh.dimension()) + 1;
1525 // const int jj = (i % mesh.dimension()) + 1;
1526 // writer.add_field("tensor_value_avg_" + std::to_string(ii) + std::to_string(jj), tvals.col(i));
1527 // }
1528 }
1529 }
1530
1531 if (opts.material_params)
1532 {
1533 const auto &params = assembler.parameters();
1534
1535 std::map<std::string, Eigen::MatrixXd> param_val;
1536 for (const auto &[p, _] : params)
1537 param_val[p] = Eigen::MatrixXd(points.rows(), 1);
1538 Eigen::MatrixXd rhos(points.rows(), 1);
1539
1540 Eigen::MatrixXd local_pts;
1541 Eigen::MatrixXi vis_faces_poly, vis_edges_poly;
1542
1543 int index = 0;
1544 const auto &sampler = ref_element_sampler;
1545 for (int e = 0; e < int(bases.size()); ++e)
1546 {
1547 const basis::ElementBases &gbs = gbases[e];
1548 const basis::ElementBases &bs = bases[e];
1549
1550 if (opts.use_sampler)
1551 {
1552 if (mesh.is_simplex(e))
1553 local_pts = sampler.simplex_points();
1554 else if (mesh.is_cube(e))
1555 local_pts = sampler.cube_points();
1556 else
1557 {
1558 if (mesh.is_volume())
1559 sampler.sample_polyhedron(polys_3d.at(e).first, polys_3d.at(e).second, local_pts, vis_faces_poly, vis_edges_poly);
1560 else
1561 sampler.sample_polygon(polys.at(e), local_pts, vis_faces_poly, vis_edges_poly);
1562 }
1563 }
1564 else
1565 {
1566 if (mesh.is_volume())
1567 {
1568 if (mesh.is_simplex(e))
1569 autogen::p_nodes_3d(disc_orders(e), local_pts);
1570 else if (mesh.is_cube(e))
1571 autogen::q_nodes_3d(disc_orders(e), local_pts);
1572 else
1573 continue;
1574 }
1575 else
1576 {
1577 if (mesh.is_simplex(e))
1578 autogen::p_nodes_2d(disc_orders(e), local_pts);
1579 else if (mesh.is_cube(e))
1580 autogen::q_nodes_2d(disc_orders(e), local_pts);
1581 else
1582 {
1583 const auto &mesh2d = static_cast<const mesh::Mesh2D &>(mesh);
1584 const int n_v = mesh2d.n_face_vertices(e);
1585 local_pts.resize(n_v, 2);
1586
1587 for (int j = 0; j < n_v; ++j)
1588 {
1589 local_pts.row(j) = mesh2d.point(mesh2d.face_vertex(e, j));
1590 }
1591 }
1592 }
1593 }
1594
1596 vals.compute(e, mesh.is_volume(), local_pts, bs, gbs);
1597
1598 for (int j = 0; j < vals.val.rows(); ++j)
1599 {
1600 for (const auto &[p, func] : params)
1601 param_val.at(p)(index) = func(local_pts.row(j), vals.val.row(j), t, e);
1602
1603 rhos(index) = density(local_pts.row(j), vals.val.row(j), t, e);
1604
1605 ++index;
1606 }
1607 }
1608
1609 assert(index == points.rows());
1610
1611 if (obstacle.n_vertices() > 0)
1612 {
1613 for (auto &[_, tmp] : param_val)
1614 {
1615 tmp.conservativeResize(tmp.size() + obstacle.n_vertices(), 1);
1616 tmp.bottomRows(obstacle.n_vertices()).setZero();
1617 }
1618
1619 rhos.conservativeResize(rhos.size() + obstacle.n_vertices(), 1);
1620 rhos.bottomRows(obstacle.n_vertices()).setZero();
1621 }
1622 for (const auto &[p, tmp] : param_val)
1623 writer.add_field(p, tmp);
1624 writer.add_field("rho", rhos);
1625 }
1626
1627 if (opts.body_ids)
1628 {
1629
1630 Eigen::MatrixXd ids(points.rows(), 1);
1631
1632 for (int i = 0; i < points.rows(); ++i)
1633 {
1634 ids(i) = mesh.get_body_id(el_id(i));
1635 }
1636
1637 if (obstacle.n_vertices() > 0)
1638 {
1639 ids.conservativeResize(ids.size() + obstacle.n_vertices(), 1);
1640 ids.bottomRows(obstacle.n_vertices()).setZero();
1641 }
1642
1643 writer.add_field("body_ids", ids);
1644 }
1645
1646 // interpolate_function(pts_index, rhs, fun, opts.boundary_only);
1647 // writer.add_field("rhs", fun);
1648
1649 if (fun.cols() != 1 && state.mixed_assembler == nullptr)
1650 {
1651 Eigen::MatrixXd traction_forces, traction_forces_fun;
1652 compute_traction_forces(state, sol, t, traction_forces, false);
1653
1655 mesh, problem.is_scalar(), bases, state.disc_orders,
1656 state.polys, state.polys_3d, ref_element_sampler,
1657 points.rows(), traction_forces, traction_forces_fun, opts.use_sampler, opts.boundary_only);
1658
1659 if (obstacle.n_vertices() > 0)
1660 {
1661 traction_forces_fun.conservativeResize(traction_forces_fun.rows() + obstacle.n_vertices(), traction_forces_fun.cols());
1662 traction_forces_fun.bottomRows(obstacle.n_vertices()).setZero();
1663 }
1664
1665 writer.add_field("traction_force", traction_forces_fun);
1666 }
1667
1668 if (fun.cols() != 1 && state.mixed_assembler == nullptr)
1669 {
1670 try
1671 {
1672 Eigen::MatrixXd potential_grad, potential_grad_fun;
1673 state.assembler->assemble_gradient(mesh.is_volume(), state.n_bases, bases, gbases, state.ass_vals_cache, t, dt, sol, sol, potential_grad);
1674
1676 mesh, problem.is_scalar(), bases, state.disc_orders,
1677 state.polys, state.polys_3d, ref_element_sampler,
1678 points.rows(), potential_grad, potential_grad_fun, opts.use_sampler, opts.boundary_only);
1679
1680 if (obstacle.n_vertices() > 0)
1681 {
1682 potential_grad_fun.conservativeResize(potential_grad_fun.rows() + obstacle.n_vertices(), potential_grad_fun.cols());
1683 potential_grad_fun.bottomRows(obstacle.n_vertices()).setZero();
1684 }
1685
1686 writer.add_field("gradient_of_potential", potential_grad_fun);
1687 }
1688 catch (std::exception &)
1689 {
1690 }
1691 }
1692
1693 // Write the solution last so it is the default for warp-by-vector
1694 if (opts.solve_export_to_file)
1695 writer.add_field("solution", fun);
1696 else
1697 solution_frames.back().solution = fun;
1698
1699 if (opts.solve_export_to_file)
1700 {
1701 if (obstacle.n_vertices() > 0)
1702 {
1703 const int orig_p = points.rows();
1704 points.conservativeResize(points.rows() + obstacle.n_vertices(), points.cols());
1705 points.bottomRows(obstacle.n_vertices()) = obstacle.v();
1706
1707 if (elements.empty())
1708 {
1709 for (int i = 0; i < tets.rows(); ++i)
1710 {
1711 elements.emplace_back();
1712 for (int j = 0; j < tets.cols(); ++j)
1713 elements.back().push_back(tets(i, j));
1714 }
1715 }
1716
1717 for (int i = 0; i < obstacle.get_face_connectivity().rows(); ++i)
1718 {
1719 elements.emplace_back();
1720 for (int j = 0; j < obstacle.get_face_connectivity().cols(); ++j)
1721 elements.back().push_back(obstacle.get_face_connectivity()(i, j) + orig_p);
1722 }
1723
1724 for (int i = 0; i < obstacle.get_edge_connectivity().rows(); ++i)
1725 {
1726 elements.emplace_back();
1727 for (int j = 0; j < obstacle.get_edge_connectivity().cols(); ++j)
1728 elements.back().push_back(obstacle.get_edge_connectivity()(i, j) + orig_p);
1729 }
1730
1731 for (int i = 0; i < obstacle.get_vertex_connectivity().size(); ++i)
1732 {
1733 elements.emplace_back();
1734 elements.back().push_back(obstacle.get_vertex_connectivity()(i) + orig_p);
1735 }
1736 }
1737
1738 if (elements.empty())
1739 writer.write_mesh(path, points, tets);
1740 else
1741 writer.write_mesh(path, points, elements, true, disc_orders.maxCoeff() == 1);
1742 }
1743 else
1744 {
1745 solution_frames.back().name = path;
1746 solution_frames.back().points = points;
1747 solution_frames.back().connectivity = tets;
1748 }
1749 }
1750
1752 const State &state,
1753 const Eigen::MatrixXd &points,
1754 const ExportOptions &opts,
1755 const std::string &name,
1756 const Eigen::VectorXd &field,
1757 paraviewo::ParaviewWriter &writer) const
1758 {
1759 Eigen::MatrixXd inerpolated_field;
1761 *state.mesh, state.problem->is_scalar(), state.bases, state.disc_orders,
1762 state.polys, state.polys_3d, ref_element_sampler,
1763 points.rows(), field, inerpolated_field, opts.use_sampler, opts.boundary_only);
1764
1765 if (state.obstacle.n_vertices() > 0)
1766 {
1767 inerpolated_field.conservativeResize(
1768 inerpolated_field.rows() + state.obstacle.n_vertices(), inerpolated_field.cols());
1769 inerpolated_field.bottomRows(state.obstacle.n_vertices()) =
1770 utils::unflatten(field.tail(state.obstacle.ndof()), inerpolated_field.cols());
1771 }
1772
1773 if (opts.solve_export_to_file)
1774 {
1775 writer.add_field(name, inerpolated_field);
1776 }
1777 // TODO: else save to solution frames
1778 }
1779
1781 const std::string &export_surface,
1782 const State &state,
1783 const Eigen::MatrixXd &sol,
1784 const Eigen::MatrixXd &pressure,
1785 const double t,
1786 const double dt_in,
1787 const ExportOptions &opts,
1788 const bool is_contact_enabled,
1789 std::vector<SolutionFrame> &solution_frames) const
1790 {
1791
1792 const Eigen::VectorXi &disc_orders = state.disc_orders;
1793 const auto &density = state.mass_matrix_assembler->density();
1794 const std::vector<basis::ElementBases> &bases = state.bases;
1795 const std::vector<basis::ElementBases> &pressure_bases = state.pressure_bases;
1796 const std::vector<basis::ElementBases> &gbases = state.geom_bases();
1797 const assembler::Assembler &assembler = *state.assembler;
1798 const assembler::Problem &problem = *state.problem;
1799 const mesh::Mesh &mesh = *state.mesh;
1800 int problem_dim = (problem.is_scalar() ? 1 : mesh.dimension());
1801
1802 Eigen::MatrixXd boundary_vis_vertices;
1803 Eigen::MatrixXd boundary_vis_local_vertices;
1804 Eigen::MatrixXi boundary_vis_elements;
1805 Eigen::MatrixXi boundary_vis_elements_ids;
1806 Eigen::MatrixXi boundary_vis_primitive_ids;
1807 Eigen::MatrixXd boundary_vis_normals;
1808 Eigen::MatrixXd displaced_boundary_vis_normals;
1809
1810 build_vis_boundary_mesh(mesh, bases, gbases, state.total_local_boundary, sol, problem_dim,
1811 boundary_vis_vertices, boundary_vis_local_vertices, boundary_vis_elements,
1812 boundary_vis_elements_ids, boundary_vis_primitive_ids, boundary_vis_normals,
1813 displaced_boundary_vis_normals);
1814
1815 Eigen::MatrixXd fun, interp_p, discr, vect, b_sidesets;
1816
1817 Eigen::MatrixXd lsol, lp, lgrad, lpgrad;
1818
1819 int actual_dim = 1;
1820 if (!problem.is_scalar())
1821 actual_dim = mesh.dimension();
1822
1823 discr.resize(boundary_vis_vertices.rows(), 1);
1824 fun.resize(boundary_vis_vertices.rows(), actual_dim);
1825 interp_p.resize(boundary_vis_vertices.rows(), 1);
1826 vect.resize(boundary_vis_vertices.rows(), mesh.dimension());
1827
1828 b_sidesets.resize(boundary_vis_vertices.rows(), 1);
1829 b_sidesets.setZero();
1830
1831 for (int i = 0; i < boundary_vis_vertices.rows(); ++i)
1832 {
1833 const auto s_id = mesh.get_boundary_id(boundary_vis_primitive_ids(i));
1834 if (s_id > 0)
1835 {
1836 b_sidesets(i) = s_id;
1837 }
1838
1839 const int el_index = boundary_vis_elements_ids(i);
1841 mesh, problem.is_scalar(), bases, gbases,
1842 el_index, boundary_vis_local_vertices.row(i), sol, lsol, lgrad);
1843 assert(lsol.size() == actual_dim);
1844 if (state.mixed_assembler != nullptr)
1845 {
1847 mesh, 1, pressure_bases, gbases,
1848 el_index, boundary_vis_local_vertices.row(i), pressure, lp, lpgrad);
1849 assert(lp.size() == 1);
1850 interp_p(i) = lp(0);
1851 }
1852
1853 discr(i) = disc_orders(el_index);
1854 for (int j = 0; j < actual_dim; ++j)
1855 {
1856 fun(i, j) = lsol(j);
1857 }
1858
1859 if (actual_dim == 1)
1860 {
1861 assert(lgrad.size() == mesh.dimension());
1862 for (int j = 0; j < mesh.dimension(); ++j)
1863 {
1864 vect(i, j) = lgrad(j);
1865 }
1866 }
1867 else
1868 {
1869 assert(lgrad.size() == actual_dim * actual_dim);
1870 std::vector<assembler::Assembler::NamedMatrix> tensor_flat;
1871 const basis::ElementBases &gbs = gbases[el_index];
1872 const basis::ElementBases &bs = bases[el_index];
1873 assembler.compute_tensor_value(assembler::OutputData(t, el_index, bs, gbs, boundary_vis_local_vertices.row(i), sol), tensor_flat);
1874 // TF computed only from cauchy stress
1875 assert(tensor_flat[0].first == "cauchy_stess");
1876 assert(tensor_flat[0].second.size() == actual_dim * actual_dim);
1877
1878 Eigen::Map<Eigen::MatrixXd> tensor(tensor_flat[0].second.data(), actual_dim, actual_dim);
1879 vect.row(i) = displaced_boundary_vis_normals.row(i) * tensor;
1880
1881 double area = 0;
1882 if (mesh.is_volume())
1883 {
1884 if (mesh.is_simplex(el_index))
1885 area = mesh.tri_area(boundary_vis_primitive_ids(i));
1886 else if (mesh.is_cube(el_index))
1887 area = mesh.quad_area(boundary_vis_primitive_ids(i));
1888 }
1889 else
1890 area = mesh.edge_length(boundary_vis_primitive_ids(i));
1891
1892 vect.row(i) *= area;
1893 }
1894 }
1895
1896 std::shared_ptr<paraviewo::ParaviewWriter> tmpw;
1897 if (opts.use_hdf5)
1898 tmpw = std::make_shared<paraviewo::HDF5VTUWriter>();
1899 else
1900 tmpw = std::make_shared<paraviewo::VTUWriter>();
1901 paraviewo::ParaviewWriter &writer = *tmpw;
1902
1903 if (opts.solve_export_to_file)
1904 {
1905
1906 writer.add_field("normals", boundary_vis_normals);
1907 writer.add_field("displaced_normals", displaced_boundary_vis_normals);
1908 if (state.mixed_assembler != nullptr)
1909 writer.add_field("pressure", interp_p);
1910 writer.add_field("discr", discr);
1911 writer.add_field("sidesets", b_sidesets);
1912
1913 if (actual_dim == 1)
1914 writer.add_field("solution_grad", vect);
1915 else
1916 {
1917 writer.add_field("traction_force", vect);
1918 }
1919 }
1920 else
1921 {
1922 if (state.mixed_assembler != nullptr)
1923 solution_frames.back().pressure = interp_p;
1924 }
1925
1926 if (opts.material_params)
1927 {
1928 const auto &params = assembler.parameters();
1929
1930 std::map<std::string, Eigen::MatrixXd> param_val;
1931 for (const auto &[p, _] : params)
1932 param_val[p] = Eigen::MatrixXd(boundary_vis_vertices.rows(), 1);
1933 Eigen::MatrixXd rhos(boundary_vis_vertices.rows(), 1);
1934
1935 for (int i = 0; i < boundary_vis_vertices.rows(); ++i)
1936 {
1937 double lambda, mu;
1938
1939 for (const auto &[p, func] : params)
1940 param_val.at(p)(i) = func(boundary_vis_local_vertices.row(i), boundary_vis_vertices.row(i), t, boundary_vis_elements_ids(i));
1941
1942 rhos(i) = density(boundary_vis_local_vertices.row(i), boundary_vis_vertices.row(i), t, boundary_vis_elements_ids(i));
1943 }
1944
1945 for (const auto &[p, tmp] : param_val)
1946 writer.add_field(p, tmp);
1947 writer.add_field("rho", rhos);
1948 }
1949
1950 if (opts.body_ids)
1951 {
1952
1953 Eigen::MatrixXd ids(boundary_vis_vertices.rows(), 1);
1954
1955 for (int i = 0; i < boundary_vis_vertices.rows(); ++i)
1956 {
1957 ids(i) = mesh.get_body_id(boundary_vis_elements_ids(i));
1958 }
1959
1960 writer.add_field("body_ids", ids);
1961 }
1962
1963 // Write the solution last so it is the default for warp-by-vector
1964 if (opts.solve_export_to_file)
1965 writer.add_field("solution", fun);
1966 else
1967 solution_frames.back().solution = fun;
1968
1969 if (opts.solve_export_to_file)
1970 writer.write_mesh(export_surface, boundary_vis_vertices, boundary_vis_elements);
1971 else
1972 {
1973 solution_frames.back().name = export_surface;
1974 solution_frames.back().points = boundary_vis_vertices;
1975 solution_frames.back().connectivity = boundary_vis_elements;
1976 }
1977 }
1978
1980 const std::string &export_surface,
1981 const State &state,
1982 const Eigen::MatrixXd &sol,
1983 const Eigen::MatrixXd &pressure,
1984 const double t,
1985 const double dt_in,
1986 const ExportOptions &opts,
1987 const bool is_contact_enabled,
1988 std::vector<SolutionFrame> &solution_frames) const
1989 {
1990 const mesh::Mesh &mesh = *state.mesh;
1991 const ipc::CollisionMesh &collision_mesh = state.collision_mesh;
1992 const double dhat = state.args["contact"]["dhat"];
1993 const double friction_coefficient = state.args["contact"]["friction_coefficient"];
1994 const double epsv = state.args["contact"]["epsv"];
1995 const std::shared_ptr<solver::ContactForm> &contact_form = state.solve_data.contact_form;
1996 const std::shared_ptr<solver::FrictionForm> &friction_form = state.solve_data.friction_form;
1997
1998 if (opts.solve_export_to_file)
1999 {
2000 std::shared_ptr<paraviewo::ParaviewWriter> tmpw;
2001 if (opts.use_hdf5)
2002 tmpw = std::make_shared<paraviewo::HDF5VTUWriter>();
2003 else
2004 tmpw = std::make_shared<paraviewo::VTUWriter>();
2005 paraviewo::ParaviewWriter &writer = *tmpw;
2006
2007 const int problem_dim = mesh.dimension();
2008 const Eigen::MatrixXd full_displacements = utils::unflatten(sol, problem_dim);
2009 const Eigen::MatrixXd surface_displacements = collision_mesh.map_displacements(full_displacements);
2010
2011 const Eigen::MatrixXd displaced_surface = collision_mesh.displace_vertices(full_displacements);
2012
2013 ipc::Collisions collision_set;
2014 collision_set.set_use_convergent_formulation(state.args["contact"]["use_convergent_formulation"]);
2015 collision_set.build(
2016 collision_mesh, displaced_surface, dhat,
2017 /*dmin=*/0, state.args["solver"]["contact"]["CCD"]["broad_phase"]);
2018
2019 ipc::BarrierPotential barrier_potential(dhat);
2020
2021 const double barrier_stiffness = contact_form != nullptr ? contact_form->barrier_stiffness() : 1;
2022
2023 if (opts.contact_forces)
2024 {
2025 Eigen::MatrixXd forces = -barrier_stiffness * barrier_potential.gradient(collision_set, collision_mesh, displaced_surface);
2026
2027 Eigen::MatrixXd forces_reshaped = utils::unflatten(forces, problem_dim);
2028
2029 assert(forces_reshaped.rows() == surface_displacements.rows());
2030 assert(forces_reshaped.cols() == surface_displacements.cols());
2031 writer.add_field("contact_forces", forces_reshaped);
2032 }
2033
2034 if (opts.friction_forces)
2035 {
2036 ipc::FrictionCollisions friction_collision_set;
2037 friction_collision_set.build(
2038 collision_mesh, displaced_surface, collision_set,
2039 barrier_potential, barrier_stiffness, friction_coefficient);
2040
2041 ipc::FrictionPotential friction_potential(epsv);
2042
2043 Eigen::MatrixXd velocities;
2044 if (state.solve_data.time_integrator != nullptr)
2045 velocities = state.solve_data.time_integrator->v_prev();
2046 else
2047 velocities = sol;
2048 velocities = collision_mesh.map_displacements(utils::unflatten(velocities, collision_mesh.dim()));
2049
2050 Eigen::MatrixXd forces = -friction_potential.gradient(
2051 friction_collision_set, collision_mesh, velocities);
2052
2053 Eigen::MatrixXd forces_reshaped = utils::unflatten(forces, problem_dim);
2054
2055 assert(forces_reshaped.rows() == surface_displacements.rows());
2056 assert(forces_reshaped.cols() == surface_displacements.cols());
2057 writer.add_field("friction_forces", forces_reshaped);
2058 }
2059
2060 assert(collision_mesh.rest_positions().rows() == surface_displacements.rows());
2061 assert(collision_mesh.rest_positions().cols() == surface_displacements.cols());
2062
2063 // Write the solution last so it is the default for warp-by-vector
2064 writer.add_field("solution", surface_displacements);
2065
2066 writer.write_mesh(
2067 export_surface.substr(0, export_surface.length() - 4) + "_contact.vtu",
2068 collision_mesh.rest_positions(),
2069 problem_dim == 3 ? collision_mesh.faces() : collision_mesh.edges());
2070 }
2071 }
2072
2074 const std::string &name,
2075 const State &state,
2076 const Eigen::MatrixXd &sol,
2077 const double t,
2078 const ExportOptions &opts,
2079 std::vector<SolutionFrame> &solution_frames) const
2080 {
2081 const std::vector<basis::ElementBases> &gbases = state.geom_bases();
2082 const mesh::Mesh &mesh = *state.mesh;
2083 const assembler::Problem &problem = *state.problem;
2084
2085 if (!opts.solve_export_to_file) // TODO?
2086 return;
2087 const auto &sampler = ref_element_sampler;
2088
2089 Eigen::MatrixXi vis_faces_poly, vis_edges_poly;
2090 Eigen::MatrixXd vis_pts_poly;
2091
2092 const auto &current_bases = gbases;
2093 int seg_total_size = 0;
2094 int pts_total_size = 0;
2095 int faces_total_size = 0;
2096
2097 for (size_t i = 0; i < current_bases.size(); ++i)
2098 {
2099 const auto &bs = current_bases[i];
2100
2101 if (mesh.is_simplex(i))
2102 {
2103 pts_total_size += sampler.simplex_points().rows();
2104 seg_total_size += sampler.simplex_edges().rows();
2105 faces_total_size += sampler.simplex_faces().rows();
2106 }
2107 else if (mesh.is_cube(i))
2108 {
2109 pts_total_size += sampler.cube_points().rows();
2110 seg_total_size += sampler.cube_edges().rows();
2111 faces_total_size += sampler.cube_faces().rows();
2112 }
2113 else
2114 {
2115 if (mesh.is_volume())
2116 sampler.sample_polyhedron(state.polys_3d.at(i).first, state.polys_3d.at(i).second, vis_pts_poly, vis_faces_poly, vis_edges_poly);
2117 else
2118 sampler.sample_polygon(state.polys.at(i), vis_pts_poly, vis_faces_poly, vis_edges_poly);
2119
2120 pts_total_size += vis_pts_poly.rows();
2121 seg_total_size += vis_edges_poly.rows();
2122 faces_total_size += vis_faces_poly.rows();
2123 }
2124 }
2125
2126 Eigen::MatrixXd points(pts_total_size, mesh.dimension());
2127 Eigen::MatrixXi edges(seg_total_size, 2);
2128 Eigen::MatrixXi faces(faces_total_size, 3);
2129 points.setZero();
2130
2131 Eigen::MatrixXd mapped, tmp;
2132 int seg_index = 0, pts_index = 0, face_index = 0;
2133 for (size_t i = 0; i < current_bases.size(); ++i)
2134 {
2135 const auto &bs = current_bases[i];
2136
2137 if (mesh.is_simplex(i))
2138 {
2139 bs.eval_geom_mapping(sampler.simplex_points(), mapped);
2140 edges.block(seg_index, 0, sampler.simplex_edges().rows(), edges.cols()) = sampler.simplex_edges().array() + pts_index;
2141 seg_index += sampler.simplex_edges().rows();
2142
2143 faces.block(face_index, 0, sampler.simplex_faces().rows(), 3) = sampler.simplex_faces().array() + pts_index;
2144 face_index += sampler.simplex_faces().rows();
2145
2146 points.block(pts_index, 0, mapped.rows(), points.cols()) = mapped;
2147 pts_index += mapped.rows();
2148 }
2149 else if (mesh.is_cube(i))
2150 {
2151 bs.eval_geom_mapping(sampler.cube_points(), mapped);
2152 edges.block(seg_index, 0, sampler.cube_edges().rows(), edges.cols()) = sampler.cube_edges().array() + pts_index;
2153 seg_index += sampler.cube_edges().rows();
2154
2155 faces.block(face_index, 0, sampler.cube_faces().rows(), 3) = sampler.cube_faces().array() + pts_index;
2156 face_index += sampler.cube_faces().rows();
2157
2158 points.block(pts_index, 0, mapped.rows(), points.cols()) = mapped;
2159 pts_index += mapped.rows();
2160 }
2161 else
2162 {
2163 if (mesh.is_volume())
2164 sampler.sample_polyhedron(state.polys_3d.at(i).first, state.polys_3d.at(i).second, vis_pts_poly, vis_faces_poly, vis_edges_poly);
2165 else
2166 sampler.sample_polygon(state.polys.at(i), vis_pts_poly, vis_faces_poly, vis_edges_poly);
2167
2168 edges.block(seg_index, 0, vis_edges_poly.rows(), edges.cols()) = vis_edges_poly.array() + pts_index;
2169 seg_index += vis_edges_poly.rows();
2170
2171 faces.block(face_index, 0, vis_faces_poly.rows(), 3) = vis_faces_poly.array() + pts_index;
2172 face_index += vis_faces_poly.rows();
2173
2174 points.block(pts_index, 0, vis_pts_poly.rows(), points.cols()) = vis_pts_poly;
2175 pts_index += vis_pts_poly.rows();
2176 }
2177 }
2178
2179 assert(pts_index == points.rows());
2180 assert(face_index == faces.rows());
2181
2182 if (mesh.is_volume())
2183 {
2184 // reverse all faces
2185 for (long i = 0; i < faces.rows(); ++i)
2186 {
2187 const int v0 = faces(i, 0);
2188 const int v1 = faces(i, 1);
2189 const int v2 = faces(i, 2);
2190
2191 int tmpc = faces(i, 2);
2192 faces(i, 2) = faces(i, 1);
2193 faces(i, 1) = tmpc;
2194 }
2195 }
2196 else
2197 {
2198 Eigen::Matrix2d mmat;
2199 for (long i = 0; i < faces.rows(); ++i)
2200 {
2201 const int v0 = faces(i, 0);
2202 const int v1 = faces(i, 1);
2203 const int v2 = faces(i, 2);
2204
2205 mmat.row(0) = points.row(v2) - points.row(v0);
2206 mmat.row(1) = points.row(v1) - points.row(v0);
2207
2208 if (mmat.determinant() > 0)
2209 {
2210 int tmpc = faces(i, 2);
2211 faces(i, 2) = faces(i, 1);
2212 faces(i, 1) = tmpc;
2213 }
2214 }
2215 }
2216
2217 Eigen::MatrixXd fun;
2219 mesh, problem.is_scalar(), state.bases, state.disc_orders,
2220 state.polys, state.polys_3d, ref_element_sampler,
2221 pts_index, sol, fun, /*use_sampler*/ true, false);
2222
2223 Eigen::MatrixXd exact_fun, err;
2224
2225 if (problem.has_exact_sol())
2226 {
2227 problem.exact(points, t, exact_fun);
2228 err = (fun - exact_fun).eval().rowwise().norm();
2229 }
2230
2231 std::shared_ptr<paraviewo::ParaviewWriter> tmpw;
2232 if (opts.use_hdf5)
2233 tmpw = std::make_shared<paraviewo::HDF5VTUWriter>();
2234 else
2235 tmpw = std::make_shared<paraviewo::VTUWriter>();
2236 paraviewo::ParaviewWriter &writer = *tmpw;
2237
2238 if (problem.has_exact_sol())
2239 {
2240 writer.add_field("exact", exact_fun);
2241 writer.add_field("error", err);
2242 }
2243
2244 if (fun.cols() != 1)
2245 {
2246 std::vector<assembler::Assembler::NamedMatrix> scalar_val;
2248 mesh, problem.is_scalar(), state.bases, gbases,
2249 state.disc_orders, state.polys, state.polys_3d,
2250 *state.assembler,
2251 ref_element_sampler, pts_index, sol, t, scalar_val, /*use_sampler*/ true, false);
2252 for (const auto &v : scalar_val)
2253 writer.add_field(v.first, v.second);
2254 }
2255 // Write the solution last so it is the default for warp-by-vector
2256 writer.add_field("solution", fun);
2257
2258 writer.write_mesh(name, points, edges);
2259 }
2260
2262 const std::string &path,
2263 const State &state,
2264 const Eigen::MatrixXd &sol,
2265 const ExportOptions &opts,
2266 std::vector<SolutionFrame> &solution_frames) const
2267 {
2268 const auto &dirichlet_nodes = state.dirichlet_nodes;
2269 const auto &dirichlet_nodes_position = state.dirichlet_nodes_position;
2270 const mesh::Mesh &mesh = *state.mesh;
2271 const assembler::Problem &problem = *state.problem;
2272
2273 int actual_dim = 1;
2274 if (!problem.is_scalar())
2275 actual_dim = mesh.dimension();
2276
2277 Eigen::MatrixXd fun(dirichlet_nodes_position.size(), actual_dim);
2278 Eigen::MatrixXd b_sidesets(dirichlet_nodes_position.size(), 1);
2279 b_sidesets.setZero();
2280 Eigen::MatrixXd points(dirichlet_nodes_position.size(), mesh.dimension());
2281 std::vector<std::vector<int>> cells(dirichlet_nodes_position.size());
2282
2283 for (int i = 0; i < dirichlet_nodes_position.size(); ++i)
2284 {
2285 const int n_id = dirichlet_nodes[i];
2286 const auto s_id = mesh.get_node_id(n_id);
2287 if (s_id > 0)
2288 {
2289 b_sidesets(i) = s_id;
2290 }
2291
2292 for (int j = 0; j < actual_dim; ++j)
2293 {
2294 fun(i, j) = sol(n_id * actual_dim + j);
2295 }
2296
2297 points.row(i) = dirichlet_nodes_position[i];
2298 cells[i].push_back(i);
2299 }
2300
2301 std::shared_ptr<paraviewo::ParaviewWriter> tmpw;
2302 if (opts.use_hdf5)
2303 tmpw = std::make_shared<paraviewo::HDF5VTUWriter>();
2304 else
2305 tmpw = std::make_shared<paraviewo::VTUWriter>();
2306 paraviewo::ParaviewWriter &writer = *tmpw;
2307
2308 if (opts.solve_export_to_file)
2309 {
2310 writer.add_field("sidesets", b_sidesets);
2311 // Write the solution last so it is the default for warp-by-vector
2312 writer.add_field("solution", fun);
2313 writer.write_mesh(path, points, cells, false, false);
2314 }
2315 }
2316
2318 const std::string &name,
2319 const std::function<std::string(int)> &vtu_names,
2320 int time_steps, double t0, double dt, int skip_frame) const
2321 {
2322 paraviewo::PVDWriter::save_pvd(name, vtu_names, time_steps, t0, dt, skip_frame);
2323 }
2324
2325 void OutGeometryData::init_sampler(const polyfem::mesh::Mesh &mesh, const double vismesh_rel_area)
2326 {
2327 ref_element_sampler.init(mesh.is_volume(), mesh.n_elements(), vismesh_rel_area);
2328 }
2329
2330 void OutGeometryData::build_grid(const polyfem::mesh::Mesh &mesh, const double spacing)
2331 {
2332 if (spacing <= 0)
2333 return;
2334
2335 RowVectorNd min, max;
2336 mesh.bounding_box(min, max);
2337 const RowVectorNd delta = max - min;
2338 const int nx = delta[0] / spacing + 1;
2339 const int ny = delta[1] / spacing + 1;
2340 const int nz = delta.cols() >= 3 ? (delta[2] / spacing + 1) : 1;
2341 const int n = nx * ny * nz;
2342
2343 grid_points.resize(n, delta.cols());
2344 int index = 0;
2345 for (int i = 0; i < nx; ++i)
2346 {
2347 const double x = (delta[0] / (nx - 1)) * i + min[0];
2348
2349 for (int j = 0; j < ny; ++j)
2350 {
2351 const double y = (delta[1] / (ny - 1)) * j + min[1];
2352
2353 if (delta.cols() <= 2)
2354 {
2355 grid_points.row(index++) << x, y;
2356 }
2357 else
2358 {
2359 for (int k = 0; k < nz; ++k)
2360 {
2361 const double z = (delta[2] / (nz - 1)) * k + min[2];
2362 grid_points.row(index++) << x, y, z;
2363 }
2364 }
2365 }
2366 }
2367
2368 assert(index == n);
2369
2370 std::vector<std::array<Eigen::Vector3d, 2>> boxes;
2371 mesh.elements_boxes(boxes);
2372
2373 SimpleBVH::BVH bvh;
2374 bvh.init(boxes);
2375
2376 const double eps = 1e-6;
2377
2378 grid_points_to_elements.resize(grid_points.rows(), 1);
2379 grid_points_to_elements.setConstant(-1);
2380
2381 grid_points_bc.resize(grid_points.rows(), mesh.is_volume() ? 4 : 3);
2382
2383 for (int i = 0; i < grid_points.rows(); ++i)
2384 {
2385 const Eigen::Vector3d min(
2386 grid_points(i, 0) - eps,
2387 grid_points(i, 1) - eps,
2388 (mesh.is_volume() ? grid_points(i, 2) : 0) - eps);
2389
2390 const Eigen::Vector3d max(
2391 grid_points(i, 0) + eps,
2392 grid_points(i, 1) + eps,
2393 (mesh.is_volume() ? grid_points(i, 2) : 0) + eps);
2394
2395 std::vector<unsigned int> candidates;
2396
2397 bvh.intersect_box(min, max, candidates);
2398
2399 for (const auto cand : candidates)
2400 {
2401 if (!mesh.is_simplex(cand))
2402 {
2403 logger().warn("Element {} is not simplex, skipping", cand);
2404 continue;
2405 }
2406
2407 Eigen::MatrixXd coords;
2408 mesh.barycentric_coords(grid_points.row(i), cand, coords);
2409
2410 for (int d = 0; d < coords.size(); ++d)
2411 {
2412 if (fabs(coords(d)) < 1e-8)
2413 coords(d) = 0;
2414 else if (fabs(coords(d) - 1) < 1e-8)
2415 coords(d) = 1;
2416 }
2417
2418 if (coords.array().minCoeff() >= 0 && coords.array().maxCoeff() <= 1)
2419 {
2420 grid_points_to_elements(i) = cand;
2421 grid_points_bc.row(i) = coords;
2422 break;
2423 }
2424 }
2425 }
2426 }
2427
2428 void OutStatsData::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)
2429 {
2430 Eigen::MatrixXd samples_simplex, samples_cube, mapped, p0, p1, p;
2431
2432 mesh_size = 0;
2433 average_edge_length = 0;
2434 min_edge_length = std::numeric_limits<double>::max();
2435
2436 if (!use_curved_mesh_size)
2437 {
2438 mesh_in.get_edges(p0, p1);
2439 p = p0 - p1;
2440 min_edge_length = p.rowwise().norm().minCoeff();
2441 average_edge_length = p.rowwise().norm().mean();
2442 mesh_size = p.rowwise().norm().maxCoeff();
2443
2444 logger().info("hmin: {}", min_edge_length);
2445 logger().info("hmax: {}", mesh_size);
2446 logger().info("havg: {}", average_edge_length);
2447
2448 return;
2449 }
2450
2451 if (mesh_in.is_volume())
2452 {
2453 utils::EdgeSampler::sample_3d_simplex(n_samples, samples_simplex);
2454 utils::EdgeSampler::sample_3d_cube(n_samples, samples_cube);
2455 }
2456 else
2457 {
2458 utils::EdgeSampler::sample_2d_simplex(n_samples, samples_simplex);
2459 utils::EdgeSampler::sample_2d_cube(n_samples, samples_cube);
2460 }
2461
2462 int n = 0;
2463 for (size_t i = 0; i < bases_in.size(); ++i)
2464 {
2465 if (mesh_in.is_polytope(i))
2466 continue;
2467 int n_edges;
2468
2469 if (mesh_in.is_simplex(i))
2470 {
2471 n_edges = mesh_in.is_volume() ? 6 : 3;
2472 bases_in[i].eval_geom_mapping(samples_simplex, mapped);
2473 }
2474 else
2475 {
2476 n_edges = mesh_in.is_volume() ? 12 : 4;
2477 bases_in[i].eval_geom_mapping(samples_cube, mapped);
2478 }
2479
2480 for (int j = 0; j < n_edges; ++j)
2481 {
2482 double current_edge = 0;
2483 for (int k = 0; k < n_samples - 1; ++k)
2484 {
2485 p0 = mapped.row(j * n_samples + k);
2486 p1 = mapped.row(j * n_samples + k + 1);
2487 p = p0 - p1;
2488
2489 current_edge += p.norm();
2490 }
2491
2492 mesh_size = std::max(current_edge, mesh_size);
2493 min_edge_length = std::min(current_edge, min_edge_length);
2494 average_edge_length += current_edge;
2495 ++n;
2496 }
2497 }
2498
2499 average_edge_length /= n;
2500
2501 logger().info("hmin: {}", min_edge_length);
2502 logger().info("hmax: {}", mesh_size);
2503 logger().info("havg: {}", average_edge_length);
2504 }
2505
2507 {
2508 sigma_avg = 0;
2509 sigma_max = 0;
2510 sigma_min = 0;
2511
2512 n_flipped = 0;
2513 }
2514
2515 void OutStatsData::count_flipped_elements(const polyfem::mesh::Mesh &mesh, const std::vector<polyfem::basis::ElementBases> &gbases)
2516 {
2517 using namespace mesh;
2518
2519 logger().info("Counting flipped elements...");
2520 const auto &els_tag = mesh.elements_tag();
2521
2522 // flipped_elements.clear();
2523 for (size_t i = 0; i < gbases.size(); ++i)
2524 {
2525 if (mesh.is_polytope(i))
2526 continue;
2527
2529 if (!vals.is_geom_mapping_positive(mesh.is_volume(), gbases[i]))
2530 {
2531 ++n_flipped;
2532
2533 static const std::vector<std::string> element_type_names{{
2534 "Simplex",
2535 "RegularInteriorCube",
2536 "RegularBoundaryCube",
2537 "SimpleSingularInteriorCube",
2538 "MultiSingularInteriorCube",
2539 "SimpleSingularBoundaryCube",
2540 "InterfaceCube",
2541 "MultiSingularBoundaryCube",
2542 "BoundaryPolytope",
2543 "InteriorPolytope",
2544 "Undefined",
2545 }};
2546
2547 log_and_throw_error("element {} is flipped, type {}", i, element_type_names[static_cast<int>(els_tag[i])]);
2548 }
2549 }
2550
2551 logger().info(" done");
2552
2553 // dynamic_cast<Mesh3D *>(mesh.get())->save({56}, 1, "mesh.HYBRID");
2554
2555 // std::sort(flipped_elements.begin(), flipped_elements.end());
2556 // auto it = std::unique(flipped_elements.begin(), flipped_elements.end());
2557 // flipped_elements.resize(std::distance(flipped_elements.begin(), it));
2558 }
2559
2561 const int n_bases,
2562 const std::vector<polyfem::basis::ElementBases> &bases,
2563 const std::vector<polyfem::basis::ElementBases> &gbases,
2564 const polyfem::mesh::Mesh &mesh,
2565 const assembler::Problem &problem,
2566 const double tend,
2567 const Eigen::MatrixXd &sol)
2568 {
2569 if (n_bases <= 0)
2570 {
2571 logger().error("Build the bases first!");
2572 return;
2573 }
2574 if (sol.size() <= 0)
2575 {
2576 logger().error("Solve the problem first!");
2577 return;
2578 }
2579
2580 int actual_dim = 1;
2581 if (!problem.is_scalar())
2582 actual_dim = mesh.dimension();
2583
2584 igl::Timer timer;
2585 timer.start();
2586 logger().info("Computing errors...");
2587 using std::max;
2588
2589 const int n_el = int(bases.size());
2590
2591 Eigen::MatrixXd v_exact, v_approx;
2592 Eigen::MatrixXd v_exact_grad(0, 0), v_approx_grad;
2593
2594 l2_err = 0;
2595 h1_err = 0;
2596 grad_max_err = 0;
2597 h1_semi_err = 0;
2598 linf_err = 0;
2599 lp_err = 0;
2600 // double pred_norm = 0;
2601
2602 static const int p = 8;
2603
2604 // Eigen::MatrixXd err_per_el(n_el, 5);
2606
2607 for (int e = 0; e < n_el; ++e)
2608 {
2609 vals.compute(e, mesh.is_volume(), bases[e], gbases[e]);
2610
2611 if (problem.has_exact_sol())
2612 {
2613 problem.exact(vals.val, tend, v_exact);
2614 problem.exact_grad(vals.val, tend, v_exact_grad);
2615 }
2616
2617 v_approx.resize(vals.val.rows(), actual_dim);
2618 v_approx.setZero();
2619
2620 v_approx_grad.resize(vals.val.rows(), mesh.dimension() * actual_dim);
2621 v_approx_grad.setZero();
2622
2623 const int n_loc_bases = int(vals.basis_values.size());
2624
2625 for (int i = 0; i < n_loc_bases; ++i)
2626 {
2627 const auto &val = vals.basis_values[i];
2628
2629 for (size_t ii = 0; ii < val.global.size(); ++ii)
2630 {
2631 for (int d = 0; d < actual_dim; ++d)
2632 {
2633 v_approx.col(d) += val.global[ii].val * sol(val.global[ii].index * actual_dim + d) * val.val;
2634 v_approx_grad.block(0, d * val.grad_t_m.cols(), v_approx_grad.rows(), val.grad_t_m.cols()) += val.global[ii].val * sol(val.global[ii].index * actual_dim + d) * val.grad_t_m;
2635 }
2636 }
2637 }
2638
2639 const auto err = problem.has_exact_sol() ? (v_exact - v_approx).eval().rowwise().norm().eval() : (v_approx).eval().rowwise().norm().eval();
2640 const auto err_grad = problem.has_exact_sol() ? (v_exact_grad - v_approx_grad).eval().rowwise().norm().eval() : (v_approx_grad).eval().rowwise().norm().eval();
2641
2642 // for(long i = 0; i < err.size(); ++i)
2643 // errors.push_back(err(i));
2644
2645 linf_err = std::max(linf_err, err.maxCoeff());
2646 grad_max_err = std::max(linf_err, err_grad.maxCoeff());
2647
2648 // {
2649 // const auto &mesh3d = *dynamic_cast<Mesh3D *>(mesh.get());
2650 // const auto v0 = mesh3d.point(mesh3d.cell_vertex(e, 0));
2651 // const auto v1 = mesh3d.point(mesh3d.cell_vertex(e, 1));
2652 // const auto v2 = mesh3d.point(mesh3d.cell_vertex(e, 2));
2653 // const auto v3 = mesh3d.point(mesh3d.cell_vertex(e, 3));
2654
2655 // Eigen::Matrix<double, 6, 3> ee;
2656 // ee.row(0) = v0 - v1;
2657 // ee.row(1) = v1 - v2;
2658 // ee.row(2) = v2 - v0;
2659
2660 // ee.row(3) = v0 - v3;
2661 // ee.row(4) = v1 - v3;
2662 // ee.row(5) = v2 - v3;
2663
2664 // Eigen::Matrix<double, 6, 1> en = ee.rowwise().norm();
2665
2666 // // Eigen::Matrix<double, 3*4, 1> alpha;
2667 // // alpha(0) = angle3(e.row(0), -e.row(1)); alpha(1) = angle3(e.row(1), -e.row(2)); alpha(2) = angle3(e.row(2), -e.row(0));
2668 // // alpha(3) = angle3(e.row(0), -e.row(4)); alpha(4) = angle3(e.row(4), e.row(3)); alpha(5) = angle3(-e.row(3), -e.row(0));
2669 // // alpha(6) = angle3(-e.row(4), -e.row(1)); alpha(7) = angle3(e.row(1), -e.row(5)); alpha(8) = angle3(e.row(5), e.row(4));
2670 // // alpha(9) = angle3(-e.row(2), -e.row(5)); alpha(10) = angle3(e.row(5), e.row(3)); alpha(11) = angle3(-e.row(3), e.row(2));
2671
2672 // const double S = (ee.row(0).cross(ee.row(1)).norm() + ee.row(0).cross(ee.row(4)).norm() + ee.row(4).cross(ee.row(1)).norm() + ee.row(2).cross(ee.row(5)).norm()) / 2;
2673 // const double V = std::abs(ee.row(3).dot(ee.row(2).cross(-ee.row(0))))/6;
2674 // const double rho = 3 * V / S;
2675 // const double hp = en.maxCoeff();
2676 // const int pp = disc_orders(e);
2677 // const int p_ref = args["space"]["discr_order"];
2678
2679 // err_per_el(e, 0) = err.mean();
2680 // err_per_el(e, 1) = err.maxCoeff();
2681 // err_per_el(e, 2) = std::pow(hp, pp+1)/(rho/hp); // /std::pow(average_edge_length, p_ref+1) * (sqrt(6)/12);
2682 // err_per_el(e, 3) = rho/hp;
2683 // err_per_el(e, 4) = (vals.det.array() * vals.quadrature.weights.array()).sum();
2684
2685 // // pred_norm += (pow(std::pow(hp, pp+1)/(rho/hp),p) * vals.det.array() * vals.quadrature.weights.array()).sum();
2686 // }
2687
2688 l2_err += (err.array() * err.array() * vals.det.array() * vals.quadrature.weights.array()).sum();
2689 h1_err += (err_grad.array() * err_grad.array() * vals.det.array() * vals.quadrature.weights.array()).sum();
2690 lp_err += (err.array().pow(p) * vals.det.array() * vals.quadrature.weights.array()).sum();
2691 }
2692
2693 h1_semi_err = sqrt(fabs(h1_err));
2694 h1_err = sqrt(fabs(l2_err) + fabs(h1_err));
2695 l2_err = sqrt(fabs(l2_err));
2696
2697 lp_err = pow(fabs(lp_err), 1. / p);
2698
2699 // pred_norm = pow(fabs(pred_norm), 1./p);
2700
2701 timer.stop();
2702 const double computing_errors_time = timer.getElapsedTime();
2703 logger().info(" took {}s", computing_errors_time);
2704
2705 logger().info("-- L2 error: {}", l2_err);
2706 logger().info("-- Lp error: {}", lp_err);
2707 logger().info("-- H1 error: {}", h1_err);
2708 logger().info("-- H1 semi error: {}", h1_semi_err);
2709 // logger().info("-- Perd norm: {}", pred_norm);
2710
2711 logger().info("-- Linf error: {}", linf_err);
2712 logger().info("-- grad max error: {}", grad_max_err);
2713
2714 // {
2715 // std::ofstream out("errs.txt");
2716 // out<<err_per_el;
2717 // out.close();
2718 // }
2719 }
2720
2722 {
2723 using namespace polyfem::mesh;
2724
2725 simplex_count = 0;
2726 regular_count = 0;
2727 regular_boundary_count = 0;
2728 simple_singular_count = 0;
2729 multi_singular_count = 0;
2730 boundary_count = 0;
2731 non_regular_boundary_count = 0;
2732 non_regular_count = 0;
2733 undefined_count = 0;
2734 multi_singular_boundary_count = 0;
2735
2736 const auto &els_tag = mesh.elements_tag();
2737
2738 for (size_t i = 0; i < els_tag.size(); ++i)
2739 {
2740 const ElementType type = els_tag[i];
2741
2742 switch (type)
2743 {
2744 case ElementType::SIMPLEX:
2745 simplex_count++;
2746 break;
2747 case ElementType::REGULAR_INTERIOR_CUBE:
2748 regular_count++;
2749 break;
2750 case ElementType::REGULAR_BOUNDARY_CUBE:
2751 regular_boundary_count++;
2752 break;
2753 case ElementType::SIMPLE_SINGULAR_INTERIOR_CUBE:
2754 simple_singular_count++;
2755 break;
2756 case ElementType::MULTI_SINGULAR_INTERIOR_CUBE:
2757 multi_singular_count++;
2758 break;
2759 case ElementType::SIMPLE_SINGULAR_BOUNDARY_CUBE:
2760 boundary_count++;
2761 break;
2762 case ElementType::INTERFACE_CUBE:
2763 case ElementType::MULTI_SINGULAR_BOUNDARY_CUBE:
2764 multi_singular_boundary_count++;
2765 break;
2766 case ElementType::BOUNDARY_POLYTOPE:
2767 non_regular_boundary_count++;
2768 break;
2769 case ElementType::INTERIOR_POLYTOPE:
2770 non_regular_count++;
2771 break;
2772 case ElementType::UNDEFINED:
2773 undefined_count++;
2774 break;
2775 default:
2776 throw std::runtime_error("Unknown element type");
2777 }
2778 }
2779
2780 logger().info("simplex_count: \t{}", simplex_count);
2781 logger().info("regular_count: \t{}", regular_count);
2782 logger().info("regular_boundary_count: \t{}", regular_boundary_count);
2783 logger().info("simple_singular_count: \t{}", simple_singular_count);
2784 logger().info("multi_singular_count: \t{}", multi_singular_count);
2785 logger().info("boundary_count: \t{}", boundary_count);
2786 logger().info("multi_singular_boundary_count: \t{}", multi_singular_boundary_count);
2787 logger().info("non_regular_count: \t{}", non_regular_count);
2788 logger().info("non_regular_boundary_count: \t{}", non_regular_boundary_count);
2789 logger().info("undefined_count: \t{}", undefined_count);
2790 logger().info("total count:\t {}", mesh.n_elements());
2791 }
2792
2794 const nlohmann::json &args,
2795 const int n_bases, const int n_pressure_bases,
2796 const Eigen::MatrixXd &sol,
2797 const mesh::Mesh &mesh,
2798 const Eigen::VectorXi &disc_orders,
2799 const assembler::Problem &problem,
2800 const OutRuntimeData &runtime,
2801 const std::string &formulation,
2802 const bool isoparametric,
2803 const int sol_at_node_id,
2804 nlohmann::json &j)
2805 {
2806
2807 j["args"] = args;
2808
2809 j["geom_order"] = mesh.orders().size() > 0 ? mesh.orders().maxCoeff() : 1;
2810 j["geom_order_min"] = mesh.orders().size() > 0 ? mesh.orders().minCoeff() : 1;
2811 j["discr_order_min"] = disc_orders.minCoeff();
2812 j["discr_order_max"] = disc_orders.maxCoeff();
2813 j["iso_parametric"] = isoparametric;
2814 j["problem"] = problem.name();
2815 j["mat_size"] = mat_size;
2816 j["num_bases"] = n_bases;
2817 j["num_pressure_bases"] = n_pressure_bases;
2818 j["num_non_zero"] = nn_zero;
2819 j["num_flipped"] = n_flipped;
2820 j["num_dofs"] = num_dofs;
2821 j["num_vertices"] = mesh.n_vertices();
2822 j["num_elements"] = mesh.n_elements();
2823
2824 j["num_p1"] = (disc_orders.array() == 1).count();
2825 j["num_p2"] = (disc_orders.array() == 2).count();
2826 j["num_p3"] = (disc_orders.array() == 3).count();
2827 j["num_p4"] = (disc_orders.array() == 4).count();
2828 j["num_p5"] = (disc_orders.array() == 5).count();
2829
2830 j["mesh_size"] = mesh_size;
2831 j["max_angle"] = max_angle;
2832
2833 j["sigma_max"] = sigma_max;
2834 j["sigma_min"] = sigma_min;
2835 j["sigma_avg"] = sigma_avg;
2836
2837 j["min_edge_length"] = min_edge_length;
2838 j["average_edge_length"] = average_edge_length;
2839
2840 j["err_l2"] = l2_err;
2841 j["err_h1"] = h1_err;
2842 j["err_h1_semi"] = h1_semi_err;
2843 j["err_linf"] = linf_err;
2844 j["err_linf_grad"] = grad_max_err;
2845 j["err_lp"] = lp_err;
2846
2847 j["spectrum"] = {spectrum(0), spectrum(1), spectrum(2), spectrum(3)};
2848 j["spectrum_condest"] = std::abs(spectrum(3)) / std::abs(spectrum(0));
2849
2850 // j["errors"] = errors;
2851
2852 j["time_building_basis"] = runtime.building_basis_time;
2853 j["time_loading_mesh"] = runtime.loading_mesh_time;
2854 j["time_computing_poly_basis"] = runtime.computing_poly_basis_time;
2855 j["time_assembling_stiffness_mat"] = runtime.assembling_stiffness_mat_time;
2856 j["time_assembling_mass_mat"] = runtime.assembling_mass_mat_time;
2857 j["time_assigning_rhs"] = runtime.assigning_rhs_time;
2858 j["time_solving"] = runtime.solving_time;
2859 // j["time_computing_errors"] = runtime.computing_errors_time;
2860
2861 j["solver_info"] = solver_info;
2862
2863 j["count_simplex"] = simplex_count;
2864 j["count_regular"] = regular_count;
2865 j["count_regular_boundary"] = regular_boundary_count;
2866 j["count_simple_singular"] = simple_singular_count;
2867 j["count_multi_singular"] = multi_singular_count;
2868 j["count_boundary"] = boundary_count;
2869 j["count_non_regular_boundary"] = non_regular_boundary_count;
2870 j["count_non_regular"] = non_regular_count;
2871 j["count_undefined"] = undefined_count;
2872 j["count_multi_singular_boundary"] = multi_singular_boundary_count;
2873
2874 j["is_simplicial"] = mesh.n_elements() == simplex_count;
2875
2876 j["peak_memory"] = getPeakRSS() / (1024 * 1024);
2877
2878 const int actual_dim = problem.is_scalar() ? 1 : mesh.dimension();
2879
2880 std::vector<double> mmin(actual_dim);
2881 std::vector<double> mmax(actual_dim);
2882
2883 for (int d = 0; d < actual_dim; ++d)
2884 {
2885 mmin[d] = std::numeric_limits<double>::max();
2886 mmax[d] = -std::numeric_limits<double>::max();
2887 }
2888
2889 for (int i = 0; i < sol.size(); i += actual_dim)
2890 {
2891 for (int d = 0; d < actual_dim; ++d)
2892 {
2893 mmin[d] = std::min(mmin[d], sol(i + d));
2894 mmax[d] = std::max(mmax[d], sol(i + d));
2895 }
2896 }
2897
2898 std::vector<double> sol_at_node(actual_dim);
2899
2900 if (sol_at_node_id >= 0)
2901 {
2902 const int node_id = sol_at_node_id;
2903
2904 for (int d = 0; d < actual_dim; ++d)
2905 {
2906 sol_at_node[d] = sol(node_id * actual_dim + d);
2907 }
2908 }
2909
2910 j["sol_at_node"] = sol_at_node;
2911 j["sol_min"] = mmin;
2912 j["sol_max"] = mmax;
2913
2914#if defined(POLYFEM_WITH_CPP_THREADS)
2915 j["num_threads"] = utils::get_n_threads();
2916#elif defined(POLYFEM_WITH_TBB)
2917 j["num_threads"] = utils::get_n_threads();
2918#else
2919 j["num_threads"] = 1;
2920#endif
2921
2922 j["formulation"] = formulation;
2923
2924 logger().info("done");
2925 }
2926
2927 EnergyCSVWriter::EnergyCSVWriter(const std::string &path, const solver::SolveData &solve_data)
2928 : file(path), solve_data(solve_data)
2929 {
2930 file << "i,";
2931 for (const auto &[name, _] : solve_data.named_forms())
2932 {
2933 file << name << ",";
2934 }
2935 file << "total_energy" << std::endl;
2936 }
2937
2939 {
2940 file.close();
2941 }
2942
2943 void EnergyCSVWriter::write(const int i, const Eigen::MatrixXd &sol)
2944 {
2945 const double s = solve_data.time_integrator
2946 ? solve_data.time_integrator->acceleration_scaling()
2947 : 1;
2948 file << i << ",";
2949 for (const auto &[_, form] : solve_data.named_forms())
2950 {
2951 // Divide by acceleration scaling to get the energy (units of J)
2952 file << ((form && form->enabled()) ? form->value(sol) : 0) / s << ",";
2953 }
2954 file << solve_data.nl_problem->value(sol) / s << "\n";
2955 file.flush();
2956 }
2957
2958 RuntimeStatsCSVWriter::RuntimeStatsCSVWriter(const std::string &path, const State &state, const double t0, const double dt)
2959 : file(path), state(state), t0(t0), dt(dt)
2960 {
2961 file << "step,time,forward,remeshing,global_relaxation,peak_mem,#V,#T" << std::endl;
2962 }
2963
2968
2969 void RuntimeStatsCSVWriter::write(const int t, const double forward, const double remeshing, const double global_relaxation, const Eigen::MatrixXd &sol)
2970 {
2971 total_forward_solve_time += forward;
2972 total_remeshing_time += remeshing;
2973 total_global_relaxation_time += global_relaxation;
2974
2975 // logger().debug(
2976 // "Forward (cur, avg, total): {} s, {} s, {} s",
2977 // forward, total_forward_solve_time / t, total_forward_solve_time);
2978 // logger().debug(
2979 // "Remeshing (cur, avg, total): {} s, {} s, {} s",
2980 // remeshing, total_remeshing_time / t, total_remeshing_time);
2981 // logger().debug(
2982 // "Global relaxation (cur, avg, total): {} s, {} s, {} s",
2983 // global_relaxation, total_global_relaxation_time / t, total_global_relaxation_time);
2984
2985 const double peak_mem = getPeakRSS() / double(1 << 30);
2986 // logger().debug("Peak mem: {} GiB", peak_mem);
2987
2988 file << fmt::format(
2989 "{},{},{},{},{},{},{},{}\n",
2990 t, t0 + dt * t, forward, remeshing, global_relaxation, peak_mem,
2991 state.n_bases, state.mesh->n_elements());
2992 file.flush();
2993 }
2994
2995} // namespace polyfem::io
double val
Definition Assembler.cpp:86
ElementAssemblyValues vals
Definition Assembler.cpp:22
std::vector< Eigen::VectorXi > faces
int y
int z
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
int n_bases
number of bases
Definition State.hpp:178
assembler::AssemblyValsCache ass_vals_cache
used to store assembly values for small problems
Definition State.hpp:196
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:223
Eigen::VectorXi in_node_to_node
Inpute nodes (including high-order) to polyfem nodes, only for isoparametric.
Definition State.hpp:450
mesh::Obstacle obstacle
Obstacles used in collisions.
Definition State.hpp:468
std::shared_ptr< assembler::Assembler > assembler
assemblers
Definition State.hpp:155
ipc::CollisionMesh collision_mesh
IPC collision mesh.
Definition State.hpp:515
std::vector< basis::ElementBases > pressure_bases
FE pressure bases for mixed elements, the size is #elements.
Definition State.hpp:173
std::shared_ptr< assembler::Mass > mass_matrix_assembler
Definition State.hpp:157
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:466
std::vector< int > dirichlet_nodes
per node dirichlet
Definition State.hpp:443
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
json args
main input arguments containing all defaults
Definition State.hpp:101
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:171
std::vector< mesh::LocalBoundary > total_local_boundary
mapping from elements to nodes for all mesh
Definition State.hpp:431
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:324
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::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
virtual std::map< std::string, ParamFunc > parameters() const =0
virtual void compute_tensor_value(const OutputData &data, std::vector< NamedMatrix > &result) const
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,...
std::vector< Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > > jac_it
const std::string & name() const
Definition Problem.hpp:23
virtual void exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
Definition Problem.hpp:46
virtual bool is_scalar() const =0
virtual bool has_exact_sol() const =0
virtual void exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
Definition Problem.hpp:45
virtual bool is_time_dependent() const
Definition Problem.hpp:50
Represents one basis function and its gradient.
Definition Basis.hpp:44
const std::vector< Local2Global > & global() const
Definition Basis.hpp:104
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d).
Eigen::VectorXi local_nodes_for_primitive(const int local_index, const mesh::Mesh &mesh) const
std::vector< Basis > bases
one basis function per node in the element
EnergyCSVWriter(const std::string &path, const solver::SolveData &solve_data)
Definition OutData.cpp:2927
const solver::SolveData & solve_data
Definition OutData.hpp:520
void write(const int i, const Eigen::MatrixXd &sol)
Definition OutData.cpp:2943
static void interpolate_at_local_vals(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const int el_index, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &fun, Eigen::MatrixXd &result, Eigen::MatrixXd &result_grad)
interpolate solution and gradient at element (calls interpolate_at_local_vals with sol)
static void compute_stress_at_quadrature_points(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXi &disc_orders, const assembler::Assembler &assembler, const Eigen::MatrixXd &fun, const double t, Eigen::MatrixXd &result, Eigen::VectorXd &von_mises)
compute von mises stress at quadrature points for the function fun, also compute the interpolated fun...
static void mark_flipped_cells(const mesh::Mesh &mesh, const std::vector< basis::ElementBases > &gbasis, const std::vector< basis::ElementBases > &basis, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const utils::RefElementSampler &sampler, const int n_points, const Eigen::MatrixXd &fun, Eigen::Vector< bool, -1 > &result, const bool use_sampler, const bool boundary_only)
static void interpolate_function(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const utils::RefElementSampler &sampler, const int n_points, const Eigen::MatrixXd &fun, Eigen::MatrixXd &result, const bool use_sampler, const bool boundary_only)
interpolate the function fun.
static void average_grad_based_function(const mesh::Mesh &mesh, const bool is_problem_scalar, const int n_bases, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const assembler::Assembler &assembler, const utils::RefElementSampler &sampler, const double t, const int n_points, const Eigen::MatrixXd &fun, std::vector< assembler::Assembler::NamedMatrix > &result_scalar, std::vector< assembler::Assembler::NamedMatrix > &result_tensor, const bool use_sampler, const bool boundary_only)
computes scalar quantity of funtion (ie von mises for elasticity and norm of velocity for fluid) the ...
static void compute_scalar_value(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const assembler::Assembler &assembler, const utils::RefElementSampler &sampler, const int n_points, const Eigen::MatrixXd &fun, const double t, std::vector< assembler::Assembler::NamedMatrix > &result, const bool use_sampler, const bool boundary_only)
computes scalar quantity of funtion (ie von mises for elasticity and norm of velocity for fluid)
static void compute_tensor_value(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const assembler::Assembler &assembler, const utils::RefElementSampler &sampler, const int n_points, const Eigen::MatrixXd &fun, const double t, std::vector< assembler::Assembler::NamedMatrix > &result, const bool use_sampler, const bool boundary_only)
compute tensor quantity (ie stress tensor or velocy)
void build_vis_mesh(const mesh::Mesh &mesh, const Eigen::VectorXi &disc_orders, const std::vector< basis::ElementBases > &gbases, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const bool boundary_only, Eigen::MatrixXd &points, Eigen::MatrixXi &tets, Eigen::MatrixXi &el_id, Eigen::MatrixXd &discr) const
builds visualzation mesh, upsampled mesh used for visualization the visualization mesh is a dense mes...
Definition OutData.cpp:645
void build_high_order_vis_mesh(const mesh::Mesh &mesh, const Eigen::VectorXi &disc_orders, const std::vector< basis::ElementBases > &bases, Eigen::MatrixXd &points, std::vector< std::vector< int > > &elements, Eigen::MatrixXi &el_id, Eigen::MatrixXd &discr) const
builds high-der visualzation mesh per element all disconnected it also retuns the mapping to element ...
Definition OutData.cpp:788
void save_points(const std::string &path, const State &state, const Eigen::MatrixXd &sol, const ExportOptions &opts, std::vector< SolutionFrame > &solution_frames) const
saves the nodal values
Definition OutData.cpp:2261
void save_volume_vector_field(const State &state, const Eigen::MatrixXd &points, const ExportOptions &opts, const std::string &name, const Eigen::VectorXd &field, paraviewo::ParaviewWriter &writer) const
Definition OutData.cpp:1751
Eigen::MatrixXd grid_points_bc
grid mesh boundaries
Definition OutData.hpp:278
Eigen::MatrixXd grid_points
grid mesh points to export solution sampled on a grid
Definition OutData.hpp:274
void save_wire(const std::string &name, const State &state, const Eigen::MatrixXd &sol, const double t, const ExportOptions &opts, std::vector< SolutionFrame > &solution_frames) const
saves the wireframe
Definition OutData.cpp:2073
void save_surface(const std::string &export_surface, const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const double t, const double dt_in, const ExportOptions &opts, const bool is_contact_enabled, std::vector< SolutionFrame > &solution_frames) const
saves the surface vtu file for for surface quantites, eg traction forces
Definition OutData.cpp:1780
void build_vis_boundary_mesh(const mesh::Mesh &mesh, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const std::vector< mesh::LocalBoundary > &total_local_boundary, const Eigen::MatrixXd &solution, const int problem_dim, Eigen::MatrixXd &boundary_vis_vertices, Eigen::MatrixXd &boundary_vis_local_vertices, Eigen::MatrixXi &boundary_vis_elements, Eigen::MatrixXi &boundary_vis_elements_ids, Eigen::MatrixXi &boundary_vis_primitive_ids, Eigen::MatrixXd &boundary_vis_normals, Eigen::MatrixXd &displaced_boundary_vis_normals) const
builds the boundary mesh for visualization
Definition OutData.cpp:397
void build_grid(const polyfem::mesh::Mesh &mesh, const double spacing)
builds the grid to export the solution
Definition OutData.cpp:2330
void save_volume(const std::string &path, const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const double t, const double dt, const ExportOptions &opts, std::vector< SolutionFrame > &solution_frames) const
saves the volume vtu file
Definition OutData.cpp:1197
static void extract_boundary_mesh(const mesh::Mesh &mesh, const int n_bases, const std::vector< basis::ElementBases > &bases, const std::vector< mesh::LocalBoundary > &total_local_boundary, Eigen::MatrixXd &node_positions, Eigen::MatrixXi &boundary_edges, Eigen::MatrixXi &boundary_triangles, std::vector< Eigen::Triplet< double > > &displacement_map_entries)
extracts the boundary mesh
Definition OutData.cpp:145
void save_pvd(const std::string &name, const std::function< std::string(int)> &vtu_names, int time_steps, double t0, double dt, int skip_frame=1) const
save a PVD of a time dependent simulation
Definition OutData.cpp:2317
void export_data(const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const bool is_time_dependent, const double tend_in, const double dt, const ExportOptions &opts, const std::string &vis_mesh_path, const std::string &nodes_path, const std::string &solution_path, const std::string &stress_path, const std::string &mises_path, const bool is_contact_enabled, std::vector< SolutionFrame > &solution_frames) const
exports everytihng, txt, vtu, etc
Definition OutData.cpp:940
void save_vtu(const std::string &path, const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const double t, const double dt, const ExportOptions &opts, const bool is_contact_enabled, std::vector< SolutionFrame > &solution_frames) const
saves the vtu file for time t
Definition OutData.cpp:1114
void save_contact_surface(const std::string &export_surface, const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const double t, const double dt_in, const ExportOptions &opts, const bool is_contact_enabled, std::vector< SolutionFrame > &solution_frames) const
saves the surface vtu file for for constact quantites, eg contact or friction forces
Definition OutData.cpp:1979
void init_sampler(const polyfem::mesh::Mesh &mesh, const double vismesh_rel_area)
unitalize the ref element sampler
Definition OutData.cpp:2325
Eigen::MatrixXi grid_points_to_elements
grid mesh mapping to fe elements
Definition OutData.hpp:276
utils::RefElementSampler ref_element_sampler
used to sample the solution
Definition OutData.hpp:271
stores all runtime data
Definition OutData.hpp:365
double loading_mesh_time
time to load the mesh
Definition OutData.hpp:370
double assembling_stiffness_mat_time
time to assembly
Definition OutData.hpp:374
double assigning_rhs_time
time to computing the rhs
Definition OutData.hpp:378
double assembling_mass_mat_time
time to assembly mass
Definition OutData.hpp:376
double building_basis_time
time to construct the basis
Definition OutData.hpp:368
double solving_time
time to solve
Definition OutData.hpp:380
double computing_poly_basis_time
time to build the polygonal/polyhedral bases
Definition OutData.hpp:372
void save_json(const nlohmann::json &args, const int n_bases, const int n_pressure_bases, const Eigen::MatrixXd &sol, const mesh::Mesh &mesh, const Eigen::VectorXi &disc_orders, const assembler::Problem &problem, const OutRuntimeData &runtime, const std::string &formulation, const bool isoparametric, const int sol_at_node_id, nlohmann::json &j)
saves the output statistic to a json object
Definition OutData.cpp:2793
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:2515
void compute_errors(const int n_bases, const std::vector< polyfem::basis::ElementBases > &bases, const std::vector< polyfem::basis::ElementBases > &gbases, const polyfem::mesh::Mesh &mesh, const assembler::Problem &problem, const double tend, const Eigen::MatrixXd &sol)
compute errors
Definition OutData.cpp:2560
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:2428
void reset()
clears all stats
Definition OutData.cpp:2506
void compute_mesh_stats(const polyfem::mesh::Mesh &mesh)
compute stats (counts els type, mesh lenght, etc), step 1 of solve
Definition OutData.cpp:2721
void write(const int t, const double forward, const double remeshing, const double global_relaxation, const Eigen::MatrixXd &sol)
Definition OutData.cpp:2969
RuntimeStatsCSVWriter(const std::string &path, const State &state, const double t0, const double dt)
Definition OutData.cpp:2958
Boundary primitive IDs for a single element.
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
int n_elements() const
utitlity to return the number of elements, cells or faces in 3d and 2d
Definition Mesh.hpp:161
virtual int n_vertices() const =0
number of vertices
virtual int get_body_id(const int primitive) const
Get the volume selection of an element (cell in 3d, face in 2d)
Definition Mesh.hpp:501
virtual double edge_length(const int gid) const
edge length
Definition Mesh.hpp:306
bool is_polytope(const int el_id) const
checks if element is polygon compatible
Definition Mesh.cpp:365
virtual void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1) const =0
Get all the edges.
virtual double tri_area(const int gid) const
area of a tri face of a tet mesh
Definition Mesh.hpp:324
bool is_simplicial() const
checks if the mesh is simplicial
Definition Mesh.hpp:577
virtual bool is_conforming() const =0
if the mesh is conforming
virtual void bounding_box(RowVectorNd &min, RowVectorNd &max) const =0
computes the bbox of the mesh
virtual void barycentric_coords(const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const =0
constructs barycentric coodiantes for a point p.
bool is_cube(const int el_id) const
checks if element is cube compatible
Definition Mesh.cpp:352
const Eigen::MatrixXi & orders() const
order of each element
Definition Mesh.hpp:283
virtual int get_boundary_id(const int primitive) const
Get the boundary selection of an element (face in 3d, edge in 2d)
Definition Mesh.hpp:475
bool is_simplex(const int el_id) const
checks if element is simples compatible
Definition Mesh.cpp:422
virtual double quad_area(const int gid) const
area of a quad face of an hex mesh
Definition Mesh.hpp:315
bool is_linear() const
check if the mesh is linear
Definition Mesh.hpp:591
virtual bool is_volume() const =0
checks if mesh is volume
bool has_poly() const
checks if the mesh has polytopes
Definition Mesh.hpp:563
int dimension() const
utily for dimension
Definition Mesh.hpp:151
virtual int n_faces() const =0
number of faces
const std::vector< ElementType > & elements_tag() const
Returns the elements types.
Definition Mesh.hpp:415
virtual int n_face_vertices(const int f_id) const =0
number of vertices of a face
virtual void elements_boxes(std::vector< std::array< Eigen::Vector3d, 2 > > &boxes) const =0
constructs a box around every element (3d cell, 2d face)
virtual bool is_boundary_element(const int element_global_id) const =0
is cell boundary
virtual int get_node_id(const int node_id) const
Get the boundary selection of a node.
Definition Mesh.hpp:484
const Eigen::MatrixXi & get_edge_connectivity() const
Definition Obstacle.hpp:47
const Eigen::MatrixXi & get_face_connectivity() const
Definition Obstacle.hpp:46
const Eigen::MatrixXd & v() const
Definition Obstacle.hpp:41
const Eigen::VectorXi & get_vertex_connectivity() const
Definition Obstacle.hpp:48
class to store time stepping data
Definition SolveData.hpp:57
std::shared_ptr< solver::FrictionForm > friction_form
std::shared_ptr< solver::NLProblem > nl_problem
std::shared_ptr< solver::ContactForm > contact_form
std::vector< std::pair< std::string, std::shared_ptr< solver::Form > > > named_forms() const
std::shared_ptr< time_integrator::ImplicitTimeIntegrator > time_integrator
static void normal_for_quad_edge(int index, Eigen::MatrixXd &normal)
static void normal_for_tri_edge(int index, Eigen::MatrixXd &normal)
static void normal_for_quad_face(int index, Eigen::MatrixXd &normal)
static void sample_parametric_tri_face(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static void normal_for_tri_face(int index, Eigen::MatrixXd &normal)
static void sample_parametric_quad_face(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static void normal_for_polygon_edge(int face_id, int edge_id, const mesh::Mesh &mesh, Eigen::MatrixXd &normal)
static void sample_parametric_quad_edge(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static void sample_polygon_edge(int face_id, int edge_id, int n_samples, const mesh::Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static bool boundary_quadrature(const mesh::LocalBoundary &local_boundary, const int order, const mesh::Mesh &mesh, const bool skip_computation, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::MatrixXd &normals, Eigen::VectorXd &weights, Eigen::VectorXi &global_primitive_ids)
static void sample_parametric_tri_edge(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static void sample_3d_simplex(const int resolution, Eigen::MatrixXd &samples)
static void sample_3d_cube(const int resolution, Eigen::MatrixXd &samples)
static void sample_2d_cube(const int resolution, Eigen::MatrixXd &samples)
static void sample_2d_simplex(const int resolution, Eigen::MatrixXd &samples)
void init(const bool is_volume, const int n_elements, const double target_rel_area)
const Eigen::MatrixXd & simplex_points() const
const Eigen::MatrixXi & simplex_volume() const
size_t getPeakRSS(void)
Returns the peak (maximum so far) resident set size (physical memory use) measured in bytes,...
Definition getRSS.c:37
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)
ElementType
Type of Element, check [Poly-Spline Finite Element Method] for a complete description.
Definition Mesh.hpp:23
size_t get_n_threads()
Definition par_for.hpp:51
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
void append_rows_of_zeros(DstMat &dst, const size_t n_zero_rows)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
nlohmann::json json
Definition Common.hpp:9
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71
std::string file_extension() const
return the extension of the output paraview files depending on use_hdf5
Definition OutData.hpp:90
ExportOptions(const json &args, const bool is_mesh_linear, const bool is_problem_scalar, const bool solve_export_to_file)
initialize the flags based on the input args
Definition OutData.cpp:1081