41 inline std::vector<unsigned int>
42 invert_numbering(
const std::vector<unsigned int> &in)
44 std::vector<unsigned int> out(in.size());
45 for (
unsigned int i = 0; i < in.size(); ++i)
61 Polynomials::Hierarchical::generate_complete_basis(degree)),
74 , face_renumber(face_fe_q_hierarchical_to_hierarchic_numbering(degree))
95 std::vector<FullMatrix<double>> dofs_cell(
100 std::vector<FullMatrix<double>> dofs_subcell(
134 std::ostringstream namebuf;
135 namebuf <<
"FE_Q_Hierarchical<" << dim <<
">(" << this->degree <<
")";
137 return namebuf.str();
143std::unique_ptr<FiniteElement<dim, dim>>
146 return std::make_unique<FE_Q_Hierarchical<dim>>(*this);
163 Assert(matrix.m() == this->n_dofs_per_cell(),
177 if (this->n_dofs_per_cell() >= source_fe->n_dofs_per_cell())
179 const std::vector<unsigned int> dof_map =
180 this->get_embedding_dofs(source_fe->degree);
181 for (
unsigned int j = 0; j < dof_map.size(); ++j)
182 matrix[dof_map[j]][j] = 1.;
187 const std::vector<unsigned int> dof_map =
188 source_fe->get_embedding_dofs(this->degree);
189 for (
unsigned int j = 0; j < dof_map.size(); ++j)
190 matrix[j][dof_map[j]] = 1.;
198 "Interpolation is supported only between FE_Q_Hierarchical"));
205 const unsigned int child,
211 "Prolongation matrices are only available for isotropic refinement!"));
215 return this->prolongation[refinement_case - 1][child];
228std::vector<std::pair<unsigned int, unsigned int>>
242 return std::vector<std::pair<unsigned int, unsigned int>>(
243 1, std::make_pair(0U, 0U));
251 return std::vector<std::pair<unsigned int, unsigned int>>();
256 return std::vector<std::pair<unsigned int, unsigned int>>();
261std::vector<std::pair<unsigned int, unsigned int>>
271 const unsigned int this_dpl = this->n_dofs_per_line();
277 std::vector<std::pair<unsigned int, unsigned int>> res;
278 for (
unsigned int i = 0; i <
std::min(this_dpl, other_dpl); ++i)
279 res.emplace_back(i, i);
289 return std::vector<std::pair<unsigned int, unsigned int>>();
294 return std::vector<std::pair<unsigned int, unsigned int>>();
299std::vector<std::pair<unsigned int, unsigned int>>
302 const unsigned int)
const
313 const unsigned int face_no = 0;
315 const unsigned int this_dpq = this->n_dofs_per_quad(face_no);
321 std::vector<std::pair<unsigned int, unsigned int>> res;
322 for (
unsigned int i = 0; i <
std::min(this_dpq, other_dpq); ++i)
323 res.emplace_back(i, i);
333 return std::vector<std::pair<unsigned int, unsigned int>>();
338 return std::vector<std::pair<unsigned int, unsigned int>>();
346 const unsigned int codim)
const
354 if (
dynamic_cast<const FE_DGQ<dim> *
>(&fe_other) !=
nullptr)
365 if (this->degree < fe_hierarchical_other->degree)
367 else if (this->degree == fe_hierarchical_other->degree)
375 if (fe_nothing->is_dominating())
400 const unsigned int dofs_1d =
401 2 * this->n_dofs_per_vertex() + this->n_dofs_per_line();
450 for (
unsigned int c = 0; c < GeometryInfo<1>::max_children_per_cell; ++c)
451 for (
unsigned int j = 0; j < dofs_1d; ++j)
452 for (
unsigned int k = 0; k < dofs_1d; ++k)
455 if ((j <= 1) && (k <= 1))
457 if (((c == 0) && (j == 0) && (k == 0)) ||
458 ((c == 1) && (j == 1) && (k == 1)))
459 dofs_cell[c](j, k) = 1.;
461 dofs_cell[c](j, k) = 0.;
463 if (((c == 0) && (j == 1)) || ((c == 1) && (j == 0)))
464 dofs_subcell[c](j, k) = .5;
465 else if (((c == 0) && (k == 0)) || ((c == 1) && (k == 1)))
466 dofs_subcell[c](j, k) = 1.;
468 dofs_subcell[c](j, k) = 0.;
471 else if ((j <= 1) && (k >= 2))
473 if (((c == 0) && (j == 1) && ((k % 2) == 0)) ||
474 ((c == 1) && (j == 0) && ((k % 2) == 0)))
475 dofs_subcell[c](j, k) = -1.;
478 else if ((j >= 2) && (k >= 2) && (j <= k))
481 for (
unsigned int i = 1; i <= j; ++i)
483 (
static_cast<double>(k - i + 1)) / (
static_cast<double>(i));
487 dofs_subcell[c](j, k) =
489 std::pow(.5,
static_cast<double>(k)) * factor :
490 -
std::pow(.5,
static_cast<double>(k)) * factor;
492 std::pow(2.,
static_cast<double>(j)) * factor;
496 dofs_subcell[c](j, k) =
497 std::pow(.5,
static_cast<double>(k)) * factor;
500 std::pow(2.,
static_cast<double>(j)) * factor :
501 -
std::pow(2.,
static_cast<double>(j)) * factor;
514 const unsigned int dofs_1d =
515 2 * this->n_dofs_per_vertex() + this->n_dofs_per_line();
517 this->interface_constraints.TableBase<2, double>::reinit(
518 this->interface_constraints_size());
531 for (
unsigned int i = 0; i < dofs_1d; ++i)
532 this->interface_constraints(0, i) = dofs_subcell[0](1, i);
534 for (
unsigned int c = 0; c < GeometryInfo<1>::max_children_per_cell;
536 for (
unsigned int i = 0; i < dofs_1d; ++i)
537 for (
unsigned int j = 2; j < dofs_1d; ++j)
538 this->interface_constraints(1 + c * (this->degree - 1) + j - 2,
539 i) = dofs_subcell[c](j, i);
545 for (
unsigned int i = 0; i < dofs_1d * dofs_1d; ++i)
548 this->interface_constraints(0, face_renumber[i]) =
549 dofs_subcell[0](1, i % dofs_1d) *
550 dofs_subcell[0](1, (i - (i % dofs_1d)) / dofs_1d);
553 this->interface_constraints(1, face_renumber[i]) =
554 dofs_subcell[0](0, i % dofs_1d) *
555 dofs_subcell[0](1, (i - (i % dofs_1d)) / dofs_1d);
556 this->interface_constraints(2, face_renumber[i]) =
557 dofs_subcell[1](1, i % dofs_1d) *
558 dofs_subcell[0](1, (i - (i % dofs_1d)) / dofs_1d);
559 this->interface_constraints(3, face_renumber[i]) =
560 dofs_subcell[0](1, i % dofs_1d) *
561 dofs_subcell[0](0, (i - (i % dofs_1d)) / dofs_1d);
562 this->interface_constraints(4, face_renumber[i]) =
563 dofs_subcell[1](0, i % dofs_1d) *
564 dofs_subcell[1](1, (i - (i % dofs_1d)) / dofs_1d);
567 for (
unsigned int j = 0; j < (this->degree - 1); ++j)
569 this->interface_constraints(5 + j, face_renumber[i]) =
570 dofs_subcell[0](1, i % dofs_1d) *
571 dofs_subcell[0](2 + j, (i - (i % dofs_1d)) / dofs_1d);
572 this->interface_constraints(5 + (this->degree - 1) + j,
574 dofs_subcell[0](1, i % dofs_1d) *
575 dofs_subcell[1](2 + j, (i - (i % dofs_1d)) / dofs_1d);
576 this->interface_constraints(5 + 2 * (this->degree - 1) + j,
578 dofs_subcell[0](2 + j, i % dofs_1d) *
579 dofs_subcell[1](0, (i - (i % dofs_1d)) / dofs_1d);
580 this->interface_constraints(5 + 3 * (this->degree - 1) + j,
582 dofs_subcell[1](2 + j, i % dofs_1d) *
583 dofs_subcell[0](1, (i - (i % dofs_1d)) / dofs_1d);
587 for (
unsigned int j = 0; j < (this->degree - 1); ++j)
590 this->interface_constraints(5 + 4 * (this->degree - 1) + j,
592 dofs_subcell[0](0, i % dofs_1d) *
593 dofs_subcell[0](2 + j, (i - (i % dofs_1d)) / dofs_1d);
594 this->interface_constraints(5 + 4 * (this->degree - 1) +
595 (this->degree - 1) + j,
597 dofs_subcell[0](0, i % dofs_1d) *
598 dofs_subcell[1](2 + j, (i - (i % dofs_1d)) / dofs_1d);
600 this->interface_constraints(5 + 4 * (this->degree - 1) +
601 2 * (this->degree - 1) + j,
603 dofs_subcell[1](1, i % dofs_1d) *
604 dofs_subcell[0](2 + j, (i - (i % dofs_1d)) / dofs_1d);
605 this->interface_constraints(5 + 4 * (this->degree - 1) +
606 3 * (this->degree - 1) + j,
608 dofs_subcell[1](1, i % dofs_1d) *
609 dofs_subcell[1](2 + j, (i - (i % dofs_1d)) / dofs_1d);
611 this->interface_constraints(5 + 4 * (this->degree - 1) +
612 4 * (this->degree - 1) + j,
614 dofs_subcell[0](2 + j, i % dofs_1d) *
615 dofs_subcell[0](0, (i - (i % dofs_1d)) / dofs_1d);
616 this->interface_constraints(5 + 4 * (this->degree - 1) +
617 5 * (this->degree - 1) + j,
619 dofs_subcell[1](2 + j, i % dofs_1d) *
620 dofs_subcell[0](0, (i - (i % dofs_1d)) / dofs_1d);
622 this->interface_constraints(5 + 4 * (this->degree - 1) +
623 6 * (this->degree - 1) + j,
625 dofs_subcell[0](2 + j, i % dofs_1d) *
626 dofs_subcell[1](1, (i - (i % dofs_1d)) / dofs_1d);
627 this->interface_constraints(5 + 4 * (this->degree - 1) +
628 7 * (this->degree - 1) + j,
630 dofs_subcell[1](2 + j, i % dofs_1d) *
631 dofs_subcell[1](1, (i - (i % dofs_1d)) / dofs_1d);
635 for (
unsigned int j = 0; j < (this->degree - 1); ++j)
636 for (
unsigned int k = 0; k < (this->degree - 1); ++k)
639 this->interface_constraints(5 + 12 * (this->degree - 1) +
640 j + k * (this->degree - 1),
642 dofs_subcell[0](2 + j, i % dofs_1d) *
643 dofs_subcell[0](2 + k, (i - (i % dofs_1d)) / dofs_1d);
645 this->interface_constraints(5 + 12 * (this->degree - 1) +
646 j + k * (this->degree - 1) +
650 dofs_subcell[1](2 + j, i % dofs_1d) *
651 dofs_subcell[0](2 + k, (i - (i % dofs_1d)) / dofs_1d);
653 this->interface_constraints(5 + 12 * (this->degree - 1) +
654 j + k * (this->degree - 1) +
655 2 * (this->degree - 1) *
658 dofs_subcell[0](2 + j, i % dofs_1d) *
659 dofs_subcell[1](2 + k, (i - (i % dofs_1d)) / dofs_1d);
661 this->interface_constraints(5 + 12 * (this->degree - 1) +
662 j + k * (this->degree - 1) +
663 3 * (this->degree - 1) *
666 dofs_subcell[1](2 + j, i % dofs_1d) *
667 dofs_subcell[1](2 + k, (i - (i % dofs_1d)) / dofs_1d);
688 const unsigned int dofs_1d =
689 2 * this->n_dofs_per_vertex() + this->n_dofs_per_line();
692 const std::vector<unsigned int> &renumber =
695 for (
unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell; ++c)
697 this->prolongation[iso][c].reinit(this->n_dofs_per_cell(),
698 this->n_dofs_per_cell());
699 this->restriction[iso][c].reinit(this->n_dofs_per_cell(),
700 this->n_dofs_per_cell());
707 for (
unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
710 this->prolongation[iso][c].fill(dofs_subcell[c]);
711 this->restriction[iso][c].fill(dofs_cell[c]);
726 for (
unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
727 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
732 for (
unsigned int c = 0;
733 c < GeometryInfo<2>::max_children_per_cell;
737 const unsigned int c0 = c % 2;
739 const unsigned int c1 = c / 2;
741 this->prolongation[iso][c](j, i) =
742 dofs_subcell[c0](renumber[j] % dofs_1d,
743 renumber[i] % dofs_1d) *
744 dofs_subcell[c1]((renumber[j] - (renumber[j] % dofs_1d)) /
746 (renumber[i] - (renumber[i] % dofs_1d)) /
749 this->restriction[iso][c](j, i) =
750 dofs_cell[c0](renumber[j] % dofs_1d,
751 renumber[i] % dofs_1d) *
752 dofs_cell[c1]((renumber[j] - (renumber[j] % dofs_1d)) /
754 (renumber[i] - (renumber[i] % dofs_1d)) /
762 for (
unsigned int c = 0;
763 c < GeometryInfo<3>::max_children_per_cell;
767 const unsigned int c0 = c % 2;
769 const unsigned int c1 = (c % 4) / 2;
771 const unsigned int c2 = c / 4;
773 this->prolongation[iso][c](j, i) =
774 dofs_subcell[c0](renumber[j] % dofs_1d,
775 renumber[i] % dofs_1d) *
777 ((renumber[j] - (renumber[j] % dofs_1d)) / dofs_1d) %
779 ((renumber[i] - (renumber[i] % dofs_1d)) / dofs_1d) %
782 ((renumber[j] - (renumber[j] % dofs_1d)) / dofs_1d -
783 (((renumber[j] - (renumber[j] % dofs_1d)) / dofs_1d) %
786 ((renumber[i] - (renumber[i] % dofs_1d)) / dofs_1d -
787 (((renumber[i] - (renumber[i] % dofs_1d)) / dofs_1d) %
791 this->restriction[iso][c](j, i) =
792 dofs_cell[c0](renumber[j] % dofs_1d,
793 renumber[i] % dofs_1d) *
795 ((renumber[j] - (renumber[j] % dofs_1d)) / dofs_1d) %
797 ((renumber[i] - (renumber[i] % dofs_1d)) / dofs_1d) %
800 ((renumber[j] - (renumber[j] % dofs_1d)) / dofs_1d -
801 (((renumber[j] - (renumber[j] % dofs_1d)) / dofs_1d) %
804 ((renumber[i] - (renumber[i] % dofs_1d)) / dofs_1d -
805 (((renumber[i] - (renumber[i] % dofs_1d)) / dofs_1d) %
824 unsigned int n = this->degree + 1;
825 for (
unsigned int i = 1; i < dim; ++i)
826 n *= this->degree + 1;
828 this->generalized_support_points.resize(n);
832 const std::vector<unsigned int> &index_map_inverse =
864 for (
unsigned int iz = 0; iz <= ((dim > 2) ? this->degree : 0); ++iz)
865 for (
unsigned int iy = 0; iy <= ((dim > 1) ? this->degree : 0); ++iy)
866 for (
unsigned int ix = 0; ix <= this->degree; ++ix)
892 this->generalized_support_points[index_map_inverse[k++]] = p;
911 const unsigned int)
const
923 const unsigned int)
const
935 const unsigned int face_no)
const
942 (
dynamic_cast<const FEQHierarchical *
>(&x_source_fe) !=
946 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
948 this->n_dofs_per_face(face_no)));
972 interpolation_matrix = 0;
984 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
985 interpolation_matrix(i, i) = 1;
992 for (
unsigned int i = 0; i < GeometryInfo<3>::vertices_per_face; ++i)
993 interpolation_matrix(i, i) = 1;
995 for (
unsigned int i = 0; i < this->degree - 1; ++i)
997 for (
unsigned int j = 0; j < GeometryInfo<3>::lines_per_face; ++j)
998 interpolation_matrix(i + j * (x_source_fe.
degree - 1) +
1000 i + j * (this->degree - 1) +
1003 for (
unsigned int j = 0; j < this->degree - 1; ++j)
1005 (x_source_fe.
degree - 1) +
1008 (this->degree - 1) +
1022 const unsigned int subface,
1024 const unsigned int face_no)
const
1031 (
dynamic_cast<const FEQHierarchical *
>(&x_source_fe) !=
1035 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
1037 this->n_dofs_per_face(face_no)));
1069 interpolation_matrix(0, 0) = 1.0;
1070 interpolation_matrix(1, 0) = 0.5;
1071 interpolation_matrix(1, 1) = 0.5;
1073 for (
unsigned int dof = 2;
1074 dof < this->n_dofs_per_face(face_no);)
1076 interpolation_matrix(1, dof) = -1.0;
1080 int factorial_i = 1;
1082 for (
unsigned int i = 2; i < this->n_dofs_per_face(face_no);
1085 interpolation_matrix(i, i) =
std::pow(0.5, i);
1087 int factorial_j = factorial_i;
1088 int factorial_ij = 1;
1090 for (
unsigned int j = i + 1;
1091 j < this->n_dofs_per_face(face_no);
1094 factorial_ij *= j - i;
1097 if (((i + j) & 1) != 0u)
1098 interpolation_matrix(i, j) =
1099 -1.0 *
std::pow(0.5, j) * factorial_j /
1100 (factorial_i * factorial_ij);
1103 interpolation_matrix(i, j) =
1105 (factorial_i * factorial_ij);
1114 interpolation_matrix(0, 0) = 0.5;
1115 interpolation_matrix(0, 1) = 0.5;
1117 for (
unsigned int dof = 2;
1118 dof < this->n_dofs_per_face(face_no);)
1120 interpolation_matrix(0, dof) = -1.0;
1124 interpolation_matrix(1, 1) = 1.0;
1126 int factorial_i = 1;
1128 for (
unsigned int i = 2; i < this->n_dofs_per_face(face_no);
1131 interpolation_matrix(i, i) =
std::pow(0.5, i);
1133 int factorial_j = factorial_i;
1134 int factorial_ij = 1;
1136 for (
unsigned int j = i + 1;
1137 j < this->n_dofs_per_face(face_no);
1140 factorial_ij *= j - i;
1142 interpolation_matrix(i, j) =
1144 (factorial_i * factorial_ij);
1159 interpolation_matrix(0, 0) = 1.0;
1160 interpolation_matrix(1, 0) = 0.5;
1161 interpolation_matrix(1, 1) = 0.5;
1162 interpolation_matrix(2, 0) = 0.5;
1163 interpolation_matrix(2, 2) = 0.5;
1165 for (
unsigned int i = 0; i < this->degree - 1;)
1167 for (
unsigned int line = 0;
1168 line < GeometryInfo<3>::lines_per_face;
1170 interpolation_matrix(3,
1171 i + line * (this->degree - 1) +
1174 for (
unsigned int j = 0; j < this->degree - 1;)
1176 interpolation_matrix(3,
1177 i + (j + 4) * this->degree - j) =
1182 interpolation_matrix(1, i + 2 * (this->degree + 1)) =
1184 interpolation_matrix(2, i + 4) = -1.0;
1188 for (
unsigned int vertex = 0;
1189 vertex < GeometryInfo<3>::vertices_per_face;
1191 interpolation_matrix(3, vertex) = 0.25;
1193 int factorial_i = 1;
1195 for (
unsigned int i = 2; i <= this->degree; ++i)
1198 interpolation_matrix(i + 2, i + 2) = tmp;
1199 interpolation_matrix(i + 2 * source_fe.
degree,
1200 i + 2 * this->degree) = tmp;
1202 interpolation_matrix(i + source_fe.
degree + 1, i + 2) =
1204 interpolation_matrix(i + source_fe.
degree + 1,
1205 i + this->degree + 1) = tmp;
1206 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1207 i + 2 * this->degree) = tmp;
1208 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1209 i + 3 * this->degree - 1) = tmp;
1212 for (
unsigned int j = 0; j < this->degree - 1;)
1214 interpolation_matrix(i + source_fe.
degree + 1,
1215 (i + 2) * this->degree + j + 2 -
1217 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1218 i + (j + 4) * this->degree - j -
1223 int factorial_k = 1;
1225 for (
unsigned int j = 2; j <= this->degree; ++j)
1227 interpolation_matrix(i + (j + 2) * source_fe.
degree -
1229 i + (j + 2) * this->degree - j) =
1232 int factorial_kl = 1;
1233 int factorial_l = factorial_k;
1235 for (
unsigned int k = j + 1; k < this->degree; ++k)
1237 factorial_kl *= k - j;
1240 if (((j + k) & 1) != 0u)
1241 interpolation_matrix(
1242 i + (j + 2) * source_fe.
degree - j,
1243 i + (k + 2) * this->degree - k) =
1244 -1.0 *
std::pow(0.5, i + k) * factorial_l /
1245 (factorial_k * factorial_kl);
1248 interpolation_matrix(
1249 i + (j + 2) * source_fe.
degree - j,
1250 i + (k + 2) * this->degree - k) =
1251 std::pow(0.5, i + k) * factorial_l /
1252 (factorial_k * factorial_kl);
1257 int factorial_j = factorial_i;
1258 int factorial_ij = 1;
1260 for (
unsigned int j = i + 1; j <= this->degree; ++j)
1262 factorial_ij *= j - i;
1265 if (((i + j) & 1) != 0u)
1267 tmp = -1.0 *
std::pow(0.5, j) * factorial_j /
1268 (factorial_i * factorial_ij);
1269 interpolation_matrix(i + 2, j + 2) = tmp;
1270 interpolation_matrix(i + 2 * source_fe.
degree,
1271 j + 2 * this->degree) = tmp;
1274 for (
unsigned int k = 2; k <= this->degree; ++k)
1276 interpolation_matrix(
1277 i + (k + 2) * source_fe.
degree - k,
1278 j + (k + 2) * this->degree - k) =
1281 int factorial_l = factorial_k;
1282 int factorial_kl = 1;
1284 for (
unsigned int l = k + 1;
1288 factorial_kl *= l - k;
1291 if (((k + l) & 1) != 0u)
1292 interpolation_matrix(
1293 i + (k + 2) * source_fe.
degree - k,
1294 j + (l + 2) * this->degree - l) =
1297 (factorial_k * factorial_kl);
1300 interpolation_matrix(
1301 i + (k + 2) * source_fe.
degree - k,
1302 j + (l + 2) * this->degree - l) =
1303 tmp *
std::pow(0.5, l) * factorial_l /
1304 (factorial_k * factorial_kl);
1309 interpolation_matrix(i + source_fe.
degree + 1,
1311 interpolation_matrix(i + source_fe.
degree + 1,
1312 j + this->degree + 1) = tmp;
1313 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1314 j + 2 * this->degree) = tmp;
1315 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1316 j + 3 * this->degree - 1) =
1320 for (
unsigned int k = 0; k < this->degree - 1;)
1322 interpolation_matrix(i + source_fe.
degree + 1,
1323 (j + 2) * this->degree +
1325 interpolation_matrix(
1326 i + 3 * source_fe.
degree - 1,
1327 j + (k + 4) * this->degree - k - 2) = tmp;
1333 tmp =
std::pow(0.5, j) * factorial_j /
1334 (factorial_i * factorial_ij);
1335 interpolation_matrix(i + 2, j + 2) = tmp;
1336 interpolation_matrix(i + 2 * source_fe.
degree,
1337 j + 2 * this->degree) = tmp;
1340 for (
unsigned int k = 2; k <= this->degree; ++k)
1342 interpolation_matrix(
1343 i + (k + 2) * source_fe.
degree - k,
1344 j + (k + 2) * this->degree - k) =
1347 int factorial_l = factorial_k;
1348 int factorial_kl = 1;
1350 for (
unsigned int l = k + 1;
1354 factorial_kl *= l - k;
1357 if (((k + l) & 1) != 0u)
1358 interpolation_matrix(
1359 i + (k + 2) * source_fe.
degree - k,
1360 j + (l + 2) * this->degree - l) =
1363 (factorial_k * factorial_kl);
1366 interpolation_matrix(
1367 i + (k + 2) * source_fe.
degree - k,
1368 j + (l + 2) * this->degree - l) =
1369 tmp *
std::pow(0.5, l) * factorial_l /
1370 (factorial_k * factorial_kl);
1375 interpolation_matrix(i + source_fe.
degree + 1,
1377 interpolation_matrix(i + source_fe.
degree + 1,
1378 j + this->degree + 1) = tmp;
1379 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1380 j + 2 * this->degree) = tmp;
1381 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1382 j + 3 * this->degree - 1) =
1386 for (
unsigned int k = 0; k < this->degree - 1;)
1388 interpolation_matrix(i + source_fe.
degree + 1,
1389 (j + 2) * this->degree +
1391 interpolation_matrix(
1392 i + 3 * source_fe.
degree - 1,
1393 j + (k + 4) * this->degree - k - 2) = tmp;
1405 interpolation_matrix(0, 0) = 0.5;
1406 interpolation_matrix(0, 1) = 0.5;
1407 interpolation_matrix(1, 1) = 1.0;
1408 interpolation_matrix(3, 1) = 0.5;
1409 interpolation_matrix(3, 3) = 0.5;
1411 for (
unsigned int i = 0; i < this->degree - 1;)
1413 for (
unsigned int line = 0;
1414 line < GeometryInfo<3>::lines_per_face;
1416 interpolation_matrix(2,
1417 i + line * (this->degree - 1) +
1420 for (
unsigned int j = 0; j < this->degree - 1;)
1422 interpolation_matrix(2,
1423 i + (j + 4) * this->degree - j) =
1428 interpolation_matrix(0, i + 2 * (this->degree + 1)) =
1430 interpolation_matrix(3, i + this->degree + 3) = -1.0;
1434 for (
unsigned int vertex = 0;
1435 vertex < GeometryInfo<3>::vertices_per_face;
1437 interpolation_matrix(2, vertex) = 0.25;
1439 int factorial_i = 1;
1441 for (
unsigned int i = 2; i <= this->degree; ++i)
1444 interpolation_matrix(i + 2, i + 2) = tmp;
1445 interpolation_matrix(i + 2, i + this->degree + 1) = tmp;
1446 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1447 i + 2 * this->degree) = tmp;
1448 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1449 i + 3 * this->degree - 1) = tmp;
1452 for (
unsigned int j = 0; j < this->degree - 1;)
1454 interpolation_matrix(i + 2,
1455 j + (i + 2) * this->degree + 2 -
1457 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1458 i + (j + 4) * this->degree - j -
1464 interpolation_matrix(i + source_fe.
degree + 1,
1465 i + this->degree + 1) = tmp;
1466 interpolation_matrix(i + 2 * source_fe.
degree,
1467 i + 2 * this->degree) = tmp;
1469 int factorial_j = factorial_i;
1470 int factorial_ij = 1;
1472 for (
unsigned int j = i + 1; j <= this->degree; ++j)
1474 factorial_ij *= j - i;
1476 tmp =
std::pow(0.5, j) * factorial_j /
1477 (factorial_i * factorial_ij);
1478 interpolation_matrix(i + 2 * source_fe.
degree,
1479 j + 2 * this->degree) = tmp;
1480 int factorial_k = 1;
1482 for (
unsigned int k = 2; k <= this->degree; ++k)
1484 interpolation_matrix(
1485 i + (k + 2) * source_fe.
degree - k,
1486 j + (k + 2) * this->degree - k) =
1489 int factorial_l = factorial_k;
1490 int factorial_kl = 1;
1492 for (
unsigned int l = k + 1; l <= this->degree;
1495 factorial_kl *= l - k;
1498 if (((k + l) & 1) != 0u)
1499 interpolation_matrix(
1500 i + (k + 2) * source_fe.
degree - k,
1501 j + (l + 2) * this->degree - l) =
1504 (factorial_k * factorial_kl);
1507 interpolation_matrix(
1508 i + (k + 2) * source_fe.
degree - k,
1509 j + (l + 2) * this->degree - l) =
1510 tmp *
std::pow(0.5, l) * factorial_l /
1511 (factorial_k * factorial_kl);
1517 for (
unsigned int k = 0; k < this->degree - 1;)
1519 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1520 j + (k + 4) * this->degree -
1526 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1527 j + 2 * this->degree) = tmp;
1528 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1529 j + 3 * this->degree - 1) = tmp;
1531 if (((i + j) & 1) != 0u)
1534 interpolation_matrix(i + 2, j + 2) = tmp;
1535 interpolation_matrix(i + 2, j + this->degree + 1) =
1537 interpolation_matrix(i + source_fe.
degree + 1,
1538 j + this->degree + 1) =
1542 for (
unsigned int k = 0; k < this->degree - 1;)
1544 interpolation_matrix(i + 2,
1545 k + (j + 2) * this->degree +
1551 int factorial_k = 1;
1553 for (
unsigned int j = 2; j <= this->degree; ++j)
1555 interpolation_matrix(i + (j + 2) * source_fe.
degree -
1557 i + (j + 2) * this->degree - j) =
1560 int factorial_l = factorial_k;
1561 int factorial_kl = 1;
1563 for (
unsigned int k = j + 1; k <= this->degree; ++k)
1565 factorial_kl *= k - j;
1568 if (((j + k) & 1) != 0u)
1569 interpolation_matrix(
1570 i + (j + 2) * source_fe.
degree - j,
1571 i + (k + 2) * this->degree - k) =
1572 -1.0 *
std::pow(0.5, i + k) * factorial_l /
1573 (factorial_k * factorial_kl);
1576 interpolation_matrix(
1577 i + (j + 2) * source_fe.
degree - j,
1578 i + (k + 2) * this->degree - k) =
1579 std::pow(0.5, i + k) * factorial_l /
1580 (factorial_k * factorial_kl);
1590 interpolation_matrix(0, 0) = 0.5;
1591 interpolation_matrix(0, 2) = 0.5;
1592 interpolation_matrix(2, 2) = 1.0;
1593 interpolation_matrix(3, 2) = 0.5;
1594 interpolation_matrix(3, 3) = 0.5;
1596 for (
unsigned int i = 0; i < this->degree - 1;)
1598 for (
unsigned int line = 0;
1599 line < GeometryInfo<3>::lines_per_face;
1601 interpolation_matrix(1,
1602 i + line * (this->degree - 1) +
1605 for (
unsigned int j = 0; j < this->degree - 1;)
1607 interpolation_matrix(1,
1608 i + (j + 4) * this->degree - j) =
1613 interpolation_matrix(0, i + 4) = -1.0;
1614 interpolation_matrix(3, i + 3 * this->degree + 1) = -1.0;
1618 for (
unsigned int vertex = 0;
1619 vertex < GeometryInfo<3>::vertices_per_face;
1621 interpolation_matrix(1, vertex) = 0.25;
1623 int factorial_i = 1;
1625 for (
unsigned int i = 2; i <= this->degree; ++i)
1628 interpolation_matrix(i + 2, i + 2) = tmp;
1629 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1630 i + 3 * this->degree - 1) = tmp;
1632 interpolation_matrix(i + source_fe.
degree + 1, i + 2) =
1634 interpolation_matrix(i + source_fe.
degree + 1,
1635 i + this->degree + 1) = tmp;
1636 interpolation_matrix(i + 2 * source_fe.
degree,
1637 i + 2 * this->degree) = tmp;
1638 interpolation_matrix(i + 2 * source_fe.
degree,
1639 i + 3 * this->degree - 1) = tmp;
1642 for (
unsigned int j = 0; j < this->degree - 1;)
1644 interpolation_matrix(i + source_fe.
degree + 1,
1645 j + (i + 2) * this->degree + 2 -
1647 interpolation_matrix(i + 2 * source_fe.
degree,
1648 i + (j + 4) * this->degree - j -
1653 int factorial_k = 1;
1655 for (
unsigned int j = 2; j <= this->degree; ++j)
1657 interpolation_matrix(i + (j + 2) * source_fe.
degree -
1659 i + (j + 2) * this->degree - j) =
1662 int factorial_kl = 1;
1663 int factorial_l = factorial_k;
1665 for (
unsigned int k = j + 1; k <= this->degree; ++k)
1667 factorial_kl *= k - j;
1669 interpolation_matrix(
1670 i + (j + 2) * source_fe.
degree - j,
1671 i + (k + 2) * this->degree - k) =
1672 std::pow(0.5, i + k) * factorial_l /
1673 (factorial_k * factorial_kl);
1678 int factorial_j = factorial_i;
1679 int factorial_ij = 1;
1681 for (
unsigned int j = i + 1; j <= this->degree; ++j)
1683 factorial_ij *= j - i;
1685 tmp =
std::pow(0.5, j) * factorial_j /
1686 (factorial_i * factorial_ij);
1687 interpolation_matrix(i + 2, j + 2) = tmp;
1690 for (
unsigned int k = 0; k < this->degree - 1;)
1692 interpolation_matrix(i + source_fe.
degree + 1,
1693 k + (j + 2) * this->degree +
1699 interpolation_matrix(i + source_fe.
degree + 1,
1701 interpolation_matrix(i + source_fe.
degree + 1,
1702 j + this->degree + 1) = tmp;
1704 if (((i + j) & 1) != 0u)
1707 interpolation_matrix(i + 2 * source_fe.
degree,
1708 j + 2 * this->degree) = tmp;
1709 interpolation_matrix(i + 2 * source_fe.
degree,
1710 j + 3 * this->degree - 1) = tmp;
1712 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1713 j + 3 * this->degree - 1) = tmp;
1714 int factorial_k = 1;
1716 for (
unsigned int k = 2; k <= this->degree; ++k)
1718 interpolation_matrix(
1719 i + (k + 2) * source_fe.
degree - k,
1720 j + (k + 2) * this->degree - k) =
1723 int factorial_l = factorial_k;
1724 int factorial_kl = 1;
1726 for (
unsigned int l = k + 1; l <= this->degree;
1729 factorial_kl *= l - k;
1731 interpolation_matrix(
1732 i + (k + 2) * source_fe.
degree - k,
1733 j + (l + 2) * this->degree - l) =
1734 tmp *
std::pow(0.5, l) * factorial_l /
1735 (factorial_k * factorial_kl);
1741 for (
unsigned int k = 0; k < this->degree - 1;)
1743 interpolation_matrix(i + 2 * source_fe.
degree,
1744 j + (k + 4) * this->degree -
1756 for (
unsigned int vertex = 0;
1757 vertex < GeometryInfo<3>::vertices_per_face;
1759 interpolation_matrix(0, vertex) = 0.25;
1761 for (
unsigned int i = 0; i < this->degree - 1;)
1763 for (
unsigned int line = 0;
1764 line < GeometryInfo<3>::lines_per_face;
1766 interpolation_matrix(0,
1767 i + line * (this->degree - 1) +
1770 for (
unsigned int j = 0; j < this->degree - 1;)
1772 interpolation_matrix(0,
1773 i + (j + 4) * this->degree - j) =
1778 interpolation_matrix(1, i + 4) = -1.0;
1779 interpolation_matrix(2, i + 3 * this->degree + 1) = -1.0;
1783 interpolation_matrix(1, 0) = 0.5;
1784 interpolation_matrix(1, 1) = 0.5;
1785 interpolation_matrix(2, 2) = 0.5;
1786 interpolation_matrix(2, 3) = 0.5;
1787 interpolation_matrix(3, 3) = 1.0;
1789 int factorial_i = 1;
1791 for (
unsigned int i = 2; i <= this->degree; ++i)
1794 interpolation_matrix(i + 2, i + 2) = tmp;
1795 interpolation_matrix(i + 2, i + this->degree + 1) = tmp;
1796 interpolation_matrix(i + 2 * source_fe.
degree,
1797 i + 2 * this->degree) = tmp;
1798 interpolation_matrix(i + 2 * source_fe.
degree,
1799 i + 3 * this->degree - 1) = tmp;
1802 for (
unsigned int j = 0; j < this->degree - 1;)
1804 interpolation_matrix(i + 2,
1805 j + (i + 2) * this->degree + 2 -
1807 interpolation_matrix(i + 2 * source_fe.
degree,
1808 i + (j + 4) * this->degree - 2) =
1814 interpolation_matrix(i + source_fe.
degree + 1,
1815 i + this->degree + 1) = tmp;
1816 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1817 i + 3 * this->degree - 1) = tmp;
1818 int factorial_k = 1;
1820 for (
unsigned int j = 2; j <= this->degree; ++j)
1822 interpolation_matrix(i + (j + 2) * source_fe.
degree -
1824 i + (j + 2) * this->degree - j) =
1827 int factorial_l = factorial_k;
1828 int factorial_kl = 1;
1830 for (
unsigned int k = j + 1; k <= this->degree; ++k)
1832 factorial_kl *= k - j;
1834 interpolation_matrix(
1835 i + (j + 2) * source_fe.
degree - j,
1836 i + (k + 2) * this->degree - k) =
1837 std::pow(0.5, i + k) * factorial_l /
1838 (factorial_k * factorial_kl);
1843 int factorial_j = factorial_i;
1844 int factorial_ij = 1;
1846 for (
unsigned int j = i + 1; j <= this->degree; ++j)
1848 factorial_ij *= j - i;
1850 tmp =
std::pow(0.5, j + 1) * factorial_j /
1851 (factorial_i * factorial_ij);
1852 interpolation_matrix(i + 2, j + 2) = tmp;
1853 interpolation_matrix(i + 2, j + this->degree + 1) =
1855 interpolation_matrix(i + 2 * source_fe.
degree,
1856 j + 2 * this->degree) = tmp;
1857 interpolation_matrix(i + 2 * source_fe.
degree,
1858 j + 3 * this->degree - 1) = tmp;
1860 interpolation_matrix(i + source_fe.
degree + 1,
1861 j + this->degree + 1) = tmp;
1862 interpolation_matrix(i + 3 * source_fe.
degree - 1,
1863 j + 3 * this->degree - 1) = tmp;
1864 int factorial_k = 1;
1866 for (
unsigned int k = 2; k <= this->degree; ++k)
1868 interpolation_matrix(
1869 i + (k + 2) * source_fe.
degree - k,
1870 j + (k + 2) * this->degree - k) =
1873 int factorial_l = factorial_k;
1874 int factorial_kl = 1;
1876 for (
unsigned int l = k + 1; l <= this->degree;
1879 factorial_kl *= l - k;
1881 interpolation_matrix(
1882 i + (k + 2) * source_fe.
degree - k,
1883 j + (l + 2) * this->degree - l) =
1884 tmp *
std::pow(0.5, l) * factorial_l /
1885 (factorial_k * factorial_kl);
1891 for (
unsigned int k = 0; k < this->degree - 1;)
1893 interpolation_matrix(i + 2,
1894 k + (j + 2) * this->degree +
1896 interpolation_matrix(i + 2 * source_fe.
degree,
1897 j + (k + 4) * this->degree -
1915 const unsigned int codim = dim - 1;
1920 const unsigned int face_no = 0;
1923 unsigned int n = this->degree + 1;
1924 for (
unsigned int i = 1; i < codim; ++i)
1925 n *= this->degree + 1;
1927 this->generalized_face_support_points[face_no].resize(n);
1932 for (
unsigned int iz = 0; iz <= ((codim > 2) ? this->degree : 0); ++iz)
1933 for (
unsigned int iy = 0; iy <= ((codim > 1) ? this->degree : 0); ++iy)
1934 for (
unsigned int ix = 0; ix <= this->degree; ++ix)
1960 this->generalized_face_support_points[face_no][face_renumber[k++]] =
1968std::vector<unsigned int>
1971 std::vector<unsigned int> dpo(dim + 1, 1U);
1972 for (
unsigned int i = 1; i < dpo.size(); ++i)
1973 dpo[i] = dpo[i - 1] * (deg - 1);
1980std::vector<unsigned int>
1991 const unsigned int n = degree + 1;
2030 unsigned int next_index = 0;
2032 h2l[next_index++] = 0;
2033 h2l[next_index++] = 1;
2034 h2l[next_index++] = n;
2035 h2l[next_index++] = n + 1;
2038 h2l[next_index++] = (2 + i) * n;
2041 h2l[next_index++] = (2 + i) * n + 1;
2044 h2l[next_index++] = 2 + i;
2047 h2l[next_index++] = n + 2 + i;
2054 h2l[next_index++] = (2 + i) * n + 2 + j;
2063 unsigned int next_index = 0;
2064 const unsigned int n2 = n * n;
2067 h2l[next_index++] = 0;
2068 h2l[next_index++] = 1;
2069 h2l[next_index++] = n;
2070 h2l[next_index++] = n + 1;
2072 h2l[next_index++] = n2;
2073 h2l[next_index++] = n2 + 1;
2074 h2l[next_index++] = n2 + n;
2075 h2l[next_index++] = n2 + n + 1;
2080 h2l[next_index++] = (2 + i) * n;
2082 h2l[next_index++] = (2 + i) * n + 1;
2084 h2l[next_index++] = 2 + i;
2086 h2l[next_index++] = n + 2 + i;
2089 h2l[next_index++] = n2 + (2 + i) * n;
2091 h2l[next_index++] = n2 + (2 + i) * n + 1;
2093 h2l[next_index++] = n2 + 2 + i;
2095 h2l[next_index++] = n2 + n + 2 + i;
2098 h2l[next_index++] = (2 + i) * n2;
2100 h2l[next_index++] = (2 + i) * n2 + 1;
2102 h2l[next_index++] = (2 + i) * n2 + n;
2104 h2l[next_index++] = (2 + i) * n2 + n + 1;
2109 const unsigned int face_no = 0;
2119 h2l[next_index++] = (2 + i) * n2 + (2 + j) * n;
2123 h2l[next_index++] = (2 + i) * n2 + (2 + j) * n + 1;
2127 h2l[next_index++] = (2 + i) * n2 + 2 + j;
2131 h2l[next_index++] = (2 + i) * n2 + n + 2 + j;
2135 h2l[next_index++] = (2 + i) * n + 2 + j;
2139 h2l[next_index++] = n2 + (2 + i) * n + 2 + j;
2148 h2l[next_index++] = (2 + i) * n2 + (2 + j) * n + 2 + k;
2163std::vector<unsigned int>
2165 const unsigned int degree)
2169 return internal::FE_Q_Hierarchical::invert_numbering(
2177std::vector<unsigned int>
2181 return std::vector<unsigned int>();
2188 const unsigned int face_index)
const
2200 return (((shape_index == 0) && (face_index == 0)) ||
2201 ((shape_index == 1) && (face_index == 1)));
2209 const unsigned int face_index)
const
2218 if (((dim == 2) && (shape_index >=
2219 this->get_first_quad_index(0 ))) ||
2220 ((dim == 3) && (shape_index >= this->get_first_hex_index())))
2225 if (shape_index < this->get_first_line_index())
2232 const unsigned int vertex_no = shape_index;
2235 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_face; ++i)
2241 else if (shape_index < this->get_first_quad_index(0))
2244 const unsigned int line_index =
2245 (shape_index - this->get_first_line_index()) / this->n_dofs_per_line();
2249 for (
unsigned int i = 0; i < GeometryInfo<dim>::lines_per_face; ++i)
2254 else if (shape_index < this->get_first_hex_index())
2257 const unsigned int quad_index =
2258 (shape_index - this->get_first_quad_index(0 )) /
2259 this->n_dofs_per_quad(face_index);
2260 Assert(
static_cast<signed int>(quad_index) <
2273 return (quad_index == face_index);
2295std::vector<unsigned int>
2298 Assert((sub_degree > 0) && (sub_degree <= this->degree),
2303 std::vector<unsigned int> embedding_dofs(sub_degree + 1);
2304 for (
unsigned int i = 0; i < (sub_degree + 1); ++i)
2305 embedding_dofs[i] = i;
2307 return embedding_dofs;
2310 if (sub_degree == 1)
2312 std::vector<unsigned int> embedding_dofs(
2315 embedding_dofs[i] = i;
2317 return embedding_dofs;
2319 else if (sub_degree == this->degree)
2321 std::vector<unsigned int> embedding_dofs(this->n_dofs_per_cell());
2322 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
2323 embedding_dofs[i] = i;
2325 return embedding_dofs;
2328 if ((dim == 2) || (dim == 3))
2330 std::vector<unsigned int> embedding_dofs(
2331 (dim == 2) ? (sub_degree + 1) * (sub_degree + 1) :
2332 (sub_degree + 1) * (sub_degree + 1) * (sub_degree + 1));
2334 for (
unsigned int i = 0;
2336 (sub_degree + 1) * (sub_degree + 1) :
2337 (sub_degree + 1) * (sub_degree + 1) * (sub_degree + 1));
2342 embedding_dofs[i] = i;
2347 const unsigned int j =
2349 const unsigned int line =
2354 line * (this->degree - 1) + j;
2362 const unsigned int j =
2366 const unsigned int k =
2371 const unsigned int face =
2374 k * (sub_degree - 1) - j) /
2375 ((sub_degree - 1) * (sub_degree - 1));
2380 face * (this->degree - 1) * (this->degree - 1) +
2381 k * (this->degree - 1) + j;
2389 (sub_degree - 1) * (sub_degree - 1))
2391 const unsigned int j =
2397 const unsigned int k =
2405 const unsigned int l =
2410 j - k * (sub_degree - 1)) /
2411 ((sub_degree - 1) * (sub_degree - 1));
2417 (this->degree - 1) +
2418 l * (this->degree - 1) * (this->degree - 1) +
2419 k * (this->degree - 1) + j;
2423 return embedding_dofs;
2428 return std::vector<unsigned int>();
2435std::pair<Table<2, bool>, std::vector<unsigned int>>
2440 constant_modes(0, i) =
true;
2442 i < this->n_dofs_per_cell();
2444 constant_modes(0, i) =
false;
2445 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
2446 constant_modes, std::vector<unsigned int>(1, 0));
2462#include "fe_q_hierarchical.inst"
const std::unique_ptr< ScalarPolynomialsBase< dim > > poly_space
void build_dofs_cell(std::vector< FullMatrix< double > > &dofs_cell, std::vector< FullMatrix< double > > &dofs_subcell) const
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual void get_interpolation_matrix(const FiniteElement< dim > &source, FullMatrix< double > &matrix) const override
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
void initialize_embedding_and_restriction(const std::vector< FullMatrix< double > > &dofs_cell, const std::vector< FullMatrix< double > > &dofs_subcell)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
friend class FE_Q_Hierarchical
std::vector< unsigned int > get_embedding_dofs(const unsigned int sub_degree) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim > &fe_other, const unsigned int face_no=0) const override
void initialize_generalized_support_points()
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim > &fe_other, const unsigned int codim=0) const override final
static std::vector< unsigned int > hierarchic_to_fe_q_hierarchical_numbering(const FiniteElementData< dim > &fe)
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
virtual void get_subface_interpolation_matrix(const FiniteElement< dim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim > &fe_other) const override
static std::vector< unsigned int > face_fe_q_hierarchical_to_hierarchic_numbering(const unsigned int degree)
void initialize_generalized_face_support_points()
void initialize_constraints(const std::vector< FullMatrix< double > > &dofs_subcell)
virtual std::string get_name() const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual std::size_t memory_consumption() const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual bool hp_constraints_are_implemented() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim > &fe_other) const override
unsigned int n_dofs_per_vertex() const
const unsigned int degree
unsigned int n_dofs_per_cell() const
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_components() const
unsigned int n_unique_faces() const
unsigned int n_dofs_per_quad(unsigned int face_no=0) const
unsigned int n_dofs_per_hex() const
virtual std::string get_name() const =0
const std::vector< unsigned int > & get_numbering() const
void set_numbering(const std::vector< unsigned int > &renumber)
const std::vector< unsigned int > & get_numbering_inverse() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
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 & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ either_element_can_dominate
@ other_element_dominates
@ neither_element_dominates
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices()