38 count_nonzeros(
const std::vector<unsigned int> &vec)
40 return std::count_if(vec.begin(), vec.end(), [](
const unsigned int i) {
53 template <
int dim,
int spacedim = dim>
56 const unsigned int base_no)
64 unsigned int out_index = 0;
65 for (
unsigned int system_index = 0; system_index < fe.
n_dofs_per_cell();
72 const unsigned int base_component =
74 const unsigned int base_index =
79 table[base_component][base_index] = out_index;
90 template <
int dim,
int spacedim = dim>
91 std::vector<typename FESystem<dim, spacedim>::BaseOffsets>
93 const unsigned int base_no)
95 std::vector<typename FESystem<dim, spacedim>::BaseOffsets> table;
98 unsigned int out_index = 0;
99 for (
unsigned int system_index = 0; system_index < fe.
n_dofs_per_cell();
104 const unsigned int base_index =
107 table.emplace_back();
109 table.back().n_nonzero_components =
111 unsigned int in_index = 0;
112 for (
unsigned int i = 0; i < base_index; ++i)
115 table.back().in_index = in_index;
116 table.back().out_index = out_index;
131 template <
int dim,
int spacedim = dim>
135 const unsigned int base_no,
136 const unsigned int n_q_points,
146 const unsigned int n_dofs_per_cell =
148 for (
unsigned int component = 0; component < n_components; ++component)
149 for (
unsigned int b = 0; b < n_dofs_per_cell; ++b)
151 const unsigned int out_index = base_to_system_table[component][b];
154 for (
unsigned int q = 0; q < n_q_points; ++q)
159 for (
unsigned int q = 0; q < n_q_points; ++q)
164 for (
unsigned int q = 0; q < n_q_points; ++q)
169 for (
unsigned int q = 0; q < n_q_points; ++q)
179 template <
int dim,
int spacedim = dim>
183 const unsigned int base_no,
184 const unsigned int n_q_points,
196 for (
const auto &offset : offsets)
199 for (
unsigned int s = 0; s < offset.n_nonzero_components; ++s)
200 for (
unsigned int q = 0; q < n_q_points; ++q)
205 for (
unsigned int s = 0; s < offset.n_nonzero_components; ++s)
206 for (
unsigned int q = 0; q < n_q_points; ++q)
211 for (
unsigned int s = 0; s < offset.n_nonzero_components; ++s)
212 for (
unsigned int q = 0; q < n_q_points; ++q)
217 for (
unsigned int s = 0; s < offset.n_nonzero_components; ++s)
218 for (
unsigned int q = 0; q < n_q_points; ++q)
227template <
int dim,
int spacedim>
229 const unsigned int n_base_elements)
230 : base_fe_datas(n_base_elements)
231 , base_fe_output_objects(n_base_elements)
236template <
int dim,
int spacedim>
245template <
int dim,
int spacedim>
248 const unsigned int base_no)
const
251 return *base_fe_datas[base_no];
256template <
int dim,
int spacedim>
259 const unsigned int base_no,
263 base_fe_datas[base_no] = std::move(ptr);
268template <
int dim,
int spacedim>
271 const unsigned int base_no)
const
274 return base_fe_output_objects[base_no];
282template <
int dim,
int spacedim>
286template <
int dim,
int spacedim>
288 const unsigned int n_elements)
296 std::vector<const FiniteElement<dim, spacedim> *> fes;
298 std::vector<unsigned int> multiplicities;
299 multiplicities.push_back(n_elements);
300 initialize(fes, multiplicities);
305template <
int dim,
int spacedim>
307 const unsigned int n1,
309 const unsigned int n2)
317 , base_elements(static_cast<
int>(n1 > 0) + static_cast<
int>(n2 > 0))
319 std::vector<const FiniteElement<dim, spacedim> *> fes;
322 std::vector<unsigned int> multiplicities;
323 multiplicities.push_back(n1);
324 multiplicities.push_back(n2);
325 initialize(fes, multiplicities);
330template <
int dim,
int spacedim>
332 const unsigned int n1,
334 const unsigned int n2,
336 const unsigned int n3)
351 , base_elements(static_cast<
int>(n1 > 0) + static_cast<
int>(n2 > 0) +
352 static_cast<
int>(n3 > 0))
354 std::vector<const FiniteElement<dim, spacedim> *> fes;
358 std::vector<unsigned int> multiplicities;
359 multiplicities.push_back(n1);
360 multiplicities.push_back(n2);
361 multiplicities.push_back(n3);
362 initialize(fes, multiplicities);
367template <
int dim,
int spacedim>
369 const unsigned int n1,
371 const unsigned int n2,
373 const unsigned int n3,
375 const unsigned int n4)
401 , base_elements(static_cast<
int>(n1 > 0) + static_cast<
int>(n2 > 0) +
402 static_cast<
int>(n3 > 0) + static_cast<
int>(n4 > 0))
404 std::vector<const FiniteElement<dim, spacedim> *> fes;
409 std::vector<unsigned int> multiplicities;
410 multiplicities.push_back(n1);
411 multiplicities.push_back(n2);
412 multiplicities.push_back(n3);
413 multiplicities.push_back(n4);
414 initialize(fes, multiplicities);
419template <
int dim,
int spacedim>
421 const unsigned int n1,
423 const unsigned int n2,
425 const unsigned int n3,
427 const unsigned int n4,
429 const unsigned int n5)
432 multiply_dof_numbers(&fe1, n1, &fe2, n2, &fe3, n3, &fe4, n4, &fe5, n5),
453 , base_elements(static_cast<
int>(n1 > 0) + static_cast<
int>(n2 > 0) +
454 static_cast<
int>(n3 > 0) + static_cast<
int>(n4 > 0) +
455 static_cast<
int>(n5 > 0))
457 std::vector<const FiniteElement<dim, spacedim> *> fes;
463 std::vector<unsigned int> multiplicities;
464 multiplicities.push_back(n1);
465 multiplicities.push_back(n2);
466 multiplicities.push_back(n3);
467 multiplicities.push_back(n4);
468 multiplicities.push_back(n5);
469 initialize(fes, multiplicities);
474template <
int dim,
int spacedim>
477 const std::vector<unsigned int> & multiplicities)
484 , base_elements(count_nonzeros(multiplicities))
486 initialize(fes, multiplicities);
491template <
int dim,
int spacedim>
502 std::ostringstream namebuf;
505 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
507 namebuf << base_element(i).get_name();
508 if (this->element_multiplicity(i) != 1)
509 namebuf <<
'^' << this->element_multiplicity(i);
510 if (i != this->n_base_elements() - 1)
515 return namebuf.str();
520template <
int dim,
int spacedim>
521std::unique_ptr<FiniteElement<dim, spacedim>>
524 std::vector<const FiniteElement<dim, spacedim> *> fes;
525 std::vector<unsigned int> multiplicities;
527 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
529 fes.push_back(&base_element(i));
530 multiplicities.push_back(this->element_multiplicity(i));
532 return std::make_unique<FESystem<dim, spacedim>>(fes, multiplicities);
537template <
int dim,
int spacedim>
540 const unsigned int first_component,
541 const unsigned int n_selected_components)
const
543 Assert(first_component + n_selected_components <= this->n_components(),
544 ExcMessage(
"Invalid arguments (not a part of this FiniteElement)."));
546 const unsigned int base_index =
547 this->component_to_base_table[first_component].first.first;
548 const unsigned int component_in_base =
549 this->component_to_base_table[first_component].first.second;
550 const unsigned int base_components =
551 this->base_element(base_index).n_components();
555 if (n_selected_components <= base_components)
556 return this->base_element(base_index)
557 .get_sub_fe(component_in_base, n_selected_components);
559 Assert(n_selected_components == this->n_components(),
560 ExcMessage(
"You can not select a part of a FiniteElement."));
566template <
int dim,
int spacedim>
572 Assert(this->is_primitive(i),
576 return (base_element(this->system_to_base_table[i].
first.first)
577 .shape_value(this->system_to_base_table[i].second, p));
582template <
int dim,
int spacedim>
585 const unsigned int i,
587 const unsigned int component)
const
594 if (this->nonzero_components[i][component] ==
false)
602 const unsigned int base = this->component_to_base_index(component).first;
603 const unsigned int component_in_base =
604 this->component_to_base_index(component).second;
612 return (base_element(base).shape_value_component(
613 this->system_to_base_table[i].
second, p, component_in_base));
618template <
int dim,
int spacedim>
624 Assert(this->is_primitive(i),
628 return (base_element(this->system_to_base_table[i].
first.first)
629 .shape_grad(this->system_to_base_table[i].second, p));
634template <
int dim,
int spacedim>
637 const unsigned int i,
639 const unsigned int component)
const
645 if (this->nonzero_components[i][component] ==
false)
650 const unsigned int base = this->component_to_base_index(component).first;
651 const unsigned int component_in_base =
652 this->component_to_base_index(component).second;
657 return (base_element(base).shape_grad_component(
658 this->system_to_base_table[i].
second, p, component_in_base));
663template <
int dim,
int spacedim>
669 Assert(this->is_primitive(i),
673 return (base_element(this->system_to_base_table[i].
first.first)
674 .shape_grad_grad(this->system_to_base_table[i].second, p));
679template <
int dim,
int spacedim>
682 const unsigned int i,
684 const unsigned int component)
const
690 if (this->nonzero_components[i][component] ==
false)
695 const unsigned int base = this->component_to_base_index(component).first;
696 const unsigned int component_in_base =
697 this->component_to_base_index(component).second;
702 return (base_element(base).shape_grad_grad_component(
703 this->system_to_base_table[i].
second, p, component_in_base));
708template <
int dim,
int spacedim>
714 Assert(this->is_primitive(i),
718 return (base_element(this->system_to_base_table[i].
first.first)
719 .shape_3rd_derivative(this->system_to_base_table[i].second, p));
724template <
int dim,
int spacedim>
727 const unsigned int i,
729 const unsigned int component)
const
735 if (this->nonzero_components[i][component] ==
false)
740 const unsigned int base = this->component_to_base_index(component).first;
741 const unsigned int component_in_base =
742 this->component_to_base_index(component).second;
747 return (base_element(base).shape_3rd_derivative_component(
748 this->system_to_base_table[i].
second, p, component_in_base));
753template <
int dim,
int spacedim>
759 Assert(this->is_primitive(i),
763 return (base_element(this->system_to_base_table[i].
first.first)
764 .shape_4th_derivative(this->system_to_base_table[i].second, p));
769template <
int dim,
int spacedim>
772 const unsigned int i,
774 const unsigned int component)
const
780 if (this->nonzero_components[i][component] ==
false)
785 const unsigned int base = this->component_to_base_index(component).first;
786 const unsigned int component_in_base =
787 this->component_to_base_index(component).second;
792 return (base_element(base).shape_4th_derivative_component(
793 this->system_to_base_table[i].
second, p, component_in_base));
798template <
int dim,
int spacedim>
808 Assert((interpolation_matrix.
m() == this->n_dofs_per_cell()) ||
811 this->n_dofs_per_cell()));
813 (this->n_dofs_per_cell() == 0),
823 (x_source_fe.
get_name().find(
"FESystem<") == 0) ||
837 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
841 spacedim>::ExcInterpolationNotImplemented()));
850 std::vector<FullMatrix<double>> base_matrices(this->n_base_elements());
851 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
853 base_matrices[i].reinit(base_element(i).n_dofs_per_cell(),
855 base_element(i).get_interpolation_matrix(source_fe.
base_element(i),
862 interpolation_matrix = 0;
863 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
865 if (this->system_to_base_table[i].
first ==
867 interpolation_matrix(i, j) =
868 (base_matrices[this->system_to_base_table[i].first.first](
869 this->system_to_base_table[i].second,
875template <
int dim,
int spacedim>
878 const unsigned int child,
885 "Restriction matrices are only available for refined cells!"));
889 if (this->restriction[refinement_case - 1][child].n() == 0)
891 std::lock_guard<std::mutex> lock(this->mutex);
894 if (this->restriction[refinement_case - 1][child].n() ==
895 this->n_dofs_per_cell())
896 return this->restriction[refinement_case - 1][child];
899 bool do_restriction =
true;
902 std::vector<const FullMatrix<double> *> base_matrices(
903 this->n_base_elements());
905 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
908 &base_element(i).get_restriction_matrix(child, refinement_case);
909 if (base_matrices[i]->n() != base_element(i).n_dofs_per_cell())
910 do_restriction =
false;
919 this->n_dofs_per_cell());
929 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
930 for (
unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
937 if (this->system_to_base_table[i].
first !=
938 this->system_to_base_table[j].
first)
942 const unsigned int base =
943 this->system_to_base_table[i].first.first;
945 const unsigned int base_index_i =
946 this->system_to_base_table[i].second,
948 this->system_to_base_table[j].second;
953 (*base_matrices[base])(base_index_i, base_index_j);
957 this->restriction[refinement_case - 1][child]));
961 return this->restriction[refinement_case - 1][child];
966template <
int dim,
int spacedim>
969 const unsigned int child,
976 "Restriction matrices are only available for refined cells!"));
981 if (this->prolongation[refinement_case - 1][child].n() == 0)
983 std::lock_guard<std::mutex> lock(this->mutex);
985 if (this->prolongation[refinement_case - 1][child].n() ==
986 this->n_dofs_per_cell())
987 return this->prolongation[refinement_case - 1][child];
989 bool do_prolongation =
true;
990 std::vector<const FullMatrix<double> *> base_matrices(
991 this->n_base_elements());
992 for (
unsigned int i = 0; i < this->n_base_elements(); ++i)
995 &base_element(i).get_prolongation_matrix(child, refinement_case);
996 if (base_matrices[i]->n() != base_element(i).n_dofs_per_cell())
997 do_prolongation =
false;
1002 if (do_prolongation)
1005 this->n_dofs_per_cell());
1007 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
1008 for (
unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
1010 if (this->system_to_base_table[i].
first !=
1011 this->system_to_base_table[j].
first)
1013 const unsigned int base =
1014 this->system_to_base_table[i].first.first;
1016 const unsigned int base_index_i =
1017 this->system_to_base_table[i].second,
1019 this->system_to_base_table[j].second;
1021 (*base_matrices[base])(base_index_i, base_index_j);
1024 this->prolongation[refinement_case - 1][child]));
1028 return this->prolongation[refinement_case - 1][child];
1032template <
int dim,
int spacedim>
1035 const unsigned int face,
1036 const bool face_orientation,
1037 const bool face_flip,
1038 const bool face_rotation)
const
1043 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
1044 face_base_index = this->face_system_to_base_index(face_dof_index, face);
1046 const unsigned int base_face_to_cell_index =
1047 this->base_element(face_base_index.first.first)
1048 .face_to_cell_index(face_base_index.second,
1059 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int> target =
1060 std::make_pair(face_base_index.first, base_face_to_cell_index);
1061 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
1062 if (this->system_to_base_index(i) == target)
1077template <
int dim,
int spacedim>
1084 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1085 out |= base_element(base_no).requires_update_flags(flags);
1091template <
int dim,
int spacedim>
1092std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1108 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1109 data_ptr = std::make_unique<InternalData>(this->n_base_elements());
1110 auto &data =
dynamic_cast<InternalData &
>(*data_ptr);
1111 data.update_each = requires_update_flags(flags);
1125 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1128 &base_fe_output_object = data.get_fe_output_object(base_no);
1131 base_element(base_no),
1132 flags | base_element(base_no).requires_update_flags(flags));
1140 auto base_fe_data = base_element(base_no).get_data(flags,
1143 base_fe_output_object);
1145 data.set_fe_data(base_no, std::move(base_fe_data));
1154template <
int dim,
int spacedim>
1155std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1171 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1172 data_ptr = std::make_unique<InternalData>(this->n_base_elements());
1173 auto &data =
dynamic_cast<InternalData &
>(*data_ptr);
1174 data.update_each = requires_update_flags(flags);
1188 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1191 &base_fe_output_object = data.get_fe_output_object(base_no);
1194 base_element(base_no),
1195 flags | base_element(base_no).requires_update_flags(flags));
1203 auto base_fe_data = base_element(base_no).get_face_data(
1204 flags, mapping, quadrature, base_fe_output_object);
1206 data.set_fe_data(base_no, std::move(base_fe_data));
1217template <
int dim,
int spacedim>
1218std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1234 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1235 data_ptr = std::make_unique<InternalData>(this->n_base_elements());
1236 auto &data =
dynamic_cast<InternalData &
>(*data_ptr);
1238 data.update_each = requires_update_flags(flags);
1252 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1255 &base_fe_output_object = data.get_fe_output_object(base_no);
1258 base_element(base_no),
1259 flags | base_element(base_no).requires_update_flags(flags));
1267 auto base_fe_data = base_element(base_no).get_subface_data(
1268 flags, mapping, quadrature, base_fe_output_object);
1270 data.set_fe_data(base_no, std::move(base_fe_data));
1278template <
int dim,
int spacedim>
1286 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
1294 compute_fill(mapping,
1296 invalid_face_number,
1297 invalid_face_number,
1308template <
int dim,
int spacedim>
1312 const unsigned int face_no,
1316 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
1324 compute_fill(mapping,
1327 invalid_face_number,
1338template <
int dim,
int spacedim>
1342 const unsigned int face_no,
1343 const unsigned int sub_no,
1347 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
1355 compute_fill(mapping,
1369template <
int dim,
int spacedim>
1370template <
class Q_or_QC>
1375 const unsigned int face_no,
1376 const unsigned int sub_no,
1377 const Q_or_QC & quadrature,
1390 Assert(
dynamic_cast<const InternalData *
>(&fe_internal) !=
nullptr,
1392 const InternalData &fe_data =
static_cast<const InternalData &
>(fe_internal);
1409 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1413 fe_data.get_fe_data(base_no);
1416 &base_data = fe_data.get_fe_output_object(base_no);
1422 const Quadrature<dim - 1> * sub_face_quadrature =
nullptr;
1426 if (face_no == invalid_face_number)
1428 const Subscriptor *quadrature_base_pointer = &quadrature;
1430 quadrature_base_pointer) !=
nullptr,
1435 n_q_points = cell_quadrature->
size();
1437 else if (sub_no == invalid_face_number)
1439 const Subscriptor *quadrature_base_pointer = &quadrature;
1441 quadrature_base_pointer) !=
nullptr,
1447 quadrature_base_pointer);
1449 (*face_quadrature)[face_quadrature->size() == 1 ? 0 : face_no]
1454 const Subscriptor *quadrature_base_pointer = &quadrature;
1456 quadrature_base_pointer) !=
nullptr,
1459 sub_face_quadrature =
1460 static_cast<const Quadrature<dim - 1
> *>(quadrature_base_pointer);
1461 n_q_points = sub_face_quadrature->
size();
1473 if (face_no == invalid_face_number)
1482 else if (sub_no == invalid_face_number)
1495 *sub_face_quadrature,
1514 primitive_offset_tables[base_no],
1525 nonprimitive_offset_tables[base_no],
1533template <
int dim,
int spacedim>
1539 if (this->n_unique_faces() != 1)
1542 const unsigned int face_no = 0;
1548 for (
unsigned int base = 0; base < this->n_base_elements(); ++base)
1549 if (base_element(base).constraints_are_implemented() ==
false)
1552 this->interface_constraints.TableBase<2, double>::reinit(
1553 this->interface_constraints_size());
1563 for (
unsigned int n = 0; n < this->interface_constraints.n(); ++n)
1564 for (
unsigned int m = 0; m < this->interface_constraints.m(); ++m)
1573 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
1574 n_index = this->face_system_to_base_table[face_no][n];
1578 std::pair<std::pair<unsigned int, unsigned int>,
unsigned int> m_index;
1595 if (m < this->n_dofs_per_vertex())
1596 m_index = this->system_to_base_table[m];
1600 const unsigned int index_in_line =
1601 (m - this->n_dofs_per_vertex()) % this->n_dofs_per_line();
1602 const unsigned int sub_line =
1603 (m - this->n_dofs_per_vertex()) / this->n_dofs_per_line();
1610 const unsigned int tmp1 =
1611 2 * this->n_dofs_per_vertex() + index_in_line;
1613 this->face_system_to_base_table[face_no][tmp1].first;
1625 this->face_system_to_base_table[face_no][tmp1].
second >=
1627 base_element(m_index.first.first).n_dofs_per_vertex(),
1629 const unsigned int tmp2 =
1630 this->face_system_to_base_table[face_no][tmp1].second -
1631 2 * base_element(m_index.first.first).n_dofs_per_vertex();
1632 Assert(tmp2 < base_element(m_index.first.first)
1636 base_element(m_index.first.first).n_dofs_per_vertex() +
1637 base_element(m_index.first.first).n_dofs_per_line() *
1651 if (m < 5 * this->n_dofs_per_vertex())
1652 m_index = this->system_to_base_table[m];
1655 if (m < 5 * this->n_dofs_per_vertex() +
1656 12 * this->n_dofs_per_line())
1659 const unsigned int index_in_line =
1660 (m - 5 * this->n_dofs_per_vertex()) %
1661 this->n_dofs_per_line();
1662 const unsigned int sub_line =
1663 (m - 5 * this->n_dofs_per_vertex()) /
1664 this->n_dofs_per_line();
1667 const unsigned int tmp1 =
1668 4 * this->n_dofs_per_vertex() + index_in_line;
1670 this->face_system_to_base_table[face_no][tmp1].first;
1673 this->face_system_to_base_table[face_no][tmp1].
second >=
1675 base_element(m_index.first.first).n_dofs_per_vertex(),
1677 const unsigned int tmp2 =
1678 this->face_system_to_base_table[face_no][tmp1].second -
1679 4 * base_element(m_index.first.first).n_dofs_per_vertex();
1680 Assert(tmp2 < base_element(m_index.first.first)
1685 base_element(m_index.first.first).n_dofs_per_vertex() +
1686 base_element(m_index.first.first).n_dofs_per_line() *
1694 const unsigned int index_in_quad =
1695 (m - 5 * this->n_dofs_per_vertex() -
1696 12 * this->n_dofs_per_line()) %
1697 this->n_dofs_per_quad(face_no);
1698 Assert(index_in_quad < this->n_dofs_per_quad(face_no),
1700 const unsigned int sub_quad =
1701 ((m - 5 * this->n_dofs_per_vertex() -
1702 12 * this->n_dofs_per_line()) /
1703 this->n_dofs_per_quad(face_no));
1706 const unsigned int tmp1 = 4 * this->n_dofs_per_vertex() +
1707 4 * this->n_dofs_per_line() +
1710 this->face_system_to_base_table[face_no].size(),
1713 this->face_system_to_base_table[face_no][tmp1].first;
1716 this->face_system_to_base_table[face_no][tmp1].
second >=
1717 4 * base_element(m_index.first.first)
1718 .n_dofs_per_vertex() +
1720 base_element(m_index.first.first).n_dofs_per_line(),
1722 const unsigned int tmp2 =
1723 this->face_system_to_base_table[face_no][tmp1].second -
1725 base_element(m_index.first.first).n_dofs_per_vertex() -
1726 4 * base_element(m_index.first.first).n_dofs_per_line();
1727 Assert(tmp2 < base_element(m_index.first.first)
1728 .n_dofs_per_quad(face_no),
1732 base_element(m_index.first.first).n_dofs_per_vertex() +
1733 12 * base_element(m_index.first.first).n_dofs_per_line() +
1734 base_element(m_index.first.first)
1735 .n_dofs_per_quad(face_no) *
1750 if (n_index.first == m_index.first)
1751 this->interface_constraints(m, n) =
1752 (base_element(n_index.first.first)
1753 .constraints()(m_index.second, n_index.second));
1759template <
int dim,
int spacedim>
1763 const std::vector<unsigned int> & multiplicities)
1765 Assert(fes.size() == multiplicities.size(),
1768 ExcMessage(
"Need to pass at least one finite element."));
1769 Assert(count_nonzeros(multiplicities) > 0,
1770 ExcMessage(
"You only passed FiniteElements with multiplicity 0."));
1775 this->base_to_block_indices.reinit(0, 0);
1777 for (
unsigned int i = 0; i < fes.size(); ++i)
1778 if (multiplicities[i] > 0)
1779 this->base_to_block_indices.push_back(multiplicities[i]);
1784 unsigned int ind = 0;
1785 for (
unsigned int i = 0; i < fes.size(); ++i)
1786 if (multiplicities[i] > 0)
1789 base_elements[ind] = {fes[i]->clone(), multiplicities[i]};
1803 this->system_to_component_table.resize(this->n_dofs_per_cell());
1806 this->system_to_component_table,
1807 this->component_to_base_table,
1810 this->face_system_to_component_table.resize(this->n_unique_faces());
1812 for (
unsigned int face_no = 0; face_no < this->n_unique_faces(); ++face_no)
1814 this->face_system_to_component_table[face_no].resize(
1815 this->n_dofs_per_face(face_no));
1818 this->face_system_to_base_table[face_no],
1819 this->face_system_to_component_table[face_no],
1840 for (
unsigned int base_el = 0; base_el < this->n_base_elements(); ++base_el)
1841 if (!base_element(base_el).has_support_points() &&
1842 base_element(base_el).n_dofs_per_cell() != 0)
1851 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
1853 const unsigned int base = this->system_to_base_table[i].first.first,
1854 base_index = this->system_to_base_table[i].second;
1859 base_element(base).unit_support_points[base_index];
1864 primitive_offset_tables.resize(this->n_base_elements());
1866 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1867 if (base_element(base_no).is_primitive())
1868 primitive_offset_tables[base_no] =
1873 nonprimitive_offset_tables.resize(this->n_base_elements());
1875 for (
unsigned int base_no = 0; base_no < this->n_base_elements(); ++base_no)
1876 if (!base_element(base_no).is_primitive())
1877 nonprimitive_offset_tables[base_no] =
1884 for (
unsigned int face_no = 0; face_no < this->n_unique_faces();
1897 bool flag_has_no_support_points =
false;
1899 for (
unsigned int base_el = 0; base_el < this->n_base_elements();
1901 if (!base_element(base_el).has_support_points() &&
1902 (base_element(base_el).n_dofs_per_face(face_no) > 0))
1904 this->unit_face_support_points[face_no].resize(0);
1905 flag_has_no_support_points =
true;
1910 if (flag_has_no_support_points)
1915 this->unit_face_support_points[face_no].resize(
1916 this->n_dofs_per_face(face_no));
1918 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
1920 const unsigned int base_i =
1921 this->face_system_to_base_table[face_no][i].first.first;
1922 const unsigned int index_in_base =
1923 this->face_system_to_base_table[face_no][i].second;
1927 base_element(base_i).unit_face_support_points[face_no].size(),
1930 this->unit_face_support_points[face_no][i] =
1931 base_element(base_i)
1932 .unit_face_support_points[face_no][index_in_base];
1945 generalized_support_points_index_table.resize(this->n_base_elements());
1947 for (
unsigned int base = 0; base < this->n_base_elements(); ++base)
1958 if (!base_element(base).has_generalized_support_points())
1961 for (
const auto &point :
1962 base_element(base).get_generalized_support_points())
1966 std::find(std::begin(this->generalized_support_points),
1967 std::end(this->generalized_support_points),
1970 if (p == std::end(this->generalized_support_points))
1973 const auto n = this->generalized_support_points.size();
1974 generalized_support_points_index_table[base].push_back(n);
1975 this->generalized_support_points.push_back(point);
1980 const auto n = p - std::begin(this->generalized_support_points);
1981 generalized_support_points_index_table[base].push_back(n);
1988 for (
unsigned int i = 0; i < base_elements.size(); ++i)
1990 if (!base_element(i).has_generalized_support_points())
1993 const auto &points =
1994 base_elements[i].first->get_generalized_support_points();
1995 for (
unsigned int j = 0; j < points.size(); ++j)
1997 const auto n = generalized_support_points_index_table[i][j];
1998 Assert(this->generalized_support_points[n] == points[j],
2008 for (
unsigned int face_no = 0; face_no < this->n_unique_faces();
2013 Assert(this->adjust_quad_dof_index_for_face_orientation_table[face_no]
2014 .n_elements() == 8 * this->n_dofs_per_quad(face_no),
2019 unsigned int index = 0;
2020 for (
unsigned int b = 0;
b < this->n_base_elements(); ++
b)
2023 this->base_element(b)
2024 .adjust_quad_dof_index_for_face_orientation_table[face_no];
2025 for (
unsigned int c = 0; c < this->element_multiplicity(b); ++c)
2027 for (
unsigned int i = 0; i < temp.size(0); ++i)
2028 for (
unsigned int j = 0; j < 8; ++j)
2029 this->adjust_quad_dof_index_for_face_orientation_table
2030 [face_no](index + i, j) = temp(i, j);
2031 index += temp.size(0);
2038 Assert(this->adjust_line_dof_index_for_line_orientation_table.size() ==
2039 this->n_dofs_per_line(),
2041 unsigned int index = 0;
2042 for (
unsigned int b = 0;
b < this->n_base_elements(); ++
b)
2044 const std::vector<int> &temp2 =
2045 this->base_element(b)
2046 .adjust_line_dof_index_for_line_orientation_table;
2047 for (
unsigned int c = 0; c < this->element_multiplicity(b); ++c)
2052 this->adjust_line_dof_index_for_line_orientation_table.begin() +
2054 index += temp2.size();
2061 init_tasks.join_all();
2066template <
int dim,
int spacedim>
2070 for (
unsigned int b = 0;
b < this->n_base_elements(); ++
b)
2071 if (base_element(b).hp_constraints_are_implemented() ==
false)
2079template <
int dim,
int spacedim>
2084 const unsigned int face_no)
const
2086 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
2088 this->n_dofs_per_face(face_no)));
2103 if (
const auto *fe_other_system =
2107 interpolation_matrix = 0;
2111 unsigned int base_index = 0, base_index_other = 0;
2112 unsigned int multiplicity = 0, multiplicity_other = 0;
2120 fe_other_system->base_element(
2127 base_to_base_interpolation.reinit(base_other.n_dofs_per_face(face_no),
2130 base_to_base_interpolation,
2136 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
2137 if (this->face_system_to_base_index(i, face_no).first ==
2138 std::make_pair(base_index, multiplicity))
2139 for (
unsigned int j = 0;
2140 j < fe_other_system->n_dofs_per_face(face_no);
2142 if (fe_other_system->face_system_to_base_index(j, face_no)
2144 std::make_pair(base_index_other, multiplicity_other))
2145 interpolation_matrix(j, i) = base_to_base_interpolation(
2146 fe_other_system->face_system_to_base_index(j, face_no)
2148 this->face_system_to_base_index(i, face_no).second);
2154 if (multiplicity == this->element_multiplicity(base_index))
2159 ++multiplicity_other;
2160 if (multiplicity_other ==
2161 fe_other_system->element_multiplicity(base_index_other))
2163 multiplicity_other = 0;
2169 if (base_index == this->n_base_elements())
2171 Assert(base_index_other == fe_other_system->n_base_elements(),
2178 Assert(base_index_other != fe_other_system->n_base_elements(),
2189 spacedim>::ExcInterpolationNotImplemented()));
2195template <
int dim,
int spacedim>
2199 const unsigned int subface,
2201 const unsigned int face_no)
const
2204 (x_source_fe.
get_name().find(
"FESystem<") == 0) ||
2208 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
2210 this->n_dofs_per_face(face_no)));
2227 if (fe_other_system !=
nullptr)
2230 interpolation_matrix = 0;
2234 unsigned int base_index = 0, base_index_other = 0;
2235 unsigned int multiplicity = 0, multiplicity_other = 0;
2250 base_to_base_interpolation.reinit(base_other.n_dofs_per_face(face_no),
2254 base_to_base_interpolation,
2260 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
2261 if (this->face_system_to_base_index(i, face_no).first ==
2262 std::make_pair(base_index, multiplicity))
2263 for (
unsigned int j = 0;
2268 std::make_pair(base_index_other, multiplicity_other))
2269 interpolation_matrix(j, i) = base_to_base_interpolation(
2272 this->face_system_to_base_index(i, face_no).second);
2278 if (multiplicity == this->element_multiplicity(base_index))
2283 ++multiplicity_other;
2284 if (multiplicity_other ==
2287 multiplicity_other = 0;
2293 if (base_index == this->n_base_elements())
2310 fe_other_system !=
nullptr,
2312 spacedim>::ExcInterpolationNotImplemented()));
2318template <
int dim,
int spacedim>
2319template <
int structdim>
2320std::vector<std::pair<unsigned int, unsigned int>>
2323 const unsigned int face_no)
const
2342 unsigned int base_index = 0, base_index_other = 0;
2343 unsigned int multiplicity = 0, multiplicity_other = 0;
2347 unsigned int dof_offset = 0, dof_offset_other = 0;
2349 std::vector<std::pair<unsigned int, unsigned int>> identities;
2363 std::vector<std::pair<unsigned int, unsigned int>> base_identities;
2380 for (
const auto &base_identity : base_identities)
2381 identities.emplace_back(base_identity.first + dof_offset,
2382 base_identity.second + dof_offset_other);
2385 dof_offset += base.template n_dofs_per_object<structdim>();
2387 base_other.template n_dofs_per_object<structdim>();
2393 if (multiplicity == this->element_multiplicity(base_index))
2398 ++multiplicity_other;
2399 if (multiplicity_other ==
2402 multiplicity_other = 0;
2408 if (base_index == this->n_base_elements())
2426 return std::vector<std::pair<unsigned int, unsigned int>>();
2432template <
int dim,
int spacedim>
2433std::vector<std::pair<unsigned int, unsigned int>>
2437 return hp_object_dof_identities<0>(fe_other);
2440template <
int dim,
int spacedim>
2441std::vector<std::pair<unsigned int, unsigned int>>
2445 return hp_object_dof_identities<1>(fe_other);
2450template <
int dim,
int spacedim>
2451std::vector<std::pair<unsigned int, unsigned int>>
2454 const unsigned int face_no)
const
2456 return hp_object_dof_identities<2>(fe_other, face_no);
2461template <
int dim,
int spacedim>
2465 const unsigned int codim)
const
2476 Assert(this->n_components() == fe_sys_other->n_components(),
2478 Assert(this->n_base_elements() == fe_sys_other->n_base_elements(),
2485 for (
unsigned int b = 0;
b < this->n_base_elements(); ++
b)
2487 Assert(this->base_element(b).n_components() ==
2488 fe_sys_other->base_element(b).n_components(),
2490 Assert(this->element_multiplicity(b) ==
2491 fe_sys_other->element_multiplicity(b),
2497 (this->base_element(b).compare_for_domination(
2498 fe_sys_other->base_element(b), codim));
2499 domination = domination & base_domination;
2511template <
int dim,
int spacedim>
2516 return *base_elements[
index].first;
2521template <
int dim,
int spacedim>
2524 const unsigned int shape_index,
2525 const unsigned int face_index)
const
2527 return (base_element(this->system_to_base_index(shape_index).
first.first)
2528 .has_support_on_face(this->system_to_base_index(shape_index).second,
2534template <
int dim,
int spacedim>
2540 (this->unit_support_points.size() == 0),
2549 return (base_element(this->system_to_base_index(index).
first.first)
2550 .unit_support_point(this->system_to_base_index(index).second));
2555template <
int dim,
int spacedim>
2558 const unsigned int index,
2559 const unsigned int face_no)
const
2563 (this->unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
2564 .size() == this->n_dofs_per_face(face_no)) ||
2565 (this->unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
2570 if (this->unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
2573 ->unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
2579 base_element(this->face_system_to_base_index(index, face_no).
first.first)
2580 .unit_face_support_point(
2581 this->face_system_to_base_index(index, face_no).second, face_no));
2586template <
int dim,
int spacedim>
2587std::pair<Table<2, bool>, std::vector<unsigned int>>
2593 Table<2, bool> constant_modes(this->n_components(), this->n_dofs_per_cell());
2594 std::vector<unsigned int> components;
2595 for (
unsigned int i = 0; i < base_elements.size(); ++i)
2597 const std::pair<Table<2, bool>, std::vector<unsigned int>> base_table =
2598 base_elements[i].first->get_constant_modes();
2600 const unsigned int element_multiplicity = this->element_multiplicity(i);
2605 const unsigned int comp = components.size();
2606 if (constant_modes.n_rows() <
2607 comp + base_table.first.n_rows() * element_multiplicity)
2609 Table<2, bool> new_constant_modes(comp + base_table.first.n_rows() *
2610 element_multiplicity,
2611 constant_modes.n_cols());
2612 for (
unsigned int r = 0; r < comp; ++r)
2613 for (
unsigned int c = 0; c < this->n_dofs_per_cell(); ++c)
2614 new_constant_modes(r, c) = constant_modes(r, c);
2615 constant_modes.swap(new_constant_modes);
2620 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
2622 std::pair<std::pair<unsigned int, unsigned int>,
unsigned int> ind =
2623 this->system_to_base_index(k);
2624 if (ind.first.first == i)
2625 for (
unsigned int c = 0; c < base_table.first.n_rows(); ++c)
2626 constant_modes(comp +
2627 ind.first.second * base_table.first.n_rows() + c,
2628 k) = base_table.first(c, ind.second);
2630 for (
unsigned int r = 0; r < element_multiplicity; ++r)
2631 for (
const unsigned int c : base_table.second)
2632 components.push_back(
2633 comp + r * this->base_elements[i].
first->n_components() + c);
2636 return std::pair<Table<2, bool>, std::vector<unsigned int>>(constant_modes,
2642template <
int dim,
int spacedim>
2646 std::vector<double> & dof_values)
const
2648 Assert(this->has_generalized_support_points(),
2649 ExcMessage(
"The FESystem does not have generalized support points"));
2652 this->get_generalized_support_points().size());
2655 std::vector<double> base_dof_values;
2656 std::vector<Vector<double>> base_point_values;
2661 unsigned int current_vector_component = 0;
2662 for (
unsigned int base = 0; base < base_elements.size(); ++base)
2667 const auto & base_element = this->base_element(base);
2668 const unsigned int multiplicity = this->element_multiplicity(base);
2669 const unsigned int n_base_dofs = base_element.n_dofs_per_cell();
2670 const unsigned int n_base_components = base_element.n_components();
2675 if (n_base_dofs == 0)
2677 current_vector_component += multiplicity * n_base_components;
2681 if (base_element.has_generalized_support_points())
2683 const unsigned int n_base_points =
2684 base_element.get_generalized_support_points().size();
2686 base_dof_values.resize(n_base_dofs);
2687 base_point_values.resize(n_base_points);
2689 for (
unsigned int m = 0; m < multiplicity;
2690 ++m, current_vector_component += n_base_components)
2694 for (
unsigned int j = 0; j < base_point_values.size(); ++j)
2696 base_point_values[j].reinit(n_base_components,
false);
2699 generalized_support_points_index_table[base][j];
2704 std::begin(point_values[n]) + current_vector_component;
2705 const auto end =
begin + n_base_components;
2706 std::copy(begin, end, std::begin(base_point_values[j]));
2710 .convert_generalized_support_point_values_to_dof_values(
2711 base_point_values, base_dof_values);
2720 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
2721 if (this->system_to_base_index(i).first ==
2722 std::make_pair(base, m))
2724 base_dof_values[this->system_to_base_index(i).second];
2736 for (
unsigned int m = 0; m < multiplicity; ++m)
2737 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
2738 if (this->system_to_base_index(i).first ==
2739 std::make_pair(base, m))
2740 dof_values[i] = std::numeric_limits<double>::signaling_NaN();
2742 current_vector_component += multiplicity * n_base_components;
2749template <
int dim,
int spacedim>
2757 sizeof(base_elements));
2758 for (
unsigned int i = 0; i < base_elements.size(); ++i)
2766#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 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 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 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
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
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 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 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
bool is_primitive() const
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 void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, 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
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
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
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
@ 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.
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)
Task< RT > new_task(const std::function< RT()> &function)
@ neither_element_dominates
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type 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 unsigned int n_q_points, 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