16 #include <deal.II/base/thread_management.h> 17 #include <deal.II/base/quadrature_lib.h> 18 #include <deal.II/base/table.h> 19 #include <deal.II/base/template_constraints.h> 20 #include <deal.II/base/utilities.h> 21 #include <deal.II/lac/sparsity_pattern.h> 22 #include <deal.II/lac/dynamic_sparsity_pattern.h> 23 #include <deal.II/lac/trilinos_sparsity_pattern.h> 24 #include <deal.II/lac/block_sparsity_pattern.h> 25 #include <deal.II/lac/vector.h> 26 #include <deal.II/lac/constraint_matrix.h> 27 #include <deal.II/grid/tria.h> 28 #include <deal.II/grid/tria_iterator.h> 29 #include <deal.II/grid/intergrid_map.h> 30 #include <deal.II/grid/filtered_iterator.h> 31 #include <deal.II/grid/grid_tools.h> 32 #include <deal.II/dofs/dof_handler.h> 33 #include <deal.II/dofs/dof_accessor.h> 34 #include <deal.II/fe/fe.h> 35 #include <deal.II/fe/fe_values.h> 36 #include <deal.II/fe/fe_tools.h> 37 #include <deal.II/hp/fe_collection.h> 38 #include <deal.II/hp/q_collection.h> 39 #include <deal.II/hp/fe_values.h> 40 #include <deal.II/dofs/dof_tools.h> 41 #include <deal.II/distributed/tria.h> 42 #include <deal.II/distributed/shared_tria.h> 47 DEAL_II_NAMESPACE_OPEN
66 template <
int dim,
int spacedim>
67 std::vector<unsigned char>
71 std::vector<unsigned char> local_component_association (fe.
dofs_per_cell,
80 local_component_association[i] =
91 const unsigned int first_comp =
96 component_mask).n_selected_components(fe.
n_components()) == 0)
97 local_component_association[i] = first_comp;
104 for (
unsigned int c=first_comp; c<fe.
n_components(); ++c)
105 if (component_mask[c] ==
true)
107 local_component_association[i] = c;
112 Assert (std::find (local_component_association.begin(),
113 local_component_association.end(),
116 local_component_association.end(),
119 return local_component_association;
139 template <
typename DoFHandlerType>
141 get_component_association (
const DoFHandlerType &dof,
143 std::vector<unsigned char> &dofs_by_component)
145 const ::hp::FECollection<DoFHandlerType::dimension,DoFHandlerType::space_dimension>
146 fe_collection (dof.get_fe());
147 Assert (fe_collection.n_components() < 256, ExcNotImplemented());
148 Assert (dofs_by_component.size() == dof.n_locally_owned_dofs(),
149 ExcDimensionMismatch(dofs_by_component.size(),
150 dof.n_locally_owned_dofs()));
159 std::vector<std::vector<unsigned char> >
160 local_component_association (fe_collection.size());
161 for (
unsigned int f=0; f<fe_collection.size(); ++f)
165 local_component_association[f]
166 = get_local_component_association (fe, component_mask);
170 std::vector<types::global_dof_index> indices;
171 for (
typename DoFHandlerType::active_cell_iterator c=dof.begin_active();
173 if (c->is_locally_owned())
175 const unsigned int fe_index = c->active_fe_index();
176 const unsigned int dofs_per_cell = c->get_fe().dofs_per_cell;
177 indices.resize(dofs_per_cell);
178 c->get_dof_indices(indices);
179 for (
unsigned int i=0; i<dofs_per_cell; ++i)
180 if (dof.locally_owned_dofs().is_element(indices[i]))
181 dofs_by_component[dof.locally_owned_dofs().index_within_set(indices[i])]
182 = local_component_association[fe_index][i];
194 template <
typename DoFHandlerType>
197 get_block_association (
const DoFHandlerType &dof,
198 std::vector<unsigned char> &dofs_by_block)
200 const ::hp::FECollection<DoFHandlerType::dimension,DoFHandlerType::space_dimension>
201 fe_collection (dof.get_fe());
202 Assert (fe_collection.n_components() < 256, ExcNotImplemented());
203 Assert (dofs_by_block.size() == dof.n_locally_owned_dofs(),
204 ExcDimensionMismatch(dofs_by_block.size(),
205 dof.n_locally_owned_dofs()));
214 std::vector<std::vector<unsigned char> > local_block_association
215 (fe_collection.size());
216 for (
unsigned int f=0; f<fe_collection.size(); ++f)
221 (
unsigned char)(-1));
225 Assert (std::find (local_block_association[f].begin(),
226 local_block_association[f].end(),
229 local_block_association[f].end(),
234 std::vector<types::global_dof_index> indices;
235 for (
typename DoFHandlerType::active_cell_iterator c=dof.begin_active();
237 if (c->is_locally_owned())
239 const unsigned int fe_index = c->active_fe_index();
240 const unsigned int dofs_per_cell = c->get_fe().dofs_per_cell;
241 indices.resize(dofs_per_cell);
242 c->get_dof_indices(indices);
243 for (
unsigned int i=0; i<dofs_per_cell; ++i)
244 if (dof.locally_owned_dofs().is_element(indices[i]))
245 dofs_by_block[dof.locally_owned_dofs().index_within_set(indices[i])]
246 = local_block_association[fe_index][i];
253 template <
typename DoFHandlerType,
typename Number>
257 const unsigned int component)
259 const unsigned int dim = DoFHandlerType::dimension;
260 const unsigned int spacedim = DoFHandlerType::space_dimension;
267 Assert (fe_is_primitive(dof_handler) ==
true,
273 const bool consider_components = (n_components(dof_handler) != 1);
276 if (consider_components ==
false)
280 std::vector<unsigned char> component_dofs (dof_handler.n_locally_owned_dofs());
281 internal::get_component_association (dof_handler,
282 dof_handler.get_fe().component_mask
286 for (
unsigned int i=0; i<dof_data.
size(); ++i)
287 if (component_dofs[i] == static_cast<unsigned char>(component))
292 std::vector<unsigned char> touch_count (dof_handler.n_dofs(), 0);
294 typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active(),
295 endc = dof_handler.
end();
296 std::vector<types::global_dof_index> dof_indices;
297 dof_indices.reserve (max_dofs_per_cell(dof_handler));
299 for (
unsigned int present_cell = 0; cell!=endc; ++cell, ++present_cell)
301 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
302 dof_indices.resize (dofs_per_cell);
303 cell->get_dof_indices (dof_indices);
305 for (
unsigned int i=0; i<dofs_per_cell; ++i)
308 if (!consider_components ||
309 (cell->get_fe().system_to_component_index(i).first == component))
312 dof_data(dof_indices[i]) += cell_data(present_cell);
315 ++touch_count[dof_indices[i]];
325 Assert (consider_components || (touch_count[i]!=0),
327 if (touch_count[i] != 0)
328 dof_data(i) /= touch_count[i];
334 template <
int dim,
int spacedim>
338 std::vector<bool> &selected_dofs)
344 ExcMessage (
"The given component mask is not sized correctly to represent the " 345 "components of the given finite element."));
369 internal::get_component_association (dof, component_mask,
373 if (component_mask[dofs_by_component[i]] ==
true)
374 selected_dofs[i] =
true;
380 template <
int dim,
int spacedim>
384 std::vector<bool> &selected_dofs)
390 ExcMessage (
"The given component mask is not sized correctly to represent the " 391 "components of the given finite element."));
393 ExcDimensionMismatch(selected_dofs.size(), dof.
n_dofs()));
399 std::fill_n (selected_dofs.begin(), dof.
n_dofs(),
false);
404 std::fill_n (selected_dofs.begin(), dof.
n_dofs(),
true);
410 std::fill_n (selected_dofs.begin(), dof.
n_dofs(),
false);
414 std::vector<unsigned char> dofs_by_component (dof.
n_dofs());
415 internal::get_component_association (dof, component_mask,
419 if (component_mask[dofs_by_component[i]] ==
true)
420 selected_dofs[i] =
true;
425 template <
int dim,
int spacedim>
429 std::vector<bool> &selected_dofs)
438 template <
int dim,
int spacedim>
442 std::vector<bool> &selected_dofs)
451 template<
typename DoFHandlerType>
454 const DoFHandlerType &dof,
456 std::vector<bool> &selected_dofs)
461 ExcMessage (
"The given component mask is not sized correctly to represent the " 462 "components of the given finite element."));
463 Assert(selected_dofs.size() == dof.n_dofs(level),
464 ExcDimensionMismatch(selected_dofs.size(), dof.n_dofs(level)));
470 std::fill_n (selected_dofs.begin(), dof.n_dofs(level),
false);
475 std::fill_n (selected_dofs.begin(), dof.n_dofs(level),
true);
480 std::fill_n (selected_dofs.begin(), dof.n_dofs(level),
false);
484 std::vector<unsigned char> local_component_asssociation
485 = internal::get_local_component_association (fe, component_mask);
488 local_selected_dofs[i] = component_mask[local_component_asssociation[i]];
491 std::vector<types::global_dof_index> indices(fe.
dofs_per_cell);
492 typename DoFHandlerType::level_cell_iterator c;
493 for (c = dof.begin(level) ; c != dof.end(level) ; ++ c)
495 c->get_mg_dof_indices(indices);
497 selected_dofs[indices[i]] = local_selected_dofs[i];
503 template<
typename DoFHandlerType>
506 const DoFHandlerType &dof,
508 std::vector<bool> &selected_dofs)
517 template <
typename DoFHandlerType>
521 std::vector<bool> &selected_dofs,
522 const std::set<types::boundary_id> &boundary_ids)
525 (&dof_handler.get_triangulation())
527 ExcMessage (
"This function can not be used with distributed triangulations." 528 "See the documentation for more information."));
532 indices, boundary_ids);
535 selected_dofs.clear ();
536 selected_dofs.resize (dof_handler.n_dofs(),
false);
539 indices.fill_binary_vector(selected_dofs);
543 template <
typename DoFHandlerType>
548 const std::set<types::boundary_id> &boundary_ids)
551 ExcMessage (
"Component mask has invalid size."));
553 ExcInvalidBoundaryIndicator());
554 const unsigned int dim=DoFHandlerType::dimension;
557 selected_dofs.
clear ();
558 selected_dofs.
set_size(dof_handler.n_dofs());
562 const bool check_boundary_id = (boundary_ids.size() != 0);
566 const bool check_vector_component
570 n_components(dof_handler)));
572 std::vector<types::global_dof_index> face_dof_indices;
573 face_dof_indices.reserve (max_dofs_per_face(dof_handler));
581 for (
typename DoFHandlerType::active_cell_iterator cell=dof_handler.begin_active();
582 cell!=dof_handler.end(); ++cell)
586 if (cell->is_artificial() ==
false)
587 for (
unsigned int face=0;
588 face<GeometryInfo<DoFHandlerType::dimension>::faces_per_cell; ++face)
589 if (cell->at_boundary(face))
590 if (! check_boundary_id ||
591 (boundary_ids.find (cell->face(face)->boundary_id())
592 != boundary_ids.end()))
595 &fe = cell->get_fe();
598 face_dof_indices.resize (dofs_per_face);
599 cell->face(face)->get_dof_indices (face_dof_indices,
600 cell->active_fe_index());
603 if (!check_vector_component)
604 selected_dofs.
add_index (face_dof_indices[i]);
612 const unsigned int cell_index
633 selected_dofs.
add_index (face_dof_indices[i]);
637 const unsigned int first_nonzero_comp
642 if (component_mask[first_nonzero_comp] ==
true)
643 selected_dofs.
add_index (face_dof_indices[i]);
651 template <
typename DoFHandlerType>
655 std::vector<bool> &selected_dofs,
656 const std::set<types::boundary_id> &boundary_ids)
659 ExcMessage (
"This component mask has the wrong size."));
661 ExcInvalidBoundaryIndicator());
665 const bool check_boundary_id = (boundary_ids.size() != 0);
669 const bool check_vector_component
673 selected_dofs.clear ();
674 selected_dofs.resize (dof_handler.n_dofs(),
false);
675 std::vector<types::global_dof_index> cell_dof_indices;
676 cell_dof_indices.reserve (max_dofs_per_cell(dof_handler));
684 for (
typename DoFHandlerType::active_cell_iterator cell=dof_handler.begin_active();
685 cell!=dof_handler.end(); ++cell)
686 for (
unsigned int face=0;
687 face<GeometryInfo<DoFHandlerType::dimension>::faces_per_cell; ++face)
688 if (cell->at_boundary(face))
689 if (! check_boundary_id ||
690 (boundary_ids.find (cell->face(face)->boundary_id())
691 != boundary_ids.end()))
697 cell_dof_indices.resize (dofs_per_cell);
698 cell->get_dof_indices (cell_dof_indices);
703 if (!check_vector_component)
704 selected_dofs[cell_dof_indices[i]] =
true;
711 selected_dofs[cell_dof_indices[i]]
716 const unsigned int first_nonzero_comp
721 selected_dofs[cell_dof_indices[i]]
722 = (component_mask[first_nonzero_comp]
737 template <
int spacedim>
739 std::vector<bool> &selected_dofs)
741 Assert(selected_dofs.size() == dof_handler.n_dofs(),
742 ExcDimensionMismatch(selected_dofs.size(), dof_handler.n_dofs()));
744 std::fill_n (selected_dofs.begin(), dof_handler.n_dofs(),
false);
750 template <
int spacedim>
752 std::vector<bool> &selected_dofs)
754 const unsigned int dim = 2;
756 Assert(selected_dofs.size() == dof_handler.n_dofs(),
757 ExcDimensionMismatch(selected_dofs.size(), dof_handler.n_dofs()));
759 std::fill_n (selected_dofs.begin(), dof_handler.n_dofs(),
false);
765 typename ::DoFHandler<dim,spacedim>::active_cell_iterator
766 cell = dof_handler.begin_active(),
767 endc = dof_handler.end();
768 for (; cell!=endc; ++cell)
769 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
770 if (cell->face(face)->has_children())
772 const typename ::DoFHandler<dim,spacedim>::line_iterator
773 line = cell->face(face);
775 for (
unsigned int dof=0; dof!=fe.dofs_per_vertex; ++dof)
776 selected_dofs[line->child(0)->vertex_dof_index(1,dof)] =
true;
778 for (
unsigned int child=0; child<2; ++child)
779 for (
unsigned int dof=0; dof!=fe.dofs_per_line; ++dof)
780 selected_dofs[line->child(child)->dof_index(dof)] =
true;
785 template <
int spacedim>
787 std::vector<bool> &selected_dofs)
789 const unsigned int dim = 3;
791 Assert(selected_dofs.size() == dof_handler.n_dofs(),
792 ExcDimensionMismatch(selected_dofs.size(), dof_handler.n_dofs()));
794 std::fill_n (selected_dofs.begin(), dof_handler.n_dofs(),
false);
801 typename ::DoFHandler<dim,spacedim>::active_cell_iterator
802 cell = dof_handler.begin_active(),
803 endc = dof_handler.end();
804 for (; cell!=endc; ++cell)
805 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
806 if (cell->face(f)->has_children())
808 const typename ::DoFHandler<dim,spacedim>::face_iterator
809 face = cell->face(f);
811 for (
unsigned int dof=0; dof!=fe.dofs_per_vertex; ++dof)
812 selected_dofs[face->child(0)->vertex_dof_index(2,dof)] =
true;
816 for (
unsigned int line=0; line<4; ++line)
817 for (
unsigned int dof=0; dof!=fe.dofs_per_vertex; ++dof)
818 selected_dofs[face->line(line)->child(0)->vertex_dof_index(1,dof)] =
true;
823 for (
unsigned int dof=0; dof<fe.dofs_per_line; ++dof)
824 selected_dofs[face->child(0)->line(1)->dof_index(dof)] =
true;
825 for (
unsigned int dof=0; dof<fe.dofs_per_line; ++dof)
826 selected_dofs[face->child(1)->line(2)->dof_index(dof)] =
true;
827 for (
unsigned int dof=0; dof<fe.dofs_per_line; ++dof)
828 selected_dofs[face->child(2)->line(3)->dof_index(dof)] =
true;
829 for (
unsigned int dof=0; dof<fe.dofs_per_line; ++dof)
830 selected_dofs[face->child(3)->line(0)->dof_index(dof)] =
true;
833 for (
unsigned int line=0; line<4; ++line)
834 for (
unsigned int child=0; child<2; ++child)
835 for (
unsigned int dof=0; dof!=fe.dofs_per_line; ++dof)
836 selected_dofs[face->line(line)->child(child)->dof_index(dof)] =
true;
839 for (
unsigned int child=0; child<4; ++child)
840 for (
unsigned int dof=0; dof!=fe.dofs_per_quad; ++dof)
841 selected_dofs[face->child(child)->dof_index(dof)] =
true;
849 template <
int dim,
int spacedim>
853 std::vector<bool> &selected_dofs)
855 internal::extract_hanging_node_dofs (dof_handler,
861 template <
typename DoFHandlerType>
865 std::vector<bool> &selected_dofs)
867 Assert(selected_dofs.size() == dof_handler.n_dofs(),
868 ExcDimensionMismatch(selected_dofs.size(), dof_handler.n_dofs()));
871 std::fill_n (selected_dofs.begin(), dof_handler.n_dofs(),
false);
873 std::vector<types::global_dof_index> local_dof_indices;
874 local_dof_indices.reserve (max_dofs_per_cell(dof_handler));
878 typename DoFHandlerType::active_cell_iterator
879 cell = dof_handler.begin_active(),
880 endc = dof_handler.end();
881 for (; cell!=endc; ++cell)
882 if (cell->subdomain_id() == subdomain_id)
884 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
885 local_dof_indices.resize (dofs_per_cell);
886 cell->get_dof_indices (local_dof_indices);
887 for (
unsigned int i=0; i<dofs_per_cell; ++i)
888 selected_dofs[local_dof_indices[i]] =
true;
894 template <
typename DoFHandlerType>
900 dof_set = dof_handler.locally_owned_dofs();
906 template <
typename DoFHandlerType>
912 dof_set = dof_handler.locally_owned_dofs();
917 std::vector<types::global_dof_index> dof_indices;
918 std::set<types::global_dof_index> global_dof_indices;
920 typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active(),
921 endc = dof_handler.
end();
922 for (; cell!=endc; ++cell)
923 if (cell->is_locally_owned())
925 dof_indices.resize(cell->get_fe().dofs_per_cell);
926 cell->get_dof_indices(dof_indices);
928 for (std::vector<types::global_dof_index>::iterator it=dof_indices.begin();
929 it!=dof_indices.end();
932 global_dof_indices.insert(*it);
935 dof_set.
add_indices(global_dof_indices.begin(), global_dof_indices.end());
942 template <
typename DoFHandlerType>
948 dof_set = dof_handler.locally_owned_dofs();
958 std::vector<types::global_dof_index> dof_indices;
959 std::vector<types::global_dof_index> dofs_on_ghosts;
961 typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active(),
962 endc = dof_handler.
end();
963 for (; cell!=endc; ++cell)
964 if (cell->is_ghost())
966 dof_indices.resize(cell->get_fe().dofs_per_cell);
967 cell->get_dof_indices(dof_indices);
968 for (
unsigned int i=0; i<dof_indices.size(); ++i)
970 dofs_on_ghosts.push_back(dof_indices[i]);
974 std::sort(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
975 dof_set.
add_indices(dofs_on_ghosts.begin(), std::unique(dofs_on_ghosts.begin(),
976 dofs_on_ghosts.end()));
981 template <
typename DoFHandlerType>
984 const unsigned int level,
988 dof_set = dof_handler.locally_owned_mg_dofs(level);
998 std::vector<types::global_dof_index> dof_indices;
999 std::vector<types::global_dof_index> dofs_on_ghosts;
1001 typename DoFHandlerType::cell_iterator cell = dof_handler.
begin(level),
1002 endc = dof_handler.end(level);
1003 for (; cell!=endc; ++cell)
1008 if (
id == dof_handler.get_triangulation().locally_owned_subdomain()
1012 dof_indices.resize(cell->get_fe().dofs_per_cell);
1013 cell->get_mg_dof_indices(dof_indices);
1014 for (
unsigned int i=0; i<dof_indices.size(); ++i)
1016 dofs_on_ghosts.push_back(dof_indices[i]);
1020 std::sort(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1021 dof_set.
add_indices(dofs_on_ghosts.begin(), std::unique(dofs_on_ghosts.begin(),
1022 dofs_on_ghosts.end()));
1027 template <
typename DoFHandlerType>
1031 std::vector<std::vector<bool> > &constant_modes)
1033 const unsigned int n_components = dof_handler.get_fe().n_components();
1035 ExcDimensionMismatch(n_components,
1036 component_mask.
size()));
1038 std::vector<unsigned char> dofs_by_component (dof_handler.n_locally_owned_dofs());
1039 internal::get_component_association (dof_handler, component_mask,
1041 unsigned int n_selected_dofs = 0;
1042 for (
unsigned int i=0; i<n_components; ++i)
1043 if (component_mask[i] ==
true)
1044 n_selected_dofs += std::count (dofs_by_component.begin(),
1045 dofs_by_component.end(), i);
1048 const IndexSet &locally_owned_dofs = dof_handler.locally_owned_dofs();
1049 std::vector<unsigned int> component_numbering(locally_owned_dofs.n_elements(),
1051 for (
unsigned int i=0, count=0; i<locally_owned_dofs.n_elements(); ++i)
1052 if (component_mask[dofs_by_component[i]])
1053 component_numbering[i] = count++;
1060 const ::hp::FECollection<DoFHandlerType::dimension,DoFHandlerType::space_dimension>
1061 fe_collection (dof_handler.get_fe());
1062 std::vector<Table<2,bool> > element_constant_modes;
1063 std::vector<std::vector<std::pair<unsigned int, unsigned int> > >
1064 constant_mode_to_component_translation(n_components);
1065 unsigned int n_constant_modes = 0;
1066 for (
unsigned int f=0; f<fe_collection.size(); ++f)
1068 std::pair<Table<2,bool>, std::vector<unsigned int> > data
1069 = fe_collection[f].get_constant_modes();
1070 element_constant_modes.push_back(data.first);
1072 for (
unsigned int i=0; i<data.second.size(); ++i)
1073 if (component_mask[data.second[i]])
1074 constant_mode_to_component_translation[data.second[i]].
1075 push_back(std::make_pair(n_constant_modes++,i));
1077 element_constant_modes[0].n_rows());
1081 constant_modes.clear ();
1082 constant_modes.resize (n_constant_modes, std::vector<bool>(n_selected_dofs,
1087 typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active(),
1088 endc = dof_handler.end();
1089 std::vector<types::global_dof_index> dof_indices;
1090 for (; cell!=endc; ++cell)
1091 if (cell->is_locally_owned())
1093 dof_indices.resize(cell->get_fe().dofs_per_cell);
1094 cell->get_dof_indices(dof_indices);
1096 for (
unsigned int i=0; i<dof_indices.size(); ++i)
1097 if (locally_owned_dofs.is_element(dof_indices[i]))
1099 const unsigned int loc_index =
1100 locally_owned_dofs.index_within_set(dof_indices[i]);
1101 const unsigned int comp = dofs_by_component[loc_index];
1102 if (component_mask[comp])
1103 for (
unsigned int j=0; j<constant_mode_to_component_translation[comp].size(); ++j)
1104 constant_modes[constant_mode_to_component_translation[comp][j].first]
1105 [component_numbering[loc_index]] =
1106 element_constant_modes[cell->active_fe_index()]
1107 (constant_mode_to_component_translation[comp][j].second,i);
1114 template <
typename DoFHandlerType>
1117 std::vector<unsigned int> &active_fe_indices)
1119 AssertDimension (active_fe_indices.size(), dof_handler.get_triangulation().n_active_cells());
1121 typename DoFHandlerType::active_cell_iterator
1122 cell = dof_handler.begin_active(),
1123 endc = dof_handler.end();
1124 for (; cell!=endc; ++cell)
1125 active_fe_indices[cell->active_cell_index()] = cell->active_fe_index();
1128 template <
typename DoFHandlerType>
1129 std::vector<IndexSet>
1134 Assert ((
dynamic_cast<const parallel::distributed::
1135 Triangulation<DoFHandlerType::dimension,DoFHandlerType::space_dimension> *
> 1136 (&dof_handler.get_triangulation()) == 0),
1137 ExcMessage (
"For parallel::distributed::Triangulation objects and " 1138 "associated DoF handler objects, asking for any information " 1139 "related to a subdomain other than the locally owned one does " 1140 "not make sense."));
1143 std::vector< ::types::subdomain_id > subdomain_association (dof_handler.n_dofs ());
1144 ::DoFTools::get_subdomain_association (dof_handler, subdomain_association);
1146 const unsigned int n_subdomains = 1 + (*std::max_element (subdomain_association.begin (),
1147 subdomain_association.end () ));
1149 std::vector<::IndexSet> index_sets (n_subdomains,::
IndexSet(dof_handler.n_dofs()));
1157 index < subdomain_association.size (); ++index)
1160 if (subdomain_association[index] != this_subdomain)
1162 index_sets[this_subdomain].add_range (i_min, index);
1164 this_subdomain = subdomain_association[index];
1169 if (i_min == subdomain_association.size () - 1)
1171 index_sets[this_subdomain].add_index (i_min);
1177 index_sets[this_subdomain].add_range (
1178 i_min, subdomain_association.size ());
1181 for (
unsigned int i = 0; i < n_subdomains; i++)
1182 index_sets[i].compress ();
1187 template <
typename DoFHandlerType>
1188 std::vector<IndexSet>
1193 Assert ((
dynamic_cast<const parallel::distributed::
1194 Triangulation<DoFHandlerType::dimension,DoFHandlerType::space_dimension> *
> 1195 (&dof_handler.get_triangulation()) == 0),
1196 ExcMessage (
"For parallel::distributed::Triangulation objects and " 1197 "associated DoF handler objects, asking for any information " 1198 "related to a subdomain other than the locally owned one does " 1199 "not make sense."));
1207 const ::types::subdomain_id n_subdomains = dof_set.size();
1213 subdomain_id < n_subdomains; ++subdomain_id)
1216 std_cxx11::function<bool (const typename DoFHandlerType::active_cell_iterator &)> predicate
1218 const std::vector<typename DoFHandlerType::active_cell_iterator>
1223 std::vector<types::global_dof_index> local_dof_indices;
1224 std::set<types::global_dof_index> subdomain_halo_global_dof_indices;
1225 for (
typename std::vector<typename DoFHandlerType::active_cell_iterator>::const_iterator
1226 it_cell = active_halo_layer.begin(); it_cell!=active_halo_layer.end(); ++it_cell)
1228 const typename DoFHandlerType::active_cell_iterator &cell = *it_cell;
1229 Assert(cell->subdomain_id() != subdomain_id,
1230 ExcMessage(
"The subdomain ID of the halo cell should not match that of the vector entry."));
1232 local_dof_indices.resize(cell->get_fe().dofs_per_cell);
1233 cell->get_dof_indices(local_dof_indices);
1235 for (std::vector<types::global_dof_index>::iterator it=local_dof_indices.begin();
1236 it!=local_dof_indices.end();
1238 subdomain_halo_global_dof_indices.insert(*it);
1241 dof_set[subdomain_id].add_indices(subdomain_halo_global_dof_indices.begin(),
1242 subdomain_halo_global_dof_indices.end());
1244 dof_set[subdomain_id].compress();
1250 template <
typename DoFHandlerType>
1253 std::vector<types::subdomain_id> &subdomain_association)
1257 Assert ((
dynamic_cast<const parallel::distributed::
1258 Triangulation<DoFHandlerType::dimension,DoFHandlerType::space_dimension> *
> 1259 (&dof_handler.get_triangulation()) == 0),
1260 ExcMessage (
"For parallel::distributed::Triangulation objects and " 1261 "associated DoF handler objects, asking for any subdomain other " 1262 "than the locally owned one does not make sense."));
1264 Assert(subdomain_association.size() == dof_handler.n_dofs(),
1265 ExcDimensionMismatch(subdomain_association.size(),
1266 dof_handler.n_dofs()));
1268 Assert(dof_handler.n_dofs() > 0,
1270 "This could happen when the function is called before NumberCache is written."));
1276 std::vector<types::subdomain_id> cell_owners (dof_handler.get_triangulation().n_active_cells());
1279 DoFHandlerType::space_dimension
>*> (&dof_handler.get_triangulation())))
1281 cell_owners = tr->get_true_subdomain_ids_of_cells();
1282 Assert (tr->get_true_subdomain_ids_of_cells().size() == tr->n_active_cells(),
1283 ExcInternalError());
1287 for (
typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active();
1288 cell!= dof_handler.end(); cell++)
1289 if (cell->is_locally_owned())
1290 cell_owners[cell->active_cell_index()] = cell->subdomain_id();
1294 std::fill_n (subdomain_association.begin(), dof_handler.n_dofs(),
1297 std::vector<types::global_dof_index> local_dof_indices;
1298 local_dof_indices.reserve (max_dofs_per_cell(dof_handler));
1302 bool coin_flip =
true;
1306 typename DoFHandlerType::active_cell_iterator
1307 cell = dof_handler.begin_active(),
1308 endc = dof_handler.end();
1309 for (; cell!=endc; ++cell)
1312 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
1313 local_dof_indices.resize (dofs_per_cell);
1314 cell->get_dof_indices (local_dof_indices);
1320 for (
unsigned int i=0; i<dofs_per_cell; ++i)
1321 if (subdomain_association[local_dof_indices[i]] ==
1323 subdomain_association[local_dof_indices[i]] = subdomain_id;
1326 if (coin_flip ==
true)
1327 subdomain_association[local_dof_indices[i]] = subdomain_id;
1328 coin_flip = !coin_flip;
1332 Assert (std::find (subdomain_association.begin(),
1333 subdomain_association.end(),
1335 == subdomain_association.end(),
1336 ExcInternalError());
1341 template <
typename DoFHandlerType>
1346 std::vector<types::subdomain_id> subdomain_association (dof_handler.n_dofs());
1349 return std::count (subdomain_association.begin(),
1350 subdomain_association.end(),
1356 template <
typename DoFHandlerType>
1367 (subdomain == dof_handler.get_triangulation().locally_owned_subdomain()),
1368 ExcMessage (
"For parallel::distributed::Triangulation objects and " 1369 "associated DoF handler objects, asking for any subdomain other " 1370 "than the locally owned one does not make sense."));
1372 IndexSet index_set (dof_handler.n_dofs());
1374 std::vector<types::global_dof_index> local_dof_indices;
1375 local_dof_indices.reserve (max_dofs_per_cell(dof_handler));
1382 std::vector<types::global_dof_index> subdomain_indices;
1384 typename DoFHandlerType::active_cell_iterator
1385 cell = dof_handler.begin_active(),
1386 endc = dof_handler.
end();
1387 for (; cell!=endc; ++cell)
1388 if ((cell->is_artificial() ==
false)
1390 (cell->subdomain_id() == subdomain))
1392 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
1393 local_dof_indices.resize (dofs_per_cell);
1394 cell->get_dof_indices (local_dof_indices);
1395 subdomain_indices.insert(subdomain_indices.end(),
1396 local_dof_indices.begin(),
1397 local_dof_indices.end());
1400 std::sort (subdomain_indices.begin(), subdomain_indices.end());
1401 subdomain_indices.erase (std::unique(subdomain_indices.begin(),
1402 subdomain_indices.end()),
1403 subdomain_indices.end());
1406 index_set.add_indices (subdomain_indices.begin(), subdomain_indices.end());
1407 index_set.compress ();
1414 template <
typename DoFHandlerType>
1418 std::vector<unsigned int> &n_dofs_on_subdomain)
1420 Assert (n_dofs_on_subdomain.size() == dof_handler.get_fe().n_components(),
1421 ExcDimensionMismatch (n_dofs_on_subdomain.size(),
1422 dof_handler.get_fe().n_components()));
1423 std::fill (n_dofs_on_subdomain.begin(), n_dofs_on_subdomain.end(), 0);
1431 DoFHandlerType::space_dimension>::active_cell_iterator
1433 cell!=dof_handler.get_triangulation().end(); ++cell)
1434 if (cell->subdomain_id() == subdomain)
1440 ExcMessage (
"There are no cells for the given subdomain!"));
1444 std::vector<types::subdomain_id> subdomain_association (dof_handler.n_dofs());
1447 std::vector<unsigned char> component_association (dof_handler.n_dofs());
1448 internal::get_component_association (dof_handler, std::vector<bool>(),
1449 component_association);
1451 for (
unsigned int c=0; c<dof_handler.get_fe().n_components(); ++c)
1454 if ((subdomain_association[i] == subdomain) &&
1455 (component_association[i] ==
static_cast<unsigned char>(c)))
1456 ++n_dofs_on_subdomain[c];
1468 template <
int dim,
int spacedim>
1471 const std::vector<unsigned char> &dofs_by_component,
1472 const std::vector<unsigned int> &target_component,
1473 const bool only_once,
1474 std::vector<types::global_dof_index> &dofs_per_component,
1475 unsigned int &component)
1486 resolve_components(base, dofs_by_component, target_component,
1487 only_once, dofs_per_component, component);
1490 for (
unsigned int dd=0; dd<d; ++dd,++component)
1491 dofs_per_component[target_component[component]]
1492 += std::count(dofs_by_component.begin(),
1493 dofs_by_component.end(),
1500 for (
unsigned int dd=1; dd<d; ++dd)
1501 dofs_per_component[target_component[component-d+dd]] =
1502 dofs_per_component[target_component[component-d]];
1509 template <
int dim,
int spacedim>
1512 const std::vector<unsigned char> &dofs_by_component,
1513 const std::vector<unsigned int> &target_component,
1514 const bool only_once,
1515 std::vector<types::global_dof_index> &dofs_per_component,
1516 unsigned int &component)
1521 for (
unsigned int fe=1; fe<fe_collection.
size(); ++fe)
1523 Assert (fe_collection[fe].n_components() == fe_collection[0].n_components(),
1524 ExcNotImplemented());
1525 Assert (fe_collection[fe].n_base_elements() == fe_collection[0].n_base_elements(),
1526 ExcNotImplemented());
1527 for (
unsigned int b=0; b<fe_collection[0].n_base_elements(); ++b)
1529 Assert (fe_collection[fe].base_element(b).n_components() == fe_collection[0].base_element(b).n_components(),
1530 ExcNotImplemented());
1531 Assert (fe_collection[fe].base_element(b).n_base_elements() == fe_collection[0].base_element(b).n_base_elements(),
1532 ExcNotImplemented());
1536 resolve_components (fe_collection[0], dofs_by_component,
1537 target_component, only_once, dofs_per_component,
1551 template <
int dim,
int spacedim>
1561 template <
int dim,
int spacedim>
1562 bool all_elements_are_primitive (const ::hp::FECollection<dim,spacedim> &fe_collection)
1564 for (
unsigned int i=0; i<fe_collection.size(); ++i)
1565 if (fe_collection[i].is_primitive() ==
false)
1573 template <
typename DoFHandlerType>
1576 std::vector<types::global_dof_index> &dofs_per_component,
1578 std::vector<unsigned int> target_component)
1580 const unsigned int n_components = dof_handler.get_fe().n_components();
1582 std::fill (dofs_per_component.begin(), dofs_per_component.end(),
1587 if (target_component.size()==0)
1589 target_component.resize(n_components);
1590 for (
unsigned int i=0; i<n_components; ++i)
1591 target_component[i] = i;
1594 Assert (target_component.size()==n_components,
1595 ExcDimensionMismatch(target_component.size(),
1599 const unsigned int max_component
1600 = *std::max_element (target_component.begin(),
1601 target_component.end());
1602 const unsigned int n_target_components = max_component + 1;
1603 (void)n_target_components;
1609 if (n_components == 1)
1611 dofs_per_component[0] = dof_handler.n_locally_owned_dofs();
1618 std::vector<unsigned char> dofs_by_component (dof_handler.n_locally_owned_dofs());
1619 internal::get_component_association (dof_handler,
ComponentMask(),
1623 unsigned int component = 0;
1624 internal::resolve_components(dof_handler.get_fe(),
1625 dofs_by_component, target_component,
1626 only_once, dofs_per_component, component);
1627 Assert (n_components == component, ExcInternalError());
1631 Assert ((internal::all_elements_are_primitive(dof_handler.get_fe()) ==
false)
1633 (std::accumulate (dofs_per_component.begin(),
1634 dofs_per_component.end(),
1636 == dof_handler.n_locally_owned_dofs()),
1637 ExcInternalError());
1640 #ifdef DEAL_II_WITH_MPI 1641 const unsigned int dim = DoFHandlerType::dimension;
1642 const unsigned int spacedim = DoFHandlerType::space_dimension;
1646 (&dof_handler.get_triangulation())))
1648 std::vector<types::global_dof_index> local_dof_count = dofs_per_component;
1650 MPI_Allreduce ( &local_dof_count[0], &dofs_per_component[0], n_target_components,
1651 DEAL_II_DOF_INDEX_MPI_TYPE,
1652 MPI_SUM, tria->get_communicator());
1659 template <
typename DoFHandlerType>
1662 std::vector<types::global_dof_index> &dofs_per_block,
1663 const std::vector<unsigned int> &target_block_)
1665 std::vector<unsigned int> target_block = target_block_;
1667 const ::hp::FECollection<DoFHandlerType::dimension,DoFHandlerType::space_dimension>
1668 fe_collection (dof_handler.get_fe());
1669 Assert (fe_collection.size() < 256, ExcNotImplemented());
1671 for (
unsigned int this_fe=0; this_fe<fe_collection.size(); ++this_fe)
1674 std::fill (dofs_per_block.begin(), dofs_per_block.end(),
1679 if (target_block.size()==0)
1681 target_block.resize(fe.
n_blocks());
1682 for (
unsigned int i=0; i<fe.
n_blocks(); ++i)
1683 target_block[i] = i;
1687 ExcDimensionMismatch(target_block.size(),
1692 const unsigned int max_block
1693 = *std::max_element (target_block.begin(),
1694 target_block.end());
1695 const unsigned int n_target_blocks = max_block + 1;
1696 (void)n_target_blocks;
1698 const unsigned int n_blocks = fe.
n_blocks();
1706 dofs_per_block[0] = dof_handler.n_dofs();
1710 std::vector<unsigned char> dofs_by_block (dof_handler.n_locally_owned_dofs());
1711 internal::get_block_association (dof_handler, dofs_by_block);
1714 for (
unsigned int block=0; block<fe.
n_blocks(); ++block)
1715 dofs_per_block[target_block[block]]
1716 += std::count(dofs_by_block.begin(), dofs_by_block.end(),
1719 #ifdef DEAL_II_WITH_MPI 1724 (&dof_handler.get_triangulation())))
1726 std::vector<types::global_dof_index> local_dof_count = dofs_per_block;
1727 MPI_Allreduce ( &local_dof_count[0], &dofs_per_block[0],
1729 DEAL_II_DOF_INDEX_MPI_TYPE,
1730 MPI_SUM, tria->get_communicator());
1738 template <
typename DoFHandlerType>
1741 std::vector<types::global_dof_index> &mapping)
1744 mapping.insert (mapping.end(), dof_handler.n_dofs(),
1745 DoFHandlerType::invalid_dof_index);
1747 std::vector<types::global_dof_index> dofs_on_face;
1748 dofs_on_face.reserve (max_dofs_per_face(dof_handler));
1757 typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active(),
1758 endc = dof_handler.end();
1759 for (; cell!=endc; ++cell)
1760 for (
unsigned int f=0; f<GeometryInfo<DoFHandlerType::dimension>::faces_per_cell; ++f)
1761 if (cell->at_boundary(f))
1763 const unsigned int dofs_per_face = cell->get_fe().dofs_per_face;
1764 dofs_on_face.resize (dofs_per_face);
1765 cell->face(f)->get_dof_indices (dofs_on_face,
1766 cell->active_fe_index());
1767 for (
unsigned int i=0; i<dofs_per_face; ++i)
1768 if (mapping[dofs_on_face[i]] == DoFHandlerType::invalid_dof_index)
1769 mapping[dofs_on_face[i]] = next_boundary_index++;
1777 template <
typename DoFHandlerType>
1779 (
const DoFHandlerType &dof_handler,
1780 const std::set<types::boundary_id> &boundary_ids,
1781 std::vector<types::global_dof_index> &mapping)
1784 ExcInvalidBoundaryIndicator());
1787 mapping.insert (mapping.end(), dof_handler.n_dofs(),
1788 DoFHandlerType::invalid_dof_index);
1791 if (boundary_ids.size() == 0)
1794 std::vector<types::global_dof_index> dofs_on_face;
1795 dofs_on_face.reserve (max_dofs_per_face(dof_handler));
1798 typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active(),
1799 endc = dof_handler.end();
1800 for (; cell!=endc; ++cell)
1801 for (
unsigned int f=0; f<GeometryInfo<DoFHandlerType::dimension>::faces_per_cell; ++f)
1802 if (boundary_ids.find (cell->face(f)->boundary_id()) !=
1805 const unsigned int dofs_per_face = cell->get_fe().dofs_per_face;
1806 dofs_on_face.resize (dofs_per_face);
1807 cell->face(f)->get_dof_indices (dofs_on_face, cell->active_fe_index());
1808 for (
unsigned int i=0; i<dofs_per_face; ++i)
1809 if (mapping[dofs_on_face[i]] == DoFHandlerType::invalid_dof_index)
1810 mapping[dofs_on_face[i]] = next_boundary_index++;
1814 dof_handler.n_boundary_dofs (boundary_ids));
1821 template <
typename DoFHandlerType>
1825 const DoFHandlerType &dof_handler,
1828 const unsigned int dim = DoFHandlerType::dimension;
1829 const unsigned int spacedim = DoFHandlerType::space_dimension;
1834 for (
unsigned int fe_index = 0; fe_index < fe_collection.size(); ++fe_index)
1837 Assert(fe_collection[fe_index].has_support_points(),
1841 fe_collection[fe_index].get_unit_support_points()));
1853 typename DoFHandlerType::active_cell_iterator cell =
1854 dof_handler.begin_active(), endc = dof_handler.end();
1856 std::vector<types::global_dof_index> local_dof_indices;
1857 for (; cell != endc; ++cell)
1859 if (cell->is_artificial() ==
false)
1861 hp_fe_values.
reinit(cell);
1864 local_dof_indices.resize(cell->get_fe().dofs_per_cell);
1865 cell->get_dof_indices(local_dof_indices);
1867 const std::vector<Point<spacedim> > &points =
1869 for (
unsigned int i = 0; i < cell->get_fe().dofs_per_cell; ++i)
1871 support_points[local_dof_indices[i]] = points[i];
1876 template <
typename DoFHandlerType>
1880 const DoFHandlerType &dof_handler,
1884 std::map<types::global_dof_index,Point<DoFHandlerType::space_dimension> > x_support_points;
1891 Assert (x_support_points.find(i) != x_support_points.end(),
1892 ExcInternalError());
1893 support_points[i] = x_support_points[i];
1899 template <
int dim,
int spacedim>
1910 ExcMessage (
"This function can not be used with distributed triangulations." 1911 "See the documentation for more information."));
1917 internal::map_dofs_to_support_points (mapping_collection,
1923 template<
int dim,
int spacedim>
1934 ExcMessage (
"This function can not be used with distributed triangulations." 1935 "See the documentation for more information."));
1939 internal::map_dofs_to_support_points (mapping,
1945 template <
int dim,
int spacedim>
1951 support_points.clear();
1957 internal::map_dofs_to_support_points (mapping_collection,
1963 template<
int dim,
int spacedim>
1969 support_points.clear();
1973 internal::map_dofs_to_support_points (mapping,
1979 template<
int dim,
int spacedim>
1987 const unsigned int nb = fe.
n_blocks();
1989 tables_by_block.resize(1);
1990 tables_by_block[0].reinit(nb, nb);
1991 tables_by_block[0].fill(
none);
1999 tables_by_block[0](ib,jb) |= table(i,j);
2005 template<
int dim,
int spacedim>
2012 tables_by_block.resize(fe_collection.
size());
2014 for (
unsigned int f=0; f<fe_collection.
size(); ++f)
2018 const unsigned int nb = fe.
n_blocks();
2019 tables_by_block[f].reinit(nb, nb);
2020 tables_by_block[f].fill(
none);
2027 tables_by_block[f](ib,jb) |= table(i,j);
2035 template <
typename DoFHandlerType,
class Sparsity>
2037 const DoFHandlerType &dof_handler,
2038 const unsigned int level,
2039 const std::vector<bool> &selected_dofs,
2042 typename DoFHandlerType::level_cell_iterator cell;
2043 typename DoFHandlerType::level_cell_iterator endc = dof_handler.end(level);
2044 std::vector<types::global_dof_index> indices;
2047 for (cell=dof_handler.begin(level); cell != endc; ++i, ++cell)
2049 indices.resize(cell->get_fe().dofs_per_cell);
2050 cell->get_mg_dof_indices(indices);
2052 if (selected_dofs.size()!=0)
2057 if (selected_dofs.size() == 0)
2058 block_list.add(i,indices[j]-offset);
2061 if (selected_dofs[j])
2062 block_list.add(i,indices[j]-offset);
2069 template <
typename DoFHandlerType>
2071 const DoFHandlerType &dof_handler,
2072 const unsigned int level,
2073 const bool interior_only)
2076 block_list.
reinit(1, dof_handler.n_dofs(level), dof_handler.n_dofs(level));
2077 typename DoFHandlerType::level_cell_iterator cell;
2078 typename DoFHandlerType::level_cell_iterator endc = dof_handler.end(level);
2080 std::vector<types::global_dof_index> indices;
2081 std::vector<bool> exclude;
2083 for (cell=dof_handler.begin(level); cell != endc; ++cell)
2085 indices.resize(cell->get_fe().dofs_per_cell);
2086 cell->get_mg_dof_indices(indices);
2092 std::fill(exclude.begin(), exclude.end(),
false);
2095 for (
unsigned int face=0; face<GeometryInfo<DoFHandlerType::dimension>::faces_per_cell; ++face)
2096 if (cell->at_boundary(face) || cell->neighbor(face)->level() != cell->level())
2097 for (
unsigned int i=0; i<dpf; ++i)
2101 block_list.
add(0, indices[j]);
2106 block_list.
add(0, indices[j]);
2112 template <
typename DoFHandlerType>
2114 const DoFHandlerType &dof_handler,
2115 const unsigned int level,
2116 const bool interior_dofs_only,
2117 const bool boundary_dofs)
2119 Assert(level > 0 && level < dof_handler.get_triangulation().n_levels(),
2120 ExcIndexRange(level, 1, dof_handler.get_triangulation().n_levels()));
2122 typename DoFHandlerType::level_cell_iterator pcell = dof_handler.begin(level-1);
2123 typename DoFHandlerType::level_cell_iterator endc = dof_handler.end(level-1);
2125 std::vector<types::global_dof_index> indices;
2126 std::vector<bool> exclude;
2128 for (
unsigned int block = 0; pcell != endc; ++pcell)
2130 if (!pcell->has_children())
2133 for (
unsigned int child=0; child<pcell->n_children(); ++child)
2135 const typename DoFHandlerType::level_cell_iterator cell = pcell->child(child);
2140 indices.resize(n_dofs);
2141 exclude.resize(n_dofs);
2142 std::fill(exclude.begin(), exclude.end(),
false);
2143 cell->get_mg_dof_indices(indices);
2145 if (interior_dofs_only)
2151 for (
unsigned int d=0; d<DoFHandlerType::dimension; ++d)
2154 for (
unsigned int i=0; i<dpf; ++i)
2161 for (
unsigned int face=0; face< GeometryInfo<DoFHandlerType::dimension>::faces_per_cell; ++face)
2162 if (cell->at_boundary(face))
2163 for (
unsigned int i=0; i<dpf; ++i)
2167 for (
unsigned int i=0; i<n_dofs; ++i)
2169 block_list.
add(block, indices[i]);
2177 template <
typename DoFHandlerType>
2179 const DoFHandlerType &dof_handler,
2180 const unsigned int level,
2181 const bool interior_only,
2182 const bool boundary_patches,
2183 const bool level_boundary_patches,
2184 const bool single_cell_patches)
2186 typename DoFHandlerType::level_cell_iterator cell;
2187 typename DoFHandlerType::level_cell_iterator endc = dof_handler.end(level);
2191 std::vector<unsigned int> vertex_cell_count(dof_handler.get_triangulation().n_vertices(), 0);
2194 std::vector<bool> vertex_boundary(dof_handler.get_triangulation().n_vertices(),
false);
2196 std::vector<unsigned int> vertex_mapping(dof_handler.get_triangulation().n_vertices(),
2200 std::vector<unsigned int> vertex_dof_count(dof_handler.get_triangulation().n_vertices(), 0);
2204 for (cell=dof_handler.begin(level); cell != endc; ++cell)
2205 for (
unsigned int v=0; v<GeometryInfo<DoFHandlerType::dimension>::vertices_per_cell; ++v)
2207 const unsigned int vg = cell->vertex_index(v);
2208 vertex_dof_count[vg] += cell->get_fe().dofs_per_cell;
2209 ++vertex_cell_count[vg];
2210 for (
unsigned int d=0; d<DoFHandlerType::dimension; ++d)
2213 if (cell->at_boundary(face))
2214 vertex_boundary[vg] =
true;
2215 else if ((!level_boundary_patches)
2216 && (cell->neighbor(face)->level() != (int) level))
2217 vertex_boundary[vg] =
true;
2223 for (
unsigned int vg=0; vg<vertex_dof_count.size(); ++vg)
2224 if ((!single_cell_patches && vertex_cell_count[vg] < 2)
2226 (!boundary_patches && vertex_boundary[vg]))
2227 vertex_dof_count[vg] = 0;
2230 unsigned int n_vertex_count=0;
2231 for (
unsigned int vg=0; vg<vertex_mapping.size(); ++vg)
2232 if (vertex_dof_count[vg] != 0)
2233 vertex_mapping[vg] = n_vertex_count++;
2236 for (
unsigned int vg=0; vg<vertex_mapping.size(); ++vg)
2237 if (vertex_dof_count[vg] != 0)
2238 vertex_dof_count[vertex_mapping[vg]] = vertex_dof_count[vg];
2242 vertex_dof_count.resize(n_vertex_count);
2246 block_list.
reinit(vertex_dof_count.size(), dof_handler.n_dofs(level), vertex_dof_count);
2248 std::vector<types::global_dof_index> indices;
2249 std::vector<bool> exclude;
2251 for (cell=dof_handler.begin(level); cell != endc; ++cell)
2255 cell->get_mg_dof_indices(indices);
2257 for (
unsigned int v=0; v<GeometryInfo<DoFHandlerType::dimension>::vertices_per_cell; ++v)
2259 const unsigned int vg = cell->vertex_index(v);
2260 const unsigned int block = vertex_mapping[vg];
2269 std::fill(exclude.begin(), exclude.end(),
false);
2272 for (
unsigned int d=0; d<DoFHandlerType::dimension; ++d)
2276 for (
unsigned int i=0; i<dpf; ++i)
2279 for (
unsigned int j=0; j<indices.size(); ++j)
2281 block_list.
add(block, indices[j]);
2285 for (
unsigned int j=0; j<indices.size(); ++j)
2286 block_list.
add(block, indices[j]);
2294 template <
typename DoFHandlerType>
2298 std::set<types::global_dof_index> dofs_on_patch;
2299 std::vector<types::global_dof_index> local_dof_indices;
2304 for (
unsigned int i=0; i<patch.size(); ++i)
2306 const typename DoFHandlerType::active_cell_iterator cell = patch[i];
2307 Assert (cell->is_artificial() ==
false,
2308 ExcMessage(
"This function can not be called with cells that are " 2309 "not either locally owned or ghost cells."));
2310 local_dof_indices.resize (cell->get_fe().dofs_per_cell);
2311 cell->get_dof_indices (local_dof_indices);
2312 dofs_on_patch.insert (local_dof_indices.begin(),
2313 local_dof_indices.end());
2317 return dofs_on_patch.size();
2322 template <
typename DoFHandlerType>
2323 std::vector<types::global_dof_index>
2326 std::set<types::global_dof_index> dofs_on_patch;
2327 std::vector<types::global_dof_index> local_dof_indices;
2332 for (
unsigned int i=0; i<patch.size(); ++i)
2334 const typename DoFHandlerType::active_cell_iterator cell = patch[i];
2335 Assert (cell->is_artificial() ==
false,
2336 ExcMessage(
"This function can not be called with cells that are " 2337 "not either locally owned or ghost cells."));
2338 local_dof_indices.resize (cell->get_fe().dofs_per_cell);
2339 cell->get_dof_indices (local_dof_indices);
2340 dofs_on_patch.insert (local_dof_indices.begin(),
2341 local_dof_indices.end());
2344 Assert (dofs_on_patch.size() == count_dofs_on_patch<DoFHandlerType>(patch),
2345 ExcInternalError());
2351 return std::vector<types::global_dof_index> (dofs_on_patch.begin(),
2352 dofs_on_patch.end());
2362 #include "dof_tools.inst" 2366 DEAL_II_NAMESPACE_CLOSE
unsigned int n_active_cells() const
static const unsigned int invalid_unsigned_int
#define AssertDimension(dim1, dim2)
ElementIterator end() const
const types::subdomain_id invalid_subdomain_id
void add_index(const size_type index)
::ExceptionBase & ExcMessage(std::string arg1)
#define AssertIndexRange(index, range)
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
active_cell_iterator begin_active(const unsigned int level=0) const
bool represents_n_components(const unsigned int n) const
bool represents_the_all_selected_mask() const
unsigned int size() const
Transformed quadrature points.
const FiniteElement< dim, spacedim > & get_fe() const
const unsigned int dofs_per_line
unsigned int n_blocks() const
const ::FEValues< dim, spacedim > & get_present_fe_values() const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const hp::FECollection< dim, spacedim > & get_fe() const
unsigned int component_to_block_index(const unsigned int component) const
unsigned int n_locally_owned_dofs() const
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType, lda > > cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
void push_back(const FiniteElement< dim, spacedim > &new_fe)
void add(const size_type i, const size_type j)
unsigned int global_dof_index
#define Assert(cond, exc)
unsigned int element_multiplicity(const unsigned int index) const
bool is_primitive(const unsigned int i) const
Abstract base class for mapping classes.
const ComponentMask & get_nonzero_components(const unsigned int i) const
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index) const
types::global_dof_index n_dofs() const
unsigned int size() 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
unsigned int subdomain_id
const unsigned int dofs_per_cell
void reinit(const size_type m, const size_type n, const unsigned int max_per_row)
void set_size(const size_type size)
const types::subdomain_id artificial_subdomain_id
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
const Triangulation< dim, spacedim > & get_triangulation() 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 first_selected_component(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
const unsigned int dofs_per_face
types::global_dof_index n_dofs() const
const Triangulation< dim, spacedim > & get_triangulation() const
const unsigned int dofs_per_vertex
bool is_element(const size_type index) const
const types::boundary_id internal_face_boundary_id
unsigned int n_base_elements() const
ElementIterator begin() const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Triangulation< dim, spacedim > & get_triangulation()