Reference documentation for deal.II version 9.6.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 - 2024 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
19
20#include <deal.II/fe/fe_dgq.h>
22#include <deal.II/fe/fe_q.h>
24#include <deal.II/fe/fe_tools.h>
25#include <deal.II/fe/mapping.h>
26
28
30
31namespace
32{
33 // TODO: replace this with QProjector once QProjector learns how to project
34 // quadrature points onto separate faces of simplices
35 template <int dim>
37 face_midpoint(const unsigned int face_no)
38 {
40 const auto face_reference_cell =
41 reference_cell.face_reference_cell(face_no);
43 Point<dim> midpoint;
44 for (const auto face_vertex_no : face_reference_cell.vertex_indices())
45 {
46 const auto vertex_no = reference_cell.face_to_cell_vertices(
47 face_no,
48 face_vertex_no,
50
51 midpoint += reference_cell.template vertex<dim>(vertex_no);
52 }
53
54 midpoint /= reference_cell.face_reference_cell(0).n_vertices();
55 return midpoint;
56 }
57
62 std::vector<unsigned int>
63 get_dpo_vector_fe_p(const unsigned int dim, const unsigned int degree)
64 {
65 switch (dim)
66 {
67 case 1:
68 switch (degree)
69 {
70 case 1:
71 return {1, 0};
72 case 2:
73 return {1, 1};
74 case 3:
75 return {1, 2};
76 default:
78 }
79 case 2:
80 switch (degree)
81 {
82 case 1:
83 return {1, 0, 0};
84 case 2:
85 return {1, 1, 0};
86 case 3:
87 return {1, 2, 1};
88 default:
90 }
91 case 3:
92 switch (degree)
93 {
94 case 1:
95 return {1, 0, 0, 0};
96 case 2:
97 return {1, 1, 0, 0};
98 case 3:
99 return {1, 2, 1, 0};
100 default:
102 }
103 }
104
106 return {};
108
109
110
115 template <int dim>
116 std::vector<Point<dim>>
117 unit_support_points_fe_p(const unsigned int degree)
118 {
119 Assert(dim != 0, ExcInternalError());
120 Assert(degree <= 3, ExcNotImplemented());
121 std::vector<Point<dim>> unit_points;
122 const auto reference_cell = ReferenceCells::get_simplex<dim>();
123
124 // Piecewise constants are a special case: use a support point at the
125 // centroid and only the centroid
126 if (degree == 0)
127 {
128 unit_points.emplace_back(reference_cell.template barycenter<dim>());
129 return unit_points;
130 }
131
132 // otherwise write everything as linear combinations of vertices
133 const auto dpo = get_dpo_vector_fe_p(dim, degree);
134 Assert(dpo.size() == dim + 1, ExcInternalError());
135 Assert(dpo[0] == 1, ExcNotImplemented());
136 // vertices:
137 for (const unsigned int d : reference_cell.vertex_indices())
138 unit_points.push_back(reference_cell.template vertex<dim>(d));
139 // lines:
140 for (const unsigned int l : reference_cell.line_indices())
141 {
142 const Point<dim> p0 =
143 unit_points[reference_cell.line_to_cell_vertices(l, 0)];
144 const Point<dim> p1 =
145 unit_points[reference_cell.line_to_cell_vertices(l, 1)];
146 for (unsigned int p = 0; p < dpo[1]; ++p)
147 unit_points.push_back((double(dpo[1] - p) / (dpo[1] + 1)) * p0 +
148 (double(p + 1) / (dpo[1] + 1)) * p1);
149 }
150 // quads:
151 if (dim > 1 && dpo[2] > 0)
152 {
153 Assert(dpo[2] == 1, ExcNotImplemented());
154 if (dim == 2)
155 unit_points.push_back(reference_cell.template barycenter<dim>());
156 if (dim == 3)
157 for (const unsigned int f : reference_cell.face_indices())
158 unit_points.push_back(face_midpoint<dim>(f));
159 }
160
161 return unit_points;
162 }
163
164 template <>
165 std::vector<Point<0>>
166 unit_support_points_fe_p(const unsigned int /*degree*/)
167 {
168 return {Point<0>()};
169 }
170
175 template <int dim>
176 std::vector<std::vector<Point<dim - 1>>>
177 unit_face_support_points_fe_p(
178 const unsigned int degree,
179 typename FiniteElementData<dim>::Conformity conformity)
180 {
181 // Discontinuous elements don't have face support points
183 return {};
184
185 // this concept doesn't exist in 1d so just return an empty vector
186 if (dim == 1)
187 return {};
188
189 std::vector<std::vector<Point<dim - 1>>> unit_face_points;
190
191 // all faces have the same support points
192 for (auto face_n : ReferenceCells::get_simplex<dim>().face_indices())
193 {
194 (void)face_n;
195 unit_face_points.emplace_back(
196 unit_support_points_fe_p<dim - 1>(degree));
197 }
198
199 return unit_face_points;
200 }
201
207 template <int dim>
209 constraints_fe_p(const unsigned int /*degree*/)
210 {
211 // no constraints in 1d
212 // constraints in 3d not implemented yet
213 return FullMatrix<double>();
214 }
215
216 template <>
218 constraints_fe_p<2>(const unsigned int degree)
219 {
220 constexpr int dim = 2;
221 Assert(degree <= 3, ExcNotImplemented());
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 if (degree == 2)
233 {
234 // midpoint on subface 0
235 constraint_points.emplace_back(0.25);
236 // midpoint on subface 1
237 constraint_points.emplace_back(0.75);
238 }
239
240 // Now construct relation between destination (child) and source (mother)
241 // dofs.
242
243 const unsigned int n_dofs_constrained = constraint_points.size();
244 unsigned int n_dofs_per_face = degree + 1;
245 FullMatrix<double> interface_constraints(n_dofs_constrained,
246 n_dofs_per_face);
247
248 const auto poly = BarycentricPolynomials<dim - 1>::get_fe_p_basis(degree);
249
250 for (unsigned int i = 0; i < n_dofs_constrained; ++i)
251 for (unsigned int j = 0; j < n_dofs_per_face; ++j)
252 {
253 interface_constraints(i, j) =
254 poly.compute_value(j, constraint_points[i]);
255
256 // if the value is small up to round-off, then simply set it to zero
257 // to avoid unwanted fill-in of the constraint matrices (which would
258 // then increase the number of other DoFs a constrained DoF would
259 // couple to)
260 if (std::fabs(interface_constraints(i, j)) < 1e-13)
261 interface_constraints(i, j) = 0;
262 }
263 return interface_constraints;
264 }
265
266
267
272 std::vector<unsigned int>
273 get_dpo_vector_fe_dgp(const unsigned int dim, const unsigned int degree)
274 {
275 // First treat the case of piecewise constant elements:
276 if (degree == 0)
277 {
278 std::vector<unsigned int> dpo(dim + 1, 0U);
279 dpo[dim] = 1;
280 return dpo;
281 }
282 else
283 {
284 // This element has the same degrees of freedom as the continuous one,
285 // but they are all counted for the interior of the cell because
286 // it is continuous. Rather than hard-code how many DoFs the element
287 // has, we just get the numbers from the continuous case and add them
288 // up
289 const auto continuous_dpo = get_dpo_vector_fe_p(dim, degree);
290
291 switch (dim)
292 {
293 case 1:
294 return {0U,
295 ReferenceCells::Line.n_vertices() * continuous_dpo[0] +
296 continuous_dpo[dim]};
297
298 case 2:
299 return {0U,
300 0U,
302 continuous_dpo[0] +
303 ReferenceCells::Triangle.n_lines() * continuous_dpo[1] +
304 continuous_dpo[dim]};
305
306 case 3:
307 return {
308 0U,
309 0U,
310 0U,
311 ReferenceCells::Tetrahedron.n_vertices() * continuous_dpo[0] +
312 ReferenceCells::Tetrahedron.n_lines() * continuous_dpo[1] +
313 ReferenceCells::Tetrahedron.n_faces() * continuous_dpo[2] +
314 continuous_dpo[dim]};
315 }
316
318 return {};
319 }
320 }
321} // namespace
322
323
324
325template <int dim, int spacedim>
327 const BarycentricPolynomials<dim> polynomials,
328 const FiniteElementData<dim> &fe_data,
329 const bool prolongation_is_additive,
330 const std::vector<Point<dim>> &unit_support_points,
331 const std::vector<std::vector<Point<dim - 1>>> unit_face_support_points,
332 const FullMatrix<double> &interface_constraints)
333 : ::FE_Poly<dim, spacedim>(
334 polynomials,
335 fe_data,
336 std::vector<bool>(fe_data.dofs_per_cell, prolongation_is_additive),
337 std::vector<ComponentMask>(fe_data.dofs_per_cell,
338 ComponentMask(std::vector<bool>(1, true))))
339{
342 this->interface_constraints = interface_constraints;
343}
344
345
346
347template <int dim, int spacedim>
348std::pair<Table<2, bool>, std::vector<unsigned int>>
350{
351 Table<2, bool> constant_modes(1, this->n_dofs_per_cell());
352 constant_modes.fill(true);
353 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
354 constant_modes, std::vector<unsigned int>(1, 0));
355}
356
357
358
359template <int dim, int spacedim>
360const FullMatrix<double> &
362 const unsigned int child,
363 const RefinementCase<dim> &refinement_case) const
364{
365 if (dim == 3)
366 Assert(RefinementCase<dim>(refinement_case) ==
368 static_cast<char>(IsotropicRefinementChoice::cut_tet_68)) ||
369 RefinementCase<dim>(refinement_case) ==
371 static_cast<char>(IsotropicRefinementChoice::cut_tet_57)) ||
372 RefinementCase<dim>(refinement_case) ==
374 static_cast<char>(IsotropicRefinementChoice::cut_tet_49)),
376 else
377 Assert(refinement_case ==
380 AssertDimension(dim, spacedim);
381
382 // initialization upon first request
383 if (this->prolongation[refinement_case - 1][child].n() == 0)
384 {
385 std::lock_guard<std::mutex> lock(prolongation_matrix_mutex);
386
387 // if matrix got updated while waiting for the lock
388 if (this->prolongation[refinement_case - 1][child].n() ==
389 this->n_dofs_per_cell())
390 return this->prolongation[refinement_case - 1][child];
391
392 // now do the work. need to get a non-const version of data in order to
393 // be able to modify them inside a const function
394 auto &this_nonconst = const_cast<FE_SimplexPoly<dim, spacedim> &>(*this);
395
396 if (dim == 2)
397 {
398 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
400 isotropic_matrices.back().resize(
401 this->reference_cell().n_children(
402 RefinementCase<dim>(refinement_case)),
403 FullMatrix<double>(this->n_dofs_per_cell(),
404 this->n_dofs_per_cell()));
405
406 FETools::compute_embedding_matrices(*this, isotropic_matrices, true);
407
408 this_nonconst.prolongation[refinement_case - 1] =
409 std::move(isotropic_matrices.back());
410 }
411 else if (dim == 3)
412 {
413 std::vector<std::vector<FullMatrix<double>>> matrices(
414 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49),
415 std::vector<FullMatrix<double>>(
416 this->reference_cell().n_children(
417 RefinementCase<dim>(refinement_case)),
418 FullMatrix<double>(this->n_dofs_per_cell(),
419 this->n_dofs_per_cell())));
420 FETools::compute_embedding_matrices(*this, matrices, true);
421 for (unsigned int refinement_direction = static_cast<unsigned int>(
423 refinement_direction <=
424 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49);
425 refinement_direction++)
426 this_nonconst.prolongation[refinement_direction - 1] =
427 std::move(matrices[refinement_direction - 1]);
428 }
429 else
431 }
432
433 // finally return the matrix
434 return this->prolongation[refinement_case - 1][child];
435}
436
437
438
439template <int dim, int spacedim>
440const FullMatrix<double> &
442 const unsigned int child,
443 const RefinementCase<dim> &refinement_case) const
444{
445 if (dim == 3)
446 Assert(RefinementCase<dim>(refinement_case) ==
448 static_cast<char>(IsotropicRefinementChoice::cut_tet_68)) ||
449 RefinementCase<dim>(refinement_case) ==
451 static_cast<char>(IsotropicRefinementChoice::cut_tet_57)) ||
452 RefinementCase<dim>(refinement_case) ==
454 static_cast<char>(IsotropicRefinementChoice::cut_tet_49)),
456 else
459 AssertDimension(dim, spacedim);
460
461 // initialization upon first request
462 if (this->restriction[refinement_case - 1][child].n() == 0)
463 {
464 std::lock_guard<std::mutex> lock(restriction_matrix_mutex);
465
466 // if matrix got updated while waiting for the lock
467 if (this->restriction[refinement_case - 1][child].n() ==
468 this->n_dofs_per_cell())
469 return this->restriction[refinement_case - 1][child];
470
471 // get the restriction matrix
472 // Refine a unit cell. As the parent cell is a unit
473 // cell, the reference cell of the children equals the parent, i.e. they
474 // have the support points at the same locations. So we just have to check
475 // if a support point of the parent is one of the interpolation points of
476 // the child. If this is not the case we find the interpolation of the
477 // point.
478
479 const double eps = 1e-12;
480 FullMatrix<double> restriction_mat(this->n_dofs_per_cell(),
481 this->n_dofs_per_cell());
482
483 // first get all support points on the reference cell
484 const std::vector<Point<dim>> unit_support_points =
485 this->get_unit_support_points();
486
487 // now create children on the reference cell
489 GridGenerator::reference_cell(tria, this->reference_cell());
490 tria.begin_active()->set_refine_flag(
492 if (dim == 3)
493 tria.begin_active()->set_refine_choice(refinement_case);
495
496 const auto &child_cell = tria.begin(0)->child(child);
497
498 // iterate over all support points and transform them to the unit cell of
499 // the child
500 for (unsigned int i = 0; i < unit_support_points.size(); i++)
501 {
502 std::vector<Point<dim>> transformed_point(1);
503 const std::vector<Point<spacedim>> unit_support_point = {
504 dim == 2 ? Point<spacedim>(unit_support_points[i][0],
505 unit_support_points[i][1]) :
506 Point<spacedim>(unit_support_points[i][0],
507 unit_support_points[i][1],
508 unit_support_points[i][2])};
509 this->reference_cell()
512 child_cell,
513 make_array_view(unit_support_point),
514 make_array_view(transformed_point));
515
516 // if point is inside the unit cell iterate over all shape functions
517 if (this->reference_cell().contains_point(transformed_point[0], eps))
518 for (unsigned int j = 0; j < this->n_dofs_per_cell(); j++)
519 restriction_mat[i][j] =
520 this->shape_value(j, transformed_point[0]);
521 }
522#ifdef DEBUG
523 for (unsigned int i = 0; i < this->n_dofs_per_cell(); i++)
524 {
525 double sum = 0.;
526
527 for (unsigned int j = 0; j < this->n_dofs_per_cell(); j++)
528 sum += restriction_mat[i][j];
529
530 Assert(std::fabs(sum - 1) < eps || std::fabs(sum) < eps,
532 "The entries in a row of the local "
533 "restriction matrix do not add to zero or one. "
534 "This typically indicates that the "
535 "polynomial interpolation is "
536 "ill-conditioned such that round-off "
537 "prevents the sum to be one."));
538 }
539#endif
540
541 // Remove small entries from the matrix
542 for (unsigned int i = 0; i < restriction_mat.m(); ++i)
543 for (unsigned int j = 0; j < restriction_mat.n(); ++j)
544 {
545 if (std::fabs(restriction_mat(i, j)) < eps)
546 restriction_mat(i, j) = 0.;
547 if (std::fabs(restriction_mat(i, j) - 1) < eps)
548 restriction_mat(i, j) = 1.;
549 }
550
551 const_cast<FullMatrix<double> &>(
552 this->restriction[refinement_case - 1][child]) =
553 std::move(restriction_mat);
554 }
555
556 // finally return the matrix
557 return this->restriction[refinement_case - 1][child];
558}
559
560
561
562template <int dim, int spacedim>
563void
565 const FiniteElement<dim, spacedim> &source_fe,
566 FullMatrix<double> &interpolation_matrix,
567 const unsigned int face_no) const
568{
569 Assert(interpolation_matrix.m() == source_fe.n_dofs_per_face(face_no),
570 ExcDimensionMismatch(interpolation_matrix.m(),
571 source_fe.n_dofs_per_face(face_no)));
572
573 // see if source is a P or Q element
574 if ((dynamic_cast<const FE_SimplexPoly<dim, spacedim> *>(&source_fe) !=
575 nullptr) ||
576 (dynamic_cast<const FE_Q_Base<dim, spacedim> *>(&source_fe) != nullptr))
577 {
578 const Quadrature<dim - 1> quad_face_support(
579 source_fe.get_unit_face_support_points(face_no));
580
581 const double eps = 2e-13 * this->degree * (dim - 1);
582
583 std::vector<Point<dim>> face_quadrature_points(quad_face_support.size());
584 QProjector<dim>::project_to_face(this->reference_cell(),
585 quad_face_support,
586 face_no,
587 face_quadrature_points);
588
589 for (unsigned int i = 0; i < source_fe.n_dofs_per_face(face_no); ++i)
590 for (unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
591 {
592 double matrix_entry =
593 this->shape_value(this->face_to_cell_index(j, 0),
594 face_quadrature_points[i]);
595
596 // Correct the interpolated value. I.e. if it is close to 1 or
597 // 0, make it exactly 1 or 0. Unfortunately, this is required to
598 // avoid problems with higher order elements.
599 if (std::fabs(matrix_entry - 1.0) < eps)
600 matrix_entry = 1.0;
601 if (std::fabs(matrix_entry) < eps)
602 matrix_entry = 0.0;
603
604 interpolation_matrix(i, j) = matrix_entry;
605 }
606
607#ifdef DEBUG
608 for (unsigned int j = 0; j < source_fe.n_dofs_per_face(face_no); ++j)
609 {
610 double sum = 0.;
611
612 for (unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
613 sum += interpolation_matrix(j, i);
614
615 Assert(std::fabs(sum - 1) < eps, ExcInternalError());
616 }
617#endif
618 }
619 else if (dynamic_cast<const FE_Nothing<dim> *>(&source_fe) != nullptr)
620 {
621 // nothing to do here, the FE_Nothing has no degrees of freedom anyway
622 }
623 else
625 false,
626 (typename FiniteElement<dim,
627 spacedim>::ExcInterpolationNotImplemented()));
628}
629
630
631
632template <int dim, int spacedim>
633void
635 const FiniteElement<dim, spacedim> &source_fe,
636 const unsigned int subface,
637 FullMatrix<double> &interpolation_matrix,
638 const unsigned int face_no) const
639{
640 Assert(interpolation_matrix.m() == source_fe.n_dofs_per_face(face_no),
641 ExcDimensionMismatch(interpolation_matrix.m(),
642 source_fe.n_dofs_per_face(face_no)));
643
644 // see if source is a P or Q element
645 if ((dynamic_cast<const FE_SimplexPoly<dim, spacedim> *>(&source_fe) !=
646 nullptr) ||
647 (dynamic_cast<const FE_Q_Base<dim, spacedim> *>(&source_fe) != nullptr))
648 {
649 const Quadrature<dim - 1> quad_face_support(
650 source_fe.get_unit_face_support_points(face_no));
651
652 const double eps = 2e-13 * this->degree * (dim - 1);
653
654 std::vector<Point<dim>> subface_quadrature_points(
655 quad_face_support.size());
656 QProjector<dim>::project_to_subface(this->reference_cell(),
657 quad_face_support,
658 face_no,
659 subface,
660 subface_quadrature_points);
661
662 for (unsigned int i = 0; i < source_fe.n_dofs_per_face(face_no); ++i)
663 for (unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
664 {
665 double matrix_entry =
666 this->shape_value(this->face_to_cell_index(j, 0),
667 subface_quadrature_points[i]);
668
669 // Correct the interpolated value. I.e. if it is close to 1 or
670 // 0, make it exactly 1 or 0. Unfortunately, this is required to
671 // avoid problems with higher order elements.
672 if (std::fabs(matrix_entry - 1.0) < eps)
673 matrix_entry = 1.0;
674 if (std::fabs(matrix_entry) < eps)
675 matrix_entry = 0.0;
676
677 interpolation_matrix(i, j) = matrix_entry;
678 }
679
680#ifdef DEBUG
681 for (unsigned int j = 0; j < source_fe.n_dofs_per_face(face_no); ++j)
682 {
683 double sum = 0.;
684
685 for (unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
686 sum += interpolation_matrix(j, i);
687
688 Assert(std::fabs(sum - 1) < eps, ExcInternalError());
689 }
690#endif
691 }
692 else if (dynamic_cast<const FE_Nothing<dim> *>(&source_fe) != nullptr)
693 {
694 // nothing to do here, the FE_Nothing has no degrees of freedom anyway
695 }
696 else
698 false,
699 (typename FiniteElement<dim,
700 spacedim>::ExcInterpolationNotImplemented()));
701}
702
703
704
705template <int dim, int spacedim>
706bool
711
712
713
714template <int dim, int spacedim>
715void
718 const std::vector<Vector<double>> &support_point_values,
719 std::vector<double> &nodal_values) const
720{
721 AssertDimension(support_point_values.size(),
722 this->get_unit_support_points().size());
723 AssertDimension(support_point_values.size(), nodal_values.size());
724 AssertDimension(this->dofs_per_cell, nodal_values.size());
725
726 for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
727 {
728 AssertDimension(support_point_values[i].size(), 1);
729
730 nodal_values[i] = support_point_values[i](0);
731 }
732}
733
734
735
736template <int dim, int spacedim>
738 : FE_SimplexPoly<dim, spacedim>(
739 BarycentricPolynomials<dim>::get_fe_p_basis(degree),
740 FiniteElementData<dim>(get_dpo_vector_fe_p(dim, degree),
741 ReferenceCells::get_simplex<dim>(),
742 1,
743 degree,
744 FiniteElementData<dim>::H1),
745 false,
746 unit_support_points_fe_p<dim>(degree),
747 unit_face_support_points_fe_p<dim>(degree, FiniteElementData<dim>::H1),
748 constraints_fe_p<dim>(degree))
749{
750 if (degree > 2)
751 for (unsigned int i = 0; i < this->n_dofs_per_line(); ++i)
753 this->n_dofs_per_line() - 1 - i - i;
754 // We do not support multiple DoFs per quad yet
756}
757
758
759
760template <int dim, int spacedim>
761std::unique_ptr<FiniteElement<dim, spacedim>>
763{
764 return std::make_unique<FE_SimplexP<dim, spacedim>>(*this);
765}
766
767
768
769template <int dim, int spacedim>
770std::string
772{
773 std::ostringstream namebuf;
774 namebuf << "FE_SimplexP<" << Utilities::dim_string(dim, spacedim) << ">("
775 << this->degree << ")";
776
777 return namebuf.str();
778}
779
780
781
782template <int dim, int spacedim>
785 const FiniteElement<dim, spacedim> &fe_other,
786 const unsigned int codim) const
787{
788 Assert(codim <= dim, ExcImpossibleInDim(dim));
789
790 // vertex/line/face domination
791 // (if fe_other is derived from FE_SimplexDGP)
792 // ------------------------------------
793 if (codim > 0)
794 if (dynamic_cast<const FE_SimplexDGP<dim, spacedim> *>(&fe_other) !=
795 nullptr)
796 // there are no requirements between continuous and discontinuous
797 // elements
799
800 // vertex/line/face domination
801 // (if fe_other is not derived from FE_SimplexDGP)
802 // & cell domination
803 // ----------------------------------------
804 if (const FE_SimplexP<dim, spacedim> *fe_p_other =
805 dynamic_cast<const FE_SimplexP<dim, spacedim> *>(&fe_other))
806 {
807 if (this->degree < fe_p_other->degree)
809 else if (this->degree == fe_p_other->degree)
811 else
813 }
814 else if (const FE_Q<dim, spacedim> *fe_q_other =
815 dynamic_cast<const FE_Q<dim, spacedim> *>(&fe_other))
816 {
817 if (this->degree < fe_q_other->degree)
819 else if (this->degree == fe_q_other->degree)
821 else
823 }
824 else if (const FE_Nothing<dim, spacedim> *fe_nothing =
825 dynamic_cast<const FE_Nothing<dim, spacedim> *>(&fe_other))
826 {
827 if (fe_nothing->is_dominating())
829 else
830 // the FE_Nothing has no degrees of freedom and it is typically used
831 // in a context where we don't require any continuity along the
832 // interface
834 }
835
838}
839
840
841
842template <int dim, int spacedim>
843std::vector<std::pair<unsigned int, unsigned int>>
845 const FiniteElement<dim, spacedim> &fe_other) const
846{
847 AssertDimension(dim, 2);
848
849 if (dynamic_cast<const FE_SimplexP<dim, spacedim> *>(&fe_other) != nullptr)
850 {
851 // there should be exactly one single DoF of each FE at a vertex, and
852 // they should have identical value
853 return {{0U, 0U}};
854 }
855 else if (dynamic_cast<const FE_Q<dim, spacedim> *>(&fe_other) != nullptr)
856 {
857 // there should be exactly one single DoF of each FE at a vertex, and
858 // they should have identical value
859 return {{0U, 0U}};
860 }
861 else if (dynamic_cast<const FE_Nothing<dim> *>(&fe_other) != nullptr)
862 {
863 // the FE_Nothing has no degrees of freedom, so there are no
864 // equivalencies to be recorded
865 return {};
866 }
867 else if (fe_other.n_unique_faces() == 1 && fe_other.n_dofs_per_face(0) == 0)
868 {
869 // if the other element has no DoFs on faces at all,
870 // then it would be impossible to enforce any kind of
871 // continuity even if we knew exactly what kind of element
872 // we have -- simply because the other element declares
873 // that it is discontinuous because it has no DoFs on
874 // its faces. in that case, just state that we have no
875 // constraints to declare
876 return {};
877 }
878 else
879 {
881 return {};
882 }
883}
884
885
886
887template <int dim, int spacedim>
888std::vector<std::pair<unsigned int, unsigned int>>
890 const FiniteElement<dim, spacedim> &fe_other) const
891{
892 AssertDimension(dim, 2);
893
894 if (const FE_SimplexP<dim, spacedim> *fe_p_other =
895 dynamic_cast<const FE_SimplexP<dim, spacedim> *>(&fe_other))
896 {
897 // dofs are located along lines, so two dofs are identical if they are
898 // located at identical positions.
899 // Therefore, read the points in unit_support_points for the
900 // first coordinate direction. For FE_SimplexP, they are currently
901 // hard-coded and we iterate over points on the first line which begin
902 // after the 3 vertex points in the complete list of unit support points
903
904 std::vector<std::pair<unsigned int, unsigned int>> identities;
905
906 for (unsigned int i = 0; i < this->degree - 1; ++i)
907 for (unsigned int j = 0; j < fe_p_other->degree - 1; ++j)
908 if (std::fabs(this->unit_support_points[i + 3][0] -
909 fe_p_other->unit_support_points[i + 3][0]) < 1e-14)
910 identities.emplace_back(i, j);
911 else
912 {
913 // If nodes are not located in the same place, we have to
914 // interpolate. This is then not handled through the
915 // current function, but via interpolation matrices that
916 // result in constraints, rather than identities. Since
917 // that happens in a different function, there is nothing
918 // for us to do here.
919 }
920
921 return identities;
922 }
923 else if (const FE_Q<dim, spacedim> *fe_q_other =
924 dynamic_cast<const FE_Q<dim, spacedim> *>(&fe_other))
925 {
926 // dofs are located along lines, so two dofs are identical if they are
927 // located at identical positions. if we had only equidistant points, we
928 // could simply check for similarity like (i+1)*q == (j+1)*p, but we
929 // might have other support points (e.g. Gauss-Lobatto
930 // points). Therefore, read the points in unit_support_points for the
931 // first coordinate direction. For FE_Q, we take the lexicographic
932 // ordering of the line support points in the first direction (i.e.,
933 // x-direction), which we access between index 1 and p-1 (index 0 and p
934 // are vertex dofs). For FE_SimplexP, they are currently hard-coded and we
935 // iterate over points on the first line which begin after the 3 vertex
936 // points in the complete list of unit support points
937
938 const std::vector<unsigned int> &index_map_inverse_q_other =
939 fe_q_other->get_poly_space_numbering_inverse();
940
941 std::vector<std::pair<unsigned int, unsigned int>> identities;
942
943 for (unsigned int i = 0; i < this->degree - 1; ++i)
944 for (unsigned int j = 0; j < fe_q_other->degree - 1; ++j)
945 if (std::fabs(this->unit_support_points[i + 3][0] -
946 fe_q_other->get_unit_support_points()
947 [index_map_inverse_q_other[j + 1]][0]) < 1e-14)
948 identities.emplace_back(i, j);
949 else
950 {
951 // If nodes are not located in the same place, we have to
952 // interpolate. This will then also
953 // capture the case where the FE_Q has a different polynomial
954 // degree than the current element. In either case, the resulting
955 // constraints are computed elsewhere, rather than via the
956 // identities this function returns: Since
957 // that happens in a different function, there is nothing
958 // for us to do here.
959 }
960
961 return identities;
962 }
963 else if (dynamic_cast<const FE_Nothing<dim> *>(&fe_other) != nullptr)
964 {
965 // The FE_Nothing has no degrees of freedom, so there are no
966 // equivalencies to be recorded. (If the FE_Nothing is dominating,
967 // then this will also leads to constraints, but we are not concerned
968 // with this here.)
969 return {};
970 }
971 else if (fe_other.n_unique_faces() == 1 && fe_other.n_dofs_per_face(0) == 0)
972 {
973 // if the other element has no elements on faces at all,
974 // then it would be impossible to enforce any kind of
975 // continuity even if we knew exactly what kind of element
976 // we have -- simply because the other element declares
977 // that it is discontinuous because it has no DoFs on
978 // its faces. in that case, just state that we have no
979 // constraints to declare
980 return {};
981 }
982 else
983 {
985 return {};
986 }
987}
988
989
990
991template <int dim, int spacedim>
993 : FE_SimplexPoly<dim, spacedim>(
994 BarycentricPolynomials<dim>::get_fe_p_basis(degree),
995 FiniteElementData<dim>(get_dpo_vector_fe_dgp(dim, degree),
996 ReferenceCells::get_simplex<dim>(),
997 1,
998 degree,
999 FiniteElementData<dim>::L2),
1000 true,
1001 unit_support_points_fe_p<dim>(degree),
1002 unit_face_support_points_fe_p<dim>(degree, FiniteElementData<dim>::L2),
1003 constraints_fe_p<dim>(degree))
1004{}
1005
1006
1007
1008template <int dim, int spacedim>
1009std::unique_ptr<FiniteElement<dim, spacedim>>
1011{
1012 return std::make_unique<FE_SimplexDGP<dim, spacedim>>(*this);
1013}
1014
1015
1016
1017template <int dim, int spacedim>
1018std::string
1020{
1021 std::ostringstream namebuf;
1022 namebuf << "FE_SimplexDGP<" << Utilities::dim_string(dim, spacedim) << ">("
1023 << this->degree << ")";
1024
1025 return namebuf.str();
1026}
1027
1028
1029template <int dim, int spacedim>
1032 const FiniteElement<dim, spacedim> &fe_other,
1033 const unsigned int codim) const
1034{
1035 Assert(codim <= dim, ExcImpossibleInDim(dim));
1036
1037 // vertex/line/face domination
1038 // ---------------------------
1039 if (codim > 0)
1040 // this is a discontinuous element, so by definition there will
1041 // be no constraints wherever this element comes together with
1042 // any other kind of element
1044
1045 // cell domination
1046 // ---------------
1047 if (const FE_SimplexDGP<dim, spacedim> *fe_dgp_other =
1048 dynamic_cast<const FE_SimplexDGP<dim, spacedim> *>(&fe_other))
1049 {
1050 if (this->degree < fe_dgp_other->degree)
1052 else if (this->degree == fe_dgp_other->degree)
1054 else
1056 }
1057 else if (const FE_DGQ<dim, spacedim> *fe_dgq_other =
1058 dynamic_cast<const FE_DGQ<dim, spacedim> *>(&fe_other))
1059 {
1060 if (this->degree < fe_dgq_other->degree)
1062 else if (this->degree == fe_dgq_other->degree)
1064 else
1066 }
1067 else if (const FE_Nothing<dim, spacedim> *fe_nothing =
1068 dynamic_cast<const FE_Nothing<dim, spacedim> *>(&fe_other))
1069 {
1070 if (fe_nothing->is_dominating())
1072 else
1073 // the FE_Nothing has no degrees of freedom and it is typically used
1074 // in a context where we don't require any continuity along the
1075 // interface
1077 }
1078
1081}
1082
1083
1084
1085template <int dim, int spacedim>
1086std::vector<std::pair<unsigned int, unsigned int>>
1088 const FiniteElement<dim, spacedim> &fe_other) const
1089{
1090 (void)fe_other;
1091
1092 return {};
1093}
1094
1095
1096
1097template <int dim, int spacedim>
1098std::vector<std::pair<unsigned int, unsigned int>>
1100 const FiniteElement<dim, spacedim> &fe_other) const
1101{
1102 (void)fe_other;
1103
1104 return {};
1105}
1106
1107
1108
1109template <int dim, int spacedim>
1110const FullMatrix<double> &
1112 const unsigned int child,
1113 const RefinementCase<dim> &refinement_case) const
1114{
1115 if (dim == 3)
1116 Assert(RefinementCase<dim>(refinement_case) ==
1118 static_cast<char>(IsotropicRefinementChoice::cut_tet_68)) ||
1119 RefinementCase<dim>(refinement_case) ==
1121 static_cast<char>(IsotropicRefinementChoice::cut_tet_57)) ||
1122 RefinementCase<dim>(refinement_case) ==
1124 static_cast<char>(IsotropicRefinementChoice::cut_tet_49)),
1126 else
1129 AssertDimension(dim, spacedim);
1130
1131 // initialization upon first request
1132 if (this->restriction[refinement_case - 1][child].n() == 0)
1133 {
1134 std::lock_guard<std::mutex> lock(this->restriction_matrix_mutex);
1135
1136 // if matrix got updated while waiting for the lock
1137 if (this->restriction[refinement_case - 1][child].n() ==
1138 this->n_dofs_per_cell())
1139 return this->restriction[refinement_case - 1][child];
1140
1141 // now do the work. need to get a non-const version of data in order to
1142 // be able to modify them inside a const function
1143 auto &this_nonconst = const_cast<FE_SimplexDGP<dim, spacedim> &>(*this);
1144
1145 if (dim == 2)
1146 {
1147 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
1149 isotropic_matrices.back().resize(
1150 this->reference_cell().n_children(
1151 RefinementCase<dim>(refinement_case)),
1152 FullMatrix<double>(this->n_dofs_per_cell(),
1153 this->n_dofs_per_cell()));
1154
1155 FETools::compute_projection_matrices(*this, isotropic_matrices, true);
1156
1157 this_nonconst.restriction[refinement_case - 1] =
1158 std::move(isotropic_matrices.back());
1159 }
1160 else if (dim == 3)
1161 {
1162 std::vector<std::vector<FullMatrix<double>>> matrices(
1163 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49),
1164 std::vector<FullMatrix<double>>(
1165 this->reference_cell().n_children(
1166 RefinementCase<dim>(refinement_case)),
1167 FullMatrix<double>(this->n_dofs_per_cell(),
1168 this->n_dofs_per_cell())));
1169 FETools::compute_projection_matrices(*this, matrices, true);
1170 for (unsigned int refinement_direction = static_cast<unsigned int>(
1172 refinement_direction <=
1173 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49);
1174 refinement_direction++)
1175 this_nonconst.restriction[refinement_direction - 1] =
1176 std::move(matrices[refinement_direction - 1]);
1177 }
1178 else
1180 }
1181
1182 // finally return the matrix
1183 return this->restriction[refinement_case - 1][child];
1184}
1185
1186// explicit instantiations
1187#include "fe_simplex_p.inst"
1188
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:949
static BarycentricPolynomials< dim > get_fe_p_basis(const unsigned int degree)
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
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source_fe, FullMatrix< double > &interpolation_matrix, const unsigned int face_no) const override
const unsigned int degree
Definition fe_data.h:452
unsigned int n_dofs_per_line() const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_unique_faces() const
std::vector< std::vector< Point< dim - 1 > > > unit_face_support_points
Definition fe.h:2458
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:2500
std::vector< Point< dim > > unit_support_points
Definition fe.h:2451
FullMatrix< double > interface_constraints
Definition fe.h:2439
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:111
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)
unsigned int n_vertices() const
static constexpr unsigned char default_combined_face_orientation()
unsigned int n_faces() const
unsigned int n_lines() const
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:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
unsigned int vertex_indices[2]
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define AssertThrow(cond, exc)
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:294
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)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
constexpr const ReferenceCell Tetrahedron
constexpr const ReferenceCell Triangle
constexpr const ReferenceCell Line
constexpr const ReferenceCell & get_simplex()
std::string dim_string(const int dim, const int spacedim)
Definition utilities.cc:555
STL namespace.