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