36#include <deal.II/multigrid/mg_transfer_component.templates.h>
79 template <
int dim,
typename number,
int spacedim>
81 reinit_vector_by_components(
84 const std::vector<bool> & sel,
85 const std::vector<unsigned int> & target_comp,
86 std::vector<std::vector<types::global_dof_index>> &ndofs)
88 std::vector<bool> selected = sel;
89 std::vector<unsigned int> target_component = target_comp;
101 if (target_component.size() == 0)
103 target_component.resize(ncomp);
104 for (
unsigned int i = 0; i < ncomp; ++i)
105 target_component[i] = i;
110 if (selected.size() == 0)
112 selected.resize(target_component.size());
113 std::fill_n(selected.begin(), ncomp,
false);
114 for (
const unsigned int component : target_component)
115 selected[component] = true;
118 Assert(selected.size() == target_component.size(),
122 const unsigned int n_selected =
123 std::accumulate(selected.begin(), selected.end(), 0u);
125 if (ndofs.size() == 0)
127 std::vector<std::vector<types::global_dof_index>> new_dofs(
129 std::vector<types::global_dof_index>(target_component.size()));
130 std::swap(ndofs, new_dofs);
134 for (
unsigned int level = v.min_level();
level <= v.max_level(); ++
level)
136 v[
level].reinit(n_selected, 0);
138 for (
unsigned int i = 0;
139 i < selected.size() && (k < v[
level].n_blocks());
146 v[
level].collect_sizes();
178 template <
int dim,
typename number,
int spacedim>
180 reinit_vector_by_components(
184 const std::vector<unsigned int> & target_component,
185 std::vector<std::vector<types::global_dof_index>> &ndofs)
188 ExcMessage(
"The component mask does not have the correct size."));
190 unsigned int selected_block = 0;
191 for (
unsigned int i = 0; i < target_component.size(); ++i)
192 if (component_mask[i])
193 selected_block = target_component[i];
195 if (ndofs.size() == 0)
197 std::vector<std::vector<types::global_dof_index>> new_dofs(
199 std::vector<types::global_dof_index>(target_component.size()));
200 std::swap(ndofs, new_dofs);
204 for (
unsigned int level = v.min_level();
level <= v.max_level(); ++
level)
206 v[
level].reinit(ndofs[
level][selected_block]);
212template <
typename number>
213template <
int dim,
class InVector,
int spacedim>
218 const InVector & src)
const
223 ExcMatricesNotBuilt());
225 reinit_vector_by_components(
226 mg_dof_handler, dst, mg_component_mask, mg_target_component, sizes);
244 using IT = std::vector<
245 std::pair<types::global_dof_index, unsigned int>>::const_iterator;
246 for (IT i = copy_to_and_from_indices[
level].begin();
247 i != copy_to_and_from_indices[
level].end();
249 dst[
level](i->second) = src(i->first);
255template <
int dim,
int spacedim>
308 const unsigned int n_components =
315 ExcMessage(
"Component mask has wrong size."));
318 sizes.resize(n_levels);
319 for (
unsigned int l = 0; l < n_levels; ++l)
320 sizes[l].resize(n_components);
334 level_component_start = k;
375 for (
unsigned int i = 0; i < n_levels - 1; ++i)
384 std::vector<types::global_dof_index> dof_indices_parent(dofs_per_cell);
385 std::vector<types::global_dof_index> dof_indices_child(dofs_per_cell);
405 for (
unsigned int i = 0; i < n_components; ++i)
406 for (
unsigned int j = 0; j < n_components; ++j)
420 if (cell->has_children())
422 cell->get_mg_dof_indices(dof_indices_parent);
424 for (
unsigned int child = 0; child < cell->n_children(); ++child)
431 child, cell->refinement_case());
433 cell->child(child)->get_mg_dof_indices(dof_indices_child);
438 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
439 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
440 if (prolongation(i, j) != 0)
442 const unsigned int icomp =
444 const unsigned int jcomp =
448 dof_indices_child[i], dof_indices_parent[j]);
460 if (cell->has_children())
462 cell->get_mg_dof_indices(dof_indices_parent);
464 for (
unsigned int child = 0; child < cell->n_children(); ++child)
471 child, cell->refinement_case());
473 cell->child(child)->get_mg_dof_indices(dof_indices_child);
477 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
478 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
479 if (prolongation(i, j) != 0)
481 const unsigned int icomp =
483 const unsigned int jcomp =
487 dof_indices_child[i],
488 dof_indices_parent[j],
501 std::vector<std::vector<types::global_dof_index>> dofs_per_component(
503 std::vector<types::global_dof_index>(n_components));
513 for (
unsigned int iblock = 0; iblock < n_components; ++iblock)
514 for (
unsigned int jblock = 0; jblock < n_components; ++jblock)
515 if (iblock == jblock)
523 ->block(iblock, jblock)
526 ->block(iblock, jblock)
528 for (; begin != end; ++begin)
535 dofs_per_component[
level]);
540 std::set<types::global_dof_index>::const_iterator
543 const bool is_boundary_index =
546 if (is_boundary_index)
549 ->block(iblock, jblock)
550 .set(i, column_number, 0);
561template <
typename number>
562template <
int dim,
int spacedim>
567 unsigned int mg_select,
568 const std::vector<unsigned int> & t_component,
569 const std::vector<unsigned int> & mg_t_component,
570 const std::vector<std::set<types::global_dof_index>> &bdry_indices)
575 target_component = t_component;
576 mg_target_component = mg_t_component;
577 boundary_indices = bdry_indices;
579 selected_component = select;
580 mg_selected_component = mg_select;
583 std::vector<bool> tmp(ncomp,
false);
584 for (
unsigned int c = 0; c < ncomp; ++c)
585 if (t_component[c] == selected_component)
591 std::vector<bool> tmp(ncomp,
false);
592 for (
unsigned int c = 0; c < ncomp; ++c)
593 if (mg_t_component[c] == mg_selected_component)
602 for (
unsigned int i = 0; i < target_component.size(); ++i)
604 if (target_component[i] == select)
606 selected_component = i;
611 for (
unsigned int i = 0; i < mg_target_component.size(); ++i)
613 if (mg_target_component[i] == mg_select)
615 mg_selected_component = i;
625 interface_dofs[l].clear();
626 interface_dofs[l].set_size(mg_dof.
n_dofs(l));
632 std::vector<types::global_dof_index> temp_copy_indices;
633 std::vector<types::global_dof_index> global_dof_indices(fe.
n_dofs_per_cell());
634 std::vector<types::global_dof_index> level_dof_indices(fe.
n_dofs_per_cell());
638 copy_to_and_from_indices[
level].clear();
644 temp_copy_indices.resize(0);
650 for (; level_cell != level_end; ++level_cell)
656 level_cell->get_dof_indices(global_dof_indices);
657 level_cell->get_mg_dof_indices(level_dof_indices);
661 const unsigned int component =
663 if (component_mask[component] &&
664 !interface_dofs[
level].is_element(level_dof_indices[i]))
667 mg_component_start[
level][mg_target_component[component]];
669 component_start[target_component[component]];
670 temp_copy_indices[level_dof_indices[i] - level_start] =
671 global_dof_indices[i] - global_start;
679 std::count_if(temp_copy_indices.begin(),
680 temp_copy_indices.end(),
682 return index != numbers::invalid_dof_index;
684 copy_to_and_from_indices[
level].resize(n_active_dofs);
688 copy_to_and_from_indices[
level][counter++] =
689 std::pair<types::global_dof_index, unsigned int>(
690 temp_copy_indices[i], i);
698#include "mg_transfer_component.inst"
size_type local_to_global(const unsigned int block, const size_type index) const
bool represents_n_components(const unsigned int n) const
cell_iterator end() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
active_cell_iterator begin_active(const unsigned int level=0) const
active_cell_iterator end_active(const unsigned int level) const
types::global_dof_index n_dofs() const
cell_iterator begin(const unsigned int level=0) const
unsigned int n_dofs_per_cell() const
unsigned int n_components() const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
std::vector< unsigned int > mg_target_component
ComponentMask mg_component_mask
std::vector< unsigned int > target_component
std::vector< types::global_dof_index > component_start
std::vector< std::shared_ptr< BlockSparseMatrix< double > > > prolongation_matrices
std::vector< std::vector< std::pair< types::global_dof_index, unsigned int > > > copy_to_and_from_indices
std::vector< std::vector< types::global_dof_index > > mg_component_start
std::vector< std::shared_ptr< BlockSparsityPattern > > prolongation_sparsities
void build(const DoFHandler< dim, spacedim > &dof_handler)
std::vector< std::vector< types::global_dof_index > > sizes
std::vector< std::set< types::global_dof_index > > boundary_indices
void do_copy_to_mg(const DoFHandler< dim, spacedim > &mg_dof, MGLevelObject< Vector< number > > &dst, const InVector &src) const
void build(const DoFHandler< dim, spacedim > &dof, unsigned int selected, unsigned int mg_selected, const std::vector< unsigned int > &target_component=std::vector< unsigned int >(), const std::vector< unsigned int > &mg_target_component=std::vector< unsigned int >(), const std::vector< std::set< types::global_dof_index > > &boundary_indices=std::vector< std::set< types::global_dof_index > >())
unsigned int n_levels() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::cell_iterator cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
const types::global_dof_index invalid_dof_index