17 #include <deal.II/base/logstream.h> 19 #include <deal.II/dofs/dof_accessor.h> 20 #include <deal.II/dofs/dof_tools.h> 22 #include <deal.II/fe/fe.h> 24 #include <deal.II/grid/tria.h> 25 #include <deal.II/grid/tria_iterator.h> 27 #include <deal.II/lac/block_sparse_matrix.h> 28 #include <deal.II/lac/block_vector.h> 29 #include <deal.II/lac/sparse_matrix.h> 30 #include <deal.II/lac/vector.h> 32 #include <deal.II/multigrid/mg_tools.h> 33 #include <deal.II/multigrid/mg_transfer_block.h> 34 #include <deal.II/multigrid/mg_transfer_block.templates.h> 41 DEAL_II_NAMESPACE_OPEN
52 template <
int dim,
typename number,
int spacedim>
54 reinit_vector_by_blocks(
55 const ::DoFHandler<dim, spacedim> & mg_dof,
57 const std::vector<bool> & sel,
58 std::vector<std::vector<types::global_dof_index>> &ndofs)
60 std::vector<bool> selected = sel;
62 const unsigned int n_selected =
63 std::accumulate(selected.begin(), selected.end(), 0u);
65 if (ndofs.size() == 0)
67 std::vector<std::vector<types::global_dof_index>> new_dofs(
68 mg_dof.get_triangulation().n_levels(),
69 std::vector<types::global_dof_index>(selected.size()));
74 for (
unsigned int level = v.min_level(); level <= v.max_level(); ++level)
76 v[level].reinit(n_selected, 0);
78 for (
unsigned int i = 0;
79 i < selected.size() && (k < v[level].n_blocks());
84 v[level].block(k++).reinit(ndofs[level][i]);
86 v[level].collect_sizes();
98 template <
int dim,
typename number,
int spacedim>
100 reinit_vector_by_blocks(
101 const ::DoFHandler<dim, spacedim> & mg_dof,
103 const unsigned int selected_block,
104 std::vector<std::vector<types::global_dof_index>> &ndofs)
106 const unsigned int n_blocks = mg_dof.get_fe().n_blocks();
107 Assert(selected_block < n_blocks,
110 std::vector<bool> selected(n_blocks,
false);
111 selected[selected_block] =
true;
113 if (ndofs.size() == 0)
115 std::vector<std::vector<types::global_dof_index>> new_dofs(
116 mg_dof.get_triangulation().n_levels(),
117 std::vector<types::global_dof_index>(selected.size()));
122 for (
unsigned int level = v.min_level(); level <= v.max_level(); ++level)
124 v[level].reinit(ndofs[level][selected_block]);
130 template <
typename number>
131 template <
int dim,
typename number2,
int spacedim>
138 reinit_vector_by_blocks(mg_dof_handler, dst, selected_block, sizes);
147 for (IT i = copy_indices[selected_block][level].begin();
148 i != copy_indices[selected_block][level].end();
150 dst[level](i->second) = src.
block(selected_block)(i->first);
156 template <
typename number>
157 template <
int dim,
typename number2,
int spacedim>
162 const Vector<number2> & src)
const 164 reinit_vector_by_blocks(mg_dof_handler, dst, selected_block, sizes);
172 for (IT i = copy_indices[selected_block][level].begin();
173 i != copy_indices[selected_block][level].end();
175 dst[level](i->second) = src(i->first);
181 template <
typename number>
182 template <
int dim,
typename number2,
int spacedim>
189 reinit_vector_by_blocks(mg_dof_handler, dst, selected, sizes);
194 for (
unsigned int block = 0; block < selected.size(); ++block)
196 for (IT i = copy_indices[block][level].begin();
197 i != copy_indices[block][level].end();
199 dst[level].block(mg_block[block])(i->second) =
200 src.
block(block)(i->first);
206 template <
int dim,
int spacedim>
212 const unsigned int n_blocks = fe.
n_blocks();
224 for (
unsigned int i = 0; i < n_blocks; ++i)
232 sizes.resize(n_levels, std::vector<types::global_dof_index>(fe.
n_blocks()));
245 level_block_start = k;
268 for (
unsigned int block = 0; block < n_blocks; ++block)
285 prolongation_sparsities.resize(0);
287 prolongation_sparsities.reserve(n_levels - 1);
289 for (
unsigned int i = 0; i < n_levels - 1; ++i)
298 std::vector<types::global_dof_index> dof_indices_parent(dofs_per_cell);
299 std::vector<types::global_dof_index> dof_indices_child(dofs_per_cell);
309 for (
unsigned int level = 0; level < n_levels - 1; ++level)
319 prolongation_sparsities[level]->reinit(n_blocks, n_blocks);
320 for (
unsigned int i = 0; i < n_blocks; ++i)
321 for (
unsigned int j = 0; j < n_blocks; ++j)
323 prolongation_sparsities[level]->block(i, j).reinit(
324 sizes[level + 1][i],
sizes[level][j], dofs_per_cell + 1);
326 prolongation_sparsities[level]->block(i, j).reinit(
329 prolongation_sparsities[level]->collect_sizes();
333 cell != mg_dof.
end(level);
335 if (cell->has_children())
337 cell->get_mg_dof_indices(dof_indices_parent);
339 Assert(cell->n_children() ==
342 for (
unsigned int child = 0; child < cell->n_children(); ++child)
348 mg_dof.
get_fe().get_prolongation_matrix(
349 child, cell->refinement_case());
351 cell->child(child)->get_mg_dof_indices(dof_indices_child);
356 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
357 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
358 if (prolongation(i, j) != 0)
360 const unsigned int icomp =
362 const unsigned int jcomp =
364 if ((icomp == jcomp) &&
selected[icomp])
365 prolongation_sparsities[level]->add(
366 dof_indices_child[i], dof_indices_parent[j]);
370 prolongation_sparsities[level]->compress();
376 cell != mg_dof.
end(level);
378 if (cell->has_children())
380 cell->get_mg_dof_indices(dof_indices_parent);
382 Assert(cell->n_children() ==
385 for (
unsigned int child = 0; child < cell->n_children(); ++child)
391 mg_dof.
get_fe().get_prolongation_matrix(
392 child, cell->refinement_case());
394 cell->child(child)->get_mg_dof_indices(dof_indices_child);
398 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
399 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
400 if (prolongation(i, j) != 0)
402 const unsigned int icomp =
404 const unsigned int jcomp =
406 if ((icomp == jcomp) &&
selected[icomp])
408 dof_indices_child[i],
409 dof_indices_parent[j],
421 std::vector<types::global_dof_index> constrain_indices;
422 std::vector<std::vector<bool>> constraints_per_block(n_blocks);
423 for (
int level = n_levels - 2; level >= 0; --level)
434 constrain_indices.resize(0);
439 for (; dof != endd; ++dof)
440 constrain_indices[*dof] = 1;
442 unsigned int index = 0;
443 for (
unsigned int block = 0; block < n_blocks; ++block)
447 constraints_per_block[block].resize(0);
448 constraints_per_block[block].resize(n_dofs,
false);
450 constraints_per_block[block][i] =
451 (constrain_indices[index] == 1);
457 ->block(block, block)
461 for (; start_row != end_row; ++start_row)
463 if (constraints_per_block[block][start_row->column()])
464 start_row->value() = 0.;
474 template <
typename number>
475 template <
int dim,
int spacedim>
483 unsigned int n_blocks = mg_dof.
get_fe().n_blocks();
485 selected_block = select;
486 selected.resize(n_blocks,
false);
487 selected[select] =
true;
491 std::vector<types::global_dof_index> temp_copy_indices;
492 std::vector<types::global_dof_index> global_dof_indices(fe.
dofs_per_cell);
493 std::vector<types::global_dof_index> level_dof_indices(fe.
dofs_per_cell);
502 temp_copy_indices.resize(0);
503 temp_copy_indices.resize(sizes[level][selected_block],
508 for (; level_cell != level_end; ++level_cell)
514 level_cell->get_dof_indices(global_dof_indices);
515 level_cell->get_mg_dof_indices(level_dof_indices);
522 if (mg_constrained_dofs !=
nullptr)
524 if (!mg_constrained_dofs->at_refinement_edge(
525 level, level_dof_indices[i]))
526 temp_copy_indices[level_dof_indices[i] -
527 mg_block_start[level][block]] =
528 global_dof_indices[i] - block_start[block];
531 temp_copy_indices[level_dof_indices[i] -
532 mg_block_start[level][block]] =
533 global_dof_indices[i] - block_start[block];
545 std::count_if(temp_copy_indices.begin(),
546 temp_copy_indices.end(),
547 std::bind(std::not_equal_to<types::global_dof_index>(),
548 std::placeholders::_1,
550 copy_indices[selected_block][level].resize(n_active_dofs);
554 copy_indices[selected_block][level][counter++] =
555 std::pair<types::global_dof_index, unsigned int>(
556 temp_copy_indices[i], i);
563 template <
typename number>
564 template <
int dim,
int spacedim>
568 const std::vector<bool> & sel)
571 unsigned int n_blocks = mg_dof.
get_fe().n_blocks();
575 Assert(sel.size() == n_blocks,
579 if (selected.size() == 0)
580 selected = std::vector<bool>(n_blocks,
true);
584 std::vector<std::vector<types::global_dof_index>> temp_copy_indices(n_blocks);
585 std::vector<types::global_dof_index> global_dof_indices(fe.
dofs_per_cell);
586 std::vector<types::global_dof_index> level_dof_indices(fe.
dofs_per_cell);
594 for (
unsigned int block = 0; block < n_blocks; ++block)
597 temp_copy_indices[block].resize(0);
598 temp_copy_indices[block].resize(sizes[level][block],
604 for (; level_cell != level_end; ++level_cell)
610 level_cell->get_dof_indices(global_dof_indices);
611 level_cell->get_mg_dof_indices(level_dof_indices);
617 temp_copy_indices[block][level_dof_indices[i] -
618 mg_block_start[level][block]] =
619 global_dof_indices[i] - block_start[block];
623 for (
unsigned int block = 0; block < n_blocks; ++block)
627 temp_copy_indices[block].begin(),
628 temp_copy_indices[block].end(),
629 std::bind(std::not_equal_to<types::global_dof_index>(),
630 std::placeholders::_1,
632 copy_indices[block][level].resize(n_active_dofs);
635 i < temp_copy_indices[block].size();
638 copy_indices[block][level][counter++] =
639 std::pair<types::global_dof_index, unsigned int>(
640 temp_copy_indices[block][i], i);
649 #include "mg_transfer_block.inst" 652 DEAL_II_NAMESPACE_CLOSE
void build_matrices(const DoFHandler< dim, spacedim > &dof, const DoFHandler< dim, spacedim > &mg_dof)
const Triangulation< dim, spacedim > & get_triangulation() const
static const unsigned int invalid_unsigned_int
ElementIterator end() const
cell_iterator begin(const unsigned int level=0) const
cell_iterator end() const
std::vector< types::global_dof_index > block_start
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
SmartPointer< const MGConstrainedDoFs, MGTransferBlockBase > mg_constrained_dofs
unsigned int n_blocks() const
active_cell_iterator begin_active(const unsigned int level=0) const
const IndexSet & get_boundary_indices(const unsigned int level) const
std::vector< std::vector< types::global_dof_index > > sizes
std::vector< std::shared_ptr< BlockSparseMatrix< double > > > prolongation_matrices
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
active_cell_iterator end_active(const unsigned int level) const
const unsigned int dofs_per_cell
std::vector< std::vector< std::vector< std::pair< unsigned int, unsigned int > > > > copy_indices
std::vector< unsigned int > mg_block
void swap(Vector< Number > &u, Vector< Number > &v)
unsigned int global_dof_index
void build_matrices(const DoFHandler< dim, spacedim > &dof, const DoFHandler< dim, spacedim > &mg_dof, unsigned int selected)
std::vector< std::vector< types::global_dof_index > > mg_block_start
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
void copy_to_mg(const DoFHandler< dim, spacedim > &mg_dof, MGLevelObject< Vector< number >> &dst, const Vector< number2 > &src) const
void build_matrices(const DoFHandler< dim, spacedim > &dof, const DoFHandler< dim, spacedim > &mg_dof, const std::vector< bool > &selected)
static ::ExceptionBase & ExcNotImplemented()
typename ActiveSelector::cell_iterator cell_iterator
void copy_to_mg(const DoFHandler< dim, spacedim > &mg_dof, MGLevelObject< BlockVector< number >> &dst, const BlockVector< number2 > &src) const
typename ActiveSelector::active_cell_iterator active_cell_iterator
bool have_boundary_indices() const
const types::global_dof_index invalid_dof_index
ElementIterator begin() const
BlockType & block(const unsigned int i)
size_type n_elements() const
std::vector< bool > selected
static ::ExceptionBase & ExcInternalError()