PolyFEM
Loading...
Searching...
No Matches
CMesh3D.cpp
Go to the documentation of this file.
5
7
8#include <igl/barycentric_coordinates.h>
9
10#include <geogram/mesh/mesh_io.h>
11#include <fstream>
12
13using namespace polyfem::utils;
14
15namespace polyfem
16{
17 namespace mesh
18 {
19 void CMesh3D::refine(const int n_refinement, const double t)
20 {
21 if (n_refinement <= 0)
22 {
23 return;
24 }
25
26 // TODO refine high order mesh!
27 orders_.resize(0, 0);
29 {
31 }
32 else
33 {
34 for (size_t i = 0; i < elements_tag().size(); ++i)
35 {
37 mesh_.elements[i].hex = false;
38 }
39
40 bool reverse_grow = false;
41 std::vector<int> parent_nodes;
42 MeshProcessing3D::refine_catmul_clark_polar(mesh_, n_refinement, reverse_grow, parent_nodes);
43 }
44
47
48 in_ordered_vertices_ = Eigen::VectorXi::LinSpaced(n_vertices(), 0, n_vertices() - 1);
49 assert(in_ordered_vertices_[0] == 0);
50 assert(in_ordered_vertices_[1] == 1);
51 assert(in_ordered_vertices_[2] == 2);
52 assert(in_ordered_vertices_[in_ordered_vertices_.size() - 1] == n_vertices() - 1);
53
54 in_ordered_edges_.resize(mesh_.edges.size(), 2);
55
56 for (int e = 0; e < (int)mesh_.edges.size(); ++e)
57 {
58 assert(mesh_.edges[e].vs.size() == 2);
59 for (int lv = 0; lv < 2; ++lv)
60 {
61 in_ordered_edges_(e, lv) = mesh_.edges[e].vs[lv];
62 }
63 }
64 assert(in_ordered_edges_.size() > 0);
65
66 in_ordered_faces_.resize(mesh_.faces.size(), mesh_.faces[0].vs.size());
67
68 for (int f = 0; f < (int)mesh_.faces.size(); ++f)
69 {
70 assert(in_ordered_faces_.cols() == mesh_.faces[f].vs.size());
71
72 for (int lv = 0; lv < in_ordered_faces_.cols(); ++lv)
73 {
74 in_ordered_faces_(f, lv) = mesh_.faces[f].vs[lv];
75 }
76 }
77 assert(in_ordered_faces_.size() > 0);
78 }
79
80 bool CMesh3D::load(const std::string &path)
81 {
82 edge_nodes_.clear();
83 face_nodes_.clear();
84 cell_nodes_.clear();
85
86 if (!StringUtils::endswith(path, ".HYBRID"))
87 {
88 GEO::Mesh M;
89 GEO::mesh_load(path, M);
90 return load(M);
91 }
92
93 FILE *f = fopen(path.data(), "rt");
94 if (!f)
95 return false;
96
97 int nv, np, nh;
98 auto _ = fscanf(f, "%d %d %d", &nv, &np, &nh);
99 nh /= 3;
100
101 mesh_.points.resize(3, nv);
102 mesh_.vertices.resize(nv);
103
104 for (int i = 0; i < nv; i++)
105 {
106 double x, y, z;
107 auto _ = fscanf(f, "%lf %lf %lf", &x, &y, &z);
108 mesh_.points(0, i) = x;
109 mesh_.points(1, i) = y;
110 mesh_.points(2, i) = z;
111 Vertex v;
112 v.id = i;
113 mesh_.vertices[i] = v;
114 }
115 mesh_.faces.resize(np);
116 for (int i = 0; i < np; i++)
117 {
118 Face &p = mesh_.faces[i];
119 p.id = i;
120 int nw;
121
122 auto _ = fscanf(f, "%d", &nw);
123 p.vs.resize(nw);
124 for (int j = 0; j < nw; j++)
125 {
126 auto _ = fscanf(f, "%u", &(p.vs[j]));
127 }
128 }
129 mesh_.elements.resize(nh);
130 for (int i = 0; i < nh; i++)
131 {
132 Element &h = mesh_.elements[i];
133 h.id = i;
134
135 int nf;
136 auto _ = fscanf(f, "%d", &nf);
137 h.fs.resize(nf);
138
139 for (int j = 0; j < nf; j++)
140 {
141 auto _ = fscanf(f, "%u", &(h.fs[j]));
142 }
143
144 for (auto fid : h.fs)
145 h.vs.insert(h.vs.end(), mesh_.faces[fid].vs.begin(), mesh_.faces[fid].vs.end());
146 sort(h.vs.begin(), h.vs.end());
147 h.vs.erase(unique(h.vs.begin(), h.vs.end()), h.vs.end());
148
149 int tmp;
150 auto __ = fscanf(f, "%d", &tmp);
151 for (int j = 0; j < nf; j++)
152 {
153 int s;
154 auto _ = fscanf(f, "%d", &s);
155 h.fs_flag.push_back(s);
156 }
157 }
158 for (int i = 0; i < nh; i++)
159 {
160 int tmp;
161 auto _ = fscanf(f, "%d", &tmp);
162 mesh_.elements[i].hex = tmp;
163 }
164
165 char s[1024], sread[1024];
166 int find = false, num = 0;
167 while (!feof(f) && !find)
168 {
169 auto _ = fgets(s, 1024, f);
170 if (sscanf(s, "%s%d", sread, &num) == 2 && (strcmp(sread, "KERNEL") == 0))
171 find = true;
172 }
173 if (find)
174 {
175 for (int i = 0; i < nh; i++)
176 {
177 double x, y, z;
178 auto _ = fscanf(f, "%lf %lf %lf", &x, &y, &z);
179 mesh_.elements[i].v_in_Kernel.push_back(x);
180 mesh_.elements[i].v_in_Kernel.push_back(y);
181 mesh_.elements[i].v_in_Kernel.push_back(z);
182 }
183 }
184
185 fclose(f);
186
187 // remove horrible kernels and replace with barycenters
188 for (int c = 0; c < n_cells(); ++c)
189 {
190 auto bary = cell_barycenter(c);
191 for (int d = 0; d < 3; ++d)
192 mesh_.elements[c].v_in_Kernel[d] = bary(d);
193 }
194
196 // if(is_simplicial())
197 // MeshProcessing3D::orient_volume_mesh(mesh_);
199 return true;
200 }
201
202 // load from a geogram surface mesh (for debugging), or volume mesh
203 // if loading a surface mesh, it assumes there is only one polyhedral cell, and the last vertex id a point in the kernel
204 bool CMesh3D::load(const GEO::Mesh &M)
205 {
206 edge_nodes_.clear();
207 face_nodes_.clear();
208 cell_nodes_.clear();
209
210 assert(M.vertices.dimension() == 3);
211
212 // Set vertices
213 const int nv = M.vertices.nb();
214 mesh_.points.resize(3, nv);
215 mesh_.vertices.resize(nv);
216 for (int i = 0; i < nv; ++i)
217 {
218 mesh_.points(0, i) = M.vertices.point(i)[0];
219 mesh_.points(1, i) = M.vertices.point(i)[1];
220 mesh_.points(2, i) = M.vertices.point(i)[2];
221 Vertex v;
222 v.id = i;
223 mesh_.vertices[i] = v;
224 }
225
226 // Set cells
227 if (M.cells.nb() == 0)
228 {
229
230 bool last_isolated = true;
231
232 // Set faces
233 mesh_.faces.resize(M.facets.nb());
234 for (int i = 0; i < (int)M.facets.nb(); ++i)
235 {
236 Face &face = mesh_.faces[i];
237 face.id = i;
238
239 face.vs.resize(M.facets.nb_vertices(i));
240 for (int j = 0; j < (int)M.facets.nb_vertices(i); ++j)
241 {
242 face.vs[j] = M.facets.vertex(i, j);
243 if ((int)face.vs[j] == nv - 1)
244 {
245 last_isolated = false;
246 }
247 }
248 }
249
250 // Assumes there is only 1 polyhedron described by a closed input surface
251 mesh_.elements.resize(1);
252 for (int i = 0; i < 1; ++i)
253 {
254 Element &cell = mesh_.elements[i];
255 cell.id = i;
256
257 int nf = M.facets.nb();
258 cell.fs.resize(nf);
259
260 for (int j = 0; j < nf; ++j)
261 {
262 cell.fs[j] = j;
263 }
264
265 for (auto fid : cell.fs)
266 {
267 cell.vs.insert(cell.vs.end(), mesh_.faces[fid].vs.begin(), mesh_.faces[fid].vs.end());
268 }
269 sort(cell.vs.begin(), cell.vs.end());
270 cell.vs.erase(unique(cell.vs.begin(), cell.vs.end()), cell.vs.end());
271
272 for (int j = 0; j < nf; ++j)
273 {
274 cell.fs_flag.push_back(1);
275 }
276
277 if (last_isolated)
278 {
279 cell.v_in_Kernel.push_back(M.vertices.point(nv - 1)[0]);
280 cell.v_in_Kernel.push_back(M.vertices.point(nv - 1)[1]);
281 cell.v_in_Kernel.push_back(M.vertices.point(nv - 1)[2]);
282 }
283 else
284 {
285 // Compute a point in the kernel (assumes the barycenter is ok)
286 Eigen::RowVector3d p(0, 0, 0);
287 for (int v : cell.vs)
288 {
289 p += mesh_.points.col(v).transpose();
290 }
291 p /= cell.vs.size();
292 cell.v_in_Kernel.push_back(p[0]);
293 cell.v_in_Kernel.push_back(p[1]);
294 cell.v_in_Kernel.push_back(p[2]);
295 }
296 }
297
298 for (int i = 0; i < 1; ++i)
299 {
300 mesh_.elements[i].hex = false;
301 } // FIME me here!
302
303 mesh_.type = M.cells.are_simplices() ? MeshType::TET : MeshType::HYB;
304 }
305 else
306 {
307
308 // Set faces
309 mesh_.faces.clear();
310 // for (int f = 0; f < (int) M.facets.nb(); ++f) {
311 // Face &face = mesh_.faces[f];
312 // face.id = -1;
313 // face.vs.clear();
314 // // face.id = f;
315
316 // // face.vs.resize(M.facets.nb_vertices(f));
317 // // for (int lv = 0; lv < (int) M.facets.nb_vertices(f); ++lv) {
318 // // face.vs[lv] = M.facets.vertex(f, lv);
319 // // }
320 // }
321
322 auto opposite_cell_facet = [&M](int c, int cf) {
323 GEO::index_t c2 = M.cell_facets.adjacent_cell(cf);
324 if (c2 == GEO::NO_FACET)
325 {
326 return -1;
327 }
328 for (int lf = 0; lf < (int)M.cells.nb_facets(c2); ++lf)
329 {
330 if (c == (int)M.cells.adjacent(c2, lf))
331 {
332 return (int)M.cells.facet(c2, lf);
333 }
334 }
335 assert(false);
336 return -1;
337 };
338
339 std::vector<int> cell_facet_to_facet(M.cell_facets.nb(), -1);
340
341 // Creates 1 hex or polyhedral element for each cell of the input mesh
342 int facet_counter = 0;
343 mesh_.elements.resize(M.cells.nb());
344 bool is_hex = true;
345 for (int c = 0; c < (int)M.cells.nb(); ++c)
346 {
347 Element &cell = mesh_.elements[c];
348 cell.id = c;
349 cell.hex = (M.cells.type(c) == GEO::MESH_HEX);
350
351 is_hex = is_hex && cell.hex;
352
353 int nf = M.cells.nb_facets(c);
354 cell.fs.resize(nf);
355
356 for (int lf = 0; lf < nf; ++lf)
357 {
358 int cf = M.cells.facet(c, lf);
359 int cf2 = opposite_cell_facet(c, cf);
360 if (cf2 < 0 || cell_facet_to_facet[cf2] < 0)
361 {
362 mesh_.faces.emplace_back();
363 Face &face = mesh_.faces[facet_counter];
364 assert(face.vs.empty());
365 face.vs.resize(M.cells.facet_nb_vertices(c, lf));
366 for (int lv = 0; lv < (int)M.cells.facet_nb_vertices(c, lf); ++lv)
367 {
368 face.vs[lv] = M.cells.facet_vertex(c, lf, lv);
369 }
370 cell.fs_flag.push_back(0);
371 cell.fs[lf] = face.id = facet_counter;
372 cell_facet_to_facet[cf] = facet_counter;
373 ++facet_counter;
374 }
375 else
376 {
377 cell.fs[lf] = cell_facet_to_facet[cf2];
378 cell.fs_flag.push_back(1);
379 }
380 }
381
382 for (auto fid : cell.fs)
383 {
384 cell.vs.insert(cell.vs.end(), mesh_.faces[fid].vs.begin(), mesh_.faces[fid].vs.end());
385 }
386 sort(cell.vs.begin(), cell.vs.end());
387 cell.vs.erase(unique(cell.vs.begin(), cell.vs.end()), cell.vs.end());
388
389 // Compute a point in the kernel (assumes the barycenter is ok)
390 Eigen::RowVector3d p(0, 0, 0);
391 for (int v : cell.vs)
392 {
393 p += mesh_.points.col(v).transpose();
394 }
395 p /= cell.vs.size();
396 mesh_.elements[c].v_in_Kernel.push_back(p[0]);
397 mesh_.elements[c].v_in_Kernel.push_back(p[1]);
398 mesh_.elements[c].v_in_Kernel.push_back(p[2]);
399 }
400 mesh_.type = is_hex ? MeshType::HEX : (M.cells.are_simplices() ? MeshType::TET : MeshType::HYB);
401 }
402
404 // if (is_simplicial()) {
405 // MeshProcessing3D::orient_volume_mesh(mesh_);
406 // }
408 return true;
409 }
410
411 bool CMesh3D::save(const std::string &path) const
412 {
413
414 if (!StringUtils::endswith(path, ".HYBRID"))
415 {
416 GEO::Mesh M;
417 to_geogram_mesh(*this, M);
418 GEO::mesh_save(M, path);
419 return true;
420 }
421
422 std::fstream f(path, std::ios::out);
423
424 f << mesh_.points.cols() << " " << mesh_.faces.size() << " " << 3 * mesh_.elements.size() << std::endl;
425 for (int i = 0; i < mesh_.points.cols(); i++)
426 f << mesh_.points(0, i) << " " << mesh_.points(1, i) << " " << mesh_.points(2, i) << std::endl;
427
428 for (auto f_ : mesh_.faces)
429 {
430 f << f_.vs.size() << " ";
431 for (auto vid : f_.vs)
432 f << vid << " ";
433 f << std::endl;
434 }
435
436 for (uint32_t i = 0; i < mesh_.elements.size(); i++)
437 {
438 f << mesh_.elements[i].fs.size() << " ";
439 for (auto fid : mesh_.elements[i].fs)
440 f << fid << " ";
441 f << std::endl;
442 f << mesh_.elements[i].fs_flag.size() << " ";
443 for (auto f_flag : mesh_.elements[i].fs_flag)
444 f << f_flag << " ";
445 f << std::endl;
446 }
447
448 for (uint32_t i = 0; i < mesh_.elements.size(); i++)
449 {
450 f << mesh_.elements[i].hex << std::endl;
451 }
452
453 f << "KERNEL"
454 << " " << mesh_.elements.size() << std::endl;
455 for (uint32_t i = 0; i < mesh_.elements.size(); i++)
456 {
457 f << mesh_.elements[i].v_in_Kernel[0] << " " << mesh_.elements[i].v_in_Kernel[1] << " " << mesh_.elements[i].v_in_Kernel[2] << std::endl;
458 }
459 f.close();
460
461 return true;
462 }
463
464 bool CMesh3D::build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
465 {
466 assert(F.cols() == 4 || F.cols() == 5 || F.cols() == 6 || F.cols() == 8);
467 edge_nodes_.clear();
468 face_nodes_.clear();
469 cell_nodes_.clear();
470
471 GEO::Mesh M;
472 M.vertices.create_vertices((int)V.rows());
473 for (int i = 0; i < (int)M.vertices.nb(); ++i)
474 {
475 GEO::vec3 &p = M.vertices.point(i);
476 p[0] = V(i, 0);
477 p[1] = V(i, 1);
478 p[2] = V(i, 2);
479 }
480
481 static const std::vector<int> permute_tet = {0, 1, 2, 3};
482 static const std::vector<int> permute_pyramid = {0, 1, 2, 3, 4};
483 static const std::vector<int> permute_prism = {0, 1, 2, 3, 4, 5};
484 // polyfem uses the msh file format for hexes ordering!
485 static const std::vector<int> permute_hex = {1, 0, 2, 3, 5, 4, 6, 7};
486
487 auto add_cell = [&](int c, GEO::MeshCellType t, int nv, const std::vector<int> &perm) {
488 GEO::index_t cid = M.cells.create_cells(1, t);
489 for (int lv = 0; lv < nv; ++lv)
490 {
491 int vi = F(c, perm[lv]);
492 assert(vi >= 0 && vi < V.rows());
493 M.cells.set_vertex(cid, lv, GEO::index_t(vi));
494 }
495 };
496
497 for (int c = 0; c < F.rows(); ++c)
498 {
499 int nV;
500 for (nV = 0; nV < F.cols(); ++nV)
501 {
502 if (F(c, nV) == -1)
503 break;
504 }
505
506#ifndef NDEBUG
507 for (int k = nV; k < F.cols(); ++k)
508 {
509 assert(F(c, k) == -1);
510 }
511#endif
512 if (nV == 4)
513 add_cell(c, GEO::MESH_TET, 4, permute_tet);
514 else if (nV == 5)
515 add_cell(c, GEO::MESH_PYRAMID, 5, permute_pyramid);
516 else if (nV == 6)
517 add_cell(c, GEO::MESH_PRISM, 6, permute_prism);
518 else if (nV == 8)
519 add_cell(c, GEO::MESH_HEX, 8, permute_hex);
520 else
521 log_and_throw_error(fmt::format("Invalid number of vertices {}", nV));
522 }
523 M.cells.connect();
524
525 return load(M);
526 }
527
528 void CMesh3D::attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector<std::vector<int>> &nodes)
529 {
530 edge_nodes_.clear();
531 face_nodes_.clear();
532 cell_nodes_.clear();
533
534 edge_nodes_.resize(n_edges());
535 face_nodes_.resize(n_faces());
536 cell_nodes_.resize(n_cells());
537
538 orders_.resize(n_cells(), 1);
539
540 const auto attach_p2 = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
541 auto &n = edge_nodes_[index.edge];
542
543 if (n.nodes.size() > 0)
544 return;
545
546 n.v1 = index.vertex;
547 n.v2 = switch_vertex(index).vertex;
548
549 const int n_v1 = index.vertex;
550 const int n_v2 = switch_vertex(index).vertex;
551
552 int node_index = 0;
553
554 if ((n_v1 == nodes_ids[0] && n_v2 == nodes_ids[1]) || (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[1]))
555 node_index = 4;
556 else if ((n_v1 == nodes_ids[1] && n_v2 == nodes_ids[2]) || (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[2]))
557 node_index = 5;
558 else if ((n_v1 == nodes_ids[2] && n_v2 == nodes_ids[3]) || (n_v2 == nodes_ids[2] && n_v1 == nodes_ids[3]))
559 node_index = 8;
560
561 else if ((n_v1 == nodes_ids[0] && n_v2 == nodes_ids[3]) || (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[3]))
562 node_index = 7;
563 else if ((n_v1 == nodes_ids[0] && n_v2 == nodes_ids[2]) || (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[2]))
564 node_index = 6;
565 else
566 node_index = 9;
567
568 n.nodes.resize(1, 3);
569 n.nodes << V(nodes_ids[node_index], 0), V(nodes_ids[node_index], 1), V(nodes_ids[node_index], 2);
570 n.nodes_ids.push_back(nodes_ids[node_index]);
571 };
572
573 const auto attach_p3 = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
574 auto &n = edge_nodes_[index.edge];
575
576 if (n.nodes.size() > 0)
577 return;
578
579 n.v1 = index.vertex;
580 n.v2 = switch_vertex(index).vertex;
581
582 const int n_v1 = index.vertex;
583 const int n_v2 = switch_vertex(index).vertex;
584
585 int node_index1 = 0;
586 int node_index2 = 0;
587 if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[1])
588 {
589 node_index1 = 4;
590 node_index2 = 5;
591 }
592 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[1])
593 {
594 node_index1 = 5;
595 node_index2 = 4;
596 }
597 else if (n_v1 == nodes_ids[1] && n_v2 == nodes_ids[2])
598 {
599 node_index1 = 6;
600 node_index2 = 7;
601 }
602 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[2])
603 {
604 node_index1 = 7;
605 node_index2 = 6;
606 }
607 else if (n_v1 == nodes_ids[2] && n_v2 == nodes_ids[3])
608 {
609 node_index1 = 13;
610 node_index2 = 12;
611 }
612 else if (n_v2 == nodes_ids[2] && n_v1 == nodes_ids[3])
613 {
614 node_index1 = 12;
615 node_index2 = 13;
616 }
617
618 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[3])
619 {
620 node_index1 = 11;
621 node_index2 = 10;
622 }
623 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[3])
624 {
625 node_index1 = 10;
626 node_index2 = 11;
627 }
628 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[2])
629 {
630 node_index1 = 9;
631 node_index2 = 8;
632 }
633 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[2])
634 {
635 node_index1 = 8;
636 node_index2 = 9;
637 }
638
639 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[3])
640 {
641 node_index1 = 14;
642 node_index2 = 15;
643 }
644 else
645 {
646 node_index1 = 15;
647 node_index2 = 14;
648 }
649
650 n.nodes.resize(2, 3);
651 n.nodes.row(0) << V(nodes_ids[node_index1], 0), V(nodes_ids[node_index1], 1), V(nodes_ids[node_index1], 2);
652 n.nodes.row(1) << V(nodes_ids[node_index2], 0), V(nodes_ids[node_index2], 1), V(nodes_ids[node_index2], 2);
653 n.nodes_ids.push_back(nodes_ids[node_index1]);
654 n.nodes_ids.push_back(nodes_ids[node_index2]);
655 };
656
657 const auto attach_p3_face = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids, int id) {
658 auto &n = face_nodes_[index.face];
659 if (n.nodes.size() <= 0)
660 {
661 n.v1 = face_vertex(index.face, 0);
662 n.v2 = face_vertex(index.face, 1);
663 n.v3 = face_vertex(index.face, 2);
664 n.nodes.resize(1, 3);
665 n.nodes << V(nodes_ids[id], 0), V(nodes_ids[id], 1), V(nodes_ids[id], 2);
666 n.nodes_ids.push_back(nodes_ids[id]);
667 }
668 };
669
670 const auto attach_p4 = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
671 auto &n = edge_nodes_[index.edge];
672
673 if (n.nodes.size() > 0)
674 return;
675
676 n.v1 = index.vertex;
677 n.v2 = switch_vertex(index).vertex;
678
679 const int n_v1 = index.vertex;
680 const int n_v2 = switch_vertex(index).vertex;
681
682 int node_index1 = 0;
683 int node_index2 = 0;
684 int node_index3 = 0;
685 if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[1])
686 {
687 node_index1 = 4;
688 node_index2 = 5;
689 node_index3 = 6;
690 }
691 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[1])
692 {
693 node_index1 = 6;
694 node_index2 = 5;
695 node_index3 = 4;
696 }
697
698 else if (n_v1 == nodes_ids[1] && n_v2 == nodes_ids[2])
699 {
700 node_index1 = 7;
701 node_index2 = 8;
702 node_index3 = 9;
703 }
704 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[2])
705 {
706 node_index1 = 9;
707 node_index2 = 8;
708 node_index3 = 7;
709 }
710
711 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[2])
712 {
713 node_index1 = 10;
714 node_index2 = 11;
715 node_index3 = 12;
716 }
717 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[2])
718 {
719 node_index1 = 12;
720 node_index2 = 11;
721 node_index3 = 10;
722 }
723
724 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[3])
725 {
726 node_index1 = 13;
727 node_index2 = 14;
728 node_index3 = 15;
729 }
730 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[3])
731 {
732 node_index1 = 15;
733 node_index2 = 14;
734 node_index3 = 13;
735 }
736
737 else if (n_v2 == nodes_ids[2] && n_v1 == nodes_ids[3])
738 {
739 node_index1 = 16;
740 node_index2 = 17;
741 node_index3 = 18;
742 }
743 else if (n_v1 == nodes_ids[2] && n_v2 == nodes_ids[3])
744 {
745 node_index1 = 18;
746 node_index2 = 17;
747 node_index3 = 16;
748 }
749
750 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[3])
751 {
752 node_index1 = 19;
753 node_index2 = 20;
754 node_index3 = 21;
755 }
756 else if (n_v1 == nodes_ids[1] && n_v2 == nodes_ids[3])
757 {
758 node_index1 = 21;
759 node_index2 = 20;
760 node_index3 = 19;
761 }
762 else
763 {
764 assert(false);
765 }
766
767 n.nodes.resize(3, 3);
768 n.nodes.row(0) << V(nodes_ids[node_index1], 0), V(nodes_ids[node_index1], 1), V(nodes_ids[node_index1], 2);
769 n.nodes.row(1) << V(nodes_ids[node_index2], 0), V(nodes_ids[node_index2], 1), V(nodes_ids[node_index2], 2);
770 n.nodes.row(2) << V(nodes_ids[node_index3], 0), V(nodes_ids[node_index3], 1), V(nodes_ids[node_index3], 2);
771 n.nodes_ids.push_back(nodes_ids[node_index1]);
772 n.nodes_ids.push_back(nodes_ids[node_index2]);
773 n.nodes_ids.push_back(nodes_ids[node_index3]);
774 };
775
776 const auto attach_p4_face = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
777 auto &n = face_nodes_[index.face];
778 if (n.nodes.size() <= 0)
779 {
780 n.v1 = face_vertex(index.face, 0);
781 n.v2 = face_vertex(index.face, 1);
782 n.v3 = face_vertex(index.face, 2);
783
784 std::array<int, 3> vid = {{n.v1, n.v2, n.v3}};
785 std::sort(vid.begin(), vid.end());
786
787 std::array<int, 3> c1 = {{nodes_ids[0], nodes_ids[1], nodes_ids[2]}}; // 22
788 std::array<int, 3> c2 = {{nodes_ids[0], nodes_ids[1], nodes_ids[3]}}; // 25
789 std::array<int, 3> c3 = {{nodes_ids[0], nodes_ids[2], nodes_ids[3]}}; // 28
790 std::array<int, 3> c4 = {{nodes_ids[1], nodes_ids[2], nodes_ids[3]}}; // 31
791
792 std::sort(c1.begin(), c1.end());
793 std::sort(c2.begin(), c2.end());
794 std::sort(c3.begin(), c3.end());
795 std::sort(c4.begin(), c4.end());
796
797 int id = 0;
798 int index0 = 0;
799 int index1 = 1;
800 int index2 = 2;
801 if (vid == c1)
802 {
803 id = 22;
804 n.v1 = nodes_ids[0];
805 n.v2 = nodes_ids[1];
806 n.v3 = nodes_ids[2];
807
808 index0 = 0;
809 index1 = 2;
810 index2 = 1; // ok
811 }
812 else if (vid == c2)
813 {
814 id = 25;
815 n.v1 = nodes_ids[0];
816 n.v2 = nodes_ids[1];
817 n.v3 = nodes_ids[3];
818
819 index0 = 0;
820 index1 = 1;
821 index2 = 2; // ok
822 }
823 else if (vid == c3)
824 {
825 id = 28;
826 n.v1 = nodes_ids[0];
827 n.v2 = nodes_ids[2];
828 n.v3 = nodes_ids[3];
829
830 index0 = 0;
831 index1 = 2;
832 index2 = 1; // ok
833 }
834 else if (vid == c4)
835 {
836 id = 31;
837 n.v1 = nodes_ids[1];
838 n.v2 = nodes_ids[2];
839 n.v3 = nodes_ids[3];
840
841 index0 = 1;
842 index1 = 2;
843 index2 = 0; // ok
844 }
845 else
846 {
847 // the face nees to be one of the 4 above
848 assert(false);
849 }
850
851 n.nodes.resize(3, 3);
852 assert(id + index0 < nodes_ids.size());
853 assert(id + index1 < nodes_ids.size());
854 assert(id + index2 < nodes_ids.size());
855 n.nodes.row(0) << V(nodes_ids[id + index0], 0), V(nodes_ids[id + index0], 1), V(nodes_ids[id + index0], 2);
856 n.nodes.row(1) << V(nodes_ids[id + index1], 0), V(nodes_ids[id + index1], 1), V(nodes_ids[id + index1], 2);
857 n.nodes.row(2) << V(nodes_ids[id + index2], 0), V(nodes_ids[id + index2], 1), V(nodes_ids[id + index2], 2);
858 n.nodes_ids.push_back(nodes_ids[id + index0]);
859 n.nodes_ids.push_back(nodes_ids[id + index1]);
860 n.nodes_ids.push_back(nodes_ids[id + index2]);
861 }
862 };
863
864 const auto attach_p4_cell = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
865 auto &n = cell_nodes_[index.element];
866 assert(nodes_ids.size() == 35);
867 assert(n.nodes.size() == 0);
868
869 if (n.nodes.size() <= 0)
870 {
871 n.v1 = cell_vertex(index.element, 0);
872 n.v2 = cell_vertex(index.element, 1);
873 n.v3 = cell_vertex(index.element, 2);
874 n.v4 = cell_vertex(index.element, 3);
875 n.nodes.resize(1, 3);
876
877 n.nodes << V(nodes_ids[34], 0), V(nodes_ids[34], 1), V(nodes_ids[34], 2);
878 n.nodes_ids.push_back(nodes_ids[34]);
879 }
880 };
881
882 assert(nodes.size() == n_cells());
883
884 for (int c = 0; c < n_cells(); ++c)
885 {
886 auto index = get_index_from_element(c);
887
888 const auto &nodes_ids = nodes[c];
889
890 if (nodes_ids.size() == 4)
891 {
892 orders_(c) = 1;
893 continue;
894 }
895 // P2
896 else if (nodes_ids.size() == 10)
897 {
898 orders_(c) = 2;
899
900 for (int le = 0; le < 3; ++le)
901 {
902 attach_p2(index, nodes_ids);
903 index = next_around_face(index);
904 }
905
906 index = switch_vertex(switch_edge(switch_face(index)));
907 attach_p2(index, nodes_ids);
908
909 index = switch_edge(index);
910 attach_p2(index, nodes_ids);
911
912 index = switch_edge(switch_face(index));
913 attach_p2(index, nodes_ids);
914 }
915 // P3
916 else if (nodes_ids.size() == 20)
917 {
918 orders_(c) = 3;
919
920 for (int le = 0; le < 3; ++le)
921 {
922 attach_p3(index, nodes_ids);
923 index = next_around_face(index);
924 }
925
926 {
927 index = switch_vertex(switch_edge(switch_face(index)));
928 attach_p3(index, nodes_ids);
929
930 index = switch_edge(index);
931 attach_p3(index, nodes_ids);
932
933 index = switch_edge(switch_face(index));
934 attach_p3(index, nodes_ids);
935 }
936
937 {
938 index = get_index_from_element(c);
939 std::array<int, 4> indices;
940
941 {
942 std::array<int, 3> f16 = {{nodes_ids[0], nodes_ids[1], nodes_ids[2]}};
943 std::array<int, 3> f17 = {{nodes_ids[3], nodes_ids[1], nodes_ids[0]}};
944 std::array<int, 3> f18 = {{nodes_ids[0], nodes_ids[2], nodes_ids[3]}};
945 std::array<int, 3> f19 = {{nodes_ids[1], nodes_ids[2], nodes_ids[3]}};
946 std::sort(f16.begin(), f16.end());
947 std::sort(f17.begin(), f17.end());
948 std::sort(f18.begin(), f18.end());
949 std::sort(f19.begin(), f19.end());
950
951 auto tmp = index;
952 std::array<int, 3> f0 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
953 tmp = switch_face(index);
954 std::array<int, 3> f1 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
955 tmp = switch_face(next_around_face(index));
956 std::array<int, 3> f2 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
958 std::array<int, 3> f3 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
959
960 std::sort(f0.begin(), f0.end());
961 std::sort(f1.begin(), f1.end());
962 std::sort(f2.begin(), f2.end());
963 std::sort(f3.begin(), f3.end());
964
965 const std::array<std::array<int, 3>, 4> faces = {{f0, f1, f2, f3}};
966 const std::array<std::array<int, 3>, 4> nodes = {{f16, f17, f18, f19}};
967 for (int i = 0; i < 4; ++i)
968 {
969 const auto &f = faces[i];
970 bool found = false;
971 for (int j = 0; j < 4; ++j)
972 {
973 if (nodes[j] == f)
974 {
975 indices[i] = j + 16;
976 found = true;
977 break;
978 }
979 }
980
981 assert(found);
982 }
983 }
984
985 attach_p3_face(index, nodes_ids, indices[0]);
986 attach_p3_face(switch_face(index), nodes_ids, indices[1]);
987 attach_p3_face(switch_face(next_around_face(index)), nodes_ids, indices[2]);
988 attach_p3_face(switch_face(next_around_face(next_around_face(index))), nodes_ids, indices[3]);
989 }
990 }
991 // P4
992 else if (nodes_ids.size() == 35)
993 {
994 orders_(c) = 4;
995 for (int le = 0; le < 3; ++le)
996 {
997 attach_p4(index, nodes_ids);
998 index = next_around_face(index);
999 }
1000
1001 {
1002 index = switch_vertex(switch_edge(switch_face(index)));
1003 attach_p4(index, nodes_ids);
1004
1005 index = switch_edge(index);
1006 attach_p4(index, nodes_ids);
1007
1008 index = switch_edge(switch_face(index));
1009 attach_p4(index, nodes_ids);
1010 }
1011
1012 {
1013 index = get_index_from_element(c);
1014
1015 attach_p4_face(index, nodes_ids);
1016 attach_p4_face(switch_face(index), nodes_ids);
1017 attach_p4_face(switch_face(next_around_face(index)), nodes_ids);
1018 attach_p4_face(switch_face(next_around_face(next_around_face(index))), nodes_ids);
1019 }
1020
1021 attach_p4_cell(get_index_from_element(c), nodes_ids);
1022 }
1023 // unsupported
1024 else
1025 {
1026 assert(false);
1027 }
1028 }
1029 }
1030
1032 {
1033 const auto &V = mesh_.points;
1034 min = V.rowwise().minCoeff().transpose();
1035 max = V.rowwise().maxCoeff().transpose();
1036 }
1037
1039 {
1040 RowVectorNd minV, maxV;
1041 bounding_box(minV, maxV);
1042 auto &V = mesh_.points;
1043 const auto shift = V.rowwise().minCoeff().eval();
1044 const double scaling = 1.0 / (V.rowwise().maxCoeff() - V.rowwise().minCoeff()).maxCoeff();
1045 V = (V.colwise() - shift) * scaling;
1046
1047 for (int i = 0; i < n_cells(); ++i)
1048 {
1049 for (int d = 0; d < 3; ++d)
1050 {
1051 auto val = mesh_.elements[i].v_in_Kernel[d];
1052 mesh_.elements[i].v_in_Kernel[d] = (val - shift(d)) * scaling;
1053 }
1054 }
1055
1056 for (auto &n : edge_nodes_)
1057 {
1058 if (n.nodes.size() > 0)
1059 n.nodes = (n.nodes.rowwise() - shift.transpose()) * scaling;
1060 }
1061 for (auto &n : face_nodes_)
1062 {
1063 if (n.nodes.size() > 0)
1064 n.nodes = (n.nodes.rowwise() - shift.transpose()) * scaling;
1065 }
1066 for (auto &n : cell_nodes_)
1067 {
1068 if (n.nodes.size() > 0)
1069 n.nodes = (n.nodes.rowwise() - shift.transpose()) * scaling;
1070 }
1071
1072 logger().debug("-- bbox before normalization:");
1073 logger().debug(" min : {}", minV);
1074 logger().debug(" max : {}", maxV);
1075 logger().debug(" extent: {}", maxV - minV);
1076 bounding_box(minV, maxV);
1077 logger().debug("-- bbox after normalization:");
1078 logger().debug(" min : {}", minV);
1079 logger().debug(" max : {}", maxV);
1080 logger().debug(" extent: {}", maxV - minV);
1081
1082 // V.row(1) /= 100.;
1083
1084 // for(int i = 0; i < n_cells(); ++i)
1085 // {
1086 // mesh_.elements[i].v_in_Kernel[1] /= 100.;
1087 // }
1088
1089 Eigen::MatrixXd p0, p1, p;
1090 get_edges(p0, p1);
1091 p = p0 - p1;
1092 logger().debug("-- edge length after normalization:");
1093 logger().debug(" min: ", p.rowwise().norm().minCoeff());
1094 logger().debug(" max: ", p.rowwise().norm().maxCoeff());
1095 logger().debug(" avg: ", p.rowwise().norm().mean());
1096 }
1097
1098 // void CMesh3D::triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector<int> &ranges) const
1099 // {
1100 // ranges.clear();
1101
1102 // std::vector<Eigen::MatrixXi> local_tris(mesh_.elements.size());
1103 // std::vector<Eigen::MatrixXd> local_pts(mesh_.elements.size());
1104 // Eigen::MatrixXi tets;
1105
1106 // int total_tris = 0;
1107 // int total_pts = 0;
1108
1109 // ranges.push_back(0);
1110
1111 // Eigen::MatrixXd face_barys;
1112 // face_barycenters(face_barys);
1113
1114 // Eigen::MatrixXd cell_barys;
1115 // cell_barycenters(cell_barys);
1116
1117 // for (std::size_t e = 0; e < mesh_.elements.size(); ++e)
1118 // {
1119 // const Element &el = mesh_.elements[e];
1120
1121 // const int n_vertices = el.vs.size();
1122 // const int n_faces = el.fs.size();
1123
1124 // Eigen::MatrixXd local_pt(n_vertices + n_faces, 3);
1125
1126 // std::map<int, int> global_to_local;
1127
1128 // for (int i = 0; i < n_vertices; ++i)
1129 // {
1130 // const int global_index = el.vs[i];
1131 // local_pt.row(i) = mesh_.points.col(global_index).transpose();
1132 // global_to_local[global_index] = i;
1133 // }
1134
1135 // int n_local_faces = 0;
1136 // for (int i = 0; i < n_faces; ++i)
1137 // {
1138 // const Face &f = mesh_.faces[el.fs[i]];
1139 // n_local_faces += f.vs.size();
1140
1141 // local_pt.row(n_vertices + i) = face_barys.row(f.id); // node_from_face(f.id);
1142 // }
1143
1144 // Eigen::MatrixXi local_faces(n_local_faces, 3);
1145
1146 // int face_index = 0;
1147 // for (int i = 0; i < n_faces; ++i)
1148 // {
1149 // const Face &f = mesh_.faces[el.fs[i]];
1150 // const int n_face_vertices = f.vs.size();
1151
1152 // const Eigen::RowVector3d e0 = (point(f.vs[0]) - local_pt.row(n_vertices + i));
1153 // const Eigen::RowVector3d e1 = (point(f.vs[1]) - local_pt.row(n_vertices + i));
1154 // const Eigen::RowVector3d normal = e0.cross(e1);
1155 // // const Eigen::RowVector3d check_dir = (node_from_element(e)-p);
1156 // const Eigen::RowVector3d check_dir = (cell_barys.row(e) - point(f.vs[1]));
1157
1158 // const bool reverse_order = normal.dot(check_dir) > 0;
1159
1160 // for (int j = 0; j < n_face_vertices; ++j)
1161 // {
1162 // const int jp = (j + 1) % n_face_vertices;
1163 // if (reverse_order)
1164 // {
1165 // local_faces(face_index, 0) = global_to_local[f.vs[jp]];
1166 // local_faces(face_index, 1) = global_to_local[f.vs[j]];
1167 // }
1168 // else
1169 // {
1170 // local_faces(face_index, 0) = global_to_local[f.vs[j]];
1171 // local_faces(face_index, 1) = global_to_local[f.vs[jp]];
1172 // }
1173 // local_faces(face_index, 2) = n_vertices + i;
1174
1175 // ++face_index;
1176 // }
1177 // }
1178
1179 // local_pts[e] = local_pt;
1180 // local_tris[e] = local_faces;
1181
1182 // total_tris += local_tris[e].rows();
1183 // total_pts += local_pts[e].rows();
1184
1185 // ranges.push_back(total_tris);
1186
1187 // assert(local_pts[e].rows() == local_pt.rows());
1188 // }
1189
1190 // tris.resize(total_tris, 3);
1191 // pts.resize(total_pts, 3);
1192
1193 // int tri_index = 0;
1194 // int pts_index = 0;
1195 // for (std::size_t i = 0; i < local_tris.size(); ++i)
1196 // {
1197 // tris.block(tri_index, 0, local_tris[i].rows(), local_tris[i].cols()) = local_tris[i].array() + pts_index;
1198 // tri_index += local_tris[i].rows();
1199
1200 // pts.block(pts_index, 0, local_pts[i].rows(), local_pts[i].cols()) = local_pts[i];
1201 // pts_index += local_pts[i].rows();
1202 // }
1203 // }
1204
1205 bool CMesh3D::is_boundary_element(const int element_global_id) const
1206 {
1207 const auto &fs = mesh_.elements[element_global_id].fs;
1208
1209 for (auto f_id : fs)
1210 {
1211 if (is_boundary_face(f_id))
1212 return true;
1213 }
1214
1215 const auto &vs = mesh_.elements[element_global_id].vs;
1216
1217 for (auto v_id : vs)
1218 {
1219 if (is_boundary_vertex(v_id))
1220 return true;
1221 }
1222
1223 return false;
1224 }
1225
1226 void CMesh3D::compute_boundary_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &marker)
1227 {
1228 boundary_ids_.resize(n_faces());
1229
1230 for (int f = 0; f < n_faces(); ++f)
1231 {
1232 const bool is_boundary = is_boundary_face(f);
1233 std::vector<int> vs(n_face_vertices(f));
1234 for (int vid = 0; vid < vs.size(); ++vid)
1235 vs[vid] = face_vertex(f, vid);
1236
1237 const auto p = face_barycenter(f);
1238
1239 std::sort(vs.begin(), vs.end());
1240 boundary_ids_[f] = marker(f, vs, p, is_boundary);
1241 }
1242 }
1243
1244 void CMesh3D::compute_body_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &)> &marker)
1245 {
1246 body_ids_.resize(n_elements());
1247 std::fill(body_ids_.begin(), body_ids_.end(), -1);
1248
1249 for (int e = 0; e < n_elements(); ++e)
1250 {
1251 const auto bary = cell_barycenter(e);
1252 body_ids_[e] = marker(e, element_vertices(e), bary);
1253 }
1254 }
1255
1256 void CMesh3D::set_point(const int global_index, const RowVectorNd &p)
1257 {
1258 mesh_.points.col(global_index) = p.transpose();
1259 if (mesh_.vertices[global_index].v.size() == 3)
1260 {
1261 mesh_.vertices[global_index].v[0] = p[0];
1262 mesh_.vertices[global_index].v[1] = p[1];
1263 mesh_.vertices[global_index].v[2] = p[2];
1264 }
1265 }
1266
1267 RowVectorNd CMesh3D::point(const int global_index) const
1268 {
1269 RowVectorNd pt = mesh_.points.col(global_index).transpose();
1270 return pt;
1271 }
1272
1273 RowVectorNd CMesh3D::kernel(const int c) const
1274 {
1275 RowVectorNd pt(3);
1276 pt << mesh_.elements[c].v_in_Kernel[0], mesh_.elements[c].v_in_Kernel[1], mesh_.elements[c].v_in_Kernel[2];
1277 return pt;
1278 }
1279
1281 {
1282 std::vector<ElementType> &ele_tag = elements_tag_;
1283 ele_tag.clear();
1284
1285 ele_tag.resize(mesh_.elements.size());
1286 for (auto &t : ele_tag)
1288
1289 // boundary flags
1290 std::vector<bool> bv_flag(mesh_.vertices.size(), false), be_flag(mesh_.edges.size(), false), bf_flag(mesh_.faces.size(), false);
1291 for (auto f : mesh_.faces)
1292 if (f.boundary)
1293 bf_flag[f.id] = true;
1294 else
1295 {
1296 for (auto nhid : f.neighbor_hs)
1297 if (!mesh_.elements[nhid].hex)
1298 bf_flag[f.id] = true;
1299 }
1300 for (uint32_t i = 0; i < mesh_.faces.size(); ++i)
1301 if (bf_flag[i])
1302 for (uint32_t j = 0; j < mesh_.faces[i].vs.size(); ++j)
1303 {
1304 uint32_t eid = mesh_.faces[i].es[j];
1305 be_flag[eid] = true;
1306 bv_flag[mesh_.faces[i].vs[j]] = true;
1307 }
1308
1309 for (auto &ele : mesh_.elements)
1310 {
1311 if (ele.hex)
1312 {
1313 bool attaching_non_hex = false, on_boundary = false;
1314 ;
1315 for (auto vid : ele.vs)
1316 {
1317 for (auto eleid : mesh_.vertices[vid].neighbor_hs)
1318 if (!mesh_.elements[eleid].hex)
1319 {
1320 attaching_non_hex = true;
1321 break;
1322 }
1323 if (mesh_.vertices[vid].boundary)
1324 {
1325 on_boundary = true;
1326 break;
1327 }
1328 if (on_boundary || attaching_non_hex)
1329 break;
1330 }
1331 if (attaching_non_hex)
1332 {
1333 ele_tag[ele.id] = ElementType::INTERFACE_CUBE;
1334 continue;
1335 }
1336
1337 if (on_boundary)
1338 {
1340 // has no boundary edge--> singular
1341 bool boundary_edge = false, boundary_edge_singular = false, interior_edge_singular = false;
1342 int n_interior_edge_singular = 0;
1343 for (auto eid : ele.es)
1344 {
1345 int en = 0;
1346 if (be_flag[eid])
1347 {
1348 boundary_edge = true;
1349 for (auto nhid : mesh_.edges[eid].neighbor_hs)
1350 if (mesh_.elements[nhid].hex)
1351 en++;
1352 if (en > 2)
1353 boundary_edge_singular = true;
1354 }
1355 else
1356 {
1357 for (auto nhid : mesh_.edges[eid].neighbor_hs)
1358 if (mesh_.elements[nhid].hex)
1359 en++;
1360 if (en != 4)
1361 {
1362 interior_edge_singular = true;
1363 n_interior_edge_singular++;
1364 }
1365 }
1366 }
1367 if (!boundary_edge || boundary_edge_singular || n_interior_edge_singular > 1)
1368 continue;
1369
1370 bool has_singular_v = false, has_iregular_v = false;
1371 int n_in_irregular_v = 0;
1372 for (auto vid : ele.vs)
1373 {
1374 int vn = 0;
1375 if (bv_flag[vid])
1376 {
1377 int nh = 0;
1378 for (auto nhid : mesh_.vertices[vid].neighbor_hs)
1379 if (mesh_.elements[nhid].hex)
1380 nh++;
1381 if (nh > 4)
1382 has_iregular_v = true;
1383 continue; // not sure the conditions
1384 }
1385 else
1386 {
1387 if (mesh_.vertices[vid].neighbor_hs.size() != 8)
1388 n_in_irregular_v++;
1389 int n_irregular_e = 0;
1390 for (auto eid : mesh_.vertices[vid].neighbor_es)
1391 {
1392 if (mesh_.edges[eid].neighbor_hs.size() != 4)
1393 n_irregular_e++;
1394 }
1395 if (n_irregular_e != 0 && n_irregular_e != 2)
1396 {
1397 has_singular_v = true;
1398 break;
1399 }
1400 }
1401 }
1402 int n_irregular_e = 0;
1403 for (auto eid : ele.es)
1404 if (!be_flag[eid] && mesh_.edges[eid].neighbor_hs.size() != 4)
1405 n_irregular_e++;
1406 if (has_singular_v)
1407 continue;
1408 if (!has_singular_v)
1409 {
1410 if (n_irregular_e == 1)
1411 {
1413 }
1414 else if (n_irregular_e == 0 && n_in_irregular_v == 0 && !has_iregular_v)
1415 ele_tag[ele.id] = ElementType::REGULAR_BOUNDARY_CUBE;
1416 else
1417 continue;
1418 }
1419 continue;
1420 }
1421
1422 // type 1
1423 bool has_irregular_v = false;
1424 for (auto vid : ele.vs)
1425 if (mesh_.vertices[vid].neighbor_hs.size() != 8)
1426 {
1427 has_irregular_v = true;
1428 break;
1429 }
1430 if (!has_irregular_v)
1431 {
1432 ele_tag[ele.id] = ElementType::REGULAR_INTERIOR_CUBE;
1433 continue;
1434 }
1435 // type 2
1436 bool has_singular_v = false;
1437 int n_irregular_v = 0;
1438 for (auto vid : ele.vs)
1439 {
1440 if (mesh_.vertices[vid].neighbor_hs.size() != 8)
1441 n_irregular_v++;
1442 int n_irregular_e = 0;
1443 for (auto eid : mesh_.vertices[vid].neighbor_es)
1444 {
1445 if (mesh_.edges[eid].neighbor_hs.size() != 4)
1446 n_irregular_e++;
1447 }
1448 if (n_irregular_e != 0 && n_irregular_e != 2)
1449 {
1450 has_singular_v = true;
1451 break;
1452 }
1453 }
1454 if (!has_singular_v && n_irregular_v == 2)
1455 {
1457 continue;
1458 }
1459
1461 }
1462 else
1463 {
1464 ele_tag[ele.id] = ElementType::INTERIOR_POLYTOPE;
1465 for (auto fid : ele.fs)
1466 if (mesh_.faces[fid].boundary)
1467 {
1468 ele_tag[ele.id] = ElementType::BOUNDARY_POLYTOPE;
1469 break;
1470 }
1471 }
1472 }
1473
1474 // TODO correct?
1475 for (auto &ele : mesh_.elements)
1476 {
1477 if (ele.vs.size() == 4)
1478 ele_tag[ele.id] = ElementType::SIMPLEX;
1479 else if (ele.vs.size() == 5)
1480 ele_tag[ele.id] = ElementType::PYRAMID;
1481 else if (ele.vs.size() == 6)
1482 ele_tag[ele.id] = ElementType::PRISM;
1483 }
1484 }
1485
1486 double CMesh3D::quad_area(const int gid) const
1487 {
1488 const int n_vertices = n_face_vertices(gid);
1489 assert(n_vertices == 4);
1490
1491 const auto &vertices = mesh_.faces[gid].vs;
1492
1493 const auto v1 = point(vertices[0]);
1494 const auto v2 = point(vertices[1]);
1495 const auto v3 = point(vertices[2]);
1496 const auto v4 = point(vertices[3]);
1497
1498 const Eigen::Vector3d e0 = (v2 - v1).transpose();
1499 const Eigen::Vector3d e1 = (v3 - v1).transpose();
1500
1501 const Eigen::Vector3d e2 = (v2 - v4).transpose();
1502 const Eigen::Vector3d e3 = (v3 - v4).transpose();
1503
1504 return e0.cross(e1).norm() / 2 + e2.cross(e3).norm() / 2;
1505 }
1506
1508 {
1509 const int v0 = mesh_.edges[e].vs[0];
1510 const int v1 = mesh_.edges[e].vs[1];
1511 return 0.5 * (point(v0) + point(v1));
1512 }
1513
1515 {
1516 const int n_vertices = n_face_vertices(f);
1517 RowVectorNd bary(3);
1518 bary.setZero();
1519
1520 const auto &vertices = mesh_.faces[f].vs;
1521 for (int lv = 0; lv < n_vertices; ++lv)
1522 {
1523 bary += point(vertices[lv]);
1524 }
1525 return bary / n_vertices;
1526 }
1527
1529 {
1530 const int n_vertices = n_cell_vertices(c);
1531 RowVectorNd bary(3);
1532 bary.setZero();
1533
1534 const auto &vertices = mesh_.elements[c].vs;
1535 for (int lv = 0; lv < n_vertices; ++lv)
1536 {
1537 bary += point(vertices[lv]);
1538 }
1539 return bary / n_vertices;
1540 }
1541
1543 {
1544 m.vertices.clear();
1545 m.edges.clear();
1546 m.faces.clear();
1547 m.vertices.resize(gm.vertices.nb());
1548 m.faces.resize(gm.facets.nb());
1549 for (uint32_t i = 0; i < m.vertices.size(); i++)
1550 {
1551 Vertex v;
1552 v.id = i;
1553 v.v.push_back(gm.vertices.point_ptr(i)[0]);
1554 v.v.push_back(gm.vertices.point_ptr(i)[1]);
1555 v.v.push_back(gm.vertices.point_ptr(i)[2]);
1556 m.vertices[i] = v;
1557 }
1558 m.points.resize(3, m.vertices.size());
1559 for (uint32_t i = 0; i < m.vertices.size(); i++)
1560 {
1561 m.points(0, i) = m.vertices[i].v[0];
1562 m.points(1, i) = m.vertices[i].v[1];
1563 m.points(2, i) = m.vertices[i].v[2];
1564 }
1565
1566 if (m.type == MeshType::TRI || m.type == MeshType::QUA || m.type == MeshType::H_SUR)
1567 {
1568 for (uint32_t i = 0; i < m.faces.size(); i++)
1569 {
1570 Face f;
1571 f.id = i;
1572 f.vs.resize(gm.facets.nb_vertices(i));
1573 for (uint32_t j = 0; j < f.vs.size(); j++)
1574 {
1575 f.vs[j] = gm.facets.vertex(i, j);
1576 }
1577 m.faces[i] = f;
1578 }
1580 }
1581 }
1582
1583 void CMesh3D::elements_boxes(std::vector<std::array<Eigen::Vector3d, 2>> &boxes) const
1584 {
1585 boxes.resize(n_elements());
1586
1587 for (int i = 0; i < n_elements(); ++i)
1588 {
1589 auto &box = boxes[i];
1590 box[0].setConstant(std::numeric_limits<double>::max());
1591 box[1].setConstant(std::numeric_limits<double>::min());
1592
1593 for (int j = 0; j < n_cell_vertices(i); ++j)
1594 {
1595 const int v_id = cell_vertex(i, j);
1596 for (int d = 0; d < 3; ++d)
1597 {
1598 box[0][d] = std::min(box[0][d], point(v_id)[d]);
1599 box[1][d] = std::max(box[1][d], point(v_id)[d]);
1600 }
1601 }
1602 }
1603 }
1604
1605 void CMesh3D::barycentric_coords(const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const
1606 {
1607 assert(is_simplex(el_id));
1608
1609 const auto indices = get_ordered_vertices_from_tet(el_id);
1610
1611 const auto A = point(indices[0]);
1612 const auto B = point(indices[1]);
1613 const auto C = point(indices[2]);
1614 const auto D = point(indices[3]);
1615
1616 igl::barycentric_coordinates(p, A, B, C, D, coord);
1617 }
1618
1619 void CMesh3D::append(const Mesh &mesh)
1620 {
1621 assert(typeid(mesh) == typeid(CMesh3D));
1622 Mesh::append(mesh);
1623
1624 const CMesh3D &mesh3d = dynamic_cast<const CMesh3D &>(mesh);
1625 mesh_.append(mesh3d.mesh_);
1626
1629 }
1630
1631 std::unique_ptr<Mesh> CMesh3D::copy() const
1632 {
1633 return std::make_unique<CMesh3D>(*this);
1634 }
1635
1636 } // namespace mesh
1637} // namespace polyfem
int V
double val
Definition Assembler.cpp:86
int y
int z
int x
RowVectorNd cell_barycenter(const int c) const override
cell barycenter
Definition CMesh3D.cpp:1528
RowVectorNd edge_barycenter(const int e) const override
edge barycenter
Definition CMesh3D.cpp:1507
bool save(const std::string &path) const override
Definition CMesh3D.cpp:411
static void geomesh_2_mesh_storage(const GEO::Mesh &gm, Mesh3DStorage &m)
Definition CMesh3D.cpp:1542
RowVectorNd point(const int vertex_id) const override
point coordinates
Definition CMesh3D.cpp:1267
void normalize() override
normalize the mesh
Definition CMesh3D.cpp:1038
void refine(const int n_refinement, const double t) override
refine the mesh
Definition CMesh3D.cpp:19
int n_edges() const override
number of edges
Definition CMesh3D.hpp:34
int face_vertex(const int f_id, const int lv_id) const override
id of the face vertex
Definition CMesh3D.hpp:44
void compute_boundary_ids(const std::function< int(const size_t, const std::vector< int > &, const RowVectorNd &, bool)> &marker) override
computes boundary selections based on a function
Definition CMesh3D.cpp:1226
RowVectorNd face_barycenter(const int f) const override
face barycenter
Definition CMesh3D.cpp:1514
int n_vertices() const override
number of vertices
Definition CMesh3D.hpp:35
Mesh3DStorage mesh_
Definition CMesh3D.hpp:130
bool is_boundary_face(const int face_global_id) const override
is face boundary
Definition CMesh3D.hpp:52
double quad_area(const int gid) const override
area of a quad face of an hex mesh
Definition CMesh3D.cpp:1486
Navigation3D::Index next_around_face(Navigation3D::Index idx) const override
Definition CMesh3D.hpp:94
bool is_boundary_vertex(const int vertex_global_id) const override
is vertex boundary
Definition CMesh3D.hpp:50
RowVectorNd kernel(const int cell_id) const override
Definition CMesh3D.cpp:1273
void elements_boxes(std::vector< std::array< Eigen::Vector3d, 2 > > &boxes) const override
constructs a box around every element (3d cell, 2d face)
Definition CMesh3D.cpp:1583
void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector< std::vector< int > > &nodes) override
attach high order nodes
Definition CMesh3D.cpp:528
void set_point(const int global_index, const RowVectorNd &p) override
Set the point.
Definition CMesh3D.cpp:1256
void bounding_box(RowVectorNd &min, RowVectorNd &max) const override
computes the bbox of the mesh
Definition CMesh3D.cpp:1031
void append(const Mesh &mesh) override
appends a new mesh to the end of this
Definition CMesh3D.cpp:1619
int n_cells() const override
number of cells
Definition CMesh3D.hpp:32
int n_faces() const override
number of faces
Definition CMesh3D.hpp:33
Navigation3D::Index switch_vertex(Navigation3D::Index idx) const override
Definition CMesh3D.hpp:87
bool load(const std::string &path) override
loads a mesh from the path
Definition CMesh3D.cpp:80
Navigation3D::Index get_index_from_element(int hi, int lf, int lv) const override
Definition CMesh3D.hpp:77
std::unique_ptr< Mesh > copy() const override
Create a copy of the mesh.
Definition CMesh3D.cpp:1631
void compute_elements_tag() override
compute element types, see ElementType
Definition CMesh3D.cpp:1280
Navigation3D::Index switch_edge(Navigation3D::Index idx) const override
Definition CMesh3D.hpp:88
int n_cell_vertices(const int c_id) const override
number of vertices of a cell
Definition CMesh3D.hpp:38
Navigation3D::Index switch_face(Navigation3D::Index idx) const override
Definition CMesh3D.hpp:89
void barycentric_coords(const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const override
constructs barycentric coodiantes for a point p.
Definition CMesh3D.cpp:1605
bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F) override
build a mesh from matrices
Definition CMesh3D.cpp:464
bool is_boundary_element(const int element_global_id) const override
is cell boundary
Definition CMesh3D.cpp:1205
int n_face_vertices(const int f_id) const override
number of vertices of a face
Definition CMesh3D.hpp:37
int cell_vertex(const int c_id, const int lv_id) const override
id of the vertex of a cell
Definition CMesh3D.hpp:41
void compute_body_ids(const std::function< int(const size_t, const std::vector< int > &, const RowVectorNd &)> &marker) override
computes boundary selections based on a function
Definition CMesh3D.cpp:1244
void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1) const override
Get all the edges.
Definition Mesh3D.cpp:33
virtual std::array< int, 4 > get_ordered_vertices_from_tet(const int element_index) const
Definition Mesh3D.cpp:373
void append(const Mesh3DStorage &other)
std::vector< Element > elements
std::vector< Vertex > vertices
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:41
int n_elements() const
utitlity to return the number of elements, cells or faces in 3d and 2d
Definition Mesh.hpp:163
Eigen::MatrixXi orders_
list of geometry orders, one per cell
Definition Mesh.hpp:705
std::vector< ElementType > elements_tag_
list of element types
Definition Mesh.hpp:697
Eigen::MatrixXi in_ordered_faces_
Order of the input faces, TODO: change to std::vector of Eigen::Vector.
Definition Mesh.hpp:723
bool is_simplex(const int el_id) const
checks if element is simplex
Definition Mesh.cpp:422
std::vector< int > boundary_ids_
list of surface labels
Definition Mesh.hpp:701
std::vector< CellNodes > cell_nodes_
high-order nodes associates to cells
Definition Mesh.hpp:714
std::vector< int > element_vertices(const int el_id) const
list of vids of an element
Definition Mesh.hpp:231
std::vector< int > body_ids_
list of volume labels
Definition Mesh.hpp:703
std::vector< FaceNodes > face_nodes_
high-order nodes associates to faces
Definition Mesh.hpp:712
const std::vector< ElementType > & elements_tag() const
Returns the elements types.
Definition Mesh.hpp:428
std::vector< EdgeNodes > edge_nodes_
high-order nodes associates to edges
Definition Mesh.hpp:710
std::vector< std::vector< int > > faces() const
list of sorted faces.
Definition Mesh.cpp:453
Eigen::MatrixXi in_ordered_edges_
Order of the input edges.
Definition Mesh.hpp:721
virtual void append(const Mesh &mesh)
appends a new mesh to the end of this
Definition Mesh.cpp:504
Eigen::VectorXi in_ordered_vertices_
Order of the input vertices.
Definition Mesh.hpp:719
void refine_red_refinement_tet(Mesh3DStorage &M, int iter)
void build_connectivity(Mesh3DStorage &hmi)
void refine_catmul_clark_polar(Mesh3DStorage &M, int iter, bool reverse, std::vector< int > &Parents)
void prepare_mesh(Mesh3DStorage &M)
int find(const Eigen::VectorXi &vec, int x)
Definition NCMesh2D.cpp:289
@ REGULAR_INTERIOR_CUBE
Triangle/tet element.
@ REGULAR_BOUNDARY_CUBE
Quad/Hex incident to more than 1 singular vertices (should not happen in 2D)
@ MULTI_SINGULAR_BOUNDARY_CUBE
Quad incident to exactly 1 singular vertex (in 2D); hex incident to exactly 1 singular interior edge,...
@ PRISM
Boundary polytope.
@ MULTI_SINGULAR_INTERIOR_CUBE
Quad/hex incident to exactly 1 singular vertex (in 2D) or edge (in 3D)
@ SIMPLE_SINGULAR_INTERIOR_CUBE
Regular quad/hex inside a 3^n patch.
@ INTERFACE_CUBE
Boundary hex that is not regular nor SimpleSingularBoundaryCube.
@ INTERIOR_POLYTOPE
Quad/hex that is at the interface with a polytope (if a cube has both external boundary and and inter...
@ SIMPLE_SINGULAR_BOUNDARY_CUBE
Boundary quad/hex, where all boundary vertices/edges are incident to at most 2 quads/hexes.
@ BOUNDARY_POLYTOPE
Interior polytope.
void to_geogram_mesh(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, GEO::Mesh &M)
Converts a triangle mesh to a Geogram mesh.
bool endswith(const std::string &str, const std::string &suffix)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:44
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::vector< double > v_in_Kernel
std::vector< uint32_t > fs
std::vector< uint32_t > vs
std::vector< bool > fs_flag
std::vector< uint32_t > vs
std::vector< double > v