deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
reference_cell.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2020 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
19
28
31#include <deal.II/grid/tria.h>
32
33#include <algorithm>
34#include <iostream>
35#include <memory>
36
38
39namespace
40{
41 namespace VTKCellType
42 {
43 // Define VTK constants for linear, quadratic and
44 // high-order Lagrange geometrices
45 enum : unsigned int
46 {
47 VTK_VERTEX = 1,
48 // Linear cells
49 VTK_LINE = 3,
50 VTK_TRIANGLE = 5,
51 VTK_QUAD = 9,
52 VTK_TETRA = 10,
53 VTK_HEXAHEDRON = 12,
54 VTK_WEDGE = 13,
55 VTK_PYRAMID = 14,
56 // Quadratic cells
57 VTK_QUADRATIC_EDGE = 21,
58 VTK_QUADRATIC_TRIANGLE = 22,
59 VTK_QUADRATIC_QUAD = 23,
60 VTK_QUADRATIC_TETRA = 24,
61 VTK_QUADRATIC_HEXAHEDRON = 25,
62 VTK_QUADRATIC_WEDGE = 26,
63 VTK_QUADRATIC_PYRAMID = 27,
64 // Lagrange cells
65 VTK_LAGRANGE_CURVE = 68,
66 VTK_LAGRANGE_TRIANGLE = 69,
67 VTK_LAGRANGE_QUADRILATERAL = 70,
68 VTK_LAGRANGE_TETRAHEDRON = 71,
69 VTK_LAGRANGE_HEXAHEDRON = 72,
70 VTK_LAGRANGE_WEDGE = 73,
71 VTK_LAGRANGE_PYRAMID = 74,
72 // Invalid code
74 };
75
76 } // namespace VTKCellType
77
78} // namespace
79
81
84
87
88std::string
90{
91 switch (this->kind)
92 {
94 return "Vertex";
96 return "Line";
98 return "Tri";
100 return "Quad";
102 return "Tet";
104 return "Pyramid";
106 return "Wedge";
108 return "Hex";
110 return "Invalid";
111 default:
113 }
114
115 return "Invalid";
116}
117
118
119
120template <int dim>
121std::pair<unsigned int, RefinementCase<dim - 1>>
123 const types::geometric_orientation combined_face_orientation,
124 const internal::SubfaceCase<dim> subface_case,
125 const unsigned int subface_no) const
126{
127 if constexpr (dim == 3)
128 {
129 static const RefinementCase<dim - 1>
130 equivalent_refine_case[internal::SubfaceCase<dim>::case_isotropic + 1]
132 // case_none. there should be only
133 // invalid values here. However, as
134 // this function is also called (in
135 // tests) for cells which have no
136 // refined faces, use isotropic
137 // refinement instead
138 {RefinementCase<dim - 1>::cut_xy,
139 RefinementCase<dim - 1>::cut_xy,
140 RefinementCase<dim - 1>::cut_xy,
141 RefinementCase<dim - 1>::cut_xy},
142 // case_x
143 {RefinementCase<dim - 1>::cut_x,
144 RefinementCase<dim - 1>::cut_x,
145 RefinementCase<dim - 1>::no_refinement,
146 RefinementCase<dim - 1>::no_refinement},
147 // case_x1y
148 {RefinementCase<dim - 1>::cut_xy,
149 RefinementCase<dim - 1>::cut_xy,
150 RefinementCase<dim - 1>::cut_x,
151 RefinementCase<dim - 1>::no_refinement},
152 // case_x2y
153 {RefinementCase<dim - 1>::cut_x,
154 RefinementCase<dim - 1>::cut_xy,
155 RefinementCase<dim - 1>::cut_xy,
156 RefinementCase<dim - 1>::no_refinement},
157 // case_x1y2y
158 {RefinementCase<dim - 1>::cut_xy,
159 RefinementCase<dim - 1>::cut_xy,
160 RefinementCase<dim - 1>::cut_xy,
161 RefinementCase<dim - 1>::cut_xy},
162 // case_y
163 {RefinementCase<dim - 1>::cut_y,
164 RefinementCase<dim - 1>::cut_y,
165 RefinementCase<dim - 1>::no_refinement,
166 RefinementCase<dim - 1>::no_refinement},
167 // case_y1x
168 {RefinementCase<dim - 1>::cut_xy,
169 RefinementCase<dim - 1>::cut_xy,
170 RefinementCase<dim - 1>::cut_y,
171 RefinementCase<dim - 1>::no_refinement},
172 // case_y2x
173 {RefinementCase<dim - 1>::cut_y,
174 RefinementCase<dim - 1>::cut_xy,
175 RefinementCase<dim - 1>::cut_xy,
176 RefinementCase<dim - 1>::no_refinement},
177 // case_y1x2x
178 {RefinementCase<dim - 1>::cut_xy,
179 RefinementCase<dim - 1>::cut_xy,
180 RefinementCase<dim - 1>::cut_xy,
181 RefinementCase<dim - 1>::cut_xy},
182 // case_xy (case_isotropic)
183 {RefinementCase<dim - 1>::cut_xy,
184 RefinementCase<dim - 1>::cut_xy,
185 RefinementCase<dim - 1>::cut_xy,
186 RefinementCase<dim - 1>::cut_xy}};
187
188 constexpr unsigned int X = numbers::invalid_unsigned_int;
189 static const unsigned int
190 equivalent_subface_number[internal::SubfaceCase<dim>::case_isotropic +
192 // case_none, see above
193 {0, 1, 2, 3},
194 // case_x
195 {0, 1, X, X},
196 // case_x1y
197 {0, 2, 1, X},
198 // case_x2y
199 {0, 1, 3, X},
200 // case_x1y2y
201 {0, 2, 1, 3},
202 // case_y
203 {0, 1, X, X},
204 // case_y1x
205 {0, 1, 1, X},
206 // case_y2x
207 {0, 2, 3, X},
208 // case_y1x2x
209 {0, 1, 2, 3},
210 // case_xy (case_isotropic)
211 {0, 1, 2, 3}};
212
213 static const RefinementCase<dim - 1> rotated_refinement_case[4] = {
214 RefinementCase<dim - 1>::no_refinement,
215 RefinementCase<dim - 1>::cut_y,
216 RefinementCase<dim - 1>::cut_x,
217 RefinementCase<dim - 1>::cut_xy};
218 const auto [face_orientation, face_rotation, face_flip] =
219 internal::split_face_orientation(combined_face_orientation);
220
221 const auto equivalent_refinement_case =
222 equivalent_refine_case[subface_case][subface_no];
223 const unsigned int equivalent_subface_no =
224 equivalent_subface_number[subface_case][subface_no];
225 // make sure, that we got a valid subface and RefineCase
228 Assert(equivalent_subface_no != X, ExcInternalError());
229 // now, finally respect non-standard faces
230 const RefinementCase<dim - 1> final_refinement_case =
231 (face_orientation == face_rotation ?
232 rotated_refinement_case[equivalent_refinement_case] :
234
235 const unsigned int final_subface_no =
237 final_refinement_case),
238 /*face_no = */ 4,
239 equivalent_subface_no,
240 face_orientation,
241 face_flip,
242 face_rotation,
244
245 return std::make_pair(final_subface_no, final_refinement_case);
246 }
247 else
248 {
249 (void)combined_face_orientation;
250 (void)subface_case;
251 (void)subface_no;
252
254 return {};
255 }
256}
257
258
259
260template <int dim, int spacedim>
261std::unique_ptr<Mapping<dim, spacedim>>
262ReferenceCell::get_default_mapping(const unsigned int degree) const
263{
265
266 if (is_hyper_cube())
267 return std::make_unique<MappingQ<dim, spacedim>>(degree);
268 else if (is_simplex())
269 return std::make_unique<MappingFE<dim, spacedim>>(
271 else if (*this == ReferenceCells::Pyramid)
272 return std::make_unique<MappingFE<dim, spacedim>>(
274 else if (*this == ReferenceCells::Wedge)
275 return std::make_unique<MappingFE<dim, spacedim>>(
277 else
278 {
280 }
281
282 return std::make_unique<MappingQ<dim, spacedim>>(degree);
283}
284
285
286
287template <int dim, int spacedim>
290{
292
293 if (is_hyper_cube())
294 {
296 }
297 else if (is_simplex())
298 {
300 return mapping;
301 }
302 else if (*this == ReferenceCells::Pyramid)
303 {
306 return mapping;
307 }
308 else if (*this == ReferenceCells::Wedge)
309 {
312 return mapping;
313 }
314 else
315 {
317 }
318
319 return StaticMappingQ1<dim, spacedim>::mapping; // never reached
320}
321
322
323
324template <int dim>
326ReferenceCell::get_gauss_type_quadrature(const unsigned n_points_1d) const
327{
329
330 if (is_hyper_cube())
331 return QGauss<dim>(n_points_1d);
332 else if (is_simplex())
333 return QGaussSimplex<dim>(n_points_1d);
334 else if (*this == ReferenceCells::Pyramid)
335 return QGaussPyramid<dim>(n_points_1d);
336 else if (*this == ReferenceCells::Wedge)
337 return QGaussWedge<dim>(n_points_1d);
338 else
340
341 return Quadrature<dim>(); // never reached
342}
343
344
345
346template <int dim>
347const Quadrature<dim> &
349{
351
352 // A function that is used to fill a quadrature object of the
353 // desired type the first time we encounter a particular
354 // reference cell
355 const auto create_quadrature = [](const ReferenceCell &reference_cell) {
356 std::vector<Point<dim>> vertices(reference_cell.n_vertices());
357 for (const unsigned int v : reference_cell.vertex_indices())
358 vertices[v] = reference_cell.vertex<dim>(v);
359
360 return Quadrature<dim>(vertices);
361 };
362
363 if (is_hyper_cube())
364 {
365 static const Quadrature<dim> quadrature = create_quadrature(*this);
366 return quadrature;
367 }
368 else if (is_simplex())
369 {
370 static const Quadrature<dim> quadrature = create_quadrature(*this);
371 return quadrature;
372 }
373 else if (*this == ReferenceCells::Pyramid)
374 {
375 static const Quadrature<dim> quadrature = create_quadrature(*this);
376 return quadrature;
377 }
378 else if (*this == ReferenceCells::Wedge)
379 {
380 static const Quadrature<dim> quadrature = create_quadrature(*this);
381 return quadrature;
382 }
383 else
385
386 static const Quadrature<dim> dummy;
387 return dummy; // never reached
388}
389
390
391
392unsigned int
393ReferenceCell::exodusii_vertex_to_deal_vertex(const unsigned int vertex_n) const
394{
395 AssertIndexRange(vertex_n, n_vertices());
396
397 switch (this->kind)
398 {
401 return vertex_n;
403 {
404 constexpr std::array<unsigned int, 4> exodus_to_deal{{0, 1, 3, 2}};
405 return exodus_to_deal[vertex_n];
406 }
408 return vertex_n;
410 {
411 constexpr std::array<unsigned int, 8> exodus_to_deal{
412 {0, 1, 3, 2, 4, 5, 7, 6}};
413 return exodus_to_deal[vertex_n];
414 }
416 {
417 constexpr std::array<unsigned int, 6> exodus_to_deal{
418 {2, 1, 0, 5, 4, 3}};
419 return exodus_to_deal[vertex_n];
420 }
422 {
423 constexpr std::array<unsigned int, 5> exodus_to_deal{{0, 1, 3, 2, 4}};
424 return exodus_to_deal[vertex_n];
425 }
426 default:
428 }
429
431}
432
433
434
435unsigned int
436ReferenceCell::exodusii_face_to_deal_face(const unsigned int face_n) const
437{
438 AssertIndexRange(face_n, n_faces());
439
440 switch (this->kind)
441 {
443 return 0;
446 return face_n;
448 {
449 constexpr std::array<unsigned int, 4> exodus_to_deal{{2, 1, 3, 0}};
450 return exodus_to_deal[face_n];
451 }
453 {
454 constexpr std::array<unsigned int, 4> exodus_to_deal{{1, 3, 2, 0}};
455 return exodus_to_deal[face_n];
456 }
458 {
459 constexpr std::array<unsigned int, 6> exodus_to_deal{
460 {2, 1, 3, 0, 4, 5}};
461 return exodus_to_deal[face_n];
462 }
464 {
465 constexpr std::array<unsigned int, 6> exodus_to_deal{{3, 4, 2, 0, 1}};
466 return exodus_to_deal[face_n];
467 }
469 {
470 constexpr std::array<unsigned int, 5> exodus_to_deal{{3, 2, 4, 1, 0}};
471 return exodus_to_deal[face_n];
472 }
473 default:
475 }
476
478}
479
480
481
482unsigned int
483ReferenceCell::unv_vertex_to_deal_vertex(const unsigned int vertex_n) const
484{
485 AssertIndexRange(vertex_n, n_vertices());
486 // Information on this file format isn't easy to find - the documents here
487 //
488 // https://www.ceas3.uc.edu/sdrluff/
489 //
490 // Don't actually explain anything about the sections we care about (2412) in
491 // any detail. For node numbering I worked backwards from what is actually in
492 // our test files (since that's supposed to work), which all use some
493 // non-standard clockwise numbering scheme which starts at the bottom right
494 // vertex.
495 if (*this == ReferenceCells::Line)
496 {
497 return vertex_n;
498 }
499 else if (*this == ReferenceCells::Quadrilateral)
500 {
501 constexpr std::array<unsigned int, 4> unv_to_deal{{1, 0, 2, 3}};
502 return unv_to_deal[vertex_n];
503 }
504 else if (*this == ReferenceCells::Hexahedron)
505 {
506 constexpr std::array<unsigned int, 8> unv_to_deal{
507 {6, 7, 5, 4, 2, 3, 1, 0}};
508 return unv_to_deal[vertex_n];
509 }
510
512
514}
515
516
517
518unsigned int
520{
521 switch (this->kind)
522 {
524 return VTKCellType::VTK_VERTEX;
526 return VTKCellType::VTK_LINE;
528 return VTKCellType::VTK_TRIANGLE;
530 return VTKCellType::VTK_QUAD;
532 return VTKCellType::VTK_TETRA;
534 return VTKCellType::VTK_PYRAMID;
536 return VTKCellType::VTK_WEDGE;
538 return VTKCellType::VTK_HEXAHEDRON;
540 return VTKCellType::VTK_INVALID;
541 default:
543 }
544
545 return VTKCellType::VTK_INVALID;
546}
547
548
549
550unsigned int
552{
553 switch (this->kind)
554 {
556 return VTKCellType::VTK_VERTEX;
558 return VTKCellType::VTK_QUADRATIC_EDGE;
560 return VTKCellType::VTK_QUADRATIC_TRIANGLE;
562 return VTKCellType::VTK_QUADRATIC_QUAD;
564 return VTKCellType::VTK_QUADRATIC_TETRA;
566 return VTKCellType::VTK_QUADRATIC_PYRAMID;
568 return VTKCellType::VTK_QUADRATIC_WEDGE;
570 return VTKCellType::VTK_QUADRATIC_HEXAHEDRON;
572 return VTKCellType::VTK_INVALID;
573 default:
575 }
576
577 return VTKCellType::VTK_INVALID;
578}
579
580
581
582unsigned int
584{
585 switch (this->kind)
586 {
588 return VTKCellType::VTK_VERTEX;
590 return VTKCellType::VTK_LAGRANGE_CURVE;
592 return VTKCellType::VTK_LAGRANGE_TRIANGLE;
594 return VTKCellType::VTK_LAGRANGE_QUADRILATERAL;
596 return VTKCellType::VTK_LAGRANGE_TETRAHEDRON;
598 return VTKCellType::VTK_LAGRANGE_PYRAMID;
600 return VTKCellType::VTK_LAGRANGE_WEDGE;
602 return VTKCellType::VTK_LAGRANGE_HEXAHEDRON;
604 return VTKCellType::VTK_INVALID;
605 default:
607 }
608
609 return VTKCellType::VTK_INVALID;
610}
611
612
613
614template <>
615unsigned int
617 const std::array<unsigned, 0> &,
618 const std::array<unsigned, 0> &,
619 const bool) const
620{
622 return 0;
623}
624
625
626
631template <>
632unsigned int
634 const std::array<unsigned, 1> &node_indices,
635 const std::array<unsigned, 1> &nodes_per_direction,
636 const bool) const
637{
638 const unsigned int i = node_indices[0];
639
640 const bool ibdy = (i == 0 || i == nodes_per_direction[0]);
641 // How many boundaries do we lie on at once?
642 const int nbdy = (ibdy ? 1 : 0);
643
644 if (nbdy == 1) // Vertex DOF
645 { // ijk is a corner node. Return the proper index (somewhere in [0,7]):
646 return i ? 1 : 0;
647 }
648
649 const int offset = 2;
650 return (i - 1) + offset;
651}
652
653
654
659template <>
660unsigned int
662 const std::array<unsigned, 2> &node_indices,
663 const std::array<unsigned, 2> &nodes_per_direction,
664 const bool) const
665{
667
668 const unsigned int i = node_indices[0];
669 const unsigned int j = node_indices[1];
670
671 const bool ibdy = (i == 0 || i == nodes_per_direction[0]);
672 const bool jbdy = (j == 0 || j == nodes_per_direction[1]);
673 // How many boundaries do we lie on at once?
674 const int nbdy = (ibdy ? 1 : 0) + (jbdy ? 1 : 0);
675
676 if (nbdy == 2) // Vertex DOF
677 { // ijk is a corner node. Return the proper index (somewhere in [0,3]):
678 return (i != 0u ? (j != 0u ? 2 : 1) : (j != 0u ? 3 : 0));
679 }
680
681 int offset = 4;
682 if (nbdy == 1) // Edge DOF
683 {
684 if (!ibdy)
685 { // On i axis
686 return (i - 1) +
687 (j != 0u ?
688 nodes_per_direction[0] - 1 + nodes_per_direction[1] - 1 :
689 0) +
690 offset;
691 }
692
693 if (!jbdy)
694 { // On j axis
695 return (j - 1) +
696 (i != 0u ? nodes_per_direction[0] - 1 :
697 2 * (nodes_per_direction[0] - 1) +
698 nodes_per_direction[1] - 1) +
699 offset;
700 }
701 }
702
703 offset += 2 * (nodes_per_direction[0] - 1 + nodes_per_direction[1] - 1);
704 // nbdy == 0: Face DOF
705 return offset + (i - 1) + (nodes_per_direction[0] - 1) * ((j - 1));
706}
707
708
709
720template <>
721unsigned int
723 const std::array<unsigned, 3> &node_indices,
724 const std::array<unsigned, 3> &nodes_per_direction,
725 const bool legacy_format) const
726{
728
729 const unsigned int i = node_indices[0];
730 const unsigned int j = node_indices[1];
731 const unsigned int k = node_indices[2];
732
733 const bool ibdy = (i == 0 || i == nodes_per_direction[0]);
734 const bool jbdy = (j == 0 || j == nodes_per_direction[1]);
735 const bool kbdy = (k == 0 || k == nodes_per_direction[2]);
736 // How many boundaries do we lie on at once?
737 const int nbdy = (ibdy ? 1 : 0) + (jbdy ? 1 : 0) + (kbdy ? 1 : 0);
738
739 if (nbdy == 3) // Vertex DOF
740 { // ijk is a corner node. Return the proper index (somewhere in [0,7]):
741 return (i != 0u ? (j != 0u ? 2 : 1) : (j != 0u ? 3 : 0)) +
742 (k != 0u ? 4 : 0);
743 }
744
745 int offset = 8;
746 if (nbdy == 2) // Edge DOF
747 {
748 if (!ibdy)
749 { // On i axis
750 return (i - 1) +
751 (j != 0u ?
752 nodes_per_direction[0] - 1 + nodes_per_direction[1] - 1 :
753 0) +
754 (k != 0u ? 2 * (nodes_per_direction[0] - 1 +
755 nodes_per_direction[1] - 1) :
756 0) +
757 offset;
758 }
759 if (!jbdy)
760 { // On j axis
761 return (j - 1) +
762 (i != 0u ? nodes_per_direction[0] - 1 :
763 2 * (nodes_per_direction[0] - 1) +
764 nodes_per_direction[1] - 1) +
765 (k != 0u ? 2 * (nodes_per_direction[0] - 1 +
766 nodes_per_direction[1] - 1) :
767 0) +
768 offset;
769 }
770 // !kbdy, On k axis
771 offset +=
772 4 * (nodes_per_direction[0] - 1) + 4 * (nodes_per_direction[1] - 1);
773 if (legacy_format)
774 return (k - 1) +
775 (nodes_per_direction[2] - 1) *
776 (i != 0u ? (j != 0u ? 3 : 1) : (j != 0u ? 2 : 0)) +
777 offset;
778 else
779 return (k - 1) +
780 (nodes_per_direction[2] - 1) *
781 (i != 0u ? (j != 0u ? 2 : 1) : (j != 0u ? 3 : 0)) +
782 offset;
783 }
784
785 offset += 4 * (nodes_per_direction[0] - 1 + nodes_per_direction[1] - 1 +
786 nodes_per_direction[2] - 1);
787 if (nbdy == 1) // Face DOF
788 {
789 if (ibdy) // On i-normal face
790 {
791 return (j - 1) + ((nodes_per_direction[1] - 1) * (k - 1)) +
792 (i != 0u ? (nodes_per_direction[1] - 1) *
793 (nodes_per_direction[2] - 1) :
794 0) +
795 offset;
796 }
797 offset += 2 * (nodes_per_direction[1] - 1) * (nodes_per_direction[2] - 1);
798 if (jbdy) // On j-normal face
799 {
800 return (i - 1) + ((nodes_per_direction[0] - 1) * (k - 1)) +
801 (j != 0u ? (nodes_per_direction[2] - 1) *
802 (nodes_per_direction[0] - 1) :
803 0) +
804 offset;
805 }
806 offset += 2 * (nodes_per_direction[2] - 1) * (nodes_per_direction[0] - 1);
807 // kbdy, On k-normal face
808 return (i - 1) + ((nodes_per_direction[0] - 1) * (j - 1)) +
809 (k != 0u ?
810 (nodes_per_direction[0] - 1) * (nodes_per_direction[1] - 1) :
811 0) +
812 offset;
813 }
814
815 // nbdy == 0: Body DOF
816 offset += 2 * ((nodes_per_direction[1] - 1) * (nodes_per_direction[2] - 1) +
817 (nodes_per_direction[2] - 1) * (nodes_per_direction[0] - 1) +
818 (nodes_per_direction[0] - 1) * (nodes_per_direction[1] - 1));
819 return offset + (i - 1) +
820 (nodes_per_direction[0] - 1) *
821 ((j - 1) + (nodes_per_direction[1] - 1) * ((k - 1)));
822}
823
824
825
826unsigned int
827ReferenceCell::vtk_vertex_to_deal_vertex(const unsigned int vertex_index) const
828{
829 AssertIndexRange(vertex_index, n_vertices());
830
831 // For some of the following, deal.II uses the same ordering as VTK
832 // and in that case, we only need to return 'vertex_index' (i.e.,
833 // use the identity mapping). For some others, we need to translate.
834 //
835 // For the ordering, see the VTK manual (for example at
836 // http://www.princeton.edu/~efeibush/viscourse/vtk.pdf, page 9).
837 switch (this->kind)
838 {
840 return vertex_index;
842 return vertex_index;
844 return vertex_index;
846 {
847 static constexpr std::array<unsigned int, 4> index_translation_table =
848 {{0, 1, 3, 2}};
849 return index_translation_table[vertex_index];
850 }
852 return vertex_index;
854 {
855 static constexpr std::array<unsigned int, 5> index_translation_table =
856 {{0, 1, 3, 2, 4}};
857 return index_translation_table[vertex_index];
858 }
860 return vertex_index;
862 {
863 static constexpr std::array<unsigned int, 8> index_translation_table =
864 {{0, 1, 3, 2, 4, 5, 7, 6}};
865 return index_translation_table[vertex_index];
866 }
868 {
871 }
872 default:
874 }
875
877}
878
879
880
881unsigned int
883{
884 /*
885 From the GMSH documentation:
886
887 elm-type
888 defines the geometrical type of the n-th element:
889
890 1
891 Line (2 nodes).
892
893 2
894 Triangle (3 nodes).
895
896 3
897 Quadrangle (4 nodes).
898
899 4
900 Tetrahedron (4 nodes).
901
902 5
903 Hexahedron (8 nodes).
904
905 6
906 Prism (6 nodes).
907
908 7
909 Pyramid (5 nodes).
910
911 8
912 Second order line (3 nodes: 2 associated with the vertices and 1 with the
913 edge).
914
915 9
916 Second order triangle (6 nodes: 3 associated with the vertices and 3 with
917 the edges).
918
919 10 Second order quadrangle (9 nodes: 4 associated with the
920 vertices, 4 with the edges and 1 with the face).
921
922 11 Second order tetrahedron (10 nodes: 4 associated with the vertices and 6
923 with the edges).
924
925 12 Second order hexahedron (27 nodes: 8 associated with the vertices, 12
926 with the edges, 6 with the faces and 1 with the volume).
927
928 13 Second order prism (18 nodes: 6 associated with the vertices, 9 with the
929 edges and 3 with the quadrangular faces).
930
931 14 Second order pyramid (14 nodes: 5 associated with the vertices, 8 with
932 the edges and 1 with the quadrangular face).
933
934 15 Point (1 node).
935 */
936
937 switch (this->kind)
938 {
940 return 15;
942 return 1;
944 return 2;
946 return 3;
948 return 4;
950 return 7;
952 return 6;
954 return 5;
956 default:
958 }
959
961}
962
963
964
965namespace
966{
967 // Compute the nearest point to @p on the line segment and the square of its
968 // distance to @p.
969 template <int dim>
970 std::pair<Point<dim>, double>
971 project_to_line(const Point<dim> &x0,
972 const Point<dim> &x1,
973 const Point<dim> &p)
974 {
975 Assert(x0 != x1, ExcInternalError());
976 // t is the convex combination coefficient (x = (1 - t) * x0 + t * x1)
977 // defining the position of the closest point on the line (not line segment)
978 // to p passing through x0 and x1. This formula is equivalent to the
979 // standard 'project a vector onto another vector', where each vector is
980 // shifted to start at x0.
981 const double t = ((x1 - x0) * (p - x0)) / ((x1 - x0).norm_square());
982
983 // Only consider points between x0 and x1
984 if (t <= 0)
985 return std::make_pair(x0, x0.distance_square(p));
986 else if (t <= 1)
987 {
988 const auto p2 = x0 + t * (x1 - x0);
989 return std::make_pair(p2, p2.distance_square(p));
990 }
991 else
992 return std::make_pair(x1, x1.distance_square(p));
993 }
994
995 // template base-case
996 template <int dim>
997 std::pair<Point<dim>, double>
998 project_to_quad(const std::array<Point<dim>, 3> & /*vertices*/,
999 const Point<dim> & /*p*/,
1000 const ReferenceCell /*reference_cell*/)
1001 {
1003 return std::make_pair(Point<dim>(),
1004 std::numeric_limits<double>::signaling_NaN());
1005 }
1006
1024 template <>
1025 std::pair<Point<3>, double>
1026 project_to_quad(const std::array<Point<3>, 3> &vertices,
1027 const Point<3> &p,
1028 const ReferenceCell face_reference_cell)
1029 {
1030 Assert(face_reference_cell == ReferenceCells::Triangle ||
1031 face_reference_cell == ReferenceCells::Quadrilateral,
1033
1034 // Make the problem slightly easier by shifting everything to avoid a point
1035 // at the origin (this way we can invert the matrix of vertices). Use 2.0 so
1036 // that the bottom left vertex of a Pyramid is now at x = 1.
1037 std::array<Point<3>, 3> shifted_vertices = vertices;
1038 const Tensor<1, 3> shift{{2.0, 2.0, 2.0}};
1039 for (Point<3> &shifted_vertex : shifted_vertices)
1040 shifted_vertex += shift;
1041 const Point<3> shifted_p = p + shift;
1042
1043 // As we are projecting onto a face of a reference cell, the vectors
1044 // describing its local coordinate system should be orthogonal. We don't
1045 // know which of the three vectors computed from p are mutually orthogonal
1046 // for triangles so that case requires an extra check.
1047 Tensor<1, 3> e0;
1048 Tensor<1, 3> e1;
1049 const Point<3> vertex = shifted_vertices[0];
1050 // Triangles are difficult because of two cases:
1051 // 1. the top face of a Tetrahedron, which does not have a right angle
1052 // 2. wedges and pyramids, whose faces do not lie on the reference simplex
1053 //
1054 // Deal with both by creating a locally orthogonal (but not necessarily
1055 // orthonormal) coordinate system and testing if the projected point is in
1056 // the triangle by expressing it as a convex combination of the vertices.
1057 if (face_reference_cell == ReferenceCells::Triangle)
1058 {
1059 e0 = shifted_vertices[1] - shifted_vertices[0];
1060 e1 = shifted_vertices[2] - shifted_vertices[0];
1061 e1 -= (e0 * e1) * e0 / (e0.norm_square());
1062 }
1063 else
1064 {
1065 e0 = shifted_vertices[1] - shifted_vertices[0];
1066 e1 = shifted_vertices[2] - shifted_vertices[0];
1067 }
1068 Assert(std::abs(e0 * e1) <= 1e-14, ExcInternalError());
1069 // the quadrilaterals on pyramids and wedges don't necessarily have edge
1070 // lengths of 1 so we cannot skip the denominator
1071 const double c0 = e0 * (shifted_p - vertex) / e0.norm_square();
1072 const double c1 = e1 * (shifted_p - vertex) / e1.norm_square();
1073 const Point<3> projected_shifted_p = vertex + c0 * e0 + c1 * e1;
1074
1075 bool in_quad = false;
1076 if (face_reference_cell == ReferenceCells::Triangle)
1077 {
1078 Tensor<2, 3> shifted_vertex_matrix;
1079 for (unsigned int i = 0; i < 3; ++i)
1080 shifted_vertex_matrix[i] = shifted_vertices[i];
1081 const Tensor<1, 3> combination_coordinates =
1082 invert(transpose(shifted_vertex_matrix)) * projected_shifted_p;
1083 bool is_convex_combination = true;
1084 for (unsigned int i = 0; i < 3; ++i)
1085 is_convex_combination = is_convex_combination &&
1086 (0.0 <= combination_coordinates[i]) &&
1087 (combination_coordinates[i] <= 1.0);
1088 in_quad = is_convex_combination;
1089 }
1090 else
1091 in_quad = (0.0 <= c0 && c0 <= 1.0 && 0.0 <= c1 && c1 <= 1.0);
1092
1093 if (in_quad)
1094 return std::make_pair(projected_shifted_p - shift,
1095 shifted_p.distance_square(projected_shifted_p));
1096 else
1097 return std::make_pair(Point<3>(), std::numeric_limits<double>::max());
1098 }
1099} // namespace
1100
1101
1102
1103template <int dim>
1106{
1108
1109 // Handle simple cases first:
1110 if (dim == 0)
1111 return p;
1112 if (contains_point(p, 0.0))
1113 return p;
1114 if (dim == 1)
1115 return project_to_line(vertex<dim>(0), vertex<dim>(1), p).first;
1116
1117 // Find the closest vertex so that we only need to check adjacent faces and
1118 // lines.
1119 Point<dim> result;
1120 unsigned int closest_vertex_no = 0;
1121 double closest_vertex_distance_square = vertex<dim>(0).distance_square(p);
1122 for (unsigned int i = 1; i < n_vertices(); ++i)
1123 {
1124 const double new_vertex_distance_square =
1125 vertex<dim>(i).distance_square(p);
1126 if (new_vertex_distance_square < closest_vertex_distance_square)
1127 {
1128 closest_vertex_no = i;
1129 closest_vertex_distance_square = new_vertex_distance_square;
1130 }
1131 }
1132
1133 double min_distance_square = std::numeric_limits<double>::max();
1134 if (dim == 2)
1135 {
1136 for (const unsigned int face_no :
1137 faces_for_given_vertex(closest_vertex_no))
1138 {
1139 const Point<dim> v0 = vertex<dim>(line_to_cell_vertices(face_no, 0));
1140 const Point<dim> v1 = vertex<dim>(line_to_cell_vertices(face_no, 1));
1141
1142 auto pair = project_to_line(v0, v1, p);
1143 if (pair.second < min_distance_square)
1144 {
1145 result = pair.first;
1146 min_distance_square = pair.second;
1147 }
1148 }
1149 }
1150 else
1151 {
1152 // Check faces and then lines.
1153 //
1154 // For reference cells with sloped faces (i.e., all 3D shapes except
1155 // Hexahedra), we might be able to do a valid normal projection to a face
1156 // with a different slope which is on the 'other side' of the reference
1157 // cell. To catch that case we have to unconditionally check lines after
1158 // checking faces.
1159 //
1160 // For pyramids the closest vertex might not be on the closest face: for
1161 // example, the origin is closest to vertex 4 which is not on the bottom
1162 // plane. Get around that by just checking all faces for pyramids.
1163 const std::array<unsigned int, 5> all_pyramid_faces{{0, 1, 2, 3, 4}};
1164 const auto &faces = *this == ReferenceCells::Pyramid ?
1165 ArrayView<const unsigned int>(all_pyramid_faces) :
1166 faces_for_given_vertex(closest_vertex_no);
1167 for (const unsigned int face_no : faces)
1168 {
1169 auto face_cell = face_reference_cell(face_no);
1170 // We only need the first three points since for quads the last point
1171 // is redundant
1172 std::array<Point<dim>, 3> vertices;
1173 for (unsigned int vertex_no = 0; vertex_no < 3; ++vertex_no)
1174 vertices[vertex_no] = vertex<dim>(face_to_cell_vertices(
1175 face_no, vertex_no, numbers::default_geometric_orientation));
1176
1177 auto pair = project_to_quad(vertices, p, face_cell);
1178 if (pair.second < min_distance_square)
1179 {
1180 result = pair.first;
1181 min_distance_square = pair.second;
1182 }
1183 }
1184
1185 for (const unsigned int face_no :
1186 faces_for_given_vertex(closest_vertex_no))
1187 {
1188 auto face_cell = face_reference_cell(face_no);
1189 for (const unsigned int face_line_no : face_cell.line_indices())
1190 {
1191 const auto cell_line_no =
1192 face_to_cell_lines(face_no,
1193 face_line_no,
1195 const auto v0 =
1196 vertex<dim>(line_to_cell_vertices(cell_line_no, 0));
1197 const auto v1 =
1198 vertex<dim>(line_to_cell_vertices(cell_line_no, 1));
1199 auto pair = project_to_line(v0, v1, p);
1200 if (pair.second < min_distance_square)
1201 {
1202 result = pair.first;
1203 min_distance_square = pair.second;
1204 }
1205 }
1206 }
1207 }
1208
1209 Assert(min_distance_square < std::numeric_limits<double>::max(),
1211
1212 // If necessary, slightly adjust the computed point so that it is closer to
1213 // being on the surface of the reference cell. Due to roundoff it is difficult
1214 // to place points on sloped surfaces (e.g., for Pyramids) so this check isn't
1215 // perfect but does improve the accuracy of the projected point.
1216 if (!contains_point(result, 0.0))
1217 {
1218 constexpr unsigned int x_index = 0;
1219 constexpr unsigned int y_index = (dim >= 2 ? 1 : 0);
1220 constexpr unsigned int z_index = (dim >= 3 ? 2 : 0);
1221 switch (this->kind)
1222 {
1225 break;
1226 // the bounds for each dimension of a hypercube are mutually
1227 // independent:
1231 for (unsigned int d = 0; d < dim; ++d)
1232 result[d] = std::clamp(result[d], 0.0, 1.0);
1233 // simplices can use the standard definition of a simplex:
1234 break;
1236 result[x_index] = std::clamp(result[x_index], 0.0, 1.0);
1237 result[y_index] =
1238 std::clamp(result[y_index], 0.0, 1.0 - result[x_index]);
1239 break;
1241 result[x_index] = std::clamp(result[x_index], 0.0, 1.0);
1242 result[y_index] =
1243 std::clamp(result[y_index], 0.0, 1.0 - result[x_index]);
1244 result[z_index] =
1245 std::clamp(result[z_index],
1246 0.0,
1247 1.0 - result[x_index] - result[y_index]);
1248 break;
1249 // wedges and pyramids are more ad-hoc:
1251 result[x_index] = std::clamp(result[x_index], 0.0, 1.0);
1252 result[y_index] =
1253 std::clamp(result[y_index], 0.0, 1.0 - result[x_index]);
1254 result[z_index] = std::clamp(result[z_index], 0.0, 1.0);
1255 break;
1257 {
1258 result[x_index] = std::clamp(result[x_index], -1.0, 1.0);
1259 result[y_index] = std::clamp(result[y_index], -1.0, 1.0);
1260 // It suffices to transform everything to the first quadrant to
1261 // adjust z:
1262 const auto x_abs = std::abs(result[x_index]);
1263 const auto y_abs = std::abs(result[y_index]);
1264
1265 if (y_abs <= x_abs)
1266 result[z_index] = std::clamp(result[z_index], 0.0, 1.0 - x_abs);
1267 else
1268 result[z_index] = std::clamp(result[z_index], 0.0, 1.0 - y_abs);
1269 }
1270 break;
1271 default:
1273 }
1274 }
1275 // We should be within 4 * eps of the cell by this point. The roundoff error
1276 // comes from, e.g., computing (1 - x) + x when moving points onto the top of
1277 // a Pyramid.
1278 Assert(contains_point(result, 4.0 * std::numeric_limits<double>::epsilon()),
1280
1281 return result;
1282}
1283
1284
1285
1286std::ostream &
1287operator<<(std::ostream &out, const ReferenceCell &reference_cell)
1288{
1289 AssertThrow(out.fail() == false, ExcIO());
1290
1291 // Output as an integer to avoid outputting it as a character with
1292 // potentially non-printing value:
1293 out << static_cast<unsigned int>(reference_cell.kind);
1294 return out;
1295}
1296
1297
1298
1299std::istream &
1300operator>>(std::istream &in, ReferenceCell &reference_cell)
1301{
1302 AssertThrow(in.fail() == false, ExcIO());
1303
1304 // Read the information as an integer and convert it to the correct type
1305 unsigned int value;
1306 in >> value;
1307 reference_cell.kind = static_cast<decltype(reference_cell.kind)>(value);
1308
1309 // Ensure that the object we read is valid
1310 Assert(
1311 (reference_cell == ReferenceCells::Vertex) ||
1312 (reference_cell == ReferenceCells::Line) ||
1313 (reference_cell == ReferenceCells::Triangle) ||
1314 (reference_cell == ReferenceCells::Quadrilateral) ||
1315 (reference_cell == ReferenceCells::Tetrahedron) ||
1316 (reference_cell == ReferenceCells::Hexahedron) ||
1317 (reference_cell == ReferenceCells::Wedge) ||
1318 (reference_cell == ReferenceCells::Pyramid) ||
1319 (reference_cell == ReferenceCells::Invalid),
1320 ExcMessage(
1321 "The reference cell kind just read does not correspond to one of the "
1322 "valid choices. There must be an error."));
1323
1324 return in;
1325}
1326
1327// explicitly instantiate dimension 0 quadrature in addition to the standard
1328// dimensions
1329template Quadrature<0>
1330ReferenceCell::get_gauss_type_quadrature(const unsigned n_points_1D) const;
1331template const Quadrature<0> &
1333
1334#include "grid/reference_cell.inst"
1335
Implementation of the classic affine transformation mapping used for simplices.
Definition mapping_p1.h:73
Abstract base class for mapping classes.
Definition mapping.h:320
Definition point.h:113
constexpr numbers::NumberTraits< Number >::real_type distance_square(const Point< dim, Number > &p) const
ArrayView< const unsigned int > faces_for_given_vertex(const unsigned int vertex_index) const
unsigned int n_vertices() const
Quadrature< dim > get_gauss_type_quadrature(const unsigned n_points_1d) const
bool is_hyper_cube() const
unsigned int vtk_linear_type() const
std::pair< unsigned int, RefinementCase< dim - 1 > > equivalent_refinement_case(const types::geometric_orientation combined_face_orientation, const internal::SubfaceCase< dim > subface_case, const unsigned int subface_no) const
Point< dim > vertex(const unsigned int v) const
const Quadrature< dim > & get_nodal_type_quadrature() const
unsigned int face_to_cell_lines(const unsigned int face, const unsigned int line, const types::geometric_orientation face_orientation) const
static constexpr ndarray< unsigned int, 2, 2 > line_vertex_permutations
unsigned int vtk_lexicographic_to_node_index(const std::array< unsigned, dim > &node_indices, const std::array< unsigned, dim > &nodes_per_direction, const bool legacy_format) const
bool contains_point(const Point< dim > &p, const double tolerance=0) const
unsigned int gmsh_element_type() const
static constexpr ndarray< unsigned int, 8, 4 > quadrilateral_vertex_permutations
unsigned int n_faces() const
unsigned int line_to_cell_vertices(const unsigned int line, const unsigned int vertex) const
std::uint8_t kind
constexpr ReferenceCell()
unsigned int exodusii_vertex_to_deal_vertex(const unsigned int vertex_n) const
unsigned int unv_vertex_to_deal_vertex(const unsigned int vertex_n) const
unsigned int vtk_vertex_to_deal_vertex(const unsigned int vertex_index) const
unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const types::geometric_orientation face_orientation) const
std::unique_ptr< Mapping< dim, spacedim > > get_default_mapping(const unsigned int degree) const
unsigned int vtk_quadratic_type() const
unsigned int vtk_lagrange_type() const
unsigned int get_dimension() const
ReferenceCell face_reference_cell(const unsigned int face_no) const
unsigned int exodusii_face_to_deal_face(const unsigned int face_n) const
const Mapping< dim, spacedim > & get_default_linear_mapping() const
std::string to_string() const
bool is_simplex() const
Point< dim > closest_point(const Point< dim > &p) const
static constexpr ndarray< unsigned int, 6, 3 > triangle_vertex_permutations
constexpr numbers::NumberTraits< Number >::real_type norm_square() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
const bool IsBlockVector< VectorType >::value
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
void shift(const Tensor< 1, spacedim > &shift_vector, Triangulation< dim, spacedim > &triangulation)
constexpr ReferenceCell Triangle
constexpr ReferenceCell Hexahedron
constexpr ReferenceCell Pyramid
constexpr ReferenceCell Wedge
constexpr ReferenceCell Vertex
constexpr ReferenceCell Invalid
constexpr ReferenceCell Quadrilateral
constexpr ReferenceCell Tetrahedron
constexpr ReferenceCell Line
std::tuple< bool, bool, bool > split_face_orientation(const types::geometric_orientation combined_face_orientation)
constexpr unsigned int invalid_unsigned_int
Definition types.h:238
constexpr types::geometric_orientation default_geometric_orientation
Definition types.h:352
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned char geometric_orientation
Definition types.h:40
typename internal::ndarray::HelperArray< T, Ns... >::type ndarray
Definition ndarray.h:107
std::istream & operator>>(std::istream &in, ReferenceCell &reference_cell)
std::ostream & operator<<(std::ostream &out, const ReferenceCell &reference_cell)
static unsigned int child_cell_on_face(const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim - 1 > &face_refinement_case=RefinementCase< dim - 1 >::isotropic_refinement)
static constexpr unsigned int max_children_per_face
static MappingQ< dim, spacedim > mapping
Definition mapping_q1.h:97
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)