17 #include <deal.II/base/logstream.h> 19 #include <deal.II/lac/vector.h> 20 #include <deal.II/lac/block_vector.h> 21 #include <deal.II/lac/sparse_matrix.h> 22 #include <deal.II/lac/block_sparse_matrix.h> 23 #include <deal.II/grid/tria.h> 24 #include <deal.II/grid/tria_iterator.h> 25 #include <deal.II/dofs/dof_tools.h> 26 #include <deal.II/fe/fe.h> 27 #include <deal.II/dofs/dof_accessor.h> 28 #include <deal.II/multigrid/mg_transfer_block.h> 29 #include <deal.II/multigrid/mg_transfer_block.templates.h> 30 #include <deal.II/multigrid/mg_tools.h> 37 DEAL_II_NAMESPACE_OPEN
48 template <
int dim,
typename number,
int spacedim>
50 reinit_vector_by_blocks (
51 const ::DoFHandler<dim,spacedim> &mg_dof,
53 const std::vector<bool> &sel,
54 std::vector<std::vector<types::global_dof_index> > &ndofs)
56 std::vector<bool> selected=sel;
58 const unsigned int n_selected
59 = std::accumulate(selected.begin(),
63 if (ndofs.size() == 0)
65 std::vector<std::vector<types::global_dof_index> >
66 new_dofs(mg_dof.get_triangulation().n_levels(),
67 std::vector<types::global_dof_index>(selected.size()));
72 for (
unsigned int level=v.min_level();
73 level<=v.max_level(); ++level)
75 v[level].reinit(n_selected, 0);
77 for (
unsigned int i=0; i<selected.size() && (k<v[level].n_blocks()); ++i)
81 v[level].block(k++).reinit(ndofs[level][i]);
83 v[level].collect_sizes();
95 template <
int dim,
typename number,
int spacedim>
97 reinit_vector_by_blocks (
98 const ::DoFHandler<dim,spacedim> &mg_dof,
100 const unsigned int selected_block,
101 std::vector<std::vector<types::global_dof_index> > &ndofs)
103 const unsigned int n_blocks = mg_dof.get_fe().n_blocks();
106 std::vector<bool> selected(n_blocks,
false);
107 selected[selected_block] =
true;
109 if (ndofs.size() == 0)
111 std::vector<std::vector<types::global_dof_index> >
112 new_dofs(mg_dof.get_triangulation().n_levels(),
113 std::vector<types::global_dof_index>(selected.size()));
118 for (
unsigned int level=v.min_level();
119 level<=v.max_level(); ++level)
121 v[level].reinit(ndofs[level][selected_block]);
127 template <
typename number>
128 template <
int dim,
typename number2,
int spacedim>
135 reinit_vector_by_blocks(mg_dof_handler, dst, selected_block, sizes);
140 for (
unsigned int level=mg_dof_handler.
get_triangulation().n_levels(); level != 0;)
143 for (IT i= copy_indices[selected_block][level].begin();
144 i != copy_indices[selected_block][level].end(); ++i)
145 dst[level](i->second) = src.
block(selected_block)(i->first);
151 template <
typename number>
152 template <
int dim,
typename number2,
int spacedim>
157 const Vector<number2> &src)
const 159 reinit_vector_by_blocks(mg_dof_handler, dst, selected_block, sizes);
163 for (
unsigned int level=mg_dof_handler.
get_triangulation().n_levels(); level != 0;)
166 for (IT i= copy_indices[selected_block][level].begin();
167 i != copy_indices[selected_block][level].end(); ++i)
168 dst[level](i->second) = src(i->first);
174 template <
typename number>
175 template <
int dim,
typename number2,
int spacedim>
182 reinit_vector_by_blocks(mg_dof_handler, dst, selected, sizes);
183 for (
unsigned int level=mg_dof_handler.
get_triangulation().n_levels(); level != 0;)
186 for (
unsigned int block=0; block<selected.size(); ++block)
188 for (IT i= copy_indices[block][level].begin();
189 i != copy_indices[block][level].end(); ++i)
190 dst[level].block(mg_block[block])(i->second) = src.
block(block)(i->first);
196 template <
int dim,
int spacedim>
202 const unsigned int n_blocks = fe.
n_blocks();
214 for (
unsigned int i=0; i<n_blocks; ++i)
222 sizes.resize(n_levels, std::vector<types::global_dof_index>(fe.
n_blocks()));
258 for (
unsigned int block=0; block<n_blocks; ++block)
275 prolongation_sparsities.resize (0);
277 prolongation_sparsities.reserve (n_levels - 1);
279 for (
unsigned int i=0; i<n_levels-1; ++i)
288 std::vector<types::global_dof_index> dof_indices_parent (dofs_per_cell);
289 std::vector<types::global_dof_index> dof_indices_child (dofs_per_cell);
299 for (
unsigned int level=0; level<n_levels-1; ++level)
309 prolongation_sparsities[level]->reinit (n_blocks, n_blocks);
310 for (
unsigned int i=0; i<n_blocks; ++i)
311 for (
unsigned int j=0; j<n_blocks; ++j)
313 prolongation_sparsities[level]->block(i,j)
314 .reinit(
sizes[level+1][i],
318 prolongation_sparsities[level]->block(i,j)
319 .reinit(
sizes[level+1][i],
323 prolongation_sparsities[level]->collect_sizes();
326 cell != mg_dof.
end(level); ++cell)
327 if (cell->has_children())
329 cell->get_mg_dof_indices (dof_indices_parent);
333 for (
unsigned int child=0; child<cell->n_children(); ++child)
339 = mg_dof.
get_fe().get_prolongation_matrix (child, cell->refinement_case());
341 cell->child(child)->get_mg_dof_indices (dof_indices_child);
346 for (
unsigned int i=0; i<dofs_per_cell; ++i)
347 for (
unsigned int j=0; j<dofs_per_cell; ++j)
348 if (prolongation(i,j) != 0)
350 const unsigned int icomp
352 const unsigned int jcomp
354 if ((icomp==jcomp) &&
selected[icomp])
355 prolongation_sparsities[level]->add(dof_indices_child[i],
356 dof_indices_parent[j]);
360 prolongation_sparsities[level]->compress ();
365 cell != mg_dof.
end(level); ++cell)
366 if (cell->has_children())
368 cell->get_mg_dof_indices (dof_indices_parent);
372 for (
unsigned int child=0; child<cell->n_children(); ++child)
378 = mg_dof.
get_fe().get_prolongation_matrix (child, cell->refinement_case());
380 cell->child(child)->get_mg_dof_indices (dof_indices_child);
384 for (
unsigned int i=0; i<dofs_per_cell; ++i)
385 for (
unsigned int j=0; j<dofs_per_cell; ++j)
386 if (prolongation(i,j) != 0)
390 if ((icomp==jcomp) &&
selected[icomp])
392 dof_indices_parent[j],
403 std::vector<types::global_dof_index> constrain_indices;
404 std::vector<std::vector<bool> > constraints_per_block (n_blocks);
405 for (
int level=n_levels-2; level>=0; --level)
415 constrain_indices.resize (0);
420 for (; dof != endd; ++dof)
421 constrain_indices[*dof] = 1;
423 unsigned int index = 0;
424 for (
unsigned int block=0; block<n_blocks; ++block)
427 constraints_per_block[block].resize(0);
428 constraints_per_block[block].resize(n_dofs,
false);
430 constraints_per_block[block][i] = (constrain_indices[index] == 1);
437 for (; start_row != end_row; ++start_row)
439 if (constraints_per_block[block][start_row->column()])
440 start_row->value() = 0.;
450 template <
typename number>
451 template <
int dim,
int spacedim>
458 unsigned int n_blocks = mg_dof.
get_fe().n_blocks();
460 selected_block = select;
461 selected.resize(n_blocks,
false);
462 selected[select] =
true;
466 std::vector<types::global_dof_index> temp_copy_indices;
467 std::vector<types::global_dof_index> global_dof_indices (fe.
dofs_per_cell);
468 std::vector<types::global_dof_index> level_dof_indices (fe.
dofs_per_cell);
477 temp_copy_indices.resize (0);
478 temp_copy_indices.resize (sizes[level][selected_block],
483 for (; level_cell!=level_end; ++level_cell)
489 level_cell->get_dof_indices(global_dof_indices);
490 level_cell->get_mg_dof_indices (level_dof_indices);
497 if (mg_constrained_dofs !=
nullptr)
499 if (!mg_constrained_dofs->at_refinement_edge(level,level_dof_indices[i]))
500 temp_copy_indices[level_dof_indices[i] - mg_block_start[level][block]]
501 = global_dof_indices[i] - block_start[block];
504 temp_copy_indices[level_dof_indices[i] - mg_block_start[level][block]]
505 = global_dof_indices[i] - block_start[block];
517 std::count_if (temp_copy_indices.begin(), temp_copy_indices.end(),
518 std::bind (std::not_equal_to<types::global_dof_index>(),
519 std::placeholders::_1,
521 copy_indices[selected_block][level].resize (n_active_dofs);
525 copy_indices[selected_block][level][counter++] =
526 std::pair<types::global_dof_index, unsigned int> (temp_copy_indices[i], i);
534 template <
typename number>
535 template <
int dim,
int spacedim>
539 const std::vector<bool> &sel)
542 unsigned int n_blocks = mg_dof.
get_fe().n_blocks();
546 Assert(sel.size() == n_blocks,
550 if (selected.size() == 0)
551 selected = std::vector<bool> (n_blocks,
true);
555 std::vector<std::vector<types::global_dof_index> > temp_copy_indices (n_blocks);
556 std::vector<types::global_dof_index> global_dof_indices (fe.
dofs_per_cell);
557 std::vector<types::global_dof_index> level_dof_indices (fe.
dofs_per_cell);
565 for (
unsigned int block=0; block<n_blocks; ++block)
568 temp_copy_indices[block].resize (0);
569 temp_copy_indices[block].resize (sizes[level][block],
575 for (; level_cell!=level_end; ++level_cell)
581 level_cell->get_dof_indices(global_dof_indices);
582 level_cell->get_mg_dof_indices (level_dof_indices);
588 temp_copy_indices[block][level_dof_indices[i] - mg_block_start[level][block]]
589 = global_dof_indices[i] - block_start[block];
593 for (
unsigned int block=0; block<n_blocks; ++block)
597 std::count_if (temp_copy_indices[block].begin(),
598 temp_copy_indices[block].end(),
599 std::bind (std::not_equal_to<types::global_dof_index>(),
600 std::placeholders::_1,
602 copy_indices[block][level].resize (n_active_dofs);
606 copy_indices[block][level][counter++] =
607 std::pair<types::global_dof_index, unsigned int>
608 (temp_copy_indices[block][i], i);
617 #include "mg_transfer_block.inst" 620 DEAL_II_NAMESPACE_CLOSE
void build_matrices(const DoFHandler< dim, spacedim > &dof, const DoFHandler< dim, spacedim > &mg_dof)
static const unsigned int invalid_unsigned_int
ElementIterator end() const
cell_iterator begin(const unsigned int level=0) const
cell_iterator end() const
ActiveSelector::active_cell_iterator active_cell_iterator
std::vector< types::global_dof_index > block_start
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
SmartPointer< const MGConstrainedDoFs, MGTransferBlockBase > mg_constrained_dofs
void copy_to_mg(const DoFHandler< dim, spacedim > &mg_dof, MGLevelObject< Vector< number > > &dst, const Vector< number2 > &src) const
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
ActiveSelector::cell_iterator cell_iterator
void copy_to_mg(const DoFHandler< dim, spacedim > &mg_dof, MGLevelObject< BlockVector< number > > &dst, const BlockVector< number2 > &src) const
std::vector< std::vector< types::global_dof_index > > sizes
unsigned int global_dof_index
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 FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) 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)
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 build_matrices(const DoFHandler< dim, spacedim > &dof, const DoFHandler< dim, spacedim > &mg_dof, const std::vector< bool > &selected)
static ::ExceptionBase & ExcNotImplemented()
const Triangulation< dim, spacedim > & get_triangulation() const
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()