36 #include <deal.II/multigrid/mg_transfer_component.templates.h>
81 template <
int dim,
typename number,
int spacedim>
83 reinit_vector_by_components(
84 const ::DoFHandler<dim, spacedim> & mg_dof,
86 const std::vector<bool> & sel,
87 const std::vector<unsigned int> & target_comp,
88 std::vector<std::vector<types::global_dof_index>> &ndofs)
90 std::vector<bool> selected = sel;
91 std::vector<unsigned int> target_component = target_comp;
92 const unsigned int ncomp = mg_dof.get_fe(0).n_components();
103 if (target_component.size() == 0)
105 target_component.resize(ncomp);
106 for (
unsigned int i = 0; i < ncomp; ++i)
107 target_component[i] = i;
112 if (selected.size() == 0)
114 selected.resize(target_component.size());
115 std::fill_n(selected.begin(), ncomp,
false);
116 for (
const unsigned int component : target_component)
117 selected[component] =
true;
120 Assert(selected.size() == target_component.size(),
124 const unsigned int n_selected =
125 std::accumulate(selected.begin(), selected.end(), 0u);
127 if (ndofs.size() == 0)
129 std::vector<std::vector<types::global_dof_index>> new_dofs(
130 mg_dof.get_triangulation().n_levels(),
131 std::vector<types::global_dof_index>(target_component.size()));
136 for (
unsigned int level = v.min_level();
level <= v.max_level(); ++
level)
138 v[
level].reinit(n_selected, 0);
140 for (
unsigned int i = 0;
148 v[
level].collect_sizes();
180 template <
int dim,
typename number,
int spacedim>
182 reinit_vector_by_components(
183 const ::DoFHandler<dim, spacedim> & mg_dof,
186 const std::vector<unsigned int> & target_component,
187 std::vector<std::vector<types::global_dof_index>> &ndofs)
190 ExcMessage(
"The component mask does not have the correct size."));
192 unsigned int selected_block = 0;
193 for (
unsigned int i = 0; i < target_component.size(); ++i)
194 if (component_mask[i])
195 selected_block = target_component[i];
197 if (ndofs.size() == 0)
199 std::vector<std::vector<types::global_dof_index>> new_dofs(
200 mg_dof.get_triangulation().n_levels(),
201 std::vector<types::global_dof_index>(target_component.size()));
206 for (
unsigned int level = v.min_level();
level <= v.max_level(); ++
level)
208 v[
level].reinit(ndofs[
level][selected_block]);
214 template <
typename number>
215 template <
int dim,
class InVector,
int spacedim>
220 const InVector & src)
const
225 ExcMatricesNotBuilt());
227 reinit_vector_by_components(
228 mg_dof_handler, dst, mg_component_mask, mg_target_component, sizes);
246 using IT = std::vector<
247 std::pair<types::global_dof_index, unsigned int>>::const_iterator;
248 for (IT i = copy_to_and_from_indices[
level].
begin();
249 i != copy_to_and_from_indices[
level].end();
251 dst[
level](i->second) = src(i->first);
257 template <
int dim,
int spacedim>
267 template <
int dim,
int spacedim>
285 mg_dof.
get_fe(0).n_components()));
305 mg_dof.
get_fe(0).n_components()));
327 ExcMessage(
"Component mask has wrong size."));
330 sizes.resize(n_levels);
331 for (
unsigned int l = 0;
l < n_levels; ++
l)
346 level_component_start = k;
387 for (
unsigned int i = 0; i < n_levels - 1; ++i)
396 std::vector<types::global_dof_index> dof_indices_parent(dofs_per_cell);
397 std::vector<types::global_dof_index> dof_indices_child(dofs_per_cell);
432 if (cell->has_children())
434 cell->get_mg_dof_indices(dof_indices_parent);
436 for (
unsigned int child = 0; child < cell->n_children(); ++child)
442 mg_dof.
get_fe().get_prolongation_matrix(
443 child, cell->refinement_case());
445 cell->child(child)->get_mg_dof_indices(dof_indices_child);
450 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
451 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
452 if (prolongation(i, j) != 0)
454 const unsigned int icomp =
456 const unsigned int jcomp =
460 dof_indices_child[i], dof_indices_parent[j]);
472 if (cell->has_children())
474 cell->get_mg_dof_indices(dof_indices_parent);
476 for (
unsigned int child = 0; child < cell->n_children(); ++child)
482 mg_dof.
get_fe().get_prolongation_matrix(
483 child, cell->refinement_case());
485 cell->child(child)->get_mg_dof_indices(dof_indices_child);
489 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
490 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
491 if (prolongation(i, j) != 0)
493 const unsigned int icomp =
495 const unsigned int jcomp =
499 dof_indices_child[i],
500 dof_indices_parent[j],
513 std::vector<std::vector<types::global_dof_index>> dofs_per_component(
525 for (
unsigned int iblock = 0; iblock <
n_components; ++iblock)
526 for (
unsigned int jblock = 0; jblock <
n_components; ++jblock)
527 if (iblock == jblock)
535 ->block(iblock, jblock)
538 ->block(iblock, jblock)
547 dofs_per_component[
level]);
552 std::set<types::global_dof_index>::const_iterator
555 const bool is_boundary_index =
558 if (is_boundary_index)
561 ->block(iblock, jblock)
562 .set(i, column_number, 0);
573 template <
typename number>
574 template <
int dim,
int spacedim>
580 unsigned int mg_select,
581 const std::vector<unsigned int> & t_component,
582 const std::vector<unsigned int> & mg_t_component,
583 const std::vector<std::set<types::global_dof_index>> &bdry_indices)
585 build(mg_dof, select, mg_select, t_component, mg_t_component, bdry_indices);
590 template <
typename number>
591 template <
int dim,
int spacedim>
596 unsigned int mg_select,
597 const std::vector<unsigned int> & t_component,
598 const std::vector<unsigned int> & mg_t_component,
599 const std::vector<std::set<types::global_dof_index>> &bdry_indices)
602 unsigned int ncomp = mg_dof.
get_fe(0).n_components();
604 target_component = t_component;
605 mg_target_component = mg_t_component;
606 boundary_indices = bdry_indices;
608 selected_component = select;
609 mg_selected_component = mg_select;
612 std::vector<bool> tmp(ncomp,
false);
613 for (
unsigned int c = 0; c < ncomp; ++c)
614 if (t_component[c] == selected_component)
620 std::vector<bool> tmp(ncomp,
false);
621 for (
unsigned int c = 0; c < ncomp; ++c)
622 if (mg_t_component[c] == mg_selected_component)
631 for (
unsigned int i = 0; i < target_component.size(); ++i)
633 if (target_component[i] == select)
635 selected_component = i;
640 for (
unsigned int i = 0; i < mg_target_component.size(); ++i)
642 if (mg_target_component[i] == mg_select)
644 mg_selected_component = i;
654 interface_dofs[
l].clear();
655 interface_dofs[
l].set_size(mg_dof.
n_dofs(
l));
661 std::vector<types::global_dof_index> temp_copy_indices;
662 std::vector<types::global_dof_index> global_dof_indices(fe.
dofs_per_cell);
663 std::vector<types::global_dof_index> level_dof_indices(fe.
dofs_per_cell);
667 copy_to_and_from_indices[
level].clear();
673 temp_copy_indices.resize(0);
679 for (; level_cell != level_end; ++level_cell)
685 level_cell->get_dof_indices(global_dof_indices);
686 level_cell->get_mg_dof_indices(level_dof_indices);
690 const unsigned int component =
692 if (component_mask[component] &&
693 !interface_dofs[
level].is_element(level_dof_indices[i]))
696 mg_component_start[
level][mg_target_component[component]];
698 component_start[target_component[component]];
699 temp_copy_indices[level_dof_indices[i] - level_start] =
700 global_dof_indices[i] - global_start;
708 std::count_if(temp_copy_indices.begin(),
709 temp_copy_indices.end(),
711 return index != numbers::invalid_dof_index;
713 copy_to_and_from_indices[
level].resize(n_active_dofs);
717 copy_to_and_from_indices[
level][counter++] =
718 std::pair<types::global_dof_index, unsigned int>(
719 temp_copy_indices[i], i);
727 #include "mg_transfer_component.inst"