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() == 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 if (F.cols() == 4)
482 M.cells.create_tets((int)F.rows());
483 else if (F.cols() == 6)
484 M.cells.create_prisms((int)F.rows());
485 else if (F.cols() == 8)
486 M.cells.create_hexes((int)F.rows());
487 else
488 {
489 throw std::runtime_error("Mesh format not supported");
490 }
491
492 static const std::vector<int> permute_tet = {0, 1, 2, 3};
493 static const std::vector<int> permute_prism = {0, 1, 2, 3, 4, 5};
494 // polyfem uses the msh file format for hexes ordering!
495 static const std::vector<int> permute_hex = {1, 0, 2, 3, 5, 4, 6, 7};
496 const std::vector<int> permute = F.cols() == 4 ? permute_tet : (F.cols() == 6 ? permute_prism : permute_hex);
497
498 for (int c = 0; c < (int)M.cells.nb(); ++c)
499 {
500 for (int lv = 0; lv < F.cols(); ++lv)
501 {
502 M.cells.set_vertex(c, lv, F(c, permute[lv]));
503 }
504 }
505 M.cells.connect();
506
507 return load(M);
508 }
509
510 void CMesh3D::attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector<std::vector<int>> &nodes)
511 {
512 edge_nodes_.clear();
513 face_nodes_.clear();
514 cell_nodes_.clear();
515
516 edge_nodes_.resize(n_edges());
517 face_nodes_.resize(n_faces());
518 cell_nodes_.resize(n_cells());
519
520 orders_.resize(n_cells(), 1);
521
522 const auto attach_p2 = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
523 auto &n = edge_nodes_[index.edge];
524
525 if (n.nodes.size() > 0)
526 return;
527
528 n.v1 = index.vertex;
529 n.v2 = switch_vertex(index).vertex;
530
531 const int n_v1 = index.vertex;
532 const int n_v2 = switch_vertex(index).vertex;
533
534 int node_index = 0;
535
536 if ((n_v1 == nodes_ids[0] && n_v2 == nodes_ids[1]) || (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[1]))
537 node_index = 4;
538 else if ((n_v1 == nodes_ids[1] && n_v2 == nodes_ids[2]) || (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[2]))
539 node_index = 5;
540 else if ((n_v1 == nodes_ids[2] && n_v2 == nodes_ids[3]) || (n_v2 == nodes_ids[2] && n_v1 == nodes_ids[3]))
541 node_index = 8;
542
543 else if ((n_v1 == nodes_ids[0] && n_v2 == nodes_ids[3]) || (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[3]))
544 node_index = 7;
545 else if ((n_v1 == nodes_ids[0] && n_v2 == nodes_ids[2]) || (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[2]))
546 node_index = 6;
547 else
548 node_index = 9;
549
550 n.nodes.resize(1, 3);
551 n.nodes << V(nodes_ids[node_index], 0), V(nodes_ids[node_index], 1), V(nodes_ids[node_index], 2);
552 n.nodes_ids.push_back(nodes_ids[node_index]);
553 };
554
555 const auto attach_p3 = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
556 auto &n = edge_nodes_[index.edge];
557
558 if (n.nodes.size() > 0)
559 return;
560
561 n.v1 = index.vertex;
562 n.v2 = switch_vertex(index).vertex;
563
564 const int n_v1 = index.vertex;
565 const int n_v2 = switch_vertex(index).vertex;
566
567 int node_index1 = 0;
568 int node_index2 = 0;
569 if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[1])
570 {
571 node_index1 = 4;
572 node_index2 = 5;
573 }
574 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[1])
575 {
576 node_index1 = 5;
577 node_index2 = 4;
578 }
579 else if (n_v1 == nodes_ids[1] && n_v2 == nodes_ids[2])
580 {
581 node_index1 = 6;
582 node_index2 = 7;
583 }
584 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[2])
585 {
586 node_index1 = 7;
587 node_index2 = 6;
588 }
589 else if (n_v1 == nodes_ids[2] && n_v2 == nodes_ids[3])
590 {
591 node_index1 = 13;
592 node_index2 = 12;
593 }
594 else if (n_v2 == nodes_ids[2] && n_v1 == nodes_ids[3])
595 {
596 node_index1 = 12;
597 node_index2 = 13;
598 }
599
600 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[3])
601 {
602 node_index1 = 11;
603 node_index2 = 10;
604 }
605 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[3])
606 {
607 node_index1 = 10;
608 node_index2 = 11;
609 }
610 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[2])
611 {
612 node_index1 = 9;
613 node_index2 = 8;
614 }
615 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[2])
616 {
617 node_index1 = 8;
618 node_index2 = 9;
619 }
620
621 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[3])
622 {
623 node_index1 = 14;
624 node_index2 = 15;
625 }
626 else
627 {
628 node_index1 = 15;
629 node_index2 = 14;
630 }
631
632 n.nodes.resize(2, 3);
633 n.nodes.row(0) << V(nodes_ids[node_index1], 0), V(nodes_ids[node_index1], 1), V(nodes_ids[node_index1], 2);
634 n.nodes.row(1) << V(nodes_ids[node_index2], 0), V(nodes_ids[node_index2], 1), V(nodes_ids[node_index2], 2);
635 n.nodes_ids.push_back(nodes_ids[node_index1]);
636 n.nodes_ids.push_back(nodes_ids[node_index2]);
637 };
638
639 const auto attach_p3_face = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids, int id) {
640 auto &n = face_nodes_[index.face];
641 if (n.nodes.size() <= 0)
642 {
643 n.v1 = face_vertex(index.face, 0);
644 n.v2 = face_vertex(index.face, 1);
645 n.v3 = face_vertex(index.face, 2);
646 n.nodes.resize(1, 3);
647 n.nodes << V(nodes_ids[id], 0), V(nodes_ids[id], 1), V(nodes_ids[id], 2);
648 n.nodes_ids.push_back(nodes_ids[id]);
649 }
650 };
651
652 const auto attach_p4 = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
653 auto &n = edge_nodes_[index.edge];
654
655 if (n.nodes.size() > 0)
656 return;
657
658 n.v1 = index.vertex;
659 n.v2 = switch_vertex(index).vertex;
660
661 const int n_v1 = index.vertex;
662 const int n_v2 = switch_vertex(index).vertex;
663
664 int node_index1 = 0;
665 int node_index2 = 0;
666 int node_index3 = 0;
667 if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[1])
668 {
669 node_index1 = 4;
670 node_index2 = 5;
671 node_index3 = 6;
672 }
673 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[1])
674 {
675 node_index1 = 6;
676 node_index2 = 5;
677 node_index3 = 4;
678 }
679
680 else if (n_v1 == nodes_ids[1] && n_v2 == nodes_ids[2])
681 {
682 node_index1 = 7;
683 node_index2 = 8;
684 node_index3 = 9;
685 }
686 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[2])
687 {
688 node_index1 = 9;
689 node_index2 = 8;
690 node_index3 = 7;
691 }
692
693 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[2])
694 {
695 node_index1 = 10;
696 node_index2 = 11;
697 node_index3 = 12;
698 }
699 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[2])
700 {
701 node_index1 = 12;
702 node_index2 = 11;
703 node_index3 = 10;
704 }
705
706 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[3])
707 {
708 node_index1 = 13;
709 node_index2 = 14;
710 node_index3 = 15;
711 }
712 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[3])
713 {
714 node_index1 = 15;
715 node_index2 = 14;
716 node_index3 = 13;
717 }
718
719 else if (n_v2 == nodes_ids[2] && n_v1 == nodes_ids[3])
720 {
721 node_index1 = 16;
722 node_index2 = 17;
723 node_index3 = 18;
724 }
725 else if (n_v1 == nodes_ids[2] && n_v2 == nodes_ids[3])
726 {
727 node_index1 = 18;
728 node_index2 = 17;
729 node_index3 = 16;
730 }
731
732 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[3])
733 {
734 node_index1 = 19;
735 node_index2 = 20;
736 node_index3 = 21;
737 }
738 else if (n_v1 == nodes_ids[1] && n_v2 == nodes_ids[3])
739 {
740 node_index1 = 21;
741 node_index2 = 20;
742 node_index3 = 19;
743 }
744 else
745 {
746 assert(false);
747 }
748
749 n.nodes.resize(3, 3);
750 n.nodes.row(0) << V(nodes_ids[node_index1], 0), V(nodes_ids[node_index1], 1), V(nodes_ids[node_index1], 2);
751 n.nodes.row(1) << V(nodes_ids[node_index2], 0), V(nodes_ids[node_index2], 1), V(nodes_ids[node_index2], 2);
752 n.nodes.row(2) << V(nodes_ids[node_index3], 0), V(nodes_ids[node_index3], 1), V(nodes_ids[node_index3], 2);
753 n.nodes_ids.push_back(nodes_ids[node_index1]);
754 n.nodes_ids.push_back(nodes_ids[node_index2]);
755 n.nodes_ids.push_back(nodes_ids[node_index3]);
756 };
757
758 const auto attach_p4_face = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
759 auto &n = face_nodes_[index.face];
760 if (n.nodes.size() <= 0)
761 {
762 n.v1 = face_vertex(index.face, 0);
763 n.v2 = face_vertex(index.face, 1);
764 n.v3 = face_vertex(index.face, 2);
765
766 std::array<int, 3> vid = {{n.v1, n.v2, n.v3}};
767 std::sort(vid.begin(), vid.end());
768
769 std::array<int, 3> c1 = {{nodes_ids[0], nodes_ids[1], nodes_ids[2]}}; // 22
770 std::array<int, 3> c2 = {{nodes_ids[0], nodes_ids[1], nodes_ids[3]}}; // 25
771 std::array<int, 3> c3 = {{nodes_ids[0], nodes_ids[2], nodes_ids[3]}}; // 28
772 std::array<int, 3> c4 = {{nodes_ids[1], nodes_ids[2], nodes_ids[3]}}; // 31
773
774 std::sort(c1.begin(), c1.end());
775 std::sort(c2.begin(), c2.end());
776 std::sort(c3.begin(), c3.end());
777 std::sort(c4.begin(), c4.end());
778
779 int id = 0;
780 int index0 = 0;
781 int index1 = 1;
782 int index2 = 2;
783 if (vid == c1)
784 {
785 id = 22;
786 n.v1 = nodes_ids[0];
787 n.v2 = nodes_ids[1];
788 n.v3 = nodes_ids[2];
789
790 index0 = 0;
791 index1 = 2;
792 index2 = 1; // ok
793 }
794 else if (vid == c2)
795 {
796 id = 25;
797 n.v1 = nodes_ids[0];
798 n.v2 = nodes_ids[1];
799 n.v3 = nodes_ids[3];
800
801 index0 = 0;
802 index1 = 1;
803 index2 = 2; // ok
804 }
805 else if (vid == c3)
806 {
807 id = 28;
808 n.v1 = nodes_ids[0];
809 n.v2 = nodes_ids[2];
810 n.v3 = nodes_ids[3];
811
812 index0 = 0;
813 index1 = 2;
814 index2 = 1; // ok
815 }
816 else if (vid == c4)
817 {
818 id = 31;
819 n.v1 = nodes_ids[1];
820 n.v2 = nodes_ids[2];
821 n.v3 = nodes_ids[3];
822
823 index0 = 1;
824 index1 = 2;
825 index2 = 0; // ok
826 }
827 else
828 {
829 // the face nees to be one of the 4 above
830 assert(false);
831 }
832
833 n.nodes.resize(3, 3);
834 assert(id + index0 < nodes_ids.size());
835 assert(id + index1 < nodes_ids.size());
836 assert(id + index2 < nodes_ids.size());
837 n.nodes.row(0) << V(nodes_ids[id + index0], 0), V(nodes_ids[id + index0], 1), V(nodes_ids[id + index0], 2);
838 n.nodes.row(1) << V(nodes_ids[id + index1], 0), V(nodes_ids[id + index1], 1), V(nodes_ids[id + index1], 2);
839 n.nodes.row(2) << V(nodes_ids[id + index2], 0), V(nodes_ids[id + index2], 1), V(nodes_ids[id + index2], 2);
840 n.nodes_ids.push_back(nodes_ids[id + index0]);
841 n.nodes_ids.push_back(nodes_ids[id + index1]);
842 n.nodes_ids.push_back(nodes_ids[id + index2]);
843 }
844 };
845
846 const auto attach_p4_cell = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
847 auto &n = cell_nodes_[index.element];
848 assert(nodes_ids.size() == 35);
849 assert(n.nodes.size() == 0);
850
851 if (n.nodes.size() <= 0)
852 {
853 n.v1 = cell_vertex(index.element, 0);
854 n.v2 = cell_vertex(index.element, 1);
855 n.v3 = cell_vertex(index.element, 2);
856 n.v4 = cell_vertex(index.element, 3);
857 n.nodes.resize(1, 3);
858
859 n.nodes << V(nodes_ids[34], 0), V(nodes_ids[34], 1), V(nodes_ids[34], 2);
860 n.nodes_ids.push_back(nodes_ids[34]);
861 }
862 };
863
864 assert(nodes.size() == n_cells());
865
866 for (int c = 0; c < n_cells(); ++c)
867 {
868 auto index = get_index_from_element(c);
869
870 const auto &nodes_ids = nodes[c];
871
872 if (nodes_ids.size() == 4)
873 {
874 orders_(c) = 1;
875 continue;
876 }
877 // P2
878 else if (nodes_ids.size() == 10)
879 {
880 orders_(c) = 2;
881
882 for (int le = 0; le < 3; ++le)
883 {
884 attach_p2(index, nodes_ids);
885 index = next_around_face(index);
886 }
887
888 index = switch_vertex(switch_edge(switch_face(index)));
889 attach_p2(index, nodes_ids);
890
891 index = switch_edge(index);
892 attach_p2(index, nodes_ids);
893
894 index = switch_edge(switch_face(index));
895 attach_p2(index, nodes_ids);
896 }
897 // P3
898 else if (nodes_ids.size() == 20)
899 {
900 orders_(c) = 3;
901
902 for (int le = 0; le < 3; ++le)
903 {
904 attach_p3(index, nodes_ids);
905 index = next_around_face(index);
906 }
907
908 {
909 index = switch_vertex(switch_edge(switch_face(index)));
910 attach_p3(index, nodes_ids);
911
912 index = switch_edge(index);
913 attach_p3(index, nodes_ids);
914
915 index = switch_edge(switch_face(index));
916 attach_p3(index, nodes_ids);
917 }
918
919 {
920 index = get_index_from_element(c);
921 std::array<int, 4> indices;
922
923 {
924 std::array<int, 3> f16 = {{nodes_ids[0], nodes_ids[1], nodes_ids[2]}};
925 std::array<int, 3> f17 = {{nodes_ids[3], nodes_ids[1], nodes_ids[0]}};
926 std::array<int, 3> f18 = {{nodes_ids[0], nodes_ids[2], nodes_ids[3]}};
927 std::array<int, 3> f19 = {{nodes_ids[1], nodes_ids[2], nodes_ids[3]}};
928 std::sort(f16.begin(), f16.end());
929 std::sort(f17.begin(), f17.end());
930 std::sort(f18.begin(), f18.end());
931 std::sort(f19.begin(), f19.end());
932
933 auto tmp = index;
934 std::array<int, 3> f0 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
935 tmp = switch_face(index);
936 std::array<int, 3> f1 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
937 tmp = switch_face(next_around_face(index));
938 std::array<int, 3> f2 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
940 std::array<int, 3> f3 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
941
942 std::sort(f0.begin(), f0.end());
943 std::sort(f1.begin(), f1.end());
944 std::sort(f2.begin(), f2.end());
945 std::sort(f3.begin(), f3.end());
946
947 const std::array<std::array<int, 3>, 4> faces = {{f0, f1, f2, f3}};
948 const std::array<std::array<int, 3>, 4> nodes = {{f16, f17, f18, f19}};
949 for (int i = 0; i < 4; ++i)
950 {
951 const auto &f = faces[i];
952 bool found = false;
953 for (int j = 0; j < 4; ++j)
954 {
955 if (nodes[j] == f)
956 {
957 indices[i] = j + 16;
958 found = true;
959 break;
960 }
961 }
962
963 assert(found);
964 }
965 }
966
967 attach_p3_face(index, nodes_ids, indices[0]);
968 attach_p3_face(switch_face(index), nodes_ids, indices[1]);
969 attach_p3_face(switch_face(next_around_face(index)), nodes_ids, indices[2]);
970 attach_p3_face(switch_face(next_around_face(next_around_face(index))), nodes_ids, indices[3]);
971 }
972 }
973 // P4
974 else if (nodes_ids.size() == 35)
975 {
976 orders_(c) = 4;
977 for (int le = 0; le < 3; ++le)
978 {
979 attach_p4(index, nodes_ids);
980 index = next_around_face(index);
981 }
982
983 {
984 index = switch_vertex(switch_edge(switch_face(index)));
985 attach_p4(index, nodes_ids);
986
987 index = switch_edge(index);
988 attach_p4(index, nodes_ids);
989
990 index = switch_edge(switch_face(index));
991 attach_p4(index, nodes_ids);
992 }
993
994 {
995 index = get_index_from_element(c);
996
997 attach_p4_face(index, nodes_ids);
998 attach_p4_face(switch_face(index), nodes_ids);
999 attach_p4_face(switch_face(next_around_face(index)), nodes_ids);
1000 attach_p4_face(switch_face(next_around_face(next_around_face(index))), nodes_ids);
1001 }
1002
1003 attach_p4_cell(get_index_from_element(c), nodes_ids);
1004 }
1005 // unsupported
1006 else
1007 {
1008 assert(false);
1009 }
1010 }
1011 }
1012
1014 {
1015 const auto &V = mesh_.points;
1016 min = V.rowwise().minCoeff().transpose();
1017 max = V.rowwise().maxCoeff().transpose();
1018 }
1019
1021 {
1022 RowVectorNd minV, maxV;
1023 bounding_box(minV, maxV);
1024 auto &V = mesh_.points;
1025 const auto shift = V.rowwise().minCoeff().eval();
1026 const double scaling = 1.0 / (V.rowwise().maxCoeff() - V.rowwise().minCoeff()).maxCoeff();
1027 V = (V.colwise() - shift) * scaling;
1028
1029 for (int i = 0; i < n_cells(); ++i)
1030 {
1031 for (int d = 0; d < 3; ++d)
1032 {
1033 auto val = mesh_.elements[i].v_in_Kernel[d];
1034 mesh_.elements[i].v_in_Kernel[d] = (val - shift(d)) * scaling;
1035 }
1036 }
1037
1038 for (auto &n : edge_nodes_)
1039 {
1040 if (n.nodes.size() > 0)
1041 n.nodes = (n.nodes.rowwise() - shift.transpose()) * scaling;
1042 }
1043 for (auto &n : face_nodes_)
1044 {
1045 if (n.nodes.size() > 0)
1046 n.nodes = (n.nodes.rowwise() - shift.transpose()) * scaling;
1047 }
1048 for (auto &n : cell_nodes_)
1049 {
1050 if (n.nodes.size() > 0)
1051 n.nodes = (n.nodes.rowwise() - shift.transpose()) * scaling;
1052 }
1053
1054 logger().debug("-- bbox before normalization:");
1055 logger().debug(" min : {}", minV);
1056 logger().debug(" max : {}", maxV);
1057 logger().debug(" extent: {}", maxV - minV);
1058 bounding_box(minV, maxV);
1059 logger().debug("-- bbox after normalization:");
1060 logger().debug(" min : {}", minV);
1061 logger().debug(" max : {}", maxV);
1062 logger().debug(" extent: {}", maxV - minV);
1063
1064 // V.row(1) /= 100.;
1065
1066 // for(int i = 0; i < n_cells(); ++i)
1067 // {
1068 // mesh_.elements[i].v_in_Kernel[1] /= 100.;
1069 // }
1070
1071 Eigen::MatrixXd p0, p1, p;
1072 get_edges(p0, p1);
1073 p = p0 - p1;
1074 logger().debug("-- edge length after normalization:");
1075 logger().debug(" min: ", p.rowwise().norm().minCoeff());
1076 logger().debug(" max: ", p.rowwise().norm().maxCoeff());
1077 logger().debug(" avg: ", p.rowwise().norm().mean());
1078 }
1079
1080 // void CMesh3D::triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector<int> &ranges) const
1081 // {
1082 // ranges.clear();
1083
1084 // std::vector<Eigen::MatrixXi> local_tris(mesh_.elements.size());
1085 // std::vector<Eigen::MatrixXd> local_pts(mesh_.elements.size());
1086 // Eigen::MatrixXi tets;
1087
1088 // int total_tris = 0;
1089 // int total_pts = 0;
1090
1091 // ranges.push_back(0);
1092
1093 // Eigen::MatrixXd face_barys;
1094 // face_barycenters(face_barys);
1095
1096 // Eigen::MatrixXd cell_barys;
1097 // cell_barycenters(cell_barys);
1098
1099 // for (std::size_t e = 0; e < mesh_.elements.size(); ++e)
1100 // {
1101 // const Element &el = mesh_.elements[e];
1102
1103 // const int n_vertices = el.vs.size();
1104 // const int n_faces = el.fs.size();
1105
1106 // Eigen::MatrixXd local_pt(n_vertices + n_faces, 3);
1107
1108 // std::map<int, int> global_to_local;
1109
1110 // for (int i = 0; i < n_vertices; ++i)
1111 // {
1112 // const int global_index = el.vs[i];
1113 // local_pt.row(i) = mesh_.points.col(global_index).transpose();
1114 // global_to_local[global_index] = i;
1115 // }
1116
1117 // int n_local_faces = 0;
1118 // for (int i = 0; i < n_faces; ++i)
1119 // {
1120 // const Face &f = mesh_.faces[el.fs[i]];
1121 // n_local_faces += f.vs.size();
1122
1123 // local_pt.row(n_vertices + i) = face_barys.row(f.id); // node_from_face(f.id);
1124 // }
1125
1126 // Eigen::MatrixXi local_faces(n_local_faces, 3);
1127
1128 // int face_index = 0;
1129 // for (int i = 0; i < n_faces; ++i)
1130 // {
1131 // const Face &f = mesh_.faces[el.fs[i]];
1132 // const int n_face_vertices = f.vs.size();
1133
1134 // const Eigen::RowVector3d e0 = (point(f.vs[0]) - local_pt.row(n_vertices + i));
1135 // const Eigen::RowVector3d e1 = (point(f.vs[1]) - local_pt.row(n_vertices + i));
1136 // const Eigen::RowVector3d normal = e0.cross(e1);
1137 // // const Eigen::RowVector3d check_dir = (node_from_element(e)-p);
1138 // const Eigen::RowVector3d check_dir = (cell_barys.row(e) - point(f.vs[1]));
1139
1140 // const bool reverse_order = normal.dot(check_dir) > 0;
1141
1142 // for (int j = 0; j < n_face_vertices; ++j)
1143 // {
1144 // const int jp = (j + 1) % n_face_vertices;
1145 // if (reverse_order)
1146 // {
1147 // local_faces(face_index, 0) = global_to_local[f.vs[jp]];
1148 // local_faces(face_index, 1) = global_to_local[f.vs[j]];
1149 // }
1150 // else
1151 // {
1152 // local_faces(face_index, 0) = global_to_local[f.vs[j]];
1153 // local_faces(face_index, 1) = global_to_local[f.vs[jp]];
1154 // }
1155 // local_faces(face_index, 2) = n_vertices + i;
1156
1157 // ++face_index;
1158 // }
1159 // }
1160
1161 // local_pts[e] = local_pt;
1162 // local_tris[e] = local_faces;
1163
1164 // total_tris += local_tris[e].rows();
1165 // total_pts += local_pts[e].rows();
1166
1167 // ranges.push_back(total_tris);
1168
1169 // assert(local_pts[e].rows() == local_pt.rows());
1170 // }
1171
1172 // tris.resize(total_tris, 3);
1173 // pts.resize(total_pts, 3);
1174
1175 // int tri_index = 0;
1176 // int pts_index = 0;
1177 // for (std::size_t i = 0; i < local_tris.size(); ++i)
1178 // {
1179 // tris.block(tri_index, 0, local_tris[i].rows(), local_tris[i].cols()) = local_tris[i].array() + pts_index;
1180 // tri_index += local_tris[i].rows();
1181
1182 // pts.block(pts_index, 0, local_pts[i].rows(), local_pts[i].cols()) = local_pts[i];
1183 // pts_index += local_pts[i].rows();
1184 // }
1185 // }
1186
1187 bool CMesh3D::is_boundary_element(const int element_global_id) const
1188 {
1189 const auto &fs = mesh_.elements[element_global_id].fs;
1190
1191 for (auto f_id : fs)
1192 {
1193 if (is_boundary_face(f_id))
1194 return true;
1195 }
1196
1197 const auto &vs = mesh_.elements[element_global_id].vs;
1198
1199 for (auto v_id : vs)
1200 {
1201 if (is_boundary_vertex(v_id))
1202 return true;
1203 }
1204
1205 return false;
1206 }
1207
1208 void CMesh3D::compute_boundary_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &marker)
1209 {
1210 boundary_ids_.resize(n_faces());
1211
1212 for (int f = 0; f < n_faces(); ++f)
1213 {
1214 const bool is_boundary = is_boundary_face(f);
1215 std::vector<int> vs(n_face_vertices(f));
1216 for (int vid = 0; vid < vs.size(); ++vid)
1217 vs[vid] = face_vertex(f, vid);
1218
1219 const auto p = face_barycenter(f);
1220
1221 std::sort(vs.begin(), vs.end());
1222 boundary_ids_[f] = marker(f, vs, p, is_boundary);
1223 }
1224 }
1225
1226 void CMesh3D::compute_body_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &)> &marker)
1227 {
1228 body_ids_.resize(n_elements());
1229 std::fill(body_ids_.begin(), body_ids_.end(), -1);
1230
1231 for (int e = 0; e < n_elements(); ++e)
1232 {
1233 const auto bary = cell_barycenter(e);
1234 body_ids_[e] = marker(e, element_vertices(e), bary);
1235 }
1236 }
1237
1238 void CMesh3D::set_point(const int global_index, const RowVectorNd &p)
1239 {
1240 mesh_.points.col(global_index) = p.transpose();
1241 if (mesh_.vertices[global_index].v.size() == 3)
1242 {
1243 mesh_.vertices[global_index].v[0] = p[0];
1244 mesh_.vertices[global_index].v[1] = p[1];
1245 mesh_.vertices[global_index].v[2] = p[2];
1246 }
1247 }
1248
1249 RowVectorNd CMesh3D::point(const int global_index) const
1250 {
1251 RowVectorNd pt = mesh_.points.col(global_index).transpose();
1252 return pt;
1253 }
1254
1255 RowVectorNd CMesh3D::kernel(const int c) const
1256 {
1257 RowVectorNd pt(3);
1258 pt << mesh_.elements[c].v_in_Kernel[0], mesh_.elements[c].v_in_Kernel[1], mesh_.elements[c].v_in_Kernel[2];
1259 return pt;
1260 }
1261
1263 {
1264 std::vector<ElementType> &ele_tag = elements_tag_;
1265 ele_tag.clear();
1266
1267 ele_tag.resize(mesh_.elements.size());
1268 for (auto &t : ele_tag)
1270
1271 // boundary flags
1272 std::vector<bool> bv_flag(mesh_.vertices.size(), false), be_flag(mesh_.edges.size(), false), bf_flag(mesh_.faces.size(), false);
1273 for (auto f : mesh_.faces)
1274 if (f.boundary)
1275 bf_flag[f.id] = true;
1276 else
1277 {
1278 for (auto nhid : f.neighbor_hs)
1279 if (!mesh_.elements[nhid].hex)
1280 bf_flag[f.id] = true;
1281 }
1282 for (uint32_t i = 0; i < mesh_.faces.size(); ++i)
1283 if (bf_flag[i])
1284 for (uint32_t j = 0; j < mesh_.faces[i].vs.size(); ++j)
1285 {
1286 uint32_t eid = mesh_.faces[i].es[j];
1287 be_flag[eid] = true;
1288 bv_flag[mesh_.faces[i].vs[j]] = true;
1289 }
1290
1291 for (auto &ele : mesh_.elements)
1292 {
1293 if (ele.hex)
1294 {
1295 bool attaching_non_hex = false, on_boundary = false;
1296 ;
1297 for (auto vid : ele.vs)
1298 {
1299 for (auto eleid : mesh_.vertices[vid].neighbor_hs)
1300 if (!mesh_.elements[eleid].hex)
1301 {
1302 attaching_non_hex = true;
1303 break;
1304 }
1305 if (mesh_.vertices[vid].boundary)
1306 {
1307 on_boundary = true;
1308 break;
1309 }
1310 if (on_boundary || attaching_non_hex)
1311 break;
1312 }
1313 if (attaching_non_hex)
1314 {
1315 ele_tag[ele.id] = ElementType::INTERFACE_CUBE;
1316 continue;
1317 }
1318
1319 if (on_boundary)
1320 {
1322 // has no boundary edge--> singular
1323 bool boundary_edge = false, boundary_edge_singular = false, interior_edge_singular = false;
1324 int n_interior_edge_singular = 0;
1325 for (auto eid : ele.es)
1326 {
1327 int en = 0;
1328 if (be_flag[eid])
1329 {
1330 boundary_edge = true;
1331 for (auto nhid : mesh_.edges[eid].neighbor_hs)
1332 if (mesh_.elements[nhid].hex)
1333 en++;
1334 if (en > 2)
1335 boundary_edge_singular = true;
1336 }
1337 else
1338 {
1339 for (auto nhid : mesh_.edges[eid].neighbor_hs)
1340 if (mesh_.elements[nhid].hex)
1341 en++;
1342 if (en != 4)
1343 {
1344 interior_edge_singular = true;
1345 n_interior_edge_singular++;
1346 }
1347 }
1348 }
1349 if (!boundary_edge || boundary_edge_singular || n_interior_edge_singular > 1)
1350 continue;
1351
1352 bool has_singular_v = false, has_iregular_v = false;
1353 int n_in_irregular_v = 0;
1354 for (auto vid : ele.vs)
1355 {
1356 int vn = 0;
1357 if (bv_flag[vid])
1358 {
1359 int nh = 0;
1360 for (auto nhid : mesh_.vertices[vid].neighbor_hs)
1361 if (mesh_.elements[nhid].hex)
1362 nh++;
1363 if (nh > 4)
1364 has_iregular_v = true;
1365 continue; // not sure the conditions
1366 }
1367 else
1368 {
1369 if (mesh_.vertices[vid].neighbor_hs.size() != 8)
1370 n_in_irregular_v++;
1371 int n_irregular_e = 0;
1372 for (auto eid : mesh_.vertices[vid].neighbor_es)
1373 {
1374 if (mesh_.edges[eid].neighbor_hs.size() != 4)
1375 n_irregular_e++;
1376 }
1377 if (n_irregular_e != 0 && n_irregular_e != 2)
1378 {
1379 has_singular_v = true;
1380 break;
1381 }
1382 }
1383 }
1384 int n_irregular_e = 0;
1385 for (auto eid : ele.es)
1386 if (!be_flag[eid] && mesh_.edges[eid].neighbor_hs.size() != 4)
1387 n_irregular_e++;
1388 if (has_singular_v)
1389 continue;
1390 if (!has_singular_v)
1391 {
1392 if (n_irregular_e == 1)
1393 {
1395 }
1396 else if (n_irregular_e == 0 && n_in_irregular_v == 0 && !has_iregular_v)
1397 ele_tag[ele.id] = ElementType::REGULAR_BOUNDARY_CUBE;
1398 else
1399 continue;
1400 }
1401 continue;
1402 }
1403
1404 // type 1
1405 bool has_irregular_v = false;
1406 for (auto vid : ele.vs)
1407 if (mesh_.vertices[vid].neighbor_hs.size() != 8)
1408 {
1409 has_irregular_v = true;
1410 break;
1411 }
1412 if (!has_irregular_v)
1413 {
1414 ele_tag[ele.id] = ElementType::REGULAR_INTERIOR_CUBE;
1415 continue;
1416 }
1417 // type 2
1418 bool has_singular_v = false;
1419 int n_irregular_v = 0;
1420 for (auto vid : ele.vs)
1421 {
1422 if (mesh_.vertices[vid].neighbor_hs.size() != 8)
1423 n_irregular_v++;
1424 int n_irregular_e = 0;
1425 for (auto eid : mesh_.vertices[vid].neighbor_es)
1426 {
1427 if (mesh_.edges[eid].neighbor_hs.size() != 4)
1428 n_irregular_e++;
1429 }
1430 if (n_irregular_e != 0 && n_irregular_e != 2)
1431 {
1432 has_singular_v = true;
1433 break;
1434 }
1435 }
1436 if (!has_singular_v && n_irregular_v == 2)
1437 {
1439 continue;
1440 }
1441
1443 }
1444 else
1445 {
1446 ele_tag[ele.id] = ElementType::INTERIOR_POLYTOPE;
1447 for (auto fid : ele.fs)
1448 if (mesh_.faces[fid].boundary)
1449 {
1450 ele_tag[ele.id] = ElementType::BOUNDARY_POLYTOPE;
1451 break;
1452 }
1453 }
1454 }
1455
1456 // TODO correct?
1457 for (auto &ele : mesh_.elements)
1458 {
1459 if (ele.vs.size() == 4)
1460 ele_tag[ele.id] = ElementType::SIMPLEX;
1461 else if (ele.vs.size() == 6)
1462 ele_tag[ele.id] = ElementType::PRISM;
1463 }
1464 }
1465
1466 double CMesh3D::quad_area(const int gid) const
1467 {
1468 const int n_vertices = n_face_vertices(gid);
1469 assert(n_vertices == 4);
1470
1471 const auto &vertices = mesh_.faces[gid].vs;
1472
1473 const auto v1 = point(vertices[0]);
1474 const auto v2 = point(vertices[1]);
1475 const auto v3 = point(vertices[2]);
1476 const auto v4 = point(vertices[3]);
1477
1478 const Eigen::Vector3d e0 = (v2 - v1).transpose();
1479 const Eigen::Vector3d e1 = (v3 - v1).transpose();
1480
1481 const Eigen::Vector3d e2 = (v2 - v4).transpose();
1482 const Eigen::Vector3d e3 = (v3 - v4).transpose();
1483
1484 return e0.cross(e1).norm() / 2 + e2.cross(e3).norm() / 2;
1485 }
1486
1488 {
1489 const int v0 = mesh_.edges[e].vs[0];
1490 const int v1 = mesh_.edges[e].vs[1];
1491 return 0.5 * (point(v0) + point(v1));
1492 }
1493
1495 {
1496 const int n_vertices = n_face_vertices(f);
1497 RowVectorNd bary(3);
1498 bary.setZero();
1499
1500 const auto &vertices = mesh_.faces[f].vs;
1501 for (int lv = 0; lv < n_vertices; ++lv)
1502 {
1503 bary += point(vertices[lv]);
1504 }
1505 return bary / n_vertices;
1506 }
1507
1509 {
1510 const int n_vertices = n_cell_vertices(c);
1511 RowVectorNd bary(3);
1512 bary.setZero();
1513
1514 const auto &vertices = mesh_.elements[c].vs;
1515 for (int lv = 0; lv < n_vertices; ++lv)
1516 {
1517 bary += point(vertices[lv]);
1518 }
1519 return bary / n_vertices;
1520 }
1521
1523 {
1524 m.vertices.clear();
1525 m.edges.clear();
1526 m.faces.clear();
1527 m.vertices.resize(gm.vertices.nb());
1528 m.faces.resize(gm.facets.nb());
1529 for (uint32_t i = 0; i < m.vertices.size(); i++)
1530 {
1531 Vertex v;
1532 v.id = i;
1533 v.v.push_back(gm.vertices.point_ptr(i)[0]);
1534 v.v.push_back(gm.vertices.point_ptr(i)[1]);
1535 v.v.push_back(gm.vertices.point_ptr(i)[2]);
1536 m.vertices[i] = v;
1537 }
1538 m.points.resize(3, m.vertices.size());
1539 for (uint32_t i = 0; i < m.vertices.size(); i++)
1540 {
1541 m.points(0, i) = m.vertices[i].v[0];
1542 m.points(1, i) = m.vertices[i].v[1];
1543 m.points(2, i) = m.vertices[i].v[2];
1544 }
1545
1546 if (m.type == MeshType::TRI || m.type == MeshType::QUA || m.type == MeshType::H_SUR)
1547 {
1548 for (uint32_t i = 0; i < m.faces.size(); i++)
1549 {
1550 Face f;
1551 f.id = i;
1552 f.vs.resize(gm.facets.nb_vertices(i));
1553 for (uint32_t j = 0; j < f.vs.size(); j++)
1554 {
1555 f.vs[j] = gm.facets.vertex(i, j);
1556 }
1557 m.faces[i] = f;
1558 }
1560 }
1561 }
1562
1563 void CMesh3D::elements_boxes(std::vector<std::array<Eigen::Vector3d, 2>> &boxes) const
1564 {
1565 boxes.resize(n_elements());
1566
1567 for (int i = 0; i < n_elements(); ++i)
1568 {
1569 auto &box = boxes[i];
1570 box[0].setConstant(std::numeric_limits<double>::max());
1571 box[1].setConstant(std::numeric_limits<double>::min());
1572
1573 for (int j = 0; j < n_cell_vertices(i); ++j)
1574 {
1575 const int v_id = cell_vertex(i, j);
1576 for (int d = 0; d < 3; ++d)
1577 {
1578 box[0][d] = std::min(box[0][d], point(v_id)[d]);
1579 box[1][d] = std::max(box[1][d], point(v_id)[d]);
1580 }
1581 }
1582 }
1583 }
1584
1585 void CMesh3D::barycentric_coords(const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const
1586 {
1587 assert(is_simplex(el_id));
1588
1589 const auto indices = get_ordered_vertices_from_tet(el_id);
1590
1591 const auto A = point(indices[0]);
1592 const auto B = point(indices[1]);
1593 const auto C = point(indices[2]);
1594 const auto D = point(indices[3]);
1595
1596 igl::barycentric_coordinates(p, A, B, C, D, coord);
1597 }
1598
1599 void CMesh3D::append(const Mesh &mesh)
1600 {
1601 assert(typeid(mesh) == typeid(CMesh3D));
1602 Mesh::append(mesh);
1603
1604 const CMesh3D &mesh3d = dynamic_cast<const CMesh3D &>(mesh);
1605 mesh_.append(mesh3d.mesh_);
1606
1609 }
1610
1611 std::unique_ptr<Mesh> CMesh3D::copy() const
1612 {
1613 return std::make_unique<CMesh3D>(*this);
1614 }
1615
1616 } // namespace mesh
1617} // 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:1508
RowVectorNd edge_barycenter(const int e) const override
edge barycenter
Definition CMesh3D.cpp:1487
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:1522
RowVectorNd point(const int vertex_id) const override
point coordinates
Definition CMesh3D.cpp:1249
void normalize() override
normalize the mesh
Definition CMesh3D.cpp:1020
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:1208
RowVectorNd face_barycenter(const int f) const override
face barycenter
Definition CMesh3D.cpp:1494
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:1466
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:1255
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:1563
void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector< std::vector< int > > &nodes) override
attach high order nodes
Definition CMesh3D.cpp:510
void set_point(const int global_index, const RowVectorNd &p) override
Set the point.
Definition CMesh3D.cpp:1238
void bounding_box(RowVectorNd &min, RowVectorNd &max) const override
computes the bbox of the mesh
Definition CMesh3D.cpp:1013
void append(const Mesh &mesh) override
appends a new mesh to the end of this
Definition CMesh3D.cpp:1599
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:1611
void compute_elements_tag() override
compute element types, see ElementType
Definition CMesh3D.cpp:1262
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:1585
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:1187
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:1226
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:350
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:40
int n_elements() const
utitlity to return the number of elements, cells or faces in 3d and 2d
Definition Mesh.hpp:162
Eigen::MatrixXi orders_
list of geometry orders, one per cell
Definition Mesh.hpp:684
std::vector< ElementType > elements_tag_
list of element types
Definition Mesh.hpp:676
Eigen::MatrixXi in_ordered_faces_
Order of the input faces, TODO: change to std::vector of Eigen::Vector.
Definition Mesh.hpp:702
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:680
std::vector< CellNodes > cell_nodes_
high-order nodes associates to cells
Definition Mesh.hpp:693
std::vector< int > element_vertices(const int el_id) const
list of vids of an element
Definition Mesh.hpp:230
std::vector< int > body_ids_
list of volume labels
Definition Mesh.hpp:682
std::vector< FaceNodes > face_nodes_
high-order nodes associates to faces
Definition Mesh.hpp:691
const std::vector< ElementType > & elements_tag() const
Returns the elements types.
Definition Mesh.hpp:421
std::vector< EdgeNodes > edge_nodes_
high-order nodes associates to edges
Definition Mesh.hpp:689
std::vector< std::vector< int > > faces() const
list of sorted faces.
Definition Mesh.cpp:448
Eigen::MatrixXi in_ordered_edges_
Order of the input edges.
Definition Mesh.hpp:700
virtual void append(const Mesh &mesh)
appends a new mesh to the end of this
Definition Mesh.cpp:499
Eigen::VectorXi in_ordered_vertices_
Order of the input vertices.
Definition Mesh.hpp:698
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
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