Reference documentation for deal.II version 9.0.0
mg_transfer_block.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 
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  for (unsigned int level=mg_dof_handler.get_triangulation().n_levels(); level != 0;)
141  {
142  --level;
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);
146  }
147 }
148 
149 
150 
151 template <typename number>
152 template <int dim, typename number2, int spacedim>
153 void
155  const DoFHandler<dim,spacedim> &mg_dof_handler,
157  const Vector<number2> &src) const
158 {
159  reinit_vector_by_blocks(mg_dof_handler, dst, selected_block, sizes);
160  // For MGTransferBlockSelect, the
161  // multilevel block is always the
162  // first, since only one block is selected.
163  for (unsigned int level=mg_dof_handler.get_triangulation().n_levels(); level != 0;)
164  {
165  --level;
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);
169  }
170 }
171 
172 
173 
174 template <typename number>
175 template <int dim, typename number2, int spacedim>
176 void
178  const DoFHandler<dim,spacedim> &mg_dof_handler,
180  const BlockVector<number2> &src) const
181 {
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;)
184  {
185  --level;
186  for (unsigned int block=0; block<selected.size(); ++block)
187  if (selected[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);
191  }
192 }
193 
194 
195 
196 template <int dim, int spacedim>
198  const DoFHandler<dim,spacedim> &,
199  const DoFHandler<dim,spacedim> &mg_dof)
200 {
201  const FiniteElement<dim> &fe = mg_dof.get_fe();
202  const unsigned int n_blocks = fe.n_blocks();
203  const unsigned int dofs_per_cell = fe.dofs_per_cell;
204  const unsigned int n_levels = mg_dof.get_triangulation().n_levels();
205 
206  Assert (selected.size() == n_blocks,
207  ExcDimensionMismatch(selected.size(), n_blocks));
208 
209  // Compute the mapping between real
210  // blocks and blocks used for
211  // multigrid computations.
212  mg_block.resize(n_blocks);
213  n_mg_blocks = 0;
214  for (unsigned int i=0; i<n_blocks; ++i)
215  if (selected[i])
216  mg_block[i] = n_mg_blocks++;
217  else
219 
220  // Compute the lengths of all blocks
221  sizes.clear ();
222  sizes.resize(n_levels, std::vector<types::global_dof_index>(fe.n_blocks()));
224 
225  // Fill some index vectors
226  // for later use.
228  // Compute start indices from sizes
229  for (unsigned int l=0; l<mg_block_start.size(); ++l)
230  {
232  for (unsigned int i=0; i<mg_block_start[l].size(); ++i)
233  {
235  mg_block_start[l][i] = k;
236  k += t;
237  }
238  }
239 
240  block_start.resize(n_blocks);
241  DoFTools::count_dofs_per_block (static_cast<const DoFHandler<dim,spacedim>&>(mg_dof),
242  block_start);
243 
245  for (unsigned int i=0; i<block_start.size(); ++i)
246  {
248  block_start[i] = k;
249  k += t;
250  }
251  // Build index vectors for
252  // copy_to_mg and
253  // copy_from_mg. These vectors must
254  // be prebuilt, since the
255  // get_dof_indices functions are
256  // too slow
257  copy_indices.resize(n_blocks);
258  for (unsigned int block=0; block<n_blocks; ++block)
259  if (selected[block])
260  copy_indices[block].resize(n_levels);
261 
262 // Building the prolongation matrices starts here!
263 
264  // reset the size of the array of
265  // matrices. call resize(0) first,
266  // in order to delete all elements
267  // and clear their memory. then
268  // repopulate these arrays
269  //
270  // note that on resize(0), the
271  // shared_ptr class takes care of
272  // deleting the object it points to
273  // by itself
274  prolongation_matrices.resize (0);
275  prolongation_sparsities.resize (0);
276  prolongation_matrices.reserve (n_levels - 1);
277  prolongation_sparsities.reserve (n_levels - 1);
278 
279  for (unsigned int i=0; i<n_levels-1; ++i)
280  {
281  prolongation_sparsities.emplace_back (new BlockSparsityPattern);
283  }
284 
285  // two fields which will store the
286  // indices of the multigrid dofs
287  // for a cell and one of its children
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);
290 
291  // for each level: first build the
292  // sparsity pattern of the matrices
293  // and then build the matrices
294  // themselves. note that we only
295  // need to take care of cells on
296  // the coarser level which have
297  // children
298 
299  for (unsigned int level=0; level<n_levels-1; ++level)
300  {
301  // reset the dimension of the
302  // structure. note that for
303  // the number of entries per
304  // row, the number of parent
305  // dofs coupling to a child dof
306  // is necessary. this, is the
307  // number of degrees of freedom
308  // per cell
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)
312  if (i==j)
313  prolongation_sparsities[level]->block(i,j)
314  .reinit(sizes[level+1][i],
315  sizes[level][j],
316  dofs_per_cell+1);
317  else
318  prolongation_sparsities[level]->block(i,j)
319  .reinit(sizes[level+1][i],
320  sizes[level][j],
321  0);
322 
323  prolongation_sparsities[level]->collect_sizes();
324 
325  for (typename DoFHandler<dim,spacedim>::cell_iterator cell=mg_dof.begin(level);
326  cell != mg_dof.end(level); ++cell)
327  if (cell->has_children())
328  {
329  cell->get_mg_dof_indices (dof_indices_parent);
330 
333  for (unsigned int child=0; child<cell->n_children(); ++child)
334  {
335  // set an alias to the
336  // prolongation matrix for
337  // this child
338  const FullMatrix<double> &prolongation
339  = mg_dof.get_fe().get_prolongation_matrix (child, cell->refinement_case());
340 
341  cell->child(child)->get_mg_dof_indices (dof_indices_child);
342 
343  // now tag the entries in the
344  // matrix which will be used
345  // for this pair of parent/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)
349  {
350  const unsigned int icomp
351  = fe.system_to_block_index(i).first;
352  const unsigned int jcomp
353  = fe.system_to_block_index(j).first;
354  if ((icomp==jcomp) && selected[icomp])
355  prolongation_sparsities[level]->add(dof_indices_child[i],
356  dof_indices_parent[j]);
357  };
358  };
359  };
360  prolongation_sparsities[level]->compress ();
361 
362  prolongation_matrices[level]->reinit (*prolongation_sparsities[level]);
363  // now actually build the matrices
364  for (typename DoFHandler<dim,spacedim>::cell_iterator cell=mg_dof.begin(level);
365  cell != mg_dof.end(level); ++cell)
366  if (cell->has_children())
367  {
368  cell->get_mg_dof_indices (dof_indices_parent);
369 
372  for (unsigned int child=0; child<cell->n_children(); ++child)
373  {
374  // set an alias to the
375  // prolongation matrix for
376  // this child
377  const FullMatrix<double> &prolongation
378  = mg_dof.get_fe().get_prolongation_matrix (child, cell->refinement_case());
379 
380  cell->child(child)->get_mg_dof_indices (dof_indices_child);
381 
382  // now set the entries in the
383  // matrix
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)
387  {
388  const unsigned int icomp = fe.system_to_block_index(i).first;
389  const unsigned int jcomp = fe.system_to_block_index(j).first;
390  if ((icomp==jcomp) && selected[icomp])
391  prolongation_matrices[level]->set(dof_indices_child[i],
392  dof_indices_parent[j],
393  prolongation(i,j));
394  }
395  }
396  }
397  }
398  // impose boundary conditions
399  // but only in the column of
400  // the prolongation matrix
402  {
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)
406  {
408  continue;
409 
410  // need to delete all the columns in the
411  // matrix that are on the boundary. to achieve
412  // this, create an array as long as there are
413  // matrix columns, and find which columns we
414  // need to filter away.
415  constrain_indices.resize (0);
416  constrain_indices.resize (prolongation_matrices[level]->n(), 0);
420  for (; dof != endd; ++dof)
421  constrain_indices[*dof] = 1;
422 
423  unsigned int index = 0;
424  for (unsigned int block=0; block<n_blocks; ++block)
425  {
426  const types::global_dof_index n_dofs = prolongation_matrices[level]->block(block, block).m();
427  constraints_per_block[block].resize(0);
428  constraints_per_block[block].resize(n_dofs, false);
429  for (types::global_dof_index i=0; i<n_dofs; ++i, ++index)
430  constraints_per_block[block][i] = (constrain_indices[index] == 1);
431 
432  for (types::global_dof_index i=0; i<n_dofs; ++i)
433  {
435  start_row = prolongation_matrices[level]->block(block, block).begin(i),
436  end_row = prolongation_matrices[level]->block(block, block).end(i);
437  for (; start_row != end_row; ++start_row)
438  {
439  if (constraints_per_block[block][start_row->column()])
440  start_row->value() = 0.;
441  }
442  }
443  }
444  }
445  }
446 }
447 
448 
449 
450 template <typename number>
451 template <int dim, int spacedim>
453  const DoFHandler<dim,spacedim> &dof,
454  const DoFHandler<dim,spacedim> &mg_dof,
455  unsigned int select)
456 {
457  const FiniteElement<dim> &fe = mg_dof.get_fe();
458  unsigned int n_blocks = mg_dof.get_fe().n_blocks();
459 
460  selected_block = select;
461  selected.resize(n_blocks, false);
462  selected[select] = true;
463 
465 
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);
469 
470  for (int level=dof.get_triangulation().n_levels()-1; level>=0; --level)
471  {
473  level_cell = mg_dof.begin_active(level);
475  level_end = mg_dof.end_active(level);
476 
477  temp_copy_indices.resize (0);
478  temp_copy_indices.resize (sizes[level][selected_block],
480 
481  // Compute coarse level right hand side
482  // by restricting from fine level.
483  for (; level_cell!=level_end; ++level_cell)
484  {
485  // get the dof numbers of
486  // this cell for the global
487  // and the level-wise
488  // numbering
489  level_cell->get_dof_indices(global_dof_indices);
490  level_cell->get_mg_dof_indices (level_dof_indices);
491 
492  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
493  {
494  const unsigned int block = fe.system_to_block_index(i).first;
495  if (selected[block])
496  {
497  if (mg_constrained_dofs != nullptr)
498  {
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];
502  }
503  else
504  temp_copy_indices[level_dof_indices[i] - mg_block_start[level][block]]
505  = global_dof_indices[i] - block_start[block];
506  }
507  }
508  }
509 
510  // now all the active dofs got a valid entry,
511  // the other ones have an invalid entry. Count
512  // the invalid entries and then resize the
513  // copy_indices object. Then, insert the pairs
514  // of global index and level index into
515  // copy_indices.
516  const types::global_dof_index n_active_dofs =
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);
522  types::global_dof_index counter = 0;
523  for (types::global_dof_index i=0; i<temp_copy_indices.size(); ++i)
524  if (temp_copy_indices[i] != numbers::invalid_dof_index)
525  copy_indices[selected_block][level][counter++] =
526  std::pair<types::global_dof_index, unsigned int> (temp_copy_indices[i], i);
527  Assert (counter == n_active_dofs, ExcInternalError());
528  }
529 }
530 
531 
532 
533 
534 template <typename number>
535 template <int dim, int spacedim>
537  const DoFHandler<dim,spacedim> &dof,
538  const DoFHandler<dim,spacedim> &mg_dof,
539  const std::vector<bool> &sel)
540 {
541  const FiniteElement<dim> &fe = mg_dof.get_fe();
542  unsigned int n_blocks = mg_dof.get_fe().n_blocks();
543 
544  if (sel.size() != 0)
545  {
546  Assert(sel.size() == n_blocks,
547  ExcDimensionMismatch(sel.size(), n_blocks));
548  selected = sel;
549  }
550  if (selected.size() == 0)
551  selected = std::vector<bool> (n_blocks, true);
552 
554 
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);
558  for (int level=dof.get_triangulation().n_levels()-1; level>=0; --level)
559  {
561  level_cell = mg_dof.begin_active(level);
563  level_end = mg_dof.end_active(level);
564 
565  for (unsigned int block=0; block<n_blocks; ++block)
566  if (selected[block])
567  {
568  temp_copy_indices[block].resize (0);
569  temp_copy_indices[block].resize (sizes[level][block],
571  }
572 
573  // Compute coarse level right hand side
574  // by restricting from fine level.
575  for (; level_cell!=level_end; ++level_cell)
576  {
577  // get the dof numbers of
578  // this cell for the global
579  // and the level-wise
580  // numbering
581  level_cell->get_dof_indices(global_dof_indices);
582  level_cell->get_mg_dof_indices (level_dof_indices);
583 
584  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
585  {
586  const unsigned int block = fe.system_to_block_index(i).first;
587  if (selected[block])
588  temp_copy_indices[block][level_dof_indices[i] - mg_block_start[level][block]]
589  = global_dof_indices[i] - block_start[block];
590  }
591  }
592 
593  for (unsigned int block=0; block<n_blocks; ++block)
594  if (selected[block])
595  {
596  const types::global_dof_index n_active_dofs =
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);
603  types::global_dof_index counter = 0;
604  for (types::global_dof_index i=0; i<temp_copy_indices[block].size(); ++i)
605  if (temp_copy_indices[block][i] != numbers::invalid_dof_index)
606  copy_indices[block][level][counter++] =
607  std::pair<types::global_dof_index, unsigned int>
608  (temp_copy_indices[block][i], i);
609  Assert (counter == n_active_dofs, ExcInternalError());
610  }
611  }
612 }
613 
614 
615 
616 // explicit instantiations
617 #include "mg_transfer_block.inst"
618 
619 
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
Definition: types.h:173
ElementIterator end() const
Definition: index_set.h:1516
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:723
cell_iterator end() const
Definition: dof_handler.cc:751
ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:229
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
Definition: dof_handler.cc:735
const IndexSet & get_boundary_indices(const unsigned int level) const
ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:257
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::shared_ptr< BlockSparseMatrix< double > > > prolongation_matrices
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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
const unsigned int dofs_per_cell
Definition: fe_base.h:295
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)
Definition: vector.h:1312
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:1836
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
Definition: fe.h:3176
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 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
Definition: types.h:187
ElementIterator begin() const
Definition: index_set.h:1453
BlockType & block(const unsigned int i)
size_type n_elements() const
Definition: index_set.h:1717
std::vector< bool > selected
static ::ExceptionBase & ExcInternalError()