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()));
68 std::swap(ndofs, new_dofs);
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();
104 Assert(selected_block < n_blocks, ExcIndexRange(selected_block, 0, 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()));
114 std::swap(ndofs, new_dofs);
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);
141 for (
unsigned int level=mg_dof_handler.
get_triangulation().n_levels(); level != 0;)
144 for (IT i= copy_indices[selected_block][level].begin();
145 i != copy_indices[selected_block][level].end(); ++i)
146 dst[level](i->second) = src.
block(selected_block)(i->first);
148 restrict_and_add (level+1, dst[level], dst[level+1]);
155 template <
typename number>
156 template <
int dim,
typename number2,
int spacedim>
163 reinit_vector_by_blocks(mg_dof_handler, dst, selected_block, sizes);
168 for (
unsigned int level=mg_dof_handler.
get_triangulation().n_levels(); level != 0;)
171 for (IT i= copy_indices[selected_block][level].begin();
172 i != copy_indices[selected_block][level].end(); ++i)
173 dst[level](i->second) = src(i->first);
175 restrict_and_add (level+1, dst[level], dst[level+1]);
182 template <
typename number>
183 template <
int dim,
typename number2,
int spacedim>
190 reinit_vector_by_blocks(mg_dof_handler, dst, selected, sizes);
192 for (
unsigned int level=mg_dof_handler.
get_triangulation().n_levels(); level != 0;)
195 for (
unsigned int block=0; block<selected.size(); ++block)
197 for (IT i= copy_indices[block][level].begin();
198 i != copy_indices[block][level].end(); ++i)
199 dst[level].block(mg_block[block])(i->second) = src.
block(block)(i->first);
201 restrict_and_add (level+1, dst[level], dst[level+1]);
208 template <
int dim,
int spacedim>
214 const unsigned int n_blocks = fe.
n_blocks();
218 Assert (selected.size() == n_blocks,
219 ExcDimensionMismatch(selected.size(), n_blocks));
226 for (
unsigned int i=0; i<n_blocks; ++i)
234 sizes.resize(n_levels, std::vector<types::global_dof_index>(fe.
n_blocks()));
270 for (
unsigned int block=0; block<n_blocks; ++block)
287 prolongation_sparsities.resize (0);
289 for (
unsigned int i=0; i<n_levels-1; ++i)
291 prolongation_sparsities
300 std::vector<types::global_dof_index> dof_indices_parent (dofs_per_cell);
301 std::vector<types::global_dof_index> dof_indices_child (dofs_per_cell);
311 for (
unsigned int level=0; level<n_levels-1; ++level)
321 prolongation_sparsities[level]->reinit (n_blocks, n_blocks);
322 for (
unsigned int i=0; i<n_blocks; ++i)
323 for (
unsigned int j=0; j<n_blocks; ++j)
325 prolongation_sparsities[level]->block(i,j)
326 .reinit(
sizes[level+1][i],
330 prolongation_sparsities[level]->block(i,j)
331 .reinit(
sizes[level+1][i],
335 prolongation_sparsities[level]->collect_sizes();
338 cell != mg_dof.
end(level); ++cell)
339 if (cell->has_children())
341 cell->get_mg_dof_indices (dof_indices_parent);
344 ExcNotImplemented());
345 for (
unsigned int child=0; child<cell->n_children(); ++child)
351 = mg_dof.
get_fe().get_prolongation_matrix (child, cell->refinement_case());
353 cell->child(child)->get_mg_dof_indices (dof_indices_child);
358 for (
unsigned int i=0; i<dofs_per_cell; ++i)
359 for (
unsigned int j=0; j<dofs_per_cell; ++j)
360 if (prolongation(i,j) != 0)
362 const unsigned int icomp
364 const unsigned int jcomp
366 if ((icomp==jcomp) && selected[icomp])
367 prolongation_sparsities[level]->add(dof_indices_child[i],
368 dof_indices_parent[j]);
372 prolongation_sparsities[level]->compress ();
377 cell != mg_dof.
end(level); ++cell)
378 if (cell->has_children())
380 cell->get_mg_dof_indices (dof_indices_parent);
383 ExcNotImplemented());
384 for (
unsigned int child=0; child<cell->n_children(); ++child)
390 = mg_dof.
get_fe().get_prolongation_matrix (child, cell->refinement_case());
392 cell->child(child)->get_mg_dof_indices (dof_indices_child);
396 for (
unsigned int i=0; i<dofs_per_cell; ++i)
397 for (
unsigned int j=0; j<dofs_per_cell; ++j)
398 if (prolongation(i,j) != 0)
402 if ((icomp==jcomp) && selected[icomp])
404 dof_indices_parent[j],
415 std::vector<types::global_dof_index> constrain_indices;
416 std::vector<std::vector<bool> > constraints_per_block (n_blocks);
417 for (
int level=n_levels-2; level>=0; --level)
427 constrain_indices.resize (0);
432 for (; dof != endd; ++dof)
433 constrain_indices[*dof] = 1;
435 unsigned int index = 0;
436 for (
unsigned int block=0; block<n_blocks; ++block)
439 constraints_per_block[block].resize(0);
440 constraints_per_block[block].resize(n_dofs, 0);
442 constraints_per_block[block][i] = (constrain_indices[index] == 1);
449 for (; start_row != end_row; ++start_row)
451 if (constraints_per_block[block][start_row->column()])
452 start_row->value() = 0.;
462 template <
typename number>
463 template <
int dim,
int spacedim>
470 unsigned int n_blocks = mg_dof.
get_fe().n_blocks();
472 selected_block = select;
473 selected.resize(n_blocks,
false);
474 selected[select] =
true;
478 std::vector<types::global_dof_index> temp_copy_indices;
479 std::vector<types::global_dof_index> global_dof_indices (fe.
dofs_per_cell);
480 std::vector<types::global_dof_index> level_dof_indices (fe.
dofs_per_cell);
489 temp_copy_indices.resize (0);
490 temp_copy_indices.resize (
sizes[level][selected_block],
495 for (; level_cell!=level_end; ++level_cell)
501 level_cell->get_dof_indices(global_dof_indices);
502 level_cell->get_mg_dof_indices (level_dof_indices);
512 temp_copy_indices[level_dof_indices[i] -
mg_block_start[level][block]]
516 temp_copy_indices[level_dof_indices[i] -
mg_block_start[level][block]]
529 std::count_if (temp_copy_indices.begin(), temp_copy_indices.end(),
530 std::bind2nd(std::not_equal_to<types::global_dof_index>(),
532 copy_indices[selected_block][level].resize (n_active_dofs);
537 std::pair<types::global_dof_index, unsigned int> (temp_copy_indices[i], i);
538 Assert (counter == n_active_dofs, ExcInternalError());
545 template <
typename number>
546 template <
int dim,
int spacedim>
550 const std::vector<bool> &sel)
553 unsigned int n_blocks = mg_dof.
get_fe().n_blocks();
557 Assert(sel.size() == n_blocks,
558 ExcDimensionMismatch(sel.size(), n_blocks));
561 if (selected.size() == 0)
562 selected = std::vector<bool> (n_blocks,
true);
566 std::vector<std::vector<types::global_dof_index> > temp_copy_indices (n_blocks);
567 std::vector<types::global_dof_index> global_dof_indices (fe.
dofs_per_cell);
568 std::vector<types::global_dof_index> level_dof_indices (fe.
dofs_per_cell);
576 for (
unsigned int block=0; block<n_blocks; ++block)
579 temp_copy_indices[block].resize (0);
580 temp_copy_indices[block].resize (
sizes[level][block],
586 for (; level_cell!=level_end; ++level_cell)
592 level_cell->get_dof_indices(global_dof_indices);
593 level_cell->get_mg_dof_indices (level_dof_indices);
599 temp_copy_indices[block][level_dof_indices[i] -
mg_block_start[level][block]]
604 for (
unsigned int block=0; block<n_blocks; ++block)
608 std::count_if (temp_copy_indices[block].begin(),
609 temp_copy_indices[block].end(),
610 std::bind2nd(std::not_equal_to<types::global_dof_index>(),
617 std::pair<types::global_dof_index, unsigned int>
618 (temp_copy_indices[block][i], i);
619 Assert (counter == n_active_dofs, ExcInternalError());
627 #include "mg_transfer_block.inst" 630 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
const FiniteElement< dim, spacedim > & get_fe() const
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_cxx11::shared_ptr< BlockSparseMatrix< double > > > prolongation_matrices
#define Assert(cond, exc)
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
bool at_refinement_edge(const unsigned int level, const types::global_dof_index index) const
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)
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