16 #include <deal.II/base/logstream.h> 17 #include <deal.II/base/mg_level_object.h> 18 #include <deal.II/base/thread_management.h> 20 #include <deal.II/dofs/dof_accessor.h> 21 #include <deal.II/dofs/dof_handler.h> 22 #include <deal.II/dofs/dof_tools.h> 24 #include <deal.II/fe/component_mask.h> 25 #include <deal.II/fe/fe.h> 27 #include <deal.II/grid/tria.h> 28 #include <deal.II/grid/tria_iterator.h> 30 #include <deal.II/lac/block_sparse_matrix.h> 31 #include <deal.II/lac/block_sparsity_pattern.h> 32 #include <deal.II/lac/block_vector.h> 33 #include <deal.II/lac/dynamic_sparsity_pattern.h> 34 #include <deal.II/lac/sparse_matrix.h> 35 #include <deal.II/lac/sparsity_pattern.h> 37 #include <deal.II/multigrid/mg_base.h> 38 #include <deal.II/multigrid/mg_constrained_dofs.h> 39 #include <deal.II/multigrid/mg_tools.h> 41 #include <deal.II/numerics/matrix_tools.h> 47 DEAL_II_NAMESPACE_OPEN
57 std::vector<unsigned int> &,
69 std::vector<unsigned int> &,
82 std::vector<unsigned int> &,
93 std::vector<unsigned int> &,
103 template <
int dim,
int spacedim>
106 const unsigned int level,
107 std::vector<unsigned int> & row_lengths,
115 std::fill(row_lengths.begin(), row_lengths.end(), 0);
118 std::vector<bool> old_flags;
131 std::vector<types::global_dof_index> cell_indices;
132 std::vector<types::global_dof_index> neighbor_indices;
141 for (cell = dofs.
begin(level); cell != end; ++cell)
145 cell->get_mg_dof_indices(cell_indices);
172 row_lengths[cell_indices[i++]] += increment;
183 increment = (dim > 1) ?
188 row_lengths[cell_indices[i++]] += increment;
191 increment = (dim > 2) ?
196 row_lengths[cell_indices[i++]] += increment;
201 row_lengths[cell_indices[i++]] += increment;
214 for (
unsigned int iface = 0; iface < GeometryInfo<dim>::faces_per_cell;
217 bool level_boundary = cell->at_boundary(iface);
221 neighbor = cell->neighbor(iface);
222 if (static_cast<unsigned int>(neighbor->level()) != level)
223 level_boundary =
true;
228 for (
unsigned int local_dof = 0; local_dof < fe.
dofs_per_cell;
248 const unsigned int dof_increment =
250 for (
unsigned int local_dof = 0; local_dof < fe.
dofs_per_cell;
252 row_lengths[cell_indices[local_dof]] += dof_increment;
257 if (face->user_flag_set())
259 face->set_user_flag();
273 neighbor->get_mg_dof_indices(neighbor_indices);
274 for (
unsigned int local_dof = 0; local_dof < fe.
dofs_per_cell;
277 for (
unsigned int local_dof = 0; local_dof < nfe.
dofs_per_cell;
279 row_lengths[neighbor_indices[local_dof]] += fe.
dofs_per_face;
287 template <
int dim,
int spacedim>
290 const unsigned int level,
291 std::vector<unsigned int> & row_lengths,
300 std::fill(row_lengths.begin(), row_lengths.end(), 0);
303 std::vector<bool> old_flags;
316 std::vector<types::global_dof_index> cell_indices;
317 std::vector<types::global_dof_index> neighbor_indices;
323 std::vector<Table<2, DoFTools::Coupling>> couple_cell;
324 std::vector<Table<2, DoFTools::Coupling>> couple_face;
338 const unsigned int fe_index = cell->active_fe_index();
352 cell->get_mg_dof_indices(cell_indices);
375 unsigned int increment;
387 row_lengths[cell_indices[i]] += increment;
412 ((dim > 1) ? (dim - 1) :
415 row_lengths[cell_indices[i]] += increment;
432 ((dim > 2) ? (dim - 2) :
435 row_lengths[cell_indices[i]] += increment;
453 row_lengths[cell_indices[i]] += increment;
469 for (
unsigned int iface = 0; iface < GeometryInfo<dim>::faces_per_cell;
472 bool level_boundary = cell->at_boundary(iface);
476 neighbor = cell->neighbor(iface);
477 if (static_cast<unsigned int>(neighbor->level()) != level)
478 level_boundary =
true;
483 for (
unsigned int local_dof = 0; local_dof < fe.
dofs_per_cell;
504 for (
unsigned int local_dof = 0; local_dof < fe.
dofs_per_cell;
506 if (couple_face[fe_index](
510 const unsigned int dof_increment =
513 row_lengths[cell_indices[local_dof]] += dof_increment;
519 if (face->user_flag_set())
521 face->set_user_flag();
543 neighbor->get_mg_dof_indices(neighbor_indices);
547 for (
unsigned int local_dof = 0; local_dof < fe.
dofs_per_cell;
549 if (couple_cell[fe_index](
552 row_lengths[cell_indices[local_dof]] +=
557 for (
unsigned int local_dof = 0; local_dof < nfe.
dofs_per_cell;
559 if (couple_cell[fe_index](
562 row_lengths[neighbor_indices[local_dof]] +=
571 template <
typename DoFHandlerType,
typename SparsityPatternType>
574 SparsityPatternType & sparsity,
575 const unsigned int level)
580 Assert(sparsity.n_rows() == n_dofs,
582 Assert(sparsity.n_cols() == n_dofs,
585 const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
586 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
587 typename DoFHandlerType::cell_iterator cell = dof.begin(level),
588 endc = dof.end(level);
589 for (; cell != endc; ++cell)
590 if (dof.get_triangulation().locally_owned_subdomain() ==
592 cell->level_subdomain_id() ==
593 dof.get_triangulation().locally_owned_subdomain())
595 cell->get_mg_dof_indices(dofs_on_this_cell);
597 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
598 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
599 sparsity.add(dofs_on_this_cell[i], dofs_on_this_cell[j]);
605 template <
int dim,
typename SparsityPatternType,
int spacedim>
608 SparsityPatternType & sparsity,
609 const unsigned int level)
614 Assert(sparsity.n_rows() == n_dofs,
616 Assert(sparsity.n_cols() == n_dofs,
619 const unsigned int dofs_per_cell = dof.
get_fe().dofs_per_cell;
620 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
621 std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
623 endc = dof.
end(level);
624 for (; cell != endc; ++cell)
626 if (!cell->is_locally_owned_on_level())
629 cell->get_mg_dof_indices(dofs_on_this_cell);
631 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
632 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
633 sparsity.add(dofs_on_this_cell[i], dofs_on_this_cell[j]);
636 for (
unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell;
639 bool use_face =
false;
640 if ((!cell->at_boundary(face)) &&
641 (static_cast<unsigned int>(cell->neighbor_level(face)) ==
644 else if (cell->has_periodic_neighbor(face) &&
645 (
static_cast<unsigned int>(
646 cell->periodic_neighbor_level(face)) == level))
652 cell->neighbor_or_periodic_neighbor(face);
653 neighbor->get_mg_dof_indices(dofs_on_other_cell);
657 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
659 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
661 sparsity.add(dofs_on_this_cell[i],
662 dofs_on_other_cell[j]);
665 if (neighbor->is_locally_owned_on_level() ==
false)
666 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
667 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
669 sparsity.add(dofs_on_other_cell[i],
670 dofs_on_other_cell[j]);
671 sparsity.add(dofs_on_other_cell[i],
672 dofs_on_this_cell[j]);
681 template <
int dim,
typename SparsityPatternType,
int spacedim>
684 SparsityPatternType & sparsity,
685 const unsigned int level)
697 Assert(sparsity.n_rows() == coarse_dofs,
699 Assert(sparsity.n_cols() == fine_dofs,
702 const unsigned int dofs_per_cell = dof.
get_fe().dofs_per_cell;
703 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
704 std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
706 endc = dof.
end(level);
707 for (; cell != endc; ++cell)
709 if (!cell->is_locally_owned_on_level())
712 cell->get_mg_dof_indices(dofs_on_this_cell);
714 for (
unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell;
718 bool use_face =
false;
719 if ((!cell->at_boundary(face)) &&
720 (static_cast<unsigned int>(cell->neighbor_level(face)) !=
723 else if (cell->has_periodic_neighbor(face) &&
724 (
static_cast<unsigned int>(
725 cell->periodic_neighbor_level(face)) != level))
731 cell->neighbor_or_periodic_neighbor(face);
732 neighbor->get_mg_dof_indices(dofs_on_other_cell);
734 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
736 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
738 sparsity.add(dofs_on_other_cell[i],
739 dofs_on_this_cell[j]);
740 sparsity.add(dofs_on_other_cell[j],
741 dofs_on_this_cell[i]);
751 template <
int dim,
typename SparsityPatternType,
int spacedim>
754 SparsityPatternType & sparsity,
755 const unsigned int level,
765 Assert(sparsity.n_rows() == n_dofs,
767 Assert(sparsity.n_cols() == n_dofs,
769 Assert(int_mask.n_rows() == n_comp,
771 Assert(int_mask.n_cols() == n_comp,
773 Assert(flux_mask.n_rows() == n_comp,
775 Assert(flux_mask.n_cols() == n_comp,
779 std::vector<types::global_dof_index> dofs_on_this_cell(total_dofs);
780 std::vector<types::global_dof_index> dofs_on_other_cell(total_dofs);
785 endc = dof.
end(level);
793 for (
unsigned int i = 0; i < total_dofs; ++i)
794 for (
unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
805 std::vector<bool> user_flags;
810 for (; cell != endc; ++cell)
812 if (!cell->is_locally_owned_on_level())
815 cell->get_mg_dof_indices(dofs_on_this_cell);
817 for (
unsigned int i = 0; i < total_dofs; ++i)
818 for (
unsigned int j = 0; j < total_dofs; ++j)
820 sparsity.add(dofs_on_this_cell[i], dofs_on_this_cell[j]);
823 for (
unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell;
828 if (cell_face->user_flag_set())
831 if (cell->at_boundary(face) && !cell->has_periodic_neighbor(face))
833 for (
unsigned int i = 0; i < total_dofs; ++i)
835 const bool i_non_zero_i = support_on_face(i, face);
836 for (
unsigned int j = 0; j < total_dofs; ++j)
838 const bool j_non_zero_i = support_on_face(j, face);
841 sparsity.add(dofs_on_this_cell[i],
842 dofs_on_this_cell[j]);
844 i_non_zero_i && j_non_zero_i)
845 sparsity.add(dofs_on_this_cell[i],
846 dofs_on_this_cell[j]);
853 cell->neighbor_or_periodic_neighbor(face);
855 if (neighbor->level() < cell->level())
858 unsigned int neighbor_face =
859 cell->has_periodic_neighbor(face) ?
860 cell->periodic_neighbor_of_periodic_neighbor(face) :
861 cell->neighbor_of_neighbor(face);
863 neighbor->get_mg_dof_indices(dofs_on_other_cell);
864 for (
unsigned int i = 0; i < total_dofs; ++i)
866 const bool i_non_zero_i = support_on_face(i, face);
867 const bool i_non_zero_e = support_on_face(i, neighbor_face);
868 for (
unsigned int j = 0; j < total_dofs; ++j)
870 const bool j_non_zero_i = support_on_face(j, face);
871 const bool j_non_zero_e =
872 support_on_face(j, neighbor_face);
875 sparsity.add(dofs_on_this_cell[i],
876 dofs_on_other_cell[j]);
877 sparsity.add(dofs_on_other_cell[i],
878 dofs_on_this_cell[j]);
879 sparsity.add(dofs_on_this_cell[i],
880 dofs_on_this_cell[j]);
881 sparsity.add(dofs_on_other_cell[i],
882 dofs_on_other_cell[j]);
886 if (i_non_zero_i && j_non_zero_e)
887 sparsity.add(dofs_on_this_cell[i],
888 dofs_on_other_cell[j]);
889 if (i_non_zero_e && j_non_zero_i)
890 sparsity.add(dofs_on_other_cell[i],
891 dofs_on_this_cell[j]);
892 if (i_non_zero_i && j_non_zero_i)
893 sparsity.add(dofs_on_this_cell[i],
894 dofs_on_this_cell[j]);
895 if (i_non_zero_e && j_non_zero_e)
896 sparsity.add(dofs_on_other_cell[i],
897 dofs_on_other_cell[j]);
902 sparsity.add(dofs_on_this_cell[j],
903 dofs_on_other_cell[i]);
904 sparsity.add(dofs_on_other_cell[j],
905 dofs_on_this_cell[i]);
906 sparsity.add(dofs_on_this_cell[j],
907 dofs_on_this_cell[i]);
908 sparsity.add(dofs_on_other_cell[j],
909 dofs_on_other_cell[i]);
913 if (j_non_zero_i && i_non_zero_e)
914 sparsity.add(dofs_on_this_cell[j],
915 dofs_on_other_cell[i]);
916 if (j_non_zero_e && i_non_zero_i)
917 sparsity.add(dofs_on_other_cell[j],
918 dofs_on_this_cell[i]);
919 if (j_non_zero_i && i_non_zero_i)
920 sparsity.add(dofs_on_this_cell[j],
921 dofs_on_this_cell[i]);
922 if (j_non_zero_e && i_non_zero_e)
923 sparsity.add(dofs_on_other_cell[j],
924 dofs_on_other_cell[i]);
928 neighbor->face(neighbor_face)->set_user_flag();
935 .load_user_flags(user_flags);
940 template <
int dim,
typename SparsityPatternType,
int spacedim>
943 SparsityPatternType & sparsity,
944 const unsigned int level,
961 Assert(sparsity.n_rows() == coarse_dofs,
963 Assert(sparsity.n_cols() == fine_dofs,
965 Assert(flux_mask.n_rows() == n_comp,
967 Assert(flux_mask.n_cols() == n_comp,
970 const unsigned int dofs_per_cell = dof.
get_fe().dofs_per_cell;
971 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
972 std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
977 endc = dof.
end(level);
982 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
983 for (
unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
986 for (; cell != endc; ++cell)
988 if (!cell->is_locally_owned_on_level())
991 cell->get_mg_dof_indices(dofs_on_this_cell);
993 for (
unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell;
997 bool use_face =
false;
998 if ((!cell->at_boundary(face)) &&
999 (static_cast<unsigned int>(cell->neighbor_level(face)) !=
1002 else if (cell->has_periodic_neighbor(face) &&
1003 (
static_cast<unsigned int>(
1004 cell->periodic_neighbor_level(face)) != level))
1010 cell->neighbor_or_periodic_neighbor(face);
1011 neighbor->get_mg_dof_indices(dofs_on_other_cell);
1013 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1015 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
1019 sparsity.add(dofs_on_other_cell[i],
1020 dofs_on_this_cell[j]);
1021 sparsity.add(dofs_on_other_cell[j],
1022 dofs_on_this_cell[i]);
1033 template <
typename DoFHandlerType,
typename SparsityPatternType>
1037 SparsityPatternType & sparsity,
1038 const unsigned int level)
1043 Assert(sparsity.n_rows() == n_dofs,
1045 Assert(sparsity.n_cols() == n_dofs,
1048 const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
1049 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
1050 typename DoFHandlerType::cell_iterator cell = dof.begin(level),
1051 endc = dof.end(level);
1052 for (; cell != endc; ++cell)
1053 if (cell->is_locally_owned_on_level())
1055 cell->get_mg_dof_indices(dofs_on_this_cell);
1056 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1057 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
1059 level, dofs_on_this_cell[i], dofs_on_this_cell[j]))
1060 sparsity.add(dofs_on_this_cell[i], dofs_on_this_cell[j]);
1066 template <
int dim,
int spacedim>
1070 std::vector<std::vector<types::global_dof_index>> &result,
1072 std::vector<unsigned int> target_component)
1076 const unsigned int nlevels =
1079 Assert(result.size() == nlevels,
1082 if (target_component.size() == 0)
1084 target_component.resize(n_components);
1085 for (
unsigned int i = 0; i < n_components; ++i)
1086 target_component[i] = i;
1089 Assert(target_component.size() == n_components,
1092 for (
unsigned int l = 0; l < nlevels; ++l)
1094 result[l].resize(n_components);
1095 std::fill(result[l].begin(), result[l].end(), 0U);
1101 if (n_components == 1)
1103 result[l][0] = dof_handler.
n_dofs(l);
1110 std::vector<std::vector<bool>> dofs_in_component(
1111 n_components, std::vector<bool>(dof_handler.
n_dofs(l),
false));
1112 std::vector<ComponentMask> component_select(n_components);
1114 for (
unsigned int i = 0; i < n_components; ++i)
1116 void (*fun_ptr)(
const unsigned int level,
1119 std::vector<bool> &) =
1122 std::vector<bool> tmp(n_components,
false);
1129 component_select[i],
1130 dofs_in_component[i]);
1135 unsigned int component = 0;
1144 for (
unsigned int dd = 0; dd < d; ++dd)
1147 result[l][target_component[component]] +=
1148 std::count(dofs_in_component[component].begin(),
1149 dofs_in_component[component].end(),
1157 std::accumulate(result[l].begin(),
1168 template <
typename DoFHandlerType>
1171 const DoFHandlerType & dof_handler,
1172 std::vector<std::vector<types::global_dof_index>> &dofs_per_block,
1173 std::vector<unsigned int> target_block)
1176 DoFHandlerType::space_dimension> &fe =
1177 dof_handler.get_fe();
1178 const unsigned int n_blocks = fe.n_blocks();
1179 const unsigned int n_levels =
1180 dof_handler.get_triangulation().n_global_levels();
1184 for (
unsigned int l = 0; l < n_levels; ++l)
1185 std::fill(dofs_per_block[l].begin(), dofs_per_block[l].end(), 0U);
1189 if (target_block.size() == 0)
1191 target_block.resize(n_blocks);
1192 for (
unsigned int i = 0; i < n_blocks; ++i)
1193 target_block[i] = i;
1195 Assert(target_block.size() == n_blocks,
1198 const unsigned int max_block =
1199 *std::max_element(target_block.begin(), target_block.end());
1200 const unsigned int n_target_blocks = max_block + 1;
1201 (void)n_target_blocks;
1203 for (
unsigned int l = 0; l < n_levels; ++l)
1212 for (
unsigned int l = 0; l < n_levels; ++l)
1213 dofs_per_block[l][0] = dof_handler.n_dofs(l);
1219 for (
unsigned int l = 0; l < n_levels; ++l)
1221 std::vector<std::vector<bool>> dofs_in_block(
1222 n_blocks, std::vector<bool>(dof_handler.n_dofs(l),
false));
1223 std::vector<BlockMask> block_select(n_blocks);
1225 for (
unsigned int i = 0; i < n_blocks; ++i)
1227 void (*fun_ptr)(
const unsigned int level,
1228 const DoFHandlerType &,
1230 std::vector<bool> &) =
1231 &DoFTools::extract_level_dofs<DoFHandlerType>;
1233 std::vector<bool> tmp(n_blocks,
false);
1235 block_select[i] = tmp;
1238 fun_ptr, l, dof_handler, block_select[i], dofs_in_block[i]);
1243 for (
unsigned int block = 0; block < fe.n_blocks(); ++block)
1244 dofs_per_block[l][target_block[block]] +=
1245 std::count(dofs_in_block[block].begin(),
1246 dofs_in_block[block].end(),
1253 template <
int dim,
int spacedim>
1259 std::vector<std::set<types::global_dof_index>> &boundary_indices,
1266 std::set<types::boundary_id> boundary_ids;
1267 for (
const auto &boundary_function : function_map)
1268 boundary_ids.insert(boundary_function.first);
1270 std::vector<IndexSet> boundary_indexset;
1273 boundary_indices[i].insert(boundary_indexset[i].begin(),
1274 boundary_indexset[i].end());
1278 template <
int dim,
int spacedim>
1283 std::vector<IndexSet> &boundary_indices,
1290 std::set<types::boundary_id> boundary_ids;
1291 for (
const auto &boundary_function : function_map)
1292 boundary_ids.insert(boundary_function.first);
1299 template <
int dim,
int spacedim>
1302 const std::set<types::boundary_id> &boundary_ids,
1303 std::vector<IndexSet> & boundary_indices,
1309 if (boundary_ids.size() == 0)
1313 if (boundary_indices[i].size() == 0)
1317 const bool fe_is_system = (n_components != 1);
1319 std::vector<types::global_dof_index> local_dofs;
1329 for (; cell != endc; ++cell)
1336 const unsigned int level = cell->level();
1339 for (
unsigned int face_no = 0;
1340 face_no < GeometryInfo<dim>::faces_per_cell;
1342 if (cell->at_boundary(face_no) ==
true)
1345 cell->face(face_no);
1348 if (boundary_ids.find(bi) != boundary_ids.end())
1350 face->get_mg_dof_indices(level, local_dofs);
1351 boundary_indices[level].add_indices(local_dofs.begin(),
1361 "It's probably worthwhile to select at least one component."));
1365 for (; cell != endc; ++cell)
1369 for (
unsigned int face_no = 0;
1370 face_no < GeometryInfo<dim>::faces_per_cell;
1373 if (cell->at_boundary(face_no) ==
false)
1377 const unsigned int level = cell->level();
1380 cell->face(face_no);
1382 face->boundary_id();
1383 if (boundary_ids.find(boundary_component) != boundary_ids.end())
1386 for (
unsigned int i = 0; i < cell->get_fe().dofs_per_cell;
1390 cell->get_fe().get_nonzero_components(i);
1394 bool selected =
false;
1395 for (
unsigned int c = 0; c < n_components; ++c)
1396 if (nonzero_component_array[c] ==
true &&
1397 component_mask[c] ==
true)
1403 for (
unsigned int c = 0; c < n_components; ++c)
1405 nonzero_component_array[c] ==
false ||
1406 component_mask[c] ==
true,
1408 "You are using a non-primitive FiniteElement " 1409 "and try to constrain just some of its components!"));
1415 face->get_mg_dof_indices(level, local_dofs);
1418 for (
unsigned int i = 0; i < local_dofs.size(); ++i)
1420 unsigned int component =
1431 cell->get_fe().get_nonzero_components(i);
1432 for (
unsigned int c = 0; c < n_components; ++c)
1433 if (nonzero_component_array[c] ==
true)
1441 if (component_mask[component] ==
true)
1442 boundary_indices[level].add_index(local_dofs[i]);
1446 boundary_indices[level].add_indices(local_dofs.begin(),
1454 template <
int dim,
int spacedim>
1458 std::vector<std::set<types::global_dof_index>> &non_interface_dofs)
1460 Assert(non_interface_dofs.size() ==
1463 non_interface_dofs.size(),
1471 std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1472 std::vector<bool> cell_dofs(dofs_per_cell,
false);
1473 std::vector<bool> cell_dofs_interface(dofs_per_cell,
false);
1476 endc = mg_dof_handler.
end();
1479 for (; cell != endc; ++cell)
1483 cell->level_subdomain_id() !=
1487 std::fill(cell_dofs.begin(), cell_dofs.end(),
false);
1488 std::fill(cell_dofs_interface.begin(),
1489 cell_dofs_interface.end(),
1492 for (
unsigned int face_nr = 0;
1493 face_nr < GeometryInfo<dim>::faces_per_cell;
1497 cell->face(face_nr);
1498 if (!face->at_boundary() || cell->has_periodic_neighbor(face_nr))
1502 cell->neighbor_or_periodic_neighbor(face_nr);
1504 if ((neighbor->level() < cell->level()))
1506 for (
unsigned int j = 0; j < dofs_per_face; ++j)
1512 for (
unsigned int j = 0; j < dofs_per_face; ++j)
1519 for (
unsigned int j = 0; j < dofs_per_face; ++j)
1524 const unsigned int level = cell->level();
1525 cell->get_mg_dof_indices(local_dof_indices);
1527 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1528 if (cell_dofs[i] && !cell_dofs_interface[i])
1529 non_interface_dofs[level].insert(local_dof_indices[i]);
1534 template <
int dim,
int spacedim>
1537 std::vector<IndexSet> & interface_dofs)
1539 Assert(interface_dofs.size() ==
1542 interface_dofs.size(),
1545 std::vector<std::vector<types::global_dof_index>> tmp_interface_dofs(
1546 interface_dofs.size());
1550 const unsigned int dofs_per_cell = fe.dofs_per_cell;
1551 const unsigned int dofs_per_face = fe.dofs_per_face;
1553 std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1555 std::vector<bool> cell_dofs(dofs_per_cell,
false);
1558 endc = mg_dof_handler.
end();
1560 for (; cell != endc; ++cell)
1569 bool has_coarser_neighbor =
false;
1571 std::fill(cell_dofs.begin(), cell_dofs.end(),
false);
1573 for (
unsigned int face_nr = 0;
1574 face_nr < GeometryInfo<dim>::faces_per_cell;
1578 cell->face(face_nr);
1579 if (!face->at_boundary() || cell->has_periodic_neighbor(face_nr))
1583 cell->neighbor_or_periodic_neighbor(face_nr);
1588 .locally_owned_subdomain() !=
1590 neighbor->level_subdomain_id() ==
1595 if (neighbor->level() < cell->level())
1597 for (
unsigned int j = 0; j < dofs_per_face; ++j)
1598 cell_dofs[fe.face_to_cell_index(j, face_nr)] =
true;
1600 has_coarser_neighbor =
true;
1605 if (has_coarser_neighbor ==
false)
1608 const unsigned int level = cell->level();
1609 cell->get_mg_dof_indices(local_dof_indices);
1611 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1614 tmp_interface_dofs[level].push_back(local_dof_indices[i]);
1618 for (
unsigned int l = 0;
1622 interface_dofs[l].clear();
1623 std::sort(tmp_interface_dofs[l].begin(), tmp_interface_dofs[l].end());
1624 interface_dofs[l].add_indices(tmp_interface_dofs[l].begin(),
1625 std::unique(tmp_interface_dofs[l].begin(),
1626 tmp_interface_dofs[l].end()));
1627 interface_dofs[l].compress();
1633 template <
int dim,
int spacedim>
1646 for (; cell != endc; ++cell)
1647 if (cell->is_locally_owned())
1649 min_level = cell->level();
1653 unsigned int global_min = min_level;
1668 #include "mg_tools.inst" 1670 DEAL_II_NAMESPACE_CLOSE
const unsigned int first_hex_index
const Triangulation< dim, spacedim > & get_triangulation() const
static const unsigned int invalid_unsigned_int
#define AssertDimension(dim1, dim2)
const types::subdomain_id invalid_subdomain_id
cell_iterator begin(const unsigned int level=0) const
Task< RT > new_task(const std::function< RT()> &function)
void load_user_flags(std::istream &in)
cell_iterator end() const
#define AssertIndexRange(index, range)
active_cell_iterator begin_active(const unsigned int level=0) const
bool is_primitive() const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
typename ActiveSelector::face_iterator face_iterator
bool is_interface_matrix_entry(const unsigned int level, const types::global_dof_index i, const types::global_dof_index j) const
cell_iterator end() const
active_cell_iterator begin_active(const unsigned int level=0) const
static ::ExceptionBase & ExcMessage(std::string arg1)
const unsigned int first_quad_index
#define Assert(cond, exc)
unsigned int element_multiplicity(const unsigned int index) const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index) const
types::global_dof_index n_dofs() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) 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
const unsigned int dofs_per_cell
types::global_dof_index first_block_of_base(const unsigned int b) const
unsigned int global_dof_index
const types::subdomain_id artificial_subdomain_id
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
unsigned int n_components() const
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
virtual unsigned int n_global_levels() const
T min(const T &t, const MPI_Comm &mpi_communicator)
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
const unsigned int dofs_per_face
const unsigned int first_line_index
void save_user_flags(std::ostream &out) const
static ::ExceptionBase & ExcNotImplemented()
typename ActiveSelector::cell_iterator cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
const types::global_dof_index invalid_dof_index
unsigned int n_base_elements() const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
static ::ExceptionBase & ExcInternalError()