PolyFEM
Loading...
Searching...
No Matches
LagrangeBasis3d.cpp
Go to the documentation of this file.
1
3#include "LagrangeBasis3d.hpp"
4
9
11
15
17
18#include <cassert>
19#include <array>
21
22using namespace polyfem;
23using namespace polyfem::assembler;
24using namespace polyfem::basis;
25using namespace polyfem::mesh;
26using namespace polyfem::quadrature;
27
28/*
29Axes:
30 z y
31 │╱
32 o──x
33Boundaries:
34X axis: left/right
35Y axis: front/back
36Z axis: bottom/top
37Corner nodes:
38v3──x──v2
39 │⋱ ╱ ╲
40 x ⋱╱ ╲
41 │ x x x
42 │ ╱ ⋱ ╲
43 │╱ ⋱╲
44v0─────x──────v1
45v0 = (0, 0, 0)
46v1 = (1, 0, 0)
47v2 = (0, 1, 0)
48v3 = (0, 0, 1)
49Edge nodes:
50e0 = (0.5, 0, 0)
51e1 = (0.5, 0.5, 0)
52e2 = ( 0, 0.5, 0)
53e3 = ( 0, 0, 0.5)
54e4 = (0.5, 0, 0.5)
55e5 = ( 0, 0.5, 0.5)
56Corner nodes:
57 v7──────x─────v6
58 ╱┆ ╱ ╱│
59 ╱ ┆ ╱ ╱ │
60 x┄┄┼┄┄┄x┄┄┄┄┄┄x │
61 ╱ x ╱ ╱┆ x
62 ╱ ┆ ╱ ╱ ┆ ╱│
63v4─────┼x─────v5 ┆╱ │
64 │ ┆┆ │ x │
65 │ v3┼┄┄┄┄┄x┼┄⌿┼┄v2
66 │ ╱ ┆ │╱ ┆ ╱
67 x┄┄┄⌿┄┄x┄┄┄┄┄┄x ┆╱
68 │ x ┆ │ x
69 │ ╱ ┆ │ ╱
70 │╱ ┆ │╱
71v0──────x─────v1
72v0 = (0, 0, 0)
73v1 = (1, 0, 0)
74v2 = (1, 1, 0)
75v3 = (0, 1, 0)
76v4 = (0, 0, 1)
77v5 = (1, 0, 1)
78v6 = (1, 1, 1)
79v7 = (0, 1, 1)
80Edge nodes:
81 x─────e10─────x
82 ╱┆ ╱ ╱│
83 ╱ ┆ ╱ ╱ │
84 e11┄┼┄┄┄x┄┄┄┄┄e9 │
85 ╱ e7 ╱ ╱┆ e6
86 ╱ ┆ ╱ ╱ ┆ ╱│
87 x─────e8──────x ┆╱ │
88 │ ┆┆ │ x │
89 │ x┼┄┄┄┄┄e2┄⌿┼┄┄x
90 │ ╱ ┆ │╱ ┆ ╱
91e4┄┄┄⌿┄┄x┄┄┄┄┄e5 ┆╱
92 │ e3 ┆ │ e1
93 │ ╱ ┆ │ ╱
94 │╱ ┆ │╱
95 x─────e0──────x
96e0 = (0.5, 0, 0)
97e1 = ( 1, 0.5, 0)
98e2 = (0.5, 1, 0)
99e3 = ( 0, 0.5, 0)
100e4 = ( 0, 0, 0.5)
101e5 = ( 1, 0, 0.5)
102e6 = ( 1, 1, 0.5)
103e7 = ( 0, 1, 0.5)
104e8 = (0.5, 0, 1)
105e9 = ( 1, 0.5, 1)
106e10 = (0.5, 1, 1)
107e11 = ( 0, 0.5, 1)
108Face nodes:
109 v7──────x─────v6
110 ╱┆ ╱ ╱│
111 ╱ ┆ ╱ ╱ │
112 x┄┄┼┄┄f5┄┄┄┄┄┄x │
113 ╱ x ╱ f3 ╱┆ x
114 ╱ ┆ ╱ ╱ ┆ ╱│
115v4─────┼x─────v5 ┆╱ │
116 │ f0 ┆┆ │ f1 │
117 │ v3┼┄┄┄┄┄x┼┄⌿┼┄v2
118 │ ╱ ┆ │╱ ┆ ╱
119 x┄┄┄⌿┄f2┄┄┄┄┄┄x ┆╱
120 │ x ┆ f4 │ x
121 │ ╱ ┆ │ ╱
122 │╱ ┆ │╱
123v0──────x─────v1
124f0 = ( 0, 0.5, 0.5)
125f1 = ( 1, 0.5, 0.5)
126f2 = (0.5, 0, 0.5)
127f3 = (0.5, 1, 0.5)
128f4 = (0.5, 0.5, 0)
129f5 = (0.5, 0.5, 1)
130*/
131
132namespace
133{
134 template <class InputIterator, class T>
135 int find_index(InputIterator first, InputIterator last, const T &val)
136 {
137 return std::distance(first, std::find(first, last, val));
138 }
139
140 Navigation3D::Index find_quad_face(const Mesh3D &mesh, int c, int v1, int v2, int v3, int v4)
141 {
142 std::array<int, 4> v = {{v1, v2, v3, v4}};
143 std::sort(v.begin(), v.end());
144 for (int lf = 0; lf < mesh.n_cell_faces(c); ++lf)
145 {
146 auto idx = mesh.get_index_from_element(c, lf, 0);
147 if (mesh.n_face_vertices(idx.face) != 4)
148 continue;
149
150 std::array<int, 4> u;
151 for (int lv = 0; lv < mesh.n_face_vertices(idx.face); ++lv)
152 {
153 u[lv] = idx.vertex;
154 idx = mesh.next_around_face(idx);
155 }
156 std::sort(u.begin(), u.end());
157 if (u == v)
158 {
159 return idx;
160 }
161 }
162 assert(false);
163 return Navigation3D::Index();
164 }
165
166 std::array<int, 4> tet_vertices_local_to_global(const Mesh3D &mesh, int c)
167 {
168 // Vertex nodes
169 assert(mesh.is_simplex(c));
170 std::array<int, 4> l2g;
171 int lv = 0;
172 for (int vi : mesh.get_ordered_vertices_from_tet(c))
173 {
174 l2g[lv++] = vi;
175 }
176
177 return l2g;
178 }
179
180 std::array<int, 8> hex_vertices_local_to_global(const Mesh3D &mesh, int c)
181 {
182 assert(mesh.is_cube(c));
183
184 // Vertex nodes
185 std::array<int, 8> l2g;
186 int lv = 0;
187 for (int vi : mesh.get_ordered_vertices_from_hex(c))
188 {
189 l2g[lv++] = vi;
190 }
191
192 return l2g;
193 }
194
195 std::array<int, 6> prism_vertices_local_to_global(const Mesh3D &mesh, int c)
196 {
197 assert(mesh.is_prism(c));
198
199 // Vertex nodes
200 std::array<int, 6> l2g;
201 int lv = 0;
202 for (int vi : mesh.get_ordered_vertices_from_prism(c))
203 {
204 l2g[lv++] = vi;
205 }
206
207 return l2g;
208 }
209
210 int lowest_order_elem_on_edge(const polyfem::mesh::NCMesh3D &mesh, const Eigen::VectorXi &discr_orders, const int eid)
211 {
212 auto elem_list = mesh.edge_neighs(eid);
213 int min = std::numeric_limits<int>::max();
214 int elem = -1;
215 for (const auto e : elem_list)
216 if (discr_orders[e] < min)
217 elem = e;
218 return elem;
219 }
220
221 void tet_local_to_global(const bool is_geom_bases, const int p, const Mesh3D &mesh, int c, const Eigen::VectorXi &discr_order, const Eigen::VectorXi &edge_orders, const Eigen::VectorXi &face_orders, std::vector<int> &res, polyfem::mesh::MeshNodes &nodes, std::vector<std::vector<int>> &edge_virtual_nodes, std::vector<std::vector<int>> &face_virtual_nodes)
222 {
223 const int n_edge_nodes = p > 1 ? ((p - 1) * 6) : 0;
224 const int nn = p > 2 ? (p - 2) : 0;
225 const int n_loc_f = (nn * (nn + 1) / 2);
226 const int n_face_nodes = n_loc_f * 4;
227 int n_cell_nodes = 0;
228 for (int pp = 4; pp <= p; ++pp)
229 n_cell_nodes += ((pp - 3) * ((pp - 3) + 1) / 2);
230
231 if (p == 0)
232 {
233 res.push_back(nodes.node_id_from_cell(c));
234 return;
235 }
236
237 // std::vector<int> res;
238 res.reserve(4 + n_edge_nodes + n_face_nodes + n_cell_nodes);
239
240 // Edge nodes
241 Eigen::Matrix<Navigation3D::Index, 4, 1> f;
242
243 auto v = tet_vertices_local_to_global(mesh, c);
244 Eigen::Matrix<int, 4, 3> fv;
245 fv.row(0) << v[0], v[1], v[2];
246 fv.row(1) << v[0], v[1], v[3];
247 fv.row(2) << v[1], v[2], v[3];
248 fv.row(3) << v[2], v[0], v[3];
249
250 for (long lf = 0; lf < fv.rows(); ++lf)
251 {
252 const auto index = mesh.get_index_from_element_face(c, fv(lf, 0), fv(lf, 1), fv(lf, 2));
253 f[lf] = index;
254 }
255
256 Eigen::Matrix<Navigation3D::Index, 6, 1> e;
257 Eigen::Matrix<int, 6, 2> ev;
258 ev.row(0) << v[0], v[1];
259 ev.row(1) << v[1], v[2];
260 ev.row(2) << v[2], v[0];
261
262 ev.row(3) << v[0], v[3];
263 ev.row(4) << v[1], v[3];
264 ev.row(5) << v[2], v[3];
265
266 for (int le = 0; le < e.rows(); ++le)
267 {
268 // const auto index = find_edge(mesh, c, ev(le, 0), ev(le, 1));
269 const auto index = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
270 e[le] = index;
271 }
272
273 // vertices
274 for (size_t lv = 0; lv < v.size(); ++lv)
275 {
276 if (!mesh.is_conforming() && !is_geom_bases)
277 {
278 const auto &ncmesh = dynamic_cast<const NCMesh3D &>(mesh);
279 // hanging vertex
280 if (ncmesh.leader_edge_of_vertex(v[lv]) >= 0 || ncmesh.leader_face_of_vertex(v[lv]) >= 0)
281 res.push_back(-lv - 1);
282 else
283 res.push_back(nodes.node_id_from_primitive(v[lv]));
284 }
285 else
286 res.push_back(nodes.node_id_from_primitive(v[lv]));
287 }
288
289 // Edges
290 for (int le = 0; le < e.rows(); ++le)
291 {
292 const auto index = e[le];
293 auto neighs = mesh.edge_neighs(index.edge);
294 int min_p = discr_order.size() > 0 ? discr_order(c) : 0;
295
296 if (is_geom_bases)
297 {
298 auto node_ids = nodes.node_ids_from_edge(index, p - 1);
299 res.insert(res.end(), node_ids.begin(), node_ids.end());
300 }
301 else
302 {
303 if (!mesh.is_conforming() && !is_geom_bases)
304 {
305 const auto &ncmesh = dynamic_cast<const NCMesh3D &>(mesh);
306 // slave edge
307 if (ncmesh.leader_edge_of_edge(index.edge) >= 0 || ncmesh.leader_face_of_edge(index.edge) >= 0)
308 {
309 for (int tmp = 0; tmp < p - 1; ++tmp)
310 res.push_back(-le - 1);
311 }
312 // master or conforming edge with constrained order
313 else if (edge_orders[index.edge] < discr_order(c))
314 {
315 for (int tmp = 0; tmp < p - 1; ++tmp)
316 res.push_back(-le - 1);
317
318 int min_order_elem = lowest_order_elem_on_edge(ncmesh, discr_order, index.edge);
319 // master edge, add extra nodes
320 if (min_order_elem == c)
321 edge_virtual_nodes[index.edge] = nodes.node_ids_from_edge(index, edge_orders[index.edge] - 1);
322 }
323 else
324 {
325 auto node_ids = nodes.node_ids_from_edge(index, p - 1);
326 res.insert(res.end(), node_ids.begin(), node_ids.end());
327 }
328 }
329 else
330 {
331 for (auto cid : neighs)
332 {
333 min_p = std::min(min_p, discr_order.size() > 0 ? discr_order(cid) : 0);
334 }
335
336 if (discr_order.size() > 0 && discr_order(c) > min_p)
337 {
338 for (int tmp = 0; tmp < p - 1; ++tmp)
339 res.push_back(-le - 10);
340 }
341 else
342 {
343 auto node_ids = nodes.node_ids_from_edge(index, p - 1);
344 res.insert(res.end(), node_ids.begin(), node_ids.end());
345 }
346 }
347 }
348 }
349
350 // faces
351 for (int lf = 0; lf < f.rows(); ++lf)
352 {
353 const auto index = f[lf];
354 const auto other_cell = mesh.switch_element(index).element;
355
356 const bool skip_other = discr_order.size() > 0 && other_cell >= 0 && discr_order(c) > discr_order(other_cell);
357
358 if (is_geom_bases)
359 {
360 auto node_ids = nodes.node_ids_from_face(index, p - 2);
361 res.insert(res.end(), node_ids.begin(), node_ids.end());
362 }
363 else
364 {
365 if (!mesh.is_conforming() && !is_geom_bases)
366 {
367 const auto &ncmesh = dynamic_cast<const NCMesh3D &>(mesh);
368 // slave face
369 if (ncmesh.leader_face_of_face(index.face) >= 0)
370 {
371 for (int tmp = 0; tmp < n_loc_f; ++tmp)
372 res.push_back(-lf - 1);
373 }
374 // master face or conforming face with constrained order
375 else if (face_orders[index.face] < discr_order[c])
376 {
377 for (int tmp = 0; tmp < n_loc_f; ++tmp)
378 res.push_back(-lf - 1);
379 // master face
380 if (ncmesh.n_follower_faces(index.face) > 0 && face_orders[index.face] > 2)
381 face_virtual_nodes[index.face] = nodes.node_ids_from_face(index, face_orders[index.face] - 2);
382 }
383 else
384 {
385 auto node_ids = nodes.node_ids_from_face(index, p - 2);
386 res.insert(res.end(), node_ids.begin(), node_ids.end());
387 }
388 }
389 else
390 {
391 if (skip_other)
392 {
393 for (int tmp = 0; tmp < n_loc_f; ++tmp)
394 res.push_back(-lf - 1);
395 }
396 else
397 {
398 auto node_ids = nodes.node_ids_from_face(index, p - 2);
399 res.insert(res.end(), node_ids.begin(), node_ids.end());
400 }
401 }
402 }
403 }
404
405 // cells
406 if (n_cell_nodes > 0)
407 {
408 const auto index = f[0];
409
410 auto node_ids = nodes.node_ids_from_cell(index, p - 3);
411 res.insert(res.end(), node_ids.begin(), node_ids.end());
412 }
413
414 assert(res.size() == size_t(4 + n_edge_nodes + n_face_nodes + n_cell_nodes));
415 }
416
417 void hex_local_to_global(const bool serendipity, const int q, const Mesh3D &mesh, int c, const Eigen::VectorXi &discr_order, std::vector<int> &res, MeshNodes &nodes)
418 {
419 assert(mesh.is_cube(c));
420
421 const int n_edge_nodes = ((q - 1) * 12);
422 const int nn = (q - 1);
423 const int n_loc_f = serendipity ? 0 : (nn * nn);
424 const int n_face_nodes = serendipity ? 0 : (n_loc_f * 6);
425 const int n_cell_nodes = serendipity ? 0 : (nn * nn * nn);
426
427 if (q == 0)
428 {
429 res.push_back(nodes.node_id_from_cell(c));
430 return;
431 }
432
433 // std::vector<int> res;
434 res.reserve(8 + n_edge_nodes + n_face_nodes + n_cell_nodes);
435
436 // Vertex nodes
437 auto v = hex_vertices_local_to_global(mesh, c);
438
439 // Edge nodes
440 Eigen::Matrix<Navigation3D::Index, 12, 1> e;
441 Eigen::Matrix<int, 12, 2> ev;
442 ev.row(0) << v[0], v[1];
443 ev.row(1) << v[1], v[2];
444 ev.row(2) << v[2], v[3];
445 ev.row(3) << v[3], v[0];
446 ev.row(4) << v[0], v[4];
447 ev.row(5) << v[1], v[5];
448 ev.row(6) << v[2], v[6];
449 ev.row(7) << v[3], v[7];
450 ev.row(8) << v[4], v[5];
451 ev.row(9) << v[5], v[6];
452 ev.row(10) << v[6], v[7];
453 ev.row(11) << v[7], v[4];
454 for (int le = 0; le < e.rows(); ++le)
455 {
456 // e[le] = find_edge(mesh, c, ev(le, 0), ev(le, 1)).edge;
457 e[le] = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
458 }
459
460 // Face nodes
461 Eigen::Matrix<Navigation3D::Index, 6, 1> f;
462 Eigen::Matrix<int, 6, 4> fv;
463 fv.row(0) << v[0], v[3], v[4], v[7];
464 fv.row(1) << v[1], v[2], v[5], v[6];
465 fv.row(2) << v[0], v[1], v[5], v[4];
466 fv.row(3) << v[3], v[2], v[6], v[7];
467 fv.row(4) << v[0], v[1], v[2], v[3];
468 fv.row(5) << v[4], v[5], v[6], v[7];
469 for (int lf = 0; lf < f.rows(); ++lf)
470 {
471 const auto index = find_quad_face(mesh, c, fv(lf, 0), fv(lf, 1), fv(lf, 2), fv(lf, 3));
472 f[lf] = index;
473 }
474
475 // vertices
476 for (size_t lv = 0; lv < v.size(); ++lv)
477 {
478 res.push_back(nodes.node_id_from_primitive(v[lv]));
479 }
480 assert(res.size() == size_t(8));
481
482 // Edges
483 for (int le = 0; le < e.rows(); ++le)
484 {
485 const auto index = e[le];
486 auto neighs = mesh.edge_neighs(index.edge);
487 int min_q = discr_order.size() > 0 ? discr_order(c) : 0;
488
489 for (auto cid : neighs)
490 {
491 min_q = std::min(min_q, discr_order.size() > 0 ? discr_order(cid) : 0);
492 }
493
494 if (discr_order.size() > 0 && discr_order(c) > min_q)
495 {
496 for (int tmp = 0; tmp < q - 1; ++tmp)
497 res.push_back(-le - 10);
498 }
499 else
500 {
501 auto node_ids = nodes.node_ids_from_edge(index, q - 1);
502 res.insert(res.end(), node_ids.begin(), node_ids.end());
503 }
504 }
505 assert(res.size() == size_t(8 + n_edge_nodes));
506
507 // faces
508 for (int lf = 0; lf < f.rows(); ++lf)
509 {
510 const auto index = f[lf];
511 const auto other_cell = mesh.switch_element(index).element;
512
513 const bool skip_other = discr_order.size() > 0 && other_cell >= 0 && discr_order(c) > discr_order(other_cell);
514
515 if (skip_other)
516 {
517 for (int tmp = 0; tmp < n_loc_f; ++tmp)
518 res.push_back(-lf - 1);
519 }
520 else
521 {
522 auto node_ids = nodes.node_ids_from_face(index, serendipity ? 0 : (q - 1));
523 assert(node_ids.size() == n_loc_f);
524 res.insert(res.end(), node_ids.begin(), node_ids.end());
525 }
526 }
527 assert(res.size() == size_t(8 + n_edge_nodes + n_face_nodes));
528
529 // cells
530 if (n_cell_nodes > 0)
531 {
532 const auto index = f[0];
533
534 auto node_ids = nodes.node_ids_from_cell(index, q - 1);
535 res.insert(res.end(), node_ids.begin(), node_ids.end());
536 }
537
538 assert(res.size() == size_t(8 + n_edge_nodes + n_face_nodes + n_cell_nodes));
539 }
540
541 void prism_local_to_global(const int p, const int q, const Mesh3D &mesh, int c, const Eigen::VectorXi &discr_order, std::vector<int> &res, MeshNodes &nodes)
542 {
543 assert(mesh.is_prism(c));
544
545 const int n_edge_nodest = p > 1 ? ((p - 1) * 3) : 0;
546 const int nnt = p > 2 ? (p - 2) : 0;
547 const int n_face_nodest = nnt * (nnt + 1) / 2;
548
549 const int nnq = q > 1 ? (q - 1) : 0;
550 const int n_face_nodesq = nnq * n_edge_nodest / 3;
551
552 const int n_edge_nodes = n_edge_nodest * 2 + nnq * 3;
553 const int n_face_nodes = n_face_nodest * 2 + n_face_nodesq * 3;
554 const int n_cell_nodes = n_face_nodest * nnq;
555
556 if (p == 0 && q == 0)
557 {
558 res.push_back(nodes.node_id_from_cell(c));
559 return;
560 }
561
562 res.reserve(6 + n_edge_nodes + n_face_nodes + n_cell_nodes);
563
564 // Vertex nodes
565 auto v = prism_vertices_local_to_global(mesh, c);
566
567 // Edge nodes
568 Eigen::Matrix<Navigation3D::Index, 9, 1> e;
569 Eigen::Matrix<int, 9, 2> ev;
570 ev.row(0) << v[0], v[1];
571 ev.row(1) << v[1], v[2];
572 ev.row(2) << v[2], v[0];
573 ev.row(3) << v[3], v[4];
574 ev.row(4) << v[4], v[5];
575 ev.row(5) << v[5], v[3];
576 ev.row(6) << v[0], v[3];
577 ev.row(7) << v[1], v[4];
578 ev.row(8) << v[2], v[5];
579
580 for (int le = 0; le < e.rows(); ++le)
581 {
582 e[le] = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
583 }
584
585 // Face nodes
586 Eigen::Matrix<Navigation3D::Index, 5, 1> f;
587 Eigen::Matrix<int, 2, 3> fvt;
588 fvt.row(0) << v[0], v[1], v[2];
589 fvt.row(1) << v[3], v[4], v[5];
590 for (int lf = 0; lf < fvt.rows(); ++lf)
591 {
592 const auto index = mesh.get_index_from_element_face(c, fvt(lf, 0), fvt(lf, 1), fvt(lf, 2));
593 f[lf] = index;
594 }
595
596 Eigen::Matrix<int, 3, 4> fvq;
597 fvq.row(0) << v[0], v[1], v[4], v[3];
598 fvq.row(1) << v[1], v[2], v[5], v[4];
599 fvq.row(2) << v[2], v[0], v[3], v[5];
600 for (int lf = 0; lf < fvq.rows(); ++lf)
601 {
602 const auto index = find_quad_face(mesh, c, fvq(lf, 0), fvq(lf, 1), fvq(lf, 2), fvq(lf, 3));
603 f[lf + 2] = index;
604 }
605
606 // vertices
607 for (size_t lv = 0; lv < v.size(); ++lv)
608 {
609 res.push_back(nodes.node_id_from_primitive(v[lv]));
610 }
611 assert(res.size() == size_t(6));
612
613 // Edges
614 for (int le = 0; le < e.rows(); ++le)
615 {
616 const auto index = e[le];
617 auto neighs = mesh.edge_neighs(index.edge);
618 // int min_q = discr_order.size() > 0 ? discr_order(c) : 0;
619
620 // TODO prism: not conforming
621 // for (auto cid : neighs)
622 // {
623 // min_q = std::min(min_q, discr_order.size() > 0 ? discr_order(cid) : 0);
624 // }
625
626 // if (discr_order.size() > 0 && discr_order(c) > min_q)
627 // {
628 // for (int tmp = 0; tmp < q - 1; ++tmp)
629 // res.push_back(-le - 10);
630 // }
631 // else
632 {
633 auto node_ids = nodes.node_ids_from_edge(index, (le < 6 ? p : q) - 1);
634 res.insert(res.end(), node_ids.begin(), node_ids.end());
635 }
636 }
637 assert(res.size() == size_t(6 + n_edge_nodes));
638
639 // faces
640 for (int lf = 0; lf < f.rows(); ++lf)
641 {
642 const auto index = f[lf];
643 const auto other_cell = mesh.switch_element(index).element;
644
645 // TODO prism: not conforming
646 // const bool skip_other = discr_order.size() > 0 && other_cell >= 0 && discr_order(c) > discr_order(other_cell);
647
648 // if (skip_other)
649 // {
650 // for (int tmp = 0; tmp < n_loc_f; ++tmp)
651 // res.push_back(-lf - 1);
652 // }
653 // else
654 {
655 // todo prism, nodes are not necessarly a square
656 auto node_ids = nodes.node_ids_from_face(index, lf < 2 ? (p - 2) : ((q - 1) * (p - 1)));
657 // assert(node_ids.size() == n_loc_f);
658 res.insert(res.end(), node_ids.begin(), node_ids.end());
659 }
660 }
661 assert(res.size() == size_t(6 + n_edge_nodes + n_face_nodes));
662
663 // cells
664 if (n_cell_nodes > 0)
665 {
666 // TODO: implement me internal prism nodes
667 log_and_throw_error("Prism internal nodes not implemented yet");
668
669 const auto index = f[0];
670
671 auto node_ids = nodes.node_ids_from_cell(index, q - 1);
672 res.insert(res.end(), node_ids.begin(), node_ids.end());
673 }
674
675 assert(res.size() == size_t(6 + n_edge_nodes + n_face_nodes + n_cell_nodes));
676 }
677 // -----------------------------------------------------------------------------
678
692 void compute_nodes(
693 const Mesh3D &mesh,
694 const Eigen::VectorXi &discr_ordersp,
695 const Eigen::VectorXi &discr_ordersq,
696 const Eigen::VectorXi &edge_orders,
697 const Eigen::VectorXi &face_orders,
698 const bool serendipity,
699 const bool has_polys,
700 const bool is_geom_bases,
701 MeshNodes &nodes,
702 std::vector<std::vector<int>> &edge_virtual_nodes,
703 std::vector<std::vector<int>> &face_virtual_nodes,
704 std::vector<std::vector<int>> &element_nodes_id,
705 std::vector<LocalBoundary> &local_boundary,
706 std::map<int, InterfaceData> &poly_face_to_data)
707 {
708 // Step 1: Assign global node ids for each quads
709 local_boundary.clear();
710 // local_boundary.resize(mesh.n_faces());
711 element_nodes_id.resize(mesh.n_faces());
712
713 if (!mesh.is_conforming())
714 {
715 const auto &ncmesh = dynamic_cast<const NCMesh3D &>(mesh);
716 edge_virtual_nodes.resize(ncmesh.n_edges());
717 face_virtual_nodes.resize(ncmesh.n_faces());
718 }
719
720 for (int c = 0; c < mesh.n_cells(); ++c)
721 {
722 const int discr_order = discr_ordersp(c);
723 const int discr_orderq = discr_ordersq(c);
724
725 if (mesh.is_cube(c))
726 {
727 hex_local_to_global(serendipity, discr_order, mesh, c, discr_ordersp, element_nodes_id[c], nodes);
728
729 auto v = hex_vertices_local_to_global(mesh, c);
730 Eigen::Matrix<int, 6, 4> fv;
731 fv.row(0) << v[0], v[3], v[4], v[7];
732 fv.row(1) << v[1], v[2], v[5], v[6];
733 fv.row(2) << v[0], v[1], v[5], v[4];
734 fv.row(3) << v[3], v[2], v[6], v[7];
735 fv.row(4) << v[0], v[1], v[2], v[3];
736 fv.row(5) << v[4], v[5], v[6], v[7];
737
738 LocalBoundary lb(c, BoundaryType::QUAD);
739 for (int i = 0; i < fv.rows(); ++i)
740 {
741 const int f = find_quad_face(mesh, c, fv(i, 0), fv(i, 1), fv(i, 2), fv(i, 3)).face;
742
743 if (mesh.is_boundary_face(f) || mesh.get_boundary_id(f) > 0)
744 {
745 lb.add_boundary_primitive(f, i);
746 }
747 }
748
749 if (!lb.empty())
750 local_boundary.emplace_back(lb);
751 }
752 else if (mesh.is_simplex(c))
753 {
754 // element_nodes_id[c] = polyfem::LagrangeBasis3d::tet_local_to_global(discr_order, mesh, c, discr_orders, nodes);
755 tet_local_to_global(is_geom_bases, discr_order, mesh, c, discr_ordersp, edge_orders, face_orders, element_nodes_id[c], nodes, edge_virtual_nodes, face_virtual_nodes);
756
757 auto v = tet_vertices_local_to_global(mesh, c);
758 Eigen::Matrix<int, 4, 3> fv;
759 fv.row(0) << v[0], v[1], v[2];
760 fv.row(1) << v[0], v[1], v[3];
761 fv.row(2) << v[1], v[2], v[3];
762 fv.row(3) << v[2], v[0], v[3];
763
764 LocalBoundary lb(c, BoundaryType::TRI);
765 for (long i = 0; i < fv.rows(); ++i)
766 {
767 const int f = mesh.get_index_from_element_face(c, fv(i, 0), fv(i, 1), fv(i, 2)).face;
768
769 if (mesh.is_boundary_face(f))
770 {
771 lb.add_boundary_primitive(f, i);
772 }
773 }
774
775 if (!lb.empty())
776 local_boundary.emplace_back(lb);
777 }
778 else if (mesh.is_prism(c))
779 {
780 // todo non conforming prisms
781 prism_local_to_global(discr_order, discr_orderq, mesh, c, discr_ordersp, element_nodes_id[c], nodes);
782
783 auto v = prism_vertices_local_to_global(mesh, c);
784 Eigen::Matrix<int, 2, 3> fvt;
785 fvt.row(0) << v[0], v[1], v[2];
786 fvt.row(1) << v[3], v[4], v[5];
787
788 LocalBoundary lb(c, BoundaryType::PRISM);
789 for (long i = 0; i < fvt.rows(); ++i)
790 {
791 const int f = mesh.get_index_from_element_face(c, fvt(i, 0), fvt(i, 1), fvt(i, 2)).face;
792
793 if (mesh.is_boundary_face(f))
794 {
795 lb.add_boundary_primitive(f, i);
796 }
797 }
798
799 Eigen::Matrix<int, 3, 4> fvq;
800 fvq.row(0) << v[0], v[1], v[4], v[3];
801 fvq.row(1) << v[1], v[2], v[5], v[4];
802 fvq.row(2) << v[2], v[0], v[3], v[5];
803
804 for (long i = 0; i < fvq.rows(); ++i)
805 {
806 const int f = find_quad_face(mesh, c, fvq(i, 0), fvq(i, 1), fvq(i, 2), fvq(i, 3)).face;
807
808 if (mesh.is_boundary_face(f))
809 {
810 lb.add_boundary_primitive(f, i + 2);
811 }
812 }
813
814 if (!lb.empty())
815 local_boundary.emplace_back(lb);
816 }
817 }
818
819 if (!has_polys)
820 return;
821
822 // Step 2: Iterate over edges of polygons and compute interface weights
823 Eigen::VectorXi indices;
824 for (int c = 0; c < mesh.n_cells(); ++c)
825 {
826 // Skip non-polytopes
827 if (!mesh.is_polytope(c))
828 {
829 continue;
830 }
831
832 for (int lf = 0; lf < mesh.n_cell_faces(c); ++lf)
833 {
834 auto index = mesh.get_index_from_element(c, lf, 0);
835 auto index2 = mesh.switch_element(index);
836 int c2 = index2.element;
837 assert(c2 >= 0);
838
839 const int discr_order = discr_ordersp(c2);
840 const int discr_orderq = discr_ordersq(c2);
841 if (mesh.is_cube(c2))
842 {
843 indices = LagrangeBasis3d::hex_face_local_nodes(serendipity, discr_order, mesh, index2);
844 }
845 else if (mesh.is_simplex(c2))
846 {
847 indices = LagrangeBasis3d::tet_face_local_nodes(discr_order, mesh, index2);
848 }
849 else if (mesh.is_prism(c2))
850 {
851 indices = LagrangeBasis3d::prism_face_local_nodes(discr_order, discr_orderq, mesh, index2);
852 }
853 else
854 continue;
855
856 InterfaceData data;
857 data.local_indices.insert(data.local_indices.begin(), indices.data(), indices.data() + indices.size());
858 assert(indices.size() == data.local_indices.size());
859 poly_face_to_data[index2.face] = data;
860 }
861 }
862 }
869 void local_to_global(const Eigen::MatrixXd &verts, const Eigen::MatrixXd &uv, Eigen::MatrixXd &pts)
870 {
871 const int dim = verts.cols();
872 const int N = uv.rows();
873 assert(dim == 3);
874 assert(uv.cols() == dim);
875 assert(verts.rows() == dim + 1);
876
877 pts.setZero(N, dim);
878 for (int i = 0; i < N; i++)
879 pts.row(i) = uv(i, 0) * verts.row(1) + uv(i, 1) * verts.row(2) + uv(i, 2) * verts.row(3) + (1.0 - uv(i, 0) - uv(i, 1) - uv(i, 2)) * verts.row(0);
880 }
881
882 void local_to_global_face(const Eigen::MatrixXd &verts, const Eigen::MatrixXd &uv, Eigen::MatrixXd &pts)
883 {
884 const int dim = verts.cols();
885 const int N = uv.rows();
886 assert(dim == 3);
887 assert(uv.cols() == 2);
888 assert(verts.rows() == 3);
889
890 pts.setZero(N, dim);
891 for (int i = 0; i < N; i++)
892 pts.row(i) = uv(i, 0) * verts.row(1) + uv(i, 1) * verts.row(2) + (1.0 - uv(i, 0) - uv(i, 1)) * verts.row(0);
893 }
894
901 void global_to_local(const Eigen::MatrixXd &verts, const Eigen::MatrixXd &pts, Eigen::MatrixXd &uv)
902 {
903 const int dim = verts.cols();
904 const int N = pts.rows();
905 assert(dim == 3);
906 assert(verts.rows() == dim + 1);
907 assert(pts.cols() == dim);
908
909 Eigen::Matrix3d J;
910 for (int i = 0; i < dim; i++)
911 J.col(i) = verts.row(i + 1) - verts.row(0);
912
913 Eigen::Matrix3d Jinv = J.inverse();
914
915 uv.setZero(N, dim);
916 polyfem::utils::maybe_parallel_for(N, [&](int start, int end, int thread_id) {
917 for (int i = start; i < end; i++)
918 {
919 auto point = pts.row(i) - verts.row(0);
920 uv.row(i) = Jinv * point.transpose();
921 }
922 });
923 }
924
925 void global_to_local_face(const Eigen::MatrixXd &verts, const Eigen::MatrixXd &pts, Eigen::MatrixXd &uv)
926 {
927 const int dim = verts.cols();
928 const int N = pts.rows();
929 assert(dim == 3);
930 assert(verts.rows() == 3);
931 assert(pts.cols() == dim);
932
933 Eigen::Matrix3d J;
934 for (int i = 0; i < 2; i++)
935 J.col(i) = verts.row(i + 1) - verts.row(0);
936
937 Eigen::Vector3d a = J.col(0);
938 Eigen::Vector3d b = J.col(1);
939 Eigen::Vector3d virtual_vert = a.cross(b);
940 J.col(2) = virtual_vert;
941
942 uv.setZero(N, 2);
943 polyfem::utils::maybe_parallel_for(N, [&](int start, int end, int thread_id) {
944 for (int i = start; i < end; i++)
945 {
946 auto point = pts.row(i) - verts.row(0);
947 Eigen::Vector3d x = J.colPivHouseholderQr().solve(point.transpose());
948 uv.row(i) = x.block(0, 0, 2, 1);
949 assert(std::abs(x(2)) < 1e-8);
950 }
951 });
952 }
953
954 void global_to_local_edge(const Eigen::MatrixXd &verts, const Eigen::MatrixXd &pts, Eigen::VectorXd &uv)
955 {
956 const int dim = verts.cols();
957 const int N = pts.rows();
958 assert(dim == 3);
959 assert(verts.rows() == 2);
960 assert(pts.cols() == dim);
961
962 auto edge = verts.row(1) - verts.row(0);
963 double squared_length = edge.squaredNorm();
964
965 uv.setZero(N);
966 polyfem::utils::maybe_parallel_for(N, [&](int start, int end, int thread_id) {
967 for (int i = start; i < end; i++)
968 {
969 auto vec = pts.row(i) - verts.row(0);
970 uv(i) = (vec.dot(edge)) / squared_length;
971 }
972 });
973 }
974
975 bool check_edge_face_orders(const polyfem::mesh::NCMesh3D &mesh, const Eigen::VectorXi &elem_orders, const Eigen::VectorXi &edge_orders, const Eigen::VectorXi &face_orders)
976 {
977 // same order for overlapping faces
978 for (int i = 0; i < mesh.n_faces(); i++)
979 if (mesh.leader_face_of_face(i) >= 0)
980 if (face_orders[mesh.leader_face_of_face(i)] != face_orders[i])
981 return false;
982
983 // face order no smaller than order of its edges
984 for (int i = 0; i < mesh.n_faces(); i++)
985 {
986 if (mesh.n_face_cells(i) == 0)
987 continue;
988 for (int j = 0; j < mesh.n_face_vertices(i); j++)
989 {
990 const int e_id = mesh.face_edge(i, j);
991 if (edge_orders[e_id] > face_orders[i])
992 return false;
993 }
994 }
995
996 // same order for overlapping edges
997 for (int i = 0; i < mesh.n_edges(); i++)
998 if (mesh.leader_edge_of_edge(i) >= 0)
999 if (edge_orders[mesh.leader_edge_of_edge(i)] != edge_orders[i])
1000 return false;
1001
1002 // face order no larger than order of interior edges
1003 for (int i = 0; i < mesh.n_edges(); i++)
1004 {
1005 if (mesh.n_edge_cells(i) == 0)
1006 continue;
1007 if (mesh.leader_face_of_edge(i) >= 0 && mesh.leader_edge_of_edge(i) < 0)
1008 if (face_orders[mesh.leader_face_of_edge(i)] > edge_orders[i])
1009 return false;
1010 }
1011 return true;
1012 }
1013
1021 void compute_edge_face_orders(const polyfem::mesh::NCMesh3D &mesh, const Eigen::VectorXi &elem_orders, Eigen::VectorXi &edge_orders, Eigen::VectorXi &face_orders)
1022 {
1023 const int max_order = elem_orders.maxCoeff();
1024 edge_orders.setConstant(mesh.n_edges(), max_order);
1025 face_orders.setConstant(mesh.n_faces(), max_order);
1026
1027 for (int i = 0; i < mesh.n_cells(); i++)
1028 for (int j = 0; j < mesh.n_cell_faces(i); j++)
1029 face_orders[mesh.cell_face(i, j)] = std::min(face_orders[mesh.cell_face(i, j)], elem_orders[i]);
1030
1031 for (int i = 0; i < mesh.n_cells(); i++)
1032 for (int j = 0; j < mesh.n_cell_edges(i); j++)
1033 edge_orders[mesh.cell_edge(i, j)] = std::min(edge_orders[mesh.cell_edge(i, j)], elem_orders[i]);
1034
1035 while (!check_edge_face_orders(mesh, elem_orders, edge_orders, face_orders))
1036 {
1037 // same order for overlapping faces
1038 for (int i = 0; i < mesh.n_faces(); i++)
1039 if (mesh.leader_face_of_face(i) >= 0)
1040 face_orders[mesh.leader_face_of_face(i)] = std::min(face_orders[mesh.leader_face_of_face(i)], face_orders[i]);
1041
1042 for (int i = 0; i < mesh.n_faces(); i++)
1043 if (mesh.leader_face_of_face(i) >= 0)
1044 face_orders[i] = std::min(face_orders[mesh.leader_face_of_face(i)], face_orders[i]);
1045
1046 // face order no smaller than order of its edges
1047 for (int i = 0; i < mesh.n_faces(); i++)
1048 {
1049 if (mesh.n_face_cells(i) == 0)
1050 continue;
1051 for (int j = 0; j < mesh.n_face_vertices(i); j++)
1052 {
1053 const int e_id = mesh.face_edge(i, j);
1054 edge_orders[e_id] = std::min(edge_orders[e_id], face_orders[i]);
1055 }
1056 }
1057
1058 // same order for overlapping edges
1059 for (int i = 0; i < mesh.n_edges(); i++)
1060 if (mesh.leader_edge_of_edge(i) >= 0)
1061 edge_orders[mesh.leader_edge_of_edge(i)] = std::min(edge_orders[mesh.leader_edge_of_edge(i)], edge_orders[i]);
1062
1063 for (int i = 0; i < mesh.n_edges(); i++)
1064 if (mesh.leader_edge_of_edge(i) >= 0)
1065 edge_orders[i] = std::min(edge_orders[mesh.leader_edge_of_edge(i)], edge_orders[i]);
1066
1067 // face order no larger than order of interior edges
1068 for (int i = 0; i < mesh.n_edges(); i++)
1069 {
1070 if (mesh.n_edge_cells(i) == 0)
1071 continue;
1072 if (mesh.leader_face_of_edge(i) >= 0 && mesh.leader_edge_of_edge(i) < 0)
1073 face_orders[mesh.leader_face_of_edge(i)] = std::min(face_orders[mesh.leader_face_of_edge(i)], edge_orders[i]);
1074 }
1075 }
1076 }
1077} // anonymous namespace
1078
1079Eigen::VectorXi LagrangeBasis3d::tet_face_local_nodes(const int p, const Mesh3D &mesh, Navigation3D::Index index)
1080{
1081 const int nn = p > 2 ? (p - 2) : 0;
1082 const int n_edge_nodes = (p - 1) * 6;
1083 const int n_face_nodes = nn * (nn + 1) / 2;
1084
1085 const int c = index.element;
1086 assert(mesh.is_simplex(c));
1087
1088 // Local to global mapping of node indices
1089 const auto l2g = tet_vertices_local_to_global(mesh, c);
1090
1091 // Extract requested interface
1092 Eigen::VectorXi result(3 + (p - 1) * 3 + n_face_nodes);
1093 result[0] = find_index(l2g.begin(), l2g.end(), index.vertex);
1094 result[1] = find_index(l2g.begin(), l2g.end(), mesh.next_around_face(index).vertex);
1095 result[2] = find_index(l2g.begin(), l2g.end(), mesh.next_around_face(mesh.next_around_face(index)).vertex);
1096
1097 Eigen::Matrix<Navigation3D::Index, 6, 1> e;
1098 Eigen::Matrix<int, 6, 2> ev;
1099 ev.row(0) << l2g[0], l2g[1];
1100 ev.row(1) << l2g[1], l2g[2];
1101 ev.row(2) << l2g[2], l2g[0];
1102
1103 ev.row(3) << l2g[0], l2g[3];
1104 ev.row(4) << l2g[1], l2g[3];
1105 ev.row(5) << l2g[2], l2g[3];
1106
1107 Navigation3D::Index tmp = index;
1108
1109 for (int le = 0; le < e.rows(); ++le)
1110 {
1111 // const auto index = find_edge(mesh, c, ev(le, 0), ev(le, 1));
1112 const auto l_index = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
1113 e[le] = l_index;
1114 }
1115
1116 int ii = 3;
1117 for (int k = 0; k < 3; ++k)
1118 {
1119 bool reverse = false;
1120 int le = 0;
1121 for (; le < ev.rows(); ++le)
1122 {
1123 // const auto l_index = find_edge(mesh, c, ev(le, 0), ev(le, 1));
1124 // const auto l_index = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
1125 const auto l_index = e[le];
1126 if (l_index.edge == tmp.edge)
1127 {
1128 if (l_index.vertex == tmp.vertex)
1129 reverse = false;
1130 else
1131 {
1132 reverse = true;
1133 if (mesh.switch_vertex(tmp).vertex != l_index.vertex)
1134 assert(false);
1135 }
1136
1137 break;
1138 }
1139 }
1140 assert(le < 6);
1141
1142 if (!reverse)
1143 {
1144
1145 for (int i = 0; i < p - 1; ++i)
1146 {
1147 result[ii++] = 4 + le * (p - 1) + i;
1148 }
1149 }
1150 else
1151 {
1152 for (int i = 0; i < p - 1; ++i)
1153 {
1154 result[ii++] = 4 + (le + 1) * (p - 1) - i - 1;
1155 }
1156 }
1157
1158 tmp = mesh.next_around_face(tmp);
1159 }
1160
1161 // faces
1162
1163 Eigen::Matrix<int, 4, 3> fv;
1164 fv.row(0) << l2g[0], l2g[1], l2g[2];
1165 fv.row(1) << l2g[0], l2g[1], l2g[3];
1166 fv.row(2) << l2g[1], l2g[2], l2g[3];
1167 fv.row(3) << l2g[2], l2g[0], l2g[3];
1168
1169 long lf = 0;
1170 for (; lf < fv.rows(); ++lf)
1171 {
1172 const auto l_index = mesh.get_index_from_element_face(c, fv(lf, 0), fv(lf, 1), fv(lf, 2));
1173 if (l_index.face == index.face)
1174 break;
1175 }
1176
1177 assert(lf < fv.rows());
1178
1179 if (n_face_nodes == 0)
1180 {
1181 }
1182 else if (n_face_nodes == 1)
1183 result[ii++] = 4 + n_edge_nodes + lf;
1184 else // if (n_face_nodes == 3)
1185 {
1186
1187 const auto get_order = [&p, &nn, &n_face_nodes](const std::array<int, 3> &corners) {
1188 int index;
1189 int start;
1190 int offset;
1191
1192 std::vector<int> order1(n_face_nodes); // A-> B
1193 for (int k = 0; k < n_face_nodes; ++k)
1194 order1[k] = k;
1195
1196 std::vector<int> order2(n_face_nodes); // B-> A
1197 index = 0;
1198 start = nn - 1;
1199 for (int k = 0; k < nn; ++k)
1200 {
1201 for (int l = 0; l < nn - k; ++l)
1202 {
1203 order2[index] = start - l;
1204 index++;
1205 }
1206 start += (nn - 1) - k;
1207 }
1208
1209 std::vector<int> order3(n_face_nodes); // A->C
1210 index = 0;
1211 for (int k = 0; k < nn; ++k)
1212 {
1213 offset = k;
1214 for (int l = 0; l < nn - k; ++l)
1215 {
1216 order3[index] = offset;
1217 offset += nn - l;
1218
1219 index++;
1220 }
1221 }
1222
1223 std::vector<int> order4(n_face_nodes); // C-> A
1224 index = 0;
1225 start = n_face_nodes - 1;
1226 for (int k = 0; k < nn; ++k)
1227 {
1228 offset = 0;
1229 for (int l = 0; l < nn - k; ++l)
1230 {
1231 order4[index] = start - offset;
1232 offset += k + 2 + l;
1233 index++;
1234 }
1235
1236 start += -k - 1;
1237 }
1238
1239 std::vector<int> order5(n_face_nodes); // B-> C
1240 index = 0;
1241 start = nn - 1;
1242 for (int k = 0; k < nn; ++k)
1243 {
1244 offset = 0;
1245 for (int l = 0; l < nn - k; ++l)
1246 {
1247 order5[index] = start + offset;
1248 offset += nn - 1 - l;
1249 index++;
1250 }
1251
1252 start--;
1253 }
1254
1255 std::vector<int> order6(n_face_nodes); // C-> B
1256 index = 0;
1257 start = n_face_nodes;
1258 for (int k = 0; k < nn; ++k)
1259 {
1260 offset = 0;
1261 start = start - k - 1;
1262 for (int l = 0; l < nn - k; ++l)
1263 {
1264 order6[index] = start - offset;
1265 offset += l + 1 + k;
1266 index++;
1267 }
1268 }
1269
1270 if (corners[0] == order1[0] && corners[1] == order1[nn - 1])
1271 {
1272 assert(corners[2] == order1[n_face_nodes - 1]);
1273 return order1;
1274 }
1275
1276 if (corners[0] == order2[0] && corners[1] == order2[nn - 1])
1277 {
1278 assert(corners[2] == order2[n_face_nodes - 1]);
1279 return order2;
1280 }
1281
1282 if (corners[0] == order3[0] && corners[1] == order3[nn - 1])
1283 {
1284 assert(corners[2] == order3[n_face_nodes - 1]);
1285 return order3;
1286 }
1287
1288 if (corners[0] == order4[0] && corners[1] == order4[nn - 1])
1289 {
1290 assert(corners[2] == order4[n_face_nodes - 1]);
1291 return order4;
1292 }
1293
1294 if (corners[0] == order5[0] && corners[1] == order5[nn - 1])
1295 {
1296 assert(corners[2] == order5[n_face_nodes - 1]);
1297 return order5;
1298 }
1299
1300 if (corners[0] == order6[0] && corners[1] == order6[nn - 1])
1301 {
1302 assert(corners[2] == order6[n_face_nodes - 1]);
1303 return order6;
1304 }
1305
1306 assert(false);
1307 return order1;
1308 };
1309
1310 Eigen::MatrixXd nodes;
1311 autogen::p_nodes_3d(p, nodes);
1312 // auto pos = LagrangeBasis3d::linear_tet_face_local_nodes_coordinates(mesh, index);
1313 // Local to global mapping of node indices
1314
1315 // Extract requested interface
1316 std::array<int, 3> idx;
1317 for (int lv = 0; lv < 3; ++lv)
1318 {
1319 idx[lv] = find_index(l2g.begin(), l2g.end(), index.vertex);
1320 index = mesh.next_around_face(index);
1321 }
1322 Eigen::Matrix3d pos(3, 3);
1323 int cnt = 0;
1324 for (int i : idx)
1325 {
1326 pos.row(cnt++) = nodes.row(i);
1327 }
1328
1329 const Eigen::RowVector3d bary = pos.colwise().mean();
1330
1331 const int offset = 4 + n_edge_nodes;
1332 bool found = false;
1333 for (int lff = 0; lff < 4; ++lff)
1334 {
1335 Eigen::MatrixXd loc_nodes = nodes.block(offset + lff * n_face_nodes, 0, n_face_nodes, 3);
1336 Eigen::RowVector3d node_bary = loc_nodes.colwise().mean();
1337
1338 if ((node_bary - bary).norm() < 1e-10)
1339 {
1340 std::array<int, 3> corners;
1341 int sum = 0;
1342 for (int m = 0; m < 3; ++m)
1343 {
1344 auto t = pos.row(m);
1345 int min_n = -1;
1346 double min_dis = 10000;
1347
1348 for (int n = 0; n < n_face_nodes; ++n)
1349 {
1350 double dis = (loc_nodes.row(n) - t).squaredNorm();
1351 if (dis < min_dis)
1352 {
1353 min_dis = dis;
1354 min_n = n;
1355 }
1356 }
1357
1358 assert(min_n >= 0);
1359 assert(min_n < n_face_nodes);
1360 corners[m] = min_n;
1361 }
1362
1363 const auto indices = get_order(corners);
1364 for (int min_n : indices)
1365 {
1366 sum += min_n;
1367 result[ii++] = 4 + n_edge_nodes + min_n + lf * n_face_nodes;
1368 }
1369
1370 assert(sum == (n_face_nodes - 1) * n_face_nodes / 2);
1371
1372 found = true;
1373 assert(lff == lf);
1374
1375 break;
1376 }
1377 }
1378
1379 assert(found);
1380 }
1381 // else
1382 // {
1383 // assert(n_face_nodes == 0);
1384 // }
1385
1386 assert(ii == result.size());
1387 return result;
1388}
1389
1390Eigen::VectorXi LagrangeBasis3d::hex_face_local_nodes(const bool serendipity, const int q, const Mesh3D &mesh, Navigation3D::Index index)
1391{
1392 const int nn = q - 1;
1393 const int n_edge_nodes = nn * 12;
1394 const int n_face_nodes = serendipity ? 0 : nn * nn;
1395
1396 const int c = index.element;
1397 assert(mesh.is_cube(c));
1398
1399 // Local to global mapping of node indices
1400 const auto l2g = hex_vertices_local_to_global(mesh, c);
1401
1402 // Extract requested interface
1403 Eigen::VectorXi result(4 + nn * 4 + n_face_nodes);
1404 result[0] = find_index(l2g.begin(), l2g.end(), index.vertex);
1405 result[1] = find_index(l2g.begin(), l2g.end(), mesh.next_around_face(index).vertex);
1406 result[2] = find_index(l2g.begin(), l2g.end(), mesh.next_around_face(mesh.next_around_face(index)).vertex);
1407 result[3] = find_index(l2g.begin(), l2g.end(), mesh.next_around_face(mesh.next_around_face(mesh.next_around_face(index))).vertex);
1408
1409 Eigen::Matrix<Navigation3D::Index, 12, 1> e;
1410 Eigen::Matrix<int, 12, 2> ev;
1411 ev.row(0) << l2g[0], l2g[1];
1412 ev.row(1) << l2g[1], l2g[2];
1413 ev.row(2) << l2g[2], l2g[3];
1414 ev.row(3) << l2g[3], l2g[0];
1415 ev.row(4) << l2g[0], l2g[4];
1416 ev.row(5) << l2g[1], l2g[5];
1417 ev.row(6) << l2g[2], l2g[6];
1418 ev.row(7) << l2g[3], l2g[7];
1419 ev.row(8) << l2g[4], l2g[5];
1420 ev.row(9) << l2g[5], l2g[6];
1421 ev.row(10) << l2g[6], l2g[7];
1422 ev.row(11) << l2g[7], l2g[4];
1423
1424 Navigation3D::Index tmp = index;
1425
1426 for (int le = 0; le < e.rows(); ++le)
1427 {
1428 const auto l_index = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
1429 e[le] = l_index;
1430 }
1431
1432 int ii = 4;
1433 for (int k = 0; k < 4; ++k)
1434 {
1435 bool reverse = false;
1436 int le = 0;
1437 for (; le < ev.rows(); ++le)
1438 {
1439 // const auto l_index = find_edge(mesh, c, ev(le, 0), ev(le, 1));
1440 // const auto l_index = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
1441 const auto l_index = e[le];
1442 if (l_index.edge == tmp.edge)
1443 {
1444 if (l_index.vertex == tmp.vertex)
1445 reverse = false;
1446 else
1447 {
1448 reverse = true;
1449 assert(mesh.switch_vertex(tmp).vertex == l_index.vertex);
1450 }
1451
1452 break;
1453 }
1454 }
1455 assert(le < 12);
1456
1457 if (!reverse)
1458 {
1459
1460 for (int i = 0; i < q - 1; ++i)
1461 {
1462 result[ii++] = 8 + le * (q - 1) + i;
1463 }
1464 }
1465 else
1466 {
1467 for (int i = 0; i < q - 1; ++i)
1468 {
1469 result[ii++] = 8 + (le + 1) * (q - 1) - i - 1;
1470 }
1471 }
1472
1473 tmp = mesh.next_around_face(tmp);
1474 }
1475
1476 // faces
1477
1478 Eigen::Matrix<int, 6, 4> fv;
1479 fv.row(0) << l2g[0], l2g[3], l2g[4], l2g[7];
1480 fv.row(1) << l2g[1], l2g[2], l2g[5], l2g[6];
1481 fv.row(2) << l2g[0], l2g[1], l2g[5], l2g[4];
1482 fv.row(3) << l2g[3], l2g[2], l2g[6], l2g[7];
1483 fv.row(4) << l2g[0], l2g[1], l2g[2], l2g[3];
1484 fv.row(5) << l2g[4], l2g[5], l2g[6], l2g[7];
1485
1486 long lf = 0;
1487 for (; lf < fv.rows(); ++lf)
1488 {
1489 const auto l_index = find_quad_face(mesh, c, fv(lf, 0), fv(lf, 1), fv(lf, 2), fv(lf, 3));
1490 if (l_index.face == index.face)
1491 break;
1492 }
1493
1494 assert(lf < fv.rows());
1495
1496 if (n_face_nodes == 1)
1497 result[ii++] = 8 + n_edge_nodes + lf;
1498 else if (n_face_nodes != 0)
1499 {
1500 Eigen::MatrixXd nodes;
1501 autogen::q_nodes_3d(q, nodes);
1502 // auto pos = LagrangeBasis3d::linear_tet_face_local_nodes_coordinates(mesh, index);
1503 // Local to global mapping of node indices
1504
1505 // Extract requested interface
1506 std::array<int, 4> idx;
1507 for (int lv = 0; lv < 4; ++lv)
1508 {
1509 idx[lv] = find_index(l2g.begin(), l2g.end(), index.vertex);
1510 index = mesh.next_around_face(index);
1511 }
1512 Eigen::Matrix<double, 4, 3> pos(4, 3);
1513 int cnt = 0;
1514 for (int i : idx)
1515 {
1516 pos.row(cnt++) = nodes.row(i);
1517 }
1518
1519 const Eigen::RowVector3d bary = pos.colwise().mean();
1520
1521 const int offset = 8 + n_edge_nodes;
1522 bool found = false;
1523 for (int lff = 0; lff < 6; ++lff)
1524 {
1525 Eigen::Matrix<double, 4, 3> loc_nodes = nodes.block<4, 3>(offset + lff * n_face_nodes, 0);
1526 Eigen::RowVector3d node_bary = loc_nodes.colwise().mean();
1527
1528 if ((node_bary - bary).norm() < 1e-10)
1529 {
1530 int sum = 0;
1531 for (int m = 0; m < 4; ++m)
1532 {
1533 auto t = pos.row(m);
1534 int min_n = -1;
1535 double min_dis = 10000;
1536
1537 for (int n = 0; n < 4; ++n)
1538 {
1539 double dis = (loc_nodes.row(n) - t).squaredNorm();
1540 if (dis < min_dis)
1541 {
1542 min_dis = dis;
1543 min_n = n;
1544 }
1545 }
1546
1547 assert(min_n >= 0);
1548 assert(min_n < 4);
1549
1550 sum += min_n;
1551
1552 result[ii++] = 8 + n_edge_nodes + min_n + lf * n_face_nodes;
1553 }
1554
1555 assert(sum == 6); // 0 + 1 + 2 + 3
1556
1557 found = true;
1558 assert(lff == lf);
1559 }
1560
1561 if (found)
1562 break;
1563 }
1564
1565 assert(found);
1566 }
1567
1568 assert(ii == result.size());
1569 return result;
1570}
1571
1572Eigen::VectorXi LagrangeBasis3d::prism_face_local_nodes(const int p, const int q, const Mesh3D &mesh, Navigation3D::Index index)
1573{
1574 const int c = index.element;
1575 assert(mesh.is_prism(c));
1576
1577 // Local to global mapping of node indices
1578 const auto l2g = prism_vertices_local_to_global(mesh, c);
1579
1580 const int global_n_edges_nodes = (p - 1) * 6 + (q - 1) * 3;
1581
1582 if (mesh.n_face_vertices(index.face) == 3)
1583 {
1584 const int nn = p > 2 ? (p - 2) : 0;
1585 const int n_edge_nodes = (p - 1) * 6;
1586 const int n_face_nodes = nn * (nn + 1) / 2;
1587
1588 // Extract requested interface
1589 Eigen::VectorXi result(3 + (p - 1) * 3 + n_face_nodes);
1590 result[0] = find_index(l2g.begin(), l2g.end(), index.vertex);
1591 result[1] = find_index(l2g.begin(), l2g.end(), mesh.next_around_face(index).vertex);
1592 result[2] = find_index(l2g.begin(), l2g.end(), mesh.next_around_face(mesh.next_around_face(index)).vertex);
1593
1594 Eigen::Matrix<Navigation3D::Index, 6, 1> e;
1595 Eigen::Matrix<int, 6, 2> ev;
1596 ev.row(0) << l2g[0], l2g[1];
1597 ev.row(1) << l2g[1], l2g[2];
1598 ev.row(2) << l2g[2], l2g[0];
1599
1600 ev.row(3) << l2g[3], l2g[4];
1601 ev.row(4) << l2g[4], l2g[5];
1602 ev.row(5) << l2g[5], l2g[3];
1603
1604 Navigation3D::Index tmp = index;
1605
1606 for (int le = 0; le < e.rows(); ++le)
1607 {
1608 // const auto index = find_edge(mesh, c, ev(le, 0), ev(le, 1));
1609 const auto l_index = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
1610 e[le] = l_index;
1611 }
1612
1613 int ii = 3;
1614 for (int k = 0; k < 3; ++k)
1615 {
1616 bool reverse = false;
1617 int le = 0;
1618 for (; le < ev.rows(); ++le)
1619 {
1620 // const auto l_index = find_edge(mesh, c, ev(le, 0), ev(le, 1));
1621 // const auto l_index = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
1622 const auto l_index = e[le];
1623 if (l_index.edge == tmp.edge)
1624 {
1625 if (l_index.vertex == tmp.vertex)
1626 reverse = false;
1627 else
1628 {
1629 reverse = true;
1630 if (mesh.switch_vertex(tmp).vertex != l_index.vertex)
1631 assert(false);
1632 }
1633
1634 break;
1635 }
1636 }
1637 assert(le < 6);
1638
1639 if (!reverse)
1640 {
1641
1642 for (int i = 0; i < p - 1; ++i)
1643 {
1644 result[ii++] = 6 + le * (p - 1) + i;
1645 }
1646 }
1647 else
1648 {
1649 for (int i = 0; i < p - 1; ++i)
1650 {
1651 result[ii++] = 6 + (le + 1) * (p - 1) - i - 1;
1652 }
1653 }
1654
1655 tmp = mesh.next_around_face(tmp);
1656 }
1657
1658 // faces
1659
1660 Eigen::Matrix<int, 2, 3> fv;
1661 fv.row(0) << l2g[0], l2g[1], l2g[2];
1662 fv.row(1) << l2g[3], l2g[4], l2g[5];
1663
1664 long lf = 0;
1665 for (; lf < fv.rows(); ++lf)
1666 {
1667 const auto l_index = mesh.get_index_from_element_face(c, fv(lf, 0), fv(lf, 1), fv(lf, 2));
1668 if (l_index.face == index.face)
1669 break;
1670 }
1671
1672 assert(lf < fv.rows());
1673
1674 if (n_face_nodes == 0)
1675 {
1676 }
1677 else if (n_face_nodes == 1)
1678 result[ii++] = 6 + global_n_edges_nodes + lf;
1679 else // if (n_face_nodes == 3)
1680 {
1681
1682 const auto get_order = [&p, &q, &nn, &n_face_nodes](const std::array<int, 3> &corners) {
1683 int index;
1684 int start;
1685 int offset;
1686
1687 std::vector<int> order1(n_face_nodes); // A-> B
1688 for (int k = 0; k < n_face_nodes; ++k)
1689 order1[k] = k;
1690
1691 std::vector<int> order2(n_face_nodes); // B-> A
1692 index = 0;
1693 start = nn - 1;
1694 for (int k = 0; k < nn; ++k)
1695 {
1696 for (int l = 0; l < nn - k; ++l)
1697 {
1698 order2[index] = start - l;
1699 index++;
1700 }
1701 start += (nn - 1) - k;
1702 }
1703
1704 std::vector<int> order3(n_face_nodes); // A->C
1705 index = 0;
1706 for (int k = 0; k < nn; ++k)
1707 {
1708 offset = k;
1709 for (int l = 0; l < nn - k; ++l)
1710 {
1711 order3[index] = offset;
1712 offset += nn - l;
1713
1714 index++;
1715 }
1716 }
1717
1718 std::vector<int> order4(n_face_nodes); // C-> A
1719 index = 0;
1720 start = n_face_nodes - 1;
1721 for (int k = 0; k < nn; ++k)
1722 {
1723 offset = 0;
1724 for (int l = 0; l < nn - k; ++l)
1725 {
1726 order4[index] = start - offset;
1727 offset += k + 2 + l;
1728 index++;
1729 }
1730
1731 start += -k - 1;
1732 }
1733
1734 std::vector<int> order5(n_face_nodes); // B-> C
1735 index = 0;
1736 start = nn - 1;
1737 for (int k = 0; k < nn; ++k)
1738 {
1739 offset = 0;
1740 for (int l = 0; l < nn - k; ++l)
1741 {
1742 order5[index] = start + offset;
1743 offset += nn - 1 - l;
1744 index++;
1745 }
1746
1747 start--;
1748 }
1749
1750 std::vector<int> order6(n_face_nodes); // C-> B
1751 index = 0;
1752 start = n_face_nodes;
1753 for (int k = 0; k < nn; ++k)
1754 {
1755 offset = 0;
1756 start = start - k - 1;
1757 for (int l = 0; l < nn - k; ++l)
1758 {
1759 order6[index] = start - offset;
1760 offset += l + 1 + k;
1761 index++;
1762 }
1763 }
1764
1765 if (corners[0] == order1[0] && corners[1] == order1[nn - 1])
1766 {
1767 assert(corners[2] == order1[n_face_nodes - 1]);
1768 return order1;
1769 }
1770
1771 if (corners[0] == order2[0] && corners[1] == order2[nn - 1])
1772 {
1773 assert(corners[2] == order2[n_face_nodes - 1]);
1774 return order2;
1775 }
1776
1777 if (corners[0] == order3[0] && corners[1] == order3[nn - 1])
1778 {
1779 assert(corners[2] == order3[n_face_nodes - 1]);
1780 return order3;
1781 }
1782
1783 if (corners[0] == order4[0] && corners[1] == order4[nn - 1])
1784 {
1785 assert(corners[2] == order4[n_face_nodes - 1]);
1786 return order4;
1787 }
1788
1789 if (corners[0] == order5[0] && corners[1] == order5[nn - 1])
1790 {
1791 assert(corners[2] == order5[n_face_nodes - 1]);
1792 return order5;
1793 }
1794
1795 if (corners[0] == order6[0] && corners[1] == order6[nn - 1])
1796 {
1797 assert(corners[2] == order6[n_face_nodes - 1]);
1798 return order6;
1799 }
1800
1801 assert(false);
1802 return order1;
1803 };
1804
1805 Eigen::MatrixXd nodes;
1806 autogen::prism_nodes_3d(p, q, nodes);
1807 // auto pos = LagrangeBasis3d::linear_tet_face_local_nodes_coordinates(mesh, index);
1808 // Local to global mapping of node indices
1809
1810 // Extract requested interface
1811 std::array<int, 3> idx;
1812 for (int lv = 0; lv < 3; ++lv)
1813 {
1814 idx[lv] = find_index(l2g.begin(), l2g.end(), index.vertex);
1815 index = mesh.next_around_face(index);
1816 }
1817 Eigen::Matrix3d pos(3, 3);
1818 int cnt = 0;
1819 for (int i : idx)
1820 {
1821 pos.row(cnt++) = nodes.row(i);
1822 }
1823
1824 const Eigen::RowVector3d bary = pos.colwise().mean();
1825
1826 const int offset = 6 + global_n_edges_nodes;
1827 bool found = false;
1828 for (int lff = 0; lff < 2; ++lff)
1829 {
1830 Eigen::MatrixXd loc_nodes = nodes.block(offset + lff * n_face_nodes, 0, n_face_nodes, 3);
1831 Eigen::RowVector3d node_bary = loc_nodes.colwise().mean();
1832
1833 if ((node_bary - bary).norm() < 1e-10)
1834 {
1835 std::array<int, 3> corners;
1836 int sum = 0;
1837 for (int m = 0; m < 3; ++m)
1838 {
1839 auto t = pos.row(m);
1840 int min_n = -1;
1841 double min_dis = 10000;
1842
1843 for (int n = 0; n < n_face_nodes; ++n)
1844 {
1845 double dis = (loc_nodes.row(n) - t).squaredNorm();
1846 if (dis < min_dis)
1847 {
1848 min_dis = dis;
1849 min_n = n;
1850 }
1851 }
1852
1853 assert(min_n >= 0);
1854 assert(min_n < n_face_nodes);
1855 corners[m] = min_n;
1856 }
1857
1858 const auto indices = get_order(corners);
1859 for (int min_n : indices)
1860 {
1861 sum += min_n;
1862 result[ii++] = 6 + global_n_edges_nodes + min_n + lf * n_face_nodes;
1863 }
1864
1865 assert(sum == (n_face_nodes - 1) * n_face_nodes / 2);
1866
1867 found = true;
1868 assert(lff == lf);
1869
1870 break;
1871 }
1872 }
1873
1874 assert(found);
1875 }
1876 // else
1877 // {
1878 // assert(n_face_nodes == 0);
1879 // }
1880
1881 assert(ii == result.size());
1882 return result;
1883 }
1884 else
1885 {
1886 const int nn = q - 1;
1887 const int n_edge_nodes = nn * 12;
1888 const int n_face_nodes = nn * nn;
1889
1890 // Extract requested interface
1891 Eigen::VectorXi result(4 + nn * 4 + n_face_nodes);
1892 result[0] = find_index(l2g.begin(), l2g.end(), index.vertex);
1893 result[1] = find_index(l2g.begin(), l2g.end(), mesh.next_around_face(index).vertex);
1894 result[2] = find_index(l2g.begin(), l2g.end(), mesh.next_around_face(mesh.next_around_face(index)).vertex);
1895 result[3] = find_index(l2g.begin(), l2g.end(), mesh.next_around_face(mesh.next_around_face(mesh.next_around_face(index))).vertex);
1896
1897 Eigen::Matrix<Navigation3D::Index, 9, 1> e;
1898 Eigen::Matrix<int, 9, 2> ev;
1899 ev.row(0) << l2g[0], l2g[1];
1900 ev.row(1) << l2g[1], l2g[2];
1901 ev.row(2) << l2g[2], l2g[0];
1902
1903 ev.row(3) << l2g[3], l2g[4];
1904 ev.row(4) << l2g[4], l2g[5];
1905 ev.row(5) << l2g[5], l2g[3];
1906
1907 ev.row(6) << l2g[0], l2g[3];
1908 ev.row(7) << l2g[1], l2g[4];
1909 ev.row(8) << l2g[2], l2g[5];
1910
1911 Navigation3D::Index tmp = index;
1912
1913 for (int le = 0; le < e.rows(); ++le)
1914 {
1915 const auto l_index = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
1916 e[le] = l_index;
1917 }
1918
1919 int ii = 4;
1920 for (int k = 0; k < 4; ++k)
1921 {
1922 bool reverse = false;
1923 int le = 0;
1924 for (; le < ev.rows(); ++le)
1925 {
1926 // const auto l_index = find_edge(mesh, c, ev(le, 0), ev(le, 1));
1927 // const auto l_index = mesh.get_index_from_element_edge(c, ev(le, 0), ev(le, 1));
1928 const auto l_index = e[le];
1929 if (l_index.edge == tmp.edge)
1930 {
1931 if (l_index.vertex == tmp.vertex)
1932 reverse = false;
1933 else
1934 {
1935 reverse = true;
1936 assert(mesh.switch_vertex(tmp).vertex == l_index.vertex);
1937 }
1938
1939 break;
1940 }
1941 }
1942 assert(le < 9);
1943
1944 if (!reverse)
1945 {
1946
1947 for (int i = 0; i < q - 1; ++i)
1948 {
1949 result[ii++] = 6 + le * (q - 1) + i;
1950 }
1951 }
1952 else
1953 {
1954 for (int i = 0; i < q - 1; ++i)
1955 {
1956 result[ii++] = 6 + (le + 1) * (q - 1) - i - 1;
1957 }
1958 }
1959
1960 tmp = mesh.next_around_face(tmp);
1961 }
1962
1963 // faces
1964
1965 Eigen::Matrix<int, 3, 4> fv;
1966 fv.row(0) << l2g[0], l2g[3], l2g[4], l2g[1];
1967 fv.row(1) << l2g[1], l2g[4], l2g[5], l2g[2];
1968 fv.row(2) << l2g[2], l2g[5], l2g[3], l2g[0];
1969
1970 long lf = 0;
1971 for (; lf < fv.rows(); ++lf)
1972 {
1973 const auto l_index = find_quad_face(mesh, c, fv(lf, 0), fv(lf, 1), fv(lf, 2), fv(lf, 3));
1974 if (l_index.face == index.face)
1975 break;
1976 }
1977
1978 assert(lf < fv.rows());
1979
1980 if (n_face_nodes == 1)
1981 result[ii++] = 6 + global_n_edges_nodes + lf;
1982 else if (n_face_nodes != 0)
1983 {
1984 Eigen::MatrixXd nodes;
1985 autogen::q_nodes_3d(q, nodes);
1986 // auto pos = LagrangeBasis3d::linear_tet_face_local_nodes_coordinates(mesh, index);
1987 // Local to global mapping of node indices
1988
1989 // Extract requested interface
1990 std::array<int, 4> idx;
1991 for (int lv = 0; lv < 4; ++lv)
1992 {
1993 idx[lv] = find_index(l2g.begin(), l2g.end(), index.vertex);
1994 index = mesh.next_around_face(index);
1995 }
1996 Eigen::Matrix<double, 4, 3> pos(4, 3);
1997 int cnt = 0;
1998 for (int i : idx)
1999 {
2000 pos.row(cnt++) = nodes.row(i);
2001 }
2002
2003 const Eigen::RowVector3d bary = pos.colwise().mean();
2004
2005 const int offset = 8 + n_edge_nodes;
2006 bool found = false;
2007 for (int lff = 0; lff < 6; ++lff)
2008 {
2009 Eigen::Matrix<double, 4, 3> loc_nodes = nodes.block<4, 3>(offset + lff * n_face_nodes, 0);
2010 Eigen::RowVector3d node_bary = loc_nodes.colwise().mean();
2011
2012 if ((node_bary - bary).norm() < 1e-10)
2013 {
2014 int sum = 0;
2015 for (int m = 0; m < 4; ++m)
2016 {
2017 auto t = pos.row(m);
2018 int min_n = -1;
2019 double min_dis = 10000;
2020
2021 for (int n = 0; n < 4; ++n)
2022 {
2023 double dis = (loc_nodes.row(n) - t).squaredNorm();
2024 if (dis < min_dis)
2025 {
2026 min_dis = dis;
2027 min_n = n;
2028 }
2029 }
2030
2031 assert(min_n >= 0);
2032 assert(min_n < 4);
2033
2034 sum += min_n;
2035
2036 result[ii++] = 8 + n_edge_nodes + min_n + lf * n_face_nodes;
2037 }
2038
2039 assert(sum == 6); // 0 + 1 + 2 + 3
2040
2041 found = true;
2042 assert(lff == lf);
2043 }
2044
2045 if (found)
2046 break;
2047 }
2048
2049 assert(found);
2050 }
2051
2052 assert(ii == result.size());
2053 return result;
2054 }
2055}
2056
2058 const Mesh3D &mesh,
2059 const std::string &assembler,
2060 const int quadrature_order,
2061 const int mass_quadrature_order,
2062 const int discr_orderp,
2063 const int discr_orderq,
2064 const bool bernstein,
2065 const bool serendipity,
2066 const bool has_polys,
2067 const bool is_geom_bases,
2068 const bool use_corner_quadrature,
2069 std::vector<ElementBases> &bases,
2070 std::vector<LocalBoundary> &local_boundary,
2071 std::map<int, InterfaceData> &poly_face_to_data,
2072 std::shared_ptr<MeshNodes> &mesh_nodes)
2073{
2074 Eigen::VectorXi discr_ordersp(mesh.n_cells());
2075 discr_ordersp.setConstant(discr_orderp);
2076
2077 Eigen::VectorXi discr_ordersq(mesh.n_cells());
2078 discr_ordersq.setConstant(discr_orderq);
2079
2080 return build_bases(mesh, assembler, quadrature_order, mass_quadrature_order, discr_ordersp, discr_ordersq, bernstein, serendipity, has_polys, is_geom_bases, use_corner_quadrature, bases, local_boundary, poly_face_to_data, mesh_nodes);
2081}
2082
2084 const Mesh3D &mesh,
2085 const std::string &assembler,
2086 const int quadrature_order,
2087 const int mass_quadrature_order,
2088 const Eigen::VectorXi &discr_ordersp,
2089 const Eigen::VectorXi &discr_ordersq,
2090 const bool bernstein,
2091 const bool serendipity,
2092 const bool has_polys,
2093 const bool is_geom_bases,
2094 const bool use_corner_quadrature,
2095 std::vector<ElementBases> &bases,
2096 std::vector<LocalBoundary> &local_boundary,
2097 std::map<int, InterfaceData> &poly_face_to_data,
2098 std::shared_ptr<MeshNodes> &mesh_nodes)
2099{
2100 assert(mesh.is_volume());
2101 assert(discr_ordersp.size() == mesh.n_cells());
2102 assert(discr_ordersq.size() == mesh.n_cells());
2103
2104 // Navigation3D::get_index_from_element_face_time = 0;
2105 // Navigation3D::switch_vertex_time = 0;
2106 // Navigation3D::switch_edge_time = 0;
2107 // Navigation3D::switch_face_time = 0;
2108 // Navigation3D::switch_element_time = 0;
2109
2110 const int max_p = discr_ordersp.maxCoeff();
2111 const int max_q = discr_ordersq.maxCoeff();
2112 const int mmax = std::max(max_p, max_q);
2113 // assert(max_p < 5); //P5 not supported
2114
2115 const int nn = mmax > 1 ? (mmax - 1) : 0;
2116 const int n_face_nodes = nn * nn;
2117 const int n_cells_nodes = nn * nn * nn;
2118
2119 Eigen::VectorXi edge_orders, face_orders;
2120 if (!mesh.is_conforming())
2121 {
2122 const auto &ncmesh = dynamic_cast<const NCMesh3D &>(mesh);
2123 compute_edge_face_orders(ncmesh, discr_ordersp, edge_orders, face_orders);
2124 }
2125
2126 mesh_nodes = std::make_shared<MeshNodes>(mesh, has_polys, !is_geom_bases, nn, n_face_nodes * (is_geom_bases ? 2 : 1), mmax == 0 ? 1 : n_cells_nodes);
2127 MeshNodes &nodes = *mesh_nodes;
2128 std::vector<std::vector<int>> element_nodes_id, edge_virtual_nodes, face_virtual_nodes;
2129 compute_nodes(mesh, discr_ordersp, discr_ordersq, edge_orders, face_orders, serendipity, has_polys, is_geom_bases, nodes, edge_virtual_nodes, face_virtual_nodes, element_nodes_id, local_boundary, poly_face_to_data);
2130 // boundary_nodes = nodes.boundary_nodes();
2131
2132 bases.resize(mesh.n_cells());
2133 std::vector<int> interface_elements;
2134 interface_elements.reserve(mesh.n_faces());
2135
2136 for (int e = 0; e < mesh.n_cells(); ++e)
2137 {
2138 ElementBases &b = bases[e];
2139 const int discr_order = discr_ordersp(e);
2140 const int discr_orderq = discr_ordersq(e);
2141 const int n_el_bases = (int)element_nodes_id[e].size();
2142 b.bases.resize(n_el_bases);
2143
2144 bool skip_interface_element = false;
2145
2146 for (int j = 0; j < n_el_bases; ++j)
2147 {
2148 const int global_index = element_nodes_id[e][j];
2149 if (global_index < 0)
2150 {
2151 skip_interface_element = true;
2152 break;
2153 }
2154 }
2155
2156 if (skip_interface_element)
2157 {
2158 interface_elements.push_back(e);
2159 }
2160
2161 if (mesh.is_cube(e))
2162 {
2163 const int real_order = quadrature_order > 0 ? quadrature_order : AssemblerUtils::quadrature_order(assembler, discr_order, AssemblerUtils::BasisType::CUBE_LAGRANGE, 3);
2164 const int real_mass_order = mass_quadrature_order > 0 ? mass_quadrature_order : AssemblerUtils::quadrature_order("Mass", discr_order, AssemblerUtils::BasisType::CUBE_LAGRANGE, 3);
2165 b.set_quadrature([real_order](Quadrature &quad) {
2166 HexQuadrature hex_quadrature;
2167 hex_quadrature.get_quadrature(real_order, quad);
2168 });
2169 b.set_mass_quadrature([real_mass_order](Quadrature &quad) {
2170 HexQuadrature hex_quadrature;
2171 hex_quadrature.get_quadrature(real_mass_order, quad);
2172 });
2173
2174 b.set_local_node_from_primitive_func([serendipity, discr_order, e](const int primitive_id, const Mesh &mesh) {
2175 const auto &mesh3d = dynamic_cast<const Mesh3D &>(mesh);
2176 Navigation3D::Index index;
2177
2178 for (int lf = 0; lf < 6; ++lf)
2179 {
2180 index = mesh3d.get_index_from_element(e, lf, 0);
2181 if (index.face == primitive_id)
2182 break;
2183 }
2184 assert(index.face == primitive_id);
2185 return hex_face_local_nodes(serendipity, discr_order, mesh3d, index);
2186 });
2187
2188 for (int j = 0; j < n_el_bases; ++j)
2189 {
2190 const int global_index = element_nodes_id[e][j];
2191
2192 b.bases[j].init(discr_order, global_index, j, nodes.node_position(global_index));
2193
2194 const int dtmp = serendipity ? -2 : discr_order;
2195
2196 b.bases[j].set_basis([dtmp, j](const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) { autogen::q_basis_value_3d(dtmp, j, uv, val); });
2197 b.bases[j].set_grad([dtmp, j](const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) { autogen::q_grad_basis_value_3d(dtmp, j, uv, val); });
2198 }
2199 }
2200 else if (mesh.is_simplex(e))
2201 {
2202 const int real_order = quadrature_order > 0 ? quadrature_order : AssemblerUtils::quadrature_order(assembler, discr_order, AssemblerUtils::BasisType::SIMPLEX_LAGRANGE, 3);
2203 const int real_mass_order = mass_quadrature_order > 0 ? mass_quadrature_order : AssemblerUtils::quadrature_order("Mass", discr_order, AssemblerUtils::BasisType::SIMPLEX_LAGRANGE, 3);
2204
2205 b.set_quadrature([real_order, use_corner_quadrature](Quadrature &quad) {
2206 TetQuadrature tet_quadrature(use_corner_quadrature);
2207 tet_quadrature.get_quadrature(real_order, quad);
2208 });
2209 b.set_mass_quadrature([real_mass_order, use_corner_quadrature](Quadrature &quad) {
2210 TetQuadrature tet_quadrature(use_corner_quadrature);
2211 tet_quadrature.get_quadrature(real_mass_order, quad);
2212 });
2213
2214 b.set_local_node_from_primitive_func([discr_order, e](const int primitive_id, const Mesh &mesh) {
2215 const auto &mesh3d = dynamic_cast<const Mesh3D &>(mesh);
2216 Navigation3D::Index index;
2217
2218 for (int lf = 0; lf < mesh3d.n_cell_faces(e); ++lf)
2219 {
2220 index = mesh3d.get_index_from_element(e, lf, 0);
2221 if (index.face == primitive_id)
2222 break;
2223 }
2224 assert(index.face == primitive_id);
2225 return tet_face_local_nodes(discr_order, mesh3d, index);
2226 });
2227
2228 const bool rational = is_geom_bases && mesh.is_rational() && !mesh.cell_weights(e).empty();
2229 assert(!rational);
2230
2231 for (int j = 0; j < n_el_bases; ++j)
2232 {
2233 const int global_index = element_nodes_id[e][j];
2234 if (!skip_interface_element)
2235 {
2236 b.bases[j].init(discr_order, global_index, j, nodes.node_position(global_index));
2237 }
2238
2239 b.bases[j].set_basis([bernstein, discr_order, j](const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) { autogen::p_basis_value_3d(bernstein, discr_order, j, uv, val); });
2240 b.bases[j].set_grad([bernstein, discr_order, j](const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) { autogen::p_grad_basis_value_3d(bernstein, discr_order, j, uv, val); });
2241 }
2242 }
2243 else if (mesh.is_prism(e))
2244 {
2245 const int orderp = quadrature_order > 0 ? quadrature_order : AssemblerUtils::quadrature_order(assembler, discr_order, AssemblerUtils::BasisType::PRISM_LAGRANGE, 2);
2246 const int orderq = quadrature_order > 0 ? quadrature_order : AssemblerUtils::quadrature_order(assembler, discr_orderq, AssemblerUtils::BasisType::PRISM_LAGRANGE, 1);
2247
2248 const int mass_orderp = mass_quadrature_order > 0 ? mass_quadrature_order : AssemblerUtils::quadrature_order("Mass", discr_order, AssemblerUtils::BasisType::PRISM_LAGRANGE, 2);
2249 const int mass_orderq = mass_quadrature_order > 0 ? mass_quadrature_order : AssemblerUtils::quadrature_order("Mass", discr_orderq, AssemblerUtils::BasisType::PRISM_LAGRANGE, 1);
2250
2251 b.set_quadrature([orderp, orderq](Quadrature &quad) {
2252 PrismQuadrature tet_quadrature;
2253 tet_quadrature.get_quadrature(orderp, orderq, quad);
2254 });
2255 b.set_mass_quadrature([mass_orderp, mass_orderq](Quadrature &quad) {
2256 PrismQuadrature tet_quadrature;
2257 tet_quadrature.get_quadrature(mass_orderp, mass_orderq, quad);
2258 });
2259
2260 b.set_local_node_from_primitive_func([discr_order, discr_orderq, e](const int primitive_id, const Mesh &mesh) {
2261 const auto &mesh3d = dynamic_cast<const Mesh3D &>(mesh);
2262 Navigation3D::Index index;
2263
2264 for (int lf = 0; lf < mesh3d.n_cell_faces(e); ++lf)
2265 {
2266 index = mesh3d.get_index_from_element(e, lf, 0);
2267 if (index.face == primitive_id)
2268 break;
2269 }
2270 assert(index.face == primitive_id);
2271 return prism_face_local_nodes(discr_order, discr_orderq, mesh3d, index);
2272 });
2273
2274 for (int j = 0; j < n_el_bases; ++j)
2275 {
2276 const int global_index = element_nodes_id[e][j];
2277 if (!skip_interface_element)
2278 {
2279 b.bases[j].init(discr_order, global_index, j, nodes.node_position(global_index));
2280 }
2281
2282 b.bases[j].set_basis([discr_order, discr_orderq, j](const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) { autogen::prism_basis_value_3d(discr_order, discr_orderq, j, uv, val); });
2283 b.bases[j].set_grad([discr_order, discr_orderq, j](const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) { autogen::prism_grad_basis_value_3d(discr_order, discr_orderq, j, uv, val); });
2284 }
2285 }
2286 else
2287 {
2288 // Polyhedra bases are built later on
2289 // assert(false);
2290 }
2291 }
2292
2293 if (!is_geom_bases)
2294 {
2295 if (!mesh.is_conforming())
2296 {
2297 const auto &ncmesh = dynamic_cast<const NCMesh3D &>(mesh);
2298
2299 std::vector<std::vector<int>> elementOrder;
2300 {
2301 const int max_order = discr_ordersp.maxCoeff(), min_order = discr_ordersp.minCoeff();
2302 int max_level = 0;
2303 for (int e = 0; e < ncmesh.n_cells(); e++)
2304 if (max_level < ncmesh.cell_ref_level(e))
2305 max_level = ncmesh.cell_ref_level(e);
2306
2307 elementOrder.resize((max_level + 1) * (max_order - min_order + 1));
2308 int N = 0;
2309 int cur_level = 0;
2310 while (cur_level <= max_level)
2311 {
2312 int order = min_order;
2313 while (order <= max_order)
2314 {
2315 int cur_bucket = (max_order - min_order + 1) * cur_level + (order - min_order);
2316 for (int i = 0; i < ncmesh.n_cells(); i++)
2317 {
2318 if (ncmesh.cell_ref_level(i) != cur_level || discr_ordersp[i] != order)
2319 continue;
2320
2321 N++;
2322 elementOrder[cur_bucket].push_back(i);
2323 }
2324 order++;
2325 }
2326 cur_level++;
2327 }
2328 }
2329
2330 for (const auto &bucket : elementOrder)
2331 {
2332 if (bucket.size() == 0)
2333 continue;
2334 polyfem::utils::maybe_parallel_for((int)bucket.size(), [&](int start, int end, int thread_id) {
2335 for (int e_aux = start; e_aux < end; e_aux++)
2336 {
2337 const int e = bucket[e_aux];
2338 ElementBases &b = bases[e];
2339 const int discr_order = discr_ordersp(e);
2340 const int n_edge_nodes = discr_order - 1;
2341 const int n_face_nodes = (discr_order - 1) * (discr_order - 2) / 2;
2342 const int n_el_bases = element_nodes_id[e].size();
2343
2344 auto v = tet_vertices_local_to_global(mesh, e);
2345
2346 Eigen::Matrix<Navigation3D::Index, 4, 1> cell_faces;
2347 Eigen::Matrix<int, 4, 3> fv;
2348 fv.row(0) << v[0], v[1], v[2];
2349 fv.row(1) << v[0], v[1], v[3];
2350 fv.row(2) << v[1], v[2], v[3];
2351 fv.row(3) << v[2], v[0], v[3];
2352
2353 for (long lf = 0; lf < fv.rows(); ++lf)
2354 {
2355 const auto index = mesh.get_index_from_element_face(e, fv(lf, 0), fv(lf, 1), fv(lf, 2));
2356 cell_faces[lf] = index;
2357 }
2358
2359 Eigen::Matrix<Navigation3D::Index, 6, 1> cell_edges;
2360 Eigen::Matrix<int, 6, 2> ev;
2361 ev.row(0) << v[0], v[1];
2362 ev.row(1) << v[1], v[2];
2363 ev.row(2) << v[2], v[0];
2364
2365 ev.row(3) << v[0], v[3];
2366 ev.row(4) << v[1], v[3];
2367 ev.row(5) << v[2], v[3];
2368
2369 for (int le = 0; le < ev.rows(); ++le)
2370 {
2371 // const auto index = find_edge(mesh, c, ev(le, 0), ev(le, 1));
2372 const auto index = mesh.get_index_from_element_edge(e, ev(le, 0), ev(le, 1));
2373 cell_edges[le] = index;
2374 }
2375
2376 Eigen::MatrixXd verts(4, 3);
2377 for (int i = 0; i < ncmesh.n_cell_vertices(e); i++)
2378 verts.row(i) = ncmesh.point(v[i]);
2379
2380 for (int j = 0; j < n_el_bases; ++j)
2381 {
2382 const int global_index = element_nodes_id[e][j];
2383
2384 if (global_index >= 0)
2385 {
2386 b.bases[j].init(discr_order, global_index, j, nodes.node_position(global_index));
2387 }
2388 else
2389 {
2390 // vertex node - hanging vertex
2391 if (j < 4)
2392 {
2393 int large_elem = -1;
2394 if (ncmesh.leader_edge_of_vertex(v[j]) >= 0)
2395 {
2396 large_elem = lowest_order_elem_on_edge(ncmesh, discr_ordersp, ncmesh.leader_edge_of_vertex(v[j]));
2397 }
2398 else if (ncmesh.leader_face_of_vertex(v[j]) >= 0)
2399 {
2400 std::vector<int> ids;
2401 ncmesh.get_face_elements_neighs(ncmesh.leader_face_of_vertex(v[j]), ids);
2402 assert(ids.size() == 1);
2403 large_elem = ids[0];
2404 }
2405 else
2406 assert(false);
2407
2408 Eigen::MatrixXd large_elem_verts(4, 3);
2409 auto v_large = tet_vertices_local_to_global(mesh, large_elem);
2410 for (int i = 0; i < ncmesh.n_cell_vertices(large_elem); i++)
2411 large_elem_verts.row(i) = ncmesh.point(v_large[i]);
2412
2413 Eigen::MatrixXd node_position;
2414 global_to_local(large_elem_verts, verts.row(j), node_position);
2415
2416 // evaluate the basis of the large element at this node
2417 const auto &other_bases = bases[large_elem];
2418 std::vector<AssemblyValues> w;
2419 other_bases.evaluate_bases(node_position, w);
2420
2421 // apply basis projection
2422 for (long i = 0; i < w.size(); ++i)
2423 {
2424 assert(w[i].val.size() == 1);
2425 if (std::abs(w[i].val(0)) < 1e-12)
2426 continue;
2427
2428 assert(other_bases.bases[i].global().size() > 0);
2429 for (size_t ii = 0; ii < other_bases.bases[i].global().size(); ++ii)
2430 {
2431 const auto &other_global = other_bases.bases[i].global()[ii];
2432 assert(other_global.index >= 0);
2433 b.bases[j].global().emplace_back(other_global.index, other_global.node, w[i].val(0) * other_global.val);
2434 }
2435 }
2436 }
2437 // edge node - slave edge / edge on face / constrained order
2438 else if (j < 4 + 6 * n_edge_nodes)
2439 {
2440 const int local_edge_id = (j - 4) / n_edge_nodes;
2441 const int edge_id = cell_edges[local_edge_id].edge;
2442 bool need_extra_fake_nodes = false;
2443 int large_elem = -1;
2444
2445 // slave edge
2446 if (ncmesh.leader_edge_of_edge(edge_id) >= 0)
2447 {
2448 std::vector<int> ids;
2449 ncmesh.get_edge_elements_neighs(ncmesh.leader_edge_of_edge(edge_id), ids);
2450 large_elem = ids[0];
2451 }
2452 // edge on face
2453 else if (ncmesh.leader_face_of_edge(edge_id) >= 0)
2454 {
2455 std::vector<int> ids;
2456 ncmesh.get_face_elements_neighs(ncmesh.leader_face_of_edge(edge_id), ids);
2457 assert(ids.size() == 1);
2458 large_elem = ids[0];
2459 }
2460 // constrained order
2461 else if (discr_order > edge_orders[edge_id])
2462 {
2463 int min_order_elem = lowest_order_elem_on_edge(ncmesh, discr_ordersp, edge_id);
2464 // if haven't built min_order_elem? directly contribute to extra nodes
2465 if (discr_ordersp[min_order_elem] < discr_order)
2466 large_elem = min_order_elem;
2467
2468 // constrained order, master edge -- need extra fake nodes
2469 if (large_elem < 0)
2470 {
2471 // assert((edge.order < 2 || edge.global_ids.size() > 0) && edge.slaves.size() > 0);
2472 need_extra_fake_nodes = true;
2473 }
2474 }
2475 else
2476 assert(false);
2477
2478 assert(large_elem >= 0 || need_extra_fake_nodes);
2479 Eigen::MatrixXd lnodes;
2480 autogen::p_nodes_3d(discr_order, lnodes);
2481 Eigen::MatrixXd local_position = lnodes.row(j);
2482 if (need_extra_fake_nodes)
2483 {
2484 Eigen::MatrixXd global_position, edge_verts(2, 3);
2485 Eigen::VectorXd point_weight;
2486
2487 edge_verts.row(0) = ncmesh.point(ncmesh.edge_vertex(edge_id, 0));
2488 edge_verts.row(1) = ncmesh.point(ncmesh.edge_vertex(edge_id, 1));
2489
2490 local_to_global(verts, local_position, global_position);
2491 global_to_local_edge(edge_verts, global_position, point_weight);
2492
2493 std::function<double(const int, const int, const double)> basis_1d = [](const int order, const int id, const double x) -> double {
2494 assert(id <= order && id >= 0);
2495 double y = 1;
2496 for (int o = 0; o <= order; o++)
2497 {
2498 if (o != id)
2499 y *= (x * order - o) / (id - o);
2500 }
2501 return y;
2502 };
2503
2504 // contribution to edge nodes
2505 for (int i = 0; i < edge_virtual_nodes[edge_id].size(); i++)
2506 {
2507 const int global_index = edge_virtual_nodes[edge_id][i];
2508 // const double weight = basis_1d(edge_orders[edge_id], i+1, edge_weight);
2509 Eigen::VectorXd node_weight;
2510 global_to_local_edge(edge_verts, nodes.node_position(global_index), node_weight);
2511 const int basis_id = std::lround(node_weight(0) * edge_orders[edge_id]);
2512 const double weight = basis_1d(edge_orders[edge_id], basis_id, point_weight(0));
2513 if (std::abs(weight) < 1e-12)
2514 continue;
2515 b.bases[j].global().emplace_back(global_index, nodes.node_position(global_index), weight);
2516 }
2517
2518 // contribution to vertex nodes
2519 for (int i = 0; i < 2; i++)
2520 {
2521 const int lv = ev(local_edge_id, i);
2522 const auto &global_ = b.bases[lv].global();
2523 Eigen::VectorXd node_weight;
2524 global_to_local_edge(edge_verts, verts.row(lv), node_weight);
2525 const int basis_id = std::lround(node_weight(0) * edge_orders[edge_id]);
2526 const double weight = basis_1d(edge_orders[edge_id], basis_id, point_weight(0));
2527 if (std::abs(weight) > 1e-12)
2528 {
2529 assert(global_.size() > 0);
2530 for (size_t ii = 0; ii < global_.size(); ++ii)
2531 b.bases[j].global().emplace_back(global_[ii].index, global_[ii].node, weight * global_[ii].val);
2532 }
2533 }
2534 }
2535 else
2536 {
2537 Eigen::MatrixXd global_position, large_elem_verts(4, 3);
2538 auto v_large = tet_vertices_local_to_global(mesh, large_elem);
2539 for (int i = 0; i < ncmesh.n_cell_vertices(large_elem); i++)
2540 large_elem_verts.row(i) = ncmesh.point(v_large[i]);
2541 local_to_global(verts, local_position, global_position);
2542 global_to_local(large_elem_verts, global_position, local_position);
2543
2544 // evaluate the basis of the large element at this node
2545 const auto &other_bases = bases[large_elem];
2546 std::vector<AssemblyValues> w;
2547 other_bases.evaluate_bases(local_position, w);
2548
2549 // apply basis projection
2550 for (long i = 0; i < w.size(); ++i)
2551 {
2552 assert(w[i].val.size() == 1);
2553 if (std::abs(w[i].val(0)) < 1e-12)
2554 continue;
2555
2556 assert(other_bases.bases[i].global().size() > 0);
2557 for (size_t ii = 0; ii < other_bases.bases[i].global().size(); ++ii)
2558 {
2559 const auto &other_global = other_bases.bases[i].global()[ii];
2560 assert(other_global.index >= 0);
2561 b.bases[j].global().emplace_back(other_global.index, other_global.node, w[i].val(0) * other_global.val);
2562 }
2563 }
2564 }
2565 }
2566 // face node - slave face / constrained order
2567 else if (j < 4 + 6 * n_edge_nodes + 4 * n_face_nodes)
2568 {
2569 const int local_face_id = (j - (4 + 6 * n_edge_nodes)) / n_face_nodes;
2570 const int face_id = cell_faces(local_face_id).face;
2571 int large_elem = -1;
2572 bool need_extra_fake_nodes = false;
2573
2574 std::vector<int> ids;
2575 ncmesh.get_face_elements_neighs(ncmesh.leader_face_of_face(face_id), ids);
2576
2577 Eigen::MatrixXd face_verts(3, 3);
2578 for (int i = 0; i < ncmesh.n_face_vertices(face_id); i++)
2579 face_verts.row(i) = ncmesh.point(ncmesh.face_vertex(face_id, i));
2580
2581 // slave face
2582 if (ncmesh.leader_face_of_face(face_id) >= 0)
2583 {
2584 assert(ids.size() == 1);
2585 large_elem = ids[0];
2586 }
2587 // constrained order, conforming face
2588
2589 else if (face_orders[face_id] < discr_order && ids.size() == 2)
2590 {
2591 large_elem = ids[0] == e ? ids[1] : ids[0];
2592 }
2593 // constrained order, master face -- need extra fake nodes
2594 else if (face_orders[face_id] < discr_order && ncmesh.n_follower_faces(face_id) > 0)
2595 {
2596 // assert(ncface.global_ids.size() > 0 || ncface.order < 3);
2597 need_extra_fake_nodes = true;
2598 }
2599 else
2600 assert(false);
2601
2602 assert(large_elem >= 0 || need_extra_fake_nodes);
2603 Eigen::MatrixXd lnodes;
2604 autogen::p_nodes_3d(discr_order, lnodes);
2605 Eigen::MatrixXd local_position = lnodes.row(j);
2606 if (need_extra_fake_nodes)
2607 {
2608 Eigen::MatrixXd global_position;
2609 local_to_global(verts, local_position, global_position);
2610
2611 Eigen::MatrixXd tmp;
2612 global_to_local_face(face_verts, global_position, tmp);
2613 Eigen::VectorXd face_weight = tmp.transpose();
2614
2615 std::function<double(const int, const int, const double)> basis_aux = [](const int order, const int id, const double x) -> double {
2616 assert(id <= order && id >= 0);
2617 double y = 1;
2618 for (int o = 0; o < id; o++)
2619 y *= (x * order - o) / (id - o);
2620 return y;
2621 };
2622
2623 std::function<double(const int, const int, const int, const Eigen::Vector2d)> basis_2d = [&basis_aux](const int order, const int i, const int j, const Eigen::Vector2d uv) -> double {
2624 assert(i + j <= order && i >= 0 && j >= 0);
2625 double u = uv(0), v = uv(1);
2626 return basis_aux(order, i, u) * basis_aux(order, j, v) * basis_aux(order, order - i - j, 1 - u - v);
2627 };
2628
2629 // contribution to face nodes
2630 for (int global_ : face_virtual_nodes[face_id])
2631 {
2632 auto low_order_node = nodes.node_position(global_);
2633 Eigen::MatrixXd low_order_node_face_weight;
2634 global_to_local_face(face_verts, low_order_node, low_order_node_face_weight);
2635 int x = round(low_order_node_face_weight(0) * face_orders[face_id]), y = round(low_order_node_face_weight(1) * face_orders[face_id]);
2636 const double weight = basis_2d(face_orders[face_id], x, y, face_weight);
2637 if (std::abs(weight) < 1e-12)
2638 continue;
2639 b.bases[j].global().emplace_back(global_, nodes.node_position(global_), weight);
2640 }
2641
2642 // contribution to vertex nodes
2643 for (int i = 0; i < 3; i++)
2644 {
2645 const auto &global_ = b.bases[fv(local_face_id, i)].global();
2646 auto low_order_node = ncmesh.point(fv(local_face_id, i));
2647 Eigen::MatrixXd low_order_node_face_weight;
2648 global_to_local_face(face_verts, low_order_node, low_order_node_face_weight);
2649 int x = round(low_order_node_face_weight(0) * face_orders[face_id]), y = round(low_order_node_face_weight(1) * face_orders[face_id]);
2650 double weight = basis_2d(face_orders[face_id], x, y, face_weight);
2651 if (std::abs(weight) > 1e-12)
2652 {
2653 assert(global_.size() > 0);
2654 for (size_t ii = 0; ii < global_.size(); ++ii)
2655 b.bases[j].global().emplace_back(global_[ii].index, global_[ii].node, weight * global_[ii].val);
2656 }
2657 }
2658
2659 // contribution to edge nodes, two steps
2660 for (int x = 0, idx = 0; x <= face_orders[face_id]; x++)
2661 {
2662 for (int y = 0; x + y <= face_orders[face_id]; y++)
2663 {
2664 const int z = face_orders[face_id] - x - y;
2665 int flag = (int)(x == 0) + (int)(y == 0) + (int)(z == 0);
2666 if (flag != 1)
2667 continue;
2668
2669 // first step
2670 const double weight = basis_2d(face_orders[face_id], x, y, face_weight);
2671 if (std::abs(weight) < 1e-12)
2672 continue;
2673 Eigen::MatrixXd face_weight(1, 2);
2674 face_weight << (double)x / face_orders[face_id], (double)y / face_orders[face_id];
2675 Eigen::MatrixXd pos, local_pos;
2676 local_to_global_face(face_verts, face_weight, pos);
2677 global_to_local(verts, pos, local_pos);
2678 Local2Global step1(idx, local_pos, weight);
2679 idx++;
2680
2681 {
2682 // evaluate the basis of the large element at this node
2683 const auto &other_bases = bases[e];
2684 std::vector<AssemblyValues> w;
2685 other_bases.evaluate_bases(local_pos, w);
2686
2687 // apply basis projection
2688 for (long i = 0; i < w.size(); ++i)
2689 {
2690 assert(w[i].val.size() == 1);
2691 if (std::abs(w[i].val(0)) < 1e-12)
2692 continue;
2693
2694 assert(other_bases.bases[i].global().size() > 0);
2695 for (size_t ii = 0; ii < other_bases.bases[i].global().size(); ++ii)
2696 {
2697 const auto &other_global = other_bases.bases[i].global()[ii];
2698 assert(other_global.index >= 0);
2699 b.bases[j].global().emplace_back(other_global.index, other_global.node, step1.val * w[i].val(0) * other_global.val);
2700 }
2701 }
2702 }
2703 }
2704 }
2705 }
2706 else
2707 {
2708 Eigen::MatrixXd global_position, large_elem_verts(4, 3);
2709 auto v_large = tet_vertices_local_to_global(mesh, large_elem);
2710 for (int i = 0; i < ncmesh.n_cell_vertices(large_elem); i++)
2711 large_elem_verts.row(i) = ncmesh.point(v_large[i]);
2712 local_to_global(verts, local_position, global_position);
2713 global_to_local(large_elem_verts, global_position, local_position);
2714
2715 // evaluate the basis of the large element at this node
2716 const auto &other_bases = bases[large_elem];
2717 std::vector<AssemblyValues> w;
2718 other_bases.evaluate_bases(local_position, w);
2719
2720 // apply basis projection
2721 for (long i = 0; i < w.size(); ++i)
2722 {
2723 assert(w[i].val.size() == 1);
2724 if (std::abs(w[i].val(0)) < 1e-12)
2725 continue;
2726
2727 assert(other_bases.bases[i].global().size() > 0);
2728 for (size_t ii = 0; ii < other_bases.bases[i].global().size(); ++ii)
2729 {
2730 const auto &other_global = other_bases.bases[i].global()[ii];
2731 assert(other_global.index >= 0);
2732 b.bases[j].global().emplace_back(other_global.index, other_global.node, w[i].val(0) * other_global.val);
2733 }
2734 }
2735 }
2736 }
2737 else
2738 assert(false);
2739
2740 auto &global_ = b.bases[j].global();
2741 if (global_.size() <= 1)
2742 continue;
2743
2744 std::map<int, Local2Global> list;
2745 for (size_t ii = 0; ii < global_.size(); ii++)
2746 {
2747 auto pair = list.insert({global_[ii].index, global_[ii]});
2748 if (!pair.second && pair.first != list.end())
2749 {
2750 assert((pair.first->second.node - global_[ii].node).norm() < 1e-12);
2751 pair.first->second.val += global_[ii].val;
2752 }
2753 }
2754
2755 global_.clear();
2756 for (auto it = list.begin(); it != list.end(); ++it)
2757 {
2758 if (std::abs(it->second.val) > 1e-12)
2759 {
2760 global_.push_back(it->second);
2761 }
2762 }
2763 }
2764 }
2765 }
2766 });
2767 }
2768 }
2769 else
2770 {
2771 for (int pp = 2; pp <= autogen::MAX_P_BASES; ++pp)
2772 {
2773 for (int e : interface_elements)
2774 {
2775 ElementBases &b = bases[e];
2776 // todo non conforming
2777 const int discr_order = discr_ordersp(e);
2778 const int n_el_bases = element_nodes_id[e].size();
2779 assert(discr_order > 1);
2780 if (discr_order != pp)
2781 continue;
2782
2783 if (mesh.is_cube(e) || mesh.is_prism(e))
2784 {
2785 // TODO
2786 assert(false);
2787 }
2788 else if (mesh.is_simplex(e))
2789 {
2790 for (int j = 0; j < n_el_bases; ++j)
2791 {
2792 const int global_index = element_nodes_id[e][j];
2793
2794 if (global_index >= 0)
2795 b.bases[j].init(discr_order, global_index, j, nodes.node_position(global_index));
2796 else
2797 {
2798 const int lnn = max_p > 2 ? (discr_order - 2) : 0;
2799 const int ln_edge_nodes = discr_order - 1;
2800 const int ln_face_nodes = lnn * (lnn + 1) / 2;
2801
2802 const auto v = tet_vertices_local_to_global(mesh, e);
2803 Navigation3D::Index index;
2804 if (global_index <= -30)
2805 {
2806 assert(false);
2807 // const auto lv = -(global_index + 30);
2808 // assert(lv>=0 && lv < 4);
2809 // assert(j < 4);
2810
2811 // if(lv == 3)
2812 // {
2813 // index = mesh.switch_element(find_edge(mesh, e, v[lv], v[0]));
2814 // if(index.element < 0)
2815 // index = mesh.switch_element(find_edge(mesh, e, v[lv], v[1]));
2816 // if(index.element < 0)
2817 // index = mesh.switch_element(find_edge(mesh, e, v[lv], v[2]));
2818 // }
2819 // else
2820 // {
2821 // index = mesh.switch_element(find_edge(mesh, e, v[lv], v[(lv+1)%3]));
2822 // if(index.element < 0)
2823 // index = mesh.switch_element(find_edge(mesh, e, v[lv], v[(lv+2)%3]));
2824 // if(index.element < 0)
2825 // index = mesh.switch_element(find_edge(mesh, e, v[lv], v[3]));
2826 // }
2827 }
2828 else if (global_index <= -10)
2829 {
2830 const auto le = -(global_index + 10);
2831 assert(le >= 0 && le < 6);
2832 assert(j >= 4 && j < 4 + 6 * ln_edge_nodes);
2833
2834 Eigen::Matrix<int, 6, 2> ev;
2835 ev.row(0) << v[0], v[1];
2836 ev.row(1) << v[1], v[2];
2837 ev.row(2) << v[2], v[0];
2838
2839 ev.row(3) << v[0], v[3];
2840 ev.row(4) << v[1], v[3];
2841 ev.row(5) << v[2], v[3];
2842
2843 // const auto edge_index = find_edge(mesh, e, ev(le, 0), ev(le, 1));
2844 const auto edge_index = mesh.get_index_from_element_edge(e, ev(le, 0), ev(le, 1));
2845 auto neighs = mesh.edge_neighs(edge_index.edge);
2846 int min_p = discr_order;
2847 int min_cell = edge_index.element;
2848
2849 for (auto cid : neighs)
2850 {
2851 if (discr_ordersp[cid] < min_p)
2852 {
2853 min_p = discr_ordersp[cid];
2854 min_cell = cid;
2855 }
2856 }
2857
2858 bool found = false;
2859 for (int lf = 0; lf < 4; ++lf)
2860 {
2861 for (int lv = 0; lv < 4; ++lv)
2862 {
2863 index = mesh.get_index_from_element(min_cell, lf, lv);
2864
2865 if (index.vertex == edge_index.vertex)
2866 {
2867 if (index.edge != edge_index.edge)
2868 {
2869 auto tmp = index;
2870 index = mesh.switch_edge(tmp);
2871
2872 if (index.edge != edge_index.edge)
2873 {
2874 index = mesh.switch_edge(mesh.switch_face(tmp));
2875 }
2876 }
2877 found = true;
2878 break;
2879 }
2880 }
2881
2882 if (found)
2883 break;
2884 }
2885
2886 assert(found);
2887 assert(index.vertex == edge_index.vertex && index.edge == edge_index.edge);
2888 assert(index.element != edge_index.element);
2889 }
2890 else
2891 {
2892 const auto lf = -(global_index + 1);
2893 assert(lf >= 0 && lf < 4);
2894 assert(j >= 4 + 6 * ln_edge_nodes && j < 4 + 6 * ln_edge_nodes + 4 * ln_face_nodes);
2895
2896 Eigen::Matrix<int, 4, 3> fv;
2897 fv.row(0) << v[0], v[1], v[2];
2898 fv.row(1) << v[0], v[1], v[3];
2899 fv.row(2) << v[1], v[2], v[3];
2900 fv.row(3) << v[2], v[0], v[3];
2901
2902 index = mesh.switch_element(mesh.get_index_from_element_face(e, fv(lf, 0), fv(lf, 1), fv(lf, 2)));
2903 }
2904
2905 const auto other_cell = index.element;
2906 assert(other_cell >= 0);
2907 assert(discr_order > discr_ordersp(other_cell));
2908
2909 auto indices = tet_face_local_nodes(discr_order, mesh, index);
2910 Eigen::MatrixXd lnodes;
2911 autogen::p_nodes_3d(discr_order, lnodes);
2912 Eigen::RowVector3d node_position; // = lnodes.row(indices(ii));
2913
2914 if (j < 4)
2915 node_position = lnodes.row(indices(0));
2916 else if (j < 4 + 6 * ln_edge_nodes)
2917 node_position = lnodes.row(indices(((j - 4) % ln_edge_nodes) + 3));
2918 else if (j < 4 + 6 * ln_edge_nodes + 4 * ln_face_nodes)
2919 {
2920 // node_position = lnodes.row(indices(((j - 4 - 6*ln_edge_nodes) % ln_face_nodes) + 3 + 3*ln_edge_nodes));
2921 auto me_indices = tet_face_local_nodes(discr_order, mesh, mesh.switch_element(index));
2922 int ii;
2923 for (ii = 0; ii < me_indices.size(); ++ii)
2924 {
2925 if (me_indices(ii) == j)
2926 break;
2927 }
2928
2929 assert(ii >= 3 + 3 * ln_edge_nodes);
2930 assert(ii < me_indices.size());
2931
2932 node_position = lnodes.row(indices(ii));
2933 }
2934 else
2935 assert(false);
2936
2937 const auto &other_bases = bases[other_cell];
2938 // Eigen::MatrixXd w;
2939 std::vector<AssemblyValues> w;
2940 other_bases.evaluate_bases(node_position, w);
2941
2942 assert(b.bases[j].global().size() == 0);
2943
2944 for (long i = 0; i < w.size(); ++i)
2945 {
2946 assert(w[i].val.size() == 1);
2947 if (std::abs(w[i].val(0)) < 1e-8)
2948 continue;
2949
2950 // assert(other_bases.bases[i].global().size() == 1);
2951 for (size_t ii = 0; ii < other_bases.bases[i].global().size(); ++ii)
2952 {
2953 const auto &other_global = other_bases.bases[i].global()[ii];
2954 // logger().trace("e {} j {} gid {}", e, j, other_global.index);
2955 b.bases[j].global().emplace_back(other_global.index, other_global.node, w[i].val(0) * other_global.val);
2956 }
2957 }
2958 }
2959 }
2960 }
2961 else
2962 {
2963 // Polygon bases are built later on
2964 }
2965 }
2966 }
2967 }
2968 }
2969
2970 return nodes.n_nodes();
2971}
Eigen::MatrixXd vec
Definition Assembler.cpp:72
double val
Definition Assembler.cpp:86
int x
static int quadrature_order(const std::string &assembler, const int basis_degree, const BasisType &b_type, const int dim)
utility for retrieving the needed quadrature order to precisely integrate the given form on the given...
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d).
void set_quadrature(const QuadratureFunction &fun)
void set_local_node_from_primitive_func(LocalNodeFromPrimitiveFunc fun)
sets mapping from local nodes to global nodes
std::vector< Basis > bases
one basis function per node in the element
void set_mass_quadrature(const QuadratureFunction &fun)
static Eigen::VectorXi hex_face_local_nodes(const bool serendipity, const int q, const mesh::Mesh3D &mesh, mesh::Navigation3D::Index index)
static Eigen::VectorXi prism_face_local_nodes(const int p, const int q, const mesh::Mesh3D &mesh, mesh::Navigation3D::Index index)
static int build_bases(const mesh::Mesh3D &mesh, const std::string &assembler, const int quadrature_order, const int mass_quadrature_order, const int discr_orderp, const int discr_orderq, const bool bernstein, const bool serendipity, const bool has_polys, const bool is_geom_bases, const bool use_corner_quadrature, std::vector< ElementBases > &bases, std::vector< mesh::LocalBoundary > &local_boundary, std::map< int, InterfaceData > &poly_face_to_data, std::shared_ptr< mesh::MeshNodes > &mesh_nodes)
Builds FE basis functions over the entire mesh (P1, P2 over tets, Q1, Q2 over hes).
static Eigen::VectorXi tet_face_local_nodes(const int p, const mesh::Mesh3D &mesh, mesh::Navigation3D::Index index)
Boundary primitive IDs for a single element.
virtual Navigation3D::Index get_index_from_element(int hi, int lf, int lv) const =0
virtual Navigation3D::Index switch_element(Navigation3D::Index idx) const =0
virtual Navigation3D::Index get_index_from_element_edge(int hi, int v0, int v1) const =0
virtual int n_cell_faces(const int c_id) const =0
virtual std::vector< uint32_t > edge_neighs(const int e_gid) const =0
bool is_volume() const override
checks if mesh is volume
Definition Mesh3D.hpp:28
virtual Navigation3D::Index next_around_face(Navigation3D::Index idx) const =0
virtual Navigation3D::Index get_index_from_element_face(int hi, int v0, int v1, int v2) const =0
virtual Navigation3D::Index switch_vertex(Navigation3D::Index idx) const =0
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:40
bool is_polytope(const int el_id) const
checks if element is polygon compatible
Definition Mesh.cpp:365
bool is_rational() const
check if curved mesh has rational polynomials elements
Definition Mesh.hpp:288
virtual bool is_conforming() const =0
if the mesh is conforming
bool is_cube(const int el_id) const
checks if element is cube compatible
Definition Mesh.cpp:352
virtual bool is_boundary_face(const int face_global_id) const =0
is face boundary
virtual int get_boundary_id(const int primitive) const
Get the boundary selection of an element (face in 3d, edge in 2d)
Definition Mesh.hpp:481
bool is_simplex(const int el_id) const
checks if element is simplex
Definition Mesh.cpp:422
bool is_prism(const int el_id) const
checks if element is a prism
Definition Mesh.cpp:427
const std::vector< double > & cell_weights(const int cell_index) const
weights for rational polynomial meshes
Definition Mesh.hpp:556
virtual int n_cells() const =0
number of cells
virtual int n_faces() const =0
number of faces
virtual int n_face_vertices(const int f_id) const =0
number of vertices of a face
int cell_ref_level(const int c_id) const
Definition NCMesh3D.hpp:227
int n_faces() const override
number of faces
Definition NCMesh3D.hpp:189
int face_edge(const int f_id, const int le_id) const
Definition NCMesh3D.cpp:18
int n_cell_faces(const int c_id) const override
Definition NCMesh3D.hpp:219
int n_edge_cells(const int e_id) const
Definition NCMesh3D.hpp:216
int n_cells() const override
number of cells
Definition NCMesh3D.hpp:188
int n_edges() const override
number of edges
Definition NCMesh3D.hpp:197
int n_face_vertices(const int f_id) const override
number of vertices of a face
Definition NCMesh3D.hpp:214
std::vector< uint32_t > edge_neighs(const int e_gid) const override
Definition NCMesh3D.cpp:446
int leader_face_of_edge(const int e_id) const
Definition NCMesh3D.hpp:298
int cell_face(const int c_id, const int lf_id) const override
Definition NCMesh3D.hpp:221
int leader_edge_of_edge(const int e_id) const
Definition NCMesh3D.hpp:286
int n_cell_edges(const int c_id) const override
Definition NCMesh3D.hpp:218
int cell_edge(const int c_id, const int le_id) const override
Definition NCMesh3D.hpp:222
int leader_face_of_face(const int f_id) const
Definition NCMesh3D.hpp:304
int n_face_cells(const int f_id) const
Definition NCMesh3D.hpp:215
void get_quadrature(const int order, Quadrature &quad)
void get_quadrature(const int order, const int order_h, Quadrature &quad)
list tmp
Definition p_bases.py:365
str nodes
Definition p_bases.py:398
list indices
Definition p_bases.py:258
Used for test only.
void q_grad_basis_value_3d(const int q, const int local_index, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val)
void prism_basis_value_3d(const int p, const int q, const int li, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val)
void prism_grad_basis_value_3d(const int p, const int q, const int li, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val)
void prism_nodes_3d(const int p, const int q, Eigen::MatrixXd &val)
void p_grad_basis_value_3d(const bool bernstein, const int p, const int local_index, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val)
void p_nodes_3d(const int p, Eigen::MatrixXd &val)
void p_basis_value_3d(const bool bernstein, const int p, const int local_index, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val)
void q_nodes_3d(const int q, Eigen::MatrixXd &val)
void q_basis_value_3d(const int q, const int local_index, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val)
void maybe_parallel_for(int size, const std::function< void(int, int, int)> &partial_for)
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:73
std::vector< int > local_indices