16 #include <deal.II/base/logstream.h> 17 #include <deal.II/base/thread_management.h> 18 #include <deal.II/lac/sparsity_pattern.h> 19 #include <deal.II/lac/block_sparsity_pattern.h> 20 #include <deal.II/lac/dynamic_sparsity_pattern.h> 21 #include <deal.II/lac/sparsity_pattern.h> 22 #include <deal.II/lac/block_vector.h> 23 #include <deal.II/lac/sparse_matrix.h> 24 #include <deal.II/lac/block_sparse_matrix.h> 25 #include <deal.II/lac/block_vector.h> 26 #include <deal.II/grid/tria.h> 27 #include <deal.II/grid/tria_iterator.h> 28 #include <deal.II/dofs/dof_accessor.h> 29 #include <deal.II/multigrid/mg_tools.h> 30 #include <deal.II/multigrid/mg_base.h> 31 #include <deal.II/base/mg_level_object.h> 32 #include <deal.II/dofs/dof_handler.h> 33 #include <deal.II/dofs/dof_tools.h> 34 #include <deal.II/fe/fe.h> 35 #include <deal.II/fe/component_mask.h> 36 #include <deal.II/numerics/matrix_tools.h> 42 DEAL_II_NAMESPACE_OPEN
55 std::vector<unsigned int> &,
68 std::vector<unsigned int> &,
82 std::vector<unsigned int> &,
94 std::vector<unsigned int> &,
104 template <
int dim,
int spacedim>
108 const unsigned int level,
109 std::vector<unsigned int> &row_lengths,
112 Assert (row_lengths.size() == dofs.n_dofs(),
117 std::fill(row_lengths.begin(), row_lengths.end(), 0);
120 std::vector<bool> old_flags;
132 std::vector<types::global_dof_index> cell_indices;
133 std::vector<types::global_dof_index> neighbor_indices;
142 for (cell = dofs.begin(level); cell != end; ++cell)
146 cell->get_mg_dof_indices(cell_indices);
172 row_lengths[cell_indices[i++]] += increment;
187 row_lengths[cell_indices[i++]] += increment;
194 row_lengths[cell_indices[i++]] += increment;
198 row_lengths[cell_indices[i++]] += increment;
211 for (
unsigned int iface=0; iface<GeometryInfo<dim>::faces_per_cell; ++iface)
213 bool level_boundary = cell->at_boundary(iface);
217 neighbor = cell->neighbor(iface);
218 if (static_cast<unsigned int>(neighbor->level()) != level)
219 level_boundary =
true;
224 for (
unsigned int local_dof=0; local_dof<fe.
dofs_per_cell; ++local_dof)
230 typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(iface);
243 for (
unsigned int local_dof=0; local_dof<fe.
dofs_per_cell; ++local_dof)
244 row_lengths[cell_indices[local_dof]] += dof_increment;
249 if (face->user_flag_set())
251 face->set_user_flag();
265 neighbor->get_mg_dof_indices(neighbor_indices);
266 for (
unsigned int local_dof=0; local_dof<fe.
dofs_per_cell; ++local_dof)
268 for (
unsigned int local_dof=0; local_dof<nfe.
dofs_per_cell; ++local_dof)
269 row_lengths[neighbor_indices[local_dof]] += fe.
dofs_per_face;
277 template <
int dim,
int spacedim>
281 const unsigned int level,
282 std::vector<unsigned int> &row_lengths,
286 Assert (row_lengths.size() == dofs.n_dofs(),
291 std::fill(row_lengths.begin(), row_lengths.end(), 0);
294 std::vector<bool> old_flags;
306 std::vector<types::global_dof_index> cell_indices;
307 std::vector<types::global_dof_index> neighbor_indices;
313 std::vector<Table<2, DoFTools::Coupling> > couple_cell;
314 std::vector<Table<2, DoFTools::Coupling> > couple_face;
325 for (cell = dofs.begin_active(); cell != end; ++cell)
328 const unsigned int fe_index = cell->active_fe_index();
340 cell->get_mg_dof_indices(cell_indices);
363 unsigned int increment;
373 row_lengths[cell_indices[i]] += increment;
399 row_lengths[cell_indices[i]] += increment;
417 row_lengths[cell_indices[i]] += increment;
433 row_lengths[cell_indices[i]] += increment;
449 for (
unsigned int iface=0; iface<GeometryInfo<dim>::faces_per_cell; ++iface)
451 bool level_boundary = cell->at_boundary(iface);
455 neighbor = cell->neighbor(iface);
456 if (static_cast<unsigned int>(neighbor->level()) != level)
457 level_boundary =
true;
462 for (
unsigned int local_dof=0; local_dof<fe.
dofs_per_cell; ++local_dof)
468 typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(iface);
480 for (
unsigned int local_dof=0; local_dof<fe.
dofs_per_cell; ++local_dof)
486 row_lengths[cell_indices[local_dof]] += dof_increment;
492 if (face->user_flag_set())
494 face->set_user_flag();
516 neighbor->get_mg_dof_indices(neighbor_indices);
519 for (
unsigned int local_dof=0; local_dof<fe.
dofs_per_cell; ++local_dof)
522 row_lengths[cell_indices[local_dof]]
526 for (
unsigned int local_dof=0; local_dof<nfe.
dofs_per_cell; ++local_dof)
529 row_lengths[neighbor_indices[local_dof]]
538 template <
typename DoFHandlerType,
typename SparsityPatternType>
540 SparsityPatternType &sparsity,
541 const unsigned int level)
546 Assert (sparsity.n_rows() == n_dofs,
548 Assert (sparsity.n_cols() == n_dofs,
551 const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
552 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
553 typename DoFHandlerType::cell_iterator cell = dof.begin(level),
554 endc = dof.end(level);
555 for (; cell!=endc; ++cell)
557 || cell->level_subdomain_id()==dof.get_triangulation().locally_owned_subdomain())
559 cell->get_mg_dof_indices (dofs_on_this_cell);
561 for (
unsigned int i=0; i<dofs_per_cell; ++i)
562 for (
unsigned int j=0; j<dofs_per_cell; ++j)
563 sparsity.add (dofs_on_this_cell[i],
564 dofs_on_this_cell[j]);
570 template <
int dim,
typename SparsityPatternType,
int spacedim>
573 SparsityPatternType &sparsity,
574 const unsigned int level)
579 Assert (sparsity.n_rows() == n_dofs,
581 Assert (sparsity.n_cols() == n_dofs,
584 const unsigned int dofs_per_cell = dof.
get_fe().dofs_per_cell;
585 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
586 std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
588 endc = dof.
end(level);
589 for (; cell!=endc; ++cell)
591 if (!cell->is_locally_owned_on_level())
continue;
593 cell->get_mg_dof_indices (dofs_on_this_cell);
595 for (
unsigned int i=0; i<dofs_per_cell; ++i)
596 for (
unsigned int j=0; j<dofs_per_cell; ++j)
597 sparsity.add (dofs_on_this_cell[i],
598 dofs_on_this_cell[j]);
601 for (
unsigned int face = 0;
602 face < GeometryInfo<dim>::faces_per_cell;
605 bool use_face =
false;
606 if ((! cell->at_boundary(face)) &&
607 (static_cast<unsigned int>(cell->neighbor_level(face)) == level))
609 else if (cell->has_periodic_neighbor(face) &&
610 (
static_cast<unsigned int>(cell->periodic_neighbor_level(face)) == level))
616 neighbor = cell->neighbor_or_periodic_neighbor(face);
617 neighbor->get_mg_dof_indices (dofs_on_other_cell);
621 for (
unsigned int i=0; i<dofs_per_cell; ++i)
623 for (
unsigned int j=0; j<dofs_per_cell; ++j)
625 sparsity.add (dofs_on_this_cell[i],
626 dofs_on_other_cell[j]);
629 if (neighbor->is_locally_owned_on_level() ==
false)
630 for (
unsigned int i=0; i<dofs_per_cell; ++i)
631 for (
unsigned int j=0; j<dofs_per_cell; ++j)
633 sparsity.add (dofs_on_other_cell[i],
634 dofs_on_other_cell[j]);
635 sparsity.add (dofs_on_other_cell[i],
636 dofs_on_this_cell[j]);
645 template <
int dim,
typename SparsityPatternType,
int spacedim>
648 SparsityPatternType &sparsity,
649 const unsigned int level)
661 Assert (sparsity.n_rows() == coarse_dofs,
663 Assert (sparsity.n_cols() == fine_dofs,
666 const unsigned int dofs_per_cell = dof.
get_fe().dofs_per_cell;
667 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
668 std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
670 endc = dof.
end(level);
671 for (; cell!=endc; ++cell)
673 if (!cell->is_locally_owned_on_level())
continue;
675 cell->get_mg_dof_indices (dofs_on_this_cell);
677 for (
unsigned int face = 0;
678 face < GeometryInfo<dim>::faces_per_cell;
682 bool use_face =
false;
683 if ((! cell->at_boundary(face)) &&
684 (static_cast<unsigned int>(cell->neighbor_level(face)) != level))
686 else if (cell->has_periodic_neighbor(face) &&
687 (
static_cast<unsigned int>(cell->periodic_neighbor_level(face)) != level))
693 neighbor = cell->neighbor_or_periodic_neighbor(face);
694 neighbor->get_mg_dof_indices (dofs_on_other_cell);
696 for (
unsigned int i=0; i<dofs_per_cell; ++i)
698 for (
unsigned int j=0; j<dofs_per_cell; ++j)
700 sparsity.add (dofs_on_other_cell[i],
701 dofs_on_this_cell[j]);
702 sparsity.add (dofs_on_other_cell[j],
703 dofs_on_this_cell[i]);
713 template <
int dim,
typename SparsityPatternType,
int spacedim>
716 SparsityPatternType &sparsity,
717 const unsigned int level,
727 Assert (sparsity.n_rows() == n_dofs,
729 Assert (sparsity.n_cols() == n_dofs,
731 Assert (int_mask.n_rows() == n_comp,
733 Assert (int_mask.n_cols() == n_comp,
735 Assert (flux_mask.n_rows() == n_comp,
737 Assert (flux_mask.n_cols() == n_comp,
741 std::vector<types::global_dof_index> dofs_on_this_cell(total_dofs);
742 std::vector<types::global_dof_index> dofs_on_other_cell(total_dofs);
746 endc = dof.
end(level);
752 for (
unsigned int i=0; i<total_dofs; ++i)
753 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
764 std::vector<bool> user_flags;
768 for (; cell!=endc; ++cell)
770 if (!cell->is_locally_owned_on_level())
continue;
772 cell->get_mg_dof_indices (dofs_on_this_cell);
774 for (
unsigned int i=0; i<total_dofs; ++i)
775 for (
unsigned int j=0; j<total_dofs; ++j)
777 sparsity.add (dofs_on_this_cell[i],
778 dofs_on_this_cell[j]);
781 for (
unsigned int face = 0;
782 face < GeometryInfo<dim>::faces_per_cell;
785 typename DoFHandler<dim,spacedim>::face_iterator cell_face = cell->face(face);
786 if (cell_face->user_flag_set ())
789 if (cell->at_boundary (face) && !cell->has_periodic_neighbor(face))
791 for (
unsigned int i=0; i<total_dofs; ++i)
793 const bool i_non_zero_i = support_on_face (i, face);
794 for (
unsigned int j=0; j<total_dofs; ++j)
796 const bool j_non_zero_i = support_on_face (j, face);
799 sparsity.add (dofs_on_this_cell[i],
800 dofs_on_this_cell[j]);
802 && i_non_zero_i && j_non_zero_i)
803 sparsity.add (dofs_on_this_cell[i],
804 dofs_on_this_cell[j]);
811 neighbor = cell->neighbor_or_periodic_neighbor(face);
813 if (neighbor->level() < cell->level())
816 unsigned int neighbor_face = cell->has_periodic_neighbor(face)?
817 cell->periodic_neighbor_of_periodic_neighbor(face):
818 cell->neighbor_of_neighbor(face);
820 neighbor->get_mg_dof_indices (dofs_on_other_cell);
821 for (
unsigned int i=0; i<total_dofs; ++i)
823 const bool i_non_zero_i = support_on_face (i, face);
824 const bool i_non_zero_e = support_on_face (i, neighbor_face);
825 for (
unsigned int j=0; j<total_dofs; ++j)
827 const bool j_non_zero_i = support_on_face (j, face);
828 const bool j_non_zero_e = support_on_face (j, neighbor_face);
831 sparsity.add (dofs_on_this_cell[i],
832 dofs_on_other_cell[j]);
833 sparsity.add (dofs_on_other_cell[i],
834 dofs_on_this_cell[j]);
835 sparsity.add (dofs_on_this_cell[i],
836 dofs_on_this_cell[j]);
837 sparsity.add (dofs_on_other_cell[i],
838 dofs_on_other_cell[j]);
842 if (i_non_zero_i && j_non_zero_e)
843 sparsity.add (dofs_on_this_cell[i],
844 dofs_on_other_cell[j]);
845 if (i_non_zero_e && j_non_zero_i)
846 sparsity.add (dofs_on_other_cell[i],
847 dofs_on_this_cell[j]);
848 if (i_non_zero_i && j_non_zero_i)
849 sparsity.add (dofs_on_this_cell[i],
850 dofs_on_this_cell[j]);
851 if (i_non_zero_e && j_non_zero_e)
852 sparsity.add (dofs_on_other_cell[i],
853 dofs_on_other_cell[j]);
858 sparsity.add (dofs_on_this_cell[j],
859 dofs_on_other_cell[i]);
860 sparsity.add (dofs_on_other_cell[j],
861 dofs_on_this_cell[i]);
862 sparsity.add (dofs_on_this_cell[j],
863 dofs_on_this_cell[i]);
864 sparsity.add (dofs_on_other_cell[j],
865 dofs_on_other_cell[i]);
869 if (j_non_zero_i && i_non_zero_e)
870 sparsity.add (dofs_on_this_cell[j],
871 dofs_on_other_cell[i]);
872 if (j_non_zero_e && i_non_zero_i)
873 sparsity.add (dofs_on_other_cell[j],
874 dofs_on_this_cell[i]);
875 if (j_non_zero_i && i_non_zero_i)
876 sparsity.add (dofs_on_this_cell[j],
877 dofs_on_this_cell[i]);
878 if (j_non_zero_e && i_non_zero_e)
879 sparsity.add (dofs_on_other_cell[j],
880 dofs_on_other_cell[i]);
884 neighbor->face(neighbor_face)->set_user_flag ();
895 template <
int dim,
typename SparsityPatternType,
int spacedim>
898 SparsityPatternType &sparsity,
899 const unsigned int level,
916 Assert (sparsity.n_rows() == coarse_dofs,
918 Assert (sparsity.n_cols() == fine_dofs,
920 Assert (flux_mask.n_rows() == n_comp,
922 Assert (flux_mask.n_cols() == n_comp,
925 const unsigned int dofs_per_cell = dof.
get_fe().dofs_per_cell;
926 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
927 std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
931 endc = dof.
end(level);
936 for (
unsigned int i=0; i<dofs_per_cell; ++i)
937 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
940 for (; cell!=endc; ++cell)
942 if (!cell->is_locally_owned_on_level())
continue;
944 cell->get_mg_dof_indices (dofs_on_this_cell);
946 for (
unsigned int face = 0;
947 face < GeometryInfo<dim>::faces_per_cell;
951 bool use_face =
false;
952 if ((! cell->at_boundary(face)) &&
953 (static_cast<unsigned int>(cell->neighbor_level(face)) != level))
955 else if (cell->has_periodic_neighbor(face) &&
956 (
static_cast<unsigned int>(cell->periodic_neighbor_level(face)) != level))
962 neighbor = cell->neighbor_or_periodic_neighbor(face);
963 neighbor->get_mg_dof_indices (dofs_on_other_cell);
965 for (
unsigned int i=0; i<dofs_per_cell; ++i)
967 for (
unsigned int j=0; j<dofs_per_cell; ++j)
971 sparsity.add (dofs_on_other_cell[i],
972 dofs_on_this_cell[j]);
973 sparsity.add (dofs_on_other_cell[j],
974 dofs_on_this_cell[i]);
985 template <
int dim,
int spacedim>
988 std::vector<std::vector<types::global_dof_index> > &result,
990 std::vector<unsigned int> target_component)
996 Assert (result.size() == nlevels,
999 if (target_component.size() == 0)
1001 target_component.resize(n_components);
1002 for (
unsigned int i=0; i<n_components; ++i)
1003 target_component[i] = i;
1006 Assert(target_component.size() == n_components,
1009 for (
unsigned int l=0; l<nlevels; ++l)
1011 result[l].resize (n_components);
1012 std::fill (result[l].begin(),result[l].end(), 0U);
1018 if (n_components == 1)
1020 result[l][0] = dof_handler.
n_dofs(l);
1027 std::vector<std::vector<bool> >
1028 dofs_in_component (n_components,
1029 std::vector<bool>(dof_handler.
n_dofs(l),
1031 std::vector<ComponentMask> component_select (n_components);
1033 for (
unsigned int i=0; i<n_components; ++i)
1035 void (*fun_ptr) (
const unsigned int level,
1038 std::vector<bool> &)
1041 std::vector<bool> tmp(n_components,
false);
1047 component_select[i],
1048 dofs_in_component[i]);
1053 unsigned int component = 0;
1062 for (
unsigned int dd=0; dd<d; ++dd)
1065 result[l][target_component[component]]
1066 += std::count(dofs_in_component[component].begin(),
1067 dofs_in_component[component].end(),
1076 std::accumulate (result[l].begin(),
1077 result[l].end(), 0U)
1087 template <
typename DoFHandlerType>
1090 (
const DoFHandlerType &dof_handler,
1091 std::vector<std::vector<types::global_dof_index> > &dofs_per_block,
1092 std::vector<unsigned int> target_block)
1095 const unsigned int n_blocks = fe.
n_blocks();
1096 const unsigned int n_levels = dof_handler.get_triangulation().n_global_levels();
1100 for (
unsigned int l=0; l<n_levels; ++l)
1101 std::fill (dofs_per_block[l].begin(), dofs_per_block[l].end(), 0U);
1105 if (target_block.size()==0)
1107 target_block.resize(n_blocks);
1108 for (
unsigned int i=0; i<n_blocks; ++i)
1109 target_block[i] = i;
1111 Assert(target_block.size()==n_blocks,
1114 const unsigned int max_block
1115 = *std::max_element (target_block.begin(),
1116 target_block.end());
1117 const unsigned int n_target_blocks = max_block + 1;
1118 (void)n_target_blocks;
1120 for (
unsigned int l=0; l<n_levels; ++l)
1129 for (
unsigned int l=0; l<n_levels; ++l)
1130 dofs_per_block[l][0] = dof_handler.n_dofs(l);
1136 for (
unsigned int l=0; l<n_levels; ++l)
1138 std::vector<std::vector<bool> >
1139 dofs_in_block (n_blocks, std::vector<bool>(dof_handler.n_dofs(l),
false));
1140 std::vector<BlockMask> block_select (n_blocks);
1142 for (
unsigned int i=0; i<n_blocks; ++i)
1144 void (*fun_ptr) (
const unsigned int level,
1145 const DoFHandlerType &,
1147 std::vector<bool> &)
1148 = &DoFTools::extract_level_dofs<DoFHandlerType>;
1150 std::vector<bool> tmp(n_blocks,
false);
1152 block_select[i] = tmp;
1155 l, dof_handler, block_select[i],
1161 for (
unsigned int block=0; block<fe.
n_blocks(); ++block)
1162 dofs_per_block[l][target_block[block]]
1163 += std::count(dofs_in_block[block].begin(),
1164 dofs_in_block[block].end(),
1171 template <
int dim,
int spacedim>
1176 std::vector<std::set<types::global_dof_index> > &boundary_indices,
1183 std::set<types::boundary_id> boundary_ids;
1184 typename std::map<types::boundary_id, const Function<dim>* >::const_iterator it
1185 = function_map.begin();
1186 for (; it!=function_map.end(); ++it)
1187 boundary_ids.insert(it->first);
1189 std::vector<IndexSet> boundary_indexset;
1192 boundary_indices[i].insert(boundary_indexset[i].begin(),
1193 boundary_indexset[i].end());
1197 template <
int dim,
int spacedim>
1201 std::vector<IndexSet> &boundary_indices,
1208 std::set<types::boundary_id> boundary_ids;
1209 typename std::map<types::boundary_id, const Function<dim>* >::const_iterator it = function_map.begin();
1210 for (; it!=function_map.end(); ++it)
1211 boundary_ids.insert(it->first);
1218 template <
int dim,
int spacedim>
1221 const std::set<types::boundary_id> &boundary_ids,
1222 std::vector<IndexSet> &boundary_indices,
1228 if (boundary_ids.size() == 0)
1232 if (boundary_indices[i].size() == 0)
1235 const unsigned int n_components = DoFTools::n_components(dof);
1236 const bool fe_is_system = (n_components != 1);
1238 std::vector<types::global_dof_index> local_dofs;
1239 local_dofs.reserve (DoFTools::max_dofs_per_face(dof));
1240 std::fill (local_dofs.begin (),
1250 for (; cell!=endc; ++cell)
1256 const unsigned int level = cell->level();
1259 for (
unsigned int face_no = 0; face_no < GeometryInfo<dim>::faces_per_cell;
1261 if (cell->at_boundary(face_no) ==
true)
1263 const typename DoFHandler<dim,spacedim>::face_iterator
1264 face = cell->face(face_no);
1267 if (boundary_ids.find(bi) != boundary_ids.end())
1269 face->get_mg_dof_indices(level, local_dofs);
1270 boundary_indices[level].add_indices(local_dofs.begin(), local_dofs.end());
1278 ExcMessage(
"It's probably worthwhile to select at least one component."));
1283 for (; cell!=endc; ++cell)
1286 for (
unsigned int face_no = 0; face_no < GeometryInfo<dim>::faces_per_cell;
1289 if (cell->at_boundary(face_no) ==
false)
1293 const unsigned int level = cell->level();
1295 typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(face_no);
1297 if (boundary_ids.find(boundary_component) != boundary_ids.end())
1301 for (
unsigned int i=0; i<cell->get_fe().dofs_per_cell; ++i)
1304 = cell->get_fe().get_nonzero_components (i);
1308 bool selected =
false;
1309 for (
unsigned int c=0; c<n_components; ++c)
1310 if (nonzero_component_array[c] ==
true 1311 && component_mask[c]==
true)
1317 for (
unsigned int c=0; c<n_components; ++c)
1318 Assert (nonzero_component_array[c] ==
false || component_mask[c] ==
true,
1319 ExcMessage (
"You are using a non-primitive FiniteElement " 1320 "and try to constrain just some of its components!"));
1325 face->get_mg_dof_indices (level, local_dofs);
1328 for (
unsigned int i=0; i<local_dofs.size(); ++i)
1339 = cell->get_fe().get_nonzero_components (i);
1340 for (
unsigned int c=0; c<n_components; ++c)
1341 if (nonzero_component_array[c] ==
true)
1348 if (component_mask[component] ==
true)
1349 boundary_indices[level].add_index(local_dofs[i]);
1353 boundary_indices[level].add_indices(local_dofs.begin(), local_dofs.end());
1360 template <
int dim,
int spacedim>
1363 std::vector<std::set<types::global_dof_index> > &non_interface_dofs)
1374 std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
1375 std::vector<bool> cell_dofs(dofs_per_cell,
false);
1376 std::vector<bool> cell_dofs_interface(dofs_per_cell,
false);
1379 endc = mg_dof_handler.
end();
1382 for (; cell!=endc; ++cell)
1385 && cell->level_subdomain_id()!=mg_dof_handler.
get_triangulation().locally_owned_subdomain())
1388 std::fill (cell_dofs.begin(), cell_dofs.end(),
false);
1389 std::fill (cell_dofs_interface.begin(), cell_dofs_interface.end(),
false);
1391 for (
unsigned int face_nr=0; face_nr<GeometryInfo<dim>::faces_per_cell; ++face_nr)
1393 const typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(face_nr);
1394 if (!face->at_boundary())
1398 neighbor = cell->neighbor(face_nr);
1400 if ((neighbor->level() < cell->level()))
1402 for (
unsigned int j=0; j<dofs_per_face; ++j)
1407 for (
unsigned int j=0; j<dofs_per_face; ++j)
1414 for (
unsigned int j=0; j<dofs_per_face; ++j)
1419 const unsigned int level = cell->level();
1420 cell->get_mg_dof_indices (local_dof_indices);
1422 for (
unsigned int i=0; i<dofs_per_cell; ++i)
1423 if (cell_dofs[i] && !cell_dofs_interface[i])
1424 non_interface_dofs[level].insert(local_dof_indices[i]);
1429 template <
int dim,
int spacedim>
1432 std::vector<IndexSet> &interface_dofs)
1438 std::vector<std::vector<types::global_dof_index> >
1439 tmp_interface_dofs(interface_dofs.size());
1443 const unsigned int dofs_per_cell = fe.dofs_per_cell;
1444 const unsigned int dofs_per_face = fe.dofs_per_face;
1446 std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
1448 std::vector<bool> cell_dofs(dofs_per_cell,
false);
1451 endc = mg_dof_handler.
end();
1453 for (; cell!=endc; ++cell)
1461 bool has_coarser_neighbor =
false;
1463 std::fill (cell_dofs.begin(), cell_dofs.end(),
false);
1465 for (
unsigned int face_nr=0; face_nr<GeometryInfo<dim>::faces_per_cell; ++face_nr)
1467 const typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(face_nr);
1468 if (!face->at_boundary())
1472 neighbor = cell->neighbor(face_nr);
1481 if (neighbor->level() < cell->level())
1483 for (
unsigned int j=0; j<dofs_per_face; ++j)
1484 cell_dofs[fe.face_to_cell_index(j,face_nr)] =
true;
1486 has_coarser_neighbor =
true;
1491 if (has_coarser_neighbor ==
false)
1494 const unsigned int level = cell->level();
1495 cell->get_mg_dof_indices (local_dof_indices);
1497 for (
unsigned int i=0; i<dofs_per_cell; ++i)
1500 tmp_interface_dofs[level].push_back(local_dof_indices[i]);
1504 for (
unsigned int l=0; l<mg_dof_handler.
get_triangulation().n_global_levels(); ++l)
1506 interface_dofs[l].clear();
1507 std::sort(tmp_interface_dofs[l].begin(), tmp_interface_dofs[l].end());
1508 interface_dofs[l].add_indices(tmp_interface_dofs[l].begin(),
1509 std::unique(tmp_interface_dofs[l].begin(),
1510 tmp_interface_dofs[l].end()));
1511 interface_dofs[l].compress();
1519 #include "mg_tools.inst" 1521 DEAL_II_NAMESPACE_CLOSE
const unsigned int first_hex_index
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
void load_user_flags(std::istream &in)
cell_iterator end() const
ActiveSelector::active_cell_iterator active_cell_iterator
bool is_primitive() const
const FiniteElement< dim, spacedim > & get_fe() const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
unsigned int n_blocks() const
static ::ExceptionBase & ExcMessage(std::string arg1)
ActiveSelector::cell_iterator cell_iterator
const unsigned int first_quad_index
unsigned int global_dof_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
std::map< types::boundary_id, const Function< dim, Number > * > type
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
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
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()
const Triangulation< dim, spacedim > & get_triangulation() const
unsigned char boundary_id
unsigned int n_base_elements() const
Task< RT > new_task(const std_cxx11::function< RT()> &function)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
static ::ExceptionBase & ExcInternalError()