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 // std::cout << M.cells.adjacent(c, lf) << std::endl;
331 // std::cout << c << ' ' << M.cells.facet(c2, lf) << std::endl;
332 if (c == (int)M.cells.adjacent(c2, lf))
333 {
334 return (int)M.cells.facet(c2, lf);
335 }
336 }
337 assert(false);
338 return -1;
339 };
340
341 std::vector<int> cell_facet_to_facet(M.cell_facets.nb(), -1);
342
343 // Creates 1 hex or polyhedral element for each cell of the input mesh
344 int facet_counter = 0;
345 mesh_.elements.resize(M.cells.nb());
346 bool is_hex = true;
347 for (int c = 0; c < (int)M.cells.nb(); ++c)
348 {
349 Element &cell = mesh_.elements[c];
350 cell.id = c;
351 cell.hex = (M.cells.type(c) == GEO::MESH_HEX);
352
353 is_hex = is_hex && cell.hex;
354
355 int nf = M.cells.nb_facets(c);
356 cell.fs.resize(nf);
357
358 for (int lf = 0; lf < nf; ++lf)
359 {
360 int cf = M.cells.facet(c, lf);
361 int cf2 = opposite_cell_facet(c, cf);
362 // std::cout << "cf2: " << cf2 << std::endl;
363 // std::cout << "face_counter: " << facet_counter << std::endl;
364 if (cf2 < 0 || cell_facet_to_facet[cf2] < 0)
365 {
366 mesh_.faces.emplace_back();
367 Face &face = mesh_.faces[facet_counter];
368 // std::cout << "fid: " << face.id << std::endl;
369 assert(face.vs.empty());
370 face.vs.resize(M.cells.facet_nb_vertices(c, lf));
371 for (int lv = 0; lv < (int)M.cells.facet_nb_vertices(c, lf); ++lv)
372 {
373 face.vs[lv] = M.cells.facet_vertex(c, lf, lv);
374 }
375 cell.fs_flag.push_back(0);
376 cell.fs[lf] = face.id = facet_counter;
377 cell_facet_to_facet[cf] = facet_counter;
378 ++facet_counter;
379 }
380 else
381 {
382 cell.fs[lf] = cell_facet_to_facet[cf2];
383 cell.fs_flag.push_back(1);
384 }
385 }
386
387 for (auto fid : cell.fs)
388 {
389 cell.vs.insert(cell.vs.end(), mesh_.faces[fid].vs.begin(), mesh_.faces[fid].vs.end());
390 }
391 sort(cell.vs.begin(), cell.vs.end());
392 cell.vs.erase(unique(cell.vs.begin(), cell.vs.end()), cell.vs.end());
393
394 // Compute a point in the kernel (assumes the barycenter is ok)
395 Eigen::RowVector3d p(0, 0, 0);
396 for (int v : cell.vs)
397 {
398 p += mesh_.points.col(v).transpose();
399 }
400 p /= cell.vs.size();
401 mesh_.elements[c].v_in_Kernel.push_back(p[0]);
402 mesh_.elements[c].v_in_Kernel.push_back(p[1]);
403 mesh_.elements[c].v_in_Kernel.push_back(p[2]);
404 }
405 mesh_.type = is_hex ? MeshType::HEX : (M.cells.are_simplices() ? MeshType::TET : MeshType::HYB);
406 }
407
409 // if (is_simplicial()) {
410 // MeshProcessing3D::orient_volume_mesh(mesh_);
411 // }
413 return true;
414 }
415
416 bool CMesh3D::save(const std::string &path) const
417 {
418
419 if (!StringUtils::endswith(path, ".HYBRID"))
420 {
421 GEO::Mesh M;
422 to_geogram_mesh(*this, M);
423 GEO::mesh_save(M, path);
424 return true;
425 }
426
427 std::fstream f(path, std::ios::out);
428
429 f << mesh_.points.cols() << " " << mesh_.faces.size() << " " << 3 * mesh_.elements.size() << std::endl;
430 for (int i = 0; i < mesh_.points.cols(); i++)
431 f << mesh_.points(0, i) << " " << mesh_.points(1, i) << " " << mesh_.points(2, i) << std::endl;
432
433 for (auto f_ : mesh_.faces)
434 {
435 f << f_.vs.size() << " ";
436 for (auto vid : f_.vs)
437 f << vid << " ";
438 f << std::endl;
439 }
440
441 for (uint32_t i = 0; i < mesh_.elements.size(); i++)
442 {
443 f << mesh_.elements[i].fs.size() << " ";
444 for (auto fid : mesh_.elements[i].fs)
445 f << fid << " ";
446 f << std::endl;
447 f << mesh_.elements[i].fs_flag.size() << " ";
448 for (auto f_flag : mesh_.elements[i].fs_flag)
449 f << f_flag << " ";
450 f << std::endl;
451 }
452
453 for (uint32_t i = 0; i < mesh_.elements.size(); i++)
454 {
455 f << mesh_.elements[i].hex << std::endl;
456 }
457
458 f << "KERNEL"
459 << " " << mesh_.elements.size() << std::endl;
460 for (uint32_t i = 0; i < mesh_.elements.size(); i++)
461 {
462 f << mesh_.elements[i].v_in_Kernel[0] << " " << mesh_.elements[i].v_in_Kernel[1] << " " << mesh_.elements[i].v_in_Kernel[2] << std::endl;
463 }
464 f.close();
465
466 return true;
467 }
468
469 bool CMesh3D::build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
470 {
471 assert(F.cols() == 4 || F.cols() == 8);
472 edge_nodes_.clear();
473 face_nodes_.clear();
474 cell_nodes_.clear();
475
476 GEO::Mesh M;
477 M.vertices.create_vertices((int)V.rows());
478 for (int i = 0; i < (int)M.vertices.nb(); ++i)
479 {
480 GEO::vec3 &p = M.vertices.point(i);
481 p[0] = V(i, 0);
482 p[1] = V(i, 1);
483 p[2] = V(i, 2);
484 }
485
486 if (F.cols() == 4)
487 M.cells.create_tets((int)F.rows());
488 else if (F.cols() == 8)
489 M.cells.create_hexes((int)F.rows());
490 else
491 {
492 throw std::runtime_error("Mesh format not supported");
493 }
494
495 static const std::vector<int> permute_tet = {0, 1, 2, 3};
496 // polyfem uses the msh file format for hexes ordering!
497 static const std::vector<int> permute_hex = {1, 0, 2, 3, 5, 4, 6, 7};
498 const std::vector<int> permute = F.cols() == 4 ? permute_tet : permute_hex;
499
500 for (int c = 0; c < (int)M.cells.nb(); ++c)
501 {
502 for (int lv = 0; lv < F.cols(); ++lv)
503 {
504 M.cells.set_vertex(c, lv, F(c, permute[lv]));
505 }
506 }
507 M.cells.connect();
508
509 return load(M);
510 }
511
512 void CMesh3D::attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector<std::vector<int>> &nodes)
513 {
514 edge_nodes_.clear();
515 face_nodes_.clear();
516 cell_nodes_.clear();
517
518 edge_nodes_.resize(n_edges());
519 face_nodes_.resize(n_faces());
520 cell_nodes_.resize(n_cells());
521
522 orders_.resize(n_cells(), 1);
523
524 const auto attach_p2 = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
525 auto &n = edge_nodes_[index.edge];
526
527 if (n.nodes.size() > 0)
528 return;
529
530 n.v1 = index.vertex;
531 n.v2 = switch_vertex(index).vertex;
532
533 const int n_v1 = index.vertex;
534 const int n_v2 = switch_vertex(index).vertex;
535
536 int node_index = 0;
537
538 if ((n_v1 == nodes_ids[0] && n_v2 == nodes_ids[1]) || (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[1]))
539 node_index = 4;
540 else if ((n_v1 == nodes_ids[1] && n_v2 == nodes_ids[2]) || (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[2]))
541 node_index = 5;
542 else if ((n_v1 == nodes_ids[2] && n_v2 == nodes_ids[3]) || (n_v2 == nodes_ids[2] && n_v1 == nodes_ids[3]))
543 node_index = 8;
544
545 else if ((n_v1 == nodes_ids[0] && n_v2 == nodes_ids[3]) || (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[3]))
546 node_index = 7;
547 else if ((n_v1 == nodes_ids[0] && n_v2 == nodes_ids[2]) || (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[2]))
548 node_index = 6;
549 else
550 node_index = 9;
551
552 n.nodes.resize(1, 3);
553 n.nodes << V(nodes_ids[node_index], 0), V(nodes_ids[node_index], 1), V(nodes_ids[node_index], 2);
554 n.nodes_ids.push_back(nodes_ids[node_index]);
555 };
556
557 const auto attach_p3 = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
558 auto &n = edge_nodes_[index.edge];
559
560 if (n.nodes.size() > 0)
561 return;
562
563 n.v1 = index.vertex;
564 n.v2 = switch_vertex(index).vertex;
565
566 const int n_v1 = index.vertex;
567 const int n_v2 = switch_vertex(index).vertex;
568
569 int node_index1 = 0;
570 int node_index2 = 0;
571 if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[1])
572 {
573 node_index1 = 4;
574 node_index2 = 5;
575 }
576 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[1])
577 {
578 node_index1 = 5;
579 node_index2 = 4;
580 }
581 else if (n_v1 == nodes_ids[1] && n_v2 == nodes_ids[2])
582 {
583 node_index1 = 6;
584 node_index2 = 7;
585 }
586 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[2])
587 {
588 node_index1 = 7;
589 node_index2 = 6;
590 }
591 else if (n_v1 == nodes_ids[2] && n_v2 == nodes_ids[3])
592 {
593 node_index1 = 13;
594 node_index2 = 12;
595 }
596 else if (n_v2 == nodes_ids[2] && n_v1 == nodes_ids[3])
597 {
598 node_index1 = 12;
599 node_index2 = 13;
600 }
601
602 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[3])
603 {
604 node_index1 = 11;
605 node_index2 = 10;
606 }
607 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[3])
608 {
609 node_index1 = 10;
610 node_index2 = 11;
611 }
612 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[2])
613 {
614 node_index1 = 9;
615 node_index2 = 8;
616 }
617 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[2])
618 {
619 node_index1 = 8;
620 node_index2 = 9;
621 }
622
623 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[3])
624 {
625 node_index1 = 14;
626 node_index2 = 15;
627 }
628 else
629 {
630 node_index1 = 15;
631 node_index2 = 14;
632 }
633
634 n.nodes.resize(2, 3);
635 n.nodes.row(0) << V(nodes_ids[node_index1], 0), V(nodes_ids[node_index1], 1), V(nodes_ids[node_index1], 2);
636 n.nodes.row(1) << V(nodes_ids[node_index2], 0), V(nodes_ids[node_index2], 1), V(nodes_ids[node_index2], 2);
637 n.nodes_ids.push_back(nodes_ids[node_index1]);
638 n.nodes_ids.push_back(nodes_ids[node_index2]);
639 };
640
641 const auto attach_p3_face = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids, int id) {
642 auto &n = face_nodes_[index.face];
643 if (n.nodes.size() <= 0)
644 {
645 n.v1 = face_vertex(index.face, 0);
646 n.v2 = face_vertex(index.face, 1);
647 n.v3 = face_vertex(index.face, 2);
648 n.nodes.resize(1, 3);
649 n.nodes << V(nodes_ids[id], 0), V(nodes_ids[id], 1), V(nodes_ids[id], 2);
650 n.nodes_ids.push_back(nodes_ids[id]);
651 }
652 };
653
654 const auto attach_p4 = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
655 auto &n = edge_nodes_[index.edge];
656
657 if (n.nodes.size() > 0)
658 return;
659
660 n.v1 = index.vertex;
661 n.v2 = switch_vertex(index).vertex;
662
663 const int n_v1 = index.vertex;
664 const int n_v2 = switch_vertex(index).vertex;
665
666 int node_index1 = 0;
667 int node_index2 = 0;
668 int node_index3 = 0;
669 if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[1])
670 {
671 node_index1 = 4;
672 node_index2 = 5;
673 node_index3 = 6;
674 }
675 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[1])
676 {
677 node_index1 = 6;
678 node_index2 = 5;
679 node_index3 = 4;
680 }
681
682 else if (n_v1 == nodes_ids[1] && n_v2 == nodes_ids[2])
683 {
684 node_index1 = 7;
685 node_index2 = 8;
686 node_index3 = 9;
687 }
688 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[2])
689 {
690 node_index1 = 9;
691 node_index2 = 8;
692 node_index3 = 7;
693 }
694
695 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[2])
696 {
697 node_index1 = 10;
698 node_index2 = 11;
699 node_index3 = 12;
700 }
701 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[2])
702 {
703 node_index1 = 12;
704 node_index2 = 11;
705 node_index3 = 10;
706 }
707
708 else if (n_v2 == nodes_ids[0] && n_v1 == nodes_ids[3])
709 {
710 node_index1 = 13;
711 node_index2 = 14;
712 node_index3 = 15;
713 }
714 else if (n_v1 == nodes_ids[0] && n_v2 == nodes_ids[3])
715 {
716 node_index1 = 15;
717 node_index2 = 14;
718 node_index3 = 13;
719 }
720
721 else if (n_v2 == nodes_ids[2] && n_v1 == nodes_ids[3])
722 {
723 node_index1 = 16;
724 node_index2 = 17;
725 node_index3 = 18;
726 }
727 else if (n_v1 == nodes_ids[2] && n_v2 == nodes_ids[3])
728 {
729 node_index1 = 18;
730 node_index2 = 17;
731 node_index3 = 16;
732 }
733
734 else if (n_v2 == nodes_ids[1] && n_v1 == nodes_ids[3])
735 {
736 node_index1 = 19;
737 node_index2 = 20;
738 node_index3 = 21;
739 }
740 else if (n_v1 == nodes_ids[1] && n_v2 == nodes_ids[3])
741 {
742 node_index1 = 21;
743 node_index2 = 20;
744 node_index3 = 19;
745 }
746 else
747 {
748 assert(false);
749 }
750
751 n.nodes.resize(3, 3);
752 n.nodes.row(0) << V(nodes_ids[node_index1], 0), V(nodes_ids[node_index1], 1), V(nodes_ids[node_index1], 2);
753 n.nodes.row(1) << V(nodes_ids[node_index2], 0), V(nodes_ids[node_index2], 1), V(nodes_ids[node_index2], 2);
754 n.nodes.row(2) << V(nodes_ids[node_index3], 0), V(nodes_ids[node_index3], 1), V(nodes_ids[node_index3], 2);
755 n.nodes_ids.push_back(nodes_ids[node_index1]);
756 n.nodes_ids.push_back(nodes_ids[node_index2]);
757 n.nodes_ids.push_back(nodes_ids[node_index3]);
758 };
759
760 const auto attach_p4_face = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
761 auto &n = face_nodes_[index.face];
762 if (n.nodes.size() <= 0)
763 {
764 n.v1 = face_vertex(index.face, 0);
765 n.v2 = face_vertex(index.face, 1);
766 n.v3 = face_vertex(index.face, 2);
767
768 std::array<int, 3> vid = {{n.v1, n.v2, n.v3}};
769 std::sort(vid.begin(), vid.end());
770
771 std::array<int, 3> c1 = {{nodes_ids[0], nodes_ids[1], nodes_ids[2]}}; // 22
772 std::array<int, 3> c2 = {{nodes_ids[0], nodes_ids[1], nodes_ids[3]}}; // 25
773 std::array<int, 3> c3 = {{nodes_ids[0], nodes_ids[2], nodes_ids[3]}}; // 28
774 std::array<int, 3> c4 = {{nodes_ids[1], nodes_ids[2], nodes_ids[3]}}; // 31
775
776 std::sort(c1.begin(), c1.end());
777 std::sort(c2.begin(), c2.end());
778 std::sort(c3.begin(), c3.end());
779 std::sort(c4.begin(), c4.end());
780
781 int id = 0;
782 int index0 = 0;
783 int index1 = 1;
784 int index2 = 2;
785 if (vid == c1)
786 {
787 id = 22;
788 n.v1 = nodes_ids[0];
789 n.v2 = nodes_ids[1];
790 n.v3 = nodes_ids[2];
791
792 index0 = 0;
793 index1 = 2;
794 index2 = 1; // ok
795 }
796 else if (vid == c2)
797 {
798 id = 25;
799 n.v1 = nodes_ids[0];
800 n.v2 = nodes_ids[1];
801 n.v3 = nodes_ids[3];
802
803 index0 = 0;
804 index1 = 1;
805 index2 = 2; // ok
806 }
807 else if (vid == c3)
808 {
809 id = 28;
810 n.v1 = nodes_ids[0];
811 n.v2 = nodes_ids[2];
812 n.v3 = nodes_ids[3];
813
814 index0 = 0;
815 index1 = 2;
816 index2 = 1; // ok
817 }
818 else if (vid == c4)
819 {
820 id = 31;
821 n.v1 = nodes_ids[1];
822 n.v2 = nodes_ids[2];
823 n.v3 = nodes_ids[3];
824
825 index0 = 1;
826 index1 = 2;
827 index2 = 0; // ok
828 }
829 else
830 {
831 // the face nees to be one of the 4 above
832 assert(false);
833 }
834
835 n.nodes.resize(3, 3);
836 assert(id + index0 < nodes_ids.size());
837 assert(id + index1 < nodes_ids.size());
838 assert(id + index2 < nodes_ids.size());
839 n.nodes.row(0) << V(nodes_ids[id + index0], 0), V(nodes_ids[id + index0], 1), V(nodes_ids[id + index0], 2);
840 n.nodes.row(1) << V(nodes_ids[id + index1], 0), V(nodes_ids[id + index1], 1), V(nodes_ids[id + index1], 2);
841 n.nodes.row(2) << V(nodes_ids[id + index2], 0), V(nodes_ids[id + index2], 1), V(nodes_ids[id + index2], 2);
842 n.nodes_ids.push_back(nodes_ids[id + index0]);
843 n.nodes_ids.push_back(nodes_ids[id + index1]);
844 n.nodes_ids.push_back(nodes_ids[id + index2]);
845 }
846 };
847
848 const auto attach_p4_cell = [&](const Navigation3D::Index &index, const std::vector<int> &nodes_ids) {
849 auto &n = cell_nodes_[index.element];
850 assert(nodes_ids.size() == 35);
851 assert(n.nodes.size() == 0);
852
853 if (n.nodes.size() <= 0)
854 {
855 n.v1 = cell_vertex(index.element, 0);
856 n.v2 = cell_vertex(index.element, 1);
857 n.v3 = cell_vertex(index.element, 2);
858 n.v4 = cell_vertex(index.element, 3);
859 n.nodes.resize(1, 3);
860
861 n.nodes << V(nodes_ids[34], 0), V(nodes_ids[34], 1), V(nodes_ids[34], 2);
862 n.nodes_ids.push_back(nodes_ids[34]);
863 }
864 };
865
866 assert(nodes.size() == n_cells());
867
868 for (int c = 0; c < n_cells(); ++c)
869 {
870 auto index = get_index_from_element(c);
871
872 const auto &nodes_ids = nodes[c];
873
874 if (nodes_ids.size() == 4)
875 {
876 orders_(c) = 1;
877 continue;
878 }
879 // P2
880 else if (nodes_ids.size() == 10)
881 {
882 orders_(c) = 2;
883
884 for (int le = 0; le < 3; ++le)
885 {
886 attach_p2(index, nodes_ids);
887 index = next_around_face(index);
888 }
889
890 index = switch_vertex(switch_edge(switch_face(index)));
891 attach_p2(index, nodes_ids);
892
893 index = switch_edge(index);
894 attach_p2(index, nodes_ids);
895
896 index = switch_edge(switch_face(index));
897 attach_p2(index, nodes_ids);
898 }
899 // P3
900 else if (nodes_ids.size() == 20)
901 {
902 orders_(c) = 3;
903
904 for (int le = 0; le < 3; ++le)
905 {
906 attach_p3(index, nodes_ids);
907 index = next_around_face(index);
908 }
909
910 {
911 index = switch_vertex(switch_edge(switch_face(index)));
912 attach_p3(index, nodes_ids);
913
914 index = switch_edge(index);
915 attach_p3(index, nodes_ids);
916
917 index = switch_edge(switch_face(index));
918 attach_p3(index, nodes_ids);
919 }
920
921 {
922 index = get_index_from_element(c);
923 std::array<int, 4> indices;
924
925 {
926 std::array<int, 3> f16 = {{nodes_ids[0], nodes_ids[1], nodes_ids[2]}};
927 std::array<int, 3> f17 = {{nodes_ids[3], nodes_ids[1], nodes_ids[0]}};
928 std::array<int, 3> f18 = {{nodes_ids[0], nodes_ids[2], nodes_ids[3]}};
929 std::array<int, 3> f19 = {{nodes_ids[1], nodes_ids[2], nodes_ids[3]}};
930 std::sort(f16.begin(), f16.end());
931 std::sort(f17.begin(), f17.end());
932 std::sort(f18.begin(), f18.end());
933 std::sort(f19.begin(), f19.end());
934
935 auto tmp = index;
936 std::array<int, 3> f0 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
937 tmp = switch_face(index);
938 std::array<int, 3> f1 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
939 tmp = switch_face(next_around_face(index));
940 std::array<int, 3> f2 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
942 std::array<int, 3> f3 = {{face_vertex(tmp.face, 0), face_vertex(tmp.face, 1), face_vertex(tmp.face, 2)}};
943
944 std::sort(f0.begin(), f0.end());
945 std::sort(f1.begin(), f1.end());
946 std::sort(f2.begin(), f2.end());
947 std::sort(f3.begin(), f3.end());
948
949 const std::array<std::array<int, 3>, 4> faces = {{f0, f1, f2, f3}};
950 const std::array<std::array<int, 3>, 4> nodes = {{f16, f17, f18, f19}};
951 for (int i = 0; i < 4; ++i)
952 {
953 const auto &f = faces[i];
954 bool found = false;
955 for (int j = 0; j < 4; ++j)
956 {
957 if (nodes[j] == f)
958 {
959 indices[i] = j + 16;
960 found = true;
961 break;
962 }
963 }
964
965 assert(found);
966 }
967 }
968
969 attach_p3_face(index, nodes_ids, indices[0]);
970 attach_p3_face(switch_face(index), nodes_ids, indices[1]);
971 attach_p3_face(switch_face(next_around_face(index)), nodes_ids, indices[2]);
972 attach_p3_face(switch_face(next_around_face(next_around_face(index))), nodes_ids, indices[3]);
973 }
974 }
975 // P4
976 else if (nodes_ids.size() == 35)
977 {
978 orders_(c) = 4;
979 for (int le = 0; le < 3; ++le)
980 {
981 attach_p4(index, nodes_ids);
982 index = next_around_face(index);
983 }
984
985 {
986 index = switch_vertex(switch_edge(switch_face(index)));
987 attach_p4(index, nodes_ids);
988
989 index = switch_edge(index);
990 attach_p4(index, nodes_ids);
991
992 index = switch_edge(switch_face(index));
993 attach_p4(index, nodes_ids);
994 }
995
996 {
997 index = get_index_from_element(c);
998
999 attach_p4_face(index, nodes_ids);
1000 attach_p4_face(switch_face(index), nodes_ids);
1001 attach_p4_face(switch_face(next_around_face(index)), nodes_ids);
1002 attach_p4_face(switch_face(next_around_face(next_around_face(index))), nodes_ids);
1003 }
1004
1005 attach_p4_cell(get_index_from_element(c), nodes_ids);
1006 }
1007 // unsupported
1008 else
1009 {
1010 assert(false);
1011 }
1012 }
1013 }
1014
1016 {
1017 const auto &V = mesh_.points;
1018 min = V.rowwise().minCoeff().transpose();
1019 max = V.rowwise().maxCoeff().transpose();
1020 }
1021
1023 {
1024 RowVectorNd minV, maxV;
1025 bounding_box(minV, maxV);
1026 auto &V = mesh_.points;
1027 const auto shift = V.rowwise().minCoeff().eval();
1028 const double scaling = 1.0 / (V.rowwise().maxCoeff() - V.rowwise().minCoeff()).maxCoeff();
1029 V = (V.colwise() - shift) * scaling;
1030
1031 for (int i = 0; i < n_cells(); ++i)
1032 {
1033 for (int d = 0; d < 3; ++d)
1034 {
1035 auto val = mesh_.elements[i].v_in_Kernel[d];
1036 mesh_.elements[i].v_in_Kernel[d] = (val - shift(d)) * scaling;
1037 }
1038 }
1039
1040 for (auto &n : edge_nodes_)
1041 {
1042 if (n.nodes.size() > 0)
1043 n.nodes = (n.nodes.rowwise() - shift.transpose()) * scaling;
1044 }
1045 for (auto &n : face_nodes_)
1046 {
1047 if (n.nodes.size() > 0)
1048 n.nodes = (n.nodes.rowwise() - shift.transpose()) * scaling;
1049 }
1050 for (auto &n : cell_nodes_)
1051 {
1052 if (n.nodes.size() > 0)
1053 n.nodes = (n.nodes.rowwise() - shift.transpose()) * scaling;
1054 }
1055
1056 logger().debug("-- bbox before normalization:");
1057 logger().debug(" min : {}", minV);
1058 logger().debug(" max : {}", maxV);
1059 logger().debug(" extent: {}", maxV - minV);
1060 bounding_box(minV, maxV);
1061 logger().debug("-- bbox after normalization:");
1062 logger().debug(" min : {}", minV);
1063 logger().debug(" max : {}", maxV);
1064 logger().debug(" extent: {}", maxV - minV);
1065
1066 // V.row(1) /= 100.;
1067
1068 // for(int i = 0; i < n_cells(); ++i)
1069 // {
1070 // mesh_.elements[i].v_in_Kernel[1] /= 100.;
1071 // }
1072
1073 Eigen::MatrixXd p0, p1, p;
1074 get_edges(p0, p1);
1075 p = p0 - p1;
1076 logger().debug("-- edge length after normalization:");
1077 logger().debug(" min: ", p.rowwise().norm().minCoeff());
1078 logger().debug(" max: ", p.rowwise().norm().maxCoeff());
1079 logger().debug(" avg: ", p.rowwise().norm().mean());
1080 }
1081
1082 // void CMesh3D::triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector<int> &ranges) const
1083 // {
1084 // ranges.clear();
1085
1086 // std::vector<Eigen::MatrixXi> local_tris(mesh_.elements.size());
1087 // std::vector<Eigen::MatrixXd> local_pts(mesh_.elements.size());
1088 // Eigen::MatrixXi tets;
1089
1090 // int total_tris = 0;
1091 // int total_pts = 0;
1092
1093 // ranges.push_back(0);
1094
1095 // Eigen::MatrixXd face_barys;
1096 // face_barycenters(face_barys);
1097
1098 // Eigen::MatrixXd cell_barys;
1099 // cell_barycenters(cell_barys);
1100
1101 // for (std::size_t e = 0; e < mesh_.elements.size(); ++e)
1102 // {
1103 // const Element &el = mesh_.elements[e];
1104
1105 // const int n_vertices = el.vs.size();
1106 // const int n_faces = el.fs.size();
1107
1108 // Eigen::MatrixXd local_pt(n_vertices + n_faces, 3);
1109
1110 // std::map<int, int> global_to_local;
1111
1112 // for (int i = 0; i < n_vertices; ++i)
1113 // {
1114 // const int global_index = el.vs[i];
1115 // local_pt.row(i) = mesh_.points.col(global_index).transpose();
1116 // global_to_local[global_index] = i;
1117 // }
1118
1119 // int n_local_faces = 0;
1120 // for (int i = 0; i < n_faces; ++i)
1121 // {
1122 // const Face &f = mesh_.faces[el.fs[i]];
1123 // n_local_faces += f.vs.size();
1124
1125 // local_pt.row(n_vertices + i) = face_barys.row(f.id); // node_from_face(f.id);
1126 // }
1127
1128 // Eigen::MatrixXi local_faces(n_local_faces, 3);
1129
1130 // int face_index = 0;
1131 // for (int i = 0; i < n_faces; ++i)
1132 // {
1133 // const Face &f = mesh_.faces[el.fs[i]];
1134 // const int n_face_vertices = f.vs.size();
1135
1136 // const Eigen::RowVector3d e0 = (point(f.vs[0]) - local_pt.row(n_vertices + i));
1137 // const Eigen::RowVector3d e1 = (point(f.vs[1]) - local_pt.row(n_vertices + i));
1138 // const Eigen::RowVector3d normal = e0.cross(e1);
1139 // // const Eigen::RowVector3d check_dir = (node_from_element(e)-p);
1140 // const Eigen::RowVector3d check_dir = (cell_barys.row(e) - point(f.vs[1]));
1141
1142 // const bool reverse_order = normal.dot(check_dir) > 0;
1143
1144 // for (int j = 0; j < n_face_vertices; ++j)
1145 // {
1146 // const int jp = (j + 1) % n_face_vertices;
1147 // if (reverse_order)
1148 // {
1149 // local_faces(face_index, 0) = global_to_local[f.vs[jp]];
1150 // local_faces(face_index, 1) = global_to_local[f.vs[j]];
1151 // }
1152 // else
1153 // {
1154 // local_faces(face_index, 0) = global_to_local[f.vs[j]];
1155 // local_faces(face_index, 1) = global_to_local[f.vs[jp]];
1156 // }
1157 // local_faces(face_index, 2) = n_vertices + i;
1158
1159 // ++face_index;
1160 // }
1161 // }
1162
1163 // local_pts[e] = local_pt;
1164 // local_tris[e] = local_faces;
1165
1166 // total_tris += local_tris[e].rows();
1167 // total_pts += local_pts[e].rows();
1168
1169 // ranges.push_back(total_tris);
1170
1171 // assert(local_pts[e].rows() == local_pt.rows());
1172 // }
1173
1174 // tris.resize(total_tris, 3);
1175 // pts.resize(total_pts, 3);
1176
1177 // int tri_index = 0;
1178 // int pts_index = 0;
1179 // for (std::size_t i = 0; i < local_tris.size(); ++i)
1180 // {
1181 // tris.block(tri_index, 0, local_tris[i].rows(), local_tris[i].cols()) = local_tris[i].array() + pts_index;
1182 // tri_index += local_tris[i].rows();
1183
1184 // pts.block(pts_index, 0, local_pts[i].rows(), local_pts[i].cols()) = local_pts[i];
1185 // pts_index += local_pts[i].rows();
1186 // }
1187 // }
1188
1189 bool CMesh3D::is_boundary_element(const int element_global_id) const
1190 {
1191 const auto &fs = mesh_.elements[element_global_id].fs;
1192
1193 for (auto f_id : fs)
1194 {
1195 if (is_boundary_face(f_id))
1196 return true;
1197 }
1198
1199 const auto &vs = mesh_.elements[element_global_id].vs;
1200
1201 for (auto v_id : vs)
1202 {
1203 if (is_boundary_vertex(v_id))
1204 return true;
1205 }
1206
1207 return false;
1208 }
1209
1210 void CMesh3D::compute_boundary_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &marker)
1211 {
1212 boundary_ids_.resize(n_faces());
1213
1214 for (int f = 0; f < n_faces(); ++f)
1215 {
1216 const bool is_boundary = is_boundary_face(f);
1217 std::vector<int> vs(n_face_vertices(f));
1218 for (int vid = 0; vid < vs.size(); ++vid)
1219 vs[vid] = face_vertex(f, vid);
1220
1221 const auto p = face_barycenter(f);
1222
1223 std::sort(vs.begin(), vs.end());
1224 boundary_ids_[f] = marker(f, vs, p, is_boundary);
1225 }
1226 }
1227
1228 void CMesh3D::compute_body_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &)> &marker)
1229 {
1230 body_ids_.resize(n_elements());
1231 std::fill(body_ids_.begin(), body_ids_.end(), -1);
1232
1233 for (int e = 0; e < n_elements(); ++e)
1234 {
1235 const auto bary = cell_barycenter(e);
1236 body_ids_[e] = marker(e, element_vertices(e), bary);
1237 }
1238 }
1239
1240 void CMesh3D::set_point(const int global_index, const RowVectorNd &p)
1241 {
1242 mesh_.points.col(global_index) = p.transpose();
1243 if (mesh_.vertices[global_index].v.size() == 3)
1244 {
1245 mesh_.vertices[global_index].v[0] = p[0];
1246 mesh_.vertices[global_index].v[1] = p[1];
1247 mesh_.vertices[global_index].v[2] = p[2];
1248 }
1249 }
1250
1251 RowVectorNd CMesh3D::point(const int global_index) const
1252 {
1253 RowVectorNd pt = mesh_.points.col(global_index).transpose();
1254 return pt;
1255 }
1256
1257 RowVectorNd CMesh3D::kernel(const int c) const
1258 {
1259 RowVectorNd pt(3);
1260 pt << mesh_.elements[c].v_in_Kernel[0], mesh_.elements[c].v_in_Kernel[1], mesh_.elements[c].v_in_Kernel[2];
1261 return pt;
1262 }
1263
1265 {
1266 std::vector<ElementType> &ele_tag = elements_tag_;
1267 ele_tag.clear();
1268
1269 ele_tag.resize(mesh_.elements.size());
1270 for (auto &t : ele_tag)
1272
1273 // boundary flags
1274 std::vector<bool> bv_flag(mesh_.vertices.size(), false), be_flag(mesh_.edges.size(), false), bf_flag(mesh_.faces.size(), false);
1275 for (auto f : mesh_.faces)
1276 if (f.boundary)
1277 bf_flag[f.id] = true;
1278 else
1279 {
1280 for (auto nhid : f.neighbor_hs)
1281 if (!mesh_.elements[nhid].hex)
1282 bf_flag[f.id] = true;
1283 }
1284 for (uint32_t i = 0; i < mesh_.faces.size(); ++i)
1285 if (bf_flag[i])
1286 for (uint32_t j = 0; j < mesh_.faces[i].vs.size(); ++j)
1287 {
1288 uint32_t eid = mesh_.faces[i].es[j];
1289 be_flag[eid] = true;
1290 bv_flag[mesh_.faces[i].vs[j]] = true;
1291 }
1292
1293 for (auto &ele : mesh_.elements)
1294 {
1295 if (ele.hex)
1296 {
1297 bool attaching_non_hex = false, on_boundary = false;
1298 ;
1299 for (auto vid : ele.vs)
1300 {
1301 for (auto eleid : mesh_.vertices[vid].neighbor_hs)
1302 if (!mesh_.elements[eleid].hex)
1303 {
1304 attaching_non_hex = true;
1305 break;
1306 }
1307 if (mesh_.vertices[vid].boundary)
1308 {
1309 on_boundary = true;
1310 break;
1311 }
1312 if (on_boundary || attaching_non_hex)
1313 break;
1314 }
1315 if (attaching_non_hex)
1316 {
1317 ele_tag[ele.id] = ElementType::INTERFACE_CUBE;
1318 continue;
1319 }
1320
1321 if (on_boundary)
1322 {
1324 // has no boundary edge--> singular
1325 bool boundary_edge = false, boundary_edge_singular = false, interior_edge_singular = false;
1326 int n_interior_edge_singular = 0;
1327 for (auto eid : ele.es)
1328 {
1329 int en = 0;
1330 if (be_flag[eid])
1331 {
1332 boundary_edge = true;
1333 for (auto nhid : mesh_.edges[eid].neighbor_hs)
1334 if (mesh_.elements[nhid].hex)
1335 en++;
1336 if (en > 2)
1337 boundary_edge_singular = true;
1338 }
1339 else
1340 {
1341 for (auto nhid : mesh_.edges[eid].neighbor_hs)
1342 if (mesh_.elements[nhid].hex)
1343 en++;
1344 if (en != 4)
1345 {
1346 interior_edge_singular = true;
1347 n_interior_edge_singular++;
1348 }
1349 }
1350 }
1351 if (!boundary_edge || boundary_edge_singular || n_interior_edge_singular > 1)
1352 continue;
1353
1354 bool has_singular_v = false, has_iregular_v = false;
1355 int n_in_irregular_v = 0;
1356 for (auto vid : ele.vs)
1357 {
1358 int vn = 0;
1359 if (bv_flag[vid])
1360 {
1361 int nh = 0;
1362 for (auto nhid : mesh_.vertices[vid].neighbor_hs)
1363 if (mesh_.elements[nhid].hex)
1364 nh++;
1365 if (nh > 4)
1366 has_iregular_v = true;
1367 continue; // not sure the conditions
1368 }
1369 else
1370 {
1371 if (mesh_.vertices[vid].neighbor_hs.size() != 8)
1372 n_in_irregular_v++;
1373 int n_irregular_e = 0;
1374 for (auto eid : mesh_.vertices[vid].neighbor_es)
1375 {
1376 if (mesh_.edges[eid].neighbor_hs.size() != 4)
1377 n_irregular_e++;
1378 }
1379 if (n_irregular_e != 0 && n_irregular_e != 2)
1380 {
1381 has_singular_v = true;
1382 break;
1383 }
1384 }
1385 }
1386 int n_irregular_e = 0;
1387 for (auto eid : ele.es)
1388 if (!be_flag[eid] && mesh_.edges[eid].neighbor_hs.size() != 4)
1389 n_irregular_e++;
1390 if (has_singular_v)
1391 continue;
1392 if (!has_singular_v)
1393 {
1394 if (n_irregular_e == 1)
1395 {
1397 }
1398 else if (n_irregular_e == 0 && n_in_irregular_v == 0 && !has_iregular_v)
1399 ele_tag[ele.id] = ElementType::REGULAR_BOUNDARY_CUBE;
1400 else
1401 continue;
1402 }
1403 continue;
1404 }
1405
1406 // type 1
1407 bool has_irregular_v = false;
1408 for (auto vid : ele.vs)
1409 if (mesh_.vertices[vid].neighbor_hs.size() != 8)
1410 {
1411 has_irregular_v = true;
1412 break;
1413 }
1414 if (!has_irregular_v)
1415 {
1416 ele_tag[ele.id] = ElementType::REGULAR_INTERIOR_CUBE;
1417 continue;
1418 }
1419 // type 2
1420 bool has_singular_v = false;
1421 int n_irregular_v = 0;
1422 for (auto vid : ele.vs)
1423 {
1424 if (mesh_.vertices[vid].neighbor_hs.size() != 8)
1425 n_irregular_v++;
1426 int n_irregular_e = 0;
1427 for (auto eid : mesh_.vertices[vid].neighbor_es)
1428 {
1429 if (mesh_.edges[eid].neighbor_hs.size() != 4)
1430 n_irregular_e++;
1431 }
1432 if (n_irregular_e != 0 && n_irregular_e != 2)
1433 {
1434 has_singular_v = true;
1435 break;
1436 }
1437 }
1438 if (!has_singular_v && n_irregular_v == 2)
1439 {
1441 continue;
1442 }
1443
1445 }
1446 else
1447 {
1448 ele_tag[ele.id] = ElementType::INTERIOR_POLYTOPE;
1449 for (auto fid : ele.fs)
1450 if (mesh_.faces[fid].boundary)
1451 {
1452 ele_tag[ele.id] = ElementType::BOUNDARY_POLYTOPE;
1453 break;
1454 }
1455 }
1456 }
1457
1458 // TODO correct?
1459 for (auto &ele : mesh_.elements)
1460 {
1461 if (ele.vs.size() == 4)
1462 ele_tag[ele.id] = ElementType::SIMPLEX;
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:416
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:1251
void normalize() override
normalize the mesh
Definition CMesh3D.cpp:1022
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:1210
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:1257
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:512
void set_point(const int global_index, const RowVectorNd &p) override
Set the point.
Definition CMesh3D.cpp:1240
void bounding_box(RowVectorNd &min, RowVectorNd &max) const override
computes the bbox of the mesh
Definition CMesh3D.cpp:1015
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:1264
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:469
bool is_boundary_element(const int element_global_id) const override
is cell boundary
Definition CMesh3D.cpp:1189
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:1228
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:346
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:39
int n_elements() const
utitlity to return the number of elements, cells or faces in 3d and 2d
Definition Mesh.hpp:161
Eigen::MatrixXi orders_
list of geometry orders, one per cell
Definition Mesh.hpp:664
std::vector< ElementType > elements_tag_
list of element types
Definition Mesh.hpp:656
Eigen::MatrixXi in_ordered_faces_
Order of the input faces, TODO: change to std::vector of Eigen::Vector.
Definition Mesh.hpp:682
bool is_simplex(const int el_id) const
checks if element is simples compatible
Definition Mesh.cpp:422
std::vector< int > boundary_ids_
list of surface labels
Definition Mesh.hpp:660
std::vector< CellNodes > cell_nodes_
high-order nodes associates to cells
Definition Mesh.hpp:673
std::vector< int > element_vertices(const int el_id) const
list of vids of an element
Definition Mesh.hpp:229
std::vector< int > body_ids_
list of volume labels
Definition Mesh.hpp:662
std::vector< FaceNodes > face_nodes_
high-order nodes associates to faces
Definition Mesh.hpp:671
const std::vector< ElementType > & elements_tag() const
Returns the elements types.
Definition Mesh.hpp:415
std::vector< EdgeNodes > edge_nodes_
high-order nodes associates to edges
Definition Mesh.hpp:669
std::vector< std::vector< int > > faces() const
list of sorted faces.
Definition Mesh.cpp:443
Eigen::MatrixXi in_ordered_edges_
Order of the input edges.
Definition Mesh.hpp:680
virtual void append(const Mesh &mesh)
appends a new mesh to the end of this
Definition Mesh.cpp:494
Eigen::VectorXi in_ordered_vertices_
Order of the input vertices.
Definition Mesh.hpp:678
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,...
@ 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:42
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