33#include <deal.II/multigrid/mg_transfer_block.templates.h>
51 template <
int dim,
typename number,
int spacedim>
53 reinit_vector_by_blocks(
56 const std::vector<bool> &sel,
57 std::vector<std::vector<types::global_dof_index>> &ndofs)
59 std::vector<bool> selected = sel;
61 const unsigned int n_selected =
62 std::accumulate(selected.begin(), selected.end(), 0u);
66 std::vector<std::vector<types::global_dof_index>> new_dofs(
68 std::vector<types::global_dof_index>(selected.size()));
69 std::swap(ndofs, new_dofs);
73 for (
unsigned int level = v.min_level();
level <= v.max_level(); ++
level)
75 v[
level].reinit(n_selected, 0);
77 for (
unsigned int i = 0;
78 i < selected.size() && (k < v[
level].n_blocks());
85 v[
level].collect_sizes();
97 template <
int dim,
typename number,
int spacedim>
99 reinit_vector_by_blocks(
102 const unsigned int selected_block,
103 std::vector<std::vector<types::global_dof_index>> &ndofs)
108 std::vector<bool> selected(n_blocks,
false);
109 selected[selected_block] =
true;
113 std::vector<std::vector<types::global_dof_index>> new_dofs(
115 std::vector<types::global_dof_index>(selected.size()));
116 std::swap(ndofs, new_dofs);
120 for (
unsigned int level = v.min_level();
level <= v.max_level(); ++
level)
122 v[
level].reinit(ndofs[
level][selected_block]);
128template <
typename number>
129template <
int dim,
typename number2,
int spacedim>
136 reinit_vector_by_blocks(dof_handler, dst, selected_block, sizes);
145 for (IT i = copy_indices[selected_block][
level].begin();
146 i != copy_indices[selected_block][
level].end();
148 dst[
level](i->second) = src.
block(selected_block)(i->first);
154template <
typename number>
155template <
int dim,
typename number2,
int spacedim>
162 reinit_vector_by_blocks(dof_handler, dst, selected_block, sizes);
170 for (IT i = copy_indices[selected_block][
level].begin();
171 i != copy_indices[selected_block][
level].end();
173 dst[
level](i->second) = src(i->first);
179template <
typename number>
180template <
int dim,
typename number2,
int spacedim>
187 reinit_vector_by_blocks(dof_handler, dst, selected, sizes);
192 for (
unsigned int block = 0; block < selected.size(); ++block)
194 for (IT i = copy_indices[block][
level].begin();
195 i != copy_indices[block][
level].end();
197 dst[
level].block(mg_block[block])(i->second) =
198 src.
block(block)(i->first);
204template <
int dim,
int spacedim>
209 const unsigned int n_blocks = fe.
n_blocks();
221 for (
unsigned int i = 0; i < n_blocks; ++i)
229 sizes.resize(n_levels, std::vector<types::global_dof_index>(fe.
n_blocks()));
242 level_block_start = k;
264 for (
unsigned int block = 0; block < n_blocks; ++block)
285 for (
unsigned int i = 0; i < n_levels - 1; ++i)
294 std::vector<types::global_dof_index> dof_indices_parent(dofs_per_cell);
295 std::vector<types::global_dof_index> dof_indices_child(dofs_per_cell);
316 for (
unsigned int i = 0; i < n_blocks; ++i)
317 for (
unsigned int j = 0; j < n_blocks; ++j)
331 if (cell->has_children())
333 cell->get_mg_dof_indices(dof_indices_parent);
335 Assert(cell->n_children() ==
338 for (
unsigned int child = 0; child < cell->n_children(); ++child)
345 child, cell->refinement_case());
347 cell->child(child)->get_mg_dof_indices(dof_indices_child);
352 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
353 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
354 if (prolongation(i, j) != 0)
356 const unsigned int icomp =
358 const unsigned int jcomp =
360 if ((icomp == jcomp) &&
selected[icomp])
362 dof_indices_child[i], dof_indices_parent[j]);
374 if (cell->has_children())
376 cell->get_mg_dof_indices(dof_indices_parent);
378 Assert(cell->n_children() ==
381 for (
unsigned int child = 0; child < cell->n_children(); ++child)
388 child, cell->refinement_case());
390 cell->child(child)->get_mg_dof_indices(dof_indices_child);
394 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
395 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
396 if (prolongation(i, j) != 0)
398 const unsigned int icomp =
400 const unsigned int jcomp =
402 if ((icomp == jcomp) &&
selected[icomp])
404 dof_indices_child[i],
405 dof_indices_parent[j],
417 std::vector<types::global_dof_index> constrain_indices;
418 std::vector<std::vector<bool>> constraints_per_block(n_blocks);
430 constrain_indices.resize(0);
435 for (; dof != endd; ++dof)
436 constrain_indices[*dof] = 1;
438 unsigned int index = 0;
439 for (
unsigned int block = 0; block < n_blocks; ++block)
443 constraints_per_block[block].resize(0);
444 constraints_per_block[block].resize(n_dofs,
false);
446 constraints_per_block[block][i] =
447 (constrain_indices[index] == 1);
453 ->block(block, block)
457 for (; start_row != end_row; ++start_row)
459 if (constraints_per_block[block][start_row->column()])
460 start_row->value() = 0.;
468template <
typename number>
469template <
int dim,
int spacedim>
478 selected_block = select;
479 selected.resize(n_blocks,
false);
480 selected[select] =
true;
484 std::vector<types::global_dof_index> temp_copy_indices;
485 std::vector<types::global_dof_index> global_dof_indices(fe.
n_dofs_per_cell());
486 std::vector<types::global_dof_index> level_dof_indices(fe.
n_dofs_per_cell());
496 temp_copy_indices.resize(0);
497 temp_copy_indices.resize(sizes[
level][selected_block],
502 for (; level_cell != level_end; ++level_cell)
508 level_cell->get_dof_indices(global_dof_indices);
509 level_cell->get_mg_dof_indices(level_dof_indices);
516 if (mg_constrained_dofs !=
nullptr)
518 if (!mg_constrained_dofs->at_refinement_edge(
519 level, level_dof_indices[i]))
520 temp_copy_indices[level_dof_indices[i] -
521 mg_block_start[
level][block]] =
522 global_dof_indices[i] - block_start[block];
525 temp_copy_indices[level_dof_indices[i] -
526 mg_block_start[
level][block]] =
527 global_dof_indices[i] - block_start[block];
539 std::count_if(temp_copy_indices.begin(),
540 temp_copy_indices.end(),
542 return index != numbers::invalid_dof_index;
544 copy_indices[selected_block][
level].resize(n_active_dofs);
548 copy_indices[selected_block][
level][counter++] =
549 std::pair<types::global_dof_index, unsigned int>(
550 temp_copy_indices[i], i);
556template <
typename number>
557template <
int dim,
int spacedim>
560 const std::vector<bool> &sel)
567 Assert(sel.size() == n_blocks,
571 if (selected.empty())
572 selected = std::vector<bool>(n_blocks,
true);
576 std::vector<std::vector<types::global_dof_index>> temp_copy_indices(n_blocks);
577 std::vector<types::global_dof_index> global_dof_indices(fe.
n_dofs_per_cell());
578 std::vector<types::global_dof_index> level_dof_indices(fe.
n_dofs_per_cell());
587 for (
unsigned int block = 0; block < n_blocks; ++block)
590 temp_copy_indices[block].resize(0);
591 temp_copy_indices[block].resize(sizes[
level][block],
597 for (; level_cell != level_end; ++level_cell)
603 level_cell->get_dof_indices(global_dof_indices);
604 level_cell->get_mg_dof_indices(level_dof_indices);
610 temp_copy_indices[block][level_dof_indices[i] -
611 mg_block_start[
level][block]] =
612 global_dof_indices[i] - block_start[block];
616 for (
unsigned int block = 0; block < n_blocks; ++block)
620 std::count_if(temp_copy_indices[block].begin(),
621 temp_copy_indices[block].end(),
625 copy_indices[block][
level].resize(n_active_dofs);
628 i < temp_copy_indices[block].size();
631 copy_indices[block][
level][counter++] =
632 std::pair<types::global_dof_index, unsigned int>(
633 temp_copy_indices[block][i], i);
642#include "mg_transfer_block.inst"
BlockType & block(const unsigned int i)
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
cell_iterator begin(const unsigned int level=0) const
unsigned int n_dofs_per_cell() const
unsigned int n_blocks() const
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
size_type n_elements() const
ElementIterator begin() const
ElementIterator end() const
bool have_boundary_indices() const
const IndexSet & get_boundary_indices(const unsigned int level) const
std::vector< unsigned int > mg_block
std::vector< std::vector< types::global_dof_index > > sizes
void build(const DoFHandler< dim, spacedim > &dof_handler)
std::vector< bool > selected
std::vector< types::global_dof_index > block_start
std::vector< std::shared_ptr< BlockSparseMatrix< double > > > prolongation_matrices
std::vector< std::shared_ptr< BlockSparsityPattern > > prolongation_sparsities
ObserverPointer< const MGConstrainedDoFs, MGTransferBlockBase > mg_constrained_dofs
std::vector< std::vector< types::global_dof_index > > mg_block_start
std::vector< std::vector< std::vector< std::pair< unsigned int, unsigned int > > > > copy_indices
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< Vector< number > > &dst, const Vector< number2 > &src) const
void build(const DoFHandler< dim, spacedim > &dof_handler, unsigned int selected)
void build(const DoFHandler< dim, spacedim > &dof_handler, const std::vector< bool > &selected)
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< BlockVector< number > > &dst, const BlockVector< number2 > &src) const
unsigned int n_levels() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
typename ActiveSelector::cell_iterator cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
std::enable_if_t< IsBlockVector< VectorType >::value, unsigned int > n_blocks(const VectorType &vector)
static const unsigned int invalid_unsigned_int
const types::global_dof_index invalid_dof_index