41#ifdef DEAL_II_WITH_MPI
61 check_primary_dof_list(
63 const std::vector<types::global_dof_index> &primary_dof_list)
65 const unsigned int N = primary_dof_list.size();
68 for (
unsigned int i = 0; i < N; ++i)
69 for (
unsigned int j = 0; j < N; ++j)
70 tmp(i, j) = face_interpolation_matrix(primary_dof_list[i], j);
81 double diagonal_sum = 0;
82 for (
unsigned int i = 0; i < N; ++i)
83 diagonal_sum += std::fabs(tmp(i, i));
84 const double typical_diagonal_element = diagonal_sum / N;
88 std::vector<unsigned int> p(N);
89 for (
unsigned int i = 0; i < N; ++i)
92 for (
unsigned int j = 0; j < N; ++j)
96 double max = std::fabs(tmp(j, j));
98 for (
unsigned int i = j + 1; i < N; ++i)
100 if (std::fabs(tmp(i, j)) > max)
102 max = std::fabs(tmp(i, j));
109 if (max < 1.e-12 * typical_diagonal_element)
115 for (
unsigned int k = 0; k < N; ++k)
116 std::swap(tmp(j, k), tmp(r, k));
118 std::swap(p[j], p[r]);
122 const double hr = 1. / tmp(j, j);
124 for (
unsigned int k = 0; k < N; ++k)
128 for (
unsigned int i = 0; i < N; ++i)
132 tmp(i, k) -= tmp(i, j) * tmp(j, k) * hr;
135 for (
unsigned int i = 0; i < N; ++i)
171 template <
int dim,
int spacedim>
173 select_primary_dofs_for_face_restriction(
177 std::vector<bool> & primary_dof_mask)
183 const unsigned int face_no = 0;
218 std::vector<types::global_dof_index> primary_dof_list;
219 unsigned int index = 0;
224 unsigned int dofs_added = 0;
233 primary_dof_list.push_back(index + i);
236 if (check_primary_dof_list(face_interpolation_matrix,
237 primary_dof_list) ==
true)
242 primary_dof_list.pop_back();
255 unsigned int dofs_added = 0;
261 primary_dof_list.push_back(index + i);
262 if (check_primary_dof_list(face_interpolation_matrix,
263 primary_dof_list) ==
true)
266 primary_dof_list.pop_back();
278 unsigned int dofs_added = 0;
284 primary_dof_list.push_back(index + i);
285 if (check_primary_dof_list(face_interpolation_matrix,
286 primary_dof_list) ==
true)
289 primary_dof_list.pop_back();
300 std::fill(primary_dof_mask.begin(), primary_dof_mask.end(),
false);
301 for (
const auto dof : primary_dof_list)
302 primary_dof_mask[dof] = true;
311 template <
int dim,
int spacedim>
313 ensure_existence_of_primary_dof_mask(
317 std::unique_ptr<std::vector<bool>> &primary_dof_mask)
323 const unsigned int face_no = 0;
325 if (primary_dof_mask ==
nullptr)
329 select_primary_dofs_for_face_restriction(fe1,
331 face_interpolation_matrix,
343 template <
int dim,
int spacedim>
345 ensure_existence_of_face_matrix(
354 const unsigned int face_no = 0;
356 if (matrix ==
nullptr)
358 matrix = std::make_unique<FullMatrix<double>>(
369 template <
int dim,
int spacedim>
371 ensure_existence_of_subface_matrix(
374 const unsigned int subface,
381 const unsigned int face_no = 0;
383 if (matrix ==
nullptr)
385 matrix = std::make_unique<FullMatrix<double>>(
402 ensure_existence_of_split_face_matrix(
404 const std::vector<bool> & primary_dof_mask,
409 Assert(std::count(primary_dof_mask.begin(),
410 primary_dof_mask.end(),
412 static_cast<signed int>(face_interpolation_matrix.
n()),
415 if (split_matrix ==
nullptr)
417 split_matrix = std::make_unique<
420 const unsigned int n_primary_dofs = face_interpolation_matrix.
n();
421 const unsigned int n_dofs = face_interpolation_matrix.
m();
427 split_matrix->first.reinit(n_primary_dofs, n_primary_dofs);
428 split_matrix->second.reinit(n_dofs - n_primary_dofs,
431 unsigned int nth_primary_dof = 0, nth_dependent_dof = 0;
433 for (
unsigned int i = 0; i < n_dofs; ++i)
434 if (primary_dof_mask[i] ==
true)
436 for (
unsigned int j = 0; j < n_primary_dofs; ++j)
437 split_matrix->first(nth_primary_dof, j) =
438 face_interpolation_matrix(i, j);
443 for (
unsigned int j = 0; j < n_primary_dofs; ++j)
444 split_matrix->second(nth_dependent_dof, j) =
445 face_interpolation_matrix(i, j);
454 split_matrix->first.gauss_jordan();
464 template <
int dim,
int spacedim>
486 template <
typename number1,
typename number2>
489 const std::vector<types::global_dof_index> &primary_dofs,
490 const std::vector<types::global_dof_index> &dependent_dofs,
494 Assert(face_constraints.
n() == primary_dofs.size(),
496 Assert(face_constraints.
m() == dependent_dofs.size(),
498 face_constraints.
m()));
500 const unsigned int n_primary_dofs = primary_dofs.size();
501 const unsigned int n_dependent_dofs = dependent_dofs.size();
505 for (
unsigned int row = 0; row != n_dependent_dofs; ++row)
508 for (
unsigned int col = 0; col != n_primary_dofs; ++col)
513 for (
unsigned int row = 0; row != n_dependent_dofs; ++row)
533 bool is_trivial_constraint =
false;
535 for (
unsigned int i = 0; i < n_primary_dofs; ++i)
536 if (face_constraints(row, i) == 1.0)
537 if (dependent_dofs[row] == primary_dofs[i])
539 is_trivial_constraint =
true;
541 for (
unsigned int ii = 0; ii < n_primary_dofs; ++ii)
543 Assert(face_constraints(row, ii) == 0.0,
549 if (is_trivial_constraint ==
true)
555 for (
unsigned int i = 0; i < n_primary_dofs; ++i)
556 abs_sum +=
std::abs(face_constraints(row, i));
564 constraints.
add_line(dependent_dofs[row]);
565 for (
unsigned int i = 0; i < n_primary_dofs; ++i)
566 if ((face_constraints(row, i) != 0) &&
567 (std::fabs(face_constraints(row, i)) >= 1e-14 * abs_sum))
568 constraints.
add_entry(dependent_dofs[row],
570 face_constraints(row, i));
578 template <
typename number>
587 template <
typename number>
591 std::integral_constant<int, 1>)
597 template <
typename number>
606 template <
typename number>
610 std::integral_constant<int, 1>)
616 template <
typename number,
int spacedim>
626 template <
typename number,
int spacedim>
631 std::integral_constant<int, 1>)
636 template <
int dim_,
int spacedim,
typename number>
641 std::integral_constant<int, 2>)
643 const unsigned int dim = 2;
645 std::vector<types::global_dof_index> dofs_on_mother;
646 std::vector<types::global_dof_index> dofs_on_children;
658 endc = dof_handler.
end();
659 for (; cell != endc; ++cell)
663 if (cell->is_artificial())
666 for (
const unsigned int face : cell->face_indices())
667 if (cell->face(face)->has_children())
673 Assert(cell->face(face)->n_active_fe_indices() == 1,
675 Assert(cell->face(face)->fe_index_is_active(
676 cell->active_fe_index()) ==
true,
678 for (
unsigned int c = 0; c < cell->face(face)->n_children();
680 if (!cell->neighbor_child_on_subface(face, c)
682 Assert(cell->face(face)->child(c)->n_active_fe_indices() ==
688 for (
unsigned int c = 0; c < cell->face(face)->n_children();
690 if (!cell->neighbor_child_on_subface(face, c)
692 Assert(cell->face(face)->child(c)->fe_index_is_active(
693 cell->active_fe_index()) ==
true,
700 const unsigned int n_dofs_on_mother =
707 dofs_on_mother.resize(n_dofs_on_mother);
710 dofs_on_children.clear();
711 dofs_on_children.reserve(n_dofs_on_children);
721 this_face = cell->face(face);
725 unsigned int next_index = 0;
726 for (
unsigned int vertex = 0; vertex < 2; ++vertex)
729 dofs_on_mother[next_index++] =
730 this_face->vertex_dof_index(vertex, dof, fe_index);
732 dofs_on_mother[next_index++] =
733 this_face->dof_index(dof, fe_index);
737 dofs_on_children.push_back(
738 this_face->child(0)->vertex_dof_index(1, dof, fe_index));
739 for (
unsigned int child = 0; child < 2; ++child)
742 if (cell->neighbor_child_on_subface(face, child)
747 dofs_on_children.push_back(
748 this_face->child(child)->dof_index(dof, fe_index));
751 Assert(dofs_on_children.size() <= n_dofs_on_children,
755 for (
unsigned int row = 0; row != dofs_on_children.size();
758 constraints.
add_line(dofs_on_children[row]);
759 for (
unsigned int i = 0; i != dofs_on_mother.size(); ++i)
760 constraints.
add_entry(dofs_on_children[row],
773 if (!cell->at_boundary(face) &&
774 !cell->neighbor(face)->is_artificial())
776 Assert(cell->face(face)->n_active_fe_indices() == 1,
778 Assert(cell->face(face)->fe_index_is_active(
779 cell->active_fe_index()) ==
true,
788 template <
int dim_,
int spacedim,
typename number>
793 std::integral_constant<int, 3>)
795 const unsigned int dim = 3;
797 std::vector<types::global_dof_index> dofs_on_mother;
798 std::vector<types::global_dof_index> dofs_on_children;
810 endc = dof_handler.
end();
811 for (; cell != endc; ++cell)
815 if (cell->is_artificial())
818 for (
const unsigned int face : cell->face_indices())
819 if (cell->face(face)->has_children())
827 Assert(cell->face(face)->refinement_case() ==
836 Assert(cell->face(face)->fe_index_is_active(
837 cell->active_fe_index()) ==
true,
839 for (
unsigned int c = 0; c < cell->face(face)->n_children();
841 if (!cell->neighbor_child_on_subface(face, c)
844 cell->face(face)->child(c)->n_active_fe_indices(), 1);
849 for (
unsigned int c = 0; c < cell->face(face)->n_children();
851 if (!cell->neighbor_child_on_subface(face, c)
854 Assert(cell->face(face)->child(c)->fe_index_is_active(
855 cell->active_fe_index()) ==
true,
857 for (
unsigned int e = 0; e < 4; ++e)
862 ->n_active_fe_indices() == 1,
867 ->fe_index_is_active(
868 cell->active_fe_index()) ==
true,
872 for (
unsigned int e = 0; e < 4; ++e)
874 Assert(cell->face(face)->line(e)->n_active_fe_indices() ==
877 Assert(cell->face(face)->line(e)->fe_index_is_active(
878 cell->active_fe_index()) ==
true,
887 const unsigned int n_dofs_on_children =
894 dofs_on_mother.resize(n_dofs_on_mother);
897 dofs_on_children.clear();
898 dofs_on_children.reserve(n_dofs_on_children);
908 this_face = cell->face(face);
912 unsigned int next_index = 0;
913 for (
unsigned int vertex = 0; vertex < 4; ++vertex)
916 dofs_on_mother[next_index++] =
917 this_face->vertex_dof_index(vertex, dof, fe_index);
918 for (
unsigned int line = 0; line < 4; ++line)
920 dofs_on_mother[next_index++] =
921 this_face->line(line)->dof_index(dof, fe_index);
924 dofs_on_mother[next_index++] =
925 this_face->dof_index(dof, fe_index);
934 ((this_face->child(0)->vertex_index(3) ==
935 this_face->child(1)->vertex_index(2)) &&
936 (this_face->child(0)->vertex_index(3) ==
937 this_face->child(2)->vertex_index(1)) &&
938 (this_face->child(0)->vertex_index(3) ==
939 this_face->child(3)->vertex_index(0))),
943 dofs_on_children.push_back(
944 this_face->child(0)->vertex_dof_index(3, dof));
947 for (
unsigned int line = 0; line < 4; ++line)
950 dofs_on_children.push_back(
951 this_face->line(line)->child(0)->vertex_dof_index(
958 dofs_on_children.push_back(
959 this_face->child(0)->line(1)->dof_index(dof, fe_index));
961 dofs_on_children.push_back(
962 this_face->child(2)->line(1)->dof_index(dof, fe_index));
964 dofs_on_children.push_back(
965 this_face->child(0)->line(3)->dof_index(dof, fe_index));
967 dofs_on_children.push_back(
968 this_face->child(1)->line(3)->dof_index(dof, fe_index));
971 for (
unsigned int line = 0; line < 4; ++line)
972 for (
unsigned int child = 0; child < 2; ++child)
976 dofs_on_children.push_back(
977 this_face->line(line)->child(child)->dof_index(
982 for (
unsigned int child = 0; child < 4; ++child)
985 if (cell->neighbor_child_on_subface(face, child)
990 dofs_on_children.push_back(
991 this_face->child(child)->dof_index(dof, fe_index));
995 Assert(dofs_on_children.size() <= n_dofs_on_children,
999 for (
unsigned int row = 0; row != dofs_on_children.size();
1002 constraints.
add_line(dofs_on_children[row]);
1003 for (
unsigned int i = 0; i != dofs_on_mother.size(); ++i)
1004 constraints.
add_entry(dofs_on_children[row],
1017 if (!cell->at_boundary(face) &&
1018 !cell->neighbor(face)->is_artificial())
1020 Assert(cell->face(face)->n_active_fe_indices() == 1,
1022 Assert(cell->face(face)->fe_index_is_active(
1023 cell->active_fe_index()) ==
true,
1032 template <
int dim,
int spacedim,
typename number>
1050 std::vector<types::global_dof_index> primary_dofs;
1051 std::vector<types::global_dof_index> dependent_dofs;
1052 std::vector<types::global_dof_index> scratch_dofs;
1058 n_finite_elements(dof_handler), n_finite_elements(dof_handler));
1060 subface_interpolation_matrices(
1061 n_finite_elements(dof_handler),
1062 n_finite_elements(dof_handler),
1072 split_face_interpolation_matrices(n_finite_elements(dof_handler),
1073 n_finite_elements(dof_handler));
1079 n_finite_elements(dof_handler), n_finite_elements(dof_handler));
1090 if (cell->is_artificial())
1093 for (
const unsigned int face : cell->face_indices())
1094 if (cell->face(face)->has_children())
1099 if (cell->get_fe().n_dofs_per_face(face) == 0)
1102 Assert(cell->face(face)->refinement_case() ==
1113 Assert(cell->face(face)->n_active_fe_indices() == 1,
1115 Assert(cell->face(face)->fe_index_is_active(
1116 cell->active_fe_index()) ==
true,
1118 for (
unsigned int c = 0; c < cell->face(face)->n_children();
1120 if (!cell->neighbor_child_on_subface(face, c)
1122 Assert(cell->face(face)->child(c)->n_active_fe_indices() ==
1138 std::set<types::fe_index> fe_ind_face_subface;
1139 fe_ind_face_subface.insert(cell->active_fe_index());
1142 for (
unsigned int c = 0;
1143 c < cell->face(face)->n_active_descendants();
1146 const auto subcell =
1147 cell->neighbor_child_on_subface(face, c);
1148 if (!subcell->is_artificial())
1150 mother_face_dominates =
1151 mother_face_dominates &
1152 (cell->get_fe().compare_for_domination(
1153 subcell->get_fe(), 1));
1154 fe_ind_face_subface.insert(
1155 subcell->active_fe_index());
1159 switch (mother_face_dominates)
1171 primary_dofs.resize(
1172 cell->get_fe().n_dofs_per_face(face));
1174 cell->face(face)->get_dof_indices(
1175 primary_dofs, cell->active_fe_index());
1181 for (
unsigned int c = 0;
1182 c < cell->face(face)->n_children();
1185 if (cell->neighbor_child_on_subface(face, c)
1190 active_face_iterator subface =
1191 cell->face(face)->child(c);
1193 Assert(subface->n_active_fe_indices() == 1,
1197 subface->nth_active_fe_index(0);
1208 if (cell->get_fe().compare_for_domination(
1209 subface->
get_fe(subface_fe_index),
1216 dependent_dofs.resize(
1217 subface->
get_fe(subface_fe_index)
1219 subface->get_dof_indices(dependent_dofs,
1225 (void)dependent_dof;
1253 ensure_existence_of_subface_matrix(
1255 subface->
get_fe(subface_fe_index),
1257 subface_interpolation_matrices
1258 [cell->active_fe_index()][subface_fe_index][c]);
1262 filter_constraints(primary_dofs,
1264 *(subface_interpolation_matrices
1265 [cell->active_fe_index()]
1266 [subface_fe_index][c]),
1286 const ::hp::FECollection<dim, spacedim>
1316 {fe_ind_face_subface.begin(),
1317 fe_ind_face_subface.end()},
1323 "Could not find a least face dominating FE."));
1326 dof_handler.
get_fe(dominating_fe_index);
1331 cell->get_fe().n_dofs_per_face(face),
1334 ensure_existence_of_face_matrix(
1337 face_interpolation_matrices[dominating_fe_index]
1338 [cell->active_fe_index()]);
1342 ensure_existence_of_primary_dof_mask(
1345 (*face_interpolation_matrices
1346 [dominating_fe_index][cell->active_fe_index()]),
1347 primary_dof_masks[dominating_fe_index]
1348 [cell->active_fe_index()]);
1350 ensure_existence_of_split_face_matrix(
1351 *face_interpolation_matrices[dominating_fe_index]
1352 [cell->active_fe_index()],
1353 (*primary_dof_masks[dominating_fe_index]
1354 [cell->active_fe_index()]),
1355 split_face_interpolation_matrices
1356 [dominating_fe_index][cell->active_fe_index()]);
1359 &restrict_mother_to_virtual_primary_inv =
1360 (split_face_interpolation_matrices
1361 [dominating_fe_index][cell->active_fe_index()]
1365 &restrict_mother_to_virtual_dependent =
1366 (split_face_interpolation_matrices
1367 [dominating_fe_index][cell->active_fe_index()]
1372 constraint_matrix.reinit(
1373 cell->get_fe().n_dofs_per_face(face) -
1376 restrict_mother_to_virtual_dependent.
mmult(
1378 restrict_mother_to_virtual_primary_inv);
1382 scratch_dofs.resize(
1383 cell->get_fe().n_dofs_per_face(face));
1384 cell->face(face)->get_dof_indices(
1385 scratch_dofs, cell->active_fe_index());
1388 primary_dofs.clear();
1389 dependent_dofs.clear();
1390 for (
unsigned int i = 0;
1391 i < cell->get_fe().n_dofs_per_face(face);
1393 if ((*primary_dof_masks[dominating_fe_index]
1395 ->active_fe_index()])[i] ==
1397 primary_dofs.push_back(scratch_dofs[i]);
1399 dependent_dofs.push_back(scratch_dofs[i]);
1404 cell->get_fe().n_dofs_per_face(face) -
1407 filter_constraints(primary_dofs,
1416 for (
unsigned int sf = 0;
1417 sf < cell->face(face)->n_children();
1422 if (cell->neighbor_child_on_subface(face, sf)
1423 ->is_artificial() ||
1424 (dim == 2 && cell->is_ghost() &&
1425 cell->neighbor_child_on_subface(face, sf)
1431 ->n_active_fe_indices() == 1,
1435 cell->face(face)->child(sf)->nth_active_fe_index(
1438 dof_handler.
get_fe(subface_fe_index);
1445 ensure_existence_of_subface_matrix(
1449 subface_interpolation_matrices
1450 [dominating_fe_index][subface_fe_index][sf]);
1453 &restrict_subface_to_virtual = *(
1454 subface_interpolation_matrices
1455 [dominating_fe_index][subface_fe_index][sf]);
1457 constraint_matrix.reinit(
1461 restrict_subface_to_virtual.
mmult(
1463 restrict_mother_to_virtual_primary_inv);
1465 dependent_dofs.resize(
1467 cell->face(face)->child(sf)->get_dof_indices(
1468 dependent_dofs, subface_fe_index);
1470 filter_constraints(primary_dofs,
1493 Assert(cell->face(face)->fe_index_is_active(
1494 cell->active_fe_index()) ==
true,
1501 if (!cell->at_boundary(face) &&
1502 cell->neighbor(face)->is_artificial())
1508 !cell->face(face)->at_boundary() &&
1509 (cell->neighbor(face)->active_fe_index() !=
1510 cell->active_fe_index()) &&
1511 (!cell->face(face)->has_children() &&
1512 !cell->neighbor_is_coarser(face)))
1515 spacedim>::level_cell_iterator
1516 neighbor = cell->neighbor(face);
1520 cell->get_fe().compare_for_domination(neighbor->get_fe(),
1527 primary_dofs.resize(
1528 cell->get_fe().n_dofs_per_face(face));
1529 cell->face(face)->get_dof_indices(
1530 primary_dofs, cell->active_fe_index());
1535 if (primary_dofs.size() == 0)
1538 dependent_dofs.resize(
1539 neighbor->get_fe().n_dofs_per_face(face));
1540 cell->face(face)->get_dof_indices(
1541 dependent_dofs, neighbor->active_fe_index());
1545 ensure_existence_of_face_matrix(
1548 face_interpolation_matrices
1549 [cell->active_fe_index()]
1550 [neighbor->active_fe_index()]);
1556 *(face_interpolation_matrices
1557 [cell->active_fe_index()]
1558 [neighbor->active_fe_index()]),
1601 if (cell < neighbor)
1616 cell->active_fe_index();
1618 neighbor->active_fe_index();
1619 std::set<types::fe_index> fes;
1620 fes.insert(this_fe_index);
1621 fes.insert(neighbor_fe_index);
1622 const ::hp::FECollection<dim, spacedim>
1628 {fes.begin(), fes.end()}, 1);
1633 "Could not find the dominating FE for " +
1634 cell->get_fe().get_name() +
" and " +
1635 neighbor->get_fe().get_name() +
1636 " inside FECollection."));
1639 fe_collection[dominating_fe_index];
1647 cell->get_fe().n_dofs_per_face(face),
1650 ensure_existence_of_face_matrix(
1653 face_interpolation_matrices
1654 [dominating_fe_index][cell->active_fe_index()]);
1658 ensure_existence_of_primary_dof_mask(
1661 (*face_interpolation_matrices
1662 [dominating_fe_index]
1663 [cell->active_fe_index()]),
1664 primary_dof_masks[dominating_fe_index]
1665 [cell->active_fe_index()]);
1667 ensure_existence_of_split_face_matrix(
1668 *face_interpolation_matrices
1669 [dominating_fe_index][cell->active_fe_index()],
1670 (*primary_dof_masks[dominating_fe_index]
1671 [cell->active_fe_index()]),
1672 split_face_interpolation_matrices
1673 [dominating_fe_index][cell->active_fe_index()]);
1676 double> &restrict_mother_to_virtual_primary_inv =
1677 (split_face_interpolation_matrices
1678 [dominating_fe_index][cell->active_fe_index()]
1682 double> &restrict_mother_to_virtual_dependent =
1683 (split_face_interpolation_matrices
1684 [dominating_fe_index][cell->active_fe_index()]
1689 constraint_matrix.reinit(
1690 cell->get_fe().n_dofs_per_face(face) -
1693 restrict_mother_to_virtual_dependent.
mmult(
1695 restrict_mother_to_virtual_primary_inv);
1699 scratch_dofs.resize(
1700 cell->get_fe().n_dofs_per_face(face));
1701 cell->face(face)->get_dof_indices(
1702 scratch_dofs, cell->active_fe_index());
1705 primary_dofs.clear();
1706 dependent_dofs.clear();
1707 for (
unsigned int i = 0;
1708 i < cell->get_fe().n_dofs_per_face(face);
1710 if ((*primary_dof_masks[dominating_fe_index]
1711 [cell->active_fe_index()])
1713 primary_dofs.push_back(scratch_dofs[i]);
1715 dependent_dofs.push_back(scratch_dofs[i]);
1721 dependent_dofs.size(),
1722 cell->get_fe().n_dofs_per_face(face) -
1725 filter_constraints(primary_dofs,
1734 neighbor->get_fe().n_dofs_per_face(face),
1737 ensure_existence_of_face_matrix(
1740 face_interpolation_matrices
1741 [dominating_fe_index]
1742 [neighbor->active_fe_index()]);
1745 &restrict_secondface_to_virtual =
1746 *(face_interpolation_matrices
1747 [dominating_fe_index]
1748 [neighbor->active_fe_index()]);
1750 constraint_matrix.reinit(
1751 neighbor->get_fe().n_dofs_per_face(face),
1754 restrict_secondface_to_virtual.
mmult(
1756 restrict_mother_to_virtual_primary_inv);
1758 dependent_dofs.resize(
1759 neighbor->get_fe().n_dofs_per_face(face));
1760 cell->face(face)->get_dof_indices(
1761 dependent_dofs, neighbor->active_fe_index());
1763 filter_constraints(primary_dofs,
1789 template <
int dim,
int spacedim,
typename number>
1796 "The given DoFHandler does not have any DoFs. Did you forget to "
1797 "call dof_handler.distribute_dofs()?"));
1807 dof_handler, constraints, std::integral_constant<int, dim>());
1839 template <
typename FaceIterator,
typename number>
1841 set_periodicity_constraints(
1842 const FaceIterator & face_1,
1847 const bool face_orientation,
1848 const bool face_flip,
1849 const bool face_rotation,
1850 const number periodicity_factor)
1852 static const int dim = FaceIterator::AccessorType::dimension;
1853 static const int spacedim = FaceIterator::AccessorType::space_dimension;
1863 face_1->get_fe(face_1->nth_active_fe_index(0)).n_unique_faces(), 1);
1865 face_2->get_fe(face_2->nth_active_fe_index(0)).n_unique_faces(), 1);
1866 const unsigned int face_no = 0;
1871 if (face_2->has_children())
1873 Assert(face_2->n_children() ==
1877 const unsigned int dofs_per_face =
1878 face_1->get_fe(face_1->nth_active_fe_index(0))
1879 .n_dofs_per_face(face_no);
1884 for (
unsigned int c = 0; c < face_2->n_children(); ++c)
1889 const auto &fe = face_1->get_fe(face_1->nth_active_fe_index(0));
1890 fe.get_subface_interpolation_matrix(fe,
1892 subface_interpolation,
1894 subface_interpolation.mmult(child_transformation, transformation);
1896 set_periodicity_constraints(face_1,
1898 child_transformation,
1904 periodicity_factor);
1916 Assert(face_1->get_fe(face_1_index) == face_2->get_fe(face_2_index),
1918 "Matching periodic cells need to use the same finite element"));
1924 "The number of components in the mask has to be either "
1925 "zero or equal to the number of components in the finite "
1930 std::vector<types::global_dof_index> dofs_1(dofs_per_face);
1931 std::vector<types::global_dof_index> dofs_2(dofs_per_face);
1933 face_1->get_dof_indices(dofs_1, face_1_index);
1934 face_2->get_dof_indices(dofs_2, face_2_index);
1952 for (
unsigned int i = 0; i < dofs_per_face; ++i)
1974 std::map<unsigned int, unsigned int> cell_to_rotated_face_index;
1977 for (
unsigned int i = 0; i < dofs_per_face; ++i)
1996 for (
unsigned int i = 0; i < dofs_per_face; ++i)
2017 bool is_identity_constrained =
false;
2019 number constraint_factor = periodicity_factor;
2021 constexpr double eps = 1.e-13;
2022 for (
unsigned int jj = 0; jj < dofs_per_face; ++jj)
2024 const auto entry = transformation(i, jj);
2027 if (is_identity_constrained)
2032 is_identity_constrained =
false;
2035 is_identity_constrained =
true;
2037 constraint_factor = entry * periodicity_factor;
2046 if (!is_identity_constrained)
2054 affine_constraints.
add_line(dofs_2[i]);
2056 for (
unsigned int jj = 0; jj < dofs_per_face; ++jj)
2060 const unsigned int j =
2062 jj, 0, face_orientation, face_flip, face_rotation)];
2064 if (
std::abs(transformation(i, jj)) > eps)
2067 transformation(i, jj));
2078 const unsigned int j =
2080 target, 0, face_orientation, face_flip, face_rotation)];
2082 auto dof_left = dofs_1[j];
2083 auto dof_right = dofs_2[i];
2089 (dof_left < dof_right &&
2092 std::swap(dof_left, dof_right);
2093 constraint_factor = 1. / constraint_factor;
2109 bool constraints_are_cyclic =
true;
2110 number cycle_constraint_factor = constraint_factor;
2112 for (
auto test_dof = dof_right; test_dof != dof_left;)
2116 constraints_are_cyclic =
false;
2120 const auto &constraint_entries =
2122 if (constraint_entries.size() == 1)
2124 test_dof = constraint_entries[0].first;
2125 cycle_constraint_factor *= constraint_entries[0].second;
2129 constraints_are_cyclic =
false;
2144 if (constraints_are_cyclic)
2146 if (
std::abs(cycle_constraint_factor - number(1.)) > eps)
2147 affine_constraints.
add_line(dof_left);
2151 affine_constraints.
add_line(dof_left);
2176 ExcMessage(
"The periodicity constraint is too large. "
2177 "The parameter periodicity_factor might "
2178 "be too large or too small."));
2190 template <
int dim,
int spacedim>
2192 compute_transformation(
2195 const std::vector<unsigned int> & first_vector_components)
2200 const unsigned int face_no = 0;
2206 if (
matrix.m() == n_dofs_per_face)
2213 if (first_vector_components.empty() &&
matrix.m() == 0)
2226 quadrature(fe.get_unit_face_support_points(face_no));
2230 using DoFTuple =
std::array<
unsigned int, spacedim>;
2235 for (
unsigned int i = 0; i < n_dofs_per_face; ++i)
2237 std::vector<unsigned int>::const_iterator comp_it =
2238 std::find(first_vector_components.begin(),
2239 first_vector_components.end(),
2241 if (comp_it != first_vector_components.end())
2243 const unsigned int first_vector_component = *comp_it;
2246 DoFTuple vector_dofs;
2248 unsigned int n_found = 1;
2253 "Error: the finite element does not have enough components "
2254 "to define rotated periodic boundaries."));
2256 for (
unsigned int k = 0; k < n_dofs_per_face; ++k)
2257 if ((k != i) && (quadrature.point(k) == quadrature.point(i)) &&
2259 first_vector_component) &&
2261 first_vector_component + spacedim))
2265 first_vector_component] = k;
2273 for (
unsigned int i = 0; i < spacedim; ++i)
2275 transformation[vector_dofs[i]][vector_dofs[i]] = 0.;
2276 for (
unsigned int j = 0; j < spacedim; ++j)
2277 transformation[vector_dofs[i]][vector_dofs[j]] =
2282 return transformation;
2290 template <
typename FaceIterator,
typename number>
2293 const FaceIterator & face_1,
2297 const bool face_orientation,
2298 const bool face_flip,
2299 const bool face_rotation,
2301 const std::vector<unsigned int> & first_vector_components,
2302 const number periodicity_factor)
2304 static const int dim = FaceIterator::AccessorType::dimension;
2305 static const int spacedim = FaceIterator::AccessorType::space_dimension;
2307 Assert((dim != 1) || (face_orientation ==
true && face_flip ==
false &&
2308 face_rotation ==
false),
2310 "(face_orientation, face_flip, face_rotation) "
2311 "is invalid for 1d"));
2313 Assert((dim != 2) || (face_orientation ==
true && face_rotation ==
false),
2315 "(face_orientation, face_flip, face_rotation) "
2316 "is invalid for 2d"));
2319 ExcMessage(
"face_1 and face_2 are equal! Cannot constrain DoFs "
2320 "on the very same face"));
2322 Assert(face_1->at_boundary() && face_2->at_boundary(),
2323 ExcMessage(
"Faces for periodicity constraints must be on the "
2326 Assert(matrix.m() == matrix.n(),
2327 ExcMessage(
"The supplied (rotation or interpolation) matrix must "
2328 "be a square matrix"));
2330 Assert(first_vector_components.empty() || matrix.m() == spacedim,
2331 ExcMessage(
"first_vector_components is nonempty, so matrix must "
2332 "be a rotation matrix exactly of size spacedim"));
2335 if (!face_1->has_children())
2340 face_1->get_fe(face_1->nth_active_fe_index(0)).n_unique_faces(), 1);
2341 const unsigned int face_no = 0;
2344 const unsigned int n_dofs_per_face =
2345 face_1->get_fe(face_1->nth_active_fe_index(0))
2346 .n_dofs_per_face(face_no);
2348 Assert(matrix.m() == 0 ||
2349 (first_vector_components.empty() &&
2350 matrix.m() == n_dofs_per_face) ||
2351 (!first_vector_components.empty() && matrix.m() == spacedim),
2353 "The matrix must have either size 0 or spacedim "
2354 "(if first_vector_components is nonempty) "
2355 "or the size must be equal to the # of DoFs on the face "
2356 "(if first_vector_components is empty)."));
2359 if (!face_2->has_children())
2364 face_2->get_fe(face_2->nth_active_fe_index(0)).n_unique_faces(), 1);
2365 const unsigned int face_no = 0;
2368 const unsigned int n_dofs_per_face =
2369 face_2->get_fe(face_2->nth_active_fe_index(0))
2370 .n_dofs_per_face(face_no);
2372 Assert(matrix.m() == 0 ||
2373 (first_vector_components.empty() &&
2374 matrix.m() == n_dofs_per_face) ||
2375 (!first_vector_components.empty() && matrix.m() == spacedim),
2377 "The matrix must have either size 0 or spacedim "
2378 "(if first_vector_components is nonempty) "
2379 "or the size must be equal to the # of DoFs on the face "
2380 "(if first_vector_components is empty)."));
2387 static const int lookup_table_2d[2][2] = {
2393 static const int lookup_table_3d[2][2][2][4] = {
2417 if (face_1->has_children() && face_2->has_children())
2422 Assert(face_1->n_children() ==
2428 for (
unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face;
2436 j = lookup_table_2d[face_flip][i];
2439 j = lookup_table_3d[face_orientation][face_flip]
2454 first_vector_components,
2455 periodicity_factor);
2465 face_1->has_children() ?
2466 face_2->get_fe(face_2->nth_active_fe_index(0)) :
2467 face_1->get_fe(face_1->nth_active_fe_index(0));
2472 const unsigned int face_no = 0;
2478 if (n_dofs_per_face == 0)
2482 compute_transformation(fe, matrix, first_vector_components);
2484 if (!face_2->has_children())
2488 if (first_vector_components.empty() && matrix.m() == 0)
2490 set_periodicity_constraints(face_2,
2498 periodicity_factor);
2503 inverse.
invert(transformation);
2505 set_periodicity_constraints(face_2,
2513 periodicity_factor);
2526 set_periodicity_constraints(face_1,
2533 face_rotation ^ face_flip :
2536 periodicity_factor);
2543 template <
int dim,
int spacedim,
typename number>
2550 const std::vector<unsigned int> &first_vector_components,
2551 const number periodicity_factor)
2554 for (
auto &pair : periodic_faces)
2557 const FaceIterator face_1 = pair.cell[0]->face(pair.face_idx[0]);
2558 const FaceIterator face_2 = pair.cell[1]->face(pair.face_idx[1]);
2560 Assert(face_1->at_boundary() && face_2->at_boundary(),
2571 pair.orientation[0],
2572 pair.orientation[1],
2573 pair.orientation[2],
2575 first_vector_components,
2576 periodicity_factor);
2584 template <
int dim,
int spacedim,
typename number>
2589 const unsigned int direction,
2592 const number periodicity_factor)
2597 ExcMessage(
"The boundary indicators b_id1 and b_id2 must be "
2598 "different to denote different boundaries."));
2606 dof_handler, b_id1, b_id2, direction, matched_faces);
2608 make_periodicity_constraints<dim, spacedim>(matched_faces,
2611 std::vector<unsigned int>(),
2612 periodicity_factor);
2617 template <
int dim,
int spacedim,
typename number>
2621 const unsigned int direction,
2624 const number periodicity_factor)
2640 make_periodicity_constraints<dim, spacedim>(matched_faces,
2643 std::vector<unsigned int>(),
2644 periodicity_factor);
2657 template <
int dim,
int spacedim>
2662#ifdef DEAL_II_WITH_MPI
2663 std::vector<::LinearAlgebra::distributed::Vector<double>>
2678 template <
int dim,
int spacedim>
2680 compute_intergrid_weights_3(
2684 const unsigned int coarse_component,
2743 for (
unsigned int local_dof = 0; local_dof < copy_data.
dofs_per_cell;
2749 const unsigned int local_parameter_dof =
2758 coarse_to_fine_grid_map[cell]->set_dof_values_by_interpolation(
2759 parameter_dofs[local_parameter_dof],
2771 template <
int dim,
int spacedim>
2773 copy_intergrid_weights_3(
2774 const Assembler::CopyData<dim, spacedim> & copy_data,
2775 const unsigned int coarse_component,
2777 const std::vector<types::global_dof_index> &weight_mapping,
2778 const bool is_called_in_parallel,
2779 std::vector<std::map<types::global_dof_index, float>> &weights)
2781 unsigned int pos = 0;
2782 for (
unsigned int local_dof = 0; local_dof < copy_data.dofs_per_cell;
2809 i < copy_data.global_parameter_representation[pos].size();
2815 if (copy_data.global_parameter_representation[pos](i) != 0)
2818 wi = copy_data.parameter_dof_indices[local_dof],
2819 wj = weight_mapping[i];
2821 copy_data.global_parameter_representation[pos](i);
2824 else if (!is_called_in_parallel)
2829 Assert(copy_data.global_parameter_representation[pos](i) ==
2845 template <
int dim,
int spacedim>
2847 compute_intergrid_weights_2(
2849 const unsigned int coarse_component,
2852 const std::vector<types::global_dof_index> & weight_mapping,
2853 std::vector<std::map<types::global_dof_index, float>> &weights)
2855 Assembler::Scratch scratch;
2856 Assembler::CopyData<dim, spacedim> copy_data;
2858 unsigned int n_interesting_dofs = 0;
2859 for (
unsigned int local_dof = 0;
2860 local_dof < coarse_grid.
get_fe().n_dofs_per_cell();
2864 ++n_interesting_dofs;
2866 copy_data.global_parameter_representation.resize(n_interesting_dofs);
2868 bool is_called_in_parallel =
false;
2869 for (std::size_t i = 0;
2870 i < copy_data.global_parameter_representation.size();
2873#ifdef DEAL_II_WITH_MPI
2874 MPI_Comm communicator = MPI_COMM_SELF;
2877 const typename ::parallel::TriangulationBase<dim,
2879 &
tria =
dynamic_cast<const typename ::parallel::
2880 TriangulationBase<dim, spacedim> &
>(
2881 coarse_to_fine_grid_map.get_destination_grid()
2882 .get_triangulation());
2884 is_called_in_parallel =
true;
2886 catch (std::bad_cast &)
2894 coarse_to_fine_grid_map.get_destination_grid(),
2895 locally_relevant_dofs);
2897 copy_data.global_parameter_representation[i].reinit(
2898 coarse_to_fine_grid_map.get_destination_grid()
2899 .locally_owned_dofs(),
2900 locally_relevant_dofs,
2904 copy_data.global_parameter_representation[i].reinit(n_fine_dofs);
2911 &coarse_to_fine_grid_map,
2915 const Assembler::Scratch & scratch_data,
2916 Assembler::CopyData<dim, spacedim> ©_data) {
2917 compute_intergrid_weights_3<dim, spacedim>(cell,
2921 coarse_grid.get_fe(),
2922 coarse_to_fine_grid_map,
2930 is_called_in_parallel,
2931 &weights](
const Assembler::CopyData<dim, spacedim> ©_data) {
2932 copy_intergrid_weights_3<dim, spacedim>(copy_data,
2934 coarse_grid.get_fe(),
2936 is_called_in_parallel,
2947#ifdef DEAL_II_WITH_MPI
2948 for (std::size_t i = 0;
2949 i < copy_data.global_parameter_representation.size();
2951 copy_data.global_parameter_representation[i].update_ghost_values();
2962 template <
int dim,
int spacedim>
2964 compute_intergrid_weights_1(
2966 const unsigned int coarse_component,
2968 const unsigned int fine_component,
2970 std::vector<std::map<types::global_dof_index, float>> &weights,
2971 std::vector<types::global_dof_index> & weight_mapping)
2975 &fine_fe = fine_grid.
get_fe();
2979 n_fine_dofs = fine_grid.
n_dofs();
2982 const unsigned int fine_dofs_per_cell = fine_fe.n_dofs_per_cell();
2986 const unsigned int coarse_dofs_per_cell_component =
3000 Assert(&coarse_to_fine_grid_map.get_source_grid() == &coarse_grid,
3002 Assert(&coarse_to_fine_grid_map.get_destination_grid() == &fine_grid,
3014 fine_fe.component_to_base_index(fine_component).first),
3020 for (
const auto &cell : coarse_grid.active_cell_iterators())
3046 std::vector<::Vector<double>> parameter_dofs(
3047 coarse_dofs_per_cell_component,
3052 for (
unsigned int local_coarse_dof = 0;
3053 local_coarse_dof < coarse_dofs_per_cell_component;
3055 for (
unsigned int fine_dof = 0; fine_dof < fine_fe.n_dofs_per_cell();
3057 if (fine_fe.system_to_component_index(fine_dof) ==
3058 std::make_pair(fine_component, local_coarse_dof))
3060 parameter_dofs[local_coarse_dof](fine_dof) = 1.;
3067 unsigned int n_parameters_on_fine_grid = 0;
3071 std::vector<bool> dof_is_interesting(fine_grid.
n_dofs(),
false);
3072 std::vector<types::global_dof_index> local_dof_indices(
3073 fine_fe.n_dofs_per_cell());
3075 for (
const auto &cell : fine_grid.active_cell_iterators() |
3078 cell->get_dof_indices(local_dof_indices);
3079 for (
unsigned int i = 0; i < fine_fe.n_dofs_per_cell(); ++i)
3080 if (fine_fe.system_to_component_index(i).first ==
3082 dof_is_interesting[local_dof_indices[i]] =
true;
3085 n_parameters_on_fine_grid = std::count(dof_is_interesting.begin(),
3086 dof_is_interesting.end(),
3093 weights.resize(n_coarse_dofs);
3095 weight_mapping.clear();
3099 std::vector<types::global_dof_index> local_dof_indices(
3100 fine_fe.n_dofs_per_cell());
3101 unsigned int next_free_index = 0;
3102 for (
const auto &cell : fine_grid.active_cell_iterators() |
3105 cell->get_dof_indices(local_dof_indices);
3106 for (
unsigned int i = 0; i < fine_fe.n_dofs_per_cell(); ++i)
3109 if ((fine_fe.system_to_component_index(i).first ==
3111 (weight_mapping[local_dof_indices[i]] ==
3114 weight_mapping[local_dof_indices[i]] = next_free_index;
3119 Assert(next_free_index == n_parameters_on_fine_grid,
3131 compute_intergrid_weights_2(coarse_grid,
3133 coarse_to_fine_grid_map,
3153 for (
unsigned int col = 0; col < n_parameters_on_fine_grid; ++col)
3157 if (weights[row].find(col) != weights[row].
end())
3158 sum += weights[row][col];
3159 Assert((std::fabs(sum - 1) < 1.e-12) ||
3166 return n_parameters_on_fine_grid;
3175 template <
int dim,
int spacedim>
3179 const unsigned int coarse_component,
3181 const unsigned int fine_component,
3187 ExcMessage(
"This function is not yet implemented for DoFHandlers "
3188 "using hp-capabilities."));
3209 std::vector<std::map<types::global_dof_index, float>> weights;
3215 std::vector<types::global_dof_index> weight_mapping;
3217 const unsigned int n_parameters_on_fine_grid =
3218 internal::compute_intergrid_weights_1(coarse_grid,
3222 coarse_to_fine_grid_map,
3225 (void)n_parameters_on_fine_grid;
3229 n_fine_dofs = fine_grid.
n_dofs();
3237 mask[coarse_component] =
true;
3239 coarse_dof_is_parameter =
3240 extract_dofs<dim, spacedim>(coarse_grid,
ComponentMask(mask));
3252 std::vector<types::global_dof_index> representants(
3255 parameter_dof < n_coarse_dofs;
3257 if (coarse_dof_is_parameter.
is_element(parameter_dof))
3264 std::map<types::global_dof_index, float>::const_iterator i =
3265 weights[parameter_dof].begin();
3266 for (; i != weights[parameter_dof].end(); ++i)
3276 for (; global_dof < weight_mapping.size(); ++global_dof)
3277 if (weight_mapping[global_dof] ==
3283 representants[parameter_dof] = global_dof;
3299 std::vector<std::pair<types::global_dof_index, double>> constraint_line;
3318 std::map<types::global_dof_index, float>::const_iterator col_entry =
3320 for (; first_used_row < n_coarse_dofs; ++first_used_row)
3322 col_entry = weights[first_used_row].find(col);
3323 if (col_entry != weights[first_used_row].end())
3327 Assert(col_entry != weights[first_used_row].end(),
3330 if ((col_entry->second == 1) &&
3331 (representants[first_used_row] == global_dof))
3341 constraint_line.clear();
3343 row < n_coarse_dofs;
3346 const std::map<types::global_dof_index, float>::const_iterator j =
3347 weights[row].find(col);
3348 if ((j != weights[row].end()) && (j->second != 0))
3349 constraint_line.emplace_back(representants[row], j->second);
3352 constraints.
add_entries(global_dof, constraint_line);
3358 template <
int dim,
int spacedim>
3362 const unsigned int coarse_component,
3364 const unsigned int fine_component,
3366 std::vector<std::map<types::global_dof_index, float>>
3367 &transfer_representation)
3371 ExcMessage(
"This function is not yet implemented for DoFHandlers "
3372 "using hp-capabilities."));
3393 std::vector<std::map<types::global_dof_index, float>> weights;
3399 std::vector<types::global_dof_index> weight_mapping;
3401 internal::compute_intergrid_weights_1(coarse_grid,
3405 coarse_to_fine_grid_map,
3411 std::count_if(weight_mapping.begin(),
3412 weight_mapping.end(),
3414 return dof != numbers::invalid_dof_index;
3418 std::vector<types::global_dof_index> inverse_weight_mapping(
3427 Assert((inverse_weight_mapping[parameter_dof] ==
3431 inverse_weight_mapping[parameter_dof] = i;
3438 transfer_representation.clear();
3439 transfer_representation.resize(n_rows);
3444 std::map<types::global_dof_index, float>::const_iterator j =
3446 for (; j != weights[i].end(); ++j)
3451 transfer_representation[p][i] = j->second;
3458 template <
int dim,
int spacedim,
typename number>
3467 ExcMessage(
"The number of components in the mask has to be either "
3468 "zero or equal to the number of components in the finite "
3477 std::vector<types::global_dof_index> face_dofs;
3480 std::vector<types::global_dof_index> cell_dofs;
3486 for (; cell != endc; ++cell)
3487 if (!cell->is_artificial() && cell->at_boundary())
3493 cell->get_dof_indices(cell_dofs);
3495 for (
const auto face_no : cell->face_indices())
3498 cell->face(face_no);
3502 if (face->at_boundary() &&
3504 (face->boundary_id() == boundary_id)))
3508 face->get_dof_indices(face_dofs, cell->active_fe_index());
3516 const std::vector<types::global_dof_index>::iterator
3517 it_index_on_cell = std::find(cell_dofs.begin(),
3520 Assert(it_index_on_cell != cell_dofs.end(),
3522 const unsigned int index_on_cell =
3523 std::distance(cell_dofs.begin(), it_index_on_cell);
3527 for (
unsigned int c = 0; c < n_components; ++c)
3528 if (nonzero_component_array[c] && component_mask[c])
3535 zero_boundary_constraints.
add_line(face_dof);
3544 template <
int dim,
int spacedim,
typename number>
3553 zero_boundary_constraints,
3564#include "dof_tools_constraints.inst"
void add_entries(const size_type constrained_dof_index, const std::vector< std::pair< size_type, number > > &col_weight_pairs)
void add_line(const size_type line_n)
void add_entry(const size_type constrained_dof_index, const size_type column, const number weight)
void set_inhomogeneity(const size_type constrained_dof_index, const number value)
bool is_constrained(const size_type line_n) const
const std::vector< std::pair< size_type, number > > * get_constraint_entries(const size_type line_n) const
bool represents_n_components(const unsigned int n) const
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
cell_iterator end() const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
bool has_active_dofs() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
active_cell_iterator begin_active(const unsigned int level=0) const
bool has_hp_capabilities() const
types::global_dof_index n_dofs() const
unsigned int n_dofs_per_vertex() const
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_line() const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_components() const
unsigned int n_unique_faces() const
unsigned int n_dofs_per_quad(unsigned int face_no=0) const
const ComponentMask & get_nonzero_components(const unsigned int i) const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
std::pair< unsigned int, unsigned int > component_to_base_index(const unsigned int component) const
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
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::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
const FullMatrix< double > & constraints(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index, const unsigned int face_no=0) const
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
void invert(const FullMatrix< number2 > &M)
bool is_element(const size_type index) const
virtual MPI_Comm get_communicator() const
bool get_anisotropic_refinement_flag() const
unsigned int n_cells() const
unsigned int size() const
unsigned int find_dominating_fe_extended(const std::set< unsigned int > &fes, const unsigned int codim=0) const
bool hp_constraints_are_implemented() const
unsigned int max_dofs_per_face() const
unsigned int n_components() const
unsigned int max_dofs_per_cell() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcInvalidIterator()
static ::ExceptionBase & ExcGridNotCoarser()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcFiniteElementsDontMatch()
static ::ExceptionBase & ExcNoComponentSelected()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcGridsDontMatch()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename ActiveSelector::cell_iterator cell_iterator
typename ActiveSelector::line_iterator line_iterator
typename ActiveSelector::face_iterator face_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
void compute_intergrid_transfer_representation(const DoFHandler< dim, spacedim > &coarse_grid, const unsigned int coarse_component, const DoFHandler< dim, spacedim > &fine_grid, const unsigned int fine_component, const InterGridMap< DoFHandler< dim, spacedim > > &coarse_to_fine_grid_map, std::vector< std::map< types::global_dof_index, float > > &transfer_representation)
void make_zero_boundary_constraints(const DoFHandler< dim, spacedim > &dof, const types::boundary_id boundary_id, AffineConstraints< number > &zero_boundary_constraints, const ComponentMask &component_mask=ComponentMask())
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void compute_intergrid_constraints(const DoFHandler< dim, spacedim > &coarse_grid, const unsigned int coarse_component, const DoFHandler< dim, spacedim > &fine_grid, const unsigned int fine_component, const InterGridMap< DoFHandler< dim, spacedim > > &coarse_to_fine_grid_map, AffineConstraints< double > &constraints)
@ either_element_can_dominate
@ other_element_dominates
@ neither_element_dominates
@ matrix
Contents is actually a matrix.
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
VectorType::value_type * end(VectorType &V)
T sum(const T &t, const MPI_Comm mpi_communicator)
void run(const std::vector< std::vector< Iterator > > &colored_iterators, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
const types::boundary_id invalid_boundary_id
const types::fe_index invalid_fe_index
static const unsigned int invalid_unsigned_int
const types::global_dof_index invalid_dof_index
typename type_identity< T >::type type_identity_t
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned short int fe_index
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
const ::Triangulation< dim, spacedim > & tria