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
fe_simplex_p.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
15#include <deal.II/base/config.h>
16
20#include <deal.II/base/types.h>
21
22#include <deal.II/fe/fe_dgq.h>
24#include <deal.II/fe/fe_q.h>
26#include <deal.II/fe/fe_tools.h>
27#include <deal.II/fe/mapping.h>
28
31
33
34namespace
35{
36 // TODO: replace this with QProjector once QProjector learns how to project
37 // quadrature points onto separate faces of simplices
38 template <int dim>
40 face_midpoint(const unsigned int face_no)
41 {
43 const auto face_reference_cell =
44 reference_cell.face_reference_cell(face_no);
45
46 Point<dim> midpoint;
47 for (const auto face_vertex_no : face_reference_cell.vertex_indices())
48 {
49 const auto vertex_no = reference_cell.face_to_cell_vertices(
50 face_no, face_vertex_no, numbers::default_geometric_orientation);
51
52 midpoint += reference_cell.template vertex<dim>(vertex_no);
53 }
54
55 midpoint /= reference_cell.face_reference_cell(0).n_vertices();
56 return midpoint;
57 }
58
63 std::vector<unsigned int>
64 get_dpo_vector_fe_p(const unsigned int dim, const unsigned int degree)
65 {
66 switch (dim)
67 {
68 case 1:
69 switch (degree)
70 {
71 case 1:
72 return {1, 0};
73 case 2:
74 return {1, 1};
75 case 3:
76 return {1, 2};
77 default:
79 }
80 case 2:
81 switch (degree)
82 {
83 case 1:
84 return {1, 0, 0};
85 case 2:
86 return {1, 1, 0};
87 case 3:
88 return {1, 2, 1};
89 default:
91 }
92 case 3:
93 switch (degree)
94 {
95 case 1:
96 return {1, 0, 0, 0};
97 case 2:
98 return {1, 1, 0, 0};
99 case 3:
100 return {1, 2, 1, 0};
101 default:
103 }
104 }
105
107 return {};
108 }
109
110
116 template <int dim>
117 std::vector<Point<dim>>
118 unit_support_points_fe_p(const unsigned int degree)
119 {
120 Assert(dim != 0, ExcInternalError());
122 std::vector<Point<dim>> unit_points;
124
125 // Piecewise constants are a special case: use a support point at the
126 // centroid and only the centroid
127 if (degree == 0)
128 {
129 unit_points.emplace_back(reference_cell.template barycenter<dim>());
130 return unit_points;
131 }
132
133 // otherwise write everything as linear combinations of vertices
134 const auto dpo = get_dpo_vector_fe_p(dim, degree);
135 Assert(dpo.size() == dim + 1, ExcInternalError());
136 Assert(dpo[0] == 1, ExcNotImplemented());
137 // vertices:
138 for (const unsigned int d : reference_cell.vertex_indices())
139 unit_points.push_back(reference_cell.template vertex<dim>(d));
140 // lines:
141 for (const unsigned int l : reference_cell.line_indices())
142 {
143 const Point<dim> p0 =
144 unit_points[reference_cell.line_to_cell_vertices(l, 0)];
145 const Point<dim> p1 =
146 unit_points[reference_cell.line_to_cell_vertices(l, 1)];
147 for (unsigned int p = 0; p < dpo[1]; ++p)
148 unit_points.push_back((double(dpo[1] - p) / (dpo[1] + 1)) * p0 +
149 (double(p + 1) / (dpo[1] + 1)) * p1);
150 }
151 // quads:
152 if (dim > 1 && dpo[2] > 0)
153 {
154 Assert(dpo[2] == 1, ExcNotImplemented());
155 if (dim == 2)
156 unit_points.push_back(reference_cell.template barycenter<dim>());
157 if (dim == 3)
158 for (const unsigned int f : reference_cell.face_indices())
159 unit_points.push_back(face_midpoint<dim>(f));
160 }
161
162 return unit_points;
163 }
164
165 template <>
166 std::vector<Point<0>>
167 unit_support_points_fe_p(const unsigned int /*degree*/)
168 {
169 return {Point<0>()};
170 }
171
176 template <int dim>
177 std::vector<std::vector<Point<dim - 1>>>
178 unit_face_support_points_fe_p(
179 const unsigned int degree,
180 typename FiniteElementData<dim>::Conformity conformity)
181 {
182 // Discontinuous elements don't have face support points
184 return {};
185
186 // this concept doesn't exist in 1d so just return an empty vector
187 if (dim == 1)
188 return {};
189
190 std::vector<std::vector<Point<dim - 1>>> unit_face_points;
191
192 // all faces have the same support points
193 for (auto face_n : ReferenceCells::get_simplex<dim>().face_indices())
194 {
195 (void)face_n;
196 unit_face_points.emplace_back(
197 unit_support_points_fe_p<dim - 1>(degree));
198 }
199
200 return unit_face_points;
201 }
202
208 template <int dim>
210 constraints_fe_p(const unsigned int /*degree*/)
211 {
212 // no constraints in 1d
213 // constraints in 3d not implemented yet
214 return FullMatrix<double>();
215 }
216
217 template <>
219 constraints_fe_p<2>(const unsigned int degree)
220 {
221 constexpr int dim = 2;
222
223 // the following implements the 2d case
224 // (the 3d case is not implemented yet)
225 //
226 // consult FE_Q_Base::Implementation::initialize_constraints()
227 // for more information
228
229 std::vector<Point<dim - 1>> constraint_points;
230 // midpoint
231 constraint_points.emplace_back(0.5);
232 // subface 0
233 for (unsigned int i = 1; i < degree; ++i)
234 constraint_points.push_back(
236 Point<dim - 1>(i / double(degree)), 0));
237 // subface 1
238 for (unsigned int i = 1; i < degree; ++i)
239 constraint_points.push_back(
241 Point<dim - 1>(i / double(degree)), 1));
242
243 // Now construct relation between destination (child) and source (mother)
244 // dofs.
245
246 const unsigned int n_dofs_constrained = constraint_points.size();
247 unsigned int n_dofs_per_face = degree + 1;
248 FullMatrix<double> interface_constraints(n_dofs_constrained,
249 n_dofs_per_face);
250
251 const auto poly = BarycentricPolynomials<dim - 1>::get_fe_p_basis(degree);
252
253 for (unsigned int i = 0; i < n_dofs_constrained; ++i)
254 for (unsigned int j = 0; j < n_dofs_per_face; ++j)
255 {
256 interface_constraints(i, j) =
257 poly.compute_value(j, constraint_points[i]);
258
259 // if the value is small up to round-off, then simply set it to zero
260 // to avoid unwanted fill-in of the constraint matrices (which would
261 // then increase the number of other DoFs a constrained DoF would
262 // couple to)
263 if (std::fabs(interface_constraints(i, j)) < 1e-13)
264 interface_constraints(i, j) = 0;
265 }
266 return interface_constraints;
267 }
268
269
270
275 std::vector<unsigned int>
276 get_dpo_vector_fe_dgp(const unsigned int dim, const unsigned int degree)
277 {
278 // First treat the case of piecewise constant elements:
279 if (degree == 0)
280 {
281 std::vector<unsigned int> dpo(dim + 1, 0U);
282 dpo[dim] = 1;
283 return dpo;
284 }
285 else
286 {
287 // This element has the same degrees of freedom as the continuous one,
288 // but they are all counted for the interior of the cell because
289 // it is continuous. Rather than hard-code how many DoFs the element
290 // has, we just get the numbers from the continuous case and add them
291 // up
292 const auto continuous_dpo = get_dpo_vector_fe_p(dim, degree);
293
294 switch (dim)
295 {
296 case 1:
297 return {0U,
298 ReferenceCells::Line.n_vertices() * continuous_dpo[0] +
299 continuous_dpo[dim]};
300
301 case 2:
302 return {0U,
303 0U,
304 ReferenceCells::Triangle.n_vertices() *
305 continuous_dpo[0] +
306 ReferenceCells::Triangle.n_lines() * continuous_dpo[1] +
307 continuous_dpo[dim]};
308
309 case 3:
310 return {
311 0U,
312 0U,
313 0U,
314 ReferenceCells::Tetrahedron.n_vertices() * continuous_dpo[0] +
315 ReferenceCells::Tetrahedron.n_lines() * continuous_dpo[1] +
316 ReferenceCells::Tetrahedron.n_faces() * continuous_dpo[2] +
317 continuous_dpo[dim]};
318 }
319
321 return {};
322 }
323 }
324} // namespace
325
326
327
328template <int dim, int spacedim>
330 const BarycentricPolynomials<dim> polynomials,
331 const FiniteElementData<dim> &fe_data,
332 const bool prolongation_is_additive,
333 const std::vector<Point<dim>> &unit_support_points,
334 const std::vector<std::vector<Point<dim - 1>>> unit_face_support_points,
336 : ::FE_Poly<dim, spacedim>(
337 polynomials,
338 fe_data,
339 std::vector<bool>(fe_data.dofs_per_cell, prolongation_is_additive),
340 std::vector<ComponentMask>(fe_data.dofs_per_cell,
341 ComponentMask(std::vector<bool>(1, true))))
342{
345 this->interface_constraints = interface_constraints;
346}
347
348
349
350template <int dim, int spacedim>
351std::pair<Table<2, bool>, std::vector<unsigned int>>
353{
354 Table<2, bool> constant_modes(1, this->n_dofs_per_cell());
355 constant_modes.fill(true);
356 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
357 constant_modes, std::vector<unsigned int>(1, 0));
358}
359
360
361
362template <int dim, int spacedim>
363const FullMatrix<double> &
365 const unsigned int child,
366 const RefinementCase<dim> &refinement_case) const
367{
368 if (dim == 3)
369 Assert(RefinementCase<dim>(refinement_case) ==
371 static_cast<char>(IsotropicRefinementChoice::cut_tet_68)) ||
372 RefinementCase<dim>(refinement_case) ==
374 static_cast<char>(IsotropicRefinementChoice::cut_tet_57)) ||
375 RefinementCase<dim>(refinement_case) ==
377 static_cast<char>(IsotropicRefinementChoice::cut_tet_49)),
379 else
380 Assert(refinement_case ==
383 AssertDimension(dim, spacedim);
384
385 // initialization upon first request
386 if (this->prolongation[refinement_case - 1][child].n() == 0)
387 {
388 std::lock_guard<std::mutex> lock(prolongation_matrix_mutex);
389
390 // if matrix got updated while waiting for the lock
391 if (this->prolongation[refinement_case - 1][child].n() ==
392 this->n_dofs_per_cell())
393 return this->prolongation[refinement_case - 1][child];
394
395 // now do the work. need to get a non-const version of data in order to
396 // be able to modify them inside a const function
397 auto &this_nonconst = const_cast<FE_SimplexPoly<dim, spacedim> &>(*this);
398
399 if (dim == 2)
400 {
401 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
403 isotropic_matrices.back().resize(
404 this->reference_cell().n_children(
405 RefinementCase<dim>(refinement_case)),
407 this->n_dofs_per_cell()));
408
409 FETools::compute_embedding_matrices(*this, isotropic_matrices, true);
410
411 this_nonconst.prolongation[refinement_case - 1] =
412 std::move(isotropic_matrices.back());
413 }
414 else if (dim == 3)
415 {
416 std::vector<std::vector<FullMatrix<double>>> matrices(
417 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49),
418 std::vector<FullMatrix<double>>(
419 this->reference_cell().n_children(
420 RefinementCase<dim>(refinement_case)),
422 this->n_dofs_per_cell())));
423 FETools::compute_embedding_matrices(*this, matrices, true);
424 for (unsigned int refinement_direction = static_cast<unsigned int>(
426 refinement_direction <=
427 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49);
428 refinement_direction++)
429 this_nonconst.prolongation[refinement_direction - 1] =
430 std::move(matrices[refinement_direction - 1]);
431 }
432 else
434 }
435
436 // finally return the matrix
437 return this->prolongation[refinement_case - 1][child];
438}
439
440
441
442template <int dim, int spacedim>
443unsigned int
445 const unsigned int face_dof_index,
446 const unsigned int face,
447 const types::geometric_orientation combined_orientation) const
448{
449 // Function largely lifted from FE_Q_Base::face_to_cell_index()
450 AssertIndexRange(face_dof_index, this->n_dofs_per_face(face));
451 AssertIndexRange(face, this->reference_cell().n_faces());
452
453 AssertIndexRange(combined_orientation,
454 this->reference_cell().n_face_orientations(face));
455
456 // we need to distinguish between DoFs on vertices, lines and in 3d quads.
457 // do so in a sequence of if-else statements
458 if (face_dof_index < this->get_first_face_line_index(face))
459 // DoF is on a vertex
460 {
461 // get the number of the vertex on the face that corresponds to this DoF,
462 // along with the number of the DoF on this vertex
463 const unsigned int face_vertex =
464 face_dof_index / this->n_dofs_per_vertex();
465 const unsigned int dof_index_on_vertex =
466 face_dof_index % this->n_dofs_per_vertex();
467
468 // then get the number of this vertex on the cell and translate
469 // this to a DoF number on the cell
470 return this->reference_cell().face_to_cell_vertices(
471 face, face_vertex, combined_orientation) *
472 this->n_dofs_per_vertex() +
473 dof_index_on_vertex;
474 }
475 else if (face_dof_index < this->get_first_face_quad_index(face))
476 // DoF is on a line
477 {
478 // do the same kind of translation as before. we need to only consider
479 // DoFs on the lines, i.e., ignoring those on the vertices
480 const unsigned int index =
481 face_dof_index - this->get_first_face_line_index(face);
482
483 const unsigned int face_line = index / this->n_dofs_per_line();
484 const unsigned int dof_index_on_line = index % this->n_dofs_per_line();
485
486 // we now also need to adjust the line index for the case of
487 // face orientation, face flips, etc
488 unsigned int adjusted_dof_index_on_line = 0;
489 switch (dim)
490 {
491 case 1:
493 break;
494
495 case 2:
496 if (combined_orientation == numbers::default_geometric_orientation)
497 adjusted_dof_index_on_line = dof_index_on_line;
498 else
499 adjusted_dof_index_on_line =
500 this->n_dofs_per_line() - dof_index_on_line - 1;
501 break;
502
503 case 3:
504 {
505 // Line orientations are not consistent between faces of
506 // tetrahedra: e.g., the cell's first line is (0, 1) on face 0 and
507 // (1, 0) on face 1. Hence we have to determine the relative line
508 // orientation to determine how to index dofs along lines.
509 //
510 // face_to_cell_line_orientation() may only be called with
511 // canonical (face, line) pairs. Extract that information and then
512 // also the new canonical face_line_index.
513 const auto cell_line_index =
514 this->reference_cell().face_to_cell_lines(face,
515 face_line,
516 combined_orientation);
517 const auto [canonical_face, canonical_line] =
518 this->reference_cell().standard_line_to_face_and_line_index(
519 cell_line_index);
520 // Here we don't take into account what the actual orientation of
521 // the line is: that's usually handled inside, e.g.,
522 // face->get_dof_indices().
523 const auto face_line_orientation =
524 this->reference_cell().face_to_cell_line_orientation(
525 canonical_line,
526 canonical_face,
527 combined_orientation,
528 /*combined_line_orientation =*/
530
531 if (face_line_orientation ==
533 adjusted_dof_index_on_line = dof_index_on_line;
534 else
535 {
536 Assert(face_line_orientation ==
539 adjusted_dof_index_on_line =
540 this->n_dofs_per_line() - dof_index_on_line - 1;
541 }
542 }
543 break;
544
545 default:
547 }
548
549 return (this->get_first_line_index() +
550 this->reference_cell().face_to_cell_lines(face,
551 face_line,
552 combined_orientation) *
553 this->n_dofs_per_line() +
554 adjusted_dof_index_on_line);
555 }
556 else
557 // DoF is on a quad
558 {
559 Assert(dim >= 3, ExcInternalError());
560
561 // ignore vertex and line dofs
562 const unsigned int index =
563 face_dof_index - this->get_first_face_quad_index(face);
564
565 // Since polynomials on simplices are typically defined in barycentric
566 // coordinates, handling higher degree elements is not that hard: if the
567 // DoFs come in triples then they simply need to be rotated in the correct
568 // way (as though they were the vertices). However, we don't implement
569 // higher than cubic elements, so that is not yet done here. For example,
570 // for a P4 element a quad will have 3 DoFs, and each of those DoFs will
571 // be on its corresponding vertex's bisector. Hence they can be rotated by
572 // simply calling
573 // ReferenceCells::Triangle::permute_by_combined_orientation().
574 Assert((this->n_dofs_per_quad(face) <= 1) ||
575 combined_orientation == numbers::default_geometric_orientation,
577 return (this->get_first_quad_index(face) + index);
578 }
579}
580
581
582
583template <int dim, int spacedim>
584const FullMatrix<double> &
586 const unsigned int child,
587 const RefinementCase<dim> &refinement_case) const
588{
589 if (dim == 3)
590 Assert(RefinementCase<dim>(refinement_case) ==
592 static_cast<char>(IsotropicRefinementChoice::cut_tet_68)) ||
593 RefinementCase<dim>(refinement_case) ==
595 static_cast<char>(IsotropicRefinementChoice::cut_tet_57)) ||
596 RefinementCase<dim>(refinement_case) ==
598 static_cast<char>(IsotropicRefinementChoice::cut_tet_49)),
600 else
603 AssertDimension(dim, spacedim);
604
605 // initialization upon first request
606 if (this->restriction[refinement_case - 1][child].n() == 0)
607 {
608 std::lock_guard<std::mutex> lock(restriction_matrix_mutex);
609
610 // if matrix got updated while waiting for the lock
611 if (this->restriction[refinement_case - 1][child].n() ==
612 this->n_dofs_per_cell())
613 return this->restriction[refinement_case - 1][child];
614
615 // get the restriction matrix
616 // Refine a unit cell. As the parent cell is a unit
617 // cell, the reference cell of the children equals the parent, i.e. they
618 // have the support points at the same locations. So we just have to check
619 // if a support point of the parent is one of the interpolation points of
620 // the child. If this is not the case we find the interpolation of the
621 // point.
622
623 const double eps = 1e-12;
624 FullMatrix<double> restriction_mat(this->n_dofs_per_cell(),
625 this->n_dofs_per_cell());
626
627 // first get all support points on the reference cell
628 const std::vector<Point<dim>> unit_support_points =
630
631 // now create children on the reference cell
634 tria.begin_active()->set_refine_flag(
636 if (dim == 3)
637 tria.begin_active()->set_refine_choice(refinement_case);
639
640 const auto &child_cell = tria.begin(0)->child(child);
641
642 // iterate over all support points and transform them to the unit cell of
643 // the child
644 for (unsigned int i = 0; i < unit_support_points.size(); i++)
645 {
646 std::vector<Point<dim>> transformed_point(1);
647 const std::vector<Point<spacedim>> unit_support_point = {
648 dim == 2 ? Point<spacedim>(unit_support_points[i][0],
649 unit_support_points[i][1]) :
652 unit_support_points[i][2])};
653 this->reference_cell()
656 child_cell,
658 make_array_view(transformed_point));
659
660 // if point is inside the unit cell iterate over all shape functions
661 if (this->reference_cell().contains_point(transformed_point[0], eps))
662 for (unsigned int j = 0; j < this->n_dofs_per_cell(); j++)
663 restriction_mat[i][j] =
664 this->shape_value(j, transformed_point[0]);
665 }
666 if constexpr (running_in_debug_mode())
667 {
668 for (unsigned int i = 0; i < this->n_dofs_per_cell(); i++)
669 {
670 double sum = 0.;
671
672 for (unsigned int j = 0; j < this->n_dofs_per_cell(); j++)
673 sum += restriction_mat[i][j];
674
675 Assert(std::fabs(sum - 1) < eps || std::fabs(sum) < eps,
677 "The entries in a row of the local "
678 "restriction matrix do not add to zero or one. "
679 "This typically indicates that the "
680 "polynomial interpolation is "
681 "ill-conditioned such that round-off "
682 "prevents the sum to be one."));
683 }
684 }
685
686 // Remove small entries from the matrix
687 for (unsigned int i = 0; i < restriction_mat.m(); ++i)
688 for (unsigned int j = 0; j < restriction_mat.n(); ++j)
689 {
690 if (std::fabs(restriction_mat(i, j)) < eps)
691 restriction_mat(i, j) = 0.;
692 if (std::fabs(restriction_mat(i, j) - 1) < eps)
693 restriction_mat(i, j) = 1.;
694 }
695
696 const_cast<FullMatrix<double> &>(
697 this->restriction[refinement_case - 1][child]) =
698 std::move(restriction_mat);
699 }
700
701 // finally return the matrix
702 return this->restriction[refinement_case - 1][child];
703}
704
705
706
707template <int dim, int spacedim>
708void
710 const FiniteElement<dim, spacedim> &source_fe,
711 FullMatrix<double> &interpolation_matrix,
712 const unsigned int face_no) const
713{
714 Assert(interpolation_matrix.m() == source_fe.n_dofs_per_face(face_no),
715 ExcDimensionMismatch(interpolation_matrix.m(),
716 source_fe.n_dofs_per_face(face_no)));
717
718 // see if source is a P or Q element
719 if ((dynamic_cast<const FE_SimplexPoly<dim, spacedim> *>(&source_fe) !=
720 nullptr) ||
721 (dynamic_cast<const FE_Q_Base<dim, spacedim> *>(&source_fe) != nullptr))
722 {
723 const Quadrature<dim - 1> quad_face_support(
724 source_fe.get_unit_face_support_points(face_no));
725
726 const double eps = 2e-13 * this->degree * (dim - 1);
727
728 const std::vector<Point<dim>> face_quadrature_points =
730 quad_face_support,
731 face_no,
733 .get_points();
734
735 for (unsigned int i = 0; i < source_fe.n_dofs_per_face(face_no); ++i)
736 for (unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
737 {
738 double matrix_entry =
739 this->shape_value(this->face_to_cell_index(j, 0),
740 face_quadrature_points[i]);
741
742 // Correct the interpolated value. I.e. if it is close to 1 or
743 // 0, make it exactly 1 or 0. Unfortunately, this is required to
744 // avoid problems with higher order elements.
745 if (std::fabs(matrix_entry - 1.0) < eps)
746 matrix_entry = 1.0;
747 if (std::fabs(matrix_entry) < eps)
748 matrix_entry = 0.0;
749
750 interpolation_matrix(i, j) = matrix_entry;
751 }
752
753 if constexpr (running_in_debug_mode())
754 {
755 for (unsigned int j = 0; j < source_fe.n_dofs_per_face(face_no); ++j)
756 {
757 double sum = 0.;
758
759 for (unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
760 sum += interpolation_matrix(j, i);
761
762 Assert(std::fabs(sum - 1) < eps, ExcInternalError());
763 }
764 }
765 }
766 else if (dynamic_cast<const FE_Nothing<dim> *>(&source_fe) != nullptr)
767 {
768 // nothing to do here, the FE_Nothing has no degrees of freedom anyway
769 }
770 else
772 false,
773 (typename FiniteElement<dim,
775}
776
777
778
779template <int dim, int spacedim>
780void
782 const FiniteElement<dim, spacedim> &source_fe,
783 const unsigned int subface,
784 FullMatrix<double> &interpolation_matrix,
785 const unsigned int face_no) const
786{
787 Assert(interpolation_matrix.m() == source_fe.n_dofs_per_face(face_no),
788 ExcDimensionMismatch(interpolation_matrix.m(),
789 source_fe.n_dofs_per_face(face_no)));
790
791 // see if source is a P or Q element
792 if ((dynamic_cast<const FE_SimplexPoly<dim, spacedim> *>(&source_fe) !=
793 nullptr) ||
794 (dynamic_cast<const FE_Q_Base<dim, spacedim> *>(&source_fe) != nullptr))
795 {
796 const Quadrature<dim - 1> quad_face_support(
797 source_fe.get_unit_face_support_points(face_no));
798
799 const double eps = 2e-13 * this->degree * (dim - 1);
800
801 const Quadrature<dim> subface_quadrature =
803 this->reference_cell(),
804 quad_face_support,
805 face_no,
806 subface,
809
810 for (unsigned int i = 0; i < source_fe.n_dofs_per_face(face_no); ++i)
811 for (unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
812 {
813 double matrix_entry =
814 this->shape_value(this->face_to_cell_index(j, 0),
815 subface_quadrature.point(i));
816
817 // Correct the interpolated value. I.e. if it is close to 1 or
818 // 0, make it exactly 1 or 0. Unfortunately, this is required to
819 // avoid problems with higher order elements.
820 if (std::fabs(matrix_entry - 1.0) < eps)
821 matrix_entry = 1.0;
822 if (std::fabs(matrix_entry) < eps)
823 matrix_entry = 0.0;
824
825 interpolation_matrix(i, j) = matrix_entry;
826 }
827
828 if constexpr (running_in_debug_mode())
829 {
830 for (unsigned int j = 0; j < source_fe.n_dofs_per_face(face_no); ++j)
831 {
832 double sum = 0.;
833
834 for (unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
835 sum += interpolation_matrix(j, i);
836
837 Assert(std::fabs(sum - 1) < eps, ExcInternalError());
838 }
839 }
840 }
841 else if (dynamic_cast<const FE_Nothing<dim> *>(&source_fe) != nullptr)
842 {
843 // nothing to do here, the FE_Nothing has no degrees of freedom anyway
844 }
845 else
847 false,
848 (typename FiniteElement<dim,
850}
851
852
853
854template <int dim, int spacedim>
855bool
860
861
862
863template <int dim, int spacedim>
864void
867 const std::vector<Vector<double>> &support_point_values,
868 std::vector<double> &nodal_values) const
869{
870 AssertDimension(support_point_values.size(),
871 this->get_unit_support_points().size());
872 AssertDimension(support_point_values.size(), nodal_values.size());
873 AssertDimension(this->dofs_per_cell, nodal_values.size());
874
875 for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
876 {
877 AssertDimension(support_point_values[i].size(), 1);
878
879 nodal_values[i] = support_point_values[i](0);
880 }
881}
882
883
884
885template <int dim, int spacedim>
887 : FE_SimplexPoly<dim, spacedim>(
888 BarycentricPolynomials<dim>::get_fe_p_basis(degree),
889 FiniteElementData<dim>(get_dpo_vector_fe_p(dim, degree),
890 ReferenceCells::get_simplex<dim>(),
891 1,
892 degree,
893 FiniteElementData<dim>::H1),
894 false,
895 unit_support_points_fe_p<dim>(degree),
896 unit_face_support_points_fe_p<dim>(degree, FiniteElementData<dim>::H1),
897 constraints_fe_p<dim>(degree))
898{
899 if (degree > 2)
900 for (unsigned int i = 0; i < this->n_dofs_per_line(); ++i)
902 this->n_dofs_per_line() - 1 - i - i;
903 // We do not support multiple DoFs per quad yet
905}
906
907
908
909template <int dim, int spacedim>
910std::unique_ptr<FiniteElement<dim, spacedim>>
912{
913 return std::make_unique<FE_SimplexP<dim, spacedim>>(*this);
914}
915
916
917
918template <int dim, int spacedim>
919std::string
921{
922 std::ostringstream namebuf;
923 namebuf << "FE_SimplexP<" << Utilities::dim_string(dim, spacedim) << ">("
924 << this->degree << ")";
925
926 return namebuf.str();
927}
928
929
930
931template <int dim, int spacedim>
934 const FiniteElement<dim, spacedim> &fe_other,
935 const unsigned int codim) const
936{
937 Assert(codim <= dim, ExcImpossibleInDim(dim));
938
939 // vertex/line/face domination
940 // (if fe_other is derived from FE_SimplexDGP)
941 // ------------------------------------
942 if (codim > 0)
943 if (dynamic_cast<const FE_SimplexDGP<dim, spacedim> *>(&fe_other) !=
944 nullptr)
945 // there are no requirements between continuous and discontinuous
946 // elements
948
949 // vertex/line/face domination
950 // (if fe_other is not derived from FE_SimplexDGP)
951 // & cell domination
952 // ----------------------------------------
953 if (const FE_SimplexP<dim, spacedim> *fe_p_other =
954 dynamic_cast<const FE_SimplexP<dim, spacedim> *>(&fe_other))
955 {
956 if (this->degree < fe_p_other->degree)
958 else if (this->degree == fe_p_other->degree)
960 else
962 }
963 else if (const FE_Q<dim, spacedim> *fe_q_other =
964 dynamic_cast<const FE_Q<dim, spacedim> *>(&fe_other))
965 {
966 if (this->degree < fe_q_other->degree)
968 else if (this->degree == fe_q_other->degree)
970 else
972 }
973 else if (const FE_Nothing<dim, spacedim> *fe_nothing =
974 dynamic_cast<const FE_Nothing<dim, spacedim> *>(&fe_other))
975 {
976 if (fe_nothing->is_dominating())
978 else
979 // the FE_Nothing has no degrees of freedom and it is typically used
980 // in a context where we don't require any continuity along the
981 // interface
983 }
984
987}
988
989
990
991template <int dim, int spacedim>
992std::vector<std::pair<unsigned int, unsigned int>>
994 const FiniteElement<dim, spacedim> &fe_other) const
995{
996 AssertDimension(dim, 2);
997
998 if (dynamic_cast<const FE_SimplexP<dim, spacedim> *>(&fe_other) != nullptr)
999 {
1000 // there should be exactly one single DoF of each FE at a vertex, and
1001 // they should have identical value
1002 return {{0U, 0U}};
1003 }
1004 else if (dynamic_cast<const FE_Q<dim, spacedim> *>(&fe_other) != nullptr)
1005 {
1006 // there should be exactly one single DoF of each FE at a vertex, and
1007 // they should have identical value
1008 return {{0U, 0U}};
1009 }
1010 else if (dynamic_cast<const FE_Nothing<dim> *>(&fe_other) != nullptr)
1011 {
1012 // the FE_Nothing has no degrees of freedom, so there are no
1013 // equivalencies to be recorded
1014 return {};
1015 }
1016 else if (fe_other.n_unique_faces() == 1 && fe_other.n_dofs_per_face(0) == 0)
1017 {
1018 // if the other element has no DoFs on faces at all,
1019 // then it would be impossible to enforce any kind of
1020 // continuity even if we knew exactly what kind of element
1021 // we have -- simply because the other element declares
1022 // that it is discontinuous because it has no DoFs on
1023 // its faces. in that case, just state that we have no
1024 // constraints to declare
1025 return {};
1026 }
1027 else
1028 {
1030 return {};
1031 }
1032}
1033
1034
1035
1036template <int dim, int spacedim>
1037std::vector<std::pair<unsigned int, unsigned int>>
1039 const FiniteElement<dim, spacedim> &fe_other) const
1040{
1041 AssertDimension(dim, 2);
1042
1043 if (const FE_SimplexP<dim, spacedim> *fe_p_other =
1044 dynamic_cast<const FE_SimplexP<dim, spacedim> *>(&fe_other))
1045 {
1046 // dofs are located along lines, so two dofs are identical if they are
1047 // located at identical positions.
1048 // Therefore, read the points in unit_support_points for the
1049 // first coordinate direction. For FE_SimplexP, they are currently
1050 // hard-coded and we iterate over points on the first line which begin
1051 // after the 3 vertex points in the complete list of unit support points
1052
1053 std::vector<std::pair<unsigned int, unsigned int>> identities;
1054
1055 for (unsigned int i = 0; i < this->degree - 1; ++i)
1056 for (unsigned int j = 0; j < fe_p_other->degree - 1; ++j)
1057 if (std::fabs(this->unit_support_points[i + 3][0] -
1058 fe_p_other->unit_support_points[i + 3][0]) < 1e-14)
1059 identities.emplace_back(i, j);
1060 else
1061 {
1062 // If nodes are not located in the same place, we have to
1063 // interpolate. This is then not handled through the
1064 // current function, but via interpolation matrices that
1065 // result in constraints, rather than identities. Since
1066 // that happens in a different function, there is nothing
1067 // for us to do here.
1068 }
1069
1070 return identities;
1071 }
1072 else if (const FE_Q<dim, spacedim> *fe_q_other =
1073 dynamic_cast<const FE_Q<dim, spacedim> *>(&fe_other))
1074 {
1075 // dofs are located along lines, so two dofs are identical if they are
1076 // located at identical positions. if we had only equidistant points, we
1077 // could simply check for similarity like (i+1)*q == (j+1)*p, but we
1078 // might have other support points (e.g. Gauss-Lobatto
1079 // points). Therefore, read the points in unit_support_points for the
1080 // first coordinate direction. For FE_Q, we take the lexicographic
1081 // ordering of the line support points in the first direction (i.e.,
1082 // x-direction), which we access between index 1 and p-1 (index 0 and p
1083 // are vertex dofs). For FE_SimplexP, they are currently hard-coded and we
1084 // iterate over points on the first line which begin after the 3 vertex
1085 // points in the complete list of unit support points
1086
1087 const std::vector<unsigned int> &index_map_inverse_q_other =
1088 fe_q_other->get_poly_space_numbering_inverse();
1089
1090 std::vector<std::pair<unsigned int, unsigned int>> identities;
1091
1092 for (unsigned int i = 0; i < this->degree - 1; ++i)
1093 for (unsigned int j = 0; j < fe_q_other->degree - 1; ++j)
1094 if (std::fabs(this->unit_support_points[i + 3][0] -
1095 fe_q_other->get_unit_support_points()
1096 [index_map_inverse_q_other[j + 1]][0]) < 1e-14)
1097 identities.emplace_back(i, j);
1098 else
1099 {
1100 // If nodes are not located in the same place, we have to
1101 // interpolate. This will then also
1102 // capture the case where the FE_Q has a different polynomial
1103 // degree than the current element. In either case, the resulting
1104 // constraints are computed elsewhere, rather than via the
1105 // identities this function returns: Since
1106 // that happens in a different function, there is nothing
1107 // for us to do here.
1108 }
1109
1110 return identities;
1111 }
1112 else if (dynamic_cast<const FE_Nothing<dim> *>(&fe_other) != nullptr)
1113 {
1114 // The FE_Nothing has no degrees of freedom, so there are no
1115 // equivalencies to be recorded. (If the FE_Nothing is dominating,
1116 // then this will also leads to constraints, but we are not concerned
1117 // with this here.)
1118 return {};
1119 }
1120 else if (fe_other.n_unique_faces() == 1 && fe_other.n_dofs_per_face(0) == 0)
1121 {
1122 // if the other element has no elements on faces at all,
1123 // then it would be impossible to enforce any kind of
1124 // continuity even if we knew exactly what kind of element
1125 // we have -- simply because the other element declares
1126 // that it is discontinuous because it has no DoFs on
1127 // its faces. in that case, just state that we have no
1128 // constraints to declare
1129 return {};
1130 }
1131 else
1132 {
1134 return {};
1135 }
1136}
1137
1138
1139
1140template <int dim, int spacedim>
1142 : FE_SimplexPoly<dim, spacedim>(
1143 BarycentricPolynomials<dim>::get_fe_p_basis(degree),
1144 FiniteElementData<dim>(get_dpo_vector_fe_dgp(dim, degree),
1145 ReferenceCells::get_simplex<dim>(),
1146 1,
1147 degree,
1148 FiniteElementData<dim>::L2),
1149 true,
1150 unit_support_points_fe_p<dim>(degree),
1151 unit_face_support_points_fe_p<dim>(degree, FiniteElementData<dim>::L2),
1152 constraints_fe_p<dim>(degree))
1153{}
1154
1155
1156
1157template <int dim, int spacedim>
1158std::unique_ptr<FiniteElement<dim, spacedim>>
1160{
1161 return std::make_unique<FE_SimplexDGP<dim, spacedim>>(*this);
1162}
1163
1164
1165
1166template <int dim, int spacedim>
1167std::string
1169{
1170 std::ostringstream namebuf;
1171 namebuf << "FE_SimplexDGP<" << Utilities::dim_string(dim, spacedim) << ">("
1172 << this->degree << ")";
1173
1174 return namebuf.str();
1175}
1176
1177
1178template <int dim, int spacedim>
1181 const FiniteElement<dim, spacedim> &fe_other,
1182 const unsigned int codim) const
1183{
1184 Assert(codim <= dim, ExcImpossibleInDim(dim));
1185
1186 // vertex/line/face domination
1187 // ---------------------------
1188 if (codim > 0)
1189 // this is a discontinuous element, so by definition there will
1190 // be no constraints wherever this element comes together with
1191 // any other kind of element
1193
1194 // cell domination
1195 // ---------------
1196 if (const FE_SimplexDGP<dim, spacedim> *fe_dgp_other =
1197 dynamic_cast<const FE_SimplexDGP<dim, spacedim> *>(&fe_other))
1198 {
1199 if (this->degree < fe_dgp_other->degree)
1201 else if (this->degree == fe_dgp_other->degree)
1203 else
1205 }
1206 else if (const FE_DGQ<dim, spacedim> *fe_dgq_other =
1207 dynamic_cast<const FE_DGQ<dim, spacedim> *>(&fe_other))
1208 {
1209 if (this->degree < fe_dgq_other->degree)
1211 else if (this->degree == fe_dgq_other->degree)
1213 else
1215 }
1216 else if (const FE_Nothing<dim, spacedim> *fe_nothing =
1217 dynamic_cast<const FE_Nothing<dim, spacedim> *>(&fe_other))
1218 {
1219 if (fe_nothing->is_dominating())
1221 else
1222 // the FE_Nothing has no degrees of freedom and it is typically used
1223 // in a context where we don't require any continuity along the
1224 // interface
1226 }
1227
1230}
1231
1232
1233
1234template <int dim, int spacedim>
1235std::vector<std::pair<unsigned int, unsigned int>>
1237 const FiniteElement<dim, spacedim> &fe_other) const
1238{
1239 (void)fe_other;
1240
1241 return {};
1242}
1243
1244
1245
1246template <int dim, int spacedim>
1247std::vector<std::pair<unsigned int, unsigned int>>
1249 const FiniteElement<dim, spacedim> &fe_other) const
1250{
1251 (void)fe_other;
1252
1253 return {};
1254}
1255
1256
1257
1258template <int dim, int spacedim>
1259const FullMatrix<double> &
1261 const unsigned int child,
1262 const RefinementCase<dim> &refinement_case) const
1263{
1264 if (dim == 3)
1265 Assert(RefinementCase<dim>(refinement_case) ==
1267 static_cast<char>(IsotropicRefinementChoice::cut_tet_68)) ||
1268 RefinementCase<dim>(refinement_case) ==
1270 static_cast<char>(IsotropicRefinementChoice::cut_tet_57)) ||
1271 RefinementCase<dim>(refinement_case) ==
1273 static_cast<char>(IsotropicRefinementChoice::cut_tet_49)),
1275 else
1278 AssertDimension(dim, spacedim);
1279
1280 // initialization upon first request
1281 if (this->restriction[refinement_case - 1][child].n() == 0)
1282 {
1283 std::lock_guard<std::mutex> lock(this->restriction_matrix_mutex);
1284
1285 // if matrix got updated while waiting for the lock
1286 if (this->restriction[refinement_case - 1][child].n() ==
1287 this->n_dofs_per_cell())
1288 return this->restriction[refinement_case - 1][child];
1289
1290 // now do the work. need to get a non-const version of data in order to
1291 // be able to modify them inside a const function
1292 auto &this_nonconst = const_cast<FE_SimplexDGP<dim, spacedim> &>(*this);
1293
1294 if (dim == 2)
1295 {
1296 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
1298 isotropic_matrices.back().resize(
1299 this->reference_cell().n_children(
1300 RefinementCase<dim>(refinement_case)),
1302 this->n_dofs_per_cell()));
1303
1304 FETools::compute_projection_matrices(*this, isotropic_matrices, true);
1305
1306 this_nonconst.restriction[refinement_case - 1] =
1307 std::move(isotropic_matrices.back());
1308 }
1309 else if (dim == 3)
1310 {
1311 std::vector<std::vector<FullMatrix<double>>> matrices(
1312 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49),
1313 std::vector<FullMatrix<double>>(
1314 this->reference_cell().n_children(
1315 RefinementCase<dim>(refinement_case)),
1317 this->n_dofs_per_cell())));
1318 FETools::compute_projection_matrices(*this, matrices, true);
1319 for (unsigned int refinement_direction = static_cast<unsigned int>(
1321 refinement_direction <=
1322 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49);
1323 refinement_direction++)
1324 this_nonconst.restriction[refinement_direction - 1] =
1325 std::move(matrices[refinement_direction - 1]);
1326 }
1327 else
1329 }
1330
1331 // finally return the matrix
1332 return this->restriction[refinement_case - 1][child];
1333}
1334
1335// explicit instantiations
1336#include "fe/fe_simplex_p.inst"
1337
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition array_view.h:954
static BarycentricPolynomials< dim > get_fe_p_basis(const unsigned int degree)
TriaIterator< CellAccessor< dim, spacedim > > child(const unsigned int i) const
FE_Poly(const ScalarPolynomialsBase< dim > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
Definition fe_q.h:554
FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim) const override
std::string get_name() const override
std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
FE_SimplexDGP(const unsigned int degree)
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
FE_SimplexP(const unsigned int degree)
std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim) const override
std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
std::string get_name() const override
std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
FE_SimplexPoly(const BarycentricPolynomials< dim > polynomials, const FiniteElementData< dim > &fe_data, const bool prolongation_is_additive, const std::vector< Point< dim > > &unit_support_points, const std::vector< std::vector< Point< dim - 1 > > > unit_face_support_points, const FullMatrix< double > &interface_constraints)
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &x_source_fe, const unsigned int subface, FullMatrix< double > &interpolation_matrix, const unsigned int face_no) const override
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const override
bool hp_constraints_are_implemented() const override
Threads::Mutex restriction_matrix_mutex
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const types::geometric_orientation combined_orientation=numbers::default_geometric_orientation) const override
Threads::Mutex prolongation_matrix_mutex
void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source_fe, FullMatrix< double > &interpolation_matrix, const unsigned int face_no) const override
unsigned int get_first_line_index() const
unsigned int n_dofs_per_vertex() const
const unsigned int degree
Definition fe_data.h:452
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_line() const
unsigned int get_first_quad_index(const unsigned int quad_no=0) const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
ReferenceCell reference_cell() const
unsigned int get_first_face_line_index(const unsigned int face_no=0) const
unsigned int get_first_face_quad_index(const unsigned int face_no=0) const
unsigned int n_unique_faces() const
unsigned int n_dofs_per_quad(unsigned int face_no=0) const
const unsigned int dofs_per_cell
Definition fe_data.h:436
FiniteElementData(const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const BlockIndices &block_indices=BlockIndices())
Definition fe_data.cc:114
virtual Point< dim > unit_support_point(const unsigned int index) const
std::vector< std::vector< Point< dim - 1 > > > unit_face_support_points
Definition fe.h:2569
std::vector< std::vector< FullMatrix< double > > > restriction
Definition fe.h:2524
const std::vector< Point< dim > > & get_unit_support_points() const
const std::vector< Point< dim - 1 > > & get_unit_face_support_points(const unsigned int face_no=0) const
std::vector< int > adjust_line_dof_index_for_line_orientation_table
Definition fe.h:2611
FiniteElement(const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
std::vector< Point< dim > > unit_support_points
Definition fe.h:2562
FullMatrix< double > interface_constraints
Definition fe.h:2550
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition fe.h:2538
size_type n() const
size_type m() const
virtual void transform_points_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim > > &real_points, const ArrayView< Point< dim > > &unit_points) const
Definition point.h:113
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 void project_to_face(const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
cell_iterator begin(const unsigned int level=0) const
virtual void execute_coarsening_and_refinement()
active_cell_iterator begin_active(const unsigned int level=0) const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
constexpr bool running_in_debug_mode()
Definition config.h:78
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcInterpolationNotImplemented()
#define AssertThrow(cond, exc)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:316
void compute_embedding_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false, const double threshold=1.e-12)
void compute_projection_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
constexpr char U
constexpr ReferenceCell Triangle
constexpr ReferenceCell Tetrahedron
constexpr const ReferenceCell & get_simplex()
constexpr ReferenceCell Line
std::string dim_string(const int dim, const int spacedim)
Definition utilities.cc:551
constexpr types::geometric_orientation reverse_line_orientation
Definition types.h:365
constexpr types::geometric_orientation default_geometric_orientation
Definition types.h:352
STL namespace.
unsigned char geometric_orientation
Definition types.h:40
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)