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
qprojector.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
18
21
22#include <boost/container/small_vector.hpp>
23
24
26
27
28namespace internal
29{
30 namespace QProjector
31 {
32 namespace
33 {
34 // Reflect points across the y = x line.
35 std::vector<Point<2>>
36 reflect(const std::vector<Point<2>> &points)
37 {
38 std::vector<Point<2>> q_points;
39 q_points.reserve(points.size());
40 for (const Point<2> &p : points)
41 q_points.emplace_back(p[1], p[0]);
42
43 return q_points;
44 }
45
46
47 // Rotate points in the plane around the positive z axis (i.e.,
48 // counter-clockwise).
49 std::vector<Point<2>>
50 rotate(const std::vector<Point<2>> &points, const unsigned int n_times)
51 {
52 std::vector<Point<2>> q_points;
53 q_points.reserve(points.size());
54 switch (n_times % 4)
55 {
56 case 0:
57 // 0 degree. the point remains as it is.
58 for (const Point<2> &p : points)
59 q_points.push_back(p);
60 break;
61 case 1:
62 // 90 degree counterclockwise
63 for (const Point<2> &p : points)
64 q_points.emplace_back(1.0 - p[1], p[0]);
65 break;
66 case 2:
67 // 180 degree counterclockwise
68 for (const Point<2> &p : points)
69 q_points.emplace_back(1.0 - p[0], 1.0 - p[1]);
70 break;
71 case 3:
72 // 270 degree counterclockwise
73 for (const Point<2> &p : points)
74 q_points.emplace_back(p[1], 1.0 - p[0]);
75 break;
76 }
77
78 return q_points;
79 }
80
87 void
88 project_to_hex_face_and_append(
89 const std::vector<Point<2>> &points,
90 const unsigned int face_no,
91 std::vector<Point<3>> &q_points,
93 const unsigned int subface_no = 0)
94 {
95 // one coordinate is at a const value. for faces 0, 2 and 4 this value
96 // is 0.0, for faces 1, 3 and 5 it is 1.0
97 const double const_value = face_no % 2;
98
99 // local 2d coordinates are xi and eta, global 3d coordinates are x, y
100 // and z. those have to be mapped. the following indices tell, which
101 // global coordinate (0->x, 1->y, 2->z) corresponds to which local one
102 const unsigned int xi_index = (1 + face_no / 2) % 3,
103 eta_index = (2 + face_no / 2) % 3,
104 const_index = face_no / 2;
105
106 // for a standard face (no refinement), we use the default values of
107 // the xi and eta scales and translations, otherwise the xi and eta
108 // values will be scaled (by factor 0.5 or factor 1.0) depending on
109 // the refinement case and translated (by 0.0 or 0.5) depending on the
110 // refinement case and subface_no
111 double xi_scale = 1.0, eta_scale = 1.0, xi_translation = 0.0,
112 eta_translation = 0.0;
113
114 // set the scale and translation parameter for individual subfaces
115 switch (ref_case)
116 {
118 break;
120 xi_scale = 0.5;
121 xi_translation = subface_no % 2 * 0.5;
122 break;
124 eta_scale = 0.5;
125 eta_translation = subface_no % 2 * 0.5;
126 break;
128 xi_scale = 0.5;
129 eta_scale = 0.5;
130 xi_translation = int(subface_no % 2) * 0.5;
131 eta_translation = int(subface_no / 2) * 0.5;
132 break;
133 default:
135 break;
136 }
137
138 // finally, compute the scaled, translated, projected quadrature
139 // points
140 for (const Point<2> &p : points)
141 {
142 Point<3> cell_point;
143 cell_point[xi_index] = xi_scale * p[0] + xi_translation;
144 cell_point[eta_index] = eta_scale * p[1] + eta_translation;
145 cell_point[const_index] = const_value;
146 q_points.push_back(cell_point);
147 }
148 }
149
150 std::vector<Point<2>>
151 mutate_points_with_offset(
152 const std::vector<Point<2>> &points,
153 const types::geometric_orientation combined_orientation)
154 {
155 // These rotations are backwards (relative to the standard notion of,
156 // e.g., what rotation index 7 means) since they are rotations about the
157 // positive z axis in 2d: i.e., they are done from the perspective of
158 // 'inside' a cell instead of the perspective of an abutting cell.
159 //
160 // For example: consider points on face 4 of a hexahedron with
161 // orientation 3. In 2d, rotating such points clockwise is the same as
162 // rotating them counter-clockwise from the perspective of the abutting
163 // face. Hence, such points must be rotated 90 degrees
164 // counter-clockwise.
165 switch (combined_orientation)
166 {
167 case 1:
168 return reflect(points);
169 case 3:
170 return rotate(reflect(points), 3);
171 case 5:
172 return rotate(reflect(points), 2);
173 case 7:
174 return rotate(reflect(points), 1);
175 case 0:
176 return points;
177 case 2:
178 return rotate(points, 1);
179 case 4:
180 return rotate(points, 2);
181 case 6:
182 return rotate(points, 3);
183 default:
185 }
186 return {};
187 }
188
189
190
217 template <int dim>
218 void
219 append_subobject_rule(
220 const ReferenceCell &face_reference_cell,
221 const Quadrature<dim - 1> &quadrature,
222 const std::vector<Point<dim>> &vertices,
223 const double measure,
224 const types::geometric_orientation combined_orientation,
225 std::vector<Point<dim>> &points,
226 std::vector<double> &weights)
227 {
228 const auto support_points =
229 face_reference_cell.permute_by_combined_orientation(
230 make_const_array_view(vertices),
231 face_reference_cell.get_inverse_combined_orientation(
232 combined_orientation));
233
234 for (unsigned int j = 0; j < quadrature.size(); ++j)
235 {
236 Point<dim> mapped_point;
237
238 // map reference quadrature point
239 for (const unsigned int vertex_no :
240 face_reference_cell.vertex_indices())
241 mapped_point +=
242 support_points[vertex_no] *
243 face_reference_cell.d_linear_shape_function(quadrature.point(j),
244 vertex_no);
245
246 points.push_back(mapped_point);
247
248 // rescale quadrature weights so that the sum of the weights on
249 // each face equals the measure of that face.
250 weights.push_back(quadrature.weight(j) * measure /
251 face_reference_cell.volume());
252 }
253 }
254 } // namespace
255 } // namespace QProjector
256} // namespace internal
257
258
259
260template <>
261void
263 const Quadrature<0> &quadrature,
264 const unsigned int face_no,
265 std::vector<Point<1>> &q_points)
266{
267 AssertDimension(quadrature.size(), q_points.size());
268 const auto face_quadrature =
269 QProjector<1>::project_to_face(reference_cell,
270 quadrature,
271 face_no,
273 q_points = face_quadrature.get_points();
274}
275
276
277
278template <>
279void
281 const Quadrature<1> &quadrature,
282 const unsigned int face_no,
283 std::vector<Point<2>> &q_points)
284{
285 AssertDimension(quadrature.size(), q_points.size());
286 const auto face_quadrature =
287 QProjector<2>::project_to_face(reference_cell,
288 quadrature,
289 face_no,
291 q_points = face_quadrature.get_points();
292}
293
294
295
296template <>
297void
299 const Quadrature<2> &quadrature,
300 const unsigned int face_no,
301 std::vector<Point<3>> &q_points)
302{
304 (void)reference_cell;
305
307 Assert(q_points.size() == quadrature.size(),
308 ExcDimensionMismatch(q_points.size(), quadrature.size()));
309 q_points.clear();
310 internal::QProjector::project_to_hex_face_and_append(quadrature.get_points(),
311 face_no,
312 q_points);
313}
314
315
316
317template <int dim>
320 const Quadrature<dim - 1> &quadrature,
321 const unsigned int face_no,
322 const bool face_orientation,
323 const bool face_flip,
324 const bool face_rotation)
325{
327 reference_cell,
328 quadrature,
329 face_no,
331 face_rotation,
332 face_flip));
333}
334
335
336
337template <int dim>
340 const ReferenceCell &reference_cell,
341 const Quadrature<dim - 1> &quadrature,
342 const unsigned int face_no,
343 const types::geometric_orientation combined_orientation)
344{
345 AssertIndexRange(face_no, reference_cell.n_faces());
346 AssertIndexRange(combined_orientation,
347 reference_cell.n_face_orientations(face_no));
348 AssertDimension(reference_cell.get_dimension(), dim);
349
350 std::vector<Point<dim>> points;
351 std::vector<double> weights;
352
353 const ReferenceCell face_reference_cell =
354 reference_cell.face_reference_cell(face_no);
355 std::vector<Point<dim>> face_vertices(face_reference_cell.n_vertices());
356 for (const unsigned int vertex_no : face_reference_cell.vertex_indices())
357 face_vertices[vertex_no] =
358 reference_cell.face_vertex_location<dim>(face_no, vertex_no);
359 internal::QProjector::append_subobject_rule(face_reference_cell,
360 quadrature,
361 face_vertices,
362 reference_cell.face_measure(
363 face_no),
364 combined_orientation,
365 points,
366 weights);
367
368 return Quadrature<dim>(std::move(points), std::move(weights));
369}
370
371
372
373template <>
374void
376 const Quadrature<0> &,
377 const unsigned int face_no,
378 const unsigned int,
379 std::vector<Point<1>> &q_points,
380 const RefinementCase<0> &)
381{
382 Assert(reference_cell == ReferenceCells::Line, ExcNotImplemented());
383 (void)reference_cell;
384
385 const unsigned int dim = 1;
387 AssertDimension(q_points.size(), 1);
388
389 q_points[0] = Point<dim>(static_cast<double>(face_no));
390}
391
392
393
394template <>
395void
397 const Quadrature<1> &quadrature,
398 const unsigned int face_no,
399 const unsigned int subface_no,
400 std::vector<Point<2>> &q_points,
401 const RefinementCase<1> &ref_case)
402{
403 AssertDimension(quadrature.size(), q_points.size());
404 const auto face_quadrature =
405 project_to_subface(reference_cell,
406 quadrature,
407 face_no,
408 subface_no,
410 ref_case);
411 q_points = face_quadrature.get_points();
412}
413
414
415
416template <>
417void
419 const Quadrature<2> &quadrature,
420 const unsigned int face_no,
421 const unsigned int subface_no,
422 std::vector<Point<3>> &q_points,
423 const RefinementCase<2> &ref_case)
424{
426 (void)reference_cell;
427 AssertDimension(quadrature.size(), q_points.size());
428 const auto face_quadrature =
429 project_to_subface(reference_cell,
430 quadrature,
431 face_no,
432 subface_no,
434 ref_case);
435 q_points = face_quadrature.get_points();
436}
437
438
439
440template <int dim>
443 const ReferenceCell &reference_cell,
444 const Quadrature<dim - 1> &quadrature,
445 const unsigned int face_no,
446 const unsigned int subface_no,
447 const bool,
448 const bool,
449 const bool,
451{
453 reference_cell,
454 quadrature,
455 face_no,
456 subface_no,
458}
459
460
461
462template <int dim>
465 const ReferenceCell &reference_cell,
466 const SubQuadrature &quadrature,
467 const unsigned int face_no,
468 const unsigned int subface_no,
469 const types::geometric_orientation combined_orientation,
470 const RefinementCase<dim - 1> &ref_case)
471{
472 AssertIndexRange(face_no, reference_cell.n_faces());
473 AssertIndexRange(combined_orientation,
474 reference_cell.n_face_orientations(face_no));
475 AssertDimension(reference_cell.get_dimension(), dim);
476 AssertIndexRange(subface_no,
477 reference_cell.face_reference_cell(face_no)
478 .template n_children<dim - 1>(ref_case));
479
480 std::vector<Point<dim>> q_points;
481 std::vector<double> q_weights = quadrature.get_weights();
482 q_points.reserve(quadrature.size());
483
484 if constexpr (dim == 1)
485 {
486 AssertDimension(quadrature.size(), 1);
487 q_points.emplace_back(static_cast<double>(face_no));
488 }
489 else if constexpr (dim == 2)
490 {
491 if (reference_cell == ReferenceCells::Triangle)
492 // use linear polynomial to map the reference quadrature points
493 // correctly on faces, i.e., BarycentricPolynomials<1>(1)
494 for (unsigned int p = 0; p < quadrature.size(); ++p)
495 {
496 if (face_no == 0)
497 {
498 if (subface_no == 0)
499 q_points.emplace_back(quadrature.point(p)[0] / 2, 0);
500 else
501 q_points.emplace_back(0.5 + quadrature.point(p)[0] / 2, 0);
502 }
503 else if (face_no == 1)
504 {
505 if (subface_no == 0)
506 q_points.emplace_back(1 - quadrature.point(p)[0] / 2,
507 quadrature.point(p)[0] / 2);
508 else
509 q_points.emplace_back(0.5 - quadrature.point(p)[0] / 2,
510 0.5 + quadrature.point(p)[0] / 2);
511 }
512 else if (face_no == 2)
513 {
514 if (subface_no == 0)
515 q_points.emplace_back(0, 1 - quadrature.point(p)[0] / 2);
516 else
517 q_points.emplace_back(0, 0.5 - quadrature.point(p)[0] / 2);
518 }
519 else
521 }
522 else if (reference_cell == ReferenceCells::Quadrilateral)
523 for (unsigned int p = 0; p < quadrature.size(); ++p)
524 {
525 if (face_no == 0)
526 {
527 if (subface_no == 0)
528 q_points.emplace_back(0, quadrature.point(p)[0] / 2);
529 else
530 q_points.emplace_back(0, quadrature.point(p)[0] / 2 + 0.5);
531 }
532 else if (face_no == 1)
533 {
534 if (subface_no == 0)
535 q_points.emplace_back(1, quadrature.point(p)[0] / 2);
536 else
537 q_points.emplace_back(1, quadrature.point(p)[0] / 2 + 0.5);
538 }
539 else if (face_no == 2)
540 {
541 if (subface_no == 0)
542 q_points.emplace_back(quadrature.point(p)[0] / 2, 0);
543 else
544 q_points.emplace_back(quadrature.point(p)[0] / 2 + 0.5, 0);
545 }
546 else if (face_no == 3)
547 {
548 if (subface_no == 0)
549 q_points.emplace_back(quadrature.point(p)[0] / 2, 1);
550 else
551 q_points.emplace_back(quadrature.point(p)[0] / 2 + 0.5, 1);
552 }
553 else
555 }
556 else
558
559 if (combined_orientation == numbers::reverse_line_orientation)
560 {
561 std::reverse(q_points.begin(), q_points.end());
562 std::reverse(q_weights.begin(), q_weights.end());
563 }
564 for (auto &w : q_weights)
565 w *= reference_cell.face_measure(face_no);
566 }
567 else if constexpr (dim == 3)
568 {
570 internal::QProjector::project_to_hex_face_and_append(
571 quadrature.get_points(), face_no, q_points, ref_case, subface_no);
572 }
573 else
574 {
576 }
577
578 return Quadrature<dim>(std::move(q_points), std::move(q_weights));
579}
580
581
582
583template <int dim>
586 const ReferenceCell &reference_cell,
587 const hp::QCollection<dim - 1> &quadrature)
588{
589 std::vector<Point<dim>> points;
590 std::vector<double> weights;
591
592 for (const unsigned int face_no : reference_cell.face_indices())
593 {
594 const ReferenceCell face_reference_cell =
595 reference_cell.face_reference_cell(face_no);
596 std::vector<Point<dim>> face_vertices(face_reference_cell.n_vertices());
597 for (const unsigned int vertex_no : face_reference_cell.vertex_indices())
598 face_vertices[vertex_no] =
599 reference_cell.face_vertex_location<dim>(face_no, vertex_no);
600
601 for (types::geometric_orientation combined_orientation = 0;
602 combined_orientation < reference_cell.n_face_orientations(face_no);
603 ++combined_orientation)
604 internal::QProjector::append_subobject_rule(
605 face_reference_cell,
606 quadrature[quadrature.size() == 1 ? 0 : face_no],
607 face_vertices,
608 reference_cell.face_measure(face_no),
609 combined_orientation,
610 points,
611 weights);
612 }
613
614 return Quadrature<dim>(std::move(points), std::move(weights));
615}
616
617
618
619template <>
622 const Quadrature<0> &quadrature)
623{
624 Assert(reference_cell == ReferenceCells::Line, ExcNotImplemented());
625 (void)reference_cell;
626
627 const unsigned int dim = 1;
628
629 const unsigned int n_points = 1, n_faces = GeometryInfo<dim>::faces_per_cell,
630 subfaces_per_face =
632
633 // first fix quadrature points
634 std::vector<Point<dim>> q_points;
635 q_points.reserve(n_points * n_faces * subfaces_per_face);
636 std::vector<Point<dim>> help(n_points);
637
638 // project to each face and copy
639 // results
640 for (unsigned int face = 0; face < n_faces; ++face)
641 for (unsigned int subface = 0; subface < subfaces_per_face; ++subface)
642 {
643 project_to_subface(reference_cell, quadrature, face, subface, help);
644 std::copy(help.begin(), help.end(), std::back_inserter(q_points));
645 }
646
647 // next copy over weights
648 std::vector<double> weights;
649 weights.reserve(n_points * n_faces * subfaces_per_face);
650 for (unsigned int face = 0; face < n_faces; ++face)
651 for (unsigned int subface = 0; subface < subfaces_per_face; ++subface)
652 std::copy(quadrature.get_weights().begin(),
653 quadrature.get_weights().end(),
654 std::back_inserter(weights));
655
656 Assert(q_points.size() == n_points * n_faces * subfaces_per_face,
658 Assert(weights.size() == n_points * n_faces * subfaces_per_face,
660
661 return Quadrature<dim>(std::move(q_points), std::move(weights));
662}
663
664
665
666template <>
669 const SubQuadrature &quadrature)
670{
671 Assert(reference_cell == ReferenceCells::Quadrilateral ||
672 reference_cell == ReferenceCells::Triangle,
674
675 const unsigned int dim = 2;
676
677 std::vector<Point<dim>> q_points;
678 std::vector<double> weights;
679
680 // project to each face and copy
681 // results
682 for (unsigned int face = 0; face < reference_cell.n_faces(); ++face)
683 for (types::geometric_orientation orientation = 0;
684 orientation < reference_cell.n_face_orientations(face);
685 ++orientation)
686 for (unsigned int subface = 0;
687 subface <
688 reference_cell.face_reference_cell(face).n_isotropic_children();
689 ++subface)
690 {
691 const unsigned int local_subface =
692 orientation == numbers::reverse_line_orientation ? 1 - subface :
693 subface;
694 const auto sub_quadrature =
695 project_to_subface(reference_cell,
696 quadrature,
697 face,
698 local_subface,
699 orientation,
701 q_points.insert(q_points.end(),
702 sub_quadrature.get_points().begin(),
703 sub_quadrature.get_points().end());
704 weights.insert(weights.end(),
705 sub_quadrature.get_weights().begin(),
706 sub_quadrature.get_weights().end());
707 }
708
709 return Quadrature<dim>(std::move(q_points), std::move(weights));
710}
711
712
713
714template <>
717 const SubQuadrature &quadrature)
718{
719 if (reference_cell == ReferenceCells::Triangle ||
720 reference_cell == ReferenceCells::Tetrahedron)
721 return Quadrature<3>(); // nothing to do
722
724
725 const unsigned int dim = 3;
726
727 const unsigned int n_points = quadrature.size(),
729 total_subfaces_per_face = 2 + 2 + 4;
730
731 // first fix quadrature points
732 std::vector<Point<dim>> q_points;
733 q_points.reserve(n_points * n_faces * total_subfaces_per_face * 8);
734
735 std::vector<double> weights;
736 weights.reserve(n_points * n_faces * total_subfaces_per_face * 8);
737
738 // do the following for all possible mutations of a face (mutation==0
739 // corresponds to a face with standard orientation, no flip and no rotation)
740 for (unsigned char offset = 0; offset < 8; ++offset)
741 {
742 const auto mutation =
743 internal::QProjector::mutate_points_with_offset(quadrature.get_points(),
744 offset);
745
746 // project to each face and copy results
747 for (unsigned int face = 0; face < n_faces; ++face)
748 for (unsigned int ref_case = RefinementCase<dim - 1>::cut_xy;
749 ref_case >= RefinementCase<dim - 1>::cut_x;
750 --ref_case)
751 for (unsigned int subface = 0;
753 RefinementCase<dim - 1>(ref_case));
754 ++subface)
755 {
756 internal::QProjector::project_to_hex_face_and_append(
757 mutation,
758 face,
759 q_points,
760 RefinementCase<dim - 1>(ref_case),
761 subface);
762
763 // next copy over weights
764 std::copy(quadrature.get_weights().begin(),
765 quadrature.get_weights().end(),
766 std::back_inserter(weights));
767 }
768 }
769
770 Assert(q_points.size() == n_points * n_faces * total_subfaces_per_face * 8,
772 Assert(weights.size() == n_points * n_faces * total_subfaces_per_face * 8,
774
775 return Quadrature<dim>(std::move(q_points), std::move(weights));
776}
777
778
779
780template <int dim>
783 const Quadrature<dim> &quadrature,
784 const unsigned int child_no)
785{
786 Assert(reference_cell == ReferenceCells::get_hypercube<dim>(),
788 (void)reference_cell;
789
791
792 const unsigned int n_q_points = quadrature.size();
793
794 std::vector<Point<dim>> q_points(n_q_points);
795 for (unsigned int i = 0; i < n_q_points; ++i)
796 q_points[i] =
798 child_no);
799
800 // for the weights, things are
801 // equally simple: copy them and
802 // scale them
803 std::vector<double> weights = quadrature.get_weights();
804 for (unsigned int i = 0; i < n_q_points; ++i)
805 weights[i] *= (1. / GeometryInfo<dim>::max_children_per_cell);
806
807 return Quadrature<dim>(q_points, weights);
808}
809
810
811
812template <int dim>
815 const Quadrature<dim> &quadrature)
816{
817 Assert(reference_cell == ReferenceCells::get_hypercube<dim>(),
819 (void)reference_cell;
820
821 const unsigned int n_points = quadrature.size(),
823
824 std::vector<Point<dim>> q_points(n_points * n_children);
825 std::vector<double> weights(n_points * n_children);
826
827 // project to each child and copy
828 // results
829 for (unsigned int child = 0; child < n_children; ++child)
830 {
831 Quadrature<dim> help =
832 project_to_child(reference_cell, quadrature, child);
833 for (unsigned int i = 0; i < n_points; ++i)
834 {
835 q_points[child * n_points + i] = help.point(i);
836 weights[child * n_points + i] = help.weight(i);
837 }
838 }
839 return Quadrature<dim>(q_points, weights);
840}
841
842
843
844template <int dim>
847 const Quadrature<1> &quadrature,
848 const Point<dim> &p1,
849 const Point<dim> &p2)
850{
851 Assert(reference_cell == ReferenceCells::get_hypercube<dim>(),
853 (void)reference_cell;
854
855 const unsigned int n = quadrature.size();
856 std::vector<Point<dim>> points(n);
857 std::vector<double> weights(n);
858 const double length = p1.distance(p2);
859
860 for (unsigned int k = 0; k < n; ++k)
861 {
862 const double alpha = quadrature.point(k)[0];
863 points[k] = alpha * p2;
864 points[k] += (1. - alpha) * p1;
865 weights[k] = length * quadrature.weight(k);
866 }
867 return Quadrature<dim>(points, weights);
868}
869
870
871
872template <int dim>
875 const unsigned int face_no,
876 const bool face_orientation,
877 const bool face_flip,
878 const bool face_rotation,
879 const unsigned int n_quadrature_points)
880{
881 return face(reference_cell,
882 face_no,
884 face_rotation,
885 face_flip),
886 n_quadrature_points);
887}
888
889
890
891template <int dim>
894 const ReferenceCell &reference_cell,
895 const unsigned int face_no,
896 const types::geometric_orientation combined_orientation,
897 const unsigned int n_quadrature_points)
898{
899 AssertIndexRange(face_no, reference_cell.n_faces());
900 AssertIndexRange(combined_orientation,
901 reference_cell.n_face_orientations(face_no));
902 AssertDimension(reference_cell.get_dimension(), dim);
903
904
905 return {(reference_cell.n_face_orientations(face_no) * face_no +
906 combined_orientation) *
907 n_quadrature_points};
908}
909
910
911
912template <int dim>
915 const ReferenceCell &reference_cell,
916 const unsigned int face_no,
917 const bool face_orientation,
918 const bool face_flip,
919 const bool face_rotation,
920 const hp::QCollection<dim - 1> &quadrature)
921{
922 return face(reference_cell,
923 face_no,
925 face_rotation,
926 face_flip),
927 quadrature);
928}
929
930
931
932template <int dim>
935 const ReferenceCell &reference_cell,
936 const unsigned int face_no,
937 const types::geometric_orientation combined_orientation,
938 const hp::QCollection<dim - 1> &quadrature)
939{
940 AssertIndexRange(face_no, reference_cell.n_faces());
941 AssertIndexRange(combined_orientation,
942 reference_cell.n_face_orientations(face_no));
943 AssertDimension(reference_cell.get_dimension(), dim);
944
945 unsigned int offset = 0;
946 if (quadrature.size() == 1)
947 offset =
948 reference_cell.n_face_orientations(0) * quadrature[0].size() * face_no;
949 else
950 for (unsigned int i = 0; i < face_no; ++i)
951 offset += reference_cell.n_face_orientations(i) * quadrature[i].size();
952
953 return {offset + combined_orientation *
954 quadrature[quadrature.size() == 1 ? 0 : face_no].size()};
955}
956
957
958
959template <int dim>
962 const ReferenceCell &reference_cell,
963 const unsigned int face_no,
964 const unsigned int subface_no,
965 const bool face_orientation,
966 const bool face_flip,
967 const bool face_rotation,
968 const unsigned int n_quadrature_points,
969 const internal::SubfaceCase<dim> ref_case)
970{
972 reference_cell,
973 face_no,
974 subface_no,
976 face_rotation,
977 face_flip),
978 n_quadrature_points,
979 ref_case);
980}
981
982
983
984template <>
987 const ReferenceCell &reference_cell,
988 const unsigned int face_no,
989 const unsigned int subface_no,
990 const types::geometric_orientation /*combined_orientation*/,
991 const unsigned int n_quadrature_points,
993{
994 Assert(reference_cell == ReferenceCells::Line, ExcNotImplemented());
995 (void)reference_cell;
996
1000
1001 return ((face_no * GeometryInfo<1>::max_children_per_face + subface_no) *
1002 n_quadrature_points);
1003}
1004
1005
1006
1007template <>
1010 const ReferenceCell &reference_cell,
1011 const unsigned int face_no,
1012 const unsigned int subface_no,
1013 const types::geometric_orientation combined_orientation,
1014 const unsigned int n_quadrature_points,
1016{
1017 Assert(reference_cell == ReferenceCells::Quadrilateral ||
1018 reference_cell == ReferenceCells::Triangle,
1020
1021 const unsigned int n_faces = reference_cell.n_faces();
1022 const unsigned int n_subfaces =
1023 reference_cell.face_reference_cell(face_no).n_isotropic_children();
1024 const unsigned int n_orientations =
1025 reference_cell.n_face_orientations(face_no);
1026
1027 AssertIndexRange(face_no, n_faces);
1028 AssertIndexRange(subface_no, n_subfaces);
1029 AssertIndexRange(combined_orientation, n_orientations);
1030
1031 return ((face_no * n_orientations * n_subfaces +
1032 combined_orientation * n_subfaces + subface_no) *
1033 n_quadrature_points);
1034}
1035
1036
1037
1038template <>
1041 const ReferenceCell &reference_cell,
1042 const unsigned int face_no,
1043 const unsigned int subface_no,
1044 const types::geometric_orientation combined_orientation,
1045 const unsigned int n_quadrature_points,
1046 const internal::SubfaceCase<3> ref_case)
1047{
1048 const unsigned int dim = 3;
1049
1051 (void)reference_cell;
1052
1056
1057 // in 3d, we have to account for faces that
1058 // have non-standard orientation. thus, we
1059 // have to store _eight_ data sets per face
1060 // or subface already for the isotropic
1061 // case. Additionally, we have three
1062 // different refinement cases, resulting in
1063 // <tt>4 + 2 + 2 = 8</tt> different subfaces
1064 // for each face.
1065 const unsigned int total_subfaces_per_face = 8;
1066
1067 // set up a table with the offsets for a
1068 // given refinement case respecting the
1069 // corresponding number of subfaces. the
1070 // index corresponds to (RefineCase::Type - 1)
1071
1072 // note, that normally we should use the
1073 // obvious offsets 0,2,6. However, prior to
1074 // the implementation of anisotropic
1075 // refinement, in many places of the library
1076 // the convention was used, that the first
1077 // dataset with offset 0 corresponds to a
1078 // standard (isotropic) face
1079 // refinement. therefore we use the offsets
1080 // 6,4,0 here to stick to that (implicit)
1081 // convention
1082 static const unsigned int ref_case_offset[3] = {
1083 6, // cut_x
1084 4, // cut_y
1085 0 // cut_xy
1086 };
1087
1088 const auto [final_subface_no, refinement_case] =
1089 reference_cell.equivalent_refinement_case(combined_orientation,
1090 ref_case,
1091 subface_no);
1092
1093 return ((face_no * total_subfaces_per_face +
1094 ref_case_offset[refinement_case - 1] + final_subface_no) +
1095 GeometryInfo<dim>::faces_per_cell * total_subfaces_per_face *
1096 combined_orientation) *
1097 n_quadrature_points;
1098}
1099
1100
1101
1102template <int dim>
1105 const SubQuadrature &quadrature,
1106 const unsigned int face_no)
1107{
1108 Assert(reference_cell == ReferenceCells::get_hypercube<dim>(),
1110 (void)reference_cell;
1111
1112 std::vector<Point<dim>> points(quadrature.size());
1113 project_to_face(reference_cell, quadrature, face_no, points);
1114 return Quadrature<dim>(points, quadrature.get_weights());
1115}
1116
1117
1118
1119template <int dim>
1122 const SubQuadrature &quadrature,
1123 const unsigned int face_no,
1124 const unsigned int subface_no,
1125 const RefinementCase<dim - 1> &ref_case)
1126{
1127 Assert(reference_cell == ReferenceCells::get_hypercube<dim>(),
1129 (void)reference_cell;
1130
1131 std::vector<Point<dim>> points(quadrature.size());
1133 reference_cell, quadrature, face_no, subface_no, points, ref_case);
1134 return Quadrature<dim>(points, quadrature.get_weights());
1135}
1136
1137
1138// explicit instantiations; note: we need them all for all dimensions
1139template class QProjector<1>;
1140template class QProjector<2>;
1141template class QProjector<3>;
1142
auto make_const_array_view(const Container &container) -> decltype(make_array_view(container))
Definition point.h:113
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
Class storing the offset index into a Quadrature rule created by project_to_all_faces() or project_to...
Definition qprojector.h:334
static DataSetDescriptor subface(const ReferenceCell &reference_cell, const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points, const internal::SubfaceCase< dim > ref_case=internal::SubfaceCase< dim >::case_isotropic)
static DataSetDescriptor face(const ReferenceCell &reference_cell, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
Class which transforms dim - 1-dimensional quadrature rules to dim-dimensional face quadratures.
Definition qprojector.h:69
static Quadrature< dim > project_to_all_subfaces(const ReferenceCell &reference_cell, const SubQuadrature &quadrature)
static Quadrature< dim > project_to_all_faces(const ReferenceCell &reference_cell, const hp::QCollection< dim - 1 > &quadrature)
static Quadrature< dim > project_to_child(const ReferenceCell &reference_cell, const Quadrature< dim > &quadrature, const unsigned int child_no)
static Quadrature< dim > project_to_oriented_subface(const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const internal::SubfaceCase< dim > ref_case)
static Quadrature< dim > project_to_oriented_face(const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation)
static Quadrature< dim > project_to_line(const ReferenceCell &reference_cell, const Quadrature< 1 > &quadrature, const Point< dim > &p1, const Point< dim > &p2)
static void project_to_subface(const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim > > &q_points, const RefinementCase< dim - 1 > &ref_case=RefinementCase< dim - 1 >::isotropic_refinement)
static Quadrature< dim > project_to_all_children(const ReferenceCell &reference_cell, const Quadrature< dim > &quadrature)
static void project_to_face(const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
Quadrature< dim - 1 > SubQuadrature
Definition qprojector.h:75
const Point< dim > & point(const unsigned int i) const
double weight(const unsigned int i) const
const std::vector< double > & get_weights() const
const std::vector< Point< dim > > & get_points() const
unsigned int size() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
unsigned int n_vertices() const
double d_linear_shape_function(const Point< dim > &xi, const unsigned int i) const
types::geometric_orientation get_inverse_combined_orientation(const types::geometric_orientation orientation) const
boost::container::small_vector< T, 8 > permute_by_combined_orientation(const ArrayView< const T > &vertices, const types::geometric_orientation orientation) const
double volume() const
unsigned int size() const
Definition collection.h:316
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ASSERT_UNREACHABLE()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
constexpr ReferenceCell Triangle
constexpr ReferenceCell Hexahedron
constexpr const ReferenceCell & get_hypercube()
constexpr ReferenceCell Quadrilateral
constexpr ReferenceCell Tetrahedron
constexpr ReferenceCell Line
types::geometric_orientation combined_face_orientation(const bool face_orientation, const bool face_rotation, const bool face_flip)
constexpr types::geometric_orientation reverse_line_orientation
Definition types.h:365
constexpr types::geometric_orientation default_geometric_orientation
Definition types.h:352
unsigned char geometric_orientation
Definition types.h:40
static constexpr unsigned int max_children_per_face
static constexpr unsigned int faces_per_cell
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
static Point< dim > child_to_cell_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)
static constexpr unsigned int max_children_per_cell