Reference documentation for deal.II version 9.0.0
mg_transfer_component.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/logstream.h>
18 #include <deal.II/base/function.h>
19 
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>
33 
34 #include <algorithm>
35 #include <numeric>
36 #include <iostream>
37 
38 DEAL_II_NAMESPACE_OPEN
39 
40 
41 namespace
42 {
77  template <int dim, typename number, int spacedim>
78  void
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)
85  {
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();
89 
90  // If the selected and
91  // target_component have size 0,
92  // they must be replaced by default
93  // values.
94  //
95  // Since we already made copies
96  // directly after this function was
97  // called, we use the arguments
98  // directly.
99  if (target_component.size() == 0)
100  {
101  target_component.resize(ncomp);
102  for (unsigned int i=0; i<ncomp; ++i)
103  target_component[i] = i;
104  }
105 
106  // If selected is an empty vector,
107  // all components are selected.
108  if (selected.size() == 0)
109  {
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;
114  }
115 
116  Assert (selected.size() == target_component.size(),
117  ExcDimensionMismatch(selected.size(), target_component.size()));
118 
119  // Compute the number of blocks needed
120  const unsigned int n_selected
121  = std::accumulate(selected.begin(),
122  selected.end(),
123  0u);
124 
125  if (ndofs.size() == 0)
126  {
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()));
130  std::swap(ndofs, new_dofs);
131  MGTools::count_dofs_per_block (mg_dof, ndofs, target_component);
132  }
133 
134  for (unsigned int level=v.min_level();
135  level<=v.max_level(); ++level)
136  {
137  v[level].reinit(n_selected, 0);
138  unsigned int k=0;
139  for (unsigned int i=0; i<selected.size() && (k<v[level].n_blocks()); ++i)
140  {
141  if (selected[i])
142  {
143  v[level].block(k++).reinit(ndofs[level][i]);
144  }
145  v[level].collect_sizes();
146  }
147  }
148  }
149 
150 
177  template <int dim, typename number, int spacedim>
178  void
179  reinit_vector_by_components (
180  const ::DoFHandler<dim,spacedim> &mg_dof,
182  const ComponentMask &component_mask,
183  const std::vector<unsigned int> &target_component,
184  std::vector<std::vector<types::global_dof_index> > &ndofs)
185  {
186  Assert (component_mask.represents_n_components(target_component.size()),
187  ExcMessage ("The component mask does not have the correct size."));
188 
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];
193 
194  if (ndofs.size() == 0)
195  {
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()));
199  std::swap(ndofs, new_dofs);
200  MGTools::count_dofs_per_block (mg_dof, ndofs,
201  target_component);
202  }
203 
204  for (unsigned int level=v.min_level();
205  level<=v.max_level(); ++level)
206  {
207  v[level].reinit(ndofs[level][selected_block]);
208  }
209  }
210 }
211 
212 
213 template <typename number>
214 template <int dim, class InVector, int spacedim>
215 void
217  const DoFHandler<dim,spacedim> &mg_dof_handler,
219  const InVector &src) const
220 {
221  dst=0;
222 
223  Assert(sizes.size()==mg_dof_handler.get_triangulation().n_levels(),
224  ExcMatricesNotBuilt());
225 
226  reinit_vector_by_components(mg_dof_handler, dst,
227  mg_component_mask,
228  mg_target_component, sizes);
229 
230  // traverse the grid top-down
231  // (i.e. starting with the most
232  // refined grid). this way, we can
233  // always get that part of one
234  // level of the output vector which
235  // corresponds to a region which is
236  // more refined, by restriction of
237  // the respective vector on the
238  // next finer level, which we then
239  // already have built.
240 
241  for (unsigned int level=mg_dof_handler.get_triangulation().n_levels(); level!=0;)
242  {
243  --level;
244 
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);
249  }
250 }
251 
252 
253 template <int dim, int spacedim>
255  const DoFHandler<dim,spacedim> &,
256  const DoFHandler<dim,spacedim> &mg_dof)
257 {
258  // Fill target component with
259  // standard values (identity) if it
260  // is empty
261  if (target_component.size() == 0)
262  {
263  target_component.resize(mg_dof.get_fe(0).n_components());
264  for (unsigned int i=0; i<target_component.size(); ++i)
265  target_component[i] = i;
266  }
267  else
268  {
269  // otherwise, check it for consistency
270  Assert (target_component.size() == mg_dof.get_fe(0).n_components(),
272  mg_dof.get_fe(0).n_components()));
273 
274  for (unsigned int i=0; i<target_component.size(); ++i)
275  {
276  Assert(i<target_component.size(),
277  ExcIndexRange(i,0,target_component.size()));
278  }
279  }
280  // Do the same for the multilevel
281  // components. These may be
282  // different.
283  if (mg_target_component.size() == 0)
284  {
285  mg_target_component.resize(mg_dof.get_fe(0).n_components());
286  for (unsigned int i=0; i<mg_target_component.size(); ++i)
288  }
289  else
290  {
291  Assert (mg_target_component.size() == mg_dof.get_fe(0).n_components(),
293  mg_dof.get_fe(0).n_components()));
294 
295  for (unsigned int i=0; i<mg_target_component.size(); ++i)
296  {
297  Assert(i<mg_target_component.size(),
298  ExcIndexRange(i,0,mg_target_component.size()));
299  }
300  }
301 
302  const FiniteElement<dim> &fe = mg_dof.get_fe();
303 
304  // Effective number of components
305  // is the maximum entry in
306  // mg_target_component. This
307  // assumes that the values in that
308  // vector don't have holes.
309  const unsigned int n_components =
310  *std::max_element(mg_target_component.begin(), mg_target_component.end()) + 1;
311  const unsigned int dofs_per_cell = fe.dofs_per_cell;
312  const unsigned int n_levels = mg_dof.get_triangulation().n_levels();
313 
315  ExcMessage ("Component mask has wrong size."));
316 
317  // Compute the lengths of all blocks
318  sizes.resize(n_levels);
319  for (unsigned int l=0; l<n_levels; ++l)
320  sizes[l].resize(n_components);
321 
323 
324  // Fill some index vectors
325  // for later use.
327  // Compute start indices from sizes
328  for (unsigned int l=0; l<mg_component_start.size(); ++l)
329  {
331  for (unsigned int i=0; i<mg_component_start[l].size(); ++i)
332  {
334  mg_component_start[l][i] = k;
335  k += t;
336  }
337  }
338 
339  component_start.resize(*std::max_element (target_component.begin(),
340  target_component.end()) + 1);
341  DoFTools::
343 
345  for (unsigned int i=0; i<component_start.size(); ++i)
346  {
348  component_start[i] = k;
349  k += t;
350  }
351 
352  // Build index vectors for
353  // copy_to_mg and
354  // copy_from_mg. These vectors must
355  // be prebuilt, since the
356  // get_dof_indices functions are
357  // too slow
358 
359  copy_to_and_from_indices.resize(n_levels);
360 
361 // Building the prolongation matrices starts here!
362 
363  // reset the size of the array of
364  // matrices. call resize(0) first,
365  // in order to delete all elements
366  // and clear their memory. then
367  // repopulate these arrays
368  //
369  // note that on resize(0), the
370  // shared_ptr class takes care of
371  // deleting the object it points to
372  // by itself
373  prolongation_matrices.resize (0);
374  prolongation_sparsities.resize (0);
375  prolongation_matrices.reserve (n_levels - 1);
376  prolongation_sparsities.reserve (n_levels - 1);
377 
378  for (unsigned int i=0; i<n_levels-1; ++i)
379  {
380  prolongation_sparsities.emplace_back (new BlockSparsityPattern);
382  }
383 
384  // two fields which will store the
385  // indices of the multigrid dofs
386  // for a cell and one of its children
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);
389 
390  // for each level: first build the
391  // sparsity pattern of the matrices
392  // and then build the matrices
393  // themselves. note that we only
394  // need to take care of cells on
395  // the coarser level which have
396  // children
397  for (unsigned int level=0; level<n_levels-1; ++level)
398  {
399  // reset the dimension of the
400  // structure. note that for
401  // the number of entries per
402  // row, the number of parent
403  // dofs coupling to a child dof
404  // is necessary. this, is the
405  // number of degrees of freedom
406  // per cell
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)
410  if (i==j)
411  prolongation_sparsities[level]->block(i,j)
412  .reinit(sizes[level+1][i],
413  sizes[level][j],
414  dofs_per_cell+1);
415  else
416  prolongation_sparsities[level]->block(i,j)
417  .reinit(sizes[level+1][i],
418  sizes[level][j],
419  0);
420 
421  prolongation_sparsities[level]->collect_sizes();
422 
423  for (typename DoFHandler<dim,spacedim>::cell_iterator cell=mg_dof.begin(level);
424  cell != mg_dof.end(level); ++cell)
425  if (cell->has_children())
426  {
427  cell->get_mg_dof_indices (dof_indices_parent);
428 
429  for (unsigned int child=0; child<cell->n_children(); ++child)
430  {
431  // set an alias to the
432  // prolongation matrix for
433  // this child
434  const FullMatrix<double> &prolongation
435  = mg_dof.get_fe().get_prolongation_matrix (child, cell->refinement_case());
436 
437  cell->child(child)->get_mg_dof_indices (dof_indices_child);
438 
439  // now tag the entries in the
440  // matrix which will be used
441  // for this pair of parent/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)
445  {
446  const unsigned int icomp
447  = fe.system_to_component_index(i).first;
448  const unsigned int jcomp
449  = fe.system_to_component_index(j).first;
450  if ((icomp==jcomp) && mg_component_mask[icomp])
451  prolongation_sparsities[level]->add(dof_indices_child[i],
452  dof_indices_parent[j]);
453  };
454  };
455  };
456  prolongation_sparsities[level]->compress ();
457 
458  prolongation_matrices[level]->reinit (*prolongation_sparsities[level]);
459  // now actually build the matrices
460  for (typename DoFHandler<dim,spacedim>::cell_iterator cell=mg_dof.begin(level);
461  cell != mg_dof.end(level); ++cell)
462  if (cell->has_children())
463  {
464  cell->get_mg_dof_indices (dof_indices_parent);
465 
466  for (unsigned int child=0; child<cell->n_children(); ++child)
467  {
468  // set an alias to the
469  // prolongation matrix for
470  // this child
471  const FullMatrix<double> &prolongation
472  = mg_dof.get_fe().get_prolongation_matrix (child, cell->refinement_case());
473 
474  cell->child(child)->get_mg_dof_indices (dof_indices_child);
475 
476  // now set the entries in the
477  // matrix
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)
481  {
482  const unsigned int icomp = fe.system_to_component_index(i).first;
483  const unsigned int jcomp = fe.system_to_component_index(j).first;
484  if ((icomp==jcomp) && mg_component_mask[icomp])
485  prolongation_matrices[level]->set(dof_indices_child[i],
486  dof_indices_parent[j],
487  prolongation(i,j));
488  }
489  }
490  }
491  }
492  // impose boundary conditions
493  // but only in the column of
494  // the prolongation matrix
495  //TODO: this way is not very efficient
496 
497  if (boundary_indices.size() != 0)
498  {
499  std::vector<std::vector<types::global_dof_index> >
500  dofs_per_component(mg_dof.get_triangulation().n_levels(),
501  std::vector<types::global_dof_index>(n_components));
502 
503  MGTools::count_dofs_per_block (mg_dof, dofs_per_component, mg_target_component);
504  for (unsigned int level=0; level<n_levels-1; ++level)
505  {
506  if (boundary_indices[level].size() == 0)
507  continue;
508 
509  for (unsigned int iblock=0; iblock<n_components; ++iblock)
510  for (unsigned int jblock=0; jblock<n_components; ++jblock)
511  if (iblock==jblock)
512  {
513  const types::global_dof_index n_dofs = prolongation_matrices[level]->block(iblock,jblock).m();
514  for (types::global_dof_index i=0; i<n_dofs; ++i)
515  {
516  SparseMatrix<double>::iterator anfang = prolongation_matrices[level]->block(iblock,jblock).begin(i),
517  ende = prolongation_matrices[level]->block(iblock,jblock).end(i);
518  for (; anfang != ende; ++anfang)
519  {
520  const types::global_dof_index column_number = anfang->column();
521 
522  //convert global indices into local ones
523  const BlockIndices block_indices_coarse (dofs_per_component[level]);
524  const types::global_dof_index global_j = block_indices_coarse.local_to_global(iblock, column_number);
525 
526  std::set<types::global_dof_index>::const_iterator found_dof =
527  boundary_indices[level].find(global_j);
528 
529  const bool is_boundary_index =
530  (found_dof != boundary_indices[level].end());
531 
532  if (is_boundary_index)
533  {
534  prolongation_matrices[level]->block(iblock,jblock)
535  .set(i,column_number,0);
536  }
537  }
538  }
539  }
540  }
541  }
542 }
543 
544 
545 template <typename number>
546 template <int dim, int spacedim>
548  const DoFHandler<dim,spacedim> &dof,
549  const DoFHandler<dim,spacedim> &mg_dof,
550  unsigned int select,
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)
555 {
556  const FiniteElement<dim> &fe = mg_dof.get_fe();
557  unsigned int ncomp = mg_dof.get_fe(0).n_components();
558 
559  target_component = t_component;
560  mg_target_component = mg_t_component;
561  boundary_indices = bdry_indices;
562 
563  selected_component = select;
564  mg_selected_component = mg_select;
565 
566  {
567  std::vector<bool> tmp(ncomp, false);
568  for (unsigned int c=0; c<ncomp; ++c)
569  if (t_component[c] == selected_component)
570  tmp[c] = true;
571  component_mask = ComponentMask(tmp);
572  }
573 
574  {
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)
578  tmp[c] = true;
579  mg_component_mask = ComponentMask(tmp);
580  }
581 
582  // If components are renumbered,
583  // find the first original
584  // component corresponding to the
585  // target component.
586  for (unsigned int i=0; i<target_component.size(); ++i)
587  {
588  if (target_component[i] == select)
589  {
590  selected_component = i;
591  break;
592  }
593  }
594 
595  for (unsigned int i=0; i<mg_target_component.size(); ++i)
596  {
597  if (mg_target_component[i] == mg_select)
598  {
599  mg_selected_component = i;
600  break;
601  }
602  }
603 
605 
606  interface_dofs.resize(mg_dof.get_triangulation().n_levels());
607  for (unsigned int l=0; l<mg_dof.get_triangulation().n_levels(); ++l)
608  {
609  interface_dofs[l].clear();
610  interface_dofs[l].set_size(mg_dof.n_dofs(l));
611  }
612  MGTools::extract_inner_interface_dofs(mg_dof, interface_dofs);
613 
614  // use a temporary vector to create the
615  // relation between global and level dofs
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);
619  for (int level=dof.get_triangulation().n_levels()-1; level>=0; --level)
620  {
621  copy_to_and_from_indices[level].clear();
623  level_cell = mg_dof.begin_active(level);
625  level_end = mg_dof.end_active(level);
626 
627  temp_copy_indices.resize (0);
628  temp_copy_indices.resize (mg_dof.n_dofs(level), numbers::invalid_dof_index);
629 
630  // Compute coarse level right hand side
631  // by restricting from fine level.
632  for (; level_cell!=level_end; ++level_cell)
633  {
634  // get the dof numbers of
635  // this cell for the global
636  // and the level-wise
637  // numbering
638  level_cell->get_dof_indices(global_dof_indices);
639  level_cell->get_mg_dof_indices (level_dof_indices);
640 
641  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
642  {
643  const unsigned int component
644  = fe.system_to_component_index(i).first;
645  if (component_mask[component] &&
646  !interface_dofs[level].is_element(level_dof_indices[i]))
647  {
648  const types::global_dof_index level_start
649  = mg_component_start[level][mg_target_component[component]];
650  const types::global_dof_index global_start
651  = component_start[target_component[component]];
652  temp_copy_indices[level_dof_indices[i]-level_start] =
653  global_dof_indices[i] - global_start;
654  }
655  }
656  }
657 
658  // write indices from vector into the map from
659  // global to level dofs
660  const types::global_dof_index n_active_dofs =
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);
666  types::global_dof_index counter = 0;
667  for (types::global_dof_index i=0; i<temp_copy_indices.size(); ++i)
668  if (temp_copy_indices[i] != numbers::invalid_dof_index)
669  copy_to_and_from_indices[level][counter++] =
670  std::pair<types::global_dof_index, unsigned int> (temp_copy_indices[i], i);
671  Assert (counter == n_active_dofs, ExcInternalError());
672  }
673 }
674 
675 
676 
677 // explicit instantiations
678 #include "mg_transfer_component.inst"
679 
680 
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
Definition: dof_handler.cc:723
std::vector< std::set< types::global_dof_index > > boundary_indices
cell_iterator end() const
Definition: dof_handler.cc:751
std::vector< unsigned int > target_component
bool represents_n_components(const unsigned int n) const
ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:229
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: dof_handler.cc:735
static ::ExceptionBase & ExcMessage(std::string arg1)
ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:257
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:1142
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
Definition: dof_handler.cc:773
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
Definition: fe_base.h:295
void swap(Vector< Number > &u, Vector< Number > &v)
Definition: vector.h:1312
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
Definition: fe.h:3032
void count_dofs_per_block(const DoFHandlerType &dof, std::vector< types::global_dof_index > &dofs_per_block, const std::vector< unsigned int > &target_block=std::vector< unsigned int >())
Definition: dof_tools.cc:1836
unsigned int n_components() const
void count_dofs_per_block(const DoFHandlerType &dof_handler, std::vector< std::vector< types::global_dof_index > > &dofs_per_block, std::vector< unsigned int > target_block=std::vector< unsigned int >())
Definition: mg_tools.cc:1125
void extract_inner_interface_dofs(const DoFHandler< dim, spacedim > &mg_dof_handler, std::vector< IndexSet > &interface_dofs)
Definition: mg_tools.cc:1466
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
Definition: types.h:187
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()