17 #include <deal.II/base/logstream.h> 18 #include <deal.II/base/function.h> 20 #include <deal.II/lac/vector.h> 21 #include <deal.II/lac/block_vector.h> 22 #include <deal.II/lac/block_indices.h> 23 #include <deal.II/lac/sparse_matrix.h> 24 #include <deal.II/lac/block_sparse_matrix.h> 25 #include <deal.II/grid/tria.h> 26 #include <deal.II/grid/tria_iterator.h> 27 #include <deal.II/dofs/dof_tools.h> 28 #include <deal.II/fe/fe.h> 29 #include <deal.II/dofs/dof_accessor.h> 30 #include <deal.II/multigrid/mg_transfer_component.h> 31 #include <deal.II/multigrid/mg_transfer_component.templates.h> 32 #include <deal.II/multigrid/mg_tools.h> 38 DEAL_II_NAMESPACE_OPEN
77 template <
int dim,
typename number,
int spacedim>
79 reinit_vector_by_components (
80 const ::DoFHandler<dim,spacedim> &mg_dof,
82 const std::vector<bool> &sel,
83 const std::vector<unsigned int> &target_comp,
84 std::vector<std::vector<types::global_dof_index> > &ndofs)
86 std::vector<bool> selected=sel;
87 std::vector<unsigned int> target_component=target_comp;
88 const unsigned int ncomp = mg_dof.get_fe(0).n_components();
99 if (target_component.size() == 0)
101 target_component.resize(ncomp);
102 for (
unsigned int i=0; i<ncomp; ++i)
103 target_component[i] = i;
108 if (selected.size() == 0)
110 selected.resize(target_component.size());
111 std::fill_n (selected.begin(), ncomp,
false);
112 for (
unsigned int i=0; i<target_component.size(); ++i)
113 selected[target_component[i]] =
true;
116 Assert (selected.size() == target_component.size(),
120 const unsigned int n_selected
121 = std::accumulate(selected.begin(),
125 if (ndofs.size() == 0)
127 std::vector<std::vector<types::global_dof_index> >
128 new_dofs(mg_dof.get_triangulation().n_levels(),
129 std::vector<types::global_dof_index>(target_component.size()));
134 for (
unsigned int level=v.min_level();
135 level<=v.max_level(); ++level)
137 v[level].reinit(n_selected, 0);
139 for (
unsigned int i=0; i<selected.size() && (k<v[level].n_blocks()); ++i)
143 v[level].block(k++).reinit(ndofs[level][i]);
145 v[level].collect_sizes();
177 template <
int dim,
typename number,
int spacedim>
179 reinit_vector_by_components (
180 const ::DoFHandler<dim,spacedim> &mg_dof,
183 const std::vector<unsigned int> &target_component,
184 std::vector<std::vector<types::global_dof_index> > &ndofs)
187 ExcMessage (
"The component mask does not have the correct size."));
189 unsigned int selected_block = 0;
190 for (
unsigned int i=0; i<target_component.size(); ++i)
191 if (component_mask[i])
192 selected_block = target_component[i];
194 if (ndofs.size() == 0)
196 std::vector<std::vector<types::global_dof_index> >
197 new_dofs(mg_dof.get_triangulation().n_levels(),
198 std::vector<types::global_dof_index>(target_component.size()));
204 for (
unsigned int level=v.min_level();
205 level<=v.max_level(); ++level)
207 v[level].reinit(ndofs[level][selected_block]);
213 template <
typename number>
214 template <
int dim,
class InVector,
int spacedim>
219 const InVector &src)
const 224 ExcMatricesNotBuilt());
226 reinit_vector_by_components(mg_dof_handler, dst,
228 mg_target_component, sizes);
241 for (
unsigned int level=mg_dof_handler.
get_triangulation().n_levels(); level!=0;)
245 typedef std::vector<std::pair<types::global_dof_index, unsigned int> >::const_iterator IT;
246 for (IT i=copy_to_and_from_indices[level].begin();
247 i != copy_to_and_from_indices[level].end(); ++i)
248 dst[level](i->second) = src(i->first);
253 template <
int dim,
int spacedim>
272 mg_dof.
get_fe(0).n_components()));
293 mg_dof.
get_fe(0).n_components()));
309 const unsigned int n_components =
315 ExcMessage (
"Component mask has wrong size."));
318 sizes.resize(n_levels);
319 for (
unsigned int l=0; l<n_levels; ++l)
320 sizes[l].resize(n_components);
374 prolongation_sparsities.resize (0);
376 prolongation_sparsities.reserve (n_levels - 1);
378 for (
unsigned int i=0; i<n_levels-1; ++i)
387 std::vector<types::global_dof_index> dof_indices_parent (dofs_per_cell);
388 std::vector<types::global_dof_index> dof_indices_child (dofs_per_cell);
397 for (
unsigned int level=0; level<n_levels-1; ++level)
407 prolongation_sparsities[level]->reinit (n_components, n_components);
408 for (
unsigned int i=0; i<n_components; ++i)
409 for (
unsigned int j=0; j<n_components; ++j)
411 prolongation_sparsities[level]->block(i,j)
412 .reinit(
sizes[level+1][i],
416 prolongation_sparsities[level]->block(i,j)
417 .reinit(
sizes[level+1][i],
421 prolongation_sparsities[level]->collect_sizes();
424 cell != mg_dof.
end(level); ++cell)
425 if (cell->has_children())
427 cell->get_mg_dof_indices (dof_indices_parent);
429 for (
unsigned int child=0; child<cell->n_children(); ++child)
435 = mg_dof.
get_fe().get_prolongation_matrix (child, cell->refinement_case());
437 cell->child(child)->get_mg_dof_indices (dof_indices_child);
442 for (
unsigned int i=0; i<dofs_per_cell; ++i)
443 for (
unsigned int j=0; j<dofs_per_cell; ++j)
444 if (prolongation(i,j) != 0)
446 const unsigned int icomp
448 const unsigned int jcomp
451 prolongation_sparsities[level]->add(dof_indices_child[i],
452 dof_indices_parent[j]);
456 prolongation_sparsities[level]->compress ();
461 cell != mg_dof.
end(level); ++cell)
462 if (cell->has_children())
464 cell->get_mg_dof_indices (dof_indices_parent);
466 for (
unsigned int child=0; child<cell->n_children(); ++child)
472 = mg_dof.
get_fe().get_prolongation_matrix (child, cell->refinement_case());
474 cell->child(child)->get_mg_dof_indices (dof_indices_child);
478 for (
unsigned int i=0; i<dofs_per_cell; ++i)
479 for (
unsigned int j=0; j<dofs_per_cell; ++j)
480 if (prolongation(i,j) != 0)
486 dof_indices_parent[j],
499 std::vector<std::vector<types::global_dof_index> >
501 std::vector<types::global_dof_index>(n_components));
504 for (
unsigned int level=0; level<n_levels-1; ++level)
509 for (
unsigned int iblock=0; iblock<n_components; ++iblock)
510 for (
unsigned int jblock=0; jblock<n_components; ++jblock)
518 for (; anfang != ende; ++anfang)
523 const BlockIndices block_indices_coarse (dofs_per_component[level]);
526 std::set<types::global_dof_index>::const_iterator found_dof =
529 const bool is_boundary_index =
532 if (is_boundary_index)
535 .set(i,column_number,0);
545 template <
typename number>
546 template <
int dim,
int spacedim>
551 unsigned int mg_select,
552 const std::vector<unsigned int> &t_component,
553 const std::vector<unsigned int> &mg_t_component,
554 const std::vector<std::set<types::global_dof_index> > &bdry_indices)
557 unsigned int ncomp = mg_dof.
get_fe(0).n_components();
559 target_component = t_component;
560 mg_target_component = mg_t_component;
561 boundary_indices = bdry_indices;
563 selected_component = select;
564 mg_selected_component = mg_select;
567 std::vector<bool> tmp(ncomp,
false);
568 for (
unsigned int c=0; c<ncomp; ++c)
569 if (t_component[c] == selected_component)
575 std::vector<bool> tmp(ncomp,
false);
576 for (
unsigned int c=0; c<ncomp; ++c)
577 if (mg_t_component[c] == mg_selected_component)
586 for (
unsigned int i=0; i<target_component.size(); ++i)
588 if (target_component[i] == select)
590 selected_component = i;
595 for (
unsigned int i=0; i<mg_target_component.size(); ++i)
597 if (mg_target_component[i] == mg_select)
599 mg_selected_component = i;
609 interface_dofs[l].clear();
610 interface_dofs[l].set_size(mg_dof.
n_dofs(l));
616 std::vector<types::global_dof_index> temp_copy_indices;
617 std::vector<types::global_dof_index> global_dof_indices (fe.
dofs_per_cell);
618 std::vector<types::global_dof_index> level_dof_indices (fe.
dofs_per_cell);
621 copy_to_and_from_indices[level].clear();
627 temp_copy_indices.resize (0);
632 for (; level_cell!=level_end; ++level_cell)
638 level_cell->get_dof_indices(global_dof_indices);
639 level_cell->get_mg_dof_indices (level_dof_indices);
643 const unsigned int component
645 if (component_mask[component] &&
646 !interface_dofs[level].is_element(level_dof_indices[i]))
649 = mg_component_start[level][mg_target_component[component]];
651 = component_start[target_component[component]];
652 temp_copy_indices[level_dof_indices[i]-level_start] =
653 global_dof_indices[i] - global_start;
661 std::count_if (temp_copy_indices.begin(), temp_copy_indices.end(),
662 std::bind (std::not_equal_to<types::global_dof_index>(),
663 std::placeholders::_1,
665 copy_to_and_from_indices[level].resize (n_active_dofs);
669 copy_to_and_from_indices[level][counter++] =
670 std::pair<types::global_dof_index, unsigned int> (temp_copy_indices[i], i);
678 #include "mg_transfer_component.inst" 681 DEAL_II_NAMESPACE_CLOSE
void build_matrices(const DoFHandler< dim, spacedim > &dof, const DoFHandler< dim, spacedim > &mg_dof, unsigned int selected, unsigned int mg_selected, const std::vector< unsigned int > &target_component=std::vector< unsigned int >(), const std::vector< unsigned int > &mg_target_component=std::vector< unsigned int >(), const std::vector< std::set< types::global_dof_index > > &boundary_indices=std::vector< std::set< types::global_dof_index > >())
std::vector< std::shared_ptr< BlockSparseMatrix< double > > > prolongation_matrices
cell_iterator begin(const unsigned int level=0) const
std::vector< std::set< types::global_dof_index > > boundary_indices
cell_iterator end() const
std::vector< unsigned int > target_component
bool represents_n_components(const unsigned int n) const
ActiveSelector::active_cell_iterator active_cell_iterator
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
active_cell_iterator begin_active(const unsigned int level=0) const
static ::ExceptionBase & ExcMessage(std::string arg1)
ActiveSelector::cell_iterator cell_iterator
ComponentMask mg_component_mask
unsigned int global_dof_index
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< types::global_dof_index > component_start
active_cell_iterator end_active(const unsigned int level) const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
std::vector< std::vector< std::pair< types::global_dof_index, unsigned int > > > copy_to_and_from_indices
types::global_dof_index n_dofs() const
size_type local_to_global(const unsigned int block, const size_type index) const
const unsigned int dofs_per_cell
void swap(Vector< Number > &u, Vector< Number > &v)
std::vector< unsigned int > mg_target_component
void do_copy_to_mg(const DoFHandler< dim, spacedim > &mg_dof, MGLevelObject< Vector< number > > &dst, const InVector &src) const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
unsigned int n_components() const
const Triangulation< dim, spacedim > & get_triangulation() const
std::vector< std::vector< types::global_dof_index > > mg_component_start
const types::global_dof_index invalid_dof_index
void build_matrices(const DoFHandler< dim, spacedim > &dof, const DoFHandler< dim, spacedim > &mg_dof)
std::vector< std::vector< types::global_dof_index > > sizes
static ::ExceptionBase & ExcInternalError()