33#include <deal.II/multigrid/mg_transfer_block.templates.h>
51 template <
int dim,
typename number,
int spacedim>
56 const std::vector<bool> &
sel,
57 std::vector<std::vector<types::global_dof_index>> &
ndofs)
59 std::vector<bool> selected =
sel;
62 std::accumulate(selected.begin(), selected.end(), 0
u);
66 std::vector<std::vector<types::global_dof_index>>
new_dofs(
67 dof_handler.get_triangulation().n_levels(),
68 std::vector<types::global_dof_index>(selected.size()));
73 for (
unsigned int level = v.min_level();
level <= v.max_level(); ++
level)
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>
102 const unsigned int selected_block,
103 std::vector<std::vector<types::global_dof_index>> &
ndofs)
105 const unsigned int n_blocks = dof_handler.get_fe().n_blocks();
108 std::vector<bool> selected(n_blocks,
false);
109 selected[selected_block] =
true;
113 std::vector<std::vector<types::global_dof_index>>
new_dofs(
114 dof_handler.get_triangulation().n_levels(),
115 std::vector<types::global_dof_index>(selected.size()));
120 for (
unsigned int level = v.min_level();
level <= v.max_level(); ++
level)
128template <
typename number>
129template <
int dim,
typename number2,
int spacedim>
141 for (
unsigned int level = dof_handler.get_triangulation().n_levels();
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>
166 for (
unsigned int level = dof_handler.get_triangulation().n_levels();
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>
188 for (
unsigned int level = dof_handler.get_triangulation().n_levels();
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();
210 const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
211 const unsigned int n_levels = dof_handler.get_triangulation().n_levels();
221 for (
unsigned int i = 0; i < n_blocks; ++i)
229 sizes.resize(n_levels, std::vector<types::global_dof_index>(fe.n_blocks()));
264 for (
unsigned int block = 0; block < n_blocks; ++block)
285 for (
unsigned int i = 0; i < n_levels - 1; ++i)
316 for (
unsigned int i = 0; i < n_blocks; ++i)
317 for (
unsigned int j = 0;
j < n_blocks; ++
j)
328 dof_handler.begin(
level);
329 cell != dof_handler.end(
level);
331 if (cell->has_children())
335 Assert(cell->n_children() ==
338 for (
unsigned int child = 0; child < cell->n_children(); ++child)
344 dof_handler.get_fe().get_prolongation_matrix(
345 child, cell->refinement_case());
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 =
357 fe.system_to_block_index(i).first;
358 const unsigned int jcomp =
359 fe.system_to_block_index(
j).first;
371 dof_handler.begin(
level);
372 cell != dof_handler.end(
level);
374 if (cell->has_children())
378 Assert(cell->n_children() ==
381 for (
unsigned int child = 0; child < cell->n_children(); ++child)
387 dof_handler.get_fe().get_prolongation_matrix(
388 child, cell->refinement_case());
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 =
399 fe.system_to_block_index(i).first;
400 const unsigned int jcomp =
401 fe.system_to_block_index(
j).first;
435 for (; dof !=
endd; ++dof)
438 unsigned int index = 0;
439 for (
unsigned int block = 0; block < n_blocks; ++block)
453 ->block(block, block)
468template <
typename number>
469template <
int dim,
int spacedim>
476 unsigned int n_blocks = dof_handler.get_fe().n_blocks();
478 selected_block = select;
479 selected.resize(n_blocks,
false);
480 selected[select] =
true;
486 std::vector<types::global_dof_index> level_dof_indices(fe.n_dofs_per_cell());
488 for (
int level = dof_handler.get_triangulation().n_levels() - 1;
level >= 0;
492 dof_handler.begin_active(
level);
494 dof_handler.end_active(
level);
502 for (; level_cell !=
level_end; ++level_cell)
509 level_cell->get_mg_dof_indices(level_dof_indices);
511 for (
unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
513 const unsigned int block = fe.system_to_block_index(i).first;
516 if (mg_constrained_dofs !=
nullptr)
518 if (!mg_constrained_dofs->at_refinement_edge(
519 level, level_dof_indices[i]))
521 mg_block_start[
level][block]] =
526 mg_block_start[
level][block]] =
542 return index != numbers::invalid_dof_index;
548 copy_indices[selected_block][
level][counter++] =
549 std::pair<types::global_dof_index, unsigned int>(
556template <
typename number>
557template <
int dim,
int spacedim>
560 const std::vector<bool> &
sel)
563 unsigned int n_blocks = dof_handler.get_fe().n_blocks();
571 if (selected.empty())
572 selected = std::vector<bool>(n_blocks,
true);
578 std::vector<types::global_dof_index> level_dof_indices(fe.n_dofs_per_cell());
579 for (
int level = dof_handler.get_triangulation().n_levels() - 1;
level >= 0;
583 dof_handler.begin_active(
level);
585 dof_handler.end_active(
level);
587 for (
unsigned int block = 0; block < n_blocks; ++block)
597 for (; level_cell !=
level_end; ++level_cell)
604 level_cell->get_mg_dof_indices(level_dof_indices);
606 for (
unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
608 const unsigned int block = fe.system_to_block_index(i).first;
611 mg_block_start[
level][block]] =
616 for (
unsigned int block = 0; block < n_blocks; ++block)
631 copy_indices[block][
level][counter++] =
632 std::pair<types::global_dof_index, unsigned int>(
642#include "mg_transfer_block.inst"
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
#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