40 count_nonzeros(
const std::vector<unsigned int> &vec)
42 return std::count_if(vec.begin(), vec.end(), [](
const unsigned int i) {
55 template <
int dim,
int spacedim = dim>
58 const unsigned int base_no)
66 unsigned int out_index = 0;
67 for (
unsigned int system_index = 0; system_index < fe.
n_dofs_per_cell();
74 const unsigned int base_component =
76 const unsigned int base_index =
81 table[base_component][base_index] = out_index;
92 template <
int dim,
int spacedim = dim>
93 std::vector<typename FESystem<dim, spacedim>::BaseOffsets>
95 const unsigned int base_no)
97 std::vector<typename FESystem<dim, spacedim>::BaseOffsets> table;
100 unsigned int out_index = 0;
101 for (
unsigned int system_index = 0; system_index < fe.
n_dofs_per_cell();
106 const unsigned int base_index =
109 table.emplace_back();
111 table.back().n_nonzero_components =
113 unsigned int in_index = 0;
114 for (
unsigned int i = 0; i < base_index; ++i)
117 table.back().in_index = in_index;
118 table.back().out_index = out_index;
133 template <
int dim,
int spacedim = dim>
137 const unsigned int base_no,
147 const unsigned int n_dofs_per_cell =
150 auto copy_row = [](
const auto row_in,
auto row_out) {
151 std::copy(row_in.begin(), row_in.end(), row_out.begin());
155 for (
unsigned int component = 0; component < n_components; ++component)
156 for (
unsigned int b = 0; b < n_dofs_per_cell; ++b)
159 output_data.
shape_values[base_to_system_table[component][b]]);
162 for (
unsigned int component = 0; component < n_components; ++component)
163 for (
unsigned int b = 0; b < n_dofs_per_cell; ++b)
169 for (
unsigned int component = 0; component < n_components; ++component)
170 for (
unsigned int b = 0; b < n_dofs_per_cell; ++b)
176 for (
unsigned int component = 0; component < n_components; ++component)
177 for (
unsigned int b = 0; b < n_dofs_per_cell; ++b)
188 template <
int dim,
int spacedim = dim>
192 const unsigned int base_no,
193 const unsigned int n_q_points,
205 for (
const auto &offset : offsets)
208 for (
unsigned int s = 0; s < offset.n_nonzero_components; ++s)
209 for (
unsigned int q = 0; q < n_q_points; ++q)
214 for (
unsigned int s = 0; s < offset.n_nonzero_components; ++s)
215 for (
unsigned int q = 0; q < n_q_points; ++q)
220 for (
unsigned int s = 0; s < offset.n_nonzero_components; ++s)
221 for (
unsigned int q = 0; q < n_q_points; ++q)
226 for (
unsigned int s = 0; s < offset.n_nonzero_components; ++s)
227 for (
unsigned int q = 0; q < n_q_points; ++q)
236template <
int dim,
int spacedim>
238 const unsigned int n_base_elements)
239 : base_fe_datas(n_base_elements)
240 , base_fe_output_objects(n_base_elements)
245template <
int dim,
int spacedim>
254template <
int dim,
int spacedim>
257 const unsigned int base_no)
const
260 return *base_fe_datas[base_no];
265template <
int dim,
int spacedim>
268 const unsigned int base_no,
272 base_fe_datas[base_no] = std::move(ptr);
277template <
int dim,
int spacedim>
280 const unsigned int base_no)
const
283 return base_fe_output_objects[base_no];
291template <
int dim,
int spacedim>
295template <
int dim,
int spacedim>
297 const unsigned int n_elements)
305 std::vector<const FiniteElement<dim, spacedim> *> fes;
307 std::vector<unsigned int> multiplicities;
308 multiplicities.push_back(n_elements);
309 initialize(fes, multiplicities);
314template <
int dim,
int spacedim>
316 const unsigned int n1,
318 const unsigned int n2)
326 , base_elements(static_cast<
int>(n1 > 0) + static_cast<
int>(n2 > 0))
328 std::vector<const FiniteElement<dim, spacedim> *> fes;
331 std::vector<unsigned int> multiplicities;
332 multiplicities.push_back(n1);
333 multiplicities.push_back(n2);
334 initialize(fes, multiplicities);
339template <
int dim,
int spacedim>
341 const unsigned int n1,
343 const unsigned int n2,
345 const unsigned int n3)
360 , base_elements(static_cast<
int>(n1 > 0) + static_cast<
int>(n2 > 0) +
361 static_cast<
int>(n3 > 0))
363 std::vector<const FiniteElement<dim, spacedim> *> fes;
367 std::vector<unsigned int> multiplicities;
368 multiplicities.push_back(n1);
369 multiplicities.push_back(n2);
370 multiplicities.push_back(n3);
371 initialize(fes, multiplicities);
376template <
int dim,
int spacedim>
378 const unsigned int n1,
380 const unsigned int n2,
382 const unsigned int n3,
384 const unsigned int n4)
410 , base_elements(static_cast<
int>(n1 > 0) + static_cast<
int>(n2 > 0) +
411 static_cast<
int>(n3 > 0) + static_cast<
int>(n4 > 0))
413 std::vector<const FiniteElement<dim, spacedim> *> fes;
418 std::vector<unsigned int> multiplicities;
419 multiplicities.push_back(n1);
420 multiplicities.push_back(n2);
421 multiplicities.push_back(n3);
422 multiplicities.push_back(n4);
423 initialize(fes, multiplicities);
428template <
int dim,
int spacedim>
430 const unsigned int n1,
432 const unsigned int n2,
434 const unsigned int n3,
436 const unsigned int n4,
438 const unsigned int n5)
441 multiply_dof_numbers(&fe1, n1, &fe2, n2, &fe3, n3, &fe4, n4, &fe5, n5),
462 , base_elements(static_cast<
int>(n1 > 0) + static_cast<
int>(n2 > 0) +
463 static_cast<
int>(n3 > 0) + static_cast<
int>(n4 > 0) +
464 static_cast<
int>(n5 > 0))
466 std::vector<const FiniteElement<dim, spacedim> *> fes;
472 std::vector<unsigned int> multiplicities;
473 multiplicities.push_back(n1);
474 multiplicities.push_back(n2);
475 multiplicities.push_back(n3);
476 multiplicities.push_back(n4);
477 multiplicities.push_back(n5);
478 initialize(fes, multiplicities);
483template <
int dim,
int spacedim>
486 const std::vector<unsigned int> & multiplicities)
493 , base_elements(count_nonzeros(multiplicities))
495 initialize(fes, multiplicities);
500template <
int dim,
int spacedim>
511 std::ostringstream namebuf;
514 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
516 namebuf << base_element(i).get_name();
517 if (this->element_multiplicity(i) != 1)
518 namebuf <<
'^' << this->element_multiplicity(i);
519 if (i != this->n_base_elements() - 1)
524 return namebuf.str();
529template <
int dim,
int spacedim>
530std::unique_ptr<FiniteElement<dim, spacedim>>
533 std::vector<const FiniteElement<dim, spacedim> *> fes;
534 std::vector<unsigned int> multiplicities;
536 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
538 fes.push_back(&base_element(i));
539 multiplicities.push_back(this->element_multiplicity(i));
541 return std::make_unique<FESystem<dim, spacedim>>(fes, multiplicities);
546template <
int dim,
int spacedim>
549 const unsigned int first_component,
550 const unsigned int n_selected_components)
const
552 Assert(first_component + n_selected_components <= this->n_components(),
553 ExcMessage(
"Invalid arguments (not a part of this FiniteElement)."));
555 const unsigned int base_index =
556 this->component_to_base_table[first_component].first.first;
557 const unsigned int component_in_base =
558 this->component_to_base_table[first_component].first.second;
559 const unsigned int base_components =
560 this->base_element(base_index).n_components();
564 if (n_selected_components <= base_components)
565 return this->base_element(base_index)
566 .get_sub_fe(component_in_base, n_selected_components);
568 Assert(n_selected_components == this->n_components(),
569 ExcMessage(
"You can not select a part of a FiniteElement."));
575template <
int dim,
int spacedim>
581 Assert(this->is_primitive(i),
585 return (base_element(this->system_to_base_table[i].
first.first)
586 .shape_value(this->system_to_base_table[i].second, p));
591template <
int dim,
int spacedim>
594 const unsigned int i,
596 const unsigned int component)
const
603 if (this->nonzero_components[i][component] ==
false)
611 const unsigned int base = this->component_to_base_index(component).first;
612 const unsigned int component_in_base =
613 this->component_to_base_index(component).second;
621 return (base_element(base).shape_value_component(
622 this->system_to_base_table[i].
second, p, component_in_base));
627template <
int dim,
int spacedim>
633 Assert(this->is_primitive(i),
637 return (base_element(this->system_to_base_table[i].
first.first)
638 .shape_grad(this->system_to_base_table[i].second, p));
643template <
int dim,
int spacedim>
646 const unsigned int i,
648 const unsigned int component)
const
654 if (this->nonzero_components[i][component] ==
false)
659 const unsigned int base = this->component_to_base_index(component).first;
660 const unsigned int component_in_base =
661 this->component_to_base_index(component).second;
666 return (base_element(base).shape_grad_component(
667 this->system_to_base_table[i].
second, p, component_in_base));
672template <
int dim,
int spacedim>
678 Assert(this->is_primitive(i),
682 return (base_element(this->system_to_base_table[i].
first.first)
683 .shape_grad_grad(this->system_to_base_table[i].second, p));
688template <
int dim,
int spacedim>
691 const unsigned int i,
693 const unsigned int component)
const
699 if (this->nonzero_components[i][component] ==
false)
704 const unsigned int base = this->component_to_base_index(component).first;
705 const unsigned int component_in_base =
706 this->component_to_base_index(component).second;
711 return (base_element(base).shape_grad_grad_component(
712 this->system_to_base_table[i].
second, p, component_in_base));
717template <
int dim,
int spacedim>
723 Assert(this->is_primitive(i),
727 return (base_element(this->system_to_base_table[i].
first.first)
728 .shape_3rd_derivative(this->system_to_base_table[i].second, p));
733template <
int dim,
int spacedim>
736 const unsigned int i,
738 const unsigned int component)
const
744 if (this->nonzero_components[i][component] ==
false)
749 const unsigned int base = this->component_to_base_index(component).first;
750 const unsigned int component_in_base =
751 this->component_to_base_index(component).second;
756 return (base_element(base).shape_3rd_derivative_component(
757 this->system_to_base_table[i].
second, p, component_in_base));
762template <
int dim,
int spacedim>
768 Assert(this->is_primitive(i),
772 return (base_element(this->system_to_base_table[i].
first.first)
773 .shape_4th_derivative(this->system_to_base_table[i].second, p));
778template <
int dim,
int spacedim>
781 const unsigned int i,
783 const unsigned int component)
const
789 if (this->nonzero_components[i][component] ==
false)
794 const unsigned int base = this->component_to_base_index(component).first;
795 const unsigned int component_in_base =
796 this->component_to_base_index(component).second;
801 return (base_element(base).shape_4th_derivative_component(
802 this->system_to_base_table[i].
second, p, component_in_base));
807template <
int dim,
int spacedim>
817 Assert((interpolation_matrix.
m() == this->n_dofs_per_cell()) ||
820 this->n_dofs_per_cell()));
822 (this->n_dofs_per_cell() == 0),
832 (x_source_fe.
get_name().find(
"FESystem<") == 0) ||
846 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
850 spacedim>::ExcInterpolationNotImplemented()));
859 std::vector<FullMatrix<double>> base_matrices(this->n_base_elements());
860 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
862 base_matrices[i].reinit(base_element(i).n_dofs_per_cell(),
864 base_element(i).get_interpolation_matrix(source_fe.
base_element(i),
871 interpolation_matrix = 0;
872 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
874 if (this->system_to_base_table[i].
first ==
876 interpolation_matrix(i, j) =
877 (base_matrices[this->system_to_base_table[i].first.first](
878 this->system_to_base_table[i].second,
884template <
int dim,
int spacedim>
887 const unsigned int child,
894 "Restriction matrices are only available for refined cells!"));
898 if (this->restriction[refinement_case - 1][child].n() == 0)
900 std::lock_guard<std::mutex> lock(this->mutex);
903 if (this->restriction[refinement_case - 1][child].n() ==
904 this->n_dofs_per_cell())
905 return this->restriction[refinement_case - 1][child];
908 bool do_restriction =
true;
911 std::vector<const FullMatrix<double> *> base_matrices(
912 this->n_base_elements());
914 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
917 &base_element(i).get_restriction_matrix(child, refinement_case);
918 if (base_matrices[i]->n() != base_element(i).n_dofs_per_cell())
919 do_restriction =
false;
928 this->n_dofs_per_cell());
938 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
939 for (
unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
946 if (this->system_to_base_table[i].
first !=
947 this->system_to_base_table[j].
first)
951 const unsigned int base =
952 this->system_to_base_table[i].first.first;
954 const unsigned int base_index_i =
955 this->system_to_base_table[i].second,
957 this->system_to_base_table[j].second;
962 (*base_matrices[base])(base_index_i, base_index_j);
966 this->restriction[refinement_case - 1][child]));
970 return this->restriction[refinement_case - 1][child];
975template <
int dim,
int spacedim>
978 const unsigned int child,
985 "Restriction matrices are only available for refined cells!"));
990 if (this->prolongation[refinement_case - 1][child].n() == 0)
992 std::lock_guard<std::mutex> lock(this->mutex);
994 if (this->prolongation[refinement_case - 1][child].n() ==
995 this->n_dofs_per_cell())
996 return this->prolongation[refinement_case - 1][child];
998 bool do_prolongation =
true;
999 std::vector<const FullMatrix<double> *> base_matrices(
1000 this->n_base_elements());
1001 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
1004 &base_element(i).get_prolongation_matrix(child, refinement_case);
1005 if (base_matrices[i]->n() != base_element(i).n_dofs_per_cell())
1006 do_prolongation =
false;
1011 if (do_prolongation)
1014 this->n_dofs_per_cell());
1016 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
1017 for (
unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
1019 if (this->system_to_base_table[i].
first !=
1020 this->system_to_base_table[j].
first)
1022 const unsigned int base =
1023 this->system_to_base_table[i].first.first;
1025 const unsigned int base_index_i =
1026 this->system_to_base_table[i].second,
1028 this->system_to_base_table[j].second;
1030 (*base_matrices[base])(base_index_i, base_index_j);
1033 this->prolongation[refinement_case - 1][child]));
1037 return this->prolongation[refinement_case - 1][child];
1041template <
int dim,
int spacedim>
1044 const unsigned int face,
1045 const bool face_orientation,
1046 const bool face_flip,
1047 const bool face_rotation)
const
1052 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
1053 face_base_index = this->face_system_to_base_index(face_dof_index, face);
1055 const unsigned int base_face_to_cell_index =
1056 this->base_element(face_base_index.first.first)
1057 .face_to_cell_index(face_base_index.second,
1068 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int> target =
1069 std::make_pair(face_base_index.first, base_face_to_cell_index);
1070 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
1071 if (this->system_to_base_index(i) == target)
1086template <
int dim,
int spacedim>
1093 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1094 out |= base_element(base_no).requires_update_flags(flags);
1100template <
int dim,
int spacedim>
1101std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1117 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1118 data_ptr = std::make_unique<InternalData>(this->n_base_elements());
1119 auto &data =
dynamic_cast<InternalData &
>(*data_ptr);
1120 data.update_each = requires_update_flags(flags);
1134 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1137 &base_fe_output_object = data.get_fe_output_object(base_no);
1140 base_element(base_no),
1141 flags | base_element(base_no).requires_update_flags(flags));
1149 auto base_fe_data = base_element(base_no).get_data(flags,
1152 base_fe_output_object);
1154 data.set_fe_data(base_no, std::move(base_fe_data));
1163template <
int dim,
int spacedim>
1164std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1180 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1181 data_ptr = std::make_unique<InternalData>(this->n_base_elements());
1182 auto &data =
dynamic_cast<InternalData &
>(*data_ptr);
1183 data.update_each = requires_update_flags(flags);
1197 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1200 &base_fe_output_object = data.get_fe_output_object(base_no);
1203 base_element(base_no),
1204 flags | base_element(base_no).requires_update_flags(flags));
1212 auto base_fe_data = base_element(base_no).get_face_data(
1213 flags, mapping, quadrature, base_fe_output_object);
1215 data.set_fe_data(base_no, std::move(base_fe_data));
1226template <
int dim,
int spacedim>
1227std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1243 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1244 data_ptr = std::make_unique<InternalData>(this->n_base_elements());
1245 auto &data =
dynamic_cast<InternalData &
>(*data_ptr);
1247 data.update_each = requires_update_flags(flags);
1261 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1264 &base_fe_output_object = data.get_fe_output_object(base_no);
1267 base_element(base_no),
1268 flags | base_element(base_no).requires_update_flags(flags));
1276 auto base_fe_data = base_element(base_no).get_subface_data(
1277 flags, mapping, quadrature, base_fe_output_object);
1279 data.set_fe_data(base_no, std::move(base_fe_data));
1287template <
int dim,
int spacedim>
1302 compute_fill(mapping,
1304 invalid_face_number,
1305 invalid_face_number,
1316template <
int dim,
int spacedim>
1320 const unsigned int face_no,
1331 compute_fill(mapping,
1334 invalid_face_number,
1345template <
int dim,
int spacedim>
1349 const unsigned int face_no,
1350 const unsigned int sub_no,
1361 compute_fill(mapping,
1375template <
int dim,
int spacedim>
1376template <
class Q_or_QC>
1381 const unsigned int face_no,
1382 const unsigned int sub_no,
1383 const Q_or_QC & quadrature,
1396 Assert(
dynamic_cast<const InternalData *
>(&fe_internal) !=
nullptr,
1398 const InternalData &fe_data =
static_cast<const InternalData &
>(fe_internal);
1415 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1419 fe_data.get_fe_data(base_no);
1422 &base_data = fe_data.get_fe_output_object(base_no);
1428 const Quadrature<dim - 1> * sub_face_quadrature =
nullptr;
1432 if (face_no == invalid_face_number)
1434 const Subscriptor *quadrature_base_pointer = &quadrature;
1436 quadrature_base_pointer) !=
nullptr,
1441 n_q_points = cell_quadrature->
size();
1443 else if (sub_no == invalid_face_number)
1445 const Subscriptor *quadrature_base_pointer = &quadrature;
1447 quadrature_base_pointer) !=
nullptr,
1453 quadrature_base_pointer);
1455 (*face_quadrature)[face_quadrature->size() == 1 ? 0 : face_no]
1460 const Subscriptor *quadrature_base_pointer = &quadrature;
1462 quadrature_base_pointer) !=
nullptr,
1465 sub_face_quadrature =
1466 static_cast<const Quadrature<dim - 1
> *>(quadrature_base_pointer);
1467 n_q_points = sub_face_quadrature->
size();
1479 if (face_no == invalid_face_number)
1488 else if (sub_no == invalid_face_number)
1501 *sub_face_quadrature,
1519 primitive_offset_tables[base_no],
1530 nonprimitive_offset_tables[base_no],
1538template <
int dim,
int spacedim>
1546 for (
unsigned int base = 0; base < this->n_base_elements(); ++base)
1547 if (base_element(base).constraints_are_implemented() ==
false)
1553 const unsigned int face_no = 0;
1555 this->interface_constraints.TableBase<2, double>::reinit(
1556 this->interface_constraints_size());
1566 for (
unsigned int n = 0; n < this->interface_constraints.n(); ++n)
1567 for (
unsigned int m = 0; m < this->interface_constraints.m(); ++m)
1576 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
1577 n_index = this->face_system_to_base_table[face_no][n];
1581 std::pair<std::pair<unsigned int, unsigned int>,
unsigned int> m_index;
1598 if (m < this->n_dofs_per_vertex())
1599 m_index = this->system_to_base_table[m];
1603 const unsigned int index_in_line =
1604 (m - this->n_dofs_per_vertex()) % this->n_dofs_per_line();
1605 const unsigned int sub_line =
1606 (m - this->n_dofs_per_vertex()) / this->n_dofs_per_line();
1613 const unsigned int tmp1 =
1614 2 * this->n_dofs_per_vertex() + index_in_line;
1616 this->face_system_to_base_table[face_no][tmp1].first;
1628 this->face_system_to_base_table[face_no][tmp1].
second >=
1630 base_element(m_index.first.first).n_dofs_per_vertex(),
1632 const unsigned int tmp2 =
1633 this->face_system_to_base_table[face_no][tmp1].second -
1634 2 * base_element(m_index.first.first).n_dofs_per_vertex();
1635 Assert(tmp2 < base_element(m_index.first.first)
1639 base_element(m_index.first.first).n_dofs_per_vertex() +
1640 base_element(m_index.first.first).n_dofs_per_line() *
1654 if (m < 5 * this->n_dofs_per_vertex())
1655 m_index = this->system_to_base_table[m];
1658 if (m < 5 * this->n_dofs_per_vertex() +
1659 12 * this->n_dofs_per_line())
1662 const unsigned int index_in_line =
1663 (m - 5 * this->n_dofs_per_vertex()) %
1664 this->n_dofs_per_line();
1665 const unsigned int sub_line =
1666 (m - 5 * this->n_dofs_per_vertex()) /
1667 this->n_dofs_per_line();
1670 const unsigned int tmp1 =
1671 4 * this->n_dofs_per_vertex() + index_in_line;
1673 this->face_system_to_base_table[face_no][tmp1].first;
1676 this->face_system_to_base_table[face_no][tmp1].
second >=
1678 base_element(m_index.first.first).n_dofs_per_vertex(),
1680 const unsigned int tmp2 =
1681 this->face_system_to_base_table[face_no][tmp1].second -
1682 4 * base_element(m_index.first.first).n_dofs_per_vertex();
1683 Assert(tmp2 < base_element(m_index.first.first)
1688 base_element(m_index.first.first).n_dofs_per_vertex() +
1689 base_element(m_index.first.first).n_dofs_per_line() *
1697 const unsigned int index_in_quad =
1698 (m - 5 * this->n_dofs_per_vertex() -
1699 12 * this->n_dofs_per_line()) %
1700 this->n_dofs_per_quad(face_no);
1701 Assert(index_in_quad < this->n_dofs_per_quad(face_no),
1703 const unsigned int sub_quad =
1704 ((m - 5 * this->n_dofs_per_vertex() -
1705 12 * this->n_dofs_per_line()) /
1706 this->n_dofs_per_quad(face_no));
1709 const unsigned int tmp1 = 4 * this->n_dofs_per_vertex() +
1710 4 * this->n_dofs_per_line() +
1713 this->face_system_to_base_table[face_no].size(),
1716 this->face_system_to_base_table[face_no][tmp1].first;
1719 this->face_system_to_base_table[face_no][tmp1].
second >=
1720 4 * base_element(m_index.first.first)
1721 .n_dofs_per_vertex() +
1723 base_element(m_index.first.first).n_dofs_per_line(),
1725 const unsigned int tmp2 =
1726 this->face_system_to_base_table[face_no][tmp1].second -
1728 base_element(m_index.first.first).n_dofs_per_vertex() -
1729 4 * base_element(m_index.first.first).n_dofs_per_line();
1730 Assert(tmp2 < base_element(m_index.first.first)
1731 .n_dofs_per_quad(face_no),
1735 base_element(m_index.first.first).n_dofs_per_vertex() +
1736 12 * base_element(m_index.first.first).n_dofs_per_line() +
1737 base_element(m_index.first.first)
1738 .n_dofs_per_quad(face_no) *
1753 if (n_index.first == m_index.first)
1754 this->interface_constraints(m, n) =
1755 (base_element(n_index.first.first)
1756 .constraints()(m_index.second, n_index.second));
1762template <
int dim,
int spacedim>
1766 const std::vector<unsigned int> & multiplicities)
1768 Assert(fes.size() == multiplicities.size(),
1771 ExcMessage(
"Need to pass at least one finite element."));
1772 Assert(count_nonzeros(multiplicities) > 0,
1773 ExcMessage(
"You only passed FiniteElements with multiplicity 0."));
1778 this->base_to_block_indices.reinit(0, 0);
1780 for (
unsigned int i = 0; i < fes.size(); ++i)
1781 if (multiplicities[i] > 0)
1782 this->base_to_block_indices.push_back(multiplicities[i]);
1787 unsigned int ind = 0;
1788 for (
unsigned int i = 0; i < fes.size(); ++i)
1789 if (multiplicities[i] > 0)
1792 base_elements[ind] = {fes[i]->clone(), multiplicities[i]};
1806 this->system_to_component_table.resize(this->n_dofs_per_cell());
1809 this->system_to_component_table,
1810 this->component_to_base_table,
1813 this->face_system_to_component_table.resize(this->n_unique_faces());
1815 for (
unsigned int face_no = 0; face_no < this->n_unique_faces(); ++face_no)
1817 this->face_system_to_component_table[face_no].resize(
1818 this->n_dofs_per_face(face_no));
1821 this->face_system_to_base_table[face_no],
1822 this->face_system_to_component_table[face_no],
1843 for (
unsigned int base_el = 0; base_el < this->n_base_elements(); ++base_el)
1844 if (!base_element(base_el).has_support_points() &&
1845 base_element(base_el).n_dofs_per_cell() != 0)
1854 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
1856 const unsigned int base = this->system_to_base_table[i].first.first,
1857 base_index = this->system_to_base_table[i].second;
1862 base_element(base).unit_support_points[base_index];
1867 primitive_offset_tables.resize(this->n_base_elements());
1869 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1870 if (base_element(base_no).is_primitive())
1871 primitive_offset_tables[base_no] =
1876 nonprimitive_offset_tables.resize(this->n_base_elements());
1878 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1879 if (!base_element(base_no).is_primitive())
1880 nonprimitive_offset_tables[base_no] =
1887 for (
unsigned int face_no = 0; face_no < this->n_unique_faces();
1900 bool flag_has_no_support_points =
false;
1902 for (
unsigned int base_el = 0; base_el < this->n_base_elements();
1904 if (!base_element(base_el).has_support_points() &&
1905 (base_element(base_el).n_dofs_per_face(face_no) > 0))
1907 this->unit_face_support_points[face_no].resize(0);
1908 flag_has_no_support_points =
true;
1913 if (flag_has_no_support_points)
1918 this->unit_face_support_points[face_no].resize(
1919 this->n_dofs_per_face(face_no));
1921 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
1923 const unsigned int base_i =
1924 this->face_system_to_base_table[face_no][i].first.first;
1925 const unsigned int index_in_base =
1926 this->face_system_to_base_table[face_no][i].second;
1930 base_element(base_i).unit_face_support_points[face_no].size(),
1933 this->unit_face_support_points[face_no][i] =
1934 base_element(base_i)
1935 .unit_face_support_points[face_no][index_in_base];
1948 generalized_support_points_index_table.resize(this->n_base_elements());
1950 for (
unsigned int base = 0; base < this->n_base_elements(); ++base)
1961 if (!base_element(base).has_generalized_support_points())
1964 for (
const auto &point :
1965 base_element(base).get_generalized_support_points())
1969 std::find(std::begin(this->generalized_support_points),
1970 std::end(this->generalized_support_points),
1973 if (p == std::end(this->generalized_support_points))
1976 const auto n = this->generalized_support_points.size();
1977 generalized_support_points_index_table[base].push_back(n);
1978 this->generalized_support_points.push_back(point);
1983 const auto n = p - std::begin(this->generalized_support_points);
1984 generalized_support_points_index_table[base].push_back(n);
1991 for (
unsigned int i = 0; i < base_elements.size(); ++i)
1993 if (!base_element(i).has_generalized_support_points())
1996 const auto &points =
1997 base_elements[i].first->get_generalized_support_points();
1998 for (
unsigned int j = 0; j < points.size(); ++j)
2000 const auto n = generalized_support_points_index_table[i][j];
2001 Assert(this->generalized_support_points[n] == points[j],
2011 for (
unsigned int face_no = 0; face_no < this->n_unique_faces();
2016 Assert(this->adjust_quad_dof_index_for_face_orientation_table[face_no]
2019 this->n_dofs_per_quad(face_no),
2024 unsigned int index = 0;
2025 for (
unsigned int b = 0;
b < this->n_base_elements(); ++
b)
2028 this->base_element(b)
2029 .adjust_quad_dof_index_for_face_orientation_table[face_no];
2030 for (
unsigned int c = 0; c < this->element_multiplicity(b); ++c)
2032 for (
unsigned int i = 0; i < temp.size(0); ++i)
2033 for (
unsigned int j = 0;
2037 this->adjust_quad_dof_index_for_face_orientation_table
2038 [face_no](index + i, j) = temp(i, j);
2039 index += temp.size(0);
2046 Assert(this->adjust_line_dof_index_for_line_orientation_table.size() ==
2047 this->n_dofs_per_line(),
2049 unsigned int index = 0;
2050 for (
unsigned int b = 0;
b < this->n_base_elements(); ++
b)
2052 const std::vector<int> &temp2 =
2053 this->base_element(b)
2054 .adjust_line_dof_index_for_line_orientation_table;
2055 for (
unsigned int c = 0; c < this->element_multiplicity(b); ++c)
2060 this->adjust_line_dof_index_for_line_orientation_table.begin() +
2062 index += temp2.size();
2069 init_tasks.join_all();
2074template <
int dim,
int spacedim>
2078 for (
unsigned int b = 0;
b < this->n_base_elements(); ++
b)
2079 if (base_element(b).hp_constraints_are_implemented() ==
false)
2087template <
int dim,
int spacedim>
2092 const unsigned int face_no)
const
2094 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
2096 this->n_dofs_per_face(face_no)));
2111 if (
const auto *fe_other_system =
2115 interpolation_matrix = 0;
2119 unsigned int base_index = 0, base_index_other = 0;
2120 unsigned int multiplicity = 0, multiplicity_other = 0;
2128 fe_other_system->base_element(
2135 base_to_base_interpolation.reinit(base_other.n_dofs_per_face(face_no),
2138 base_to_base_interpolation,
2144 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
2145 if (this->face_system_to_base_index(i, face_no).first ==
2146 std::make_pair(base_index, multiplicity))
2147 for (
unsigned int j = 0;
2148 j < fe_other_system->n_dofs_per_face(face_no);
2150 if (fe_other_system->face_system_to_base_index(j, face_no)
2152 std::make_pair(base_index_other, multiplicity_other))
2153 interpolation_matrix(j, i) = base_to_base_interpolation(
2154 fe_other_system->face_system_to_base_index(j, face_no)
2156 this->face_system_to_base_index(i, face_no).second);
2162 if (multiplicity == this->element_multiplicity(base_index))
2167 ++multiplicity_other;
2168 if (multiplicity_other ==
2169 fe_other_system->element_multiplicity(base_index_other))
2171 multiplicity_other = 0;
2177 if (base_index == this->n_base_elements())
2179 Assert(base_index_other == fe_other_system->n_base_elements(),
2186 Assert(base_index_other != fe_other_system->n_base_elements(),
2197 spacedim>::ExcInterpolationNotImplemented()));
2203template <
int dim,
int spacedim>
2207 const unsigned int subface,
2209 const unsigned int face_no)
const
2212 (x_source_fe.
get_name().find(
"FESystem<") == 0) ||
2216 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
2218 this->n_dofs_per_face(face_no)));
2235 if (fe_other_system !=
nullptr)
2238 interpolation_matrix = 0;
2242 unsigned int base_index = 0, base_index_other = 0;
2243 unsigned int multiplicity = 0, multiplicity_other = 0;
2258 base_to_base_interpolation.reinit(base_other.n_dofs_per_face(face_no),
2262 base_to_base_interpolation,
2268 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
2269 if (this->face_system_to_base_index(i, face_no).first ==
2270 std::make_pair(base_index, multiplicity))
2271 for (
unsigned int j = 0;
2276 std::make_pair(base_index_other, multiplicity_other))
2277 interpolation_matrix(j, i) = base_to_base_interpolation(
2280 this->face_system_to_base_index(i, face_no).second);
2286 if (multiplicity == this->element_multiplicity(base_index))
2291 ++multiplicity_other;
2292 if (multiplicity_other ==
2295 multiplicity_other = 0;
2301 if (base_index == this->n_base_elements())
2318 fe_other_system !=
nullptr,
2320 spacedim>::ExcInterpolationNotImplemented()));
2326template <
int dim,
int spacedim>
2327template <
int structdim>
2328std::vector<std::pair<unsigned int, unsigned int>>
2331 const unsigned int face_no)
const
2350 unsigned int base_index = 0, base_index_other = 0;
2351 unsigned int multiplicity = 0, multiplicity_other = 0;
2355 unsigned int dof_offset = 0, dof_offset_other = 0;
2357 std::vector<std::pair<unsigned int, unsigned int>> identities;
2371 std::vector<std::pair<unsigned int, unsigned int>> base_identities;
2388 for (
const auto &base_identity : base_identities)
2389 identities.emplace_back(base_identity.
first + dof_offset,
2390 base_identity.
second + dof_offset_other);
2393 dof_offset += base.template n_dofs_per_object<structdim>();
2395 base_other.template n_dofs_per_object<structdim>();
2401 if (multiplicity == this->element_multiplicity(base_index))
2406 ++multiplicity_other;
2407 if (multiplicity_other ==
2410 multiplicity_other = 0;
2416 if (base_index == this->n_base_elements())
2434 return std::vector<std::pair<unsigned int, unsigned int>>();
2440template <
int dim,
int spacedim>
2441std::vector<std::pair<unsigned int, unsigned int>>
2445 return hp_object_dof_identities<0>(fe_other);
2448template <
int dim,
int spacedim>
2449std::vector<std::pair<unsigned int, unsigned int>>
2453 return hp_object_dof_identities<1>(fe_other);
2458template <
int dim,
int spacedim>
2459std::vector<std::pair<unsigned int, unsigned int>>
2462 const unsigned int face_no)
const
2464 return hp_object_dof_identities<2>(fe_other, face_no);
2469template <
int dim,
int spacedim>
2473 const unsigned int codim)
const
2484 Assert(this->n_components() == fe_sys_other->n_components(),
2486 Assert(this->n_base_elements() == fe_sys_other->n_base_elements(),
2493 for (
unsigned int b = 0;
b < this->n_base_elements(); ++
b)
2495 Assert(this->base_element(b).n_components() ==
2496 fe_sys_other->base_element(b).n_components(),
2498 Assert(this->element_multiplicity(b) ==
2499 fe_sys_other->element_multiplicity(b),
2505 (this->base_element(b).compare_for_domination(
2506 fe_sys_other->base_element(b), codim));
2507 domination = domination & base_domination;
2519template <
int dim,
int spacedim>
2524 return *base_elements[
index].first;
2529template <
int dim,
int spacedim>
2532 const unsigned int shape_index,
2533 const unsigned int face_index)
const
2535 return (base_element(this->system_to_base_index(shape_index).
first.first)
2536 .has_support_on_face(this->system_to_base_index(shape_index).second,
2542template <
int dim,
int spacedim>
2548 (this->unit_support_points.size() == 0),
2557 return (base_element(this->system_to_base_index(index).
first.first)
2558 .unit_support_point(this->system_to_base_index(index).second));
2563template <
int dim,
int spacedim>
2566 const unsigned int index,
2567 const unsigned int face_no)
const
2571 (this->unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
2572 .size() == this->n_dofs_per_face(face_no)) ||
2573 (this->unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
2575 (typename
FiniteElement<dim, spacedim>::ExcFEHasNoSupportPoints()));
2578 if (this->unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
2581 ->unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
2587 base_element(this->face_system_to_base_index(index, face_no).
first.first)
2588 .unit_face_support_point(
2589 this->face_system_to_base_index(index, face_no).second, face_no));
2594template <
int dim,
int spacedim>
2595std::pair<Table<2, bool>, std::vector<unsigned int>>
2601 Table<2, bool> constant_modes(this->n_components(), this->n_dofs_per_cell());
2602 std::vector<unsigned int> components;
2603 for (
unsigned int i = 0; i < base_elements.size(); ++i)
2605 const std::pair<Table<2, bool>, std::vector<unsigned int>> base_table =
2606 base_elements[i].first->get_constant_modes();
2608 const unsigned int element_multiplicity = this->element_multiplicity(i);
2613 const unsigned int comp = components.size();
2614 if (constant_modes.n_rows() <
2615 comp + base_table.first.n_rows() * element_multiplicity)
2617 Table<2, bool> new_constant_modes(comp + base_table.first.n_rows() *
2618 element_multiplicity,
2619 constant_modes.n_cols());
2620 for (
unsigned int r = 0; r < comp; ++r)
2621 for (
unsigned int c = 0; c < this->n_dofs_per_cell(); ++c)
2622 new_constant_modes(r, c) = constant_modes(r, c);
2623 constant_modes.swap(new_constant_modes);
2628 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
2630 std::pair<std::pair<unsigned int, unsigned int>,
unsigned int> ind =
2631 this->system_to_base_index(k);
2632 if (ind.first.first == i)
2633 for (
unsigned int c = 0; c < base_table.first.n_rows(); ++c)
2634 constant_modes(comp +
2635 ind.first.second * base_table.first.n_rows() + c,
2636 k) = base_table.first(c, ind.second);
2638 for (
unsigned int r = 0; r < element_multiplicity; ++r)
2639 for (
const unsigned int c : base_table.
second)
2640 components.push_back(
2641 comp + r * this->base_elements[i].
first->n_components() + c);
2644 return std::pair<Table<2, bool>, std::vector<unsigned int>>(constant_modes,
2650template <
int dim,
int spacedim>
2654 std::vector<double> & dof_values)
const
2656 Assert(this->has_generalized_support_points(),
2657 ExcMessage(
"The FESystem does not have generalized support points"));
2660 this->get_generalized_support_points().size());
2663 std::vector<double> base_dof_values;
2664 std::vector<Vector<double>> base_point_values;
2669 unsigned int current_vector_component = 0;
2670 for (
unsigned int base = 0; base < base_elements.size(); ++base)
2675 const auto & base_element = this->base_element(base);
2676 const unsigned int multiplicity = this->element_multiplicity(base);
2677 const unsigned int n_base_dofs = base_element.n_dofs_per_cell();
2678 const unsigned int n_base_components = base_element.n_components();
2683 if (n_base_dofs == 0)
2685 current_vector_component += multiplicity * n_base_components;
2689 if (base_element.has_generalized_support_points())
2691 const unsigned int n_base_points =
2692 base_element.get_generalized_support_points().size();
2694 base_dof_values.resize(n_base_dofs);
2695 base_point_values.resize(n_base_points);
2697 for (
unsigned int m = 0; m < multiplicity;
2698 ++m, current_vector_component += n_base_components)
2702 for (
unsigned int j = 0; j < base_point_values.size(); ++j)
2704 base_point_values[j].reinit(n_base_components,
false);
2707 generalized_support_points_index_table[base][j];
2712 std::begin(point_values[n]) + current_vector_component;
2713 const auto end =
begin + n_base_components;
2714 std::copy(begin, end, std::begin(base_point_values[j]));
2718 .convert_generalized_support_point_values_to_dof_values(
2719 base_point_values, base_dof_values);
2728 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
2729 if (this->system_to_base_index(i).first ==
2730 std::make_pair(base, m))
2732 base_dof_values[this->system_to_base_index(i).second];
2744 for (
unsigned int m = 0; m < multiplicity; ++m)
2745 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
2746 if (this->system_to_base_index(i).first ==
2747 std::make_pair(base, m))
2748 dof_values[i] = std::numeric_limits<double>::signaling_NaN();
2750 current_vector_component += multiplicity * n_base_components;
2757template <
int dim,
int spacedim>
2765 sizeof(base_elements));
2766 for (
unsigned int i = 0; i < base_elements.size(); ++i)
2774#include "fe_system.inst"
void set_fe_data(const unsigned int base_no, std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase >)
InternalData(const unsigned int n_base_elements)
std::vector< std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > > base_fe_datas
FiniteElement< dim, spacedim >::InternalDataBase & get_fe_data(const unsigned int base_no) const
internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > & get_fe_output_object(const unsigned int base_no) const
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const override
virtual Point< dim > unit_support_point(const unsigned int index) const override
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const override
std::vector< std::pair< std::unique_ptr< const FiniteElement< dim, spacedim > >, unsigned int > > base_elements
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
virtual std::string get_name() const override
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
virtual const FiniteElement< dim, spacedim > & get_sub_fe(const unsigned int first_component, const unsigned int n_selected_components) const override
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const override
void compute_fill(const Mapping< dim, spacedim > &mapping, const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Q_or_QC &quadrature, const CellSimilarity::Similarity cell_similarity, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) 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 bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const hp::QCollection< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual bool hp_constraints_are_implemented() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
void build_interface_constraints()
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual Point< dim - 1 > unit_face_support_point(const unsigned int index, const unsigned int face_no=0) const override
void initialize(const std::vector< const FiniteElement< dim, spacedim > * > &fes, const std::vector< unsigned int > &multiplicities)
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &dof_values) const override
std::vector< std::pair< unsigned int, unsigned int > > hp_object_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual std::size_t memory_consumption() const override
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const override
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_components() const
virtual std::string get_name() const =0
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
bool is_primitive() const
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > system_to_base_index(const unsigned int index) const
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > face_system_to_base_index(const unsigned int index, const unsigned int face_no=0) const
unsigned int element_multiplicity(const unsigned int index) const
unsigned int n_nonzero_components(const unsigned int i) const
unsigned int n_base_elements() const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const
virtual std::size_t memory_consumption() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Abstract base class for mapping classes.
unsigned int size() const
unsigned int max_n_quadrature_points() 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 & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_gradients
Shape function gradients.
@ update_default
No update.
Task< RT > new_task(const std::function< RT()> &function)
@ neither_element_dominates
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
std::enable_if_t< std::is_fundamental< T >::value, std::size_t > memory_consumption(const T &t)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
std::string dim_string(const int dim, const int spacedim)
std::vector< Point< dim > > unit_support_points(const std::vector< Point< 1 > > &line_support_points, const std::vector< unsigned int > &renumbering)
std::vector< typename FESystem< dim, spacedim >::BaseOffsets > setup_nonprimitive_offset_table(const FESystem< dim, spacedim > &fe, const unsigned int base_no)
Table< 2, unsigned int > setup_primitive_offset_table(const FESystem< dim, spacedim > &fe, const unsigned int base_no)
void copy_nonprimitive_base_element_values(const FESystem< dim, spacedim > &fe, const unsigned int base_no, const unsigned int n_q_points, const UpdateFlags base_flags, const std::vector< typename FESystem< dim, spacedim >::BaseOffsets > &offsets, const FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &base_data, FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data)
void copy_primitive_base_element_values(const FESystem< dim, spacedim > &fe, const unsigned int base_no, const UpdateFlags base_flags, const Table< 2, unsigned int > &base_to_system_table, const FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &base_data, FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data)
static const unsigned int invalid_unsigned_int