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