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