Reference documentation for deal.II version 8.4.2
mg_transfer_block.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2015 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 
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>
31 
32 #include <algorithm>
33 #include <numeric>
34 #include <iostream>
35 #include <utility>
36 
37 DEAL_II_NAMESPACE_OPEN
38 
39 namespace
40 {
48  template <int dim, typename number, int spacedim>
49  void
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)
55  {
56  std::vector<bool> selected=sel;
57  // Compute the number of blocks needed
58  const unsigned int n_selected
59  = std::accumulate(selected.begin(),
60  selected.end(),
61  0U);
62 
63  if (ndofs.size() == 0)
64  {
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);
69  MGTools::count_dofs_per_block (mg_dof, ndofs);
70  }
71 
72  for (unsigned int level=v.min_level();
73  level<=v.max_level(); ++level)
74  {
75  v[level].reinit(n_selected, 0);
76  unsigned int k=0;
77  for (unsigned int i=0; i<selected.size() && (k<v[level].n_blocks()); ++i)
78  {
79  if (selected[i])
80  {
81  v[level].block(k++).reinit(ndofs[level][i]);
82  }
83  v[level].collect_sizes();
84  }
85  }
86  }
87 
88 
95  template <int dim, typename number, int spacedim>
96  void
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)
102  {
103  const unsigned int n_blocks = mg_dof.get_fe().n_blocks();
104  Assert(selected_block < n_blocks, ExcIndexRange(selected_block, 0, n_blocks));
105 
106  std::vector<bool> selected(n_blocks, false);
107  selected[selected_block] = true;
108 
109  if (ndofs.size() == 0)
110  {
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);
115  MGTools::count_dofs_per_block (mg_dof, ndofs);
116  }
117 
118  for (unsigned int level=v.min_level();
119  level<=v.max_level(); ++level)
120  {
121  v[level].reinit(ndofs[level][selected_block]);
122  }
123  }
124 }
125 
126 
127 template <typename number>
128 template <int dim, typename number2, int spacedim>
129 void
131  const DoFHandler<dim,spacedim> &mg_dof_handler,
133  const BlockVector<number2> &src) const
134 {
135  reinit_vector_by_blocks(mg_dof_handler, dst, selected_block, sizes);
136  // For MGTransferBlockSelect, the
137  // multilevel block is always the
138  // first, since only one block is
139  // selected.
140  bool first = true;
141  for (unsigned int level=mg_dof_handler.get_triangulation().n_levels(); level != 0;)
142  {
143  --level;
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);
147  if (!first)
148  restrict_and_add (level+1, dst[level], dst[level+1]);
149  first = false;
150  }
151 }
152 
153 
154 
155 template <typename number>
156 template <int dim, typename number2, int spacedim>
157 void
159  const DoFHandler<dim,spacedim> &mg_dof_handler,
161  const Vector<number2> &src) const
162 {
163  reinit_vector_by_blocks(mg_dof_handler, dst, selected_block, sizes);
164  // For MGTransferBlockSelect, the
165  // multilevel block is always the
166  // first, since only one block is selected.
167  bool first = true;
168  for (unsigned int level=mg_dof_handler.get_triangulation().n_levels(); level != 0;)
169  {
170  --level;
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);
174  if (!first)
175  restrict_and_add (level+1, dst[level], dst[level+1]);
176  first = false;
177  }
178 }
179 
180 
181 
182 template <typename number>
183 template <int dim, typename number2, int spacedim>
184 void
186  const DoFHandler<dim,spacedim> &mg_dof_handler,
188  const BlockVector<number2> &src) const
189 {
190  reinit_vector_by_blocks(mg_dof_handler, dst, selected, sizes);
191  bool first = true;
192  for (unsigned int level=mg_dof_handler.get_triangulation().n_levels(); level != 0;)
193  {
194  --level;
195  for (unsigned int block=0; block<selected.size(); ++block)
196  if (selected[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);
200  if (!first)
201  restrict_and_add (level+1, dst[level], dst[level+1]);
202  first = false;
203  }
204 }
205 
206 
207 
208 template <int dim, int spacedim>
210  const DoFHandler<dim,spacedim> &,
211  const DoFHandler<dim,spacedim> &mg_dof)
212 {
213  const FiniteElement<dim> &fe = mg_dof.get_fe();
214  const unsigned int n_blocks = fe.n_blocks();
215  const unsigned int dofs_per_cell = fe.dofs_per_cell;
216  const unsigned int n_levels = mg_dof.get_triangulation().n_levels();
217 
218  Assert (selected.size() == n_blocks,
219  ExcDimensionMismatch(selected.size(), n_blocks));
220 
221  // Compute the mapping between real
222  // blocks and blocks used for
223  // multigrid computations.
224  mg_block.resize(n_blocks);
225  n_mg_blocks = 0;
226  for (unsigned int i=0; i<n_blocks; ++i)
227  if (selected[i])
228  mg_block[i] = n_mg_blocks++;
229  else
231 
232  // Compute the lengths of all blocks
233  sizes.clear ();
234  sizes.resize(n_levels, std::vector<types::global_dof_index>(fe.n_blocks()));
236 
237  // Fill some index vectors
238  // for later use.
240  // Compute start indices from sizes
241  for (unsigned int l=0; l<mg_block_start.size(); ++l)
242  {
244  for (unsigned int i=0; i<mg_block_start[l].size(); ++i)
245  {
247  mg_block_start[l][i] = k;
248  k += t;
249  }
250  }
251 
252  block_start.resize(n_blocks);
253  DoFTools::count_dofs_per_block (static_cast<const DoFHandler<dim,spacedim>&>(mg_dof),
254  block_start);
255 
257  for (unsigned int i=0; i<block_start.size(); ++i)
258  {
260  block_start[i] = k;
261  k += t;
262  }
263  // Build index vectors for
264  // copy_to_mg and
265  // copy_from_mg. These vectors must
266  // be prebuilt, since the
267  // get_dof_indices functions are
268  // too slow
269  copy_indices.resize(n_blocks);
270  for (unsigned int block=0; block<n_blocks; ++block)
271  if (selected[block])
272  copy_indices[block].resize(n_levels);
273 
274 // Building the prolongation matrices starts here!
275 
276  // reset the size of the array of
277  // matrices. call resize(0) first,
278  // in order to delete all elements
279  // and clear their memory. then
280  // repopulate these arrays
281  //
282  // note that on resize(0), the
283  // shared_ptr class takes care of
284  // deleting the object it points to
285  // by itself
286  prolongation_matrices.resize (0);
287  prolongation_sparsities.resize (0);
288 
289  for (unsigned int i=0; i<n_levels-1; ++i)
290  {
291  prolongation_sparsities
292  .push_back (std_cxx11::shared_ptr<BlockSparsityPattern> (new BlockSparsityPattern));
294  .push_back (std_cxx11::shared_ptr<BlockSparseMatrix<double> > (new BlockSparseMatrix<double>));
295  }
296 
297  // two fields which will store the
298  // indices of the multigrid dofs
299  // for a cell and one of its children
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);
302 
303  // for each level: first build the
304  // sparsity pattern of the matrices
305  // and then build the matrices
306  // themselves. note that we only
307  // need to take care of cells on
308  // the coarser level which have
309  // children
310 
311  for (unsigned int level=0; level<n_levels-1; ++level)
312  {
313  // reset the dimension of the
314  // structure. note that for
315  // the number of entries per
316  // row, the number of parent
317  // dofs coupling to a child dof
318  // is necessary. this, is the
319  // number of degrees of freedom
320  // per cell
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)
324  if (i==j)
325  prolongation_sparsities[level]->block(i,j)
326  .reinit(sizes[level+1][i],
327  sizes[level][j],
328  dofs_per_cell+1);
329  else
330  prolongation_sparsities[level]->block(i,j)
331  .reinit(sizes[level+1][i],
332  sizes[level][j],
333  0);
334 
335  prolongation_sparsities[level]->collect_sizes();
336 
337  for (typename DoFHandler<dim,spacedim>::cell_iterator cell=mg_dof.begin(level);
338  cell != mg_dof.end(level); ++cell)
339  if (cell->has_children())
340  {
341  cell->get_mg_dof_indices (dof_indices_parent);
342 
344  ExcNotImplemented());
345  for (unsigned int child=0; child<cell->n_children(); ++child)
346  {
347  // set an alias to the
348  // prolongation matrix for
349  // this child
350  const FullMatrix<double> &prolongation
351  = mg_dof.get_fe().get_prolongation_matrix (child, cell->refinement_case());
352 
353  cell->child(child)->get_mg_dof_indices (dof_indices_child);
354 
355  // now tag the entries in the
356  // matrix which will be used
357  // for this pair of parent/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)
361  {
362  const unsigned int icomp
363  = fe.system_to_block_index(i).first;
364  const unsigned int jcomp
365  = fe.system_to_block_index(j).first;
366  if ((icomp==jcomp) && selected[icomp])
367  prolongation_sparsities[level]->add(dof_indices_child[i],
368  dof_indices_parent[j]);
369  };
370  };
371  };
372  prolongation_sparsities[level]->compress ();
373 
374  prolongation_matrices[level]->reinit (*prolongation_sparsities[level]);
375  // now actually build the matrices
376  for (typename DoFHandler<dim,spacedim>::cell_iterator cell=mg_dof.begin(level);
377  cell != mg_dof.end(level); ++cell)
378  if (cell->has_children())
379  {
380  cell->get_mg_dof_indices (dof_indices_parent);
381 
383  ExcNotImplemented());
384  for (unsigned int child=0; child<cell->n_children(); ++child)
385  {
386  // set an alias to the
387  // prolongation matrix for
388  // this child
389  const FullMatrix<double> &prolongation
390  = mg_dof.get_fe().get_prolongation_matrix (child, cell->refinement_case());
391 
392  cell->child(child)->get_mg_dof_indices (dof_indices_child);
393 
394  // now set the entries in the
395  // matrix
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)
399  {
400  const unsigned int icomp = fe.system_to_block_index(i).first;
401  const unsigned int jcomp = fe.system_to_block_index(j).first;
402  if ((icomp==jcomp) && selected[icomp])
403  prolongation_matrices[level]->set(dof_indices_child[i],
404  dof_indices_parent[j],
405  prolongation(i,j));
406  }
407  }
408  }
409  }
410  // impose boundary conditions
411  // but only in the column of
412  // the prolongation matrix
414  {
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)
418  {
420  continue;
421 
422  // need to delete all the columns in the
423  // matrix that are on the boundary. to achieve
424  // this, create an array as long as there are
425  // matrix columns, and find which columns we
426  // need to filter away.
427  constrain_indices.resize (0);
428  constrain_indices.resize (prolongation_matrices[level]->n(), 0);
432  for (; dof != endd; ++dof)
433  constrain_indices[*dof] = 1;
434 
435  unsigned int index = 0;
436  for (unsigned int block=0; block<n_blocks; ++block)
437  {
438  const types::global_dof_index n_dofs = prolongation_matrices[level]->block(block, block).m();
439  constraints_per_block[block].resize(0);
440  constraints_per_block[block].resize(n_dofs, 0);
441  for (types::global_dof_index i=0; i<n_dofs; ++i, ++index)
442  constraints_per_block[block][i] = (constrain_indices[index] == 1);
443 
444  for (types::global_dof_index i=0; i<n_dofs; ++i)
445  {
447  start_row = prolongation_matrices[level]->block(block, block).begin(i),
448  end_row = prolongation_matrices[level]->block(block, block).end(i);
449  for (; start_row != end_row; ++start_row)
450  {
451  if (constraints_per_block[block][start_row->column()])
452  start_row->value() = 0.;
453  }
454  }
455  }
456  }
457  }
458 }
459 
460 
461 
462 template <typename number>
463 template <int dim, int spacedim>
465  const DoFHandler<dim,spacedim> &dof,
466  const DoFHandler<dim,spacedim> &mg_dof,
467  unsigned int select)
468 {
469  const FiniteElement<dim> &fe = mg_dof.get_fe();
470  unsigned int n_blocks = mg_dof.get_fe().n_blocks();
471 
472  selected_block = select;
473  selected.resize(n_blocks, false);
474  selected[select] = true;
475 
477 
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);
481 
482  for (int level=dof.get_triangulation().n_levels()-1; level>=0; --level)
483  {
485  level_cell = mg_dof.begin_active(level);
487  level_end = mg_dof.end_active(level);
488 
489  temp_copy_indices.resize (0);
490  temp_copy_indices.resize (sizes[level][selected_block],
492 
493  // Compute coarse level right hand side
494  // by restricting from fine level.
495  for (; level_cell!=level_end; ++level_cell)
496  {
497  // get the dof numbers of
498  // this cell for the global
499  // and the level-wise
500  // numbering
501  level_cell->get_dof_indices(global_dof_indices);
502  level_cell->get_mg_dof_indices (level_dof_indices);
503 
504  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
505  {
506  const unsigned int block = fe.system_to_block_index(i).first;
507  if (selected[block])
508  {
509  if (mg_constrained_dofs != 0)
510  {
511  if (!mg_constrained_dofs->at_refinement_edge(level,level_dof_indices[i]))
512  temp_copy_indices[level_dof_indices[i] - mg_block_start[level][block]]
513  = global_dof_indices[i] - block_start[block];
514  }
515  else
516  temp_copy_indices[level_dof_indices[i] - mg_block_start[level][block]]
517  = global_dof_indices[i] - block_start[block];
518  }
519  }
520  }
521 
522  // now all the active dofs got a valid entry,
523  // the other ones have an invalid entry. Count
524  // the invalid entries and then resize the
525  // copy_indices object. Then, insert the pairs
526  // of global index and level index into
527  // copy_indices.
528  const types::global_dof_index n_active_dofs =
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);
533  types::global_dof_index counter = 0;
534  for (types::global_dof_index i=0; i<temp_copy_indices.size(); ++i)
535  if (temp_copy_indices[i] != numbers::invalid_dof_index)
536  copy_indices[selected_block][level][counter++] =
537  std::pair<types::global_dof_index, unsigned int> (temp_copy_indices[i], i);
538  Assert (counter == n_active_dofs, ExcInternalError());
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  const std::vector<bool> &sel)
551 {
552  const FiniteElement<dim> &fe = mg_dof.get_fe();
553  unsigned int n_blocks = mg_dof.get_fe().n_blocks();
554 
555  if (sel.size() != 0)
556  {
557  Assert(sel.size() == n_blocks,
558  ExcDimensionMismatch(sel.size(), n_blocks));
559  selected = sel;
560  }
561  if (selected.size() == 0)
562  selected = std::vector<bool> (n_blocks, true);
563 
565 
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);
569  for (int level=dof.get_triangulation().n_levels()-1; level>=0; --level)
570  {
572  level_cell = mg_dof.begin_active(level);
574  level_end = mg_dof.end_active(level);
575 
576  for (unsigned int block=0; block<n_blocks; ++block)
577  if (selected[block])
578  {
579  temp_copy_indices[block].resize (0);
580  temp_copy_indices[block].resize (sizes[level][block],
582  }
583 
584  // Compute coarse level right hand side
585  // by restricting from fine level.
586  for (; level_cell!=level_end; ++level_cell)
587  {
588  // get the dof numbers of
589  // this cell for the global
590  // and the level-wise
591  // numbering
592  level_cell->get_dof_indices(global_dof_indices);
593  level_cell->get_mg_dof_indices (level_dof_indices);
594 
595  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
596  {
597  const unsigned int block = fe.system_to_block_index(i).first;
598  if (selected[block])
599  temp_copy_indices[block][level_dof_indices[i] - mg_block_start[level][block]]
600  = global_dof_indices[i] - block_start[block];
601  }
602  }
603 
604  for (unsigned int block=0; block<n_blocks; ++block)
605  if (selected[block])
606  {
607  const types::global_dof_index n_active_dofs =
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>(),
612  copy_indices[block][level].resize (n_active_dofs);
613  types::global_dof_index counter = 0;
614  for (types::global_dof_index i=0; i<temp_copy_indices[block].size(); ++i)
615  if (temp_copy_indices[block][i] != numbers::invalid_dof_index)
616  copy_indices[block][level][counter++] =
617  std::pair<types::global_dof_index, unsigned int>
618  (temp_copy_indices[block][i], i);
619  Assert (counter == n_active_dofs, ExcInternalError());
620  }
621  }
622 }
623 
624 
625 
626 // explicit instantiations
627 #include "mg_transfer_block.inst"
628 
629 
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
Definition: types.h:164
ElementIterator end() const
Definition: index_set.h:1175
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:833
cell_iterator end() const
Definition: dof_handler.cc:861
ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:221
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
Definition: dof_handler.cc:845
const IndexSet & get_boundary_indices(const unsigned int level) const
ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:249
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
Definition: types.h:88
std::vector< std_cxx11::shared_ptr< BlockSparseMatrix< double > > > prolongation_matrices
#define Assert(cond, exc)
Definition: exceptions.h:294
active_cell_iterator end_active(const unsigned int level) const
Definition: dof_handler.cc:883
const unsigned int dofs_per_cell
Definition: fe_base.h:283
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
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:1661
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
Definition: fe.h:2885
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:1069
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
Definition: types.h:178
ElementIterator begin() const
Definition: index_set.h:1165
BlockType & block(const unsigned int i)
size_type n_elements() const
Definition: index_set.h:1380