deal.II version GIT relicensing-2287-g6548a49e0a 2024-12-20 18:30:00+00:00
\(\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 {
39 const auto reference_cell = ReferenceCells::get_simplex<dim>();
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
222 // the following implements the 2d case
223 // (the 3d case is not implemented yet)
224 //
225 // consult FE_Q_Base::Implementation::initialize_constraints()
226 // for more information
227
228 std::vector<Point<dim - 1>> constraint_points;
229 // midpoint
230 constraint_points.emplace_back(0.5);
231 // subface 0
232 for (unsigned int i = 1; i < degree; ++i)
233 constraint_points.push_back(
235 Point<dim - 1>(i / double(degree)), 0));
236 // subface 1
237 for (unsigned int i = 1; i < degree; ++i)
238 constraint_points.push_back(
240 Point<dim - 1>(i / double(degree)), 1));
241
242 // Now construct relation between destination (child) and source (mother)
243 // dofs.
244
245 const unsigned int n_dofs_constrained = constraint_points.size();
246 unsigned int n_dofs_per_face = degree + 1;
247 FullMatrix<double> interface_constraints(n_dofs_constrained,
248 n_dofs_per_face);
249
250 const auto poly = BarycentricPolynomials<dim - 1>::get_fe_p_basis(degree);
251
252 for (unsigned int i = 0; i < n_dofs_constrained; ++i)
253 for (unsigned int j = 0; j < n_dofs_per_face; ++j)
254 {
255 interface_constraints(i, j) =
256 poly.compute_value(j, constraint_points[i]);
257
258 // if the value is small up to round-off, then simply set it to zero
259 // to avoid unwanted fill-in of the constraint matrices (which would
260 // then increase the number of other DoFs a constrained DoF would
261 // couple to)
262 if (std::fabs(interface_constraints(i, j)) < 1e-13)
263 interface_constraints(i, j) = 0;
264 }
265 return interface_constraints;
266 }
267
268
269
274 std::vector<unsigned int>
275 get_dpo_vector_fe_dgp(const unsigned int dim, const unsigned int degree)
276 {
277 // First treat the case of piecewise constant elements:
278 if (degree == 0)
279 {
280 std::vector<unsigned int> dpo(dim + 1, 0U);
281 dpo[dim] = 1;
282 return dpo;
283 }
284 else
285 {
286 // This element has the same degrees of freedom as the continuous one,
287 // but they are all counted for the interior of the cell because
288 // it is continuous. Rather than hard-code how many DoFs the element
289 // has, we just get the numbers from the continuous case and add them
290 // up
291 const auto continuous_dpo = get_dpo_vector_fe_p(dim, degree);
292
293 switch (dim)
294 {
295 case 1:
296 return {0U,
297 ReferenceCells::Line.n_vertices() * continuous_dpo[0] +
298 continuous_dpo[dim]};
299
300 case 2:
301 return {0U,
302 0U,
304 continuous_dpo[0] +
305 ReferenceCells::Triangle.n_lines() * continuous_dpo[1] +
306 continuous_dpo[dim]};
307
308 case 3:
309 return {
310 0U,
311 0U,
312 0U,
313 ReferenceCells::Tetrahedron.n_vertices() * continuous_dpo[0] +
314 ReferenceCells::Tetrahedron.n_lines() * continuous_dpo[1] +
315 ReferenceCells::Tetrahedron.n_faces() * continuous_dpo[2] +
316 continuous_dpo[dim]};
317 }
318
320 return {};
321 }
322 }
323} // namespace
324
325
326
327template <int dim, int spacedim>
329 const BarycentricPolynomials<dim> polynomials,
330 const FiniteElementData<dim> &fe_data,
331 const bool prolongation_is_additive,
332 const std::vector<Point<dim>> &unit_support_points,
333 const std::vector<std::vector<Point<dim - 1>>> unit_face_support_points,
334 const FullMatrix<double> &interface_constraints)
335 : ::FE_Poly<dim, spacedim>(
336 polynomials,
337 fe_data,
338 std::vector<bool>(fe_data.dofs_per_cell, prolongation_is_additive),
339 std::vector<ComponentMask>(fe_data.dofs_per_cell,
340 ComponentMask(std::vector<bool>(1, true))))
341{
344 this->interface_constraints = interface_constraints;
345}
346
347
348
349template <int dim, int spacedim>
350std::pair<Table<2, bool>, std::vector<unsigned int>>
352{
353 Table<2, bool> constant_modes(1, this->n_dofs_per_cell());
354 constant_modes.fill(true);
355 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
356 constant_modes, std::vector<unsigned int>(1, 0));
357}
358
359
360
361template <int dim, int spacedim>
362const FullMatrix<double> &
364 const unsigned int child,
365 const RefinementCase<dim> &refinement_case) const
366{
367 if (dim == 3)
368 Assert(RefinementCase<dim>(refinement_case) ==
370 static_cast<char>(IsotropicRefinementChoice::cut_tet_68)) ||
371 RefinementCase<dim>(refinement_case) ==
373 static_cast<char>(IsotropicRefinementChoice::cut_tet_57)) ||
374 RefinementCase<dim>(refinement_case) ==
376 static_cast<char>(IsotropicRefinementChoice::cut_tet_49)),
378 else
379 Assert(refinement_case ==
382 AssertDimension(dim, spacedim);
383
384 // initialization upon first request
385 if (this->prolongation[refinement_case - 1][child].n() == 0)
386 {
387 std::lock_guard<std::mutex> lock(prolongation_matrix_mutex);
388
389 // if matrix got updated while waiting for the lock
390 if (this->prolongation[refinement_case - 1][child].n() ==
391 this->n_dofs_per_cell())
392 return this->prolongation[refinement_case - 1][child];
393
394 // now do the work. need to get a non-const version of data in order to
395 // be able to modify them inside a const function
396 auto &this_nonconst = const_cast<FE_SimplexPoly<dim, spacedim> &>(*this);
397
398 if (dim == 2)
399 {
400 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
402 isotropic_matrices.back().resize(
403 this->reference_cell().n_children(
404 RefinementCase<dim>(refinement_case)),
405 FullMatrix<double>(this->n_dofs_per_cell(),
406 this->n_dofs_per_cell()));
407
408 FETools::compute_embedding_matrices(*this, isotropic_matrices, true);
409
410 this_nonconst.prolongation[refinement_case - 1] =
411 std::move(isotropic_matrices.back());
412 }
413 else if (dim == 3)
414 {
415 std::vector<std::vector<FullMatrix<double>>> matrices(
416 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49),
417 std::vector<FullMatrix<double>>(
418 this->reference_cell().n_children(
419 RefinementCase<dim>(refinement_case)),
420 FullMatrix<double>(this->n_dofs_per_cell(),
421 this->n_dofs_per_cell())));
422 FETools::compute_embedding_matrices(*this, matrices, true);
423 for (unsigned int refinement_direction = static_cast<unsigned int>(
425 refinement_direction <=
426 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49);
427 refinement_direction++)
428 this_nonconst.prolongation[refinement_direction - 1] =
429 std::move(matrices[refinement_direction - 1]);
430 }
431 else
433 }
434
435 // finally return the matrix
436 return this->prolongation[refinement_case - 1][child];
437}
438
439
440
441template <int dim, int spacedim>
442const FullMatrix<double> &
444 const unsigned int child,
445 const RefinementCase<dim> &refinement_case) const
446{
447 if (dim == 3)
448 Assert(RefinementCase<dim>(refinement_case) ==
450 static_cast<char>(IsotropicRefinementChoice::cut_tet_68)) ||
451 RefinementCase<dim>(refinement_case) ==
453 static_cast<char>(IsotropicRefinementChoice::cut_tet_57)) ||
454 RefinementCase<dim>(refinement_case) ==
456 static_cast<char>(IsotropicRefinementChoice::cut_tet_49)),
458 else
461 AssertDimension(dim, spacedim);
462
463 // initialization upon first request
464 if (this->restriction[refinement_case - 1][child].n() == 0)
465 {
466 std::lock_guard<std::mutex> lock(restriction_matrix_mutex);
467
468 // if matrix got updated while waiting for the lock
469 if (this->restriction[refinement_case - 1][child].n() ==
470 this->n_dofs_per_cell())
471 return this->restriction[refinement_case - 1][child];
472
473 // get the restriction matrix
474 // Refine a unit cell. As the parent cell is a unit
475 // cell, the reference cell of the children equals the parent, i.e. they
476 // have the support points at the same locations. So we just have to check
477 // if a support point of the parent is one of the interpolation points of
478 // the child. If this is not the case we find the interpolation of the
479 // point.
480
481 const double eps = 1e-12;
482 FullMatrix<double> restriction_mat(this->n_dofs_per_cell(),
483 this->n_dofs_per_cell());
484
485 // first get all support points on the reference cell
486 const std::vector<Point<dim>> unit_support_points =
487 this->get_unit_support_points();
488
489 // now create children on the reference cell
491 GridGenerator::reference_cell(tria, this->reference_cell());
492 tria.begin_active()->set_refine_flag(
494 if (dim == 3)
495 tria.begin_active()->set_refine_choice(refinement_case);
497
498 const auto &child_cell = tria.begin(0)->child(child);
499
500 // iterate over all support points and transform them to the unit cell of
501 // the child
502 for (unsigned int i = 0; i < unit_support_points.size(); i++)
503 {
504 std::vector<Point<dim>> transformed_point(1);
505 const std::vector<Point<spacedim>> unit_support_point = {
506 dim == 2 ? Point<spacedim>(unit_support_points[i][0],
507 unit_support_points[i][1]) :
508 Point<spacedim>(unit_support_points[i][0],
509 unit_support_points[i][1],
510 unit_support_points[i][2])};
511 this->reference_cell()
512 .template get_default_linear_mapping<dim, spacedim>()
513 .transform_points_real_to_unit_cell(
514 child_cell,
515 make_array_view(unit_support_point),
516 make_array_view(transformed_point));
517
518 // if point is inside the unit cell iterate over all shape functions
519 if (this->reference_cell().contains_point(transformed_point[0], eps))
520 for (unsigned int j = 0; j < this->n_dofs_per_cell(); j++)
521 restriction_mat[i][j] =
522 this->shape_value(j, transformed_point[0]);
523 }
524#ifdef DEBUG
525 for (unsigned int i = 0; i < this->n_dofs_per_cell(); i++)
526 {
527 double sum = 0.;
528
529 for (unsigned int j = 0; j < this->n_dofs_per_cell(); j++)
530 sum += restriction_mat[i][j];
531
532 Assert(std::fabs(sum - 1) < eps || std::fabs(sum) < eps,
534 "The entries in a row of the local "
535 "restriction matrix do not add to zero or one. "
536 "This typically indicates that the "
537 "polynomial interpolation is "
538 "ill-conditioned such that round-off "
539 "prevents the sum to be one."));
540 }
541#endif
542
543 // Remove small entries from the matrix
544 for (unsigned int i = 0; i < restriction_mat.m(); ++i)
545 for (unsigned int j = 0; j < restriction_mat.n(); ++j)
546 {
547 if (std::fabs(restriction_mat(i, j)) < eps)
548 restriction_mat(i, j) = 0.;
549 if (std::fabs(restriction_mat(i, j) - 1) < eps)
550 restriction_mat(i, j) = 1.;
551 }
552
553 const_cast<FullMatrix<double> &>(
554 this->restriction[refinement_case - 1][child]) =
555 std::move(restriction_mat);
556 }
557
558 // finally return the matrix
559 return this->restriction[refinement_case - 1][child];
560}
561
562
563
564template <int dim, int spacedim>
565void
567 const FiniteElement<dim, spacedim> &source_fe,
568 FullMatrix<double> &interpolation_matrix,
569 const unsigned int face_no) const
570{
571 Assert(interpolation_matrix.m() == source_fe.n_dofs_per_face(face_no),
572 ExcDimensionMismatch(interpolation_matrix.m(),
573 source_fe.n_dofs_per_face(face_no)));
574
575 // see if source is a P or Q element
576 if ((dynamic_cast<const FE_SimplexPoly<dim, spacedim> *>(&source_fe) !=
577 nullptr) ||
578 (dynamic_cast<const FE_Q_Base<dim, spacedim> *>(&source_fe) != nullptr))
579 {
580 const Quadrature<dim - 1> quad_face_support(
581 source_fe.get_unit_face_support_points(face_no));
582
583 const double eps = 2e-13 * this->degree * (dim - 1);
584
585 std::vector<Point<dim>> face_quadrature_points(quad_face_support.size());
586 QProjector<dim>::project_to_face(this->reference_cell(),
587 quad_face_support,
588 face_no,
589 face_quadrature_points);
590
591 for (unsigned int i = 0; i < source_fe.n_dofs_per_face(face_no); ++i)
592 for (unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
593 {
594 double matrix_entry =
595 this->shape_value(this->face_to_cell_index(j, 0),
596 face_quadrature_points[i]);
597
598 // Correct the interpolated value. I.e. if it is close to 1 or
599 // 0, make it exactly 1 or 0. Unfortunately, this is required to
600 // avoid problems with higher order elements.
601 if (std::fabs(matrix_entry - 1.0) < eps)
602 matrix_entry = 1.0;
603 if (std::fabs(matrix_entry) < eps)
604 matrix_entry = 0.0;
605
606 interpolation_matrix(i, j) = matrix_entry;
607 }
608
609#ifdef DEBUG
610 for (unsigned int j = 0; j < source_fe.n_dofs_per_face(face_no); ++j)
611 {
612 double sum = 0.;
613
614 for (unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
615 sum += interpolation_matrix(j, i);
616
617 Assert(std::fabs(sum - 1) < eps, ExcInternalError());
618 }
619#endif
620 }
621 else if (dynamic_cast<const FE_Nothing<dim> *>(&source_fe) != nullptr)
622 {
623 // nothing to do here, the FE_Nothing has no degrees of freedom anyway
624 }
625 else
627 false,
628 (typename FiniteElement<dim,
629 spacedim>::ExcInterpolationNotImplemented()));
630}
631
632
633
634template <int dim, int spacedim>
635void
637 const FiniteElement<dim, spacedim> &source_fe,
638 const unsigned int subface,
639 FullMatrix<double> &interpolation_matrix,
640 const unsigned int face_no) const
641{
642 Assert(interpolation_matrix.m() == source_fe.n_dofs_per_face(face_no),
643 ExcDimensionMismatch(interpolation_matrix.m(),
644 source_fe.n_dofs_per_face(face_no)));
645
646 // see if source is a P or Q element
647 if ((dynamic_cast<const FE_SimplexPoly<dim, spacedim> *>(&source_fe) !=
648 nullptr) ||
649 (dynamic_cast<const FE_Q_Base<dim, spacedim> *>(&source_fe) != nullptr))
650 {
651 const Quadrature<dim - 1> quad_face_support(
652 source_fe.get_unit_face_support_points(face_no));
653
654 const double eps = 2e-13 * this->degree * (dim - 1);
655
656 std::vector<Point<dim>> subface_quadrature_points(
657 quad_face_support.size());
658 QProjector<dim>::project_to_subface(this->reference_cell(),
659 quad_face_support,
660 face_no,
661 subface,
662 subface_quadrature_points);
663
664 for (unsigned int i = 0; i < source_fe.n_dofs_per_face(face_no); ++i)
665 for (unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
666 {
667 double matrix_entry =
668 this->shape_value(this->face_to_cell_index(j, 0),
669 subface_quadrature_points[i]);
670
671 // Correct the interpolated value. I.e. if it is close to 1 or
672 // 0, make it exactly 1 or 0. Unfortunately, this is required to
673 // avoid problems with higher order elements.
674 if (std::fabs(matrix_entry - 1.0) < eps)
675 matrix_entry = 1.0;
676 if (std::fabs(matrix_entry) < eps)
677 matrix_entry = 0.0;
678
679 interpolation_matrix(i, j) = matrix_entry;
680 }
681
682#ifdef DEBUG
683 for (unsigned int j = 0; j < source_fe.n_dofs_per_face(face_no); ++j)
684 {
685 double sum = 0.;
686
687 for (unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
688 sum += interpolation_matrix(j, i);
689
690 Assert(std::fabs(sum - 1) < eps, ExcInternalError());
691 }
692#endif
693 }
694 else if (dynamic_cast<const FE_Nothing<dim> *>(&source_fe) != nullptr)
695 {
696 // nothing to do here, the FE_Nothing has no degrees of freedom anyway
697 }
698 else
700 false,
701 (typename FiniteElement<dim,
702 spacedim>::ExcInterpolationNotImplemented()));
703}
704
705
706
707template <int dim, int spacedim>
708bool
713
714
715
716template <int dim, int spacedim>
717void
720 const std::vector<Vector<double>> &support_point_values,
721 std::vector<double> &nodal_values) const
722{
723 AssertDimension(support_point_values.size(),
724 this->get_unit_support_points().size());
725 AssertDimension(support_point_values.size(), nodal_values.size());
726 AssertDimension(this->dofs_per_cell, nodal_values.size());
727
728 for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
729 {
730 AssertDimension(support_point_values[i].size(), 1);
731
732 nodal_values[i] = support_point_values[i](0);
733 }
734}
735
736
737
738template <int dim, int spacedim>
740 : FE_SimplexPoly<dim, spacedim>(
741 BarycentricPolynomials<dim>::get_fe_p_basis(degree),
742 FiniteElementData<dim>(get_dpo_vector_fe_p(dim, degree),
743 ReferenceCells::get_simplex<dim>(),
744 1,
745 degree,
746 FiniteElementData<dim>::H1),
747 false,
748 unit_support_points_fe_p<dim>(degree),
749 unit_face_support_points_fe_p<dim>(degree, FiniteElementData<dim>::H1),
750 constraints_fe_p<dim>(degree))
751{
752 if (degree > 2)
753 for (unsigned int i = 0; i < this->n_dofs_per_line(); ++i)
755 this->n_dofs_per_line() - 1 - i - i;
756 // We do not support multiple DoFs per quad yet
758}
759
760
761
762template <int dim, int spacedim>
763std::unique_ptr<FiniteElement<dim, spacedim>>
765{
766 return std::make_unique<FE_SimplexP<dim, spacedim>>(*this);
767}
768
769
770
771template <int dim, int spacedim>
772std::string
774{
775 std::ostringstream namebuf;
776 namebuf << "FE_SimplexP<" << Utilities::dim_string(dim, spacedim) << ">("
777 << this->degree << ")";
778
779 return namebuf.str();
780}
781
782
783
784template <int dim, int spacedim>
787 const FiniteElement<dim, spacedim> &fe_other,
788 const unsigned int codim) const
789{
790 Assert(codim <= dim, ExcImpossibleInDim(dim));
791
792 // vertex/line/face domination
793 // (if fe_other is derived from FE_SimplexDGP)
794 // ------------------------------------
795 if (codim > 0)
796 if (dynamic_cast<const FE_SimplexDGP<dim, spacedim> *>(&fe_other) !=
797 nullptr)
798 // there are no requirements between continuous and discontinuous
799 // elements
801
802 // vertex/line/face domination
803 // (if fe_other is not derived from FE_SimplexDGP)
804 // & cell domination
805 // ----------------------------------------
806 if (const FE_SimplexP<dim, spacedim> *fe_p_other =
807 dynamic_cast<const FE_SimplexP<dim, spacedim> *>(&fe_other))
808 {
809 if (this->degree < fe_p_other->degree)
811 else if (this->degree == fe_p_other->degree)
813 else
815 }
816 else if (const FE_Q<dim, spacedim> *fe_q_other =
817 dynamic_cast<const FE_Q<dim, spacedim> *>(&fe_other))
818 {
819 if (this->degree < fe_q_other->degree)
821 else if (this->degree == fe_q_other->degree)
823 else
825 }
826 else if (const FE_Nothing<dim, spacedim> *fe_nothing =
827 dynamic_cast<const FE_Nothing<dim, spacedim> *>(&fe_other))
828 {
829 if (fe_nothing->is_dominating())
831 else
832 // the FE_Nothing has no degrees of freedom and it is typically used
833 // in a context where we don't require any continuity along the
834 // interface
836 }
837
840}
841
842
843
844template <int dim, int spacedim>
845std::vector<std::pair<unsigned int, unsigned int>>
847 const FiniteElement<dim, spacedim> &fe_other) const
848{
849 AssertDimension(dim, 2);
850
851 if (dynamic_cast<const FE_SimplexP<dim, spacedim> *>(&fe_other) != nullptr)
852 {
853 // there should be exactly one single DoF of each FE at a vertex, and
854 // they should have identical value
855 return {{0U, 0U}};
856 }
857 else if (dynamic_cast<const FE_Q<dim, spacedim> *>(&fe_other) != nullptr)
858 {
859 // there should be exactly one single DoF of each FE at a vertex, and
860 // they should have identical value
861 return {{0U, 0U}};
862 }
863 else if (dynamic_cast<const FE_Nothing<dim> *>(&fe_other) != nullptr)
864 {
865 // the FE_Nothing has no degrees of freedom, so there are no
866 // equivalencies to be recorded
867 return {};
868 }
869 else if (fe_other.n_unique_faces() == 1 && fe_other.n_dofs_per_face(0) == 0)
870 {
871 // if the other element has no DoFs on faces at all,
872 // then it would be impossible to enforce any kind of
873 // continuity even if we knew exactly what kind of element
874 // we have -- simply because the other element declares
875 // that it is discontinuous because it has no DoFs on
876 // its faces. in that case, just state that we have no
877 // constraints to declare
878 return {};
879 }
880 else
881 {
883 return {};
884 }
885}
886
887
888
889template <int dim, int spacedim>
890std::vector<std::pair<unsigned int, unsigned int>>
892 const FiniteElement<dim, spacedim> &fe_other) const
893{
894 AssertDimension(dim, 2);
895
896 if (const FE_SimplexP<dim, spacedim> *fe_p_other =
897 dynamic_cast<const FE_SimplexP<dim, spacedim> *>(&fe_other))
898 {
899 // dofs are located along lines, so two dofs are identical if they are
900 // located at identical positions.
901 // Therefore, read the points in unit_support_points for the
902 // first coordinate direction. For FE_SimplexP, they are currently
903 // hard-coded and we iterate over points on the first line which begin
904 // after the 3 vertex points in the complete list of unit support points
905
906 std::vector<std::pair<unsigned int, unsigned int>> identities;
907
908 for (unsigned int i = 0; i < this->degree - 1; ++i)
909 for (unsigned int j = 0; j < fe_p_other->degree - 1; ++j)
910 if (std::fabs(this->unit_support_points[i + 3][0] -
911 fe_p_other->unit_support_points[i + 3][0]) < 1e-14)
912 identities.emplace_back(i, j);
913 else
914 {
915 // If nodes are not located in the same place, we have to
916 // interpolate. This is then not handled through the
917 // current function, but via interpolation matrices that
918 // result in constraints, rather than identities. Since
919 // that happens in a different function, there is nothing
920 // for us to do here.
921 }
922
923 return identities;
924 }
925 else if (const FE_Q<dim, spacedim> *fe_q_other =
926 dynamic_cast<const FE_Q<dim, spacedim> *>(&fe_other))
927 {
928 // dofs are located along lines, so two dofs are identical if they are
929 // located at identical positions. if we had only equidistant points, we
930 // could simply check for similarity like (i+1)*q == (j+1)*p, but we
931 // might have other support points (e.g. Gauss-Lobatto
932 // points). Therefore, read the points in unit_support_points for the
933 // first coordinate direction. For FE_Q, we take the lexicographic
934 // ordering of the line support points in the first direction (i.e.,
935 // x-direction), which we access between index 1 and p-1 (index 0 and p
936 // are vertex dofs). For FE_SimplexP, they are currently hard-coded and we
937 // iterate over points on the first line which begin after the 3 vertex
938 // points in the complete list of unit support points
939
940 const std::vector<unsigned int> &index_map_inverse_q_other =
941 fe_q_other->get_poly_space_numbering_inverse();
942
943 std::vector<std::pair<unsigned int, unsigned int>> identities;
944
945 for (unsigned int i = 0; i < this->degree - 1; ++i)
946 for (unsigned int j = 0; j < fe_q_other->degree - 1; ++j)
947 if (std::fabs(this->unit_support_points[i + 3][0] -
948 fe_q_other->get_unit_support_points()
949 [index_map_inverse_q_other[j + 1]][0]) < 1e-14)
950 identities.emplace_back(i, j);
951 else
952 {
953 // If nodes are not located in the same place, we have to
954 // interpolate. This will then also
955 // capture the case where the FE_Q has a different polynomial
956 // degree than the current element. In either case, the resulting
957 // constraints are computed elsewhere, rather than via the
958 // identities this function returns: Since
959 // that happens in a different function, there is nothing
960 // for us to do here.
961 }
962
963 return identities;
964 }
965 else if (dynamic_cast<const FE_Nothing<dim> *>(&fe_other) != nullptr)
966 {
967 // The FE_Nothing has no degrees of freedom, so there are no
968 // equivalencies to be recorded. (If the FE_Nothing is dominating,
969 // then this will also leads to constraints, but we are not concerned
970 // with this here.)
971 return {};
972 }
973 else if (fe_other.n_unique_faces() == 1 && fe_other.n_dofs_per_face(0) == 0)
974 {
975 // if the other element has no elements on faces at all,
976 // then it would be impossible to enforce any kind of
977 // continuity even if we knew exactly what kind of element
978 // we have -- simply because the other element declares
979 // that it is discontinuous because it has no DoFs on
980 // its faces. in that case, just state that we have no
981 // constraints to declare
982 return {};
983 }
984 else
985 {
987 return {};
988 }
989}
990
991
992
993template <int dim, int spacedim>
995 : FE_SimplexPoly<dim, spacedim>(
996 BarycentricPolynomials<dim>::get_fe_p_basis(degree),
997 FiniteElementData<dim>(get_dpo_vector_fe_dgp(dim, degree),
998 ReferenceCells::get_simplex<dim>(),
999 1,
1000 degree,
1001 FiniteElementData<dim>::L2),
1002 true,
1003 unit_support_points_fe_p<dim>(degree),
1004 unit_face_support_points_fe_p<dim>(degree, FiniteElementData<dim>::L2),
1005 constraints_fe_p<dim>(degree))
1006{}
1007
1008
1009
1010template <int dim, int spacedim>
1011std::unique_ptr<FiniteElement<dim, spacedim>>
1013{
1014 return std::make_unique<FE_SimplexDGP<dim, spacedim>>(*this);
1015}
1016
1017
1018
1019template <int dim, int spacedim>
1020std::string
1022{
1023 std::ostringstream namebuf;
1024 namebuf << "FE_SimplexDGP<" << Utilities::dim_string(dim, spacedim) << ">("
1025 << this->degree << ")";
1026
1027 return namebuf.str();
1028}
1029
1030
1031template <int dim, int spacedim>
1034 const FiniteElement<dim, spacedim> &fe_other,
1035 const unsigned int codim) const
1036{
1037 Assert(codim <= dim, ExcImpossibleInDim(dim));
1038
1039 // vertex/line/face domination
1040 // ---------------------------
1041 if (codim > 0)
1042 // this is a discontinuous element, so by definition there will
1043 // be no constraints wherever this element comes together with
1044 // any other kind of element
1046
1047 // cell domination
1048 // ---------------
1049 if (const FE_SimplexDGP<dim, spacedim> *fe_dgp_other =
1050 dynamic_cast<const FE_SimplexDGP<dim, spacedim> *>(&fe_other))
1051 {
1052 if (this->degree < fe_dgp_other->degree)
1054 else if (this->degree == fe_dgp_other->degree)
1056 else
1058 }
1059 else if (const FE_DGQ<dim, spacedim> *fe_dgq_other =
1060 dynamic_cast<const FE_DGQ<dim, spacedim> *>(&fe_other))
1061 {
1062 if (this->degree < fe_dgq_other->degree)
1064 else if (this->degree == fe_dgq_other->degree)
1066 else
1068 }
1069 else if (const FE_Nothing<dim, spacedim> *fe_nothing =
1070 dynamic_cast<const FE_Nothing<dim, spacedim> *>(&fe_other))
1071 {
1072 if (fe_nothing->is_dominating())
1074 else
1075 // the FE_Nothing has no degrees of freedom and it is typically used
1076 // in a context where we don't require any continuity along the
1077 // interface
1079 }
1080
1083}
1084
1085
1086
1087template <int dim, int spacedim>
1088std::vector<std::pair<unsigned int, unsigned int>>
1090 const FiniteElement<dim, spacedim> &fe_other) const
1091{
1092 (void)fe_other;
1093
1094 return {};
1095}
1096
1097
1098
1099template <int dim, int spacedim>
1100std::vector<std::pair<unsigned int, unsigned int>>
1102 const FiniteElement<dim, spacedim> &fe_other) const
1103{
1104 (void)fe_other;
1105
1106 return {};
1107}
1108
1109
1110
1111template <int dim, int spacedim>
1112const FullMatrix<double> &
1114 const unsigned int child,
1115 const RefinementCase<dim> &refinement_case) const
1116{
1117 if (dim == 3)
1118 Assert(RefinementCase<dim>(refinement_case) ==
1120 static_cast<char>(IsotropicRefinementChoice::cut_tet_68)) ||
1121 RefinementCase<dim>(refinement_case) ==
1123 static_cast<char>(IsotropicRefinementChoice::cut_tet_57)) ||
1124 RefinementCase<dim>(refinement_case) ==
1126 static_cast<char>(IsotropicRefinementChoice::cut_tet_49)),
1128 else
1131 AssertDimension(dim, spacedim);
1132
1133 // initialization upon first request
1134 if (this->restriction[refinement_case - 1][child].n() == 0)
1135 {
1136 std::lock_guard<std::mutex> lock(this->restriction_matrix_mutex);
1137
1138 // if matrix got updated while waiting for the lock
1139 if (this->restriction[refinement_case - 1][child].n() ==
1140 this->n_dofs_per_cell())
1141 return this->restriction[refinement_case - 1][child];
1142
1143 // now do the work. need to get a non-const version of data in order to
1144 // be able to modify them inside a const function
1145 auto &this_nonconst = const_cast<FE_SimplexDGP<dim, spacedim> &>(*this);
1146
1147 if (dim == 2)
1148 {
1149 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
1151 isotropic_matrices.back().resize(
1152 this->reference_cell().n_children(
1153 RefinementCase<dim>(refinement_case)),
1154 FullMatrix<double>(this->n_dofs_per_cell(),
1155 this->n_dofs_per_cell()));
1156
1157 FETools::compute_projection_matrices(*this, isotropic_matrices, true);
1158
1159 this_nonconst.restriction[refinement_case - 1] =
1160 std::move(isotropic_matrices.back());
1161 }
1162 else if (dim == 3)
1163 {
1164 std::vector<std::vector<FullMatrix<double>>> matrices(
1165 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49),
1166 std::vector<FullMatrix<double>>(
1167 this->reference_cell().n_children(
1168 RefinementCase<dim>(refinement_case)),
1169 FullMatrix<double>(this->n_dofs_per_cell(),
1170 this->n_dofs_per_cell())));
1171 FETools::compute_projection_matrices(*this, matrices, true);
1172 for (unsigned int refinement_direction = static_cast<unsigned int>(
1174 refinement_direction <=
1175 static_cast<unsigned int>(IsotropicRefinementChoice::cut_tet_49);
1176 refinement_direction++)
1177 this_nonconst.restriction[refinement_direction - 1] =
1178 std::move(matrices[refinement_direction - 1]);
1179 }
1180 else
1182 }
1183
1184 // finally return the matrix
1185 return this->restriction[refinement_case - 1][child];
1186}
1187
1188// explicit instantiations
1189#include "fe_simplex_p.inst"
1190
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:2495
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:2537
std::vector< Point< dim > > unit_support_points
Definition fe.h:2488
FullMatrix< double > interface_constraints
Definition fe.h:2476
size_type n() const
size_type m() 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:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
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()
std::size_t size
Definition mpi.cc:734
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.